v0.15.0
Loading...
Searching...
No Matches
PostProcBrokenMeshInMoabBase.hpp
Go to the documentation of this file.
1/**
2 * @file PostProcBrokenMeshInMoabBase.hpp
3 * @author your name (you@domain.com)
4 * @brief
5 * @version 0.1
6 * @date 2022-09-14
7 *
8 * @copyright Copyright (c) 2022
9 *
10 */
11
12#ifndef __POSTPROCBROKENMESHINMOABBASE_HPP
13#define __POSTPROCBROKENMESHINMOABBASE_HPP
14
15namespace MoFEM {
16
17/**
18 * Each element is subdivided on smaller elements, i.e. a reference mesh on
19 * single element is created. Nodes of such reference mesh are used as
20 * integration points at which field values are calculated and to
21 * each node a "moab" tag is attached to store those values.
22 */
24
25 std::vector<MatrixDouble>
26 levelShapeFunctions; //< values of shape functions on refences element
27 // refinement levels
28 std::vector<MatrixDouble>
29 levelGaussPtsOnRefMesh; //< gauss points at refinement levels
30 std::vector<ublas::matrix<int>> levelRef; //< connectivity at refinement level
31
32 EntityHandle startingVertEleHandle; //< starting handle of first vertex
33 std::vector<double *> verticesOnEleArrays; //< array of vertices coordinate
34 EntityHandle startingEleHandle; //< starting handle of first element
35 EntityHandle *eleConn; //< array of elements connectivity
36
37 int countEle; //< count elements
38 int countVertEle; //< count vertices on the mesh
39
40 int nbVertices; //< numer of vertices on the mesh
41 int nbEles; //< number of elements on the mesh
42
44 virtual ~PostProcGenerateRefMeshBase() = default;
45
46 virtual MoFEMErrorCode getOptions(std::string prefix); //< get options from command
48
49 PetscBool hoNodes; //< if true mid nodes are added
50 int defMaxLevel; //< default max number of refinement levels
51};
52
54 boost::shared_ptr<PostProcGenerateRefMeshBase>;
55
56/**
57 * @brief Element for postprocessing. Uses MoAB to generate post-processing
58 * mesh.
59 *
60 * @tparam T Finite Element Implementation
61 */
62template <EntityType T> struct PostProcGenerateRefMesh;
63
64template <>
69
70template <>
75
76template <>
81
82template <>
87
88template <>
93
94template <typename E> struct PostProcBrokenMeshInMoabBase : public E {
95
97 std::string opts_prefix = "");
99 boost::shared_ptr<moab::Core> core_mesh_ptr,
100 std::string opts_prefix = "");
101
103
104 /**
105 * @brief Get vector of vectors associated to integration points
106 *
107 * @return std::vector<EntityHandle>&
108 */
109 inline auto &getMapGaussPts();
110
111 /**
112 * @brief Get postprocessing mesh
113 *
114 * @return moab::Interface&
115 */
116 inline auto &getPostProcMesh();
117
118 /**
119 * @brief Get postprocessing elements
120 *
121 * @return auto&
122 */
123 inline auto &getPostProcElements();
124
125 /**
126 * \brief wrote results in (MOAB) format, use "file_name.h5m"
127 * @param file_name file name (should always end with .h5m)
128 * @return error code
129 * \ingroup mofem_fs_post_proc
130 */
131 MoFEMErrorCode writeFile(const std::string file_name);
132
133 /**
134 * @brief Set tags to be transferred to post-processing mesh
135 *
136 * @param tags_to_transfer
137 * @return MoFEMErrorCode
138 */
139 MoFEMErrorCode setTagsToTransfer(std::vector<Tag> tags_to_transfer);
140
141protected:
142
144
145 /**
146 * @brief Generate vertices and elements
147 *
148 * @return MoFEMErrorCode
149 */
151
153
154 int getRule(int order);
155
156 /**
157 * @brief Determine refinement level based on fields approx ordre.
158 *
159 * level = (order - 1) / 2
160 *
161 * @return int
162 */
163 virtual int getMaxLevel() const;
164
165 boost::shared_ptr<moab::Core> coreMeshPtr = boost::make_shared<moab::Core>();
166
167 std::vector<EntityHandle> mapGaussPts;
169
170 std::map<EntityType, PostProcGenerateRefMeshPtr>
171 refElementsMap; ///< Storing data about element types, and data for
172 ///< ReadUtilIface to create element entities on
173 ///< post-process mesh.
174 std::vector<Tag> tagsToTransfer; ///< Set of tags on mesh to transfer to
175 ///< postprocessing mesh
176 std::string optionsPrefix = ""; ///< Prefix for options
177
180
181 virtual MoFEMErrorCode transferTags(); ///< transfer tags from mesh to
182 ///< post-process mesh
183};
184
185template <typename E>
187 std::vector<Tag> tags_to_transfer) {
189 tagsToTransfer.swap(tags_to_transfer);
191}
192
193template <typename E> int PostProcBrokenMeshInMoabBase<E>::getMaxLevel() const {
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};
206
207template <typename E>
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};
283
284template <typename E>
286 MoFEM::Interface &m_field, std::string opts_prefix)
287 : E(m_field), optionsPrefix(opts_prefix) {}
288
289template <typename E>
291 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
292 std::string opts_prefix)
293 : PostProcBrokenMeshInMoabBase(m_field, opts_prefix) {
294 coreMeshPtr = core_mesh_ptr;
295}
296
297template <typename E>
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}
304
305template <typename E> int PostProcBrokenMeshInMoabBase<E>::getRule(int order) {
306 return -1;
307};
308
309template <typename E>
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);
417 CHKERR transferTags();
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};
435
436template <typename E>
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}
550
551template <typename E>
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}
622
623template <typename E>
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
635 CHKERR preProcPostProc();
636
638}
639
640template <typename E>
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
649 CHKERR postProcPostProc();
650 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
651
653}
654
656 return mapGaussPts;
657}
658
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}
665
666template <typename E>
668 return postProcElements;
669}
670
671template <typename E>
673PostProcBrokenMeshInMoabBase<E>::writeFile(const std::string file_name) {
675 ParallelComm *pcomm_post_proc_mesh =
676 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
677 if (pcomm_post_proc_mesh == NULL)
678 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
679 "ParallelComm not allocated");
680 CHKERR getPostProcMesh().write_file(file_name.c_str(), "MOAB",
681 "PARALLEL=WRITE_PART");
683};
684
685template <typename E> struct PostProcBrokenMeshInMoab;
686
687template <>
689 : public PostProcBrokenMeshInMoabBase<VolumeElementForcesAndSourcesCore> {
691 VolumeElementForcesAndSourcesCore>::PostProcBrokenMeshInMoabBase;
692};
693
694template <>
696 : public PostProcBrokenMeshInMoabBase<FaceElementForcesAndSourcesCore> {
698 FaceElementForcesAndSourcesCore>::PostProcBrokenMeshInMoabBase;
699};
700
701template <>
703 : public PostProcBrokenMeshInMoabBase<EdgeElementForcesAndSourcesCore> {
705 EdgeElementForcesAndSourcesCore>::PostProcBrokenMeshInMoabBase;
706};
707
708/**
709 * @brief Post post-proc data at points from hash maps
710 *
711 * @tparam DIM1 dimension of vector in data_map_vec and first column of
712 * data_map_may
713 * @tparam DIM2 dimension of second column in data_map_mat
714 */
715template <int DIM1, int DIM2,
717struct OpPostProcMapInMoab : public O {
718
719 using DataMapVec = std::map<std::string, boost::shared_ptr<VectorDouble>>;
720 using DataMapMat = std::map<std::string, boost::shared_ptr<MatrixDouble>>;
721
722 /**
723 * @brief Construct a new OpPostProcMapInMoab object
724 *
725 * @param post_proc_mesh postprocessing mesh
726 * @param map_gauss_pts map of gauss points to nodes of postprocessing mesh
727 * @param data_map_scalar hash map of scalar values (string is name of the
728 * tag)
729 * @param data_map_vec hash map of vector values
730 * @param data_map_mat hash map of second order tensor values
731 * @param data_symm_map_mat hash map of symmetric second order tensor values
732 */
733 OpPostProcMapInMoab(moab::Interface &post_proc_mesh,
734 std::vector<EntityHandle> &map_gauss_pts,
735 DataMapVec data_map_scalar, DataMapMat data_map_vec,
736 DataMapMat data_map_mat, DataMapMat data_symm_map_mat)
737 : O(NOSPACE, O::OPSPACE), postProcMesh(post_proc_mesh),
738 mapGaussPts(map_gauss_pts), dataMapScalar(data_map_scalar),
739 dataMapVec(data_map_vec), dataMapMat(data_map_mat),
740 dataMapSymmMat(data_symm_map_mat) {
741 // Operator is only executed for vertices
742 std::fill(&O::doEntities[MBEDGE], &O::doEntities[MBMAXTYPE], false);
743 }
744 MoFEMErrorCode doWork(int side, EntityType type,
746
747private:
748 moab::Interface &postProcMesh;
749 std::vector<EntityHandle> &mapGaussPts;
754};
755
756template <int DIM1, int DIM2, typename O>
761
762 std::array<double, 9> def;
763 std::fill(def.begin(), def.end(), 0);
764
765 auto get_tag = [&](const std::string name, size_t size) {
766 Tag th;
767 CHKERR postProcMesh.tag_get_handle(name.c_str(), size, MB_TYPE_DOUBLE, th,
768 MB_TAG_CREAT | MB_TAG_SPARSE,
769 def.data());
770 return th;
771 };
772
773 MatrixDouble3by3 mat(3, 3);
774
775 auto set_vector_3d = [&](auto &t) -> MatrixDouble3by3 & {
776 mat.clear();
777 for (size_t r = 0; r != DIM1; ++r)
778 mat(0, r) = t(r);
779 return mat;
780 };
781
782 auto set_matrix_3d = [&](auto &t) -> MatrixDouble3by3 & {
783 mat.clear();
784 for (size_t r = 0; r != DIM1; ++r)
785 for (size_t c = 0; c != DIM2; ++c)
786 mat(r, c) = t(r, c);
787 return mat;
788 };
789
790 auto set_matrix_symm_3d = [&](auto &t) -> MatrixDouble3by3 & {
791 mat.clear();
792 for (size_t r = 0; r != DIM1; ++r)
793 for (size_t c = 0; c != DIM1; ++c)
794 mat(r, c) = t(r, c);
795 return mat;
796 };
797
798 auto set_scalar = [&](auto t) -> MatrixDouble3by3 & {
799 mat.clear();
800 mat(0, 0) = t;
801 return mat;
802 };
803
804 auto set_float_precision = [](const double x) {
805 if (std::abs(x) < std::numeric_limits<float>::min())
806 return static_cast<float>(0.);
807 else
808 return static_cast<float>(x);
809 };
810
811 auto set_float_max = [](const double x) {
812 if (x > std::numeric_limits<float>::max())
813 return std::numeric_limits<float>::max();
814 else
815 return static_cast<float>(x);
816 };
817
818 auto set_float_lowest = [](const double x) {
819 if (x < std::numeric_limits<float>::lowest())
820 return std::numeric_limits<float>::lowest();
821 else
822 return static_cast<float>(x);
823 };
824
825 auto set_tag = [&](auto th, auto gg, MatrixDouble3by3 &mat) {
826 for (auto &v : mat.data())
827 v = set_float_max(set_float_lowest(set_float_precision(v)));
828 return postProcMesh.tag_set_data(th, &mapGaussPts[gg], 1,
829 &*mat.data().begin());
830 };
831
832 for (auto &m : dataMapScalar) {
833 if (m.second) {
834 auto th = get_tag(m.first, 1);
835 auto t_scl = getFTensor0FromVec(*m.second);
836 auto nb_integration_pts = O::getGaussPts().size2();
837 for (int gg = 0; gg != nb_integration_pts; ++gg) {
838 CHKERR set_tag(th, gg, set_scalar(t_scl));
839 ++t_scl;
840 }
841 }
842 }
843
844 for (auto &m : dataMapVec) {
845 if (m.second) {
846 auto th = get_tag(m.first, 3);
847 auto t_vec = getFTensor1FromMat<DIM1>(*m.second);
848 auto nb_integration_pts = O::getGaussPts().size2();
849 for (int gg = 0; gg != nb_integration_pts; ++gg) {
850 CHKERR set_tag(th, gg, set_vector_3d(t_vec));
851 ++t_vec;
852 }
853 }
854 }
855
856 for (auto &m : dataMapMat) {
857 if (m.second) {
858 auto th = get_tag(m.first, 9);
859 auto t_mat = getFTensor2FromMat<DIM1, DIM2>(*m.second);
860 auto nb_integration_pts = O::getGaussPts().size2();
861 for (int gg = 0; gg != nb_integration_pts; ++gg) {
862 CHKERR set_tag(th, gg, set_matrix_3d(t_mat));
863 ++t_mat;
864 }
865 }
866 }
867
868 for (auto &m : dataMapSymmMat) {
869 if (m.second) {
870 auto th = get_tag(m.first, 9);
871 auto t_mat = getFTensor2SymmetricFromMat<DIM1>(*m.second);
872 auto nb_integration_pts = O::getGaussPts().size2();
873 for (int gg = 0; gg != nb_integration_pts; ++gg) {
874 CHKERR set_tag(th, gg, set_matrix_symm_3d(t_mat));
875 ++t_mat;
876 }
877 }
878 }
879
881}
882
883template <typename PostProcEle> struct PostProcBrokenMeshInMoabBaseBeginImpl;
884
885template <typename E>
887 : protected PostProcBrokenMeshInMoabBase<E> {
888
890 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
891 std::string opts_prefix = "")
892 : PostProcBrokenMeshInMoabBase<E>(m_field, core_mesh_ptr, opts_prefix) {}
893
896 ParallelComm *pcomm_post_proc_mesh =
897 ParallelComm::get_pcomm(&this->getPostProcMesh(), MYPCOMM_INDEX);
898 if (pcomm_post_proc_mesh != NULL)
899 delete pcomm_post_proc_mesh;
900 CHKERR this->getPostProcMesh().delete_mesh();
901 pcomm_post_proc_mesh =
902 new ParallelComm(&(this->getPostProcMesh()), PETSC_COMM_WORLD);
904 }
905
906 inline FEMethod *getFEMethod() { return this; }
907
910};
911
912template <typename PostProcEle> struct PostProcBrokenMeshInMoabBaseContImpl;
913
914template <typename E>
916 : public PostProcBrokenMeshInMoabBase<E> {
917
919 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
920 std::string opts_prefix = "")
921 : PostProcBrokenMeshInMoabBase<E>(m_field, core_mesh_ptr, opts_prefix) {}
922
925
926protected:
928};
929
930template <typename PostProcEle> struct PostProcBrokenMeshInMoabBaseEndImpl;
931
932template <typename E>
934 : protected PostProcBrokenMeshInMoabBase<E> {
935
937 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
938 std::string opts_prefix = "")
939 : PostProcBrokenMeshInMoabBase<E>(m_field, core_mesh_ptr, opts_prefix) {}
940
945 ParallelComm *pcomm_post_proc_mesh =
946 ParallelComm::get_pcomm(&this->getPostProcMesh(), MYPCOMM_INDEX);
947 if (pcomm_post_proc_mesh == nullptr)
948 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
949 "PComm not allocated");
950 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
952 }
953
954 inline FEMethod *getFEMethod() { return this; }
955
957};
958
959/**
960 * @brief Enable to run stack of post-processing elements. Use this to begin stack.
961 *
962 * See @ref scalar_check_approximation.cpp
963 *
964 */
967
968/**
969 * @brief Enable to run stack of post-processing elements.
970 *
971 * See @ref scalar_check_approximation.cpp
972 *
973 * @tparam E
974 */
975template <typename E>
978
979/**
980 * @brief Enable to run stack of post-processing elements. Use this to end stack.
981 *
982 */
985
986
987} // namespace MoFEM
988
989#endif //__POSTPROCBROKENMESHINMOABBASE_HPP
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#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 CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(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 ...
constexpr int order
Order displacement.
MoFEMErrorCode writeFile(const std::string file_name)
wrote results in (MOAB) format, use "file_name.h5m"
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
constexpr int DIM2
Definition level_set.cpp:22
constexpr int DIM1
Definition level_set.cpp:21
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto type_from_handle(const EntityHandle h)
get type from entity handle
boost::shared_ptr< PostProcGenerateRefMeshBase > PostProcGenerateRefMeshPtr
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
constexpr double t
plate stiffness
Definition plate.cpp:58
FTensor::Index< 'm', 3 > m
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Structure for user loop methods on finite elements.
@ OPSPACE
operator do Work is execute on space data
Post post-proc data at points from hash maps.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
std::vector< EntityHandle > & mapGaussPts
std::map< std::string, boost::shared_ptr< MatrixDouble > > DataMapMat
std::map< std::string, boost::shared_ptr< VectorDouble > > DataMapVec
OpPostProcMapInMoab(moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, DataMapVec data_map_scalar, DataMapMat data_map_vec, DataMapMat data_map_mat, DataMapMat data_symm_map_mat)
Construct a new OpPostProcMapInMoab object.
PostProcBrokenMeshInMoabBaseBeginImpl(MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, std::string opts_prefix="")
PostProcBrokenMeshInMoabBaseContImpl(MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, std::string opts_prefix="")
PostProcBrokenMeshInMoabBaseEndImpl(MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, std::string opts_prefix="")
auto & getPostProcElements()
Get postprocessing elements.
MoFEMErrorCode preProcess()
Generate vertices and elements.
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, std::string opts_prefix="")
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field, std::string opts_prefix="")
auto & getMapGaussPts()
Get vector of vectors associated to integration points.
auto & getPostProcMesh()
Get postprocessing mesh.
virtual int getMaxLevel() const
Determine refinement level based on fields approx ordre.
MoFEMErrorCode setTagsToTransfer(std::vector< Tag > tags_to_transfer)
Set tags to be transferred to post-processing mesh.
std::map< EntityType, PostProcGenerateRefMeshPtr > refElementsMap
virtual MoFEMErrorCode generateReferenceElementMesh()=0
virtual ~PostProcGenerateRefMeshBase()=default
std::vector< ublas::matrix< int > > levelRef
virtual MoFEMErrorCode getOptions(std::string prefix)
Element for postprocessing. Uses MoAB to generate post-processing mesh.