v0.15.0
Loading...
Searching...
No Matches
MoFEM::ForcesAndSourcesCore Struct Reference

structure to get information form mofem into EntitiesFieldData More...

#include "src/finite_elements/ForcesAndSourcesCore.hpp"

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

Classes

struct  UserDataOperator
 

Public Types

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_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , 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
}
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 

Public Member Functions

 ForcesAndSourcesCore (Interface &m_field)
 
boost::ptr_deque< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator.
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object.
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object.
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop
 
int getMaxDataOrder () const
 Get max order of approximation for data fields.
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows.
 
int getMaxColOrder () const
 Get max order of approximation for field in columns.
 
auto & getEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data.
 
auto & getDataOnElementBySpaceArray ()
 Get data on entities and space.
 
auto & getDerivedDataOnElementBySpaceArray ()
 Get derived data on entities and space.
 
- 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
 
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
 
int getLoopSize () const
 get loop size
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities.
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements.
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements.
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object.
 
- 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
 
- 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.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
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.
 
- 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.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TaoMethod ()
 
virtual ~TaoMethod ()=default
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data.
 

Public Attributes

InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule.
 
GaussHookFun setRuleHook
 Set function to calculate integration rule.
 
MatrixDouble gaussPts
 Matrix of integration points.
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method
 
int loopSize
 local number oe methods to process
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities.
 
int rAnk
 processor rank
 
int sIze
 number of processors in communicator
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities
 
const ProblemproblemPtr
 raw pointer to problem
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context.
 
KSP ksp
 KSP solver.
 
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 dx
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver
 
Vec & snes_x
 state vector
 
Vec & snes_dx
 solution update
 
Vec & snes_f
 residual
 
Mat & snes_A
 jacobian matrix
 
Mat & snes_B
 preconditioner of jacobian matrix
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver)
 
PetscReal ts_aa
 shift for U_tt shift for U_tt
 
PetscReal ts_t
 time
 
PetscReal ts_dt
 time step size
 
Vec & ts_u
 state vector
 
Vec & ts_u_t
 time derivative of state vector
 
Vec & ts_u_tt
 second time derivative of state vector
 
Vec & ts_F
 residual vector
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 
Tao tao
 tao solver
 
Vec & tao_x
 
Vec & tao_f
 state vector
 
Mat & tao_A
 gradient vector
 
Mat & tao_B
 hessian matrix
 

Protected Member Functions

MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 get sense (orientation) of entity
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 Get the entity data order.
 
template<EntityType type>
MoFEMErrorCode getEntitySense (EntitiesFieldData &data) const
 Get the entity sense (orientation)
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const
 Get the entity data order for given space.
 
MoFEMErrorCode getFaceNodes (EntitiesFieldData &data) const
 Get nodes on faces.
 
MoFEMErrorCode getSpacesAndBaseOnEntities (EntitiesFieldData &data) const
 Get field approximation space and base on entities.
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions.
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement ()
 Calculate base functions.
 
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement ()
 Calculate Bernstein-Bezier base.
 
MoFEMErrorCode createDataOnElement (EntityType type)
 Create a entity data on element object.
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators.
 
Indices
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
 
MoFEMErrorCode getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get row node indices from FENumeredDofEntity_multiIndex
 
MoFEMErrorCode getColNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get col node indices from FENumeredDofEntity_multiIndex
 
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 getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices
 
MoFEMErrorCode getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices
 
Data
MoFEMErrorCode getBitRefLevelOnData ()
 
MoFEMErrorCode getNodesFieldData (EntitiesFieldData &data, const int bit_number) const
 Get data on nodes.
 
MoFEMErrorCode getEntityFieldData (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
template<typename EXTRACTOR >
MoFEMErrorCode getNoFieldEntityFieldData (EntitiesFieldData &data, const int bit_number, EXTRACTOR &&extractor) const
 Get field data on entities where no field is defined.
 
Data form NumeredDofEntity_multiIndex
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
 
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
 
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
 
Deprecated (do not use)
virtual int getRule (int order)
 
virtual MoFEMErrorCode setGaussPts (int order)
 

Protected Attributes

const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEdataOnElement
 Entity data on element entity rows fields.
 
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields.
 
EntitiesFieldDatadataNoField
 
EntitiesFieldDatadataH1
 
EntitiesFieldDatadataHcurl
 
EntitiesFieldDatadataHdiv
 
EntitiesFieldDatadataL2
 
boost::ptr_deque< UserDataOperatoropPtrVector
 Vector of finite element users data operators.
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity.
 
MatrixDouble coordsAtGaussPts
 coordinated at gauss points
 
double elementMeasure
 

Private Member Functions

MoFEMErrorCode setSideFEPtr (const ForcesAndSourcesCore *side_fe_ptr)
 Set the pointer to face element on the side.
 
MoFEMErrorCode setRefineFEPtr (const ForcesAndSourcesCore *refine_fe_ptr)
 Set the pointer to face element refined.
 

Private Attributes

boost::shared_ptr< BaseFunctionelementPolynomialBasePtr
 Pointer to entity polynomial base.
 
boost::shared_ptr< BaseFunctionuserPolynomialBasePtr
 Pointer to user polynomial base.
 
ForcesAndSourcesCoresidePtrFE
 Element to integrate on the sides.
 
ForcesAndSourcesCorerefinePtrFE
 Element to integrate parent or child.
 

Friends

class UserDataOperator
 
class VolumeElementForcesAndSourcesCoreOnSide
 
class FaceElementForcesAndSourcesCoreOnSide
 
class FaceElementForcesAndSourcesCoreOnChildParent
 
class EdgeElementForcesAndSourcesCoreOnChildParent
 
class VolumeElementForcesAndSourcesCoreOnContactPrismSide
 
template<int DIM>
struct OpCopyGeomDataToE
 
template<typename E >
struct OpBrokenLoopSide
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- 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 CtxSetDX = PetscData::Switches(CTX_SET_DX)
 
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

Member Typedef Documentation

◆ GaussHookFun

boost::function<MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> MoFEM::ForcesAndSourcesCore::GaussHookFun

Definition at line 33 of file ForcesAndSourcesCore.hpp.

◆ RuleHookFun

boost::function<int(int order_row, int order_col, int order_data)> MoFEM::ForcesAndSourcesCore::RuleHookFun

Definition at line 28 of file ForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ ForcesAndSourcesCore()

MoFEM::ForcesAndSourcesCore::ForcesAndSourcesCore ( Interface & m_field)

Definition at line 40 of file ForcesAndSourcesCore.cpp.

41 :
42
43 mField(m_field), getRuleHook(0), setRuleHook(0),
45
46 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // NOSPACE,
47 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // NOFIELD
48 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // H1
49 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // HCURL
50 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // HDIV
51 boost::make_shared<EntitiesFieldData>(MBENTITYSET) // L2
52
53 },
55
56 nullptr,
57 boost::make_shared<DerivedEntitiesFieldData>(
58 dataOnElement[NOFIELD]), // NOFIELD
59 boost::make_shared<DerivedEntitiesFieldData>(dataOnElement[H1]), // H1
60 boost::make_shared<DerivedEntitiesFieldData>(
61 dataOnElement[HCURL]), // HCURL
62 boost::make_shared<DerivedEntitiesFieldData>(
63 dataOnElement[HDIV]), // HDIV
64 boost::make_shared<DerivedEntitiesFieldData>(dataOnElement[L2]) // L2
65
66 },
70 lastEvaluatedElementEntityType(MBMAXTYPE), sidePtrFE(nullptr),
71 refinePtrFE(nullptr) {
72
73 dataOnElement[NOSPACE]->dataOnEntities[MBENTITYSET].push_back(
74 new EntitiesFieldData::EntData());
75
76 dataOnElement[NOFIELD]->dataOnEntities[MBENTITYSET].push_back(
77 new EntitiesFieldData::EntData());
78 derivedDataOnElement[NOFIELD]->dataOnEntities[MBENTITYSET].push_back(
79 new EntitiesFieldData::EntData());
80}
@ L2
field with C-1 continuity
Definition definitions.h:88
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition definitions.h:84
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
GaussHookFun setRuleHook
Set function to calculate integration rule.
ForcesAndSourcesCore * refinePtrFE
Element to integrate parent or child.
RuleHookFun getRuleHook
Hook to get rule.

Member Function Documentation

◆ calBernsteinBezierBaseFunctionsOnElement()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement ( )
protected

Calculate Bernstein-Bezier base.

Returns
MoFEMErrorCode

Definition at line 1131 of file ForcesAndSourcesCore.cpp.

1131 {
1133
1134 const auto ele_type = numeredEntFiniteElementPtr->getEntType();
1135
1136 auto get_nodal_base_data = [&](EntitiesFieldData &data, auto field_ptr) {
1138 auto &space = data.dataOnEntities[MBVERTEX][0].getSpace();
1139 auto &base = data.dataOnEntities[MBVERTEX][0].getBase();
1140 auto &bb_node_order = data.dataOnEntities[MBVERTEX][0].getBBNodeOrder();
1141
1142 auto &field_ents = getDataFieldEnts();
1143 auto bit_number = field_ptr->getBitNumber();
1144 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
1145 bit_number, get_id_for_min_type<MBVERTEX>());
1146 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
1147 cmp_uid_lo);
1148 if (lo != field_ents.end()) {
1149 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
1150 bit_number, get_id_for_max_type<MBVERTEX>());
1151 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
1152 if (lo != hi) {
1153
1154 for (auto it = lo; it != hi; ++it)
1155 if (auto first_e = it->lock()) {
1156 space = first_e->getSpace();
1157 base = first_e->getApproxBase();
1158 const int num_nodes = getNumberOfNodes();
1159 bb_node_order.resize(num_nodes, false);
1160 bb_node_order.clear();
1161
1162 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
1163
1164 for (; it != hi; ++it) {
1165 if (auto e = it->lock()) {
1166 const auto &sn = e->getSideNumberPtr();
1167 const int side_number = sn->side_number;
1168 const int brother_side_number = sn->brother_side_number;
1169 if (brother_side_number != -1)
1170 brother_ents_vec.emplace_back(e);
1171 bb_node_order[side_number] = e->getMaxOrder();
1172 }
1173 }
1174
1175 for (auto &it : brother_ents_vec) {
1176 if (const auto e = it.lock()) {
1177 const auto &sn = e->getSideNumberPtr();
1178 const int side_number = sn->side_number;
1179 const int brother_side_number = sn->brother_side_number;
1180 bb_node_order[brother_side_number] = bb_node_order[side_number];
1181 }
1182 }
1183
1184 break;
1185 }
1186 }
1187 }
1188
1190 };
1191
1192 auto get_entity_base_data = [&](EntitiesFieldData &data, auto field_ptr,
1193 const EntityType type_lo,
1194 const EntityType type_hi) {
1196 for (EntityType t = MBEDGE; t != MBPOLYHEDRON; ++t) {
1197 for (auto &dat : data.dataOnEntities[t]) {
1198 dat.getOrder() = 0;
1199 dat.getBase() = NOBASE;
1200 dat.getSpace() = NOSPACE;
1201 dat.getFieldData().resize(0, false);
1202 dat.getFieldDofs().resize(0, false);
1203 }
1204 }
1205
1206 auto &field_ents = getDataFieldEnts();
1207 auto bit_number = field_ptr->getBitNumber();
1208 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
1209 bit_number, get_id_for_min_type(type_lo));
1210 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
1211 cmp_uid_lo);
1212 if (lo != field_ents.end()) {
1213 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
1214 bit_number, get_id_for_max_type(type_hi));
1215 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
1216
1217 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
1218 for (; lo != hi; ++lo) {
1219 if (auto e = lo->lock()) {
1220 if (auto cache = e->entityCacheDataDofs.lock()) {
1221 if (cache->loHi[0] != cache->loHi[1]) {
1222 if (const auto side = e->getSideNumberPtr()->side_number;
1223 side >= 0) {
1224 const EntityType type = e->getEntType();
1225 auto &dat = data.dataOnEntities[type][side];
1226 const int brother_side =
1227 e->getSideNumberPtr()->brother_side_number;
1228 if (brother_side != -1)
1229 brother_ents_vec.emplace_back(e);
1230 dat.getBase() = e->getApproxBase();
1231 dat.getSpace() = e->getSpace();
1232 const auto ent_order = e->getMaxOrder();
1233 dat.getOrder() =
1234 dat.getOrder() > ent_order ? dat.getOrder() : ent_order;
1235 }
1236 }
1237 }
1238 }
1239 }
1240
1241 for (auto &ent_ptr : brother_ents_vec) {
1242 if (auto e = ent_ptr.lock()) {
1243 const EntityType type = e->getEntType();
1244 const int side = e->getSideNumberPtr()->side_number;
1245 const int brother_side = e->getSideNumberPtr()->brother_side_number;
1246 auto &dat = data.dataOnEntities[type][side];
1247 auto &dat_brother = data.dataOnEntities[type][brother_side];
1248 dat_brother.getBase() = dat.getBase();
1249 dat_brother.getSpace() = dat.getSpace();
1250 dat_brother.getOrder() = dat.getOrder();
1251 }
1252 }
1253 }
1255 };
1256
1257 for (auto &e : getDataFieldEnts()) {
1258 if (auto ent_data_ptr = e.lock()) {
1259 if (ent_data_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
1260 auto space = ent_data_ptr->getSpace();
1261 for (EntityType t = MBVERTEX; t != MBPOLYHEDRON; ++t) {
1262 for (auto &dat : (*dataOnElement[space]).dataOnEntities[t]) {
1263 for (auto &ptr : dat.getBBAlphaIndicesByOrderArray())
1264 ptr.reset();
1265 for (auto &ptr : dat.getBBNByOrderArray())
1266 ptr.reset();
1267 for (auto &ptr : dat.getBBDiffNByOrderArray())
1268 ptr.reset();
1269 }
1270 }
1271 }
1272 }
1273 }
1274
1275 auto check_space = [&](const auto space) {
1276 switch (space) {
1277 case H1:
1278 for (auto t = MBVERTEX; t <= ele_type; ++t) {
1279 if (dataOnElement[H1]->spacesOnEntities[t].test(H1))
1280 return true;
1281 }
1282 return false;
1283 case HCURL:
1284 for (auto t = MBEDGE; t <= ele_type; ++t) {
1285 if (dataOnElement[HCURL]->spacesOnEntities[t].test(HCURL))
1286 return true;
1287 }
1288 return false;
1289 case HDIV:
1290 for (auto t = MBTRI; t <= ele_type; ++t) {
1291 if (dataOnElement[HDIV]->spacesOnEntities[t].test(HDIV))
1292 return true;
1293 }
1294 return false;
1295 case L2:
1296 return dataOnElement[L2]->spacesOnEntities[ele_type].test(L2);
1297 break;
1298 default:
1299 THROW_MESSAGE("Not implemented");
1300 }
1301 };
1302
1303 std::set<string> fields_list;
1304 for (auto &e : getDataFieldEnts()) {
1305 if (auto ent_data_ptr = e.lock()) {
1306 if (ent_data_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
1307 auto field_name = ent_data_ptr->getName();
1308 if (fields_list.find(field_name) == fields_list.end()) {
1309 auto field_ptr = ent_data_ptr->getFieldRawPtr();
1310 auto space = ent_data_ptr->getSpace();
1311 CHKERR get_nodal_base_data(*dataOnElement[space], field_ptr);
1312 CHKERR get_entity_base_data(*dataOnElement[space], field_ptr, MBEDGE,
1313 ele_type);
1314 if (check_space(space)) {
1315#ifndef NDEBUG
1316 if (ent_data_ptr->getContinuity() != CONTINUOUS)
1317 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1318 "Broken space not implemented");
1319#endif // NDEBUG
1321 -> getValue(gaussPts,
1322 boost::make_shared<EntPolynomialBaseCtx>(
1323 *dataOnElement[space], field_name,
1324 static_cast<FieldSpace>(space), CONTINUOUS,
1326 fields_list.insert(field_name);
1327 }
1328 }
1329 }
1330 }
1331 }
1332
1334};
@ NOBASE
Definition definitions.h:59
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition definitions.h:64
FieldSpace
approximation spaces
Definition definitions.h:82
@ CONTINUOUS
Regular field.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
static auto cmp_uid_lo(const boost::weak_ptr< FieldEntity > &a, const UId &b)
EntityHandle get_id_for_max_type()
static auto cmp_uid_hi(const UId &b, const boost::weak_ptr< FieldEntity > &a)
EntityHandle get_id_for_min_type()
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
const FieldEntity_vector_view & getDataFieldEnts() const
auto getNumberOfNodes() const
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
MatrixDouble gaussPts
Matrix of integration points.

◆ calHierarchicalBaseFunctionsOnElement() [1/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calHierarchicalBaseFunctionsOnElement ( )
protected

Calculate base functions.

Returns
Error code

Use the some node base. Node base is usually used for construction other bases.

Definition at line 1111 of file ForcesAndSourcesCore.cpp.

1111 {
1113 /// Use the some node base. Node base is usually used for construction other
1114 /// bases.
1115 for (int space = HCURL; space != LASTSPACE; ++space) {
1116 dataOnElement[space]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
1117 dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
1118 dataOnElement[space]->dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(
1119 NOBASE) =
1120 dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(
1121 NOBASE);
1122 }
1123 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
1125 static_cast<FieldApproximationBase>(b));
1126 }
1128}
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.

◆ calHierarchicalBaseFunctionsOnElement() [2/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calHierarchicalBaseFunctionsOnElement ( const FieldApproximationBase b)
protected

Calculate base functions.

Returns
Error code

Definition at line 1048 of file ForcesAndSourcesCore.cpp.

1049 {
1051 if (dataOnElement[H1]->bAse.test(b)) {
1052
1053 switch (static_cast<FieldApproximationBase>(b)) {
1054 case NOBASE:
1056 // BERNSTEIN_BEZIER_BASE is not hierarchical base
1058 }
1059 default:
1060 break;
1061 }
1062
1063 auto get_elem_base = [&](auto base) {
1064 if (b != USER_BASE) {
1068 "Functions generating approximation base not defined");
1069 return getElementPolynomialBase();
1070 } else {
1071 if (!getUserPolynomialBase())
1074 "Functions generating user approximation base not defined");
1075 return getUserPolynomialBase();
1076 }
1077 };
1078
1079 auto get_ctx = [&](auto space, auto base, auto continuity) {
1080 return boost::make_shared<EntPolynomialBaseCtx>(
1081 *dataOnElement[space], static_cast<FieldSpace>(space), continuity,
1082 static_cast<FieldApproximationBase>(base), NOBASE);
1083 };
1084
1085 auto get_value = [&](auto elem_base, auto &gaussPts, auto &&ctx) {
1086 return elem_base->getValue(gaussPts, ctx);
1087 };
1088
1089 for (int space = H1; space != LASTSPACE; ++space) {
1090
1091#ifndef NDEBUG
1092 if (dataOnElement[H1]->basesOnSpaces[space].test(b) &&
1093 dataOnElement[H1]->brokenBasesOnSpaces[space].test(b)) {
1094 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1095 "Discontinuous and continuous bases on the same space");
1096 }
1097#endif // NDEBUG
1098
1099 if (dataOnElement[H1]->basesOnSpaces[space].test(b))
1100 CHKERR get_value(get_elem_base(b), gaussPts,
1101 get_ctx(space, b, CONTINUOUS));
1102 else if (dataOnElement[H1]->brokenBasesOnSpaces[space].test(b))
1103 CHKERR get_value(get_elem_base(b), gaussPts,
1104 get_ctx(space, b, DISCONTINUOUS));
1105 }
1106
1107 }
1109}
@ USER_BASE
user implemented approximation base
Definition definitions.h:68
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ DISCONTINUOUS
Broken continuity (No effect on L2 space)
static MoFEMErrorCode get_value(MatrixDouble &pts_x, MatrixDouble &pts_t, TYPE *ctx)
auto & getUserPolynomialBase()
Get the User Polynomial Base object.

◆ createDataOnElement()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::createDataOnElement ( EntityType type)
protected

Create a entity data on element object.

Returns
MoFEMErrorCode

Definition at line 1336 of file ForcesAndSourcesCore.cpp.

1336 {
1338
1339 // Data on elements for proper spaces
1340 for (int space = H1; space != LASTSPACE; ++space) {
1341 dataOnElement[space]->setElementType(type);
1342 derivedDataOnElement[space]->setElementType(type);
1343 }
1344
1346}

◆ getBitRefLevelOnData()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getBitRefLevelOnData ( )
protected

Get bit ref level in entities, and set it to data

Definition at line 549 of file ForcesAndSourcesCore.cpp.

549 {
551
552 for (auto &data : dataOnElement) {
553 if (data) {
554 for (auto &dat : data->dataOnEntities) {
555 for (auto &ent_dat : dat) {
556 ent_dat.getEntDataBitRefLevel().clear();
557 }
558 }
559 }
560 }
561
562 auto &field_ents = getDataFieldEnts();
563 for (auto it : field_ents) {
564 if (auto e = it.lock()) {
565 const FieldSpace space = e->getSpace();
566 if (space > NOFIELD) {
567 const EntityType type = e->getEntType();
568 const signed char side = e->getSideNumberPtr()->side_number;
569 if (side >= 0) {
570 if (auto &data = dataOnElement[space]) {
571 if (type == MBVERTEX) {
572 auto &dat = data->dataOnEntities[type][0];
573 dat.getEntDataBitRefLevel().resize(getNumberOfNodes(), false);
574 dat.getEntDataBitRefLevel()[side] = e->getBitRefLevel();
575 } else {
576 auto &dat = data->dataOnEntities[type][side];
577 dat.getEntDataBitRefLevel().resize(1,
578 false); // one entity per sie
579 dat.getEntDataBitRefLevel()[0] = e->getBitRefLevel();
580 }
581 }
582 }
583 }
584 }
585 }
586
588};

◆ getColNodesIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getColNodesIndices ( EntitiesFieldData & data,
const int bit_number ) const
protected

get col node indices from FENumeredDofEntity_multiIndex

Definition at line 328 of file ForcesAndSourcesCore.cpp.

329 {
330
331 struct Extractor {
332 boost::weak_ptr<EntityCacheNumeredDofs>
333 operator()(boost::shared_ptr<FieldEntity> &e) {
334 return e->entityCacheColDofs;
335 }
336 };
337
338 return getNodesIndices(bit_number, getColFieldEnts(),
339 data.dataOnEntities[MBVERTEX][0].getIndices(),
340 data.dataOnEntities[MBVERTEX][0].getLocalIndices(),
341 Extractor());
342}
auto & getColFieldEnts() const
MoFEMErrorCode getNodesIndices(const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
get node indices

◆ getDataOnElementBySpaceArray()

auto & MoFEM::ForcesAndSourcesCore::getDataOnElementBySpaceArray ( )
inline

Get data on entities and space.

Entities data are stored by space, by entity type, and entity side.

Returns
std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>

Definition at line 156 of file ForcesAndSourcesCore.hpp.

156{ return dataOnElement; }

◆ getDerivedDataOnElementBySpaceArray()

auto & MoFEM::ForcesAndSourcesCore::getDerivedDataOnElementBySpaceArray ( )
inline

Get derived data on entities and space.

Entities data are stored by space, by entity type, and entity side. Derived data is used to store data on columns, so it shares information about shape functions wih rows.

Returns
std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>

Definition at line 167 of file ForcesAndSourcesCore.hpp.

167{ return derivedDataOnElement; }

◆ getElementPolynomialBase()

auto & MoFEM::ForcesAndSourcesCore::getElementPolynomialBase ( )
inline

Get the Entity Polynomial Base object.

Returns
boost::shared_ptr<BaseFunction>&&

Definition at line 90 of file ForcesAndSourcesCore.hpp.

boost::shared_ptr< BaseFunction > elementPolynomialBasePtr
Pointer to entity polynomial base.

◆ getEntData()

auto & MoFEM::ForcesAndSourcesCore::getEntData ( const FieldSpace space,
const EntityType type,
const int side )
inline

Get the entity data.

Parameters
space
type
side
Returns
EntitiesFieldData::EntData&

Definition at line 144 of file ForcesAndSourcesCore.hpp.

145 {
146 return dataOnElement[space]->dataOnEntities[type][side];
147 }

◆ getEntityColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityColIndices ( EntitiesFieldData & data,
const int bit_number,
const EntityType type_lo = MBVERTEX,
const EntityType type_hi = MBPOLYHEDRON ) const
protected

Definition at line 427 of file ForcesAndSourcesCore.cpp.

429 {
430
431 struct Extractor {
432 boost::weak_ptr<EntityCacheNumeredDofs>
433 operator()(boost::shared_ptr<FieldEntity> &e) {
434 return e->entityCacheColDofs;
435 }
436 };
437
438 return getEntityIndices(data, bit_number, getColFieldEnts(), type_lo, type_hi,
439 Extractor());
440}
MoFEMErrorCode getEntityIndices(EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const

◆ getEntityDataOrder() [1/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityDataOrder ( const EntityType type,
const FieldSpace space,
boost::ptr_vector< EntitiesFieldData::EntData > & data ) const
protected

Get the entity data order.

Parameters
type
space
data
Returns
MoFEMErrorCode

Definition at line 139 of file ForcesAndSourcesCore.cpp.

141 {
143
144 auto set_order = [&]() {
146 auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable();
147
148 for (unsigned int s = 0; s != data.size(); ++s)
149 data[s].getOrder() = 0;
150
151 const FieldEntity_vector_view &data_field_ent = getDataFieldEnts();
152
153 for (auto r = fieldsPtr->get<BitFieldId_space_mi_tag>().equal_range(space);
154 r.first != r.second; ++r.first) {
155
156 const auto field_bit_number = (*r.first)->getBitNumber();
157 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
158 field_bit_number, get_id_for_min_type(type));
159 auto lo = std::lower_bound(data_field_ent.begin(), data_field_ent.end(),
160 lo_uid, cmp_uid_lo);
161 if (lo != data_field_ent.end()) {
162 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
163 field_bit_number, get_id_for_max_type(type));
164 auto hi =
165 std::upper_bound(lo, data_field_ent.end(), hi_uid, cmp_uid_hi);
166 for (; lo != hi; ++lo) {
167
168 if (auto ptr = lo->lock()) {
169
170 auto &e = *ptr;
171 auto sit = side_table.find(e.getEnt());
172 if (sit != side_table.end()) {
173 auto &side = *sit;
174 if (const auto side_number = side->side_number;
175 side_number >= 0) {
176 ApproximationOrder ent_order = e.getMaxOrder();
177 auto &dat = data[side_number];
178 dat.getOrder() =
179 dat.getOrder() > ent_order ? dat.getOrder() : ent_order;
180 }
181 } else
182 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
183 "Entity on side of the element not found");
184 }
185 }
186 }
187 }
188
190 };
191
192 auto set_order_on_brother = [&]() {
194 auto &side_table =
195 numeredEntFiniteElementPtr->getSideNumberTable().get<0>();
196 auto sit = side_table.lower_bound(get_id_for_min_type(type));
197 if (sit != side_table.end()) {
198 auto hi_sit = side_table.upper_bound(get_id_for_max_type(type));
199 for (; sit != hi_sit; ++sit) {
200 const int brother_side_number = (*sit)->brother_side_number;
201 if (brother_side_number != -1) {
202 const int side_number = (*sit)->side_number;
203 data[brother_side_number].getOrder() = data[side_number].getOrder();
204 }
205 }
206 }
208 };
209
210 CHKERR set_order();
211 CHKERR set_order_on_brother();
212
214}
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
int ApproximationOrder
Approximation on the entity.
Definition Types.hpp:26
std::vector< boost::weak_ptr< FieldEntity > > FieldEntity_vector_view
int r
Definition sdf.py:8
const Field_multiIndex * fieldsPtr
raw pointer to fields container

◆ getEntityDataOrder() [2/2]

template<EntityType type>
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityDataOrder ( EntitiesFieldData & data,
const FieldSpace space ) const
inlineprotected

Get the entity data order for given space.

Template Parameters
type
Parameters
data
space
Returns
MoFEMErrorCode

Definition at line 213 of file ForcesAndSourcesCore.hpp.

214 {
215 return getEntityDataOrder(type, space, data.dataOnEntities[type]);
216 }
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.

◆ getEntityFieldData()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityFieldData ( EntitiesFieldData & data,
const int bit_number,
const EntityType type_lo = MBVERTEX,
const EntityType type_hi = MBPOLYHEDRON ) const
protected

Definition at line 719 of file ForcesAndSourcesCore.cpp.

721 {
723 for (EntityType t = type_lo; t != type_hi; ++t) {
724 for (auto &dat : data.dataOnEntities[t]) {
725 dat.getOrder() = 0;
726 dat.getBase() = NOBASE;
727 dat.getSpace() = NOSPACE;
728 dat.getFieldData().resize(0, false);
729 dat.getFieldDofs().resize(0, false);
730 dat.getFieldEntities().resize(0, false);
731 }
732 }
733
734 auto &field_ents = getDataFieldEnts();
735 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
736 bit_number, get_id_for_min_type(type_lo));
737 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
738 cmp_uid_lo);
739 if (lo != field_ents.end()) {
740 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
741 bit_number, get_id_for_max_type(type_hi));
742 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
743 if (lo != hi) {
744
745 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
746
747 for (auto it = lo; it != hi; ++it)
748 if (auto e = it->lock()) {
749 auto side_ptr = e->getSideNumberPtr();
750 if (const auto side = side_ptr->side_number; side >= 0) {
751 const EntityType type = e->getEntType();
752 auto &dat = data.dataOnEntities[type][side];
753 auto &ent_field = dat.getFieldEntities();
754 auto &ent_field_dofs = dat.getFieldDofs();
755 auto &ent_field_data = dat.getFieldData();
756
757 const int brother_side = side_ptr->brother_side_number;
758 if (brother_side != -1)
759 brother_ents_vec.emplace_back(e);
760
761 dat.getBase() = e->getApproxBase();
762 dat.getSpace() = e->getSpace();
763 const int ent_order = e->getMaxOrder();
764 dat.getOrder() =
765 dat.getOrder() > ent_order ? dat.getOrder() : ent_order;
766
767 auto ent_data = e->getEntFieldData();
768 ent_field_data.resize(ent_data.size(), false);
769 noalias(ent_field_data) = ent_data;
770 ent_field_dofs.resize(ent_data.size(), false);
771 std::fill(ent_field_dofs.begin(), ent_field_dofs.end(), nullptr);
772 ent_field.resize(1, false);
773 ent_field[0] = e.get();
774 if (auto cache = e->entityCacheDataDofs.lock()) {
775 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
776 ent_field_dofs[(*dit)->getEntDofIdx()] =
777 reinterpret_cast<FEDofEntity *>((*dit).get());
778 }
779 }
780 }
781 }
782
783 for (auto &it : brother_ents_vec) {
784 if (const auto e = it.lock()) {
785 const EntityType type = e->getEntType();
786 const int side = e->getSideNumberPtr()->side_number;
787 const int brother_side = e->getSideNumberPtr()->brother_side_number;
788 auto &dat = data.dataOnEntities[type][side];
789 auto &dat_brother = data.dataOnEntities[type][brother_side];
790 dat_brother.getBase() = dat.getBase();
791 dat_brother.getSpace() = dat.getSpace();
792 dat_brother.getOrder() = dat.getOrder();
793 dat_brother.getFieldData() = dat.getFieldData();
794 dat_brother.getFieldDofs() = dat.getFieldDofs();
795 dat_brother.getFieldEntities() = dat.getFieldEntities();
796 }
797 }
798 }
799 }
800
802}

◆ getEntityIndices()

template<typename EXTRACTOR >
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityIndices ( EntitiesFieldData & data,
const int bit_number,
FieldEntity_vector_view & ents_field,
const EntityType type_lo,
const EntityType type_hi,
EXTRACTOR && extractor ) const
protected

Definition at line 345 of file ForcesAndSourcesCore.cpp.

348 {
350
351 for (EntityType t = type_lo; t != type_hi; ++t) {
352 for (auto &dat : data.dataOnEntities[t]) {
353 dat.getIndices().resize(0, false);
354 dat.getLocalIndices().resize(0, false);
355 }
356 }
357
358 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
359 bit_number, get_id_for_min_type(type_lo));
360 auto lo = std::lower_bound(ents_field.begin(), ents_field.end(), lo_uid,
361 cmp_uid_lo);
362 if (lo != ents_field.end()) {
363 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
364 bit_number, get_id_for_max_type(type_hi));
365 auto hi = std::upper_bound(lo, ents_field.end(), hi_uid, cmp_uid_hi);
366 if (lo != hi) {
367
368 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
369
370 for (auto it = lo; it != hi; ++it)
371 if (auto e = it->lock()) {
372
373 const EntityType type = e->getEntType();
374 auto side_ptr = e->getSideNumberPtr();
375 if (const auto side = side_ptr->side_number; side >= 0) {
376 const auto nb_dofs_on_ent = e->getNbDofsOnEnt();
377 const auto brother_side = side_ptr->brother_side_number;
378 auto &dat = data.dataOnEntities[type][side];
379 auto &ent_field_indices = dat.getIndices();
380 auto &ent_field_local_indices = dat.getLocalIndices();
381
382 ent_field_indices.resize(nb_dofs_on_ent, false);
383 ent_field_local_indices.resize(nb_dofs_on_ent, false);
384 std::fill(ent_field_indices.data().begin(),
385 ent_field_indices.data().end(), -1);
386 std::fill(ent_field_local_indices.data().begin(),
387 ent_field_local_indices.data().end(), -1);
388
389 if (auto cache = extractor(e).lock()) {
390 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
391 const int idx = (*dit)->getEntDofIdx();
392 ent_field_indices[idx] = (*dit)->getPetscGlobalDofIdx();
393 ent_field_local_indices[idx] = (*dit)->getPetscLocalDofIdx();
394 }
395 }
396
397 if (brother_side != -1) {
398 auto &dat_brother = data.dataOnEntities[type][brother_side];
399 dat_brother.getIndices().resize(nb_dofs_on_ent, false);
400 dat_brother.getLocalIndices().resize(nb_dofs_on_ent, false);
401 noalias(dat_brother.getIndices()) = dat.getIndices();
402 noalias(dat_brother.getLocalIndices()) = dat.getLocalIndices();
403 }
404 }
405 }
406 }
407 }
408
410}

◆ getEntityRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityRowIndices ( EntitiesFieldData & data,
const int bit_number,
const EntityType type_lo = MBVERTEX,
const EntityType type_hi = MBPOLYHEDRON ) const
protected

Definition at line 412 of file ForcesAndSourcesCore.cpp.

414 {
415
416 struct Extractor {
417 boost::weak_ptr<EntityCacheNumeredDofs>
418 operator()(boost::shared_ptr<FieldEntity> &e) {
419 return e->entityCacheRowDofs;
420 }
421 };
422
423 return getEntityIndices(data, bit_number, getRowFieldEnts(), type_lo, type_hi,
424 Extractor());
425}
auto & getRowFieldEnts() const

◆ getEntitySense() [1/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntitySense ( const EntityType type,
boost::ptr_vector< EntitiesFieldData::EntData > & data ) const
protected

get sense (orientation) of entity

Parameters
typetype of entity
dataentity data
Returns
error code

Definition at line 84 of file ForcesAndSourcesCore.cpp.

86 {
88
89 auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable().get<0>();
90 auto sit = side_table.lower_bound(get_id_for_min_type(type));
91 if (sit != side_table.end()) {
92 auto hi_sit = side_table.upper_bound(get_id_for_max_type(type));
93 for (; sit != hi_sit; ++sit) {
94 if (const auto side_number = (*sit)->side_number; side_number >= 0) {
95 const int brother_side_number = (*sit)->brother_side_number;
96 const int sense = (*sit)->sense;
97
98 data[side_number].getSense() = sense;
99 if (brother_side_number != -1)
100 data[brother_side_number].getSense() = sense;
101 }
102 }
103 }
105}

◆ getEntitySense() [2/2]

template<EntityType type>
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntitySense ( EntitiesFieldData & data) const
inlineprotected

Get the entity sense (orientation)

Template Parameters
type
Parameters
data
Returns
MoFEMErrorCode

Definition at line 200 of file ForcesAndSourcesCore.hpp.

200 {
201 return getEntitySense(type, data.dataOnEntities[type]);
202 }
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity

◆ getFaceNodes()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getFaceNodes ( EntitiesFieldData & data) const
protected

Get nodes on faces.

Definition at line 918 of file ForcesAndSourcesCore.cpp.

918 {
920 auto &face_nodes = data.facesNodes;
921 auto &face_nodes_order = data.facesNodesOrder;
922 auto &side_table = const_cast<SideNumber_multiIndex &>(
923 numeredEntFiniteElementPtr->getSideNumberTable());
924 const auto ent = numeredEntFiniteElementPtr->getEnt();
925 const auto type = numeredEntFiniteElementPtr->getEntType();
926 const auto nb_faces = CN::NumSubEntities(type, 2);
927 const EntityHandle *conn_ele;
928 int num_nodes_ele;
929 CHKERR mField.get_moab().get_connectivity(ent, conn_ele, num_nodes_ele, true);
930 auto side_ptr_it = side_table.get<1>().lower_bound(
931 boost::make_tuple(CN::TypeDimensionMap[2].first, 0));
932 auto hi_side_ptr_it = side_table.get<1>().upper_bound(boost::make_tuple(
933 CN::TypeDimensionMap[2].second, std::numeric_limits<signed char>::max()));
934
935 for (; side_ptr_it != hi_side_ptr_it; ++side_ptr_it) {
936 const auto side = (*side_ptr_it)->side_number;
937 if (side < 0)
938 continue;
939 const auto sense = (*side_ptr_it)->sense;
940 const auto offset = (*side_ptr_it)->offset;
941
942 EntityType face_type;
943 int nb_nodes_face;
944 auto face_indices =
945 CN::SubEntityVertexIndices(type, 2, side, face_type, nb_nodes_face);
946 face_nodes.resize(nb_faces, nb_nodes_face);
947 face_nodes_order.resize(nb_faces, nb_nodes_face);
948
949 if (sense == 1)
950 for (int n = 0; n != nb_nodes_face; ++n)
951 face_nodes_order(side, n) = (n + offset) % nb_nodes_face;
952 else
953 for (int n = 0; n != nb_nodes_face; ++n)
954 face_nodes_order(side, n) =
955 (nb_nodes_face - (n - offset) % nb_nodes_face) % nb_nodes_face;
956
957 for (int n = 0; n != nb_nodes_face; ++n)
958 face_nodes(side, n) = face_indices[face_nodes_order(side, n)];
959
960#ifndef NDEBUG
961 const auto face_ent = (*side_ptr_it)->ent;
962 auto check = [&]() {
964 const EntityHandle *conn_face;
965 // int nb_nodes_face;
966 CHKERR mField.get_moab().get_connectivity(face_ent, conn_face,
967 nb_nodes_face, true);
968 face_nodes.resize(nb_faces, nb_nodes_face);
969 for (int nn = 0; nn != nb_nodes_face; ++nn) {
970 if (face_nodes(side, nn) !=
971 std::distance(
972 conn_ele,
973 std::find(conn_ele, &conn_ele[num_nodes_ele], conn_face[nn]))) {
974 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
975 "Wrong face numeration");
976 }
977 }
979 };
980 CHKERR check();
981#endif
982 }
983
985}
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
const double n
refractive index of diffusive medium
virtual moab::Interface & get_moab()=0

◆ getMaxColOrder()

int MoFEM::ForcesAndSourcesCore::getMaxColOrder ( ) const

Get max order of approximation for field in columns.

Definition at line 135 of file ForcesAndSourcesCore.cpp.

135 {
137}
static int getMaxOrder(const ENTMULTIINDEX &multi_index)

◆ getMaxDataOrder()

int MoFEM::ForcesAndSourcesCore::getMaxDataOrder ( ) const

Get max order of approximation for data fields.

getMeasure getMaxDataOrder () return maximal order on entities, for all data on the element. So for example if finite element is triangle, and triangle base function have order 4 and on edges base function have order 2, this function return 4.

If finite element has for example 2 or more approximated fields, for example Pressure (order 3) and displacement field (order 5), this function returns 5.

Definition at line 120 of file ForcesAndSourcesCore.cpp.

120 {
121 int max_order = 0;
122 for (auto e : getDataFieldEnts()) {
123 if (auto ptr = e.lock()) {
124 const int order = ptr->getMaxOrder();
125 max_order = (max_order < order) ? order : max_order;
126 }
127 }
128 return max_order;
129}
constexpr int order

◆ getMaxRowOrder()

int MoFEM::ForcesAndSourcesCore::getMaxRowOrder ( ) const

Get max order of approximation for field in rows.

Definition at line 131 of file ForcesAndSourcesCore.cpp.

131 {
133}

◆ getNodesFieldData()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNodesFieldData ( EntitiesFieldData & data,
const int bit_number ) const
protected

Get data on nodes.

Parameters
dataData structure
field_nameField name
Returns
Error code

Definition at line 591 of file ForcesAndSourcesCore.cpp.

592 {
593
594 auto get_nodes_field_data = [&](VectorDouble &nodes_data,
595 VectorFieldEntities &field_entities,
596 VectorDofs &nodes_dofs, FieldSpace &space,
598 VectorInt &bb_node_order) {
600
601 nodes_data.resize(0, false);
602 nodes_dofs.resize(0, false);
603 field_entities.resize(0, false);
604
605 auto field_it =
606 fieldsPtr->get<BitFieldId_mi_tag>().find(BitFEId().set(bit_number - 1));
607 if (field_it != fieldsPtr->get<BitFieldId_mi_tag>().end()) {
608
609#ifndef NDEBUG
610 if ((*field_it)->getBitNumber() != bit_number)
611 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong bit number");
612#endif
613 const int nb_dofs_on_vert = (*field_it)->getNbOfCoeffs();
614 space = (*field_it)->getSpace();
615 base = (*field_it)->getApproxBase();
616
617 auto &field_ents = getDataFieldEnts();
618 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
619 bit_number, get_id_for_min_type<MBVERTEX>());
620 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
621 cmp_uid_lo);
622 if (lo != field_ents.end()) {
623 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
624 bit_number, get_id_for_max_type<MBVERTEX>());
625 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
626 if (lo != hi) {
627
628 int nb_dofs = 0;
629 for (auto it = lo; it != hi; ++it) {
630 if (auto e = it->lock()) {
631 if (auto cache = e->entityCacheDataDofs.lock()) {
632 if (cache->loHi[0] != cache->loHi[1]) {
633 nb_dofs += std::distance(cache->loHi[0], cache->loHi[1]);
634 break;
635 }
636 }
637 }
638 }
639
640 if (nb_dofs) {
641
642 const int num_nodes = getNumberOfNodes();
643 bb_node_order.resize(num_nodes, false);
644 bb_node_order.clear();
645 const int max_nb_dofs = nb_dofs_on_vert * num_nodes;
646 nodes_data.resize(max_nb_dofs, false);
647 nodes_dofs.resize(max_nb_dofs, false);
648 field_entities.resize(num_nodes, false);
649 std::fill(nodes_data.begin(), nodes_data.end(), 0);
650 std::fill(nodes_dofs.begin(), nodes_dofs.end(), nullptr);
651 std::fill(field_entities.begin(), field_entities.end(), nullptr);
652
653 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
654
655 for (auto it = lo; it != hi; ++it) {
656 if (auto e = it->lock()) {
657 const auto &sn = e->getSideNumberPtr();
658 // Some field entities on skeleton can have negative side
659 // number
660 if (const auto side_number = sn->side_number;
661 side_number >= 0) {
662 const int brother_side_number = sn->brother_side_number;
663
664 field_entities[side_number] = e.get();
665 if (brother_side_number != -1) {
666 brother_ents_vec.emplace_back(e);
667 field_entities[side_number] = field_entities[side_number];
668 }
669
670 bb_node_order[side_number] = e->getMaxOrder();
671 int pos = side_number * nb_dofs_on_vert;
672 auto ent_filed_data_vec = e->getEntFieldData();
673 if (auto cache = e->entityCacheDataDofs.lock()) {
674 for (auto dit = cache->loHi[0]; dit != cache->loHi[1];
675 ++dit) {
676 const auto dof_idx = (*dit)->getEntDofIdx();
677 nodes_data[pos + dof_idx] = ent_filed_data_vec[dof_idx];
678 nodes_dofs[pos + dof_idx] =
679 reinterpret_cast<FEDofEntity *>((*dit).get());
680 }
681 }
682 }
683 }
684 }
685
686 for (auto &it : brother_ents_vec) {
687 if (const auto e = it.lock()) {
688 const auto &sn = e->getSideNumberPtr();
689 const int side_number = sn->side_number;
690 const int brother_side_number = sn->brother_side_number;
691 bb_node_order[brother_side_number] = bb_node_order[side_number];
692 int pos = side_number * nb_dofs_on_vert;
693 int brother_pos = brother_side_number * nb_dofs_on_vert;
694 for (int ii = 0; ii != nb_dofs_on_vert; ++ii) {
695 nodes_data[brother_pos] = nodes_data[pos];
696 nodes_dofs[brother_pos] = nodes_dofs[pos];
697 ++pos;
698 ++brother_pos;
699 }
700 }
701 }
702 }
703 }
704 }
705 }
706
708 };
709
710 return get_nodes_field_data(
711 data.dataOnEntities[MBVERTEX][0].getFieldData(),
712 data.dataOnEntities[MBVERTEX][0].getFieldEntities(),
713 data.dataOnEntities[MBVERTEX][0].getFieldDofs(),
714 data.dataOnEntities[MBVERTEX][0].getSpace(),
715 data.dataOnEntities[MBVERTEX][0].getBase(),
716 data.dataOnEntities[MBVERTEX][0].getBBNodeOrder());
717}
std::bitset< BITFEID_SIZE > BitFEId
Finite element Id.
Definition Types.hpp:43
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
ublas::vector< FieldEntity *, FieldEntAllocator > VectorFieldEntities
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs

◆ getNodesIndices()

template<typename EXTRACTOR >
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNodesIndices ( const int bit_number,
FieldEntity_vector_view & ents_field,
VectorInt & nodes_indices,
VectorInt & local_nodes_indices,
EXTRACTOR && extractor ) const
protected

get node indices

Definition at line 219 of file ForcesAndSourcesCore.cpp.

222 {
224
225 auto field_it = fieldsPtr->get<BitFieldId_mi_tag>().find(
226 BitFieldId().set(bit_number - 1));
227 if (field_it != fieldsPtr->get<BitFieldId_mi_tag>().end()) {
228
229#ifndef NDEBUG
230 if ((*field_it)->getBitNumber() != bit_number)
231 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong bit number");
232#endif
233 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
234 bit_number, get_id_for_min_type<MBVERTEX>());
235 auto lo = std::lower_bound(ents_field.begin(), ents_field.end(), lo_uid,
236 cmp_uid_lo);
237 if (lo != ents_field.end()) {
238 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
239 bit_number, get_id_for_max_type<MBVERTEX>());
240 auto hi = std::upper_bound(lo, ents_field.end(), hi_uid, cmp_uid_hi);
241
242 const int num_nodes = getNumberOfNodes();
243 const int nb_dofs_on_vert = (*field_it)->getNbOfCoeffs();
244 const int max_nb_dofs = nb_dofs_on_vert * num_nodes;
245
246 int nb_dofs = 0;
247 for (auto it = lo; it != hi; ++it) {
248 if (auto e = it->lock()) {
249 if (auto cache = extractor(e).lock()) {
250 if (cache->loHi[0] != cache->loHi[1]) {
251 nb_dofs += std::distance(cache->loHi[0], cache->loHi[1]);
252 break;
253 }
254 }
255 }
256 }
257
258 if (nb_dofs) {
259 nodes_indices.resize(max_nb_dofs, false);
260 local_nodes_indices.resize(max_nb_dofs, false);
261 } else {
262 nodes_indices.resize(0, false);
263 local_nodes_indices.resize(0, false);
264 }
265
266 if (nb_dofs != max_nb_dofs) {
267 std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
268 std::fill(local_nodes_indices.begin(), local_nodes_indices.end(), -1);
269 }
270
271 for (auto it = lo; it != hi; ++it) {
272 if (auto e = it->lock()) {
273 auto side_ptr = e->getSideNumberPtr();
274 if (const auto side_number = side_ptr->side_number;
275 side_number >= 0) {
276 const auto brother_side_number = side_ptr->brother_side_number;
277 if (auto cache = extractor(e).lock()) {
278 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
279 auto &dof = **dit;
280 const int idx = dof.getPetscGlobalDofIdx();
281 const int local_idx = dof.getPetscLocalDofIdx();
282 const int pos =
283 side_number * nb_dofs_on_vert + dof.getDofCoeffIdx();
284 nodes_indices[pos] = idx;
285 local_nodes_indices[pos] = local_idx;
286 if (brother_side_number != -1) {
287 const int elem_idx = brother_side_number * nb_dofs_on_vert +
288 (*dit)->getDofCoeffIdx();
289 nodes_indices[elem_idx] = idx;
290 local_nodes_indices[elem_idx] = local_idx;
291 }
292 }
293 }
294 }
295 }
296 }
297 } else {
298 nodes_indices.resize(0, false);
299 local_nodes_indices.resize(0, false);
300 }
301
302 } else {
303 nodes_indices.resize(0, false);
304 local_nodes_indices.resize(0, false);
305 }
306
308}
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition Types.hpp:42

◆ getNoFieldColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldColIndices ( EntitiesFieldData & data,
const int bit_number ) const
protected

get col NoField indices

◆ getNoFieldEntityFieldData()

template<typename EXTRACTOR >
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldEntityFieldData ( EntitiesFieldData & data,
const int bit_number,
EXTRACTOR && extractor ) const
protected

Get field data on entities where no field is defined.

This is not standards, since op_data, dataOnElement[NOFIELD], has to be set for rows and columns for each NOFIELD, also size of sides is not determined.

Parameters
data

param bit_number

Returns
MoFEMErrorCode

Definition at line 806 of file ForcesAndSourcesCore.cpp.

808 {
810
811 // resize data on entities
812 std::array<FieldEntity_vector_view, MBMAXTYPE> field_ents_by_type;
813 for (auto it : getDataFieldEnts())
814 if (auto e = it.lock())
815 if (e->getBitNumber() == bit_number) {
816 field_ents_by_type[e->getEntType()].emplace_back(e);
817 }
818
819 auto set_entity_data = [&](auto &data_on_element) {
821
822 // resize data containers
823 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
824 auto &data = data_on_element.dataOnEntities[t];
825 data.resize(field_ents_by_type[t].size());
826 }
827
828 // reset data
829 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
830 auto &data = data_on_element.dataOnEntities[t];
831 for (auto &d : data) {
832 d.getEntDataBitRefLevel().clear();
833 d.getOrder() = 0;
834 d.getBase() = NOBASE;
835 d.getSpace() = NOSPACE;
836 d.getFieldData().resize(0, false);
837 d.getFieldDofs().resize(0, false);
838 d.getFieldEntities().resize(0, false);
839 d.getIndices().resize(0, false);
840 d.getLocalIndices().resize(0, false);
841 }
842 }
843
844 // set bit ref level on entities
845 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
846 int side = 0;
847 for (auto it : field_ents_by_type[t]) {
848 if (auto e = it.lock()) {
849 auto &data = data_on_element.dataOnEntities[t][side];
850 data.getEntDataBitRefLevel().resize(1, false);
851 data.getEntDataBitRefLevel()[0] = e->getBitRefLevel();
852 ++side;
853 }
854 }
855 }
856
857 // set field data on entities
858 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
859 int side = 0;
860 for (auto it : field_ents_by_type[t]) {
861 if (auto e = it.lock()) {
862 auto &data = data_on_element.dataOnEntities[t][side];
863 data.getBase() = e->getApproxBase();
864 data.getSpace() = e->getSpace();
865 data.getFieldData() = e->getEntFieldData();
866 data.getFieldEntities().resize(1, false);
867 data.getFieldEntities()[0] = e.get();
868 auto &ent_field_dofs = data.getFieldDofs();
869 ent_field_dofs.resize(data.getFieldData().size(), false);
870 std::fill(ent_field_dofs.begin(), ent_field_dofs.end(), nullptr);
871 if (auto cache = e->entityCacheDataDofs.lock()) {
872 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
873 ent_field_dofs[(*dit)->getEntDofIdx()] =
874 reinterpret_cast<FEDofEntity *>((*dit).get());
875 }
876 }
877 ++side;
878 }
879 }
880 }
881
882 // set indices on entities
883 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
884 int side = 0;
885 for (auto it : field_ents_by_type[t]) {
886 if (auto e = it.lock()) {
887 auto &data = data_on_element.dataOnEntities[t][side];
888 auto &indices = data.getIndices();
889 auto &local_indices = data.getLocalIndices();
890 indices.resize(data.getFieldData().size(), false);
891 local_indices.resize(data.getFieldData().size(), false);
892 std::fill(indices.begin(), indices.end(), -1);
893 std::fill(local_indices.begin(), local_indices.end(), -1);
894
895 if (auto cache = extractor(e).lock()) {
896 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
897 indices[(*dit)->getEntDofIdx()] = (*dit)->getPetscGlobalDofIdx();
898 local_indices[(*dit)->getEntDofIdx()] =
899 (*dit)->getPetscLocalDofIdx();
900 }
901 }
902
903 ++side;
904 }
905 }
906 }
908 };
909
910 CHKERR set_entity_data(data);
911
913}

◆ getNoFieldRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldRowIndices ( EntitiesFieldData & data,
const int bit_number ) const
protected

get col NoField indices

◆ getOpPtrVector()

boost::ptr_deque< UserDataOperator > & MoFEM::ForcesAndSourcesCore::getOpPtrVector ( )
inline

◆ getProblemNodesColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getProblemNodesColIndices ( const std::string & field_name,
VectorInt & nodes_indices ) const
protected

Definition at line 533 of file ForcesAndSourcesCore.cpp.

534 {
536 nodes_indices);
537}
const Problem * problemPtr
raw pointer to problem
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
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem

◆ getProblemNodesIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getProblemNodesIndices ( const std::string & field_name,
const NumeredDofEntity_multiIndex & dofs,
VectorInt & nodes_indices ) const
protected

get indices of nodal indices which are declared for problem but not this particular element

Definition at line 444 of file ForcesAndSourcesCore.cpp.

446 {
448
449 const Field *field_struture = mField.get_field_structure(field_name);
450 if (field_struture->getSpace() == H1) {
451
452 const int num_nodes = getNumberOfNodes();
453 nodes_indices.resize(field_struture->getNbOfCoeffs() * num_nodes, false);
454 std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
455
456 auto &side_table = const_cast<SideNumber_multiIndex &>(
457 numeredEntFiniteElementPtr->getSideNumberTable());
458 auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 0));
459 auto hi_siit =
460 side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 10000));
461
462 int nn = 0;
463 for (; siit != hi_siit; siit++, nn++) {
464 if (siit->get()->side_number >= 0) {
465 auto bit_number = mField.get_field_bit_number(field_name);
466 const EntityHandle ent = siit->get()->ent;
467 auto dit = dofs.get<Unique_mi_tag>().lower_bound(
469 auto hi_dit = dofs.get<Unique_mi_tag>().upper_bound(
471 for (; dit != hi_dit; dit++) {
472 nodes_indices[siit->get()->side_number * (*dit)->getNbOfCoeffs() +
473 (*dit)->getDofCoeffIdx()] =
474 (*dit)->getPetscGlobalDofIdx();
475 }
476 }
477 }
478 } else {
479 nodes_indices.resize(0, false);
480 }
481
483}
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
static UId getLoLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
static UId getHiLocalEntityBitNumber(const char bit_number, const EntityHandle ent)

◆ getProblemNodesRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getProblemNodesRowIndices ( const std::string & field_name,
VectorInt & nodes_indices ) const
protected

Definition at line 519 of file ForcesAndSourcesCore.cpp.

520 {
522 nodes_indices);
523}
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem

◆ getProblemTypeColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getProblemTypeColIndices ( const std::string & field_name,
EntityType type,
int side_number,
VectorInt & indices ) const
protected

Definition at line 540 of file ForcesAndSourcesCore.cpp.

542 {
544 type, side_number, indices);
545}
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

◆ getProblemTypeIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getProblemTypeIndices ( const std::string & field_name,
const NumeredDofEntity_multiIndex & dofs,
EntityType type,
int side_number,
VectorInt & indices ) const
protected

get indices by type (generic function) which are declared for problem but not this particular element

Definition at line 485 of file ForcesAndSourcesCore.cpp.

487 {
489
490 indices.resize(0);
491
492 auto &side_table = const_cast<SideNumber_multiIndex &>(
493 numeredEntFiniteElementPtr->getSideNumberTable());
494 auto siit =
495 side_table.get<1>().lower_bound(boost::make_tuple(type, side_number));
496 auto hi_siit =
497 side_table.get<1>().upper_bound(boost::make_tuple(type, side_number));
498
499 for (; siit != hi_siit; siit++) {
500 if (siit->get()->side_number >= 0) {
501
502 const EntityHandle ent = siit->get()->ent;
503 auto bit_number = mField.get_field_bit_number(field_name);
504 auto dit = dofs.get<Unique_mi_tag>().lower_bound(
506 auto hi_dit = dofs.get<Unique_mi_tag>().upper_bound(
508 indices.resize(std::distance(dit, hi_dit));
509 for (; dit != hi_dit; dit++) {
510
511 indices[(*dit)->getEntDofIdx()] = (*dit)->getPetscGlobalDofIdx();
512 }
513 }
514 }
515
517}

◆ getProblemTypeRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getProblemTypeRowIndices ( const std::string & field_name,
EntityType type,
int side_number,
VectorInt & indices ) const
protected

Definition at line 526 of file ForcesAndSourcesCore.cpp.

528 {
530 type, side_number, indices);
531}

◆ getRowNodesIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getRowNodesIndices ( EntitiesFieldData & data,
const int bit_number ) const
protected

get row node indices from FENumeredDofEntity_multiIndex

Definition at line 311 of file ForcesAndSourcesCore.cpp.

312 {
313
314 struct Extractor {
315 boost::weak_ptr<EntityCacheNumeredDofs>
316 operator()(boost::shared_ptr<FieldEntity> &e) {
317 return e->entityCacheRowDofs;
318 }
319 };
320
321 return getNodesIndices(bit_number, getRowFieldEnts(),
322 data.dataOnEntities[MBVERTEX][0].getIndices(),
323 data.dataOnEntities[MBVERTEX][0].getLocalIndices(),
324 Extractor());
325}

◆ getRule() [1/2]

int MoFEM::ForcesAndSourcesCore::getRule ( int order)
protectedvirtual
Deprecated
Use getRule(int row_order, int col_order, int data order)

Reimplemented in AnalyticalDirichletBC::ApproxField::MyTriFE, BodyForceConstantField::MyVolumeFE, BoneRemodeling::DensityMapFe, BoneRemodeling::Remodeling::Fe, CohesiveElement::CohesiveInterfaceElement::MyPrism, ConvectiveMassElement::MyVolumeFE, EdgeForce::MyFE, EdgeSlidingConstrains::MyEdgeFE, EdgeSlidingConstrains::MyEdgeFE, FluidPressure::MyTriangleFE, KelvinVoigtDamper::DamperFE, MagneticElement::TriFE, MagneticElement::VolumeFE, MixTransport::MixTransportElement::MyTriFE, MixTransport::MixTransportElement::MyVolumeFE, MoFEM::EdgeElementForcesAndSourcesCoreOnChildParent, MoFEM::FaceElementForcesAndSourcesCoreOnChildParent, MoFEM::FaceElementForcesAndSourcesCoreOnSide, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, MoFEM::VolumeElementForcesAndSourcesCoreOnContactPrismSide, MoFEM::VolumeElementForcesAndSourcesCoreOnSide, NeumannForcesSurface::MyTriangleFE, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, NonlinearElasticElement::MyVolumeFE, PostProcEdgeOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, SurfaceSlidingConstrains::MyTriangleFE, SurfaceSlidingConstrains::MyTriangleFE, ThermalElement::MyTriFE, ThermalElement::MyVolumeFE, and ThermalStressElement::MyVolumeFE.

Definition at line 2115 of file ForcesAndSourcesCore.cpp.

2115{ return 2 * order; }

◆ getRule() [2/2]

int MoFEM::ForcesAndSourcesCore::getRule ( int order_row,
int order_col,
int order_data )
protectedvirtual

another variant of getRule

Parameters
order_roworder of base function on row
order_colorder of base function on columns
order_dataorder of base function approximating data
Returns
integration rule

This function is overloaded by the user. The integration rule is set such that specific operator implemented by the user is integrated accurately. For example if user implement bilinear operator

\[ b(u,v) = \int_\mathcal{T} \frac{\partial u_i}{\partial x_j}\frac{\partial v_i}{\partial x_j} \textrm{d}\mathcal{T} \]

then if \(u\) and \(v\) are polynomial of given order, then exact integral would be

int getRule(int order) { return
2*(order-1); };
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule

The integration points and weights are set appropriately for given entity type and integration rule from quad.c

Method ForcesAndSourcesCore::getRule takes at argument takes maximal polynomial order set on the element on all fields defined on the element. If a user likes to have more control, another variant of this function can be called which distinguishing between field orders on rows, columns and data, the i.e. first argument of a bilinear form, the second argument of bilinear form and field coefficients on the element.

Note
If user set rule to -1 or any other negative integer, then method ForcesAndSourcesCore::setGaussPts is called. In that method user can implement own (specific) integration method.
Bug
this function should be const

Reimplemented in EdgeFE, QuadFE, and TriFE.

Definition at line 2103 of file ForcesAndSourcesCore.cpp.

2104 {
2105 return getRuleHook ? getRuleHook(order_row, order_col, order_data)
2106 : getRule(order_data);
2107}

◆ getSpacesAndBaseOnEntities()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities ( EntitiesFieldData & data) const
protected

Get field approximation space and base on entities.

Definition at line 989 of file ForcesAndSourcesCore.cpp.

990 {
992
993 if (nInTheLoop == 0) {
994 data.bAse.reset();
995 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
996 data.spacesOnEntities[t].reset();
997 data.basesOnEntities[t].reset();
998 }
999 for (int s = 0; s != LASTSPACE; ++s) {
1000 data.basesOnSpaces[s].reset();
1001 data.brokenBasesOnSpaces[s].reset();
1002 }
1003 }
1004
1005 auto fe_type = numeredEntFiniteElementPtr->getEntType();
1006
1007 if (getDataFieldEntsPtr())
1008 for (auto e : getDataFieldEnts()) {
1009 if (auto ptr = e.lock()) {
1010 // get data from entity
1011 const EntityType type = ptr->getEntType();
1012 if (DefaultElementAdjacency::getDefTypeMap(fe_type, type)) {
1013 const FieldSpace space = ptr->getSpace();
1014 const FieldContinuity continuity = ptr->getContinuity();
1015 const FieldApproximationBase approx = ptr->getApproxBase();
1016
1017 // set data
1018 switch (continuity) {
1019 case CONTINUOUS:
1020 data.basesOnSpaces[space].set(approx);
1021 break;
1022 case DISCONTINUOUS:
1023 data.brokenBasesOnSpaces[space].set(approx);
1024 break;
1025 default:
1026 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1027 "Continuity not defined");
1028 }
1029
1030 data.bAse.set(approx);
1031 data.spacesOnEntities[type].set(space);
1032 data.basesOnEntities[type].set(approx);
1033 }
1034 }
1035 }
1036 else
1037 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1038 "data fields ents not allocated on element");
1039
1040 for (auto space = 0; space != LASTSPACE; ++space)
1041 if ((data.brokenBasesOnSpaces[space] & data.basesOnSpaces[space]).any())
1042 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1043 "Discontinuous and continuous bases on the same space");
1044
1046}
FieldContinuity
Field continuity.
Definition definitions.h:99
int nInTheLoop
number currently of processed method
static bool getDefTypeMap(const EntityType fe_type, const EntityType ent_type)
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr() const

◆ getUserPolynomialBase()

auto & MoFEM::ForcesAndSourcesCore::getUserPolynomialBase ( )
inline

Get the User Polynomial Base object.

Returns
boost::shared_ptr<BaseFunction>&

Definition at line 97 of file ForcesAndSourcesCore.hpp.

97{ return userPolynomialBasePtr; }
boost::shared_ptr< BaseFunction > userPolynomialBasePtr
Pointer to user polynomial base.

◆ loopOverOperators()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::loopOverOperators ( )
protected

Iterate user data operators.

Returns
MoFEMErrorCode

Definition at line 1389 of file ForcesAndSourcesCore.cpp.

1389 {
1391
1392 using UDO = UserDataOperator;
1393
1394 std::array<std::string, 2> field_name;
1395 std::array<const Field *, 3> field_struture;
1396 std::array<int, 2>
1397 field_id; // note the this is field bit number (nor field bit)
1398 std::array<FieldSpace, 2> space;
1399 std::array<FieldApproximationBase, 2> base;
1400
1401 constexpr std::array<UDO::OpType, 2> types = {UDO::OPROW, UDO::OPCOL};
1402 std::array<int, 2> last_eval_field_id = {0, 0};
1403
1404 std::array<boost::shared_ptr<EntitiesFieldData>, 2> op_data;
1405
1406 auto swap_bases = [&](auto &op) {
1408 for (size_t ss = 0; ss != 2; ++ss) {
1409 if (op.getOpType() & types[ss] || op.getOpType() & UDO::OPROWCOL) {
1410 switch (base[ss]) {
1412 CHKERR op_data[ss]->baseSwap(field_name[ss],
1414 default:
1415 break;
1416 }
1417 }
1418 }
1419
1421 };
1422
1423 const EntityType type = numeredEntFiniteElementPtr->getEntType();
1424
1425 // evaluate entity data only, no field specific data provided or known
1426 auto evaluate_op_space = [&](auto &op) {
1428
1429 // reseat all data which all field dependent
1430 dataOnElement[op.sPace]->resetFieldDependentData();
1431 std::fill(last_eval_field_id.begin(), last_eval_field_id.end(), 0);
1432
1433 switch (op.sPace) {
1434 case NOSPACE:
1435 LOG_OP(op);
1436 try {
1437 CHKERR op.doWork(
1438 0, MBENTITYSET,
1439 dataOnElement[NOSPACE]->dataOnEntities[MBENTITYSET][0]);
1440 }
1441 CATCH_OP_ERRORS(op);
1442 break;
1443 case NOFIELD:
1444 case H1:
1445 case HCURL:
1446 case HDIV:
1447 case L2:
1448 LOG_OP(op);
1449 try {
1450
1451 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
1452 for (auto &e : dataOnElement[op.sPace]->dataOnEntities[t]) {
1453 e.getSpace() = op.sPace;
1454 }
1455 }
1456
1457 CHKERR op.opRhs(*dataOnElement[op.sPace], false);
1458 }
1459 CATCH_OP_ERRORS(op);
1460 break;
1461 default:
1462 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1463 "Not implemented for this space <%s>", FieldSpaceNames[op.sPace]);
1464 }
1465
1467 };
1468
1469 // set entity data
1470 auto set_op_entities_data = [&](auto ss, auto &op) {
1472
1473#ifndef NDEBUG
1474 if ((op.getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
1476 .none()) {
1477 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1478 "no data field < %s > on finite element < %s >",
1479 field_name[ss].c_str(), getFEName().c_str());
1480 }
1481#endif
1482
1483 op_data[ss] =
1484 !ss ? dataOnElement[space[ss]] : derivedDataOnElement[space[ss]];
1485
1486 for (auto &data : op_data[ss]->dataOnEntities[MBENTITYSET]) {
1487 CHKERR data.resetFieldDependentData();
1488 }
1489
1490 auto get_data_for_nodes = [&]() {
1492 if (!ss)
1493 CHKERR getRowNodesIndices(*op_data[ss], field_id[ss]);
1494 else
1495 CHKERR getColNodesIndices(*op_data[ss], field_id[ss]);
1496 CHKERR getNodesFieldData(*op_data[ss], field_id[ss]);
1498 };
1499
1500 // get data on entities but not nodes
1501 auto get_data_for_entities = [&]() {
1503 CHKERR getEntityFieldData(*op_data[ss], field_id[ss], MBEDGE);
1504 if (!ss)
1505 CHKERR getEntityRowIndices(*op_data[ss], field_id[ss], MBEDGE);
1506 else
1507 CHKERR getEntityColIndices(*op_data[ss], field_id[ss], MBEDGE);
1509 };
1510
1511 auto get_data_for_meshset = [&]() {
1513
1514 if (!ss) {
1515
1516 struct Extractor {
1517 boost::weak_ptr<EntityCacheNumeredDofs>
1518 operator()(boost::shared_ptr<FieldEntity> &e) {
1519 return e->entityCacheRowDofs;
1520 }
1521 };
1522
1523 CHKERR getNoFieldEntityFieldData(*op_data[ss], field_id[ss],
1524 Extractor());
1525 } else {
1526
1527 struct Extractor {
1528 boost::weak_ptr<EntityCacheNumeredDofs>
1529 operator()(boost::shared_ptr<FieldEntity> &e) {
1530 return e->entityCacheColDofs;
1531 }
1532 };
1533
1534 CHKERR getNoFieldEntityFieldData(*op_data[ss], field_id[ss],
1535 Extractor());
1536 }
1537
1539 };
1540
1541 switch (space[ss]) {
1542 case H1:
1543 CHKERR get_data_for_nodes();
1544 case HCURL:
1545 case HDIV:
1546 CHKERR get_data_for_entities();
1547 break;
1548 case L2:
1549 switch (type) {
1550 case MBVERTEX:
1551 CHKERR get_data_for_nodes();
1552 break;
1553 default:
1554 CHKERR get_data_for_entities();
1555 }
1556 break;
1557 case NOFIELD:
1558 CHKERR get_data_for_meshset();
1559 break;
1560 default:
1561 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1562 "not implemented for this space < %s >",
1563 FieldSpaceNames[space[ss]]);
1564 }
1566 };
1567
1568 // evaluate operators with field data
1569 auto evaluate_op_for_fields = [&](auto &op) {
1571
1572 if (op.getOpType() & UDO::OPROW) {
1573 LOG_OP(op);
1574 try {
1575 CHKERR op.opRhs(*op_data[0], false);
1576 }
1577 CATCH_OP_ERRORS(op);
1578 }
1579
1580 if (op.getOpType() & UDO::OPCOL) {
1581 LOG_OP(op);
1582 try {
1583 CHKERR op.opRhs(*op_data[1], false);
1584 }
1585 CATCH_OP_ERRORS(op);
1586 }
1587
1588 if (op.getOpType() & UDO::OPROWCOL) {
1589 LOG_OP(op);
1590 try {
1591 CHKERR op.opLhs(*op_data[0], *op_data[1]);
1592 }
1593 CATCH_OP_ERRORS(op);
1594 }
1596 };
1597
1598 // collect bit ref level on all entities, fields and spaces
1600
1601 auto oit = opPtrVector.begin();
1602 auto hi_oit = opPtrVector.end();
1603
1604 // iterate over all operators on element
1605 for (; oit != hi_oit; oit++) {
1606
1607 LOG_OP(*oit);
1608 try {
1609
1610 CHKERR oit->setPtrFE(this);
1611
1612 if ((oit->opType & UDO::OPSPACE) == UDO::OPSPACE) {
1613
1614 // run operator for space but not specific field, thus inices of fild
1615 // values are not resolved, howvere bases for given space are accessible
1616 CHKERR evaluate_op_space(*oit);
1617
1618 } else if (
1619
1620 (oit->opType & (UDO::OPROW | UDO::OPCOL | UDO::OPROWCOL)) ==
1621 oit->opType
1622
1623 ) {
1624
1625 field_name[0] = oit->rowFieldName;
1626 field_name[1] = oit->colFieldName;
1627
1628 for (size_t ss = 0; ss != 2; ss++) {
1629
1630#ifndef NDEBUG
1631 if (field_name[ss].empty()) {
1632 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1633 "Not set Field name in operator %ld (0-row, 1-column) in "
1634 "operator %s",
1635 ss,
1636 (boost::typeindex::type_id_runtime(*oit).pretty_name())
1637 .c_str());
1638 }
1639#endif
1640
1641 field_struture[ss] = mField.get_field_structure(field_name[ss]);
1642 field_id[ss] = field_struture[ss]->getBitNumber();
1643 space[ss] = field_struture[ss]->getSpace();
1644 base[ss] = field_struture[ss]->getApproxBase();
1645 }
1646
1647 // not that if field name do not change between operators, entity field
1648 // data are nor rebuild
1649 for (size_t ss = 0; ss != 2; ss++) {
1650
1651 if (oit->getOpType() & types[ss] ||
1652 oit->getOpType() & UDO::OPROWCOL) {
1653 if (last_eval_field_id[ss] != field_id[ss]) {
1654 CHKERR set_op_entities_data(ss, *oit);
1655 last_eval_field_id[ss] = field_id[ss];
1656 }
1657 }
1658 }
1659
1660 CHKERR swap_bases(*oit);
1661
1662 // run operator for given field or row, column or both
1663 CHKERR evaluate_op_for_fields(*oit);
1664
1665 CHKERR swap_bases(*oit);
1666
1667 } else {
1668
1669 for (int i = 0; i != 5; ++i)
1670 if (oit->opType & (1 << i))
1671 MOFEM_LOG("SELF", Sev::error) << UDO::OpTypeNames[i];
1672 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1673 "Impossible operator type");
1674 }
1675 }
1676 CATCH_OP_ERRORS(*oit);
1677 }
1679}
#define CATCH_OP_ERRORS(OP)
#define LOG_OP(OP)
static const char *const FieldSpaceNames[]
Definition definitions.h:92
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
virtual BitFieldId get_field_id(const std::string &name) const =0
Get field Id.
auto getFEName() const
get finite element name
FieldBitNumber getBitNumber() const
Get number of set bit in Field ID. Each field has uid, get getBitNumber get number of bit set for giv...
MoFEMErrorCode getEntityRowIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getRowNodesIndices(EntitiesFieldData &data, const int bit_number) const
get row node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getNoFieldEntityFieldData(EntitiesFieldData &data, const int bit_number, EXTRACTOR &&extractor) const
Get field data on entities where no field is defined.
MoFEMErrorCode getColNodesIndices(EntitiesFieldData &data, const int bit_number) const
get col node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
MoFEMErrorCode getEntityColIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const

◆ operator()()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::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::BasicMethod.

Reimplemented in MoFEM::ContactPrismElementForcesAndSourcesCore, MoFEM::EdgeElementForcesAndSourcesCore, MoFEM::FaceElementForcesAndSourcesCore, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::FlatPrismElementForcesAndSourcesCore, MoFEM::PipelineManager::MeshsetFE, MoFEM::VertexElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCoreOnSide, and NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE.

Definition at line 2149 of file ForcesAndSourcesCore.cpp.

2149 {
2151 if (operatorHook) {
2152 ierr = operatorHook();
2153 CHKERRG(ierr);
2154 } else {
2155#ifndef NDEBUG
2156 MOFEM_LOG("SELF", Sev::warning)
2157 << "No method operator() overloaded on element entity on finite "
2158 "element <"
2159 << boost::typeindex::type_id_runtime(*this).pretty_name() << ">";
2160#endif
2161 }
2163}
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
boost::function< MoFEMErrorCode()> operatorHook
Hook function for operator.

◆ postProcess()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::postProcess ( )
virtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::BasicMethod.

Reimplemented in ConvectiveMassElement::MyVolumeFE, KelvinVoigtDamper::DamperFE, MoFEM::PipelineManager::MeshsetFE, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, NonlinearElasticElement::MyVolumeFE, PostProcEdgeOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcFatPrismOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, Smoother::MyVolumeFE, and Smoother::MyVolumeFE.

Definition at line 2164 of file ForcesAndSourcesCore.cpp.

2164 {
2166 if (postProcessHook) {
2168 CHKERRG(ierr);
2169 }
2171}
boost::function< MoFEMErrorCode()> postProcessHook
Hook function for post-processing.

◆ preProcess()

◆ setGaussPts() [1/2]

◆ setGaussPts() [2/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::setGaussPts ( int order_row,
int order_col,
int order_data )
protectedvirtual

set user specific integration rule

This function allows for user defined integration rule. The key is to called matrix gaussPts, which is used by other MoFEM procedures. Matrix has number of rows equal to problem dimension plus one, where last index is used to store weight values. Number of columns is equal to number of integration points.

Note
This function is called if method ForcesAndSourcesCore::getRule is returning integer -1 or any other negative integer.

User sets

where

gaussPts.resize(dim+1,nb_gauss_pts);

number rows represents local coordinates of integration points in reference element, where last index in row is for integration weight.

Reimplemented in EdgeFE, QuadFE, and TriFE.

Definition at line 2109 of file ForcesAndSourcesCore.cpp.

2110 {
2111 return setRuleHook ? setRuleHook(this, order_row, order_col, order_data)
2112 : setGaussPts(order_data);
2113}
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule

◆ setRefineFEPtr()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::setRefineFEPtr ( const ForcesAndSourcesCore * refine_fe_ptr)
private

Set the pointer to face element refined.

Parameters
refine_fe_ptr
Returns
MoFEMErrorCode

Definition at line 1725 of file ForcesAndSourcesCore.cpp.

1726 {
1728 refinePtrFE = const_cast<ForcesAndSourcesCore *>(refine_fe_ptr);
1730}
ForcesAndSourcesCore(Interface &m_field)

◆ setSideFEPtr()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::setSideFEPtr ( const ForcesAndSourcesCore * side_fe_ptr)
private

Set the pointer to face element on the side.

Note
Function is is used by face element, while it iterates over elements on the side
Parameters
side_fe_ptr
Returns
MoFEMErrorCode

Definition at line 1719 of file ForcesAndSourcesCore.cpp.

1719 {
1721 sidePtrFE = const_cast<ForcesAndSourcesCore *>(side_fe_ptr);
1723}

Friends And Related Symbol Documentation

◆ EdgeElementForcesAndSourcesCoreOnChildParent

Definition at line 534 of file ForcesAndSourcesCore.hpp.

◆ FaceElementForcesAndSourcesCoreOnChildParent

Definition at line 533 of file ForcesAndSourcesCore.hpp.

◆ FaceElementForcesAndSourcesCoreOnSide

Definition at line 532 of file ForcesAndSourcesCore.hpp.

◆ OpBrokenLoopSide

template<typename E >
friend struct OpBrokenLoopSide
friend

Definition at line 538 of file ForcesAndSourcesCore.hpp.

◆ OpCopyGeomDataToE

template<int DIM>
friend struct OpCopyGeomDataToE
friend

Definition at line 537 of file ForcesAndSourcesCore.hpp.

◆ UserDataOperator

◆ VolumeElementForcesAndSourcesCoreOnContactPrismSide

Definition at line 535 of file ForcesAndSourcesCore.hpp.

◆ VolumeElementForcesAndSourcesCoreOnSide

Definition at line 531 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ coordsAtGaussPts

MatrixDouble MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
protected

coordinated at gauss points

Definition at line 541 of file ForcesAndSourcesCore.hpp.

◆ dataH1

EntitiesFieldData& MoFEM::ForcesAndSourcesCore::dataH1
protected

Definition at line 468 of file ForcesAndSourcesCore.hpp.

◆ dataHcurl

EntitiesFieldData& MoFEM::ForcesAndSourcesCore::dataHcurl
protected

Definition at line 469 of file ForcesAndSourcesCore.hpp.

◆ dataHdiv

EntitiesFieldData& MoFEM::ForcesAndSourcesCore::dataHdiv
protected

Definition at line 470 of file ForcesAndSourcesCore.hpp.

◆ dataL2

EntitiesFieldData& MoFEM::ForcesAndSourcesCore::dataL2
protected

Definition at line 471 of file ForcesAndSourcesCore.hpp.

◆ dataNoField

EntitiesFieldData& MoFEM::ForcesAndSourcesCore::dataNoField
protected

Definition at line 467 of file ForcesAndSourcesCore.hpp.

◆ dataOnElement

const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE> MoFEM::ForcesAndSourcesCore::dataOnElement
protected

Entity data on element entity rows fields.

Definition at line 458 of file ForcesAndSourcesCore.hpp.

◆ derivedDataOnElement

const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE> MoFEM::ForcesAndSourcesCore::derivedDataOnElement
protected

Entity data on element entity columns fields.

Definition at line 465 of file ForcesAndSourcesCore.hpp.

◆ elementMeasure

double MoFEM::ForcesAndSourcesCore::elementMeasure
protected

Depending on dimension of elements, stores length, area or volume of element.

Definition at line 542 of file ForcesAndSourcesCore.hpp.

◆ elementPolynomialBasePtr

boost::shared_ptr<BaseFunction> MoFEM::ForcesAndSourcesCore::elementPolynomialBasePtr
private

Pointer to entity polynomial base.

Definition at line 493 of file ForcesAndSourcesCore.hpp.

◆ gaussPts

MatrixDouble MoFEM::ForcesAndSourcesCore::gaussPts

Matrix of integration points.

Columns is equal to number of integration points, numver of rows depends on dimension of finite element entity, for example for tetrahedron rows are x,y,z,weight. Last row is integration weight.

FIXME: that should be moved to private class data and acessed only by member function

Examples
EshelbianPlasticity.cpp, hanging_node_approx.cpp, and prism_elements_from_surface.cpp.

Definition at line 109 of file ForcesAndSourcesCore.hpp.

◆ getRuleHook

RuleHookFun MoFEM::ForcesAndSourcesCore::getRuleHook

Hook to get rule.

Todo
check preferred format how works with gcc and clang, see http://www.boost.org/doc/libs/1_64_0/doc/html/function/tutorial.html#idp247873024
Examples
lorentz_force.cpp.

Definition at line 42 of file ForcesAndSourcesCore.hpp.

◆ lastEvaluatedElementEntityType

EntityType MoFEM::ForcesAndSourcesCore::lastEvaluatedElementEntityType
protected

Last evaluated type of element entity.

Definition at line 486 of file ForcesAndSourcesCore.hpp.

◆ mField

Interface& MoFEM::ForcesAndSourcesCore::mField

◆ opPtrVector

boost::ptr_deque<UserDataOperator> MoFEM::ForcesAndSourcesCore::opPtrVector
protected

Vector of finite element users data operators.

Definition at line 477 of file ForcesAndSourcesCore.hpp.

◆ refinePtrFE

ForcesAndSourcesCore* MoFEM::ForcesAndSourcesCore::refinePtrFE
private

Element to integrate parent or child.

u

Definition at line 521 of file ForcesAndSourcesCore.hpp.

◆ setRuleHook

GaussHookFun MoFEM::ForcesAndSourcesCore::setRuleHook

Set function to calculate integration rule.

Examples
hanging_node_approx.cpp.

Definition at line 48 of file ForcesAndSourcesCore.hpp.

◆ sidePtrFE

ForcesAndSourcesCore* MoFEM::ForcesAndSourcesCore::sidePtrFE
private

Element to integrate on the sides.

Definition at line 504 of file ForcesAndSourcesCore.hpp.

◆ userPolynomialBasePtr

boost::shared_ptr<BaseFunction> MoFEM::ForcesAndSourcesCore::userPolynomialBasePtr
private

Pointer to user polynomial base.

Definition at line 498 of file ForcesAndSourcesCore.hpp.


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