DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Mesh.hpp
Go to the documentation of this file.
1#pragma once
4#include "Geom/Elements.hpp"
5#include "DNDS/Array.hpp"
9#include "DNDS/ArrayPair.hpp"
10#include "Geom/PeriodicInfo.hpp"
12#include "DNDS/ObjectUtils.hpp"
14#include "AdjIndexInfo.hpp"
15#include "MeshConnectivity.hpp"
16#include "ReorderPlan.hpp"
17
18namespace DNDS::Direct
19{
21 struct DirectPrecControl;
22}
23#include "Mesh_DeviceView.hpp"
24
25namespace DNDS::Geom
26{
27 enum class MeshLoc : uint8_t
28 {
29 Unknown = 0,
30 Node = 1,
31 Face = 3,
32 Cell = 4
33 };
34
35 struct PartitionOptions; // forward declaration; defined after UnstructuredMesh
36
37 struct UnstructuredMesh : public DeviceTransferable<UnstructuredMesh>
38 {
41 int dim;
42 bool isPeriodic{false};
43 // state of: cell2node, cell2cell, bnd2node (only for non-bnd mesh object), bnd2cell (only for non-bnd mesh object)
45 // state of: face2cell, face2node, face2bnd
47 // state of: cell2face, bnd2face
49 // state of: node2cell, node2bnd
51 // state of: cell2cellFace
56 /// reader
66 tAdj1Pair bnd2bndOrig; // no device
67 /// periodic only, after reader
70
84
85 /// inverse relations
88
90 {
91 return std::make_tuple(
94 }
95
96 /// interpolated
103
104 std::vector<index> bnd2faceV; // no device
105 std::unordered_map<index, index> face2bndM; // no device
106 /// periodic only, after interpolated
108
118
120 {
121 return std::make_tuple(
124 }
125
126 /// constructed on demand
128
129 /// parent built
130 std::vector<index> node2parentNode; // from local-appended iNode to local-appended iNode in parent
131 std::vector<index> node2bndNode; // from local-appended iNode to local-appended iNode in bnd
132 std::vector<index> cell2parentCell;
133
134 /// for parallel out
135 std::vector<index> vtkCell2nodeOffsets;
136 std::vector<uint8_t> vtkCellType;
137 std::vector<index> vtkCell2node;
145 std::vector<index> nodeRecreated2nodeLocal;
146
148 {
149 size_t chunkSize = 128;
151 bool coll_on_data = false;
152 bool coll_on_meta = true;
154
155 /// only elevation
168
169 // tAdj1Pair bndFaces; // no comm needed for now
170
171 /// for cell local factorization
173
174 std::vector<index> localPartitionStarts;
175
176 /// wall dist:
178
181
182 int getDim() const { return dim; }
183
185 {
186 std::vector<Geom::NodeIndexPBI> ret;
187 ret.reserve(cell2node[iCell].size());
188 for (int ic2n = 0; ic2n < cell2node[iCell].size(); ic2n++)
189 if (isPeriodic)
191 else
193 return ret;
194 }
195
197 {
198 std::vector<Geom::NodeIndexPBI> ret;
199 ret.reserve(bnd2node[iBnd].size());
200 for (int ib2n = 0; ib2n < bnd2node[iBnd].size(); ib2n++)
201 if (isPeriodic)
203 else
205 return ret;
206 }
207
208 // =================================================================
209 // Generic index conversion templates
210 // =================================================================
211 // The 4 templates below replace 12 structurally identical methods
212 // (3 entity types x 4 variants). The old named methods are kept
213 // as thin inline wrappers so all existing call sites compile
214 // unchanged.
215
216 /**
217 * \brief Ensure a pair has a ghost mapping on its transformer.
218 *
219 * If `trans.pLGhostMapping` is already set, does nothing.
220 * Otherwise, creates a father-only ghost mapping (empty ghost set)
221 * so that `IndexGlobal2Local` / `IndexLocal2Global` work even
222 * before ghost layers are built.
223 *
224 * \pre `trans.pLGlobalMapping` must be set (call
225 * `createFatherGlobalMapping` first).
226 * \warning Collective — calls MPI_Alltoall internally when creating
227 * the mapping.
228 */
229 template <class TPair>
230 void EnsureGhostMapping(TPair &pair)
231 {
232 if (pair.trans.pLGhostMapping)
233 return;
234 DNDS_assert_info(pair.trans.pLGlobalMapping,
235 "EnsureGhostMapping: pLGlobalMapping must be set first");
236 pair.trans.pLGhostMapping = AdjIndexInfo::makeFatherOnlyMapping(
237 pair.trans.pLGlobalMapping,
238 pair.father->Size(), mpi);
239 }
240
241 /**
242 * \brief Global-to-local conversion using father+son ghost mapping.
243 * \return local index, or (-1 - iGlobal) when not found in the pair.
244 * UnInitIndex passes through unchanged.
245 */
246 template <class TPair>
248 {
249 DNDS_assert(pair.trans.pLGhostMapping);
250 if (iGlobal == UnInitIndex)
251 return iGlobal;
254 auto result = pair.trans.pLGhostMapping->search_indexAppend(iGlobal, rank, val);
255 if (result)
256 return val;
257 else
258 return -1 - iGlobal; // mapping to un-found in father-son
259 }
260
261 /**
262 * \brief Local-to-global conversion using father+son ghost mapping.
263 * \return global index, or the decoded global when the local is
264 * a negative "not-found" encoding. UnInitIndex passes through.
265 */
266 template <class TPair>
268 {
269 DNDS_assert(pair.trans.pLGhostMapping);
270 if (iLocal == UnInitIndex)
271 return iLocal;
272 if (iLocal < 0) // mapping to un-found in father-son
273 return -1 - iLocal;
274 else
275 return pair.trans.pLGhostMapping->operator()(-1, iLocal);
276 }
277
278 /**
279 * \brief Local-to-global conversion using only the father's global
280 * mapping (no son / ghost layer).
281 * \return global index. Asserts if iLocal exceeds father size.
282 * Negative "not-found" encoding and UnInitIndex pass through.
283 */
284 template <class TPair>
286 {
287 DNDS_assert(pair.father->pLGlobalMapping);
288 if (iLocal == UnInitIndex)
289 return UnInitIndex;
290 if (iLocal < 0)
291 return -1 - iLocal;
292 if (iLocal >= pair.father->Size())
293 DNDS_assert_info(false, "local idx not right: " + std::to_string(iLocal));
294 return pair.father->pLGlobalMapping->operator()(mpi.rank, iLocal);
295 }
296
297 /**
298 * \brief Global-to-local conversion using only the father's global
299 * mapping (no son / ghost layer).
300 * \return local index if the global maps to this rank, otherwise
301 * (-1 - iGlobal). UnInitIndex passes through.
302 */
303 template <class TPair>
305 {
306 DNDS_assert(pair.father->pLGlobalMapping);
307 if (iGlobal == UnInitIndex)
308 return UnInitIndex;
309 auto [ret, rank, val] = pair.father->pLGlobalMapping->search(iGlobal);
310 DNDS_assert_info(ret, "search failed with input: " + std::to_string(iGlobal));
311 if (rank == mpi.rank)
312 return val;
313 else
314 return -1 - iGlobal;
315 }
316
317 // =================================================================
318 // Named wrappers — Node
319 // =================================================================
324
325 // =================================================================
326 // Named wrappers — Cell
327 // =================================================================
332
333 // =================================================================
334 // Named wrappers — Bnd
335 // =================================================================
340
341 // =================================================================
342 // Named wrappers — Face
343 // =================================================================
346
347 // =================================================================
348 // Adjacency bulk-conversion helpers
349 // =================================================================
350
351 /**
352 * \brief Apply a conversion function to every entry of an adjacency
353 * array's first \p nRows rows.
354 *
355 * Replaces the recurring nested-loop pattern:
356 * \code
357 * for (index i = 0; i < nRows; i++)
358 * for (rowsize j = 0; j < adj.RowSize(i); j++)
359 * adj(i, j) = fn(adj(i, j));
360 * \endcode
361 *
362 * \tparam TAdj Any type with RowSize(index) and operator()(index, rowsize).
363 * \tparam TFn Callable index(index).
364 */
365 template <class TAdj, class TFn>
367 {
368 for (index i = 0; i < nRows; i++)
369 for (rowsize j = 0; j < adj.RowSize(i); j++)
370 adj(i, j) = fn(adj(i, j));
371 }
372
373 /**
374 * \brief OpenMP-parallelized variant of ConvertAdjEntries.
375 *
376 * Same semantics as ConvertAdjEntries but with `#pragma omp parallel for`
377 * over the outer row loop. Use when the adjacency array is large and
378 * the conversion function is thread-safe.
379 */
380 template <class TAdj, class TFn>
382 {
383#ifdef DNDS_USE_OMP
384# pragma omp parallel for
385#endif
386 for (index i = 0; i < nRows; i++)
387 for (rowsize j = 0; j < adj.RowSize(i); j++)
388 adj(i, j) = fn(adj(i, j));
389 }
390
391 // =================================================================
392 // Row-permutation helper
393 // =================================================================
394
395 /**
396 * \brief Permute the father rows of an ArrayPair according to a
397 * mapping function.
398 *
399 * For each old row index \p i in [0, nRows), the row is moved to
400 * position old2new(i). Works with both CSR (variable-row-size)
401 * and fixed-row-size array pairs: CSR pairs are Decompressed before
402 * and Compressed after the permutation; fixed-size pairs are
403 * copied directly.
404 *
405 * \tparam TPair An ArrayPair-like type with father/son, IsCSR().
406 * \tparam TFn Callable index(index) returning the new row index.
407 */
408 template <class TPair, class TFn>
409 static void PermuteRows(TPair &pair, index nRows, TFn &&old2new)
410 {
411 using TArr = typename decltype(pair.father)::element_type;
412 auto tmp = std::make_shared<TArr>(*pair.father); // deep copy via copy ctor
413 if constexpr (TPair::IsCSR())
414 pair.father->Decompress();
415 for (index i = 0; i < nRows; i++)
416 {
417 index iNew = old2new(i);
418 if constexpr (TPair::IsCSR())
419 pair.father->ResizeRow(iNew, tmp->RowSize(i));
420 for (rowsize j = 0; j < tmp->RowSize(i); j++)
421 pair(iNew, j) = (*tmp)(i, j);
422 }
423 if constexpr (TPair::IsCSR())
424 pair.father->Compress();
425 }
426
428 const tPoint &translation1,
429 const tPoint &rotationCenter1,
430 const tPoint &eulerAngles1,
431 const tPoint &translation2,
432 const tPoint &rotationCenter2,
433 const tPoint &eulerAngles2,
434 const tPoint &translation3,
435 const tPoint &rotationCenter3,
436 const tPoint &eulerAngles3);
437
438 /**
439 * \brief only requires father part of cell2node, bnd2node and coords
440 * generates node2cell and node2bnd (father part)
441 */
444
445 /**
446 * \brief needs to use RecoverNode2CellAndNode2Bnd before doing this.
447 * Requires node2cell.father and builds a version of its son.
448 */
451
452 /**
453 * @brief building ghost (son) from primary (currently only cell2cell)
454 * @details
455 * the face and bnd parts are currently only local (no comm available)
456 * only builds comm data of cell and node
457 * cells: current-father and cell2cell neighbor (face or node neighbor)
458 * nodes: needed by all cells
459 * faces/bnds: needed by all father cells
460 *
461 */
462 void BuildGhostPrimary(int nGhostLayers = 1);
466 // ForBnd: reduction of primary version, only on cell2node
469
472 void AdjGlobal2LocalC2F();
473 void AdjLocal2GlobalC2F();
474
475 void BuildGhostN2CB();
476 void AdjGlobal2LocalN2CB();
477 void AdjLocal2GlobalN2CB();
478
479 void AssertOnN2CB();
480
481 void BuildCell2CellFace();
484
485 void InterpolateFace();
486 void BuildGhostFace();
487 void MatchFaceBoundary();
489 void AssertOnFaces();
490
492
493 // =================================================================
494 // Registry
495 // =================================================================
496
497 /**
498 * \brief Populate a MeshConnectivity registry from this mesh's
499 * currently-built adjacencies.
500 *
501 * Registers all adjacency arrays whose father is non-null.
502 * For each entity kind that appears as a source (.from) of any
503 * registered adjacency, finds a pLGlobalMapping from any adj
504 * array for that entity kind.
505 *
506 * \pre All entity kinds that have registered adjacencies must
507 * have at least one adj array with a valid pLGlobalMapping
508 * on its father. Throws if this is not satisfied.
509 *
510 * \param dag MeshConnectivity to populate (meshDim is set).
511 */
512 void fillRegistry(MeshConnectivity &dag) const;
513
514 /// \overload Overload with an explicit skip set.
515 /// AdjKinds in \p skip are excluded from registration.
516 void fillRegistry(
518 const std::unordered_set<AdjKind, AdjKindHash> &skip) const;
519
520 // =================================================================
521 // Reorder (distributed entity reordering framework)
522 // =================================================================
523
524 /**
525 * \brief Build a ReorderRegistry containing all mesh members.
526 *
527 * Registers all built adj arrays (as type-erased callbacks) and all
528 * companion arrays (coords, cellElemInfo, pbi, etc.).
529 * Skips adjacencies involving destroyKinds.
530 *
531 * External code may extend the returned registry with its own arrays
532 * before passing to ReorderPlan::build.
533 */
535 const std::unordered_set<EntityKind> &destroyKinds = {});
536
537 /**
538 * \brief Reorder entities using the general framework.
539 *
540 * Builds a ReorderRegistry, computes follow maps, builds a
541 * ReorderPlan, applies it, rebuilds global mappings, and updates
542 * idx states.
543 *
544 * \pre All adjacencies in Adj_PointToGlobal state.
545 * \post All (non-destroyed) adjacencies in Adj_PointToGlobal.
546 * Ghost mappings stale (caller must rebuild ghosts).
547 * Global mappings fresh on reordered entities.
548 *
549 * \warning Collective.
550 */
551 void ReorderEntities(const ReorderInput &input);
552
553 /**
554 * \brief Build a ReorderPlan without applying it.
555 *
556 * Useful for external code to obtain the plan and apply it to
557 * its own arrays after the mesh reorder.
558 */
559 ReorderPlan buildReorderPlan(const ReorderInput &input);
560
561 /**
562 * \brief Augment a ReorderInput with default follows and compute follow maps.
563 *
564 * When `input.follows` is empty and Cell is explicitly reordered,
565 * default follows (Node/Bnd follow Cell) are added automatically.
566 * Returns a finalised ReorderInput with all follows resolved to
567 * explicit maps (ready for ReorderPlan::build).
568 *
569 * \param input Original input (may have empty follows).
570 * \param reg Registry with global mappings and adj data.
571 * \return ReorderInput with all entity kinds as explicit maps.
572 */
573 ReorderInput resolveFollows(
574 const ReorderInput &input,
575 const ReorderRegistry &reg);
576
577 // void ReorderCellLocal();
578
579 /**
580 * \return
581 * cell2cell for local mesh, which do not contain
582 * the diagonal part; should be a diag-less symmetric adjacency matrix
583 */
585
586 void ObtainLocalFactFillOrdering(Direct::SerialSymLUStructure &symLU, Direct::DirectPrecControl control); // 1 uses metis, 2 uses MMD, //TODO 10 uses geometric based searching
587 void ObtainSymmetricSymbolicFactorization(Direct::SerialSymLUStructure &symLU, Direct::DirectPrecControl control) const; // -1 use full LU, 0-3 use ilu(code),
588 /**
589 * \warning RecreatePeriodicNodes and BuildVTKConnectivity results are invalid after this;
590 * \warning bnd mesh's cell2parentCell is invalid after this
591 */
592 void ReorderLocalCells(int nParts = 1, int nPartsInner = 1);
593
594 /// Legacy implementation preserved for reference/fallback.
595 void ReorderLocalCellsLegacy(int nParts = 1, int nPartsInner = 1);
596
597 int NLocalParts() const { return !localPartitionStarts.empty() ? localPartitionStarts.size() - 1 : 1; }
599 index LocalPartEnd(int iPart) const { return !localPartitionStarts.empty() ? localPartitionStarts.at(iPart + 1) : this->NumCell(); }
600
602 {
603 DNDS_check_throw_info(coords.father, "coords not initialized");
604 return coords.father->Size();
605 }
607 {
608 DNDS_check_throw_info(cell2node.father, "cell2node not initialized");
609 return cell2node.father->Size();
610 }
612 {
613 DNDS_check_throw_info(face2node.father, "face2node not initialized");
614 return face2node.father->Size();
615 }
616 index NumBnd() const
617 {
618 DNDS_check_throw_info(bnd2node.father, "bnd2node not initialized");
619 return bnd2node.father->Size();
620 }
621
623 {
624 DNDS_check_throw_info(coords.son, "coords not initialized");
625 return coords.son->Size();
626 }
628 {
629 DNDS_check_throw_info(cell2node.son, "cell2node not initialized");
630 return cell2node.son->Size();
631 }
633 {
634 DNDS_check_throw_info(face2node.son, "face2node not initialized");
635 return face2node.son->Size();
636 }
638 {
639 DNDS_check_throw_info(bnd2node.son, "bnd2node not initialized");
640 return bnd2node.son->Size();
641 }
642
644 {
645 DNDS_check_throw_info(coords.father && coords.son, "coords not initialized");
646 return coords.Size();
647 }
649 {
650 DNDS_check_throw_info(cell2node.father && cell2node.son, "cell2node not initialized");
651 return cell2node.Size();
652 }
654 {
655 DNDS_check_throw_info(face2node.father && face2node.son, "face2node not initialized");
656 return face2node.Size();
657 }
659 {
660 DNDS_check_throw_info(bnd2node.father && bnd2node.son, "bnd2node not initialized");
661 return bnd2node.Size();
662 }
663
664 /// @warning must collectively call
666 {
667 DNDS_check_throw_info(cell2node.father, "cell2node not initialized");
668 return cell2node.father->globalSize();
669 }
670 /// @warning must collectively call
672 {
673 DNDS_check_throw_info(coords.father, "coords not initialized");
674 return coords.father->globalSize();
675 }
676 /// @warning must collectively call
678 {
679 DNDS_check_throw_info(face2node.father, "face2node not initialized");
680 return face2node.father->globalSize();
681 }
682 /// @warning must collectively call
684 {
685 DNDS_check_throw_info(bnd2node.father, "bnd2node not initialized");
686 return bnd2node.father->globalSize();
687 }
688
692
695 t_index GetBndZone(index iB) { return bndElemInfo(iB, 0).zone; }
696
697 MPIInfo &getMPI() { return mpi; }
698
700 void ElevatedNodesGetBoundarySmooth(const std::function<bool(t_index)> &FiFBndIdNeedSmooth);
705
707
708 bool IsO1();
709 bool IsO2();
710
711 /**
712 * @brief directly load coords; gets faulty if isPeriodic!
713 */
714 template <class tC2n>
716 {
717 cs.resize(Eigen::NoChange, c2n.size());
718 for (rowsize i = 0; i < c2n.size(); i++)
719 {
720 index iNode = c2n[i];
722 {
724 iNode = NodeIndexGlobal2Local(iNode), DNDS_assert_info(iNode >= 0, "iNode not found in main/ghost pair");
725 }
726 cs(EigenAll, i) = coords[iNode];
727 }
728 }
729
730 /**
731 * @brief directly load coords; gets faulty if isPeriodic!
732 */
733 template <class tC2n, class tCoordExt>
735 {
736 cs.resize(Eigen::NoChange, c2n.size());
737 for (rowsize i = 0; i < c2n.size(); i++)
738 {
739 index iNode = c2n[i];
741 {
743 iNode = NodeIndexGlobal2Local(iNode), DNDS_assert_info(iNode >= 0, "iNode not found in main/ghost pair");
744 }
745 cs(EigenAll, i) = coo[iNode];
746 }
747 }
748
749 /**
750 * @brief specially for periodicity
751 */
752 template <class tC2n, class tC2nPbi>
754 {
755 cs.resize(Eigen::NoChange, c2n.size());
756 for (rowsize i = 0; i < c2n.size(); i++)
757 {
758 index iNode = c2n[i];
760 {
762 iNode = NodeIndexGlobal2Local(iNode), DNDS_assert_info(iNode >= 0, "iNode not found in main/ghost pair");
763 }
765 }
766 }
767
768 /**
769 * @brief specially for periodicity
770 */
771 template <class tC2n, class tC2nPbi, class tCoordExt>
773 {
774 cs.resize(Eigen::NoChange, c2n.size());
775 for (rowsize i = 0; i < c2n.size(); i++)
776 {
777 index iNode = c2n[i];
779 {
781 iNode = NodeIndexGlobal2Local(iNode), DNDS_assert_info(iNode >= 0, "iNode not found in main/ghost pair");
782 }
783 cs(EigenAll, i) = periodicInfo.GetCoordByBits(coo[iNode], c2nPbi[i]);
784 }
785 }
786
794
802
810
817
824
831
838
840 {
842 return face2cell(iFace, 0) == iCell;
843 }
844
846 {
847 return CellIsFaceBack(iCell, iFace)
848 ? face2cell(iFace, 1)
849 : face2cell(iFace, 0);
850 }
851
852 /**
853 * @brief fA executes when if2c points to the donor side; fB the main side
854 *
855 * @tparam FA
856 * @tparam FB
857 * @tparam F0
858 * @param iFace
859 * @param if2c
860 * @param fA
861 * @param fB
862 * @param f0
863 * @return auto
864 */
865 template <class FA, class FB, class F0 = std::function<void(void)>>
867 index iFace, rowsize if2c, FA &&fA, FB &&fB, F0 &&f0 = []() {})
868 {
869 if (!this->isPeriodic)
870 return f0();
871 auto faceID = this->GetFaceZone(iFace);
873 return f0();
874 if ((if2c == 1 && Geom::FaceIDIsPeriodicMain(faceID)) ||
875 (if2c == 0 && Geom::FaceIDIsPeriodicDonor(faceID))) // I am donor
876 return fA();
877 if ((if2c == 1 && Geom::FaceIDIsPeriodicDonor(faceID)) ||
878 (if2c == 0 && Geom::FaceIDIsPeriodicMain(faceID))) // I am main
879 return fB();
880 }
881
882 void WriteSerialize(Serializer::SerializerBaseSSP serializerP, const std::string &name);
883 void ReadSerialize(Serializer::SerializerBaseSSP serializerP, const std::string &name);
884
885 /**
886 * @brief Reads mesh from an H5 serializer using even-split distribution, then
887 * repartitions via ParMetis for locality.
888 *
889 * @details
890 * This method enables reading a mesh written with any number of MPI ranks
891 * into any number of MPI ranks, without requiring the original partition.
892 *
893 * Flow:
894 * 1. Even-split read of all primary arrays (coords, cell2node, cellElemInfo,
895 * bnd2node, bndElemInfo, cell2nodePbi, bnd2nodePbi, origIndex arrays).
896 * 2. Build distributed node2cell and cell2cell (node-neighbor) using existing
897 * RecoverNode2CellAndNode2Bnd() + RecoverCell2CellAndBnd2Cell().
898 * 3. Filter cell2cell to facial neighbors only (shared O1 vertices >= dim).
899 * 4. Call ParMETIS_V3_PartKway on the distributed facial graph.
900 * 5. Redistribute cells, nodes, and boundaries to the new partition
901 * via ArrayTransformer push-based transfer.
902 * 6. Set adjPrimaryState = Adj_PointToGlobal.
903 *
904 * After return, the caller should proceed with the standard rebuild sequence:
905 * RecoverNode2CellAndNode2Bnd -> RecoverCell2CellAndBnd2Cell ->
906 * BuildGhostPrimary -> AdjGlobal2LocalPrimary -> ...
907 *
908 * @param serializerP H5 serializer (must be collective / non-per-rank)
909 * @param name Path name in the serializer hierarchy
910 * @param partitionOptions ParMetis partitioning options
911 */
914 const std::string &name,
915 const PartitionOptions &partitionOptions);
916
917 private:
918 // --- ReadSerializeAndDistribute sub-steps ---
919
920 /// Read scalar metadata and all primary arrays with H5 even-split.
921 /// Sets adjPrimaryState = Adj_PointToGlobal.
922 void ReadDistributed_EvenSplitRead(
924
925 /// Build facial cell2cell from the even-split data.
926 ssp<tAdj::element_type> ReadDistributed_BuildFacialCell2Cell();
927
928 /// Run ParMetis on the facial cell2cell graph.
929 std::vector<MPI_int> ReadDistributed_PartitionParMetis(
931 const PartitionOptions &partitionOptions);
932
933 /// Derive node and bnd partitions from cell partition.
934 struct EntityPartitions
935 {
936 std::vector<MPI_int> cellPartition;
937 std::vector<MPI_int> nodePartition;
938 std::vector<MPI_int> bndPartition;
939 };
940 EntityPartitions ReadDistributed_DeriveEntityPartitions(
941 std::vector<MPI_int> cellPartition);
942
943 /// Redistribute all primary arrays to the new partition.
944 void ReadDistributed_Redistribute(
945 const EntityPartitions &partitions);
946
947 /// Legacy implementation preserved for reference/fallback.
948 void ReadDistributed_RedistributeLegacy(
949 const EntityPartitions &partitions);
950
951 public:
952 template <class TFTrans>
953 void TransformCoords(TFTrans &&FTrans)
954 {
955 for (index iNode = 0; iNode < coords.Size(); iNode++)
957 }
958
960
962
964 std::string fname, std::string seriesName,
966 const tFGetName &names,
967 const tFGetData &data,
968 const tFGetName &vectorNames,
970 const tFGetName &namesPoint,
971 const tFGetData &dataPoint,
974 double t, MPI_Comm commDup);
975
976 void SetHDF5OutSetting(size_t chunkSiz, int deflateLevel, bool coll_on_data, bool coll_on_meta)
977 {
978 hdf5OutSetting.chunkSize = chunkSiz;
979 hdf5OutSetting.deflateLevel = deflateLevel;
980 hdf5OutSetting.coll_on_data = coll_on_data;
981 hdf5OutSetting.coll_on_meta = coll_on_meta;
982 }
983
984 void PrintMeshCGNS(std::string fname, const t_FBCID_2_Name &fbcid2name, const std::vector<std::string> &allNames);
985
987 {
989 int method = 0;
992 int verbose = 0;
993 // NOLINTNEXTLINE(modernize-use-equals-default): nested class default-init quirk with = default
995
997 {
998 // clang-format off
999 DNDS_FIELD(subdivide_quad, "Subdivide quads for wall distance computation",
1001 DNDS_FIELD(method, "Wall distance computation method (0: brute, 1: tree)",
1003 DNDS_FIELD(wallDistExecution, "MPI concurrency (0: all parallel, 1: serial, >1: batched N ranks)",
1005 DNDS_FIELD(minWallDist, "Minimum wall distance clamp",
1006 DNDS::Config::range(0.0));
1007 DNDS_FIELD(verbose, "Verbosity level for wall distance computation",
1009 // clang-format on
1010 }
1011 };
1012 void BuildNodeWallDist(const std::function<bool(Geom::t_index)> &fBndIsWall, WallDistOptions options = WallDistOptions{});
1013
1014 template <class F>
1025
1026 template <typename F>
1028 {
1029 op_on_device_arrays(std::forward<F>(f));
1030 }
1031
1032 // to_device(), to_host(), device() are provided by
1033 // DeviceTransferable<UnstructuredMesh> via for_each_device_member().
1034
1035 template <DeviceBackend B>
1037
1038 template <DeviceBackend B>
1040 {
1041 return t_deviceView<B>(*this, UnInitIndex);
1042 }
1043
1045 {
1046 index bytes = 0;
1047 auto acuumulate_bytes_arr = [&](auto &v)
1048 {
1049 if (v.ref.father)
1050 bytes += v.ref.father->FullSizeBytes();
1051 if (v.ref.son)
1052 bytes += v.ref.son->FullSizeBytes();
1053 };
1061 this->device_array_list_C2F(),
1064 this->device_array_list_N2CB(),
1067 return bytes;
1068 }
1069 };
1070
1071}
1072namespace DNDS::Geom
1073{
1074 using tFDataFieldName = std::function<std::string(int)>;
1076
1083
1085 {
1086 std::string metisType = "KWAY";
1088 int metisSeed = 0;
1089 int edgeWeightMethod = 0; // 0:no weight 1: faceSize weight
1090 int metisNcuts = 3;
1091
1093 {
1094 // clang-format off
1095 DNDS_FIELD(metisType, "METIS partitioning method",
1096 DNDS::Config::enum_values({"KWAY", "RB"}));
1097 DNDS_FIELD(metisUfactor, "METIS imbalance factor (ufactor)",
1099 DNDS_FIELD(metisSeed, "METIS random seed");
1100 DNDS_FIELD(edgeWeightMethod, "Edge weight method (0: none, 1: face size)",
1101 DNDS::Config::range(0, 1));
1102 DNDS_FIELD(metisNcuts, "Number of cuts for METIS to try",
1104 // clang-format on
1105 }
1106 };
1107
1109 {
1110 using PartitionOptions = Geom::PartitionOptions; // backward compatibility alias
1111
1112 private:
1113 int ascii_precision{16};
1114 std::string vtuFloatEncodeMode = "ascii";
1115
1116 public:
1118
1120
1121 bool dataIsSerialOut = false;
1122 bool dataIsSerialIn = false;
1123
1124 tCoord coordSerial; // created through reading
1125 tAdj cell2nodeSerial; // created through reading
1126 tAdj bnd2nodeSerial; // created through reading
1127 tElemInfoArray cellElemInfoSerial; // created through reading
1128 tElemInfoArray bndElemInfoSerial; // created through reading
1129 tAdj2 bnd2cellSerial; // created through reading
1130 tAdj1 cell2cellOrigSerial; // created through reading
1131 tAdj1 node2nodeOrigSerial; // created through reading
1132 tAdj1 bnd2bndOrigSerial; // created through reading
1133 tPbi cell2nodePbiSerial; // created through reading-Deduplicate
1134 tPbi bnd2nodePbiSerial; // created through reading-Deduplicate
1135 /***************************************************************/
1136 // Current Method: R/W don't manage actually used interpolation,
1137 // but manually get cell2cell or node2node
1138 // because: currently only support node based or cell based
1139 /***************************************************************/
1140
1141 // tAdj face2nodeSerial; // created through InterpolateTopology
1142 // tAdj2 face2cellSerial; // created through InterpolateTopology
1143 // tAdj cell2faceSerial; // created through InterpolateTopology
1144 // tElemInfoArray faceElemInfoSerial; // created through InterpolateTopology
1145
1146 tAdj cell2cellSerial; // optionally created with GetCell2Cell()
1147 tAdj node2nodeSerial; // optionally created with GetNode2Node()
1148
1149 tAdj cell2cellSerialFacial; // optionally created with GetCell2Cell()
1150
1151 tAdj node2cellSerial; // not used for now
1152 tAdj node2faceSerial; // not used for now
1153 tAdj node2edgeSerial; // not used for now
1154
1155 tAdj cell2faceSerial; // not used for now
1156 tAdj cell2edgeSerial; // not used for now
1157
1158 tAdj face2nodeSerial; // not used for now
1159 tAdj face2faceSerial; // not used for now
1160 tAdj face2edgeSerial; // not used for now
1161 tAdj face2cellSerial; // not used for now
1162
1163 tAdj edge2nodeSerial; // not used for now
1164 tAdj edge2cellSerial; // not used for now
1165 tAdj edge2edgeSerial; // not used for now
1166 tAdj edge2faceSerial; // not used for now
1167
1174
1175 std::vector<DNDS::MPI_int> cellPartition;
1176 std::vector<DNDS::MPI_int> nodePartition;
1177 std::vector<DNDS::MPI_int> bndPartition;
1178
1180
1181 /**
1182 * @brief directly load coords; gets faulty if isPeriodic!, analog of mesh's method
1183 */
1184 template <class tC2n, class tCoordExt>
1185 void _detail_GetCoordsSerial(const tC2n &c2n, tSmallCoords &cs, tCoordExt &coo)
1186 {
1187 cs.resize(Eigen::NoChange, c2n.size());
1188 for (rowsize i = 0; i < c2n.size(); i++)
1189 {
1190 index iNode = c2n[i];
1191 cs(EigenAll, i) = coo->operator[](iNode);
1192 }
1193 }
1194
1195 /**
1196 * @brief specially for periodicity, analog of mesh's method
1197 */
1198 template <class tC2n, class tC2nPbi, class tCoordExt>
1199 void _detail_GetCoordsOnElemSerial(const tC2n &c2n, const tC2nPbi &c2nPbi, tSmallCoords &cs, tCoordExt &coo)
1200 {
1201 cs.resize(Eigen::NoChange, c2n.size());
1202 for (rowsize i = 0; i < c2n.size(); i++)
1203 {
1204 index iNode = c2n[i];
1205 cs(EigenAll, i) = mesh->periodicInfo.GetCoordByBits(coo->operator[](iNode), c2nPbi[i]);
1206 }
1207 }
1208 /**
1209 * @brief analog of mesh's method
1210 */
1212 {
1213 if (!mesh->isPeriodic)
1214 _detail_GetCoordsSerial((*cell2nodeSerial)[iCell], cs, coo);
1215 else
1217 }
1218
1220 : mesh(std::move(n_mesh)), mRank(n_mRank) {}
1221
1222 /// @brief reads a cgns file as serial input
1223 /**
1224 * @details
1225 * the file MUST consist of one CGNS base node, and
1226 * multiple zones within.
1227 * All the zones are treated as a whole unstructured grid
1228 * Proofed on .cgns files generated from Pointwise
1229 * @warning //!Pointwise Options: with "Include Donor Information", "Treat Structured as Unstructured", "Unstructured Interfaces = Node-to-Node"
1230 */
1231 /// @todo //TODO Add some multi thread here!
1232 /// @param fName file name of .cgns file
1233 void ReadFromCGNSSerial(const std::string &fName, const t_FBCName_2_ID &FBCName_2_ID);
1234
1235 auto ReadFromCGNSSerial(const std::string &fName)
1236 {
1237 AutoAppendName2ID appended_name2id;
1238 this->ReadFromCGNSSerial(fName, appended_name2id);
1239 return appended_name2id;
1240 }
1241
1242 void ReadFromOpenFOAMAndConvertSerial(const std::string &fName, const std::map<std::string, std::string> &nameMapping, const t_FBCName_2_ID &FBCName_2_ID = FBC_Name_2_ID_Default);
1243
1244 void Deduplicate1to1Periodic(real search_eps = 1e-8);
1245
1246 // void InterpolateTopology();
1247
1248 /**
1249 * \brief build cell2cell topology, with node-neighbors included
1250 * \todo add support for only face-neighbors
1251 */
1252 void BuildCell2Cell(); // For cell based purpose
1253
1254 void BuildNode2Node(); // For node based purpose //!not yet implemented
1255
1256 void MeshPartitionCell2Cell(const PartitionOptions &options);
1257
1259
1261 {
1262 coordSerial.reset();
1263 cell2nodeSerial.reset();
1264 cell2cellSerial.reset();
1265 cellElemInfoSerial.reset();
1266 bnd2nodeSerial.reset();
1267 bnd2cellSerial.reset();
1268 bndElemInfoSerial.reset();
1269 cell2cellOrigSerial.reset();
1270 node2nodeOrigSerial.reset();
1271 bnd2bndOrigSerial.reset();
1272 mode = UnknownMode;
1274 }
1275
1276 /**
1277 * @brief should be called to build data for serial out
1278 */
1279 void BuildSerialOut();
1280
1281 void GetCurrentOutputArrays(int flag,
1282 tCoordPair &coordOut,
1283 tAdjPair &cell2nodeOut,
1284 tPbiPair &cell2nodePbiOut,
1285 tElemInfoArrayPair &cellElemInfoOut,
1286 index &nCell, index &nNode);
1287
1288 // void WriteToCGNSSerial(const std::string &fName);
1289
1290 // void WriteSolutionToTecBinary(const std::string &fName,
1291 // int nField, const tFDataFieldName &names, const tFDataFieldQuery &data,
1292 // int nFieldBnd, const tFDataFieldName &namesBnd, const tFDataFieldQuery &dataBnd);
1293
1294 /**
1295 * @brief names(idata) data(idata, ivolume)
1296 * https://tecplot.azureedge.net/products/360/current/360_data_format_guide.pdf
1297 * @todo //TODO add support for bnd export
1298 */
1300 std::string fname,
1301 int arraySiz, int arraySizPoint,
1302 const tFGetName &names,
1303 const tFGetData &data,
1304 const tFGetName &namesPoint,
1305 const tFGetData &dataPoint,
1306 double t, int flag);
1307
1308 /**
1309 * @brief names(idata) data(idata, ivolume)
1310 * @todo //TODO add support for bnd export
1311 */
1313 std::string fname, std::string seriesName,
1314 int arraySiz, int vecArraySiz, int arraySizPoint, int vecArraySizPoint,
1315 const tFGetName &names,
1316 const tFGetData &data,
1317 const tFGetName &vectorNames,
1318 const tFGetVecData &vectorData,
1319 const tFGetName &namesPoint,
1320 const tFGetData &dataPoint,
1321 const tFGetName &vectorNamesPoint,
1322 const tFGetVecData &vectorDataPoint,
1323 double t, int flag = 0);
1324
1325 void SetASCIIPrecision(int n) { ascii_precision = n; }
1326 void SetVTKFloatEncodeMode(const std::string &v)
1327 {
1328 vtuFloatEncodeMode = v;
1329 DNDS_assert(v == "ascii" || v == "binary");
1330 }
1331 };
1332} // namespace geom
Adjacency array (CSR-like index storage) built on ParArray.
Eigen-vector array: each row is an Eigen::Map over contiguous real storage.
Father-son array pairs with device views and ghost communication.
Core 2D variable-length array container with five data layouts.
pybind11-style configuration registration with macro-based field declaration and namespace-scoped tag...
#define DNDS_FIELD(name_, desc_,...)
Register a field inside a DNDS_DECLARE_CONFIG body.
Device memory abstraction layer with backend-specific storage and factory creation.
CRTP mixin that provides uniform to_device/to_host/device/getDeviceArrayBytes for classes that enumer...
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:117
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:112
#define DNDS_check_throw_info(expr, info)
Same as DNDS_check_throw but attaches a user-supplied info message to the thrown std::runtime_error.
Definition Errors.hpp:100
Layered DAG of mesh adjacency relations with composable DSL operations.
Tiny reflection-style helpers (MemberRef, MemberPtr) and for_each_member_* visitors used by config / ...
#define DNDS_MAKE_1_MEMBER_REF(x)
Construct a MemberRef capturing x and its stringified name.
Distributed entity reordering framework: ReorderRegistry, ReorderPlan, ReorderInput.
Ghost-communication engine for a father / son ParArray pair.
CRTP mixin giving a class uniform to_device / to_host / device / getDeviceArrayBytes methods.
EnumValuesTag enum_values(std::vector< std::string > vals)
Create an enum allowed-values tag.
RangeTag range(double min)
Create a minimum-only range constraint.
decltype(tAdj1Pair::father) tAdj1
std::function< std::string(int)> tFGetName
tFGetData tFDataFieldQuery
Definition Mesh.hpp:1075
int32_t t_index
Definition Geometric.hpp:6
@ SerialOutput
Definition Mesh.hpp:1081
@ UnknownMode
Definition Mesh.hpp:1079
@ SerialReadAndDistribute
Definition Mesh.hpp:1080
std::function< std::string(int)> tFDataFieldName
Definition Mesh.hpp:1074
decltype(tAdj2Pair::father) tAdj2
decltype(tAdjPair::father) tAdj
std::function< DNDS::real(int, DNDS::index)> tFGetData
std::function< DNDS::real(int, DNDS::index, DNDS::rowsize)> tFGetVecData
Eigen::Vector3d tPoint
Definition Geometric.hpp:9
DNDS::ssp< DNDS::ParArray< ElemInfo > > tElemInfoArray
std::function< std::string(t_index)> t_FBCID_2_Name
Eigen::Matrix< real, 3, Eigen::Dynamic > tSmallCoords
Definition Geometric.hpp:50
DNDS_DEVICE_CALLABLE bool FaceIDIsPeriodic(t_index id)
DNDS_DEVICE_CALLABLE bool FaceIDIsPeriodicDonor(t_index id)
DNDS_DEVICE_CALLABLE bool FaceIDIsPeriodicMain(t_index id)
std::function< t_index(const std::string &)> t_FBCName_2_ID
decltype(tPbiPair::father) tPbi
decltype(tCoordPair::father) tCoord
std::vector< std::vector< index > > tLocalMatStruct
Definition Geometric.hpp:66
void AllreduceOneIndex(index &v, MPI_Op op, const MPIInfo &mpi)
Single-scalar Allreduce helper for indices (in-place, count = 1).
Definition MPI.hpp:727
ssp< SerializerBase > SerializerBaseSSP
DNDS_CONSTANT const index UnInitIndex
Sentinel "not initialised" index value (= INT64_MIN).
Definition Defines.hpp:181
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:114
constexpr MPI_int UnInitMPIInt
Sentinel "not initialised" MPI_int value (= -1, invalid rank).
Definition MPI.hpp:66
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
Definition Defines.hpp:143
void for_each_member_list(TList &&obj_member_list, F &&f)
Invoke f(member) for every element of a std::tuple of members.
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:110
int MPI_int
MPI counterpart type for MPI_int (= C int). Used for counts and ranks in MPI calls.
Definition MPI.hpp:54
ssp< TArray > father
Owned-side array (must be resized before ghost setup).
index Size() const
Combined row count (father->Size() + son->Size()).
ssp< TArray > son
Ghost-side array (sized automatically by createMPITypes / BorrowAndPull).
static t_pLGhostMapping makeFatherOnlyMapping(const ssp< GlobalOffsetsMapping > &globalMapping, index fatherSize, const MPIInfo &mpi)
Create a ghost mapping with no ghost entries (father-only).
DNDS_DECLARE_CONFIG(PartitionOptions)
Definition Mesh.hpp:1092
DNDS_DEVICE_CALLABLE auto GetVectorByBits(const Eigen::Matrix< real, dim, nVec > &v, const NodePeriodicBits &bits)
DNDS_DEVICE_CALLABLE tPoint GetCoordByBits(const tPoint &c, const NodePeriodicBits &bits) const
void MeshPartitionCell2Cell(const PartitionOptions &options)
void GetCoordsOnCellSerial(index iCell, tSmallCoords &cs, tCoord &coo)
analog of mesh's method
Definition Mesh.hpp:1211
DNDS::ArrayTransformerType< tCoord::element_type >::Type coordSerialOutTrans
Definition Mesh.hpp:1168
void GetCurrentOutputArrays(int flag, tCoordPair &coordOut, tAdjPair &cell2nodeOut, tPbiPair &cell2nodePbiOut, tElemInfoArrayPair &cellElemInfoOut, index &nCell, index &nNode)
Definition Mesh_Plts.cpp:21
void ReadFromOpenFOAMAndConvertSerial(const std::string &fName, const std::map< std::string, std::string > &nameMapping, const t_FBCName_2_ID &FBCName_2_ID=FBC_Name_2_ID_Default)
DNDS::ArrayTransformerType< tElemInfoArray::element_type >::Type cellElemInfoSerialOutTrans
Definition Mesh.hpp:1172
void ReadFromCGNSSerial(const std::string &fName, const t_FBCName_2_ID &FBCName_2_ID)
reads a cgns file as serial input
void BuildSerialOut()
should be called to build data for serial out
Definition Mesh.cpp:193
void _detail_GetCoordsOnElemSerial(const tC2n &c2n, const tC2nPbi &c2nPbi, tSmallCoords &cs, tCoordExt &coo)
specially for periodicity, analog of mesh's method
Definition Mesh.hpp:1199
DNDS::ArrayTransformerType< tAdj::element_type >::Type cell2nodeSerialOutTrans
Definition Mesh.hpp:1169
void SetVTKFloatEncodeMode(const std::string &v)
Definition Mesh.hpp:1326
DNDS::ssp< UnstructuredMesh > mesh
Definition Mesh.hpp:1117
void PrintSerialPartVTKDataArray(std::string fname, std::string seriesName, int arraySiz, int vecArraySiz, int arraySizPoint, int vecArraySizPoint, const tFGetName &names, const tFGetData &data, const tFGetName &vectorNames, const tFGetVecData &vectorData, const tFGetName &namesPoint, const tFGetData &dataPoint, const tFGetName &vectorNamesPoint, const tFGetVecData &vectorDataPoint, double t, int flag=0)
names(idata) data(idata, ivolume)
DNDS::ArrayTransformerType< tPbi::element_type >::Type cell2nodePbiSerialOutTrans
Definition Mesh.hpp:1170
DNDS::ArrayTransformerType< tAdj::element_type >::Type bnd2nodeSerialOutTrans
Definition Mesh.hpp:1171
DNDS::ArrayTransformerType< tElemInfoArray::element_type >::Type bndElemInfoSerialOutTrans
Definition Mesh.hpp:1173
auto ReadFromCGNSSerial(const std::string &fName)
Definition Mesh.hpp:1235
void PrintSerialPartPltBinaryDataArray(std::string fname, int arraySiz, int arraySizPoint, const tFGetName &names, const tFGetData &data, const tFGetName &namesPoint, const tFGetData &dataPoint, double t, int flag)
names(idata) data(idata, ivolume) https://tecplot.azureedge.net/products/360/current/360_data_format_...
void BuildCell2Cell()
build cell2cell topology, with node-neighbors included
std::vector< DNDS::MPI_int > cellPartition
Definition Mesh.hpp:1175
std::vector< DNDS::MPI_int > bndPartition
Definition Mesh.hpp:1177
void _detail_GetCoordsSerial(const tC2n &c2n, tSmallCoords &cs, tCoordExt &coo)
directly load coords; gets faulty if isPeriodic!, analog of mesh's method
Definition Mesh.hpp:1185
UnstructuredMeshSerialRW(decltype(mesh) n_mesh, DNDS::MPI_int n_mRank)
Definition Mesh.hpp:1219
std::vector< DNDS::MPI_int > nodePartition
Definition Mesh.hpp:1176
index NumFaceGhost() const
Definition Mesh.hpp:632
tCoordPair nodeWallDist
wall dist:
Definition Mesh.hpp:177
MeshAdjState adjC2FState
Definition Mesh.hpp:48
UnstructuredMeshDeviceView< B > t_deviceView
Definition Mesh.hpp:1036
index BndIndexGlobal2Local_NoSon(index i)
Definition Mesh.hpp:339
index NodeIndexGlobal2Local_NoSon(index i)
Definition Mesh.hpp:323
tCoordPair coordsPeriodicRecreated
Definition Mesh.hpp:144
Elem::Element GetCellElement(index iC)
Definition Mesh.hpp:689
tPbiPair cell2nodePbi
periodic only, after reader
Definition Mesh.hpp:68
AdjPairTracked< tAdjPair > cell2cellFace
constructed on demand
Definition Mesh.hpp:127
UnstructuredMesh(const DNDS::MPIInfo &n_mpi, int n_dim)
Definition Mesh.hpp:179
tCoordPair coordsElevDisp
only elevation
Definition Mesh.hpp:156
void fillRegistry(MeshConnectivity &dag) const
Populate a MeshConnectivity registry from this mesh's currently-built adjacencies.
Definition Mesh.cpp:1743
std::unordered_map< index, index > face2bndM
Definition Mesh.hpp:105
void ReorderLocalCellsLegacy(int nParts=1, int nPartsInner=1)
Legacy implementation preserved for reference/fallback.
Definition Mesh.cpp:1846
void GetCoordsOnCell(index iCell, tSmallCoords &cs, tCoordPair &coo)
Definition Mesh.hpp:795
index NumNodeProc() const
Definition Mesh.hpp:643
index CellIndexLocal2Global_NoSon(index i)
Definition Mesh.hpp:330
tPbiPair face2nodePbi
periodic only, after interpolated
Definition Mesh.hpp:107
AdjPairTracked< tAdj1Pair > face2bnd
Definition Mesh.hpp:102
static void ConvertAdjEntriesOMP(TAdj &adj, index nRows, TFn &&fn)
OpenMP-parallelized variant of ConvertAdjEntries.
Definition Mesh.hpp:381
std::vector< index > node2bndNode
Definition Mesh.hpp:131
tPoint GetCoordNodeOnFace(index iFace, rowsize if2n)
Definition Mesh.hpp:818
auto CellOtherCellPeriodicHandle(index iFace, rowsize if2c, FA &&fA, FB &&fB, F0 &&f0=[]() {})
fA executes when if2c points to the donor side; fB the main side
Definition Mesh.hpp:866
index NumNodeGhost() const
Definition Mesh.hpp:622
std::vector< index > cell2parentCell
Definition Mesh.hpp:132
void GetCoordsOnCell(index iCell, tSmallCoords &cs)
Definition Mesh.hpp:787
std::vector< uint8_t > vtkCellType
Definition Mesh.hpp:136
MeshAdjState adjC2CFaceState
Definition Mesh.hpp:52
std::vector< index > node2parentNode
parent built
Definition Mesh.hpp:130
auto getBnd2NodeIndexPbiRow(index iBnd)
Definition Mesh.hpp:196
tAdjPair cell2nodePeriodicRecreated
Definition Mesh.hpp:143
AdjPairTracked< tAdj2Pair > face2cell
Definition Mesh.hpp:99
void _detail_GetCoords(const tC2n &c2n, tSmallCoords &cs, tCoordExt &coo)
directly load coords; gets faulty if isPeriodic!
Definition Mesh.hpp:734
index LocalPartStart(int iPart) const
Definition Mesh.hpp:598
index NumCellProc() const
Definition Mesh.hpp:648
index NodeIndexLocal2Global_NoSon(index i)
Definition Mesh.hpp:322
index FaceIndexGlobal2Local(DNDS::index i)
Definition Mesh.hpp:344
index NodeIndexGlobal2Local(DNDS::index i)
Definition Mesh.hpp:320
tLocalMatStruct cell2cellFaceVLocalParts
for cell local factorization
Definition Mesh.hpp:172
void _detail_GetCoordsOnElem(const tC2n &c2n, const tC2nPbi &c2nPbi, tSmallCoords &cs)
specially for periodicity
Definition Mesh.hpp:753
void ReorderEntities(const ReorderInput &input)
Reorder entities using the general framework.
std::vector< index > vtkCell2nodeOffsets
for parallel out
Definition Mesh.hpp:135
void ElevatedNodesGetBoundarySmooth(const std::function< bool(t_index)> &FiFBndIdNeedSmooth)
ReorderInput resolveFollows(const ReorderInput &input, const ReorderRegistry &reg)
Augment a ReorderInput with default follows and compute follow maps.
index BndIndexLocal2Global_NoSon(index i)
Definition Mesh.hpp:338
tCoordPair coords
reader
Definition Mesh.hpp:57
index CellIndexGlobal2Local(DNDS::index i)
Definition Mesh.hpp:328
tLocalMatStruct GetCell2CellFaceVLocal(bool onLocalPartition=false)
index IndexGlobal2Local(TPair &pair, DNDS::index iGlobal)
Global-to-local conversion using father+son ghost mapping.
Definition Mesh.hpp:247
void RecoverNode2CellAndNode2Bnd()
only requires father part of cell2node, bnd2node and coords generates node2cell and node2bnd (father ...
Definition Mesh.cpp:379
MeshAdjState adjN2CBState
Definition Mesh.hpp:50
void for_each_device_member(F &&f)
Definition Mesh.hpp:1027
void ReorderLocalCells(int nParts=1, int nPartsInner=1)
index IndexLocal2Global_NoSon(TPair &pair, index iLocal)
Local-to-global conversion using only the father's global mapping (no son / ghost layer).
Definition Mesh.hpp:285
AdjPairTracked< tAdjPair > cell2face
interpolated
Definition Mesh.hpp:97
void op_on_device_arrays(F &&f)
Definition Mesh.hpp:1015
index NumFaceProc() const
Definition Mesh.hpp:653
struct DNDS::Geom::UnstructuredMesh::HDF5OutSetting hdf5OutSetting
void BuildBisectO1FormO2(UnstructuredMesh &meshO2)
void BuildGhostPrimary(int nGhostLayers=1)
building ghost (son) from primary (currently only cell2cell)
Definition Mesh.cpp:667
void GetCoordsOnFace(index iFace, tSmallCoords &cs)
Definition Mesh.hpp:803
void SetHDF5OutSetting(size_t chunkSiz, int deflateLevel, bool coll_on_data, bool coll_on_meta)
Definition Mesh.hpp:976
ReorderRegistry buildReorderRegistry(const std::unordered_set< EntityKind > &destroyKinds={})
Build a ReorderRegistry containing all mesh members.
tPoint GetCoordNodeOnCell(index iCell, rowsize ic2n)
Definition Mesh.hpp:811
void ReadSerialize(Serializer::SerializerBaseSSP serializerP, const std::string &name)
Definition Mesh.cpp:1467
void PrintParallelVTKHDFDataArray(std::string fname, std::string seriesName, int arraySiz, int vecArraySiz, int arraySizPoint, int vecArraySizPoint, const tFGetName &names, const tFGetData &data, const tFGetName &vectorNames, const tFGetVecData &vectorData, const tFGetName &namesPoint, const tFGetData &dataPoint, const tFGetName &vectorNamesPoint, const tFGetVecData &vectorDataPoint, double t, MPI_Comm commDup)
auto getCell2NodeIndexPbiRow(index iCell)
Definition Mesh.hpp:184
index CellIndexGlobal2Local_NoSon(index i)
Definition Mesh.hpp:331
bool CellIsFaceBack(index iCell, index iFace) const
Definition Mesh.hpp:839
t_index GetFaceZone(index iF)
Definition Mesh.hpp:694
Elem::Element GetBndElement(index iB)
Definition Mesh.hpp:691
index BndIndexLocal2Global(DNDS::index i)
Definition Mesh.hpp:337
index FaceIndexLocal2Global(DNDS::index i)
Definition Mesh.hpp:345
AdjPairTracked< tAdjPair > cell2cell
Definition Mesh.hpp:61
std::vector< index > bnd2faceV
Definition Mesh.hpp:104
Elem::Element GetFaceElement(index iF)
Definition Mesh.hpp:690
void _detail_GetCoordsOnElem(const tC2n &c2n, const tC2nPbi &c2nPbi, tSmallCoords &cs, tCoordExt &coo)
specially for periodicity
Definition Mesh.hpp:772
index CellIndexLocal2Global(DNDS::index i)
Definition Mesh.hpp:329
index IndexLocal2Global(TPair &pair, DNDS::index iLocal)
Local-to-global conversion using father+son ghost mapping.
Definition Mesh.hpp:267
void _detail_GetCoords(const tC2n &c2n, tSmallCoords &cs)
directly load coords; gets faulty if isPeriodic!
Definition Mesh.hpp:715
AdjPairTracked< tAdjPair > face2node
Definition Mesh.hpp:98
MeshAdjState adjPrimaryState
Definition Mesh.hpp:44
void ObtainSymmetricSymbolicFactorization(Direct::SerialSymLUStructure &symLU, Direct::DirectPrecControl control) const
Definition Mesh.cpp:1734
AdjPairTracked< tAdj1Pair > bnd2face
Definition Mesh.hpp:101
index CellFaceOther(index iCell, index iFace) const
Definition Mesh.hpp:845
tElemInfoArrayPair faceElemInfo
Definition Mesh.hpp:100
void SetPeriodicGeometry(const tPoint &translation1, const tPoint &rotationCenter1, const tPoint &eulerAngles1, const tPoint &translation2, const tPoint &rotationCenter2, const tPoint &eulerAngles2, const tPoint &translation3, const tPoint &rotationCenter3, const tPoint &eulerAngles3)
Definition Mesh.cpp:347
struct DNDS::Geom::UnstructuredMesh::ElevationInfo elevationInfo
std::vector< index > localPartitionStarts
Definition Mesh.hpp:174
index NumBndProc() const
Definition Mesh.hpp:658
MeshAdjState adjFacialState
Definition Mesh.hpp:46
index BndIndexGlobal2Local(DNDS::index i)
Definition Mesh.hpp:336
void BuildO2FromO1Elevation(UnstructuredMesh &meshO1)
AdjPairTracked< tAdjPair > cell2node
Definition Mesh.hpp:58
void RecoverCell2CellAndBnd2Cell()
needs to use RecoverNode2CellAndNode2Bnd before doing this. Requires node2cell.father and builds a ve...
Definition Mesh.cpp:422
void WriteSerialize(Serializer::SerializerBaseSSP serializerP, const std::string &name)
Definition Mesh.cpp:1429
tElemInfoArrayPair cellElemInfo
Definition Mesh.hpp:62
void ObtainLocalFactFillOrdering(Direct::SerialSymLUStructure &symLU, Direct::DirectPrecControl control)
void PrintMeshCGNS(std::string fname, const t_FBCID_2_Name &fbcid2name, const std::vector< std::string > &allNames)
static void ConvertAdjEntries(TAdj &adj, index nRows, TFn &&fn)
Apply a conversion function to every entry of an adjacency array's first nRows rows.
Definition Mesh.hpp:366
AdjPairTracked< tAdjPair > node2bnd
Definition Mesh.hpp:87
index NumCellGhost() const
Definition Mesh.hpp:627
index NumBndGhost() const
Definition Mesh.hpp:637
void ReadSerializeAndDistribute(Serializer::SerializerBaseSSP serializerP, const std::string &name, const PartitionOptions &partitionOptions)
Reads mesh from an H5 serializer using even-split distribution, then repartitions via ParMetis for lo...
static void PermuteRows(TPair &pair, index nRows, TFn &&old2new)
Permute the father rows of an ArrayPair according to a mapping function.
Definition Mesh.hpp:409
AdjPairTracked< tAdjPair > bnd2node
Definition Mesh.hpp:59
AdjPairTracked< tAdjPair > node2cell
inverse relations
Definition Mesh.hpp:86
index LocalPartEnd(int iPart) const
Definition Mesh.hpp:599
tPoint GetCoordWallDistOnFace(index iFace, rowsize if2n)
Definition Mesh.hpp:832
ReorderPlan buildReorderPlan(const ReorderInput &input)
Build a ReorderPlan without applying it.
t_index GetCellZone(index iC)
Definition Mesh.hpp:693
index IndexGlobal2Local_NoSon(TPair &pair, index iGlobal)
Global-to-local conversion using only the father's global mapping (no son / ghost layer).
Definition Mesh.hpp:304
void ConstructBndMesh(UnstructuredMesh &bMesh)
Definition Mesh.cpp:1560
tElemInfoArrayPair bndElemInfo
Definition Mesh.hpp:63
tPoint GetCoordWallDistOnCell(index iCell, rowsize ic2n)
Definition Mesh.hpp:825
t_index GetBndZone(index iB)
Definition Mesh.hpp:695
void EnsureGhostMapping(TPair &pair)
Ensure a pair has a ghost mapping on its transformer.
Definition Mesh.hpp:230
MeshElevationState elevState
Definition Mesh.hpp:55
void BuildNodeWallDist(const std::function< bool(Geom::t_index)> &fBndIsWall, WallDistOptions options=WallDistOptions{})
index NodeIndexLocal2Global(DNDS::index i)
Definition Mesh.hpp:321
std::vector< index > vtkCell2node
Definition Mesh.hpp:137
std::vector< index > nodeRecreated2nodeLocal
Definition Mesh.hpp:145
AdjPairTracked< tAdj2Pair > bnd2cell
Definition Mesh.hpp:60
void TransformCoords(TFTrans &&FTrans)
Definition Mesh.hpp:953
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:231
int rank
This rank's 0-based index within comm (-1 until initialised).
Definition MPI.hpp:235
Eigen::Matrix< real, 5, 1 > v
auto adj
auto result
auto reg
std::vector< MPI_int > nodePartition(nNode, mpi.rank)
std::vector< MPI_int > cellPartition(nCellBefore)
const DNDS::index nCell
const DNDS::index nNode
ReorderInput input
Eigen::Vector3d n(1.0, 0.0, 0.0)