v0.15.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
MoFEM::OpBaseImpl< A, EleOp > Struct Template Reference

#include "src/finite_elements/FormsIntegrators.hpp"

Inheritance diagram for MoFEM::OpBaseImpl< A, EleOp >:
[legend]
Collaboration diagram for MoFEM::OpBaseImpl< A, EleOp >:
[legend]

Public Types

using OpType = typename EleOp::OpType
 
using EntData = EntitiesFieldData::EntData
 
using MatSetValuesHook = boost::function< MoFEMErrorCode(ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)>
 

Public Member Functions

 OpBaseImpl (const std::string row_field_name, const std::string col_field_name, const OpType type, boost::shared_ptr< Range > ents_ptr=nullptr)
 Constructor for base operator implementation.
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Do calculations for the left hand side.
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntData &row_data)
 Do calculations for the right hand side.
 

Public Attributes

TimeFun timeScalingFun
 assumes that time variable is set
 
FEFun feScalingFun
 assumes that time variable is set
 
boost::shared_ptr< RangeentsPtr
 Entities on which element is run.
 

Static Public Attributes

static MatSetValuesHook matSetValuesHook
 

Protected Member Functions

template<int DIM>
FTensor::Tensor1< FTensor::PackPtr< double *, DIM >, DIM > getNf ()
 Get local vector tensor for assembly.
 
template<int DIM>
FTensor::Tensor2< FTensor::PackPtr< double *, DIM >, DIM, DIM > getLocMat (const int rr)
 Get local matrix tensor for assembly.
 
virtual MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Integrate grad-grad operator.
 
virtual MoFEMErrorCode aSsemble (EntData &row_data, EntData &col_data, const bool trans)
 Assemble local matrix into global matrix.
 
virtual MoFEMErrorCode iNtegrate (EntData &data)
 Class dedicated to integrate operator.
 
virtual MoFEMErrorCode aSsemble (EntData &data)
 Assemble local vector into global vector.
 
virtual size_t getNbOfBaseFunctions (EntitiesFieldData::EntData &data)
 Get number of base functions.
 

Protected Attributes

int nbRows
 number of dofs on rows
 
int nbCols
 number if dof on column
 
int nbIntegrationPts
 number of integration points
 
int nbRowBaseFunctions
 number or row base functions
 
int rowSide
 row side number
 
int colSide
 column side number
 
EntityType rowType
 row type
 
EntityType colType
 column type
 
bool assembleTranspose
 
bool onlyTranspose
 
MatrixDouble locMat
 local entity block matrix
 
MatrixDouble locMatTranspose
 local entity block matrix
 
VectorDouble locF
 local entity vector
 

Detailed Description

template<AssemblyType A, typename EleOp>
struct MoFEM::OpBaseImpl< A, EleOp >
Examples
ThermoElasticOps.hpp, mofem/atom_tests/schur_test_diag_mat.cpp, mofem/tutorials/adv-2/src/ThermoElasticOps.hpp, mofem/tutorials/adv-2/thermo_elastic.cpp, mofem/tutorials/scl-8/radiation.cpp, thermo_elastic.cpp, and thermoplastic.cpp.

Definition at line 284 of file FormsIntegrators.hpp.

Member Typedef Documentation

◆ EntData

Definition at line 286 of file FormsIntegrators.hpp.

◆ MatSetValuesHook

template<AssemblyType A, typename EleOp >
using MoFEM::OpBaseImpl< A, EleOp >::MatSetValuesHook = boost::function<MoFEMErrorCode( ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)>

Definition at line 328 of file FormsIntegrators.hpp.

◆ OpType

template<AssemblyType A, typename EleOp >
using MoFEM::OpBaseImpl< A, EleOp >::OpType = typename EleOp::OpType

Definition at line 285 of file FormsIntegrators.hpp.

Constructor & Destructor Documentation

◆ OpBaseImpl()

template<AssemblyType A, typename EleOp >
MoFEM::OpBaseImpl< A, EleOp >::OpBaseImpl ( const std::string  row_field_name,
const std::string  col_field_name,
const OpType  type,
boost::shared_ptr< Range ents_ptr = nullptr 
)
inline

Constructor for base operator implementation.

Parameters
row_field_nameName of field for row entities
col_field_nameName of field for column entities
typeOperator type
ents_ptrOptional pointer to range of entities to operate on

Definition at line 295 of file FormsIntegrators.hpp.

297 : EleOp(row_field_name, col_field_name, type, false),
298 assembleTranspose(false), onlyTranspose(false), entsPtr(ents_ptr) {}
Ele::UserDataOperator EleOp
boost::shared_ptr< Range > entsPtr
Entities on which element is run.

Member Function Documentation

◆ aSsemble() [1/2]

template<AssemblyType A, typename EleOp >
MoFEMErrorCode OpBaseImpl::aSsemble ( EntData data)
protectedvirtual

Assemble local vector into global vector.

Parameters
dataEntity data
Returns
Error code

Definition at line 640 of file FormsIntegrators.hpp.

640 {
641 if (!this->timeScalingFun.empty())
642 this->locF *= this->timeScalingFun(this->getFEMethod()->ts_t);
643 if (!this->feScalingFun.empty())
644 this->locF *= this->feScalingFun(this->getFEMethod());
645
646 return VecSetValues<AssemblyTypeSelector<A>>(this->getKSPf(), data,
647 this->locF, ADD_VALUES);
648}
VectorDouble locF
local entity vector
TimeFun timeScalingFun
assumes that time variable is set
FEFun feScalingFun
assumes that time variable is set

◆ aSsemble() [2/2]

template<AssemblyType A, typename EleOp >
MoFEMErrorCode OpBaseImpl::aSsemble ( EntData row_data,
EntData col_data,
const bool  trans 
)
protectedvirtual

Assemble local matrix into global matrix.

Parameters
row_dataRow entity data
col_dataColumn entity data
transWhether to assemble transpose
Returns
Error code

Definition at line 613 of file FormsIntegrators.hpp.

615 {
617
618 if (!this->timeScalingFun.empty())
619 this->locMat *= this->timeScalingFun(this->getFEMethod()->ts_t);
620 if (!this->feScalingFun.empty())
621 this->locMat *= this->feScalingFun(this->getFEMethod());
622
623 // Assemble transpose
624 if (transpose) {
625 this->locMatTranspose.resize(this->locMat.size2(), this->locMat.size1(),
626 false);
627 noalias(this->locMatTranspose) = trans(this->locMat);
628 CHKERR matSetValuesHook(this, col_data, row_data, this->locMatTranspose);
629 }
630
631 if (!this->onlyTranspose) {
632 // assemble local matrix
633 CHKERR matSetValuesHook(this, row_data, col_data, this->locMat);
634 }
635
637}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MatrixDouble locMatTranspose
local entity block matrix
static MatSetValuesHook matSetValuesHook
MatrixDouble locMat
local entity block matrix

◆ doWork() [1/2]

template<AssemblyType A, typename EleOp >
MoFEMErrorCode OpBaseImpl::doWork ( int  row_side,
EntityType  row_type,
EntData row_data 
)

Do calculations for the right hand side.

Parameters
row_side
row_type
row_data
Returns
MoFEMErrorCode

Definition at line 582 of file FormsIntegrators.hpp.

583 {
585
586 if (entsPtr) {
587 if (entsPtr->find(this->getFEEntityHandle()) == entsPtr->end())
589 }
590
591 // get number of dofs on row
592 nbRows = row_data.getIndices().size();
593 rowSide = row_side;
594 rowType = row_type;
595
596 if (!nbRows)
598 // get number of integration points
599 nbIntegrationPts = EleOp::getGaussPts().size2();
600 // get row base functions
602 // resize and clear the right hand side vector
603 locF.resize(nbRows);
604 locF.clear();
605 // integrate local vector
606 CHKERR this->iNtegrate(row_data);
607 // assemble local vector
608 CHKERR this->aSsemble(row_data);
610}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
int rowSide
row side number
int nbRows
number of dofs on rows
int nbIntegrationPts
number of integration points
virtual MoFEMErrorCode aSsemble(EntData &row_data, EntData &col_data, const bool trans)
Assemble local matrix into global matrix.
virtual MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
Integrate grad-grad operator.
virtual size_t getNbOfBaseFunctions(EntitiesFieldData::EntData &data)
Get number of base functions.
int nbRowBaseFunctions
number or row base functions
EntityType rowType
row type

◆ doWork() [2/2]

template<AssemblyType A, typename EleOp >
MoFEMErrorCode OpBaseImpl::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)

Do calculations for the left hand side.

Parameters
row_siderow side number (local number) of entity on element
col_sidecolumn side number (local number) of entity on element
row_typetype of row entity
col_typetype of column entity
row_datadata for row
col_datadata for column
Returns
error code

Definition at line 527 of file FormsIntegrators.hpp.

530 {
532
533 if (entsPtr) {
534 if (entsPtr->find(this->getFEEntityHandle()) == entsPtr->end())
536 }
537
538 // get number of dofs on row
539 nbRows = row_data.getIndices().size();
540 // if no dofs on row, exit that work, nothing to do here
541 if (!nbRows)
543 rowSide = row_side;
544 rowType = row_type;
545 // get number of dofs on column
546 nbCols = col_data.getIndices().size();
547 // if no dofs on column, exit nothing to do here
548 if (!nbCols)
550 colSide = col_side;
551 colType = col_type;
552 // get number of integration points
553 nbIntegrationPts = EleOp::getGaussPts().size2();
554 // get row base functions
556
557 // set size of local entity bock
558 locMat.resize(nbRows, nbCols, false);
559 // clear matrix
560 locMat.clear();
561 // integrate local matrix for entity block
562 CHKERR this->iNtegrate(row_data, col_data);
563
564 // assemble local matrix
565 auto check_if_assemble_transpose = [&] {
566 if (this->sYmm) {
567 if (row_side != col_side || row_type != col_type)
568 return true;
569 else
570 return false;
571 } else if (assembleTranspose) {
572 return true;
573 }
574
575 return false;
576 };
577 CHKERR aSsemble(row_data, col_data, check_if_assemble_transpose());
579}
int colSide
column side number
EntityType colType
column type
int nbCols
number if dof on column

◆ getLocMat()

template<AssemblyType A, typename EleOp >
template<int DIM>
FTensor::Tensor2< FTensor::PackPtr< double *, DIM >, DIM, DIM > MoFEM::OpBaseImpl< A, EleOp >::getLocMat ( const int  rr)
inlineprotected

Get local matrix tensor for assembly.

Template Parameters
DIMDimension of the tensor
Parameters
rrRow offset in the local matrix
Returns
FTensor wrapper around local matrix

Definition at line 356 of file FormsIntegrators.hpp.

356 {
357 return getFTensor2FromArray<DIM, DIM, DIM>(locMat, rr);
358 }

◆ getNbOfBaseFunctions()

template<AssemblyType A, typename EleOp >
size_t OpBaseImpl::getNbOfBaseFunctions ( EntitiesFieldData::EntData data)
protectedvirtual

Get number of base functions.

Parameters
data
Returns
number of base functions

Definition at line 487 of file FormsIntegrators.hpp.

487 {
488 if (!data.getNSharedPtr(data.getBase())) {
489 #ifndef NDEBUG
490 MOFEM_LOG("SELF", Sev::warning)
491 << "Ptr to base " << ApproximationBaseNames[data.getBase()]
492 << " functions is null, returning 0";
493 #endif // NDEBUG
494 return 0;
495 }
496 auto nb_base_functions = data.getN().size2();
497 if (data.getBase() != USER_BASE) {
498 switch (data.getSpace()) {
499 case NOSPACE:
500 break;
501 case NOFIELD:
502 break;
503 case H1:
504 break;
505 case HCURL:
506 case HDIV:
507 nb_base_functions /= 3;
508#ifndef NDEBUG
509 if (data.getN().size2() % 3) {
510 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
511 "Number of base functions is not divisible by 3");
512 }
513#endif
514 break;
515 case L2:
516 break;
517 default:
518 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
519 "Space %s not implemented", FieldSpaceNames[data.getSpace()]);
520 }
521 }
522 return nb_base_functions;
523}
@ USER_BASE
user implemented approximation base
Definition definitions.h:68
@ L2
field with C-1 continuity
Definition definitions.h:88
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition definitions.h:84
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
static const char *const FieldSpaceNames[]
Definition definitions.h:92
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
static const char *const ApproximationBaseNames[]
Definition definitions.h:72
#define MOFEM_LOG(channel, severity)
Log.

◆ getNf()

template<AssemblyType A, typename EleOp >
template<int DIM>
FTensor::Tensor1< FTensor::PackPtr< double *, DIM >, DIM > MoFEM::OpBaseImpl< A, EleOp >::getNf ( )
inlineprotected

Get local vector tensor for assembly.

Template Parameters
DIMDimension of the tensor
Returns
FTensor wrapper around local vector

Definition at line 344 of file FormsIntegrators.hpp.

344 {
345 return getFTensor1FromArray<DIM, DIM>(locF);
346 }

◆ iNtegrate() [1/2]

template<AssemblyType A, typename EleOp >
virtual MoFEMErrorCode MoFEM::OpBaseImpl< A, EleOp >::iNtegrate ( EntData data)
inlineprotectedvirtual

Class dedicated to integrate operator.

Parameters
dataentity data on element row
Returns
error code

Reimplemented in Example::OpRadiationRhs, Example::OpFluxRhs, MoFEM::OpSourceImpl< 1, 1, GAUSS, SourceFunctionSpecialization::S< OpBase > >, MoFEM::OpSourceImpl< 1, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >, MoFEM::OpSourceImpl< 3, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >, MoFEM::OpBaseTimesScalarImpl< 1, S, GAUSS, OpBase >, MoFEM::OpBaseTimesVectorImpl< 1, FIELD_DIM, S, GAUSS, OpBase >, MoFEM::OpBaseTimesVectorImpl< 3, FIELD_DIM, S, GAUSS, OpBase >, MoFEM::OpGradTimesTensorImpl< 1, 1, SPACE_DIM, S, GAUSS, OpBase >, MoFEM::OpGradTimesTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >, MoFEM::OpGradTimesSymTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >, MoFEM::OpMixDivTimesUImpl< 3, FIELD_DIM, SPACE_DIM, GAUSS, OpBase, CoordSys >, MoFEM::OpMixDivTimesUImpl< 3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys >, MoFEM::OpMixDivTimesUImpl< 1, FIELD_DIM, FIELD_DIM, GAUSS, OpBase, CoordSys >, MoFEM::OpMixVecTimesDivLambdaImpl< SPACE_DIM, GAUSS, OpBase >, MoFEM::OpMixTensorTimesGradUImpl< SPACE_DIM, GAUSS, OpBase >, MoFEM::OpNormalMixVecTimesScalarImpl< 3, GAUSS, OpBase >, MoFEM::OpNormalMixVecTimesScalarImpl< 2, GAUSS, OpBase >, MoFEM::OpNormalMixVecTimesVectorFieldImpl< SPACE_DIM, GAUSS, OpBase >, MoFEM::OpConvectiveTermRhsImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpConvectiveTermRhsImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >, and MoFEM::OpBrokenSpaceConstrainDHybridImpl< FIELD_DIM, GAUSS, OpBase >.

Definition at line 402 of file FormsIntegrators.hpp.

402 {
404 }

◆ iNtegrate() [2/2]

template<AssemblyType A, typename EleOp >
virtual MoFEMErrorCode MoFEM::OpBaseImpl< A, EleOp >::iNtegrate ( EntData row_data,
EntData col_data 
)
inlineprotectedvirtual

Integrate grad-grad operator.

Parameters
row_datarow data (consist base functions on row entity)
col_datacolumn data (consist base functions on column entity)
Returns
error code

Reimplemented in Example::OpRadiationLhs, MoFEM::OpGradGradImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpGradGradImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpMassImpl< 1, 1, GAUSS, OpBase >, MoFEM::OpMassImpl< 1, FIELD_DIM, GAUSS, OpBase >, MoFEM::OpMassImpl< 3, FIELD_DIM, GAUSS, OpBase >, MoFEM::OpMassImpl< 3, 4, GAUSS, OpBase >, MoFEM::OpMassImpl< 3, 9, GAUSS, OpBase >, MoFEM::OpGradSymTensorGradImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >, MoFEM::OpGradTensorGradImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >, MoFEM::OpGradGradSymTensorGradGradImpl< 1, 1, SPACE_DIM, S, GAUSS, OpBase >, MoFEM::OpMixDivTimesScalarImpl< SPACE_DIM, GAUSS, OpBase >, MoFEM::OpMixDivTimesVecImpl< SPACE_DIM, GAUSS, OpBase, CoordSys >, MoFEM::OpMixScalarTimesDivImpl< SPACE_DIM, GAUSS, OpBase, COORDINATE_SYSTEM >, MoFEM::OpMixVectorTimesGradImpl< 3, SPACE_DIM, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpMixVectorTimesGradImpl< 1, SPACE_DIM, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpMixTensorTimesGradImpl< SPACE_DIM, GAUSS, OpBase >, MoFEM::OpConvectiveTermLhsDuImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpConvectiveTermLhsDyImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >, MoFEM::OpConvectiveTermLhsDuImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >, and MoFEM::OpConvectiveTermLhsDyImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >.

Definition at line 383 of file FormsIntegrators.hpp.

383 {
385 }

Member Data Documentation

◆ assembleTranspose

template<AssemblyType A, typename EleOp >
bool MoFEM::OpBaseImpl< A, EleOp >::assembleTranspose
protected

Definition at line 370 of file FormsIntegrators.hpp.

◆ colSide

template<AssemblyType A, typename EleOp >
int MoFEM::OpBaseImpl< A, EleOp >::colSide
protected

column side number

Definition at line 366 of file FormsIntegrators.hpp.

◆ colType

template<AssemblyType A, typename EleOp >
EntityType MoFEM::OpBaseImpl< A, EleOp >::colType
protected

column type

Definition at line 368 of file FormsIntegrators.hpp.

◆ entsPtr

template<AssemblyType A, typename EleOp >
boost::shared_ptr<Range> MoFEM::OpBaseImpl< A, EleOp >::entsPtr

Entities on which element is run.

Definition at line 326 of file FormsIntegrators.hpp.

◆ feScalingFun

template<AssemblyType A, typename EleOp >
FEFun MoFEM::OpBaseImpl< A, EleOp >::feScalingFun

assumes that time variable is set

Definition at line 325 of file FormsIntegrators.hpp.

◆ locF

template<AssemblyType A, typename EleOp >
VectorDouble MoFEM::OpBaseImpl< A, EleOp >::locF
protected

local entity vector

Definition at line 375 of file FormsIntegrators.hpp.

◆ locMat

template<AssemblyType A, typename EleOp >
MatrixDouble MoFEM::OpBaseImpl< A, EleOp >::locMat
protected

local entity block matrix

Examples
mofem/tutorials/scl-8/radiation.cpp.

Definition at line 373 of file FormsIntegrators.hpp.

◆ locMatTranspose

template<AssemblyType A, typename EleOp >
MatrixDouble MoFEM::OpBaseImpl< A, EleOp >::locMatTranspose
protected

local entity block matrix

Definition at line 374 of file FormsIntegrators.hpp.

◆ matSetValuesHook

template<AssemblyType A, typename EleOp >
OpBaseImpl< A, EleOp >::MatSetValuesHook OpBaseImpl::matSetValuesHook
static
Initial value:
=
[](ForcesAndSourcesCore::UserDataOperator *op_ptr,
const EntitiesFieldData::EntData &row_data,
const EntitiesFieldData::EntData &col_data, MatrixDouble &m) {
return MatSetValues<AssemblyTypeSelector<A>>(
op_ptr->getKSPB(), row_data, col_data, m, ADD_VALUES);
}
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
FTensor::Index< 'm', 3 > m

Definition at line 333 of file FormsIntegrators.hpp.

◆ nbCols

template<AssemblyType A, typename EleOp >
int MoFEM::OpBaseImpl< A, EleOp >::nbCols
protected

number if dof on column

Examples
mofem/tutorials/scl-8/radiation.cpp.

Definition at line 361 of file FormsIntegrators.hpp.

◆ nbIntegrationPts

template<AssemblyType A, typename EleOp >
int MoFEM::OpBaseImpl< A, EleOp >::nbIntegrationPts
protected

number of integration points

Examples
mofem/tutorials/scl-8/radiation.cpp.

Definition at line 362 of file FormsIntegrators.hpp.

◆ nbRowBaseFunctions

template<AssemblyType A, typename EleOp >
int MoFEM::OpBaseImpl< A, EleOp >::nbRowBaseFunctions
protected

number or row base functions

Definition at line 363 of file FormsIntegrators.hpp.

◆ nbRows

template<AssemblyType A, typename EleOp >
int MoFEM::OpBaseImpl< A, EleOp >::nbRows
protected

number of dofs on rows

Examples
mofem/tutorials/scl-8/radiation.cpp.

Definition at line 360 of file FormsIntegrators.hpp.

◆ onlyTranspose

template<AssemblyType A, typename EleOp >
bool MoFEM::OpBaseImpl< A, EleOp >::onlyTranspose
protected

Definition at line 371 of file FormsIntegrators.hpp.

◆ rowSide

template<AssemblyType A, typename EleOp >
int MoFEM::OpBaseImpl< A, EleOp >::rowSide
protected

row side number

Definition at line 365 of file FormsIntegrators.hpp.

◆ rowType

template<AssemblyType A, typename EleOp >
EntityType MoFEM::OpBaseImpl< A, EleOp >::rowType
protected

row type

Definition at line 367 of file FormsIntegrators.hpp.

◆ timeScalingFun

template<AssemblyType A, typename EleOp >
TimeFun MoFEM::OpBaseImpl< A, EleOp >::timeScalingFun

assumes that time variable is set

Definition at line 324 of file FormsIntegrators.hpp.


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