v0.8.23
Classes | Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
MoFEM::VolumeElementForcesAndSourcesCoreBase Struct Reference

Volume finite element baseUser is implementing own operator at Gauss point level, by class derived from VolumeElementForcesAndSourcesCore::UserDataOperator. Arbitrary number of operator can be added by pushing objects to OpPtrVector. More...

#include <src/finite_elements/VolumeElementForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::VolumeElementForcesAndSourcesCoreBase:
[legend]
Collaboration diagram for MoFEM::VolumeElementForcesAndSourcesCoreBase:
[legend]

Classes

struct  UserDataOperator
 default operator for TET element More...
 

Public Types

enum  Switches { NO_HO_GEOMETRY = 1 << 0 | 1 << 2, NO_TRANSFORM = 1 << 1 | 1 << 2, NO_HO_TRANSFORM = 1 << 2 }
 
- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION, CTX_OPERATORS, CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION, CTX_SNESSETJACOBIAN, CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION, CTX_TSSETRHSJACOBIAN, CTX_TSSETIFUNCTION, CTX_TSSETIJACOBIAN,
  CTX_TSTSMONITORSET, CTX_TSNONE
}
 

Public Member Functions

 VolumeElementForcesAndSourcesCoreBase (Interface &m_field, const EntityType type=MBTET)
 
template<int SWITCH>
MoFEMErrorCode OpSwitch ()
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
boost::ptr_vector< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator. More...
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object. More...
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
int getMaxDataOrder () const
 Get max order of approximation for data fields. More...
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows. More...
 
int getMaxColOrder () const
 Get max order of approximation for field in columns. More...
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. More...
 
const DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side) const
 Get the entity data. More...
 
DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
RuleHookFun setRuleHook
 Set function to calculate integration rule. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 
Vec snes_x
 
Vec snes_f
 
Mat snes_A
 
Mat snes_B
 
- Public Attributes inherited from MoFEM::TSMethod
TSContext ts_ctx
 
TS ts
 
Vec ts_u
 
Vec ts_u_t
 
Vec ts_F
 
Mat ts_A
 
Mat ts_B
 
PetscInt ts_step
 
PetscReal ts_a
 
PetscReal ts_t
 

Protected Member Functions

virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
virtual MoFEMErrorCode calculateVolumeAndJacobian ()
 Calculate element volume and Jacobian. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode transformBaseFunctions ()
 Transform base functions based on geometric element Jacobian. More...
 
virtual MoFEMErrorCode calculateHoJacobian ()
 Calculate Jacobian for HO geometry. More...
 
virtual MoFEMErrorCode transformHoBaseFunctions ()
 Transform base functions based on ho-geometry element Jacobian. More...
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (DataForcesAndSourcesCore &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &data) const
 Get field approximation space and base on entities. More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate base functions. More...
 
MoFEMErrorCode createDataOnElement ()
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
MoFEMErrorCode getNodesIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEntityIndices (DataForcesAndSourcesCore &data, const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNodesFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 
MoFEMErrorCode getNoFieldFieldData (DataForcesAndSourcesCore &data, const boost::string_ref field_name) const
 
MoFEMErrorCode getNodesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getProblemNodesIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get indices of nodal indices which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemTypeIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const
 get indices by type (generic function) which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemNodesRowIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeRowIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
MoFEMErrorCode getProblemNodesColIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeColIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
virtual int getRule (int order)
 
virtual MoFEMErrorCode setGaussPts (int order)
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Protected Attributes

VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble hoGaussPtsJac
 
MatrixDouble hoGaussPtsInvJac
 
VectorDouble hoGaussPtsDetJac
 
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
 higher order geometry data at Gauss pts More...
 
OpSetHoInvJacH1 opSetHoInvJacH1
 
OpSetHoContravariantPiolaTransform opHoContravariantTransform
 
OpSetHoCovariantPiolaTransform opHoCovariantTransform
 
OpSetHoInvJacHdivAndHcurl opSetHoInvJacHdivAndHcurl
 
double vOlume
 
int num_nodes
 
const EntityHandleconn
 
FTensor::Tensor2< double *, 3, 3 > tJac
 
FTensor::Tensor2< double *, 3, 3 > tInvJac
 
MatrixDouble coordsAtGaussPts
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const boost::shared_ptr< DataForcesAndSourcesCoredataOnElement [LASTSPACE]
 Entity data on element entity rows fields. More...
 
const boost::shared_ptr< DataForcesAndSourcesCorederivedDataOnElement [LASTSPACE]
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 

Friends

class UserDataOperator
 

Detailed Description

Volume finite element base

User is implementing own operator at Gauss point level, by class derived from VolumeElementForcesAndSourcesCore::UserDataOperator. Arbitrary number of operator can be added by pushing objects to OpPtrVector.

Definition at line 38 of file VolumeElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ Switches

Constructor & Destructor Documentation

◆ VolumeElementForcesAndSourcesCoreBase()

MoFEM::VolumeElementForcesAndSourcesCoreBase::VolumeElementForcesAndSourcesCoreBase ( Interface m_field,
const EntityType  type = MBTET 
)

Definition at line 23 of file VolumeElementForcesAndSourcesCore.cpp.

25  : ForcesAndSourcesCore(m_field),
26  meshPositionsFieldName("MESH_NODE_POSITIONS"), coords(12), jAc(3, 3),
27  invJac(3, 3), opSetInvJacH1(invJac),
35  tJac(&jAc(0, 0), &jAc(0, 1), &jAc(0, 2), &jAc(1, 0), &jAc(1, 1),
36  &jAc(1, 2), &jAc(2, 0), &jAc(2, 1), &jAc(2, 2)),
37  tInvJac(&invJac(0, 0), &invJac(0, 1), &invJac(0, 2), &invJac(1, 0),
38  &invJac(1, 1), &invJac(1, 2), &invJac(2, 0), &invJac(2, 1),
39  &invJac(2, 2)) {
41  boost::shared_ptr<BaseFunction>(new TetPolynomialBase());
42 }
ForcesAndSourcesCore(Interface &m_field)
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
higher order geometry data at Gauss pts
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

Member Function Documentation

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::calculateCoordinatesAtGaussPts ( )
protectedvirtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 119 of file VolumeElementForcesAndSourcesCore.cpp.

119  {
121  // Get coords at Gauss points
123 
124  double *shape_functions_ptr =
125  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
126  const size_t nb_gauss_pts = gaussPts.size2();
127  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
128  coordsAtGaussPts.clear();
129  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_coords_at_gauss_ptr(
130  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
131  &coordsAtGaussPts(0, 2));
133  shape_functions_ptr);
134  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
136  &coords[0], &coords[1], &coords[2]);
137  for (int bb = 0; bb != 4; ++bb) {
138  t_coords_at_gauss_ptr(i) += t_coords(i) * t_shape_functions;
139  ++t_coords;
140  ++t_shape_functions;
141  };
142  ++t_coords_at_gauss_ptr;
143  }
145 }
MatrixDouble gaussPts
Matrix of integration points.
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ calculateHoJacobian()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::calculateHoJacobian ( )
protectedvirtual

Calculate Jacobian for HO geometry.

MoFEM use hierarchical approximate base to describe geometry of the body. This function transform derivatives of base functions when HO geometry is set and calculate Jacobian, inverse of Jacobian and determinant of transformation.

Definition at line 218 of file VolumeElementForcesAndSourcesCore.cpp.

218  {
220  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
221  dataPtr->get<FieldName_mi_tag>().end()) {
222  const Field *field_struture =
224  BitFieldId id = field_struture->getId();
225  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none()) {
226  SETERRQ(mField.get_comm(), MOFEM_NOT_FOUND,
227  "no MESH_NODE_POSITIONS in element data");
228  }
229 
231  if (dataH1.dataOnEntities[MBVERTEX][0].getFieldData().size() != 12) {
232  SETERRQ(mField.get_comm(), MOFEM_NOT_FOUND,
233  "no MESH_NODE_POSITIONS in element data or field has wrong "
234  "number of coefficients");
235  }
238  hoGaussPtsInvJac.resize(hoGaussPtsJac.size1(), hoGaussPtsJac.size2(),
239  false);
240  // Express Jacobian as tensor
242  &hoGaussPtsJac(0, 0), &hoGaussPtsJac(0, 1), &hoGaussPtsJac(0, 2),
243  &hoGaussPtsJac(0, 3), &hoGaussPtsJac(0, 4), &hoGaussPtsJac(0, 5),
244  &hoGaussPtsJac(0, 6), &hoGaussPtsJac(0, 7), &hoGaussPtsJac(0, 8));
246  &hoGaussPtsInvJac(0, 0), &hoGaussPtsInvJac(0, 1),
247  &hoGaussPtsInvJac(0, 2), &hoGaussPtsInvJac(0, 3),
248  &hoGaussPtsInvJac(0, 4), &hoGaussPtsInvJac(0, 5),
249  &hoGaussPtsInvJac(0, 6), &hoGaussPtsInvJac(0, 7),
250  &hoGaussPtsInvJac(0, 8));
251 
252  const size_t nb_gauss_pts = gaussPts.size2();
253  hoGaussPtsDetJac.resize(nb_gauss_pts, false);
255  // Calculate inverse and determinant
256  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
257  CHKERR determinantTensor3by3(jac, det);
258  // if(det<0) {
259  // SETERRQ(mField.get_comm(),MOFEM_DATA_INCONSISTENCY,"Negative
260  // volume");
261  // }
262  CHKERR invertTensor3by3(jac, det, inv_jac);
263  ++jac;
264  ++inv_jac;
265  ++det;
266  }
267  } else {
268  hoCoordsAtGaussPts.resize(0, 0, false);
269  hoGaussPtsInvJac.resize(0, 0, false);
270  hoGaussPtsDetJac.resize(0, false);
271  }
273 }
MatrixDouble gaussPts
Matrix of integration points.
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant 3 by 3.
Definition: Templates.hpp:418
DataForcesAndSourcesCore & dataH1
MoFEMErrorCode getNodesFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
Get field data on nodes.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
Pointer to finite element data dofs.
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
higher order geometry data at Gauss pts
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
Definition: Templates.hpp:399
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:53
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MPI_Comm & get_comm() const =0
virtual const Field * get_field_structure(const std::string &name)=0
get field structure

◆ calculateVolumeAndJacobian()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::calculateVolumeAndJacobian ( )
protectedvirtual

Calculate element volume and Jacobian.

Note that at that point is assumed that geometry is exclusively defined by corner nodes.

Returns
Error code

Definition at line 94 of file VolumeElementForcesAndSourcesCore.cpp.

94  {
97  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
98  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
103  &coords[0], &coords[1], &coords[2]);
106  jAc.clear();
107  for (auto n : {0, 1, 2, 3}) {
108  tJac(i, j) += t_coords(i) * t_diff_n(j);
109  ++t_coords;
110  ++t_diff_n;
111  }
114  vOlume *= G_TET_W1[0] / 6.;
116 }
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant 3 by 3.
Definition: Templates.hpp:418
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
static constexpr std::array< double, 12 > diffShapeFunMBTET
Definition: Tools.hpp:116
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
Definition: Templates.hpp:399
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static const double G_TET_W1[]
Definition: fem_tools.h:495
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::getSpaceBaseAndOrderOnElement ( )
protectedvirtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 148 of file VolumeElementForcesAndSourcesCore.cpp.

148  {
150 
153  // H1
154  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1)) {
155  CHKERR getEntitySense<MBEDGE>(dataH1);
156  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
157  }
158  if ((dataH1.spacesOnEntities[MBTRI]).test(H1)) {
159  CHKERR getEntitySense<MBTRI>(dataH1);
160  CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
161  }
162  if ((dataH1.spacesOnEntities[MBTET]).test(H1)) {
163  CHKERR getEntityDataOrder<MBTET>(dataH1, H1);
164  }
165  // Hcurl
166  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL)) {
167  CHKERR getEntitySense<MBEDGE>(dataHcurl);
168  CHKERR getEntityDataOrder<MBEDGE>(dataHcurl, HCURL);
169  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
170  }
171  if ((dataH1.spacesOnEntities[MBTRI]).test(HCURL)) {
172  CHKERR getEntitySense<MBTRI>(dataHcurl);
174  CHKERR getEntityDataOrder<MBTRI>(dataHcurl, HCURL);
175  dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
176  }
177  if ((dataH1.spacesOnEntities[MBTET]).test(HCURL)) {
178  CHKERR getEntityDataOrder<MBTET>(dataHcurl, HCURL);
179  dataHcurl.spacesOnEntities[MBTET].set(HCURL);
180  }
181  // Hdiv
182  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV)) {
183  CHKERR getEntitySense<MBTRI>(dataHdiv);
185  CHKERR getEntityDataOrder<MBTRI>(dataHdiv, HDIV);
186  dataHdiv.spacesOnEntities[MBTRI].set(HDIV);
187  }
188  if ((dataH1.spacesOnEntities[MBTET]).test(HDIV)) {
189  CHKERR getEntityDataOrder<MBTET>(dataHdiv, HDIV);
190  dataHdiv.spacesOnEntities[MBTET].set(HDIV);
191  }
192  // L2
193  if ((dataH1.spacesOnEntities[MBTET]).test(L2)) {
194  CHKERR getEntityDataOrder<MBTET>(dataL2, L2);
195  dataL2.spacesOnEntities[MBTET].set(L2);
196  }
198 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:173
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MoFEMErrorCode getFaceTriNodes(DataForcesAndSourcesCore &data) const
Get nodes on triangles.
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
DataForcesAndSourcesCore & dataHdiv
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
continuous field
Definition: definitions.h:171
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
field with C-1 continuity
Definition: definitions.h:174

◆ OpSwitch()

template<int SWITCH>
MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::OpSwitch ( )

Definition at line 282 of file VolumeElementForcesAndSourcesCore.hpp.

282  {
284 
285  if (numeredEntFiniteElementPtr->getEntType() != MBTET)
288 
292  if (gaussPts.size2() == 0)
296 
297  if (!(NO_TRANSFORM & SWITCH))
299 
300  try {
301  MatrixDouble new_diff_n;
302  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
304  FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
306  dataH1.dataOnEntities[MBVERTEX][0];
307  if ((data.getDiffN(base).size1() == 4) &&
308  (data.getDiffN(base).size2() == 3)) {
309  const size_t nb_gauss_pts = gaussPts.size2();
310  const size_t nb_base_functions = 4;
311  new_diff_n.resize(nb_gauss_pts, 3 * nb_base_functions, false);
312  double *new_diff_n_ptr = &*new_diff_n.data().begin();
314  new_diff_n_ptr, &new_diff_n_ptr[1], &new_diff_n_ptr[2]);
315  double *t_diff_n_ptr = &*data.getDiffN(base).data().begin();
316  for (unsigned int gg = 0; gg != nb_gauss_pts; gg++) {
318  t_diff_n_ptr, &t_diff_n_ptr[1], &t_diff_n_ptr[2]);
319  for (unsigned int bb = 0; bb != nb_base_functions; bb++) {
320  t_new_diff_n(i) = t_diff_n(i);
321  ++t_new_diff_n;
322  ++t_diff_n;
323  }
324  }
325  data.getDiffN(base).resize(new_diff_n.size1(), new_diff_n.size2(),
326  false);
327  data.getDiffN(base).data().swap(new_diff_n.data());
328  }
329  }
330  }
331  CATCH_ERRORS;
332 
333  if (!(NO_HO_GEOMETRY & SWITCH))
335 
336  if (!(NO_HO_TRANSFORM & SWITCH))
338 
339  // Iterate over operators
341 
343  }
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
virtual MoFEMErrorCode calculateVolumeAndJacobian()
Calculate element volume and Jacobian.
MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
MatrixDouble gaussPts
Matrix of integration points.
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
DataForcesAndSourcesCore & dataH1
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:77
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
virtual MoFEMErrorCode transformBaseFunctions()
Transform base functions based on geometric element Jacobian.
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
virtual MoFEMErrorCode calculateHoJacobian()
Calculate Jacobian for HO geometry.
FieldApproximationBase
approximation base
Definition: definitions.h:143
virtual MoFEMErrorCode transformHoBaseFunctions()
Transform base functions based on ho-geometry element Jacobian.
DataForcesAndSourcesCore::EntData EntData
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:145
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:433

◆ setIntegrationPts()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::setIntegrationPts ( )
protectedvirtual

Set integration points.

Returns
Error code

Definition at line 44 of file VolumeElementForcesAndSourcesCore.cpp.

44  {
46  int order_data = getMaxDataOrder();
47  int order_row = getMaxRowOrder();
48  int order_col = getMaxColOrder();
49  int rule = getRule(order_row, order_col, order_data);
50 
51  if (rule >= 0) {
52  if (rule < QUAD_3D_TABLE_SIZE) {
53  if (QUAD_3D_TABLE[rule]->dim != 3) {
54  SETERRQ(mField.get_comm(), MOFEM_DATA_INCONSISTENCY, "wrong dimension");
55  }
56  if (QUAD_3D_TABLE[rule]->order < rule) {
58  "wrong order %d != %d", QUAD_3D_TABLE[rule]->order, rule);
59  }
60  size_t nb_gauss_pts = QUAD_3D_TABLE[rule]->npoints;
61  gaussPts.resize(4, nb_gauss_pts, false);
62  cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[1], 4,
63  &gaussPts(0, 0), 1);
64  cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[2], 4,
65  &gaussPts(1, 0), 1);
66  cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[3], 4,
67  &gaussPts(2, 0), 1);
68  cblas_dcopy(nb_gauss_pts, QUAD_3D_TABLE[rule]->weights, 1, &gaussPts(3, 0),
69  1);
70  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 4,
71  false);
72  double *shape_ptr =
73  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
74  cblas_dcopy(4 * nb_gauss_pts, QUAD_3D_TABLE[rule]->points, 1, shape_ptr, 1);
75  } else {
77  "rule > quadrature order %d < %d", rule, QUAD_3D_TABLE_SIZE);
78  }
79  } else {
80  CHKERR setGaussPts(order_row, order_col, order_data);
81  const size_t nb_gauss_pts = gaussPts.size2();
82  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 4,
83  false);
84  if (nb_gauss_pts > 0) {
86  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
87  &gaussPts(0, 0), &gaussPts(1, 0), &gaussPts(2, 0), nb_gauss_pts);
88  }
89  }
91 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
#define QUAD_3D_TABLE_SIZE
Definition: quad.h:186
MatrixDouble gaussPts
Matrix of integration points.
DataForcesAndSourcesCore & dataH1
int npoints
Definition: quad.h:29
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
#define CHKERR
Inline error check.
Definition: definitions.h:596
static MoFEMErrorCode shapeFunMBTET(double *shape, const double *ksi, const double *eta, const double *zeta, const double nb)
Calculate shape functions on tetrahedron.
Definition: Tools.hpp:176
int order
Definition: quad.h:28
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MPI_Comm & get_comm() const =0
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
int getMaxColOrder() const
Get max order of approximation for field in columns.
static QUAD *const QUAD_3D_TABLE[]
Definition: quad.h:187
int getMaxDataOrder() const
Get max order of approximation for data fields.
const int order
approximation order

◆ transformBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::transformBaseFunctions ( )
protectedvirtual

Transform base functions based on geometric element Jacobian.

This function apply transformation to base functions and its derivatives. For example when base functions for H-div are present the Piola-Transformarion is applied to base functions and their derivatives.

Returns
Error code

Definition at line 200 of file VolumeElementForcesAndSourcesCore.cpp.

200  {
202 
204  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
207  }
208  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
211  }
212  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
214  }
216 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:173
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
DataForcesAndSourcesCore & dataHdiv
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
field with C-1 continuity
Definition: definitions.h:174

◆ transformHoBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::transformHoBaseFunctions ( )
protectedvirtual

Transform base functions based on ho-geometry element Jacobian.

This function apply transformation to base functions and its derivatives. For example when base functions for H-div are present the Piola-Transformarion is applied to base functions and their derivatives.

Returns
Error code

Definition at line 276 of file VolumeElementForcesAndSourcesCore.cpp.

276  {
278  if (hoCoordsAtGaussPts.size1() > 0) {
279  // Transform derivatives of base functions and apply Piola transformation
280  // if needed.
281 
283  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
285  }
286  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
289  }
290  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
293  }
294  }
296 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:173
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
DataForcesAndSourcesCore & dataHdiv
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
field with C-1 continuity
Definition: definitions.h:174

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

Definition at line 439 of file VolumeElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ conn

const EntityHandle* MoFEM::VolumeElementForcesAndSourcesCoreBase::conn
protected

Definition at line 433 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::VolumeElementForcesAndSourcesCoreBase::coords
protected

Definition at line 409 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCoreBase::coordsAtGaussPts
protected

Definition at line 437 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCoreBase::hoCoordsAtGaussPts
protected

Definition at line 418 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsDetJac

VectorDouble MoFEM::VolumeElementForcesAndSourcesCoreBase::hoGaussPtsDetJac
protected

Definition at line 421 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsInvJac

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCoreBase::hoGaussPtsInvJac
protected

Definition at line 420 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsJac

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCoreBase::hoGaussPtsJac
protected

Definition at line 419 of file VolumeElementForcesAndSourcesCore.hpp.

◆ invJac

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCoreBase::invJac
protected

Definition at line 411 of file VolumeElementForcesAndSourcesCore.hpp.

◆ jAc

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCoreBase::jAc
protected

Definition at line 410 of file VolumeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::VolumeElementForcesAndSourcesCoreBase::meshPositionsFieldName

Definition at line 40 of file VolumeElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::VolumeElementForcesAndSourcesCoreBase::num_nodes
protected

Definition at line 432 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opContravariantPiolaTransform

OpSetContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCoreBase::opContravariantPiolaTransform
protected

Definition at line 414 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opCovariantPiolaTransform

OpSetCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCoreBase::opCovariantPiolaTransform
protected

Definition at line 415 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHOatGaussPoints

OpGetDataAndGradient<3, 3> MoFEM::VolumeElementForcesAndSourcesCoreBase::opHOatGaussPoints
protected

higher order geometry data at Gauss pts

Definition at line 424 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHoContravariantTransform

OpSetHoContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCoreBase::opHoContravariantTransform
protected

Definition at line 426 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHoCovariantTransform

OpSetHoCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCoreBase::opHoCovariantTransform
protected

Definition at line 427 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetHoInvJacH1

OpSetHoInvJacH1 MoFEM::VolumeElementForcesAndSourcesCoreBase::opSetHoInvJacH1
protected

Definition at line 425 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetHoInvJacHdivAndHcurl

OpSetHoInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCoreBase::opSetHoInvJacHdivAndHcurl
protected

Definition at line 428 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacH1

OpSetInvJacH1 MoFEM::VolumeElementForcesAndSourcesCoreBase::opSetInvJacH1
protected

Definition at line 413 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacHdivAndHcurl

OpSetInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCoreBase::opSetInvJacHdivAndHcurl
protected

Definition at line 416 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tInvJac

FTensor::Tensor2<double *, 3, 3> MoFEM::VolumeElementForcesAndSourcesCoreBase::tInvJac
protected

Definition at line 435 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tJac

FTensor::Tensor2<double *, 3, 3> MoFEM::VolumeElementForcesAndSourcesCoreBase::tJac
protected

Definition at line 434 of file VolumeElementForcesAndSourcesCore.hpp.

◆ vOlume

double MoFEM::VolumeElementForcesAndSourcesCoreBase::vOlume
protected

Definition at line 430 of file VolumeElementForcesAndSourcesCore.hpp.


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