v0.15.5
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, 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 ()
 
auto getPostProcMeshPcommPtr ()
 

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 >
Examples
mofem/tools/mesh_smoothing.cpp, mofem/tutorials/vec-2_nonlinear_elasticity/nonlinear_elastic.cpp, mofem/tutorials/vec-9_arc_length/arc_length.cpp, mofem/users_modules/adolc-plasticity/adolc_plasticity.cpp, and nonlinear_elastic.cpp.

Definition at line 94 of file PostProcBrokenMeshInMoabBase.hpp.

Constructor & Destructor Documentation

◆ PostProcBrokenMeshInMoabBase() [1/2]

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

Definition at line 290 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ PostProcBrokenMeshInMoabBase() [2/2]

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

Definition at line 295 of file PostProcBrokenMeshInMoabBase.hpp.

298 : PostProcBrokenMeshInMoabBase(m_field, opts_prefix) {
299 coreMeshPtr = core_mesh_ptr;
300}
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field, std::string opts_prefix="")

◆ ~PostProcBrokenMeshInMoabBase()

Definition at line 303 of file PostProcBrokenMeshInMoabBase.hpp.

303 {
305}

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

652 {
653 return mapGaussPts;
654}

◆ getMaxLevel()

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

Determine refinement level based on fields approx ordre.

level = (order - 1) / 2

Returns
int

Reimplemented in EshelbianPlasticity::ContactTree.

Definition at line 198 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ getPostProcElements()

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

Get postprocessing elements.

Returns
auto&

Definition at line 664 of file PostProcBrokenMeshInMoabBase.hpp.

◆ getPostProcMesh()

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

Get postprocessing mesh.

Returns
moab::Interface&

Definition at line 656 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ getPostProcMeshPcommPtr()

template<typename E >
auto MoFEM::PostProcBrokenMeshInMoabBase< E >::getPostProcMeshPcommPtr ( )
inlineprotected

Definition at line 184 of file PostProcBrokenMeshInMoabBase.hpp.

184 {
186 }
static auto getOrCreate(moab::Core *core_mesh_ptr)

◆ getRule()

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

Definition at line 307 of file PostProcBrokenMeshInMoabBase.hpp.

307 {
308 return -1;
309};

◆ postProcess()

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

Definition at line 638 of file PostProcBrokenMeshInMoabBase.hpp.

638 {
640
641 auto pcomm_post_proc_mesh = getPostProcMeshPcommPtr();
642 if (!pcomm_post_proc_mesh)
643 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
644 "ParallelComm not allocated");
645
647 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
648
650}
#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 MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcPostProc ( )
protectedvirtual

Definition at line 554 of file PostProcBrokenMeshInMoabBase.hpp.

554 {
556
557 auto update_elements = [&]() {
558 ReadUtilIface *iface;
559 CHKERR getPostProcMesh().query_interface(iface);
561
562 Range ents;
563 for (auto &m : refElementsMap) {
564 if (m.second->nbEles) {
565 MOFEM_TAG_AND_LOG("SELF", Sev::noisy, "PostProc")
566 << "Update < " << moab::CN::EntityTypeName(m.first)
567 << " number of processed " << m.second->countEle;
568 CHKERR iface->update_adjacencies(
569 m.second->startingEleHandle, m.second->countEle,
570 m.second->levelRef[0].size2(), m.second->eleConn);
571 ents.merge(Range(m.second->startingEleHandle,
572 m.second->startingEleHandle + m.second->countEle - 1));
573 }
574 }
575
577 };
578
579 auto remove_obsolete_entities = [&]() {
581 Range ents, adj;
582 for (auto &m : refElementsMap) {
583 if (m.second->nbEles) {
584 ents.merge(Range(m.second->startingEleHandle,
585 m.second->startingEleHandle + m.second->countEle - 1));
586 const int dim = moab::CN::Dimension(m.first);
587 for (auto d = 1; d != dim; ++d) {
588 CHKERR getPostProcMesh().get_adjacencies(ents, d, false, adj,
589 moab::Interface::UNION);
590 }
591 }
592 }
593 CHKERR getPostProcMesh().delete_entities(adj);
595 };
596
597 auto set_proc_tags = [&]() {
599 ParallelComm *pcomm_post_proc_mesh =
600 ParallelComm::get_pcomm(&(getPostProcMesh()), MYPCOMM_INDEX);
601 if (pcomm_post_proc_mesh) {
602 Range ents;
603 for (auto &m : refElementsMap) {
604 if (m.second->nbEles) {
605 ents.merge(
606 Range(m.second->startingEleHandle,
607 m.second->startingEleHandle + m.second->countEle - 1));
608 }
609 }
610
611 int rank = E::mField.get_comm_rank();
612 CHKERR getPostProcMesh().tag_clear_data(pcomm_post_proc_mesh->part_tag(),
613 ents, &rank);
614 }
616 };
617
618 CHKERR update_elements();
619 CHKERR remove_obsolete_entities();
620 CHKERR set_proc_tags();
621
623}
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define MYPCOMM_INDEX
default communicator number PCOMM
#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
auto & getPostProcMesh()
Get postprocessing mesh.
std::map< EntityType, PostProcGenerateRefMeshPtr > refElementsMap

◆ preProcess()

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

Generate vertices and elements.

Returns
MoFEMErrorCode

Definition at line 626 of file PostProcBrokenMeshInMoabBase.hpp.

626 {
628
629 CHKERR getPostProcMesh().delete_mesh();
630
631
633
635}

◆ preProcPostProc()

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

Definition at line 439 of file PostProcBrokenMeshInMoabBase.hpp.

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

Definition at line 312 of file PostProcBrokenMeshInMoabBase.hpp.

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

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

192 {
194 tagsToTransfer.swap(tags_to_transfer);
196}

◆ transferTags()

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

transfer tags from mesh to post-process mesh

Definition at line 213 of file PostProcBrokenMeshInMoabBase.hpp.

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