v0.13.1
Classes | Public Member Functions | Protected Attributes | Friends | List of all members
MoFEM::FatPrismElementForcesAndSourcesCore Struct Reference

FatPrism finite element. More...

#include <src/finite_elements/FatPrismElementForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 default operator for Flat Prism element More...
 

Public Member Functions

 FatPrismElementForcesAndSourcesCore (Interface &m_field)
 
virtual int getRuleTrianglesOnly (int order)
 
virtual int getRuleThroughThickness (int order)
 
virtual MoFEMErrorCode setGaussPtsTrianglesOnly (int order_triangles_only)
 
virtual MoFEMErrorCode setGaussPtsThroughThickness (int order_thickness)
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore
 VolumeElementForcesAndSourcesCore (Interface &m_field, const EntityType type=MBTET)
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
boost::ptr_vector< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator. More...
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object. More...
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
int getMaxDataOrder () const
 Get max order of approximation for data fields. More...
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows. More...
 
int getMaxColOrder () const
 Get max order of approximation for field in columns. More...
 
auto & getEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
auto & getDataOnElementBySpaceArray ()
 Get data on entities and space. More...
 
auto & getDerivedDataOnElementBySpaceArray ()
 Get derived data on entities and space. More...
 
- 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 More...
 
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
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- 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 More...
 
- Public Member Functions inherited from MoFEM::PetscData
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
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. More...
 
- 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. More...
 

Protected Attributes

double aRea [2]
 
VectorDouble normal
 
MatrixDouble gaussPtsTrianglesOnly
 
MatrixDouble coordsAtGaussPtsTrianglesOnly
 
MatrixDouble gaussPtsThroughThickness
 
EntitiesFieldData dataH1TrianglesOnly
 
EntitiesFieldData dataH1TroughThickness
 
MatrixDouble hoCoordsAtGaussPtsF3
 
MatrixDouble nOrmals_at_GaussPtF3
 
MatrixDouble tAngent1_at_GaussPtF3
 
MatrixDouble tAngent2_at_GaussPtF3
 
MatrixDouble hoCoordsAtGaussPtsF4
 
MatrixDouble nOrmals_at_GaussPtF4
 
MatrixDouble tAngent1_at_GaussPtF4
 
MatrixDouble tAngent2_at_GaussPtF4
 
OpGetCoordsAndNormalsOnPrism opHOCoordsAndNormals
 
- Protected Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCore
VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
doublevOlume
 
int num_nodes
 
const EntityHandle * conn
 
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. More...
 
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
EntitiesFieldDatadataNoField
 
EntitiesFieldDatadataH1
 
EntitiesFieldDatadataHcurl
 
EntitiesFieldDatadataHdiv
 
EntitiesFieldDatadataL2
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 
MatrixDouble coordsAtGaussPts
 coordinated at gauss points More...
 
double elementMeasure
 

Friends

class UserDataOperator
 

Additional Inherited Members

- Public Types inherited from MoFEM::VolumeElementForcesAndSourcesCore
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_X_T = 1 << 4 , 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
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Public Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCore
std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
GaussHookFun setRuleHook
 Set function to calculate integration rule. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
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 x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 
- 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 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)
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
virtual MoFEMErrorCode calculateVolumeAndJacobian ()
 Calculate element volume and Jacobian. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode transformBaseFunctions ()
 Transform base functions based on geometric element Jacobian. More...
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (EntitiesFieldData &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceNodes (EntitiesFieldData &data) const
 Get nodes on faces. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (EntitiesFieldData &data) const
 Get field approximation space and base on entities. More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement ()
 Calculate base functions. More...
 
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement ()
 Calculate Bernstein-Bezier base. More...
 
MoFEMErrorCode createDataOnElement (EntityType type)
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
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 More...
 
MoFEMErrorCode getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
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 More...
 
MoFEMErrorCode getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
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. More...
 
MoFEMErrorCode getNoFieldFieldData (EntitiesFieldData &data, const int bit_number) const
 
MoFEMErrorCode getNodesFieldData (EntitiesFieldData &data, const int bit_number) const
 Get data on nodes. More...
 
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 More...
 
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 More...
 
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)
 

Detailed Description

FatPrism finite element.

User is implementing own operator at Gauss points level, by own object derived from FatPrismElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to rowOpPtrVector and rowColOpPtrVector.

Todo:
Need to implement operators that will make this element work as Volume element
Examples
HookeElement.cpp, elasticity.cpp, prism_elements_from_surface.cpp, and prism_polynomial_approximation.cpp.

Definition at line 43 of file FatPrismElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FatPrismElementForcesAndSourcesCore()

MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore ( Interface m_field)
Examples
elasticity.cpp, and prism_polynomial_approximation.cpp.

Definition at line 23 of file FatPrismElementForcesAndSourcesCore.cpp.

25 : VolumeElementForcesAndSourcesCore(m_field, MBPRISM),
32 "Problem with creation data on element");
33}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:608
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
VolumeElementForcesAndSourcesCore(Interface &m_field, const EntityType type=MBTET)

Member Function Documentation

◆ getRuleThroughThickness()

virtual int MoFEM::FatPrismElementForcesAndSourcesCore::getRuleThroughThickness ( int  order)
virtual

◆ getRuleTrianglesOnly()

virtual int MoFEM::FatPrismElementForcesAndSourcesCore::getRuleTrianglesOnly ( int  order)
virtual

◆ operator()()

MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::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.

Definition at line 35 of file FatPrismElementForcesAndSourcesCore.cpp.

35 {
37
38 if (numeredEntFiniteElementPtr->getEntType() != MBPRISM)
40
41 auto get_fe_coordinates = [&]() {
43 EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
44 int num_nodes;
45 const EntityHandle *conn;
46 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
47 coords.resize(num_nodes * 3, false);
48 CHKERR mField.get_moab().get_coords(conn, num_nodes,
49 &*coords.data().begin());
51 };
52
53 auto calculate_area_of_triangles = [&] {
55 normal.resize(6, false);
58 aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
59 aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
61 };
62
63 CHKERR get_fe_coordinates();
64 CHKERR calculate_area_of_triangles();
65
69
70 auto get_h1_base_data = [&](auto &dataH1) {
72 CHKERR getEntitySense<MBEDGE>(dataH1);
73 CHKERR getEntitySense<MBTRI>(dataH1);
74 CHKERR getEntitySense<MBQUAD>(dataH1);
75 CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
76 CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
77 CHKERR getEntityDataOrder<MBQUAD>(dataH1, H1);
78 CHKERR getEntityDataOrder<MBPRISM>(dataH1, H1);
79 // Triangles only
80 CHKERR getEntitySense<MBEDGE>(dataH1TrianglesOnly);
81 CHKERR getEntitySense<MBTRI>(dataH1TrianglesOnly);
82 CHKERR getEntityDataOrder<MBEDGE>(dataH1TrianglesOnly, H1);
83 CHKERR getEntityDataOrder<MBTRI>(dataH1TrianglesOnly, H1);
84 // Through thickness
85 CHKERR getEntitySense<MBEDGE>(dataH1TroughThickness);
86 CHKERR getEntityDataOrder<MBEDGE>(dataH1TroughThickness, H1);
88 };
89
90 // H1
91 if ((dataH1.spacesOnEntities[MBEDGE]).test(H1))
92 CHKERR get_h1_base_data(dataH1);
93
94 // Hdiv
95 if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV))
96 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
97
98 // Hcurl
99 if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL))
100 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
101
102 // L2
103 if ((dataH1.spacesOnEntities[MBPRISM]).test(L2))
104 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
105
106 // get approx. on triangles, i.e. faces 3 and 4
107 auto set_gauss_points_on_triangle = [&](int &nb_gauss_pts_on_faces) {
109 int order_triangles_only = 1;
110 int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
111 for (unsigned int ee = 0; ee < 9; ee++) {
112 if (!valid_edges[ee])
113 continue;
114 order_triangles_only = std::max(
115 order_triangles_only,
116 dataH1TrianglesOnly.dataOnEntities[MBEDGE][ee].getOrder());
117 }
118 for (unsigned int ff = 3; ff <= 4; ff++) {
119 order_triangles_only = std::max(
120 order_triangles_only,
121 dataH1TrianglesOnly.dataOnEntities[MBTRI][ff].getOrder());
122 }
123 for (unsigned int qq = 0; qq < 3; qq++) {
124 order_triangles_only = std::max(
125 order_triangles_only,
126 dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getOrder());
127 }
128 order_triangles_only = std::max(
129 order_triangles_only,
130 dataH1TroughThickness.dataOnEntities[MBPRISM][0].getOrder());
131
132 // integration pts on the triangles surfaces
133 nb_gauss_pts_on_faces = 0;
134 int rule = getRuleTrianglesOnly(order_triangles_only);
135 if (rule >= 0) {
136 if (rule < QUAD_2D_TABLE_SIZE) {
137 if (QUAD_2D_TABLE[rule]->dim != 2) {
138 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
139 }
140 if (QUAD_2D_TABLE[rule]->order < rule) {
141 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
142 "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
143 }
144 nb_gauss_pts_on_faces = QUAD_2D_TABLE[rule]->npoints;
145 gaussPtsTrianglesOnly.resize(3, nb_gauss_pts_on_faces, false);
146 cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[1], 3,
147 &gaussPtsTrianglesOnly(0, 0), 1);
148 cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[2], 3,
149 &gaussPtsTrianglesOnly(1, 0), 1);
150 cblas_dcopy(nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->weights, 1,
151 &gaussPtsTrianglesOnly(2, 0), 1);
152 dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
153 nb_gauss_pts_on_faces, 3, false);
154 double *shape_ptr = &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
155 .getN(NOBASE)
156 .data()
157 .begin();
158 cblas_dcopy(3 * nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->points, 1,
159 shape_ptr, 1);
160 dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
161 1, 6, false);
162 std::copy(Tools::diffShapeFunMBTRI.begin(),
165 .getDiffN(NOBASE)
166 .data()
167 .begin());
168 } else
169 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
170 "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
171
172 } else {
173 CHKERR setGaussPtsTrianglesOnly(order_triangles_only);
174 nb_gauss_pts_on_faces = gaussPtsTrianglesOnly.size2();
175 if (nb_gauss_pts_on_faces == 0)
177 dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
178 nb_gauss_pts_on_faces, 3, false);
179 if (nb_gauss_pts_on_faces) {
181 .getN(NOBASE)
182 .data()
183 .begin(),
185 &gaussPtsTrianglesOnly(1, 0), nb_gauss_pts_on_faces);
186 dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
187 1, 6, false);
188 std::copy(Tools::diffShapeFunMBTRI.begin(),
191 .getDiffN(NOBASE)
192 .data()
193 .begin());
194 }
195 }
197 };
198
199 // approx. trough prism thickness
200 auto set_gauss_points_through_thickness =
201 [&](int &nb_gauss_pts_through_thickness) {
203 nb_gauss_pts_through_thickness = 0;
204 int order_thickness = 1;
205 for (unsigned int ee = 3; ee <= 5; ee++) {
206 order_thickness = std::max(
207 order_thickness,
208 dataH1TroughThickness.dataOnEntities[MBEDGE][ee].getOrder());
209 }
210 for (unsigned int qq = 0; qq < 3; qq++) {
211 order_thickness = std::max(
212 order_thickness,
213 dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getOrder());
214 }
215 order_thickness = std::max(
216 order_thickness,
217 dataH1TroughThickness.dataOnEntities[MBPRISM][0].getOrder());
218 // integration points
219 int rule = getRuleThroughThickness(order_thickness);
220 if (rule >= 0) {
221 if (rule < QUAD_1D_TABLE_SIZE) {
222 if (QUAD_1D_TABLE[rule]->dim != 1) {
223 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
224 "wrong dimension");
225 }
226 if (QUAD_1D_TABLE[rule]->order < rule) {
227 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
228 "wrong order %d != %d", QUAD_1D_TABLE[rule]->order,
229 rule);
230 }
231 nb_gauss_pts_through_thickness = QUAD_1D_TABLE[rule]->npoints;
232 gaussPtsThroughThickness.resize(2, nb_gauss_pts_through_thickness,
233 false);
234 cblas_dcopy(nb_gauss_pts_through_thickness,
235 &QUAD_1D_TABLE[rule]->points[1], 2,
236 &gaussPtsThroughThickness(0, 0), 1);
237 cblas_dcopy(nb_gauss_pts_through_thickness,
238 QUAD_1D_TABLE[rule]->weights, 1,
239 &gaussPtsThroughThickness(1, 0), 1);
240 } else {
241 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
242 "rule > quadrature order %d < %d", rule,
244 nb_gauss_pts_through_thickness = 0;
245 }
246 } else {
247 CHKERR setGaussPtsThroughThickness(order_thickness);
248 nb_gauss_pts_through_thickness = gaussPtsThroughThickness.size2();
249 }
251 };
252
253 // Generate integration pts.
254 auto set_gauss_points_in_volume = [&](int nb_gauss_pts_on_faces,
255 int nb_gauss_pts_through_thickness,
256 int &nb_gauss_pts) {
258 nb_gauss_pts = nb_gauss_pts_on_faces * nb_gauss_pts_through_thickness;
259 gaussPts.resize(4, nb_gauss_pts, false);
260 int gg = 0;
261 for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
262 for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
263 gaussPts(0, gg) = gaussPtsTrianglesOnly(0, ggf);
264 gaussPts(1, gg) = gaussPtsTrianglesOnly(1, ggf);
265 gaussPts(2, gg) = gaussPtsThroughThickness(0, ggt);
266 gaussPts(3, gg) =
268 }
269 }
271 };
272
273 int nb_gauss_pts, nb_gauss_pts_through_thickness, nb_gauss_pts_on_faces;
274 CHKERR set_gauss_points_on_triangle(nb_gauss_pts_on_faces);
275 if (!nb_gauss_pts_on_faces)
277 CHKERR set_gauss_points_through_thickness(nb_gauss_pts_through_thickness);
278 if (!nb_gauss_pts_through_thickness)
280 CHKERR set_gauss_points_in_volume(
281 nb_gauss_pts_on_faces, nb_gauss_pts_through_thickness, nb_gauss_pts);
282
283 auto calc_coordinates_at_triangles = [&]() {
285 coordsAtGaussPtsTrianglesOnly.resize(nb_gauss_pts_on_faces, 6, false);
286 for (int gg = 0; gg < nb_gauss_pts_on_faces; gg++) {
287 for (int dd = 0; dd < 3; dd++) {
289 cblas_ddot(3,
290 &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
291 NOBASE)(gg, 0),
292 1, &coords[dd], 3);
294 cblas_ddot(3,
295 &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
296 NOBASE)(gg, 0),
297 1, &coords[9 + dd], 3);
298 }
299 }
301 };
302
303 auto calc_vertex_base_on_prism = [&]() {
305 // Calculate "nobase" base functions on prism, this is cartesian product
306 // of base functions on triangles with base functions through thickness
307 // FIXME: This could be effectively implemented with tensors
308 dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 6,
309 false);
310 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(nb_gauss_pts, 18,
311 false);
312 for (int dd = 0; dd != 6; dd++) {
313 int gg = 0;
314 for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
315 int ddd = dd > 2 ? dd - 3 : dd;
316 double tri_n = dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
317 NOBASE)(ggf, ddd);
318 double dksi_tri_n =
319 dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
320 0, 2 * ddd + 0);
321 double deta_tri_n =
322 dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
323 0, 2 * ddd + 1);
324 for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
325 double zeta = gaussPtsThroughThickness(0, ggt);
326 double dzeta, edge_shape;
327 if (dd < 3) {
328 dzeta = diffN_MBEDGE0;
329 edge_shape = N_MBEDGE0(zeta);
330 } else {
331 dzeta = diffN_MBEDGE1;
332 edge_shape = N_MBEDGE1(zeta);
333 }
334 dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, dd) =
335 tri_n * edge_shape;
336 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 0) =
337 dksi_tri_n * edge_shape;
338 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 1) =
339 deta_tri_n * edge_shape;
340 dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 2) =
341 tri_n * dzeta;
342 }
343 }
344 }
346 };
347
348 auto calc_base_on_prism = [&]() {
350 // Calculate base functions on prism
351 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
352 if (dataH1.bAse.test(b)) {
353 switch (static_cast<FieldApproximationBase>(b)) {
356 if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
357 CHKERR FatPrismPolynomialBase().getValue(
358 gaussPts,
359 boost::shared_ptr<BaseFunctionCtx>(
360 new FatPrismPolynomialBaseCtx(
364 static_cast<FieldApproximationBase>(b), NOBASE)));
365 }
366 if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
367 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
368 "Not yet implemented");
369 }
370 if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
371 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
372 "Not yet implemented");
373 }
374 if (dataH1.spacesOnEntities[MBTET].test(L2)) {
375 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
376 "Not yet implemented");
377 }
378 break;
379 default:
380 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
381 "Not yet implemented");
382 }
383 }
384 }
386 };
387
388 auto calc_coordinate_on_prism = [&]() {
390 /// Calculate coordinates at integration points
391 coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
392 for (int gg = 0; gg < nb_gauss_pts; gg++) {
393 for (int dd = 0; dd < 3; dd++) {
394 coordsAtGaussPts(gg, dd) = cblas_ddot(
395 6, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
396 &coords[dd], 3);
397 }
398 }
400 };
401
402 auto calculate_volume = [&]() {
403 auto get_t_w = [&] {
405 };
406
407 auto get_t_coords = [&]() {
409 &coords[0], &coords[1], &coords[2]);
410 };
411
415
416 const size_t nb_gauss_pts = gaussPts.size2();
417 auto t_diff_n =
418 dataH1.dataOnEntities[MBVERTEX][0].getFTensor1DiffN<3>(NOBASE);
419
420 double vol = 0;
421 auto t_w = get_t_w();
422 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
423
424 auto t_coords = get_t_coords();
425 t_jac(i, j) = 0;
426 for (size_t n = 0; n != 6; ++n) {
427 t_jac(i, j) += t_coords(i) * t_diff_n(j);
428 ++t_diff_n;
429 ++t_coords;
430 }
431
432 double det;
433 CHKERR determinantTensor3by3(t_jac, det);
434 vol += det * t_w / 2;
435
436 ++t_w;
437 }
438
439 return vol;
440 };
441
442 auto calc_ho_triangle_face_normals = [&]() {
444
445 auto check_field = [&]() {
446 auto field_it =
447 fieldsPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName);
448 if (field_it != fieldsPtr->get<FieldName_mi_tag>().end())
449 if ((numeredEntFiniteElementPtr->getBitFieldIdData() &
450 (*field_it)->getId())
451 .any())
452 return true;
453 return false;
454 };
455
456 // Check if field meshPositionsFieldName exist
457 if (check_field()) {
458 hoCoordsAtGaussPtsF3.resize(nb_gauss_pts_on_faces, 3, false);
459 nOrmals_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
460 tAngent1_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
461 tAngent2_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
462 hoCoordsAtGaussPtsF4.resize(nb_gauss_pts_on_faces, 3, false);
463 nOrmals_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
464 tAngent1_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
465 tAngent2_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
466 const auto bit_number =
469 CHKERR getEntityFieldData(dataH1TrianglesOnly, bit_number, MBEDGE);
470 CHKERR getEntityFieldData(dataH1TrianglesOnly, bit_number, MBEDGE);
473 } else {
474 hoCoordsAtGaussPtsF3.resize(0, 0, false);
475 nOrmals_at_GaussPtF3.resize(0, 0, false);
476 tAngent1_at_GaussPtF3.resize(0, 0, false);
477 tAngent2_at_GaussPtF3.resize(0, 0, false);
478 hoCoordsAtGaussPtsF4.resize(0, 0, false);
479 nOrmals_at_GaussPtF4.resize(0, 0, false);
480 tAngent1_at_GaussPtF4.resize(0, 0, false);
481 tAngent2_at_GaussPtF4.resize(0, 0, false);
482 }
484 };
485
486 CHKERR calc_coordinates_at_triangles();
487 CHKERR calc_vertex_base_on_prism();
488 CHKERR calc_base_on_prism();
489 CHKERR calc_coordinate_on_prism();
490 vOlume = calculate_volume();
491 CHKERR calc_ho_triangle_face_normals();
492
493 // Iterate over operators
495
497}
static Index< 'n', 3 > n
FieldApproximationBase
approximation base
Definition: definitions.h:71
@ LASTBASE
Definition: definitions.h:82
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:73
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:75
@ NOBASE
Definition: definitions.h:72
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
@ L2
field with C-1 continuity
Definition: definitions.h:101
@ H1
continuous field
Definition: definitions.h:98
@ HCURL
field with continuous tangents
Definition: definitions.h:99
@ HDIV
field with continuous normal traction
Definition: definitions.h:100
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:44
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:45
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
const int dim
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:117
#define diffN_MBEDGE1
derivative of edge shape function
Definition: fem_tools.h:118
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:115
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:116
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:194
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
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
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1207
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163
const Field_multiIndex * fieldsPtr
raw pointer to fields container
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
std::bitset< LASTBASE > bAse
bases on element
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
virtual MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MatrixDouble coordsAtGaussPts
coordinated at gauss points
MatrixDouble gaussPts
Matrix of integration points.
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:113
static MoFEMErrorCode getTriNormal(const double *coords, double *normal)
Get the Tri Normal objectGet triangle normal.
Definition: Tools.cpp:363
int npoints
Definition: quad.h:29

◆ setGaussPtsThroughThickness()

virtual MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::setGaussPtsThroughThickness ( int  order_thickness)
virtual

Reimplemented in PrismFE, PostProcFatPrismOnRefinedMesh, and SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh.

Definition at line 57 of file FatPrismElementForcesAndSourcesCore.hpp.

57 {
59 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
61 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453

◆ setGaussPtsTrianglesOnly()

virtual MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::setGaussPtsTrianglesOnly ( int  order_triangles_only)
virtual

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

Definition at line 218 of file FatPrismElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

double MoFEM::FatPrismElementForcesAndSourcesCore::aRea[2]
protected

Definition at line 198 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly
protected

Definition at line 202 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TrianglesOnly

EntitiesFieldData MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly
protected

Definition at line 205 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TroughThickness

EntitiesFieldData MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness
protected

Definition at line 206 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsThroughThickness

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
protected

◆ gaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
protected

◆ hoCoordsAtGaussPtsF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3
protected

Definition at line 208 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4
protected

Definition at line 212 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ normal

VectorDouble MoFEM::FatPrismElementForcesAndSourcesCore::normal
protected

Definition at line 199 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3
protected

Definition at line 209 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4
protected

Definition at line 213 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnPrism MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals
protected

Definition at line 216 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3
protected

Definition at line 210 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4
protected

Definition at line 214 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3
protected

Definition at line 211 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4
protected

Definition at line 215 of file FatPrismElementForcesAndSourcesCore.hpp.


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