v0.8.13
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 calculateBaseFunctionsOnElement (const int b)
 Calculate base functions. More...
 
virtual MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate base functions. More...
 
virtual MoFEMErrorCode transformBaseFunctions ()
 Transform base functions based on geometric element Jacobian. More...
 
virtual MoFEMErrorCode calculateHoJacobian ()
 Calculate Jacobian for HO geometry. More...
 
virtual MoFEMErrorCode transformHoBaseFunctions ()
 Transform base functions based on ho-geometry element Jacobian. More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
virtual ~ForcesAndSourcesCore ()
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 
int getMaxDataOrder () const
 Get max order of approximation for data fields. More...
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows. More...
 
int getMaxColOrder () const
 Get max order of approximation for field in columns. More...
 
MoFEMErrorCode getSense (EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order of approximation on the entity More...
 
MoFEMErrorCode getDataOrderSpaceAndBase (const std::string &field_name, const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order on entity More...
 
MoFEMErrorCode getEdgesSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getTrisSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getQuadSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getEdgesDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getTrisDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getQuadDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getTetDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getPrismDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getEdgesDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTrisDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getQuadDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTetDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getPrismDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getNodesIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
 get node indices More...
 
MoFEMErrorCode getTypeIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices, VectorInt &local_indices) const
 get indices by type (generic function) More...
 
MoFEMErrorCode getTypeIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get indices by type (generic function) 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 getEdgesRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Edges row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEdgesColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Edges col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTrisRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tris row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTrisColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tris col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTetsRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tets row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTetsColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tets col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getQuadRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Quad row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getQuadColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Quad col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getPrismRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Prism row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getPrismColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Prism col indices from FENumeredDofEntity_multiIndex More...
 
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 getTypeFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 Get field data on entities. More...
 
MoFEMErrorCode getTypeFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 
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 getEdgesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTrisFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getQuadFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTetsFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getPrismFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) 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...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

double aRea [2]
 
VectorDouble normal
 
MatrixDouble gaussPtsTrianglesOnly
 
MatrixDouble coordsAtGaussPtsTrianglesOnly
 
MatrixDouble gaussPtsThroughThickness
 
DataForcesAndSourcesCore dataH1TrianglesOnly
 
DataForcesAndSourcesCore dataH1TroughThickness
 
MatrixDouble hoCoordsAtGaussPtsF3
 
MatrixDouble nOrmals_at_GaussPtF3
 
MatrixDouble tAngent1_at_GaussPtF3
 
MatrixDouble tAngent2_at_GaussPtF3
 
MatrixDouble hoCoordsAtGaussPtsF4
 
MatrixDouble nOrmals_at_GaussPtF4
 
MatrixDouble tAngent1_at_GaussPtF4
 
MatrixDouble tAngent2_at_GaussPtF4
 
OpGetCoordsAndNormalsOnPrism opHOCoordsAndNormals
 
- Public Attributes inherited from MoFEM::VolumeElementForcesAndSourcesCore
VectorDouble coords
 
MatrixDouble3by3 jAc
 
MatrixDouble3by3 invJac
 
DataForcesAndSourcesCore dataH1
 
DerivedDataForcesAndSourcesCore derivedDataH1
 
DataForcesAndSourcesCore dataL2
 
DerivedDataForcesAndSourcesCore derivedDataL2
 
DataForcesAndSourcesCore dataHdiv
 
DerivedDataForcesAndSourcesCore derivedDataHdiv
 
DataForcesAndSourcesCore dataHcurl
 
DerivedDataForcesAndSourcesCore derivedDataHcurl
 
DataForcesAndSourcesCore dataNoField
 
DataForcesAndSourcesCore dataNoFieldCol
 
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 gaussPts
 
MatrixDouble coordsAtGaussPts
 
unsigned int nbGaussPts
 Number of integration points. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 
int loopSize
 
int rAnk
 
int sIze
 
const RefEntity_multiIndexrefinedEntitiesPtr
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 
const ProblemproblemPtr
 
const Field_multiIndexfieldsPtr
 
const FieldEntity_multiIndexentitiesPtr
 
const DofEntity_multiIndexdofsPtr
 
const FiniteElement_multiIndexfiniteElementsPtr
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 
boost::function< MoFEMErrorCode()> postProcessHook
 
boost::function< MoFEMErrorCode()> operatorHook
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 
Vec snes_x
 
Vec snes_f
 
Mat snes_A
 
Mat snes_B
 
- Public Attributes inherited from MoFEM::TSMethod
TSContext ts_ctx
 
TS ts
 
Vec ts_u
 
Vec ts_u_t
 
Vec ts_F
 
Mat ts_A
 
Mat ts_B
 
PetscInt ts_step
 
PetscReal ts_a
 
PetscReal ts_t
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION, CTX_OPERATORS, CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION, CTX_SNESSETJACOBIAN, CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION, CTX_TSSETRHSJACOBIAN, CTX_TSSETIFUNCTION, CTX_TSSETIJACOBIAN,
  CTX_TSTSMONITORSET, CTX_TSNONE
}
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

FatPrism finite element

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

Todo:
Need to implement operators that will make this element work as Volume element

Definition at line 43 of file FatPrismElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FatPrismElementForcesAndSourcesCore()

MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore ( Interface m_field)

Definition at line 66 of file FatPrismElementForcesAndSourcesCore.hpp.

67  : VolumeElementForcesAndSourcesCore(m_field, MBPRISM),
68  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 78 of file FatPrismElementForcesAndSourcesCore.cpp.

78  {
80 
81  if (numeredEntFiniteElementPtr->getEntType() != MBPRISM)
83 
85  int num_nodes;
86  const EntityHandle *conn;
87  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
88  {
89  coords.resize(num_nodes * 3, false);
90  CHKERR mField.get_moab().get_coords(conn, num_nodes,
91  &*coords.data().begin());
92  double diff_n[6];
93  CHKERR ShapeDiffMBTRI(diff_n);
94  normal.resize(6, false);
95  CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[0], &normal[0]);
96  CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[9], &normal[3]);
97  aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
98  aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
99  }
100 
104  // H1
105  if ((dataH1.spacesOnEntities[MBEDGE]).test(H1)) {
113  // Triangles only
118  // Through thickness
121  }
122  // Hdiv
123  if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV)) {
124  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
125  }
126  // Hcurl
127  if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL)) {
128  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
129  }
130  // L2
131  if ((dataH1.spacesOnEntities[MBPRISM]).test(L2)) {
132  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
133  }
134 
135  // get approx. on triangles, i.e. faces 3 and 4
136  int nb_gauss_pts_on_faces;
137  try {
138  int order_triangles_only = 1;
139  int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
140  for (unsigned int ee = 0; ee < 9; ee++) {
141  if (!valid_edges[ee])
142  continue;
143  order_triangles_only = std::max(
144  order_triangles_only,
145  dataH1TrianglesOnly.dataOnEntities[MBEDGE][ee].getDataOrder());
146  }
147  for (unsigned int ff = 3; ff <= 4; ff++) {
148  order_triangles_only = std::max(
149  order_triangles_only,
150  dataH1TrianglesOnly.dataOnEntities[MBTRI][ff].getDataOrder());
151  }
152  for (unsigned int qq = 0; qq < 3; qq++) {
153  order_triangles_only = std::max(
154  order_triangles_only,
155  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
156  }
157  order_triangles_only = std::max(
158  order_triangles_only,
159  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
160  // integration pts on the triangles surfaces
161  int rule = getRuleTrianglesOnly(order_triangles_only);
162  if (rule >= 0) {
163  if (rule < QUAD_2D_TABLE_SIZE) {
164  if (QUAD_2D_TABLE[rule]->dim != 2) {
165  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
166  }
167  if (QUAD_2D_TABLE[rule]->order < rule) {
168  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
169  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
170  }
171  nb_gauss_pts_on_faces = QUAD_2D_TABLE[rule]->npoints;
172  gaussPtsTrianglesOnly.resize(3, nb_gauss_pts_on_faces, false);
173  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[1], 3,
174  &gaussPtsTrianglesOnly(0, 0), 1);
175  cblas_dcopy(nb_gauss_pts_on_faces, &QUAD_2D_TABLE[rule]->points[2], 3,
176  &gaussPtsTrianglesOnly(1, 0), 1);
177  cblas_dcopy(nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->weights, 1,
178  &gaussPtsTrianglesOnly(2, 0), 1);
179  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
180  nb_gauss_pts_on_faces, 3, false);
181  double *shape_ptr = &*dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0]
182  .getN(NOBASE)
183  .data()
184  .begin();
185  cblas_dcopy(3 * nb_gauss_pts_on_faces, QUAD_2D_TABLE[rule]->points, 1,
186  shape_ptr, 1);
187  } else {
188  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
189  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
190  nb_gauss_pts_on_faces = 0;
191  }
192  } else {
193  CHKERR setGaussPtsTrianglesOnly(order_triangles_only);
194  nb_gauss_pts_on_faces = gaussPtsTrianglesOnly.size2();
195  if (nb_gauss_pts_on_faces == 0)
197  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(
198  nb_gauss_pts_on_faces, 3, false);
199  if (nb_gauss_pts_on_faces) {
201  .getN(NOBASE)
202  .data()
203  .begin(),
204  &gaussPtsTrianglesOnly(0, 0),
205  &gaussPtsTrianglesOnly(1, 0), nb_gauss_pts_on_faces);
206  }
207  }
208  }
209  CATCH_ERRORS;
210 
211  // approx. trough prism thickness
212  int nb_gauss_pts_through_thickness;
213  try {
214  int order_thickness = 1;
215  for (unsigned int ee = 3; ee <= 5; ee++) {
216  order_thickness = std::max(
217  order_thickness,
218  dataH1TroughThickness.dataOnEntities[MBEDGE][ee].getDataOrder());
219  }
220  for (unsigned int qq = 0; qq < 3; qq++) {
221  order_thickness = std::max(
222  order_thickness,
223  dataH1TroughThickness.dataOnEntities[MBQUAD][qq].getDataOrder());
224  }
225  order_thickness = std::max(
226  order_thickness,
227  dataH1TroughThickness.dataOnEntities[MBPRISM][0].getDataOrder());
228  // integration points
229  int rule = getRuleThroughThickness(order_thickness);
230  if (rule >= 0) {
231  if (rule < QUAD_1D_TABLE_SIZE) {
232  if (QUAD_1D_TABLE[rule]->dim != 1) {
233  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
234  }
235  if (QUAD_1D_TABLE[rule]->order < rule) {
236  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
237  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order, rule);
238  }
239  nb_gauss_pts_through_thickness = QUAD_1D_TABLE[rule]->npoints;
240  gaussPtsThroughThickness.resize(2, nb_gauss_pts_through_thickness,
241  false);
242  cblas_dcopy(nb_gauss_pts_through_thickness,
243  &QUAD_1D_TABLE[rule]->points[1], 2,
244  &gaussPtsThroughThickness(0, 0), 1);
245  cblas_dcopy(nb_gauss_pts_through_thickness,
246  QUAD_1D_TABLE[rule]->weights, 1,
247  &gaussPtsThroughThickness(1, 0), 1);
248  } else {
249  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
250  "rule > quadrature order %d < %d", rule, QUAD_1D_TABLE_SIZE);
251  nb_gauss_pts_through_thickness = 0;
252  }
253  } else {
254  CHKERR setGaussPtsThroughThickness(order_thickness);
255  nb_gauss_pts_through_thickness = gaussPtsThroughThickness.size2();
256  }
257  if (nb_gauss_pts_through_thickness == 0)
259  }
260  CATCH_ERRORS;
261 
262  // Generate integration pts.
263  int nb_gauss_pts = nb_gauss_pts_on_faces * nb_gauss_pts_through_thickness;
264  gaussPts.resize(4, nb_gauss_pts, false);
265  {
266  int gg = 0;
267  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
268  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
269  gaussPts(0, gg) = gaussPtsTrianglesOnly(0, ggf);
270  gaussPts(1, gg) = gaussPtsTrianglesOnly(1, ggf);
271  gaussPts(2, gg) = gaussPtsThroughThickness(0, ggt);
272  gaussPts(3, gg) =
274  }
275  }
276  }
277 
278  {
279  coordsAtGaussPtsTrianglesOnly.resize(nb_gauss_pts_on_faces, 6, false);
280  for (int gg = 0; gg < nb_gauss_pts_on_faces; gg++) {
281  for (int dd = 0; dd < 3; dd++) {
283  cblas_ddot(3,
284  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
285  NOBASE)(gg, 0),
286  1, &coords[dd], 3);
288  cblas_ddot(3,
289  &dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
290  NOBASE)(gg, 0),
291  1, &coords[9 + dd], 3);
292  }
293  }
294  // linear for xi,eta and zeta
295  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(
296  1, 6, false);
298  .getDiffN(NOBASE)
299  .data()
300  .begin());
301 
302  // Calculate "nobase" base functions on prism, this is cartesian product
303  // of base functions on triangles with base functions through thickness
304  // FIXME: This could be effectively implemented with tensors
305  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 6,
306  false);
307  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(nb_gauss_pts, 18,
308  false);
309  for (int dd = 0; dd != 6; dd++) {
310  int gg = 0;
311  for (int ggf = 0; ggf < nb_gauss_pts_on_faces; ggf++) {
312  int ddd = dd > 2 ? dd - 3 : dd;
313  double tri_n = dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getN(
314  NOBASE)(ggf, ddd);
315  double dksi_tri_n =
316  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
317  0, 2 * ddd + 0);
318  double deta_tri_n =
319  dataH1TrianglesOnly.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(
320  0, 2 * ddd + 1);
321  for (int ggt = 0; ggt < nb_gauss_pts_through_thickness; ggt++, gg++) {
322  double zeta = gaussPtsThroughThickness(0, ggt);
323  double dzeta, edge_shape;
324  if (dd < 3) {
325  dzeta = diffN_MBEDGE0;
326  edge_shape = N_MBEDGE0(zeta);
327  } else {
328  dzeta = diffN_MBEDGE1;
329  edge_shape = N_MBEDGE1(zeta);
330  }
331  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, dd) =
332  tri_n * edge_shape;
333  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 0) =
334  dksi_tri_n * edge_shape;
335  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 1) =
336  deta_tri_n * edge_shape;
337  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE)(gg, 3 * dd + 2) =
338  tri_n * dzeta;
339  }
340  }
341  }
342 
343  /// Calculate coordinates at integration points
344  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
345  for (int gg = 0; gg < nb_gauss_pts; gg++) {
346  for (int dd = 0; dd < 3; dd++) {
348  6, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
349  &coords[dd], 3);
350  }
351  }
352  }
353 
354  // Calculate base functions on prism
355  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
356  if (dataH1.bAse.test(b)) {
357  switch (ApproximationBaseArray[b]) {
360  if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
361  CHKERR FatPrismPolynomialBase().getValue(
362  gaussPts,
363  boost::shared_ptr<BaseFunctionCtx>(new FatPrismPolynomialBaseCtx(
367  ApproximationBaseArray[b], NOBASE)));
368  }
369  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
370  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
371  "Not yet implemented");
372  }
373  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
374  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
375  "Not yet implemented");
376  }
377  if (dataH1.spacesOnEntities[MBTET].test(L2)) {
378  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
379  "Not yet implemented");
380  }
381  break;
382  default:
383  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
384  "Not yet implemented");
385  }
386  }
387  }
388 
389  // Calculate ho-geometry tangents and normals
390 
391  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
392  dataPtr->get<FieldName_mi_tag>().end()) {
393  hoCoordsAtGaussPtsF3.resize(nb_gauss_pts_on_faces, 3, false);
394  nOrmals_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
395  tAngent1_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
396  tAngent2_at_GaussPtF3.resize(nb_gauss_pts_on_faces, 3, false);
397  hoCoordsAtGaussPtsF4.resize(nb_gauss_pts_on_faces, 3, false);
398  nOrmals_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
399  tAngent1_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
400  tAngent2_at_GaussPtF4.resize(nb_gauss_pts_on_faces, 3, false);
410  } else {
411  hoCoordsAtGaussPtsF3.resize(0, 0, false);
412  nOrmals_at_GaussPtF3.resize(0, 0, false);
413  tAngent1_at_GaussPtF3.resize(0, 0, false);
414  tAngent2_at_GaussPtF3.resize(0, 0, false);
415  hoCoordsAtGaussPtsF4.resize(0, 0, false);
416  nOrmals_at_GaussPtF4.resize(0, 0, false);
417  tAngent1_at_GaussPtF4.resize(0, 0, false);
418  tAngent2_at_GaussPtF4.resize(0, 0, false);
419  }
420 
423  std::vector<std::string> last_eval_field_name(2);
424  DataForcesAndSourcesCore *op_data[2];
425  FieldSpace space[2];
426  FieldApproximationBase base[2];
427 
428  boost::ptr_vector<UserDataOperator>::iterator oit, hi_oit;
429  oit = opPtrVector.begin();
430  hi_oit = opPtrVector.end();
431 
432  // Run element operators
433  for (; oit != hi_oit; oit++) {
434 
435  oit->setPtrFE(this);
436 
437  if (oit->sPace != LASTSPACE) {
438 
439  // Set field
440  switch (oit->sPace) {
441  case NOSPACE:
442  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
443  case H1:
444  op_data[0] = &dataH1;
445  break;
446  case HCURL:
447  op_data[0] = &dataHcurl;
448  break;
449  case HDIV:
450  op_data[0] = &dataHdiv;
451  break;
452  case L2:
453  op_data[0] = &dataL2;
454  break;
455  case NOFIELD:
456  op_data[0] = &dataNoField;
457  break;
458  case LASTSPACE:
459  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
460  }
461 
462  // Reseat all data which all field dependent
463  op_data[0]->resetFieldDependentData();
464 
465  // Run operator
466  CHKERR oit->opRhs(*op_data[0], oit->doVertices, oit->doEdges,
467  oit->doQuads, oit->doTris, false, false);
468 
469  } else {
470 
471  for (int ss = 0; ss != 2; ss++) {
472 
473  std::string field_name = !ss ? oit->rowFieldName : oit->colFieldName;
474  const Field *field_struture = mField.get_field_structure(field_name);
475  BitFieldId data_id = field_struture->getId();
476 
477  if ((oit->getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
478  data_id)
479  .none()) {
480  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
481  "no data field < %s > on finite element < %s >",
482  field_name.c_str(), feName.c_str());
483  }
484 
485  if (oit->getOpType() & types[ss] ||
486  oit->getOpType() & UserDataOperator::OPROWCOL) {
487 
488  space[ss] = field_struture->getSpace();
489  switch (space[ss]) {
490  case NOSPACE:
491  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
492  case H1:
493  op_data[ss] = !ss ? &dataH1 : &derivedDataH1;
494  break;
495  case HCURL:
496  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
497  "not implemented yet");
498  break;
499  case HDIV:
500  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
501  "not implemented yet");
502  break;
503  case L2:
504  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
505  "not implemented yet");
506  break;
507  case NOFIELD:
508  op_data[ss] = !ss ? &dataNoField : &dataNoFieldCol;
509  break;
510  case LASTSPACE:
511  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
512  break;
513  }
514 
515  base[ss] = field_struture->getApproxBase();
516  switch (base[ss]) {
519  break;
520  default:
521  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
522  "unknown or not implemented base");
523  break;
524  }
525 
526  if (last_eval_field_name[ss] != field_name) {
527 
528  switch (space[ss]) {
529  case NOSPACE:
530  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
531  "unknown space");
532  case H1:
533  if (!ss) {
534  CHKERR getRowNodesIndices(*op_data[ss], field_name);
535  } else {
536  CHKERR getColNodesIndices(*op_data[ss], field_name);
537  }
538  CHKERR getNodesFieldData(*op_data[ss], field_name);
539  case HCURL:
540  if (!ss) {
541  CHKERR getEdgesRowIndices(*op_data[ss], field_name);
542  } else {
543  CHKERR getEdgesColIndices(*op_data[ss], field_name);
544  }
545  CHKERR getEdgesDataOrderSpaceAndBase(*op_data[ss], field_name);
546  CHKERR getEdgesFieldData(*op_data[ss], field_name);
547  case HDIV:
548  if (!ss) {
549  CHKERR getTrisRowIndices(*op_data[ss], field_name);
550  } else {
551  CHKERR getTrisColIndices(*op_data[ss], field_name);
552  }
553  CHKERR getTrisDataOrderSpaceAndBase(*op_data[ss], field_name);
554  CHKERR getTrisFieldData(*op_data[ss], field_name);
555  if (!ss) {
556  CHKERR getQuadRowIndices(*op_data[ss], field_name);
557  } else {
558  CHKERR getQuadColIndices(*op_data[ss], field_name);
559  }
560  CHKERR getQuadDataOrderSpaceAndBase(*op_data[ss], field_name);
561  CHKERR getQuadFieldData(*op_data[ss], field_name);
562  case L2:
563  if (!ss) {
564  CHKERR getPrismRowIndices(*op_data[ss], field_name);
565  } else {
566  CHKERR getPrismColIndices(*op_data[ss], field_name);
567  }
568  CHKERR getPrismDataOrderSpaceAndBase(*op_data[ss], field_name);
569  CHKERR getPrismFieldData(*op_data[ss], field_name);
570  break;
571  case NOFIELD:
572  if (!getNinTheLoop()) {
573  // NOFIELD data are the same for each element, can be retrieved
574  // only once
575  if (!ss) {
576  CHKERR getNoFieldRowIndices(*op_data[ss], field_name);
577  } else {
578  CHKERR getNoFieldColIndices(*op_data[ss], field_name);
579  }
580  CHKERR getNoFieldFieldData(*op_data[ss], field_name);
581  }
582  break;
583  case LASTSPACE:
584  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
585  "unknown space");
586  break;
587  }
588  last_eval_field_name[ss] = field_name;
589  }
590  }
591  }
592 
593  if (oit->getOpType() & UserDataOperator::OPROW) {
594  try {
595  CHKERR oit->opRhs(*op_data[0], oit->doVertices, oit->doEdges,
596  oit->doQuads, oit->doTris, false, oit->doPrisms);
597  } catch (std::exception &ex) {
598  std::ostringstream ss;
599  ss << "thorw in method: " << ex.what() << " at line " << __LINE__
600  << " in file " << __FILE__;
601  ss << " operator on row field name " << oit->rowFieldName;
602  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
603  }
604  }
605 
606  if (oit->getOpType() & UserDataOperator::OPCOL) {
607  try {
608  CHKERR oit->opRhs(*op_data[1], oit->doVertices, oit->doEdges,
609  oit->doQuads, oit->doTris, false, oit->doPrisms);
610  } catch (std::exception &ex) {
611  std::ostringstream ss;
612  ss << "thorw in method: " << ex.what() << " at line " << __LINE__
613  << " in file " << __FILE__;
614  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
615  }
616  }
617 
618  if (oit->getOpType() & UserDataOperator::OPROWCOL) {
619  try {
620  CHKERR oit->opLhs(*op_data[0], *op_data[1], oit->sYmm);
621  } catch (std::exception &ex) {
622  std::ostringstream ss;
623  ss << "thorw in method: " << ex.what() << " at line " << __LINE__
624  << " in file " << __FILE__;
625  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
626  }
627  }
628  }
629  }
630 
632 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
MoFEMErrorCode getTrisSense(DataForcesAndSourcesCore &data) const
boost::ptr_vector< UserDataOperator > opPtrVector
field with continuous normal traction
Definition: definitions.h:181
MoFEMErrorCode getNoFieldColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
MoFEMErrorCode getPrismRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Prism row indices from FENumeredDofEntity_multiIndex
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
MoFEMErrorCode getTrisDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
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 getPrismDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
MoFEMErrorCode getEdgesRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Edges row indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getColNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getPrismColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Prism col indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
MoFEMErrorCode getTrisRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Tris row indices from FENumeredDofEntity_multiIndex
scalar or vector of scalars describe (no true field)
Definition: definitions.h:178
MoFEMErrorCode getPrismFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
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 getTrisDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
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:495
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:80
std::string feName
std::bitset< BITFIELDID_SIZE > BitFieldId
Definition: Common.hpp:149
int getNinTheLoop() const
get number of evaluated element in the loop
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
MoFEMErrorCode getNoFieldRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
MoFEMErrorCode getQuadDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
MoFEMErrorCode getPrismDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
OpType
Controls loop over entities on element.
MoFEMErrorCode getNoFieldFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:183
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
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MoFEMErrorCode getQuadFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getQuadRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Quad row indices from FENumeredDofEntity_multiIndex
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:79
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:186
FieldApproximationBase
approximation base
Definition: definitions.h:140
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:81
MoFEMErrorCode getEdgesSense(DataForcesAndSourcesCore &data) const
MoFEMErrorCode getQuadColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Quad col indices from FENumeredDofEntity_multiIndex
field with continuous tangents
Definition: definitions.h:180
FieldSpace
approximation spaces
Definition: definitions.h:176
virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
#define CHKERR
Inline error check.
Definition: definitions.h:614
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 getTrisFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getEdgesColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Edges col indices from FENumeredDofEntity_multiIndex
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:175
MoFEMErrorCode getTrisColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Tris col indices from FENumeredDofEntity_multiIndex
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:439
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:465
continuous field
Definition: definitions.h:179
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getRowNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get row node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getQuadSense(DataForcesAndSourcesCore &data) const
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
std::bitset< LASTBASE > bAse
bases on element
MoFEMErrorCode getQuadDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
field with C-1 continuity
Definition: definitions.h:182
#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 354 of file FatPrismElementForcesAndSourcesCore.hpp.

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

◆ 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 349 of file FatPrismElementForcesAndSourcesCore.hpp.

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

◆ setGaussPtsThroughThickness()

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

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

Definition at line 83 of file FatPrismElementForcesAndSourcesCore.hpp.

83  {
85  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
87  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:519
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526

◆ setGaussPtsTrianglesOnly()

virtual MoFEMErrorCode MoFEM::FatPrismElementForcesAndSourcesCore::setGaussPtsTrianglesOnly ( int  order_triangles_only)
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:519
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526

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: