v0.9.0
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::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 int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
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
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
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 ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- 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 ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

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
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- 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
TSContext ts_ctx
 
TS ts
 time solver 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...
 
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...
 

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...
 

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:173
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
DataForcesAndSourcesCore & dataH1
scalar or vector of scalars describe (no true field)
Definition: definitions.h:170
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:172
DataForcesAndSourcesCore & dataHdiv
continuous field
Definition: definitions.h:171
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:174
DataForcesAndSourcesCore & dataNoField

Member Function Documentation

◆ calBernsteinBezierBaseFunctionsOnElement()

MoFEMErrorCode MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement ( )
protected

Calculate Bernstein-Bezier base.

Returns
MoFEMErrorCode

Definition at line 870 of file ForcesAndSourcesCore.cpp.

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

801  {
803  if (dataOnElement[H1]->bAse.test(b)) {
804  switch (static_cast<FieldApproximationBase>(b)) {
805  case NOBASE:
806  break;
808  // BERNSTEIN_BEZIER_BASE is not hierarchical base
809  break;
814  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
815  "Functions genrating approximation base not defined");
816 
817  for (int space = H1; space != LASTSPACE; ++space) {
818  if (dataOnElement[H1]->sPace.test(space) &&
819  dataOnElement[H1]->bAse.test(b) &&
820  dataOnElement[H1]->basesOnSpaces[space].test(b)) {
821  CHKERR getElementPolynomialBase()->getValue(
822  gaussPts,
823  boost::make_shared<EntPolynomialBaseCtx>(
824  *dataOnElement[space], static_cast<FieldSpace>(space),
825  static_cast<FieldApproximationBase>(b), NOBASE));
826  }
827  }
828  break;
829  case USER_BASE:
830  if (!getUserPolynomialBase())
831  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
832  "Functions genrating user approximation base not defined");
833 
834  for (int space = H1; space != LASTSPACE; ++space)
835  if (dataOnElement[H1]->sPace.test(space) &&
836  dataOnElement[H1]->bAse.test(b) &&
837  dataOnElement[H1]->basesOnSpaces[space].test(b)) {
838  CHKERR getUserPolynomialBase()->getValue(
839  gaussPts,
840  boost::make_shared<EntPolynomialBaseCtx>(
841  *dataOnElement[space], static_cast<FieldSpace>(space),
842  static_cast<FieldApproximationBase>(b), NOBASE));
843  }
844  break;
845  default:
847  "Base <%s> not yet implemented",
848  ApproximationBaseNames[static_cast<FieldApproximationBase>(b)]);
849  }
850  }
852 }
MatrixDouble gaussPts
Matrix of integration points.
user implemented approximation base
Definition: definitions.h:154
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:477
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:175
static const char *const ApproximationBaseNames[]
Definition: definitions.h:158
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:146
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MPI_Comm & get_comm() const =0
continuous field
Definition: definitions.h:171
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 854 of file ForcesAndSourcesCore.cpp.

854  {
856  /// Use the some node base. Node base is usually used for construction other
857  /// bases.
858  for (int space = HCURL; space != LASTSPACE; ++space) {
859  dataOnElement[space]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
860  dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
861  }
862  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
864  static_cast<FieldApproximationBase>(b));
865  }
867 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:175
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:146
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407
continuous field
Definition: definitions.h:171
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 1024 of file ForcesAndSourcesCore.cpp.

1024  {
1026 
1027  const EntityType type = numeredEntFiniteElementPtr->getEntType();
1028  if (type == lastEvaluatedElementEntityType)
1030 
1031  // Data on elements for proper spaces
1032  for (int space = H1; space != LASTSPACE; ++space) {
1033  dataOnElement[space]->setElementType(type);
1034  derivedDataOnElement[space]->setElementType(type);
1035  }
1036 
1038 
1040 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:175
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:407
continuous field
Definition: definitions.h:171
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 256 of file ForcesAndSourcesCore.cpp.

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

389 { 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 446 of file ForcesAndSourcesCore.hpp.

448  {
449  return dataOnElement[space]->dataOnEntities[type][side];
450  }
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 461 of file ForcesAndSourcesCore.hpp.

461  {
462  return dataOnElement[space]->dataOnEntities[type][side];
463  }
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 829 of file ForcesAndSourcesCore.hpp.

831  {
832  return getEntityIndices(data, field_name,
833  const_cast<FENumeredDofEntity_multiIndex &>(
834  numeredEntFiniteElementPtr->getColsDofs()),
835  type_lo, type_hi);
836 }
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 148 of file ForcesAndSourcesCore.cpp.

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

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

510  {
511  return getEntityDataOrder(type, space, data.dataOnEntities[type]);
512  }
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 560 of file ForcesAndSourcesCore.cpp.

562  {
564  for (EntityType t = type_lo; t != type_hi; ++t) {
565  for (auto &dat : data.dataOnEntities[t]) {
566  dat.getDataOrder() = 0;
567  dat.getBase() = NOBASE;
568  dat.getSpace() = NOSPACE;
569  dat.getFieldData().resize(0, false);
570  dat.getFieldDofs().resize(0, false);
571  }
572  }
573 
574  auto &dofs = const_cast<FEDofEntity_multiIndex &>(
575  numeredEntFiniteElementPtr->getDataDofs());
576  auto &dofs_by_type = dofs.get<Composite_Name_And_Type_mi_tag>();
577  auto dit = dofs_by_type.lower_bound(boost::make_tuple(field_name, type_lo));
578  if (dit == dofs_by_type.end())
580  auto hi_dit =
581  dofs_by_type.lower_bound(boost::make_tuple(field_name, type_hi));
582  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
583  for (; dit != hi_dit;) {
584 
585  auto &dof = **dit;
586  const int nb_dofs_on_ent = dof.getNbDofsOnEnt();
587  if (nb_dofs_on_ent) {
588 
589  const EntityType type = dof.getEntType();
590  const int side = dof.sideNumberPtr->side_number;
591  auto &dat = data.dataOnEntities[type][side];
592 
593  auto &ent_field_dofs = dat.getFieldDofs();
594  auto &ent_field_data = dat.getFieldData();
595  const int brother_side = dof.sideNumberPtr->brother_side_number;
596  if (brother_side != -1)
597  brother_dofs_vec.emplace_back(*dit);
598 
599  if (ent_field_data.empty()) {
600 
601  dat.getBase() = dof.getApproxBase();
602  dat.getSpace() = dof.getSpace();
603  const int ent_order = dof.getMaxOrder();
604  dat.getDataOrder() =
605  dat.getDataOrder() > ent_order ? dat.getDataOrder() : ent_order;
606  ent_field_data.resize(nb_dofs_on_ent, false);
607  noalias(ent_field_data) = dof.getEntFieldData();
608  ent_field_dofs.resize(nb_dofs_on_ent, false);
609  for (int ii = 0; ii != nb_dofs_on_ent; ++ii) {
610  ent_field_dofs[ii] = *dit;
611  ++dit;
612  }
613  }
614  }
615  }
616 
617  for (auto &dof_ptr : brother_dofs_vec) {
618  if (auto d = dof_ptr.lock()) {
619  const EntityType type = d->getEntType();
620  const int side = d->sideNumberPtr->side_number;
621  const int brother_side = d->sideNumberPtr->brother_side_number;
622  auto &dat = data.dataOnEntities[type][side];
623  auto &dat_brother = data.dataOnEntities[type][brother_side];
624  dat_brother.getBase() = dat.getBase();
625  dat_brother.getSpace() = dat.getSpace();
626  dat_brother.getDataOrder() = dat.getDataOrder();
627  dat_brother.getFieldData() = dat.getFieldData();
628  dat_brother.getFieldDofs() = dat.getFieldDofs();
629  }
630  }
631 
633 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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:407

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

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

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

822  {
823  return getEntityIndices(data, field_name,
824  const_cast<FENumeredDofEntity_multiIndex &>(
825  numeredEntFiniteElementPtr->getRowsDofs()),
826  type_lo, type_hi);
827 }
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 
107  const int side_number = (*r.first)->side_number;
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  }
116 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:407

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

496  {
497  return getEntitySense(type, data.dataOnEntities[type]);
498  }
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 670 of file ForcesAndSourcesCore.cpp.

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

◆ getMaxColOrder()

int MoFEM::ForcesAndSourcesCore::getMaxColOrder ( ) const

Get max order of approximation for field in columns.

Definition at line 144 of file ForcesAndSourcesCore.cpp.

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

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

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

473  {
474 
475  auto get_nodes_field_data = [&](FEDofEntity_multiIndex &dofs,
476  VectorDouble &nodes_data,
477  VectorDofs &nodes_dofs, FieldSpace &space,
479  VectorInt &bb_node_order) {
481 
482  auto &dofs_by_name_and_type = dofs.get<Composite_Name_And_Type_mi_tag>();
483  auto tuple = boost::make_tuple(field_name, MBVERTEX);
484  auto dit = dofs_by_name_and_type.lower_bound(tuple);
485  if (dit == dofs_by_name_and_type.end())
486  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
487  "No nodal dofs on element");
488  auto hi_dit = dofs.get<Composite_Name_And_Type_mi_tag>().upper_bound(tuple);
489 
490  if (dit != hi_dit) {
491  auto &first_dof = **dit;
492  space = first_dof.getSpace();
493  base = first_dof.getApproxBase();
494  int num_nodes;
495  CHKERR getNumberOfNodes(num_nodes);
496  bb_node_order.resize(num_nodes, false);
497  bb_node_order.clear();
498  const int nb_dof_idx = first_dof.getNbOfCoeffs();
499  const int max_nb_dofs = nb_dof_idx * num_nodes;
500  nodes_data.resize(max_nb_dofs, false);
501  nodes_dofs.resize(max_nb_dofs, false);
502  nodes_data.clear();
503 
504  std::vector<boost::weak_ptr<FEDofEntity>> brother_dofs_vec;
505  for (; dit != hi_dit;) {
506  const auto &dof_ptr = *dit;
507  const auto &dof = *dof_ptr;
508  const auto &sn = *dof.sideNumberPtr;
509  const int side_number = sn.side_number;
510  const int brother_side_number = sn.brother_side_number;
511  if (brother_side_number != -1)
512  brother_dofs_vec.emplace_back(dof_ptr);
513 
514  bb_node_order[side_number] = dof.getMaxOrder();
515  int pos = side_number * nb_dof_idx;
516  auto ent_filed_data_vec = dof.getEntFieldData();
517  for (int ii = 0; ii != nb_dof_idx; ++ii) {
518  nodes_data[pos] = ent_filed_data_vec[ii];
519  nodes_dofs[pos] = *dit;
520  ++pos;
521  ++dit;
522  }
523  }
524 
525  for (auto &dof_ptr : brother_dofs_vec) {
526  if (const auto d = dof_ptr.lock()) {
527  const auto &sn = d->sideNumberPtr;
528  const int side_number = sn->side_number;
529  const int brother_side_number = sn->brother_side_number;
530  bb_node_order[brother_side_number] = bb_node_order[side_number];
531  int pos = side_number * nb_dof_idx;
532  int brother_pos = brother_side_number * nb_dof_idx;
533  for (int ii = 0; ii != nb_dof_idx; ++ii) {
534  nodes_data[brother_pos] = nodes_data[pos];
535  nodes_dofs[brother_pos] = nodes_dofs[pos];
536  ++pos;
537  ++brother_pos;
538  }
539  }
540  }
541 
542  } else {
543  nodes_data.resize(0, false);
544  nodes_dofs.resize(0, false);
545  }
546 
548  };
549 
550  return get_nodes_field_data(
551  const_cast<FEDofEntity_multiIndex &>(
552  numeredEntFiniteElementPtr->getDataDofs()),
553  data.dataOnEntities[MBVERTEX][0].getFieldData(),
554  data.dataOnEntities[MBVERTEX][0].getFieldDofs(),
555  data.dataOnEntities[MBVERTEX][0].getSpace(),
556  data.dataOnEntities[MBVERTEX][0].getBase(),
557  data.dataOnEntities[MBVERTEX][0].getBBNodeOrder());
558 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:144
FieldSpace
approximation spaces
Definition: definitions.h:168
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407
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 196 of file ForcesAndSourcesCore.cpp.

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

◆ getNoFieldColIndices()

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

get col NoField indices

Definition at line 347 of file ForcesAndSourcesCore.cpp.

348  {
350  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
351  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
352  }
353  CHKERR getNoFieldIndices(field_name,
354  const_cast<FENumeredDofEntity_multiIndex &>(
355  numeredEntFiniteElementPtr->getColsDofs()),
356  data.dataOnEntities[MBENTITYSET][0].getIndices());
358 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407

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

637  {
639  auto dit = dofs.get<FieldName_mi_tag>().lower_bound(field_name);
640  auto hi_dit = dofs.get<FieldName_mi_tag>().upper_bound(field_name);
641  int size = std::distance(dit, hi_dit);
642  ent_field_data.resize(size, false);
643  ent_field_dofs.resize(size, false);
644  for (; dit != hi_dit; dit++) {
645  int idx = (*dit)->getDofCoeffIdx();
646  ent_field_data[idx] = (*dit)->getFieldData();
647  ent_field_dofs[idx] = *dit;
648  }
650 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508

◆ getNoFieldFieldData() [2/2]

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

Definition at line 652 of file ForcesAndSourcesCore.cpp.

653  {
655  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
656  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
657  }
659  field_name,
660  const_cast<FEDofEntity_multiIndex &>(
661  numeredEntFiniteElementPtr->getDataDofs()),
662  data.dataOnEntities[MBENTITYSET][0].getFieldData(),
663  data.dataOnEntities[MBENTITYSET][0].getFieldDofs());
665 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getNoFieldIndices()

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

get NoField indices

Definition at line 320 of file ForcesAndSourcesCore.cpp.

322  {
324  auto dit = dofs.get<FieldName_mi_tag>().lower_bound(field_name);
325  auto hi_dit = dofs.get<FieldName_mi_tag>().upper_bound(field_name);
326  indices.resize(std::distance(dit, hi_dit));
327  for (; dit != hi_dit; dit++) {
328  int idx = (*dit)->getPetscGlobalDofIdx();
329  indices[(*dit)->getDofCoeffIdx()] = idx;
330  }
332 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508

◆ getNoFieldRowIndices()

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

get col NoField indices

Definition at line 334 of file ForcesAndSourcesCore.cpp.

335  {
337  if (data.dataOnEntities[MBENTITYSET].size() == 0) {
338  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
339  }
340  CHKERR getNoFieldIndices(field_name,
341  const_cast<FENumeredDofEntity_multiIndex &>(
342  numeredEntFiniteElementPtr->getRowsDofs()),
343  data.dataOnEntities[MBENTITYSET][0].getIndices());
345 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407

◆ 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:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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 382 of file ForcesAndSourcesCore.hpp.

382 { 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 455 of file ForcesAndSourcesCore.cpp.

456  {
457  return getProblemNodesIndices(field_name, *(problemPtr->numeredDofsCols),
458  nodes_indices);
459 }
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 362 of file ForcesAndSourcesCore.cpp.

364  {
366  nodes_indices.resize(0);
367  auto &side_table = const_cast<SideNumber_multiIndex &>(
368  numeredEntFiniteElementPtr->getSideNumberTable());
369  auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 0));
370  auto hi_siit =
371  side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 10000));
372 
373  int nn = 0;
374  for (; siit != hi_siit; siit++, nn++) {
375 
376  if (siit->get()->side_number == -1)
377  continue;
378 
379  const EntityHandle ent = siit->get()->ent;
380  auto dit =
381  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().lower_bound(
382  boost::make_tuple(field_name, ent, 0));
383  auto hi_dit =
384  dofs.get<Composite_Name_And_Ent_And_EntDofIdx_mi_tag>().upper_bound(
385  boost::make_tuple(field_name, ent, 10000)); /// very large number
386  if (dit != hi_dit) {
387  if (!nn) {
388  nodes_indices.resize((*dit)->getNbOfCoeffs() *
389  std::distance(siit, hi_siit));
390  }
391  for (; dit != hi_dit; dit++) {
392  nodes_indices[siit->get()->side_number * (*dit)->getNbOfCoeffs() +
393  (*dit)->getDofCoeffIdx()] =
394  (*dit)->getPetscGlobalDofIdx();
395  }
396  }
397  }
398 
400 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getProblemNodesRowIndices()

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

Definition at line 441 of file ForcesAndSourcesCore.cpp.

442  {
443  return getProblemNodesIndices(field_name, *(problemPtr->numeredDofsRows),
444  nodes_indices);
445 }
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 462 of file ForcesAndSourcesCore.cpp.

464  {
465  return getProblemTypeIndices(field_name, *(problemPtr->numeredDofsCols), type,
466  side_number, indices);
467 }
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 402 of file ForcesAndSourcesCore.cpp.

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

450  {
451  return getProblemTypeIndices(field_name, *(problemPtr->numeredDofsRows), type,
452  side_number, indices);
453 }
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 246 of file ForcesAndSourcesCore.cpp.

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

1365  {
1366  return getRuleHook ? getRuleHook(order_row, order_col, order_data)
1367  : getRule(order_data);
1368 }
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 763 of file ForcesAndSourcesCore.cpp.

764  {
766 
767  if (nInTheLoop == 0) {
768  data.sPace.reset();
769  data.bAse.reset();
770  for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
771  data.spacesOnEntities[t].reset();
772  data.basesOnEntities[t].reset();
773  }
774  for (int s = 0; s != LASTSPACE; ++s) {
775  data.basesOnSpaces[s].reset();
776  }
777  }
778 
779  if (dataFieldEntsPtr)
780  for (auto e : *dataFieldEntsPtr) {
781  // get data from entity
782  const EntityType type = e->getEntType();
783  const FieldSpace space = e->getSpace();
784  const FieldApproximationBase approx = e->getApproxBase();
785 
786  // set data
787  data.sPace.set(space);
788  data.bAse.set(approx);
789  data.spacesOnEntities[type].set(space);
790  data.basesOnEntities[type].set(approx);
791  data.basesOnSpaces[space].set(approx);
792  }
793  else
794  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
795  "data fields ents not allocated on element");
796 
798 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:175
FieldApproximationBase
approximation base
Definition: definitions.h:144
FieldSpace
approximation spaces
Definition: definitions.h:168
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 396 of file ForcesAndSourcesCore.hpp.

396 { 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 1072 of file ForcesAndSourcesCore.cpp.

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

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

Definition at line 1410 of file ForcesAndSourcesCore.cpp.

1410  {
1412  if (operatorHook) {
1413  ierr = operatorHook();
1414  CHKERRG(ierr);
1415  }
1417 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:544
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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 1418 of file ForcesAndSourcesCore.cpp.

1418  {
1420  if (postProcessHook) {
1421  ierr = postProcessHook();
1422  CHKERRG(ierr);
1423  }
1425 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:544
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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, KelvinVoigtDamper::DamperFE, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, Smoother::MyVolumeFE, and NonlinearElasticElement::MyVolumeFE.

Definition at line 1402 of file ForcesAndSourcesCore.cpp.

1402  {
1404  if (preProcessHook) {
1405  ierr = preProcessHook();
1406  CHKERRG(ierr);
1407  }
1409 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:544
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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 1370 of file ForcesAndSourcesCore.cpp.

1371  {
1372  return setRuleHook ? setRuleHook(order_row, order_col, order_data)
1373  : setGaussPts(order_data);
1374 }
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 1380 of file ForcesAndSourcesCore.cpp.

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

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

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

Friends And Related Function Documentation

◆ FaceElementForcesAndSourcesCoreOnSideBase

Definition at line 814 of file ForcesAndSourcesCore.hpp.

◆ UserDataOperator

friend class UserDataOperator
friend

◆ VolumeElementForcesAndSourcesCoreOnSideBase

Definition at line 813 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ dataH1

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataH1
protected

Definition at line 764 of file ForcesAndSourcesCore.hpp.

◆ dataHcurl

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataHcurl
protected

Definition at line 765 of file ForcesAndSourcesCore.hpp.

◆ dataHdiv

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataHdiv
protected

Definition at line 766 of file ForcesAndSourcesCore.hpp.

◆ dataL2

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataL2
protected

Definition at line 767 of file ForcesAndSourcesCore.hpp.

◆ dataNoField

DataForcesAndSourcesCore& MoFEM::ForcesAndSourcesCore::dataNoField
protected

Definition at line 763 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 753 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 761 of file ForcesAndSourcesCore.hpp.

◆ elementPolynomialBasePtr

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

Pointer to entity polynomial base.

Definition at line 789 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 408 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 782 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 773 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 800 of file ForcesAndSourcesCore.hpp.

◆ userPolynomialBasePtr

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

Pointer to user polynomail base.

Definition at line 794 of file ForcesAndSourcesCore.hpp.


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