v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
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)
 
virtual ~PostProcBrokenMeshInMoabBase ()
 
auto & getMapGaussPts ()
 Get vector of vectors associated to integration points. More...
 
auto & getPostProcMesh ()
 Get postprocessing mesh. More...
 
auto & getPostProcElements ()
 Get postprocessing elements. More...
 
MoFEMErrorCode writeFile (const std::string file_name)
 wrote results in (MOAB) format, use "file_name.h5m" More...
 
MoFEMErrorCode setTagsToTransfer (std::vector< Tag > tags_to_transfer)
 Set tags to be transferred to post-processing mesh. More...
 

Protected Member Functions

MoFEMErrorCode setGaussPts (int order)
 
MoFEMErrorCode preProcess ()
 Generate vertices and elements. More...
 
MoFEMErrorCode postProcess ()
 
int getRule (int order)
 
virtual int getMaxLevel () const
 Determine refinement level based on fields approx ordre. More...
 
 PostProcBrokenMeshInMoabBase (MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr)
 
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< Tag > tagsToTransfer
 

Detailed Description

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

Definition at line 95 of file PostProcBrokenMeshInMoabBase.hpp.

Constructor & Destructor Documentation

◆ PostProcBrokenMeshInMoabBase() [1/2]

template<typename E >
MoFEM::PostProcBrokenMeshInMoabBase< E >::PostProcBrokenMeshInMoabBase ( MoFEM::Interface m_field)

Definition at line 266 of file PostProcBrokenMeshInMoabBase.hpp.

268 : E(m_field) {}

◆ ~PostProcBrokenMeshInMoabBase()

Definition at line 278 of file PostProcBrokenMeshInMoabBase.hpp.

278 {
279 ParallelComm *pcomm_post_proc_mesh =
280 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
281 if (pcomm_post_proc_mesh != NULL)
282 delete pcomm_post_proc_mesh;
283}
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
auto & getPostProcMesh()
Get postprocessing mesh.

◆ PostProcBrokenMeshInMoabBase() [2/2]

template<typename E >
MoFEM::PostProcBrokenMeshInMoabBase< E >::PostProcBrokenMeshInMoabBase ( MoFEM::Interface m_field,
boost::shared_ptr< moab::Core >  core_mesh_ptr 
)
protected

Definition at line 271 of file PostProcBrokenMeshInMoabBase.hpp.

274 coreMeshPtr = core_mesh_ptr;
275}
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field)
boost::shared_ptr< moab::Core > coreMeshPtr

Member Function Documentation

◆ getMapGaussPts()

template<typename E >
auto & MoFEM::PostProcBrokenMeshInMoabBase< E >::getMapGaussPts
inline

Get vector of vectors associated to integration points.

Returns
std::vector<EntityHandle>&

Definition at line 636 of file PostProcBrokenMeshInMoabBase.hpp.

636 {
637 return mapGaussPts;
638}
std::vector< EntityHandle > mapGaussPts

◆ getMaxLevel()

template<typename E >
int MoFEM::PostProcBrokenMeshInMoabBase< E >::getMaxLevel
protectedvirtual

Determine refinement level based on fields approx ordre.

level = (order - 1) / 2

Returns
int

Definition at line 192 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ getPostProcElements()

template<typename E >
auto & MoFEM::PostProcBrokenMeshInMoabBase< E >::getPostProcElements
inline

Get postprocessing elements.

Returns
auto&

Definition at line 648 of file PostProcBrokenMeshInMoabBase.hpp.

◆ getPostProcMesh()

template<typename E >
auto & MoFEM::PostProcBrokenMeshInMoabBase< E >::getPostProcMesh
inline

Get postprocessing mesh.

Returns
moab::Interface&

Definition at line 640 of file PostProcBrokenMeshInMoabBase.hpp.

640 {
641 if (!coreMeshPtr)
642 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Core mesh not set");
643 moab::Interface &post_proc_mesh_interface = *coreMeshPtr;
644 return post_proc_mesh_interface;
645}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:595
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31

◆ getRule()

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

Definition at line 285 of file PostProcBrokenMeshInMoabBase.hpp.

285 {
286 return -1;
287};

◆ postProcess()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcess
protected

Definition at line 622 of file PostProcBrokenMeshInMoabBase.hpp.

622 {
624
625 ParallelComm *pcomm_post_proc_mesh =
626 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
627 if(pcomm_post_proc_mesh == nullptr)
628 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
629
631 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
632
634}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
virtual MoFEMErrorCode postProcPostProc()

◆ postProcPostProc()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcPostProc
protectedvirtual

Definition at line 533 of file PostProcBrokenMeshInMoabBase.hpp.

533 {
535
536 auto update_elements = [&]() {
537 ReadUtilIface *iface;
538 CHKERR getPostProcMesh().query_interface(iface);
540
541 Range ents;
542 for (auto &m : refElementsMap) {
543 if (m.second->nbEles) {
544 MOFEM_TAG_AND_LOG("SELF", Sev::noisy, "PostProc")
545 << "Update < " << moab::CN::EntityTypeName(m.first)
546 << " number of processed " << m.second->countEle;
547 CHKERR iface->update_adjacencies(
548 m.second->startingEleHandle, m.second->countEle,
549 m.second->levelRef[0].size2(), m.second->eleConn);
550 ents.merge(Range(m.second->startingEleHandle,
551 m.second->startingEleHandle + m.second->countEle - 1));
552 }
553 }
554
556 };
557
558 auto remove_obsolete_entities = [&]() {
560 Range ents, adj;
561 for (auto &m : refElementsMap) {
562 if (m.second->nbEles) {
563 ents.merge(Range(m.second->startingEleHandle,
564 m.second->startingEleHandle + m.second->countEle - 1));
565 const int dim = moab::CN::Dimension(m.first);
566 for (auto d = 1; d != dim; ++d) {
567 CHKERR getPostProcMesh().get_adjacencies(ents, d, false, adj,
568 moab::Interface::UNION);
569 }
570 }
571 }
572 CHKERR getPostProcMesh().delete_entities(adj);
574 };
575
576 auto set_proc_tags = [&]() {
578 ParallelComm *pcomm_post_proc_mesh =
579 ParallelComm::get_pcomm(&(getPostProcMesh()), MYPCOMM_INDEX);
580 if (pcomm_post_proc_mesh) {
581 Range ents;
582 for (auto &m : refElementsMap) {
583 if (m.second->nbEles) {
584 ents.merge(
585 Range(m.second->startingEleHandle,
586 m.second->startingEleHandle + m.second->countEle - 1));
587 }
588 }
589
590 int rank = E::mField.get_comm_rank();
591 CHKERR getPostProcMesh().tag_clear_data(pcomm_post_proc_mesh->part_tag(),
592 ents, &rank);
593 }
595 };
596
597 CHKERR update_elements();
598 CHKERR remove_obsolete_entities();
599 CHKERR set_proc_tags();
600
602}
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
#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
const int dim
FTensor::Index< 'm', SPACE_DIM > m
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
std::map< EntityType, PostProcGenerateRefMeshPtr > refElementsMap

◆ preProcess()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::preProcess
protected

Generate vertices and elements.

Returns
MoFEMErrorCode

Definition at line 605 of file PostProcBrokenMeshInMoabBase.hpp.

605 {
607
608 ParallelComm *pcomm_post_proc_mesh =
609 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
610 if (pcomm_post_proc_mesh != NULL)
611 delete pcomm_post_proc_mesh;
612 CHKERR getPostProcMesh().delete_mesh();
613 pcomm_post_proc_mesh =
614 new ParallelComm(&(getPostProcMesh()), PETSC_COMM_WORLD);
615
617
619}
virtual MoFEMErrorCode preProcPostProc()

◆ preProcPostProc()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::preProcPostProc
protectedvirtual

Definition at line 417 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ setGaussPts()

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

Definition at line 290 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ setTagsToTransfer()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::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 185 of file PostProcBrokenMeshInMoabBase.hpp.

186 {
188 tagsToTransfer.swap(tags_to_transfer);
190}
std::vector< Tag > tagsToTransfer

◆ transferTags()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::transferTags
protectedvirtual

transfer tags from mesh to post-process mesh

Definition at line 207 of file PostProcBrokenMeshInMoabBase.hpp.

207 {
209
210 auto &calc_mesh = this->mField.get_moab();
211
212 auto name = [&](auto tag) {
213 std::string name;
214 CHK_MOAB_THROW(calc_mesh.tag_get_name(tag, name), "get name");
215 return name;
216 };
217
218 auto data_type = [&](auto tag) {
219 moab::DataType data_type;
220 CHK_MOAB_THROW(calc_mesh.tag_get_data_type(tag, data_type),
221 "get data type");
222 return data_type;
223 };
224
225 auto type = [&](auto tag) {
226 moab::TagType type;
227 CHK_MOAB_THROW(calc_mesh.tag_get_type(tag, type), "get tag type");
228 return type;
229 };
230
231 auto length = [&](auto tag) {
232 int length;
233 CHK_MOAB_THROW(calc_mesh.tag_get_length(tag, length), "get length ");
234 return length;
235 };
236
237 auto default_value = [&](auto tag) {
238 const void *def_val;
239 int size;
240 CHK_MOAB_THROW(calc_mesh.tag_get_default_value(tag, def_val, size),
241 "get default tag value");
242 return def_val;
243 };
244
245 auto data_ptr = [&](auto tag) {
246 const void *tag_data;
247 auto core_ent = this->getFEEntityHandle();
248 CHK_MOAB_THROW(calc_mesh.tag_get_by_ptr(tag, &core_ent, 1, &tag_data),
249 "get tag data");
250 return tag_data;
251 };
252
253 for (auto tag : tagsToTransfer) {
254 Tag tag_postproc;
255 CHKERR getPostProcMesh().tag_get_handle(
256 name(tag).c_str(), length(tag), data_type(tag), tag_postproc,
257 type(tag) | MB_TAG_CREAT, default_value(tag));
258 CHKERR getPostProcMesh().tag_clear_data(tag_postproc, postProcElements,
259 data_ptr(tag));
260 }
261
263};

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 162 of file PostProcBrokenMeshInMoabBase.hpp.

◆ mapGaussPts

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

Definition at line 167 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: