v0.8.15
Classes | Public Member Functions | Public Attributes | 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 preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore
 VolumeElementForcesAndSourcesCore (Interface &m_field, const EntityType type=MBTET)
 
virtual ~VolumeElementForcesAndSourcesCore ()
 
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
virtual DEPRECATED MoFEMErrorCode setIntegartionPts ()
 
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...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
virtual ~ForcesAndSourcesCore ()
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 
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 getSense (EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order of approximation on the entity More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (DataForcesAndSourcesCore &data) const
 
template<EntityType type>
MoFEMErrorCode getEntityFieldDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
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 getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &data) const
 Get field approximation space and base on entities. More...
 
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)
 set integration rule for finite element More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order)
 It will be removed in the future use other variant. More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
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...
 
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...
 
- 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)
 
- 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...
 

Public 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
 
- Public Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCore
VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
std::string meshPositionsFieldName
 
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
 
unsigned int nbGaussPts
 Number of integration points. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
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
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 
int loopSize
 
int rAnk
 
int sIze
 
const RefEntity_multiIndexrefinedEntitiesPtr
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 
const ProblemproblemPtr
 
const Field_multiIndexfieldsPtr
 
const FieldEntity_multiIndexentitiesPtr
 
const DofEntity_multiIndexdofsPtr
 
const FiniteElement_multiIndexfiniteElementsPtr
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 
boost::function< MoFEMErrorCode()> postProcessHook
 
boost::function< MoFEMErrorCode()> operatorHook
 
- 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
 

Additional Inherited Members

- 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
}
 
- 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

Definition at line 43 of file FatPrismElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FatPrismElementForcesAndSourcesCore()

MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore ( Interface m_field)

Definition at line 78 of file FatPrismElementForcesAndSourcesCore.cpp.

80  : VolumeElementForcesAndSourcesCore(m_field, MBPRISM),
81  dataH1TrianglesOnly(MBPRISM), dataH1TroughThickness(MBPRISM),
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 87 of file FatPrismElementForcesAndSourcesCore.cpp.

87  {
89 
90  if (numeredEntFiniteElementPtr->getEntType() != MBPRISM)
93 
95  int num_nodes;
96  const EntityHandle *conn;
97  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
98  {
99  coords.resize(num_nodes * 3, false);
100  CHKERR mField.get_moab().get_coords(conn, num_nodes,
101  &*coords.data().begin());
102  double diff_n[6];
103  CHKERR ShapeDiffMBTRI(diff_n);
104  normal.resize(6, false);
105  CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[0], &normal[0]);
106  CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[9], &normal[3]);
107  aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
108  aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
109  }
110 
114 
115  // H1
116  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1)) {
117  CHKERR getEntitySense<MBEDGE>(dataH1);
118  CHKERR getEntitySense<MBTRI>(dataH1);
119  CHKERR getEntitySense<MBQUAD>(dataH1);
120  CHKERR getEntityFieldDataOrder<MBEDGE>(dataH1, H1);
121  CHKERR getEntityFieldDataOrder<MBTRI>(dataH1, H1);
122  CHKERR getEntityFieldDataOrder<MBQUAD>(dataH1, H1);
123  CHKERR getEntityFieldDataOrder<MBPRISM>(dataH1, H1);
124  // Triangles only
125  CHKERR getEntitySense<MBEDGE>(dataH1TrianglesOnly);
126  CHKERR getEntitySense<MBTRI>(dataH1TrianglesOnly);
127  CHKERR getEntityFieldDataOrder<MBEDGE>(dataH1TrianglesOnly, H1);
128  CHKERR getEntityFieldDataOrder<MBTRI>(dataH1TrianglesOnly, H1);
129  // Through thickness
130  CHKERR getEntitySense<MBEDGE>(dataH1TroughThickness);
131  CHKERR getEntityFieldDataOrder<MBEDGE>(dataH1TroughThickness, H1);
132  }
133  // Hdiv
134  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV)) {
135  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
136  }
137  // Hcurl
138  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL)) {
139  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
140  }
141  // L2
142  if ((dataH1.spacesOnEntities[MBPRISM]).test(L2)) {
143  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
144  }
145 
146  // get approx. on triangles, i.e. faces 3 and 4
147  int nb_gauss_pts_on_faces;
148  try {
149  int order_triangles_only = 1;
150  int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
151  for (unsigned int ee = 0; ee < 9; ee++) {
152  if (!valid_edges[ee])
153  continue;
154  order_triangles_only = std::max(
155  order_triangles_only,
156  dataH1TrianglesOnly.dataOnEntities[MBEDGE][ee].getDataOrder());
157  }
158  for (unsigned int ff = 3; ff <= 4; ff++) {
159  order_triangles_only = std::max(
160  order_triangles_only,
161  dataH1TrianglesOnly.dataOnEntities[MBTRI][ff].getDataOrder());
162  }
163  for (unsigned int qq = 0; qq < 3; qq++) {
164  order_triangles_only = std::max(
165  order_triangles_only,
166  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
167  }
168  order_triangles_only = std::max(
169  order_triangles_only,
170  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
171  // integration pts on the triangles surfaces
172  int rule = getRuleTrianglesOnly(order_triangles_only);
173  if (rule >= 0) {
174  if (rule < QUAD_2D_TABLE_SIZE) {
175  if (QUAD_2D_TABLE[rule]->dim != 2) {
176  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
177  }
178  if (QUAD_2D_TABLE[rule]->order < rule) {
179  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
180  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
181  }
182  nb_gauss_pts_on_faces = QUAD_2D_TABLE[rule]->npoints;
183  gaussPtsTrianglesOnly.resize(3, nb_gauss_pts_on_faces, false);
184  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[1], 3,
185  &gaussPtsTrianglesOnly(0, 0), 1);
186  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[2], 3,
187  &gaussPtsTrianglesOnly(1, 0), 1);
188  cblas_dcopy(nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->weights, 1,
189  &gaussPtsTrianglesOnly(2, 0), 1);
190  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
191  nb_gauss_pts_on_faces, 3, false);
192  double *shape_ptr = &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
193  .getN(NOBASE)
194  .data()
195  .begin();
196  cblas_dcopy(3 * nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->points, 1,
197  shape_ptr, 1);
198  } else {
199  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
200  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
201  nb_gauss_pts_on_faces = 0;
202  }
203  } else {
204  CHKERR setGaussPtsTrianglesOnly(order_triangles_only);
205  nb_gauss_pts_on_faces = gaussPtsTrianglesOnly.size2();
206  if (nb_gauss_pts_on_faces == 0)
208  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
209  nb_gauss_pts_on_faces, 3, false);
210  if (nb_gauss_pts_on_faces) {
212  .getN(NOBASE)
213  .data()
214  .begin(),
215  &gaussPtsTrianglesOnly(0, 0),
216  &gaussPtsTrianglesOnly(1, 0), nb_gauss_pts_on_faces);
217  }
218  }
219  }
220  CATCH_ERRORS;
221 
222  // approx. trough prism thickness
223  int nb_gauss_pts_through_thickness;
224  try {
225  int order_thickness = 1;
226  for (unsigned int ee = 3; ee <= 5; ee++) {
227  order_thickness = std::max(
228  order_thickness,
229  dataH1TroughThickness.dataOnEntities[MBEDGE][ee].getDataOrder());
230  }
231  for (unsigned int qq = 0; qq < 3; qq++) {
232  order_thickness = std::max(
233  order_thickness,
234  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
235  }
236  order_thickness = std::max(
237  order_thickness,
238  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
239  // integration points
240  int rule = getRuleThroughThickness(order_thickness);
241  if (rule >= 0) {
242  if (rule < QUAD_1D_TABLE_SIZE) {
243  if (QUAD_1D_TABLE[rule]->dim != 1) {
244  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
245  }
246  if (QUAD_1D_TABLE[rule]->order < rule) {
247  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
248  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order, rule);
249  }
250  nb_gauss_pts_through_thickness = QUAD_1D_TABLE[rule]->npoints;
251  gaussPtsThroughThickness.resize(2, nb_gauss_pts_through_thickness,
252  false);
253  cblas_dcopy(nb_gauss_pts_through_thickness,
254  &QUAD_1D_TABLE[rule]->points[1], 2,
255  &gaussPtsThroughThickness(0, 0), 1);
256  cblas_dcopy(nb_gauss_pts_through_thickness,
257  QUAD_1D_TABLE[rule]->weights, 1,
258  &gaussPtsThroughThickness(1, 0), 1);
259  } else {
260  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
261  "rule > quadrature order %d < %d", rule, QUAD_1D_TABLE_SIZE);
262  nb_gauss_pts_through_thickness = 0;
263  }
264  } else {
265  CHKERR setGaussPtsThroughThickness(order_thickness);
266  nb_gauss_pts_through_thickness = gaussPtsThroughThickness.size2();
267  }
268  if (nb_gauss_pts_through_thickness == 0)
270  }
271  CATCH_ERRORS;
272 
273  // Generate integration pts.
274  int nb_gauss_pts = nb_gauss_pts_on_faces * nb_gauss_pts_through_thickness;
275  gaussPts.resize(4, nb_gauss_pts, false);
276  {
277  int gg = 0;
278  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
279  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
280  gaussPts(0, gg) = gaussPtsTrianglesOnly(0, ggf);
281  gaussPts(1, gg) = gaussPtsTrianglesOnly(1, ggf);
282  gaussPts(2, gg) = gaussPtsThroughThickness(0, ggt);
283  gaussPts(3, gg) =
285  }
286  }
287  }
288 
289  {
290  coordsAtGaussPtsTrianglesOnly.resize(nb_gauss_pts_on_faces, 6, false);
291  for (int gg = 0; gg < nb_gauss_pts_on_faces; gg++) {
292  for (int dd = 0; dd < 3; dd++) {
294  cblas_ddot(3,
295  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
296  NOBASE)(gg, 0),
297  1, &coords[dd], 3);
299  cblas_ddot(3,
300  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
301  NOBASE)(gg, 0),
302  1, &coords[9 + dd], 3);
303  }
304  }
305  // linear for xi,eta and zeta
306  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
307  1, 6, false);
309  .getDiffN(NOBASE)
310  .data()
311  .begin());
312 
313  // Calculate "nobase" base functions on prism, this is cartesian product
314  // of base functions on triangles with base functions through thickness
315  // FIXME: This could be effectively implemented with tensors
316  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 6,
317  false);
318  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(nb_gauss_pts, 18,
319  false);
320  for (int dd = 0; dd != 6; dd++) {
321  int gg = 0;
322  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
323  int ddd = dd > 2 ? dd - 3 : dd;
324  double tri_n = dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
325  NOBASE)(ggf, ddd);
326  double dksi_tri_n =
327  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
328  0, 2 * ddd + 0);
329  double deta_tri_n =
330  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
331  0, 2 * ddd + 1);
332  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
333  double zeta = gaussPtsThroughThickness(0, ggt);
334  double dzeta, edge_shape;
335  if (dd < 3) {
336  dzeta = diffN_MBEDGE0;
337  edge_shape = N_MBEDGE0(zeta);
338  } else {
339  dzeta = diffN_MBEDGE1;
340  edge_shape = N_MBEDGE1(zeta);
341  }
342  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, dd) =
343  tri_n * edge_shape;
344  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 0) =
345  dksi_tri_n * edge_shape;
346  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 1) =
347  deta_tri_n * edge_shape;
348  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 2) =
349  tri_n * dzeta;
350  }
351  }
352  }
353 
354  /// Calculate coordinates at integration points
355  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
356  for (int gg = 0; gg < nb_gauss_pts; gg++) {
357  for (int dd = 0; dd < 3; dd++) {
359  6, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
360  &coords[dd], 3);
361  }
362  }
363  }
364 
365  // Calculate base functions on prism
366  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
367  if (dataH1.bAse.test(b)) {
368  switch (static_cast<FieldApproximationBase>(b)) {
371  if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
372  CHKERR FatPrismPolynomialBase().getValue(
373  gaussPts,
374  boost::shared_ptr<BaseFunctionCtx>(new FatPrismPolynomialBaseCtx(
378  static_cast<FieldApproximationBase>(b), NOBASE)));
379  }
380  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
381  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
382  "Not yet implemented");
383  }
384  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
385  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
386  "Not yet implemented");
387  }
388  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
389  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
390  "Not yet implemented");
391  }
392  break;
393  default:
394  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
395  "Not yet implemented");
396  }
397  }
398  }
399 
400  // Calculate ho-geometry tangents and normals
401 
402  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
403  dataPtr->get<FieldName_mi_tag>().end()) {
404  hoCoordsAtGaussPtsF3.resize(nb_gauss_pts_on_faces, 3, false);
405  nOrmals_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
406  tAngent1_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
407  tAngent2_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
408  hoCoordsAtGaussPtsF4.resize(nb_gauss_pts_on_faces, 3, false);
409  nOrmals_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
410  tAngent1_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
411  tAngent2_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
417  } else {
418  hoCoordsAtGaussPtsF3.resize(0, 0, false);
419  nOrmals_at_GaussPtF3.resize(0, 0, false);
420  tAngent1_at_GaussPtF3.resize(0, 0, false);
421  tAngent2_at_GaussPtF3.resize(0, 0, false);
422  hoCoordsAtGaussPtsF4.resize(0, 0, false);
423  nOrmals_at_GaussPtF4.resize(0, 0, false);
424  tAngent1_at_GaussPtF4.resize(0, 0, false);
425  tAngent2_at_GaussPtF4.resize(0, 0, false);
426  }
427 
428  // Iterate over operators
430 
432 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
MatrixDouble gaussPts
Matrix of integration points.
field with continuous normal traction
Definition: definitions.h:170
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
virtual moab::Interface & get_moab()=0
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
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:219
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:459
#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:490
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
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:186
FieldApproximationBase
approximation base
Definition: definitions.h:140
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:81
field with continuous tangents
Definition: definitions.h:169
virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
#define CHKERR
Inline error check.
Definition: definitions.h:578
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:175
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:403
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:429
continuous field
Definition: definitions.h:168
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:171
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163

◆ postProcess()

MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::postProcess ( )
virtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::ForcesAndSourcesCore.

Reimplemented in PostProcFatPrismOnRefinedMesh, SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh, SolidShellModule::SolidShellPrismElement::SolidShellError, and SolidShellModule::SolidShellPrismElement::SolidShell.

Definition at line 348 of file FatPrismElementForcesAndSourcesCore.hpp.

348  {
351  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490

◆ preProcess()

MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::preProcess ( )
virtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from MoFEM::ForcesAndSourcesCore.

Reimplemented in PostProcFatPrismOnRefinedMesh, SolidShellModule::SolidShellPrismElement::SolidShellError, and SolidShellModule::SolidShellPrismElement::SolidShell.

Definition at line 343 of file FatPrismElementForcesAndSourcesCore.hpp.

343  {
346  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490

◆ setGaussPtsThroughThickness()

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

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

Definition at line 77 of file FatPrismElementForcesAndSourcesCore.hpp.

77  {
79  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
81  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490

◆ setGaussPtsTrianglesOnly()

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

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

Definition at line 71 of file FatPrismElementForcesAndSourcesCore.hpp.

71  {
73  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
75  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490

Member Data Documentation

◆ aRea

double MoFEM::FatPrismElementForcesAndSourcesCore::aRea[2]

Definition at line 46 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly

Definition at line 50 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TrianglesOnly

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly

Definition at line 53 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TroughThickness

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness

Definition at line 54 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsThroughThickness

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness

Definition at line 51 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly

Definition at line 49 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3

Definition at line 56 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4

Definition at line 60 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ normal

VectorDouble MoFEM::FatPrismElementForcesAndSourcesCore::normal

Definition at line 47 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3

Definition at line 57 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4

Definition at line 61 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnPrism MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals

Definition at line 64 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3

Definition at line 58 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4

Definition at line 62 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3

Definition at line 59 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4

Definition at line 63 of file FatPrismElementForcesAndSourcesCore.hpp.


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