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...
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. More...
 
const DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side) const
 Get the entity data. More...
 
DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
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 33 of file ForcesAndSourcesCore.cpp.

34  :
35 
36  mField(m_field), getRuleHook(0), setRuleHook(0),
38 
39  nullptr,
40  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // NOFIELD
41  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // H1
42  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // HCURL
43  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET), // HDIV
44  boost::make_shared<DataForcesAndSourcesCore>(MBENTITYSET) // L2
45 
46  },
48 
49  nullptr,
50  boost::make_shared<DerivedDataForcesAndSourcesCore>(
51  dataOnElement[NOFIELD]), // NOFIELD
52  boost::make_shared<DerivedDataForcesAndSourcesCore>(
53  dataOnElement[H1]), // H1
54  boost::make_shared<DerivedDataForcesAndSourcesCore>(
55  dataOnElement[HCURL]), // HCURL
56  boost::make_shared<DerivedDataForcesAndSourcesCore>(
57  dataOnElement[HDIV]), // HDIV
58  boost::make_shared<DerivedDataForcesAndSourcesCore>(
59  dataOnElement[L2]) // L2
60 
61  },
64  dataHdiv(*dataOnElement[HDIV].get()), dataL2(*dataOnElement[L2].get()),
65  lastEvaluatedElementEntityType(MBMAXTYPE), sidePtrFE(nullptr) {}
DataForcesAndSourcesCore & dataL2
field with continuous normal traction
Definition: definitions.h:178
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
DataForcesAndSourcesCore & dataH1
scalar or vector of scalars describe (no true field)
Definition: definitions.h:175
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:177
DataForcesAndSourcesCore & dataHdiv
continuous field
Definition: definitions.h:176
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:179
DataForcesAndSourcesCore & dataNoField

Member Function Documentation

◆ calBernsteinBezierBaseFunctionsOnElement()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement ( )
protected

Calculate Bernstein-Bezier base.

Returns
MoFEMErrorCode

Definition at line 891 of file ForcesAndSourcesCore.cpp.

891  {
893 
894  auto get_nodal_base_data = [&](DataForcesAndSourcesCore &data,
895  const std::string &field_name) {
897  auto &space = data.dataOnEntities[MBVERTEX][0].getSpace();
898  auto &base = data.dataOnEntities[MBVERTEX][0].getBase();
899  auto &bb_node_order = data.dataOnEntities[MBVERTEX][0].getBBNodeOrder();
900 
901  auto &dofs_by_name_and_type =
902  dataPtr->get<Composite_Name_And_Type_mi_tag>();
903  auto tuple = boost::make_tuple(field_name, MBVERTEX);
904  auto dit = dofs_by_name_and_type.lower_bound(tuple);
905  if (dit == dofs_by_name_and_type.end())
906  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
907  "No nodal dofs on element");
908  auto hi_dit =
909  dataPtr->get<Composite_Name_And_Type_mi_tag>().upper_bound(tuple);
910 
911  if (dit != hi_dit) {
912  auto &first_dof = **dit;
913  space = first_dof.getSpace();
914  base = first_dof.getApproxBase();
915  int num_nodes;
916  CHKERR getNumberOfNodes(num_nodes);
917  bb_node_order.resize(num_nodes, false);
918  bb_node_order.clear();
919  const int nb_dof_idx = first_dof.getNbOfCoeffs();
920 
921  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
922  for (; dit != hi_dit;) {
923  const auto &dof_ptr = *dit;
924  const auto &dof = *dof_ptr;
925  const auto &sn = *dof.sideNumberPtr;
926  const int side_number = sn.side_number;
927  const int brother_side_number = sn.brother_side_number;
928  if (brother_side_number != -1)
929  brother_dofs_vec.emplace_back(dof_ptr);
930  bb_node_order[side_number] = dof.getMaxOrder();
931  for (int ii = 0; ii != nb_dof_idx; ++ii)
932  ++dit;
933  }
934 
935  for (auto &dof_ptr : brother_dofs_vec) {
936  if (const auto d = dof_ptr.lock()) {
937  const auto &sn = d->sideNumberPtr;
938  const int side_number = sn->side_number;
939  const int brother_side_number = sn->brother_side_number;
940  bb_node_order[brother_side_number] = bb_node_order[side_number];
941  }
942  }
943  }
945  };
946 
947  auto get_entity_base_data = [&](DataForcesAndSourcesCore &data,
948  const std::string &field_name,
949  const EntityType type_lo,
950  const EntityType type_hi) {
952  for (EntityType t = type_lo; t != type_hi; ++t) {
953  for (auto &dat : data.dataOnEntities[t]) {
954  dat.getDataOrder() = 0;
955  dat.getBase() = NOBASE;
956  dat.getSpace() = NOSPACE;
957  dat.getFieldData().resize(0, false);
958  dat.getFieldDofs().resize(0, false);
959  }
960  }
961 
962  auto &dofs = const_cast<FEDofEntity_multiIndex &>(
963  numeredEntFiniteElementPtr->getDataDofs());
964  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
965  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
966  if (dit == dofs_by_type.end())
968  auto hi_dit =
969  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
970  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
971  for (; dit != hi_dit;) {
972 
973  auto &dof = **dit;
974  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
975  if (nb_dofs_on_ent) {
976  const EntityType type = dof.getEntType();
977  const int side = dof.sideNumberPtr->side_number;
978  auto &dat = data.dataOnEntities[type][side];
979 
980  const int brother_side = dof.sideNumberPtr->brother_side_number;
981  if (brother_side != -1)
982  brother_dofs_vec.emplace_back(*dit);
983 
984  dat.getBase() = dof.getApproxBase();
985  dat.getSpace() = dof.getSpace();
986  const int ent_order = dof.getMaxOrder();
987  dat.getDataOrder() =
988  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
989  for (int ii = 0; ii != nb_dofs_on_ent; ++ii) {
990  ++dit;
991  }
992  }
993  }
994 
995  for (auto &dof_ptr : brother_dofs_vec) {
996  if (auto d = dof_ptr.lock()) {
997  const EntityType type = d->getEntType();
998  const int side = d->sideNumberPtr->side_number;
999  const int brother_side = d->sideNumberPtr->brother_side_number;
1000  auto &dat = data.dataOnEntities[type][side];
1001  auto &dat_brother = data.dataOnEntities[type][brother_side];
1002  dat_brother.getBase() = dat.getBase();
1003  dat_brother.getSpace() = dat.getSpace();
1004  dat_brother.getDataOrder() = dat.getDataOrder();
1005  }
1006  }
1008  };
1009 
1010  auto &ents_data = *dataFieldEntsPtr;
1011  for (auto &e : ents_data) {
1012  if (e->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
1013  auto space = e->getSpace();
1014  for (EntityType t = MBVERTEX; t != MBPOLYHEDRON; ++t) {
1015  for (auto &dat : (*dataOnElement[space]).dataOnEntities[t]) {
1016  for (auto &ptr : dat.getBBAlphaIndicesByOrderArray())
1017  ptr.reset();
1018  for (auto &ptr : dat.getBBNByOrderArray())
1019  ptr.reset();
1020  for (auto &ptr : dat.getBBDiffNByOrderArray())
1021  ptr.reset();
1022  }
1023  }
1024  }
1025  }
1026 
1027  for (auto &e : ents_data) {
1028  if (e->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
1029  auto field_name = e->getName();
1030  auto space = e->getSpace();
1031  CHKERR get_nodal_base_data(*dataOnElement[space], field_name);
1032  CHKERR get_entity_base_data(*dataOnElement[space], field_name, MBEDGE,
1033  MBPOLYHEDRON);
1034  CHKERR getElementPolynomialBase()->getValue(
1035  gaussPts,
1036  boost::make_shared<EntPolynomialBaseCtx>(
1037  *dataOnElement[space], field_name, static_cast<FieldSpace>(space),
1039  }
1040  }
1041 
1043 };
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:482
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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:601
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.
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:412
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() [1/2]

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

Calculate base functions.

Returns
Error code

Definition at line 821 of file ForcesAndSourcesCore.cpp.

822  {
824  if (dataOnElement[H1]->bAse.test(b)) {
825  switch (static_cast<FieldApproximationBase>(b)) {
826  case NOBASE:
827  break;
829  // BERNSTEIN_BEZIER_BASE is not hierarchical base
830  break;
835  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
836  "Functions genrating approximation base not defined");
837 
838  for (int space = H1; space != LASTSPACE; ++space) {
839  if (dataOnElement[H1]->sPace.test(space) &&
840  dataOnElement[H1]->bAse.test(b) &&
841  dataOnElement[H1]->basesOnSpaces[space].test(b)) {
842  CHKERR getElementPolynomialBase()->getValue(
843  gaussPts,
844  boost::make_shared<EntPolynomialBaseCtx>(
845  *dataOnElement[space], static_cast<FieldSpace>(space),
846  static_cast<FieldApproximationBase>(b), NOBASE));
847  }
848  }
849  break;
850  case USER_BASE:
851  if (!getUserPolynomialBase())
852  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
853  "Functions genrating user approximation base not defined");
854 
855  for (int space = H1; space != LASTSPACE; ++space)
856  if (dataOnElement[H1]->sPace.test(space) &&
857  dataOnElement[H1]->bAse.test(b) &&
858  dataOnElement[H1]->basesOnSpaces[space].test(b)) {
859  CHKERR getUserPolynomialBase()->getValue(
860  gaussPts,
861  boost::make_shared<EntPolynomialBaseCtx>(
862  *dataOnElement[space], static_cast<FieldSpace>(space),
863  static_cast<FieldApproximationBase>(b), NOBASE));
864  }
865  break;
866  default:
868  "Base <%s> not yet implemented",
869  ApproximationBaseNames[static_cast<FieldApproximationBase>(b)]);
870  }
871  }
873 }
MatrixDouble gaussPts
Matrix of integration points.
user implemented approximation base
Definition: definitions.h:159
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:482
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:180
static const char *const ApproximationBaseNames[]
Definition: definitions.h:163
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:151
#define CHKERR
Inline error check.
Definition: definitions.h:601
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412
virtual MPI_Comm & get_comm() const =0
continuous field
Definition: definitions.h:176
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 875 of file ForcesAndSourcesCore.cpp.

875  {
877  /// Use the some node base. Node base is usually used for construction other
878  /// bases.
879  for (int space = HCURL; space != LASTSPACE; ++space) {
880  dataOnElement[space]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
881  dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
882  }
883  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
885  static_cast<FieldApproximationBase>(b));
886  }
888 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:180
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:151
field with continuous tangents
Definition: definitions.h:177
#define CHKERR
Inline error check.
Definition: definitions.h:601
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:412
continuous field
Definition: definitions.h:176
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 1045 of file ForcesAndSourcesCore.cpp.

1045  {
1047 
1048  const EntityType type = numeredEntFiniteElementPtr->getEntType();
1051 
1052  // Data on elements for proper spaces
1053  for (int space = H1; space != LASTSPACE; ++space) {
1054  dataOnElement[space]->setElementType(type);
1055  derivedDataOnElement[space]->setElementType(type);
1056  }
1057 
1059 
1061 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:180
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:412
continuous field
Definition: definitions.h:176
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 258 of file ForcesAndSourcesCore.cpp.

259  {
260  return getNodesIndices(field_name,
261  const_cast<FENumeredDofEntity_multiIndex &>(
262  numeredEntFiniteElementPtr->getColsDofs()),
263  data.dataOnEntities[MBVERTEX][0].getIndices(),
264  data.dataOnEntities[MBVERTEX][0].getLocalIndices());
265 }
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 484 of file ForcesAndSourcesCore.hpp.

486  {
487  return dataOnElement[space]->dataOnEntities[type][side];
488  }
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 499 of file ForcesAndSourcesCore.hpp.

499  {
500  return dataOnElement[space]->dataOnEntities[type][side];
501  }
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 865 of file ForcesAndSourcesCore.hpp.

867  {
868  return getEntityIndices(data, field_name,
869  const_cast<FENumeredDofEntity_multiIndex &>(
870  numeredEntFiniteElementPtr->getColsDofs()),
871  type_lo, type_hi);
872 }
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 149 of file ForcesAndSourcesCore.cpp.

151  {
153 
154  auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable();
155 
156  for (unsigned int s = 0; s != data.size(); ++s)
157  data[s].getDataOrder() = 0;
158 
159  auto &fields_ents =
160  dataFieldEntsPtr->get<Composite_EntType_and_Space_mi_tag>();
161 
162  for (auto r = fields_ents.equal_range(boost::make_tuple(type, space));
163  r.first != r.second; ++r.first) {
164 
165  auto &e = **r.first;
166 
167  auto sit = side_table.find(e.getEnt());
168  if (sit != side_table.end()) {
169 
170  auto &side = *sit;
171  const int side_number = side->side_number;
172  if (side_number >= 0) {
173  ApproximationOrder ent_order = e.getMaxOrder();
174  auto &dat = data[side_number];
175  dat.getDataOrder() =
176  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
177  }
178  } else
179  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
180  "Entity on side of the element not found");
181  }
182 
183  for (auto r = side_table.get<2>().equal_range(type); r.first != r.second;
184  ++r.first) {
185  const int brother_side_number = (*r.first)->brother_side_number;
186  if (brother_side_number != -1) {
187  const int side_number = (*r.first)->side_number;
188  data[brother_side_number].getDataOrder() =
189  data[side_number].getDataOrder();
190  }
191  }
192 
194 }
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:482
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:412

◆ 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 547 of file ForcesAndSourcesCore.hpp.

548  {
549  return getEntityDataOrder(type, space, data.dataOnEntities[type]);
550  }
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 572 of file ForcesAndSourcesCore.cpp.

574  {
576  for (EntityType t = type_lo; t != type_hi; ++t) {
577  for (auto &dat : data.dataOnEntities[t]) {
578  dat.getDataOrder() = 0;
579  dat.getBase() = NOBASE;
580  dat.getSpace() = NOSPACE;
581  dat.getFieldData().resize(0, false);
582  dat.getFieldDofs().resize(0, false);
583  }
584  }
585 
586  auto &dofs = const_cast<FEDofEntity_multiIndex &>(
587  numeredEntFiniteElementPtr->getDataDofs());
588  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
589  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
590  if (dit == dofs_by_type.end())
592  auto hi_dit =
593  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
594  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
595  for (; dit != hi_dit;) {
596 
597  auto &dof = **dit;
598  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
599  if (nb_dofs_on_ent) {
600 
601  const EntityType type = dof.getEntType();
602  const int side = dof.sideNumberPtr->side_number;
603  if (side >= 0) {
604 
605  auto &dat = data.dataOnEntities[type][side];
606  auto &ent_field_dofs = dat.getFieldDofs();
607  auto &ent_field_data = dat.getFieldData();
608  const int brother_side = dof.sideNumberPtr->brother_side_number;
609  if (brother_side != -1)
610  brother_dofs_vec.emplace_back(*dit);
611 
612  if (ent_field_data.empty()) {
613 
614  dat.getBase() = dof.getApproxBase();
615  dat.getSpace() = dof.getSpace();
616  const int ent_order = dof.getMaxOrder();
617  dat.getDataOrder() =
618  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
619  ent_field_data.resize(nb_dofs_on_ent, false);
620  noalias(ent_field_data) = dof.getEntFieldData();
621  ent_field_dofs.resize(nb_dofs_on_ent, false);
622  for (int ii = 0; ii != nb_dofs_on_ent; ++ii) {
623  ent_field_dofs[ii] = *dit;
624  ++dit;
625  }
626  }
627 
628  } else {
629 
630  for (int ii = 0; ii != nb_dofs_on_ent; ++ii)
631  ++dit;
632 
633  }
634 
635  }
636  }
637 
638  for (auto &dof_ptr : brother_dofs_vec) {
639  if (auto d = dof_ptr.lock()) {
640  const EntityType type = d->getEntType();
641  const int side = d->sideNumberPtr->side_number;
642  const int brother_side = d->sideNumberPtr->brother_side_number;
643  auto &dat = data.dataOnEntities[type][side];
644  auto &dat_brother = data.dataOnEntities[type][brother_side];
645  dat_brother.getBase() = dat.getBase();
646  dat_brother.getSpace() = dat.getSpace();
647  dat_brother.getDataOrder() = dat.getDataOrder();
648  dat_brother.getFieldData() = dat.getFieldData();
649  dat_brother.getFieldDofs() = dat.getFieldDofs();
650  }
651  }
652 
654 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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:412

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

270  {
272 
273  for (EntityType t = type_lo; t != type_hi; ++t) {
274  for (auto &dat : data.dataOnEntities[t]) {
275  dat.getIndices().resize(0, false);
276  dat.getLocalIndices().resize(0, false);
277  }
278  }
279 
280  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
281  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
282  if (dit == dofs_by_type.end())
284  auto hi_dit =
285  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
286  for (; dit != hi_dit; ++dit) {
287  auto &dof = **dit;
288  const EntityType type = dof.getEntType();
289  const int side = dof.sideNumberPtr->side_number;
290 
291  if (side >= 0) {
292 
293  auto &dat = data.dataOnEntities[type][side];
294  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
295  if (nb_dofs_on_ent) {
296  const int brother_side = dof.sideNumberPtr->brother_side_number;
297  auto &ent_field_indices = dat.getIndices();
298  auto &ent_field_local_indices = dat.getLocalIndices();
299  if (ent_field_indices.empty()) {
300  ent_field_indices.resize(nb_dofs_on_ent, false);
301  ent_field_local_indices.resize(nb_dofs_on_ent, false);
302  std::fill(ent_field_indices.data().begin(),
303  ent_field_indices.data().end(), -1);
304  std::fill(ent_field_local_indices.data().begin(),
305  ent_field_local_indices.data().end(), -1);
306  }
307  const int idx = dof.getEntDofIdx();
308  ent_field_indices[idx] = dof.getPetscGlobalDofIdx();
309  ent_field_local_indices[idx] = dof.getPetscLocalDofIdx();
310  if (brother_side != -1) {
311  auto &dat_brother = data.dataOnEntities[type][brother_side];
312  dat_brother.getIndices().resize(nb_dofs_on_ent, false);
313  dat_brother.getLocalIndices().resize(nb_dofs_on_ent, false);
314  dat_brother.getIndices()[idx] = dat.getIndices()[idx];
315  dat_brother.getLocalIndices()[idx] = dat.getLocalIndices()[idx];
316  }
317  }
318 
319  }
320  }
321 
323 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412

◆ 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 856 of file ForcesAndSourcesCore.hpp.

858  {
859  return getEntityIndices(data, field_name,
860  const_cast<FENumeredDofEntity_multiIndex &>(
861  numeredEntFiniteElementPtr->getRowsDofs()),
862  type_lo, type_hi);
863 }
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 99 of file ForcesAndSourcesCore.cpp.

101  {
103 
104  auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable().get<2>();
105  for (auto r = side_table.equal_range(type); r.first != r.second; ++r.first) {
106  const int side_number = (*r.first)->side_number;
107  if (side_number >= 0) {
108  const int brother_side_number = (*r.first)->brother_side_number;
109  const int sense = (*r.first)->sense;
110 
111  data[side_number].getSense() = sense;
112  if (brother_side_number != -1)
113  data[brother_side_number].getSense() = sense;
114  }
115  }
117 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
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:412

◆ 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 534 of file ForcesAndSourcesCore.hpp.

534  {
535  return getEntitySense(type, data.dataOnEntities[type]);
536  }
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 691 of file ForcesAndSourcesCore.cpp.

691  {
693  // PetscAttachDebugger();
694  data.facesNodes.resize(4, 3, false);
695  auto &side_table = const_cast<SideNumber_multiIndex &>(
696  numeredEntFiniteElementPtr->getSideNumberTable());
697  auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBTRI, 0));
698  auto hi_siit = side_table.get<1>().upper_bound(boost::make_tuple(MBTRI, 4));
699  if (std::distance(siit, hi_siit) != 4) {
700  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
701  "Should be 4 triangles on tet, side_table not initialized");
702  }
703  const int canonical_face_sense_p1[4][3] = {
704  {0, 1, 3},
705  {1, 2, 3},
706  {0, 3, 2} /**/,
707  {0, 2, 1} /**/}; // second index is offset (positive sense)
708  const int canonical_face_sense_m1[4][3] = {
709  {0, 3, 1},
710  {1, 3, 2},
711  {0, 2, 3},
712  {0, 1, 2}}; // second index is offset (negative sense
713  for (; siit != hi_siit; siit++) {
714  const boost::shared_ptr<SideNumber> side = *siit;
715  int face_conn[3] = {-1, -1, -1};
716  if (side->offset == 0) {
717  face_conn[0] = side->sense == 1
718  ? canonical_face_sense_p1[(int)side->side_number][0]
719  : canonical_face_sense_m1[(int)side->side_number][0];
720  face_conn[1] = side->sense == 1
721  ? canonical_face_sense_p1[(int)side->side_number][1]
722  : canonical_face_sense_m1[(int)side->side_number][1];
723  face_conn[2] = side->sense == 1
724  ? canonical_face_sense_p1[(int)side->side_number][2]
725  : canonical_face_sense_m1[(int)side->side_number][2];
726  }
727  if (side->offset == 1) {
728  face_conn[0] =
729  side->sense == 1
730  ? canonical_face_sense_p1[(int)side->side_number][1]
731  : canonical_face_sense_m1[(int)side->side_number][2] /**/;
732  face_conn[1] = side->sense == 1
733  ? canonical_face_sense_p1[(int)side->side_number][2]
734  : canonical_face_sense_m1[(int)side->side_number][0];
735  face_conn[2] = side->sense == 1
736  ? canonical_face_sense_p1[(int)side->side_number][0]
737  : canonical_face_sense_m1[(int)side->side_number][1];
738  }
739  if (side->offset == 2) {
740  face_conn[0] =
741  side->sense == 1
742  ? canonical_face_sense_p1[(int)side->side_number][2]
743  : canonical_face_sense_m1[(int)side->side_number][1] /**/;
744  face_conn[1] = side->sense == 1
745  ? canonical_face_sense_p1[(int)side->side_number][0]
746  : canonical_face_sense_m1[(int)side->side_number][2];
747  face_conn[2] = side->sense == 1
748  ? canonical_face_sense_p1[(int)side->side_number][1]
749  : canonical_face_sense_m1[(int)side->side_number][0];
750  }
751  for (int nn = 0; nn < 3; nn++)
752  data.facesNodes(side->side_number, nn) = face_conn[nn];
753  {
754  const EntityHandle *conn_tet;
755  int num_nodes_tet;
757  CHKERR mField.get_moab().get_connectivity(ent, conn_tet, num_nodes_tet,
758  true);
759  if (num_nodes_tet != 4)
760  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
761  "data inconsistency");
762  int num_nodes_face;
763  const EntityHandle *conn_face;
764  CHKERR mField.get_moab().get_connectivity(side->ent, conn_face,
765  num_nodes_face, true);
766  if (num_nodes_face != 3)
767  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
768  if (conn_face[0] != conn_tet[data.facesNodes(side->side_number, 0)])
769  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
770  "data inconsistency");
771  if (conn_face[1] != conn_tet[data.facesNodes(side->side_number, 1)])
772  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
773  "data inconsistency");
774  if (conn_face[2] != conn_tet[data.facesNodes(side->side_number, 2)])
775  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
776  "data inconsistency");
777  }
778  }
780 }
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:482
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:601
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412

◆ getMaxColOrder()

int MoFEM::ForcesAndSourcesCore::getMaxColOrder ( ) const

Get max order of approximation for field in columns.

Definition at line 145 of file ForcesAndSourcesCore.cpp.

145  {
146  return getMaxOrder(*colFieldEntsPtr);
147 }
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 132 of file ForcesAndSourcesCore.cpp.

132  {
133  int max_order = 0;
134  for (auto e : *dataFieldEntsPtr) {
135  const int order = e->getMaxOrder();
136  max_order = (max_order < order) ? order : max_order;
137  }
138  return max_order;
139 }
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 141 of file ForcesAndSourcesCore.cpp.

141  {
142  return getMaxOrder(*rowFieldEntsPtr);
143 }
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 484 of file ForcesAndSourcesCore.cpp.

485  {
486 
487  auto get_nodes_field_data = [&](FEDofEntity_multiIndex &dofs,
488  VectorDouble &nodes_data,
489  VectorDofs &nodes_dofs, FieldSpace &space,
491  VectorInt &bb_node_order) {
493 
494  auto &dofs_by_name_and_type = dofs.get<Composite_Name_And_Type_mi_tag>();
495  auto tuple = boost::make_tuple(field_name, MBVERTEX);
496  auto dit = dofs_by_name_and_type.lower_bound(tuple);
497  if (dit == dofs_by_name_and_type.end())
498  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
499  "No nodal dofs on element");
500  auto hi_dit = dofs.get<Composite_Name_And_Type_mi_tag>().upper_bound(tuple);
501 
502  if (dit != hi_dit) {
503  auto &first_dof = **dit;
504  space = first_dof.getSpace();
505  base = first_dof.getApproxBase();
506  int num_nodes;
507  CHKERR getNumberOfNodes(num_nodes);
508  bb_node_order.resize(num_nodes, false);
509  bb_node_order.clear();
510  const int nb_dof_idx = first_dof.getNbOfCoeffs();
511  const int max_nb_dofs = nb_dof_idx * num_nodes;
512  nodes_data.resize(max_nb_dofs, false);
513  nodes_dofs.resize(max_nb_dofs, false);
514  nodes_data.clear();
515 
516  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
517  for (; dit != hi_dit;) {
518  const auto &dof_ptr = *dit;
519  const auto &dof = *dof_ptr;
520  const auto &sn = *dof.sideNumberPtr;
521  const int side_number = sn.side_number;
522  const int brother_side_number = sn.brother_side_number;
523  if (brother_side_number != -1)
524  brother_dofs_vec.emplace_back(dof_ptr);
525 
526  bb_node_order[side_number] = dof.getMaxOrder();
527  int pos = side_number * nb_dof_idx;
528  auto ent_filed_data_vec = dof.getEntFieldData();
529  for (int ii = 0; ii != nb_dof_idx; ++ii) {
530  nodes_data[pos] = ent_filed_data_vec[ii];
531  nodes_dofs[pos] = *dit;
532  ++pos;
533  ++dit;
534  }
535  }
536 
537  for (auto &dof_ptr : brother_dofs_vec) {
538  if (const auto d = dof_ptr.lock()) {
539  const auto &sn = d->sideNumberPtr;
540  const int side_number = sn->side_number;
541  const int brother_side_number = sn->brother_side_number;
542  bb_node_order[brother_side_number] = bb_node_order[side_number];
543  int pos = side_number * nb_dof_idx;
544  int brother_pos = brother_side_number * nb_dof_idx;
545  for (int ii = 0; ii != nb_dof_idx; ++ii) {
546  nodes_data[brother_pos] = nodes_data[pos];
547  nodes_dofs[brother_pos] = nodes_dofs[pos];
548  ++pos;
549  ++brother_pos;
550  }
551  }
552  }
553 
554  } else {
555  nodes_data.resize(0, false);
556  nodes_dofs.resize(0, false);
557  }
558 
560  };
561 
562  return get_nodes_field_data(
563  const_cast<FEDofEntity_multiIndex &>(
564  numeredEntFiniteElementPtr->getDataDofs()),
565  data.dataOnEntities[MBVERTEX][0].getFieldData(),
566  data.dataOnEntities[MBVERTEX][0].getFieldDofs(),
567  data.dataOnEntities[MBVERTEX][0].getSpace(),
568  data.dataOnEntities[MBVERTEX][0].getBase(),
569  data.dataOnEntities[MBVERTEX][0].getBBNodeOrder());
570 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
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:149
FieldSpace
approximation spaces
Definition: definitions.h:173
#define CHKERR
Inline error check.
Definition: definitions.h:601
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:71
ublas::vector< boost::shared_ptr< const FEDofEntity >, DofsAllocator > VectorDofs
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:72

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

200  {
202  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
203  auto tuple = boost::make_tuple(field_name, MBVERTEX);
204  auto dit = dofs_by_type.lower_bound(tuple);
205  auto hi_dit = dofs_by_type.upper_bound(tuple);
206 
207  if (dit != hi_dit) {
208 
209  int num_nodes;
210  CHKERR getNumberOfNodes(num_nodes);
211  int max_nb_dofs = 0;
212  const int nb_dofs_on_vert = (*dit)->getNbOfCoeffs();
213  max_nb_dofs = nb_dofs_on_vert * num_nodes;
214  nodes_indices.resize(max_nb_dofs, false);
215  local_nodes_indices.resize(max_nb_dofs, false);
216  if (std::distance(dit, hi_dit) != max_nb_dofs) {
217  std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
218  std::fill(local_nodes_indices.begin(), local_nodes_indices.end(), -1);
219  }
220 
221  for (; dit != hi_dit; dit++) {
222  auto &dof = **dit;
223  const int idx = dof.getPetscGlobalDofIdx();
224  const int local_idx = dof.getPetscLocalDofIdx();
225  const int side_number = dof.sideNumberPtr->side_number;
226  const int pos = side_number * nb_dofs_on_vert + dof.getDofCoeffIdx();
227  nodes_indices[pos] = idx;
228  local_nodes_indices[pos] = local_idx;
229  const int brother_side_number =
230  (*dit)->sideNumberPtr->brother_side_number;
231  if (brother_side_number != -1) {
232  const int elem_idx =
233  brother_side_number * nb_dofs_on_vert + (*dit)->getDofCoeffIdx();
234  nodes_indices[elem_idx] = idx;
235  local_nodes_indices[elem_idx] = local_idx;
236  }
237  }
238 
239  } else {
240  nodes_indices.resize(0, false);
241  local_nodes_indices.resize(0, false);
242  }
243 
245 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
#define CHKERR
Inline error check.
Definition: definitions.h:601
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412

◆ getNoFieldColIndices()

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

get col NoField indices

Definition at line 353 of file ForcesAndSourcesCore.cpp.

354  {
356  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
357  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
358  }
359  CHKERR getNoFieldIndices(field_name,
360  const_cast<FENumeredDofEntity_multiIndex &>(
361  numeredEntFiniteElementPtr->getColsDofs()),
362  data.dataOnEntities[MBENTITYSET][0].getIndices());
364 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:601
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:412

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

658  {
660  auto dit = dofs.get<FieldName_mi_tag>().lower_bound(field_name);
661  auto hi_dit = dofs.get<FieldName_mi_tag>().upper_bound(field_name);
662  int size = std::distance(dit, hi_dit);
663  ent_field_data.resize(size, false);
664  ent_field_dofs.resize(size, false);
665  for (; dit != hi_dit; dit++) {
666  int idx = (*dit)->getDofCoeffIdx();
667  ent_field_data[idx] = (*dit)->getFieldData();
668  ent_field_dofs[idx] = *dit;
669  }
671 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513

◆ getNoFieldFieldData() [2/2]

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

Definition at line 673 of file ForcesAndSourcesCore.cpp.

674  {
676  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
677  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
678  }
680  field_name,
681  const_cast<FEDofEntity_multiIndex &>(
682  numeredEntFiniteElementPtr->getDataDofs()),
683  data.dataOnEntities[MBENTITYSET][0].getFieldData(),
684  data.dataOnEntities[MBENTITYSET][0].getFieldDofs());
686 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
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:601
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412

◆ getNoFieldIndices()

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

get NoField indices

Definition at line 326 of file ForcesAndSourcesCore.cpp.

328  {
330  auto dit = dofs.get<FieldName_mi_tag>().lower_bound(field_name);
331  auto hi_dit = dofs.get<FieldName_mi_tag>().upper_bound(field_name);
332  indices.resize(std::distance(dit, hi_dit));
333  for (; dit != hi_dit; dit++) {
334  int idx = (*dit)->getPetscGlobalDofIdx();
335  indices[(*dit)->getDofCoeffIdx()] = idx;
336  }
338 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513

◆ getNoFieldRowIndices()

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

get col NoField indices

Definition at line 340 of file ForcesAndSourcesCore.cpp.

341  {
343  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
344  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
345  }
346  CHKERR getNoFieldIndices(field_name,
347  const_cast<FENumeredDofEntity_multiIndex &>(
348  numeredEntFiniteElementPtr->getRowsDofs()),
349  data.dataOnEntities[MBENTITYSET][0].getIndices());
351 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:601
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:412

◆ getNumberOfNodes()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::getNumberOfNodes ( int &  num_nodes) const

Get number of DOFs on element.

Definition at line 67 of file ForcesAndSourcesCore.cpp.

67  {
69 
71  switch (mField.get_moab().type_from_handle(ent)) {
72  case MBVERTEX:
73  num_nodes = 1;
74  break;
75  case MBEDGE:
76  num_nodes = 2;
77  break;
78  case MBTRI:
79  num_nodes = 3;
80  break;
81  case MBQUAD:
82  num_nodes = 4;
83  break;
84  case MBTET:
85  num_nodes = 4;
86  break;
87  case MBPRISM:
88  num_nodes = 6;
89  break;
90  default:
91  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
92  }
93 
95 }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ 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, and MagneticElement.hpp.

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

468  {
469  return getProblemNodesIndices(field_name, *(problemPtr->numeredDofsCols),
470  nodes_indices);
471 }
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 368 of file ForcesAndSourcesCore.cpp.

370  {
372 
373  const Field *field_struture = mField.get_field_structure(field_name);
374  if (field_struture->getSpace() == H1) {
375 
376  int num_nodes;
377  CHKERR getNumberOfNodes(num_nodes);
378  nodes_indices.resize(field_struture->getNbOfCoeffs() * num_nodes, false);
379  std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
380 
381  auto &side_table = const_cast<SideNumber_multiIndex &>(
382  numeredEntFiniteElementPtr->getSideNumberTable());
383  auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 0));
384  auto hi_siit =
385  side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 10000));
386 
387  int nn = 0;
388  for (; siit != hi_siit; siit++, nn++) {
389 
390  if (siit->get()->side_number == -1)
391  continue;
392 
393  const EntityHandle ent = siit->get()->ent;
394  auto dit =
395  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().lower_bound(
396  boost::make_tuple(field_name, ent, 0));
397  auto hi_dit =
398  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().upper_bound(
399  boost::make_tuple(field_name, ent, 10000)); /// very large number
400  for (; dit != hi_dit; dit++) {
401  nodes_indices[siit->get()->side_number * (*dit)->getNbOfCoeffs() +
402  (*dit)->getDofCoeffIdx()] =
403  (*dit)->getPetscGlobalDofIdx();
404  }
405  }
406  } else {
407  nodes_indices.resize(0, false);
408  }
409 
410 
412 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:601
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.
continuous field
Definition: definitions.h:176
virtual const Field * get_field_structure(const std::string &name)=0
get field structure

◆ getProblemNodesRowIndices()

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

Definition at line 453 of file ForcesAndSourcesCore.cpp.

454  {
455  return getProblemNodesIndices(field_name, *(problemPtr->numeredDofsRows),
456  nodes_indices);
457 }
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 474 of file ForcesAndSourcesCore.cpp.

476  {
477  return getProblemTypeIndices(field_name, *(problemPtr->numeredDofsCols), type,
478  side_number, indices);
479 }
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 414 of file ForcesAndSourcesCore.cpp.

416  {
418 
419  indices.resize(0);
420 
421  auto &side_table = const_cast<SideNumber_multiIndex &>(
422  numeredEntFiniteElementPtr->getSideNumberTable());
423  auto siit =
424  side_table.get<1>().lower_bound(boost::make_tuple(type, side_number));
425  auto hi_siit =
426  side_table.get<1>().upper_bound(boost::make_tuple(type, side_number));
427 
428  for (; siit != hi_siit; siit++) {
429 
430  if (siit->get()->side_number == -1)
431  continue;
432 
433  const EntityHandle ent = siit->get()->ent;
434  NumeredDofEntity_multiIndex::index<
435  Composite_Name_And_Ent_And_EntDofIdx_mi_tag>::type::iterator dit,
436  hi_dit;
437  dit = dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().lower_bound(
438  boost::make_tuple(field_name, ent, 0));
439  hi_dit =
440  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().upper_bound(
441  boost::make_tuple(field_name, ent, 10000)); /// very large number
442 
443  indices.resize(std::distance(dit, hi_dit));
444  for (; dit != hi_dit; dit++) {
445 
446  indices[(*dit)->getEntDofIdx()] = (*dit)->getPetscGlobalDofIdx();
447  }
448  }
449 
451 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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 460 of file ForcesAndSourcesCore.cpp.

462  {
463  return getProblemTypeIndices(field_name, *(problemPtr->numeredDofsRows), type,
464  side_number, indices);
465 }
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 248 of file ForcesAndSourcesCore.cpp.

249  {
250  return getNodesIndices(field_name,
251  const_cast<FENumeredDofEntity_multiIndex &>(
252  numeredEntFiniteElementPtr->getRowsDofs()),
253  data.dataOnEntities[MBVERTEX][0].getIndices(),
254  data.dataOnEntities[MBVERTEX][0].getLocalIndices());
255 }
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 1383 of file ForcesAndSourcesCore.cpp.

1384  {
1385  return getRuleHook ? getRuleHook(order_row, order_col, order_data)
1386  : getRule(order_data);
1387 }
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 784 of file ForcesAndSourcesCore.cpp.

785  {
787 
788  if (nInTheLoop == 0) {
789  data.sPace.reset();
790  data.bAse.reset();
791  for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
792  data.spacesOnEntities[t].reset();
793  data.basesOnEntities[t].reset();
794  }
795  for (int s = 0; s != LASTSPACE; ++s) {
796  data.basesOnSpaces[s].reset();
797  }
798  }
799 
800  if (dataFieldEntsPtr)
801  for (auto e : *dataFieldEntsPtr) {
802  // get data from entity
803  const EntityType type = e->getEntType();
804  const FieldSpace space = e->getSpace();
805  const FieldApproximationBase approx = e->getApproxBase();
806 
807  // set data
808  data.sPace.set(space);
809  data.bAse.set(approx);
810  data.spacesOnEntities[type].set(space);
811  data.basesOnEntities[type].set(approx);
812  data.basesOnSpaces[space].set(approx);
813  }
814  else
815  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
816  "data fields ents not allocated on element");
817 
819 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:180
FieldApproximationBase
approximation base
Definition: definitions.h:149
FieldSpace
approximation spaces
Definition: definitions.h:173
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 1093 of file ForcesAndSourcesCore.cpp.

1093  {
1095 
1096  const EntityType type = numeredEntFiniteElementPtr->getEntType();
1099  std::vector<std::string> last_eval_field_name(2);
1100 
1101  boost::ptr_vector<UserDataOperator>::iterator oit, hi_oit;
1102  oit = opPtrVector.begin();
1103  hi_oit = opPtrVector.end();
1104 
1105  for (; oit != hi_oit; oit++) {
1106 
1107  try {
1108 
1109  CHKERR oit->setPtrFE(this);
1110 
1111  if (oit->opType == UserDataOperator::OPLAST) {
1112 
1113  // Set field
1114  switch (oit->sPace) {
1115  case NOSPACE:
1116  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Unknown space");
1117  case NOFIELD:
1118  case H1:
1119  case HCURL:
1120  case HDIV:
1121  case L2:
1122  break;
1123  default:
1124  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1125  "Not implemented for this space", oit->sPace);
1126  }
1127 
1128  // Reseat all data which all field dependent
1129  dataOnElement[oit->sPace]->resetFieldDependentData();
1130  last_eval_field_name[0] = "";
1131 
1132  // Run operator
1133  try {
1134  CHKERR oit->opRhs(*dataOnElement[oit->sPace], false);
1135  }
1136  CATCH_OP_ERRORS(*oit);
1137 
1138  } else {
1139 
1140  boost::shared_ptr<DataForcesAndSourcesCore> op_data[2];
1141  std::array<bool, 2> base_swap;
1142  std::array<std::pair<std::string, FieldApproximationBase>, 2>
1143  base_swap_data;
1144  auto swap_bases = [&]() {
1146  for (size_t ss = 0; ss != 2; ++ss)
1147  if (base_swap[ss])
1148  CHKERR op_data[ss]->baseSwap(base_swap_data[ss].first,
1149  base_swap_data[ss].second);
1151  };
1152 
1153  for (size_t ss = 0; ss != 2; ss++) {
1154 
1155  const std::string field_name =
1156  !ss ? oit->rowFieldName : oit->colFieldName;
1157  if (field_name.empty()) {
1158  SETERRQ2(
1159  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1160  "No field name in operator %d (0-row, 1-column) in operator %s",
1161  ss,
1162  (boost::typeindex::type_id_runtime(*oit).pretty_name())
1163  .c_str());
1164  }
1165  const Field *field_struture = mField.get_field_structure(field_name);
1166  const BitFieldId data_id = field_struture->getId();
1167  const FieldSpace space = field_struture->getSpace();
1168  const FieldApproximationBase base = field_struture->getApproxBase();
1169  op_data[ss] =
1170  !ss ? dataOnElement[space] : derivedDataOnElement[space];
1171 
1172  switch (base) {
1174  base_swap_data[ss] = std::pair<std::string, FieldApproximationBase>(
1175  field_name, AINSWORTH_BERNSTEIN_BEZIER_BASE);
1176  base_swap[ss] = true;
1177  break;
1178  default:
1179  base_swap[ss] = false;
1180  };
1181 
1182  if ((oit->getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
1183  data_id)
1184  .none()) {
1185  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1186  "no data field < %s > on finite element < %s >",
1187  field_name.c_str(), feName.c_str());
1188  }
1189 
1190  if (oit->getOpType() & types[ss] ||
1191  oit->getOpType() & UserDataOperator::OPROWCOL) {
1192 
1193  switch (space) {
1194  case NOSPACE:
1195  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1196  "unknown space");
1197  break;
1198  case NOFIELD:
1199  case H1:
1200  case HCURL:
1201  case HDIV:
1202  case L2:
1203  break;
1204  default:
1205  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1206  "Not implemented for this space", space);
1207  }
1208 
1209  if (last_eval_field_name[ss] != field_name) {
1210 
1211  CHKERR getEntityFieldData(*op_data[ss], field_name, MBEDGE);
1212  if (!ss)
1213  CHKERR getEntityRowIndices(*op_data[ss], field_name, MBEDGE);
1214  else
1215  CHKERR getEntityColIndices(*op_data[ss], field_name, MBEDGE);
1216 
1217  switch (space) {
1218  case NOSPACE:
1219  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1220  "unknown space");
1221  break;
1222  case H1:
1223  if (!ss)
1224  CHKERR getRowNodesIndices(*op_data[ss], field_name);
1225  else
1226  CHKERR getColNodesIndices(*op_data[ss], field_name);
1227  CHKERR getNodesFieldData(*op_data[ss], field_name);
1228  break;
1229  case HCURL:
1230  case HDIV:
1231  break;
1232  case L2:
1233  switch (type) {
1234  case MBVERTEX:
1235  CHKERR getNodesFieldData(*op_data[ss], field_name);
1236  break;
1237  default:
1238  break;
1239  }
1240  break;
1241  case NOFIELD:
1242  if (!getNinTheLoop()) {
1243  // NOFIELD data are the same for each element, can be
1244  // retrieved only once
1245  if (!ss) {
1246  CHKERR getNoFieldRowIndices(*op_data[ss], field_name);
1247  } else {
1248  CHKERR getNoFieldColIndices(*op_data[ss], field_name);
1249  }
1250  CHKERR getNoFieldFieldData(*op_data[ss], field_name);
1251  }
1252  break;
1253  default:
1254  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1255  "Not implemented for this space", space);
1256  }
1257  last_eval_field_name[ss] = field_name;
1258  }
1259  }
1260  }
1261 
1262  CHKERR swap_bases();
1263 
1264  if (oit->getOpType() & UserDataOperator::OPROW) {
1265  try {
1266  CHKERR oit->opRhs(*op_data[0], false);
1267  }
1268  CATCH_OP_ERRORS(*oit);
1269  }
1270 
1271  if (oit->getOpType() & UserDataOperator::OPCOL) {
1272  try {
1273  CHKERR oit->opRhs(*op_data[1], false);
1274  }
1275  CATCH_OP_ERRORS(*oit);
1276  }
1277 
1278  if (oit->getOpType() & UserDataOperator::OPROWCOL) {
1279  try {
1280  CHKERR oit->opLhs(*op_data[0], *op_data[1]);
1281  }
1282  CATCH_OP_ERRORS(*oit);
1283  }
1284 
1285  CHKERR swap_bases();
1286  }
1287  }
1288  CATCH_OP_ERRORS(*oit);
1289  }
1291 }
boost::ptr_vector< UserDataOperator > opPtrVector
Vector of finite element users data operators.
field with continuous normal traction
Definition: definitions.h:178
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:506
scalar or vector of scalars describe (no true field)
Definition: definitions.h:175
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
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:513
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:149
field with continuous tangents
Definition: definitions.h:177
FieldSpace
approximation spaces
Definition: definitions.h:173
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:53
#define CHKERR
Inline error check.
Definition: definitions.h:601
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:412
continuous field
Definition: definitions.h:176
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:179

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

Definition at line 1429 of file ForcesAndSourcesCore.cpp.

1429  {
1431  if (operatorHook) {
1432  ierr = operatorHook();
1433  CHKERRG(ierr);
1434  }
1436 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:549
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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 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 1437 of file ForcesAndSourcesCore.cpp.

1437  {
1439  if (postProcessHook) {
1440  ierr = postProcessHook();
1441  CHKERRG(ierr);
1442  }
1444 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:549
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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, 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 1421 of file ForcesAndSourcesCore.cpp.

1421  {
1423  if (preProcessHook) {
1424  ierr = preProcessHook();
1425  CHKERRG(ierr);
1426  }
1428 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:549
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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 1389 of file ForcesAndSourcesCore.cpp.

1390  {
1391  return setRuleHook ? setRuleHook(order_row, order_col, order_data)
1392  : setGaussPts(order_data);
1393 }
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 PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, NitscheMethod::MyVolumeFE, MoFEM::VolumeElementForcesAndSourcesCoreOnSideBase, and MoFEM::FaceElementForcesAndSourcesCoreOnSideBase.

Definition at line 1399 of file ForcesAndSourcesCore.cpp.

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

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

1328  {
1330  sidePtrFE = const_cast<ForcesAndSourcesCore *>(side_fe_ptr);
1332 }
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:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513

Friends And Related Function Documentation

◆ FaceElementForcesAndSourcesCoreOnSideBase

Definition at line 850 of file ForcesAndSourcesCore.hpp.

◆ UserDataOperator

friend class UserDataOperator
friend

◆ VolumeElementForcesAndSourcesCoreOnSideBase

Definition at line 849 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ dataH1

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataH1
protected

Definition at line 800 of file ForcesAndSourcesCore.hpp.

◆ dataHcurl

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataHcurl
protected

Definition at line 801 of file ForcesAndSourcesCore.hpp.

◆ dataHdiv

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataHdiv
protected

Definition at line 802 of file ForcesAndSourcesCore.hpp.

◆ dataL2

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataL2
protected

Definition at line 803 of file ForcesAndSourcesCore.hpp.

◆ dataNoField

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataNoField
protected

Definition at line 799 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 790 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 797 of file ForcesAndSourcesCore.hpp.

◆ elementPolynomialBasePtr

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

Pointer to entity polynomial base.

Definition at line 825 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

Definition at line 50 of file ForcesAndSourcesCore.hpp.

◆ lastEvaluatedElementEntityType

EntityType MoFEM::ForcesAndSourcesCore::lastEvaluatedElementEntityType
protected

Last evaluated type of element entity.

Definition at line 818 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 809 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 836 of file ForcesAndSourcesCore.hpp.

◆ userPolynomialBasePtr

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

Pointer to user polynomail base.

Definition at line 830 of file ForcesAndSourcesCore.hpp.


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