v0.15.0
Loading...
Searching...
No Matches
MoFEM::PostProcBrokenMeshInMoabBase< E > Struct Template Reference

#include "src/post_proc/PostProcBrokenMeshInMoabBase.hpp"

Inheritance diagram for MoFEM::PostProcBrokenMeshInMoabBase< E >:
[legend]
Collaboration diagram for MoFEM::PostProcBrokenMeshInMoabBase< E >:
[legend]

Public Member Functions

 PostProcBrokenMeshInMoabBase (MoFEM::Interface &m_field, std::string opts_prefix="")
 
 PostProcBrokenMeshInMoabBase (MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, std::string opts_prefix="")
 
virtual ~PostProcBrokenMeshInMoabBase ()
 
auto & getMapGaussPts ()
 Get vector of vectors associated to integration points.
 
auto & getPostProcMesh ()
 Get postprocessing mesh.
 
auto & getPostProcElements ()
 Get postprocessing elements.
 
MoFEMErrorCode writeFile (const std::string file_name)
 wrote results in (MOAB) format, use "file_name.h5m"
 
MoFEMErrorCode setTagsToTransfer (std::vector< Tag > tags_to_transfer)
 Set tags to be transferred to post-processing mesh.
 

Protected Member Functions

MoFEMErrorCode setGaussPts (int order)
 
MoFEMErrorCode preProcess ()
 Generate vertices and elements.
 
MoFEMErrorCode postProcess ()
 
int getRule (int order)
 
virtual int getMaxLevel () const
 Determine refinement level based on fields approx ordre.
 
virtual MoFEMErrorCode preProcPostProc ()
 
virtual MoFEMErrorCode postProcPostProc ()
 
virtual MoFEMErrorCode transferTags ()
 

Protected Attributes

boost::shared_ptr< moab::Core > coreMeshPtr = boost::make_shared<moab::Core>()
 
std::vector< EntityHandlemapGaussPts
 
Range postProcElements
 
std::map< EntityType, PostProcGenerateRefMeshPtrrefElementsMap
 
std::vector< TagtagsToTransfer
 
std::string optionsPrefix = ""
 Prefix for options.
 

Detailed Description

template<typename E>
struct MoFEM::PostProcBrokenMeshInMoabBase< E >

Definition at line 94 of file PostProcBrokenMeshInMoabBase.hpp.

Constructor & Destructor Documentation

◆ PostProcBrokenMeshInMoabBase() [1/2]

template<typename E >
PostProcBrokenMeshInMoabBase::PostProcBrokenMeshInMoabBase ( MoFEM::Interface & m_field,
std::string opts_prefix = "" )

Definition at line 277 of file PostProcBrokenMeshInMoabBase.hpp.

279 : E(m_field), optionsPrefix(opts_prefix) {}

◆ PostProcBrokenMeshInMoabBase() [2/2]

template<typename E >
PostProcBrokenMeshInMoabBase::PostProcBrokenMeshInMoabBase ( MoFEM::Interface & m_field,
boost::shared_ptr< moab::Core > core_mesh_ptr,
std::string opts_prefix = "" )

Definition at line 282 of file PostProcBrokenMeshInMoabBase.hpp.

285 : PostProcBrokenMeshInMoabBase(m_field, opts_prefix) {
286 coreMeshPtr = core_mesh_ptr;
287}
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field, std::string opts_prefix="")

◆ ~PostProcBrokenMeshInMoabBase()

template<typename E >
PostProcBrokenMeshInMoabBase::~PostProcBrokenMeshInMoabBase ( )
virtual

Definition at line 290 of file PostProcBrokenMeshInMoabBase.hpp.

290 {
291 ParallelComm *pcomm_post_proc_mesh =
292 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
293 if (pcomm_post_proc_mesh != NULL)
294 delete pcomm_post_proc_mesh;
295}
#define MYPCOMM_INDEX
default communicator number PCOMM
auto & getPostProcMesh()
Get postprocessing mesh.

Member Function Documentation

◆ getMapGaussPts()

template<typename E >
auto & PostProcBrokenMeshInMoabBase::getMapGaussPts ( )
inline

Get vector of vectors associated to integration points.

Returns
std::vector<EntityHandle>&

Definition at line 647 of file PostProcBrokenMeshInMoabBase.hpp.

647 {
648 return mapGaussPts;
649}

◆ getMaxLevel()

template<typename E >
int PostProcBrokenMeshInMoabBase::getMaxLevel ( ) const
protectedvirtual

Determine refinement level based on fields approx ordre.

level = (order - 1) / 2

Returns
int

Definition at line 193 of file PostProcBrokenMeshInMoabBase.hpp.

193 {
194 auto get_element_max_dofs_order = [&]() {
195 int max_order = 0;
196 auto dofs_vec = E::getDataVectorDofsPtr();
197 for (auto &dof : *dofs_vec) {
198 const int dof_order = dof->getDofOrder();
199 max_order = (max_order < dof_order) ? dof_order : max_order;
200 };
201 return max_order;
202 };
203 const auto dof_max_order = get_element_max_dofs_order();
204 return (dof_max_order > 0) ? (dof_max_order - 1) / 2 : 0;
205};

◆ getPostProcElements()

template<typename E >
auto & PostProcBrokenMeshInMoabBase::getPostProcElements ( )
inline

Get postprocessing elements.

Returns
auto&

Definition at line 659 of file PostProcBrokenMeshInMoabBase.hpp.

◆ getPostProcMesh()

template<typename E >
auto & PostProcBrokenMeshInMoabBase::getPostProcMesh ( )
inline

Get postprocessing mesh.

Returns
moab::Interface&

Definition at line 651 of file PostProcBrokenMeshInMoabBase.hpp.

651 {
652 if (!coreMeshPtr)
653 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Core mesh not set");
654 moab::Interface &post_proc_mesh_interface = *coreMeshPtr;
655 return post_proc_mesh_interface;
656}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31

◆ getRule()

template<typename E >
int PostProcBrokenMeshInMoabBase::getRule ( int order)
protected

Definition at line 297 of file PostProcBrokenMeshInMoabBase.hpp.

297 {
298 return -1;
299};

◆ postProcess()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::postProcess ( )
protected

Definition at line 633 of file PostProcBrokenMeshInMoabBase.hpp.

633 {
635
636 ParallelComm *pcomm_post_proc_mesh =
637 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
638 if(pcomm_post_proc_mesh == nullptr)
639 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
640
642 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
643
645}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

◆ postProcPostProc()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::postProcPostProc ( )
protectedvirtual

Definition at line 544 of file PostProcBrokenMeshInMoabBase.hpp.

544 {
546
547 auto update_elements = [&]() {
548 ReadUtilIface *iface;
549 CHKERR getPostProcMesh().query_interface(iface);
551
552 Range ents;
553 for (auto &m : refElementsMap) {
554 if (m.second->nbEles) {
555 MOFEM_TAG_AND_LOG("SELF", Sev::noisy, "PostProc")
556 << "Update < " << moab::CN::EntityTypeName(m.first)
557 << " number of processed " << m.second->countEle;
558 CHKERR iface->update_adjacencies(
559 m.second->startingEleHandle, m.second->countEle,
560 m.second->levelRef[0].size2(), m.second->eleConn);
561 ents.merge(Range(m.second->startingEleHandle,
562 m.second->startingEleHandle + m.second->countEle - 1));
563 }
564 }
565
567 };
568
569 auto remove_obsolete_entities = [&]() {
571 Range ents, adj;
572 for (auto &m : refElementsMap) {
573 if (m.second->nbEles) {
574 ents.merge(Range(m.second->startingEleHandle,
575 m.second->startingEleHandle + m.second->countEle - 1));
576 const int dim = moab::CN::Dimension(m.first);
577 for (auto d = 1; d != dim; ++d) {
578 CHKERR getPostProcMesh().get_adjacencies(ents, d, false, adj,
579 moab::Interface::UNION);
580 }
581 }
582 }
583 CHKERR getPostProcMesh().delete_entities(adj);
585 };
586
587 auto set_proc_tags = [&]() {
589 ParallelComm *pcomm_post_proc_mesh =
590 ParallelComm::get_pcomm(&(getPostProcMesh()), MYPCOMM_INDEX);
591 if (pcomm_post_proc_mesh) {
592 Range ents;
593 for (auto &m : refElementsMap) {
594 if (m.second->nbEles) {
595 ents.merge(
596 Range(m.second->startingEleHandle,
597 m.second->startingEleHandle + m.second->countEle - 1));
598 }
599 }
600
601 int rank = E::mField.get_comm_rank();
602 CHKERR getPostProcMesh().tag_clear_data(pcomm_post_proc_mesh->part_tag(),
603 ents, &rank);
604 }
606 };
607
608 CHKERR update_elements();
609 CHKERR remove_obsolete_entities();
610 CHKERR set_proc_tags();
611
613}
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'm', 3 > m
std::map< EntityType, PostProcGenerateRefMeshPtr > refElementsMap

◆ preProcess()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::preProcess ( )
protected

Generate vertices and elements.

Returns
MoFEMErrorCode

Definition at line 616 of file PostProcBrokenMeshInMoabBase.hpp.

616 {
618
619 ParallelComm *pcomm_post_proc_mesh =
620 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
621 if (pcomm_post_proc_mesh != NULL)
622 delete pcomm_post_proc_mesh;
623 CHKERR getPostProcMesh().delete_mesh();
624 pcomm_post_proc_mesh =
625 new ParallelComm(&(getPostProcMesh()), PETSC_COMM_WORLD);
626
628
630}

◆ preProcPostProc()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::preProcPostProc ( )
protectedvirtual

Definition at line 429 of file PostProcBrokenMeshInMoabBase.hpp.

429 {
431
432auto get_ref_ele = [&](const EntityType type) {
433 PostProcGenerateRefMeshPtr ref_ele_ptr;
434
435 auto it = refElementsMap.find(type);
436 if (it != refElementsMap.end()) {
437 ref_ele_ptr = it->second;
438 } else {
439 switch (type) {
440 case MBTET:
441 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTET>>();
442 break;
443 case MBHEX:
444 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBHEX>>();
445 break;
446 case MBTRI:
447 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
448 break;
449 case MBQUAD:
450 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBQUAD>>();
451 break;
452 case MBEDGE:
453 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBEDGE>>();
454 break;
455 default:
456 MOFEM_LOG("SELF", Sev::error)
457 << "Generation of reference elements for type < "
458 << moab::CN::EntityTypeName(type) << " > is not implemented";
459 CHK_THROW_MESSAGE(MOFEM_NOT_IMPLEMENTED, "Element not implemented");
460 }
461
462 CHK_THROW_MESSAGE(ref_ele_ptr->getOptions(optionsPrefix), "getOptions");
463 CHK_THROW_MESSAGE(ref_ele_ptr->generateReferenceElementMesh(),
464 "Error when generating reference element");
465
466 refElementsMap[type] = ref_ele_ptr;
467 }
468
469 return ref_ele_ptr;
470};
471
472 auto fe_ptr = this->problemPtr->numeredFiniteElementsPtr;
473
474 auto miit =
475 fe_ptr->template get<Composite_Name_And_Part_mi_tag>().lower_bound(
476 boost::make_tuple(this->getFEName(), this->getLoFERank()));
477 auto hi_miit =
478 fe_ptr->template get<Composite_Name_And_Part_mi_tag>().upper_bound(
479 boost::make_tuple(this->getFEName(), this->getHiFERank()));
480
481 const int number_of_ents_in_the_loop = this->getLoopSize();
482 if (std::distance(miit, hi_miit) != number_of_ents_in_the_loop) {
483 SETERRQ(E::mField.get_comm(), MOFEM_DATA_INCONSISTENCY,
484 "Wrong size of indicies. Inconsistent size number of iterated "
485 "elements iterated by problem and from range.");
486 }
487
488 for (auto &m : refElementsMap) {
489 m.second->nbVertices = 0;
490 m.second->nbEles = 0;
491 m.second->countEle = 0;
492 m.second->countVertEle = 0;
493 }
494
495 for (; miit != hi_miit; ++miit) {
496 auto type = (*miit)->getEntType();
497 auto ref_ele = get_ref_ele(type);
498
499 // Set pointer to element. So that getDataVectorDofsPtr in getMaxLevel
500 // can work
501 E::numeredEntFiniteElementPtr = *miit;
502 bool add = true;
503 if (E::exeTestHook) {
504 add = E::exeTestHook(this);
505 }
506
507 if (add) {
508 size_t level = getMaxLevel();
509 level = std::min(level, ref_ele->levelGaussPtsOnRefMesh.size() - 1);
510 ref_ele->nbVertices += ref_ele->levelGaussPtsOnRefMesh[level].size2();
511 ref_ele->nbEles += ref_ele->levelRef[level].size1();
512 }
513 }
514
515 auto alloc_vertices_and_elements_on_post_proc_mesh = [&]() {
517
518 ReadUtilIface *iface;
519 CHKERR getPostProcMesh().query_interface(iface);
520
521 for (auto &m : refElementsMap) {
522 if (m.second->nbEles) {
523 CHKERR iface->get_node_coords(3, m.second->nbVertices, 0,
524 m.second->startingVertEleHandle,
525 m.second->verticesOnEleArrays);
526 CHKERR iface->get_element_connect(
527 m.second->nbEles, m.second->levelRef[0].size2(), m.first, 0,
528 m.second->startingEleHandle, m.second->eleConn);
529
530 m.second->countEle = 0;
531 m.second->countVertEle = 0;
532 }
533 }
534
536 };
537
538 CHKERR alloc_vertices_and_elements_on_post_proc_mesh();
539
541}
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MOFEM_LOG(channel, severity)
Log.
boost::shared_ptr< PostProcGenerateRefMeshBase > PostProcGenerateRefMeshPtr
virtual int getMaxLevel() const
Determine refinement level based on fields approx ordre.

◆ setGaussPts()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::setGaussPts ( int order)
protected

Definition at line 302 of file PostProcBrokenMeshInMoabBase.hpp.

302 {
304
305 auto type = type_from_handle(this->getFEEntityHandle());
306
308
309 try {
310 ref_ele = refElementsMap.at(type);
311 } catch (const out_of_range &e) {
312 SETERRQ(
313 PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
314 "Generation of reference elements for type <%s> is not implemented",
315 moab::CN::EntityTypeName(type));
316 }
317
318 auto set_gauss_pts = [&](auto &level_gauss_pts_on_ref_mesh, auto &level_ref,
319 auto &level_shape_functions,
320
321 auto start_vert_handle, auto start_ele_handle,
322 auto &verts_array, auto &conn, auto &ver_count,
323 auto &ele_count
324
325 ) {
327
328 size_t level = getMaxLevel();
329 level = std::min(level, level_gauss_pts_on_ref_mesh.size() - 1);
330
331 auto &level_ref_gauss_pts = level_gauss_pts_on_ref_mesh[level];
332 auto &level_ref_ele = level_ref[level];
333 auto &shape_functions = level_shape_functions[level];
334 E::gaussPts.resize(level_ref_gauss_pts.size1(), level_ref_gauss_pts.size2(),
335 false);
336 noalias(E::gaussPts) = level_ref_gauss_pts;
337
338 const auto fe_ent = E::numeredEntFiniteElementPtr->getEnt();
339 auto get_fe_coords = [&]() {
340 const EntityHandle *conn;
341 int num_nodes;
343 E::mField.get_moab().get_connectivity(fe_ent, conn, num_nodes, true),
344 "error get connectivity");
345 VectorDouble coords(num_nodes * 3);
347 E::mField.get_moab().get_coords(conn, num_nodes, &*coords.begin()),
348 "error get coordinates");
349 return coords;
350 };
351
352 auto coords = get_fe_coords();
353
354 const int num_nodes = level_ref_gauss_pts.size2();
355 mapGaussPts.resize(level_ref_gauss_pts.size2());
356
357 FTensor::Index<'i', 3> i;
359 &*shape_functions.data().begin());
361 &verts_array[0][ver_count], &verts_array[1][ver_count],
362 &verts_array[2][ver_count]);
363 for (int gg = 0; gg != num_nodes; ++gg, ++ver_count) {
364
365 mapGaussPts[gg] = start_vert_handle + ver_count;
366
367 auto set_float_precision = [](const double x) {
368 if (std::abs(x) < std::numeric_limits<float>::epsilon())
369 return 0.;
370 else
371 return x;
372 };
373
374 t_coords(i) = 0;
375 auto t_ele_coords = getFTensor1FromArray<3, 3>(coords);
376 for (int nn = 0; nn != CN::VerticesPerEntity(type); ++nn) {
377 t_coords(i) += t_n * t_ele_coords(i);
378 ++t_ele_coords;
379 ++t_n;
380 }
381
382 for (auto ii : {0, 1, 2})
383 t_coords(ii) = set_float_precision(t_coords(ii));
384
385 ++t_coords;
386 }
387
388 Tag th;
389 int def_in_the_loop = -1;
390 CHKERR getPostProcMesh().tag_get_handle(
391 "NB_IN_THE_LOOP", 1, MB_TYPE_INTEGER, th, MB_TAG_CREAT | MB_TAG_SPARSE,
392 &def_in_the_loop);
393
394 postProcElements.clear();
395 const int num_el = level_ref_ele.size1();
396 const int num_nodes_on_ele = level_ref_ele.size2();
397 auto start_e = start_ele_handle + ele_count;
398 postProcElements = Range(start_e, start_e + num_el - 1);
399 for (auto tt = 0; tt != level_ref_ele.size1(); ++tt, ++ele_count) {
400 for (int nn = 0; nn != num_nodes_on_ele; ++nn) {
401 conn[num_nodes_on_ele * ele_count + nn] =
402 mapGaussPts[level_ref_ele(tt, nn)];
403 }
404 }
405
406 const int n_in_the_loop = E::nInTheLoop;
407 CHKERR getPostProcMesh().tag_clear_data(th, postProcElements,
408 &n_in_the_loop);
410
412 };
413
414 CHKERR set_gauss_pts(
415
416 ref_ele->levelGaussPtsOnRefMesh, ref_ele->levelRef,
417 ref_ele->levelShapeFunctions,
418
419 ref_ele->startingVertEleHandle, ref_ele->startingEleHandle,
420 ref_ele->verticesOnEleArrays, ref_ele->eleConn, ref_ele->countVertEle,
421 ref_ele->countEle
422
423 );
424
426};
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
FTensor::Index< 'i', SPACE_DIM > i
auto type_from_handle(const EntityHandle h)
get type from entity handle
auto getFTensor1FromArray(VectorDouble &data)
Get FTensor1 from array.

◆ setTagsToTransfer()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::setTagsToTransfer ( std::vector< Tag > tags_to_transfer)

Set tags to be transferred to post-processing mesh.

Parameters
tags_to_transfer
Returns
MoFEMErrorCode

Definition at line 186 of file PostProcBrokenMeshInMoabBase.hpp.

187 {
189 tagsToTransfer.swap(tags_to_transfer);
191}

◆ transferTags()

template<typename E >
MoFEMErrorCode PostProcBrokenMeshInMoabBase::transferTags ( )
protectedvirtual

transfer tags from mesh to post-process mesh

Definition at line 208 of file PostProcBrokenMeshInMoabBase.hpp.

208 {
210
211 auto &calc_mesh = this->mField.get_moab();
212
213 auto name = [&](auto tag) {
214 std::string name;
215 CHK_MOAB_THROW(calc_mesh.tag_get_name(tag, name), "get name");
216 return name;
217 };
218
219 auto data_type = [&](auto tag) {
220 moab::DataType data_type;
221 CHK_MOAB_THROW(calc_mesh.tag_get_data_type(tag, data_type),
222 "get data type");
223 return data_type;
224 };
225
226 auto type = [&](auto tag) {
227 moab::TagType type;
228 CHK_MOAB_THROW(calc_mesh.tag_get_type(tag, type), "get tag type");
229 return type;
230 };
231
232 auto length = [&](auto tag) {
233 int length;
234 CHK_MOAB_THROW(calc_mesh.tag_get_length(tag, length), "get length ");
235 return length;
236 };
237
238 auto default_value = [&](auto tag) {
239 const void *def_val;
240 int size;
241 CHK_MOAB_THROW(calc_mesh.tag_get_default_value(tag, def_val, size),
242 "get default tag value");
243 return def_val;
244 };
245
246 std::vector<double> tag_data_vec;
247 auto data_ptr = [&](auto tag) {
248 const void *tag_data;
249 auto core_ent = this->getFEEntityHandle();
250 CHK_MOAB_THROW(calc_mesh.tag_get_by_ptr(tag, &core_ent, 1, &tag_data),
251 "get tag data");
252 if (data_type(tag) == MB_TYPE_DOUBLE) {
253 tag_data_vec.resize(length(tag));
254 std::copy(static_cast<const double *>(tag_data),
255 static_cast<const double *>(tag_data) + length(tag),
256 tag_data_vec.begin());
257 for (auto &v : tag_data_vec)
258 v = std::abs(v) < std::numeric_limits<float>::epsilon() ? 0. : v;
259 return static_cast<const void *>(tag_data_vec.data());
260 }
261 return tag_data;
262 };
263
264 for (auto tag : tagsToTransfer) {
265 Tag tag_postproc;
266 CHKERR getPostProcMesh().tag_get_handle(
267 name(tag).c_str(), length(tag), data_type(tag), tag_postproc,
268 type(tag) | MB_TAG_CREAT, default_value(tag));
269 CHKERR getPostProcMesh().tag_clear_data(tag_postproc, postProcElements,
270 data_ptr(tag));
271 }
272
274};
const double v
phase velocity of light in medium (cm/ns)

Member Data Documentation

◆ coreMeshPtr

template<typename E >
boost::shared_ptr<moab::Core> MoFEM::PostProcBrokenMeshInMoabBase< E >::coreMeshPtr = boost::make_shared<moab::Core>()
protected

Definition at line 165 of file PostProcBrokenMeshInMoabBase.hpp.

◆ mapGaussPts

template<typename E >
std::vector<EntityHandle> MoFEM::PostProcBrokenMeshInMoabBase< E >::mapGaussPts
protected

Definition at line 167 of file PostProcBrokenMeshInMoabBase.hpp.

◆ optionsPrefix

template<typename E >
std::string MoFEM::PostProcBrokenMeshInMoabBase< E >::optionsPrefix = ""
protected

Prefix for options.

Definition at line 176 of file PostProcBrokenMeshInMoabBase.hpp.

◆ postProcElements

template<typename E >
Range MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcElements
protected

Definition at line 168 of file PostProcBrokenMeshInMoabBase.hpp.

◆ refElementsMap

template<typename E >
std::map<EntityType, PostProcGenerateRefMeshPtr> MoFEM::PostProcBrokenMeshInMoabBase< E >::refElementsMap
protected

Storing data about element types, and data for ReadUtilIface to create element entities on post-process mesh.

Definition at line 171 of file PostProcBrokenMeshInMoabBase.hpp.

◆ tagsToTransfer

template<typename E >
std::vector<Tag> MoFEM::PostProcBrokenMeshInMoabBase< E >::tagsToTransfer
protected

Set of tags on mesh to transfer to postprocessing mesh

Definition at line 174 of file PostProcBrokenMeshInMoabBase.hpp.


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