12 #ifndef __POSTPROCBROKENMESHINMOABBASE_HPP
13 #define __POSTPROCBROKENMESHINMOABBASE_HPP
25 std::vector<MatrixDouble>
28 std::vector<MatrixDouble>
54 boost::shared_ptr<PostProcGenerateRefMeshBase>;
97 std::string opts_prefix =
"");
99 boost::shared_ptr<moab::Core> core_mesh_ptr,
100 std::string opts_prefix =
"");
165 boost::shared_ptr<moab::Core>
coreMeshPtr = boost::make_shared<moab::Core>();
170 std::map<EntityType, PostProcGenerateRefMeshPtr>
185 template <
typename E>
187 std::vector<Tag> tags_to_transfer) {
189 tagsToTransfer.swap(tags_to_transfer);
194 auto get_element_max_dofs_order = [&]() {
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;
203 const auto dof_max_order = get_element_max_dofs_order();
204 return (dof_max_order > 0) ? (dof_max_order - 1) / 2 : 0;
207 template <
typename E>
211 auto &calc_mesh = this->mField.get_moab();
213 auto name = [&](
auto tag) {
219 auto data_type = [&](
auto tag) {
220 moab::DataType data_type;
226 auto type = [&](
auto tag) {
232 auto length = [&](
auto tag) {
234 CHK_MOAB_THROW(calc_mesh.tag_get_length(tag, length),
"get length ");
238 auto default_value = [&](
auto tag) {
241 CHK_MOAB_THROW(calc_mesh.tag_get_default_value(tag, def_val, size),
242 "get default tag value");
246 auto data_ptr = [&](
auto tag) {
247 const void *tag_data;
248 auto core_ent = this->getFEEntityHandle();
249 CHK_MOAB_THROW(calc_mesh.tag_get_by_ptr(tag, &core_ent, 1, &tag_data),
254 for (
auto tag : tagsToTransfer) {
256 CHKERR getPostProcMesh().tag_get_handle(
257 name(tag).c_str(), length(tag), data_type(tag), tag_postproc,
258 type(tag) | MB_TAG_CREAT, default_value(tag));
259 CHKERR getPostProcMesh().tag_clear_data(tag_postproc, postProcElements,
266 template <
typename E>
269 :
E(m_field), optionsPrefix(opts_prefix) {}
271 template <
typename E>
274 std::string opts_prefix)
276 coreMeshPtr = core_mesh_ptr;
279 template <
typename E>
281 ParallelComm *pcomm_post_proc_mesh =
283 if (pcomm_post_proc_mesh != NULL)
284 delete pcomm_post_proc_mesh;
291 template <
typename E>
300 ref_ele = refElementsMap.at(
type);
301 }
catch (
const out_of_range &e) {
304 "Generation of reference elements for type <%s> is not implemented",
305 moab::CN::EntityTypeName(
type));
308 auto set_gauss_pts = [&](
auto &level_gauss_pts_on_ref_mesh,
auto &level_ref,
309 auto &level_shape_functions,
311 auto start_vert_handle,
auto start_ele_handle,
312 auto &verts_array,
auto &conn,
auto &ver_count,
318 size_t level = getMaxLevel();
319 level = std::min(level, level_gauss_pts_on_ref_mesh.size() - 1);
321 auto &level_ref_gauss_pts = level_gauss_pts_on_ref_mesh[level];
322 auto &level_ref_ele = level_ref[level];
323 auto &shape_functions = level_shape_functions[level];
324 E::gaussPts.resize(level_ref_gauss_pts.size1(), level_ref_gauss_pts.size2(),
326 noalias(E::gaussPts) = level_ref_gauss_pts;
328 const auto fe_ent = E::numeredEntFiniteElementPtr->getEnt();
329 auto get_fe_coords = [&]() {
333 E::mField.get_moab().get_connectivity(fe_ent, conn, num_nodes,
true),
334 "error get connectivity");
337 E::mField.get_moab().get_coords(conn, num_nodes, &*coords.begin()),
338 "error get coordinates");
342 auto coords = get_fe_coords();
344 const int num_nodes = level_ref_gauss_pts.size2();
345 mapGaussPts.resize(level_ref_gauss_pts.size2());
349 &*shape_functions.data().begin());
351 &verts_array[0][ver_count], &verts_array[1][ver_count],
352 &verts_array[2][ver_count]);
353 for (
int gg = 0; gg != num_nodes; ++gg, ++ver_count) {
355 mapGaussPts[gg] = start_vert_handle + ver_count;
357 auto set_float_precision = [](
const double x) {
358 if (std::abs(x) < std::numeric_limits<float>::epsilon())
365 auto t_ele_coords = getFTensor1FromArray<3, 3>(coords);
366 for (
int nn = 0; nn != CN::VerticesPerEntity(
type); ++nn) {
367 t_coords(
i) += t_n * t_ele_coords(
i);
372 for (
auto ii : {0, 1, 2})
373 t_coords(ii) = set_float_precision(t_coords(ii));
379 int def_in_the_loop = -1;
380 CHKERR getPostProcMesh().tag_get_handle(
381 "NB_IN_THE_LOOP", 1, MB_TYPE_INTEGER,
th, MB_TAG_CREAT | MB_TAG_SPARSE,
384 postProcElements.clear();
385 const int num_el = level_ref_ele.size1();
386 const int num_nodes_on_ele = level_ref_ele.size2();
387 auto start_e = start_ele_handle + ele_count;
388 postProcElements =
Range(start_e, start_e + num_el - 1);
389 for (
auto tt = 0; tt != level_ref_ele.size1(); ++tt, ++ele_count) {
390 for (
int nn = 0; nn != num_nodes_on_ele; ++nn) {
391 conn[num_nodes_on_ele * ele_count + nn] =
392 mapGaussPts[level_ref_ele(tt, nn)];
396 const int n_in_the_loop = E::nInTheLoop;
397 CHKERR getPostProcMesh().tag_clear_data(
th, postProcElements,
406 ref_ele->levelGaussPtsOnRefMesh, ref_ele->levelRef,
407 ref_ele->levelShapeFunctions,
409 ref_ele->startingVertEleHandle, ref_ele->startingEleHandle,
410 ref_ele->verticesOnEleArrays, ref_ele->eleConn, ref_ele->countVertEle,
418 template <
typename E>
422 auto get_ref_ele = [&](
const EntityType
type) {
425 auto it = refElementsMap.find(
type);
426 if (it != refElementsMap.end()) {
427 ref_ele_ptr = it->second;
431 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTET>>();
434 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBHEX>>();
437 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
440 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBQUAD>>();
443 ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBEDGE>>();
447 <<
"Generation of reference elements for type < "
448 << moab::CN::EntityTypeName(
type) <<
" > is not implemented";
454 "Error when generating reference element");
456 refElementsMap[
type] = ref_ele_ptr;
462 auto fe_ptr = this->problemPtr->numeredFiniteElementsPtr;
465 fe_ptr->template get<Composite_Name_And_Part_mi_tag>().lower_bound(
466 boost::make_tuple(this->getFEName(), this->getLoFERank()));
468 fe_ptr->template get<Composite_Name_And_Part_mi_tag>().upper_bound(
469 boost::make_tuple(this->getFEName(), this->getHiFERank()));
471 const int number_of_ents_in_the_loop = this->getLoopSize();
472 if (std::distance(miit, hi_miit) != number_of_ents_in_the_loop) {
474 "Wrong size of indicies. Inconsistent size number of iterated "
475 "elements iterated by problem and from range.");
478 for (
auto &
m : refElementsMap) {
479 m.second->nbVertices = 0;
480 m.second->nbEles = 0;
481 m.second->countEle = 0;
482 m.second->countVertEle = 0;
485 for (; miit != hi_miit; ++miit) {
486 auto type = (*miit)->getEntType();
487 auto ref_ele = get_ref_ele(
type);
491 E::numeredEntFiniteElementPtr = *miit;
493 if (E::exeTestHook) {
494 add = E::exeTestHook(
this);
498 size_t level = getMaxLevel();
499 level = std::min(level, ref_ele->levelGaussPtsOnRefMesh.size() - 1);
500 ref_ele->nbVertices += ref_ele->levelGaussPtsOnRefMesh[level].size2();
501 ref_ele->nbEles += ref_ele->levelRef[level].size1();
505 auto alloc_vertices_and_elements_on_post_proc_mesh = [&]() {
508 ReadUtilIface *iface;
509 CHKERR getPostProcMesh().query_interface(iface);
511 for (
auto &
m : refElementsMap) {
512 if (
m.second->nbEles) {
513 CHKERR iface->get_node_coords(3,
m.second->nbVertices, 0,
514 m.second->startingVertEleHandle,
515 m.second->verticesOnEleArrays);
516 CHKERR iface->get_element_connect(
517 m.second->nbEles,
m.second->levelRef[0].size2(),
m.first, 0,
518 m.second->startingEleHandle,
m.second->eleConn);
520 m.second->countEle = 0;
521 m.second->countVertEle = 0;
528 CHKERR alloc_vertices_and_elements_on_post_proc_mesh();
533 template <
typename E>
537 auto update_elements = [&]() {
538 ReadUtilIface *iface;
539 CHKERR getPostProcMesh().query_interface(iface);
543 for (
auto &
m : refElementsMap) {
544 if (
m.second->nbEles) {
546 <<
"Update < " << moab::CN::EntityTypeName(
m.first)
547 <<
" number of processed " <<
m.second->countEle;
548 CHKERR iface->update_adjacencies(
549 m.second->startingEleHandle,
m.second->countEle,
550 m.second->levelRef[0].size2(),
m.second->eleConn);
551 ents.merge(
Range(
m.second->startingEleHandle,
552 m.second->startingEleHandle +
m.second->countEle - 1));
559 auto remove_obsolete_entities = [&]() {
562 for (
auto &
m : refElementsMap) {
563 if (
m.second->nbEles) {
564 ents.merge(
Range(
m.second->startingEleHandle,
565 m.second->startingEleHandle +
m.second->countEle - 1));
566 const int dim = moab::CN::Dimension(
m.first);
567 for (
auto d = 1;
d != dim; ++
d) {
568 CHKERR getPostProcMesh().get_adjacencies(ents,
d,
false, adj,
569 moab::Interface::UNION);
573 CHKERR getPostProcMesh().delete_entities(adj);
577 auto set_proc_tags = [&]() {
579 ParallelComm *pcomm_post_proc_mesh =
580 ParallelComm::get_pcomm(&(getPostProcMesh()),
MYPCOMM_INDEX);
581 if (pcomm_post_proc_mesh) {
583 for (
auto &
m : refElementsMap) {
584 if (
m.second->nbEles) {
586 Range(
m.second->startingEleHandle,
587 m.second->startingEleHandle +
m.second->countEle - 1));
591 int rank = E::mField.get_comm_rank();
592 CHKERR getPostProcMesh().tag_clear_data(pcomm_post_proc_mesh->part_tag(),
599 CHKERR remove_obsolete_entities();
605 template <
typename E>
609 ParallelComm *pcomm_post_proc_mesh =
611 if (pcomm_post_proc_mesh != NULL)
612 delete pcomm_post_proc_mesh;
613 CHKERR getPostProcMesh().delete_mesh();
614 pcomm_post_proc_mesh =
615 new ParallelComm(&(getPostProcMesh()), PETSC_COMM_WORLD);
622 template <
typename E>
626 ParallelComm *pcomm_post_proc_mesh =
628 if(pcomm_post_proc_mesh ==
nullptr)
631 CHKERR postProcPostProc();
632 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
645 return post_proc_mesh_interface;
648 template <
typename E>
650 return postProcElements;
653 template <
typename E>
657 ParallelComm *pcomm_post_proc_mesh =
659 if (pcomm_post_proc_mesh == NULL)
661 "ParallelComm not allocated");
662 CHKERR getPostProcMesh().write_file(file_name.c_str(),
"MOAB",
663 "PARALLEL=WRITE_PART");
697 template <
int DIM1,
int DIM2>
700 using DataMapVec = std::map<std::string, boost::shared_ptr<VectorDouble>>;
701 using DataMapMat = std::map<std::string, boost::shared_ptr<MatrixDouble>>;
715 std::vector<EntityHandle> &map_gauss_pts,
720 postProcMesh(post_proc_mesh), mapGaussPts(map_gauss_pts),
721 dataMapScalar(data_map_scalar), dataMapVec(data_map_vec),
722 dataMapMat(data_map_mat), dataMapSymmMat(data_symm_map_mat) {
724 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE],
false);
738 template <
int DIM1,
int DIM2>
744 std::array<double, 9> def;
745 std::fill(def.begin(), def.end(), 0);
747 auto get_tag = [&](
const std::string name,
size_t size) {
749 CHKERR postProcMesh.tag_get_handle(name.c_str(), size, MB_TYPE_DOUBLE,
th,
750 MB_TAG_CREAT | MB_TAG_SPARSE,
759 for (
size_t r = 0;
r !=
DIM1; ++
r)
766 for (
size_t r = 0;
r !=
DIM1; ++
r)
767 for (
size_t c = 0;
c !=
DIM2; ++
c)
774 for (
size_t r = 0;
r !=
DIM1; ++
r)
775 for (
size_t c = 0;
c !=
DIM1; ++
c)
786 auto set_float_precision = [](
const double x) {
787 if (std::abs(x) < std::numeric_limits<float>::epsilon())
794 for (
auto &
v : mat.data())
795 v = set_float_precision(
v);
796 return postProcMesh.tag_set_data(
th, &mapGaussPts[gg], 1,
797 &*mat.data().begin());
800 for (
auto &
m : dataMapScalar) {
802 auto th = get_tag(
m.first, 1);
804 auto nb_integration_pts = getGaussPts().size2();
805 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
806 CHKERR set_tag(
th, gg, set_scalar(t_scl));
812 for (
auto &
m : dataMapVec) {
814 auto th = get_tag(
m.first, 3);
815 auto t_vec = getFTensor1FromMat<DIM1>(*
m.second);
816 auto nb_integration_pts = getGaussPts().size2();
817 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
818 CHKERR set_tag(
th, gg, set_vector_3d(t_vec));
824 for (
auto &
m : dataMapMat) {
826 auto th = get_tag(
m.first, 9);
827 auto t_mat = getFTensor2FromMat<DIM1, DIM2>(*
m.second);
828 auto nb_integration_pts = getGaussPts().size2();
829 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
830 CHKERR set_tag(
th, gg, set_matrix_3d(t_mat));
836 for (
auto &
m : dataMapSymmMat) {
838 auto th = get_tag(
m.first, 9);
839 auto t_mat = getFTensor2SymmetricFromMat<DIM1>(*
m.second);
840 auto nb_integration_pts = getGaussPts().size2();
841 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
842 CHKERR set_tag(
th, gg, set_matrix_symm_3d(t_mat));
853 template <
typename E>
859 std::string opts_prefix =
"")
864 ParallelComm *pcomm_post_proc_mesh =
865 ParallelComm::get_pcomm(&this->getPostProcMesh(),
MYPCOMM_INDEX);
866 if (pcomm_post_proc_mesh != NULL)
867 delete pcomm_post_proc_mesh;
868 CHKERR this->getPostProcMesh().delete_mesh();
869 pcomm_post_proc_mesh =
870 new ParallelComm(&(this->getPostProcMesh()), PETSC_COMM_WORLD);
882 template <
typename E>
888 std::string opts_prefix =
"")
900 template <
typename E>
906 std::string opts_prefix =
"")
913 ParallelComm *pcomm_post_proc_mesh =
914 ParallelComm::get_pcomm(&this->getPostProcMesh(),
MYPCOMM_INDEX);
915 if (pcomm_post_proc_mesh ==
nullptr)
917 "PComm not allocated");
918 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
943 template <
typename E>
957 #endif //__POSTPROCBROKENMESHINMOABBASE_HPP