v0.15.5
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
187
188};
189
190template <typename E>
192 std::vector<Tag> tags_to_transfer) {
194 tagsToTransfer.swap(tags_to_transfer);
196}
197
198template <typename E> int PostProcBrokenMeshInMoabBase<E>::getMaxLevel() const {
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};
211
212template <typename E>
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};
288
289template <typename E>
291 MoFEM::Interface &m_field, std::string opts_prefix)
292 : E(m_field), optionsPrefix(opts_prefix) {}
293
294template <typename E>
296 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
297 std::string opts_prefix)
298 : PostProcBrokenMeshInMoabBase(m_field, opts_prefix) {
299 coreMeshPtr = core_mesh_ptr;
300}
301
302template <typename E>
306
307template <typename E> int PostProcBrokenMeshInMoabBase<E>::getRule(int order) {
308 return -1;
309};
310
311template <typename E>
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);
419 CHKERR transferTags();
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};
437
438template <typename E>
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}
552
553template <typename E>
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}
624
625template <typename E>
628
629 CHKERR getPostProcMesh().delete_mesh();
630
631
632 CHKERR preProcPostProc();
633
635}
636
637template <typename E>
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
646 CHKERR postProcPostProc();
647 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
648
650}
651
653 return mapGaussPts;
654}
655
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}
662
663template <typename E>
665 return postProcElements;
666}
667
668template <typename E>
670PostProcBrokenMeshInMoabBase<E>::writeFile(const std::string file_name) {
672 auto pcomm_post_proc_mesh = getPostProcMeshPcommPtr();
673 if (!pcomm_post_proc_mesh)
674 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
675 "ParallelComm not allocated");
676 CHKERR getPostProcMesh().write_file(file_name.c_str(), "MOAB",
677 "PARALLEL=WRITE_PART");
679};
680
681template <typename E> struct PostProcBrokenMeshInMoab;
682
683template <>
685 : public PostProcBrokenMeshInMoabBase<VolumeElementForcesAndSourcesCore> {
687 VolumeElementForcesAndSourcesCore>::PostProcBrokenMeshInMoabBase;
688};
689
690template <>
692 : public PostProcBrokenMeshInMoabBase<FaceElementForcesAndSourcesCore> {
694 FaceElementForcesAndSourcesCore>::PostProcBrokenMeshInMoabBase;
695};
696
697template <>
699 : public PostProcBrokenMeshInMoabBase<EdgeElementForcesAndSourcesCore> {
701 EdgeElementForcesAndSourcesCore>::PostProcBrokenMeshInMoabBase;
702};
703
704/**
705 * @brief Post post-proc data at points from hash maps
706 *
707 * @tparam DIM1 dimension of vector in data_map_vec and first column of
708 * data_map_may
709 * @tparam DIM2 dimension of second column in data_map_mat
710 */
711template <int DIM1, int DIM2,
713struct OpPostProcMapInMoab : public O {
714
715 using DataMapVec = std::map<std::string, boost::shared_ptr<VectorDouble>>;
716 using DataMapMat = std::map<std::string, boost::shared_ptr<MatrixDouble>>;
717
718 /**
719 * @brief Construct a new OpPostProcMapInMoab object
720 *
721 * @param post_proc_mesh postprocessing mesh
722 * @param map_gauss_pts map of gauss points to nodes of postprocessing mesh
723 * @param data_map_scalar hash map of scalar values (string is name of the
724 * tag)
725 * @param data_map_vec hash map of vector values
726 * @param data_map_mat hash map of second order tensor values
727 * @param data_symm_map_mat hash map of symmetric second order tensor values
728 */
729 OpPostProcMapInMoab(moab::Interface &post_proc_mesh,
730 std::vector<EntityHandle> &map_gauss_pts,
731 DataMapVec data_map_scalar, DataMapMat data_map_vec,
732 DataMapMat data_map_mat, DataMapMat data_symm_map_mat)
733 : O(NOSPACE, O::OPSPACE), postProcMesh(post_proc_mesh),
734 mapGaussPts(map_gauss_pts), dataMapScalar(data_map_scalar),
735 dataMapVec(data_map_vec), dataMapMat(data_map_mat),
736 dataMapSymmMat(data_symm_map_mat) {
737 // Operator is only executed for vertices
738 std::fill(&O::doEntities[MBEDGE], &O::doEntities[MBMAXTYPE], false);
739 }
740 MoFEMErrorCode doWork(int side, EntityType type,
742
743private:
744 moab::Interface &postProcMesh;
745 std::vector<EntityHandle> &mapGaussPts;
750};
751
752template <int DIM1, int DIM2, typename O>
757
758 std::array<double, 9> def;
759 std::fill(def.begin(), def.end(), 0);
760
761 auto get_tag = [&](const std::string name, size_t size) {
762 Tag th;
763 CHKERR postProcMesh.tag_get_handle(name.c_str(), size, MB_TYPE_DOUBLE, th,
764 MB_TAG_CREAT | MB_TAG_SPARSE,
765 def.data());
766 return th;
767 };
768
769 MatrixDouble3by3 mat(3, 3);
770
771 auto set_vector_3d = [&](auto &t) -> MatrixDouble3by3 & {
772 mat.clear();
773 for (size_t r = 0; r != DIM1; ++r)
774 mat(0, r) = t(r);
775 return mat;
776 };
777
778 auto set_matrix_3d = [&](auto &t) -> MatrixDouble3by3 & {
779 mat.clear();
780 for (size_t r = 0; r != DIM1; ++r)
781 for (size_t c = 0; c != DIM2; ++c)
782 mat(r, c) = t(r, c);
783 return mat;
784 };
785
786 auto set_matrix_symm_3d = [&](auto &t) -> MatrixDouble3by3 & {
787 mat.clear();
788 for (size_t r = 0; r != DIM1; ++r)
789 for (size_t c = 0; c != DIM1; ++c)
790 mat(r, c) = t(r, c);
791 return mat;
792 };
793
794 auto set_scalar = [&](auto t) -> MatrixDouble3by3 & {
795 mat.clear();
796 mat(0, 0) = t;
797 return mat;
798 };
799
800 auto set_float_precision = [](const double x) {
801 if (std::abs(x) < std::numeric_limits<float>::min())
802 return static_cast<float>(0.);
803 else
804 return static_cast<float>(x);
805 };
806
807 auto set_float_max = [](const double x) {
808 if (x > std::numeric_limits<float>::max())
809 return std::numeric_limits<float>::max();
810 else
811 return static_cast<float>(x);
812 };
813
814 auto set_float_lowest = [](const double x) {
815 if (x < std::numeric_limits<float>::lowest())
816 return std::numeric_limits<float>::lowest();
817 else
818 return static_cast<float>(x);
819 };
820
821 auto set_tag = [&](auto th, auto gg, MatrixDouble3by3 &mat) {
822 for (auto &v : mat.data())
823 v = set_float_max(set_float_lowest(set_float_precision(v)));
824 return postProcMesh.tag_set_data(th, &mapGaussPts[gg], 1,
825 &*mat.data().begin());
826 };
827
828 for (auto &m : dataMapScalar) {
829 if (m.second) {
830 auto th = get_tag(m.first, 1);
831 auto t_scl = getFTensor0FromVec(*m.second);
832 auto nb_integration_pts = O::getGaussPts().size2();
833 for (int gg = 0; gg != nb_integration_pts; ++gg) {
834 CHKERR set_tag(th, gg, set_scalar(t_scl));
835 ++t_scl;
836 }
837 }
838 }
839
840 for (auto &m : dataMapVec) {
841 if (m.second) {
842 auto th = get_tag(m.first, 3);
843 auto t_vec = getFTensor1FromMat<DIM1>(*m.second);
844 auto nb_integration_pts = O::getGaussPts().size2();
845 for (int gg = 0; gg != nb_integration_pts; ++gg) {
846 CHKERR set_tag(th, gg, set_vector_3d(t_vec));
847 ++t_vec;
848 }
849 }
850 }
851
852 for (auto &m : dataMapMat) {
853 if (m.second) {
854 auto th = get_tag(m.first, 9);
855 auto t_mat = getFTensor2FromMat<DIM1, DIM2>(*m.second);
856 auto nb_integration_pts = O::getGaussPts().size2();
857 for (int gg = 0; gg != nb_integration_pts; ++gg) {
858 CHKERR set_tag(th, gg, set_matrix_3d(t_mat));
859 ++t_mat;
860 }
861 }
862 }
863
864 for (auto &m : dataMapSymmMat) {
865 if (m.second) {
866 auto th = get_tag(m.first, 9);
867 auto t_mat = getFTensor2SymmetricFromMat<DIM1>(*m.second);
868 auto nb_integration_pts = O::getGaussPts().size2();
869 for (int gg = 0; gg != nb_integration_pts; ++gg) {
870 CHKERR set_tag(th, gg, set_matrix_symm_3d(t_mat));
871 ++t_mat;
872 }
873 }
874 }
875
877}
878
879template <typename PostProcEle> struct PostProcBrokenMeshInMoabBaseBeginImpl;
880
881template <typename E>
883 : protected PostProcBrokenMeshInMoabBase<E> {
884
886 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
887 std::string opts_prefix = "")
888 : PostProcBrokenMeshInMoabBase<E>(m_field, core_mesh_ptr, opts_prefix) {}
889
895
896 inline FEMethod *getFEMethod() { return this; }
897
900};
901
902template <typename PostProcEle> struct PostProcBrokenMeshInMoabBaseContImpl;
903
904template <typename E>
906 : public PostProcBrokenMeshInMoabBase<E> {
907
909 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
910 std::string opts_prefix = "")
911 : PostProcBrokenMeshInMoabBase<E>(m_field, core_mesh_ptr, opts_prefix) {}
912
915
916protected:
918};
919
920template <typename PostProcEle> struct PostProcBrokenMeshInMoabBaseEndImpl;
921
922template <typename E>
924 : protected PostProcBrokenMeshInMoabBase<E> {
925
927 MoFEM::Interface &m_field, boost::shared_ptr<moab::Core> core_mesh_ptr,
928 std::string opts_prefix = "")
929 : PostProcBrokenMeshInMoabBase<E>(m_field, core_mesh_ptr, opts_prefix) {}
930
935 auto pcomm_post_proc_mesh = this->getPostProcMeshPcommPtr();
936 if (!pcomm_post_proc_mesh)
937 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
938 "ParallelComm not allocated");
939 CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
941 }
942
943 inline FEMethod *getFEMethod() { return this; }
944
946};
947
948/**
949 * @brief Enable to run stack of post-processing elements. Use this to begin stack.
950 *
951 * See @ref scalar_check_approximation.cpp
952 *
953 */
956
957/**
958 * @brief Enable to run stack of post-processing elements.
959 *
960 * See @ref scalar_check_approximation.cpp
961 *
962 * @tparam E
963 */
964template <typename E>
967
968/**
969 * @brief Enable to run stack of post-processing elements. Use this to end stack.
970 *
971 */
974
975
976} // namespace MoFEM
977
978#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
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
static auto getOrCreate(moab::Core *core_mesh_ptr)
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.