v0.13.2
Loading...
Searching...
No Matches
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
 
MoFEMErrorCode getOptions ()
 
virtual MoFEMErrorCode generateReferenceElementMesh ()=0
 

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
 
std::string optPrefix
 

Detailed Description

Definition at line 78 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
265 const int max_level = defMaxLevel;
266
267 moab::Core core_ref;
268 moab::Interface &moab_ref = core_ref;
269
270 auto create_reference_element = [&moab_ref]() {
272 constexpr double base_coords[] = {
273
274 0, 0,
275 0,
276
277 1, 0,
278 0,
279
280 0, 1,
281 0
282
283 };
284 EntityHandle nodes[3];
285 for (int nn = 0; nn != 3; ++nn)
286 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
287 EntityHandle tri;
288 CHKERR moab_ref.create_element(MBTRI, nodes, 3, tri);
289
290 Range edges;
291 CHKERR moab_ref.get_adjacencies(&tri, 1, 1, true, edges,
292 moab::Interface::UNION);
293
295 };
296
297 CHKERR create_reference_element();
298
299 MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
300 MoFEM::Interface &m_field_ref = m_core_ref;
301
302 auto refine_ref_triangles = [this, &m_field_ref, max_level]() {
304 // seed ref mofem database by setting bit ref level to reference
305 // tetrahedron
306 CHKERR
307 m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
308 0, 2, BitRefLevel().set(0));
309
310 for (int ll = 0; ll != max_level; ++ll) {
311 MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
312 ll);
313 Range edges;
314 CHKERR m_field_ref.getInterface<BitRefManager>()
315 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
316 BitRefLevel().set(), MBEDGE, edges);
317 Range tris;
318 CHKERR m_field_ref.getInterface<BitRefManager>()
319 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
320 BitRefLevel(ll).set(), MBTRI, tris);
321 // refine mesh
322 auto m_ref = m_field_ref.getInterface<MeshRefinement>();
323 CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
324 BitRefLevel().set(ll + 1));
325 CHKERR m_ref->refineTris(tris, BitRefLevel().set(ll + 1));
326 }
328 };
329
330 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
332 Range faces;
333 CHKERR moab_ref.get_entities_by_type(0, MBTRI, faces, true);
334 Range faces_nodes;
335 CHKERR moab_ref.get_connectivity(faces, faces_nodes, false);
336 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
337 gauss_pts.resize(faces_nodes.size(), 4, false);
338 size_t gg = 0;
339 for (auto node : faces_nodes) {
340 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
341 little_map[node] = gg;
342 ++gg;
343 }
344 gauss_pts = trans(gauss_pts);
346 };
347
348 auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
349 &m_field_ref]() {
351 for (int ll = 0; ll != max_level + 1; ++ll) {
352
353 Range tris;
354 CHKERR
355 m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
356 BitRefLevel().set(ll), BitRefLevel().set(ll), MBTRI, tris);
357
358 if (hoNodes) {
359 EntityHandle meshset;
360 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
361 CHKERR moab_ref.add_entities(meshset, tris);
362 CHKERR moab_ref.convert_entities(meshset, true, false, false);
363 CHKERR moab_ref.delete_entities(&meshset, 1);
364 }
365
366 Range elem_nodes;
367 CHKERR moab_ref.get_connectivity(tris, elem_nodes, false);
368
369 auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
370 gauss_pts.resize(elem_nodes.size(), 3, false);
371 std::map<EntityHandle, int> little_map;
372 Range::iterator nit = elem_nodes.begin();
373 for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
374 CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
375 little_map[*nit] = gg;
376 }
377 gauss_pts = trans(gauss_pts);
378
379 auto &ref_tris = levelRef[ll];
380 Range::iterator tit = tris.begin();
381 for (int tt = 0; tit != tris.end(); ++tit, ++tt) {
382 const EntityHandle *conn;
383 int num_nodes;
384 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
385 if (tt == 0) {
386 ref_tris.resize(tris.size(), num_nodes);
387 }
388 for (int nn = 0; nn != num_nodes; ++nn) {
389 ref_tris(tt, nn) = little_map[conn[nn]];
390 }
391 }
392
393 auto &shape_functions = levelShapeFunctions[ll];
394 shape_functions.resize(elem_nodes.size(), 3);
395 CHKERR ShapeMBTRI(&*shape_functions.data().begin(), &gauss_pts(0, 0),
396 &gauss_pts(1, 0), elem_nodes.size());
397 }
399 };
400
401 levelRef.resize(max_level + 1);
402 levelGaussPtsOnRefMesh.resize(max_level + 1);
403 levelShapeFunctions.resize(max_level + 1);
404
405 CHKERR refine_ref_triangles();
406 CHKERR get_ref_gauss_pts_and_shape_functions();
407
409}
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
Definition: LogManager.hpp:367
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
Deprecated interface functions.
std::vector< ublas::matrix< int > > levelRef
PetscBool hoNodes
int defMaxLevel
std::vector< MatrixDouble > levelShapeFunctions
std::vector< MatrixDouble > levelGaussPtsOnRefMesh
MoFEMErrorCode getOptions()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase ( )

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