v0.9.0
Classes | Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
MoFEM::FaceElementForcesAndSourcesCoreBase 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::FaceElementForcesAndSourcesCoreBase:
[legend]
Collaboration diagram for MoFEM::FaceElementForcesAndSourcesCoreBase:
[legend]

Classes

struct  UserDataOperator
 default operator for TRI element More...
 

Public Types

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

Public Member Functions

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

Public Attributes

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

Protected Member Functions

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

Protected Attributes

double aRea
 
int num_nodes
 
const EntityHandleconn
 
VectorDouble nOrmal
 
VectorDouble tangentOne
 
VectorDouble tangentTwo
 
VectorDouble coords
 
MatrixDouble coordsAtGaussPts
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble normalsAtGaussPts
 
MatrixDouble tangentOneAtGaussPts
 
MatrixDouble tangentTwoAtGaussPts
 
OpGetCoordsAndNormalsOnFace opHOCoordsAndNormals
 
OpSetContravariantPiolaTransformOnFace opContravariantTransform
 
OpSetCovariantPiolaTransformOnFace opCovariantTransform
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 

Friends

class UserDataOperator
 
class VolumeElementForcesAndSourcesCoreOnSideBase
 

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.

Definition at line 37 of file FaceElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ Switches

Constructor & Destructor Documentation

◆ FaceElementForcesAndSourcesCoreBase()

MoFEM::FaceElementForcesAndSourcesCoreBase::FaceElementForcesAndSourcesCoreBase ( Interface m_field)

Definition at line 23 of file FaceElementForcesAndSourcesCore.cpp.

25  : ForcesAndSourcesCore(m_field),
26  meshPositionsFieldName("MESH_NODE_POSITIONS"),
33 }
ForcesAndSourcesCore(Interface &m_field)
std::string meshPositionsFieldName
Name of the field with geometry.
OpSetContravariantPiolaTransformOnFace opContravariantTransform

Member Function Documentation

◆ calculateAreaAndNormal()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateAreaAndNormal ( )
protectedvirtual

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 104 of file FaceElementForcesAndSourcesCore.cpp.

104  {
106 
108  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
109  coords.resize(num_nodes * 3, false);
110  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
111  nOrmal.resize(3, false);
112  tangentOne.resize(3, false);
113  tangentTwo.resize(3, false);
114 
115  auto calc_normal = [&](const double *diff_ptr) {
118  &coords[0], &coords[1], &coords[2]);
120  &nOrmal[0], &nOrmal[1], &nOrmal[2]);
122  &tangentOne[0], &tangentOne[1], &tangentOne[2]);
124  &tangentTwo[0], &tangentTwo[1], &tangentTwo[2]);
126  diff_ptr, &diff_ptr[1]);
127 
131 
134  t_t1(i) = 0;
135  t_t2(i) = 0;
136 
137  for (int nn = 0; nn != num_nodes; ++nn) {
138  t_t1(i) += t_coords(i) * t_diff(N0);
139  t_t2(i) += t_coords(i) * t_diff(N1);
140  ++t_coords;
141  ++t_diff;
142  }
143  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
144  aRea = sqrt(t_normal(i) * t_normal(i));
146  };
147 
148  const double *diff_ptr;
149  switch (numeredEntFiniteElementPtr->getEntType()) {
150  case MBTRI:
151  diff_ptr = Tools::diffShapeFunMBTRI.data();
152  CHKERR calc_normal(diff_ptr);
153  //FIXME: Normal should be divided not the area for triangle!!
154  aRea /= 2;
155  break;
156  case MBQUAD:
157  diff_ptr = Tools::diffShapeFunMBQUADAtCenter.data();
158  CHKERR calc_normal(diff_ptr);
159  break;
160  default:
161  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Element type not implemented");
162  }
163 
165 }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static constexpr std::array< double, 8 > diffShapeFunMBQUADAtCenter
Definition: Tools.hpp:112
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:87
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use

◆ calculateAreaAndNormalAtIntegrationPts()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateAreaAndNormalAtIntegrationPts ( )
protectedvirtual

Calculate element area and normal of the face at integration points.

Returns
Error code

Definition at line 36 of file FaceElementForcesAndSourcesCore.cpp.

36  {
38 
39  auto type = numeredEntFiniteElementPtr->getEntType();
40 
41  if (type == MBQUAD) {
42 
44  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
45  coords.resize(num_nodes * 3, false);
46  CHKERR mField.get_moab().get_coords(conn, num_nodes,
47  &*coords.data().begin());
48 
49  const size_t nb_gauss_pts = gaussPts.size2();
50  normalsAtGaussPts.resize(nb_gauss_pts, 3);
51  tangentOneAtGaussPts.resize(nb_gauss_pts, 3);
52  tangentTwoAtGaussPts.resize(nb_gauss_pts, 3);
53  normalsAtGaussPts.clear();
54  tangentOneAtGaussPts.clear();
55  tangentTwoAtGaussPts.clear();
56 
57  auto get_ftensor_from_mat_3d = [](MatrixDouble &m) {
59  &m(0, 0), &m(0, 1), &m(0, 2));
60  };
61  auto get_ftensor_from_mat_2d = [](MatrixDouble &m) {
63  &m(0, 1));
64  };
65  auto get_ftensor_from_vec_3d = [](VectorDouble &v) {
66  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(&v[0], &v[1],
67  &v[2]);
68  };
69 
70  auto t_t1 = get_ftensor_from_mat_3d(tangentOneAtGaussPts);
71  auto t_t2 = get_ftensor_from_mat_3d(tangentTwoAtGaussPts);
72  auto t_normal = get_ftensor_from_mat_3d(normalsAtGaussPts);
73 
74  auto t_diff = get_ftensor_from_mat_2d(
75  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE));
76 
80 
83 
84  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
85  auto t_coords = get_ftensor_from_vec_3d(coords);
86 
87  for (int nn = 0; nn != num_nodes; ++nn) {
88  t_t1(i) += t_coords(i) * t_diff(N0);
89  t_t2(i) += t_coords(i) * t_diff(N1);
90  ++t_diff;
91  ++t_coords;
92  }
93  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
94 
95  ++t_t1;
96  ++t_t2;
97  ++t_normal;
98  }
99  }
100 
102 }
MatrixDouble gaussPts
Matrix of integration points.
virtual moab::Interface & get_moab()=0
DataForcesAndSourcesCore & dataH1
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:74
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:72
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateCoordinatesAtGaussPts ( )
protectedvirtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 334 of file FaceElementForcesAndSourcesCore.cpp.

334  {
336 
337  const size_t nb_nodes =
338  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).size2();
339  double *shape_functions =
340  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
341  const size_t nb_gauss_pts = gaussPts.size2();
342  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
343  for (int gg = 0; gg != nb_gauss_pts; ++gg)
344  for (int dd = 0; dd != 3; ++dd)
346  nb_nodes, &shape_functions[nb_nodes * gg], 1, &coords[dd], 3);
347 
349 }
MatrixDouble gaussPts
Matrix of integration points.
DataForcesAndSourcesCore & dataH1
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12
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

◆ calculateHoNormal()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateHoNormal ( )
protectedvirtual

Calculate normal on curved elements.

Geometry is given by other field.

Returns
error code

Definition at line 351 of file FaceElementForcesAndSourcesCore.cpp.

351  {
353  // Check if field for high-order geometry is set and if it is set calculate
354  // higher-order normals and face tangent vectors.
355  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
356  dataPtr->get<FieldName_mi_tag>().end()) {
357 
358  const Field *field_struture =
360  BitFieldId id = field_struture->getId();
361 
362  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none())
363  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
364  "no MESH_NODE_POSITIONS in element data");
365 
366  // Calculate normal for high-order geometry
372 
373  } else if(numeredEntFiniteElementPtr->getEntType() == MBTRI) {
374  hoCoordsAtGaussPts.resize(0, 0, false);
375  normalsAtGaussPts.resize(0, 0, false);
376  tangentOneAtGaussPts.resize(0, 0, false);
377  tangentTwoAtGaussPts.resize(0, 0, false);
378  } else {
379  hoCoordsAtGaussPts.resize(coordsAtGaussPts.size1(),
380  coordsAtGaussPts.size2(), false);
382  }
383 
385 }
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
std::string meshPositionsFieldName
Name of the field with geometry.
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
Pointer to finite element data dofs.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getNodesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
Get data on nodes.
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:53
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual const Field * get_field_structure(const std::string &name)=0
get field structure

◆ getNumberOfNodes()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::getNumberOfNodes ( int num_nodes) const
protected

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::getSpaceBaseAndOrderOnElement ( )
protectedvirtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 284 of file FaceElementForcesAndSourcesCore.cpp.

284  {
286  // Get spaces order/base and sense of entities.
287 
289 
290  // H1
291  if (dataH1.spacesOnEntities[MBEDGE].test(H1)) {
292  CHKERR getEntitySense<MBEDGE>(dataH1);
293  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
294  }
295  if (dataH1.spacesOnEntities[MBTRI].test(H1)) {
296  CHKERR getEntitySense<MBTRI>(dataH1);
297  CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
298  }
299  if (dataH1.spacesOnEntities[MBQUAD].test(H1)) {
300  CHKERR getEntitySense<MBQUAD>(dataH1);
301  CHKERR getEntityDataOrder<MBQUAD>(dataH1, H1);
302  }
303 
304  // Hcurl
305  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
306  CHKERR getEntitySense<MBEDGE>(dataHcurl);
307  CHKERR getEntityDataOrder<MBEDGE>(dataHcurl, HCURL);
308  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
309  }
310  if (dataH1.spacesOnEntities[MBTRI].test(HCURL)) {
311  CHKERR getEntitySense<MBTRI>(dataHcurl);
312  CHKERR getEntityDataOrder<MBTRI>(dataHcurl, HCURL);
313  dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
314  }
315 
316  // Hdiv
317  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
318  CHKERR getEntitySense<MBTRI>(dataHdiv);
319  CHKERR getEntityDataOrder<MBTRI>(dataHdiv, HDIV);
320  dataHdiv.spacesOnEntities[MBTRI].set(HDIV);
321  }
322 
323  // L2
324  if (dataH1.spacesOnEntities[MBTRI].test(L2)) {
325  CHKERR getEntitySense<MBTRI>(dataL2);
326  CHKERR getEntityDataOrder<MBTRI>(dataL2, L2);
327  dataL2.spacesOnEntities[MBTRI].set(L2);
328  }
329 
331 }
DataForcesAndSourcesCore & dataL2
field with continuous normal traction
Definition: definitions.h:173
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
DataForcesAndSourcesCore & dataHdiv
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
continuous field
Definition: definitions.h:171
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
field with C-1 continuity
Definition: definitions.h:174

◆ OpSwitch()

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

Definition at line 335 of file FaceElementForcesAndSourcesCore.hpp.

335  {
337 
338  const EntityType type = numeredEntFiniteElementPtr->getEntType();
339  if (type != lastEvaluatedElementEntityType) {
340  switch (type) {
341  case MBTRI:
343  boost::shared_ptr<BaseFunction>(new TriPolynomialBase());
344  break;
345  case MBQUAD:
347  boost::shared_ptr<BaseFunction>(new QuadPolynomialBase());
348  break;
349  default:
351  }
353  }
354 
355  // Calculate normal and tangent vectors for face geometry
358 
360  if (gaussPts.size2() == 0)
362 
363  DataForcesAndSourcesCore &data_curl = *dataOnElement[HCURL];
364  DataForcesAndSourcesCore &data_div = *dataOnElement[HDIV];
365 
369 
370  switch (numeredEntFiniteElementPtr->getEntType()) {
371  case MBTRI:
372  break;
373  case MBQUAD:
375  break;
376  default:
377  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
378  "Element type not implemented");
379  }
380 
381  if (!(NO_HO_GEOMETRY & SWITCH)) {
383  }
384 
385  // Apply Piola transform to HDiv and HCurl spaces, uses previously
386  // calculated faces normal and tangent vectors.
387  if (!(NO_CONTRAVARIANT_TRANSFORM_HDIV & SWITCH))
388  if (dataH1.spacesOnEntities[MBTRI].test(HDIV))
390 
391  if (!(NO_COVARIANT_TRANSFORM_HCURL & SWITCH))
392  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL))
393  CHKERR opCovariantTransform.opRhs(data_curl);
394 
395  // Iterate over operators
397 
399 }
MatrixDouble gaussPts
Matrix of integration points.
field with continuous normal traction
Definition: definitions.h:173
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
DataForcesAndSourcesCore & dataH1
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
OpSetContravariantPiolaTransformOnFace opContravariantTransform
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
virtual MoFEMErrorCode calculateAreaAndNormalAtIntegrationPts()
Calculate element area and normal of the face at integration points.
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
virtual MoFEMErrorCode calculateHoNormal()
Calculate normal on curved elements.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

◆ setIntegrationPts()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::setIntegrationPts ( )
protectedvirtual

Set integration points.

Returns
Error code

Definition at line 167 of file FaceElementForcesAndSourcesCore.cpp.

167  {
169  // Set integration points
170  int order_data = getMaxDataOrder();
171  int order_row = getMaxRowOrder();
172  int order_col = getMaxColOrder();
173  int rule = getRule(order_row, order_col, order_data);
174 
175  auto set_integration_pts_for_tri = [&]() {
177  if (rule < QUAD_2D_TABLE_SIZE) {
178  if (QUAD_2D_TABLE[rule]->dim != 2) {
179  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
180  }
181  if (QUAD_2D_TABLE[rule]->order < rule) {
182  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
183  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
184  }
185  const size_t nb_gauss_pts = QUAD_2D_TABLE[rule]->npoints;
186  gaussPts.resize(3, nb_gauss_pts, false);
187  cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[1], 3,
188  &gaussPts(0, 0), 1);
189  cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[2], 3,
190  &gaussPts(1, 0), 1);
191  cblas_dcopy(nb_gauss_pts, QUAD_2D_TABLE[rule]->weights, 1,
192  &gaussPts(2, 0), 1);
193  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
194  false);
195  double *shape_ptr =
196  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
197  cblas_dcopy(3 * nb_gauss_pts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr,
198  1);
199  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
200  std::copy(
202  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
203  }
204  else {
205  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
206  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
207  }
209  };
210 
211  auto calc_base_for_quad = [&]() {
213  const size_t nb_gauss_pts = gaussPts.size2();
214  auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
215  auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
216  base.resize(nb_gauss_pts, 4, false);
217  diff_base.resize(nb_gauss_pts, 8, false);
218  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
219  const double ksi = gaussPts(0, gg);
220  const double zeta = gaussPts(1, gg);
221  base(gg, 0) = N_MBQUAD0(ksi, zeta);
222  base(gg, 1) = N_MBQUAD1(ksi, zeta);
223  base(gg, 2) = N_MBQUAD2(ksi, zeta);
224  base(gg, 3) = N_MBQUAD3(ksi, zeta);
225  diff_base(gg, 0) = diffN_MBQUAD0x(zeta);
226  diff_base(gg, 1) = diffN_MBQUAD0y(ksi);
227  diff_base(gg, 2) = diffN_MBQUAD1x(zeta);
228  diff_base(gg, 3) = diffN_MBQUAD1y(ksi);
229  diff_base(gg, 4) = diffN_MBQUAD2x(zeta);
230  diff_base(gg, 5) = diffN_MBQUAD2y(ksi);
231  diff_base(gg, 6) = diffN_MBQUAD3x(zeta);
232  diff_base(gg, 7) = diffN_MBQUAD3y(ksi);
233  }
235  };
236 
237  const auto type = numeredEntFiniteElementPtr->getEntType();
238 
239  if (rule >= 0) {
240  switch (type) {
241  case MBTRI:
242  CHKERR set_integration_pts_for_tri();
243  break;
244  case MBQUAD:
246  rule);
247  CHKERR calc_base_for_quad();
248  break;
249  default:
250  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
251  "Element type not implemented: %d", type);
252  }
253 
254  } else {
255  // If rule is negative, set user defined integration points
256  CHKERR setGaussPts(order_row, order_col, order_data);
257  const size_t nb_gauss_pts = gaussPts.size2();
258  auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
259  auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
260  if (nb_gauss_pts) {
261  switch (type) {
262  case MBTRI:
263  base.resize(nb_gauss_pts, 3, false);
264  diff_base.resize(3, 2, false);
265  CHKERR ShapeMBTRI(&*base.data().begin(), &gaussPts(0, 0),
266  &gaussPts(1, 0), nb_gauss_pts);
267  std::copy(
269  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
270  break;
271  case MBQUAD:
272  CHKERR calc_base_for_quad();
273  break;
274  default:
275  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
276  "Element type not implemented: %d", type);
277  }
278  }
279  }
281 }
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:67
int getMaxRowOrder() const
Get max order of approximation for field in rows.
MatrixDouble gaussPts
Matrix of integration points.
#define diffN_MBQUAD3x(y)
Definition: fem_tools.h:75
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
DataForcesAndSourcesCore & dataH1
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
int npoints
Definition: quad.h:29
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define diffN_MBQUAD0x(y)
Definition: fem_tools.h:69
#define diffN_MBQUAD0y(x)
Definition: fem_tools.h:70
#define diffN_MBQUAD2y(x)
Definition: fem_tools.h:74
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define diffN_MBQUAD2x(y)
Definition: fem_tools.h:73
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:68
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:65
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:87
#define CHKERR
Inline error check.
Definition: definitions.h:596
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:194
static MoFEMErrorCode outerProductOfEdgeIntegrationPtsForQuad(MatrixDouble &pts, const int edge0, const int edge1)
Definition: Tools.cpp:486
#define diffN_MBQUAD1x(y)
Definition: fem_tools.h:71
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:66
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
#define diffN_MBQUAD3y(x)
Definition: fem_tools.h:76
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.
#define diffN_MBQUAD1y(x)
Definition: fem_tools.h:72

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

Definition at line 306 of file FaceElementForcesAndSourcesCore.hpp.

◆ VolumeElementForcesAndSourcesCoreOnSideBase

Definition at line 307 of file FaceElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

double MoFEM::FaceElementForcesAndSourcesCoreBase::aRea
protected

Definition at line 291 of file FaceElementForcesAndSourcesCore.hpp.

◆ conn

const EntityHandle* MoFEM::FaceElementForcesAndSourcesCoreBase::conn
protected

Definition at line 293 of file FaceElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::coords
protected

Definition at line 295 of file FaceElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::coordsAtGaussPts
protected

Definition at line 296 of file FaceElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::hoCoordsAtGaussPts
protected

Definition at line 298 of file FaceElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::FaceElementForcesAndSourcesCoreBase::meshPositionsFieldName

Name of the field with geometry.

Definition at line 39 of file FaceElementForcesAndSourcesCore.hpp.

◆ nOrmal

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::nOrmal
protected

Definition at line 294 of file FaceElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::normalsAtGaussPts
protected

Definition at line 299 of file FaceElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::FaceElementForcesAndSourcesCoreBase::num_nodes
protected

Definition at line 292 of file FaceElementForcesAndSourcesCore.hpp.

◆ opContravariantTransform

OpSetContravariantPiolaTransformOnFace MoFEM::FaceElementForcesAndSourcesCoreBase::opContravariantTransform
protected

Definition at line 303 of file FaceElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnFace MoFEM::FaceElementForcesAndSourcesCoreBase::opCovariantTransform
protected

Definition at line 304 of file FaceElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnFace MoFEM::FaceElementForcesAndSourcesCoreBase::opHOCoordsAndNormals
protected

Definition at line 302 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOne

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentOne
protected

Definition at line 294 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOneAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentOneAtGaussPts
protected

Definition at line 300 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwo

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentTwo
protected

Definition at line 294 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwoAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentTwoAtGaussPts
protected

Definition at line 301 of file FaceElementForcesAndSourcesCore.hpp.


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