v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM > Struct Template Reference

Calculate gradient of tensor field. More...

#include "src/finite_elements/UserDataOperators.hpp"

Inheritance diagram for MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >:
[legend]
Collaboration diagram for MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >:
[legend]

Public Member Functions

 OpCalculateHVecTensorGradient (const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side.
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 Constructor for operators working on finite element spaces.
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 Constructor for operators working on a single field.
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 Constructor for operators working on two fields (bilinear forms)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side.
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Private Attributes

boost::shared_ptr< MatrixDoubledataPtr
 
const EntityHandle zeroType
 
const int zeroSide
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
struct MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >

Calculate gradient of tensor field.

Template Parameters
BASE_DIM
FIELD_DIM
SPACE_DIM
Examples
mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

Definition at line 2616 of file UserDataOperators.hpp.

Constructor & Destructor Documentation

◆ OpCalculateHVecTensorGradient()

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >::OpCalculateHVecTensorGradient ( const std::string  field_name,
boost::shared_ptr< MatrixDouble data_ptr,
const EntityType  zero_type = MBEDGE,
const int  zero_side = 0 
)
inline

Definition at line 2619 of file UserDataOperators.hpp.

2623 : ForcesAndSourcesCore::UserDataOperator(
2624 field_name, ForcesAndSourcesCore::UserDataOperator::OPROW),
2625 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2626 if (!dataPtr)
2627 THROW_MESSAGE("Pointer is not set");
2628 }
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
constexpr auto field_name
boost::shared_ptr< MatrixDouble > dataPtr

Member Function Documentation

◆ doWork()

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
MoFEMErrorCode MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
inlinevirtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 2630 of file UserDataOperators.hpp.

2631 {
2633 const size_t nb_integration_points = getGaussPts().size2();
2634 if (type == zeroType && side == zeroSide) {
2635 dataPtr->resize(BASE_DIM * SPACE_DIM * FIELD_DIM, nb_integration_points,
2636 false);
2637 dataPtr->clear();
2638 }
2639
2640 #ifndef NDEBUG
2641 if (data.getFieldData().size() % FIELD_DIM) {
2642 SETERRQ(
2643 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2644 "Data inconsistency, nb_dofs %% FIELD_DIM != 0, that is %ld %% %d "
2645 "!= 0",
2646 data.getFieldData().size(), FIELD_DIM);
2647 }
2648 #endif
2649
2650 const auto nb_dofs = data.getFieldData().size() / FIELD_DIM;
2651 if (!nb_dofs)
2653
2654 const size_t nb_base_functions = data.getN().size2() / BASE_DIM;
2658
2659 auto t_base_diff = data.getFTensor2DiffN<BASE_DIM, SPACE_DIM>();
2660 auto t_data = getFTensor3FromMat<BASE_DIM, FIELD_DIM, SPACE_DIM>(*dataPtr);
2661 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2662 auto t_dof = data.getFTensor1FieldData<FIELD_DIM>();
2663 int bb = 0;
2664 for (; bb != nb_dofs; ++bb) {
2665 t_data(k, i, j) += t_base_diff(i, j) * t_dof(k);
2666 ++t_base_diff;
2667 ++t_dof;
2668 }
2669 for (; bb != nb_base_functions; ++bb)
2670 ++t_base_diff;
2671 ++t_data;
2672 }
2674 }
constexpr int SPACE_DIM
constexpr int FIELD_DIM
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
constexpr int BASE_DIM
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

Member Data Documentation

◆ dataPtr

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
boost::shared_ptr<MatrixDouble> MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >::dataPtr
private

Definition at line 2677 of file UserDataOperators.hpp.

◆ zeroSide

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
const int MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >::zeroSide
private

Definition at line 2679 of file UserDataOperators.hpp.

◆ zeroType

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
const EntityHandle MoFEM::OpCalculateHVecTensorGradient< BASE_DIM, FIELD_DIM, SPACE_DIM >::zeroType
private

Definition at line 2678 of file UserDataOperators.hpp.


The documentation for this struct was generated from the following file: