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

#include <src/post_proc/PostProcBrokenMeshInMoabBase.hpp>

Inheritance diagram for MoFEM::PostProcGenerateRefMesh< MBHEX >:
[legend]
Collaboration diagram for MoFEM::PostProcGenerateRefMesh< MBHEX >:
[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 71 of file PostProcBrokenMeshInMoabBase.hpp.

Member Function Documentation

◆ generateReferenceElementMesh()

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

Implements MoFEM::PostProcGenerateRefMeshBase.

Definition at line 149 of file PostProcBrokenMeshInMoabBase.cpp.

149  {
151 
152 #ifndef NDEBUG
153  if (defMaxLevel > 0)
154  MOFEM_LOG("WORLD", Sev::warning)
155  << "Refinement for hexes is not implemented";
156 #endif
157 
158  moab::Core core_ref;
159  moab::Interface &moab_ref = core_ref;
160 
161  auto create_reference_element = [&moab_ref]() {
163  constexpr double base_coords[] = {0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1, 0,
164  0, 0, 1, 1, 0, 1, 1, 1, 1, 0, 1, 1};
165  EntityHandle nodes[8];
166  for (int nn = 0; nn < 8; nn++)
167  CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
168  EntityHandle hex;
169  CHKERR moab_ref.create_element(MBHEX, nodes, 8, hex);
171  };
172 
173  auto add_ho_nodes = [&]() {
175  Range hexes;
176  CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
177  EntityHandle meshset;
178  CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
179  CHKERR moab_ref.add_entities(meshset, hexes);
180  CHKERR moab_ref.convert_entities(meshset, true, true, true);
181  CHKERR moab_ref.delete_entities(&meshset, 1);
183  };
184 
185  auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
187  Range hexes;
188  CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
189  Range hexes_nodes;
190  CHKERR moab_ref.get_connectivity(hexes, hexes_nodes, false);
191  auto &gauss_pts = levelGaussPtsOnRefMesh[0];
192  gauss_pts.resize(hexes_nodes.size(), 4, false);
193  size_t gg = 0;
194  for (auto node : hexes_nodes) {
195  CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
196  little_map[node] = gg;
197  ++gg;
198  }
199  gauss_pts = trans(gauss_pts);
201  };
202 
203  auto set_ref_hexes = [&](std::map<EntityHandle, int> &little_map) {
205  Range hexes;
206  CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
207  size_t hh = 0;
208  auto &ref_hexes = levelRef[0];
209  for (auto hex : hexes) {
210  const EntityHandle *conn;
211  int num_nodes;
212  CHKERR moab_ref.get_connectivity(hex, conn, num_nodes, false);
213  if (ref_hexes.size2() != num_nodes) {
214  ref_hexes.resize(hexes.size(), num_nodes);
215  }
216  for (int nn = 0; nn != num_nodes; ++nn) {
217  ref_hexes(hh, nn) = little_map[conn[nn]];
218  }
219  ++hh;
220  }
222  };
223 
224  auto set_shape_functions = [&]() {
226  auto &gauss_pts = levelGaussPtsOnRefMesh[0];
227  auto &shape_functions = levelShapeFunctions[0];
228  const auto nb_gauss_pts = gauss_pts.size2();
229  shape_functions.resize(nb_gauss_pts, 8);
230  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
231  const double ksi = gauss_pts(0, gg);
232  const double zeta = gauss_pts(1, gg);
233  const double eta = gauss_pts(2, gg);
234  shape_functions(gg, 0) = N_MBHEX0(ksi, zeta, eta);
235  shape_functions(gg, 1) = N_MBHEX1(ksi, zeta, eta);
236  shape_functions(gg, 2) = N_MBHEX2(ksi, zeta, eta);
237  shape_functions(gg, 3) = N_MBHEX3(ksi, zeta, eta);
238  shape_functions(gg, 4) = N_MBHEX4(ksi, zeta, eta);
239  shape_functions(gg, 5) = N_MBHEX5(ksi, zeta, eta);
240  shape_functions(gg, 6) = N_MBHEX6(ksi, zeta, eta);
241  shape_functions(gg, 7) = N_MBHEX7(ksi, zeta, eta);
242  }
244  };
245 
246  levelRef.resize(1);
247  levelGaussPtsOnRefMesh.resize(1);
248  levelShapeFunctions.resize(1);
249 
250  CHKERR create_reference_element();
251  if (hoNodes)
252  CHKERR add_ho_nodes();
253  std::map<EntityHandle, int> little_map;
254  CHKERR set_gauss_pts(little_map);
255  CHKERR set_ref_hexes(little_map);
256  CHKERR set_shape_functions();
257 
259 }

◆ 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:
EntityHandle
N_MBHEX5
#define N_MBHEX5(x, y, z)
Definition: fem_tools.h:76
N_MBHEX0
#define N_MBHEX0(x, y, z)
Definition: fem_tools.h:71
MoFEM::PostProcGenerateRefMeshBase::nbVertices
int nbVertices
Definition: PostProcBrokenMeshInMoabBase.hpp:40
N_MBHEX7
#define N_MBHEX7(x, y, z)
Definition: fem_tools.h:78
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:177
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
N_MBHEX3
#define N_MBHEX3(x, y, z)
Definition: fem_tools.h:74
MoFEM::PostProcGenerateRefMeshBase::levelShapeFunctions
std::vector< MatrixDouble > levelShapeFunctions
Definition: PostProcBrokenMeshInMoabBase.hpp:26
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::PostProcGenerateRefMeshBase::hoNodes
PetscBool hoNodes
Definition: PostProcBrokenMeshInMoabBase.hpp:49
eta
double eta
Definition: free_surface.cpp:170
N_MBHEX1
#define N_MBHEX1(x, y, z)
Definition: fem_tools.h:72
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
Range
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
N_MBHEX2
#define N_MBHEX2(x, y, z)
Definition: fem_tools.h:73
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
MoFEM::PostProcGenerateRefMeshBase::countEle
int countEle
Definition: PostProcBrokenMeshInMoabBase.hpp:37
N_MBHEX4
#define N_MBHEX4(x, y, z)
Definition: fem_tools.h:75
N_MBHEX6
#define N_MBHEX6(x, y, z)
Definition: fem_tools.h:77
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346