![]() |
v0.15.0 |
Volume finite element base. More...
#include "src/finite_elements/VolumeElementForcesAndSourcesCore.hpp"
Classes | |
struct | UserDataOperator |
Public Types | |
enum | Switches |
![]() | |
typedef boost::function< int(int order_row, int order_col, int order_data)> | RuleHookFun |
typedef boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> | GaussHookFun |
![]() | |
enum | KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE } |
pass information about context of KSP/DM for with finite element is computed More... | |
![]() | |
enum | DataContext { CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 , CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 , CTX_SET_TIME = 1 << 7 } |
using | Switches = std::bitset<8> |
![]() | |
enum | SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE } |
![]() | |
enum | TSContext { CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN , CTX_TSTSMONITORSET , CTX_TSNONE } |
![]() | |
enum | TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE } |
Public Member Functions | |
VolumeElementForcesAndSourcesCore (Interface &m_field) | |
DEPRECATED | VolumeElementForcesAndSourcesCore (Interface &m_field, const EntityType type) |
MoFEMErrorCode | operator() () |
function is run for every finite element | |
![]() | |
ForcesAndSourcesCore (Interface &m_field) | |
boost::ptr_deque< UserDataOperator > & | getOpPtrVector () |
Use to push back operator for row operator. | |
auto & | getElementPolynomialBase () |
Get the Entity Polynomial Base object. | |
auto & | getUserPolynomialBase () |
Get the User Polynomial Base object. | |
virtual MoFEMErrorCode | preProcess () |
function is run at the beginning of loop | |
virtual MoFEMErrorCode | postProcess () |
function is run at the end of loop | |
int | getMaxDataOrder () const |
Get max order of approximation for data fields. | |
int | getMaxRowOrder () const |
Get max order of approximation for field in rows. | |
int | getMaxColOrder () const |
Get max order of approximation for field in columns. | |
auto & | getEntData (const FieldSpace space, const EntityType type, const int side) |
Get the entity data. | |
auto & | getDataOnElementBySpaceArray () |
Get data on entities and space. | |
auto & | getDerivedDataOnElementBySpaceArray () |
Get derived data on entities and space. | |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
FEMethod ()=default | |
auto | getFEName () const |
get finite element name | |
auto | getDataDofsPtr () const |
auto | getDataVectorDofsPtr () const |
const FieldEntity_vector_view & | getDataFieldEnts () const |
boost::shared_ptr< FieldEntity_vector_view > & | getDataFieldEntsPtr () const |
auto & | getRowFieldEnts () const |
auto & | getRowFieldEntsPtr () const |
auto & | getColFieldEnts () const |
auto & | getColFieldEntsPtr () const |
auto | getRowDofsPtr () const |
auto | getColDofsPtr () const |
auto | getNumberOfNodes () const |
EntityHandle | getFEEntityHandle () const |
MoFEMErrorCode | getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true) |
![]() | |
BasicMethod () | |
virtual | ~BasicMethod ()=default |
int | getNinTheLoop () const |
get number of evaluated element in the loop | |
int | getLoopSize () const |
get loop size | |
auto | getLoHiFERank () const |
Get lo and hi processor rank of iterated entities. | |
auto | getLoFERank () const |
Get upper rank in loop for iterating elements. | |
auto | getHiFERank () const |
Get upper rank in loop for iterating elements. | |
unsigned int | getFieldBitNumber (std::string field_name) const |
MoFEMErrorCode | copyBasicMethod (const BasicMethod &basic) |
Copy data from other base method to this base method. | |
boost::weak_ptr< CacheTuple > | getCacheWeakPtr () const |
Get the cache weak ptr object. | |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
KspMethod () | |
virtual | ~KspMethod ()=default |
MoFEMErrorCode | copyKsp (const KspMethod &ksp) |
copy data form another method | |
![]() | |
PetscData () | |
virtual | ~PetscData ()=default |
MoFEMErrorCode | copyPetscData (const PetscData &petsc_data) |
![]() | |
template<class IFACE > | |
MoFEMErrorCode | registerInterface (bool error_if_registration_failed=true) |
Register interface. | |
template<class IFACE > | |
MoFEMErrorCode | getInterface (IFACE *&iface) const |
Get interface reference to pointer of interface. | |
template<class IFACE > | |
MoFEMErrorCode | getInterface (IFACE **const iface) const |
Get interface pointer to pointer of interface. | |
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0> | |
IFACE | getInterface () const |
Get interface pointer to pointer of interface. | |
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0> | |
IFACE | getInterface () const |
Get reference to interface. | |
template<class IFACE > | |
IFACE * | getInterface () const |
Function returning pointer to interface. | |
virtual | ~UnknownInterface ()=default |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
SnesMethod () | |
virtual | ~SnesMethod ()=default |
MoFEMErrorCode | copySnes (const SnesMethod &snes) |
Copy snes data. | |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
TSMethod () | |
virtual | ~TSMethod ()=default |
MoFEMErrorCode | copyTs (const TSMethod &ts) |
Copy TS solver data. | |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
TaoMethod () | |
virtual | ~TaoMethod ()=default |
MoFEMErrorCode | copyTao (const TaoMethod &tao) |
Copy TAO data. | |
Public Attributes | |
std::string | meshPositionsFieldName |
![]() | |
Interface & | mField |
RuleHookFun | getRuleHook |
Hook to get rule. | |
GaussHookFun | setRuleHook |
Set function to calculate integration rule. | |
MatrixDouble | gaussPts |
Matrix of integration points. | |
![]() | |
std::string | feName |
Name of finite element. | |
boost::shared_ptr< const NumeredEntFiniteElement > | numeredEntFiniteElementPtr |
boost::function< bool(FEMethod *fe_method_ptr)> | exeTestHook |
Tet if element to skip element. | |
![]() | |
int | nInTheLoop |
number currently of processed method | |
int | loopSize |
local number oe methods to process | |
std::pair< int, int > | loHiFERank |
Llo and hi processor rank of iterated entities. | |
int | rAnk |
processor rank | |
int | sIze |
number of processors in communicator | |
const RefEntity_multiIndex * | refinedEntitiesPtr |
container of mofem dof entities | |
const RefElement_multiIndex * | refinedFiniteElementsPtr |
container of mofem finite element entities | |
const Problem * | problemPtr |
raw pointer to problem | |
const Field_multiIndex * | fieldsPtr |
raw pointer to fields container | |
const FieldEntity_multiIndex * | entitiesPtr |
raw pointer to container of field entities | |
const DofEntity_multiIndex * | dofsPtr |
raw pointer container of dofs | |
const FiniteElement_multiIndex * | finiteElementsPtr |
raw pointer to container finite elements | |
const EntFiniteElement_multiIndex * | finiteElementsEntitiesPtr |
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * | adjacenciesPtr |
boost::function< MoFEMErrorCode()> | preProcessHook |
Hook function for pre-processing. | |
boost::function< MoFEMErrorCode()> | postProcessHook |
Hook function for post-processing. | |
boost::function< MoFEMErrorCode()> | operatorHook |
Hook function for operator. | |
boost::movelib::unique_ptr< bool > | vecAssembleSwitch |
boost::movelib::unique_ptr< bool > | matAssembleSwitch |
boost::weak_ptr< CacheTuple > | cacheWeakPtr |
![]() | |
KSPContext | ksp_ctx |
Context. | |
KSP | ksp |
KSP solver. | |
Vec & | ksp_f |
Mat & | ksp_A |
Mat & | ksp_B |
![]() | |
Switches | data_ctx |
Vec | f |
Mat | A |
Mat | B |
Vec | x |
Vec | dx |
Vec | x_t |
Vec | x_tt |
![]() | |
SNESContext | snes_ctx |
SNES | snes |
snes solver | |
Vec & | snes_x |
state vector | |
Vec & | snes_dx |
solution update | |
Vec & | snes_f |
residual | |
Mat & | snes_A |
jacobian matrix | |
Mat & | snes_B |
preconditioner of jacobian matrix | |
![]() | |
TS | ts |
time solver | |
TSContext | ts_ctx |
PetscInt | ts_step |
time step number | |
PetscReal | ts_a |
shift for U_t (see PETSc Time Solver) | |
PetscReal | ts_aa |
shift for U_tt shift for U_tt | |
PetscReal | ts_t |
time | |
PetscReal | ts_dt |
time step size | |
Vec & | ts_u |
state vector | |
Vec & | ts_u_t |
time derivative of state vector | |
Vec & | ts_u_tt |
second time derivative of state vector | |
Vec & | ts_F |
residual vector | |
Mat & | ts_A |
Mat & | ts_B |
Preconditioner for ts_A. | |
![]() | |
TAOContext | tao_ctx |
Tao | tao |
tao solver | |
Vec & | tao_x |
Vec & | tao_f |
state vector | |
Mat & | tao_A |
gradient vector | |
Mat & | tao_B |
hessian matrix | |
Protected Member Functions | |
virtual MoFEMErrorCode | setIntegrationPts () |
Set integration points. | |
virtual MoFEMErrorCode | calculateVolumeAndJacobian () |
Calculate element volume and Jacobian. | |
virtual MoFEMErrorCode | calculateCoordinatesAtGaussPts () |
Calculate coordinate at integration points. | |
virtual MoFEMErrorCode | getSpaceBaseAndOrderOnElement () |
Determine approximation space and order of base functions. | |
virtual MoFEMErrorCode | transformBaseFunctions () |
Transform base functions based on geometric element Jacobian. | |
![]() | |
MoFEMErrorCode | getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const |
get sense (orientation) of entity | |
MoFEMErrorCode | getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const |
Get the entity data order. | |
template<EntityType type> | |
MoFEMErrorCode | getEntitySense (EntitiesFieldData &data) const |
Get the entity sense (orientation) | |
template<EntityType type> | |
MoFEMErrorCode | getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const |
Get the entity data order for given space. | |
MoFEMErrorCode | getFaceNodes (EntitiesFieldData &data) const |
Get nodes on faces. | |
MoFEMErrorCode | getSpacesAndBaseOnEntities (EntitiesFieldData &data) const |
Get field approximation space and base on entities. | |
virtual int | getRule (int order_row, int order_col, int order_data) |
another variant of getRule | |
virtual MoFEMErrorCode | setGaussPts (int order_row, int order_col, int order_data) |
set user specific integration rule | |
MoFEMErrorCode | calHierarchicalBaseFunctionsOnElement (const FieldApproximationBase b) |
Calculate base functions. | |
MoFEMErrorCode | calHierarchicalBaseFunctionsOnElement () |
Calculate base functions. | |
MoFEMErrorCode | calBernsteinBezierBaseFunctionsOnElement () |
Calculate Bernstein-Bezier base. | |
MoFEMErrorCode | createDataOnElement (EntityType type) |
Create a entity data on element object. | |
MoFEMErrorCode | loopOverOperators () |
Iterate user data operators. | |
template<typename EXTRACTOR > | |
MoFEMErrorCode | getNodesIndices (const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const |
get node indices | |
MoFEMErrorCode | getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const |
get row node indices from FENumeredDofEntity_multiIndex | |
MoFEMErrorCode | getColNodesIndices (EntitiesFieldData &data, const int bit_number) const |
get col node indices from FENumeredDofEntity_multiIndex | |
template<typename EXTRACTOR > | |
MoFEMErrorCode | getEntityIndices (EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const |
MoFEMErrorCode | getEntityRowIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const |
MoFEMErrorCode | getEntityColIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const |
MoFEMErrorCode | getNoFieldIndices (const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const |
get NoField indices | |
MoFEMErrorCode | getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const |
get col NoField indices | |
MoFEMErrorCode | getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const |
get col NoField indices | |
MoFEMErrorCode | getBitRefLevelOnData () |
MoFEMErrorCode | getNoFieldFieldData (const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const |
Get field data on nodes. | |
MoFEMErrorCode | getNoFieldFieldData (EntitiesFieldData &data, const int bit_number) const |
MoFEMErrorCode | getNodesFieldData (EntitiesFieldData &data, const int bit_number) const |
Get data on nodes. | |
MoFEMErrorCode | getEntityFieldData (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const |
MoFEMErrorCode | getProblemNodesIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const |
get indices of nodal indices which are declared for problem but not this particular element | |
MoFEMErrorCode | getProblemTypeIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const |
get indices by type (generic function) which are declared for problem but not this particular element | |
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 |
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 |
virtual int | getRule (int order) |
virtual MoFEMErrorCode | setGaussPts (int order) |
Friends | |
class | UserDataOperator |
Additional Inherited Members | |
![]() | |
static MoFEMErrorCode | getLibVersion (Version &version) |
Get library version. | |
static MoFEMErrorCode | getFileVersion (moab::Interface &moab, Version &version) |
Get database major version. | |
static MoFEMErrorCode | setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD)) |
Get database major version. | |
static MoFEMErrorCode | getInterfaceVersion (Version &version) |
Get database major version. | |
![]() | |
static constexpr Switches | CtxSetNone = PetscData::Switches(CTX_SET_NONE) |
static constexpr Switches | CtxSetF = PetscData::Switches(CTX_SET_F) |
static constexpr Switches | CtxSetA = PetscData::Switches(CTX_SET_A) |
static constexpr Switches | CtxSetB = PetscData::Switches(CTX_SET_B) |
static constexpr Switches | CtxSetX = PetscData::Switches(CTX_SET_X) |
static constexpr Switches | CtxSetDX = PetscData::Switches(CTX_SET_DX) |
static constexpr Switches | CtxSetX_T = PetscData::Switches(CTX_SET_X_T) |
static constexpr Switches | CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT) |
static constexpr Switches | CtxSetTime = PetscData::Switches(CTX_SET_TIME) |
Volume finite element base.
User is implementing own operator at Gauss point level, by class derived from VolumeElementForcesAndSourcesCore::UserDataOperator. Arbitrary number of operator can be added by pushing objects to OpPtrVector
Definition at line 26 of file VolumeElementForcesAndSourcesCore.hpp.
Definition at line 40 of file VolumeElementForcesAndSourcesCore.hpp.
MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore | ( | Interface & | m_field | ) |
Definition at line 9 of file VolumeElementForcesAndSourcesCore.cpp.
MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore | ( | Interface & | m_field, |
const EntityType | type ) |
Definition at line 22 of file VolumeElementForcesAndSourcesCore.cpp.
|
protectedvirtual |
Calculate coordinate at integration points.
Definition at line 257 of file VolumeElementForcesAndSourcesCore.cpp.
|
protectedvirtual |
Calculate element volume and Jacobian.
Note that at that point is assumed that geometry is exclusively defined by corner nodes.
Definition at line 200 of file VolumeElementForcesAndSourcesCore.cpp.
|
protectedvirtual |
Determine approximation space and order of base functions.
Definition at line 287 of file VolumeElementForcesAndSourcesCore.cpp.
|
virtual |
function is run for every finite element
It is used to calculate element local matrices and assembly. It can be used for post-processing.
Reimplemented from MoFEM::ForcesAndSourcesCore.
Reimplemented in MoFEM::VolumeElementForcesAndSourcesCoreOnSide.
Definition at line 418 of file VolumeElementForcesAndSourcesCore.cpp.
|
protectedvirtual |
Set integration points.
Definition at line 26 of file VolumeElementForcesAndSourcesCore.cpp.
|
protectedvirtual |
Transform base functions based on geometric element Jacobian.
This function apply transformation to base functions and its derivatives. For example when base functions for H-div are present the Piola-Transformation is applied to base functions and their derivatives.
Definition at line 324 of file VolumeElementForcesAndSourcesCore.cpp.
|
friend |
Definition at line 106 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 102 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 90 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 92 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 91 of file VolumeElementForcesAndSourcesCore.hpp.
std::string MoFEM::VolumeElementForcesAndSourcesCore::meshPositionsFieldName |
Definition at line 33 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 101 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 95 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 96 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 94 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 97 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 104 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 103 of file VolumeElementForcesAndSourcesCore.hpp.
|
protected |
Definition at line 99 of file VolumeElementForcesAndSourcesCore.hpp.