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

FatPrism finite elementUser 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. 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::VolumeElementForcesAndSourcesCoreSwitch< SWITCH >
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase
 VolumeElementForcesAndSourcesCoreBase (Interface &m_field, const EntityType type=MBTET)
 
template<int SWITCH>
MoFEMErrorCode OpSwitch ()
 
- 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 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...
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. More...
 
const DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side) const
 Get the entity data. More...
 
DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of 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
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Protected Attributes

double aRea [2]
 
VectorDouble normal
 
MatrixDouble gaussPtsTrianglesOnly
 
MatrixDouble coordsAtGaussPtsTrianglesOnly
 
MatrixDouble gaussPtsThroughThickness
 
DataForcesAndSourcesCore dataH1TrianglesOnly
 
DataForcesAndSourcesCore 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::VolumeElementForcesAndSourcesCoreBase
VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble hoGaussPtsJac
 
MatrixDouble hoGaussPtsInvJac
 
VectorDouble hoGaussPtsDetJac
 
OpGetDataAndGradient< 3, 3 > opHOatGaussPoints
 higher order geometry data at Gauss pts More...
 
OpSetHoInvJacH1 opSetHoInvJacH1
 
OpSetHoContravariantPiolaTransform opHoContravariantTransform
 
OpSetHoCovariantPiolaTransform opHoCovariantTransform
 
OpSetHoInvJacHdivAndHcurl opSetHoInvJacHdivAndHcurl
 
double vOlume
 
int num_nodes
 
const EntityHandleconn
 
FTensor::Tensor2< double *, 3, 3 > tJac
 
FTensor::Tensor2< double *, 3, 3 > tInvJac
 
MatrixDouble coordsAtGaussPts
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 

Friends

class UserDataOperator
 

Additional Inherited Members

- Public Types inherited from MoFEM::VolumeElementForcesAndSourcesCoreSwitch< SWITCH >
using UserDataOperator = VolumeElementForcesAndSourcesCoreBase::UserDataOperator
 
- Public Types inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase
enum  Switches { NO_HO_GEOMETRY = 1 << 0 | 1 << 2, NO_TRANSFORM = 1 << 1 | 1 << 2, NO_HO_TRANSFORM = 1 << 2 }
 
- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
- 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::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 Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase
std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
RuleHookFun 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::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process 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
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- 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
TSContext ts_ctx
 
TS ts
 time solver 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...
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_tt (see PETSc Time Solver) More...
 
PetscReal ts_v
 shift for U_t shift for U_t More...
 
PetscReal ts_t
 time More...
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase
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...
 
virtual MoFEMErrorCode calculateHoJacobian ()
 Calculate Jacobian for HO geometry. More...
 
virtual MoFEMErrorCode transformHoBaseFunctions ()
 Transform base functions based on ho-geometry element Jacobian. More...
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (DataForcesAndSourcesCore &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &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 ()
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
MoFEMErrorCode getNodesIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEntityIndices (DataForcesAndSourcesCore &data, const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (DataForcesAndSourcesCore &data, const boost::string_ref field_name) const
 
MoFEMErrorCode getNodesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (DataForcesAndSourcesCore &data, const std::string &field_name, 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)
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

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
elasticity.cpp, HookeElement.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),
26  dataH1TrianglesOnly(MBPRISM), dataH1TroughThickness(MBPRISM),
VolumeElementForcesAndSourcesCoreSwitch< 0 > VolumeElementForcesAndSourcesCore
Volume finite element default.

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 32 of file FatPrismElementForcesAndSourcesCore.cpp.

32  {
34 
35  if (numeredEntFiniteElementPtr->getEntType() != MBPRISM)
38 
39  auto get_fe_coordinates = [&]() {
42  int num_nodes;
43  const EntityHandle *conn;
44  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
45  coords.resize(num_nodes * 3, false);
46  CHKERR mField.get_moab().get_coords(conn, num_nodes,
47  &*coords.data().begin());
49  };
50 
51  auto calculate_area_of_triangles = [&] {
53  normal.resize(6, false);
56  aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
57  aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
59  };
60 
61  CHKERR get_fe_coordinates();
62  CHKERR calculate_area_of_triangles();
63 
67 
68  auto get_h1_base_data = [&](auto &dataH1) {
70  CHKERR getEntitySense<MBEDGE>(dataH1);
71  CHKERR getEntitySense<MBTRI>(dataH1);
72  CHKERR getEntitySense<MBQUAD>(dataH1);
73  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
74  CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
75  CHKERR getEntityDataOrder<MBQUAD>(dataH1, H1);
76  CHKERR getEntityDataOrder<MBPRISM>(dataH1, H1);
77  // Triangles only
78  CHKERR getEntitySense<MBEDGE>(dataH1TrianglesOnly);
79  CHKERR getEntitySense<MBTRI>(dataH1TrianglesOnly);
80  CHKERR getEntityDataOrder<MBEDGE>(dataH1TrianglesOnly, H1);
81  CHKERR getEntityDataOrder<MBTRI>(dataH1TrianglesOnly, H1);
82  // Through thickness
83  CHKERR getEntitySense<MBEDGE>(dataH1TroughThickness);
84  CHKERR getEntityDataOrder<MBEDGE>(dataH1TroughThickness, H1);
86  };
87 
88  // H1
89  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1))
90  CHKERR get_h1_base_data(dataH1);
91 
92  // Hdiv
93  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV))
94  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
95 
96  // Hcurl
97  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL))
98  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
99 
100  // L2
101  if ((dataH1.spacesOnEntities[MBPRISM]).test(L2))
102  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
103 
104  // get approx. on triangles, i.e. faces 3 and 4
105  auto set_gauss_points_on_triangle = [&](int &nb_gauss_pts_on_faces) {
107  int order_triangles_only = 1;
108  int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
109  for (unsigned int ee = 0; ee < 9; ee++) {
110  if (!valid_edges[ee])
111  continue;
112  order_triangles_only = std::max(
113  order_triangles_only,
114  dataH1TrianglesOnly.dataOnEntities[MBEDGE][ee].getDataOrder());
115  }
116  for (unsigned int ff = 3; ff <= 4; ff++) {
117  order_triangles_only = std::max(
118  order_triangles_only,
119  dataH1TrianglesOnly.dataOnEntities[MBTRI][ff].getDataOrder());
120  }
121  for (unsigned int qq = 0; qq < 3; qq++) {
122  order_triangles_only = std::max(
123  order_triangles_only,
124  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
125  }
126  order_triangles_only = std::max(
127  order_triangles_only,
128  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
129 
130  // integration pts on the triangles surfaces
131  nb_gauss_pts_on_faces = 0;
132  int rule = getRuleTrianglesOnly(order_triangles_only);
133  if (rule >= 0) {
134  if (rule < QUAD_2D_TABLE_SIZE) {
135  if (QUAD_2D_TABLE[rule]->dim != 2) {
136  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
137  }
138  if (QUAD_2D_TABLE[rule]->order < rule) {
139  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
140  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
141  }
142  nb_gauss_pts_on_faces = QUAD_2D_TABLE[rule]->npoints;
143  gaussPtsTrianglesOnly.resize(3, nb_gauss_pts_on_faces, false);
144  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[1], 3,
145  &gaussPtsTrianglesOnly(0, 0), 1);
146  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[2], 3,
147  &gaussPtsTrianglesOnly(1, 0), 1);
148  cblas_dcopy(nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->weights, 1,
149  &gaussPtsTrianglesOnly(2, 0), 1);
150  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
151  nb_gauss_pts_on_faces, 3, false);
152  double *shape_ptr = &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
153  .getN(NOBASE)
154  .data()
155  .begin();
156  cblas_dcopy(3 * nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->points, 1,
157  shape_ptr, 1);
158  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
159  1, 6, false);
160  std::copy(Tools::diffShapeFunMBTRI.begin(),
162  &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
163  .getDiffN(NOBASE)
164  .data()
165  .begin());
166  } else
167  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
168  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
169 
170  } else {
171  CHKERR setGaussPtsTrianglesOnly(order_triangles_only);
172  nb_gauss_pts_on_faces = gaussPtsTrianglesOnly.size2();
173  if (nb_gauss_pts_on_faces == 0)
175  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
176  nb_gauss_pts_on_faces, 3, false);
177  if (nb_gauss_pts_on_faces) {
179  .getN(NOBASE)
180  .data()
181  .begin(),
182  &gaussPtsTrianglesOnly(0, 0),
183  &gaussPtsTrianglesOnly(1, 0), nb_gauss_pts_on_faces);
184  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
185  1, 6, false);
186  std::copy(Tools::diffShapeFunMBTRI.begin(),
188  &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
189  .getDiffN(NOBASE)
190  .data()
191  .begin());
192  }
193  }
195  };
196 
197  // approx. trough prism thickness
198  auto set_gauss_points_through_thickness =
199  [&](int &nb_gauss_pts_through_thickness) {
201  nb_gauss_pts_through_thickness = 0;
202  int order_thickness = 1;
203  for (unsigned int ee = 3; ee <= 5; ee++) {
204  order_thickness = std::max(
205  order_thickness,
206  dataH1TroughThickness.dataOnEntities[MBEDGE][ee].getDataOrder());
207  }
208  for (unsigned int qq = 0; qq < 3; qq++) {
209  order_thickness = std::max(
210  order_thickness,
211  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
212  }
213  order_thickness = std::max(
214  order_thickness,
215  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
216  // integration points
217  int rule = getRuleThroughThickness(order_thickness);
218  if (rule >= 0) {
219  if (rule < QUAD_1D_TABLE_SIZE) {
220  if (QUAD_1D_TABLE[rule]->dim != 1) {
221  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
222  "wrong dimension");
223  }
224  if (QUAD_1D_TABLE[rule]->order < rule) {
225  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
226  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order,
227  rule);
228  }
229  nb_gauss_pts_through_thickness = QUAD_1D_TABLE[rule]->npoints;
230  gaussPtsThroughThickness.resize(2, nb_gauss_pts_through_thickness,
231  false);
232  cblas_dcopy(nb_gauss_pts_through_thickness,
233  &QUAD_1D_TABLE[rule]->points[1], 2,
234  &gaussPtsThroughThickness(0, 0), 1);
235  cblas_dcopy(nb_gauss_pts_through_thickness,
236  QUAD_1D_TABLE[rule]->weights, 1,
237  &gaussPtsThroughThickness(1, 0), 1);
238  } else {
239  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
240  "rule > quadrature order %d < %d", rule,
242  nb_gauss_pts_through_thickness = 0;
243  }
244  } else {
245  CHKERR setGaussPtsThroughThickness(order_thickness);
246  nb_gauss_pts_through_thickness = gaussPtsThroughThickness.size2();
247  }
249  };
250 
251  // Generate integration pts.
252  auto set_gauss_points_in_volume = [&](int nb_gauss_pts_on_faces,
253  int nb_gauss_pts_through_thickness,
254  int &nb_gauss_pts) {
256  nb_gauss_pts = nb_gauss_pts_on_faces * nb_gauss_pts_through_thickness;
257  gaussPts.resize(4, nb_gauss_pts, false);
258  int gg = 0;
259  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
260  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
261  gaussPts(0, gg) = gaussPtsTrianglesOnly(0, ggf);
262  gaussPts(1, gg) = gaussPtsTrianglesOnly(1, ggf);
263  gaussPts(2, gg) = gaussPtsThroughThickness(0, ggt);
264  gaussPts(3, gg) =
266  }
267  }
269  };
270 
271  int nb_gauss_pts, nb_gauss_pts_through_thickness, nb_gauss_pts_on_faces;
272  CHKERR set_gauss_points_on_triangle(nb_gauss_pts_on_faces);
273  if (!nb_gauss_pts_on_faces)
275  CHKERR set_gauss_points_through_thickness(nb_gauss_pts_through_thickness);
276  if (!nb_gauss_pts_through_thickness)
278  CHKERR set_gauss_points_in_volume(
279  nb_gauss_pts_on_faces, nb_gauss_pts_through_thickness, nb_gauss_pts);
280 
281  auto calc_coordinates_at_triangles = [&]() {
283  coordsAtGaussPtsTrianglesOnly.resize(nb_gauss_pts_on_faces, 6, false);
284  for (int gg = 0; gg < nb_gauss_pts_on_faces; gg++) {
285  for (int dd = 0; dd < 3; dd++) {
287  cblas_ddot(3,
288  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
289  NOBASE)(gg, 0),
290  1, &coords[dd], 3);
292  cblas_ddot(3,
293  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
294  NOBASE)(gg, 0),
295  1, &coords[9 + dd], 3);
296  }
297  }
299  };
300 
301  auto calc_vertex_base_on_prism = [&]() {
303  // Calculate "nobase" base functions on prism, this is cartesian product
304  // of base functions on triangles with base functions through thickness
305  // FIXME: This could be effectively implemented with tensors
306  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 6,
307  false);
308  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(nb_gauss_pts, 18,
309  false);
310  for (int dd = 0; dd != 6; dd++) {
311  int gg = 0;
312  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
313  int ddd = dd > 2 ? dd - 3 : dd;
314  double tri_n = dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
315  NOBASE)(ggf, ddd);
316  double dksi_tri_n =
317  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
318  0, 2 * ddd + 0);
319  double deta_tri_n =
320  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
321  0, 2 * ddd + 1);
322  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
323  double zeta = gaussPtsThroughThickness(0, ggt);
324  double dzeta, edge_shape;
325  if (dd < 3) {
326  dzeta = diffN_MBEDGE0;
327  edge_shape = N_MBEDGE0(zeta);
328  } else {
329  dzeta = diffN_MBEDGE1;
330  edge_shape = N_MBEDGE1(zeta);
331  }
332  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, dd) =
333  tri_n * edge_shape;
334  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 0) =
335  dksi_tri_n * edge_shape;
336  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 1) =
337  deta_tri_n * edge_shape;
338  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 2) =
339  tri_n * dzeta;
340  }
341  }
342  }
344  };
345 
346  auto calc_base_on_prism = [&]() {
348  // Calculate base functions on prism
349  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
350  if (dataH1.bAse.test(b)) {
351  switch (static_cast<FieldApproximationBase>(b)) {
354  if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
355  CHKERR FatPrismPolynomialBase().getValue(
356  gaussPts,
357  boost::shared_ptr<BaseFunctionCtx>(
358  new FatPrismPolynomialBaseCtx(
362  static_cast<FieldApproximationBase>(b), NOBASE)));
363  }
364  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
365  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
366  "Not yet implemented");
367  }
368  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
369  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
370  "Not yet implemented");
371  }
372  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
373  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
374  "Not yet implemented");
375  }
376  break;
377  default:
378  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
379  "Not yet implemented");
380  }
381  }
382  }
384  };
385 
386  auto calc_coordinate_on_prism = [&]() {
388  /// Calculate coordinates at integration points
389  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
390  for (int gg = 0; gg < nb_gauss_pts; gg++) {
391  for (int dd = 0; dd < 3; dd++) {
393  6, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
394  &coords[dd], 3);
395  }
396  }
398  };
399 
400  auto calc_ho_triangle_face_normals = [&]() {
402  // Calculate ho-geometry tangents and normals
403  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
404  dataPtr->get<FieldName_mi_tag>().end()) {
405  hoCoordsAtGaussPtsF3.resize(nb_gauss_pts_on_faces, 3, false);
406  nOrmals_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
407  tAngent1_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
408  tAngent2_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
409  hoCoordsAtGaussPtsF4.resize(nb_gauss_pts_on_faces, 3, false);
410  nOrmals_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
411  tAngent1_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
412  tAngent2_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
415  MBEDGE);
417  MBEDGE);
420  } else {
421  hoCoordsAtGaussPtsF3.resize(0, 0, false);
422  nOrmals_at_GaussPtF3.resize(0, 0, false);
423  tAngent1_at_GaussPtF3.resize(0, 0, false);
424  tAngent2_at_GaussPtF3.resize(0, 0, false);
425  hoCoordsAtGaussPtsF4.resize(0, 0, false);
426  nOrmals_at_GaussPtF4.resize(0, 0, false);
427  tAngent1_at_GaussPtF4.resize(0, 0, false);
428  tAngent2_at_GaussPtF4.resize(0, 0, false);
429  }
431  };
432 
433  CHKERR calc_coordinates_at_triangles();
434  CHKERR calc_vertex_base_on_prism();
435  CHKERR calc_base_on_prism();
436  CHKERR calc_coordinate_on_prism();
437  CHKERR calc_ho_triangle_face_normals();
438 
439  // Iterate over operators
441 
443 }
MatrixDouble gaussPts
Matrix of integration points.
field with continuous normal traction
Definition: definitions.h:173
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
virtual moab::Interface & get_moab()=0
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
DataForcesAndSourcesCore & dataH1
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
int npoints
Definition: quad.h:29
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:80
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
Pointer to finite element data dofs.
double cblas_dnrm2(const int N, const double *X, const int incX)
Definition: cblas_dnrm2.c:12
#define diffN_MBEDGE1
derivative of edge shape function
Definition: fem_tools.h:82
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getNodesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
Get data on nodes.
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:79
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
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:146
static MoFEMErrorCode getTriNormal(const double *coords, double *normal)
Get the Tri Normal objectGet triangle normal.
Definition: Tools.cpp:258
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:81
field with continuous tangents
Definition: definitions.h:172
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:87
virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:194
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
virtual MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
continuous field
Definition: definitions.h:171
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
std::bitset< LASTBASE > bAse
bases on element
field with C-1 continuity
Definition: definitions.h:174
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163

◆ setGaussPtsThroughThickness()

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

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

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:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508

◆ setGaussPtsTrianglesOnly()

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

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

Definition at line 51 of file FatPrismElementForcesAndSourcesCore.hpp.

51  {
53  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
55  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

Definition at line 214 of file FatPrismElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

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

Definition at line 194 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly
protected

Definition at line 198 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TrianglesOnly

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly
protected

Definition at line 201 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TroughThickness

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness
protected

Definition at line 202 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 204 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4
protected

Definition at line 208 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ normal

VectorDouble MoFEM::FatPrismElementForcesAndSourcesCore::normal
protected

Definition at line 195 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3
protected

Definition at line 205 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4
protected

Definition at line 209 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnPrism MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals
protected

Definition at line 212 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3
protected

Definition at line 206 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4
protected

Definition at line 210 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3
protected

Definition at line 207 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4
protected

Definition at line 211 of file FatPrismElementForcesAndSourcesCore.hpp.


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