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

Face finite element. 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

MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
 FaceElementForcesAndSourcesCore (Interface &m_field)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
boost::ptr_deque< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator. More...
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object. More...
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
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...
 
auto & getEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
auto & getDataOnElementBySpaceArray ()
 Get data on entities and space. More...
 
auto & getDerivedDataOnElementBySpaceArray ()
 Get derived data on entities and space. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference 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
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Public Attributes

std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
GaussHookFun 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::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. 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
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- 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
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
PetscReal ts_dt
 time step size 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...
 

Protected Member Functions

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...
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (EntitiesFieldData &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceNodes (EntitiesFieldData &data) const
 Get nodes on faces. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (EntitiesFieldData &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 (EntityType type)
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getNodesIndices (const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getEntityIndices (EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
 
MoFEMErrorCode getEntityRowIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getBitRefLevelOnData ()
 
MoFEMErrorCode getNoFieldFieldData (const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (EntitiesFieldData &data, const int bit_number) const
 
MoFEMErrorCode getNodesFieldData (EntitiesFieldData &data, const int bit_number) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (EntitiesFieldData &data, const int bit_number, 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 Attributes

doubleaRea
 
int num_nodes
 
const EntityHandleconn
 
VectorDouble nOrmal
 
VectorDouble tangentOne
 
VectorDouble tangentTwo
 
VectorDouble coords
 
MatrixDouble normalsAtGaussPts
 
MatrixDouble tangentOneAtGaussPts
 
MatrixDouble tangentTwoAtGaussPts
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
EntitiesFieldDatadataNoField
 
EntitiesFieldDatadataH1
 
EntitiesFieldDatadataHcurl
 
EntitiesFieldDatadataHdiv
 
EntitiesFieldDatadataL2
 
boost::ptr_deque< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 
MatrixDouble coordsAtGaussPts
 coordinated at gauss points More...
 
double elementMeasure
 

Friends

class UserDataOperator
 
class VolumeElementForcesAndSourcesCoreOnSide
 
template<int DIM>
struct OpCopyGeomDataToE
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
typedef boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> GaussHookFun
 
- 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::PetscData
enum  DataContext {
  CTX_SET_NONE = 0, CTX_SET_F = 1 << 0, CTX_SET_A = 1 << 1, CTX_SET_B = 1 << 2,
  CTX_SET_X = 1 << 3, CTX_SET_X_T = 1 << 4, CTX_SET_X_TT = 1 << 6, CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- 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
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

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
approx_sphere.cpp, boundary_marker.cpp, continuity_check_on_skeleton_3d.cpp, eigen_elastic.cpp, elasticity.cpp, EshelbianPlasticity.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, hello_world.cpp, helmholtz.cpp, MagneticElement.hpp, navier_stokes.cpp, phase.cpp, photon_diffusion.cpp, plot_base.cpp, PoissonOperators.hpp, prism_elements_from_surface.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, shallow_wave.cpp, simple_contact.cpp, simple_contact_thermal.cpp, simple_elasticity.cpp, simple_interface.cpp, and test_cache_on_entities.cpp.

Definition at line 23 of file FaceElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FaceElementForcesAndSourcesCore()

FaceElementForcesAndSourcesCore::FaceElementForcesAndSourcesCore ( Interface m_field)

Definition at line 9 of file FaceElementForcesAndSourcesCore.cpp.

11  : ForcesAndSourcesCore(m_field),
12  meshPositionsFieldName("MESH_NODE_POSITIONS"), aRea(elementMeasure) {}

Member Function Documentation

◆ calculateAreaAndNormal()

MoFEMErrorCode FaceElementForcesAndSourcesCore::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 113 of file FaceElementForcesAndSourcesCore.cpp.

113  {
115 
117  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
118  coords.resize(num_nodes * 3, false);
119  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
120  nOrmal.resize(3, false);
121  tangentOne.resize(3, false);
122  tangentTwo.resize(3, false);
123 
124  auto calc_normal = [&](const double *diff_ptr) {
127  &coords[0], &coords[1], &coords[2]);
129  &nOrmal[0], &nOrmal[1], &nOrmal[2]);
131  &tangentOne[0], &tangentOne[1], &tangentOne[2]);
133  &tangentTwo[0], &tangentTwo[1], &tangentTwo[2]);
135  diff_ptr, &diff_ptr[1]);
136 
140 
143  t_t1(i) = 0;
144  t_t2(i) = 0;
145 
146  for (int nn = 0; nn != num_nodes; ++nn) {
147  t_t1(i) += t_coords(i) * t_diff(N0);
148  t_t2(i) += t_coords(i) * t_diff(N1);
149  ++t_coords;
150  ++t_diff;
151  }
152  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
153  aRea = sqrt(t_normal(i) * t_normal(i));
155  };
156 
157  const double *diff_ptr;
158  switch (numeredEntFiniteElementPtr->getEntType()) {
159  case MBTRI:
160  diff_ptr = Tools::diffShapeFunMBTRI.data();
161  CHKERR calc_normal(diff_ptr);
162  // FIXME: Normal should be divided not the area for triangle!!
163  aRea /= 2;
164  break;
165  case MBQUAD:
166  diff_ptr = Tools::diffShapeFunMBQUADAtCenter.data();
167  CHKERR calc_normal(diff_ptr);
168  break;
169  default:
170  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
171  "Element type not implemented");
172  }
173 
175 }

◆ calculateAreaAndNormalAtIntegrationPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateAreaAndNormalAtIntegrationPts ( )
protectedvirtual

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

Returns
Error code

Definition at line 15 of file FaceElementForcesAndSourcesCore.cpp.

15  {
17 
18  auto type = numeredEntFiniteElementPtr->getEntType();
19 
23 
24  auto get_ftensor_from_vec_3d = [](VectorDouble &v) {
26  &v[2]);
27  };
28 
29  auto get_ftensor_n_diff = [&]() {
30  const auto &m = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
32  &m(0, 1));
33  };
34 
35  auto get_ftensor_from_mat_3d = [](MatrixDouble &m) {
37  &m(0, 0), &m(0, 1), &m(0, 2));
38  };
39 
40  if (type == MBTRI) {
41 
42  const size_t nb_gauss_pts = gaussPts.size2();
43  normalsAtGaussPts.resize(nb_gauss_pts, 3);
44  tangentOneAtGaussPts.resize(nb_gauss_pts, 3);
45  tangentTwoAtGaussPts.resize(nb_gauss_pts, 3);
46 
47  auto t_tan1 = get_ftensor_from_mat_3d(tangentOneAtGaussPts);
48  auto t_tan2 = get_ftensor_from_mat_3d(tangentTwoAtGaussPts);
49  auto t_normal = get_ftensor_from_mat_3d(normalsAtGaussPts);
50 
51  auto t_n =
54  &tangentOne[2]);
56  &tangentTwo[2]);
57 
58  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
59  t_normal(i) = t_n(i);
60  t_tan1(i) = t_t1(i);
61  t_tan2(i) = t_t2(i);
62  ++t_tan1;
63  ++t_tan2;
64  ++t_normal;
65  }
66 
67  } else if (type == MBQUAD) {
68 
70  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
71  coords.resize(num_nodes * 3, false);
72  CHKERR mField.get_moab().get_coords(conn, num_nodes,
73  &*coords.data().begin());
74 
75  const size_t nb_gauss_pts = gaussPts.size2();
76  normalsAtGaussPts.resize(nb_gauss_pts, 3);
77  tangentOneAtGaussPts.resize(nb_gauss_pts, 3);
78  tangentTwoAtGaussPts.resize(nb_gauss_pts, 3);
79  normalsAtGaussPts.clear();
80  tangentOneAtGaussPts.clear();
81  tangentTwoAtGaussPts.clear();
82 
83  auto t_t1 = get_ftensor_from_mat_3d(tangentOneAtGaussPts);
84  auto t_t2 = get_ftensor_from_mat_3d(tangentTwoAtGaussPts);
85  auto t_normal = get_ftensor_from_mat_3d(normalsAtGaussPts);
86 
89 
90  auto t_diff = get_ftensor_n_diff();
91  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
92  auto t_coords = get_ftensor_from_vec_3d(coords);
93  for (int nn = 0; nn != num_nodes; ++nn) {
94  t_t1(i) += t_coords(i) * t_diff(N0);
95  t_t2(i) += t_coords(i) * t_diff(N1);
96  ++t_diff;
97  ++t_coords;
98  }
99  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
100 
101  ++t_t1;
102  ++t_t2;
103  ++t_normal;
104  }
105  } else {
106  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
107  "Element type not implemented");
108  }
109 
111 }

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
protectedvirtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 346 of file FaceElementForcesAndSourcesCore.cpp.

346  {
348 
349  const size_t nb_nodes =
350  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).size2();
351  double *shape_functions =
352  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
353  const size_t nb_gauss_pts = gaussPts.size2();
354  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
355  for (int gg = 0; gg != nb_gauss_pts; ++gg)
356  for (int dd = 0; dd != 3; ++dd)
357  coordsAtGaussPts(gg, dd) = cblas_ddot(
358  nb_nodes, &shape_functions[nb_nodes * gg], 1, &coords[dd], 3);
359 
361 }

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode FaceElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
protectedvirtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 310 of file FaceElementForcesAndSourcesCore.cpp.

310  {
312  // Get spaces order/base and sense of entities.
313 
315 
316  auto type = numeredEntFiniteElementPtr->getEntType();
317  auto dim_type = CN::Dimension(type);
318 
319  auto get_data_on_ents = [&](auto lower_dim, auto space) {
321  auto data = dataOnElement[space];
322  for (auto dd = dim_type; dd >= lower_dim; --dd) {
323  int nb_ents = moab::CN::NumSubEntities(type, dd);
324  for (int ii = 0; ii != nb_ents; ++ii) {
325  auto sub_ent_type = moab::CN::SubEntityType(type, dd, ii);
326  if ((dataH1.spacesOnEntities[sub_ent_type]).test(space)) {
327  auto &data_on_ent = data->dataOnEntities[sub_ent_type];
328  CHKERR getEntitySense(sub_ent_type, data_on_ent);
329  CHKERR getEntityDataOrder(sub_ent_type, space, data_on_ent);
330  data->spacesOnEntities[sub_ent_type].set(space);
331  }
332  }
333  }
335  };
336 
337  CHKERR get_data_on_ents(1, H1); // H1
338  CHKERR get_data_on_ents(1, HCURL); // Hcurl
339  CHKERR get_data_on_ents(2, HDIV); // Hdiv
340  CHKERR get_data_on_ents(2, L2); // L2
341 
343 }

◆ 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 NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE.

Definition at line 372 of file FaceElementForcesAndSourcesCore.cpp.

372  {
374 
375  const auto type = numeredEntFiniteElementPtr->getEntType();
377  switch (type) {
378  case MBTRI:
380  boost::shared_ptr<BaseFunction>(new TriPolynomialBase());
381  break;
382  case MBQUAD:
384  boost::shared_ptr<BaseFunction>(new QuadPolynomialBase());
385  break;
386  default:
388  }
391  }
392 
393  // Calculate normal and tangent vectors for face geometry
396 
398  if (gaussPts.size2() == 0)
400 
405 
406  // Iterate over operators
408 
410 }

◆ setIntegrationPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::setIntegrationPts ( )
protectedvirtual

Set integration points.

Returns
Error code

Definition at line 177 of file FaceElementForcesAndSourcesCore.cpp.

177  {
179  // Set integration points
180  int order_data = getMaxDataOrder();
181  int order_row = getMaxRowOrder();
182  int order_col = getMaxColOrder();
183 
184  const auto type = numeredEntFiniteElementPtr->getEntType();
185 
186  auto get_rule_by_type = [&]() {
187  switch (type) {
188  case MBQUAD:
189  return getRule(order_row + 1, order_col + 1, order_data + 1);
190  default:
191  return getRule(order_row, order_col, order_data);
192  }
193  };
194 
195  const int rule = get_rule_by_type();
196 
197  auto set_integration_pts_for_tri = [&]() {
199  if (rule < QUAD_2D_TABLE_SIZE) {
200  if (QUAD_2D_TABLE[rule]->dim != 2) {
201  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
202  }
203  if (QUAD_2D_TABLE[rule]->order < rule) {
204  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
205  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
206  }
207  const size_t nb_gauss_pts = QUAD_2D_TABLE[rule]->npoints;
208  gaussPts.resize(3, nb_gauss_pts, false);
209  cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[1], 3,
210  &gaussPts(0, 0), 1);
211  cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[2], 3,
212  &gaussPts(1, 0), 1);
213  cblas_dcopy(nb_gauss_pts, QUAD_2D_TABLE[rule]->weights, 1,
214  &gaussPts(2, 0), 1);
215  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
216  false);
217  double *shape_ptr =
218  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
219  cblas_dcopy(3 * nb_gauss_pts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr,
220  1);
221  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
222  std::copy(
224  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
225  } else {
226  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
227  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
228  }
230  };
231 
232  auto calc_base_for_tri = [&]() {
234  const size_t nb_gauss_pts = gaussPts.size2();
235  auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
236  auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
237  base.resize(nb_gauss_pts, 3, false);
238  diff_base.resize(3, 2, false);
239  CHKERR ShapeMBTRI(&*base.data().begin(), &gaussPts(0, 0), &gaussPts(1, 0),
240  nb_gauss_pts);
241  std::copy(
243  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
245  };
246 
247  auto calc_base_for_quad = [&]() {
249  const size_t nb_gauss_pts = gaussPts.size2();
250  auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
251  auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
252  base.resize(nb_gauss_pts, 4, false);
253  diff_base.resize(nb_gauss_pts, 8, false);
254  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
255  const double ksi = gaussPts(0, gg);
256  const double zeta = gaussPts(1, gg);
257  base(gg, 0) = N_MBQUAD0(ksi, zeta);
258  base(gg, 1) = N_MBQUAD1(ksi, zeta);
259  base(gg, 2) = N_MBQUAD2(ksi, zeta);
260  base(gg, 3) = N_MBQUAD3(ksi, zeta);
261  diff_base(gg, 0) = diffN_MBQUAD0x(zeta);
262  diff_base(gg, 1) = diffN_MBQUAD0y(ksi);
263  diff_base(gg, 2) = diffN_MBQUAD1x(zeta);
264  diff_base(gg, 3) = diffN_MBQUAD1y(ksi);
265  diff_base(gg, 4) = diffN_MBQUAD2x(zeta);
266  diff_base(gg, 5) = diffN_MBQUAD2y(ksi);
267  diff_base(gg, 6) = diffN_MBQUAD3x(zeta);
268  diff_base(gg, 7) = diffN_MBQUAD3y(ksi);
269  }
271  };
272 
273  if (rule >= 0) {
274  switch (type) {
275  case MBTRI:
276  CHKERR set_integration_pts_for_tri();
277  break;
278  case MBQUAD:
280  rule);
281  CHKERR calc_base_for_quad();
282  break;
283  default:
284  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
285  "Element type not implemented: %d", type);
286  }
287 
288  } else {
289  // If rule is negative, set user defined integration points
290  CHKERR setGaussPts(order_row, order_col, order_data);
291  const size_t nb_gauss_pts = gaussPts.size2();
292  if (nb_gauss_pts) {
293  switch (type) {
294  case MBTRI:
295  CHKERR calc_base_for_tri();
296  break;
297  case MBQUAD:
298  CHKERR calc_base_for_quad();
299  break;
300  default:
301  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
302  "Element type not implemented: %d", type);
303  }
304  }
305  }
307 }

Friends And Related Function Documentation

◆ OpCopyGeomDataToE

template<int DIM>
friend struct OpCopyGeomDataToE
friend

Definition at line 88 of file FaceElementForcesAndSourcesCore.hpp.

◆ UserDataOperator

friend class UserDataOperator
friend

◆ VolumeElementForcesAndSourcesCoreOnSide

Definition at line 87 of file FaceElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

double& MoFEM::FaceElementForcesAndSourcesCore::aRea
protected

Definition at line 76 of file FaceElementForcesAndSourcesCore.hpp.

◆ conn

const EntityHandle* MoFEM::FaceElementForcesAndSourcesCore::conn
protected

Definition at line 78 of file FaceElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::coords
protected

Definition at line 80 of file FaceElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::FaceElementForcesAndSourcesCore::meshPositionsFieldName
Deprecated:
not used anymore, will be removed in next versions
Examples
MagneticElement.hpp.

Definition at line 29 of file FaceElementForcesAndSourcesCore.hpp.

◆ nOrmal

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::nOrmal
protected

Definition at line 79 of file FaceElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPts
protected

Definition at line 82 of file FaceElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::FaceElementForcesAndSourcesCore::num_nodes
protected

Definition at line 77 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOne

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOne
protected

Definition at line 79 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOneAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOneAtGaussPts
protected

Definition at line 83 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwo

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwo
protected

Definition at line 79 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwoAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwoAtGaussPts
protected

Definition at line 84 of file FaceElementForcesAndSourcesCore.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::FaceElementForcesAndSourcesCore::tangentOne
VectorDouble tangentOne
Definition: FaceElementForcesAndSourcesCore.hpp:79
MoFEM::FaceElementForcesAndSourcesCore::tangentTwo
VectorDouble tangentTwo
Definition: FaceElementForcesAndSourcesCore.hpp:79
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
Definition: ForcesAndSourcesCore.cpp:1132
diffN_MBQUAD1x
#define diffN_MBQUAD1x(y)
Definition: fem_tools.h:63
MoFEM::ForcesAndSourcesCore::getMaxRowOrder
int getMaxRowOrder() const
Get max order of approximation for field in rows.
Definition: ForcesAndSourcesCore.cpp:131
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
MoFEM::FaceElementForcesAndSourcesCore::meshPositionsFieldName
std::string meshPositionsFieldName
Definition: FaceElementForcesAndSourcesCore.hpp:29
MoFEM::ForcesAndSourcesCore::lastEvaluatedElementEntityType
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
Definition: ForcesAndSourcesCore.hpp:489
NOBASE
@ NOBASE
Definition: definitions.h:59
MoFEM::FaceElementForcesAndSourcesCore::tangentTwoAtGaussPts
MatrixDouble tangentTwoAtGaussPts
Definition: FaceElementForcesAndSourcesCore.hpp:84
diffN_MBQUAD1y
#define diffN_MBQUAD1y(x)
Definition: fem_tools.h:64
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::FaceElementForcesAndSourcesCore::calculateAreaAndNormal
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.
Definition: FaceElementForcesAndSourcesCore.cpp:113
QUAD_2D_TABLE
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::ForcesAndSourcesCore::setGaussPts
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: ForcesAndSourcesCore.cpp:1961
MoFEM::ForcesAndSourcesCore::getMaxDataOrder
int getMaxDataOrder() const
Get max order of approximation for data fields.
Definition: ForcesAndSourcesCore.cpp:120
diffN_MBQUAD2x
#define diffN_MBQUAD2x(y)
Definition: fem_tools.h:65
MoFEM::ForcesAndSourcesCore::getEntityDataOrder
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.
Definition: ForcesAndSourcesCore.cpp:139
QUAD_2D_TABLE_SIZE
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
MoFEM::ForcesAndSourcesCore::getElementPolynomialBase
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
Definition: ForcesAndSourcesCore.hpp:90
MoFEM::ForcesAndSourcesCore::getMaxColOrder
int getMaxColOrder() const
Get max order of approximation for field in columns.
Definition: ForcesAndSourcesCore.cpp:135
MoFEM::ForcesAndSourcesCore::calHierarchicalBaseFunctionsOnElement
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
Definition: ForcesAndSourcesCore.cpp:1112
FTensor::levi_civita
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
Definition: Levi_Civita.hpp:617
MoFEM::FaceElementForcesAndSourcesCore::aRea
double & aRea
Definition: FaceElementForcesAndSourcesCore.hpp:76
N_MBQUAD2
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:59
MoFEM::EntitiesFieldData::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: EntitiesFieldData.hpp:49
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:126
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::ForcesAndSourcesCore::elementMeasure
double elementMeasure
Definition: ForcesAndSourcesCore.hpp:545
FTensor::Number< 0 >
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::FaceElementForcesAndSourcesCore::conn
const EntityHandle * conn
Definition: FaceElementForcesAndSourcesCore.hpp:78
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
diffN_MBQUAD3x
#define diffN_MBQUAD3x(y)
Definition: fem_tools.h:67
MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
coordinated at gauss points
Definition: ForcesAndSourcesCore.hpp:544
MoFEM::FaceElementForcesAndSourcesCore::num_nodes
int num_nodes
Definition: FaceElementForcesAndSourcesCore.hpp:77
convert.type
type
Definition: convert.py:64
MoFEM::ForcesAndSourcesCore::dataOnElement
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
Definition: ForcesAndSourcesCore.hpp:461
MoFEM::ForcesAndSourcesCore::createDataOnElement
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
Definition: ForcesAndSourcesCore.cpp:1337
MoFEM::ForcesAndSourcesCore::ForcesAndSourcesCore
ForcesAndSourcesCore(Interface &m_field)
Definition: ForcesAndSourcesCore.cpp:40
QUAD_::npoints
int npoints
Definition: quad.h:29
MoFEM::ForcesAndSourcesCore::getRule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: ForcesAndSourcesCore.cpp:1955
N_MBQUAD1
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:58
diffN_MBQUAD2y
#define diffN_MBQUAD2y(x)
Definition: fem_tools.h:66
MoFEM::Tools::diffShapeFunMBQUADAtCenter
static constexpr std::array< double, 8 > diffShapeFunMBQUADAtCenter
Definition: Tools.hpp:212
MoFEM::Tools::diffShapeFunMBTRI
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:104
diffN_MBQUAD0x
#define diffN_MBQUAD0x(y)
Definition: fem_tools.h:61
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
Definition: ForcesAndSourcesCore.cpp:990
diffN_MBQUAD0y
#define diffN_MBQUAD0y(x)
Definition: fem_tools.h:62
FTensor::Index< 'i', 3 >
diffN_MBQUAD3y
#define diffN_MBQUAD3y(x)
Definition: fem_tools.h:68
MoFEM::ForcesAndSourcesCore::dataH1
EntitiesFieldData & dataH1
Definition: ForcesAndSourcesCore.hpp:471
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
FTensor::dd
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
MoFEM::FaceElementForcesAndSourcesCore::tangentOneAtGaussPts
MatrixDouble tangentOneAtGaussPts
Definition: FaceElementForcesAndSourcesCore.hpp:83
MoFEM::FaceElementForcesAndSourcesCore::nOrmal
VectorDouble nOrmal
Definition: FaceElementForcesAndSourcesCore.hpp:79
MoFEM::FaceElementForcesAndSourcesCore::coords
VectorDouble coords
Definition: FaceElementForcesAndSourcesCore.hpp:80
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
N_MBQUAD0
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:57
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
MoFEM::ForcesAndSourcesCore::loopOverOperators
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
Definition: ForcesAndSourcesCore.cpp:1379
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
N_MBQUAD3
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:60
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEM::FaceElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
Definition: FaceElementForcesAndSourcesCore.cpp:346
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::Tools::outerProductOfEdgeIntegrationPtsForQuad
static MoFEMErrorCode outerProductOfEdgeIntegrationPtsForQuad(MatrixDouble &pts, const int edge0, const int edge1)
Definition: Tools.cpp:610
MoFEM::FaceElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
Definition: FaceElementForcesAndSourcesCore.cpp:310
MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPts
MatrixDouble normalsAtGaussPts
Definition: FaceElementForcesAndSourcesCore.hpp:82
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEM::FaceElementForcesAndSourcesCore::calculateAreaAndNormalAtIntegrationPts
virtual MoFEMErrorCode calculateAreaAndNormalAtIntegrationPts()
Calculate element area and normal of the face at integration points.
Definition: FaceElementForcesAndSourcesCore.cpp:15
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::FaceElementForcesAndSourcesCore::setIntegrationPts
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
Definition: FaceElementForcesAndSourcesCore.cpp:177
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:24
MoFEM::ForcesAndSourcesCore::getEntitySense
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity
Definition: ForcesAndSourcesCore.cpp:84
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182