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

FatPrism finite element. More...

#include <src/finite_elements/FatPrismElementForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 default operator for Flat Prism element More...
 

Public Member Functions

 FatPrismElementForcesAndSourcesCore (Interface &m_field)
 
virtual int getRuleTrianglesOnly (int order)
 
virtual int getRuleThroughThickness (int order)
 
virtual MoFEMErrorCode setGaussPtsTrianglesOnly (int order_triangles_only)
 
virtual MoFEMErrorCode setGaussPtsThroughThickness (int order_thickness)
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore
 VolumeElementForcesAndSourcesCore (Interface &m_field, const EntityType type=MBTET)
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
boost::ptr_deque< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator. More...
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object. More...
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
int getMaxDataOrder () const
 Get max order of approximation for data fields. More...
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows. More...
 
int getMaxColOrder () const
 Get max order of approximation for field in columns. More...
 
auto & getEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
auto & getDataOnElementBySpaceArray ()
 Get data on entities and space. More...
 
auto & getDerivedDataOnElementBySpaceArray ()
 Get derived data on entities and space. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Protected Attributes

double aRea [2]
 
VectorDouble normal
 
MatrixDouble gaussPtsTrianglesOnly
 
MatrixDouble coordsAtGaussPtsTrianglesOnly
 
MatrixDouble gaussPtsThroughThickness
 
EntitiesFieldData dataH1TrianglesOnly
 
EntitiesFieldData dataH1TroughThickness
 
MatrixDouble hoCoordsAtGaussPtsF3
 
MatrixDouble nOrmals_at_GaussPtF3
 
MatrixDouble tAngent1_at_GaussPtF3
 
MatrixDouble tAngent2_at_GaussPtF3
 
MatrixDouble hoCoordsAtGaussPtsF4
 
MatrixDouble nOrmals_at_GaussPtF4
 
MatrixDouble tAngent1_at_GaussPtF4
 
MatrixDouble tAngent2_at_GaussPtF4
 
OpGetCoordsAndNormalsOnPrism opHOCoordsAndNormals
 
- Protected Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCore
VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
OpSetInvJacH1 opSetInvJacH1
 
OpSetContravariantPiolaTransform opContravariantPiolaTransform
 
OpSetCovariantPiolaTransform opCovariantPiolaTransform
 
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
 
doublevOlume
 
int num_nodes
 
const EntityHandleconn
 
FTensor::Tensor2< double *, 3, 3 > tJac
 
FTensor::Tensor2< double *, 3, 3 > tInvJac
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
EntitiesFieldDatadataNoField
 
EntitiesFieldDatadataH1
 
EntitiesFieldDatadataHcurl
 
EntitiesFieldDatadataHdiv
 
EntitiesFieldDatadataL2
 
boost::ptr_deque< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 
MatrixDouble coordsAtGaussPts
 coordinated at gauss points More...
 
double elementMeasure
 

Friends

class UserDataOperator
 

Additional Inherited Members

- Public Types inherited from MoFEM::VolumeElementForcesAndSourcesCore
enum  Switches
 
- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
typedef boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> GaussHookFun
 
- 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::PetscData
enum  DataContext {
  CTX_SET_NONE = 0, CTX_SET_F = 1 << 0, CTX_SET_A = 1 << 1, CTX_SET_B = 1 << 2,
  CTX_SET_X = 1 << 3, CTX_SET_X_T = 1 << 4, CTX_SET_X_TT = 1 << 6, CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- 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
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Public Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCore
std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
GaussHookFun setRuleHook
 Set function to calculate integration rule. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
PetscReal ts_dt
 time step size More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
virtual MoFEMErrorCode calculateVolumeAndJacobian ()
 Calculate element volume and Jacobian. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode transformBaseFunctions ()
 Transform base functions based on geometric element Jacobian. More...
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (EntitiesFieldData &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceNodes (EntitiesFieldData &data) const
 Get nodes on faces. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (EntitiesFieldData &data) const
 Get field approximation space and base on entities. More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement ()
 Calculate base functions. More...
 
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement ()
 Calculate Bernstein-Bezier base. More...
 
MoFEMErrorCode createDataOnElement (EntityType type)
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getNodesIndices (const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getEntityIndices (EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
 
MoFEMErrorCode getEntityRowIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getBitRefLevelOnData ()
 
MoFEMErrorCode getNoFieldFieldData (const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (EntitiesFieldData &data, const int bit_number) const
 
MoFEMErrorCode getNodesFieldData (EntitiesFieldData &data, const int bit_number) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getProblemNodesIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get indices of nodal indices which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemTypeIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const
 get indices by type (generic function) which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemNodesRowIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeRowIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
MoFEMErrorCode getProblemNodesColIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeColIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
virtual int getRule (int order)
 
virtual MoFEMErrorCode setGaussPts (int order)
 

Detailed Description

FatPrism finite element.

User is implementing own operator at Gauss points level, by own object derived from FatPrismElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to rowOpPtrVector and rowColOpPtrVector.

Todo:
Need to implement operators that will make this element work as Volume element
Examples
elasticity.cpp, HookeElement.cpp, prism_elements_from_surface.cpp, and prism_polynomial_approximation.cpp.

Definition at line 31 of file FatPrismElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FatPrismElementForcesAndSourcesCore()

MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore ( Interface m_field)

Member Function Documentation

◆ getRuleThroughThickness()

virtual int MoFEM::FatPrismElementForcesAndSourcesCore::getRuleThroughThickness ( int  order)
inlinevirtual

Reimplemented in PostProcFatPrismOnRefinedMesh, PrismFE, PrismFE, PrismFE, and CellEngineering::FatPrism.

Definition at line 37 of file FatPrismElementForcesAndSourcesCore.hpp.

37 { return 2 * order; };

◆ getRuleTrianglesOnly()

virtual int MoFEM::FatPrismElementForcesAndSourcesCore::getRuleTrianglesOnly ( int  order)
inlinevirtual

Reimplemented in PostProcFatPrismOnRefinedMesh, PrismFE, PrismFE, PrismFE, and CellEngineering::FatPrism.

Definition at line 36 of file FatPrismElementForcesAndSourcesCore.hpp.

36 { return 2 * order; };

◆ 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.

Calculate coordinates at integration points

Reimplemented from MoFEM::ForcesAndSourcesCore.

Definition at line 23 of file FatPrismElementForcesAndSourcesCore.cpp.

23  {
25 
26  if (numeredEntFiniteElementPtr->getEntType() != MBPRISM)
28 
29  auto get_fe_coordinates = [&]() {
32  int num_nodes;
33  const EntityHandle *conn;
34  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
35  coords.resize(num_nodes * 3, false);
36  CHKERR mField.get_moab().get_coords(conn, num_nodes,
37  &*coords.data().begin());
39  };
40 
41  auto calculate_area_of_triangles = [&] {
43  normal.resize(6, false);
46  aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
47  aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
49  };
50 
51  CHKERR get_fe_coordinates();
52  CHKERR calculate_area_of_triangles();
53 
57 
58  auto get_h1_base_data = [&](auto &dataH1) {
60  CHKERR getEntitySense<MBEDGE>(dataH1);
61  CHKERR getEntitySense<MBTRI>(dataH1);
62  CHKERR getEntitySense<MBQUAD>(dataH1);
63  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
64  CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
65  CHKERR getEntityDataOrder<MBQUAD>(dataH1, H1);
66  CHKERR getEntityDataOrder<MBPRISM>(dataH1, H1);
67  // Triangles only
68  CHKERR getEntitySense<MBEDGE>(dataH1TrianglesOnly);
69  CHKERR getEntitySense<MBTRI>(dataH1TrianglesOnly);
70  CHKERR getEntityDataOrder<MBEDGE>(dataH1TrianglesOnly, H1);
71  CHKERR getEntityDataOrder<MBTRI>(dataH1TrianglesOnly, H1);
72  // Through thickness
73  CHKERR getEntitySense<MBEDGE>(dataH1TroughThickness);
74  CHKERR getEntityDataOrder<MBEDGE>(dataH1TroughThickness, H1);
76  };
77 
78  // H1
79  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1))
80  CHKERR get_h1_base_data(dataH1);
81 
82  // Hdiv
83  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV))
84  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
85 
86  // Hcurl
87  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL))
88  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
89 
90  // L2
91  if ((dataH1.spacesOnEntities[MBPRISM]).test(L2))
92  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
93 
94  // get approx. on triangles, i.e. faces 3 and 4
95  auto set_gauss_points_on_triangle = [&](int &nb_gauss_pts_on_faces) {
97  int order_triangles_only = 1;
98  int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
99  for (unsigned int ee = 0; ee < 9; ee++) {
100  if (!valid_edges[ee])
101  continue;
102  order_triangles_only = std::max(
103  order_triangles_only,
104  dataH1TrianglesOnly.dataOnEntities[MBEDGE][ee].getOrder());
105  }
106  for (unsigned int ff = 3; ff <= 4; ff++) {
107  order_triangles_only = std::max(
108  order_triangles_only,
109  dataH1TrianglesOnly.dataOnEntities[MBTRI][ff].getOrder());
110  }
111  for (unsigned int qq = 0; qq < 3; qq++) {
112  order_triangles_only = std::max(
113  order_triangles_only,
114  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getOrder());
115  }
116  order_triangles_only = std::max(
117  order_triangles_only,
118  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getOrder());
119 
120  // integration pts on the triangles surfaces
121  nb_gauss_pts_on_faces = 0;
122  int rule = getRuleTrianglesOnly(order_triangles_only);
123  if (rule >= 0) {
124  if (rule < QUAD_2D_TABLE_SIZE) {
125  if (QUAD_2D_TABLE[rule]->dim != 2) {
126  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
127  }
128  if (QUAD_2D_TABLE[rule]->order < rule) {
129  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
130  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
131  }
132  nb_gauss_pts_on_faces = QUAD_2D_TABLE[rule]->npoints;
133  gaussPtsTrianglesOnly.resize(3, nb_gauss_pts_on_faces, false);
134  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[1], 3,
135  &gaussPtsTrianglesOnly(0, 0), 1);
136  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[2], 3,
137  &gaussPtsTrianglesOnly(1, 0), 1);
138  cblas_dcopy(nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->weights, 1,
139  &gaussPtsTrianglesOnly(2, 0), 1);
140  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
141  nb_gauss_pts_on_faces, 3, false);
142  double *shape_ptr = &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
143  .getN(NOBASE)
144  .data()
145  .begin();
146  cblas_dcopy(3 * nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->points, 1,
147  shape_ptr, 1);
148  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
149  1, 6, false);
150  std::copy(Tools::diffShapeFunMBTRI.begin(),
152  &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
153  .getDiffN(NOBASE)
154  .data()
155  .begin());
156  } else
157  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
158  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
159 
160  } else {
161  CHKERR setGaussPtsTrianglesOnly(order_triangles_only);
162  nb_gauss_pts_on_faces = gaussPtsTrianglesOnly.size2();
163  if (nb_gauss_pts_on_faces == 0)
165  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
166  nb_gauss_pts_on_faces, 3, false);
167  if (nb_gauss_pts_on_faces) {
169  .getN(NOBASE)
170  .data()
171  .begin(),
172  &gaussPtsTrianglesOnly(0, 0),
173  &gaussPtsTrianglesOnly(1, 0), nb_gauss_pts_on_faces);
174  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
175  1, 6, false);
176  std::copy(Tools::diffShapeFunMBTRI.begin(),
178  &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
179  .getDiffN(NOBASE)
180  .data()
181  .begin());
182  }
183  }
185  };
186 
187  // approx. trough prism thickness
188  auto set_gauss_points_through_thickness =
189  [&](int &nb_gauss_pts_through_thickness) {
191  nb_gauss_pts_through_thickness = 0;
192  int order_thickness = 1;
193  for (unsigned int ee = 3; ee <= 5; ee++) {
194  order_thickness = std::max(
195  order_thickness,
196  dataH1TroughThickness.dataOnEntities[MBEDGE][ee].getOrder());
197  }
198  for (unsigned int qq = 0; qq < 3; qq++) {
199  order_thickness = std::max(
200  order_thickness,
201  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getOrder());
202  }
203  order_thickness = std::max(
204  order_thickness,
205  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getOrder());
206  // integration points
207  int rule = getRuleThroughThickness(order_thickness);
208  if (rule >= 0) {
209  if (rule < QUAD_1D_TABLE_SIZE) {
210  if (QUAD_1D_TABLE[rule]->dim != 1) {
211  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
212  "wrong dimension");
213  }
214  if (QUAD_1D_TABLE[rule]->order < rule) {
215  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
216  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order,
217  rule);
218  }
219  nb_gauss_pts_through_thickness = QUAD_1D_TABLE[rule]->npoints;
220  gaussPtsThroughThickness.resize(2, nb_gauss_pts_through_thickness,
221  false);
222  cblas_dcopy(nb_gauss_pts_through_thickness,
223  &QUAD_1D_TABLE[rule]->points[1], 2,
224  &gaussPtsThroughThickness(0, 0), 1);
225  cblas_dcopy(nb_gauss_pts_through_thickness,
226  QUAD_1D_TABLE[rule]->weights, 1,
227  &gaussPtsThroughThickness(1, 0), 1);
228  } else {
229  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
230  "rule > quadrature order %d < %d", rule,
232  nb_gauss_pts_through_thickness = 0;
233  }
234  } else {
235  CHKERR setGaussPtsThroughThickness(order_thickness);
236  nb_gauss_pts_through_thickness = gaussPtsThroughThickness.size2();
237  }
239  };
240 
241  // Generate integration pts.
242  auto set_gauss_points_in_volume = [&](int nb_gauss_pts_on_faces,
243  int nb_gauss_pts_through_thickness,
244  int &nb_gauss_pts) {
246  nb_gauss_pts = nb_gauss_pts_on_faces * nb_gauss_pts_through_thickness;
247  gaussPts.resize(4, nb_gauss_pts, false);
248  int gg = 0;
249  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
250  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
251  gaussPts(0, gg) = gaussPtsTrianglesOnly(0, ggf);
252  gaussPts(1, gg) = gaussPtsTrianglesOnly(1, ggf);
253  gaussPts(2, gg) = gaussPtsThroughThickness(0, ggt);
254  gaussPts(3, gg) =
256  }
257  }
259  };
260 
261  int nb_gauss_pts, nb_gauss_pts_through_thickness, nb_gauss_pts_on_faces;
262  CHKERR set_gauss_points_on_triangle(nb_gauss_pts_on_faces);
263  if (!nb_gauss_pts_on_faces)
265  CHKERR set_gauss_points_through_thickness(nb_gauss_pts_through_thickness);
266  if (!nb_gauss_pts_through_thickness)
268  CHKERR set_gauss_points_in_volume(
269  nb_gauss_pts_on_faces, nb_gauss_pts_through_thickness, nb_gauss_pts);
270 
271  auto calc_coordinates_at_triangles = [&]() {
273  coordsAtGaussPtsTrianglesOnly.resize(nb_gauss_pts_on_faces, 6, false);
274  for (int gg = 0; gg < nb_gauss_pts_on_faces; gg++) {
275  for (int dd = 0; dd < 3; dd++) {
277  cblas_ddot(3,
278  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
279  NOBASE)(gg, 0),
280  1, &coords[dd], 3);
282  cblas_ddot(3,
283  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
284  NOBASE)(gg, 0),
285  1, &coords[9 + dd], 3);
286  }
287  }
289  };
290 
291  auto calc_vertex_base_on_prism = [&]() {
293  // Calculate "nobase" base functions on prism, this is cartesian product
294  // of base functions on triangles with base functions through thickness
295  // FIXME: This could be effectively implemented with tensors
296  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 6,
297  false);
298  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(nb_gauss_pts, 18,
299  false);
300  for (int dd = 0; dd != 6; dd++) {
301  int gg = 0;
302  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
303  int ddd = dd > 2 ? dd - 3 : dd;
304  double tri_n = dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
305  NOBASE)(ggf, ddd);
306  double dksi_tri_n =
307  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
308  0, 2 * ddd + 0);
309  double deta_tri_n =
310  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
311  0, 2 * ddd + 1);
312  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
313  double zeta = gaussPtsThroughThickness(0, ggt);
314  double dzeta, edge_shape;
315  if (dd < 3) {
316  dzeta = diffN_MBEDGE0;
317  edge_shape = N_MBEDGE0(zeta);
318  } else {
319  dzeta = diffN_MBEDGE1;
320  edge_shape = N_MBEDGE1(zeta);
321  }
322  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, dd) =
323  tri_n * edge_shape;
324  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 0) =
325  dksi_tri_n * edge_shape;
326  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 1) =
327  deta_tri_n * edge_shape;
328  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 2) =
329  tri_n * dzeta;
330  }
331  }
332  }
334  };
335 
336  auto calc_base_on_prism = [&]() {
338  // Calculate base functions on prism
339  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
340  if (dataH1.bAse.test(b)) {
341  switch (static_cast<FieldApproximationBase>(b)) {
344  if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
345  CHKERR FatPrismPolynomialBase().getValue(
346  gaussPts,
347  boost::shared_ptr<BaseFunctionCtx>(
348  new FatPrismPolynomialBaseCtx(
352  static_cast<FieldApproximationBase>(b), NOBASE)));
353  }
354  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
355  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
356  "Not yet implemented");
357  }
358  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
359  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
360  "Not yet implemented");
361  }
362  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
363  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
364  "Not yet implemented");
365  }
366  break;
367  default:
368  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
369  "Not yet implemented");
370  }
371  }
372  }
374  };
375 
376  auto calc_coordinate_on_prism = [&]() {
378  /// Calculate coordinates at integration points
379  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
380  for (int gg = 0; gg < nb_gauss_pts; gg++) {
381  for (int dd = 0; dd < 3; dd++) {
382  coordsAtGaussPts(gg, dd) = cblas_ddot(
383  6, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
384  &coords[dd], 3);
385  }
386  }
388  };
389 
390  auto calculate_volume = [&]() {
391  auto get_t_w = [&] {
393  };
394 
395  auto get_t_coords = [&]() {
397  &coords[0], &coords[1], &coords[2]);
398  };
399 
403 
404  const size_t nb_gauss_pts = gaussPts.size2();
405  auto t_diff_n =
406  dataH1.dataOnEntities[MBVERTEX][0].getFTensor1DiffN<3>(NOBASE);
407 
408  double vol = 0;
409  auto t_w = get_t_w();
410  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
411 
412  auto t_coords = get_t_coords();
413  t_jac(i, j) = 0;
414  for (size_t n = 0; n != 6; ++n) {
415  t_jac(i, j) += t_coords(i) * t_diff_n(j);
416  ++t_diff_n;
417  ++t_coords;
418  }
419 
420  double det;
421  CHKERR determinantTensor3by3(t_jac, det);
422  vol += det * t_w / 2;
423 
424  ++t_w;
425  }
426 
427  return vol;
428  };
429 
430  auto calc_ho_triangle_face_normals = [&]() {
432 
433  auto check_field = [&]() {
434  auto field_it =
435  fieldsPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName);
436  if (field_it != fieldsPtr->get<FieldName_mi_tag>().end())
437  if ((numeredEntFiniteElementPtr->getBitFieldIdData() &
438  (*field_it)->getId())
439  .any())
440  return true;
441  return false;
442  };
443 
444  // Check if field meshPositionsFieldName exist
445  if (check_field()) {
446  hoCoordsAtGaussPtsF3.resize(nb_gauss_pts_on_faces, 3, false);
447  nOrmals_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
448  tAngent1_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
449  tAngent2_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
450  hoCoordsAtGaussPtsF4.resize(nb_gauss_pts_on_faces, 3, false);
451  nOrmals_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
452  tAngent1_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
453  tAngent2_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
454  const auto bit_number =
457  CHKERR getEntityFieldData(dataH1TrianglesOnly, bit_number, MBEDGE);
458  CHKERR getEntityFieldData(dataH1TrianglesOnly, bit_number, MBEDGE);
461  } else {
462  hoCoordsAtGaussPtsF3.resize(0, 0, false);
463  nOrmals_at_GaussPtF3.resize(0, 0, false);
464  tAngent1_at_GaussPtF3.resize(0, 0, false);
465  tAngent2_at_GaussPtF3.resize(0, 0, false);
466  hoCoordsAtGaussPtsF4.resize(0, 0, false);
467  nOrmals_at_GaussPtF4.resize(0, 0, false);
468  tAngent1_at_GaussPtF4.resize(0, 0, false);
469  tAngent2_at_GaussPtF4.resize(0, 0, false);
470  }
472  };
473 
474  CHKERR calc_coordinates_at_triangles();
475  CHKERR calc_vertex_base_on_prism();
476  CHKERR calc_base_on_prism();
477  CHKERR calc_coordinate_on_prism();
478  vOlume = calculate_volume();
479  CHKERR calc_ho_triangle_face_normals();
480 
481  // Iterate over operators
483 
485 }

◆ setGaussPtsThroughThickness()

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

Reimplemented in PostProcFatPrismOnRefinedMesh, and PrismFE.

Definition at line 45 of file FatPrismElementForcesAndSourcesCore.hpp.

45  {
47  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
49  }

◆ setGaussPtsTrianglesOnly()

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

Reimplemented in PostProcFatPrismOnRefinedMesh, and PrismFE.

Definition at line 39 of file FatPrismElementForcesAndSourcesCore.hpp.

39  {
41  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
43  }

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

Definition at line 206 of file FatPrismElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

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

Definition at line 186 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly
protected

Definition at line 190 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TrianglesOnly

EntitiesFieldData MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly
protected

Definition at line 193 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TroughThickness

EntitiesFieldData MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness
protected

Definition at line 194 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsThroughThickness

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
protected

◆ gaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
protected

◆ hoCoordsAtGaussPtsF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3
protected

Definition at line 196 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4
protected

Definition at line 200 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ normal

VectorDouble MoFEM::FatPrismElementForcesAndSourcesCore::normal
protected

Definition at line 187 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3
protected

Definition at line 197 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4
protected

Definition at line 201 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnPrism MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals
protected

Definition at line 204 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3
protected

Definition at line 198 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4
protected

Definition at line 202 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3
protected

Definition at line 199 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4
protected

Definition at line 203 of file FatPrismElementForcesAndSourcesCore.hpp.


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
diffN_MBEDGE0
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:107
MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly
EntitiesFieldData dataH1TrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:193
MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3
MatrixDouble hoCoordsAtGaussPtsF3
Definition: FatPrismElementForcesAndSourcesCore.hpp:196
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::VolumeElementForcesAndSourcesCore::meshPositionsFieldName
std::string meshPositionsFieldName
Definition: VolumeElementForcesAndSourcesCore.hpp:32
MoFEM::DataOperator::opRhs
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)
Definition: DataOperators.cpp:140
LASTBASE
@ LASTBASE
Definition: definitions.h:69
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4
MatrixDouble hoCoordsAtGaussPtsF4
Definition: FatPrismElementForcesAndSourcesCore.hpp:200
NOBASE
@ NOBASE
Definition: definitions.h:59
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:596
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
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
QUAD_2D_TABLE
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
MoFEM::FatPrismElementForcesAndSourcesCore::aRea
double aRea[2]
Definition: FatPrismElementForcesAndSourcesCore.hpp:186
MoFEM::CoreInterface::get_field_bit_number
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4
MatrixDouble tAngent2_at_GaussPtF4
Definition: FatPrismElementForcesAndSourcesCore.hpp:203
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
MatrixDouble gaussPtsTrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:189
N_MBEDGE0
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:105
QUAD_2D_TABLE_SIZE
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
MoFEM::EntitiesFieldData::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: EntitiesFieldData.hpp:50
MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly
MatrixDouble coordsAtGaussPtsTrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:190
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:177
FTensor::Tensor2< double, 3, 3 >
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::ForcesAndSourcesCore::getEntityFieldData
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
Definition: ForcesAndSourcesCore.cpp:770
MoFEM::VolumeElementForcesAndSourcesCore::coords
VectorDouble coords
Definition: VolumeElementForcesAndSourcesCore.hpp:89
MoFEM::BasicMethod::fieldsPtr
const Field_multiIndex * fieldsPtr
raw pointer to fields container
Definition: LoopMethods.hpp:252
MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore
VolumeElementForcesAndSourcesCore(Interface &m_field, const EntityType type=MBTET)
Definition: VolumeElementForcesAndSourcesCore.cpp:9
QUAD_1D_TABLE
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
MoFEM::FatPrismElementForcesAndSourcesCore::setGaussPtsThroughThickness
virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
Definition: FatPrismElementForcesAndSourcesCore.hpp:45
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::VolumeElementForcesAndSourcesCore::num_nodes
int num_nodes
Definition: VolumeElementForcesAndSourcesCore.hpp:100
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
N_MBEDGE1
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:106
MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
coordinated at gauss points
Definition: ForcesAndSourcesCore.hpp:544
MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness
EntitiesFieldData dataH1TroughThickness
Definition: FatPrismElementForcesAndSourcesCore.hpp:194
MoFEM::VolumeElementForcesAndSourcesCore::conn
const EntityHandle * conn
Definition: VolumeElementForcesAndSourcesCore.hpp:101
MoFEM::ForcesAndSourcesCore::createDataOnElement
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
Definition: ForcesAndSourcesCore.cpp:1305
QUAD_::npoints
int npoints
Definition: quad.h:29
MoFEM::FatPrismElementForcesAndSourcesCore::getRuleTrianglesOnly
virtual int getRuleTrianglesOnly(int order)
Definition: FatPrismElementForcesAndSourcesCore.hpp:36
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
MoFEM::Tools::diffShapeFunMBTRI
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:104
MoFEM::FatPrismElementForcesAndSourcesCore::setGaussPtsTrianglesOnly
virtual MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
Definition: FatPrismElementForcesAndSourcesCore.hpp:39
MoFEM::EntitiesFieldData::bAse
std::bitset< LASTBASE > bAse
bases on element
Definition: EntitiesFieldData.hpp:45
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
Definition: ForcesAndSourcesCore.cpp:990
FTensor::Index< 'i', 3 >
diffN_MBEDGE1
#define diffN_MBEDGE1
derivative of edge shape function
Definition: fem_tools.h:108
convert.n
n
Definition: convert.py:82
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
MatrixDouble gaussPtsThroughThickness
Definition: FatPrismElementForcesAndSourcesCore.hpp:191
MoFEM::VolumeElementForcesAndSourcesCore::vOlume
double & vOlume
Definition: VolumeElementForcesAndSourcesCore.hpp:98
MoFEM::ForcesAndSourcesCore::dataH1
EntitiesFieldData & dataH1
Definition: ForcesAndSourcesCore.hpp:471
MoFEM::determinantTensor3by3
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1511
QUAD_1D_TABLE_SIZE
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163
FTensor::dd
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
FTensor::Tensor0
Definition: Tensor0.hpp:16
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEM::FatPrismElementForcesAndSourcesCore::getRuleThroughThickness
virtual int getRuleThroughThickness(int order)
Definition: FatPrismElementForcesAndSourcesCore.hpp:37
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
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
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:56
MoFEM::ForcesAndSourcesCore::loopOverOperators
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
Definition: ForcesAndSourcesCore.cpp:1347
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::field_it
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
Definition: Projection10NodeCoordsOnField.cpp:123
MoFEM::ForcesAndSourcesCore::getNodesFieldData
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
Definition: ForcesAndSourcesCore.cpp:642
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
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MoFEM::Tools::getTriNormal
static MoFEMErrorCode getTriNormal(const double *coords, double *normal, double *d_normal=nullptr)
Get the Tri Normal objectGet triangle normal.
Definition: Tools.cpp:353
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::FatPrismElementForcesAndSourcesCore::normal
VectorDouble normal
Definition: FatPrismElementForcesAndSourcesCore.hpp:187
MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals
OpGetCoordsAndNormalsOnPrism opHOCoordsAndNormals
Definition: FatPrismElementForcesAndSourcesCore.hpp:204
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:24
MoFEM::OpGetCoordsAndNormalsOnPrism::calculateNormals
MoFEMErrorCode calculateNormals()
Definition: DataOperators.cpp:467
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182