v0.13.1
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_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...
 
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
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 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...
 
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...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register 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
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 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 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...
 
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

double aRea
 
int num_nodes
 
const EntityHandle * conn
 
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_vector< 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
 

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
EshelbianPlasticity.cpp, MagneticElement.hpp, approx_sphere.cpp, boundary_marker.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_divergence_operator_2d.cpp, poisson_2d_homogeneous.hpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, and test_cache_on_entities.cpp.

Definition at line 25 of file FaceElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FaceElementForcesAndSourcesCore()

FaceElementForcesAndSourcesCore::FaceElementForcesAndSourcesCore ( Interface m_field)

Definition at line 11 of file FaceElementForcesAndSourcesCore.cpp.

13 : ForcesAndSourcesCore(m_field),
14 meshPositionsFieldName("MESH_NODE_POSITIONS") {}
ForcesAndSourcesCore(Interface &m_field)

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

115 {
117
118 EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
119 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
120 coords.resize(num_nodes * 3, false);
121 CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
122 nOrmal.resize(3, false);
123 tangentOne.resize(3, false);
124 tangentTwo.resize(3, false);
125
126 auto calc_normal = [&](const double *diff_ptr) {
129 &coords[0], &coords[1], &coords[2]);
131 &nOrmal[0], &nOrmal[1], &nOrmal[2]);
133 &tangentOne[0], &tangentOne[1], &tangentOne[2]);
135 &tangentTwo[0], &tangentTwo[1], &tangentTwo[2]);
137 diff_ptr, &diff_ptr[1]);
138
142
145 t_t1(i) = 0;
146 t_t2(i) = 0;
147
148 for (int nn = 0; nn != num_nodes; ++nn) {
149 t_t1(i) += t_coords(i) * t_diff(N0);
150 t_t2(i) += t_coords(i) * t_diff(N1);
151 ++t_coords;
152 ++t_diff;
153 }
154 t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
155 aRea = sqrt(t_normal(i) * t_normal(i));
157 };
158
159 const double *diff_ptr;
160 switch (numeredEntFiniteElementPtr->getEntType()) {
161 case MBTRI:
162 diff_ptr = Tools::diffShapeFunMBTRI.data();
163 CHKERR calc_normal(diff_ptr);
164 // FIXME: Normal should be divided not the area for triangle!!
165 aRea /= 2;
166 break;
167 case MBQUAD:
168 diff_ptr = Tools::diffShapeFunMBQUADAtCenter.data();
169 CHKERR calc_normal(diff_ptr);
170 break;
171 default:
172 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
173 "Element type not implemented");
174 }
175
177}
static Number< 1 > N1
static Number< 0 > N0
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
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
virtual moab::Interface & get_moab()=0
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static constexpr std::array< double, 8 > diffShapeFunMBQUADAtCenter
Definition: Tools.hpp:212
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:104

◆ calculateAreaAndNormalAtIntegrationPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateAreaAndNormalAtIntegrationPts ( )
protectedvirtual

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

Returns
Error code

Definition at line 17 of file FaceElementForcesAndSourcesCore.cpp.

17 {
19
20 auto type = numeredEntFiniteElementPtr->getEntType();
21
25
26 auto get_ftensor_from_vec_3d = [](VectorDouble &v) {
28 &v[2]);
29 };
30
31 auto get_ftensor_n_diff = [&]() {
32 const auto &m = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
34 &m(0, 1));
35 };
36
37 auto get_ftensor_from_mat_3d = [](MatrixDouble &m) {
39 &m(0, 0), &m(0, 1), &m(0, 2));
40 };
41
42 if (type == MBTRI) {
43
44 const size_t nb_gauss_pts = gaussPts.size2();
45 normalsAtGaussPts.resize(nb_gauss_pts, 3);
46 tangentOneAtGaussPts.resize(nb_gauss_pts, 3);
47 tangentTwoAtGaussPts.resize(nb_gauss_pts, 3);
48
49 auto t_tan1 = get_ftensor_from_mat_3d(tangentOneAtGaussPts);
50 auto t_tan2 = get_ftensor_from_mat_3d(tangentTwoAtGaussPts);
51 auto t_normal = get_ftensor_from_mat_3d(normalsAtGaussPts);
52
53 auto t_n =
56 &tangentOne[2]);
58 &tangentTwo[2]);
59
60 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
61 t_normal(i) = t_n(i);
62 t_tan1(i) = t_t1(i);
63 t_tan2(i) = t_t2(i);
64 ++t_tan1;
65 ++t_tan2;
66 ++t_normal;
67 }
68
69 } else if (type == MBQUAD) {
70
71 EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
72 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
73 coords.resize(num_nodes * 3, false);
74 CHKERR mField.get_moab().get_coords(conn, num_nodes,
75 &*coords.data().begin());
76
77 const size_t nb_gauss_pts = gaussPts.size2();
78 normalsAtGaussPts.resize(nb_gauss_pts, 3);
79 tangentOneAtGaussPts.resize(nb_gauss_pts, 3);
80 tangentTwoAtGaussPts.resize(nb_gauss_pts, 3);
81 normalsAtGaussPts.clear();
84
85 auto t_t1 = get_ftensor_from_mat_3d(tangentOneAtGaussPts);
86 auto t_t2 = get_ftensor_from_mat_3d(tangentTwoAtGaussPts);
87 auto t_normal = get_ftensor_from_mat_3d(normalsAtGaussPts);
88
91
92 auto t_diff = get_ftensor_n_diff();
93 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
94 auto t_coords = get_ftensor_from_vec_3d(coords);
95 for (int nn = 0; nn != num_nodes; ++nn) {
96 t_t1(i) += t_coords(i) * t_diff(N0);
97 t_t2(i) += t_coords(i) * t_diff(N1);
98 ++t_diff;
99 ++t_coords;
100 }
101 t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
102
103 ++t_t1;
104 ++t_t2;
105 ++t_normal;
106 }
107 } else {
108 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
109 "Element type not implemented");
110 }
111
113}
@ NOBASE
Definition: definitions.h:59
FTensor::Index< 'm', SPACE_DIM > m
const double v
phase velocity of light in medium (cm/ns)
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
MatrixDouble gaussPts
Matrix of integration points.

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
protectedvirtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 360 of file FaceElementForcesAndSourcesCore.cpp.

360 {
362
363 const size_t nb_nodes =
364 dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).size2();
365 double *shape_functions =
366 &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
367 const size_t nb_gauss_pts = gaussPts.size2();
368 coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
369 for (int gg = 0; gg != nb_gauss_pts; ++gg)
370 for (int dd = 0; dd != 3; ++dd)
371 coordsAtGaussPts(gg, dd) = cblas_ddot(
372 nb_nodes, &shape_functions[nb_nodes * gg], 1, &coords[dd], 3);
373
375}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
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
MatrixDouble coordsAtGaussPts
coordinated at gauss points

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode FaceElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
protectedvirtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 295 of file FaceElementForcesAndSourcesCore.cpp.

295 {
297 // Get spaces order/base and sense of entities.
298
300
301 // H1
302 if (dataH1.spacesOnEntities[MBEDGE].test(H1)) {
303 CHKERR getEntitySense<MBEDGE>(dataH1);
304 CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
305 }
306 if (dataH1.spacesOnEntities[MBTRI].test(H1)) {
307 CHKERR getEntitySense<MBTRI>(dataH1);
308 CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
309 }
310 if (dataH1.spacesOnEntities[MBQUAD].test(H1)) {
311 CHKERR getEntitySense<MBQUAD>(dataH1);
312 CHKERR getEntityDataOrder<MBQUAD>(dataH1, H1);
313 }
314
315 // Hcurl
316 if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
317 CHKERR getEntitySense<MBEDGE>(dataHcurl);
318 CHKERR getEntityDataOrder<MBEDGE>(dataHcurl, HCURL);
319 dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
320 }
321 if (dataH1.spacesOnEntities[MBTRI].test(HCURL)) {
322 CHKERR getEntitySense<MBTRI>(dataHcurl);
323 CHKERR getEntityDataOrder<MBTRI>(dataHcurl, HCURL);
324 dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
325 }
326 if (dataH1.spacesOnEntities[MBQUAD].test(HCURL)) {
327 CHKERR getEntitySense<MBQUAD>(dataHcurl);
328 CHKERR getEntityDataOrder<MBQUAD>(dataHcurl, HCURL);
329 dataHcurl.spacesOnEntities[MBQUAD].set(HCURL);
330 }
331
332 // Hdiv
333 if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
334 CHKERR getEntitySense<MBTRI>(dataHdiv);
335 CHKERR getEntityDataOrder<MBTRI>(dataHdiv, HDIV);
336 dataHdiv.spacesOnEntities[MBTRI].set(HDIV);
337 }
338 if (dataH1.spacesOnEntities[MBQUAD].test(HDIV)) {
339 CHKERR getEntitySense<MBQUAD>(dataHdiv);
340 CHKERR getEntityDataOrder<MBQUAD>(dataHdiv, HDIV);
341 dataHdiv.spacesOnEntities[MBQUAD].set(HDIV);
342 }
343
344 // L2
345 if (dataH1.spacesOnEntities[MBTRI].test(L2)) {
346 CHKERR getEntitySense<MBTRI>(dataL2);
347 CHKERR getEntityDataOrder<MBTRI>(dataL2, L2);
348 dataL2.spacesOnEntities[MBTRI].set(L2);
349 }
350 if (dataH1.spacesOnEntities[MBQUAD].test(L2)) {
351 CHKERR getEntitySense<MBQUAD>(dataL2);
352 CHKERR getEntityDataOrder<MBQUAD>(dataL2, L2);
353 dataL2.spacesOnEntities[MBQUAD].set(L2);
354 }
355
357}
@ L2
field with C-1 continuity
Definition: definitions.h:88
@ H1
continuous field
Definition: definitions.h:85
@ HCURL
field with continuous tangents
Definition: definitions.h:86
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.

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

386 {
388
389 const auto type = numeredEntFiniteElementPtr->getEntType();
391 switch (type) {
392 case MBTRI:
394 boost::shared_ptr<BaseFunction>(new TriPolynomialBase());
395 break;
396 case MBQUAD:
398 boost::shared_ptr<BaseFunction>(new QuadPolynomialBase());
399 break;
400 default:
402 }
405 }
406
407 // Calculate normal and tangent vectors for face geometry
410
412 if (gaussPts.size2() == 0)
414
419
420 // Iterate over operators
422
424}
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
virtual MoFEMErrorCode calculateAreaAndNormalAtIntegrationPts()
Calculate element area and normal of the face at integration points.
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.

◆ setIntegrationPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::setIntegrationPts ( )
protectedvirtual

Set integration points.

Returns
Error code

Definition at line 179 of file FaceElementForcesAndSourcesCore.cpp.

179 {
181 // Set integration points
182 int order_data = getMaxDataOrder();
183 int order_row = getMaxRowOrder();
184 int order_col = getMaxColOrder();
185 int rule = getRule(order_row, order_col, order_data);
186
187 auto set_integration_pts_for_tri = [&]() {
189 if (rule < QUAD_2D_TABLE_SIZE) {
190 if (QUAD_2D_TABLE[rule]->dim != 2) {
191 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
192 }
193 if (QUAD_2D_TABLE[rule]->order < rule) {
194 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
195 "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
196 }
197 const size_t nb_gauss_pts = QUAD_2D_TABLE[rule]->npoints;
198 gaussPts.resize(3, nb_gauss_pts, false);
199 cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[1], 3,
200 &gaussPts(0, 0), 1);
201 cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[2], 3,
202 &gaussPts(1, 0), 1);
203 cblas_dcopy(nb_gauss_pts, QUAD_2D_TABLE[rule]->weights, 1,
204 &gaussPts(2, 0), 1);
205 dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
206 false);
207 double *shape_ptr =
208 &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
209 cblas_dcopy(3 * nb_gauss_pts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr,
210 1);
211 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
212 std::copy(
214 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
215 } else {
216 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
217 "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
218 }
220 };
221
222 auto calc_base_for_quad = [&]() {
224 const size_t nb_gauss_pts = gaussPts.size2();
225 auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
226 auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
227 base.resize(nb_gauss_pts, 4, false);
228 diff_base.resize(nb_gauss_pts, 8, false);
229 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
230 const double ksi = gaussPts(0, gg);
231 const double zeta = gaussPts(1, gg);
232 base(gg, 0) = N_MBQUAD0(ksi, zeta);
233 base(gg, 1) = N_MBQUAD1(ksi, zeta);
234 base(gg, 2) = N_MBQUAD2(ksi, zeta);
235 base(gg, 3) = N_MBQUAD3(ksi, zeta);
236 diff_base(gg, 0) = diffN_MBQUAD0x(zeta);
237 diff_base(gg, 1) = diffN_MBQUAD0y(ksi);
238 diff_base(gg, 2) = diffN_MBQUAD1x(zeta);
239 diff_base(gg, 3) = diffN_MBQUAD1y(ksi);
240 diff_base(gg, 4) = diffN_MBQUAD2x(zeta);
241 diff_base(gg, 5) = diffN_MBQUAD2y(ksi);
242 diff_base(gg, 6) = diffN_MBQUAD3x(zeta);
243 diff_base(gg, 7) = diffN_MBQUAD3y(ksi);
244 }
246 };
247
248 const auto type = numeredEntFiniteElementPtr->getEntType();
249
250 if (rule >= 0) {
251 switch (type) {
252 case MBTRI:
253 CHKERR set_integration_pts_for_tri();
254 break;
255 case MBQUAD:
257 rule);
258 CHKERR calc_base_for_quad();
259 break;
260 default:
261 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
262 "Element type not implemented: %d", type);
263 }
264
265 } else {
266 // If rule is negative, set user defined integration points
267 CHKERR setGaussPts(order_row, order_col, order_data);
268 const size_t nb_gauss_pts = gaussPts.size2();
269 auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
270 auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
271 if (nb_gauss_pts) {
272 switch (type) {
273 case MBTRI:
274 base.resize(nb_gauss_pts, 3, false);
275 diff_base.resize(3, 2, false);
276 CHKERR ShapeMBTRI(&*base.data().begin(), &gaussPts(0, 0),
277 &gaussPts(1, 0), nb_gauss_pts);
278 std::copy(
280 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
281 break;
282 case MBQUAD:
283 CHKERR calc_base_for_quad();
284 break;
285 default:
286 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
287 "Element type not implemented: %d", type);
288 }
289 }
290 }
292}
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
const int dim
#define diffN_MBQUAD2y(x)
Definition: fem_tools.h:66
#define diffN_MBQUAD1x(y)
Definition: fem_tools.h:63
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:60
#define diffN_MBQUAD0x(y)
Definition: fem_tools.h:61
#define diffN_MBQUAD1y(x)
Definition: fem_tools.h:64
#define diffN_MBQUAD3y(x)
Definition: fem_tools.h:68
#define diffN_MBQUAD0y(x)
Definition: fem_tools.h:62
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:57
#define diffN_MBQUAD3x(y)
Definition: fem_tools.h:67
#define diffN_MBQUAD2x(y)
Definition: fem_tools.h:65
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:59
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:58
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
int getMaxRowOrder() const
Get max order of approximation for field in rows.
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.
static MoFEMErrorCode outerProductOfEdgeIntegrationPtsForQuad(MatrixDouble &pts, const int edge0, const int edge1)
Definition: Tools.cpp:580
int npoints
Definition: quad.h:29

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

◆ VolumeElementForcesAndSourcesCoreOnSide

Definition at line 89 of file FaceElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

double MoFEM::FaceElementForcesAndSourcesCore::aRea
protected

Definition at line 78 of file FaceElementForcesAndSourcesCore.hpp.

◆ conn

const EntityHandle* MoFEM::FaceElementForcesAndSourcesCore::conn
protected

Definition at line 80 of file FaceElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::coords
protected

Definition at line 82 of file FaceElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

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

Definition at line 31 of file FaceElementForcesAndSourcesCore.hpp.

◆ nOrmal

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::nOrmal
protected

Definition at line 81 of file FaceElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPts
protected

Definition at line 84 of file FaceElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::FaceElementForcesAndSourcesCore::num_nodes
protected

Definition at line 79 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOne

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOne
protected

Definition at line 81 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOneAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOneAtGaussPts
protected

Definition at line 85 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwo

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwo
protected

Definition at line 81 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwoAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwoAtGaussPts
protected

Definition at line 86 of file FaceElementForcesAndSourcesCore.hpp.


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