v0.13.2
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)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
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 funtions and integration points
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
Coordinates and access to internal data
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
Measures (area, volume, length, etc.)
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
Loops
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)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
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 operator to do calculations. More...
 
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 operator to do calculations. More...
 
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 operator to do calculations. More...
 
- 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. More...
 
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. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

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. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

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
 

Detailed Description

Examples
hanging_node_approx.cpp, hello_world.cpp, level_set.cpp, and mesh_smoothing.cpp.

Definition at line 546 of file ForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ OpType

enum MoFEM::ForcesAndSourcesCore::UserDataOperator::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 563 of file ForcesAndSourcesCore.hpp.

563 {
564 OPROW = 1 << 0, ///< operator doWork function is executed on FE rows
565 OPCOL = 1 << 1, ///< operator doWork function is executed on FE columns
566 OPROWCOL = 1 << 2, ///< operator doWork is executed on FE rows &columns
567 OPSPACE = 1 << 3, ///< operator do Work is execute on space data
568 OPLAST = 1 << 3 ///< @deprecated would be removed
569 };
@ 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 
)

This Constructor is used typically when some modification base shape functions on some approx. space is applied. Operator is run for all data on space.

User has no access to field data from this operator.

Definition at line 1888 of file ForcesAndSourcesCore.cpp.

◆ UserDataOperator() [2/3]

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

Definition at line 1893 of file ForcesAndSourcesCore.cpp.

1895 : 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 
)

Definition at line 1898 of file ForcesAndSourcesCore.cpp.

1901 : DataOperator(symm), opType(type), rowFieldName(row_field_name),
1902 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 1063 of file ForcesAndSourcesCore.hpp.

1063 {
1064 opType |= type;
1065}

◆ getCoordsAtGaussPts()

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

◆ getDataCtx()

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

Definition at line 1080 of file ForcesAndSourcesCore.hpp.

1080 {
1081 return getFEMethod()->data_ctx;
1082}
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ getFEDim()

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

Get dimension of finite element.

Returns
int

Definition at line 1019 of file ForcesAndSourcesCore.hpp.

1019 {
1021};
auto dimension_from_handle(const EntityHandle h)
get entity dimension form handle
Definition: Templates.hpp:1602
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
HookeInternalStressElement.hpp, UnsaturatedFlow.hpp, and mixed_poisson.cpp.

Definition at line 1015 of file ForcesAndSourcesCore.hpp.

1015 {
1016 return getNumeredEntFiniteElementPtr()->getEnt();
1017}
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 1075 of file ForcesAndSourcesCore.hpp.

1075 {
1076 return getFEMethod()->getFEName();
1077}
auto getFEName() const
get finite element name

◆ getFEType()

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

Get dimension of finite element.

Returns
int

Definition at line 1023 of file ForcesAndSourcesCore.hpp.

1023 {
1025};
auto type_from_handle(const EntityHandle h)
get type from entity handle
Definition: Templates.hpp:1594

◆ 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
PoissonDiscontinousGalerkin.hpp, UnsaturatedFlow.hpp, child_and_parent.cpp, hanging_node_approx.cpp, higher_derivatives.cpp, mixed_poisson.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, and simple_interface.cpp.

Definition at line 1280 of file ForcesAndSourcesCore.hpp.

1280 {
1282 &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
1283 &getCoordsAtGaussPts()(0, 2));
1284}
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
MagneticElement.hpp, PlasticOpsGeneric.hpp, PoissonDiscontinousGalerkin.hpp, ThermoElasticOps.hpp, UnsaturatedFlow.hpp, approx_sphere.cpp, child_and_parent.cpp, continuity_check_on_skeleton_3d.cpp, hanging_node_approx.cpp, higher_derivatives.cpp, level_set.cpp, mixed_poisson.cpp, phase.cpp, photon_diffusion.cpp, plate.cpp, poisson_2d_homogeneous.hpp, poisson_3d_homogeneous.hpp, reaction_diffusion.cpp, and simple_interface.cpp.

Definition at line 1247 of file ForcesAndSourcesCore.hpp.

1247 {
1248 return static_cast<ForcesAndSourcesCore *>(ptrFE)->gaussPts;
1249}
MatrixDouble gaussPts
Matrix of integration points.

◆ getKSPA()

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

Definition at line 1107 of file ForcesAndSourcesCore.hpp.

1107 {
1108#ifndef NDEBUG
1109 if (getFEMethod()->ksp_A == PETSC_NULL)
1110 THROW_MESSAGE("KSP not set A vector");
1111#endif
1112 return getFEMethod()->ksp_A;
1113}
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:561

◆ getKSPB()

Mat ForcesAndSourcesCore::UserDataOperator::getKSPB ( ) const
inline
Examples
PoissonDiscontinousGalerkin.hpp, plate.cpp, poisson_2d_homogeneous.hpp, and poisson_3d_homogeneous.hpp.

Definition at line 1115 of file ForcesAndSourcesCore.hpp.

1115 {
1116#ifndef NDEBUG
1117 if (getFEMethod()->ksp_B == PETSC_NULL)
1118 THROW_MESSAGE("KSP not set B vector");
1119#endif
1120 return getFEMethod()->ksp_B;
1121}

◆ getKSPCtx()

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

Definition at line 1085 of file ForcesAndSourcesCore.hpp.

1085 {
1086 return getFEMethod()->ksp_ctx;
1087}
KSPContext ksp_ctx
Context.
Definition: LoopMethods.hpp:86

◆ getKSPf()

Vec ForcesAndSourcesCore::UserDataOperator::getKSPf ( ) const
inline
Examples
PoissonDiscontinousGalerkin.hpp, poisson_2d_homogeneous.hpp, and poisson_3d_homogeneous.hpp.

Definition at line 1099 of file ForcesAndSourcesCore.hpp.

1099 {
1100#ifndef NDEBUG
1101 if (getFEMethod()->ksp_f == PETSC_NULL)
1102 THROW_MESSAGE("KSP not set F vector");
1103#endif
1104 return getFEMethod()->ksp_f;
1105}

◆ getLoopSize()

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

get size of elements in the loop

Returns
loop size

Definition at line 1071 of file ForcesAndSourcesCore.hpp.

1071 {
1072 return getFEMethod()->getLoopSize();
1073}
int getLoopSize() const
get loop size

◆ getMeasure() [1/2]

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

get measure of element

Returns
volume

Definition at line 1290 of file ForcesAndSourcesCore.hpp.

1290 {
1291 return static_cast<ForcesAndSourcesCore *>(ptrFE)->elementMeasure;
1292}

◆ getMeasure() [2/2]

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

get measure of element

Returns
volume
Examples
child_and_parent.cpp, hanging_node_approx.cpp, higher_derivatives.cpp, mixed_poisson.cpp, and poisson_3d_homogeneous.hpp.

Definition at line 1286 of file ForcesAndSourcesCore.hpp.

1286 {
1287 return static_cast<ForcesAndSourcesCore *>(ptrFE)->elementMeasure;
1288}

◆ getNinTheLoop()

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

get number of finite element in the loop

Returns
number of finite element

Definition at line 1067 of file ForcesAndSourcesCore.hpp.

1067 {
1068 return getFEMethod()->getNinTheLoop();
1069}
int getNinTheLoop() const
get number of evaluated element in the loop

◆ getNumberOfNodesOnElement()

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

Get the number of nodes on finite element.

Returns
int
Examples
boundary_marker.cpp.

Definition at line 1049 of file ForcesAndSourcesCore.hpp.

1049 {
1050 return ptrFE->getNumberOfNodes();
1051}
auto getNumberOfNodes() const

◆ getNumeredEntFiniteElementPtr()

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

Return raw pointer to NumeredEntFiniteElement.

Examples
ElasticityMixedFormulation.hpp, NavierStokesElement.cpp, and UnsaturatedFlow.hpp.

Definition at line 1011 of file ForcesAndSourcesCore.hpp.

1011 {
1013};
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getOpType()

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

Get operator types.

Returns
Return operator type

Definition at line 1057 of file ForcesAndSourcesCore.hpp.

1057{ 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 1630 of file ForcesAndSourcesCore.cpp.

1632 {
1634 if (ptrFE == NULL)
1635 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1636
1637 switch (type) {
1638 case MBVERTEX:
1640 break;
1641 default:
1642 CHKERR ptrFE->getProblemTypeColIndices(field_name, type, side, indices);
1643 }
1645}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEMErrorCode getProblemTypeColIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
MoFEMErrorCode getProblemNodesColIndices(const std::string &field_name, VectorInt &nodes_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 1613 of file ForcesAndSourcesCore.cpp.

1615 {
1617 if (ptrFE == NULL)
1618 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1619
1620 switch (type) {
1621 case MBVERTEX:
1623 break;
1624 default:
1625 CHKERR ptrFE->getProblemTypeRowIndices(field_name, type, side, indices);
1626 }
1628}
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

◆ getPtrFE()

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

Definition at line 1262 of file ForcesAndSourcesCore.hpp.

1262 {
1263 return ptrFE;
1264}

◆ getRefinePtrFE()

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

Definition at line 1272 of file ForcesAndSourcesCore.hpp.

1272 {
1273 return ptrFE->refinePtrFE;
1274}
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 expection. Side basses in argument of function doWork is zero. For other entity types side can be used as argument of this function.
if (type == MBVERTEX) {
for (int n = 0; n != number_of_nodes; ++n)
// Do somthing
} else {
EntityHandle ent = getSideEntity(side, type);
// Do somthing
}
}
FTensor::Index< 'n', SPACE_DIM > n
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
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
boundary_marker.cpp.

Definition at line 1041 of file ForcesAndSourcesCore.hpp.

1042 {
1043 if (auto side_ptr = getSideNumberPtr(side_number, type).lock())
1044 return side_ptr->ent;
1045 else
1046 return 0;
1047}
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 expection. 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 1028 of file ForcesAndSourcesCore.hpp.

1029 {
1030 auto &side_table_by_side_and_type =
1031 ptrFE->numeredEntFiniteElementPtr->getSideNumberTable().get<1>();
1032 auto side_it =
1033 side_table_by_side_and_type.find(boost::make_tuple(type, side_number));
1034 if (side_it != side_table_by_side_and_type.end())
1035 return *side_it;
1036 else
1037 return boost::weak_ptr<SideNumber>();
1038}

◆ getSidePtrFE()

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

Definition at line 1267 of file ForcesAndSourcesCore.hpp.

1267 {
1268 return ptrFE->sidePtrFE;
1269}
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.

◆ getSNESA()

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

Definition at line 1139 of file ForcesAndSourcesCore.hpp.

1139 {
1140#ifndef NDEBUG
1141 if (getFEMethod()->snes_A == PETSC_NULL)
1142 THROW_MESSAGE("SNES not set A vector");
1143#endif
1144 return getFEMethod()->snes_A;
1145}
Mat & snes_A
jacobian matrix

◆ getSNESB()

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

Definition at line 1147 of file ForcesAndSourcesCore.hpp.

1147 {
1148#ifndef NDEBUG
1149 if (getFEMethod()->snes_B == PETSC_NULL)
1150 THROW_MESSAGE("SNES not set A matrix");
1151#endif
1152 return getFEMethod()->snes_B;
1153}
Mat & snes_B
preconditioner of jacobian matrix

◆ getSNESCtx()

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

Definition at line 1090 of file ForcesAndSourcesCore.hpp.

1090 {
1091 return getFEMethod()->snes_ctx;
1092}
SNESContext snes_ctx

◆ getSNESf()

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

Definition at line 1123 of file ForcesAndSourcesCore.hpp.

1123 {
1124#ifndef NDEBUG
1125 if (getFEMethod()->snes_f == PETSC_NULL)
1126 THROW_MESSAGE("SNES not set F vector");
1127#endif
1128 return getFEMethod()->snes_f;
1129}
Vec & snes_f
residual

◆ getSNESx()

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

Definition at line 1131 of file ForcesAndSourcesCore.hpp.

1131 {
1132#ifndef NDEBUG
1133 if (getFEMethod()->snes_x == PETSC_NULL)
1134 THROW_MESSAGE("SNESnot set X vector");
1135#endif
1136 return getFEMethod()->snes_x;
1137}
Vec & snes_x
state vector

◆ getTSA()

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

Definition at line 1187 of file ForcesAndSourcesCore.hpp.

1187 {
1188#ifndef NDEBUG
1189 if (getFEMethod()->ts_A == PETSC_NULL)
1190 THROW_MESSAGE("TS not set A matrix");
1191#endif
1192 return getFEMethod()->ts_A;
1193}

◆ getTSa()

double ForcesAndSourcesCore::UserDataOperator::getTSa ( ) const
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 1227 of file ForcesAndSourcesCore.hpp.

1227 {
1228#ifndef NDEBUG
1230 .none() ||
1232 THROW_MESSAGE("TS not set B matrix");
1233#endif
1234 return getFEMethod()->ts_a;
1235}
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:40
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38
PetscReal ts_a
shift for U_t (see PETSc Time Solver)

◆ getTSaa()

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

Definition at line 1237 of file ForcesAndSourcesCore.hpp.

1237 {
1238#ifndef NDEBUG
1240 .none() ||
1242 THROW_MESSAGE("TS not set B matrix");
1243#endif
1244 return getFEMethod()->ts_aa;
1245}
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:41
PetscReal ts_aa
shift for U_tt shift for U_tt

◆ getTSB()

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

Definition at line 1195 of file ForcesAndSourcesCore.hpp.

1195 {
1196#ifndef NDEBUG
1197 if (getFEMethod()->ts_B == PETSC_NULL)
1198 THROW_MESSAGE("TS not set B matrix");
1199#endif
1200 return getFEMethod()->ts_B;
1201}
Mat & ts_B
Preconditioner for ts_A.

◆ getTSCtx()

const TSMethod::TSContext ForcesAndSourcesCore::UserDataOperator::getTSCtx ( ) const
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 1095 of file ForcesAndSourcesCore.hpp.

1095 {
1096 return getFEMethod()->ts_ctx;
1097}
TSContext ts_ctx

◆ getTSf()

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

Definition at line 1179 of file ForcesAndSourcesCore.hpp.

1179 {
1180#ifndef NDEBUG
1181 if (getFEMethod()->ts_F == PETSC_NULL)
1182 THROW_MESSAGE("TS not set F vector");
1183#endif
1184 return getFEMethod()->ts_F;
1185}
Vec & ts_F
residual vector

◆ getTSstep()

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

Definition at line 1203 of file ForcesAndSourcesCore.hpp.

1203 {
1204#ifndef NDEBUG
1205 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1206 THROW_MESSAGE("TS not set time");
1207#endif
1208 return getFEMethod()->ts_step;
1209}
PetscInt ts_step
time step number

◆ getTStime()

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

Definition at line 1211 of file ForcesAndSourcesCore.hpp.

1211 {
1212#ifndef NDEBUG
1213 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1214 THROW_MESSAGE("TS not set time");
1215#endif
1216 return getFEMethod()->ts_t;
1217}
PetscReal ts_t
time

◆ getTStimeStep()

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

Definition at line 1219 of file ForcesAndSourcesCore.hpp.

1219 {
1220#ifndef NDEBUG
1221 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1222 THROW_MESSAGE("TS not set time");
1223#endif
1224 return getFEMethod()->ts_dt;
1225}
PetscReal ts_dt
time step size

◆ getTSu()

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

Definition at line 1155 of file ForcesAndSourcesCore.hpp.

1155 {
1156#ifndef NDEBUG
1157 if (getFEMethod()->ts_u == PETSC_NULL)
1158 THROW_MESSAGE("TS not set U vector");
1159#endif
1160 return getFEMethod()->ts_u;
1161}
Vec & ts_u
state vector

◆ getTSu_t()

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

Definition at line 1163 of file ForcesAndSourcesCore.hpp.

1163 {
1164#ifndef NDEBUG
1165 if (getFEMethod()->ts_u_t == PETSC_NULL)
1166 THROW_MESSAGE("TS not set U_t vector");
1167#endif
1168 return getFEMethod()->ts_u_t;
1169}
Vec & ts_u_t
time derivative of state vector

◆ getTSu_tt()

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

Definition at line 1171 of file ForcesAndSourcesCore.hpp.

1171 {
1172#ifndef NDEBUG
1173 if (getFEMethod()->ts_u_tt == PETSC_NULL)
1174 THROW_MESSAGE("TS not set U_tt vector");
1175#endif
1176 return getFEMethod()->ts_u_tt;
1177}
Vec & ts_u_tt
second time derivative of state vector

◆ loopChildren()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::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 operator to do calculations.

Parameters
fe_name
child_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1787 of file ForcesAndSourcesCore.cpp.

1789 {
1791
1792 const auto *problem_ptr = getFEMethod()->problemPtr;
1793 auto &ref_ents = *getPtrFE()->mField.get_ref_ents();
1794 auto &numered_fe =
1795 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1796
1797 const auto parent_ent = getNumeredEntFiniteElementPtr()->getEnt();
1798 const auto parent_type = getNumeredEntFiniteElementPtr()->getEntType();
1799 auto range =
1800 ref_ents.get<Composite_ParentEnt_And_EntType_mi_tag>().equal_range(
1801 boost::make_tuple(parent_type, parent_ent));
1802
1803 if (auto size = std::distance(range.first, range.second)) {
1804
1805 std::vector<EntityHandle> childs_vec;
1806 childs_vec.reserve(size);
1807 for (; range.first != range.second; ++range.first)
1808 childs_vec.emplace_back((*range.first)->getEnt());
1809
1810 Range childs;
1811
1812 if ((childs_vec.back() - childs_vec.front() + 1) == size)
1813 childs = Range(childs_vec.front(), childs_vec.back());
1814 else
1815 childs.insert_list(childs_vec.begin(), childs_vec.end());
1816
1817 child_fe->feName = fe_name;
1818
1819 CHKERR child_fe->setRefineFEPtr(ptrFE);
1820 CHKERR child_fe->copyBasicMethod(*getFEMethod());
1821 CHKERR child_fe->copyPetscData(*getFEMethod());
1822 CHKERR child_fe->copyKsp(*getFEMethod());
1823 CHKERR child_fe->copySnes(*getFEMethod());
1824 CHKERR child_fe->copyTs(*getFEMethod());
1825
1826 CHKERR child_fe->preProcess();
1827
1828 int nn = 0;
1829 child_fe->loopSize = size;
1830
1831 auto fe_miit = ptrFE->mField.get_finite_elements()
1832 ->get<FiniteElement_name_mi_tag>()
1833 .find(fe_name);
1834 if (fe_miit != ptrFE->mField.get_finite_elements()
1835 ->get<FiniteElement_name_mi_tag>()
1836 .end()) {
1837
1838 for (auto p = childs.pair_begin(); p != childs.pair_end(); ++p) {
1839
1840 auto miit =
1841 numered_fe.lower_bound(EntFiniteElement::getLocalUniqueIdCalculate(
1842 p->first, (*fe_miit)->getFEUId()));
1843 auto hi_miit =
1844 numered_fe.upper_bound(EntFiniteElement::getLocalUniqueIdCalculate(
1845 p->second, (*fe_miit)->getFEUId()));
1846
1847 for (; miit != hi_miit; ++miit) {
1848
1849 if (verb >= VERBOSE)
1850 MOFEM_LOG("SELF", sev) << "Child finite element "
1851 << "(" << nn << "): " << **miit;
1852
1853 child_fe->nInTheLoop = nn++;
1854 child_fe->numeredEntFiniteElementPtr = *miit;
1855 CHKERR (*child_fe)();
1856 }
1857 }
1858 }
1859
1860 CHKERR child_fe->postProcess();
1861 };
1862
1864}
static Index< 'p', 3 > p
@ VERBOSE
Definition: definitions.h:209
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.
Definition: LogManager.hpp:301
const Problem * problemPtr
raw pointer to problem
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 call this function to loop over parent elements. This function calls finite element with is operator to do calculations.

Parameters
fe_name
parent_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1744 of file ForcesAndSourcesCore.cpp.

1746 {
1748
1749 const auto *problem_ptr = getFEMethod()->problemPtr;
1750 auto &numered_fe =
1751 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1752
1753 parent_fe->feName = fe_name;
1754 CHKERR parent_fe->setRefineFEPtr(ptrFE);
1755 CHKERR parent_fe->copyBasicMethod(*getFEMethod());
1756 CHKERR parent_fe->copyPetscData(*getFEMethod());
1757 CHKERR parent_fe->copyKsp(*getFEMethod());
1758 CHKERR parent_fe->copySnes(*getFEMethod());
1759 CHKERR parent_fe->copyTs(*getFEMethod());
1760
1761 CHKERR parent_fe->preProcess();
1762
1763 auto fe_miit = ptrFE->mField.get_finite_elements()
1764 ->get<FiniteElement_name_mi_tag>()
1765 .find(fe_name);
1766 if (fe_miit != ptrFE->mField.get_finite_elements()
1767 ->get<FiniteElement_name_mi_tag>()
1768 .end()) {
1769 const auto parent_ent = getNumeredEntFiniteElementPtr()->getParentEnt();
1770 auto miit = numered_fe.find(EntFiniteElement::getLocalUniqueIdCalculate(
1771 parent_ent, (*fe_miit)->getFEUId()));
1772 if (miit != numered_fe.end()) {
1773 if (verb >= VERBOSE)
1774 MOFEM_LOG("SELF", sev) << "Parent finite element: " << **miit;
1775 parent_fe->loopSize = 1;
1776 parent_fe->nInTheLoop = 0;
1777 parent_fe->numeredEntFiniteElementPtr = *miit;
1778 CHKERR (*parent_fe)();
1779 }
1780 }
1781
1782 CHKERR parent_fe->postProcess();
1783
1785}
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements

◆ loopSide()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::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 
)

User call this function to loop over elements on the side of face. This function calls finite element with is 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
Returns
MoFEMErrorCode

Definition at line 1661 of file ForcesAndSourcesCore.cpp.

1664 {
1666 const auto ent = ent_for_side ? ent_for_side : getFEEntityHandle();
1667 const auto *problem_ptr = getFEMethod()->problemPtr;
1668
1669 side_fe->feName = fe_name;
1670
1671 CHKERR side_fe->setSideFEPtr(ptrFE);
1672 CHKERR side_fe->copyBasicMethod(*getFEMethod());
1673 CHKERR side_fe->copyPetscData(*getFEMethod());
1674 CHKERR side_fe->copyKsp(*getFEMethod());
1675 CHKERR side_fe->copySnes(*getFEMethod());
1676 CHKERR side_fe->copyTs(*getFEMethod());
1677
1678 CHKERR side_fe->preProcess();
1679
1680 Range adjacent_ents;
1681 CHKERR ptrFE->mField.getInterface<BitRefManager>()->getAdjacenciesAny(
1682 ent, side_dim, adjacent_ents);
1683 auto &numered_fe =
1684 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1685
1686 auto fe_miit = ptrFE->mField.get_finite_elements()
1687 ->get<FiniteElement_name_mi_tag>()
1688 .find(fe_name);
1689 if (fe_miit != ptrFE->mField.get_finite_elements()
1690 ->get<FiniteElement_name_mi_tag>()
1691 .end()) {
1692 int nn = 0;
1693 side_fe->loopSize = adjacent_ents.size();
1694 for (auto fe_ent : adjacent_ents) {
1695 auto miit = numered_fe.find(EntFiniteElement::getLocalUniqueIdCalculate(
1696 fe_ent, (*fe_miit)->getFEUId()));
1697 if (miit != numered_fe.end()) {
1698 if (verb >= VERBOSE)
1699 MOFEM_LOG("SELF", sev) << "Side finite element "
1700 << "(" << nn << "): " << **miit;
1701 side_fe->nInTheLoop = nn++;
1702 side_fe->numeredEntFiniteElementPtr = *miit;
1703 CHKERR (*side_fe)();
1704 }
1705 }
1706 }
1707
1708 CHKERR side_fe->postProcess();
1709
1711}
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ loopThis()

MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::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 operator to do calculations.

Parameters
fe_name
parent_fe
verb
sev
Returns
MoFEMErrorCode

Definition at line 1713 of file ForcesAndSourcesCore.cpp.

1715 {
1717
1718 if (verb >= VERBOSE)
1719 MOFEM_LOG("SELF", sev) << "This finite element: "
1721
1722 const auto *problem_ptr = getFEMethod()->problemPtr;
1723 this_fe->feName = fe_name;
1724
1725 CHKERR this_fe->setRefineFEPtr(ptrFE);
1726 CHKERR this_fe->copyBasicMethod(*getFEMethod());
1727 CHKERR this_fe->copyPetscData(*getFEMethod());
1728 CHKERR this_fe->copyKsp(*getFEMethod());
1729 CHKERR this_fe->copySnes(*getFEMethod());
1730 CHKERR this_fe->copyTs(*getFEMethod());
1731
1732 CHKERR this_fe->preProcess();
1733
1734 this_fe->nInTheLoop = getNinTheLoop();
1735 this_fe->numeredEntFiniteElementPtr = getNumeredEntFiniteElementPtr();
1736
1737 CHKERR (*this_fe)();
1738
1739 CHKERR this_fe->postProcess();
1740
1742}
int getNinTheLoop() const
get number of finite element in the loop

◆ setOpType()

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

Set operator type.

Parameters
Operatortype

Definition at line 1059 of file ForcesAndSourcesCore.hpp.

1059 {
1060 opType = type;
1061}

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

1937 {
1939 if (!(ptrFE = dynamic_cast<ForcesAndSourcesCore *>(ptr)))
1940 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1941 "User operator and finite element do not work together");
1943}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440

Friends And Related Function Documentation

◆ ContactPrismElementForcesAndSourcesCore

Definition at line 1004 of file ForcesAndSourcesCore.hpp.

◆ EdgeElementForcesAndSourcesCore

friend class EdgeElementForcesAndSourcesCore
friend

Definition at line 1002 of file ForcesAndSourcesCore.hpp.

◆ FaceElementForcesAndSourcesCore

friend class FaceElementForcesAndSourcesCore
friend

Definition at line 1003 of file ForcesAndSourcesCore.hpp.

◆ ForcesAndSourcesCore

friend class ForcesAndSourcesCore
friend

Definition at line 1001 of file ForcesAndSourcesCore.hpp.

Member Data Documentation

◆ colFieldName

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

Definition at line 575 of file ForcesAndSourcesCore.hpp.

◆ opType

char MoFEM::ForcesAndSourcesCore::UserDataOperator::opType

Definition at line 573 of file ForcesAndSourcesCore.hpp.

◆ OpTypeNames

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

Definition at line 571 of file ForcesAndSourcesCore.hpp.

◆ ptrFE

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

Definition at line 996 of file ForcesAndSourcesCore.hpp.

◆ rowFieldName

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

Definition at line 574 of file ForcesAndSourcesCore.hpp.

◆ sPace

FieldSpace MoFEM::ForcesAndSourcesCore::UserDataOperator::sPace

Definition at line 576 of file ForcesAndSourcesCore.hpp.


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