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

Face finite elementUser is implementing own operator at Gauss point level, by own object derived from FaceElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to OpPtrVector. More...

#include <src/finite_elements/FaceElementForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 default operator for TRI element More...
 

Public Member Functions

 FaceElementForcesAndSourcesCore (Interface &m_field)
 
virtual MoFEMErrorCode calculateAreaAndNormal ()
 Calculate element area and normal of the face. More...
 
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
DEPRECATED MoFEMErrorCode setIntegartionPts ()
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate base functions. More...
 
virtual MoFEMErrorCode calculateHoNormal ()
 Calculate normal on curved elements. 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

double aRea
 
int num_nodes
 
const EntityHandleconn
 
VectorDouble nOrmal
 
VectorDouble tangentOne
 
VectorDouble tangentTwo
 
VectorDouble coords
 
MatrixDouble gaussPts
 
MatrixDouble coordsAtGaussPts
 
DataForcesAndSourcesCore dataH1
 
DerivedDataForcesAndSourcesCore derivedDataH1
 
DataForcesAndSourcesCore dataHdiv
 
DerivedDataForcesAndSourcesCore derivedDataHdiv
 
DataForcesAndSourcesCore dataHcurl
 
DerivedDataForcesAndSourcesCore derivedDataHcurl
 
DataForcesAndSourcesCore dataL2
 
DerivedDataForcesAndSourcesCore derivedDataL2
 
DataForcesAndSourcesCore dataNoField
 
DataForcesAndSourcesCore dataNoFieldCol
 
std::string meshPositionsFieldName
 Name of the field with geometry. More...
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble normalsAtGaussPt
 
MatrixDouble tangentOneAtGaussPt
 
MatrixDouble tangentTwoAtGaussPt
 
OpGetCoordsAndNormalsOnFace opHOCoordsAndNormals
 
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform
 
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
 
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

Face finite element

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

Examples:
hello_world.cpp, MagneticElement.hpp, simple_elasticity.cpp, and simple_interface.cpp.

Definition at line 39 of file FaceElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FaceElementForcesAndSourcesCore()

MoFEM::FaceElementForcesAndSourcesCore::FaceElementForcesAndSourcesCore ( Interface m_field)

Definition at line 70 of file FaceElementForcesAndSourcesCore.hpp.

71  : ForcesAndSourcesCore(m_field), dataH1(MBTRI), derivedDataH1(dataH1),
72  dataHdiv(MBTRI), derivedDataHdiv(dataHdiv), dataHcurl(MBTRI),
74  dataNoField(MBTRI), dataNoFieldCol(MBTRI),
75  meshPositionsFieldName("MESH_NODE_POSITIONS"),
ForcesAndSourcesCore(Interface &m_field)
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
std::string meshPositionsFieldName
Name of the field with geometry.
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform

Member Function Documentation

◆ calculateAreaAndNormal()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateAreaAndNormal ( )
virtual

Calculate element area and normal of the face.

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

Returns
Error code

Definition at line 161 of file FaceElementForcesAndSourcesCore.cpp.

161  {
164  rval = mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
166  coords.resize(num_nodes * 3, false);
167  rval = mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
169  double diff_n[6];
170  ierr = ShapeDiffMBTRI(diff_n);
171  CHKERRG(ierr);
172  nOrmal.resize(3, false);
173  ierr = ShapeFaceNormalMBTRI(diff_n, &*coords.data().begin(),
174  &*nOrmal.data().begin());
175  CHKERRG(ierr);
176  aRea = cblas_dnrm2(3, &*nOrmal.data().begin(), 1) * 0.5;
177  tangentOne.resize(3, false);
178  tangentTwo.resize(3, false);
179  for (int dd = 0; dd != 3; dd++) {
180  tangentOne[dd] = cblas_ddot(3, &diff_n[0], 2, &coords[dd], 3);
181  tangentTwo[dd] = cblas_ddot(3, &diff_n[1], 2, &coords[dd], 3);
182  }
184 }
#define CHKERRQ_MOAB(a)
check error code of MoAB function
Definition: definitions.h:536
virtual moab::Interface & get_moab()=0
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T *> &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
PetscErrorCode ShapeFaceNormalMBTRI(double *diffN, const double *coords, double *normal)
Definition: fem_tools.c:219
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
double cblas_dnrm2(const int N, const double *X, const int incX)
Definition: cblas_dnrm2.c:12
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:186
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Common.hpp:78

◆ calculateBaseFunctionsOnElement()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateBaseFunctionsOnElement ( )
virtual

Calculate base functions.

Returns
Error code

Definition at line 299 of file FaceElementForcesAndSourcesCore.cpp.

299  {
301  // Calculate base base functions for faces.
302 
303  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
304  if (dataH1.bAse.test(b)) {
305  switch (ApproximationBaseArray[b]) {
306  case NOBASE:
307  break;
310  if (dataH1.spacesOnEntities[MBVERTEX].test(H1) &&
311  dataH1.basesOnEntities[MBVERTEX].test(b) &&
312  dataH1.basesOnSpaces[H1].test(b)) {
313  CHKERR TriPolynomialBase().getValue(
314  gaussPts,
315  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
316  dataH1, H1, ApproximationBaseArray[b], NOBASE)));
317  }
318  if (dataH1.spacesOnEntities[MBTRI].test(L2) &&
319  dataH1.basesOnEntities[MBTRI].test(b) &&
320  dataH1.basesOnSpaces[L2].test(b)) {
321  CHKERR TriPolynomialBase().getValue(
322  gaussPts,
323  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
324  dataL2, L2, ApproximationBaseArray[b], NOBASE)));
325  }
327  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL) &&
328  dataH1.basesOnEntities[MBEDGE].test(b) &&
329  dataH1.basesOnSpaces[HCURL].test(b)) {
330  CHKERR TriPolynomialBase().getValue(
331  gaussPts,
332  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
333  dataHcurl, HCURL, ApproximationBaseArray[b], NOBASE)));
334  }
335  if (dataH1.spacesOnEntities[MBTRI].test(HDIV) &&
336  dataH1.basesOnEntities[MBTRI].test(b) &&
337  dataH1.basesOnSpaces[HDIV].test(b)) {
338  CHKERR TriPolynomialBase().getValue(
339  gaussPts,
340  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
341  dataHdiv, HDIV, ApproximationBaseArray[b], NOBASE)));
342  }
343  break;
344  default:
345  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
346  "Base <%s> not yet implemented",
347  ApproximationBaseNames[ApproximationBaseArray[b]]);
348  }
349  }
350  }
352 }
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
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

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
virtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 284 of file FaceElementForcesAndSourcesCore.cpp.

284  {
286  double *shape_functions =
287  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
288  coordsAtGaussPts.resize(nbGaussPts, 3, false);
289  for (int gg = 0; gg < nbGaussPts; gg++) {
290  for (int dd = 0; dd < 3; dd++) {
291  coordsAtGaussPts(gg, dd) =
292  cblas_ddot(3, &shape_functions[3 * gg], 1, &coords[dd], 3);
293  }
294  }
296 }
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T *> &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12

◆ calculateHoNormal()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateHoNormal ( )
virtual

Calculate normal on curved elements.

Geometry is given by other field.

Returns
error code

Definition at line 354 of file FaceElementForcesAndSourcesCore.cpp.

354  {
356  // Check if field for high-order geometry is set and if it is set calculate
357  // higher-order normals and face tangent vectors.
358  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
359  dataPtr->get<FieldName_mi_tag>().end()) {
360  const Field *field_struture =
362  BitFieldId id = field_struture->getId();
363 
364  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none()) {
365  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
366  "no MESH_NODE_POSITIONS in element data");
367  }
368 
369  // Calculate normal for high-order geometry
371  CHKERRG(ierr);
373  CHKERRG(ierr);
375  CHKERRG(ierr);
377  CHKERRG(ierr);
379  CHKERRG(ierr);
380  try {
382  CHKERRG(ierr);
384  CHKERRG(ierr);
385  } catch (std::exception &ex) {
386  std::ostringstream ss;
387  ss << "thorw in method: " << ex.what() << " at line " << __LINE__
388  << " in file " << __FILE__;
389  SETERRQ(PETSC_COMM_SELF, 1, ss.str().c_str());
390  }
391  } else {
392  hoCoordsAtGaussPts.resize(0, 0, false);
393  normalsAtGaussPt.resize(0, 0, false);
394  tangentOneAtGaussPt.resize(0, 0, false);
395  tangentTwoAtGaussPt.resize(0, 0, false);
396  }
398 }
MoFEMErrorCode getTrisDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getEdgesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
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 CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
std::bitset< BITFIELDID_SIZE > BitFieldId
Definition: Common.hpp:149
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
std::string meshPositionsFieldName
Name of the field with geometry.
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
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getEdgesDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode FaceElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
virtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 238 of file FaceElementForcesAndSourcesCore.cpp.

238  {
240  // Get spaces order/base and sense of entities.
241 
243 
244  // H1
245  if (dataH1.spacesOnEntities[MBEDGE].test(H1)) {
248  }
249  if (dataH1.spacesOnEntities[MBTRI].test(H1)) {
252  }
253 
254  // Hcurl
255  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
258  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
259  }
260  if (dataH1.spacesOnEntities[MBTRI].test(HCURL)) {
263  dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
264  }
265 
266  // Hdiv
267  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
270  dataHcurl.spacesOnEntities[MBTRI].set(HDIV);
271  }
272 
273  // L2
274  if (dataH1.spacesOnEntities[MBTRI].test(L2)) {
277  dataHcurl.spacesOnEntities[MBTRI].set(L2);
278  }
279 
281 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
MoFEMErrorCode getTrisSense(DataForcesAndSourcesCore &data) const
field with continuous normal traction
Definition: definitions.h:179
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 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 FaceElementForcesAndSourcesCore::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.

Reimplemented in NeummanForcesSurfaceComplexForLazy::MyTriangleSpatialFE.

Definition at line 400 of file FaceElementForcesAndSourcesCore.cpp.

400  {
402 
403  if (numeredEntFiniteElementPtr->getEntType() != MBTRI)
405 
406  // Calculate normal and tangent vectors for face geometry given by 3 nodes.
409 
411  if (nbGaussPts == 0)
413 
414  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
416  &*dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
417 
418  /// Use the some node base
420 
421  // Share base shape functions between spaces
422  dataHdiv.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
423  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
424  dataHcurl.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
425  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
426  dataL2.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
427  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
428  dataHdiv.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(NOBASE) =
429  dataH1.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(NOBASE);
430  dataHcurl.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(NOBASE) =
431  dataH1.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(NOBASE);
432  dataL2.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(NOBASE) =
433  dataH1.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(NOBASE);
434 
437 
438  // Apply Piola transform to HDiv and HCurl spaces, uses previously calculated
439  // faces normal and tangent vectors.
440  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
442  }
443  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
444  // cerr << dataHcurl.dataOnEntities[MBEDGE][0].getN(AINSWORTH_LEGENDRE_BASE)
445  // << endl;
447  }
448 
451  std::vector<std::string> last_eval_field_name(2);
452  DataForcesAndSourcesCore *op_data[2];
453  FieldSpace space[2];
454  FieldApproximationBase base[2];
455 
456  boost::ptr_vector<UserDataOperator>::iterator oit, hi_oit;
457 
458  oit = opPtrVector.begin();
459  hi_oit = opPtrVector.end();
460 
461  for (; oit != hi_oit; oit++) {
462 
463  oit->setPtrFE(this);
464 
465  if (oit->opType == UserDataOperator::OPLAST) {
466 
467  // Set field
468  switch (oit->sPace) {
469  case NOSPACE:
470  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
471  case H1:
472  op_data[0] = &dataH1;
473  break;
474  case HCURL:
475  op_data[0] = &dataHcurl;
476  break;
477  case HDIV:
478  op_data[0] = &dataHdiv;
479  break;
480  case L2:
481  op_data[0] = &dataL2;
482  break;
483  case NOFIELD:
484  op_data[0] = &dataNoField;
485  break;
486  case LASTSPACE:
487  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
488  }
489 
490  // Reseat all data which all field dependent
491  op_data[0]->resetFieldDependentData();
492  last_eval_field_name[0] = "";
493  last_eval_field_name[1] = "";
494 
495  // Run operator
496  CHKERR oit->opRhs(*op_data[0], oit->doVertices, oit->doEdges,
497  oit->doQuads, oit->doTris, false, false);
498 
499  } else {
500 
501  for (int ss = 0; ss != 2; ss++) {
502 
503  std::string field_name = !ss ? oit->rowFieldName : oit->colFieldName;
504  const Field *field_struture = mField.get_field_structure(field_name);
505  BitFieldId data_id = field_struture->getId();
506 
507  if ((oit->getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
508  data_id)
509  .none()) {
510  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
511  "no data field < %s > on finite element < %s >",
512  field_name.c_str(), feName.c_str());
513  }
514 
515  if (oit->getOpType() & types[ss] ||
516  oit->getOpType() & UserDataOperator::OPROWCOL) {
517 
518  space[ss] = field_struture->getSpace();
519  switch (space[ss]) {
520  case NOSPACE:
521  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
522  case H1:
523  op_data[ss] = !ss ? &dataH1 : &derivedDataH1;
524  break;
525  case HCURL:
526  op_data[ss] = !ss ? &dataHcurl : &derivedDataHcurl;
527  break;
528  case HDIV:
529  op_data[ss] = !ss ? &dataHdiv : &derivedDataHdiv;
530  break;
531  case L2:
532  op_data[ss] = !ss ? &dataL2 : &derivedDataL2;
533  break;
534  case NOFIELD:
535  op_data[ss] = !ss ? &dataNoField : &dataNoFieldCol;
536  break;
537  case LASTSPACE:
538  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
539  }
540 
541  base[ss] = field_struture->getApproxBase();
542  switch (base[ss]) {
543  case NOBASE:
547  break;
548  default:
549  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
550  "unknown or not implemented base");
551  break;
552  }
553 
554  if (last_eval_field_name[ss] != field_name) {
555 
556  switch (space[ss]) {
557  case NOSPACE:
558  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
559  "unknown space");
560  case H1:
561  if (!ss) {
562  CHKERR getRowNodesIndices(*op_data[ss], field_name);
563  } else {
564  CHKERR getColNodesIndices(*op_data[ss], field_name);
565  }
566  CHKERR getNodesFieldData(*op_data[ss], field_name);
567  case HCURL:
568  if (!ss) {
569  CHKERR getEdgesRowIndices(*op_data[ss], field_name);
570  } else {
571  CHKERR getEdgesColIndices(*op_data[ss], field_name);
572  }
573  CHKERR getEdgesDataOrderSpaceAndBase(*op_data[ss], field_name);
574  CHKERR getEdgesFieldData(*op_data[ss], field_name);
575  case HDIV:
576  case L2:
577  if (!ss) {
578  CHKERR getTrisRowIndices(*op_data[ss], field_name);
579  } else {
580  CHKERR getTrisColIndices(*op_data[ss], field_name);
581  }
582  CHKERR getTrisDataOrderSpaceAndBase(*op_data[ss], field_name);
583  CHKERR getTrisFieldData(*op_data[ss], field_name);
584  break;
585  case NOFIELD:
586  if (!getNinTheLoop()) {
587  // NOFIELD data are the same for each element, can be retrieved
588  // only once
589  if (!ss) {
590  CHKERR getNoFieldRowIndices(*op_data[ss], field_name);
591  } else {
592  CHKERR getNoFieldColIndices(*op_data[ss], field_name);
593  }
594  CHKERR getNoFieldFieldData(*op_data[ss], field_name);
595  }
596  break;
597  case LASTSPACE:
598  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
599  "unknown space");
600  }
601  last_eval_field_name[ss] = field_name;
602  }
603  }
604  }
605 
606  if (oit->getOpType() & UserDataOperator::OPROW) {
607  try {
608  ierr = oit->opRhs(*op_data[0], oit->doVertices, oit->doEdges,
609  oit->doQuads, oit->doTris, false, false);
610  CHKERRG(ierr);
611  } catch (std::exception &ex) {
612  std::ostringstream ss;
613  ss << "Operator "
614  << boost::typeindex::type_id_runtime(*oit).pretty_name()
615  << " operator number "
616  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
617  opPtrVector.begin(), oit)
618  << " throw in method: " << ex.what() << " at line " << __LINE__
619  << " in file " << __FILE__;
620  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
621  }
622  }
623 
624  if (oit->getOpType() & UserDataOperator::OPCOL) {
625  try {
626  ierr = oit->opRhs(*op_data[1], oit->doVertices, oit->doEdges,
627  oit->doQuads, oit->doTris, false, false);
628  CHKERRG(ierr);
629  } catch (std::exception &ex) {
630  std::ostringstream ss;
631  ss << "Operator "
632  << boost::typeindex::type_id_runtime(*oit).pretty_name()
633  << " operator number "
634  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
635  opPtrVector.begin(), oit)
636  << " throw in method: " << ex.what() << " at line " << __LINE__
637  << " in file " << __FILE__;
638  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
639  }
640  }
641 
642  if (oit->getOpType() & UserDataOperator::OPROWCOL) {
643  try {
644  ierr = oit->opLhs(*op_data[0], *op_data[1], oit->sYmm);
645  CHKERRG(ierr);
646  } catch (std::exception &ex) {
647  std::ostringstream ss;
648  ss << "Operator "
649  << boost::typeindex::type_id_runtime(*oit).pretty_name()
650  << " operator number "
651  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
652  opPtrVector.begin(), oit)
653  << " throw in method: " << ex.what() << " at line " << __LINE__
654  << " in file " << __FILE__;
655  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
656  }
657  }
658  }
659  }
660 
662 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
boost::ptr_vector< UserDataOperator > opPtrVector
field with continuous normal traction
Definition: definitions.h:179
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 setIntegrationPts()
Set integration points.
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
virtual MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:186
FieldApproximationBase
approximation base
Definition: definitions.h:140
virtual MoFEMErrorCode calculateHoNormal()
Calculate normal on curved elements.
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
field with continuous tangents
Definition: definitions.h:178
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 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)
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
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.
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
field with C-1 continuity
Definition: definitions.h:180

◆ setIntegartionPts()

DEPRECATED MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCore::setIntegartionPts ( )
Deprecated:
method with spelling mistake, use setIntegrationPts

Definition at line 381 of file FaceElementForcesAndSourcesCore.hpp.

381 { return setIntegrationPts(); }
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.

◆ setIntegrationPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::setIntegrationPts ( )
virtual

Set integration points.

Returns
Error code

Definition at line 186 of file FaceElementForcesAndSourcesCore.cpp.

186  {
188  // Set integration points
189  int order_data = getMaxDataOrder();
190  int order_row = getMaxRowOrder();
191  int order_col = getMaxColOrder();
192  int rule = getRule(order_row, order_col, order_data);
193  if (rule >= 0) {
194  if (rule < QUAD_2D_TABLE_SIZE) {
195  if (QUAD_2D_TABLE[rule]->dim != 2) {
196  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
197  }
198  if (QUAD_2D_TABLE[rule]->order < rule) {
199  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
200  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
201  }
203  gaussPts.resize(3, nbGaussPts, false);
204  cblas_dcopy(nbGaussPts, &QUAD_2D_TABLE[rule]->points[1], 3,
205  &gaussPts(0, 0), 1);
206  cblas_dcopy(nbGaussPts, &QUAD_2D_TABLE[rule]->points[2], 3,
207  &gaussPts(1, 0), 1);
208  cblas_dcopy(nbGaussPts, QUAD_2D_TABLE[rule]->weights, 1, &gaussPts(2, 0),
209  1);
210  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 3,
211  false);
212  double *shape_ptr =
213  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
214  cblas_dcopy(3 * nbGaussPts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr, 1);
215  } else {
216  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
217  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
218  nbGaussPts = 0;
219  }
220  } else {
221  // If rule is negative, set user defined integration points
222  ierr = setGaussPts(order_row, order_col, order_data);
223  CHKERRG(ierr);
224  nbGaussPts = gaussPts.size2();
225  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 3,
226  false);
227  if (nbGaussPts) {
228  ierr = ShapeMBTRI(
229  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
230  &gaussPts(0, 0), &gaussPts(1, 0), nbGaussPts);
231  CHKERRG(ierr);
232  }
233  }
235 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
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 CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
virtual MoFEMErrorCode setGaussPts(int order)
It will be removed in the future use other variant.
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:175
virtual int getRule(int order)
set integration rule for finite element
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.

Member Data Documentation

◆ aRea

double MoFEM::FaceElementForcesAndSourcesCore::aRea

Definition at line 41 of file FaceElementForcesAndSourcesCore.hpp.

◆ conn

const EntityHandle* MoFEM::FaceElementForcesAndSourcesCore::conn

Definition at line 44 of file FaceElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::coords

Definition at line 46 of file FaceElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::coordsAtGaussPts

Definition at line 48 of file FaceElementForcesAndSourcesCore.hpp.

◆ dataH1

DataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::dataH1

Definition at line 50 of file FaceElementForcesAndSourcesCore.hpp.

◆ dataHcurl

DataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::dataHcurl

Definition at line 54 of file FaceElementForcesAndSourcesCore.hpp.

◆ dataHdiv

DataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::dataHdiv

Definition at line 52 of file FaceElementForcesAndSourcesCore.hpp.

◆ dataL2

DataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::dataL2

Definition at line 56 of file FaceElementForcesAndSourcesCore.hpp.

◆ dataNoField

DataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::dataNoField

Definition at line 58 of file FaceElementForcesAndSourcesCore.hpp.

◆ dataNoFieldCol

DataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::dataNoFieldCol

Definition at line 58 of file FaceElementForcesAndSourcesCore.hpp.

◆ derivedDataH1

DerivedDataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::derivedDataH1

Definition at line 51 of file FaceElementForcesAndSourcesCore.hpp.

◆ derivedDataHcurl

DerivedDataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::derivedDataHcurl

Definition at line 55 of file FaceElementForcesAndSourcesCore.hpp.

◆ derivedDataHdiv

DerivedDataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::derivedDataHdiv

Definition at line 53 of file FaceElementForcesAndSourcesCore.hpp.

◆ derivedDataL2

DerivedDataForcesAndSourcesCore MoFEM::FaceElementForcesAndSourcesCore::derivedDataL2

Definition at line 57 of file FaceElementForcesAndSourcesCore.hpp.

◆ gaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::gaussPts

Definition at line 47 of file FaceElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::hoCoordsAtGaussPts

Definition at line 62 of file FaceElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::FaceElementForcesAndSourcesCore::meshPositionsFieldName

Name of the field with geometry.

Definition at line 60 of file FaceElementForcesAndSourcesCore.hpp.

◆ nbGaussPts

int MoFEM::FaceElementForcesAndSourcesCore::nbGaussPts

Number of integration points.

Definition at line 372 of file FaceElementForcesAndSourcesCore.hpp.

◆ nOrmal

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::nOrmal

Definition at line 45 of file FaceElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPt

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPt

Definition at line 63 of file FaceElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::FaceElementForcesAndSourcesCore::num_nodes

Definition at line 42 of file FaceElementForcesAndSourcesCore.hpp.

◆ opContravariantTransform

OpSetContravariantPiolaTransformOnTriangle MoFEM::FaceElementForcesAndSourcesCore::opContravariantTransform

Definition at line 67 of file FaceElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnTriangle MoFEM::FaceElementForcesAndSourcesCore::opCovariantTransform

Definition at line 68 of file FaceElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnFace MoFEM::FaceElementForcesAndSourcesCore::opHOCoordsAndNormals

Definition at line 66 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOne

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOne

Definition at line 45 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOneAtGaussPt

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOneAtGaussPt

Definition at line 64 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwo

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwo

Definition at line 45 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwoAtGaussPt

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwoAtGaussPt

Definition at line 65 of file FaceElementForcesAndSourcesCore.hpp.


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