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

Volume finite element base. More...

#include "src/finite_elements/VolumeElementForcesAndSourcesCore.hpp"

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

Classes

struct  UserDataOperator
 

Public Types

enum  Switches
 
- Public Types inherited from MoFEM::ForcesAndSourcesCore
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
 
- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::PetscData
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>
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 
- Public Types inherited from MoFEM::TaoMethod
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
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 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.
 
- Public Member Functions inherited from MoFEM::FEMethod
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_viewgetDataFieldEnts () 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)
 
- Public Member Functions inherited from MoFEM::BasicMethod
 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< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object.
 
- Public Member Functions inherited from MoFEM::KspMethod
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
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
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
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data.
 
- Public Member Functions inherited from MoFEM::TSMethod
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.
 
- Public Member Functions inherited from MoFEM::TaoMethod
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
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule.
 
GaussHookFun setRuleHook
 Set function to calculate integration rule.
 
MatrixDouble gaussPts
 Matrix of integration points.
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element.
 
- Public Attributes inherited from MoFEM::BasicMethod
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_multiIndexrefinedEntitiesPtr
 container of mofem dof entities
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities
 
const ProblemproblemPtr
 raw pointer to problem
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
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< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context.
 
KSP ksp
 KSP solver.
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec dx
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
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
 
- Public Attributes inherited from MoFEM::TSMethod
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.
 
- Public Attributes inherited from MoFEM::TaoMethod
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.
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
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)
 

Protected Attributes

VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
doublevOlume
 
int num_nodes
 
const EntityHandleconn
 
FTensor::Tensor2< double *, 3, 3 > tJac
 
FTensor::Tensor2< double *, 3, 3 > tInvJac
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEdataOnElement
 Entity data on element entity rows fields.
 
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields.
 
EntitiesFieldDatadataNoField
 
EntitiesFieldDatadataH1
 
EntitiesFieldDatadataHcurl
 
EntitiesFieldDatadataHdiv
 
EntitiesFieldDatadataL2
 
boost::ptr_deque< UserDataOperatoropPtrVector
 Vector of finite element users data operators.
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity.
 
MatrixDouble coordsAtGaussPts
 coordinated at gauss points
 
double elementMeasure
 

Friends

class UserDataOperator
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
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 Public Attributes inherited from MoFEM::PetscData
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)
 

Detailed Description

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

Examples
MagneticElement.hpp, Remodeling.hpp, field_evaluator.cpp, lorentz_force.cpp, mesh_smoothing.cpp, plastic.cpp, and test_cache_on_entities.cpp.

Definition at line 26 of file VolumeElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ Switches

Constructor & Destructor Documentation

◆ VolumeElementForcesAndSourcesCore() [1/2]

MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore ( Interface & m_field)

Definition at line 9 of file VolumeElementForcesAndSourcesCore.cpp.

12 meshPositionsFieldName("MESH_NODE_POSITIONS"), coords(24), jAc(3, 3),
16 tJac(&jAc(0, 0), &jAc(0, 1), &jAc(0, 2), &jAc(1, 0), &jAc(1, 1),
17 &jAc(1, 2), &jAc(2, 0), &jAc(2, 1), &jAc(2, 2)),
18 tInvJac(&invJac(0, 0), &invJac(0, 1), &invJac(0, 2), &invJac(1, 0),
19 &invJac(1, 1), &invJac(1, 2), &invJac(2, 0), &invJac(2, 1),
20 &invJac(2, 2)) {}
ForcesAndSourcesCore(Interface &m_field)

◆ VolumeElementForcesAndSourcesCore() [2/2]

MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore ( Interface & m_field,
const EntityType type )

Definition at line 22 of file VolumeElementForcesAndSourcesCore.cpp.

24 : VolumeElementForcesAndSourcesCore(m_field) {}

Member Function Documentation

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
protectedvirtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 257 of file VolumeElementForcesAndSourcesCore.cpp.

257 {
259 // Get coords at Gauss points
260 FTensor::Index<'i', 3> i;
261
262 auto &shape_functions = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
263 double *shape_functions_ptr = &*shape_functions.data().begin();
264 const size_t nb_base_functions = shape_functions.size2();
265 const size_t nb_gauss_pts = gaussPts.size2();
266 coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
267 coordsAtGaussPts.clear();
268 FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_coords_at_gauss_ptr(
269 &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
270 &coordsAtGaussPts(0, 2));
272 shape_functions_ptr);
273 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
275 &coords[0], &coords[1], &coords[2]);
276 for (int bb = 0; bb != nb_base_functions; ++bb) {
277 t_coords_at_gauss_ptr(i) += t_coords(i) * t_shape_functions;
278 ++t_coords;
279 ++t_shape_functions;
280 };
281 ++t_coords_at_gauss_ptr;
282 }
284}
@ NOBASE
Definition definitions.h:59
#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()
FTensor::Index< 'i', SPACE_DIM > i
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
MatrixDouble coordsAtGaussPts
coordinated at gauss points
MatrixDouble gaussPts
Matrix of integration points.

◆ calculateVolumeAndJacobian()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::calculateVolumeAndJacobian ( )
protectedvirtual

Calculate element volume and Jacobian.

Note that at that point is assumed that geometry is exclusively defined by corner nodes.

Returns
Error code

Definition at line 200 of file VolumeElementForcesAndSourcesCore.cpp.

200 {
202 const auto ent = numeredEntFiniteElementPtr->getEnt();
203 const auto type = numeredEntFiniteElementPtr->getEntType();
204 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
205 coords.resize(3 * num_nodes, false);
206 CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
207
208 auto get_tet_t_diff_n = [&]() {
212 };
213
214 auto get_hex_t_diff_n = [&]() {
219 };
220
221 auto get_t_diff_n = [&]() {
222 if (type == MBTET)
223 return get_tet_t_diff_n();
224 return get_hex_t_diff_n();
225 };
226
227 auto t_diff_n = get_t_diff_n();
228
230 &coords[0], &coords[1], &coords[2]);
231 FTensor::Index<'i', 3> i;
232 FTensor::Index<'j', 3> j;
233 jAc.clear();
234
235 for (size_t n = 0; n != num_nodes; ++n) {
236 tJac(i, j) += t_coords(i) * t_diff_n(j);
237 ++t_coords;
238 ++t_diff_n;
239 }
242
243 if (type == MBTET)
244 vOlume *= G_TET_W1[0] / 6.;
245
246#ifndef NDEBUG
247 if (vOlume <= std::numeric_limits<double>::epsilon() || vOlume != vOlume) {
248 MOFEM_LOG("SELF", Sev::warning) << "Volume is zero " << vOlume;
249 MOFEM_LOG("SELF", Sev::warning) << "Coords " << coords << endl;
250 }
251#endif
252
254}
#define CHKERR
Inline error check.
static const double G_TET_W1[]
Definition fem_tools.h:1115
#define MOFEM_LOG(channel, severity)
Log.
const double n
refractive index of diffusive medium
FTensor::Index< 'j', 3 > j
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
virtual moab::Interface & get_moab()=0
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static constexpr std::array< double, 24 > diffShapeFunMBHEXAtCenter
Definition Tools.hpp:218
static constexpr std::array< double, 12 > diffShapeFunMBTET
Definition Tools.hpp:271

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
protectedvirtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 287 of file VolumeElementForcesAndSourcesCore.cpp.

287 {
289
292
293 auto type = numeredEntFiniteElementPtr->getEntType();
294 auto dim_type = CN::Dimension(type);
295
296 auto get_data_on_ents = [&](auto lower_dim, auto space) {
298 auto data = dataOnElement[space];
299 data->facesNodes = dataH1.facesNodes;
300 data->facesNodesOrder = dataH1.facesNodesOrder;
301 for (auto dd = dim_type; dd >= lower_dim; --dd) {
302 int nb_ents = moab::CN::NumSubEntities(type, dd);
303 for (int ii = 0; ii != nb_ents; ++ii) {
304 auto sub_ent_type = moab::CN::SubEntityType(type, dd, ii);
305 if ((dataH1.spacesOnEntities[sub_ent_type]).test(space)) {
306 auto &data_on_ent = data->dataOnEntities[sub_ent_type];
307 CHKERR getEntitySense(sub_ent_type, data_on_ent);
308 CHKERR getEntityDataOrder(sub_ent_type, space, data_on_ent);
309 data->spacesOnEntities[sub_ent_type].set(space);
310 }
311 }
312 }
314 };
315
316 CHKERR get_data_on_ents(1, H1); // H1
317 CHKERR get_data_on_ents(1, HCURL); // Hcurl
318 CHKERR get_data_on_ents(2, HDIV); // Hdiv
319 CHKERR get_data_on_ents(3, L2); // L2
320
322}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
MatrixInt facesNodes
nodes on finite element faces
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
MatrixInt facesNodesOrder
order of face nodes on element
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
MoFEMErrorCode getFaceNodes(EntitiesFieldData &data) const
Get nodes on faces.
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.

◆ operator()()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::operator() ( )
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.

418 {
420
421 const auto type = numeredEntFiniteElementPtr->getEntType();
422 if (type != lastEvaluatedElementEntityType) {
423 switch (type) {
424 case MBTET:
426 boost::shared_ptr<BaseFunction>(new TetPolynomialBase(this));
427 break;
428 case MBHEX:
430 boost::shared_ptr<BaseFunction>(new HexPolynomialBase());
431 break;
432 default:
434 }
437 }
438
442 if (gaussPts.size2() == 0)
447
449
450 // Iterate over operators
452
454}
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
virtual MoFEMErrorCode transformBaseFunctions()
Transform base functions based on geometric element Jacobian.
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
virtual MoFEMErrorCode calculateVolumeAndJacobian()
Calculate element volume and Jacobian.
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.

◆ setIntegrationPts()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::setIntegrationPts ( )
protectedvirtual

Set integration points.

Returns
Error code

Definition at line 26 of file VolumeElementForcesAndSourcesCore.cpp.

26 {
28
29 int order_data = getMaxDataOrder();
30 int order_row = getMaxRowOrder();
31 int order_col = getMaxColOrder();
32 const auto type = numeredEntFiniteElementPtr->getEntType();
33
34 auto get_rule_by_type = [&]() {
35 switch (type) {
36 case MBHEX:
37 return getRule(order_row + 1, order_col + 1, order_data + 1);
38 default:
39 return getRule(order_row, order_col, order_data);
40 }
41 };
42
43 const int rule = get_rule_by_type();
44
45 auto calc_base_for_tet = [&]() {
47 const size_t nb_gauss_pts = gaussPts.size2();
48 auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
49 auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
50 base.resize(nb_gauss_pts, 4, false);
51 diff_base.resize(nb_gauss_pts, 12, false);
52 CHKERR Tools::shapeFunMBTET(&*base.data().begin(), &gaussPts(0, 0),
53 &gaussPts(1, 0), &gaussPts(2, 0), nb_gauss_pts);
54 double *diff_shape_ptr = &*diff_base.data().begin();
55 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
56 for (int nn = 0; nn != 4; ++nn) {
57 for (int dd = 0; dd != 3; ++dd, ++diff_shape_ptr) {
58 *diff_shape_ptr = Tools::diffShapeFunMBTET[3 * nn + dd];
59 }
60 }
61 }
63 };
64
65 auto calc_base_for_hex = [&]() {
67 const size_t nb_gauss_pts = gaussPts.size2();
68 auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
69 auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
70 base.resize(nb_gauss_pts, 8, false);
71 diff_base.resize(nb_gauss_pts, 24, false);
72 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
73 const double ksi = gaussPts(0, gg);
74 const double zeta = gaussPts(1, gg);
75 const double eta = gaussPts(2, gg);
76 base(gg, 0) = N_MBHEX0(ksi, zeta, eta);
77 base(gg, 1) = N_MBHEX1(ksi, zeta, eta);
78 base(gg, 2) = N_MBHEX2(ksi, zeta, eta);
79 base(gg, 3) = N_MBHEX3(ksi, zeta, eta);
80 base(gg, 4) = N_MBHEX4(ksi, zeta, eta);
81 base(gg, 5) = N_MBHEX5(ksi, zeta, eta);
82 base(gg, 6) = N_MBHEX6(ksi, zeta, eta);
83 base(gg, 7) = N_MBHEX7(ksi, zeta, eta);
84 diff_base(gg, 0 * 3 + 0) = diffN_MBHEX0x(zeta, eta);
85 diff_base(gg, 1 * 3 + 0) = diffN_MBHEX1x(zeta, eta);
86 diff_base(gg, 2 * 3 + 0) = diffN_MBHEX2x(zeta, eta);
87 diff_base(gg, 3 * 3 + 0) = diffN_MBHEX3x(zeta, eta);
88 diff_base(gg, 4 * 3 + 0) = diffN_MBHEX4x(zeta, eta);
89 diff_base(gg, 5 * 3 + 0) = diffN_MBHEX5x(zeta, eta);
90 diff_base(gg, 6 * 3 + 0) = diffN_MBHEX6x(zeta, eta);
91 diff_base(gg, 7 * 3 + 0) = diffN_MBHEX7x(zeta, eta);
92 diff_base(gg, 0 * 3 + 1) = diffN_MBHEX0y(ksi, eta);
93 diff_base(gg, 1 * 3 + 1) = diffN_MBHEX1y(ksi, eta);
94 diff_base(gg, 2 * 3 + 1) = diffN_MBHEX2y(ksi, eta);
95 diff_base(gg, 3 * 3 + 1) = diffN_MBHEX3y(ksi, eta);
96 diff_base(gg, 4 * 3 + 1) = diffN_MBHEX4y(ksi, eta);
97 diff_base(gg, 5 * 3 + 1) = diffN_MBHEX5y(ksi, eta);
98 diff_base(gg, 6 * 3 + 1) = diffN_MBHEX6y(ksi, eta);
99 diff_base(gg, 7 * 3 + 1) = diffN_MBHEX7y(ksi, eta);
100 diff_base(gg, 0 * 3 + 2) = diffN_MBHEX0z(ksi, zeta);
101 diff_base(gg, 1 * 3 + 2) = diffN_MBHEX1z(ksi, zeta);
102 diff_base(gg, 2 * 3 + 2) = diffN_MBHEX2z(ksi, zeta);
103 diff_base(gg, 3 * 3 + 2) = diffN_MBHEX3z(ksi, zeta);
104 diff_base(gg, 4 * 3 + 2) = diffN_MBHEX4z(ksi, zeta);
105 diff_base(gg, 5 * 3 + 2) = diffN_MBHEX5z(ksi, zeta);
106 diff_base(gg, 6 * 3 + 2) = diffN_MBHEX6z(ksi, zeta);
107 diff_base(gg, 7 * 3 + 2) = diffN_MBHEX7z(ksi, zeta);
108 }
110 };
111
112 auto set_integration_pts_for_hex = [&]() {
115 rule);
117 };
118
119 auto set_integration_pts_for_tet = [&]() {
121 if (rule < QUAD_3D_TABLE_SIZE) {
122 if (QUAD_3D_TABLE[rule]->dim != 3) {
123 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
124 }
125 if (QUAD_3D_TABLE[rule]->order < rule) {
127 "wrong order %d != %d", QUAD_3D_TABLE[rule]->order, rule);
128 }
129 size_t nb_gauss_pts = QUAD_3D_TABLE[rule]->npoints;
130 gaussPts.resize(4, nb_gauss_pts, false);
131 cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[1], 4,
132 &gaussPts(0, 0), 1);
133 cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[2], 4,
134 &gaussPts(1, 0), 1);
135 cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[3], 4,
136 &gaussPts(2, 0), 1);
137 cblas_dcopy(nb_gauss_pts, QUAD_3D_TABLE[rule]->weights, 1,
138 &gaussPts(3, 0), 1);
139
140 // CHKERR calc_base_for_tet();
141
142 auto &base = dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE);
143 auto &diff_base = dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE);
144 base.resize(nb_gauss_pts, 4, false);
145 diff_base.resize(nb_gauss_pts, 12, false);
146 double *shape_ptr = &*base.data().begin();
147 cblas_dcopy(4 * nb_gauss_pts, QUAD_3D_TABLE[rule]->points, 1, shape_ptr,
148 1);
149 double *diff_shape_ptr = &*diff_base.data().begin();
150 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
151 for (int nn = 0; nn != 4; ++nn) {
152 for (int dd = 0; dd != 3; ++dd, ++diff_shape_ptr) {
153 *diff_shape_ptr = Tools::diffShapeFunMBTET[3 * nn + dd];
154 }
155 }
156 }
157
158 } else {
159 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
160 "rule > quadrature order %d < %d", rule, QUAD_3D_TABLE_SIZE);
161 }
163 };
164
165 if (rule >= 0) {
166 switch (type) {
167 case MBTET:
168 CHKERR set_integration_pts_for_tet();
169 break;
170 case MBHEX:
171 CHKERR set_integration_pts_for_hex();
172 CHKERR calc_base_for_hex();
173 break;
174 default:
175 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
176 "Element type not implemented: %d", type);
177 }
178 } else {
179 // If rule is negative, set user defined integration points
180 CHKERR setGaussPts(order_row, order_col, order_data);
181 const size_t nb_gauss_pts = gaussPts.size2();
182 if (nb_gauss_pts) {
183 switch (type) {
184 case MBTET: {
185 CHKERR calc_base_for_tet();
186 } break;
187 case MBHEX:
188 CHKERR calc_base_for_hex();
189 break;
190 default:
191 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
192 "Element type not implemented: %d", type);
193 }
194 }
195 }
196
198}
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
constexpr int order
#define diffN_MBHEX3z(x, y)
Definition fem_tools.h:98
#define diffN_MBHEX7z(x, y)
Definition fem_tools.h:102
#define diffN_MBHEX7y(x, z)
Definition fem_tools.h:94
#define diffN_MBHEX6x(y, z)
Definition fem_tools.h:85
#define N_MBHEX7(x, y, z)
Definition fem_tools.h:78
#define N_MBHEX3(x, y, z)
Definition fem_tools.h:74
#define N_MBHEX5(x, y, z)
Definition fem_tools.h:76
#define diffN_MBHEX1y(x, z)
Definition fem_tools.h:88
#define diffN_MBHEX3y(x, z)
Definition fem_tools.h:90
#define diffN_MBHEX4y(x, z)
Definition fem_tools.h:91
#define diffN_MBHEX5x(y, z)
Definition fem_tools.h:84
#define diffN_MBHEX1z(x, y)
Definition fem_tools.h:96
#define diffN_MBHEX4z(x, y)
Definition fem_tools.h:99
#define diffN_MBHEX6y(x, z)
Definition fem_tools.h:93
#define diffN_MBHEX6z(x, y)
Definition fem_tools.h:101
#define N_MBHEX4(x, y, z)
Definition fem_tools.h:75
#define N_MBHEX0(x, y, z)
Definition fem_tools.h:71
#define N_MBHEX6(x, y, z)
Definition fem_tools.h:77
#define diffN_MBHEX5z(x, y)
Definition fem_tools.h:100
#define N_MBHEX2(x, y, z)
Definition fem_tools.h:73
#define diffN_MBHEX0x(y, z)
Definition fem_tools.h:79
#define diffN_MBHEX1x(y, z)
Definition fem_tools.h:80
#define diffN_MBHEX2z(x, y)
Definition fem_tools.h:97
#define diffN_MBHEX5y(x, z)
Definition fem_tools.h:92
#define diffN_MBHEX4x(y, z)
Definition fem_tools.h:83
#define N_MBHEX1(x, y, z)
Definition fem_tools.h:72
#define diffN_MBHEX3x(y, z)
Definition fem_tools.h:82
#define diffN_MBHEX2y(x, z)
Definition fem_tools.h:89
#define diffN_MBHEX0y(x, z)
Definition fem_tools.h:87
#define diffN_MBHEX0z(x, y)
Definition fem_tools.h:95
#define diffN_MBHEX2x(y, z)
Definition fem_tools.h:81
#define diffN_MBHEX7x(y, z)
Definition fem_tools.h:86
double eta
double zeta
Viscous hardening.
Definition plastic.cpp:130
#define QUAD_3D_TABLE_SIZE
Definition quad.h:186
static QUAD *const QUAD_3D_TABLE[]
Definition quad.h:187
virtual MPI_Comm & get_comm() const =0
int getMaxRowOrder() const
Get max order of approximation for field in rows.
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.
static MoFEMErrorCode outerProductOfEdgeIntegrationPtsForHex(MatrixDouble &pts, const int edge0, const int edge1, const int edge2)
Definition Tools.cpp:662
static MoFEMErrorCode shapeFunMBTET(double *shape, const double *ksi, const double *eta, const double *zeta, const double nb)
Calculate shape functions on tetrahedron.
Definition Tools.hpp:747
int order
Definition quad.h:28
int npoints
Definition quad.h:29

◆ transformBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCore::transformBaseFunctions ( )
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.

Returns
Error code

Definition at line 324 of file VolumeElementForcesAndSourcesCore.cpp.

324 {
326
327 if (numeredEntFiniteElementPtr->getEntType() == MBTET) {
328 // OK that code is not nice.
329 MatrixDouble new_diff_n;
330 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
331 FTensor::Index<'i', 3> i;
332 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
333 EntitiesFieldData::EntData &data = dataH1.dataOnEntities[MBVERTEX][0];
334 if ((data.getDiffN(base).size1() == 4) &&
335 (data.getDiffN(base).size2() == 3)) {
336 const size_t nb_gauss_pts = gaussPts.size2();
337 const size_t nb_base_functions = 4;
338 new_diff_n.resize(nb_gauss_pts, 3 * nb_base_functions, false);
339 double *new_diff_n_ptr = &*new_diff_n.data().begin();
341 new_diff_n_ptr, &new_diff_n_ptr[1], &new_diff_n_ptr[2]);
342 double *t_diff_n_ptr = &*data.getDiffN(base).data().begin();
343 for (unsigned int gg = 0; gg != nb_gauss_pts; gg++) {
345 t_diff_n_ptr, &t_diff_n_ptr[1], &t_diff_n_ptr[2]);
346 for (unsigned int bb = 0; bb != nb_base_functions; bb++) {
347 t_new_diff_n(i) = t_diff_n(i);
348 ++t_new_diff_n;
349 ++t_diff_n;
350 }
351 }
352 data.getDiffN(base).resize(new_diff_n.size1(), new_diff_n.size2(),
353 false);
354 data.getDiffN(base).swap(new_diff_n);
355 }
356 }
357 }
358
359 if (dataH1.spacesOnEntities[MBVERTEX].test(H1))
361
362 std::array<std::bitset<LASTSPACE>, 4> spaces_by_dim{
363
364 std::bitset<LASTSPACE>(0), //
365 std::bitset<LASTSPACE>(0), //
366 std::bitset<LASTSPACE>(0), //
367 std::bitset<LASTSPACE>(0)
368
369 };
370 for (auto type = MBEDGE; type != MBENTITYSET; ++type) {
371 spaces_by_dim[CN::Dimension(type)] |= dataH1.spacesOnEntities[type];
372 }
373
374 if (
375
376 spaces_by_dim[1].test(HCURL) || //
377 spaces_by_dim[2].test(HCURL) || //
378 spaces_by_dim[3].test(HCURL)
379
380 ) {
383 }
384
385 if (spaces_by_dim[2].test(HDIV) || spaces_by_dim[3].test(HDIV)) {
387 // Fix for tetrahedrons
388 if (numeredEntFiniteElementPtr->getEntType() == MBTET) {
389 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
390 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
391 for (auto t : {MBTRI, MBTET}) {
392 for (auto &d : dataHdiv.dataOnEntities[t]) {
393 d.getN(base) /= 6;
394 d.getDiffN(base) /= 6;
395 }
396 }
397 }
398 }
400 }
401
402 if (spaces_by_dim[3].test(L2)) {
404 }
405
407}
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
constexpr double t
plate stiffness
Definition plate.cpp:58
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)

Friends And Related Symbol Documentation

◆ UserDataOperator

Member Data Documentation

◆ conn

const EntityHandle* MoFEM::VolumeElementForcesAndSourcesCore::conn
protected

Definition at line 102 of file VolumeElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::VolumeElementForcesAndSourcesCore::coords
protected

Definition at line 90 of file VolumeElementForcesAndSourcesCore.hpp.

◆ invJac

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCore::invJac
protected

Definition at line 92 of file VolumeElementForcesAndSourcesCore.hpp.

◆ jAc

MatrixDouble3by3 MoFEM::VolumeElementForcesAndSourcesCore::jAc
protected

Definition at line 91 of file VolumeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::VolumeElementForcesAndSourcesCore::meshPositionsFieldName
Deprecated
DO NOT USE!

Definition at line 33 of file VolumeElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::VolumeElementForcesAndSourcesCore::num_nodes
protected

Definition at line 101 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opContravariantPiolaTransform

OpSetContravariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opContravariantPiolaTransform
protected

Definition at line 95 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opCovariantPiolaTransform

OpSetCovariantPiolaTransform MoFEM::VolumeElementForcesAndSourcesCore::opCovariantPiolaTransform
protected

Definition at line 96 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacH1

OpSetInvJacH1 MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacH1
protected

Definition at line 94 of file VolumeElementForcesAndSourcesCore.hpp.

◆ opSetInvJacHdivAndHcurl

OpSetInvJacHdivAndHcurl MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacHdivAndHcurl
protected

Definition at line 97 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tInvJac

FTensor::Tensor2<double *, 3, 3> MoFEM::VolumeElementForcesAndSourcesCore::tInvJac
protected

Definition at line 104 of file VolumeElementForcesAndSourcesCore.hpp.

◆ tJac

FTensor::Tensor2<double *, 3, 3> MoFEM::VolumeElementForcesAndSourcesCore::tJac
protected

Definition at line 103 of file VolumeElementForcesAndSourcesCore.hpp.

◆ vOlume

double& MoFEM::VolumeElementForcesAndSourcesCore::vOlume
protected

Definition at line 99 of file VolumeElementForcesAndSourcesCore.hpp.


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