v0.8.23
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 ()
 
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 boost::shared_ptr< DataForcesAndSourcesCoredataOnElement [LASTSPACE]
 Entity data on element entity rows fields. More...
 
const boost::shared_ptr< DataForcesAndSourcesCorederivedDataOnElement [LASTSPACE]
 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...
 

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
 
Vec snes_x
 
Vec snes_f
 
Mat snes_A
 
Mat snes_B
 
- Public Attributes inherited from MoFEM::TSMethod
TSContext ts_ctx
 
TS ts
 
Vec ts_u
 
Vec ts_u_t
 
Vec ts_F
 
Mat ts_A
 
Mat ts_B
 
PetscInt ts_step
 
PetscReal ts_a
 
PetscReal ts_t
 
- 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 calculateBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate base functions. 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 getNodesFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 
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
prisms_elements_from_surface.cpp.

Definition at line 43 of file FatPrismElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FatPrismElementForcesAndSourcesCore()

MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore ( Interface m_field)

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 
40  int num_nodes;
41  const EntityHandle *conn;
42  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
43  {
44  coords.resize(num_nodes * 3, false);
45  CHKERR mField.get_moab().get_coords(conn, num_nodes,
46  &*coords.data().begin());
47  double diff_n[6];
48  CHKERR ShapeDiffMBTRI(diff_n);
49  normal.resize(6, false);
50  CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[0], &normal[0]);
51  CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[9], &normal[3]);
52  aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
53  aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
54  }
55 
59 
60  // H1
61  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1)) {
62  CHKERR getEntitySense<MBEDGE>(dataH1);
63  CHKERR getEntitySense<MBTRI>(dataH1);
64  CHKERR getEntitySense<MBQUAD>(dataH1);
65  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
66  CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
67  CHKERR getEntityDataOrder<MBQUAD>(dataH1, H1);
68  CHKERR getEntityDataOrder<MBPRISM>(dataH1, H1);
69  // Triangles only
70  CHKERR getEntitySense<MBEDGE>(dataH1TrianglesOnly);
71  CHKERR getEntitySense<MBTRI>(dataH1TrianglesOnly);
72  CHKERR getEntityDataOrder<MBEDGE>(dataH1TrianglesOnly, H1);
73  CHKERR getEntityDataOrder<MBTRI>(dataH1TrianglesOnly, H1);
74  // Through thickness
75  CHKERR getEntitySense<MBEDGE>(dataH1TroughThickness);
76  CHKERR getEntityDataOrder<MBEDGE>(dataH1TroughThickness, H1);
77  }
78  // Hdiv
79  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV)) {
80  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
81  }
82  // Hcurl
83  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL)) {
84  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
85  }
86  // L2
87  if ((dataH1.spacesOnEntities[MBPRISM]).test(L2)) {
88  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
89  }
90 
91  // get approx. on triangles, i.e. faces 3 and 4
92  int nb_gauss_pts_on_faces;
93  try {
94  int order_triangles_only = 1;
95  int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
96  for (unsigned int ee = 0; ee < 9; ee++) {
97  if (!valid_edges[ee])
98  continue;
99  order_triangles_only = std::max(
100  order_triangles_only,
101  dataH1TrianglesOnly.dataOnEntities[MBEDGE][ee].getDataOrder());
102  }
103  for (unsigned int ff = 3; ff <= 4; ff++) {
104  order_triangles_only = std::max(
105  order_triangles_only,
106  dataH1TrianglesOnly.dataOnEntities[MBTRI][ff].getDataOrder());
107  }
108  for (unsigned int qq = 0; qq < 3; qq++) {
109  order_triangles_only = std::max(
110  order_triangles_only,
111  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
112  }
113  order_triangles_only = std::max(
114  order_triangles_only,
115  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
116  // integration pts on the triangles surfaces
117  int rule = getRuleTrianglesOnly(order_triangles_only);
118  if (rule >= 0) {
119  if (rule < QUAD_2D_TABLE_SIZE) {
120  if (QUAD_2D_TABLE[rule]->dim != 2) {
121  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
122  }
123  if (QUAD_2D_TABLE[rule]->order < rule) {
124  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
125  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
126  }
127  nb_gauss_pts_on_faces = QUAD_2D_TABLE[rule]->npoints;
128  gaussPtsTrianglesOnly.resize(3, nb_gauss_pts_on_faces, false);
129  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[1], 3,
130  &gaussPtsTrianglesOnly(0, 0), 1);
131  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[2], 3,
132  &gaussPtsTrianglesOnly(1, 0), 1);
133  cblas_dcopy(nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->weights, 1,
134  &gaussPtsTrianglesOnly(2, 0), 1);
135  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
136  nb_gauss_pts_on_faces, 3, false);
137  double *shape_ptr = &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
138  .getN(NOBASE)
139  .data()
140  .begin();
141  cblas_dcopy(3 * nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->points, 1,
142  shape_ptr, 1);
143  } else {
144  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
145  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
146  nb_gauss_pts_on_faces = 0;
147  }
148  } else {
149  CHKERR setGaussPtsTrianglesOnly(order_triangles_only);
150  nb_gauss_pts_on_faces = gaussPtsTrianglesOnly.size2();
151  if (nb_gauss_pts_on_faces == 0)
153  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
154  nb_gauss_pts_on_faces, 3, false);
155  if (nb_gauss_pts_on_faces) {
157  .getN(NOBASE)
158  .data()
159  .begin(),
160  &gaussPtsTrianglesOnly(0, 0),
161  &gaussPtsTrianglesOnly(1, 0), nb_gauss_pts_on_faces);
162  }
163  }
164  }
165  CATCH_ERRORS;
166 
167  // approx. trough prism thickness
168  int nb_gauss_pts_through_thickness;
169  try {
170  int order_thickness = 1;
171  for (unsigned int ee = 3; ee <= 5; ee++) {
172  order_thickness = std::max(
173  order_thickness,
174  dataH1TroughThickness.dataOnEntities[MBEDGE][ee].getDataOrder());
175  }
176  for (unsigned int qq = 0; qq < 3; qq++) {
177  order_thickness = std::max(
178  order_thickness,
179  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
180  }
181  order_thickness = std::max(
182  order_thickness,
183  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
184  // integration points
185  int rule = getRuleThroughThickness(order_thickness);
186  if (rule >= 0) {
187  if (rule < QUAD_1D_TABLE_SIZE) {
188  if (QUAD_1D_TABLE[rule]->dim != 1) {
189  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
190  }
191  if (QUAD_1D_TABLE[rule]->order < rule) {
192  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
193  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order, rule);
194  }
195  nb_gauss_pts_through_thickness = QUAD_1D_TABLE[rule]->npoints;
196  gaussPtsThroughThickness.resize(2, nb_gauss_pts_through_thickness,
197  false);
198  cblas_dcopy(nb_gauss_pts_through_thickness,
199  &QUAD_1D_TABLE[rule]->points[1], 2,
200  &gaussPtsThroughThickness(0, 0), 1);
201  cblas_dcopy(nb_gauss_pts_through_thickness,
202  QUAD_1D_TABLE[rule]->weights, 1,
203  &gaussPtsThroughThickness(1, 0), 1);
204  } else {
205  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
206  "rule > quadrature order %d < %d", rule, QUAD_1D_TABLE_SIZE);
207  nb_gauss_pts_through_thickness = 0;
208  }
209  } else {
210  CHKERR setGaussPtsThroughThickness(order_thickness);
211  nb_gauss_pts_through_thickness = gaussPtsThroughThickness.size2();
212  }
213  if (nb_gauss_pts_through_thickness == 0)
215  }
216  CATCH_ERRORS;
217 
218  // Generate integration pts.
219  int nb_gauss_pts = nb_gauss_pts_on_faces * nb_gauss_pts_through_thickness;
220  gaussPts.resize(4, nb_gauss_pts, false);
221  {
222  int gg = 0;
223  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
224  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
225  gaussPts(0, gg) = gaussPtsTrianglesOnly(0, ggf);
226  gaussPts(1, gg) = gaussPtsTrianglesOnly(1, ggf);
227  gaussPts(2, gg) = gaussPtsThroughThickness(0, ggt);
228  gaussPts(3, gg) =
230  }
231  }
232  }
233 
234  {
235  coordsAtGaussPtsTrianglesOnly.resize(nb_gauss_pts_on_faces, 6, false);
236  for (int gg = 0; gg < nb_gauss_pts_on_faces; gg++) {
237  for (int dd = 0; dd < 3; dd++) {
239  cblas_ddot(3,
240  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
241  NOBASE)(gg, 0),
242  1, &coords[dd], 3);
244  cblas_ddot(3,
245  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
246  NOBASE)(gg, 0),
247  1, &coords[9 + dd], 3);
248  }
249  }
250  // linear for xi,eta and zeta
251  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
252  1, 6, false);
254  .getDiffN(NOBASE)
255  .data()
256  .begin());
257 
258  // Calculate "nobase" base functions on prism, this is cartesian product
259  // of base functions on triangles with base functions through thickness
260  // FIXME: This could be effectively implemented with tensors
261  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 6,
262  false);
263  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(nb_gauss_pts, 18,
264  false);
265  for (int dd = 0; dd != 6; dd++) {
266  int gg = 0;
267  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
268  int ddd = dd > 2 ? dd - 3 : dd;
269  double tri_n = dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
270  NOBASE)(ggf, ddd);
271  double dksi_tri_n =
272  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
273  0, 2 * ddd + 0);
274  double deta_tri_n =
275  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
276  0, 2 * ddd + 1);
277  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
278  double zeta = gaussPtsThroughThickness(0, ggt);
279  double dzeta, edge_shape;
280  if (dd < 3) {
281  dzeta = diffN_MBEDGE0;
282  edge_shape = N_MBEDGE0(zeta);
283  } else {
284  dzeta = diffN_MBEDGE1;
285  edge_shape = N_MBEDGE1(zeta);
286  }
287  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, dd) =
288  tri_n * edge_shape;
289  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 0) =
290  dksi_tri_n * edge_shape;
291  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 1) =
292  deta_tri_n * edge_shape;
293  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 2) =
294  tri_n * dzeta;
295  }
296  }
297  }
298 
299  /// Calculate coordinates at integration points
300  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
301  for (int gg = 0; gg < nb_gauss_pts; gg++) {
302  for (int dd = 0; dd < 3; dd++) {
304  6, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
305  &coords[dd], 3);
306  }
307  }
308  }
309 
310  // Calculate base functions on prism
311  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
312  if (dataH1.bAse.test(b)) {
313  switch (static_cast<FieldApproximationBase>(b)) {
316  if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
317  CHKERR FatPrismPolynomialBase().getValue(
318  gaussPts,
319  boost::shared_ptr<BaseFunctionCtx>(new FatPrismPolynomialBaseCtx(
323  static_cast<FieldApproximationBase>(b), NOBASE)));
324  }
325  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
326  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
327  "Not yet implemented");
328  }
329  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
330  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
331  "Not yet implemented");
332  }
333  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
334  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
335  "Not yet implemented");
336  }
337  break;
338  default:
339  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
340  "Not yet implemented");
341  }
342  }
343  }
344 
345  // Calculate ho-geometry tangents and normals
346 
347  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
348  dataPtr->get<FieldName_mi_tag>().end()) {
349  hoCoordsAtGaussPtsF3.resize(nb_gauss_pts_on_faces, 3, false);
350  nOrmals_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
351  tAngent1_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
352  tAngent2_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
353  hoCoordsAtGaussPtsF4.resize(nb_gauss_pts_on_faces, 3, false);
354  nOrmals_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
355  tAngent1_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
356  tAngent2_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
362  } else {
363  hoCoordsAtGaussPtsF3.resize(0, 0, false);
364  nOrmals_at_GaussPtF3.resize(0, 0, false);
365  tAngent1_at_GaussPtF3.resize(0, 0, false);
366  tAngent2_at_GaussPtF3.resize(0, 0, false);
367  hoCoordsAtGaussPtsF4.resize(0, 0, false);
368  nOrmals_at_GaussPtF4.resize(0, 0, false);
369  tAngent1_at_GaussPtF4.resize(0, 0, false);
370  tAngent2_at_GaussPtF4.resize(0, 0, false);
371  }
372 
373  // Iterate over operators
375 
377 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
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
int npoints
Definition: quad.h:29
PetscErrorCode ShapeFaceNormalMBTRI(double *diffN, const double *coords, double *normal)
Definition: fem_tools.c:241
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
MoFEMErrorCode getNodesFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
Get field data on nodes.
#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.
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
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
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:79
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:206
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:145
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:81
field with continuous tangents
Definition: definitions.h:172
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
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
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:433
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
const int order
approximation order
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, 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: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, and SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh.

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

Member Data Documentation

◆ aRea

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

Definition at line 266 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly
protected

Definition at line 270 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TrianglesOnly

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly
protected

Definition at line 273 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TroughThickness

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness
protected

Definition at line 274 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsThroughThickness

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
protected

Definition at line 271 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
protected

Definition at line 269 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3
protected

Definition at line 276 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4
protected

Definition at line 280 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ normal

VectorDouble MoFEM::FatPrismElementForcesAndSourcesCore::normal
protected

Definition at line 267 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3
protected

Definition at line 277 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4
protected

Definition at line 281 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnPrism MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals
protected

Definition at line 284 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3
protected

Definition at line 278 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4
protected

Definition at line 282 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3
protected

Definition at line 279 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4
protected

Definition at line 283 of file FatPrismElementForcesAndSourcesCore.hpp.


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