v0.6.10
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 setIntegartionPts ()
 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 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>
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 >
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type >
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...
 
DEPRECATED MoFEMErrorCode set_snes (SNES snes)
 
DEPRECATED MoFEMErrorCode set_snes_ctx (const SNESContext &ctx)
 
DEPRECATED MoFEMErrorCode copy_snes (const SnesMethod &snes)
 
- 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...
 
DEPRECATED MoFEMErrorCode set_ts_ctx (const TSContext &ctx)
 
DEPRECATED MoFEMErrorCode copy_ts (const TSMethod &ts)
 
DEPRECATED MoFEMErrorCode set_ts (TS _ts)
 

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

Constructor & Destructor Documentation

◆ FatPrismElementForcesAndSourcesCore()

MoFEM::FatPrismElementForcesAndSourcesCore::FatPrismElementForcesAndSourcesCore ( Interface m_field)

Definition at line 65 of file FatPrismElementForcesAndSourcesCore.hpp.

65  :
66  VolumeElementForcesAndSourcesCore(m_field,MBPRISM),
67  dataH1TrianglesOnly(MBPRISM),
68  dataH1TroughThickness(MBPRISM),
72  ) {
73  }
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 79 of file FatPrismElementForcesAndSourcesCore.cpp.

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

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

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

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

◆ setGaussPtsThroughThickness()

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

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

Definition at line 84 of file FatPrismElementForcesAndSourcesCore.hpp.

84  {
86  SETERRQ(PETSC_COMM_SELF,MOFEM_NOT_IMPLEMENTED,"not implemented");
88  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:468
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:474

◆ setGaussPtsTrianglesOnly()

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

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

Definition at line 78 of file FatPrismElementForcesAndSourcesCore.hpp.

78  {
80  SETERRQ(PETSC_COMM_SELF,MOFEM_NOT_IMPLEMENTED,"not implemented");
82  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:468
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:474

Member Data Documentation

◆ aRea

double MoFEM::FatPrismElementForcesAndSourcesCore::aRea[2]

Definition at line 45 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::coordsAtGaussPtsTrianglesOnly

Definition at line 49 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TrianglesOnly

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TrianglesOnly

Definition at line 52 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ dataH1TroughThickness

DataForcesAndSourcesCore MoFEM::FatPrismElementForcesAndSourcesCore::dataH1TroughThickness

Definition at line 53 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsThroughThickness

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness

Definition at line 50 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ gaussPtsTrianglesOnly

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly

Definition at line 48 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF3

Definition at line 55 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPtsF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::hoCoordsAtGaussPtsF4

Definition at line 59 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ normal

VectorDouble MoFEM::FatPrismElementForcesAndSourcesCore::normal

Definition at line 46 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF3

Definition at line 56 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::nOrmals_at_GaussPtF4

Definition at line 60 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnPrism MoFEM::FatPrismElementForcesAndSourcesCore::opHOCoordsAndNormals

Definition at line 63 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF3

Definition at line 57 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent1_at_GaussPtF4

Definition at line 61 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF3

Definition at line 58 of file FatPrismElementForcesAndSourcesCore.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble MoFEM::FatPrismElementForcesAndSourcesCore::tAngent2_at_GaussPtF4

Definition at line 62 of file FatPrismElementForcesAndSourcesCore.hpp.


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