v0.14.0
Public Member Functions | Protected Member Functions | List of all members
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator Struct Reference

default operator for Flat Prism element More...

#include <src/finite_elements/FatPrismElementForcesAndSourcesCore.hpp>

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

Public Member Functions

double getArea (const int dd)
 get face aRea More...
 
double getAreaF3 ()
 
double getAreaF4 ()
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorAdaptor getNormalF3 ()
 
VectorAdaptor getNormalF4 ()
 
MatrixDoublegetGaussPts ()
 get Gauss pts. in the prism More...
 
MatrixDoublegetGaussPtsTrianglesOnly ()
 get Gauss pts. on triangles More...
 
MatrixDoublegetGaussPtsThroughThickness ()
 get Gauss pts. through thickness More...
 
MatrixDoublegetCoordsAtGaussPts ()
 get coordinates at Gauss pts. More...
 
MatrixDoublegetCoordsAtGaussPtsTrianglesOnly ()
 get coordinates at Gauss pts. More...
 
MatrixDoublegetHOCoordsAtGaussPtsF3 ()
 coordinate at Gauss points on face 3 (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHOCoordsAtGaussPtsF4 ()
 coordinate at Gauss points on face 4 (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetNormalsAtGaussPtsF3 ()
 if higher order geometry return normals at face F3 at Gauss pts. More...
 
MatrixDoublegetNormalsAtGaussPtsF4 ()
 if higher order geometry return normals at face F4 at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPtsF3 (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPtsF4 (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPtF3 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPtF3 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPtF4 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPtF4 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
EntitiesFieldDatagetTrianglesOnlyDataStructure ()
 
EntitiesFieldDatagetTroughThicknessDataStructure ()
 
const FatPrismElementForcesAndSourcesCoregetPrismFE ()
 return pointer to fat prism finite element More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Protected Member Functions

MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPSPACE = 1 << 3,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

default operator for Flat Prism element

Examples
prism_elements_from_surface.cpp, and prism_polynomial_approximation.cpp.

Definition at line 54 of file FatPrismElementForcesAndSourcesCore.hpp.

Member Function Documentation

◆ getArea()

double MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getArea ( const int  dd)
inline

get face aRea

Parameters
ddif dd == 0 it is for face F3 if dd == 1 is for face F4

Definition at line 210 of file FatPrismElementForcesAndSourcesCore.hpp.

210  {
211  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
212 }

◆ getAreaF3()

double MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getAreaF3 ( )
inline

Definition at line 215 of file FatPrismElementForcesAndSourcesCore.hpp.

215  {
216  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
217 }

◆ getAreaF4()

double MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getAreaF4 ( )
inline

Definition at line 219 of file FatPrismElementForcesAndSourcesCore.hpp.

219  {
220  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[1];
221 }

◆ getCoordsAtGaussPts()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts ( )
inline

get coordinates at Gauss pts.

Matrix has size (nb integration points)x(coordinates on F3 and F4 = 6), i.e. coordinates on face F3 and F4

Examples
prism_elements_from_surface.cpp.

Definition at line 260 of file FatPrismElementForcesAndSourcesCore.hpp.

260  {
261  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
263 }

◆ getCoordsAtGaussPtsTrianglesOnly()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPtsTrianglesOnly ( )
inline

get coordinates at Gauss pts.

Matrix has size (nb integration points)x(coordinates on F3 and F4 = 6), i.e. coordinates on face F3 and F4

Definition at line 266 of file FatPrismElementForcesAndSourcesCore.hpp.

266  {
267  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
269 }

◆ getGaussPts()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPts ( )
inline

get Gauss pts. in the prism

Examples
prism_elements_from_surface.cpp.

Definition at line 243 of file FatPrismElementForcesAndSourcesCore.hpp.

243  {
244  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->gaussPts;
245 }

◆ getGaussPtsThroughThickness()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPtsThroughThickness ( )
inline

get Gauss pts. through thickness

Definition at line 254 of file FatPrismElementForcesAndSourcesCore.hpp.

254  {
255  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
257 }

◆ getGaussPtsTrianglesOnly()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPtsTrianglesOnly ( )
inline

get Gauss pts. on triangles

Definition at line 248 of file FatPrismElementForcesAndSourcesCore.hpp.

248  {
249  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
251 }

◆ getHOCoordsAtGaussPtsF3()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getHOCoordsAtGaussPtsF3 ( )
inline

coordinate at Gauss points on face 3 (if hierarchical approximation of element geometry)

Definition at line 272 of file FatPrismElementForcesAndSourcesCore.hpp.

272  {
273  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
275 }

◆ getHOCoordsAtGaussPtsF4()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getHOCoordsAtGaussPtsF4 ( )
inline

coordinate at Gauss points on face 4 (if hierarchical approximation of element geometry)

Definition at line 278 of file FatPrismElementForcesAndSourcesCore.hpp.

278  {
279  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
281 }

◆ getNormal()

VectorDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormal ( )
inline

get triangle normal

Normal has 6 elements, first 3 are for face F3 another three for face F4

Definition at line 224 of file FatPrismElementForcesAndSourcesCore.hpp.

224  {
225  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->normal;
226 }

◆ getNormalF3()

VectorAdaptor MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalF3 ( )
inline

Definition at line 229 of file FatPrismElementForcesAndSourcesCore.hpp.

229  {
230  double *data =
231  &(static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[0]);
232  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
233 }

◆ getNormalF4()

VectorAdaptor MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalF4 ( )
inline

Definition at line 236 of file FatPrismElementForcesAndSourcesCore.hpp.

236  {
237  double *data =
238  &(static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[3]);
239  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
240 }

◆ getNormalsAtGaussPtsF3() [1/2]

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF3 ( )
inline

if higher order geometry return normals at face F3 at Gauss pts.

Face 3 is top face in canonical triangle numeration, see [69]

Definition at line 284 of file FatPrismElementForcesAndSourcesCore.hpp.

284  {
285  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
287 }

◆ getNormalsAtGaussPtsF3() [2/2]

ublas::matrix_row< MatrixDouble > MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF3 ( const int  gg)
inline

if higher order geometry return normals at Gauss pts.

Face 3 is top face in canonical triangle numeration, see [69]

Parameters
gggauss point number

Definition at line 296 of file FatPrismElementForcesAndSourcesCore.hpp.

297  {
298  return ublas::matrix_row<MatrixDouble>(
301  gg);
302 }

◆ getNormalsAtGaussPtsF4() [1/2]

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF4 ( )
inline

if higher order geometry return normals at face F4 at Gauss pts.

Face 4 is top face in canonical triangle numeration, see [69]

Definition at line 290 of file FatPrismElementForcesAndSourcesCore.hpp.

290  {
291  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
293 }

◆ getNormalsAtGaussPtsF4() [2/2]

ublas::matrix_row< MatrixDouble > MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF4 ( const int  gg)
inline

if higher order geometry return normals at Gauss pts.

Face 3 is top face in canonical triangle numeration, see [69]

Parameters
gggauss point number

Definition at line 305 of file FatPrismElementForcesAndSourcesCore.hpp.

306  {
307  return ublas::matrix_row<MatrixDouble>(
310  gg);
311 }

◆ getPrismFE()

const FatPrismElementForcesAndSourcesCore * MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getPrismFE ( )
inline

return pointer to fat prism finite element

Definition at line 350 of file FatPrismElementForcesAndSourcesCore.hpp.

350  {
351  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE);
352 }

◆ getTangent1AtGaussPtF3()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent1AtGaussPtF3 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 314 of file FatPrismElementForcesAndSourcesCore.hpp.

314  {
315  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
317 }

◆ getTangent1AtGaussPtF4()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent1AtGaussPtF4 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 326 of file FatPrismElementForcesAndSourcesCore.hpp.

326  {
327  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
329 }

◆ getTangent2AtGaussPtF3()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent2AtGaussPtF3 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 320 of file FatPrismElementForcesAndSourcesCore.hpp.

320  {
321  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
323 }

◆ getTangent2AtGaussPtF4()

MatrixDouble & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent2AtGaussPtF4 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 332 of file FatPrismElementForcesAndSourcesCore.hpp.

332  {
333  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
335 }

◆ getTrianglesOnlyDataStructure()

EntitiesFieldData & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getTrianglesOnlyDataStructure ( )
inline

Definition at line 338 of file FatPrismElementForcesAndSourcesCore.hpp.

338  {
339  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
341 }

◆ getTroughThicknessDataStructure()

EntitiesFieldData & MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getTroughThicknessDataStructure ( )
inline

Definition at line 344 of file FatPrismElementForcesAndSourcesCore.hpp.

344  {
345  return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
347 }

◆ setPtrFE()

MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Reimplemented from MoFEM::ForcesAndSourcesCore::UserDataOperator.

Definition at line 487 of file FatPrismElementForcesAndSourcesCore.cpp.

488  {
490  if (!(ptrFE = dynamic_cast<FatPrismElementForcesAndSourcesCore *>(ptr)))
491  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
492  "User operator and finite element do not work together");
494 }

The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly
EntitiesFieldData dataH1TrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:193
MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3
MatrixDouble hoCoordsAtGaussPtsF3
Definition: FatPrismElementForcesAndSourcesCore.hpp:196
MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4
MatrixDouble hoCoordsAtGaussPtsF4
Definition: FatPrismElementForcesAndSourcesCore.hpp:200
MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4
MatrixDouble tAngent1_at_GaussPtF4
Definition: FatPrismElementForcesAndSourcesCore.hpp:202
MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3
MatrixDouble tAngent1_at_GaussPtF3
Definition: FatPrismElementForcesAndSourcesCore.hpp:198
MoFEM::FatPrismElementForcesAndSourcesCore::aRea
double aRea[2]
Definition: FatPrismElementForcesAndSourcesCore.hpp:186
MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4
MatrixDouble tAngent2_at_GaussPtF4
Definition: FatPrismElementForcesAndSourcesCore.hpp:203
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
MatrixDouble gaussPtsTrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:189
MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly
MatrixDouble coordsAtGaussPtsTrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:190
MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore
FatPrismElementForcesAndSourcesCore(Interface &m_field)
Definition: FatPrismElementForcesAndSourcesCore.cpp:11
MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
coordinated at gauss points
Definition: ForcesAndSourcesCore.hpp:544
MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness
EntitiesFieldData dataH1TroughThickness
Definition: FatPrismElementForcesAndSourcesCore.hpp:194
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
MatrixDouble gaussPtsThroughThickness
Definition: FatPrismElementForcesAndSourcesCore.hpp:191
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4
MatrixDouble nOrmals_at_GaussPtF4
Definition: FatPrismElementForcesAndSourcesCore.hpp:201
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3
MatrixDouble tAngent2_at_GaussPtF3
Definition: FatPrismElementForcesAndSourcesCore.hpp:199
MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3
MatrixDouble nOrmals_at_GaussPtF3
Definition: FatPrismElementForcesAndSourcesCore.hpp:197
MoFEM::FatPrismElementForcesAndSourcesCore::normal
VectorDouble normal
Definition: FatPrismElementForcesAndSourcesCore.hpp:187
MoFEM::ForcesAndSourcesCore::UserDataOperator::ptrFE
ForcesAndSourcesCore * ptrFE
Definition: ForcesAndSourcesCore.hpp:984