Loading [MathJax]/extensions/tex2jax.js
OGS
NodeWiseMeshPartitioner.cpp
Go to the documentation of this file.
1
16
17#include <limits>
18#include <numeric>
19#include <range/v3/algorithm/transform.hpp>
20#include <range/v3/range/conversion.hpp>
21#include <unordered_map>
22
23#include "BaseLib/Error.h"
24#include "BaseLib/FileTools.h"
25#include "BaseLib/Logging.h"
26#include "BaseLib/RunTime.h"
28#include "MeshLib/IO/NodeData.h"
30#include "MeshLib/MeshEnums.h"
33
34namespace ApplicationUtils
35{
37 MeshLib::MeshItemType const item_type) const
38{
39 switch (item_type)
40 {
42 return nodes.size();
44 return regular_elements.size() + ghost_elements.size();
47 default:
48 OGS_FATAL("Unsupported MeshItemType {:s}.",
49 MeshLib::toString(item_type));
50 }
51}
52
54 std::ostream& os, std::vector<std::size_t> const& global_node_ids) const
55{
56 std::vector<MeshLib::IO::NodeData> nodes_buffer;
57 nodes_buffer.reserve(nodes.size());
58
59 for (const auto* node : nodes)
60 {
61 double const* coords = node->data();
62 nodes_buffer.emplace_back(global_node_ids[node->getID()], coords[0],
63 coords[1], coords[2]);
64 }
65 return os.write(reinterpret_cast<const char*>(nodes_buffer.data()),
66 sizeof(MeshLib::IO::NodeData) * nodes_buffer.size());
67}
68
74 std::vector<const MeshLib::Element*> const& elements)
75{
76 return 3 * elements.size() +
77 std::accumulate(begin(elements), end(elements), 0,
78 [](auto const nnodes, auto const* e)
79 { return nnodes + e->getNumberOfNodes(); });
80}
81
82std::ostream& Partition::writeConfig(std::ostream& os) const
83{
84 long const data[] = {
85 static_cast<long>(nodes.size()),
86 static_cast<long>(number_of_base_nodes),
87 static_cast<long>(regular_elements.size()),
88 static_cast<long>(ghost_elements.size()),
89 static_cast<long>(number_of_regular_base_nodes),
90 static_cast<long>(number_of_regular_nodes),
91 static_cast<long>(number_of_mesh_base_nodes),
92 static_cast<long>(number_of_mesh_all_nodes),
93 static_cast<long>(
95 static_cast<long>(
97 };
98
99 return os.write(reinterpret_cast<const char*>(data), sizeof(data));
100}
101
102std::size_t partitionLookup(std::size_t const& node_id,
103 std::vector<std::size_t> const& partition_ids,
104 std::span<std::size_t const> const node_id_mapping)
105{
106 return partition_ids[node_id_mapping[node_id]];
107}
108
109std::pair<std::vector<MeshLib::Node const*>, std::vector<MeshLib::Node const*>>
110splitIntoBaseAndHigherOrderNodes(std::vector<MeshLib::Node const*> const& nodes,
111 MeshLib::Mesh const& mesh)
112{
113 // Space for resulting vectors.
114 std::vector<MeshLib::Node const*> base_nodes;
115 // if linear mesh, then one reallocation, no realloc for higher order
116 // elements meshes.
117 base_nodes.reserve(nodes.size() / 2);
118 std::vector<MeshLib::Node const*> higher_order_nodes;
119 // if linear mesh, then wasted space, good estimate for quadratic
120 // order mesh, and realloc needed for higher order element meshes.
121 higher_order_nodes.reserve(nodes.size() / 2);
122
123 // Split the nodes into base nodes and extra nodes.
124 std::partition_copy(
125 begin(nodes), end(nodes), std::back_inserter(base_nodes),
126 std::back_inserter(higher_order_nodes),
127 [&](MeshLib::Node const* const n)
128 { return isBaseNode(*n, mesh.getElementsConnectedToNode(*n)); });
129
130 return {base_nodes, higher_order_nodes};
131}
132
136std::tuple<std::vector<MeshLib::Node*>, std::vector<MeshLib::Node*>>
138 std::size_t const part_id,
139 std::vector<MeshLib::Node*> const& nodes,
140 std::vector<MeshLib::Element const*> const& ghost_elements,
141 std::vector<std::size_t> const& partition_ids,
142 MeshLib::Mesh const& mesh,
143 std::span<std::size_t const> const node_id_mapping)
144{
145 std::vector<MeshLib::Node*> base_ghost_nodes;
146 std::vector<MeshLib::Node*> higher_order_ghost_nodes;
147
148 std::vector<bool> is_ghost_node(nodes.size(), false);
149 for (const auto* ghost_elem : ghost_elements)
150 {
151 for (unsigned i = 0; i < ghost_elem->getNumberOfNodes(); i++)
152 {
153 auto const& n = ghost_elem->getNode(i);
154 auto const node_id = n->getID();
155 if (is_ghost_node[node_id])
156 {
157 continue;
158 }
159
160 if (partitionLookup(node_id, partition_ids, node_id_mapping) !=
161 part_id)
162 {
163 if (isBaseNode(*n, mesh.getElementsConnectedToNode(*n)))
164 {
165 base_ghost_nodes.push_back(nodes[node_id]);
166 }
167 else
168 {
169 higher_order_ghost_nodes.push_back(nodes[node_id]);
170 }
171 is_ghost_node[node_id] = true;
172 }
173 }
174 }
175 return std::tuple<std::vector<MeshLib::Node*>, std::vector<MeshLib::Node*>>{
176 base_ghost_nodes, higher_order_ghost_nodes};
177}
178
181template <typename T>
183 Partition const& p,
184 std::size_t const offset,
186 MeshLib::PropertyVector<T>& partitioned_pv)
187{
188 auto const& nodes = p.nodes;
189 auto const nnodes = nodes.size();
190 auto const n_components = pv.getNumberOfGlobalComponents();
191 for (std::size_t i = 0; i < nnodes; ++i)
192 {
193 const auto global_id = nodes[i]->getID();
194 std::copy_n(&pv[n_components * global_id], n_components,
195 &partitioned_pv[offset + n_components * i]);
196 }
197 return n_components * nnodes;
198}
199
203template <typename T>
205 Partition const& p,
206 std::size_t const offset,
208 MeshLib::PropertyVector<T>& partitioned_pv)
209{
210 std::size_t const n_regular(p.regular_elements.size());
211 auto const n_components = pv.getNumberOfGlobalComponents();
212 for (std::size_t i = 0; i < n_regular; ++i)
213 {
214 const auto id = p.regular_elements[i]->getID();
215 std::copy_n(&pv[n_components * id], n_components,
216 &partitioned_pv[offset + n_components * i]);
217 }
218
219 std::size_t const n_ghost(p.ghost_elements.size());
220 for (std::size_t i = 0; i < n_ghost; ++i)
221 {
222 const auto id = p.ghost_elements[i]->getID();
223 std::copy_n(&pv[n_components * id], n_components,
224 &partitioned_pv[offset + n_components * (n_regular + i)]);
225 }
226 return n_components * (n_regular + n_ghost);
227}
228
233template <typename T>
235 MeshLib::Properties const& properties,
236 Partition const& p,
237 std::size_t const id_offset_partition,
238 std::vector<std::size_t> const& element_ip_data_offsets,
240 MeshLib::PropertyVector<T>& partitioned_pv)
241{
242 // Special field data such as OGS_VERSION, IntegrationPointMetaData,
243 // etc., which are not "real" integration points, are copied "as is"
244 // (i.e. fully) for every partition.
245 if (pv.getPropertyName().find("_ip") == std::string::npos)
246 {
247 std::copy_n(&pv[0], pv.size(), &partitioned_pv[id_offset_partition]);
248 return pv.size();
249 }
250
251 auto const n_components = pv.getNumberOfGlobalComponents();
252
253 std::size_t id_offset = 0;
254
255 auto const ip_meta_data = MeshLib::getIntegrationPointMetaDataSingleField(
257 auto copyFieldData =
258 [&](std::vector<const MeshLib::Element*> const& elements)
259 {
260 for (auto const element : elements)
261 {
262 int const number_of_element_field_data =
264 *element) *
265 n_components;
266 // The original element ID is not changed.
267 auto const element_id = element->getID();
268 int const begin_pos = element_ip_data_offsets[element_id];
269 int const end_pos = element_ip_data_offsets[element_id + 1];
270
271 std::copy(pv.begin() + begin_pos, pv.begin() + end_pos,
272 &partitioned_pv[id_offset + id_offset_partition]);
273 id_offset += number_of_element_field_data;
274 }
275 };
276
277 copyFieldData(p.regular_elements);
278 copyFieldData(p.ghost_elements);
279
280 return id_offset;
281}
282
284 std::vector<Partition>& partitions)
285{
286 auto const& opt_ip_meta_data_all =
288 for (auto const& [name, property] : properties)
289 {
290 auto const item_type = property->getMeshItemType();
291
293 {
294 continue;
295 }
296
297 // For special field data such as OGS_VERSION, IntegrationPointMetaData,
298 // etc., which are not "real" integration points:
299 if (property->getPropertyName().find("_ip") == std::string::npos)
300 {
301 continue;
302 }
303
304 auto const& ip_meta_data =
306 opt_ip_meta_data_all, property->getPropertyName());
307 auto countIntegrationPoints =
308 [&](std::vector<const MeshLib::Element*> const& elements)
309 {
310 std::size_t counter = 0;
311 for (auto const element : elements)
312 {
313 int const number_of_integration_points =
315 ip_meta_data, *element);
316 counter += number_of_integration_points;
317 }
318 return counter;
319 };
320
321 for (auto& p : partitions)
322 {
323 p.number_of_integration_points =
324 countIntegrationPoints(p.regular_elements) +
325 countIntegrationPoints(p.ghost_elements);
326 }
327 return;
328 }
329}
330
331template <typename T>
333 std::vector<MeshLib::Element*> const& global_mesh_elements,
334 MeshLib::Properties& partitioned_properties,
335 MeshLib::Properties const& properties,
336 std::vector<Partition> const& partitions,
337 MeshLib::PropertyVector<T> const* const pv,
338 std::map<MeshLib::MeshItemType, std::size_t> const& total_number_of_tuples)
339{
340 if (pv == nullptr)
341 {
342 return false;
343 }
344 auto const item_type = pv->getMeshItemType();
345
346 std::size_t partitioned_pv_size = total_number_of_tuples.at(item_type) *
348
349 std::vector<std::size_t> element_ip_data_offsets;
351 {
352 // Special field data such as OGS_VERSION, IntegrationPointMetaData,
353 // etc., which are not "real" integration points, are copied "as is"
354 // (i.e. fully) for every partition.
355 if (pv->getPropertyName().find("_ip") == std::string::npos)
356 {
357 partitioned_pv_size = pv->size() * partitions.size();
358 }
359
360 element_ip_data_offsets =
362 global_mesh_elements, *pv, properties);
363 }
364
365 auto partitioned_pv = partitioned_properties.createNewPropertyVector<T>(
366 pv->getPropertyName(), pv->getMeshItemType(),
368 if (partitioned_pv == nullptr)
369 {
370 OGS_FATAL(
371 "Could not create partitioned property vector {:s} for {} data "
372 "array.",
374 }
375 partitioned_pv->resize(partitioned_pv_size);
376
377 auto copy_property_vector_values =
378 [&](Partition const& p, std::size_t offset)
379 {
381 {
382 return copyFieldPropertyDataToPartitions(properties, p, offset,
383 element_ip_data_offsets,
384 *pv, *partitioned_pv);
385 }
386
387 if (item_type == MeshLib::MeshItemType::Node)
388 {
389 return copyNodePropertyVectorValues(p, offset, *pv,
390 *partitioned_pv);
391 }
392 if (item_type == MeshLib::MeshItemType::Cell)
393 {
394 return copyCellPropertyVectorValues(p, offset, *pv,
395 *partitioned_pv);
396 }
397
398 OGS_FATAL(
399 "Copying of property vector values for mesh item type {:s} is not "
400 "implemented.",
401 toString(item_type));
402 };
403
404 std::size_t position_offset(0);
405 for (auto p : partitions)
406 {
407 position_offset += copy_property_vector_values(p, position_offset);
408 }
409 return true;
410}
411
412void addVtkGhostTypeProperty(MeshLib::Properties& partitioned_properties,
413 std::vector<Partition> const& partitions,
414 std::size_t const total_number_of_cells)
415{
416 auto* vtk_ghost_type =
417 partitioned_properties.createNewPropertyVector<unsigned char>(
418 "vtkGhostType", MeshLib::MeshItemType::Cell, total_number_of_cells,
419 1);
420 if (vtk_ghost_type == nullptr)
421 {
422 OGS_FATAL("Could not create vtkGhostType cell data array.");
423 }
424
425 assert(vtk_ghost_type->size() == total_number_of_cells);
426 std::size_t offset = 0;
427 for (auto const& partition : partitions)
428 {
429 offset += partition.regular_elements.size();
430 for (std::size_t i = 0; i < partition.ghost_elements.size(); ++i)
431 {
432 if (partition.duplicate_ghost_cell[i])
433 {
434 (*vtk_ghost_type)[offset + i] |=
435 vtkDataSetAttributes::DUPLICATECELL;
436 }
437 }
438 offset += partition.ghost_elements.size();
439 }
440}
441
444 std::unique_ptr<MeshLib::Mesh> const& mesh,
445 std::vector<Partition>& partitions)
446{
447 using namespace MeshLib;
448
449 MeshLib::Properties const& properties = mesh->getProperties();
450
451 // Count the number of integration point data of all partitions:
452 setIntegrationPointNumberOfPartition(properties, partitions);
453
454 Properties partitioned_properties;
455 auto count_tuples = [&](MeshItemType const mesh_item_type)
456 {
457 return std::accumulate(
458 begin(partitions), end(partitions), 0,
459 [&](std::size_t const sum, Partition const& p)
460 { return sum + p.numberOfMeshItems(mesh_item_type); });
461 };
462
463 std::map<MeshItemType, std::size_t> const total_number_of_tuples = {
464 {MeshItemType::Cell, count_tuples(MeshItemType::Cell)},
465 {MeshItemType::Node, count_tuples(MeshItemType::Node)},
466 {MeshItemType::IntegrationPoint,
467 count_tuples(MeshItemType::IntegrationPoint)}};
468
469 DBUG(
470 "total number of tuples after partitioning defined for cells is {:d} "
471 "and for nodes {:d} and for integration points {:d}.",
472 total_number_of_tuples.at(MeshItemType::Cell),
473 total_number_of_tuples.at(MeshItemType::Node),
474 total_number_of_tuples.at(MeshItemType::IntegrationPoint));
475
476 // 1 create new PV
477 // 2 resize the PV with total_number_of_tuples
478 // 3 copy the values according to the partition info
480 properties,
481 [&](auto type, auto const property)
482 {
484 mesh->getElements(), partitioned_properties, properties,
485 partitions,
486 dynamic_cast<PropertyVector<decltype(type)> const*>(property),
487 total_number_of_tuples);
488 });
489
490 addVtkGhostTypeProperty(partitioned_properties,
491 partitions,
492 total_number_of_tuples.at(MeshItemType::Cell));
493
494 return partitioned_properties;
495}
496
498 std::vector<Partition>& partitions)
499{
500 std::vector<bool> cell_visited(mesh.getElements().size(), false);
501
502 for (auto& partition : partitions)
503 {
504 partition.duplicate_ghost_cell.resize(partition.ghost_elements.size(),
505 true);
506
507 for (std::size_t i = 0; i < partition.ghost_elements.size(); i++)
508 {
509 const auto& ghost_element = *partition.ghost_elements[i];
510 if (!cell_visited[ghost_element.getID()])
511 {
512 cell_visited[ghost_element.getID()] = true;
513 partition.duplicate_ghost_cell[i] = false;
514 }
515 }
516 }
517}
518
520 std::vector<MeshLib::Element*> const& global_mesh_elements,
521 MeshLib::Properties const& properties)
522{
523 auto const& opt_ip_meta_data_all =
525 for (auto const& [name, property] : properties)
526 {
527 auto const item_type = property->getMeshItemType();
528
530 {
531 continue;
532 }
533
534 // For special field data such as OGS_VERSION, IntegrationPointMetaData,
535 // etc., which are not "real" integration points:
536 if (property->getPropertyName().find("_ip") == std::string::npos)
537 {
538 continue;
539 }
540
541 std::size_t number_of_total_integration_points = 0;
542 auto const ip_meta_data =
544 opt_ip_meta_data_all, property->getPropertyName());
545 for (auto const element : global_mesh_elements)
546 {
547 int const number_of_integration_points =
549 *element);
550 number_of_total_integration_points += number_of_integration_points;
551 }
552
553 const auto pv =
554 dynamic_cast<MeshLib::PropertyVector<double> const*>(property);
555 std::size_t const component_number = pv->getNumberOfGlobalComponents();
556 if (pv->size() != number_of_total_integration_points * component_number)
557 {
558 OGS_FATAL(
559 "The property vector's size {:d} for integration point data "
560 "{:s} does not match its actual size {:d}. The field data in "
561 "the vtu file are wrong.",
562 pv->size(), name,
563 number_of_total_integration_points * component_number);
564 }
565 }
566}
567
568std::vector<std::vector<std::size_t>> computePartitionIDPerElement(
569 std::vector<std::size_t> const& node_partition_map,
570 std::vector<MeshLib::Element*> const& elements,
571 std::span<std::size_t const> const bulk_node_ids)
572{
573 auto node_partition_ids = ranges::views::transform(
574 [&](MeshLib::Element const* const element)
575 {
576 auto node_lookup = ranges::views::transform(
577 [&](std::size_t const i)
578 { return node_partition_map[bulk_node_ids[i]]; });
579
580 return element->nodes() | MeshLib::views::ids | node_lookup |
581 ranges::to<std::vector>;
582 });
583
584 return elements | node_partition_ids | ranges::to<std::vector>;
585}
586
588 std::vector<Partition>& partitions,
589 std::vector<std::size_t> const& nodes_partition_ids,
590 std::vector<MeshLib::Node*> const& nodes,
591 std::span<std::size_t const> const bulk_node_ids)
592{
593 for (auto const* const node : nodes)
594 {
595 partitions[nodes_partition_ids[bulk_node_ids[node->getID()]]]
596 .nodes.push_back(node);
597 }
598}
599
601 MeshLib::Mesh const& mesh)
602{
603 std::vector<MeshLib::Node const*> higher_order_nodes;
604 // after splitIntoBaseAndHigherOrderNodes() partition.nodes contains only
605 // base nodes
606 std::tie(partition.nodes, higher_order_nodes) =
608 partition.number_of_regular_base_nodes = partition.nodes.size();
609 std::copy(begin(higher_order_nodes), end(higher_order_nodes),
610 std::back_inserter(partition.nodes));
611 partition.number_of_regular_nodes = partition.nodes.size();
612}
613
615 std::vector<Partition>& partitions, MeshLib::Mesh const& mesh)
616{
617 for (auto& partition : partitions)
618 {
620 }
621}
622
623void setNumberOfNodesInPartitions(std::vector<Partition>& partitions,
624 MeshLib::Mesh const& mesh)
625{
626 auto const number_of_mesh_base_nodes = mesh.computeNumberOfBaseNodes();
627 auto const number_of_mesh_all_nodes = mesh.getNumberOfNodes();
628 for (auto& partition : partitions)
629 {
630 partition.number_of_regular_nodes = partition.nodes.size();
631 partition.number_of_mesh_base_nodes = number_of_mesh_base_nodes;
632 partition.number_of_mesh_all_nodes = number_of_mesh_all_nodes;
633 }
634}
635
637 std::vector<Partition>& partitions,
638 MeshLib::Mesh const& mesh,
639 std::vector<std::vector<std::size_t>> const& partition_ids_per_element)
640{
641 for (auto const& element : mesh.getElements())
642 {
643 auto const element_id = element->getID();
644 auto node_partition_ids = partition_ids_per_element[element_id];
645 // make partition ids unique
646 std::sort(node_partition_ids.begin(), node_partition_ids.end());
647 auto last =
648 std::unique(node_partition_ids.begin(), node_partition_ids.end());
649 node_partition_ids.erase(last, node_partition_ids.end());
650
651 // all element nodes belong to the same partition => regular element
652 if (node_partition_ids.size() == 1)
653 {
654 partitions[node_partition_ids[0]].regular_elements.push_back(
655 element);
656 }
657 else
658 {
659 for (auto const partition_id : node_partition_ids)
660 {
661 partitions[partition_id].ghost_elements.push_back(element);
662 }
663 }
664 }
665}
666
667// determine and append ghost nodes to partition.nodes in the following order
668// [base nodes, higher order nodes, base ghost nodes, higher order ghost
669// nodes]
671 std::vector<Partition>& partitions, MeshLib::Mesh const& mesh,
672 std::vector<std::size_t> const& nodes_partition_ids,
673 std::span<std::size_t const> const node_id_mapping)
674{
675 for (std::size_t part_id = 0; part_id < partitions.size(); part_id++)
676 {
677 auto& partition = partitions[part_id];
678 std::vector<MeshLib::Node*> base_ghost_nodes;
679 std::vector<MeshLib::Node*> higher_order_ghost_nodes;
680 std::tie(base_ghost_nodes, higher_order_ghost_nodes) =
682 part_id, mesh.getNodes(), partition.ghost_elements,
683 nodes_partition_ids, mesh, node_id_mapping);
684
685 std::copy(begin(base_ghost_nodes), end(base_ghost_nodes),
686 std::back_inserter(partition.nodes));
687
688 partition.number_of_base_nodes =
689 partition.number_of_regular_base_nodes + base_ghost_nodes.size();
690
691 std::copy(begin(higher_order_ghost_nodes),
692 end(higher_order_ghost_nodes),
693 std::back_inserter(partition.nodes));
694 }
695}
696
697void partitionMesh(std::vector<Partition>& partitions,
698 MeshLib::Mesh const& mesh,
699 std::vector<std::size_t> const& nodes_partition_ids,
700 std::span<std::size_t const> const bulk_node_ids)
701{
702 BaseLib::RunTime run_timer;
703 run_timer.start();
704 auto const partition_ids_per_element = computePartitionIDPerElement(
705 nodes_partition_ids, mesh.getElements(), bulk_node_ids);
706 INFO("partitionMesh(): Partition IDs per element computed in {:g} s",
707 run_timer.elapsed());
708
709 run_timer.start();
710 distributeNodesToPartitions(partitions, nodes_partition_ids,
711 mesh.getNodes(), bulk_node_ids);
712 INFO("partitionMesh(): distribute nodes to partitions took {:g} s",
713 run_timer.elapsed());
714
715 run_timer.start();
717 INFO(
718 "partitionMesh(): sorting [base nodes | higher order nodes] took {:g} "
719 "s",
720 run_timer.elapsed());
721
722 run_timer.start();
723 setNumberOfNodesInPartitions(partitions, mesh);
724 INFO(
725 "partitionMesh(): setting number of nodes and of all mesh base nodes "
726 "took {:g} s",
727 run_timer.elapsed());
728
729 run_timer.start();
730 distributeElementsIntoPartitions(partitions, mesh,
731 partition_ids_per_element);
732 INFO("partitionMesh(): distribute elements into partitions took {:g} s",
733 run_timer.elapsed());
734
735 run_timer.start();
737 partitions, mesh, nodes_partition_ids, bulk_node_ids);
738 INFO("partitionMesh(): determine / append ghost nodes took {:g} s",
739 run_timer.elapsed());
740
741 run_timer.start();
742 markDuplicateGhostCells(mesh, partitions);
743 INFO("partitionMesh(): markDuplicateGhostCells took {:g} s",
744 run_timer.elapsed());
745}
746
748{
749 std::vector<std::size_t> bulk_node_ids(_mesh->getNumberOfNodes());
750 std::iota(bulk_node_ids.begin(), bulk_node_ids.end(), 0);
751
753
755
756 // In case the field data in the vtu file are manually added, e.g. by using
757 // some tools, the size of the field property vector has to be checked.
758 checkFieldPropertyVectorSize(_mesh->getElements(), _mesh->getProperties());
759
761
763}
764
766 std::vector<Partition> const& partitions,
767 MeshLib::Properties& partitioned_properties)
768{
769 auto const bulk_node_ids_string =
771 if (partitioned_properties.hasPropertyVector(bulk_node_ids_string))
772 {
774 partitioned_properties.getPropertyVector<std::size_t>(
775 bulk_node_ids_string, MeshLib::MeshItemType::Node, 1),
776 partitions);
777 }
778 auto const bulk_element_ids_string =
780 if (partitioned_properties.hasPropertyVector<std::size_t>(
781 static_cast<std::string>(bulk_element_ids_string),
783 {
785 partitioned_properties.getPropertyVector<std::size_t>(
786 bulk_element_ids_string, MeshLib::MeshItemType::Cell, 1),
787 partitions);
788 }
789}
790
792 MeshLib::PropertyVector<std::size_t>* const bulk_node_ids_pv,
793 std::vector<Partition> const& local_partitions) const
794{
795 if (bulk_node_ids_pv == nullptr)
796 {
797 return;
798 }
799
800 auto& bulk_node_ids = *bulk_node_ids_pv;
801
802 std::size_t offset = 0; // offset in property vector for current partition
803
804 assert(_partitions.size() == local_partitions.size());
805 int const n_partitions = static_cast<int>(_partitions.size());
806 for (int partition_id = 0; partition_id < n_partitions; ++partition_id)
807 {
808 auto const& bulk_partition = _partitions[partition_id];
809 auto const& local_partition = local_partitions[partition_id];
810
811 // Create global-to-local node id mapping for the bulk partition.
812 auto const& bulk_nodes = bulk_partition.nodes;
813 auto const n_bulk_nodes = bulk_nodes.size();
814 std::map<std::size_t, std::size_t> global_to_local;
815 for (std::size_t local_node_id = 0; local_node_id < n_bulk_nodes;
816 ++local_node_id)
817 {
818 global_to_local[bulk_nodes[local_node_id]->getID()] = local_node_id;
819 }
820
821 auto const& local_nodes = local_partition.nodes;
822 auto const n_local_nodes = local_nodes.size();
823 for (std::size_t local_node_id = 0; local_node_id < n_local_nodes;
824 ++local_node_id)
825 {
826 bulk_node_ids[offset + local_node_id] =
827 global_to_local[bulk_node_ids[offset + local_node_id]];
828 }
829 offset += n_local_nodes;
830 }
831}
832
834 MeshLib::PropertyVector<std::size_t>* const bulk_element_ids_pv,
835 std::vector<Partition> const& local_partitions) const
836{
837 if (bulk_element_ids_pv == nullptr)
838 {
839 return;
840 }
841
842 auto& bulk_element_ids = *bulk_element_ids_pv;
843
844 std::size_t offset = 0; // offset in property vector for current partition
845
846 assert(_partitions.size() == local_partitions.size());
847 int const n_partitions = static_cast<int>(_partitions.size());
848 for (int partition_id = 0; partition_id < n_partitions; ++partition_id)
849 {
850 auto const& bulk_partition = _partitions[partition_id];
851 auto const& local_partition = local_partitions[partition_id];
852
853 // Create global-to-local element id mapping for the bulk partition.
854 std::map<std::size_t, std::size_t> global_to_local;
855 auto map_elements =
856 [&global_to_local](
857 std::vector<MeshLib::Element const*> const& elements,
858 std::size_t const offset)
859 {
860 auto const n_elements = elements.size();
861 for (std::size_t e = 0; e < n_elements; ++e)
862 {
863 global_to_local[elements[e]->getID()] = offset + e;
864 }
865 };
866
867 map_elements(bulk_partition.regular_elements, 0);
868 map_elements(bulk_partition.ghost_elements,
869 bulk_partition.regular_elements.size());
870
871 // Renumber the local bulk_element_ids map.
872 auto renumber_elements =
873 [&bulk_element_ids, &global_to_local](
874 std::vector<MeshLib::Element const*> const& elements,
875 std::size_t const offset)
876 {
877 auto const n_elements = elements.size();
878 for (std::size_t e = 0; e < n_elements; ++e)
879 {
880 bulk_element_ids[offset + e] =
881 global_to_local[bulk_element_ids[offset + e]];
882 }
883 return n_elements;
884 };
885
886 offset += renumber_elements(local_partition.regular_elements, offset);
887 offset += renumber_elements(local_partition.ghost_elements, offset);
888 }
889}
890
892 MeshLib::Mesh const& mesh) const
893{
894 auto const bulk_node_ids_string =
896 auto const& bulk_node_ids =
897 mesh.getProperties().getPropertyVector<std::size_t>(
898 bulk_node_ids_string, MeshLib::MeshItemType::Node, 1);
899
900 std::vector<Partition> partitions(_partitions.size());
901
902 partitionMesh(partitions, mesh, _nodes_partition_ids, *bulk_node_ids);
903
904 return partitions;
905}
906
908{
909 std::size_t node_global_id_offset = 0;
910 // Renumber the global indices.
911 for (auto& partition : _partitions)
912 {
913 for (std::size_t i = 0; i < partition.number_of_regular_nodes; i++)
914 {
915 _nodes_global_ids[partition.nodes[i]->getID()] =
916 node_global_id_offset++;
917 }
918 }
919}
920
921template <typename T>
922void writePropertyVectorValues(std::ostream& os,
924{
925 os.write(reinterpret_cast<const char*>(pv.data()), pv.size() * sizeof(T));
926}
927
928template <typename T>
930 MeshLib::MeshItemType const mesh_item_type,
931 std::ostream& out_val, std::ostream& out_meta)
932{
933 if (pv == nullptr)
934 {
935 return false;
936 }
937 // skip property of different mesh item type. Return true, because this
938 // operation was successful.
939 if (pv->getMeshItemType() != mesh_item_type)
940 {
941 return true;
942 }
943
945 pvmd.property_name = pv->getPropertyName();
949 writePropertyVectorValues(out_val, *pv);
951 return true;
952}
953
954void writeProperties(const std::string& file_name_base,
955 MeshLib::Properties const& partitioned_properties,
956 std::vector<Partition> const& partitions,
957 MeshLib::MeshItemType const mesh_item_type)
958{
959 auto const number_of_properties =
960 partitioned_properties.size(mesh_item_type);
961 if (number_of_properties == 0)
962 {
963 return;
964 }
965
966 auto const file_name_infix = toString(mesh_item_type);
967
968 auto const file_name_cfg = file_name_base + "_partitioned_" +
969 file_name_infix + "_properties_cfg" +
970 std::to_string(partitions.size()) + ".bin";
971 std::ofstream out(file_name_cfg, std::ios::binary);
972 if (!out)
973 {
974 OGS_FATAL("Could not open file '{:s}' for output.", file_name_cfg);
975 }
976
977 auto const file_name_val = file_name_base + "_partitioned_" +
978 file_name_infix + "_properties_val" +
979 std::to_string(partitions.size()) + ".bin";
980 std::ofstream out_val(file_name_val, std::ios::binary);
981 if (!out_val)
982 {
983 OGS_FATAL("Could not open file '{:s}' for output.", file_name_val);
984 }
985
986 BaseLib::writeValueBinary(out, number_of_properties);
987
989 partitioned_properties,
990 [&](auto type, auto const& property)
991 {
993 dynamic_cast<MeshLib::PropertyVector<decltype(type)> const*>(
994 property),
995 mesh_item_type, out_val, out);
996 });
997
998 unsigned long offset = 0;
999 for (const auto& partition : partitions)
1000 {
1002 offset, static_cast<unsigned long>(
1003 partition.numberOfMeshItems(mesh_item_type))};
1004 DBUG(
1005 "Write meta data for node-based PropertyVector: global offset "
1006 "{:d}, number of tuples {:d}",
1007 pvpmd.offset, pvpmd.number_of_tuples);
1009 offset += pvpmd.number_of_tuples;
1010 }
1011}
1012
1014{
1018
1019 std::ostream& writeConfig(std::ostream& os) const;
1020};
1021
1022std::ostream& ConfigOffsets::writeConfig(std::ostream& os) const
1023{
1024 os.write(reinterpret_cast<const char*>(this), sizeof(ConfigOffsets));
1025
1026 static long reserved = 0; // Value reserved in the binary format, not used
1027 // in the partitioning process.
1028 return os.write(reinterpret_cast<const char*>(&reserved), sizeof(long));
1029}
1030
1037
1039{
1040 return {static_cast<long>(partition.nodes.size()),
1041 static_cast<long>(partition.regular_elements.size() +
1043 partition.regular_elements)),
1044 static_cast<long>(partition.ghost_elements.size() +
1046 partition.ghost_elements))};
1047}
1048
1050 PartitionOffsets const& offsets)
1051{
1052 return {
1053 static_cast<long>(oldConfig.node_rank_offset +
1054 offsets.node * sizeof(MeshLib::IO::NodeData)),
1055 // Offset the ending entry of the element integer variables of
1056 // the non-ghost elements of this partition in the vector of elem_info.
1057 static_cast<long>(oldConfig.element_rank_offset +
1058 offsets.regular_elements * sizeof(long)),
1059
1060 // Offset the ending entry of the element integer variables of
1061 // the ghost elements of this partition in the vector of elem_info.
1062 static_cast<long>(oldConfig.ghost_element_rank_offset +
1063 offsets.ghost_elements * sizeof(long))};
1064}
1065
1071std::tuple<std::vector<long>, std::vector<long>> writeConfigData(
1072 const std::string& file_name_base, std::vector<Partition> const& partitions)
1073{
1074 auto const file_name_cfg = file_name_base + "_partitioned_msh_cfg" +
1075 std::to_string(partitions.size()) + ".bin";
1076 std::ofstream of_bin_cfg(file_name_cfg, std::ios::binary);
1077 if (!of_bin_cfg)
1078 {
1079 OGS_FATAL("Could not open file '{:s}' for output.", file_name_cfg);
1080 }
1081
1082 std::vector<long> partitions_element_offsets;
1083 partitions_element_offsets.reserve(partitions.size());
1084 std::vector<long> partitions_ghost_element_offsets;
1085 partitions_ghost_element_offsets.reserve(partitions.size());
1086
1087 ConfigOffsets config_offsets = {0, 0, 0}; // 0 for first partition.
1088 for (const auto& partition : partitions)
1089 {
1090 partition.writeConfig(of_bin_cfg);
1091
1092 config_offsets.writeConfig(of_bin_cfg);
1093 auto const& partition_offsets = computePartitionOffsets(partition);
1094 config_offsets =
1095 incrementConfigOffsets(config_offsets, partition_offsets);
1096
1097 partitions_element_offsets.push_back(
1098 partition_offsets.regular_elements);
1099 partitions_ghost_element_offsets.push_back(
1100 partition_offsets.ghost_elements);
1101 }
1102
1103 return std::make_tuple(partitions_element_offsets,
1104 partitions_ghost_element_offsets);
1105}
1106
1115 const MeshLib::Element& elem,
1116 const std::unordered_map<std::size_t, long>& local_node_ids,
1117 std::vector<long>& elem_info,
1118 long& counter)
1119{
1120 constexpr unsigned mat_id =
1121 0; // TODO: Material ID to be set from the mesh data
1122 const long nn = elem.getNumberOfNodes();
1123 elem_info[counter++] = mat_id;
1124 elem_info[counter++] = static_cast<long>(elem.getCellType());
1125 elem_info[counter++] = nn;
1126
1127 for (long i = 0; i < nn; i++)
1128 {
1129 auto const& n = *elem.getNode(i);
1130 elem_info[counter++] = local_node_ids.at(n.getID());
1131 }
1132}
1133
1135std::unordered_map<std::size_t, long> enumerateLocalNodeIds(
1136 std::vector<MeshLib::Node const*> const& nodes)
1137{
1138 std::unordered_map<std::size_t, long> local_ids;
1139 local_ids.reserve(nodes.size());
1140
1141 long local_node_id = 0;
1142 for (const auto* node : nodes)
1143 {
1144 local_ids[node->getID()] = local_node_id++;
1145 }
1146 return local_ids;
1147}
1148
1155void writeElements(std::string const& file_name_base,
1156 std::vector<Partition> const& partitions,
1157 std::vector<long> const& regular_element_offsets,
1158 std::vector<long> const& ghost_element_offsets)
1159{
1160 const std::string npartitions_str = std::to_string(partitions.size());
1161
1162 auto const file_name_ele =
1163 file_name_base + "_partitioned_msh_ele" + npartitions_str + ".bin";
1164 std::ofstream element_info_os(file_name_ele, std::ios::binary);
1165 if (!element_info_os)
1166 {
1167 OGS_FATAL("Could not open file '{:s}' for output.", file_name_ele);
1168 }
1169
1170 auto const file_name_ele_g =
1171 file_name_base + "_partitioned_msh_ele_g" + npartitions_str + ".bin";
1172 std::ofstream ghost_element_info_os(file_name_ele_g, std::ios::binary);
1173 if (!ghost_element_info_os)
1174 {
1175 OGS_FATAL("Could not open file '{:s}' for output.", file_name_ele_g);
1176 }
1177
1178 for (std::size_t i = 0; i < partitions.size(); i++)
1179 {
1180 const auto& partition = partitions[i];
1181 auto const local_node_ids = enumerateLocalNodeIds(partition.nodes);
1182
1183 // Vector containing the offsets of the regular elements of this
1184 // partition
1185 std::vector<long> ele_info(regular_element_offsets[i]);
1186
1187 auto writeElementData =
1188 [&local_node_ids](
1189 std::vector<MeshLib::Element const*> const& elements,
1190 long const element_offsets,
1191 std::ofstream& output_stream)
1192 {
1193 long counter = elements.size();
1194 std::vector<long> ele_info(element_offsets);
1195
1196 for (std::size_t j = 0; j < elements.size(); j++)
1197 {
1198 const auto* elem = elements[j];
1199 ele_info[j] = counter;
1200 getElementIntegerVariables(*elem, local_node_ids, ele_info,
1201 counter);
1202 }
1203 // Write vector data of regular elements
1204 output_stream.write(reinterpret_cast<const char*>(ele_info.data()),
1205 ele_info.size() * sizeof(long));
1206 };
1207
1208 // regular elements.
1209 writeElementData(partition.regular_elements, regular_element_offsets[i],
1210 element_info_os);
1211 // Ghost elements
1212 writeElementData(partition.ghost_elements, ghost_element_offsets[i],
1213 ghost_element_info_os);
1214 }
1215}
1216
1221void writeNodes(const std::string& file_name_base,
1222 std::vector<Partition> const& partitions,
1223 std::vector<std::size_t> const& global_node_ids)
1224{
1225 auto const file_name = file_name_base + "_partitioned_msh_nod" +
1226 std::to_string(partitions.size()) + ".bin";
1227 std::ofstream os(file_name, std::ios::binary);
1228 if (!os)
1229 {
1230 OGS_FATAL("Could not open file '{:s}' for output.", file_name);
1231 }
1232
1233 for (const auto& partition : partitions)
1234 {
1235 partition.writeNodes(os, global_node_ids);
1236 }
1237}
1238
1239void NodeWiseMeshPartitioner::write(const std::string& file_name_base)
1240{
1247
1248 const auto elements_offsets = writeConfigData(file_name_base, _partitions);
1249
1250 const std::vector<IntegerType>& regular_element_offsets =
1251 std::get<0>(elements_offsets);
1252 const std::vector<IntegerType>& ghost_element_offsets =
1253 std::get<1>(elements_offsets);
1254 writeElements(file_name_base, _partitions, regular_element_offsets,
1255 ghost_element_offsets);
1256
1257 writeNodes(file_name_base, _partitions, _nodes_global_ids);
1258}
1259
1261 std::string const& output_filename_base,
1262 std::vector<Partition> const& partitions,
1263 MeshLib::Properties const& partitioned_properties) const
1264{
1265 writeNodes(output_filename_base, partitions, _nodes_global_ids);
1266
1267 const auto elem_integers =
1268 writeConfigData(output_filename_base, partitions);
1269
1270 const std::vector<IntegerType>& num_elem_integers =
1271 std::get<0>(elem_integers);
1272 const std::vector<IntegerType>& num_g_elem_integers =
1273 std::get<1>(elem_integers);
1274 writeElements(output_filename_base, partitions, num_elem_integers,
1275 num_g_elem_integers);
1276
1277 writeProperties(output_filename_base, partitioned_properties, partitions,
1279 writeProperties(output_filename_base, partitioned_properties, partitions,
1281}
1282} // namespace ApplicationUtils
#define OGS_FATAL(...)
Definition Error.h:26
Filename manipulation routines.
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of mesh-related Enumerations.
Declare a class to perform node wise mesh partitioning.
void applyToPropertyVectors(Properties const &properties, Function f)
Definition of the RunTime class.
Implementation of the VtuInterface class.
void renumberBulkElementIdsProperty(MeshLib::PropertyVector< std::size_t > *const bulk_element_ids_pv, std::vector< Partition > const &local_partitions) const
std::vector< Partition > partitionOtherMesh(MeshLib::Mesh const &mesh) const
void write(const std::string &file_name_base)
void renumberBulkIdsProperty(std::vector< Partition > const &partitions, MeshLib::Properties &partitioned_properties)
void renumberBulkNodeIdsProperty(MeshLib::PropertyVector< std::size_t > *const bulk_node_ids, std::vector< Partition > const &local_partitions) const
std::unique_ptr< MeshLib::Mesh > _mesh
Pointer to a mesh object.
MeshLib::Properties _partitioned_properties
Properties where values at ghost nodes and extra nodes are inserted.
std::vector< std::size_t > _nodes_global_ids
Global IDs of all nodes after partitioning.
std::vector< Partition > _partitions
Data for all partitions.
std::vector< std::size_t > _nodes_partition_ids
Partition IDs of each nodes.
void writeOtherMesh(std::string const &output_filename_base, std::vector< Partition > const &partitions, MeshLib::Properties const &partitioned_properties) const
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
virtual CellType getCellType() const =0
virtual unsigned getNumberOfNodes() const =0
virtual const Node * getNode(unsigned idx) const =0
constexpr std::span< Node *const > nodes() const
Span of element's nodes, their pointers actually.
Definition Element.h:74
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:108
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:111
std::size_t computeNumberOfBaseNodes() const
Get the number of base nodes.
Definition Mesh.cpp:239
Properties & getProperties()
Definition Mesh.h:136
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:102
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition Mesh.cpp:257
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
Definition Properties.h:33
bool hasPropertyVector(std::string_view name) const
std::map< std::string, PropertyVectorBase * >::size_type size() const
PropertyVector< T > * createNewPropertyVector(std::string_view name, MeshItemType mesh_item_type, std::size_t n_components=1)
Definition Properties.h:17
PropertyVector< T > const * getPropertyVector(std::string_view name) const
MeshItemType getMeshItemType() const
int getNumberOfGlobalComponents() const
std::string const & getPropertyName() const
constexpr std::size_t size() const
constexpr std::size_t getNumberOfTuples() const
constexpr const PROP_VAL_TYPE * data() const
constexpr PROP_VAL_TYPE * begin()
std::size_t copyCellPropertyVectorValues(Partition const &p, std::size_t const offset, MeshLib::PropertyVector< T > const &pv, MeshLib::PropertyVector< T > &partitioned_pv)
std::unordered_map< std::size_t, long > enumerateLocalNodeIds(std::vector< MeshLib::Node const * > const &nodes)
Generates a mapping of given node ids to a new local (renumbered) node ids.
void addVtkGhostTypeProperty(MeshLib::Properties &partitioned_properties, std::vector< Partition > const &partitions, std::size_t const total_number_of_cells)
bool writePropertyVector(MeshLib::PropertyVector< T > const *const pv, MeshLib::MeshItemType const mesh_item_type, std::ostream &out_val, std::ostream &out_meta)
void checkFieldPropertyVectorSize(std::vector< MeshLib::Element * > const &global_mesh_elements, MeshLib::Properties const &properties)
void writeNodes(const std::string &file_name_base, std::vector< Partition > const &partitions, std::vector< std::size_t > const &global_node_ids)
ConfigOffsets incrementConfigOffsets(ConfigOffsets const &oldConfig, PartitionOffsets const &offsets)
std::vector< std::vector< std::size_t > > computePartitionIDPerElement(std::vector< std::size_t > const &node_partition_map, std::vector< MeshLib::Element * > const &elements, std::span< std::size_t const > const bulk_node_ids)
std::tuple< std::vector< MeshLib::Node * >, std::vector< MeshLib::Node * > > findGhostNodesInPartition(std::size_t const part_id, std::vector< MeshLib::Node * > const &nodes, std::vector< MeshLib::Element const * > const &ghost_elements, std::vector< std::size_t > const &partition_ids, MeshLib::Mesh const &mesh, std::span< std::size_t const > const node_id_mapping)
void setIntegrationPointNumberOfPartition(MeshLib::Properties const &properties, std::vector< Partition > &partitions)
void reorderNodesIntoBaseAndHigherOrderNodes(Partition &partition, MeshLib::Mesh const &mesh)
void reorderNodesIntoBaseAndHigherOrderNodesPerPartition(std::vector< Partition > &partitions, MeshLib::Mesh const &mesh)
MeshLib::Properties partitionProperties(std::unique_ptr< MeshLib::Mesh > const &mesh, std::vector< Partition > &partitions)
Partition existing properties and add vtkGhostType cell data array property.
void setNumberOfNodesInPartitions(std::vector< Partition > &partitions, MeshLib::Mesh const &mesh)
std::pair< std::vector< MeshLib::Node const * >, std::vector< MeshLib::Node const * > > splitIntoBaseAndHigherOrderNodes(std::vector< MeshLib::Node const * > const &nodes, MeshLib::Mesh const &mesh)
void distributeNodesToPartitions(std::vector< Partition > &partitions, std::vector< std::size_t > const &nodes_partition_ids, std::vector< MeshLib::Node * > const &nodes, std::span< std::size_t const > const bulk_node_ids)
PartitionOffsets computePartitionOffsets(Partition const &partition)
std::size_t partitionLookup(std::size_t const &node_id, std::vector< std::size_t > const &partition_ids, std::span< std::size_t const > const node_id_mapping)
NodeWiseMeshPartitioner::IntegerType getNumberOfIntegerVariablesOfElements(std::vector< const MeshLib::Element * > const &elements)
void partitionMesh(std::vector< Partition > &partitions, MeshLib::Mesh const &mesh, std::vector< std::size_t > const &nodes_partition_ids, std::span< std::size_t const > const bulk_node_ids)
std::tuple< std::vector< long >, std::vector< long > > writeConfigData(const std::string &file_name_base, std::vector< Partition > const &partitions)
std::size_t copyFieldPropertyDataToPartitions(MeshLib::Properties const &properties, Partition const &p, std::size_t const id_offset_partition, std::vector< std::size_t > const &element_ip_data_offsets, MeshLib::PropertyVector< T > const &pv, MeshLib::PropertyVector< T > &partitioned_pv)
void writePropertyVectorValues(std::ostream &os, MeshLib::PropertyVector< T > const &pv)
void getElementIntegerVariables(const MeshLib::Element &elem, const std::unordered_map< std::size_t, long > &local_node_ids, std::vector< long > &elem_info, long &counter)
void markDuplicateGhostCells(MeshLib::Mesh const &mesh, std::vector< Partition > &partitions)
void distributeElementsIntoPartitions(std::vector< Partition > &partitions, MeshLib::Mesh const &mesh, std::vector< std::vector< std::size_t > > const &partition_ids_per_element)
std::size_t copyNodePropertyVectorValues(Partition const &p, std::size_t const offset, MeshLib::PropertyVector< T > const &pv, MeshLib::PropertyVector< T > &partitioned_pv)
void writeProperties(const std::string &file_name_base, MeshLib::Properties const &partitioned_properties, std::vector< Partition > const &partitions, MeshLib::MeshItemType const mesh_item_type)
void writeElements(std::string const &file_name_base, std::vector< Partition > const &partitions, std::vector< long > const &regular_element_offsets, std::vector< long > const &ghost_element_offsets)
bool copyPropertyVector(std::vector< MeshLib::Element * > const &global_mesh_elements, MeshLib::Properties &partitioned_properties, MeshLib::Properties const &properties, std::vector< Partition > const &partitions, MeshLib::PropertyVector< T > const *const pv, std::map< MeshLib::MeshItemType, std::size_t > const &total_number_of_tuples)
void determineAndAppendGhostNodesToPartitions(std::vector< Partition > &partitions, MeshLib::Mesh const &mesh, std::vector< std::size_t > const &nodes_partition_ids, std::span< std::size_t const > const node_id_mapping)
void writeValueBinary(std::ostream &out, T const &val)
write value as binary into the given output stream
void writePropertyVectorPartitionMetaData(std::ostream &os, PropertyVectorPartitionMetaData const &pvpmd)
void writePropertyVectorMetaData(std::ostream &os, PropertyVectorMetaData const &pvmd)
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:227
std::optional< IntegrationPointMetaData > getIntegrationPointMetaData(MeshLib::Properties const &properties)
constexpr std::string_view getBulkIDString(MeshItemType mesh_item_type)
Definition Properties.h:191
static constexpr char const * toString(const MeshItemType t)
Returns a char array for a specific MeshItemType.
Definition MeshEnums.h:36
IntegrationPointMetaDataSingleField getIntegrationPointMetaDataSingleField(std::optional< IntegrationPointMetaData > const &ip_meta_data, std::string const &field_name)
std::vector< std::size_t > getIntegrationPointDataOffsetsOfMeshElements(std::vector< MeshLib::Element * > const &mesh_elements, MeshLib::PropertyVectorBase const &pv, MeshLib::Properties const &properties)
int getNumberOfElementIntegrationPoints(MeshLib::IntegrationPointMetaDataSingleField const &ip_meta_data, MeshLib::Element const &e)
std::ostream & writeConfig(std::ostream &os) const
std::ostream & writeConfig(std::ostream &os) const
std::size_t numberOfMeshItems(MeshLib::MeshItemType const item_type) const
std::vector< const MeshLib::Element * > regular_elements
Non ghost elements.
std::vector< const MeshLib::Element * > ghost_elements
std::ostream & writeNodes(std::ostream &os, std::vector< std::size_t > const &global_node_ids) const
std::vector< MeshLib::Node const * > nodes
nodes.
struct NodeData used for parallel reading and also partitioning
Definition NodeData.h:18