v0.14.0
Public Member Functions | List of all members
MoFEM::PostProcGenerateRefMesh< MBEDGE > Struct Reference

#include <src/post_proc/PostProcBrokenMeshInMoabBase.hpp>

Inheritance diagram for MoFEM::PostProcGenerateRefMesh< MBEDGE >:
[legend]
Collaboration diagram for MoFEM::PostProcGenerateRefMesh< MBEDGE >:
[legend]

Public Member Functions

MoFEMErrorCode generateReferenceElementMesh ()
 
 PostProcGenerateRefMeshBase ()
 
- Public Member Functions inherited from MoFEM::PostProcGenerateRefMeshBase
 PostProcGenerateRefMeshBase ()
 
virtual ~PostProcGenerateRefMeshBase ()=default
 
virtual MoFEMErrorCode getOptions (std::string prefix)
 

Additional Inherited Members

- Public Attributes inherited from MoFEM::PostProcGenerateRefMeshBase
std::vector< MatrixDoublelevelShapeFunctions
 
std::vector< MatrixDoublelevelGaussPtsOnRefMesh
 
std::vector< ublas::matrix< int > > levelRef
 
EntityHandle startingVertEleHandle
 
std::vector< double * > verticesOnEleArrays
 
EntityHandle startingEleHandle
 
EntityHandleeleConn
 
int countEle
 
int countVertEle
 
int nbVertices
 
int nbEles
 
PetscBool hoNodes
 
int defMaxLevel
 

Detailed Description

Definition at line 89 of file PostProcBrokenMeshInMoabBase.hpp.

Member Function Documentation

◆ generateReferenceElementMesh()

MoFEMErrorCode MoFEM::PostProcGenerateRefMesh< MBEDGE >::generateReferenceElementMesh ( )
virtual

Implements MoFEM::PostProcGenerateRefMeshBase.

Definition at line 530 of file PostProcBrokenMeshInMoabBase.cpp.

530  {
532 
533 #ifndef NDEBUG
534  if (defMaxLevel > 0)
535  MOFEM_LOG("WORLD", Sev::warning)
536  << "Refinement for edges is not implemented";
537 #endif
538 
539  auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
541 
542  int nb_nodes = 2;
543  if (hoNodes)
544  nb_nodes = 3;
545 
546  auto &gauss_pts = levelGaussPtsOnRefMesh[0];
547  gauss_pts.resize(2, nb_nodes, false);
548  gauss_pts.clear();
549 
550  int nn = 0;
551  for (; nn != 2; ++nn) {
552  gauss_pts(0, nn) = static_cast<double>(nn);
553  little_map[nn] = nn;
554  }
555 
556  if (nn < nb_nodes) {
557  gauss_pts(0, nn) = 0.5;
558  little_map[nn] = 2;
559  }
560 
562  };
563 
564  auto set_ref_edges = [&](std::map<EntityHandle, int> &little_map) {
566 
567  int level = 0;
568  int nb_edges = level + 1;
569 
570  int nb_nodes = 2;
571  if (hoNodes)
572  nb_nodes = 3;
573 
574  auto &ref_edges = levelRef[level];
575  ref_edges.resize(nb_edges, nb_nodes, false);
576 
577  for (int ee = 0; ee != nb_edges; ++ee) {
578  int nn = 0;
579  for (; nn != 2; ++nn) {
580  ref_edges(ee, nn) = nb_nodes * ee + nn;
581  }
582  if (nn < nb_nodes) {
583  ref_edges(ee, nn) = nb_nodes * ee + 2;
584  }
585  }
586 
588  };
589 
590  auto set_shape_functions = [&]() {
592  auto &gauss_pts = levelGaussPtsOnRefMesh[0];
593  auto &shape_functions = levelShapeFunctions[0];
594  const auto nb_gauss_pts = gauss_pts.size2();
595  shape_functions.resize(nb_gauss_pts, 2);
596  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
597  const double ksi = gauss_pts(0, gg);
598  shape_functions(gg, 0) = N_MBEDGE0(ksi);
599  shape_functions(gg, 1) = N_MBEDGE1(ksi);
600  }
602  };
603 
604  levelRef.resize(1);
605  levelGaussPtsOnRefMesh.resize(1);
606  levelShapeFunctions.resize(1);
607 
608  std::map<EntityHandle, int> little_map;
609  CHKERR set_gauss_pts(little_map);
610  CHKERR set_ref_edges(little_map);
611  CHKERR set_shape_functions();
612 
614 }

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase

Definition at line 13 of file PostProcBrokenMeshInMoabBase.cpp.

14  : hoNodes(PETSC_TRUE), defMaxLevel(0), countEle(0), countVertEle(0),
15  nbVertices(0) {}

The documentation for this struct was generated from the following files:
MoFEM::PostProcGenerateRefMeshBase::nbVertices
int nbVertices
Definition: PostProcBrokenMeshInMoabBase.hpp:40
N_MBEDGE0
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:105
MoFEM::PostProcGenerateRefMeshBase::levelShapeFunctions
std::vector< MatrixDouble > levelShapeFunctions
Definition: PostProcBrokenMeshInMoabBase.hpp:26
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
N_MBEDGE1
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:106
MoFEM::PostProcGenerateRefMeshBase::hoNodes
PetscBool hoNodes
Definition: PostProcBrokenMeshInMoabBase.hpp:49
MoFEM::PostProcGenerateRefMeshBase::levelGaussPtsOnRefMesh
std::vector< MatrixDouble > levelGaussPtsOnRefMesh
Definition: PostProcBrokenMeshInMoabBase.hpp:29
MoFEM::PostProcGenerateRefMeshBase::countVertEle
int countVertEle
Definition: PostProcBrokenMeshInMoabBase.hpp:38
MoFEM::PostProcGenerateRefMeshBase::defMaxLevel
int defMaxLevel
Definition: PostProcBrokenMeshInMoabBase.hpp:50
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::PostProcGenerateRefMeshBase::levelRef
std::vector< ublas::matrix< int > > levelRef
Definition: PostProcBrokenMeshInMoabBase.hpp:30
MoFEM::PostProcGenerateRefMeshBase::countEle
int countEle
Definition: PostProcBrokenMeshInMoabBase.hpp:37
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359