v0.15.0
Loading...
Searching...
No Matches
MoFEM::DataOperator Struct Reference

base operator to do operations at Gauss Pt. level More...

#include "src/finite_elements/DataOperators.hpp"

Inheritance diagram for MoFEM::DataOperator:
[legend]
Collaboration diagram for MoFEM::DataOperator:
[legend]

Public Types

using DoWorkLhsHookFunType
 
using DoWorkRhsHookFunType
 

Public Member Functions

 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

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
 

Private Member Functions

template<bool Symm>
MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
template<bool ErrorIfNoBase>
MoFEMErrorCode opRhs (EntitiesFieldData &data, const std::array< bool, MBMAXTYPE > &do_entities)
 

Detailed Description

Member Typedef Documentation

◆ DoWorkLhsHookFunType

Initial value:
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)>
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
DataOperator(const bool symm=true)

Definition at line 30 of file DataOperators.hpp.

◆ DoWorkRhsHookFunType

Initial value:
boost::function<MoFEMErrorCode(
DataOperator *op_ptr, int side, EntityType type,
EntitiesFieldData::EntData &data)>

Definition at line 58 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ DataOperator()

MoFEM::DataOperator::DataOperator ( const bool symm = true)

This not yet implemented, switch off.

Definition at line 21 of file DataOperators.cpp.

22 :
23
24 sYmm(symm),
25
26 doEntities{true, true, true, true, true, true,
27 true, true, true, true, true, true},
28
29 doVertices(doEntities[MBVERTEX]), doEdges(doEntities[MBEDGE]),
30 doQuads(doEntities[MBQUAD]), doTris(doEntities[MBTRI]),
31 doTets(doEntities[MBTET]), doPrisms(doEntities[MBPRISM]) {
32
33 /// This not yet implemented, switch off.
34 doEntities[MBPOLYGON] = false;
35 doEntities[MBPYRAMID] = false;
36 doEntities[MBKNIFE] = false;
37 doEntities[MBPOLYHEDRON] = false;
38}
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
bool & doTris
\deprectaed
bool & doPrisms
\deprectaed
bool sYmm
If true assume that matrix is symmetric structure.
bool & doVertices
\deprectaed If false skip vertices
bool & doTets
\deprectaed
bool & doEdges
\deprectaed If false skip edges
bool & doQuads
\deprectaed

◆ ~DataOperator()

virtual MoFEM::DataOperator::~DataOperator ( )
virtualdefault

Member Function Documentation

◆ doWork() [1/2]

virtual MoFEMErrorCode MoFEM::DataOperator::doWork ( int row_side,
int col_side,
EntityType row_type,
EntityType col_type,
EntitiesFieldData::EntData & row_data,
EntitiesFieldData::EntData & col_data )
inlinevirtual

Operator for bi-linear form, usually to calculate values on left hand side.

Reimplemented in AnalyticalDirichletBC::ApproxField::OpLhs, EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >, EshelbianPlasticity::OpAssembleVolumePositiveDefine, MetaSpringBC::OpSpringALEMaterialLhs, MetaSpringBC::OpSpringALEMaterialLhs_dX_dX, MetaSpringBC::OpSpringALEMaterialLhs_dX_dx, MetaSpringBC::OpSpringKs, MetaSpringBC::OpSpringKs_dX, MetaSpringBC::SpringALEMaterialVolOnSideLhs, MixTransport::MixTransportElement::OpTauDotSigma_HdivHdiv, MixTransport::MixTransportElement::OpVDivSigma_L2Hdiv, MyOp< OP >, NavierStokesElement::OpAssembleLhs, OpAssemble, OpAssembleBasic< VolUserDataOperator >, OpAssembleVolumePositiveDefine, OpFace, OpSimpleRodK, Poisson2DLagrangeMultiplierOperators::OpBoundaryLhsC, Poisson2DLagrangeMultiplierOperators::OpDomainLhsK, Poisson2DNonhomogeneousOperators::OpBoundaryLhs, Poisson2DNonhomogeneousOperators::OpDomainLhs, PrismOpLhs, QuadOpLhs, ReactionDiffusionEquation::OpAssembleMass, and ReactionDiffusionEquation::OpAssembleStiffLhs< DIM >.

Examples
MagneticElement.hpp, and lorentz_force.cpp.

Definition at line 40 of file DataOperators.hpp.

43 {
45 if (doWorkLhsHook) {
46 CHKERR doWorkLhsHook(this, row_side, col_side, row_type, col_type,
47 row_data, col_data);
48 } else {
49 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
50 "doWork function is not implemented for this operator");
51 }
53 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
DoWorkLhsHookFunType doWorkLhsHook

◆ doWork() [2/2]

virtual MoFEMErrorCode MoFEM::DataOperator::doWork ( int side,
EntityType type,
EntitiesFieldData::EntData & data )
inlinevirtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented in AnalyticalDirichletBC::ApproxField::OpRhs< FUNEVAL >, EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >, EshelbianPlasticity::OpAssembleVolumePositiveDefine, EshelbianPlasticity::OpAssembleVolumeStabilize, EshelbianPlasticity::OpGetScale, Example::BoundaryOp, Example::OpCalcSurfaceAverageTemperature, Example::OpFirst, Example::OpSecond, Example::OpZero, MetaSpringBC::OpGetNormalSpEle, MetaSpringBC::OpGetTangentSpEle, MetaSpringBC::OpSpringFs, MetaSpringBC::OpSpringFsMaterial, MixTransport::MixTransportElement::OpDivTauU_HdivL2, MixTransport::MixTransportElement::OpError, MixTransport::MixTransportElement::OpEvaluateBcOnFluxes, MixTransport::MixTransportElement::OpFluxDivergenceAtGaussPts, MixTransport::MixTransportElement::OpL2Source, MixTransport::MixTransportElement::OpPostProc, MixTransport::MixTransportElement::OpRhsBcOnValues, MixTransport::MixTransportElement::OpSkeleton, MixTransport::MixTransportElement::OpTauDotSigma_HdivHdiv, MixTransport::MixTransportElement::OpValuesAtGaussPts, MixTransport::MixTransportElement::OpValuesGradientAtGaussPts, MixTransport::MixTransportElement::OpVDivSigma_L2Hdiv, MoFEM::OpBaseDerivativesMass< 1 >, MoFEM::OpBaseDerivativesNext< 1 >, MoFEM::OpBaseDerivativesNext< 3 >, MoFEM::OpBaseDerivativesSetHOInvJacobian< 2 >, MoFEM::OpBrokenLoopSide< E >, MoFEM::OpCalcNormL2Tensor0, MoFEM::OpCalcNormL2Tensor1< DIM >, MoFEM::OpCalcNormL2Tensor2< DIM_1, DIM_2 >, MoFEM::OpCalculateBrokenHVecTensorDivergence< Tensor_Dim0, Tensor_Dim1, CoordSys >, MoFEM::OpCalculateBrokenHVecTensorField< Tensor_Dim0, Tensor_Dim1 >, MoFEM::OpCalculateDivergenceVectorFieldValues< Tensor_Dim, COORDINATE_SYSTEM >, MoFEM::OpCalculateHcurlVectorCurl< 1, 2 >, MoFEM::OpCalculateHcurlVectorCurl< 1, 3 >, MoFEM::OpCalculateHcurlVectorCurl< 3, 3 >, MoFEM::OpCalculateHdivVectorDivergence< BASE_DIM, SPACE_DIM >, MoFEM::OpCalculateHdivVectorDivergenceDot< Tensor_Dim1, Tensor_Dim2 >, MoFEM::OpCalculateHOCoords< FIELD_DIM >, MoFEM::OpCalculateHOJacForFaceImpl< 2 >, MoFEM::OpCalculateHOJacForFaceImpl< 3 >, MoFEM::OpCalculateHOJacForVolume, MoFEM::OpCalculateHTensorTensorField< Tensor_Dim0, Tensor_Dim1 >, MoFEM::OpCalculateHVecTensorDivergence< Tensor_Dim0, Tensor_Dim1, CoordSys >, MoFEM::OpCalculateHVecTensorField< Tensor_Dim0, Tensor_Dim1 >, MoFEM::OpCalculateHVecVectorField_General< 3, Field_Dim, double, ublas::row_major, DoubleAllocator >, MoFEM::OpCalculateHVecVectorFieldDot< 3, Field_Dim >, MoFEM::OpCalculateHVecVectorGradient< BASE_DIM, SPACE_DIM >, MoFEM::OpCalculateHVecVectorHessian< BASE_DIM, SPACE_DIM >, MoFEM::OpCalculateInvJacForFatPrism, MoFEM::OpCalculateInvJacForFlatPrism, MoFEM::OpCalculateScalarFieldGradient_General< Tensor_Dim, double, ublas::row_major, DoubleAllocator >, MoFEM::OpCalculateScalarFieldHessian< Tensor_Dim >, MoFEM::OpCalculateScalarFieldValues, MoFEM::OpCalculateScalarFieldValues_General< T, A >, MoFEM::OpCalculateScalarFieldValues_General< double, DoubleAllocator >, MoFEM::OpCalculateScalarFieldValuesFromPetscVecImpl< CTX >, MoFEM::OpCalculateTensor2FieldValues_General< Tensor_Dim0, Tensor_Dim1, T, L, A >, MoFEM::OpCalculateTensor2FieldValues_General< Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator >, MoFEM::OpCalculateTensor2FieldValuesDot< Tensor_Dim0, Tensor_Dim1 >, MoFEM::OpCalculateTensor2SymmetricFieldGradient_General< Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator >, MoFEM::OpCalculateTensor2SymmetricFieldValues< Tensor_Dim >, MoFEM::OpCalculateTensor2SymmetricFieldValuesDot< Tensor_Dim >, MoFEM::OpCalculateTraceFromMat< DIM >, MoFEM::OpCalculateTraceFromSymmMat< DIM >, MoFEM::OpCalculateVectorFieldGradient_General< Tensor_Dim0, Tensor_Dim1, S, double, ublas::row_major, DoubleAllocator >, MoFEM::OpCalculateVectorFieldGradientDot< Tensor_Dim0, Tensor_Dim1 >, MoFEM::OpCalculateVectorFieldHessian< Tensor_Dim0, Tensor_Dim1 >, MoFEM::OpCalculateVectorFieldValues_General< Tensor_Dim, T, L, A >, MoFEM::OpCalculateVectorFieldValues_General< Tensor_Dim, double, ublas::row_major, DoubleAllocator >, MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl< Tensor_Dim, CTX >, MoFEM::OpCopyGeomDataToE< 2 >, MoFEM::OpDGProjectionCoefficients, MoFEM::OpDGProjectionEvaluation, MoFEM::OpDGProjectionMassMatrix, MoFEM::OpGetCoordsAndNormalsOnPrism, MoFEM::OpGetDataAndGradient< RANK, DIM >, MoFEM::OpGetHONormalsOnFace< FIELD_DIM >, MoFEM::OpGetHOTangentOnEdge, MoFEM::OpGetHOTangentsOnEdge< FIELD_DIM >, MoFEM::OpGetTensor0fromFunc, MoFEM::OpGetTensor1fromFunc< SPACE_DIM, BASE_DIM >, MoFEM::OpHOSetContravariantPiolaTransformOnEdge3D, MoFEM::OpHOSetContravariantPiolaTransformOnFace3D, MoFEM::OpHOSetCovariantPiolaTransformOnFace3D, MoFEM::OpInvertMatrix< DIM >, MoFEM::OpLoopRange< E >, MoFEM::OpLoopSide< E >, MoFEM::OpLoopThis< E >, MoFEM::OpMakeHdivFromHcurl, MoFEM::OpMultiplyDeterminantOfJacobianAndWeightsForFatPrisms, MoFEM::OpPostProcMapInMoab< DIM1, DIM2, O >, MoFEM::OpRunParent, MoFEM::OpScaleBaseBySpaceInverseOfMeasure, MoFEM::OpScaleMatrix, MoFEM::OpSchurAssembleBegin, MoFEM::OpSchurAssembleEnd< SchurDGESV >, MoFEM::OpSchurAssembleEnd< SchurDSYSV >, MoFEM::OpSchurZeroRowsAndCols, MoFEM::OpSetBc, MoFEM::OpSetContravariantPiolaTransform, MoFEM::OpSetContravariantPiolaTransformOnEdge2D, MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 2 >, MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 3 >, MoFEM::OpSetContravariantPiolaTransformOnFace, MoFEM::OpSetCovariantPiolaTransform, MoFEM::OpSetCovariantPiolaTransformOnEdge, MoFEM::OpSetCovariantPiolaTransformOnFace2DImpl< 2 >, MoFEM::OpSetCovariantPiolaTransformOnFace, MoFEM::OpSetHOContravariantPiolaTransform, MoFEM::OpSetHOCovariantPiolaTransform, MoFEM::OpSetHOInvJacToScalarBasesImpl, MoFEM::OpSetHOInvJacVectorBase, MoFEM::OpSetHOWeights, MoFEM::OpSetHOWeightsOnEdge, MoFEM::OpSetHOWeightsOnFace, MoFEM::OpSetInvJacH1, MoFEM::OpSetInvJacH1ForFatPrism, MoFEM::OpSetInvJacH1ForFlatPrism, MoFEM::OpSetInvJacHcurlFaceImpl< 2 >, MoFEM::OpSetInvJacHcurlFaceImpl< 3 >, MoFEM::OpSetInvJacHdivAndHcurl, MoFEM::OpSetInvJacSpaceForFaceImpl< 2, 1 >, MoFEM::OpSetInvJacSpaceForFaceImpl< 2, 2 >, MoFEM::OpSetInvJacSpaceForFaceImpl< 3, 1 >, MoFEM::OpSymmetrizeTensor< DIM >, MoFEM::OpTensorTimesSymmetricTensor< DIM_01, DIM_23, S >, MoFEM::OpUnSetBc, MyOp< OP >, MyOp< OP >, NavierStokesElement::OpAssembleRhs, NavierStokesElement::OpCalcDragForce, NavierStokesElement::OpCalcDragTraction, NavierStokesElement::OpCalcVolumeFlux, NavierStokesElement::OpPostProcDrag, NavierStokesElement::OpPostProcVorticity, OpAleLhsPre_dX_dx< S >, OpAssemble, OpAssembleBasic< VolUserDataOperator >, OpAssembleVolumePositiveDefine, OpAssembleVolumeStabilize, OpCalculateEnergy, OpCalculateEshelbyStress, OpCalculateEshelbyStress, OpCalculateHomogeneousStiffness< S >, OpCalculateRotationAndSpatialGradient, OpCalculateStiffnessScaledByDensityField, OpCalculateStrain< D >, OpCalculateStrainAle, OpCalculateStress< S >, OpCheck, OpCheckValsDiffVals, OpCheckValsDiffVals, OpCheckValsDiffVals, OpDivergence, OpDivergence, OpFace, OpFace, OpFace, OpFlux, OpFlux, OpPressure, OpSimpleRodPreStress, OpVolumeAssemble, OpVolumeCalculation, OpVolumeSet, OpVolumeTest, Poisson2DLagrangeMultiplierOperators::OpBoundaryRhsG, Poisson2DLagrangeMultiplierOperators::OpDomainRhsF, Poisson2DNonhomogeneousOperators::OpBoundaryRhs, Poisson2DNonhomogeneousOperators::OpDomainRhs, PrismOp, PrismOpCheck, PrismOpRhs, QuadOpCheck, QuadOpRhs, ReactionDiffusionEquation::OpAssembleSlowRhs, ReactionDiffusionEquation::OpAssembleStiffRhs< DIM >, SkeletonFE, SkeletonFE, SkeletonFE, SkeletonFE::OpFaceSide, SkeletonFE::OpFaceSide, and SkeletonFE::OpFaceSide.

Definition at line 67 of file DataOperators.hpp.

68 {
70 if (doWorkRhsHook) {
71 CHKERR doWorkRhsHook(this, side, type, data);
72 } else {
73 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
74 "doWork function is not implemented for this operator");
75 }
77 }
DoWorkRhsHookFunType doWorkRhsHook

◆ getSymm()

bool MoFEM::DataOperator::getSymm ( ) const
inline

Get if operator uses symmetry of DOFs or not.

If symmetry is used, only not repeating combinations of entities are looped. For an example pair of (Vertex, Edge_0) and (Edge_0, Vertex) will calculate the same matrices only transposed. Implementing that this can be exploited by integrating only one pair.

Returns
true if symmetry

Definition at line 107 of file DataOperators.hpp.

107{ return sYmm; }

◆ opLhs() [1/2]

MoFEMErrorCode MoFEM::DataOperator::opLhs ( EntitiesFieldData & row_data,
EntitiesFieldData & col_data )
virtual

Definition at line 90 of file DataOperators.cpp.

91 {
92 if (getSymm())
93 return opLhs<true>(row_data, col_data);
94 else
95 return opLhs<false>(row_data, col_data);
96}
virtual MoFEMErrorCode opLhs(EntitiesFieldData &row_data, EntitiesFieldData &col_data)
bool getSymm() const
Get if operator uses symmetry of DOFs or not.

◆ opLhs() [2/2]

template<bool Symm>
MoFEMErrorCode MoFEM::DataOperator::opLhs ( EntitiesFieldData & row_data,
EntitiesFieldData & col_data )
inlineprivate

Definition at line 41 of file DataOperators.cpp.

42 {
44
45 auto do_col_entity =
46 [&](boost::ptr_vector<EntitiesFieldData::EntData> &row_ent_data,
47 const int ss, const EntityType row_type, const EntityType low_type,
48 const EntityType hi_type) {
50 for (EntityType col_type = low_type; col_type != hi_type; ++col_type) {
51 auto &col_ent_data = col_data.dataOnEntities[col_type];
52 for (size_t SS = 0; SS != col_ent_data.size(); SS++) {
53 if (col_ent_data[SS].getFieldData().size())
54 CHKERR doWork(ss, SS, row_type, col_type, row_ent_data[ss],
55 col_ent_data[SS]);
56 }
57 }
59 };
60
61 auto do_row_entity = [&](const EntityType type) {
63 auto &row_ent_data = row_data.dataOnEntities[type];
64 for (size_t ss = 0; ss != row_ent_data.size(); ++ss) {
65 if constexpr (!Symm)
66 CHKERR do_col_entity(row_ent_data, ss, type, MBVERTEX, type);
67 size_t SS = 0;
68 if constexpr (Symm)
69 SS = ss;
70 for (; SS < col_data.dataOnEntities[type].size(); ++SS) {
71 CHKERR doWork(ss, SS, type, type, row_ent_data[ss],
72 col_data.dataOnEntities[type][SS]);
73 }
74 CHKERR do_col_entity(row_ent_data, ss, type,
75 static_cast<EntityType>(type + 1), MBMAXTYPE);
76 }
78 };
79
80 for (EntityType row_type = MBVERTEX; row_type != MBMAXTYPE; ++row_type) {
81 if (doEntities[row_type]) {
82 CHKERR do_row_entity(row_type);
83 }
84 }
85
86
88}
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.

◆ opRhs() [1/2]

MoFEMErrorCode MoFEM::DataOperator::opRhs ( EntitiesFieldData & data,
const bool error_if_no_base = false )
virtual

Reimplemented in MoFEM::OpAddParentEntData.

Definition at line 140 of file DataOperators.cpp.

141 {
142 if (error_if_no_base)
143 return opRhs<true>(data, doEntities);
144 else
145 return opRhs<false>(data, doEntities);
146}
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)

◆ opRhs() [2/2]

template<bool ErrorIfNoBase>
MoFEMErrorCode MoFEM::DataOperator::opRhs ( EntitiesFieldData & data,
const std::array< bool, MBMAXTYPE > & do_entities )
inlineprivate

Definition at line 100 of file DataOperators.cpp.

101 {
103
104 auto do_entity = [&](auto type) {
106
107 auto &ent_data = data.dataOnEntities[type];
108 const size_t size = ent_data.size();
109 for (size_t ss = 0; ss != size; ++ss) {
110
111 auto &side_data = ent_data[ss];
112
113 if (ErrorIfNoBase) {
114 if (side_data.getFieldData().size() &&
115 (side_data.getBase() == NOBASE ||
116 side_data.getBase() == LASTBASE)) {
117 for (VectorDofs::iterator it = side_data.getFieldDofs().begin();
118 it != side_data.getFieldDofs().end(); it++)
119 if ((*it) && (*it)->getActive())
120 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "No base on");
121 }
122 }
123
124 CHKERR doWork(ss, type, side_data);
125 }
126
128 };
129
130 for (EntityType row_type = MBVERTEX; row_type != MBMAXTYPE; ++row_type) {
131 if (do_entities[row_type]) {
132 CHKERR do_entity(row_type);
133 }
134 }
135
136
138}
@ LASTBASE
Definition definitions.h:69
@ NOBASE
Definition definitions.h:59
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31

◆ setSymm()

void MoFEM::DataOperator::setSymm ( )
inline

set if operator is executed taking in account symmetry

Definition at line 110 of file DataOperators.hpp.

110{ sYmm = true; }

◆ unSetSymm()

void MoFEM::DataOperator::unSetSymm ( )
inline

unset if operator is executed for non symmetric problem

Definition at line 113 of file DataOperators.hpp.

113{ sYmm = false; }

Member Data Documentation

◆ doEdges

bool& MoFEM::DataOperator::doEdges

\deprectaed If false skip edges

Examples
ElasticityMixedFormulation.hpp, and NavierStokesElement.hpp.

Definition at line 91 of file DataOperators.hpp.

◆ doEntities

std::array<bool, MBMAXTYPE> MoFEM::DataOperator::doEntities

If true operator is executed for entity.

Examples
scalar_check_approximation.cpp.

Definition at line 85 of file DataOperators.hpp.

◆ doPrisms

bool& MoFEM::DataOperator::doPrisms

\deprectaed

Examples
ElasticityMixedFormulation.hpp, and NavierStokesElement.hpp.

Definition at line 95 of file DataOperators.hpp.

◆ doQuads

bool& MoFEM::DataOperator::doQuads

\deprectaed

Examples
ElasticityMixedFormulation.hpp, and NavierStokesElement.hpp.

Definition at line 92 of file DataOperators.hpp.

◆ doTets

bool& MoFEM::DataOperator::doTets

\deprectaed

Examples
ElasticityMixedFormulation.hpp, and NavierStokesElement.hpp.

Definition at line 94 of file DataOperators.hpp.

◆ doTris

bool& MoFEM::DataOperator::doTris

\deprectaed

Examples
ElasticityMixedFormulation.hpp, and NavierStokesElement.hpp.

Definition at line 93 of file DataOperators.hpp.

◆ doVertices

bool& MoFEM::DataOperator::doVertices

\deprectaed If false skip vertices

Examples
ElasticityMixedFormulation.hpp, and NavierStokesElement.hpp.

Definition at line 90 of file DataOperators.hpp.

◆ doWorkLhsHook

DoWorkLhsHookFunType MoFEM::DataOperator::doWorkLhsHook

Definition at line 35 of file DataOperators.hpp.

◆ doWorkRhsHook

DoWorkRhsHookFunType MoFEM::DataOperator::doWorkRhsHook

Definition at line 62 of file DataOperators.hpp.

◆ sYmm

bool MoFEM::DataOperator::sYmm

If true assume that matrix is symmetric structure.

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

Definition at line 82 of file DataOperators.hpp.


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