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

#include <src/post_proc/PostProcBrokenMeshInMoabBase.hpp>

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

Member Function Documentation

◆ generateReferenceElementMesh()

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

Implements MoFEM::PostProcGenerateRefMeshBase.

Definition at line 261 of file PostProcBrokenMeshInMoabBase.cpp.

261  {
263 
264  const int max_level = defMaxLevel;
265 
266  moab::Core core_ref;
267  moab::Interface &moab_ref = core_ref;
268 
269  auto create_reference_element = [&moab_ref]() {
271  constexpr double base_coords[] = {
272 
273  0, 0,
274  0,
275 
276  1, 0,
277  0,
278 
279  0, 1,
280  0
281 
282  };
283  EntityHandle nodes[3];
284  for (int nn = 0; nn != 3; ++nn)
285  CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
286  EntityHandle tri;
287  CHKERR moab_ref.create_element(MBTRI, nodes, 3, tri);
288 
289  Range edges;
290  CHKERR moab_ref.get_adjacencies(&tri, 1, 1, true, edges,
291  moab::Interface::UNION);
292 
294  };
295 
296  CHKERR create_reference_element();
297 
298  MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
299  MoFEM::Interface &m_field_ref = m_core_ref;
300 
301  auto refine_ref_triangles = [this, &m_field_ref, max_level]() {
303  // seed ref mofem database by setting bit ref level to reference
304  // tetrahedron
305  CHKERR
306  m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
307  0, 2, BitRefLevel().set(0));
308 
309  for (int ll = 0; ll != max_level; ++ll) {
310  MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
311  ll);
312  Range edges;
313  CHKERR m_field_ref.getInterface<BitRefManager>()
314  ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
315  BitRefLevel().set(), MBEDGE, edges);
316  Range tris;
317  CHKERR m_field_ref.getInterface<BitRefManager>()
318  ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
319  BitRefLevel(ll).set(), MBTRI, tris);
320  // refine mesh
321  auto m_ref = m_field_ref.getInterface<MeshRefinement>();
322  CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
323  BitRefLevel().set(ll + 1));
324  CHKERR m_ref->refineTris(tris, BitRefLevel().set(ll + 1));
325  }
327  };
328 
329  // auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
330  // MoFEMFunctionBegin;
331  // Range faces;
332  // CHKERR moab_ref.get_entities_by_type(0, MBTRI, faces, true);
333  // Range faces_nodes;
334  // CHKERR moab_ref.get_connectivity(faces, faces_nodes, false);
335  // auto &gauss_pts = levelGaussPtsOnRefMesh[0];
336  // gauss_pts.resize(faces_nodes.size(), 4, false);
337  // size_t gg = 0;
338  // for (auto node : faces_nodes) {
339  // CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
340  // little_map[node] = gg;
341  // ++gg;
342  // }
343  // gauss_pts = trans(gauss_pts);
344  // MoFEMFunctionReturn(0);
345  // };
346 
347  auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
348  &m_field_ref]() {
350  for (int ll = 0; ll != max_level + 1; ++ll) {
351 
352  Range tris;
353  CHKERR
354  m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
355  BitRefLevel().set(ll), BitRefLevel().set(ll), MBTRI, tris);
356 
357  if (hoNodes) {
358  EntityHandle meshset;
359  CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
360  CHKERR moab_ref.add_entities(meshset, tris);
361  CHKERR moab_ref.convert_entities(meshset, true, false, false);
362  CHKERR moab_ref.delete_entities(&meshset, 1);
363  }
364 
365  Range elem_nodes;
366  CHKERR moab_ref.get_connectivity(tris, elem_nodes, false);
367 
368  auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
369  gauss_pts.resize(elem_nodes.size(), 3, false);
370  std::map<EntityHandle, int> little_map;
371  Range::iterator nit = elem_nodes.begin();
372  for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
373  CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
374  little_map[*nit] = gg;
375  }
376  gauss_pts = trans(gauss_pts);
377 
378  auto &ref_tris = levelRef[ll];
379  Range::iterator tit = tris.begin();
380  for (int tt = 0; tit != tris.end(); ++tit, ++tt) {
381  const EntityHandle *conn;
382  int num_nodes;
383  CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
384  if (tt == 0) {
385  ref_tris.resize(tris.size(), num_nodes);
386  }
387  for (int nn = 0; nn != num_nodes; ++nn) {
388  ref_tris(tt, nn) = little_map[conn[nn]];
389  }
390  }
391 
392  auto &shape_functions = levelShapeFunctions[ll];
393  shape_functions.resize(elem_nodes.size(), 3);
394  CHKERR ShapeMBTRI(&*shape_functions.data().begin(), &gauss_pts(0, 0),
395  &gauss_pts(1, 0), elem_nodes.size());
396  }
398  };
399 
400  levelRef.resize(max_level + 1);
401  levelGaussPtsOnRefMesh.resize(max_level + 1);
402  levelShapeFunctions.resize(max_level + 1);
403 
404  CHKERR refine_ref_triangles();
405  CHKERR get_ref_gauss_pts_and_shape_functions();
406 
408 }

◆ 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::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
EntityHandle
MoFEM::PostProcGenerateRefMeshBase::nbVertices
int nbVertices
Definition: PostProcBrokenMeshInMoabBase.hpp:40
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::PostProcGenerateRefMeshBase::levelShapeFunctions
std::vector< MatrixDouble > levelShapeFunctions
Definition: PostProcBrokenMeshInMoabBase.hpp:26
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreTmp
Definition: Core.hpp:36
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
Range
MoFEM::PostProcGenerateRefMeshBase::levelRef
std::vector< ublas::matrix< int > > levelRef
Definition: PostProcBrokenMeshInMoabBase.hpp:30
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
MOFEM_TAG_AND_LOG_C
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
Definition: LogManager.hpp:370
MoFEM::PostProcGenerateRefMeshBase::countEle
int countEle
Definition: PostProcBrokenMeshInMoabBase.hpp:37
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
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
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182