v0.8.4
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 calculateBaseFunctionsOnElement (const int b)
 Calculate base functions. More...
 
virtual MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate 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 getSense (EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order of approximation on the entity More...
 
MoFEMErrorCode getDataOrderSpaceAndBase (const std::string &field_name, const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order on entity More...
 
MoFEMErrorCode getEdgesSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getTrisSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getQuadSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getEdgesDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getTrisDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getQuadDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getTetDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getPrismDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getEdgesDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTrisDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getQuadDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTetDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getPrismDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getNodesIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
 get node indices More...
 
MoFEMErrorCode getTypeIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices, VectorInt &local_indices) const
 get indices by type (generic function) More...
 
MoFEMErrorCode getTypeIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get indices by type (generic function) 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 getEdgesRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Edges row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEdgesColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Edges col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTrisRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tris row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTrisColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tris col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTetsRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tets row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTetsColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tets col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getQuadRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Quad row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getQuadColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Quad col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getPrismRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Prism row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getPrismColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Prism col indices from FENumeredDofEntity_multiIndex More...
 
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 getTypeFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 Get field data on entities. More...
 
MoFEMErrorCode getTypeFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 
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 getEdgesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTrisFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getQuadFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTetsFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getPrismFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) 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...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop 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)
 
- 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
 
DataForcesAndSourcesCore dataH1
 
DerivedDataForcesAndSourcesCore derivedDataH1
 
DataForcesAndSourcesCore dataL2
 
DerivedDataForcesAndSourcesCore derivedDataL2
 
DataForcesAndSourcesCore dataHdiv
 
DerivedDataForcesAndSourcesCore derivedDataHdiv
 
DataForcesAndSourcesCore dataHcurl
 
DerivedDataForcesAndSourcesCore derivedDataHcurl
 
DataForcesAndSourcesCore dataNoField
 
DataForcesAndSourcesCore dataNoFieldCol
 
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 gaussPts
 
MatrixDouble coordsAtGaussPts
 
unsigned int nbGaussPts
 Number of integration points. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 
int loopSize
 
int rAnk
 
int sIze
 
const RefEntity_multiIndexrefinedEntitiesPtr
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 
const ProblemproblemPtr
 
const Field_multiIndexfieldsPtr
 
const FieldEntity_multiIndexentitiesPtr
 
const DofEntity_multiIndexdofsPtr
 
const FiniteElement_multiIndexfiniteElementsPtr
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 
boost::function< MoFEMErrorCode()> postProcessHook
 
boost::function< MoFEMErrorCode()> operatorHook
 
- 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:
hello_world.cpp, MagneticElement.hpp, mesh_smoothing.cpp, Remodeling.hpp, simple_elasticity.cpp, and simple_interface.cpp.

Definition at line 39 of file VolumeElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ VolumeElementForcesAndSourcesCore()

VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore ( Interface m_field,
const EntityType  type = MBTET 
)

Definition at line 79 of file VolumeElementForcesAndSourcesCore.cpp.

81  : ForcesAndSourcesCore(m_field), coords(12), jAc(3, 3), invJac(3, 3),
87  meshPositionsFieldName("MESH_NODE_POSITIONS"),
93  tJac(&jAc(0, 0), &jAc(0, 1), &jAc(0, 2), &jAc(1, 0), &jAc(1, 1),
94  &jAc(1, 2), &jAc(2, 0), &jAc(2, 1), &jAc(2, 2)),
95  tInvJac(&invJac(0, 0), &invJac(0, 1), &invJac(0, 2), &invJac(1, 0),
96  &invJac(1, 1), &invJac(1, 2), &invJac(2, 0), &invJac(2, 1),
97  &invJac(2, 2)) {}
ForcesAndSourcesCore(Interface &m_field)
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
higher order geometry data at Gauss pts
OpSetHoContravariantPiolaTransform opHoContravariantTransform
OpSetContravariantPiolaTransform opContravariantPiolaTransform

◆ ~VolumeElementForcesAndSourcesCore()

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

Definition at line 76 of file VolumeElementForcesAndSourcesCore.hpp.

76 {}

Member Function Documentation

◆ calculateBaseFunctionsOnElement() [1/2]

MoFEMErrorCode VolumeElementForcesAndSourcesCore::calculateBaseFunctionsOnElement ( const int  b)
virtual

Calculate base functions.

Returns
Error code

Definition at line 257 of file VolumeElementForcesAndSourcesCore.cpp.

258  {
260  if (dataH1.bAse.test(b)) {
261  switch (ApproximationBaseArray[b]) {
262  case NOBASE:
263  break;
266  if (dataH1.spacesOnEntities[MBTET].test(L2) &&
267  dataH1.basesOnEntities[MBTET].test(b) &&
268  dataH1.basesOnSpaces[L2].test(b)) {
269  CHKERR TetPolynomialBase().getValue(
270  gaussPts,
271  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
272  dataL2, L2, ApproximationBaseArray[b], NOBASE)));
273  }
274  if (dataH1.spacesOnEntities[MBVERTEX].test(H1) &&
275  dataH1.basesOnEntities[MBVERTEX].test(b) &&
276  dataH1.basesOnSpaces[H1].test(b)) {
277  CHKERR TetPolynomialBase().getValue(
278  gaussPts,
279  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
280  dataH1, H1, ApproximationBaseArray[b], NOBASE)));
281  }
283  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL) &&
284  dataH1.basesOnEntities[MBEDGE].test(b) &&
285  dataH1.basesOnSpaces[HCURL].test(b)) {
286  CHKERR TetPolynomialBase().getValue(
287  gaussPts,
288  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
289  dataHcurl, HCURL, ApproximationBaseArray[b], NOBASE)));
290  }
291  if (dataH1.spacesOnEntities[MBTRI].test(HDIV) &&
292  dataH1.basesOnEntities[MBTRI].test(b) &&
293  dataH1.basesOnSpaces[HDIV].test(b)) {
294  CHKERR TetPolynomialBase().getValue(
295  gaussPts,
296  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
297  dataHdiv, HDIV, ApproximationBaseArray[b], NOBASE)));
298  }
299  break;
300  default:
302  "Base <%s> not yet implemented",
303  ApproximationBaseNames[ApproximationBaseArray[b]]);
304  }
305  }
307 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:179
std::bitset< LASTBASE > basesOnEntities[MBMAXTYPE]
bases on entity types
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:498
static const char *const ApproximationBaseNames[]
Definition: definitions.h:153
std::bitset< LASTBASE > basesOnSpaces[LASTSPACE]
base on spaces
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
field with continuous tangents
Definition: definitions.h:178
#define CHKERR
Inline error check.
Definition: definitions.h:617
Like AINSWORTH_LEGENDRE_BASE but with Lobatto base instead Legendre .
Definition: definitions.h:143
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
virtual MPI_Comm & get_comm() const =0
continuous field
Definition: definitions.h:177
std::bitset< LASTBASE > bAse
bases on element
Construction of base is by Demkowicz .
Definition: definitions.h:145
field with C-1 continuity
Definition: definitions.h:180

◆ calculateBaseFunctionsOnElement() [2/2]

MoFEMErrorCode VolumeElementForcesAndSourcesCore::calculateBaseFunctionsOnElement ( )
virtual

Calculate base functions.

Returns
Error code

Definition at line 310 of file VolumeElementForcesAndSourcesCore.cpp.

310  {
312  /// Use the some node base
313  dataHdiv.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
314  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
315  dataHcurl.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
316  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
317  dataL2.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
318  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
319  std::vector<FieldApproximationBase> shape_functions_for_bases;
320  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
322  }
324 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:498
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
#define CHKERR
Inline error check.
Definition: definitions.h:617
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
virtual MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode VolumeElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
virtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 178 of file VolumeElementForcesAndSourcesCore.cpp.

178  {
180  // Get coords at Gauss points
182  double *shape_functions_ptr =
183  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
184  coordsAtGaussPts.resize(nbGaussPts, 3, false);
185  coordsAtGaussPts.clear();
186  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_coords_at_gauss_ptr(
187  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
188  &coordsAtGaussPts(0, 2));
189  FTensor::Tensor0<double *> t_shape_functions(shape_functions_ptr);
190  for (unsigned int gg = 0; gg < nbGaussPts; gg++) {
192  &coords[0], &coords[1], &coords[2]);
193  for (int bb = 0; bb < 4; bb++) {
194  t_coords_at_gauss_ptr(i) += t_coords(i) * t_shape_functions;
195  ++t_coords;
196  ++t_shape_functions;
197  };
198  ++t_coords_at_gauss_ptr;
199  }
200 
202  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:498
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:443

◆ calculateHoJacobian()

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

343  {
345  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
346  dataPtr->get<FieldName_mi_tag>().end()) {
347  const Field *field_struture =
349  BitFieldId id = field_struture->getId();
350  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none()) {
351  SETERRQ(mField.get_comm(), MOFEM_NOT_FOUND,
352  "no MESH_NODE_POSITIONS in element data");
353  }
358  if (dataH1.dataOnEntities[MBVERTEX][0].getFieldData().size() != 12) {
359  SETERRQ(mField.get_comm(), MOFEM_NOT_FOUND,
360  "no MESH_NODE_POSITIONS in element data or field has wrong "
361  "number of coefficients");
362  }
367  hoGaussPtsInvJac.resize(hoGaussPtsJac.size1(), hoGaussPtsJac.size2(),
368  false);
369  // Express Jacobian as tensor
371  &hoGaussPtsJac(0, 0), &hoGaussPtsJac(0, 1), &hoGaussPtsJac(0, 2),
372  &hoGaussPtsJac(0, 3), &hoGaussPtsJac(0, 4), &hoGaussPtsJac(0, 5),
373  &hoGaussPtsJac(0, 6), &hoGaussPtsJac(0, 7), &hoGaussPtsJac(0, 8));
375  &hoGaussPtsInvJac(0, 0), &hoGaussPtsInvJac(0, 1),
376  &hoGaussPtsInvJac(0, 2), &hoGaussPtsInvJac(0, 3),
377  &hoGaussPtsInvJac(0, 4), &hoGaussPtsInvJac(0, 5),
378  &hoGaussPtsInvJac(0, 6), &hoGaussPtsInvJac(0, 7),
379  &hoGaussPtsInvJac(0, 8));
380  hoGaussPtsDetJac.resize(nbGaussPts, false);
382  // Calculate inverse and determinant
383  for (unsigned int gg = 0; gg != nbGaussPts; gg++) {
384  CHKERR determinantTensor3by3(jac, det);
385  // if(det<0) {
386  // SETERRQ(mField.get_comm(),MOFEM_DATA_INCONSISTENCY,"Negative
387  // volume");
388  // }
389  CHKERR invertTensor3by3(jac, det, inv_jac);
390  ++jac;
391  ++inv_jac;
392  ++det;
393  }
394  } else {
395  hoCoordsAtGaussPts.resize(0, 0, false);
396  hoGaussPtsInvJac.resize(0, 0, false);
397  hoGaussPtsDetJac.resize(0, false);
398  }
400 }
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.
MoFEMErrorCode getTrisDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
higher order geometry data at Gauss pts
MoFEMErrorCode getEdgesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
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:498
std::bitset< BITFIELDID_SIZE > BitFieldId
Definition: Common.hpp:149
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
unsigned int nbGaussPts
Number of integration points.
#define CHKERR
Inline error check.
Definition: definitions.h:617
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 getTrisFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getTetsFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
virtual MPI_Comm & get_comm() const =0
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getEdgesDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getTetDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const

◆ calculateVolumeAndJacobian()

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

150  {
153  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
155  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
157  double diff_n[12];
158  CHKERR ShapeDiffMBTET(diff_n);
160  &diff_n[0], &diff_n[1], &diff_n[2]);
162  &coords[0], &coords[1], &coords[2]);
165  jAc.clear();
166  for (int nn = 0; nn != 4; nn++) {
167  tJac(i, j) += t_coords(i) * t_diff_n(j);
168  ++t_coords;
169  ++t_diff_n;
170  }
173  vOlume *= G_TET_W1[0] / 6.;
175 }
#define CHKERRQ_MOAB(a)
check error code of MoAB function
Definition: definitions.h:536
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:498
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:303
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Common.hpp:78
static const double G_TET_W1[]
Definition: fem_tools.h:495
#define CHKERR
Inline error check.
Definition: definitions.h:617
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode VolumeElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
virtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 205 of file VolumeElementForcesAndSourcesCore.cpp.

205  {
209  // H1
210  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1)) {
213  }
214  if ((dataH1.spacesOnEntities[MBTRI]).test(H1)) {
217  }
218  if ((dataH1.spacesOnEntities[MBTET]).test(H1)) {
220  }
221  // Hcurl
222  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL)) {
225  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
226  }
227  if ((dataH1.spacesOnEntities[MBTRI]).test(HCURL)) {
231  dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
232  }
233  if ((dataH1.spacesOnEntities[MBTET]).test(HCURL)) {
235  dataHcurl.spacesOnEntities[MBTET].set(HCURL);
236  }
237  // Hdiv
238  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV)) {
242  dataHdiv.spacesOnEntities[MBTRI].set(HDIV);
243  }
244  if ((dataH1.spacesOnEntities[MBTET]).test(HDIV)) {
246  dataHdiv.spacesOnEntities[MBTET].set(HDIV);
247  }
248  // L2
249  if ((dataH1.spacesOnEntities[MBTET]).test(L2)) {
251  dataL2.spacesOnEntities[MBTET].set(L2);
252  }
254 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
MoFEMErrorCode getTrisSense(DataForcesAndSourcesCore &data) const
field with continuous normal traction
Definition: definitions.h:179
MoFEMErrorCode getTetDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
MoFEMErrorCode getEdgesDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
MoFEMErrorCode getTrisDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:498
MoFEMErrorCode getFaceTriNodes(DataForcesAndSourcesCore &data) const
Get nodes on triangles.
MoFEMErrorCode getEdgesSense(DataForcesAndSourcesCore &data) const
field with continuous tangents
Definition: definitions.h:178
#define CHKERR
Inline error check.
Definition: definitions.h:617
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
continuous field
Definition: definitions.h:177
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
field with C-1 continuity
Definition: definitions.h:180

◆ operator()()

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

Examples:
Remodeling.hpp.

Definition at line 424 of file VolumeElementForcesAndSourcesCore.cpp.

424  {
426 
427  if (numeredEntFiniteElementPtr->getEntType() != MBTET)
429 
433  if (nbGaussPts == 0)
438 
439  try {
440  MatrixDouble new_diff_n;
441  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
443  FieldApproximationBase base = ApproximationBaseArray[b];
445  dataH1.dataOnEntities[MBVERTEX][0];
446  if ((data.getDiffN(base).size1() == 4) &&
447  (data.getDiffN(base).size2() == 3)) {
448  const unsigned int nb_base_functions = 4;
449  new_diff_n.resize(nbGaussPts, 3 * nb_base_functions, false);
450  double *new_diff_n_ptr = &*new_diff_n.data().begin();
452  new_diff_n_ptr, &new_diff_n_ptr[1], &new_diff_n_ptr[2]);
453  double *t_diff_n_ptr = &*data.getDiffN(base).data().begin();
454  for (unsigned int gg = 0; gg != nbGaussPts; gg++) {
456  t_diff_n_ptr, &t_diff_n_ptr[1], &t_diff_n_ptr[2]);
457  for (unsigned int bb = 0; bb != nb_base_functions; bb++) {
458  t_new_diff_n(i) = t_diff_n(i);
459  ++t_new_diff_n;
460  ++t_diff_n;
461  }
462  }
463  data.getDiffN(base).resize(new_diff_n.size1(), new_diff_n.size2(),
464  false);
465  data.getDiffN(base).data().swap(new_diff_n.data());
466  }
467  }
468  } catch (std::exception &ex) {
469  std::ostringstream ss;
470  ss << "thorw in method: " << ex.what() << " at line " << __LINE__
471  << " in file " << __FILE__;
472  SETERRQ(mField.get_comm(), MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
473  }
474 
477 
480  std::vector<std::string> last_eval_field_name(2);
481  DataForcesAndSourcesCore *op_data[2];
482  FieldSpace space[2];
483  FieldApproximationBase base[2];
484 
485  boost::ptr_vector<UserDataOperator>::iterator oit, hi_oit;
486  oit = opPtrVector.begin();
487  hi_oit = opPtrVector.end();
488 
489  for (; oit != hi_oit; oit++) {
490 
491  oit->setPtrFE(this);
492 
493  if (oit->opType == UserDataOperator::OPLAST) {
494 
495  // Set field
496  switch (oit->sPace) {
497  case NOSPACE:
498  SETERRQ(mField.get_comm(), MOFEM_DATA_INCONSISTENCY, "unknown space");
499  case H1 :
500  op_data[0] = &dataH1;
501  break;
502  case HCURL:
503  op_data[0] = &dataHcurl;
504  break;
505  case HDIV:
506  op_data[0] = &dataHdiv;
507  break;
508  case L2:
509  op_data[0] = &dataL2;
510  break;
511  case NOFIELD:
512  op_data[0] = &dataNoField;
513  break;
514  case LASTSPACE:
515  SETERRQ(mField.get_comm(), MOFEM_DATA_INCONSISTENCY, "unknown space");
516  break;
517  }
518 
519  // Reseat all data which all field dependent
520  op_data[0]->resetFieldDependentData();
521  last_eval_field_name[0] = "";
522  last_eval_field_name[1] = "";
523 
524  // Run operator
525  CHKERR oit->opRhs(*op_data[0], oit->doVertices, oit->doEdges,
526  oit->doQuads, oit->doTris, oit->doTets, false);
527 
528  } else {
529 
530  for (int ss = 0; ss != 2; ss++) {
531 
532  std::string field_name = !ss ? oit->rowFieldName : oit->colFieldName;
533  const Field *field_struture = mField.get_field_structure(field_name);
534  BitFieldId data_id = field_struture->getId();
535 
536  if ((oit->getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
537  data_id)
538  .none()) {
540  "no data field < %s > on finite element < %s >",
541  field_name.c_str(), feName.c_str());
542  }
543 
544  if (oit->getOpType() & types[ss] ||
545  oit->getOpType() & UserDataOperator::OPROWCOL) {
546 
547  space[ss] = field_struture->getSpace();
548  switch (space[ss]) {
549  case NOSPACE:
551  "unknown space");
552  break;
553  case H1:
554  op_data[ss] = !ss ? &dataH1 : &derivedDataH1;
555  break;
556  case HCURL:
557  op_data[ss] = !ss ? &dataHcurl : &derivedDataHcurl;
558  break;
559  case HDIV:
560  op_data[ss] = !ss ? &dataHdiv : &derivedDataHdiv;
561  break;
562  case L2:
563  op_data[ss] = !ss ? &dataL2 : &derivedDataL2;
564  break;
565  case NOFIELD:
566  op_data[ss] = !ss ? &dataNoField : &dataNoFieldCol;
567  break;
568  case LASTSPACE:
570  "unknown space");
571  break;
572  }
573 
574  base[ss] = field_struture->getApproxBase();
575  switch (base[ss]) {
576  case NOBASE:
580  break;
581  default:
583  "unknown or not implemented base");
584  break;
585  }
586 
587  if (last_eval_field_name[ss] != field_name) {
588 
589  switch (space[ss]) {
590  case NOSPACE:
592  "unknown space");
593  break;
594  case H1:
595  if (!ss) {
596  CHKERR getRowNodesIndices(*op_data[ss], field_name);
597  } else {
598  CHKERR getColNodesIndices(*op_data[ss], field_name);
599  }
600  CHKERR getNodesFieldData(*op_data[ss], field_name);
601  case HCURL:
602  if (!ss) {
603  CHKERR getEdgesRowIndices(*op_data[ss], field_name);
604  } else {
605  CHKERR getEdgesColIndices(*op_data[ss], field_name);
606  }
607  CHKERR getEdgesDataOrderSpaceAndBase(*op_data[ss], field_name);
608  CHKERR getEdgesFieldData(*op_data[ss], field_name);
609  case HDIV:
610  if (!ss) {
611  CHKERR getTrisRowIndices(*op_data[ss], field_name);
612  } else {
613  CHKERR getTrisColIndices(*op_data[ss], field_name);
614  }
615  CHKERR getTrisDataOrderSpaceAndBase(*op_data[ss], field_name);
616  CHKERR getTrisFieldData(*op_data[ss], field_name);
617  case L2:
618  if (!ss) {
619  CHKERR getTetsRowIndices(*op_data[ss], field_name);
620  } else {
621  CHKERR getTetsColIndices(*op_data[ss], field_name);
622  }
623  CHKERR getTetDataOrderSpaceAndBase(*op_data[ss], field_name);
624  CHKERR getTetsFieldData(*op_data[ss], field_name);
625  break;
626  case NOFIELD:
627  if (!getNinTheLoop()) {
628  // NOFIELD data are the same for each element, can be
629  // retrieved only once
630  if (!ss) {
631  CHKERR getNoFieldRowIndices(*op_data[ss], field_name);
632  } else {
633  CHKERR getNoFieldColIndices(*op_data[ss], field_name);
634  }
635  CHKERR getNoFieldFieldData(*op_data[ss], field_name);
636  }
637  break;
638  case LASTSPACE:
640  "unknown space");
641  break;
642  }
643  last_eval_field_name[ss] = field_name;
644  }
645  }
646  }
647 
648  if (oit->getOpType() & UserDataOperator::OPROW) {
649  try {
650  ierr = oit->opRhs(*op_data[0], false);
651  CHKERRG(ierr);
652  } catch (std::exception &ex) {
653  std::ostringstream ss;
654  ss << "Operator "
655  << boost::typeindex::type_id_runtime(*oit).pretty_name()
656  << " operator number "
657  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
658  opPtrVector.begin(), oit)
659  << " thorw in method: " << ex.what() << " at line " << __LINE__
660  << " in file " << __FILE__;
662  ss.str().c_str());
663  }
664  }
665 
666  if (oit->getOpType() & UserDataOperator::OPCOL) {
667  try {
668  ierr = oit->opRhs(*op_data[1], false);
669  CHKERRG(ierr);
670  } catch (std::exception &ex) {
671  std::ostringstream ss;
672  ss << "Operator "
673  << boost::typeindex::type_id_runtime(*oit).pretty_name()
674  << " operator number "
675  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
676  opPtrVector.begin(), oit)
677  << " thorw in method: " << ex.what() << " at line " << __LINE__
678  << " in file " << __FILE__;
680  ss.str().c_str());
681  }
682  }
683 
684  if (oit->getOpType() & UserDataOperator::OPROWCOL) {
685  try {
686  ierr = oit->opLhs(*op_data[0], *op_data[1], oit->sYmm);
687  CHKERRG(ierr);
688  } catch (std::exception &ex) {
689  std::ostringstream ss;
690  ss << "Operator "
691  << boost::typeindex::type_id_runtime(*oit).pretty_name()
692  << " operator number "
693  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
694  opPtrVector.begin(), oit)
695  << " thorw in method: " << ex.what() << " at line " << __LINE__
696  << " in file " << __FILE__;
698  ss.str().c_str());
699  }
700  }
701  }
702  }
703 
705 }
boost::ptr_vector< UserDataOperator > opPtrVector
field with continuous normal traction
Definition: definitions.h:179
virtual MoFEMErrorCode transformHoBaseFunctions()
Transform base functions based on ho-geometry element Jacobian.
MoFEMErrorCode getNoFieldColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
MoFEMErrorCode getTrisDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getEdgesRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Edges row indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getColNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col node indices from FENumeredDofEntity_multiIndex
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
MoFEMErrorCode getTrisRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Tris row indices from FENumeredDofEntity_multiIndex
scalar or vector of scalars describe (no true field)
Definition: definitions.h:176
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:498
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
std::string feName
std::bitset< BITFIELDID_SIZE > BitFieldId
Definition: Common.hpp:149
int getNinTheLoop() const
get number of evaluated element in the loop
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
MoFEMErrorCode getNoFieldRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
OpType
Controls loop over entities on element.
MoFEMErrorCode getNoFieldFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getTetsRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Tets row indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getTetsColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Tets col indices from FENumeredDofEntity_multiIndex
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
unsigned int nbGaussPts
Number of integration points.
FieldApproximationBase
approximation base
Definition: definitions.h:140
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Common.hpp:212
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
field with continuous tangents
Definition: definitions.h:178
virtual MoFEMErrorCode transformBaseFunctions()
Transform base functions based on geometric element Jacobian.
DataForcesAndSourcesCore::EntData EntData
FieldSpace
approximation spaces
Definition: definitions.h:174
#define CHKERR
Inline error check.
Definition: definitions.h:617
MoFEMErrorCode getTrisFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
Like AINSWORTH_LEGENDRE_BASE but with Lobatto base instead Legendre .
Definition: definitions.h:143
MoFEMErrorCode getEdgesColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Edges col indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getTrisColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Tris col indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getTetsFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
virtual MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
virtual MPI_Comm & get_comm() const =0
continuous field
Definition: definitions.h:177
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getRowNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get row node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
Construction of base is by Demkowicz .
Definition: definitions.h:145
virtual MoFEMErrorCode calculateHoJacobian()
Calculate Jacobian for HO geometry.
virtual MoFEMErrorCode calculateVolumeAndJacobian()
Calculate element volume and Jacobian.
MoFEMErrorCode getTetDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
field with C-1 continuity
Definition: definitions.h:180

◆ setIntegartionPts()

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

Definition at line 317 of file VolumeElementForcesAndSourcesCore.hpp.

317  {
318  return setIntegrationPts();
319  }
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.

◆ setIntegrationPts()

MoFEMErrorCode VolumeElementForcesAndSourcesCore::setIntegrationPts ( )
virtual

Set integration points.

Returns
Error code

Definition at line 99 of file VolumeElementForcesAndSourcesCore.cpp.

99  {
101  int order_data = getMaxDataOrder();
102  int order_row = getMaxRowOrder();
103  int order_col = getMaxColOrder();
104  int rule = getRule(order_row, order_col, order_data);
105  // std::cerr << order_data << " " << order_row << " " << order_col << " " <<
106  // rule << std::endl;
107  if (rule >= 0) {
108  if (rule < QUAD_3D_TABLE_SIZE) {
109  if (QUAD_3D_TABLE[rule]->dim != 3) {
110  SETERRQ(mField.get_comm(), MOFEM_DATA_INCONSISTENCY, "wrong dimension");
111  }
112  if (QUAD_3D_TABLE[rule]->order < rule) {
114  "wrong order %d != %d", QUAD_3D_TABLE[rule]->order, rule);
115  }
117  gaussPts.resize(4, nbGaussPts, false);
118  cblas_dcopy(nbGaussPts, &QUAD_3D_TABLE[rule]->points[1], 4,
119  &gaussPts(0, 0), 1);
120  cblas_dcopy(nbGaussPts, &QUAD_3D_TABLE[rule]->points[2], 4,
121  &gaussPts(1, 0), 1);
122  cblas_dcopy(nbGaussPts, &QUAD_3D_TABLE[rule]->points[3], 4,
123  &gaussPts(2, 0), 1);
124  cblas_dcopy(nbGaussPts, QUAD_3D_TABLE[rule]->weights, 1, &gaussPts(3, 0),
125  1);
126  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 4,
127  false);
128  double *shape_ptr =
129  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
130  cblas_dcopy(4 * nbGaussPts, QUAD_3D_TABLE[rule]->points, 1, shape_ptr, 1);
131  } else {
133  "rule > quadrature order %d < %d", rule, QUAD_3D_TABLE_SIZE);
134  nbGaussPts = 0;
135  }
136  } else {
137  CHKERR setGaussPts(order_row, order_col, order_data);
138  nbGaussPts = gaussPts.size2();
139  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 4,
140  false);
141  if (nbGaussPts > 0) {
143  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
144  &gaussPts(0, 0), &gaussPts(1, 0), &gaussPts(2, 0), nbGaussPts);
145  }
146  }
148 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
#define QUAD_3D_TABLE_SIZE
Definition: quad.h:186
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
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:498
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:617
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:443
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 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 326 of file VolumeElementForcesAndSourcesCore.cpp.

326  {
329  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
332  }
333  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
336  }
337  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
339  }
341  }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:179
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:498
field with continuous tangents
Definition: definitions.h:178
#define CHKERR
Inline error check.
Definition: definitions.h:617
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)
OpSetContravariantPiolaTransform opContravariantPiolaTransform
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
field with C-1 continuity
Definition: definitions.h:180

◆ transformHoBaseFunctions()

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

402  {
404 
405  if (hoCoordsAtGaussPts.size1() > 0) {
406  // Transform derivatives of base functions and apply Piola transformation
407  // if needed.
409  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
411  }
412  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
415  }
416  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
419  }
420  }
422 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:179
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:498
OpSetHoContravariantPiolaTransform opHoContravariantTransform
field with continuous tangents
Definition: definitions.h:178
#define CHKERR
Inline error check.
Definition: definitions.h:617
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)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
field with C-1 continuity
Definition: definitions.h:180

Member Data Documentation

◆ conn

const EntityHandle* MoFEM::VolumeElementForcesAndSourcesCore::conn

Definition at line 81 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::VolumeElementForcesAndSourcesCore::coords

Definition at line 41 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::coordsAtGaussPts

Definition at line 86 of file VolumeElementForcesAndSourcesCore.hpp.

◆ dataH1

DataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::dataH1

Definition at line 45 of file VolumeElementForcesAndSourcesCore.hpp.

◆ dataHcurl

DataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::dataHcurl

Definition at line 51 of file VolumeElementForcesAndSourcesCore.hpp.

◆ dataHdiv

DataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::dataHdiv

Definition at line 49 of file VolumeElementForcesAndSourcesCore.hpp.

◆ dataL2

DataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::dataL2

Definition at line 47 of file VolumeElementForcesAndSourcesCore.hpp.

◆ dataNoField

DataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::dataNoField

Definition at line 53 of file VolumeElementForcesAndSourcesCore.hpp.

◆ dataNoFieldCol

DataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::dataNoFieldCol

Definition at line 54 of file VolumeElementForcesAndSourcesCore.hpp.

◆ derivedDataH1

DerivedDataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::derivedDataH1

Definition at line 46 of file VolumeElementForcesAndSourcesCore.hpp.

◆ derivedDataHcurl

DerivedDataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::derivedDataHcurl

Definition at line 52 of file VolumeElementForcesAndSourcesCore.hpp.

◆ derivedDataHdiv

DerivedDataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::derivedDataHdiv

Definition at line 50 of file VolumeElementForcesAndSourcesCore.hpp.

◆ derivedDataL2

DerivedDataForcesAndSourcesCore MoFEM::VolumeElementForcesAndSourcesCore::derivedDataL2

Definition at line 48 of file VolumeElementForcesAndSourcesCore.hpp.

◆ gaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::gaussPts

Definition at line 85 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::hoCoordsAtGaussPts

Definition at line 62 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsDetJac

VectorDouble MoFEM::VolumeElementForcesAndSourcesCore::hoGaussPtsDetJac

Definition at line 65 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsInvJac

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::hoGaussPtsInvJac

Definition at line 64 of file VolumeElementForcesAndSourcesCore.hpp.

◆ hoGaussPtsJac

MatrixDouble MoFEM::VolumeElementForcesAndSourcesCore::hoGaussPtsJac

Definition at line 63 of file VolumeElementForcesAndSourcesCore.hpp.

◆ invJac

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCore::invJac

Definition at line 43 of file VolumeElementForcesAndSourcesCore.hpp.

◆ jAc

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCore::jAc

Definition at line 42 of file VolumeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::VolumeElementForcesAndSourcesCore::meshPositionsFieldName

Definition at line 61 of file VolumeElementForcesAndSourcesCore.hpp.

◆ nbGaussPts

unsigned int MoFEM::VolumeElementForcesAndSourcesCore::nbGaussPts

Number of integration points.

Definition at line 305 of file VolumeElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::VolumeElementForcesAndSourcesCore::num_nodes

Definition at line 80 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opContravariantPiolaTransform

OpSetContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opContravariantPiolaTransform

Definition at line 57 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opCovariantPiolaTransform

OpSetCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opCovariantPiolaTransform

Definition at line 58 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHOatGaussPoints

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

higher order geometry data at Gauss pts

Definition at line 68 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHoContravariantTransform

OpSetHoContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opHoContravariantTransform

Definition at line 70 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opHoCovariantTransform

OpSetHoCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opHoCovariantTransform

Definition at line 71 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetHoInvJacH1

OpSetHoInvJacH1 MoFEM::VolumeElementForcesAndSourcesCore::opSetHoInvJacH1

Definition at line 69 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetHoInvJacHdivAndHcurl

OpSetHoInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCore::opSetHoInvJacHdivAndHcurl

Definition at line 72 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacH1

OpSetInvJacH1 MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacH1

Definition at line 56 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacHdivAndHcurl

OpSetInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacHdivAndHcurl

Definition at line 59 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tInvJac

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

Definition at line 83 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tJac

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

Definition at line 82 of file VolumeElementForcesAndSourcesCore.hpp.

◆ vOlume

double MoFEM::VolumeElementForcesAndSourcesCore::vOlume

Definition at line 78 of file VolumeElementForcesAndSourcesCore.hpp.


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