v0.15.0
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 ()
 

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/tutorials/vec-2/nonlinear_elastic.cpp, mofem/tutorials/vec-9/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 285 of file PostProcBrokenMeshInMoabBase.hpp.

287 : 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 290 of file PostProcBrokenMeshInMoabBase.hpp.

293 : PostProcBrokenMeshInMoabBase(m_field, opts_prefix) {
294 coreMeshPtr = core_mesh_ptr;
295}
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field, std::string opts_prefix="")

◆ ~PostProcBrokenMeshInMoabBase()

Definition at line 298 of file PostProcBrokenMeshInMoabBase.hpp.

298 {
299 ParallelComm *pcomm_post_proc_mesh =
300 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
301 if (pcomm_post_proc_mesh != NULL)
302 delete pcomm_post_proc_mesh;
303}
#define MYPCOMM_INDEX
default communicator number PCOMM
auto & getPostProcMesh()
Get postprocessing mesh.

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

655 {
656 return mapGaussPts;
657}

◆ 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 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 & MoFEM::PostProcBrokenMeshInMoabBase< E >::getPostProcElements ( )
inline

Get postprocessing elements.

Returns
auto&

Definition at line 667 of file PostProcBrokenMeshInMoabBase.hpp.

◆ getPostProcMesh()

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

Get postprocessing mesh.

Returns
moab::Interface&

Definition at line 659 of file PostProcBrokenMeshInMoabBase.hpp.

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

◆ getRule()

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

Definition at line 305 of file PostProcBrokenMeshInMoabBase.hpp.

305 {
306 return -1;
307};

◆ postProcess()

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

Definition at line 641 of file PostProcBrokenMeshInMoabBase.hpp.

641 {
643
644 ParallelComm *pcomm_post_proc_mesh =
645 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
646 if(pcomm_post_proc_mesh == nullptr)
647 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
648
650 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
651
653}
#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 552 of file PostProcBrokenMeshInMoabBase.hpp.

552 {
554
555 auto update_elements = [&]() {
556 ReadUtilIface *iface;
557 CHKERR getPostProcMesh().query_interface(iface);
559
560 Range ents;
561 for (auto &m : refElementsMap) {
562 if (m.second->nbEles) {
563 MOFEM_TAG_AND_LOG("SELF", Sev::noisy, "PostProc")
564 << "Update < " << moab::CN::EntityTypeName(m.first)
565 << " number of processed " << m.second->countEle;
566 CHKERR iface->update_adjacencies(
567 m.second->startingEleHandle, m.second->countEle,
568 m.second->levelRef[0].size2(), m.second->eleConn);
569 ents.merge(Range(m.second->startingEleHandle,
570 m.second->startingEleHandle + m.second->countEle - 1));
571 }
572 }
573
575 };
576
577 auto remove_obsolete_entities = [&]() {
579 Range ents, adj;
580 for (auto &m : refElementsMap) {
581 if (m.second->nbEles) {
582 ents.merge(Range(m.second->startingEleHandle,
583 m.second->startingEleHandle + m.second->countEle - 1));
584 const int dim = moab::CN::Dimension(m.first);
585 for (auto d = 1; d != dim; ++d) {
586 CHKERR getPostProcMesh().get_adjacencies(ents, d, false, adj,
587 moab::Interface::UNION);
588 }
589 }
590 }
591 CHKERR getPostProcMesh().delete_entities(adj);
593 };
594
595 auto set_proc_tags = [&]() {
597 ParallelComm *pcomm_post_proc_mesh =
598 ParallelComm::get_pcomm(&(getPostProcMesh()), MYPCOMM_INDEX);
599 if (pcomm_post_proc_mesh) {
600 Range ents;
601 for (auto &m : refElementsMap) {
602 if (m.second->nbEles) {
603 ents.merge(
604 Range(m.second->startingEleHandle,
605 m.second->startingEleHandle + m.second->countEle - 1));
606 }
607 }
608
609 int rank = E::mField.get_comm_rank();
610 CHKERR getPostProcMesh().tag_clear_data(pcomm_post_proc_mesh->part_tag(),
611 ents, &rank);
612 }
614 };
615
616 CHKERR update_elements();
617 CHKERR remove_obsolete_entities();
618 CHKERR set_proc_tags();
619
621}
#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 MoFEM::PostProcBrokenMeshInMoabBase< E >::preProcess ( )
protected

Generate vertices and elements.

Returns
MoFEMErrorCode

Definition at line 624 of file PostProcBrokenMeshInMoabBase.hpp.

624 {
626
627 ParallelComm *pcomm_post_proc_mesh =
628 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
629 if (pcomm_post_proc_mesh != NULL)
630 delete pcomm_post_proc_mesh;
631 CHKERR getPostProcMesh().delete_mesh();
632 pcomm_post_proc_mesh =
633 new ParallelComm(&(getPostProcMesh()), PETSC_COMM_WORLD);
634
636
638}

◆ preProcPostProc()

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

Definition at line 437 of file PostProcBrokenMeshInMoabBase.hpp.

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

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

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

◆ transferTags()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::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>::min() ? 0. : v;
259 for (auto &v : tag_data_vec)
260 v = v > std::numeric_limits<float>::max()
261 ? std::numeric_limits<float>::max()
262 : v;
263 for (auto &v : tag_data_vec)
264 v = v < std::numeric_limits<float>::lowest()
265 ? std::numeric_limits<float>::lowest()
266 : v;
267 return static_cast<const void *>(tag_data_vec.data());
268 }
269 return tag_data;
270 };
271
272 for (auto tag : tagsToTransfer) {
273 Tag tag_postproc;
274 CHKERR getPostProcMesh().tag_get_handle(
275 name(tag).c_str(), length(tag), data_type(tag), tag_postproc,
276 type(tag) | MB_TAG_CREAT, default_value(tag));
277 CHKERR getPostProcMesh().tag_clear_data(tag_postproc, postProcElements,
278 data_ptr(tag));
279 }
280
282};
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: