v0.8.20
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 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 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)
 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...
 

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

36  : VolumeElementForcesAndSourcesCore(m_field, MBPRISM),
37  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 43 of file FatPrismElementForcesAndSourcesCore.cpp.

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

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

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

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

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: