v0.9.1
Classes | Public Types | Public Member Functions | Public Attributes | Protected Member Functions | List of all members
MoFEM::ForcesAndSourcesCore Struct Reference

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

#include <src/finite_elements/ForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 Data operator to do calculations at integration points.Is inherited and implemented by user to do calculations. It can be used in many different ways but typically is used to integrate matrices (f.e. stiffness matrix) and the right hand vector (f.e. force vector). More...
 

Public Types

typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION, CTX_OPERATORS, CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0, CTX_SET_F = 1 << 0, CTX_SET_A = 1 << 1, CTX_SET_B = 1 << 2,
  CTX_SET_X = 1 << 3, CTX_SET_X_T = 1 << 4, CTX_SET_X_TT = 1 << 6, CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION, CTX_SNESSETJACOBIAN, CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION, CTX_TSSETRHSJACOBIAN, CTX_TSSETIFUNCTION, CTX_TSSETIJACOBIAN,
  CTX_TSTSMONITORSET, CTX_TSNONE
}
 

Public Member Functions

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

Public Attributes

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

Protected Member Functions

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

Deprecated (do not use)

const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 
class UserDataOperator
 
class VolumeElementForcesAndSourcesCoreOnSideBase
 
class FaceElementForcesAndSourcesCoreOnSideBase
 
boost::shared_ptr< BaseFunctionelementPolynomialBasePtr
 Pointer to entity polynomial base. More...
 
boost::shared_ptr< BaseFunctionuserPolynomialBasePtr
 Pointer to user polynomail base. More...
 
ForcesAndSourcesCoresidePtrFE
 Element to integrate on the sides. More...
 
virtual int getRule (int order)
 
virtual MoFEMErrorCode setGaussPts (int order)
 
MoFEMErrorCode setSideFEPtr (const ForcesAndSourcesCore *side_fe_ptr)
 Set the pointer to face element on the side. More...
 

Additional Inherited Members

- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

structure to get information form mofem into DataForcesAndSourcesCore

Examples
EshelbianPlasticity.cpp, hello_world.cpp, prism_elements_from_surface.cpp, prism_polynomial_approximation.cpp, and quad_polynomial_approximation.cpp.

Definition at line 34 of file ForcesAndSourcesCore.hpp.

Member Typedef Documentation

◆ RuleHookFun

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

Definition at line 41 of file ForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ ForcesAndSourcesCore()

MoFEM::ForcesAndSourcesCore::ForcesAndSourcesCore ( Interface m_field)

Definition at line 32 of file ForcesAndSourcesCore.cpp.

33  :
34 
35  mField(m_field), getRuleHook(0), setRuleHook(0),
37 
38  nullptr,
39  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // NOFIELD
40  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // H1
41  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // HCURL
42  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // HDIV
43  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET) // L2
44 
45  },
47 
48  nullptr,
49  boost::make_shared<DerivedDataForcesAndSourcesCore>(
50  dataOnElement[NOFIELD]), // NOFIELD
51  boost::make_shared<DerivedDataForcesAndSourcesCore>(
52  dataOnElement[H1]), // H1
53  boost::make_shared<DerivedDataForcesAndSourcesCore>(
54  dataOnElement[HCURL]), // HCURL
55  boost::make_shared<DerivedDataForcesAndSourcesCore>(
56  dataOnElement[HDIV]), // HDIV
57  boost::make_shared<DerivedDataForcesAndSourcesCore>(
58  dataOnElement[L2]) // L2
59 
60  },
63  dataHdiv(*dataOnElement[HDIV].get()), dataL2(*dataOnElement[L2].get()),
64  lastEvaluatedElementEntityType(MBMAXTYPE), sidePtrFE(nullptr) {}
DataForcesAndSourcesCore & dataL2
field with continuous normal traction
Definition: definitions.h:179
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
DataForcesAndSourcesCore & dataH1
scalar or vector of scalars describe (no true field)
Definition: definitions.h:176
RuleHookFun getRuleHook
Hook to get rule.
DataForcesAndSourcesCore & dataHcurl
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
field with continuous tangents
Definition: definitions.h:178
DataForcesAndSourcesCore & dataHdiv
continuous field
Definition: definitions.h:177
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
RuleHookFun setRuleHook
Set function to calculate integration rule.
field with C-1 continuity
Definition: definitions.h:180
DataForcesAndSourcesCore & dataNoField

Member Function Documentation

◆ calBernsteinBezierBaseFunctionsOnElement()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement ( )
protected

Calculate Bernstein-Bezier base.

Returns
MoFEMErrorCode

Definition at line 860 of file ForcesAndSourcesCore.cpp.

860  {
862 
863  auto get_nodal_base_data = [&](DataForcesAndSourcesCore &data,
864  const std::string &field_name) {
866  auto &space = data.dataOnEntities[MBVERTEX][0].getSpace();
867  auto &base = data.dataOnEntities[MBVERTEX][0].getBase();
868  auto &bb_node_order = data.dataOnEntities[MBVERTEX][0].getBBNodeOrder();
869 
870  auto &dofs_by_name_and_type =
871  dataPtr->get<Composite_Name_And_Type_mi_tag>();
872  auto tuple = boost::make_tuple(field_name, MBVERTEX);
873  auto dit = dofs_by_name_and_type.lower_bound(tuple);
874  if (dit == dofs_by_name_and_type.end())
875  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
876  "No nodal dofs on element");
877  auto hi_dit =
878  dataPtr->get<Composite_Name_And_Type_mi_tag>().upper_bound(tuple);
879 
880  if (dit != hi_dit) {
881  auto &first_dof = **dit;
882  space = first_dof.getSpace();
883  base = first_dof.getApproxBase();
884  int num_nodes;
885  CHKERR getNumberOfNodes(num_nodes);
886  bb_node_order.resize(num_nodes, false);
887  bb_node_order.clear();
888  const int nb_dof_idx = first_dof.getNbOfCoeffs();
889 
890  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
891  for (; dit != hi_dit;) {
892  const auto &dof_ptr = *dit;
893  const auto &dof = *dof_ptr;
894  const auto &sn = *dof.sideNumberPtr;
895  const int side_number = sn.side_number;
896  const int brother_side_number = sn.brother_side_number;
897  if (brother_side_number != -1)
898  brother_dofs_vec.emplace_back(dof_ptr);
899  bb_node_order[side_number] = dof.getMaxOrder();
900  for (int ii = 0; ii != nb_dof_idx; ++ii)
901  ++dit;
902  }
903 
904  for (auto &dof_ptr : brother_dofs_vec) {
905  if (const auto d = dof_ptr.lock()) {
906  const auto &sn = d->sideNumberPtr;
907  const int side_number = sn->side_number;
908  const int brother_side_number = sn->brother_side_number;
909  bb_node_order[brother_side_number] = bb_node_order[side_number];
910  }
911  }
912  }
914  };
915 
916  auto get_entity_base_data = [&](DataForcesAndSourcesCore &data,
917  const std::string &field_name,
918  const EntityType type_lo,
919  const EntityType type_hi) {
921  for (EntityType t = type_lo; t != type_hi; ++t) {
922  for (auto &dat : data.dataOnEntities[t]) {
923  dat.getDataOrder() = 0;
924  dat.getBase() = NOBASE;
925  dat.getSpace() = NOSPACE;
926  dat.getFieldData().resize(0, false);
927  dat.getFieldDofs().resize(0, false);
928  }
929  }
930 
931  auto &dofs = const_cast<FEDofEntity_multiIndex &>(
932  numeredEntFiniteElementPtr->getDataDofs());
933  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
934  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
935  if (dit == dofs_by_type.end())
937  auto hi_dit =
938  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
939  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
940  for (; dit != hi_dit;) {
941 
942  auto &dof = **dit;
943  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
944  if (nb_dofs_on_ent) {
945  const EntityType type = dof.getEntType();
946  const int side = dof.sideNumberPtr->side_number;
947  auto &dat = data.dataOnEntities[type][side];
948 
949  const int brother_side = dof.sideNumberPtr->brother_side_number;
950  if (brother_side != -1)
951  brother_dofs_vec.emplace_back(*dit);
952 
953  dat.getBase() = dof.getApproxBase();
954  dat.getSpace() = dof.getSpace();
955  const int ent_order = dof.getMaxOrder();
956  dat.getDataOrder() =
957  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
958  for (int ii = 0; ii != nb_dofs_on_ent; ++ii) {
959  ++dit;
960  }
961  }
962  }
963 
964  for (auto &dof_ptr : brother_dofs_vec) {
965  if (auto d = dof_ptr.lock()) {
966  const EntityType type = d->getEntType();
967  const int side = d->sideNumberPtr->side_number;
968  const int brother_side = d->sideNumberPtr->brother_side_number;
969  auto &dat = data.dataOnEntities[type][side];
970  auto &dat_brother = data.dataOnEntities[type][brother_side];
971  dat_brother.getBase() = dat.getBase();
972  dat_brother.getSpace() = dat.getSpace();
973  dat_brother.getDataOrder() = dat.getDataOrder();
974  }
975  }
977  };
978 
979  auto &ents_data = *dataFieldEntsPtr;
980  for (auto &e : ents_data) {
981  if (e->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
982  auto space = e->getSpace();
983  for (EntityType t = MBVERTEX; t != MBPOLYHEDRON; ++t) {
984  for (auto &dat : (*dataOnElement[space]).dataOnEntities[t]) {
985  for (auto &ptr : dat.getBBAlphaIndicesByOrderArray())
986  ptr.reset();
987  for (auto &ptr : dat.getBBNByOrderArray())
988  ptr.reset();
989  for (auto &ptr : dat.getBBDiffNByOrderArray())
990  ptr.reset();
991  }
992  }
993  }
994  }
995 
996  for (auto &e : ents_data) {
997  if (e->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
998  auto field_name = e->getName();
999  auto space = e->getSpace();
1000  CHKERR get_nodal_base_data(*dataOnElement[space], field_name);
1001  CHKERR get_entity_base_data(*dataOnElement[space], field_name, MBEDGE,
1002  MBPOLYHEDRON);
1003  CHKERR getElementPolynomialBase()->getValue(
1004  gaussPts,
1005  boost::make_shared<EntPolynomialBaseCtx>(
1006  *dataOnElement[space], field_name, static_cast<FieldSpace>(space),
1008  }
1009  }
1010 
1012 };
MatrixDouble gaussPts
Matrix of integration points.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
Pointer to finite element data dofs.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
#define CHKERR
Inline error check.
Definition: definitions.h:602
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_view > dataFieldEntsPtr
Pointer to finite element field entities data view.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

◆ calHierarchicalBaseFunctionsOnElement() [1/2]

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

Calculate base functions.

Returns
Error code

Definition at line 790 of file ForcesAndSourcesCore.cpp.

791  {
793  if (dataOnElement[H1]->bAse.test(b)) {
794  switch (static_cast<FieldApproximationBase>(b)) {
795  case NOBASE:
796  break;
798  // BERNSTEIN_BEZIER_BASE is not hierarchical base
799  break;
804  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
805  "Functions genrating approximation base not defined");
806 
807  for (int space = H1; space != LASTSPACE; ++space) {
808  if (dataOnElement[H1]->sPace.test(space) &&
809  dataOnElement[H1]->bAse.test(b) &&
810  dataOnElement[H1]->basesOnSpaces[space].test(b)) {
811  CHKERR getElementPolynomialBase()->getValue(
812  gaussPts,
813  boost::make_shared<EntPolynomialBaseCtx>(
814  *dataOnElement[space], static_cast<FieldSpace>(space),
815  static_cast<FieldApproximationBase>(b), NOBASE));
816  }
817  }
818  break;
819  case USER_BASE:
820  if (!getUserPolynomialBase())
821  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
822  "Functions genrating user approximation base not defined");
823 
824  for (int space = H1; space != LASTSPACE; ++space)
825  if (dataOnElement[H1]->sPace.test(space) &&
826  dataOnElement[H1]->bAse.test(b) &&
827  dataOnElement[H1]->basesOnSpaces[space].test(b)) {
828  CHKERR getUserPolynomialBase()->getValue(
829  gaussPts,
830  boost::make_shared<EntPolynomialBaseCtx>(
831  *dataOnElement[space], static_cast<FieldSpace>(space),
832  static_cast<FieldApproximationBase>(b), NOBASE));
833  }
834  break;
835  default:
837  "Base <%s> not yet implemented",
838  ApproximationBaseNames[static_cast<FieldApproximationBase>(b)]);
839  }
840  }
842 }
MatrixDouble gaussPts
Matrix of integration points.
user implemented approximation base
Definition: definitions.h:160
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
static const char *const ApproximationBaseNames[]
Definition: definitions.h:164
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
virtual MPI_Comm & get_comm() const =0
continuous field
Definition: definitions.h:177
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

◆ calHierarchicalBaseFunctionsOnElement() [2/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calHierarchicalBaseFunctionsOnElement ( )
protected

Calculate base functions.

Returns
Error code

Definition at line 844 of file ForcesAndSourcesCore.cpp.

844  {
846  /// Use the some node base. Node base is usually used for construction other
847  /// bases.
848  for (int space = HCURL; space != LASTSPACE; ++space) {
849  dataOnElement[space]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
850  dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
851  }
852  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
854  static_cast<FieldApproximationBase>(b));
855  }
857 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
field with continuous tangents
Definition: definitions.h:178
#define CHKERR
Inline error check.
Definition: definitions.h:602
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
continuous field
Definition: definitions.h:177
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.

◆ createDataOnElement()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::createDataOnElement ( )
protected

Create a entity data on element object.

Returns
MoFEMErrorCode

Definition at line 1014 of file ForcesAndSourcesCore.cpp.

1014  {
1016 
1017  const EntityType type = numeredEntFiniteElementPtr->getEntType();
1020 
1021  // Data on elements for proper spaces
1022  for (int space = H1; space != LASTSPACE; ++space) {
1023  dataOnElement[space]->setElementType(type);
1024  derivedDataOnElement[space]->setElementType(type);
1025  }
1026 
1028 
1030 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
continuous field
Definition: definitions.h:177
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.

◆ getColNodesIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getColNodesIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name 
) const
protected

get col node indices from FENumeredDofEntity_multiIndex

Definition at line 227 of file ForcesAndSourcesCore.cpp.

228  {
229  return getNodesIndices(field_name,
230  const_cast<FENumeredDofEntity_multiIndex &>(
231  numeredEntFiniteElementPtr->getColsDofs()),
232  data.dataOnEntities[MBVERTEX][0].getIndices(),
233  data.dataOnEntities[MBVERTEX][0].getLocalIndices());
234 }
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getNodesIndices(const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
get node indices

◆ getElementPolynomialBase()

auto& MoFEM::ForcesAndSourcesCore::getElementPolynomialBase ( )

Get the Entity Polynomial Base object.

Returns
boost::shared_ptr<BaseFunction>&&

Definition at line 427 of file ForcesAndSourcesCore.hpp.

427 { return elementPolynomialBasePtr; }
boost::shared_ptr< BaseFunction > elementPolynomialBasePtr
Pointer to entity polynomial base.

◆ getEntData() [1/2]

const DataForcesAndSourcesCore::EntData& MoFEM::ForcesAndSourcesCore::getEntData ( const FieldSpace  space,
const EntityType  type,
const int  side 
) const

Get the entity data.

Parameters
space
type
side
Returns
const DataForcesAndSourcesCore::EntData&

Definition at line 481 of file ForcesAndSourcesCore.hpp.

483  {
484  return dataOnElement[space]->dataOnEntities[type][side];
485  }
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.

◆ getEntData() [2/2]

DataForcesAndSourcesCore::EntData& MoFEM::ForcesAndSourcesCore::getEntData ( const FieldSpace  space,
const EntityType  type,
const int  side 
)

Get the entity data.

Parameters
space
type
side
Returns
DataForcesAndSourcesCore::EntData&

Definition at line 496 of file ForcesAndSourcesCore.hpp.

496  {
497  return dataOnElement[space]->dataOnEntities[type][side];
498  }
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.

◆ getEntityColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityColIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name,
const EntityType  type_lo = MBVERTEX,
const EntityType  type_hi = MBPOLYHEDRON 
) const
protected

Definition at line 862 of file ForcesAndSourcesCore.hpp.

864  {
865  return getEntityIndices(data, field_name,
866  const_cast<FENumeredDofEntity_multiIndex &>(
867  numeredEntFiniteElementPtr->getColsDofs()),
868  type_lo, type_hi);
869 }
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getEntityIndices(DataForcesAndSourcesCore &data, const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const

◆ getEntityDataOrder() [1/2]

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

Get the entity data order.

Parameters
type
space
data
Returns
MoFEMErrorCode

Definition at line 118 of file ForcesAndSourcesCore.cpp.

120  {
122 
123  auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable();
124 
125  for (unsigned int s = 0; s != data.size(); ++s)
126  data[s].getDataOrder() = 0;
127 
128  auto &fields_ents =
129  dataFieldEntsPtr->get<Composite_EntType_and_Space_mi_tag>();
130 
131  for (auto r = fields_ents.equal_range(boost::make_tuple(type, space));
132  r.first != r.second; ++r.first) {
133 
134  auto &e = **r.first;
135 
136  auto sit = side_table.find(e.getEnt());
137  if (sit != side_table.end()) {
138 
139  auto &side = *sit;
140  const int side_number = side->side_number;
141  if (side_number >= 0) {
142  ApproximationOrder ent_order = e.getMaxOrder();
143  auto &dat = data[side_number];
144  dat.getDataOrder() =
145  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
146  }
147  } else
148  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
149  "Entity on side of the element not found");
150  }
151 
152  for (auto r = side_table.get<2>().equal_range(type); r.first != r.second;
153  ++r.first) {
154  const int brother_side_number = (*r.first)->brother_side_number;
155  if (brother_side_number != -1) {
156  const int side_number = (*r.first)->side_number;
157  data[brother_side_number].getDataOrder() =
158  data[side_number].getDataOrder();
159  }
160  }
161 
163 }
int ApproximationOrder
Approximation on the entity.
Definition: Types.hpp:37
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_view > dataFieldEntsPtr
Pointer to finite element field entities data view.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getEntityDataOrder() [2/2]

template<EntityType type>
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityDataOrder ( DataForcesAndSourcesCore data,
const FieldSpace  space 
) const
protected

Get the entity data order for given space.

Template Parameters
type
Parameters
data
space
Returns
MoFEMErrorCode

Definition at line 544 of file ForcesAndSourcesCore.hpp.

545  {
546  return getEntityDataOrder(type, space, data.dataOnEntities[type]);
547  }
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
Get the entity data order.

◆ getEntityFieldData()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityFieldData ( DataForcesAndSourcesCore data,
const std::string &  field_name,
const EntityType  type_lo = MBVERTEX,
const EntityType  type_hi = MBPOLYHEDRON 
) const
protected

Definition at line 541 of file ForcesAndSourcesCore.cpp.

543  {
545  for (EntityType t = type_lo; t != type_hi; ++t) {
546  for (auto &dat : data.dataOnEntities[t]) {
547  dat.getDataOrder() = 0;
548  dat.getBase() = NOBASE;
549  dat.getSpace() = NOSPACE;
550  dat.getFieldData().resize(0, false);
551  dat.getFieldDofs().resize(0, false);
552  }
553  }
554 
555  auto &dofs = const_cast<FEDofEntity_multiIndex &>(
556  numeredEntFiniteElementPtr->getDataDofs());
557  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
558  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
559  if (dit == dofs_by_type.end())
561  auto hi_dit =
562  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
563  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
564  for (; dit != hi_dit;) {
565 
566  auto &dof = **dit;
567  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
568  if (nb_dofs_on_ent) {
569 
570  const EntityType type = dof.getEntType();
571  const int side = dof.sideNumberPtr->side_number;
572  if (side >= 0) {
573 
574  auto &dat = data.dataOnEntities[type][side];
575  auto &ent_field_dofs = dat.getFieldDofs();
576  auto &ent_field_data = dat.getFieldData();
577  const int brother_side = dof.sideNumberPtr->brother_side_number;
578  if (brother_side != -1)
579  brother_dofs_vec.emplace_back(*dit);
580 
581  if (ent_field_data.empty()) {
582 
583  dat.getBase() = dof.getApproxBase();
584  dat.getSpace() = dof.getSpace();
585  const int ent_order = dof.getMaxOrder();
586  dat.getDataOrder() =
587  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
588  ent_field_data.resize(nb_dofs_on_ent, false);
589  noalias(ent_field_data) = dof.getEntFieldData();
590  ent_field_dofs.resize(nb_dofs_on_ent, false);
591  for (int ii = 0; ii != nb_dofs_on_ent; ++ii) {
592  ent_field_dofs[ii] = (*dit).get();
593  ++dit;
594  }
595  }
596 
597  } else {
598 
599  for (int ii = 0; ii != nb_dofs_on_ent; ++ii)
600  ++dit;
601 
602  }
603 
604  }
605  }
606 
607  for (auto &dof_ptr : brother_dofs_vec) {
608  if (auto d = dof_ptr.lock()) {
609  const EntityType type = d->getEntType();
610  const int side = d->sideNumberPtr->side_number;
611  const int brother_side = d->sideNumberPtr->brother_side_number;
612  auto &dat = data.dataOnEntities[type][side];
613  auto &dat_brother = data.dataOnEntities[type][brother_side];
614  dat_brother.getBase() = dat.getBase();
615  dat_brother.getSpace() = dat.getSpace();
616  dat_brother.getDataOrder() = dat.getDataOrder();
617  dat_brother.getFieldData() = dat.getFieldData();
618  dat_brother.getFieldDofs() = dat.getFieldDofs();
619  }
620  }
621 
623 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getEntityIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name,
FENumeredDofEntity_multiIndex dofs,
const EntityType  type_lo = MBVERTEX,
const EntityType  type_hi = MBPOLYHEDRON 
) const
protected

Definition at line 236 of file ForcesAndSourcesCore.cpp.

239  {
241 
242  for (EntityType t = type_lo; t != type_hi; ++t) {
243  for (auto &dat : data.dataOnEntities[t]) {
244  dat.getIndices().resize(0, false);
245  dat.getLocalIndices().resize(0, false);
246  }
247  }
248 
249  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
250  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
251  if (dit == dofs_by_type.end())
253  auto hi_dit =
254  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
255  for (; dit != hi_dit; ++dit) {
256  auto &dof = **dit;
257  const EntityType type = dof.getEntType();
258  const int side = dof.sideNumberPtr->side_number;
259 
260  if (side >= 0) {
261 
262  auto &dat = data.dataOnEntities[type][side];
263  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
264  if (nb_dofs_on_ent) {
265  const int brother_side = dof.sideNumberPtr->brother_side_number;
266  auto &ent_field_indices = dat.getIndices();
267  auto &ent_field_local_indices = dat.getLocalIndices();
268  if (ent_field_indices.empty()) {
269  ent_field_indices.resize(nb_dofs_on_ent, false);
270  ent_field_local_indices.resize(nb_dofs_on_ent, false);
271  std::fill(ent_field_indices.data().begin(),
272  ent_field_indices.data().end(), -1);
273  std::fill(ent_field_local_indices.data().begin(),
274  ent_field_local_indices.data().end(), -1);
275  }
276  const int idx = dof.getEntDofIdx();
277  ent_field_indices[idx] = dof.getPetscGlobalDofIdx();
278  ent_field_local_indices[idx] = dof.getPetscLocalDofIdx();
279  if (brother_side != -1) {
280  auto &dat_brother = data.dataOnEntities[type][brother_side];
281  dat_brother.getIndices().resize(nb_dofs_on_ent, false);
282  dat_brother.getLocalIndices().resize(nb_dofs_on_ent, false);
283  dat_brother.getIndices()[idx] = dat.getIndices()[idx];
284  dat_brother.getLocalIndices()[idx] = dat.getLocalIndices()[idx];
285  }
286  }
287 
288  }
289  }
290 
292 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getEntityRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntityRowIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name,
const EntityType  type_lo = MBVERTEX,
const EntityType  type_hi = MBPOLYHEDRON 
) const
protected

Definition at line 853 of file ForcesAndSourcesCore.hpp.

855  {
856  return getEntityIndices(data, field_name,
857  const_cast<FENumeredDofEntity_multiIndex &>(
858  numeredEntFiniteElementPtr->getRowsDofs()),
859  type_lo, type_hi);
860 }
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getEntityIndices(DataForcesAndSourcesCore &data, const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const

◆ getEntitySense() [1/2]

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

get sense (orientation) of entity

Parameters
typetype of entity
dataentity data
Returns
error code

Definition at line 68 of file ForcesAndSourcesCore.cpp.

70  {
72 
73  auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable().get<2>();
74  for (auto r = side_table.equal_range(type); r.first != r.second; ++r.first) {
75  const int side_number = (*r.first)->side_number;
76  if (side_number >= 0) {
77  const int brother_side_number = (*r.first)->brother_side_number;
78  const int sense = (*r.first)->sense;
79 
80  data[side_number].getSense() = sense;
81  if (brother_side_number != -1)
82  data[brother_side_number].getSense() = sense;
83  }
84  }
86 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getEntitySense() [2/2]

template<EntityType type>
MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getEntitySense ( DataForcesAndSourcesCore data) const
protected

Get the entity sense (orientation)

Template Parameters
type
Parameters
data
Returns
MoFEMErrorCode

Definition at line 531 of file ForcesAndSourcesCore.hpp.

531  {
532  return getEntitySense(type, data.dataOnEntities[type]);
533  }
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
get sense (orientation) of entity

◆ getFaceTriNodes()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getFaceTriNodes ( DataForcesAndSourcesCore data) const
protected

Get nodes on triangles.

Definition at line 660 of file ForcesAndSourcesCore.cpp.

660  {
662  // PetscAttachDebugger();
663  data.facesNodes.resize(4, 3, false);
664  auto &side_table = const_cast<SideNumber_multiIndex &>(
665  numeredEntFiniteElementPtr->getSideNumberTable());
666  auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBTRI, 0));
667  auto hi_siit = side_table.get<1>().upper_bound(boost::make_tuple(MBTRI, 4));
668  if (std::distance(siit, hi_siit) != 4) {
669  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
670  "Should be 4 triangles on tet, side_table not initialized");
671  }
672  const int canonical_face_sense_p1[4][3] = {
673  {0, 1, 3},
674  {1, 2, 3},
675  {0, 3, 2} /**/,
676  {0, 2, 1} /**/}; // second index is offset (positive sense)
677  const int canonical_face_sense_m1[4][3] = {
678  {0, 3, 1},
679  {1, 3, 2},
680  {0, 2, 3},
681  {0, 1, 2}}; // second index is offset (negative sense
682  for (; siit != hi_siit; siit++) {
683  const boost::shared_ptr<SideNumber> side = *siit;
684  int face_conn[3] = {-1, -1, -1};
685  if (side->offset == 0) {
686  face_conn[0] = side->sense == 1
687  ? canonical_face_sense_p1[(int)side->side_number][0]
688  : canonical_face_sense_m1[(int)side->side_number][0];
689  face_conn[1] = side->sense == 1
690  ? canonical_face_sense_p1[(int)side->side_number][1]
691  : canonical_face_sense_m1[(int)side->side_number][1];
692  face_conn[2] = side->sense == 1
693  ? canonical_face_sense_p1[(int)side->side_number][2]
694  : canonical_face_sense_m1[(int)side->side_number][2];
695  }
696  if (side->offset == 1) {
697  face_conn[0] =
698  side->sense == 1
699  ? canonical_face_sense_p1[(int)side->side_number][1]
700  : canonical_face_sense_m1[(int)side->side_number][2] /**/;
701  face_conn[1] = side->sense == 1
702  ? canonical_face_sense_p1[(int)side->side_number][2]
703  : canonical_face_sense_m1[(int)side->side_number][0];
704  face_conn[2] = side->sense == 1
705  ? canonical_face_sense_p1[(int)side->side_number][0]
706  : canonical_face_sense_m1[(int)side->side_number][1];
707  }
708  if (side->offset == 2) {
709  face_conn[0] =
710  side->sense == 1
711  ? canonical_face_sense_p1[(int)side->side_number][2]
712  : canonical_face_sense_m1[(int)side->side_number][1] /**/;
713  face_conn[1] = side->sense == 1
714  ? canonical_face_sense_p1[(int)side->side_number][0]
715  : canonical_face_sense_m1[(int)side->side_number][2];
716  face_conn[2] = side->sense == 1
717  ? canonical_face_sense_p1[(int)side->side_number][1]
718  : canonical_face_sense_m1[(int)side->side_number][0];
719  }
720  for (int nn = 0; nn < 3; nn++)
721  data.facesNodes(side->side_number, nn) = face_conn[nn];
722  {
723  const EntityHandle *conn_tet;
724  int num_nodes_tet;
726  CHKERR mField.get_moab().get_connectivity(ent, conn_tet, num_nodes_tet,
727  true);
728  if (num_nodes_tet != 4)
729  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
730  "data inconsistency");
731  int num_nodes_face;
732  const EntityHandle *conn_face;
733  CHKERR mField.get_moab().get_connectivity(side->ent, conn_face,
734  num_nodes_face, true);
735  if (num_nodes_face != 3)
736  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
737  if (conn_face[0] != conn_tet[data.facesNodes(side->side_number, 0)])
738  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
739  "data inconsistency");
740  if (conn_face[1] != conn_tet[data.facesNodes(side->side_number, 1)])
741  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
742  "data inconsistency");
743  if (conn_face[2] != conn_tet[data.facesNodes(side->side_number, 2)])
744  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
745  "data inconsistency");
746  }
747  }
749 }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getMaxColOrder()

int MoFEM::ForcesAndSourcesCore::getMaxColOrder ( ) const

Get max order of approximation for field in columns.

Definition at line 114 of file ForcesAndSourcesCore.cpp.

114  {
115  return getMaxOrder(*colFieldEntsPtr);
116 }
boost::shared_ptr< const FieldEntity_vector_view > colFieldEntsPtr
Pointer to finite element field entities column view.
static int getMaxOrder(const ENTMULTIINDEX &multi_index)

◆ getMaxDataOrder()

int MoFEM::ForcesAndSourcesCore::getMaxDataOrder ( ) const

Get max order of approximation for data fields.

Method 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 101 of file ForcesAndSourcesCore.cpp.

101  {
102  int max_order = 0;
103  for (auto e : *dataFieldEntsPtr) {
104  const int order = e->getMaxOrder();
105  max_order = (max_order < order) ? order : max_order;
106  }
107  return max_order;
108 }
constexpr int order
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_view > dataFieldEntsPtr
Pointer to finite element field entities data view.

◆ getMaxRowOrder()

int MoFEM::ForcesAndSourcesCore::getMaxRowOrder ( ) const

Get max order of approximation for field in rows.

Definition at line 110 of file ForcesAndSourcesCore.cpp.

110  {
111  return getMaxOrder(*rowFieldEntsPtr);
112 }
boost::shared_ptr< const FieldEntity_vector_view > rowFieldEntsPtr
Pointer to finite element field entities row view.
static int getMaxOrder(const ENTMULTIINDEX &multi_index)

◆ getNodesFieldData()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNodesFieldData ( DataForcesAndSourcesCore data,
const std::string &  field_name 
) const
protected

Get data on nodes.

Parameters
dataData structure
field_nameField name
Returns
Error code

Definition at line 453 of file ForcesAndSourcesCore.cpp.

454  {
455 
456  auto get_nodes_field_data = [&](FEDofEntity_multiIndex &dofs,
457  VectorDouble &nodes_data,
458  VectorDofs &nodes_dofs, FieldSpace &space,
460  VectorInt &bb_node_order) {
462 
463  auto &dofs_by_name_and_type = dofs.get<Composite_Name_And_Type_mi_tag>();
464  auto tuple = boost::make_tuple(field_name, MBVERTEX);
465  auto dit = dofs_by_name_and_type.lower_bound(tuple);
466  if (dit == dofs_by_name_and_type.end())
467  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
468  "No nodal dofs on element");
469  auto hi_dit = dofs.get<Composite_Name_And_Type_mi_tag>().upper_bound(tuple);
470 
471  if (dit != hi_dit) {
472  auto &first_dof = **dit;
473  space = first_dof.getSpace();
474  base = first_dof.getApproxBase();
475  int num_nodes;
476  CHKERR getNumberOfNodes(num_nodes);
477  bb_node_order.resize(num_nodes, false);
478  bb_node_order.clear();
479  const int nb_dof_idx = first_dof.getNbOfCoeffs();
480  const int max_nb_dofs = nb_dof_idx * num_nodes;
481  nodes_data.resize(max_nb_dofs, false);
482  nodes_dofs.resize(max_nb_dofs, false);
483  nodes_data.clear();
484 
485  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
486  for (; dit != hi_dit;) {
487  const auto &dof_ptr = *dit;
488  const auto &dof = *dof_ptr;
489  const auto &sn = *dof.sideNumberPtr;
490  const int side_number = sn.side_number;
491  const int brother_side_number = sn.brother_side_number;
492  if (brother_side_number != -1)
493  brother_dofs_vec.emplace_back(dof_ptr);
494 
495  bb_node_order[side_number] = dof.getMaxOrder();
496  int pos = side_number * nb_dof_idx;
497  auto ent_filed_data_vec = dof.getEntFieldData();
498  for (int ii = 0; ii != nb_dof_idx; ++ii) {
499  nodes_data[pos] = ent_filed_data_vec[ii];
500  nodes_dofs[pos] = (*dit).get();
501  ++pos;
502  ++dit;
503  }
504  }
505 
506  for (auto &dof_ptr : brother_dofs_vec) {
507  if (const auto d = dof_ptr.lock()) {
508  const auto &sn = d->sideNumberPtr;
509  const int side_number = sn->side_number;
510  const int brother_side_number = sn->brother_side_number;
511  bb_node_order[brother_side_number] = bb_node_order[side_number];
512  int pos = side_number * nb_dof_idx;
513  int brother_pos = brother_side_number * nb_dof_idx;
514  for (int ii = 0; ii != nb_dof_idx; ++ii) {
515  nodes_data[brother_pos] = nodes_data[pos];
516  nodes_dofs[brother_pos] = nodes_dofs[pos];
517  ++pos;
518  ++brother_pos;
519  }
520  }
521  }
522 
523  } else {
524  nodes_data.resize(0, false);
525  nodes_dofs.resize(0, false);
526  }
527 
529  };
530 
531  return get_nodes_field_data(
532  const_cast<FEDofEntity_multiIndex &>(
533  numeredEntFiniteElementPtr->getDataDofs()),
534  data.dataOnEntities[MBVERTEX][0].getFieldData(),
535  data.dataOnEntities[MBVERTEX][0].getFieldDofs(),
536  data.dataOnEntities[MBVERTEX][0].getSpace(),
537  data.dataOnEntities[MBVERTEX][0].getBase(),
538  data.dataOnEntities[MBVERTEX][0].getBBNodeOrder());
539 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
FieldApproximationBase
approximation base
Definition: definitions.h:150
FieldSpace
approximation spaces
Definition: definitions.h:174
#define CHKERR
Inline error check.
Definition: definitions.h:602
multi_index_container< boost::shared_ptr< FEDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< FEDofEntity::interface_type_DofEntity, const UId &, &FEDofEntity::getGlobalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FEDofEntity::interface_type_DofEntity, EntityHandle, &FEDofEntity::getEnt > >, ordered_non_unique< tag< FieldName_mi_tag >, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef > >, ordered_non_unique< tag< EntType_mi_tag >, const_mem_fun< FEDofEntity::interface_type_RefEntity, EntityType, &FEDofEntity::getEntType > >, ordered_non_unique< tag< Composite_Name_And_Type_mi_tag >, composite_key< FEDofEntity, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef >, const_mem_fun< FEDofEntity::interface_type_RefEntity, EntityType, &FEDofEntity::getEntType > > >, ordered_non_unique< tag< Composite_Name_And_Ent_mi_tag >, composite_key< FEDofEntity, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef >, const_mem_fun< FEDofEntity::interface_type_DofEntity, EntityHandle, &FEDofEntity::getEnt > > >, ordered_non_unique< tag< Composite_Name_Type_And_Side_Number_mi_tag >, composite_key< FEDofEntity, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef >, const_mem_fun< FEDofEntity::interface_type_RefEntity, EntityType, &FEDofEntity::getEntType >, KeyFromKey< member< SideNumber, char, &SideNumber::side_number >, member< FEDofEntity::BaseFEEntity, boost::shared_ptr< SideNumber >, &FEDofEntity::sideNumberPtr > > > > > > FEDofEntity_multiIndex
MultiIndex container keeps FEDofEntity.
ublas::vector< int, IntAllocator > VectorInt
Definition: Types.hpp:72
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:73
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.

◆ getNodesIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNodesIndices ( const boost::string_ref  field_name,
FENumeredDofEntity_multiIndex dofs,
VectorInt nodes_indices,
VectorInt local_nodes_indices 
) const
protected

get node indices

Definition at line 167 of file ForcesAndSourcesCore.cpp.

169  {
171  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
172  auto tuple = boost::make_tuple(field_name, MBVERTEX);
173  auto dit = dofs_by_type.lower_bound(tuple);
174  auto hi_dit = dofs_by_type.upper_bound(tuple);
175 
176  if (dit != hi_dit) {
177 
178  int num_nodes;
179  CHKERR getNumberOfNodes(num_nodes);
180  int max_nb_dofs = 0;
181  const int nb_dofs_on_vert = (*dit)->getNbOfCoeffs();
182  max_nb_dofs = nb_dofs_on_vert * num_nodes;
183  nodes_indices.resize(max_nb_dofs, false);
184  local_nodes_indices.resize(max_nb_dofs, false);
185  if (std::distance(dit, hi_dit) != max_nb_dofs) {
186  std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
187  std::fill(local_nodes_indices.begin(), local_nodes_indices.end(), -1);
188  }
189 
190  for (; dit != hi_dit; dit++) {
191  auto &dof = **dit;
192  const int idx = dof.getPetscGlobalDofIdx();
193  const int local_idx = dof.getPetscLocalDofIdx();
194  const int side_number = dof.sideNumberPtr->side_number;
195  const int pos = side_number * nb_dofs_on_vert + dof.getDofCoeffIdx();
196  nodes_indices[pos] = idx;
197  local_nodes_indices[pos] = local_idx;
198  const int brother_side_number =
199  (*dit)->sideNumberPtr->brother_side_number;
200  if (brother_side_number != -1) {
201  const int elem_idx =
202  brother_side_number * nb_dofs_on_vert + (*dit)->getDofCoeffIdx();
203  nodes_indices[elem_idx] = idx;
204  local_nodes_indices[elem_idx] = local_idx;
205  }
206  }
207 
208  } else {
209  nodes_indices.resize(0, false);
210  local_nodes_indices.resize(0, false);
211  }
212 
214 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.

◆ getNoFieldColIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldColIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name 
) const
protected

get col NoField indices

Definition at line 322 of file ForcesAndSourcesCore.cpp.

323  {
325  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
326  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
327  }
328  CHKERR getNoFieldIndices(field_name,
329  const_cast<FENumeredDofEntity_multiIndex &>(
330  numeredEntFiniteElementPtr->getColsDofs()),
331  data.dataOnEntities[MBENTITYSET][0].getIndices());
333 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:602
MoFEMErrorCode getNoFieldIndices(const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
get NoField indices
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getNoFieldFieldData() [1/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldFieldData ( const boost::string_ref  field_name,
FEDofEntity_multiIndex dofs,
VectorDouble ent_field_data,
VectorDofs ent_field_dofs 
) const
protected

Get field data on nodes.

Definition at line 625 of file ForcesAndSourcesCore.cpp.

627  {
629  auto dit = dofs.get<FieldName_mi_tag>().lower_bound(field_name);
630  auto hi_dit = dofs.get<FieldName_mi_tag>().upper_bound(field_name);
631  int size = std::distance(dit, hi_dit);
632  ent_field_data.resize(size, false);
633  ent_field_dofs.resize(size, false);
634  for (; dit != hi_dit; dit++) {
635  int idx = (*dit)->getDofCoeffIdx();
636  ent_field_data[idx] = (*dit)->getFieldData();
637  ent_field_dofs[idx] = (*dit).get();
638  }
640 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514

◆ getNoFieldFieldData() [2/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldFieldData ( DataForcesAndSourcesCore data,
const boost::string_ref  field_name 
) const
protected

Definition at line 642 of file ForcesAndSourcesCore.cpp.

643  {
645  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
646  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
647  }
649  field_name,
650  const_cast<FEDofEntity_multiIndex &>(
651  numeredEntFiniteElementPtr->getDataDofs()),
652  data.dataOnEntities[MBENTITYSET][0].getFieldData(),
653  data.dataOnEntities[MBENTITYSET][0].getFieldDofs());
655 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
MoFEMErrorCode getNoFieldFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
Get field data on nodes.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:602
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getNoFieldIndices()

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

get NoField indices

Definition at line 295 of file ForcesAndSourcesCore.cpp.

297  {
299  auto dit = dofs.get<FieldName_mi_tag>().lower_bound(field_name);
300  auto hi_dit = dofs.get<FieldName_mi_tag>().upper_bound(field_name);
301  indices.resize(std::distance(dit, hi_dit));
302  for (; dit != hi_dit; dit++) {
303  int idx = (*dit)->getPetscGlobalDofIdx();
304  indices[(*dit)->getDofCoeffIdx()] = idx;
305  }
307 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514

◆ getNoFieldRowIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNoFieldRowIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name 
) const
protected

get col NoField indices

Definition at line 309 of file ForcesAndSourcesCore.cpp.

310  {
312  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
313  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
314  }
315  CHKERR getNoFieldIndices(field_name,
316  const_cast<FENumeredDofEntity_multiIndex &>(
317  numeredEntFiniteElementPtr->getRowsDofs()),
318  data.dataOnEntities[MBENTITYSET][0].getIndices());
320 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:602
MoFEMErrorCode getNoFieldIndices(const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
get NoField indices
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

◆ getOpPtrVector()

boost::ptr_vector<UserDataOperator>& MoFEM::ForcesAndSourcesCore::getOpPtrVector ( )

Use to push back operator for row operator.

It can be used to calculate nodal forces or other quantities on the mesh.

Examples
continuity_check_on_skeleton_3d.cpp, continuity_check_on_skeleton_with_simple_2d.cpp, elasticity_mixed_formulation.cpp, EshelbianPlasticity.cpp, forces_and_sources_testing_edge_element.cpp, forces_and_sources_testing_flat_prism_element.cpp, MagneticElement.hpp, and simple_contact.cpp.

Definition at line 420 of file ForcesAndSourcesCore.hpp.

420 { return opPtrVector; }
boost::ptr_vector< UserDataOperator > opPtrVector
Vector of finite element users data operators.

◆ getProblemNodesColIndices()

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

Definition at line 436 of file ForcesAndSourcesCore.cpp.

437  {
438  return getProblemNodesIndices(field_name, *(problemPtr->numeredDofsCols),
439  nodes_indices);
440 }
const Problem * problemPtr
raw pointer to problem
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredDofsCols
store DOFs on columns for this 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

◆ 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 337 of file ForcesAndSourcesCore.cpp.

339  {
341 
342  const Field *field_struture = mField.get_field_structure(field_name);
343  if (field_struture->getSpace() == H1) {
344 
345  int num_nodes;
346  CHKERR getNumberOfNodes(num_nodes);
347  nodes_indices.resize(field_struture->getNbOfCoeffs() * num_nodes, false);
348  std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
349 
350  auto &side_table = const_cast<SideNumber_multiIndex &>(
351  numeredEntFiniteElementPtr->getSideNumberTable());
352  auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 0));
353  auto hi_siit =
354  side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 10000));
355 
356  int nn = 0;
357  for (; siit != hi_siit; siit++, nn++) {
358 
359  if (siit->get()->side_number == -1)
360  continue;
361 
362  const EntityHandle ent = siit->get()->ent;
363  auto dit =
364  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().lower_bound(
365  boost::make_tuple(field_name, ent, 0));
366  auto hi_dit =
367  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().upper_bound(
368  boost::make_tuple(field_name, ent, 10000)); /// very large number
369  for (; dit != hi_dit; dit++) {
370  nodes_indices[siit->get()->side_number * (*dit)->getNbOfCoeffs() +
371  (*dit)->getDofCoeffIdx()] =
372  (*dit)->getPetscGlobalDofIdx();
373  }
374  }
375  } else {
376  nodes_indices.resize(0, false);
377  }
378 
379 
381 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:602
continuous field
Definition: definitions.h:177
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.

◆ getProblemNodesRowIndices()

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

Definition at line 422 of file ForcesAndSourcesCore.cpp.

423  {
424  return getProblemNodesIndices(field_name, *(problemPtr->numeredDofsRows),
425  nodes_indices);
426 }
const Problem * problemPtr
raw pointer to problem
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredDofsRows
store DOFs on rows for this 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

◆ getProblemTypeColIndices()

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

Definition at line 443 of file ForcesAndSourcesCore.cpp.

445  {
446  return getProblemTypeIndices(field_name, *(problemPtr->numeredDofsCols), type,
447  side_number, indices);
448 }
const Problem * problemPtr
raw pointer to problem
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
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredDofsCols
store DOFs on columns for this problem

◆ 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 383 of file ForcesAndSourcesCore.cpp.

385  {
387 
388  indices.resize(0);
389 
390  auto &side_table = const_cast<SideNumber_multiIndex &>(
391  numeredEntFiniteElementPtr->getSideNumberTable());
392  auto siit =
393  side_table.get<1>().lower_bound(boost::make_tuple(type, side_number));
394  auto hi_siit =
395  side_table.get<1>().upper_bound(boost::make_tuple(type, side_number));
396 
397  for (; siit != hi_siit; siit++) {
398 
399  if (siit->get()->side_number == -1)
400  continue;
401 
402  const EntityHandle ent = siit->get()->ent;
403  NumeredDofEntity_multiIndex::index<
404  Composite_Name_And_Ent_And_EntDofIdx_mi_tag>::type::iterator dit,
405  hi_dit;
406  dit = dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().lower_bound(
407  boost::make_tuple(field_name, ent, 0));
408  hi_dit =
409  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().upper_bound(
410  boost::make_tuple(field_name, ent, 10000)); /// very large number
411 
412  indices.resize(std::distance(dit, hi_dit));
413  for (; dit != hi_dit; dit++) {
414 
415  indices[(*dit)->getEntDofIdx()] = (*dit)->getPetscGlobalDofIdx();
416  }
417  }
418 
420 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getProblemTypeRowIndices()

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

Definition at line 429 of file ForcesAndSourcesCore.cpp.

431  {
432  return getProblemTypeIndices(field_name, *(problemPtr->numeredDofsRows), type,
433  side_number, indices);
434 }
const Problem * problemPtr
raw pointer to problem
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredDofsRows
store DOFs on rows for this problem
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

◆ getRowNodesIndices()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getRowNodesIndices ( DataForcesAndSourcesCore data,
const std::string &  field_name 
) const
protected

get row node indices from FENumeredDofEntity_multiIndex

Definition at line 217 of file ForcesAndSourcesCore.cpp.

218  {
219  return getNodesIndices(field_name,
220  const_cast<FENumeredDofEntity_multiIndex &>(
221  numeredEntFiniteElementPtr->getRowsDofs()),
222  data.dataOnEntities[MBVERTEX][0].getIndices(),
223  data.dataOnEntities[MBVERTEX][0].getLocalIndices());
224 }
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getNodesIndices(const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
get node indices

◆ getRule() [1/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); };

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.

Definition at line 1353 of file ForcesAndSourcesCore.cpp.

1354  {
1355  return getRuleHook ? getRuleHook(order_row, order_col, order_data)
1356  : getRule(order_data);
1357 }
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
RuleHookFun getRuleHook
Hook to get rule.

◆ getRule() [2/2]

int MoFEM::ForcesAndSourcesCore::getRule ( int  order)
protectedvirtual

◆ getSpacesAndBaseOnEntities()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities ( DataForcesAndSourcesCore data) const
protected

Get field approximation space and base on entities.

Definition at line 753 of file ForcesAndSourcesCore.cpp.

754  {
756 
757  if (nInTheLoop == 0) {
758  data.sPace.reset();
759  data.bAse.reset();
760  for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
761  data.spacesOnEntities[t].reset();
762  data.basesOnEntities[t].reset();
763  }
764  for (int s = 0; s != LASTSPACE; ++s) {
765  data.basesOnSpaces[s].reset();
766  }
767  }
768 
769  if (dataFieldEntsPtr)
770  for (auto e : *dataFieldEntsPtr) {
771  // get data from entity
772  const EntityType type = e->getEntType();
773  const FieldSpace space = e->getSpace();
774  const FieldApproximationBase approx = e->getApproxBase();
775 
776  // set data
777  data.sPace.set(space);
778  data.bAse.set(approx);
779  data.spacesOnEntities[type].set(space);
780  data.basesOnEntities[type].set(approx);
781  data.basesOnSpaces[space].set(approx);
782  }
783  else
784  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
785  "data fields ents not allocated on element");
786 
788 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
FieldApproximationBase
approximation base
Definition: definitions.h:150
FieldSpace
approximation spaces
Definition: definitions.h:174
int nInTheLoop
number currently of processed method
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_view > dataFieldEntsPtr
Pointer to finite element field entities data view.

◆ getUserPolynomialBase()

auto& MoFEM::ForcesAndSourcesCore::getUserPolynomialBase ( )

Get the User Polynomial Base object.

Returns
boost::shared_ptr<BaseFunction>&
Examples
EshelbianPlasticity.cpp.

Definition at line 434 of file ForcesAndSourcesCore.hpp.

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

◆ loopOverOperators()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::loopOverOperators ( )
protected

Iterate user data operators.

Returns
MoFEMErrorCode

Definition at line 1062 of file ForcesAndSourcesCore.cpp.

1062  {
1064 
1065  const EntityType type = numeredEntFiniteElementPtr->getEntType();
1068  std::vector<std::string> last_eval_field_name(2);
1069 
1070  boost::ptr_vector<UserDataOperator>::iterator oit, hi_oit;
1071  oit = opPtrVector.begin();
1072  hi_oit = opPtrVector.end();
1073 
1074  for (; oit != hi_oit; oit++) {
1075 
1076  try {
1077 
1078  CHKERR oit->setPtrFE(this);
1079 
1080  if (oit->opType == UserDataOperator::OPLAST) {
1081 
1082  // Set field
1083  switch (oit->sPace) {
1084  case NOSPACE:
1085  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Unknown space");
1086  case NOFIELD:
1087  case H1:
1088  case HCURL:
1089  case HDIV:
1090  case L2:
1091  break;
1092  default:
1093  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1094  "Not implemented for this space", oit->sPace);
1095  }
1096 
1097  // Reseat all data which all field dependent
1098  dataOnElement[oit->sPace]->resetFieldDependentData();
1099  last_eval_field_name[0] = "";
1100 
1101  // Run operator
1102  try {
1103  CHKERR oit->opRhs(*dataOnElement[oit->sPace], false);
1104  }
1105  CATCH_OP_ERRORS(*oit);
1106 
1107  } else {
1108 
1109  boost::shared_ptr<DataForcesAndSourcesCore> op_data[2];
1110  std::array<bool, 2> base_swap;
1111  std::array<std::pair<std::string, FieldApproximationBase>, 2>
1112  base_swap_data;
1113  auto swap_bases = [&]() {
1115  for (size_t ss = 0; ss != 2; ++ss)
1116  if (base_swap[ss])
1117  CHKERR op_data[ss]->baseSwap(base_swap_data[ss].first,
1118  base_swap_data[ss].second);
1120  };
1121 
1122  for (size_t ss = 0; ss != 2; ss++) {
1123 
1124  const std::string field_name =
1125  !ss ? oit->rowFieldName : oit->colFieldName;
1126  if (field_name.empty()) {
1127  SETERRQ2(
1128  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1129  "No field name in operator %d (0-row, 1-column) in operator %s",
1130  ss,
1131  (boost::typeindex::type_id_runtime(*oit).pretty_name())
1132  .c_str());
1133  }
1134  const Field *field_struture = mField.get_field_structure(field_name);
1135  const BitFieldId data_id = field_struture->getId();
1136  const FieldSpace space = field_struture->getSpace();
1137  const FieldApproximationBase base = field_struture->getApproxBase();
1138  op_data[ss] =
1139  !ss ? dataOnElement[space] : derivedDataOnElement[space];
1140 
1141  switch (base) {
1143  base_swap_data[ss] = std::pair<std::string, FieldApproximationBase>(
1144  field_name, AINSWORTH_BERNSTEIN_BEZIER_BASE);
1145  base_swap[ss] = true;
1146  break;
1147  default:
1148  base_swap[ss] = false;
1149  };
1150 
1151  if ((oit->getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
1152  data_id)
1153  .none()) {
1154  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1155  "no data field < %s > on finite element < %s >",
1156  field_name.c_str(), feName.c_str());
1157  }
1158 
1159  if (oit->getOpType() & types[ss] ||
1160  oit->getOpType() & UserDataOperator::OPROWCOL) {
1161 
1162  switch (space) {
1163  case NOSPACE:
1164  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1165  "unknown space");
1166  break;
1167  case NOFIELD:
1168  case H1:
1169  case HCURL:
1170  case HDIV:
1171  case L2:
1172  break;
1173  default:
1174  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1175  "Not implemented for this space", space);
1176  }
1177 
1178  if (last_eval_field_name[ss] != field_name) {
1179 
1180  CHKERR getEntityFieldData(*op_data[ss], field_name, MBEDGE);
1181  if (!ss)
1182  CHKERR getEntityRowIndices(*op_data[ss], field_name, MBEDGE);
1183  else
1184  CHKERR getEntityColIndices(*op_data[ss], field_name, MBEDGE);
1185 
1186  switch (space) {
1187  case NOSPACE:
1188  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1189  "unknown space");
1190  break;
1191  case H1:
1192  if (!ss)
1193  CHKERR getRowNodesIndices(*op_data[ss], field_name);
1194  else
1195  CHKERR getColNodesIndices(*op_data[ss], field_name);
1196  CHKERR getNodesFieldData(*op_data[ss], field_name);
1197  break;
1198  case HCURL:
1199  case HDIV:
1200  break;
1201  case L2:
1202  switch (type) {
1203  case MBVERTEX:
1204  CHKERR getNodesFieldData(*op_data[ss], field_name);
1205  break;
1206  default:
1207  break;
1208  }
1209  break;
1210  case NOFIELD:
1211  if (!getNinTheLoop()) {
1212  // NOFIELD data are the same for each element, can be
1213  // retrieved only once
1214  if (!ss) {
1215  CHKERR getNoFieldRowIndices(*op_data[ss], field_name);
1216  } else {
1217  CHKERR getNoFieldColIndices(*op_data[ss], field_name);
1218  }
1219  CHKERR getNoFieldFieldData(*op_data[ss], field_name);
1220  }
1221  break;
1222  default:
1223  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1224  "Not implemented for this space", space);
1225  }
1226  last_eval_field_name[ss] = field_name;
1227  }
1228  }
1229  }
1230 
1231  CHKERR swap_bases();
1232 
1233  if (oit->getOpType() & UserDataOperator::OPROW) {
1234  try {
1235  CHKERR oit->opRhs(*op_data[0], false);
1236  }
1237  CATCH_OP_ERRORS(*oit);
1238  }
1239 
1240  if (oit->getOpType() & UserDataOperator::OPCOL) {
1241  try {
1242  CHKERR oit->opRhs(*op_data[1], false);
1243  }
1244  CATCH_OP_ERRORS(*oit);
1245  }
1246 
1247  if (oit->getOpType() & UserDataOperator::OPROWCOL) {
1248  try {
1249  CHKERR oit->opLhs(*op_data[0], *op_data[1]);
1250  }
1251  CATCH_OP_ERRORS(*oit);
1252  }
1253 
1254  CHKERR swap_bases();
1255 
1256  }
1257  }
1258  CATCH_OP_ERRORS(*oit);
1259  }
1261 }
boost::ptr_vector< UserDataOperator > opPtrVector
Vector of finite element users data operators.
field with continuous normal traction
Definition: definitions.h:179
MoFEMErrorCode getNoFieldColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
MoFEMErrorCode getColNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col node indices from FENumeredDofEntity_multiIndex
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
scalar or vector of scalars describe (no true field)
Definition: definitions.h:176
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
std::string feName
Name of finite element.
int getNinTheLoop() const
get number of evaluated element in the loop
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
MoFEMErrorCode getNoFieldRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
OpType
Controls loop over entities on element.
MoFEMErrorCode getNoFieldFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
Get field data on nodes.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getNodesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
Get data on nodes.
#define CATCH_OP_ERRORS(OP)
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
FieldApproximationBase
approximation base
Definition: definitions.h:150
field with continuous tangents
Definition: definitions.h:178
FieldSpace
approximation spaces
Definition: definitions.h:174
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:52
#define CHKERR
Inline error check.
Definition: definitions.h:602
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getEntityColIndices(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
continuous field
Definition: definitions.h:177
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getRowNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get row node indices from FENumeredDofEntity_multiIndex
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
MoFEMErrorCode getEntityRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
field with C-1 continuity
Definition: definitions.h:180

◆ 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::VolumeElementForcesAndSourcesCoreSwitch< SWITCH >, MoFEM::FaceElementForcesAndSourcesCoreSwitch< SWITCH >, MoFEM::VolumeElementForcesAndSourcesCoreOnSideSwitch< SWITCH >, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::ContactPrismElementForcesAndSourcesCore, MoFEM::EdgeElementForcesAndSourcesCoreSwitch< SWITCH >, MoFEM::FaceElementForcesAndSourcesCoreOnSideSwitch< SWITCH >, MoFEM::FaceElementForcesAndSourcesCoreOnSideSwitch< FaceElementForcesAndSourcesCore::NO_HO_GEOMETRY >, MoFEM::FlatPrismElementForcesAndSourcesCore, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, and MoFEM::VertexElementForcesAndSourcesCore.

Definition at line 1399 of file ForcesAndSourcesCore.cpp.

1399  {
1401  if (operatorHook) {
1402  ierr = operatorHook();
1403  CHKERRG(ierr);
1404  }
1406 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:550
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
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 PostProcEdgeOnRefinedMesh, PostProcFatPrismOnRefinedMesh, SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, GelModule::Gel::GelFE, KelvinVoigtDamper::DamperFE, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, Smoother::MyVolumeFE, and NonlinearElasticElement::MyVolumeFE.

Definition at line 1407 of file ForcesAndSourcesCore.cpp.

1407  {
1409  if (postProcessHook) {
1410  ierr = postProcessHook();
1411  CHKERRG(ierr);
1412  }
1414 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:550
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::function< MoFEMErrorCode()> postProcessHook
Hook function for post-processing.

◆ preProcess()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::preProcess ( )
virtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from MoFEM::BasicMethod.

Reimplemented in EdgeSlidingConstrains::MyEdgeFE, PostProcEdgeOnRefinedMesh, PostProcFatPrismOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, GelModule::Gel::GelFE, SurfaceSlidingConstrains::MyTriangleFE, KelvinVoigtDamper::DamperFE, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, Smoother::MyVolumeFE, and NonlinearElasticElement::MyVolumeFE.

Definition at line 1391 of file ForcesAndSourcesCore.cpp.

1391  {
1393  if (preProcessHook) {
1394  ierr = preProcessHook();
1395  CHKERRG(ierr);
1396  }
1398 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:550
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::function< MoFEMErrorCode()> preProcessHook
Hook function for pre-processing.

◆ setGaussPts() [1/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.

Definition at line 1359 of file ForcesAndSourcesCore.cpp.

1360  {
1361  return setRuleHook ? setRuleHook(order_row, order_col, order_data)
1362  : setGaussPts(order_data);
1363 }
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
RuleHookFun setRuleHook
Set function to calculate integration rule.

◆ setGaussPts() [2/2]

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::setGaussPts ( int  order)
protectedvirtual
Deprecated:
setGaussPts(int row_order, int col_order, int data order);
Deprecated:
setGaussPts(int row_order, int col_order, int data order);

Reimplemented in PostProcEdgeOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, NitscheMethod::MyVolumeFE, MoFEM::VolumeElementForcesAndSourcesCoreOnSideBase, and MoFEM::FaceElementForcesAndSourcesCoreOnSideBase.

Definition at line 1369 of file ForcesAndSourcesCore.cpp.

1369  {
1371  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Sorry, not implemented");
1373 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514

◆ 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 1298 of file ForcesAndSourcesCore.cpp.

1298  {
1300  sidePtrFE = const_cast<ForcesAndSourcesCore *>(side_fe_ptr);
1302 }
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514

Friends And Related Function Documentation

◆ FaceElementForcesAndSourcesCoreOnSideBase

Definition at line 847 of file ForcesAndSourcesCore.hpp.

◆ UserDataOperator

friend class UserDataOperator
friend

◆ VolumeElementForcesAndSourcesCoreOnSideBase

Definition at line 846 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ dataH1

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataH1
protected

Definition at line 797 of file ForcesAndSourcesCore.hpp.

◆ dataHcurl

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataHcurl
protected

Definition at line 798 of file ForcesAndSourcesCore.hpp.

◆ dataHdiv

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataHdiv
protected

Definition at line 799 of file ForcesAndSourcesCore.hpp.

◆ dataL2

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataL2
protected

Definition at line 800 of file ForcesAndSourcesCore.hpp.

◆ dataNoField

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataNoField
protected

Definition at line 796 of file ForcesAndSourcesCore.hpp.

◆ dataOnElement

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

Entity data on element entity rows fields.

@/}

Definition at line 787 of file ForcesAndSourcesCore.hpp.

◆ derivedDataOnElement

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

Entity data on element entity columns fields.

Definition at line 794 of file ForcesAndSourcesCore.hpp.

◆ elementPolynomialBasePtr

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

Pointer to entity polynomial base.

Definition at line 822 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
prism_elements_from_surface.cpp.

Definition at line 446 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
field_evaluator.cpp.

Definition at line 50 of file ForcesAndSourcesCore.hpp.

◆ lastEvaluatedElementEntityType

EntityType MoFEM::ForcesAndSourcesCore::lastEvaluatedElementEntityType
protected

Last evaluated type of element entity.

Definition at line 815 of file ForcesAndSourcesCore.hpp.

◆ mField

Interface& MoFEM::ForcesAndSourcesCore::mField
Examples
prism_elements_from_surface.cpp.

Definition at line 36 of file ForcesAndSourcesCore.hpp.

◆ opPtrVector

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

Vector of finite element users data operators.

Definition at line 806 of file ForcesAndSourcesCore.hpp.

◆ setRuleHook

RuleHookFun MoFEM::ForcesAndSourcesCore::setRuleHook

Set function to calculate integration rule.

Definition at line 56 of file ForcesAndSourcesCore.hpp.

◆ sidePtrFE

ForcesAndSourcesCore* MoFEM::ForcesAndSourcesCore::sidePtrFE
private

Element to integrate on the sides.

Definition at line 833 of file ForcesAndSourcesCore.hpp.

◆ userPolynomialBasePtr

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

Pointer to user polynomail base.

Definition at line 827 of file ForcesAndSourcesCore.hpp.


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