v0.15.5
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
MoFEM::ForcesAndSourcesCore::UserDataOperator Struct Reference

#include "src/finite_elements/ForcesAndSourcesCore.hpp"

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

Public Types

enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 

Public Member Functions

 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 Constructor for operators working on finite element spaces.
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 Constructor for operators working on a single field.
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 Constructor for operators working on two fields (bilinear forms)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
Accessing KSP
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Accessing SNES
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Accessing TS
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
Base functions and integration points
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
Coordinates and access to internal data
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
Measures (area, volume, length, etc.)
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
Database
MoFEM::InterfacegetMField ()
 
moab::Interface & getMoab ()
 
Pipelines
virtual boost::weak_ptr< ForcesAndSourcesCoregetSubPipelinePtr () const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side.
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side.
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Public Attributes

char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 

Static Public Attributes

static const char *const OpTypeNames []
 

Protected Member Functions

virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Protected Attributes

ForcesAndSourcesCoreptrFE
 

Friends

class ForcesAndSourcesCore
 
class EdgeElementForcesAndSourcesCore
 
class FaceElementForcesAndSourcesCore
 
class ContactPrismElementForcesAndSourcesCore
 

Loops

}

using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 

Detailed Description

Examples
mofem/atom_tests/hanging_node_approx.cpp, mofem/atom_tests/schur_test_diag_mat.cpp, mofem/atom_tests/simple_interface.cpp, mofem/atom_tests/simple_l2_only.cpp, mofem/tutorials/adv-0_plasticity/plastic.cpp, mofem/tutorials/adv-3_level_set/level_set.cpp, mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp, mofem/tutorials/fun-0_hello_world/hello_world.cpp, mofem/tutorials/vec-5_free_surface/free_surface.cpp, mofem/tutorials/vec-7_shape_optimisation/adjoint.cpp, mofem/users_modules/adolc-plasticity/src/ADOLCPlasticityLargeStrain.hpp, plastic.cpp, and thermoplastic.cpp.

Definition at line 584 of file ForcesAndSourcesCore.hpp.

Member Typedef Documentation

◆ AdjCache

Definition at line 981 of file ForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ OpType

Controls loop over entities on element.

  • OPROW is used if row vector is assembled
  • OPCOL is usually used if column vector is assembled
  • OPROWCOL is usually used for assemble matrices.
  • OPSPACE no field is defined for such operator. Is usually used to modify base

For typical problem like Bubnov-Galerkin OPROW and OPCOL are the same. In more general case for example for non-square matrices columns and rows could have different numeration and/or different set of fields.

Enumerator
OPROW 

operator doWork function is executed on FE rows

OPCOL 

operator doWork function is executed on FE columns

OPROWCOL 

operator doWork is executed on FE rows &columns

OPSPACE 

operator do Work is execute on space data

OPLAST 
Deprecated:
would be removed

Definition at line 601 of file ForcesAndSourcesCore.hpp.

601 {
602 OPROW = 1 << 0, ///< operator doWork function is executed on FE rows
603 OPCOL = 1 << 1, ///< operator doWork function is executed on FE columns
604 OPROWCOL = 1 << 2, ///< operator doWork is executed on FE rows &columns
605 OPSPACE = 1 << 3, ///< operator do Work is execute on space data
606 OPLAST = 1 << 3 ///< @deprecated would be removed
607 };
@ 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

Constructor & Destructor Documentation

◆ UserDataOperator() [1/3]

ForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const FieldSpace  space,
const char  type = OPSPACE,
const bool  symm = true 
)

Constructor for operators working on finite element spaces.

This constructor is used when operators modify basis functions on approximation spaces. The operator runs for all data on the specified space without access to specific field data. Commonly used for:

  • Coordinate transformations (Piola transforms)
  • Integration weight modifications
  • Basis function scaling or rotation
Parameters
spaceField space to operate on (H1, HDIV, HCURL, L2, etc.)
typeOperator type (OPSPACE, OPROW, OPCOL, OPROWCOL)
symmWhether the operator maintains matrix symmetry (default: true)

Definition at line 2167 of file ForcesAndSourcesCore.cpp.

◆ UserDataOperator() [2/3]

ForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string  field_name,
const char  type,
const bool  symm = true 
)

Constructor for operators working on a single field.

This constructor creates operators that work with data from a specific field. Used for operations like:

  • Calculating field values at integration points
  • Computing field gradients or derivatives
  • Applying boundary conditions to specific fields
  • Post-processing field data
Parameters
field_nameName of the field to operate on
typeOperator type:
  • OPROW: Works with test functions (RHS assembly)
  • OPCOL: Works with trial functions (geometry setup)
  • OPROWCOL: Works with both (matrix assembly)
symmWhether the operator maintains matrix symmetry (default: true)

Definition at line 2172 of file ForcesAndSourcesCore.cpp.

2174 : DataOperator(symm), opType(type), rowFieldName(field_name),
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
constexpr auto field_name

◆ UserDataOperator() [3/3]

ForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string  row_field_name,
const std::string  col_field_name,
const char  type,
const bool  symm = true 
)

Constructor for operators working on two fields (bilinear forms)

This constructor creates operators that couple two different fields, essential for:

  • Mixed formulations (pressure-velocity coupling)
  • Multi-physics problems (thermal-mechanical coupling)
  • Contact mechanics (displacement-Lagrange multiplier)
  • Interface problems between different field types
Parameters
row_field_nameName of the field for test functions (rows)
col_field_nameName of the field for trial functions (columns)
typeOperator type (typically OPROWCOL for coupling)
symmWhether the operator maintains matrix symmetry (default: true) Note: Coupling between different fields is often non-symmetric

Definition at line 2177 of file ForcesAndSourcesCore.cpp.

2180 : DataOperator(symm), opType(type), rowFieldName(row_field_name),
2181 colFieldName(col_field_name), sPace(LASTSPACE), ptrFE(nullptr) {}

Member Function Documentation

◆ addOpType()

void ForcesAndSourcesCore::UserDataOperator::addOpType ( const OpType  type)
inline

Add operator type.

Definition at line 1166 of file ForcesAndSourcesCore.hpp.

1166 {
1167 opType |= type;
1168}

◆ getCoordsAtGaussPts()

MatrixDouble & ForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts ( )
inline

◆ getDataCtx()

const PetscData::Switches & ForcesAndSourcesCore::UserDataOperator::getDataCtx ( ) const
inline

Definition at line 1183 of file ForcesAndSourcesCore.hpp.

1183 {
1184 return getFEMethod()->data_ctx;
1185}
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Switches data_ctx
Current data context switches.

◆ getFEDim()

int ForcesAndSourcesCore::UserDataOperator::getFEDim ( ) const
inline

Get dimension of finite element.

Returns
int

Definition at line 1122 of file ForcesAndSourcesCore.hpp.

1122 {
1124};
auto dimension_from_handle(const EntityHandle h)
get entity dimension form handle
EntityHandle getFEEntityHandle() const
Return finite element entity handle.

◆ getFEEntityHandle()

EntityHandle ForcesAndSourcesCore::UserDataOperator::getFEEntityHandle ( ) const
inline

Return finite element entity handle.

Returns
Finite element entity handle
Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1118 of file ForcesAndSourcesCore.hpp.

1118 {
1119 return getNumeredEntFiniteElementPtr()->getEnt();
1120}
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.

◆ getFEMethod()

const FEMethod * ForcesAndSourcesCore::UserDataOperator::getFEMethod ( ) const
inline

◆ getFEName()

std::string ForcesAndSourcesCore::UserDataOperator::getFEName ( ) const
inline

Get name of the element.

Definition at line 1178 of file ForcesAndSourcesCore.hpp.

1178 {
1179 return getFEMethod()->getFEName();
1180}
auto getFEName() const
Get the name of the current finite element.

◆ getFEType()

EntityType ForcesAndSourcesCore::UserDataOperator::getFEType ( ) const
inline

Get dimension of finite element.

Returns
int

Definition at line 1126 of file ForcesAndSourcesCore.hpp.

1126 {
1128};
auto type_from_handle(const EntityHandle h)
get type from entity handle

◆ getFTensor0IntegrationWeight()

auto ForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight ( )
inline

◆ getFTensor1CoordsAtGaussPts()

auto ForcesAndSourcesCore::UserDataOperator::getFTensor1CoordsAtGaussPts ( )
inline

Get coordinates at integration points assuming linear geometry.

auto t_coords = getFTensor1CoordsAtGaussPts();
for(int gg = 0;gg!=nb_int_ptrs;gg++) {
// do something
++t_coords;
}
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
Examples
mofem/atom_tests/prism_polynomial_approximation.cpp, mofem/atom_tests/quad_polynomial_approximation.cpp, mofem/atom_tests/simple_interface.cpp, and mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1377 of file ForcesAndSourcesCore.hpp.

1377 {
1379 &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
1380 &getCoordsAtGaussPts()(0, 2));
1381}
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)

◆ getGaussPts()

MatrixDouble & ForcesAndSourcesCore::UserDataOperator::getGaussPts ( )
inline

matrix of integration (Gauss) points for Volume Element

For triangle: columns 0,1 are x,y coordinates respectively and column 2 is a weight value for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular volume element

For tetrahedron: columns 0,1,2 are x,y,z coordinates respectively and column 3 is a weight value for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular volume element

Examples
mofem/atom_tests/continuity_check_on_skeleton_3d.cpp, mofem/atom_tests/hdiv_check_approx_in_3d.cpp, mofem/atom_tests/simple_interface.cpp, mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp, mofem/tutorials/cor-9_reaction_diffusion/reaction_diffusion.cpp, mofem/tutorials/fun-1_integration/integration.cpp, mofem/tutorials/max-0_magnetostatics/src/MagneticElement.hpp, and mofem/tutorials/mix-1_light_intensity_equation/phase.cpp.

Definition at line 1350 of file ForcesAndSourcesCore.hpp.

1350 {
1351 return static_cast<ForcesAndSourcesCore *>(ptrFE)->gaussPts;
1352}
MatrixDouble gaussPts
Matrix of integration points.

◆ getKSPA()

Mat ForcesAndSourcesCore::UserDataOperator::getKSPA ( ) const
inline

Definition at line 1210 of file ForcesAndSourcesCore.hpp.

1210 {
1211#ifndef NDEBUG
1212 if (getFEMethod()->ksp_A == PETSC_NULLPTR)
1213 THROW_MESSAGE("KSP not set A vector");
1214#endif
1215 return getFEMethod()->ksp_A;
1216}
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Mat & ksp_A
Reference to system matrix in KSP context.

◆ getKSPB()

Mat ForcesAndSourcesCore::UserDataOperator::getKSPB ( ) const
inline

Definition at line 1218 of file ForcesAndSourcesCore.hpp.

1218 {
1219#ifndef NDEBUG
1220 if (getFEMethod()->ksp_B == PETSC_NULLPTR)
1221 THROW_MESSAGE("KSP not set B vector");
1222#endif
1223 return getFEMethod()->ksp_B;
1224}
Mat & ksp_B
Reference to preconditioner matrix in KSP context.

◆ getKSPCtx()

const KspMethod::KSPContext ForcesAndSourcesCore::UserDataOperator::getKSPCtx ( ) const
inline

Definition at line 1188 of file ForcesAndSourcesCore.hpp.

1188 {
1189 return getFEMethod()->ksp_ctx;
1190}
KSPContext ksp_ctx
Current KSP computation context.

◆ getKSPf()

Vec ForcesAndSourcesCore::UserDataOperator::getKSPf ( ) const
inline

Definition at line 1202 of file ForcesAndSourcesCore.hpp.

1202 {
1203#ifndef NDEBUG
1204 if (getFEMethod()->ksp_f == PETSC_NULLPTR)
1205 THROW_MESSAGE("KSP not set F vector");
1206#endif
1207 return getFEMethod()->ksp_f;
1208}
Vec & ksp_f
Reference to residual vector in KSP context.

◆ getLoopSize()

int ForcesAndSourcesCore::UserDataOperator::getLoopSize ( ) const
inline

get size of elements in the loop

Returns
loop size

Definition at line 1174 of file ForcesAndSourcesCore.hpp.

1174 {
1175 return getFEMethod()->getLoopSize();
1176}
int getLoopSize() const
Get total loop size.

◆ getMeasure() [1/2]

double & ForcesAndSourcesCore::UserDataOperator::getMeasure ( )
inline

get measure of element

Returns
volume

Definition at line 1387 of file ForcesAndSourcesCore.hpp.

1387 {
1388 return static_cast<ForcesAndSourcesCore *>(ptrFE)->elementMeasure;
1389}

◆ getMeasure() [2/2]

double ForcesAndSourcesCore::UserDataOperator::getMeasure ( ) const
inline

◆ getMField()

MoFEM::Interface & MoFEM::ForcesAndSourcesCore::UserDataOperator::getMField ( )
inline

Definition at line 1083 of file ForcesAndSourcesCore.hpp.

1083{ return ptrFE->mField; }

◆ getMoab()

moab::Interface & MoFEM::ForcesAndSourcesCore::UserDataOperator::getMoab ( )
inline

Definition at line 1084 of file ForcesAndSourcesCore.hpp.

1084{ return getMField().get_moab(); }
virtual moab::Interface & get_moab()=0

◆ getNinTheLoop()

int ForcesAndSourcesCore::UserDataOperator::getNinTheLoop ( ) const
inline

get number of finite element in the loop

Returns
number of finite element

Definition at line 1170 of file ForcesAndSourcesCore.hpp.

1170 {
1171 return getFEMethod()->getNinTheLoop();
1172}
int getNinTheLoop() const
Get current loop iteration index.

◆ getNumberOfNodesOnElement()

int ForcesAndSourcesCore::UserDataOperator::getNumberOfNodesOnElement ( ) const
inline

Get the number of nodes on finite element.

Returns
int
Examples
mofem/atom_tests/boundary_marker.cpp.

Definition at line 1152 of file ForcesAndSourcesCore.hpp.

1152 {
1153 return ptrFE->getNumberOfNodes();
1154}
auto getNumberOfNodes() const
Get the number of nodes in the current finite element.

◆ getNumeredEntFiniteElementPtr()

boost::shared_ptr< const NumeredEntFiniteElement > ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr ( ) const
inline

◆ getOpType()

int ForcesAndSourcesCore::UserDataOperator::getOpType ( ) const
inline

Get operator types.

Returns
Return operator type

Definition at line 1160 of file ForcesAndSourcesCore.hpp.

1160{ return opType; }

◆ getProblemColIndices()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::getProblemColIndices ( const std::string  filed_name,
const EntityType  type,
const int  side,
VectorInt indices 
) const

Get col indices.

Field could be or not declared for this element but is declared for problem

Parameters
field_name
typeentity type
sideside number, any number if type is MBVERTEX
Returns
indices

NOTE: Using those indices to assemble matrix will result in error if new non-zero values need to be created.

Definition at line 1701 of file ForcesAndSourcesCore.cpp.

1703 {
1705 if (ptrFE == NULL)
1706 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1707
1708 switch (type) {
1709 case MBVERTEX:
1711 break;
1712 default:
1713 CHKERR ptrFE->getProblemTypeColIndices(field_name, type, side, indices);
1714 }
1716}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
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

◆ getProblemRowIndices()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::getProblemRowIndices ( const std::string  filed_name,
const EntityType  type,
const int  side,
VectorInt indices 
) const

Get row indices.

Field could be or not declared for this element but is declared for problem

Parameters
field_name
typeentity type
sideside number, any number if type is MBVERTEX
Returns
indices

NOTE: Using those indices to assemble matrix will result in error if new non-zero values need to be created.

Definition at line 1684 of file ForcesAndSourcesCore.cpp.

1686 {
1688 if (ptrFE == NULL)
1689 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1690
1691 switch (type) {
1692 case MBVERTEX:
1694 break;
1695 default:
1696 CHKERR ptrFE->getProblemTypeRowIndices(field_name, type, side, indices);
1697 }
1699}
MoFEMErrorCode getProblemTypeRowIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
MoFEMErrorCode getProblemNodesRowIndices(const std::string &field_name, VectorInt &nodes_indices) const

◆ getPtrFE()

ForcesAndSourcesCore * ForcesAndSourcesCore::UserDataOperator::getPtrFE ( ) const
inline

Definition at line 1359 of file ForcesAndSourcesCore.hpp.

1359 {
1360 return ptrFE;
1361}

◆ getRefinePtrFE()

ForcesAndSourcesCore * ForcesAndSourcesCore::UserDataOperator::getRefinePtrFE ( ) const
inline

Definition at line 1369 of file ForcesAndSourcesCore.hpp.

1369 {
1370 return ptrFE->refinePtrFE;
1371}
ForcesAndSourcesCore * refinePtrFE
Element to integrate parent or child.

◆ getSideEntity()

EntityHandle ForcesAndSourcesCore::UserDataOperator::getSideEntity ( const int  side_number,
const EntityType  type 
)
inline

Get the side entity.

Note
For vertex is expectation. Side basses in argument of function doWork is zero. For other entity types side can be used as argument of this function.
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBVERTEX) {
for (int n = 0; n != number_of_nodes; ++n)
// Do somthing
} else {
EntityHandle ent = getSideEntity(side, type);
// Do somthing
}
}
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
Data on single entity (This is passed as argument to DataOperator::doWork)
EntityHandle getSideEntity(const int side_number, const EntityType type)
Get the side entity.
Parameters
side_number
type
Examples
mofem/atom_tests/boundary_marker.cpp.

Definition at line 1144 of file ForcesAndSourcesCore.hpp.

1145 {
1146 if (auto side_ptr = getSideNumberPtr(side_number, type).lock())
1147 return side_ptr->ent;
1148 else
1149 return 0;
1150}
boost::weak_ptr< SideNumber > getSideNumberPtr(const int side_number, const EntityType type)
Get the side number pointer.

◆ getSideNumberPtr()

boost::weak_ptr< SideNumber > ForcesAndSourcesCore::UserDataOperator::getSideNumberPtr ( const int  side_number,
const EntityType  type 
)
inline

Get the side number pointer.

Note
For vertex is expectation. Side basses in argument of function doWork is zero. For other entity types side can be used as argument of this function.
Parameters
side_number
type
Returns
boost::weak_ptr<SideNumber>

Definition at line 1131 of file ForcesAndSourcesCore.hpp.

1132 {
1133 auto &side_table_by_side_and_type =
1134 ptrFE->numeredEntFiniteElementPtr->getSideNumberTable().get<1>();
1135 auto side_it =
1136 side_table_by_side_and_type.find(boost::make_tuple(type, side_number));
1137 if (side_it != side_table_by_side_and_type.end())
1138 return *side_it;
1139 else
1140 return boost::weak_ptr<SideNumber>();
1141}

◆ getSidePtrFE()

ForcesAndSourcesCore * ForcesAndSourcesCore::UserDataOperator::getSidePtrFE ( ) const
inline

Definition at line 1364 of file ForcesAndSourcesCore.hpp.

1364 {
1365 return ptrFE->sidePtrFE;
1366}
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.

◆ getSNESA()

Mat ForcesAndSourcesCore::UserDataOperator::getSNESA ( ) const
inline

Definition at line 1242 of file ForcesAndSourcesCore.hpp.

1242 {
1243#ifndef NDEBUG
1244 if (getFEMethod()->snes_A == PETSC_NULLPTR)
1245 THROW_MESSAGE("SNES not set A vector");
1246#endif
1247 return getFEMethod()->snes_A;
1248}
Mat & snes_A
Reference to Jacobian matrix.

◆ getSNESB()

Mat ForcesAndSourcesCore::UserDataOperator::getSNESB ( ) const
inline

Definition at line 1250 of file ForcesAndSourcesCore.hpp.

1250 {
1251#ifndef NDEBUG
1252 if (getFEMethod()->snes_B == PETSC_NULLPTR)
1253 THROW_MESSAGE("SNES not set A matrix");
1254#endif
1255 return getFEMethod()->snes_B;
1256}
Mat & snes_B
Reference to preconditioner of Jacobian matrix.

◆ getSNESCtx()

const SnesMethod::SNESContext ForcesAndSourcesCore::UserDataOperator::getSNESCtx ( ) const
inline

Definition at line 1193 of file ForcesAndSourcesCore.hpp.

1193 {
1194 return getFEMethod()->snes_ctx;
1195}
SNESContext snes_ctx
Current SNES computation context.

◆ getSNESf()

Vec ForcesAndSourcesCore::UserDataOperator::getSNESf ( ) const
inline

Definition at line 1226 of file ForcesAndSourcesCore.hpp.

1226 {
1227#ifndef NDEBUG
1228 if (getFEMethod()->snes_f == PETSC_NULLPTR)
1229 THROW_MESSAGE("SNES not set F vector");
1230#endif
1231 return getFEMethod()->snes_f;
1232}
Vec & snes_f
Reference to residual vector.

◆ getSNESx()

Vec ForcesAndSourcesCore::UserDataOperator::getSNESx ( ) const
inline

Definition at line 1234 of file ForcesAndSourcesCore.hpp.

1234 {
1235#ifndef NDEBUG
1236 if (getFEMethod()->snes_x == PETSC_NULLPTR)
1237 THROW_MESSAGE("SNESnot set X vector");
1238#endif
1239 return getFEMethod()->snes_x;
1240}
Vec & snes_x
Reference to current solution state vector.

◆ getSubPipelinePtr()

virtual boost::weak_ptr< ForcesAndSourcesCore > MoFEM::ForcesAndSourcesCore::UserDataOperator::getSubPipelinePtr ( ) const
inlinevirtual

Reimplemented in MoFEM::OpLoopThis< E >, MoFEM::OpLoopSide< E >, and MoFEM::OpLoopRange< E >.

Definition at line 1092 of file ForcesAndSourcesCore.hpp.

1092 {
1093 return boost::weak_ptr<ForcesAndSourcesCore>();
1094 };

◆ getTSA()

Mat ForcesAndSourcesCore::UserDataOperator::getTSA ( ) const
inline

Definition at line 1290 of file ForcesAndSourcesCore.hpp.

1290 {
1291#ifndef NDEBUG
1292 if (getFEMethod()->ts_A == PETSC_NULLPTR)
1293 THROW_MESSAGE("TS not set A matrix");
1294#endif
1295 return getFEMethod()->ts_A;
1296}
Mat & ts_A
Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.

◆ getTSa()

double ForcesAndSourcesCore::UserDataOperator::getTSa ( ) const
inline

Definition at line 1330 of file ForcesAndSourcesCore.hpp.

1330 {
1331#ifndef NDEBUG
1333 .none() ||
1335 THROW_MESSAGE("TS not set B matrix");
1336#endif
1337 return getFEMethod()->ts_a;
1338}
static constexpr Switches CtxSetA
Jacobian matrix switch.
static constexpr Switches CtxSetX_T
First time derivative switch.
static constexpr Switches CtxSetB
Preconditioner matrix switch.
PetscReal ts_a
Shift parameter for U_t (see PETSc Time Solver documentation)

◆ getTSaa()

double ForcesAndSourcesCore::UserDataOperator::getTSaa ( ) const
inline

Definition at line 1340 of file ForcesAndSourcesCore.hpp.

1340 {
1341#ifndef NDEBUG
1343 .none() ||
1345 THROW_MESSAGE("TS not set B matrix");
1346#endif
1347 return getFEMethod()->ts_aa;
1348}
static constexpr Switches CtxSetX_TT
Second time derivative switch.
PetscReal ts_aa
Shift parameter for U_tt (second time derivative)

◆ getTSB()

Mat ForcesAndSourcesCore::UserDataOperator::getTSB ( ) const
inline

Definition at line 1298 of file ForcesAndSourcesCore.hpp.

1298 {
1299#ifndef NDEBUG
1300 if (getFEMethod()->ts_B == PETSC_NULLPTR)
1301 THROW_MESSAGE("TS not set B matrix");
1302#endif
1303 return getFEMethod()->ts_B;
1304}
Mat & ts_B
Reference to preconditioner matrix for ts_A.

◆ getTSCtx()

const TSMethod::TSContext ForcesAndSourcesCore::UserDataOperator::getTSCtx ( ) const
inline

Definition at line 1198 of file ForcesAndSourcesCore.hpp.

1198 {
1199 return getFEMethod()->ts_ctx;
1200}
TSContext ts_ctx
Current TS computation context.

◆ getTSf()

Vec ForcesAndSourcesCore::UserDataOperator::getTSf ( ) const
inline

Definition at line 1282 of file ForcesAndSourcesCore.hpp.

1282 {
1283#ifndef NDEBUG
1284 if (getFEMethod()->ts_F == PETSC_NULLPTR)
1285 THROW_MESSAGE("TS not set F vector");
1286#endif
1287 return getFEMethod()->ts_F;
1288}
Vec & ts_F
Reference to residual vector F(t,U,U_t,U_tt)

◆ getTSstep()

int ForcesAndSourcesCore::UserDataOperator::getTSstep ( ) const
inline

Definition at line 1306 of file ForcesAndSourcesCore.hpp.

1306 {
1307#ifndef NDEBUG
1308 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1309 THROW_MESSAGE("TS not set time");
1310#endif
1311 return getFEMethod()->ts_step;
1312}
PetscInt ts_step
Current time step number.

◆ getTStime()

double ForcesAndSourcesCore::UserDataOperator::getTStime ( ) const
inline

Definition at line 1314 of file ForcesAndSourcesCore.hpp.

1314 {
1315#ifndef NDEBUG
1316 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1317 THROW_MESSAGE("TS not set time");
1318#endif
1319 return getFEMethod()->ts_t;
1320}
PetscReal ts_t
Current time value.

◆ getTStimeStep()

double ForcesAndSourcesCore::UserDataOperator::getTStimeStep ( ) const
inline

Definition at line 1322 of file ForcesAndSourcesCore.hpp.

1322 {
1323#ifndef NDEBUG
1324 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1325 THROW_MESSAGE("TS not set time");
1326#endif
1327 return getFEMethod()->ts_dt;
1328}
PetscReal ts_dt
Current time step size.

◆ getTSu()

Vec ForcesAndSourcesCore::UserDataOperator::getTSu ( ) const
inline

Definition at line 1258 of file ForcesAndSourcesCore.hpp.

1258 {
1259#ifndef NDEBUG
1260 if (getFEMethod()->ts_u == PETSC_NULLPTR)
1261 THROW_MESSAGE("TS not set U vector");
1262#endif
1263 return getFEMethod()->ts_u;
1264}
Vec & ts_u
Reference to current state vector U(t)

◆ getTSu_t()

Vec ForcesAndSourcesCore::UserDataOperator::getTSu_t ( ) const
inline

Definition at line 1266 of file ForcesAndSourcesCore.hpp.

1266 {
1267#ifndef NDEBUG
1268 if (getFEMethod()->ts_u_t == PETSC_NULLPTR)
1269 THROW_MESSAGE("TS not set U_t vector");
1270#endif
1271 return getFEMethod()->ts_u_t;
1272}
Vec & ts_u_t
Reference to first time derivative of state vector dU/dt.

◆ getTSu_tt()

Vec ForcesAndSourcesCore::UserDataOperator::getTSu_tt ( ) const
inline

Definition at line 1274 of file ForcesAndSourcesCore.hpp.

1274 {
1275#ifndef NDEBUG
1276 if (getFEMethod()->ts_u_tt == PETSC_NULLPTR)
1277 THROW_MESSAGE("TS not set U_tt vector");
1278#endif
1279 return getFEMethod()->ts_u_tt;
1280}
Vec & ts_u_tt
Reference to second time derivative of state vector d²U/dt²

◆ loopChildren()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopChildren ( const string &  fe_name,
ForcesAndSourcesCore child_fe,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.

Parameters
fe_name
child_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1986 of file ForcesAndSourcesCore.cpp.

1988 {
1990
1991 auto fe_miit = ptrFE->mField.get_finite_elements()
1992 ->get<FiniteElement_name_mi_tag>()
1993 .find(fe_name);
1994 if (fe_miit != ptrFE->mField.get_finite_elements()
1995 ->get<FiniteElement_name_mi_tag>()
1996 .end()) {
1997
1998 const auto *problem_ptr = getFEMethod()->problemPtr;
1999 auto &ref_ents = *getPtrFE()->mField.get_ref_ents();
2000 auto &numered_fe =
2001 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
2002
2003 const auto parent_ent = getNumeredEntFiniteElementPtr()->getEnt();
2004 const auto parent_type = getNumeredEntFiniteElementPtr()->getEntType();
2005 auto range =
2006 ref_ents.get<Composite_ParentEnt_And_EntType_mi_tag>().equal_range(
2007 boost::make_tuple(parent_type, parent_ent));
2008
2009 if (auto size = std::distance(range.first, range.second)) {
2010
2011 std::vector<EntityHandle> childs_vec;
2012 childs_vec.reserve(size);
2013 for (; range.first != range.second; ++range.first)
2014 childs_vec.emplace_back((*range.first)->getEnt());
2015
2016 Range childs;
2017
2018 if ((childs_vec.back() - childs_vec.front() + 1) == size)
2019 childs = Range(childs_vec.front(), childs_vec.back());
2020 else
2021 childs.insert_list(childs_vec.begin(), childs_vec.end());
2022
2023 child_fe->feName = fe_name;
2024
2025 CHKERR child_fe->setRefineFEPtr(ptrFE);
2026 CHKERR child_fe->copyBasicMethod(*getFEMethod());
2027 CHKERR child_fe->copyPetscData(*getFEMethod());
2028 CHKERR child_fe->copyKsp(*getFEMethod());
2029 CHKERR child_fe->copySnes(*getFEMethod());
2030 CHKERR child_fe->copyTs(*getFEMethod());
2031
2032 child_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
2033
2034 CHKERR child_fe->preProcess();
2035
2036 int nn = 0;
2037 child_fe->loopSize = size;
2038
2039 for (auto p = childs.pair_begin(); p != childs.pair_end(); ++p) {
2040
2041 auto miit =
2042 numered_fe.lower_bound(EntFiniteElement::getLocalUniqueIdCalculate(
2043 p->first, (*fe_miit)->getFEUId()));
2044 auto hi_miit =
2045 numered_fe.upper_bound(EntFiniteElement::getLocalUniqueIdCalculate(
2046 p->second, (*fe_miit)->getFEUId()));
2047
2048 for (; miit != hi_miit; ++miit) {
2049
2050 if (verb >= VERBOSE)
2051 MOFEM_LOG("SELF", sev)
2052 << "Child finite element " << "(" << nn << "): " << **miit;
2053
2054 child_fe->nInTheLoop = nn++;
2055 child_fe->numeredEntFiniteElementPtr = *miit;
2056
2057 if (child_fe->exeTestHook) {
2058 if (child_fe->exeTestHook(child_fe)) {
2059 CHKERR (*child_fe)();
2060 }
2061 } else {
2062 CHKERR (*child_fe)();
2063 }
2064
2065 }
2066 }
2067 }
2068
2069 CHKERR child_fe->postProcess();
2070 }
2071
2073}
@ VERBOSE
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
virtual const RefEntity_multiIndex * get_ref_ents() const =0
Get the ref ents object.
#define MOFEM_LOG(channel, severity)
Log.
const Problem * problemPtr
Raw pointer to current MoFEM problem instance.
boost::weak_ptr< CacheTuple > cacheWeakPtr
Weak pointer to cached entity data.
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.

◆ loopParent()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopParent ( const string &  fe_name,
ForcesAndSourcesCore parent_fe,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.

Parameters
fe_name
parent_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1929 of file ForcesAndSourcesCore.cpp.

1931 {
1933
1934 auto &fes =
1935 ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
1936 auto fe_miit = fes.find(fe_name);
1937 if (fe_miit != fes.end()) {
1938
1939 const auto *problem_ptr = getFEMethod()->problemPtr;
1940 auto &numered_fe =
1941 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1942
1943 parent_fe->feName = fe_name;
1944 CHKERR parent_fe->setRefineFEPtr(ptrFE);
1945 CHKERR parent_fe->copyBasicMethod(*getFEMethod());
1946 CHKERR parent_fe->copyPetscData(*getFEMethod());
1947 CHKERR parent_fe->copyKsp(*getFEMethod());
1948 CHKERR parent_fe->copySnes(*getFEMethod());
1949 CHKERR parent_fe->copyTs(*getFEMethod());
1950
1951 parent_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1952
1953 const auto parent_ent = getNumeredEntFiniteElementPtr()->getParentEnt();
1954 auto miit = numered_fe.find(EntFiniteElement::getLocalUniqueIdCalculate(
1955 parent_ent, (*fe_miit)->getFEUId()));
1956 if (miit != numered_fe.end()) {
1957 if (verb >= VERBOSE)
1958 MOFEM_LOG("SELF", sev) << "Parent finite element: " << **miit;
1959 parent_fe->loopSize = 1;
1960 parent_fe->nInTheLoop = 0;
1961 parent_fe->numeredEntFiniteElementPtr = *miit;
1962 CHKERR parent_fe->preProcess();
1963
1964 if (parent_fe->exeTestHook) {
1965 if (parent_fe->exeTestHook(parent_fe)) {
1966 CHKERR (*parent_fe)();
1967 }
1968 } else {
1969 CHKERR (*parent_fe)();
1970 }
1971
1972 CHKERR parent_fe->postProcess();
1973 } else {
1974 if (verb >= VERBOSE)
1975 MOFEM_LOG("SELF", sev) << "Parent finite element: no parent";
1976 parent_fe->loopSize = 0;
1977 parent_fe->nInTheLoop = 0;
1978 CHKERR parent_fe->preProcess();
1979 CHKERR parent_fe->postProcess();
1980 }
1981 }
1982
1984}
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements

◆ loopRange()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopRange ( const string &  fe_name,
ForcesAndSourcesCore range_fe,
boost::shared_ptr< Range fe_range,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

Iterate over range of elements.

Parameters
fe_name
range_fe
fe_range
verb
sev
Returns
* MoFEMErrorCode

Definition at line 2075 of file ForcesAndSourcesCore.cpp.

2078 {
2080
2081 auto &fes =
2082 ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
2083 auto fe_miit = fes.find(fe_name);
2084 if (fe_miit != fes.end()) {
2085
2086 const auto *problem_ptr = getFEMethod()->problemPtr;
2087 auto &numered_fe =
2088 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
2089
2090 range_fe->feName = fe_name;
2091 CHKERR range_fe->setRefineFEPtr(ptrFE);
2092 CHKERR range_fe->copyBasicMethod(*getFEMethod());
2093 CHKERR range_fe->copyPetscData(*getFEMethod());
2094 CHKERR range_fe->copyKsp(*getFEMethod());
2095 CHKERR range_fe->copySnes(*getFEMethod());
2096 CHKERR range_fe->copyTs(*getFEMethod());
2097
2098 range_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
2099
2100 auto get_numered_fe_ptr = [&](auto &fe_uid, auto fe_range, auto execute) {
2102 if (fe_range) {
2103 int nn = 0;
2104 for (auto p = fe_range->pair_begin(); p != fe_range->pair_end(); ++p) {
2105 auto first = p->first;
2106 auto second = p->second;
2107 auto lo = numered_fe.lower_bound(
2109 auto hi = numered_fe.upper_bound(
2111 for (; lo != hi; ++lo) {
2112 CHKERR execute(lo, nn++);
2113 }
2114 }
2115 }
2117 };
2118
2119 auto execute = [&](auto lo, auto nn) {
2121 if (verb >= VERBOSE)
2122 MOFEM_LOG("SELF", sev) << "Range finite element: " << **lo;
2123 range_fe->nInTheLoop = nn;
2124 range_fe->numeredEntFiniteElementPtr = *lo;
2125 CHKERR range_fe->preProcess();
2126
2127 if (range_fe->exeTestHook) {
2128 if (range_fe->exeTestHook(range_fe)) {
2129 CHKERR (*range_fe)();
2130 }
2131 } else {
2132 CHKERR (*range_fe)();
2133 }
2134
2135 CHKERR range_fe->postProcess();
2137 };
2138
2139 CHKERR get_numered_fe_ptr((*fe_miit)->getFEUId(), fe_range, execute);
2140 }
2141
2143}

◆ loopSide()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopSide ( const string &  fe_name,
ForcesAndSourcesCore side_fe,
const size_t  dim,
const EntityHandle  ent_for_side = 0,
boost::shared_ptr< Range fe_range = nullptr,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy,
AdjCache adj_cache = nullptr 
)

User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.

Parameters
fe_namename of the side element
side_fepointer to the side element instance
dimdimension the of side element
ent_for_sideentity handle for which adjacent volume or face will be accessed
verb
sev
adj_cache
Returns
MoFEMErrorCode

Definition at line 1732 of file ForcesAndSourcesCore.cpp.

1735 {
1737
1738 const auto *problem_ptr = getFEMethod()->problemPtr;
1739 auto &numered_fe =
1740 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1741
1742 auto fe_miit = ptrFE->mField.get_finite_elements()
1743 ->get<FiniteElement_name_mi_tag>()
1744 .find(fe_name);
1745 if (fe_miit != ptrFE->mField.get_finite_elements()
1746 ->get<FiniteElement_name_mi_tag>()
1747 .end()) {
1748
1749 const auto ent = ent_for_side ? ent_for_side : getFEEntityHandle();
1750
1751 side_fe->feName = fe_name;
1752
1753 CHKERR side_fe->setSideFEPtr(ptrFE);
1754 CHKERR side_fe->copyBasicMethod(*getFEMethod());
1755 CHKERR side_fe->copyPetscData(*getFEMethod());
1756 CHKERR side_fe->copyKsp(*getFEMethod());
1757 CHKERR side_fe->copySnes(*getFEMethod());
1758 CHKERR side_fe->copyTs(*getFEMethod());
1759
1760 side_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1761
1762 CHKERR side_fe->preProcess();
1763
1764 std::vector<boost::weak_ptr<NumeredEntFiniteElement>> fe_vec;
1765 auto get_numered_fe_ptr = [&](auto &fe_uid, Range &&adjacent_ents)
1766 -> std::vector<boost::weak_ptr<NumeredEntFiniteElement>> & {
1767 fe_vec.reserve(adjacent_ents.size());
1768 for (auto fe_ent : adjacent_ents) {
1769 auto miit = numered_fe.find(
1771 if (miit != numered_fe.end()) {
1772 fe_vec.emplace_back(*miit);
1773 }
1774 }
1775 return fe_vec;
1776 };
1777
1778 auto get_bit_entity_adjacency = [&]() {
1779 Range adjacent_ents;
1781 ptrFE->mField.getInterface<BitRefManager>()->getAdjacenciesAny(
1782 ent, side_dim, adjacent_ents),
1783 "getAdjacenciesAny failed");
1784 return adjacent_ents;
1785 };
1786
1787 auto get_bit_meshset_entities = [&]() {
1788 if (type_from_handle(ent) != MBENTITYSET) {
1789 return intersect(get_bit_entity_adjacency(), *fe_range);
1790 } else {
1791 auto &bit = getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
1792 Range ents = *fe_range;
1794 ptrFE->mField.getInterface<BitRefManager>()
1795 ->filterEntitiesByRefLevel(bit, BitRefLevel().set(), ents),
1796 "filterEntitiesByRefLevel failed");
1797 return ents;
1798 }
1799 };
1800
1801 auto get_adj = [&](auto &fe_uid, auto get_adj_fun)
1802 -> std::vector<boost::weak_ptr<NumeredEntFiniteElement>> & {
1803 if (adj_cache) {
1804 try {
1805 return (*adj_cache).at(ent);
1806 } catch (const std::out_of_range &) {
1807 return (*adj_cache)[ent] = get_numered_fe_ptr(fe_uid, get_adj_fun());
1808 }
1809 } else {
1810 return get_numered_fe_ptr(fe_uid, get_adj_fun());
1811 }
1812 };
1813
1814 if (type_from_handle(ent) == MBENTITYSET) {
1815 if (!fe_range)
1816 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1817 "No range of finite elements");
1818 }
1819 auto adj = (!fe_range)
1820 ? get_adj((*fe_miit)->getFEUId(), get_bit_entity_adjacency)
1821 : get_adj((*fe_miit)->getFEUId(), get_bit_meshset_entities);
1822
1823 if (verb >= VERBOSE && !adj.empty())
1824 MOFEM_LOG("SELF", sev) << "Number of side finite elements " << adj.size();
1825
1826 int nn = 0;
1827 side_fe->loopSize = adj.size();
1828 for (auto fe_weak_ptr : adj) {
1829 if (auto fe_ptr = fe_weak_ptr.lock()) {
1830 if (verb >= VERBOSE)
1831 MOFEM_LOG("SELF", sev)
1832 << "Side finite element " << "(" << nn << "): " << *fe_ptr;
1833 side_fe->nInTheLoop = nn;
1834 side_fe->numeredEntFiniteElementPtr = fe_ptr;
1835 if (side_fe->exeTestHook) {
1836 if (side_fe->exeTestHook(side_fe)) {
1837 CHKERR (*side_fe)();
1838 }
1839 } else {
1840 CHKERR (*side_fe)();
1841 }
1842 }
1843 ++nn;
1844 }
1845
1846 CHKERR side_fe->postProcess();
1847 }
1848
1850}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
auto bit
set bit
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ loopThis()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopThis ( const string &  fe_name,
ForcesAndSourcesCore this_fe,
const int  verb = QUIET,
const LogManager::SeverityLevel  sev = Sev::noisy 
)

User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.

Parameters
fe_name
this_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1852 of file ForcesAndSourcesCore.cpp.

1854 {
1856
1857
1858
1859 auto &fes =
1860 ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
1861 auto fe_miit = fes.find(fe_name);
1862 if (fe_miit != fes.end()) {
1863
1864 this_fe->feName = fe_name;
1865
1866 CHKERR this_fe->setRefineFEPtr(ptrFE);
1867 CHKERR this_fe->copyBasicMethod(*getFEMethod());
1868 CHKERR this_fe->copyPetscData(*getFEMethod());
1869 CHKERR this_fe->copyKsp(*getFEMethod());
1870 CHKERR this_fe->copySnes(*getFEMethod());
1871 CHKERR this_fe->copyTs(*getFEMethod());
1872
1873 this_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1874
1875 this_fe->nInTheLoop = getNinTheLoop();
1876 this_fe->loopSize = getLoopSize();
1877
1878 CHKERR this_fe->preProcess();
1879
1880 if (fe_name == getFEName()) {
1881
1882 if (verb >= VERBOSE)
1883 MOFEM_LOG("SELF", sev)
1884 << "This finite element: " << *getNumeredEntFiniteElementPtr();
1885
1886 this_fe->numeredEntFiniteElementPtr = getNumeredEntFiniteElementPtr();
1887 CHKERR (*this_fe)();
1888 } else {
1889
1890 auto get_numered_fe_ptr = [&](auto &fe_uid, auto fe_ent) {
1891 auto &numered_fe =
1892 getFEMethod()
1893 ->problemPtr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1894 auto it = numered_fe.find(
1896 boost::shared_ptr<const NumeredEntFiniteElement> this_numered_fe_ptr;
1897 if (it != numered_fe.end()) {
1898 this_numered_fe_ptr = *it;
1899 }
1900 return this_numered_fe_ptr;
1901 };
1902
1903 auto this_numered_fe_ptr = get_numered_fe_ptr(
1904 (*fe_miit)->getFEUId(), getNumeredEntFiniteElementPtr()->getEnt());
1905 if (this_numered_fe_ptr) {
1906
1907 if (verb >= VERBOSE)
1908 MOFEM_LOG("SELF", sev)
1909 << "This finite element: " << *this_numered_fe_ptr;
1910
1911 this_fe->numeredEntFiniteElementPtr = this_numered_fe_ptr;
1912
1913 if (this_fe->exeTestHook) {
1914 if (this_fe->exeTestHook(this_fe)) {
1915 CHKERR (*this_fe)();
1916 }
1917 } else {
1918 CHKERR (*this_fe)();
1919 }
1920
1921 }
1922 }
1923 CHKERR this_fe->postProcess();
1924 }
1925
1927}
int getLoopSize() const
get size of elements in the loop
int getNinTheLoop() const
get number of finite element in the loop
std::string getFEName() const
Get name of the element.

◆ setOpType()

void ForcesAndSourcesCore::UserDataOperator::setOpType ( const OpType  type)
inline

Set operator type.

Parameters
Operatortype

Definition at line 1162 of file ForcesAndSourcesCore.hpp.

1162 {
1163 opType = type;
1164}

◆ setPtrFE()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Reimplemented in MoFEM::EdgeElementForcesAndSourcesCore::UserDataOperator, MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator, MoFEM::FaceElementForcesAndSourcesCoreOnSide::UserDataOperator, MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator, MoFEM::FlatPrismElementForcesAndSourcesCore::UserDataOperator, MoFEM::VertexElementForcesAndSourcesCore::UserDataOperator, MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator, and MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator.

Definition at line 2216 of file ForcesAndSourcesCore.cpp.

2216 {
2218 if (!(ptrFE = dynamic_cast<ForcesAndSourcesCore *>(ptr)))
2219 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2220 "User operator and finite element do not work together");
2222}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

Friends And Related Symbol Documentation

◆ ContactPrismElementForcesAndSourcesCore

Definition at line 1107 of file ForcesAndSourcesCore.hpp.

◆ EdgeElementForcesAndSourcesCore

friend class EdgeElementForcesAndSourcesCore
friend

Definition at line 1105 of file ForcesAndSourcesCore.hpp.

◆ FaceElementForcesAndSourcesCore

friend class FaceElementForcesAndSourcesCore
friend

Definition at line 1106 of file ForcesAndSourcesCore.hpp.

◆ ForcesAndSourcesCore

friend class ForcesAndSourcesCore
friend

Definition at line 1104 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ colFieldName

std::string MoFEM::ForcesAndSourcesCore::UserDataOperator::colFieldName

Definition at line 613 of file ForcesAndSourcesCore.hpp.

◆ opType

char MoFEM::ForcesAndSourcesCore::UserDataOperator::opType

Definition at line 611 of file ForcesAndSourcesCore.hpp.

◆ OpTypeNames

const char *const ForcesAndSourcesCore::UserDataOperator::OpTypeNames
static
Initial value:
= {
"OPROW", " OPCOL", "OPROWCOL", "OPSPACE", "OPLAST"}

Definition at line 609 of file ForcesAndSourcesCore.hpp.

◆ ptrFE

ForcesAndSourcesCore* MoFEM::ForcesAndSourcesCore::UserDataOperator::ptrFE
protected

Definition at line 1099 of file ForcesAndSourcesCore.hpp.

◆ rowFieldName

std::string MoFEM::ForcesAndSourcesCore::UserDataOperator::rowFieldName

Definition at line 612 of file ForcesAndSourcesCore.hpp.

◆ sPace

FieldSpace MoFEM::ForcesAndSourcesCore::UserDataOperator::sPace

Definition at line 614 of file ForcesAndSourcesCore.hpp.


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