11#ifndef __FORCES_AND_SOURCES_CORE__HPP__
12#define __FORCES_AND_SOURCES_CORE__HPP__
14using namespace boost::numeric;
27 typedef boost::function<int(
int order_row,
int order_col,
int order_data)>
31 int order_row,
int order_col,
178 boost::ptr_vector<EntitiesFieldData::EntData> &data)
const;
190 boost::ptr_vector<EntitiesFieldData::EntData> &data)
const;
199 template <EntityType type>
212 template <EntityType type>
223 template <
typename EXTRACTOR>
227 VectorInt &local_nodes_indices, EXTRACTOR &&extractor)
const;
231 const int bit_number)
const;
235 const int bit_number)
const;
237 template <
typename EXTRACTOR>
242 EXTRACTOR &&extractor)
const;
247 const EntityType type_hi = MBPOLYHEDRON)
const;
252 const EntityType type_hi = MBPOLYHEDRON)
const;
257 boost::shared_ptr<FENumeredDofEntity_multiIndex> dofs,
262 const int bit_number)
const;
266 const int bit_number)
const;
287 const int bit_number)
const;
296 const int bit_number)
const;
301 const EntityType type_hi = MBPOLYHEDRON)
const;
379 virtual int getRule(
int order_row,
int order_col,
int order_data);
460 const std::array<boost::shared_ptr<EntitiesFieldData>,
LASTSPACE>
467 const std::array<boost::shared_ptr<EntitiesFieldData>,
LASTSPACE>
586 const bool symm =
true);
589 const bool symm =
true);
592 const std::string col_field_name,
const char type,
593 const bool symm =
true);
597 inline boost::shared_ptr<const NumeredEntFiniteElement>
791 inline Vec
getTSu()
const;
797 inline Vec
getTSf()
const;
799 inline Mat
getTSA()
const;
801 inline Mat
getTSB()
const;
809 inline double getTSa()
const;
901 std::vector<boost::weak_ptr<NumeredEntFiniteElement>>>;
919 const int verb =
QUIET,
937 const int verb =
QUIET,
952 const int verb =
QUIET,
967 const int verb =
QUIET,
994boost::shared_ptr<const NumeredEntFiniteElement>
995ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr()
const {
999EntityHandle ForcesAndSourcesCore::UserDataOperator::getFEEntityHandle()
const {
1000 return getNumeredEntFiniteElementPtr()->getEnt();
1003int ForcesAndSourcesCore::UserDataOperator::getFEDim()
const {
1007EntityType ForcesAndSourcesCore::UserDataOperator::getFEType()
const {
1011boost::weak_ptr<SideNumber>
1012ForcesAndSourcesCore::UserDataOperator::getSideNumberPtr(
1013 const int side_number,
const EntityType type) {
1014 auto &side_table_by_side_and_type =
1015 ptrFE->numeredEntFiniteElementPtr->getSideNumberTable().get<1>();
1017 side_table_by_side_and_type.find(boost::make_tuple(type, side_number));
1018 if (side_it != side_table_by_side_and_type.end())
1021 return boost::weak_ptr<SideNumber>();
1025ForcesAndSourcesCore::UserDataOperator::getSideEntity(
const int side_number,
1027 if (
auto side_ptr = getSideNumberPtr(side_number, type).lock())
1028 return side_ptr->ent;
1033int ForcesAndSourcesCore::UserDataOperator::getNumberOfNodesOnElement()
const {
1034 return ptrFE->getNumberOfNodes();
1037const FEMethod *ForcesAndSourcesCore::UserDataOperator::getFEMethod()
const {
1041int ForcesAndSourcesCore::UserDataOperator::getOpType()
const {
return opType; }
1043void ForcesAndSourcesCore::UserDataOperator::setOpType(
const OpType type) {
1047void ForcesAndSourcesCore::UserDataOperator::addOpType(
const OpType type) {
1051int ForcesAndSourcesCore::UserDataOperator::getNinTheLoop()
const {
1052 return getFEMethod()->getNinTheLoop();
1055int ForcesAndSourcesCore::UserDataOperator::getLoopSize()
const {
1056 return getFEMethod()->getLoopSize();
1059std::string ForcesAndSourcesCore::UserDataOperator::getFEName()
const {
1060 return getFEMethod()->getFEName();
1064ForcesAndSourcesCore::UserDataOperator::getDataCtx()
const {
1065 return getFEMethod()->data_ctx;
1069ForcesAndSourcesCore::UserDataOperator::getKSPCtx()
const {
1070 return getFEMethod()->ksp_ctx;
1074ForcesAndSourcesCore::UserDataOperator::getSNESCtx()
const {
1075 return getFEMethod()->snes_ctx;
1079ForcesAndSourcesCore::UserDataOperator::getTSCtx()
const {
1080 return getFEMethod()->ts_ctx;
1083Vec ForcesAndSourcesCore::UserDataOperator::getKSPf()
const {
1085 if (getFEMethod()->
ksp_f == PETSC_NULL)
1088 return getFEMethod()->ksp_f;
1091Mat ForcesAndSourcesCore::UserDataOperator::getKSPA()
const {
1093 if (getFEMethod()->
ksp_A == PETSC_NULL)
1096 return getFEMethod()->ksp_A;
1099Mat ForcesAndSourcesCore::UserDataOperator::getKSPB()
const {
1101 if (getFEMethod()->
ksp_B == PETSC_NULL)
1104 return getFEMethod()->ksp_B;
1107Vec ForcesAndSourcesCore::UserDataOperator::getSNESf()
const {
1109 if (getFEMethod()->
snes_f == PETSC_NULL)
1112 return getFEMethod()->snes_f;
1115Vec ForcesAndSourcesCore::UserDataOperator::getSNESx()
const {
1117 if (getFEMethod()->
snes_x == PETSC_NULL)
1120 return getFEMethod()->snes_x;
1123Mat ForcesAndSourcesCore::UserDataOperator::getSNESA()
const {
1125 if (getFEMethod()->
snes_A == PETSC_NULL)
1128 return getFEMethod()->snes_A;
1131Mat ForcesAndSourcesCore::UserDataOperator::getSNESB()
const {
1133 if (getFEMethod()->
snes_B == PETSC_NULL)
1136 return getFEMethod()->snes_B;
1139Vec ForcesAndSourcesCore::UserDataOperator::getTSu()
const {
1141 if (getFEMethod()->
ts_u == PETSC_NULL)
1144 return getFEMethod()->ts_u;
1147Vec ForcesAndSourcesCore::UserDataOperator::getTSu_t()
const {
1149 if (getFEMethod()->
ts_u_t == PETSC_NULL)
1152 return getFEMethod()->ts_u_t;
1155Vec ForcesAndSourcesCore::UserDataOperator::getTSu_tt()
const {
1157 if (getFEMethod()->
ts_u_tt == PETSC_NULL)
1160 return getFEMethod()->ts_u_tt;
1163Vec ForcesAndSourcesCore::UserDataOperator::getTSf()
const {
1165 if (getFEMethod()->
ts_F == PETSC_NULL)
1168 return getFEMethod()->ts_F;
1171Mat ForcesAndSourcesCore::UserDataOperator::getTSA()
const {
1173 if (getFEMethod()->
ts_A == PETSC_NULL)
1176 return getFEMethod()->ts_A;
1179Mat ForcesAndSourcesCore::UserDataOperator::getTSB()
const {
1181 if (getFEMethod()->
ts_B == PETSC_NULL)
1184 return getFEMethod()->ts_B;
1187int ForcesAndSourcesCore::UserDataOperator::getTSstep()
const {
1189 if ((getFEMethod()->
data_ctx & PetscData::PetscData::CtxSetTime).none())
1192 return getFEMethod()->ts_step;
1195double ForcesAndSourcesCore::UserDataOperator::getTStime()
const {
1197 if ((getFEMethod()->
data_ctx & PetscData::PetscData::CtxSetTime).none())
1200 return getFEMethod()->ts_t;
1203double ForcesAndSourcesCore::UserDataOperator::getTStimeStep()
const {
1205 if ((getFEMethod()->
data_ctx & PetscData::PetscData::CtxSetTime).none())
1208 return getFEMethod()->ts_dt;
1211double ForcesAndSourcesCore::UserDataOperator::getTSa()
const {
1218 return getFEMethod()->ts_a;
1221double ForcesAndSourcesCore::UserDataOperator::getTSaa()
const {
1228 return getFEMethod()->ts_aa;
1235auto ForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight() {
1237 &(getGaussPts()(getGaussPts().size1() - 1, 0)));
1251ForcesAndSourcesCore::UserDataOperator::getSidePtrFE()
const {
1256ForcesAndSourcesCore::UserDataOperator::getRefinePtrFE()
const {
1260MatrixDouble &ForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts() {
1264auto ForcesAndSourcesCore::UserDataOperator::getFTensor1CoordsAtGaussPts() {
1266 &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
1267 &getCoordsAtGaussPts()(0, 2));
1270double ForcesAndSourcesCore::UserDataOperator::getMeasure()
const {
1274double &ForcesAndSourcesCore::UserDataOperator::getMeasure() {
1284template <
typename E>
1299 boost::shared_ptr<AdjCache> adj_cache =
nullptr)
1301 fieldName(fe_name), sideDim(side_dim), sevLevel(sev),
1302 adjCache(adj_cache) {}
1307 CHKERR loopSide(fieldName, sideFEPtr.get(), sideDim, 0,
VERBOSE, sevLevel,
1313 return sideFEPtr->getOpPtrVector();
FieldApproximationBase
approximation base
FieldSpace
approximation spaces
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
multi_index_container< boost::shared_ptr< NumeredDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, UId, &NumeredDofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredDofEntity, unsigned int, &NumeredDofEntity::pArt > >, ordered_non_unique< tag< Idx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::dofIdx > >, ordered_non_unique< tag< PetscGlobalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscGloablDofIdx > >, ordered_non_unique< tag< PetscLocalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscLocalDofIdx > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, EntityHandle, &NumeredDofEntity::getEnt > > > > NumeredDofEntity_multiIndex
MultiIndex container keeps NumeredDofEntity.
SeverityLevel
Severity levels.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasVector< int > VectorInt
implementation of Data Operators for Forces and Sources
auto type_from_handle(const EntityHandle h)
get type from entity handle
ublas::vector< FieldEntity *, FieldEntAllocator > VectorFieldEntities
DEPRECATED typedef ForcesAndSourcesCore ForcesAndSurcesCore
auto dimension_from_handle(const EntityHandle h)
get entity dimension form handle
std::vector< boost::weak_ptr< FieldEntity > > FieldEntity_vector_view
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
constexpr auto field_name
base operator to do operations at Gauss Pt. level
Deprecated interface functions.
Base face element used to integrate on skeleton.
Data on single entity (This is passed as argument to DataOperator::doWork)
data structure for finite element entity
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
structure for User Loop Methods on finite elements
EntityHandle getFEEntityHandle() const
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Base face element used to integrate on skeleton.
Base face element used to integrate on skeleton.
void setOpType(const OpType type)
Set operator type.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
MoFEMErrorCode getProblemColIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get col indices.
int getLoopSize() const
get size of elements in the loop
EntityType getFEType() const
Get dimension of finite element.
const KspMethod::KSPContext getKSPCtx() const
static const char *const OpTypeNames[]
int getNinTheLoop() const
get number of finite element in the loop
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
MoFEMErrorCode loopThis(const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
User call this function to loop over parent elements. This function calls finite element with is oper...
int getNumberOfNodesOnElement() const
Get the number of nodes on finite element.
boost::weak_ptr< SideNumber > getSideNumberPtr(const int side_number, const EntityType type)
Get the side number pointer.
const SnesMethod::SNESContext getSNESCtx() const
MoFEMErrorCode loopParent(const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
User call this function to loop over parent elements. This function calls finite element with is oper...
const PetscData::Switches & getDataCtx() const
ForcesAndSourcesCore * ptrFE
EntityHandle getSideEntity(const int side_number, const EntityType type)
Get the side entity.
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
auto getFTensor0IntegrationWeight()
Get integration weights.
const TSMethod::TSContext getTSCtx() const
double getMeasure() const
get measure of element
std::string getFEName() const
Get name of the element.
ForcesAndSourcesCore * getRefinePtrFE() const
MoFEMErrorCode getProblemRowIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get row indices.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > > AdjCache
double getTStimeStep() const
void addOpType(const OpType type)
Add operator type.
ForcesAndSourcesCore * getPtrFE() const
int getOpType() const
Get operator types.
OpType
Controls loop over entities on element.
@ OPCOL
operator doWork function is executed on FE columns
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
@ OPSPACE
operator do Work is execute on space data
MoFEMErrorCode loopChildren(const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
User call this function to loop over parent elements. This function calls finite element with is oper...
virtual MoFEMErrorCode setPtrFE(ForcesAndSourcesCore *ptr)
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
ForcesAndSourcesCore * getSidePtrFE() const
MoFEMErrorCode loopSide(const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
User call this function to loop over elements on the side of face. This function calls finite element...
int getFEDim() const
Get dimension of finite element.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
structure to get information form mofem into EntitiesFieldData
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
MoFEMErrorCode getProblemNodesRowIndices(const std::string &field_name, VectorInt &nodes_indices) const
EntitiesFieldData & dataNoField
auto & getDataOnElementBySpaceArray()
Get data on entities and space.
MoFEMErrorCode getEntityIndices(EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
MoFEMErrorCode getEntityRowIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getRowNodesIndices(EntitiesFieldData &data, const int bit_number) const
get row node indices from FENumeredDofEntity_multiIndex
int getMaxRowOrder() const
Get max order of approximation for field in rows.
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
MoFEMErrorCode getEntitySense(EntitiesFieldData &data) const
Get the entity sense (orientation)
EntitiesFieldData & dataHdiv
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
MoFEMErrorCode getProblemTypeColIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
MoFEMErrorCode getNoFieldIndices(const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const
get NoField indices
virtual MoFEMErrorCode operator()()
function is run for every finite element
MoFEMErrorCode getFaceNodes(EntitiesFieldData &data) const
Get nodes on faces.
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
MoFEMErrorCode getNoFieldRowIndices(EntitiesFieldData &data, const int bit_number) const
get col NoField indices
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
MoFEMErrorCode getNodesIndices(const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
get node indices
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
EntitiesFieldData & dataHcurl
MoFEMErrorCode setRefineFEPtr(const ForcesAndSourcesCore *refine_fe_ptr)
Set the pointer to face element refined.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
EntitiesFieldData & dataH1
MatrixDouble coordsAtGaussPts
coordinated at gauss points
boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> GaussHookFun
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
boost::ptr_deque< UserDataOperator > opPtrVector
Vector of finite element users data operators.
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
auto & getDerivedDataOnElementBySpaceArray()
Get derived data on entities and space.
MoFEMErrorCode getProblemNodesColIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
MatrixDouble gaussPts
Matrix of integration points.
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
auto & getEntData(const FieldSpace space, const EntityType type, const int side)
Get the entity data.
MoFEMErrorCode getEntityDataOrder(EntitiesFieldData &data, const FieldSpace space) const
Get the entity data order for given space.
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
MoFEMErrorCode getColNodesIndices(EntitiesFieldData &data, const int bit_number) const
get col node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode setSideFEPtr(const ForcesAndSourcesCore *side_fe_ptr)
Set the pointer to face element on the side.
int getMaxColOrder() const
Get max order of approximation for field in columns.
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
MoFEMErrorCode getProblemNodesIndices(const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
get indices of nodal indices which are declared for problem but not this particular element
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
GaussHookFun setRuleHook
Set function to calculate integration rule.
MoFEMErrorCode getBitRefLevelOnData()
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
MoFEMErrorCode getEntityColIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
boost::shared_ptr< BaseFunction > elementPolynomialBasePtr
Pointer to entity polynomial base.
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.
ForcesAndSourcesCore * refinePtrFE
Element to integrate parent or child.
MoFEMErrorCode getNoFieldColIndices(EntitiesFieldData &data, const int bit_number) const
get col NoField indices
int getMaxDataOrder() const
Get max order of approximation for data fields.
MoFEMErrorCode getNoFieldFieldData(const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const
Get field data on nodes.
boost::shared_ptr< BaseFunction > userPolynomialBasePtr
Pointer to user polynomail base.
MoFEMErrorCode getProblemTypeRowIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
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
RuleHookFun getRuleHook
Hook to get rule.
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
EntitiesFieldData & dataL2
KSPContext
pass information about context of KSP/DM for with finite element is computed
Element used to execute operators on side of the element.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< E > & getSideFEPtr()
boost::shared_ptr< AdjCache > adjCache
boost::shared_ptr< E > sideFEPtr
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
const std::string fieldName
const LogManager::SeverityLevel sevLevel
OpLoopSide(MoFEM::Interface &m_field, const std::string fe_name, const int side_dim, const LogManager::SeverityLevel sev=Sev::noisy, boost::shared_ptr< AdjCache > adj_cache=nullptr)
Construct a new Op Loop Side object.
static constexpr Switches CtxSetA
static constexpr Switches CtxSetX_TT
std::bitset< 8 > Switches
static constexpr Switches CtxSetX_T
static constexpr Switches CtxSetB
Mat & snes_B
preconditioner of jacobian matrix
Mat & snes_A
jacobian matrix
Vec & ts_F
residual vector
Vec & ts_u_tt
second time derivative of state vector
Vec & ts_u_t
time derivative of state vector
Mat & ts_B
Preconditioner for ts_A.
Base volume element used to integrate on skeleton.