v0.8.23
Classes | Public Member Functions | Public Attributes | List of all members
PostProcFaceOnRefinedMesh Struct Reference

Postprocess on face. More...

#include <users_modules/basic_finite_elements/src/PostProcOnRefMesh.hpp>

Inheritance diagram for PostProcFaceOnRefinedMesh:
[legend]
Collaboration diagram for PostProcFaceOnRefinedMesh:
[legend]

Classes

struct  CommonData
 

Public Member Functions

 PostProcFaceOnRefinedMesh (MoFEM::Interface &m_field, bool six_node_post_proc_tris=true)
 
virtual ~PostProcFaceOnRefinedMesh ()
 
int getRule (int order)
 
MoFEMErrorCode generateReferenceElementMesh ()
 
MoFEMErrorCode setGaussPts (int order)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
virtual PostProcCommonOnRefMesh::CommonDatagetCommonData ()
 
- Public Member Functions inherited from PostProcTemplateOnRefineMesh< MoFEM::FaceElementForcesAndSourcesCore >
 PostProcTemplateOnRefineMesh (MoFEM::Interface &m_field)
 
MoFEMErrorCode addFieldValuesPostProc (const std::string field_name, Vec v=PETSC_NULL)
 Add operator to post-process L2, H1, Hdiv, Hcurl field value. More...
 
MoFEMErrorCode addFieldValuesPostProc (const std::string field_name, const std::string tag_name, Vec v=PETSC_NULL)
 Add operator to post-process L2 or H1 field value. More...
 
MoFEMErrorCode addFieldValuesGradientPostProc (const std::string field_name, Vec v=PETSC_NULL)
 Add operator to post-process L2 or H1 field gradient. More...
 
MoFEMErrorCode addFieldValuesGradientPostProc (const std::string field_name, const std::string tag_name, Vec v=PETSC_NULL)
 Add operator to post-process L2 or H1 field gradient. More...
 
MoFEMErrorCode writeFile (const std::string file_name)
 wrote results in (MOAB) format, use "file_name.h5m" More...
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore
 FaceElementForcesAndSourcesCore (Interface &m_field)
 
virtual MoFEMErrorCode calculateAreaAndNormal ()
 Calculate element area and normal of the face. More...
 
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
DEPRECATED MoFEMErrorCode setIntegartionPts ()
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode calculateHoNormal ()
 Calculate normal on curved elements. More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
virtual ~ForcesAndSourcesCore ()
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 
int getMaxDataOrder () const
 Get max order of approximation for data fields. More...
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows. More...
 
int getMaxColOrder () const
 Get max order of approximation for field in columns. More...
 
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (DataForcesAndSourcesCore &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getNodesIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEntityIndices (DataForcesAndSourcesCore &data, const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNodesFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 
MoFEMErrorCode getNoFieldFieldData (DataForcesAndSourcesCore &data, const boost::string_ref field_name) const
 
MoFEMErrorCode getNodesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &data) const
 Get field approximation space and base on entities. More...
 
MoFEMErrorCode getProblemNodesIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get indices of nodal indices which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemTypeIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const
 get indices by type (generic function) which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemNodesRowIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeRowIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
MoFEMErrorCode getProblemNodesColIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeColIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
boost::ptr_vector< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator. More...
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object. More...
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object. More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate base functions. More...
 
MoFEMErrorCode createDataOnElement ()
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

bool sixNodePostProcTris
 
std::map< EntityHandle, EntityHandleelementsMap
 
CommonData commonData
 
- Public Attributes inherited from PostProcTemplateOnRefineMesh< MoFEM::FaceElementForcesAndSourcesCore >
moab::Core coreMesh
 
moab::Interface & postProcMesh
 
std::vector< EntityHandlemapGaussPts
 
- Public Attributes inherited from MoFEM::FaceElementForcesAndSourcesCore
double aRea
 
int num_nodes
 
const EntityHandleconn
 
VectorDouble nOrmal
 
VectorDouble tangentOne
 
VectorDouble tangentTwo
 
VectorDouble coords
 
MatrixDouble coordsAtGaussPts
 
std::string meshPositionsFieldName
 Name of the field with geometry. More...
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble normalsAtGaussPts
 
MatrixDouble tangentOneAtGaussPts
 
MatrixDouble tangentTwoAtGaussPts
 
OpGetCoordsAndNormalsOnFace opHOCoordsAndNormals
 
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform
 
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
 
int nbGaussPts
 Number of integration points. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
const boost::shared_ptr< DataForcesAndSourcesCoredataOnElement [LASTSPACE]
 Entity data on element entity rows fields. More...
 
const boost::shared_ptr< DataForcesAndSourcesCorederivedDataOnElement [LASTSPACE]
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
RuleHookFun setRuleHook
 Set function to calculate integration rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
- 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

Postprocess on face.

Examples
cell_forces.cpp, minimal_surface_area.cpp, and reaction_diffusion_equation.cpp.

Definition at line 737 of file PostProcOnRefMesh.hpp.

Constructor & Destructor Documentation

◆ PostProcFaceOnRefinedMesh()

PostProcFaceOnRefinedMesh::PostProcFaceOnRefinedMesh ( MoFEM::Interface m_field,
bool  six_node_post_proc_tris = true 
)

◆ ~PostProcFaceOnRefinedMesh()

virtual PostProcFaceOnRefinedMesh::~PostProcFaceOnRefinedMesh ( )
virtual

Definition at line 748 of file PostProcOnRefMesh.hpp.

748  {
749  ParallelComm *pcomm_post_proc_mesh =
750  ParallelComm::get_pcomm(&postProcMesh, MYPCOMM_INDEX);
751  if (pcomm_post_proc_mesh != NULL) {
752  delete pcomm_post_proc_mesh;
753  }
754  }
moab::Interface & postProcMesh
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:284

Member Function Documentation

◆ generateReferenceElementMesh()

MoFEMErrorCode PostProcFaceOnRefinedMesh::generateReferenceElementMesh ( )
Examples
minimal_surface_area.cpp.

Definition at line 532 of file PostProcOnRefMesh.cpp.

532  {
534 
535  gaussPts.resize(3, 3, false);
536  gaussPts.clear();
537  gaussPts(0, 0) = 0;
538  gaussPts(1, 0) = 0;
539  gaussPts(0, 1) = 1;
540  gaussPts(1, 1) = 0;
541  gaussPts(0, 2) = 0;
542  gaussPts(1, 2) = 1;
543  mapGaussPts.resize(gaussPts.size2());
544 
545  moab::Core core_ref;
546  moab::Interface &moab_ref = core_ref;
547  const EntityHandle *conn;
548  int num_nodes;
549  EntityHandle tri_conn[3];
550  MatrixDouble coords(6, 3);
551  for (int gg = 0; gg != 3; gg++) {
552  coords(gg, 0) = gaussPts(0, gg);
553  coords(gg, 1) = gaussPts(1, gg);
554  coords(gg, 2) = 0;
555  CHKERR moab_ref.create_vertex(&coords(gg, 0), tri_conn[gg]);
556  }
557 
558  EntityHandle tri;
559  CHKERR moab_ref.create_element(MBTRI, tri_conn, 3, tri);
560  Range edges;
561  CHKERR moab_ref.get_adjacencies(&tri, 1, 1, true, edges);
562  EntityHandle meshset;
563  CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
564  CHKERR moab_ref.add_entities(meshset, &tri, 1);
565  CHKERR moab_ref.add_entities(meshset, edges);
566  if (sixNodePostProcTris) {
567  CHKERR moab_ref.convert_entities(meshset, true, false, false);
568  }
569  CHKERR moab_ref.get_connectivity(tri, conn, num_nodes, false);
570  CHKERR moab_ref.get_coords(conn, num_nodes, &coords(0, 0));
571 
572  gaussPts.resize(3, num_nodes, false);
573  gaussPts.clear();
574  for (int nn = 0; nn < 3; nn++) {
575  gaussPts(0, nn) = coords(nn, 0);
576  gaussPts(1, nn) = coords(nn, 1);
577  gaussPts(0, 3 + nn) = coords(3 + nn, 0);
578  gaussPts(1, 3 + nn) = coords(3 + nn, 1);
579  }
580 
582 }
MatrixDouble gaussPts
Matrix of integration points.
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:77
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
bool sixNodePostProcTris
#define CHKERR
Inline error check.
Definition: definitions.h:595
std::vector< EntityHandle > mapGaussPts
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ getCommonData()

virtual PostProcCommonOnRefMesh::CommonData& PostProcFaceOnRefinedMesh::getCommonData ( )
virtual

Reimplemented from PostProcTemplateOnRefineMesh< MoFEM::FaceElementForcesAndSourcesCore >.

Definition at line 770 of file PostProcOnRefMesh.hpp.

770  {
771  return commonData;
772  }
CommonData commonData

◆ getRule()

int PostProcFaceOnRefinedMesh::getRule ( int  order)
virtual
Deprecated:
Use getRule(int row_order, int col_order, int data order)

Reimplemented from MoFEM::ForcesAndSourcesCore.

Definition at line 757 of file PostProcOnRefMesh.hpp.

757 { return -1; };

◆ postProcess()

MoFEMErrorCode PostProcFaceOnRefinedMesh::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.

Definition at line 657 of file PostProcOnRefMesh.cpp.

657  {
659  ParallelComm *pcomm =
660  ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
661  ParallelComm *pcomm_post_proc_mesh =
662  ParallelComm::get_pcomm(&postProcMesh, MYPCOMM_INDEX);
663  if (pcomm_post_proc_mesh == NULL) {
664  pcomm_post_proc_mesh = new ParallelComm(&postProcMesh, mField.get_comm());
665  }
666  Range tris;
667  CHKERR postProcMesh.get_entities_by_type(0, MBTRI, tris, false);
668  int rank = pcomm->rank();
669  Range::iterator pit = tris.begin();
670  for (; pit != tris.end(); pit++) {
671  CHKERR postProcMesh.tag_set_data(pcomm_post_proc_mesh->part_tag(), &*pit, 1,
672  &rank);
673  }
674  CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
676 }
virtual moab::Interface & get_moab()=0
moab::Interface & postProcMesh
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:284
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
virtual MPI_Comm & get_comm() const =0

◆ preProcess()

MoFEMErrorCode PostProcFaceOnRefinedMesh::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.

Definition at line 647 of file PostProcOnRefMesh.cpp.

647  {
649  ParallelComm *pcomm_post_proc_mesh =
650  ParallelComm::get_pcomm(&postProcMesh, MYPCOMM_INDEX);
651  if (pcomm_post_proc_mesh != NULL) {
652  delete pcomm_post_proc_mesh;
653  }
655 }
moab::Interface & postProcMesh
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:284

◆ setGaussPts()

MoFEMErrorCode PostProcFaceOnRefinedMesh::setGaussPts ( int  order)
virtual
Deprecated:
setGaussPts(int row_order, int col_order, int data order);

Reimplemented from MoFEM::ForcesAndSourcesCore.

Definition at line 584 of file PostProcOnRefMesh.cpp.

584  {
586  if (gaussPts.size1() == 0) {
587  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
588  "post-process mesh not generated");
589  }
590 
591  const EntityHandle *conn;
592  int num_nodes;
593  EntityHandle tri;
594 
595  if (elementsMap.find(numeredEntFiniteElementPtr->getEnt()) !=
596  elementsMap.end()) {
597  tri = elementsMap[numeredEntFiniteElementPtr->getEnt()];
598  } else {
599  ublas::vector<EntityHandle> tri_conn(3);
600  MatrixDouble coords_tri(3, 3);
601  VectorDouble coords(3);
602  CHKERR mField.get_moab().get_connectivity(
603  numeredEntFiniteElementPtr->getEnt(), conn, num_nodes, true);
604  CHKERR mField.get_moab().get_coords(conn, num_nodes, &coords_tri(0, 0));
605  for (int gg = 0; gg != 3; gg++) {
606  double ksi = gaussPts(0, gg);
607  double eta = gaussPts(1, gg);
608  double n0 = N_MBTRI0(ksi, eta);
609  double n1 = N_MBTRI1(ksi, eta);
610  double n2 = N_MBTRI2(ksi, eta);
611  double x =
612  n0 * coords_tri(0, 0) + n1 * coords_tri(1, 0) + n2 * coords_tri(2, 0);
613  double y =
614  n0 * coords_tri(0, 1) + n1 * coords_tri(1, 1) + n2 * coords_tri(2, 1);
615  coords[0] = x;
616  coords[1] = y;
617  coords[2] = 0;
618  CHKERR postProcMesh.create_vertex(&coords[0], tri_conn[gg]);
619  }
620  CHKERR postProcMesh.create_element(MBTRI, &tri_conn[0], 3, tri);
621  elementsMap[numeredEntFiniteElementPtr->getEnt()] = tri;
622  Range edges;
623  CHKERR postProcMesh.get_adjacencies(&tri, 1, 1, true, edges);
624  EntityHandle meshset;
625  CHKERR postProcMesh.create_meshset(MESHSET_SET, meshset);
626  CHKERR postProcMesh.add_entities(meshset, &tri, 1);
627  CHKERR postProcMesh.add_entities(meshset, edges);
628  if (sixNodePostProcTris) {
629  CHKERR postProcMesh.convert_entities(meshset, true, false, false);
630  }
631  CHKERR postProcMesh.delete_entities(&meshset, 1);
632  CHKERR postProcMesh.delete_entities(edges);
633  }
634 
635  // Set values which map nodes with integration points on the prism
636  {
637  CHKERR postProcMesh.get_connectivity(tri, conn, num_nodes, false);
638  mapGaussPts.resize(num_nodes);
639  for (int nn = 0; nn < num_nodes; nn++) {
640  mapGaussPts[nn] = conn[nn];
641  }
642  }
643 
645 }
MatrixDouble gaussPts
Matrix of integration points.
virtual moab::Interface & get_moab()=0
moab::Interface & postProcMesh
std::map< EntityHandle, EntityHandle > elementsMap
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:77
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
bool sixNodePostProcTris
#define CHKERR
Inline error check.
Definition: definitions.h:595
std::vector< EntityHandle > mapGaussPts
#define N_MBTRI2(x, y)
triangle shape function
Definition: fem_tools.h:56
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:76
#define N_MBTRI1(x, y)
triangle shape function
Definition: fem_tools.h:55
#define N_MBTRI0(x, y)
triangle shape function
Definition: fem_tools.h:54

Member Data Documentation

◆ commonData

CommonData PostProcFaceOnRefinedMesh::commonData

Definition at line 768 of file PostProcOnRefMesh.hpp.

◆ elementsMap

std::map<EntityHandle, EntityHandle> PostProcFaceOnRefinedMesh::elementsMap

Definition at line 762 of file PostProcOnRefMesh.hpp.

◆ sixNodePostProcTris

bool PostProcFaceOnRefinedMesh::sixNodePostProcTris

Definition at line 740 of file PostProcOnRefMesh.hpp.


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