v0.8.20
Classes | Public Member Functions | Public Attributes | List of all members
MoFEM::VolumeElementForcesAndSourcesCore Struct Reference

Volume finite elementUser 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::VolumeElementForcesAndSourcesCore:
[legend]
Collaboration diagram for MoFEM::VolumeElementForcesAndSourcesCore:
[legend]

Classes

struct  UserDataOperator
 default operator for TET element More...
 

Public Member Functions

 VolumeElementForcesAndSourcesCore (Interface &m_field, const EntityType type=MBTET)
 
virtual ~VolumeElementForcesAndSourcesCore ()
 
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
virtual DEPRECATED MoFEMErrorCode setIntegartionPts ()
 
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...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
virtual ~ForcesAndSourcesCore ()
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 
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 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 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 getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &data) const
 Get field approximation space and base on entities. More...
 
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)
 set integration rule for finite element More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order)
 It will be removed in the future use other variant. More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
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 postProcess ()
 function is run at the end of loop 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...
 
- 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

VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
std::string meshPositionsFieldName
 
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
 
unsigned int nbGaussPts
 Number of integration points. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
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
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. 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...
 
- 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
 

Additional Inherited Members

- 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
}
 
- 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...
 

Detailed Description

Volume finite element

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.

Examples:
elasticity.cpp, elasticity_mixed_formulation.cpp, forces_and_sources_testing_users_base.cpp, hello_world.cpp, HookeElement.cpp, MagneticElement.hpp, mesh_smoothing.cpp, Remodeling.hpp, simple_elasticity.cpp, simple_interface.cpp, SmallStrainPlasticity.hpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 38 of file VolumeElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ VolumeElementForcesAndSourcesCore()

MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore ( Interface m_field,
const EntityType  type = MBTET 
)
Examples:
SmallStrainPlasticity.hpp.

Definition at line 34 of file VolumeElementForcesAndSourcesCore.cpp.

36  : ForcesAndSourcesCore(m_field), coords(12), jAc(3, 3), invJac(3, 3),
39  meshPositionsFieldName("MESH_NODE_POSITIONS"),
45  tJac(&jAc(0, 0), &jAc(0, 1), &jAc(0, 2), &jAc(1, 0), &jAc(1, 1),
46  &jAc(1, 2), &jAc(2, 0), &jAc(2, 1), &jAc(2, 2)),
47  tInvJac(&invJac(0, 0), &invJac(0, 1), &invJac(0, 2), &invJac(1, 0),
48  &invJac(1, 1), &invJac(1, 2), &invJac(2, 0), &invJac(2, 1),
49  &invJac(2, 2)) {
51  boost::shared_ptr<BaseFunction>(new TetPolynomialBase());
52 }
ForcesAndSourcesCore(Interface &m_field)
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
higher order geometry data at Gauss pts
OpSetHoContravariantPiolaTransform opHoContravariantTransform
OpSetContravariantPiolaTransform opContravariantPiolaTransform
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

◆ ~VolumeElementForcesAndSourcesCore()

virtual MoFEM::VolumeElementForcesAndSourcesCore::~VolumeElementForcesAndSourcesCore ( )
virtual

Definition at line 64 of file VolumeElementForcesAndSourcesCore.hpp.

64 {}

Member Function Documentation

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
virtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 131 of file VolumeElementForcesAndSourcesCore.cpp.

131  {
133  // Get coords at Gauss points
135 
136  double *shape_functions_ptr =
137  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
138  coordsAtGaussPts.resize(nbGaussPts, 3, false);
139  coordsAtGaussPts.clear();
140  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_coords_at_gauss_ptr(
141  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
142  &coordsAtGaussPts(0, 2));
144  shape_functions_ptr);
145  for (unsigned int gg = 0; gg < nbGaussPts; gg++) {
147  &coords[0], &coords[1], &coords[2]);
148  for (int bb = 0; bb < 4; bb++) {
149  t_coords_at_gauss_ptr(i) += t_coords(i) * t_shape_functions;
150  ++t_coords;
151  ++t_shape_functions;
152  };
153  ++t_coords_at_gauss_ptr;
154  }
156 }
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
unsigned int nbGaussPts
Number of integration points.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

◆ calculateHoJacobian()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::calculateHoJacobian ( )
virtual

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 229 of file VolumeElementForcesAndSourcesCore.cpp.

229  {
231  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
232  dataPtr->get<FieldName_mi_tag>().end()) {
233  const Field *field_struture =
235  BitFieldId id = field_struture->getId();
236  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none()) {
237  SETERRQ(mField.get_comm(), MOFEM_NOT_FOUND,
238  "no MESH_NODE_POSITIONS in element data");
239  }
240 
242  if (dataH1.dataOnEntities[MBVERTEX][0].getFieldData().size() != 12) {
243  SETERRQ(mField.get_comm(), MOFEM_NOT_FOUND,
244  "no MESH_NODE_POSITIONS in element data or field has wrong "
245  "number of coefficients");
246  }
249  hoGaussPtsInvJac.resize(hoGaussPtsJac.size1(), hoGaussPtsJac.size2(),
250  false);
251  // Express Jacobian as tensor
253  &hoGaussPtsJac(0, 0), &hoGaussPtsJac(0, 1), &hoGaussPtsJac(0, 2),
254  &hoGaussPtsJac(0, 3), &hoGaussPtsJac(0, 4), &hoGaussPtsJac(0, 5),
255  &hoGaussPtsJac(0, 6), &hoGaussPtsJac(0, 7), &hoGaussPtsJac(0, 8));
257  &hoGaussPtsInvJac(0, 0), &hoGaussPtsInvJac(0, 1),
258  &hoGaussPtsInvJac(0, 2), &hoGaussPtsInvJac(0, 3),
259  &hoGaussPtsInvJac(0, 4), &hoGaussPtsInvJac(0, 5),
260  &hoGaussPtsInvJac(0, 6), &hoGaussPtsInvJac(0, 7),
261  &hoGaussPtsInvJac(0, 8));
262  hoGaussPtsDetJac.resize(nbGaussPts, false);
264  // Calculate inverse and determinant
265  for (unsigned int gg = 0; gg != nbGaussPts; gg++) {
266  CHKERR determinantTensor3by3(jac, det);
267  // if(det<0) {
268  // SETERRQ(mField.get_comm(),MOFEM_DATA_INCONSISTENCY,"Negative
269  // volume");
270  // }
271  CHKERR invertTensor3by3(jac, det, inv_jac);
272  ++jac;
273  ++inv_jac;
274  ++det;
275  }
276  } else {
277  hoCoordsAtGaussPts.resize(0, 0, false);
278  hoGaussPtsInvJac.resize(0, 0, false);
279  hoGaussPtsDetJac.resize(0, false);
280  }
282 }
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.
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
higher order geometry data at Gauss pts
DataForcesAndSourcesCore & dataH1
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant.
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:475
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
Pointer to finite element data dofs.
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
unsigned int nbGaussPts
Number of integration points.
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:53
#define CHKERR
Inline error check.
Definition: definitions.h:594
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:405
virtual MPI_Comm & get_comm() const =0
virtual const Field * get_field_structure(const std::string &name)=0
get field structure

◆ calculateVolumeAndJacobian()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::calculateVolumeAndJacobian ( )
virtual

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 105 of file VolumeElementForcesAndSourcesCore.cpp.

105  {
108  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
109  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
110  double diff_n[12];
111  CHKERR ShapeDiffMBTET(diff_n);
113  &diff_n[0], &diff_n[1], &diff_n[2]);
115  &coords[0], &coords[1], &coords[2]);
118  jAc.clear();
119  for (int nn = 0; nn != 4; nn++) {
120  tJac(i, j) += t_coords(i) * t_diff_n(j);
121  ++t_coords;
122  ++t_diff_n;
123  }
126  vOlume *= G_TET_W1[0] / 6.;
128 }
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.
virtual moab::Interface & get_moab()=0
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:303
static const double G_TET_W1[]
Definition: fem_tools.h:495
#define CHKERR
Inline error check.
Definition: definitions.h:594
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
virtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 159 of file VolumeElementForcesAndSourcesCore.cpp.

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

◆ operator()()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::operator() ( )
virtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from MoFEM::ForcesAndSourcesCore.

Definition at line 306 of file VolumeElementForcesAndSourcesCore.cpp.

306  {
308 
309  if (numeredEntFiniteElementPtr->getEntType() != MBTET)
312 
316  if (nbGaussPts == 0)
321 
322  try {
323  MatrixDouble new_diff_n;
324  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
326  FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
328  dataH1.dataOnEntities[MBVERTEX][0];
329  if ((data.getDiffN(base).size1() == 4) &&
330  (data.getDiffN(base).size2() == 3)) {
331  const unsigned int nb_base_functions = 4;
332  new_diff_n.resize(nbGaussPts, 3 * nb_base_functions, false);
333  double *new_diff_n_ptr = &*new_diff_n.data().begin();
335  new_diff_n_ptr, &new_diff_n_ptr[1], &new_diff_n_ptr[2]);
336  double *t_diff_n_ptr = &*data.getDiffN(base).data().begin();
337  for (unsigned int gg = 0; gg != nbGaussPts; gg++) {
339  t_diff_n_ptr, &t_diff_n_ptr[1], &t_diff_n_ptr[2]);
340  for (unsigned int bb = 0; bb != nb_base_functions; bb++) {
341  t_new_diff_n(i) = t_diff_n(i);
342  ++t_new_diff_n;
343  ++t_diff_n;
344  }
345  }
346  data.getDiffN(base).resize(new_diff_n.size1(), new_diff_n.size2(),
347  false);
348  data.getDiffN(base).data().swap(new_diff_n.data());
349  }
350  }
351  }
352  CATCH_ERRORS;
353 
356 
357  // Iterate over operators
359 
361 }
MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
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:475
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
unsigned int nbGaussPts
Number of integration points.
virtual MoFEMErrorCode calculateHoJacobian()
Calculate Jacobian for HO geometry.
FieldApproximationBase
approximation base
Definition: definitions.h:141
virtual MoFEMErrorCode calculateVolumeAndJacobian()
Calculate element volume and Jacobian.
virtual MoFEMErrorCode transformBaseFunctions()
Transform base functions based on geometric element Jacobian.
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:143
DataForcesAndSourcesCore::EntData EntData
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
#define CHKERR
Inline error check.
Definition: definitions.h:594
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:431
virtual MoFEMErrorCode transformHoBaseFunctions()
Transform base functions based on ho-geometry element Jacobian.

◆ setIntegartionPts()

virtual DEPRECATED MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::setIntegartionPts ( )
virtual
Deprecated:
function with spelling mistake, use setIntegrationPts

Definition at line 322 of file VolumeElementForcesAndSourcesCore.hpp.

322  {
323  return setIntegrationPts();
324  }
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.

◆ setIntegrationPts()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::setIntegrationPts ( )
virtual

Set integration points.

Returns
Error code

Definition at line 54 of file VolumeElementForcesAndSourcesCore.cpp.

54  {
56  int order_data = getMaxDataOrder();
57  int order_row = getMaxRowOrder();
58  int order_col = getMaxColOrder();
59  int rule = getRule(order_row, order_col, order_data);
60 
61  // rule << std::endl;
62  if (rule >= 0) {
63  if (rule < QUAD_3D_TABLE_SIZE) {
64  if (QUAD_3D_TABLE[rule]->dim != 3) {
65  SETERRQ(mField.get_comm(), MOFEM_DATA_INCONSISTENCY, "wrong dimension");
66  }
67  if (QUAD_3D_TABLE[rule]->order < rule) {
69  "wrong order %d != %d", QUAD_3D_TABLE[rule]->order, rule);
70  }
72  gaussPts.resize(4, nbGaussPts, false);
73  cblas_dcopy(nbGaussPts, &QUAD_3D_TABLE[rule]->points[1], 4,
74  &gaussPts(0, 0), 1);
75  cblas_dcopy(nbGaussPts, &QUAD_3D_TABLE[rule]->points[2], 4,
76  &gaussPts(1, 0), 1);
77  cblas_dcopy(nbGaussPts, &QUAD_3D_TABLE[rule]->points[3], 4,
78  &gaussPts(2, 0), 1);
79  cblas_dcopy(nbGaussPts, QUAD_3D_TABLE[rule]->weights, 1, &gaussPts(3, 0),
80  1);
81  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 4,
82  false);
83  double *shape_ptr =
84  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
85  cblas_dcopy(4 * nbGaussPts, QUAD_3D_TABLE[rule]->points, 1, shape_ptr, 1);
86  } else {
88  "rule > quadrature order %d < %d", rule, QUAD_3D_TABLE_SIZE);
89  nbGaussPts = 0;
90  }
91  } else {
92  CHKERR setGaussPts(order_row, order_col, order_data);
93  nbGaussPts = gaussPts.size2();
94  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 4,
95  false);
96  if (nbGaussPts > 0) {
98  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
99  &gaussPts(0, 0), &gaussPts(1, 0), &gaussPts(2, 0), nbGaussPts);
100  }
101  }
103 }
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.
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition: fem_tools.c:291
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:475
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
virtual MoFEMErrorCode setGaussPts(int order)
It will be removed in the future use other variant.
unsigned int nbGaussPts
Number of integration points.
#define CHKERR
Inline error check.
Definition: definitions.h:594
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:405
virtual MPI_Comm & get_comm() const =0
virtual int getRule(int order)
set integration rule for finite element
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.

◆ transformBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::transformBaseFunctions ( )
virtual

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 211 of file VolumeElementForcesAndSourcesCore.cpp.

211  {
213 
215  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
218  }
219  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
222  }
223  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
225  }
227 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:171
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:170
#define CHKERR
Inline error check.
Definition: definitions.h:594
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
OpSetContravariantPiolaTransform opContravariantPiolaTransform
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405
field with C-1 continuity
Definition: definitions.h:172

◆ transformHoBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::transformHoBaseFunctions ( )
virtual

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 284 of file VolumeElementForcesAndSourcesCore.cpp.

284  {
286  if (hoCoordsAtGaussPts.size1() > 0) {
287  // Transform derivatives of base functions and apply Piola transformation
288  // if needed.
289 
291  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
293  }
294  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
297  }
298  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
301  }
302  }
304 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:171
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
OpSetHoContravariantPiolaTransform opHoContravariantTransform
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:170
#define CHKERR
Inline error check.
Definition: definitions.h:594
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:405
field with C-1 continuity
Definition: definitions.h:172

Member Data Documentation

◆ conn

const EntityHandle* MoFEM::VolumeElementForcesAndSourcesCore::conn

Definition at line 69 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::VolumeElementForcesAndSourcesCore::coords

Definition at line 40 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::coordsAtGaussPts

Definition at line 73 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::hoCoordsAtGaussPts

Definition at line 50 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsDetJac

VectorDouble MoFEM::VolumeElementForcesAndSourcesCore::hoGaussPtsDetJac

Definition at line 53 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsInvJac

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::hoGaussPtsInvJac

Definition at line 52 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsJac

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::hoGaussPtsJac

Definition at line 51 of file VolumeElementForcesAndSourcesCore.hpp.

◆ invJac

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCore::invJac

Definition at line 42 of file VolumeElementForcesAndSourcesCore.hpp.

◆ jAc

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCore::jAc

Definition at line 41 of file VolumeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::VolumeElementForcesAndSourcesCore::meshPositionsFieldName

Definition at line 49 of file VolumeElementForcesAndSourcesCore.hpp.

◆ nbGaussPts

unsigned int MoFEM::VolumeElementForcesAndSourcesCore::nbGaussPts

Number of integration points.

Definition at line 310 of file VolumeElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::VolumeElementForcesAndSourcesCore::num_nodes

Definition at line 68 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opContravariantPiolaTransform

OpSetContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opContravariantPiolaTransform

Definition at line 45 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opCovariantPiolaTransform

OpSetCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opCovariantPiolaTransform

Definition at line 46 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHOatGaussPoints

OpGetDataAndGradient<3, 3> MoFEM::VolumeElementForcesAndSourcesCore::opHOatGaussPoints

higher order geometry data at Gauss pts

Definition at line 56 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHoContravariantTransform

OpSetHoContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opHoContravariantTransform

Definition at line 58 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHoCovariantTransform

OpSetHoCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opHoCovariantTransform

Definition at line 59 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetHoInvJacH1

OpSetHoInvJacH1 MoFEM::VolumeElementForcesAndSourcesCore::opSetHoInvJacH1

Definition at line 57 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetHoInvJacHdivAndHcurl

OpSetHoInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCore::opSetHoInvJacHdivAndHcurl

Definition at line 60 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacH1

OpSetInvJacH1 MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacH1

Definition at line 44 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacHdivAndHcurl

OpSetInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacHdivAndHcurl

Definition at line 47 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tInvJac

FTensor::Tensor2<double *, 3, 3> MoFEM::VolumeElementForcesAndSourcesCore::tInvJac

Definition at line 71 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tJac

FTensor::Tensor2<double *, 3, 3> MoFEM::VolumeElementForcesAndSourcesCore::tJac

Definition at line 70 of file VolumeElementForcesAndSourcesCore.hpp.

◆ vOlume

double MoFEM::VolumeElementForcesAndSourcesCore::vOlume

Definition at line 66 of file VolumeElementForcesAndSourcesCore.hpp.


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