Loading...
Searching...
No Matches
patch_kernel.hpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23\*---------------------------------------------------------------------------*/
24
25#ifndef __BITPIT_PATCH_KERNEL_HPP__
26#define __BITPIT_PATCH_KERNEL_HPP__
27
28#include <cstddef>
29#include <deque>
30#include <iostream>
31#include <memory>
32#if BITPIT_ENABLE_MPI==1
33# include <mpi.h>
34#endif
35#include <set>
36#include <string>
37#include <vector>
38#include <unordered_map>
39
40#include "bitpit_IO.hpp"
41#if BITPIT_ENABLE_MPI==1
42# include "bitpit_communications.hpp"
43#endif
44#include "bitpit_containers.hpp"
45
46#include "adaption.hpp"
47#include "cell.hpp"
48#include "compiler.hpp"
49#include "interface.hpp"
50#include "vertex.hpp"
51
52namespace bitpit {
53
55
56friend class PatchInfo;
57friend class PatchNumberingInfo;
58#if BITPIT_ENABLE_MPI==1
59friend class PatchGlobalInfo;
60#endif
61friend class PatchManager;
62
63public:
64 typedef PiercedVector<Vertex>::iterator VertexIterator;
65 typedef PiercedVector<Cell>::iterator CellIterator;
66 typedef PiercedVector<Interface>::iterator InterfaceIterator;
67
68 typedef PiercedVector<Vertex>::const_iterator VertexConstIterator;
69 typedef PiercedVector<Cell>::const_iterator CellConstIterator;
70 typedef PiercedVector<Interface>::const_iterator InterfaceConstIterator;
71
72 typedef PiercedVector<Vertex>::range VertexRange;
73 typedef PiercedVector<Cell>::range CellRange;
74 typedef PiercedVector<Interface>::range InterfaceRange;
75
76 typedef PiercedVector<Vertex>::const_range VertexConstRange;
77 typedef PiercedVector<Cell>::const_range CellConstRange;
78 typedef PiercedVector<Interface>::const_range InterfaceConstRange;
79
80 enum WriteTarget {
81 WRITE_TARGET_CELLS_ALL
82#if BITPIT_ENABLE_MPI
83 , WRITE_TARGET_CELLS_INTERNAL
84#endif
85 };
86
94 struct PointPositionLess
95 {
96 PointPositionLess(double tolerance)
97 : m_tolerance(tolerance)
98 {
99 }
100
101 virtual ~PointPositionLess() = default;
102
103 bool operator()(std::array<double, 3> point_1, std::array<double, 3> point_2) const
104 {
105 // Check if the points are coincident
106 //
107 // To check if two points are coincident, the specified tolerance
108 // will be used.
109 bool coincident = true;
110 for (int k = 0; k < 3; ++k) {
111 if (!utils::DoubleFloatingEqual()(point_1[k], point_2[k], m_tolerance)) {
112 coincident = false;
113 break;
114 }
115 }
116
117 if (coincident) {
118 return false;
119 }
120
121 // Compare the position of the points
122 //
123 // If the points are not coincident, we loop over the coordinates
124 // and we compare the first non-coincident coordinate.
125 //
126 // In order to guarantee a consistent ordering, the test to check
127 // if two coordinates are coincident should be performed with no
128 // tolerance.
129 for (int k = 0; k < 3; ++k) {
130 if (point_1[k] == point_2[k]) {
131 continue;
132 }
133
134 return (point_1[k] < point_2[k]);
135 }
136
137 return false;
138 }
139
140 const double m_tolerance;
141 };
142
148 struct VertexPositionLess : private PointPositionLess
149 {
150 VertexPositionLess(const PatchKernel &patch)
151 : PointPositionLess(patch.getTol()), m_patch(patch)
152 {
153 }
154
155 virtual ~VertexPositionLess() = default;
156
157 bool operator()(long id_1, long id_2) const
158 {
159 if (id_1 == id_2) {
160 return false;
161 }
162
163 const std::array<double, 3> &coords_1 = m_patch.getVertexCoords(id_1);
164 const std::array<double, 3> &coords_2 = m_patch.getVertexCoords(id_2);
165
166 return PointPositionLess::operator()(coords_1, coords_2);
167 }
168
169 const PatchKernel &m_patch;
170 };
171
177 struct VertexPositionGreater : private VertexPositionLess
178 {
179 VertexPositionGreater(const PatchKernel &patch)
180 : VertexPositionLess(patch)
181 {
182 }
183
184 bool operator()(long id_1, long id_2) const
185 {
186 return !VertexPositionLess::operator()(id_1, id_2);
187 }
188 };
189
195 struct CellPositionLess : private PointPositionLess
196 {
197 CellPositionLess(const PatchKernel &patch, bool native = true)
198 : PointPositionLess(patch.getTol()), m_patch(patch), m_native(native)
199 {
200 }
201
202 virtual ~CellPositionLess() = default;
203
204 bool operator()(long id_1, long id_2) const
205 {
206 if (id_1 == id_2) {
207 return false;
208 }
209
210 std::array<double, 3> centroid_1;
211 std::array<double, 3> centroid_2;
212 if (m_native) {
213 centroid_1 = m_patch.evalCellCentroid(id_1);
214 centroid_2 = m_patch.evalCellCentroid(id_2);
215 } else {
216 centroid_1 = m_patch.PatchKernel::evalCellCentroid(id_1);
217 centroid_2 = m_patch.PatchKernel::evalCellCentroid(id_2);
218 }
219
220 return PointPositionLess::operator()(centroid_1, centroid_2);
221 }
222
223 const PatchKernel &m_patch;
224 bool m_native;
225 };
226
232 struct CellPositionGreater : private CellPositionLess
233 {
234 CellPositionGreater(const PatchKernel &patch, bool native = true)
235 : CellPositionLess(patch, native)
236 {
237 }
238
239 bool operator()(long id_1, long id_2) const
240 {
241 return !CellPositionLess::operator()(id_1, id_2);
242 }
243 };
244
253 struct CellFuzzyPositionLess : private PointPositionLess
254 {
255 CellFuzzyPositionLess(PatchKernel &patch)
256 : PointPositionLess(patch.getTol()), m_patch(patch)
257 {
258 }
259
260 virtual ~CellFuzzyPositionLess() = default;
261
262 bool operator()(long id_1, long id_2) const
263 {
264 if (id_1 == id_2) {
265 return false;
266 }
267
268 // Select the first vertex of the first cell
269 long vertexId_1 = m_patch.getCell(id_1).getVertexId(0);
270
271 // The vertex of the second cell is choosen as the first vertex on
272 // that cell not equal to the selected vertex of the first cell.
273 long vertexId_2 = Vertex::NULL_ID;
274 for (long candidateVertexId_2 : m_patch.getCell(id_2).getVertexIds()) {
275 if (vertexId_1 != candidateVertexId_2) {
276 vertexId_2 = candidateVertexId_2;
277 break;
278 }
279 }
280
281 if (vertexId_2 == Vertex::NULL_ID) {
282 return false;
283 }
284
285 // Compare the two vertices
286 const std::array<double, 3> &vertexCoords_1 = m_patch.getVertex(vertexId_1).getCoords();
287 const std::array<double, 3> &vertexCoords_2 = m_patch.getVertex(vertexId_2).getCoords();
288
289 return PointPositionLess::operator()(vertexCoords_1, vertexCoords_2);
290 }
291
292 PatchKernel &m_patch;
293 };
294
303 struct CellFuzzyPositionGreater : private CellFuzzyPositionLess
304 {
305 CellFuzzyPositionGreater(PatchKernel &patch)
306 : CellFuzzyPositionLess(patch)
307 {
308 }
309
310 bool operator()(long id_1, long id_2) const
311 {
312 return !CellFuzzyPositionLess::operator()(id_1, id_2);
313 }
314 };
315
320 ADJACENCIES_NONE = -1,
321 ADJACENCIES_AUTOMATIC
322 };
323
328 INTERFACES_NONE = -1,
329 INTERFACES_AUTOMATIC
330 };
331
336 ADAPTION_DISABLED = -1,
343 };
344
349 ADAPTION_CLEAN,
350 ADAPTION_DIRTY,
351 ADAPTION_PREPARED,
352 ADAPTION_ALTERED
353 };
354
359 PARTITIONING_DISABLED = -1,
361 };
362
367 PARTITIONING_CLEAN,
368 PARTITIONING_PREPARED,
369 PARTITIONING_ALTERED
370 };
371
372 PatchKernel(PatchKernel &&other);
374
375 virtual ~PatchKernel();
376
377 template<typename patch_t>
378 static std::unique_ptr<patch_t> clone(const patch_t *original);
379
380 virtual std::unique_ptr<PatchKernel> clone() const = 0;
381
382 void setId(int id);
383
384 virtual void reset();
385 virtual void resetVertices();
386 virtual void resetCells();
387 virtual void resetInterfaces();
388
389 bool reserveVertices(size_t nVertices);
390 bool reserveCells(size_t nCells);
391 bool reserveInterfaces(size_t nInterfaces);
392
393 std::vector<adaption::Info> update(bool trackAdaption = true, bool squeezeStorage = false);
394
395 virtual void simulateCellUpdate(const long id, adaption::Marker marker, std::vector<Cell> *virtualCells, PiercedVector<Vertex, long> *virtualVertices) const;
396
397 bool isAdaptionSupported() const;
399 AdaptionStatus getAdaptionStatus(bool global = false) const;
400 std::vector<adaption::Info> adaption(bool trackAdaption = true, bool squeezeStorage = false);
401 std::vector<adaption::Info> adaptionPrepare(bool trackAdaption = true);
402 std::vector<adaption::Info> adaptionAlter(bool trackAdaption = true, bool squeezeStorage = false);
403 void adaptionCleanup();
404
405 virtual void settleAdaptionMarkers();
406
407 void markCellForRefinement(long id);
408 void markCellForCoarsening(long id);
409 void resetCellAdaptionMarker(long id);
410 adaption::Marker getCellAdaptionMarker(long id);
411 void enableCellBalancing(long id, bool enabled);
412
413 bool isDirty(bool global = false) const;
414 BITPIT_DEPRECATED(bool isExpert() const);
415
416 int getId() const;
417 int getDimension() const;
418 virtual void setDimension(int dimension);
419 bool isThreeDimensional() const;
420
421 virtual int getVolumeCodimension() const = 0;
422 virtual int getSurfaceCodimension() const = 0;
423 virtual int getLineCodimension() const = 0;
424 virtual int getPointCodimension() const = 0;
425
426 bool empty(bool global = true) const;
427
428 bool isVertexAutoIndexingEnabled() const;
429 void setVertexAutoIndexing(bool enabled);
430
431 virtual long getVertexCount() const;
432 long getInternalVertexCount() const;
433#if BITPIT_ENABLE_MPI==1
434 long getGhostVertexCount() const;
435#endif
437 const PiercedVector<Vertex> &getVertices() const;
438 Vertex &getVertex(long id);
439 const Vertex & getVertex(long id) const;
441 const Vertex &getLastInternalVertex() const;
442#if BITPIT_ENABLE_MPI==1
444 const Vertex &getFirstGhostVertex() const;
445#endif
446 const std::array<double, 3> & getVertexCoords(long id) const;
447 void getVertexCoords(std::size_t nVertices, const long *ids, std::unique_ptr<std::array<double, 3>[]> *coordinates) const;
448 void getVertexCoords(std::size_t nVertices, const long *ids, std::array<double, 3> *coordinates) const;
449 VertexIterator addVertex(const Vertex &source, long id = Vertex::NULL_ID);
450 VertexIterator addVertex(Vertex &&source, long id = Vertex::NULL_ID);
451 VertexIterator addVertex(const std::array<double, 3> &coords, long id = Vertex::NULL_ID);
453 long countBorderVertices() const;
454 long countOrphanVertices() const;
455 std::vector<long> findOrphanVertices();
457 std::vector<long> collapseCoincidentVertices();
459
460 VertexIterator getVertexIterator(long id);
461 VertexIterator vertexBegin();
462 VertexIterator vertexEnd();
463 VertexIterator internalVertexBegin();
464 VertexIterator internalVertexEnd();
465#if BITPIT_ENABLE_MPI==1
466 VertexIterator ghostVertexBegin();
467 VertexIterator ghostVertexEnd();
468#endif
469
470 VertexConstIterator getVertexConstIterator(long id) const;
471 VertexConstIterator vertexConstBegin() const;
472 VertexConstIterator vertexConstEnd() const;
473 VertexConstIterator internalVertexConstBegin() const;
474 VertexConstIterator internalVertexConstEnd() const;
475#if BITPIT_ENABLE_MPI==1
476 VertexConstIterator ghostVertexConstBegin() const;
477 VertexConstIterator ghostVertexConstEnd() const;
478#endif
479
480 bool isCellAutoIndexingEnabled() const;
481 void setCellAutoIndexing(bool enabled);
482
483 virtual long getCellCount() const;
484 long getInternalCellCount() const;
486#if BITPIT_ENABLE_MPI==1
487 long getGhostCellCount() const;
488 BITPIT_DEPRECATED(long getGhostCount() const);
489#endif
491 const PiercedVector<Cell> &getCells() const;
492 Cell &getCell(long id);
493 const Cell &getCell(long id) const;
494 virtual ElementType getCellType(long id) const;
496 const Cell &getLastInternalCell() const;
497#if BITPIT_ENABLE_MPI==1
500 const Cell & getFirstGhostCell() const;
501 BITPIT_DEPRECATED(const Cell & getFirstGhost() const);
502#endif
503 CellIterator addCell(const Cell &source, long id = Element::NULL_ID);
504 CellIterator addCell(Cell &&source, long id = Element::NULL_ID);
505 CellIterator addCell(ElementType type, long id = Element::NULL_ID);
506 CellIterator addCell(ElementType type, const std::vector<long> &connectivity, long id = Element::NULL_ID);
507 CellIterator addCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, long id = Element::NULL_ID);
508#if BITPIT_ENABLE_MPI==1
509 CellIterator addCell(const Cell &source, int owner, long id = Element::NULL_ID);
510 CellIterator addCell(const Cell &source, int owner, int haloLayer, long id = Element::NULL_ID);
511 CellIterator addCell(Cell &&source, int owner, long id = Element::NULL_ID);
512 CellIterator addCell(Cell &&source, int owner, int haloLayer, long id = Element::NULL_ID);
513 CellIterator addCell(ElementType type, int owner, long id = Element::NULL_ID);
514 CellIterator addCell(ElementType type, int owner, int haloLayer, long id = Element::NULL_ID);
515 CellIterator addCell(ElementType type, const std::vector<long> &connectivity, int owner, long id = Element::NULL_ID);
516 CellIterator addCell(ElementType type, const std::vector<long> &connectivity, int owner, int haloLayer, long id = Element::NULL_ID);
517 CellIterator addCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, int owner, long id = Element::NULL_ID);
518 CellIterator addCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, int owner, int haloLayer, long id = Element::NULL_ID);
519#endif
520 bool deleteCell(long id);
521 template<typename IdStorage>
522 bool deleteCells(const IdStorage &ids);
523#if BITPIT_ENABLE_MPI==1
524 CellIterator ghostCell2InternalCell(long id);
525 CellIterator internalCell2GhostCell(long id, int owner, int haloLayer);
526#endif
527 virtual double evalCellSize(long id) const = 0;
528 BITPIT_DEPRECATED(long countFreeCells() const);
529 long countBorderCells() const;
530 long countOrphanCells() const;
531 std::vector<long> findOrphanCells() const;
532 long countDuplicateCells() const;
533 std::vector<long> findDuplicateCells() const;
534
535 virtual std::array<double, 3> evalCellCentroid(long id) const;
536 virtual void evalCellBoundingBox(long id, std::array<double,3> *minPoint, std::array<double,3> *maxPoint) const;
537 BITPIT_DEPRECATED(ConstProxyVector<std::array<double BITPIT_COMMA 3>> getCellVertexCoordinates(long id) const);
538 void getCellVertexCoordinates(long id, std::unique_ptr<std::array<double, 3>[]> *coordinates) const;
539 void getCellVertexCoordinates(long id, std::array<double, 3> *coordinates) const;
540 std::vector<long> findCellNeighs(long id) const;
541 void findCellNeighs(long id, std::vector<long> *neighs) const;
542 std::vector<long> findCellNeighs(long id, int codimension, bool complete = true) const;
543 void findCellNeighs(long id, int codimension, bool complete, std::vector<long> *neighs) const;
544 std::vector<long> findCellFaceNeighs(long id) const;
545 void findCellFaceNeighs(long id, std::vector<long> *neighs) const;
546 std::vector<long> findCellFaceNeighs(long id, int face) const;
547 void findCellFaceNeighs(long id, int face, std::vector<long> *neighs) const;
548 std::vector<long> findCellEdgeNeighs(long id, bool complete = true) const;
549 void findCellEdgeNeighs(long id, bool complete, std::vector<long> *neighs) const;
550 std::vector<long> findCellEdgeNeighs(long id, int edge) const;
551 void findCellEdgeNeighs(long id, int edge, std::vector<long> *neighs) const;
552 std::vector<long> findCellVertexNeighs(long id, bool complete = true) const;
553 void findCellVertexNeighs(long id, bool complete, std::vector<long> *neighs) const;
554 std::vector<long> findCellVertexNeighs(long id, int vertex) const;
555 void findCellVertexNeighs(long id, int vertex, std::vector<long> *neighs) const;
556 std::vector<long> findCellVertexOneRing(long id, int vertex) const;
557 void findCellVertexOneRing(long id, int vertex, std::vector<long> *neighs) const;
558 bool findFaceNeighCell(long cellId, long neighId, int *cellFace, int *cellAdjacencyId) const;
559
560 std::set<int> getInternalCellPIDs();
561 std::vector<long> getInternalCellsByPID(int pid);
562
563 std::vector<long> findVertexOneRing(long vertexId) const;
564 void findVertexOneRing(long vertexId, std::vector<long> *ring) const;
565
566 CellIterator getCellIterator(long id);
567 CellIterator cellBegin();
568 CellIterator cellEnd();
569 CellIterator internalCellBegin();
570 BITPIT_DEPRECATED(CellIterator internalBegin());
571 CellIterator internalCellEnd();
572 BITPIT_DEPRECATED(CellIterator internalEnd());
573#if BITPIT_ENABLE_MPI==1
574 CellIterator ghostCellBegin();
575 BITPIT_DEPRECATED(CellIterator ghostBegin());
576 CellIterator ghostCellEnd();
577 BITPIT_DEPRECATED(CellIterator ghostEnd());
578#endif
579
580 CellConstIterator getCellConstIterator(long id) const;
581 CellConstIterator cellConstBegin() const;
582 CellConstIterator cellConstEnd() const;
583 CellConstIterator internalCellConstBegin() const;
584 BITPIT_DEPRECATED(CellConstIterator internalConstBegin() const);
585 CellConstIterator internalCellConstEnd() const;
586 BITPIT_DEPRECATED(CellConstIterator internalConstEnd() const);
587#if BITPIT_ENABLE_MPI==1
588 CellConstIterator ghostCellConstBegin() const;
589 BITPIT_DEPRECATED(CellConstIterator ghostConstBegin() const);
590 CellConstIterator ghostCellConstEnd() const;
591 BITPIT_DEPRECATED(CellConstIterator ghostConstEnd() const);
592#endif
593
595 void setInterfaceAutoIndexing(bool enabled);
596
597 virtual long getInterfaceCount() const;
600 Interface &getInterface(long id);
601 const Interface &getInterface(long id) const;
602 virtual ElementType getInterfaceType(long id) const;
603 InterfaceIterator addInterface(const Interface &source, long id = Element::NULL_ID);
604 InterfaceIterator addInterface(Interface &&source, long id = Element::NULL_ID);
605 InterfaceIterator addInterface(ElementType type, long id = Element::NULL_ID);
606 InterfaceIterator addInterface(ElementType type, const std::vector<long> &connectivity, long id = Element::NULL_ID);
607 InterfaceIterator addInterface(ElementType type, std::unique_ptr<long[]> &&connectStorage, long id = Element::NULL_ID);
608 bool deleteInterface(long id);
609 template<typename IdStorage>
610 bool deleteInterfaces(const IdStorage &ids);
612 long countBorderInterfaces() const;
613 long countOrphanInterfaces() const;
614 std::vector<long> findOrphanInterfaces() const;
616 bool isInterfaceOrphan(long id) const;
617 virtual std::array<double, 3> evalInterfaceCentroid(long id) const;
618 virtual void evalInterfaceBoundingBox(long id, std::array<double,3> *minPoint, std::array<double,3> *maxPoint) const;
619 BITPIT_DEPRECATED(ConstProxyVector<std::array<double BITPIT_COMMA 3>> getInterfaceVertexCoordinates(long id) const);
620 void getInterfaceVertexCoordinates(long id, std::unique_ptr<std::array<double, 3>[]> *coordinates) const;
621 void getInterfaceVertexCoordinates(long id, std::array<double, 3> *coordinates) const;
622 bool intersectInterfacePlane(long id, std::array<double, 3> P, std::array<double, 3> nP, std::array<std::array<double, 3>, 2> *intersection, std::vector<std::array<double, 3>> *polygon) const;
623
624 InterfaceIterator getInterfaceIterator(long id);
625 InterfaceIterator interfaceBegin();
626 InterfaceIterator interfaceEnd();
627
628 InterfaceConstIterator getInterfaceConstIterator(long id) const;
629 InterfaceConstIterator interfaceConstBegin() const;
630 InterfaceConstIterator interfaceConstEnd() const;
631
632 long countFaces() const;
633 long countBorderFaces() const;
634 BITPIT_DEPRECATED(long countFreeFaces() const);
635
636 bool sort();
637 bool sortVertices();
638 bool sortCells();
639 bool sortInterfaces();
640
641 bool squeeze();
642 bool squeezeVertices();
643 bool squeezeCells();
644 bool squeezeInterfaces();
645
646 long locatePoint(double x, double y, double z) const;
647 virtual long locatePoint(const std::array<double, 3> &point) const = 0;
648
650 bool areAdjacenciesDirty(bool global = false) const;
652 void initializeAdjacencies(AdjacenciesBuildStrategy strategy = ADJACENCIES_AUTOMATIC);
653 void updateAdjacencies(bool forcedUpdated = false);
654 void destroyAdjacencies();
655
657 bool areInterfacesDirty(bool global = false) const;
659 void initializeInterfaces(InterfacesBuildStrategy strategy = INTERFACES_AUTOMATIC);
660 void updateInterfaces(bool forcedUpdated = false);
661 void destroyInterfaces();
662
663 void getBoundingBox(std::array<double, 3> &minPoint, std::array<double, 3> &maxPoint) const;
664 void getBoundingBox(bool global, std::array<double, 3> &minPoint, std::array<double, 3> &maxPoint) const;
665 bool isBoundingBoxDirty(bool global = false) const;
666 void updateBoundingBox(bool forcedUpdated = false);
667
668 virtual void translate(const std::array<double, 3> &translation);
669 void translate(double sx, double sy, double sz);
670 virtual void rotate(const std::array<double, 3> &n0, const std::array<double, 3> &n1, double angle);
671 void rotate(double n0x, double n0y, double n0z, double n1x, double n1y, double n1z, double angle);
672 void scale(const std::array<double, 3> &scaling);
673 virtual void scale(const std::array<double, 3> &scaling, const std::array<double, 3> &origin);
674 void scale(double scaling);
675 void scale(double scaling, const std::array<double, 3> &origin);
676 void scale(double sx, double sy, double sz);
677 void scale(double sx, double sy, double sz, const std::array<double, 3> &origin);
678
679 void setTol(double tolerance);
680 double getTol() const;
681 void resetTol();
682 bool isTolCustomized() const;
683
684 void displayTopologyStats(std::ostream &out, unsigned int padding = 0) const;
685 void displayVertices(std::ostream &out, unsigned int padding = 0) const;
686 void displayCells(std::ostream &out, unsigned int padding = 0) const;
687 void displayInterfaces(std::ostream &out, unsigned int padding = 0) const;
688
690 const VTKUnstructuredGrid & getVTK() const;
691 WriteTarget getVTKWriteTarget() const;
692 void setVTKWriteTarget(WriteTarget targetCells);
693 const CellConstRange getVTKCellWriteRange() const;
694 void write(VTKWriteMode mode = VTKWriteMode::DEFAULT);
695 void write(VTKWriteMode mode, double time);
696 void write(const std::string &name, VTKWriteMode mode = VTKWriteMode::DEFAULT);
697 void write(const std::string &name, VTKWriteMode mode, double time);
698
699 void flushData(std::fstream &stream, const std::string &name, VTKFormat format) override;
700
701 int getDumpVersion() const;
702 bool dump(std::ostream &stream);
703 bool dump(std::ostream &stream) const;
704 void restore(std::istream &stream, bool reregister = false);
705
706 void consecutiveRenumberVertices(long offset = 0);
707 void consecutiveRenumberCells(long offset = 0);
708 void consecutiveRenumberInterfaces(long offset = 0);
709 void consecutiveRenumber(long offsetVertices, long offsetCells, long offsetInterfaces);
710
711#if BITPIT_ENABLE_MPI==1
712 const MPI_Comm & getCommunicator() const;
713 int getRank() const;
714 int getProcessorCount() const;
715
716 bool isDistributed(bool allowDirty = false) const;
717 int getOwner(bool allowDirty = false) const;
718
719 void setHaloSize(std::size_t haloSize);
720 std::size_t getHaloSize() const;
721
722 int getCellOwner(long id) const;
723 BITPIT_DEPRECATED(int getCellRank(long id) const);
724 int getCellHaloLayer(long id) const;
725
726 int getVertexOwner(long id) const;
727 BITPIT_DEPRECATED(int getVertexRank(long id) const);
728
729 bool isRankNeighbour(int rank);
730 std::vector<int> getNeighbourRanks();
731
732 const std::unordered_map<int, std::vector<long>> & getGhostVertexExchangeTargets() const;
733 const std::vector<long> & getGhostVertexExchangeTargets(int rank) const;
734 const std::unordered_map<int, std::vector<long>> & getGhostVertexExchangeSources() const;
735 const std::vector<long> & getGhostVertexExchangeSources(int rank) const;
736
737 const std::unordered_map<int, std::vector<long>> & getGhostCellExchangeTargets() const;
738 BITPIT_DEPRECATED(const std::unordered_map<int BITPIT_COMMA std::vector<long>> & getGhostExchangeTargets() const);
739 const std::vector<long> & getGhostCellExchangeTargets(int rank) const;
740 BITPIT_DEPRECATED(const std::vector<long> & getGhostExchangeTargets(int rank) const);
741 const std::unordered_map<int, std::vector<long>> & getGhostCellExchangeSources() const;
742 BITPIT_DEPRECATED(const std::unordered_map<int BITPIT_COMMA std::vector<long>> & getGhostExchangeSources() const);
743 const std::vector<long> & getGhostCellExchangeSources(int rank) const;
744 BITPIT_DEPRECATED(const std::vector<long> & getGhostExchangeSources(int rank) const);
745
746 bool isPartitioned() const;
747 bool isPartitioningSupported() const;
748 bool arePartitioningInfoDirty(bool global = true) const;
750 PartitioningStatus getPartitioningStatus(bool global = false) const;
751 double evalPartitioningUnbalance() const;
752 double evalPartitioningUnbalance(const std::unordered_map<long, double> &cellWeights) const;
753 BITPIT_DEPRECATED(std::vector<adaption::Info> partition(MPI_Comm communicator, const std::unordered_map<long, int> &cellRanks, bool trackPartitioning, bool squeezeStorage = false, std::size_t haloSize = 1));
754 std::vector<adaption::Info> partition(const std::unordered_map<long, int> &cellRanks, bool trackPartitioning, bool squeezeStorage = false);
755 BITPIT_DEPRECATED(std::vector<adaption::Info> partition(MPI_Comm communicator, const std::unordered_map<long, double> &cellWeights, bool trackPartitioning, bool squeezeStorage = false, std::size_t haloSize = 1));
756 std::vector<adaption::Info> partition(const std::unordered_map<long, double> &cellWeights, bool trackPartitioning, bool squeezeStorage = false);
757 BITPIT_DEPRECATED(std::vector<adaption::Info> partition(MPI_Comm communicator, bool trackPartitioning, bool squeezeStorage = false, std::size_t haloSize = 1));
758 std::vector<adaption::Info> partition(bool trackPartitioning, bool squeezeStorage = false);
759 BITPIT_DEPRECATED(std::vector<adaption::Info> partitioningPrepare(MPI_Comm communicator, const std::unordered_map<long, int> &cellRanks, bool trackPartitioning, std::size_t haloSize = 1));
760 std::vector<adaption::Info> partitioningPrepare(const std::unordered_map<long, int> &cellRanks, bool trackPartitioning);
761 BITPIT_DEPRECATED(std::vector<adaption::Info> partitioningPrepare(MPI_Comm communicator, const std::unordered_map<long, double> &cellWeights, bool trackPartitioning, std::size_t haloSize = 1));
762 std::vector<adaption::Info> partitioningPrepare(const std::unordered_map<long, double> &cellWeights, bool trackPartitioning);
763 BITPIT_DEPRECATED(std::vector<adaption::Info> partitioningPrepare(MPI_Comm communicator, bool trackPartitioning, std::size_t haloSize = 1));
764 std::vector<adaption::Info> partitioningPrepare(bool trackPartitioning);
765 std::vector<adaption::Info> partitioningAlter(bool trackPartitioning = true, bool squeezeStorage = false);
766 void partitioningCleanup();
767#endif
768
769 template<typename Function>
770 void processCellNeighbours(long seedId, int nLayers, Function function) const;
771 template<typename Selector, typename Function>
772 void processCellNeighbours(long seedId, int nLayers, Selector isSelected, Function function) const;
773 template<typename Function, typename SeedContainer>
774 void processCellsNeighbours(const SeedContainer &seedIds, int nLayers, Function function) const;
775 template<typename Selector, typename Function, typename SeedContainer>
776 void processCellsNeighbours(const SeedContainer &seedIds, int nLayers, Selector isSelected, Function function) const;
777 template<typename Function>
778 void processCellFaceNeighbours(long seedId, int nLayers, Function function) const;
779 template<typename Selector, typename Function>
780 void processCellFaceNeighbours(long seedId, int nLayers, Selector isSelected, Function function) const;
781 template<typename Function, typename SeedContainer>
782 void processCellsFaceNeighbours(const SeedContainer &seedIds, int nLayers, Function function) const;
783 template<typename Selector, typename Function, typename SeedContainer>
784 void processCellsFaceNeighbours(const SeedContainer &seedIds, int nLayers, Selector isSelected, Function function) const;
785
786 std::array<double, 3> evalElementCentroid(const Element &element) const;
787 void evalElementBoundingBox(const Element &element, std::array<double,3> *minPoint, std::array<double,3> *maxPoint) const;
788 BITPIT_DEPRECATED(ConstProxyVector<std::array<double BITPIT_COMMA 3>> getElementVertexCoordinates(const Element &element) const);
789 void getElementVertexCoordinates(const Element &element, std::unique_ptr<std::array<double, 3>[]> *coordinates) const;
790 void getElementVertexCoordinates(const Element &element, std::array<double, 3> *coordinates) const;
791
792protected:
793 typedef uint16_t AlterationFlags;
794 typedef std::unordered_map<long, AlterationFlags> AlterationFlagsStorage;
795
796#if BITPIT_ENABLE_MPI==1
797 const static double DEFAULT_PARTITIONING_WEIGTH;
798#endif
799
800 const static AlterationFlags FLAG_NONE = 0x0;
801 const static AlterationFlags FLAG_DELETED = (1u << 0);
802 const static AlterationFlags FLAG_ADJACENCIES_DIRTY = (1u << 1);
803 const static AlterationFlags FLAG_INTERFACES_DIRTY = (1u << 2);
804 const static AlterationFlags FLAG_DANGLING = (1u << 3);
805
806 PiercedVector<Vertex> m_vertices;
807 PiercedVector<Cell> m_cells;
808 PiercedVector<Interface> m_interfaces;
809
810 AlterationFlagsStorage m_alteredCells;
811 AlterationFlagsStorage m_alteredInterfaces;
812
813#if BITPIT_ENABLE_MPI==1
814 PatchKernel(MPI_Comm communicator, std::size_t haloSize, AdaptionMode adaptionMode, PartitioningMode partitioningMode);
815 PatchKernel(int dimension, MPI_Comm communicator, std::size_t haloSize, AdaptionMode adaptionMode, PartitioningMode partitioningMode);
816 PatchKernel(int id, int dimension, MPI_Comm communicator, std::size_t haloSize, AdaptionMode adaptionMode, PartitioningMode partitioningMode);
817#else
818 PatchKernel(AdaptionMode adaptionMode);
819 PatchKernel(int dimension, AdaptionMode adaptionMode);
820 PatchKernel(int id, int dimension, AdaptionMode adaptionMode);
821#endif
822 PatchKernel(const PatchKernel &other);
823 PatchKernel & operator=(const PatchKernel &other) = delete;
824
825 void clearBoundingBox();
826 bool isBoundingBoxFrozen() const;
827 void setBoundingBoxFrozen(bool frozen);
828 void setBoundingBoxDirty(bool dirty);
829 void setBoundingBox(const std::array<double, 3> &minPoint, const std::array<double, 3> &maxPoint);
830
831#if BITPIT_ENABLE_MPI==1
832 bool isCommunicatorSet() const;
833 virtual void setCommunicator(MPI_Comm communicator);
834#endif
835
836#if BITPIT_ENABLE_MPI==1
837 CellIterator restoreCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, int owner, int haloLayer, long id);
838#else
839 CellIterator restoreCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, long id);
840#endif
841
842 InterfaceIterator restoreInterface(ElementType type, std::unique_ptr<long[]> &&connectStorage, long id);
843
844#if BITPIT_ENABLE_MPI==1
845 VertexIterator restoreVertex(const std::array<double, 3> &coords, int owner, long id);
846#else
847 VertexIterator restoreVertex(const std::array<double, 3> &coords, long id);
848#endif
849
850 bool deleteVertex(long id);
851 template<typename IdStorage>
852 bool deleteVertices(const IdStorage &ids);
853#if BITPIT_ENABLE_MPI==1
854 VertexIterator ghostVertex2InternalVertex(long id);
855 VertexIterator internalVertex2GhostVertex(long id, int owner);
856#endif
857
858 void dumpVertices(std::ostream &stream) const;
859 void restoreVertices(std::istream &stream);
860
861 void dumpCells(std::ostream &stream) const;
862 void restoreCells(std::istream &stream);
863
864 void dumpInterfaces(std::ostream &stream) const;
865 void restoreInterfaces(std::istream &stream);
866
868#if BITPIT_ENABLE_MPI==1
870#endif
871
873#if BITPIT_ENABLE_MPI==1
875#endif
876
877 std::unordered_map<long, std::vector<long>> binGroupVertices(const PiercedVector<Vertex> &vertices, int nBins);
878 std::unordered_map<long, std::vector<long>> binGroupVertices(int nBins);
879
881 void resetAdjacencies();
883 virtual void _resetAdjacencies(bool release);
884 virtual void _updateAdjacencies();
885
888 virtual void _resetInterfaces(bool release);
889 virtual void _updateInterfaces();
890
891 bool testCellAlterationFlags(long id, AlterationFlags flags) const;
892 AlterationFlags getCellAlterationFlags(long id) const;
893 void resetCellAlterationFlags(long id, AlterationFlags flags = FLAG_NONE);
894 void setCellAlterationFlags(AlterationFlags flags);
895 void setCellAlterationFlags(long id, AlterationFlags flags);
896 void unsetCellAlterationFlags(AlterationFlags flags);
897 void unsetCellAlterationFlags(long id, AlterationFlags flags);
898
899 bool testInterfaceAlterationFlags(long id, AlterationFlags flags) const;
900 AlterationFlags getInterfaceAlterationFlags(long id) const;
901 void resetInterfaceAlterationFlags(long id, AlterationFlags flags = FLAG_NONE);
902 void setInterfaceAlterationFlags(AlterationFlags flags);
903 void setInterfaceAlterationFlags(long id, AlterationFlags flags);
904 void unsetInterfaceAlterationFlags(AlterationFlags flags);
905 void unsetInterfaceAlterationFlags(long id, AlterationFlags flags);
906
907 bool testAlterationFlags(AlterationFlags availableFlags, AlterationFlags requestedFlags) const;
908
909 void setAdaptionMode(AdaptionMode mode);
911 virtual std::vector<adaption::Info> _adaptionPrepare(bool trackAdaption);
912 virtual std::vector<adaption::Info> _adaptionAlter(bool trackAdaption);
913 virtual void _adaptionCleanup();
914 virtual bool _markCellForRefinement(long id);
915 virtual bool _markCellForCoarsening(long id);
916 virtual bool _resetCellAdaptionMarker(long id);
917 virtual adaption::Marker _getCellAdaptionMarker(long id);
918 virtual bool _enableCellBalancing(long id, bool enabled);
919
920 virtual void _setTol(double tolerance);
921 virtual void _resetTol();
922
923 virtual int _getDumpVersion() const = 0;
924 virtual void _dump(std::ostream &stream) const = 0;
925 virtual void _restore(std::istream &stream) = 0;
926
927 virtual long _getCellNativeIndex(long id) const;
928
929 virtual void _findCellNeighs(long id, const std::vector<long> *blackList, std::vector<long> *neighs) const;
930 virtual void _findCellFaceNeighs(long id, int face, const std::vector<long> *blackList, std::vector<long> *neighs) const;
931 virtual void _findCellEdgeNeighs(long id, int edge, const std::vector<long> *blackList, std::vector<long> *neighs) const;
932 virtual void _findCellVertexNeighs(long id, int vertex, const std::vector<long> *blackList, std::vector<long> *neighs) const;
933
934 virtual void _writePrepare();
935 virtual void _writeFinalize();
936
937 BITPIT_DEPRECATED(void setExpert(bool expert));
938
939 void extractEnvelope(PatchKernel &envelope) const;
940
941 void addPointToBoundingBox(const std::array<double, 3> &point);
942 void removePointFromBoundingBox(const std::array<double, 3> &point);
943#if BITPIT_ENABLE_MPI==1
944 virtual std::size_t _getMaxHaloSize();
945 virtual void _setHaloSize(std::size_t haloSize);
946
947 void setPartitioned(bool partitioned);
950 virtual std::vector<adaption::Info> _partitioningPrepare(const std::unordered_map<long, double> &cellWeights, double defaultWeight, bool trackPartitioning);
951 virtual std::vector<adaption::Info> _partitioningPrepare(const std::unordered_map<long, int> &cellRanks, bool trackPartitioning);
952 virtual std::vector<adaption::Info> _partitioningAlter(bool trackPartitioning);
953 virtual void _partitioningCleanup();
954
955 virtual std::vector<long> _findGhostCellExchangeSources(int rank);
956#endif
957
958 template<typename item_t, typename id_t = long>
959 std::unordered_map<id_t, id_t> consecutiveItemRenumbering(PiercedVector<item_t, id_t> &container, long offset);
960
961 template<typename item_t, typename id_t = long>
962 void mappedItemRenumbering(PiercedVector<item_t, id_t> &container, const std::unordered_map<id_t, id_t> &renumberMap);
963
964 virtual int findAdjoinNeighFace(const Cell &cell, int cellFace, const Cell &neigh) const;
965 virtual bool isSameFace(const Cell &cell_A, int face_A, const Cell &cell_B, int face_B) const;
966
967private:
968 struct GhostVertexInfo {
969 int owner;
970 };
971
972 struct GhostCellInfo {
973 int owner;
974 int haloLayer;
975 };
976
977 std::unique_ptr<IndexGenerator<long>> m_vertexIdGenerator;
978 std::unique_ptr<IndexGenerator<long>> m_interfaceIdGenerator;
979 std::unique_ptr<IndexGenerator<long>> m_cellIdGenerator;
980
981 long m_nInternalVertices;
982#if BITPIT_ENABLE_MPI==1
983 long m_nGhostVertices;
984#endif
985
986 long m_lastInternalVertexId;
987#if BITPIT_ENABLE_MPI==1
988 long m_firstGhostVertexId;
989#endif
990
991 long m_nInternalCells;
992#if BITPIT_ENABLE_MPI==1
993 long m_nGhostCells;
994#endif
995
996 long m_lastInternalCellId;
997#if BITPIT_ENABLE_MPI==1
998 long m_firstGhostCellId;
999#endif
1000
1001 VTKUnstructuredGrid m_vtk ;
1002 WriteTarget m_vtkWriteTarget;
1003 PiercedStorage<long, long> m_vtkVertexMap;
1004
1005 bool m_boxFrozen;
1006 bool m_boxDirty;
1007 std::array<double, 3> m_boxMinPoint;
1008 std::array<double, 3> m_boxMaxPoint;
1009 std::array<int, 3> m_boxMinCounter;
1010 std::array<int, 3> m_boxMaxCounter;
1011
1012 AdjacenciesBuildStrategy m_adjacenciesBuildStrategy;
1013
1014 InterfacesBuildStrategy m_interfacesBuildStrategy;
1015
1016 AdaptionMode m_adaptionMode;
1017 AdaptionStatus m_adaptionStatus;
1018
1019 int m_id;
1020 int m_dimension;
1021
1022 bool m_toleranceCustom;
1023 double m_tolerance;
1024
1025 int m_rank;
1026 int m_nProcessors;
1027#if BITPIT_ENABLE_MPI==1
1028 MPI_Comm m_communicator;
1029 PartitioningMode m_partitioningMode;
1030 PartitioningStatus m_partitioningStatus;
1031
1032 int m_owner;
1033
1034 std::size_t m_haloSize;
1035
1036 int m_partitioningCellsTag;
1037 int m_partitioningVerticesTag;
1038 bool m_partitioningSerialization;
1039 std::unordered_map<long, int> m_partitioningOutgoings;
1040 std::vector<std::pair<int, int>> m_partitioningGlobalExchanges;
1041
1042 bool m_partitioningInfoDirty;
1043
1044 std::unordered_map<long, GhostVertexInfo> m_ghostVertexInfo;
1045 std::unordered_map<int, std::vector<long>> m_ghostVertexExchangeTargets;
1046 std::unordered_map<int, std::vector<long>> m_ghostVertexExchangeSources;
1047
1048 std::unordered_map<long, GhostCellInfo> m_ghostCellInfo;
1049 std::unordered_map<int, std::vector<long>> m_ghostCellExchangeTargets;
1050 std::unordered_map<int, std::vector<long>> m_ghostCellExchangeSources;
1051
1052 void setGhostVertexInfo(long id, int owner);
1053 void unsetGhostVertexInfo(long id);
1054 void clearGhostVerticesInfo();
1055
1056 void setGhostCellInfo(long id, int owner, int haloLayer);
1057 void unsetGhostCellInfo(long id);
1058 void clearGhostCellsInfo();
1059
1060 void computeCellHaloLayer(int id);
1061
1062 void _partitioningAlter_deleteGhosts();
1063
1064 std::unordered_map<long, int> _partitioningAlter_evalGhostCellOwnershipChanges();
1065 void _partitioningAlter_applyGhostCellOwnershipChanges(int sendRank, std::unordered_map<long, int> *ghostCellOwnershipChanges);
1066
1067 std::vector<adaption::Info> _partitioningAlter_sendCells(const std::unordered_set<int> &recvRanks, bool trackPartitioning, std::unordered_map<long, int> *ghostCellOwnershipChanges);
1068 std::vector<adaption::Info> _partitioningAlter_receiveCells(const std::unordered_set<int> &sendRanks, bool trackPartitioning, std::unordered_map<long, int> *ghostCellOwnershipChanges);
1069
1070 void setPartitioningInfoDirty(bool dirty);
1071
1072 void updatePartitioningInfo(bool forcedUpdated = false);
1073
1074 void updateGhostCellExchangeInfo();
1075
1076 void updateGhostVertexOwners();
1077 void updateGhostVertexExchangeInfo();
1078
1079 void updateOwner();
1080 int evalOwner() const;
1081
1082 std::unordered_map<long, int> evaluateExchangeVertexOwners() const;
1083
1084 template<typename ExcludeList>
1085 bool confirmCellHaloLayer(const Cell &cell, int haloLayer, const ExcludeList &excludeList = ExcludeList()) const;
1086#endif
1087
1088#if BITPIT_ENABLE_MPI==1
1089 void initialize(MPI_Comm communicator, std::size_t haloSize);
1090 void initializeHaloSize(std::size_t haloSize);
1091 void initializeCommunicator(MPI_Comm communicator);
1092
1093 void freeCommunicator();
1094#else
1095 void initialize();
1096 void initializeSerialCommunicator();
1097#endif
1098
1099 void finalizeAlterations(bool squeezeStorage = false);
1100
1101 InterfaceIterator buildCellInterface(Cell *cell_1, int face_1, Cell *cell_2, int face_2, long interfaceId = Element::NULL_ID);
1102
1103 void _setId(int id);
1104
1105 bool testElementAlterationFlags(long id, AlterationFlags flags, const AlterationFlagsStorage &flagsStorage) const;
1106 AlterationFlags getElementAlterationFlags(long id, const AlterationFlagsStorage &flagsStorage) const;
1107 void resetElementAlterationFlags(long id, AlterationFlags flags, AlterationFlagsStorage *flagsStorage) const;
1108 void setElementAlterationFlags(long id, AlterationFlags flags, AlterationFlagsStorage *flagsStorage) const;
1109 void unsetElementAlterationFlags(AlterationFlags flags, AlterationFlagsStorage *flagsStorage) const;
1110 void unsetElementAlterationFlags(long id, AlterationFlags flags, AlterationFlagsStorage *flagsStorage) const;
1111
1112 void mergeAdaptionInfo(std::vector<adaption::Info> &&source, std::vector<adaption::Info> &destination);
1113
1114 void setRestoredCellAlterationFlags(long id);
1115 void setAddedCellAlterationFlags(long id);
1116 void setDeletedCellAlterationFlags(long id);
1117
1118 void setRestoredInterfaceAlterationFlags(long id);
1119 void setAddedInterfaceAlterationFlags(long id);
1120 void setDeletedInterfaceAlterationFlags(long id);
1121
1122 void dumpVertexAutoIndexing(std::ostream &stream) const;
1123 void restoreVertexAutoIndexing(std::istream &stream);
1124 void createVertexIndexGenerator(bool populate);
1125 void importVertexIndexGenerator(const PatchKernel &source);
1126
1127 void dumpCellAutoIndexing(std::ostream &stream) const;
1128 void restoreCellAutoIndexing(std::istream &stream);
1129 void createCellIndexGenerator(bool populate);
1130 void importCellIndexGenerator(const PatchKernel &source);
1131
1132 void dumpInterfaceAutoIndexing(std::ostream &stream) const;
1133 void restoreInterfaceAutoIndexing(std::istream &stream);
1134 void createInterfaceIndexGenerator(bool populate);
1135 void importInterfaceIndexGenerator(const PatchKernel &source);
1136
1137 VertexIterator _addInternalVertex(const std::array<double, 3> &coords, long id);
1138
1139 void _restoreInternalVertex(const VertexIterator &iterator, const std::array<double, 3> &coords);
1140#if BITPIT_ENABLE_MPI==1
1141 void _restoreGhostVertex(const VertexIterator &iterator, const std::array<double, 3> &coords, int owner);
1142#endif
1143
1144 void _deleteInternalVertex(long id);
1145#if BITPIT_ENABLE_MPI==1
1146 void _deleteGhostVertex(long id);
1147#endif
1148
1149 CellIterator _addInternalCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, long id);
1150#if BITPIT_ENABLE_MPI==1
1151 CellIterator _addGhostCell(ElementType type, std::unique_ptr<long[]> &&connectStorage, int owner, int haloLayer, long id);
1152#endif
1153
1154 void _restoreInternalCell(const CellIterator &iterator, ElementType type, std::unique_ptr<long[]> &&connectStorage);
1155#if BITPIT_ENABLE_MPI==1
1156 void _restoreGhostCell(const CellIterator &iterator, ElementType type, std::unique_ptr<long[]> &&connectStorage, int owner, int haloLayer);
1157#endif
1158
1159 void _deleteInternalCell(long id);
1160#if BITPIT_ENABLE_MPI==1
1161 void _deleteGhostCell(long id);
1162#endif
1163
1164 InterfaceIterator _addInterface(ElementType type, std::unique_ptr<long[]> &&connectStorage, long id);
1165
1166 void _restoreInterface(const InterfaceIterator &iterator, ElementType type, std::unique_ptr<long[]> &&connectStorage);
1167
1168 void _deleteInterface(long id);
1169
1170 void replaceVTKStreamer(const VTKBaseStreamer *original, VTKBaseStreamer *updated);
1171
1172};
1173
1174}
1175
1176// Template implementation
1177#include "patch_kernel.tpp"
1178#include "patch_kernel_parallel.tpp"
1179
1180#endif
The Cell class defines the cells.
Definition cell.hpp:42
The Element class provides an interface for defining elements.
Definition element.hpp:46
The Interface class defines the interfaces among cells.
Definition interface.hpp:37
The PatchKernel class provides an interface for defining patches.
void initializeAdjacencies(AdjacenciesBuildStrategy strategy=ADJACENCIES_AUTOMATIC)
VertexIterator vertexEnd()
void unsetCellAlterationFlags(AlterationFlags flags)
CellConstIterator cellConstBegin() const
CellIterator cellBegin()
void processCellNeighbours(long seedId, int nLayers, Function function) const
Vertex & getLastInternalVertex()
CellIterator getCellIterator(long id)
void mappedItemRenumbering(PiercedVector< item_t, id_t > &container, const std::unordered_map< id_t, id_t > &renumberMap)
void dumpInterfaces(std::ostream &stream) const
CellIterator internalCellBegin()
InterfacesBuildStrategy getInterfacesBuildStrategy() const
CellConstIterator internalConstBegin() const
virtual bool _enableCellBalancing(long id, bool enabled)
CellConstIterator internalCellConstEnd() const
CellIterator internalBegin()
AdjacenciesBuildStrategy getAdjacenciesBuildStrategy() const
virtual bool _markCellForRefinement(long id)
std::vector< long > getInternalCellsByPID(int pid)
void displayVertices(std::ostream &out, unsigned int padding=0) const
int getVertexOwner(long id) const
InterfaceIterator restoreInterface(ElementType type, std::unique_ptr< long[]> &&connectStorage, long id)
const std::unordered_map< int, std::vector< long > > & getGhostCellExchangeSources() const
std::unordered_map< id_t, id_t > consecutiveItemRenumbering(PiercedVector< item_t, id_t > &container, long offset)
virtual void _findCellEdgeNeighs(long id, int edge, const std::vector< long > *blackList, std::vector< long > *neighs) const
void markCellForRefinement(long id)
std::vector< long > findVertexOneRing(long vertexId) const
virtual std::array< double, 3 > evalInterfaceCentroid(long id) const
void setCellAutoIndexing(bool enabled)
virtual std::vector< long > _findGhostCellExchangeSources(int rank)
virtual long getInterfaceCount() const
std::vector< adaption::Info > update(bool trackAdaption=true, bool squeezeStorage=false)
bool intersectInterfacePlane(long id, std::array< double, 3 > P, std::array< double, 3 > nP, std::array< std::array< double, 3 >, 2 > *intersection, std::vector< std::array< double, 3 > > *polygon) const
InterfaceIterator interfaceBegin()
PiercedVector< Cell > & getCells()
std::vector< long > findOrphanInterfaces() const
@ ADAPTION_AUTOMATIC
No adaption can be performed.
bool testAlterationFlags(AlterationFlags availableFlags, AlterationFlags requestedFlags) const
VertexConstIterator vertexConstEnd() const
void dumpVertices(std::ostream &stream) const
virtual void _writeFinalize()
InterfaceIterator interfaceEnd()
Interface & getInterface(long id)
PatchKernel(PatchKernel &&other)
virtual void _resetInterfaces(bool release)
VertexIterator restoreVertex(const std::array< double, 3 > &coords, int owner, long id)
AdaptionMode getAdaptionMode() const
CellIterator restoreCell(ElementType type, std::unique_ptr< long[]> &&connectStorage, int owner, int haloLayer, long id)
void setExpert(bool expert)
long countBorderInterfaces() const
InterfaceIterator addInterface(const Interface &source, long id=Element::NULL_ID)
void dumpCells(std::ostream &stream) const
void setBoundingBoxDirty(bool dirty)
CellConstIterator internalCellConstBegin() const
bool isVertexAutoIndexingEnabled() const
virtual void _writePrepare()
VertexIterator vertexBegin()
InterfaceConstIterator interfaceConstEnd() const
InterfaceConstIterator interfaceConstBegin() const
std::vector< long > collapseCoincidentVertices()
void consecutiveRenumberInterfaces(long offset=0)
bool isCellAutoIndexingEnabled() const
VertexConstIterator vertexConstBegin() const
virtual std::vector< adaption::Info > _partitioningPrepare(const std::unordered_map< long, double > &cellWeights, double defaultWeight, bool trackPartitioning)
WriteTarget getVTKWriteTarget() const
virtual long getVertexCount() const
const std::unordered_map< int, std::vector< long > > & getGhostCellExchangeTargets() const
long countOrphanInterfaces() const
virtual void evalCellBoundingBox(long id, std::array< double, 3 > *minPoint, std::array< double, 3 > *maxPoint) const
void restoreCells(std::istream &stream)
virtual void _findCellFaceNeighs(long id, int face, const std::vector< long > *blackList, std::vector< long > *neighs) const
void addPointToBoundingBox(const std::array< double, 3 > &point)
std::vector< long > findDuplicateCells() const
virtual ElementType getInterfaceType(long id) const
VertexIterator addVertex(const Vertex &source, long id=Vertex::NULL_ID)
bool empty(bool global=true) const
std::vector< long > findCellNeighs(long id) const
bool deleteVertex(long id)
void getBoundingBox(std::array< double, 3 > &minPoint, std::array< double, 3 > &maxPoint) const
PartitioningMode getPartitioningMode() const
std::vector< adaption::Info > adaptionPrepare(bool trackAdaption=true)
adaption::Marker getCellAdaptionMarker(long id)
void processCellsNeighbours(const SeedContainer &seedIds, int nLayers, Function function) const
VTKUnstructuredGrid & getVTK()
std::unordered_map< long, std::vector< long > > binGroupVertices(const PiercedVector< Vertex > &vertices, int nBins)
void setBoundingBoxFrozen(bool frozen)
VertexIterator internalVertex2GhostVertex(long id, int owner)
virtual bool isSameFace(const Cell &cell_A, int face_A, const Cell &cell_B, int face_B) const
virtual std::vector< adaption::Info > _partitioningAlter(bool trackPartitioning)
Cell & getCell(long id)
const CellConstRange getVTKCellWriteRange() const
virtual bool _resetCellAdaptionMarker(long id)
void consecutiveRenumberCells(long offset=0)
virtual void simulateCellUpdate(const long id, adaption::Marker marker, std::vector< Cell > *virtualCells, PiercedVector< Vertex, long > *virtualVertices) const
void setVertexAutoIndexing(bool enabled)
VertexIterator getVertexIterator(long id)
VertexConstIterator getVertexConstIterator(long id) const
PiercedVector< Vertex > & getVertices()
virtual bool _markCellForCoarsening(long id)
std::vector< long > findCellFaceNeighs(long id) const
virtual void _setTol(double tolerance)
bool deleteVertices(const IdStorage &ids)
long locatePoint(double x, double y, double z) const
InterfaceConstIterator getInterfaceConstIterator(long id) const
void flushData(std::fstream &stream, const std::string &name, VTKFormat format) override
virtual adaption::Marker _getCellAdaptionMarker(long id)
int getVertexRank(long id) const
std::vector< adaption::Info > partition(MPI_Comm communicator, const std::unordered_map< long, int > &cellRanks, bool trackPartitioning, bool squeezeStorage=false, std::size_t haloSize=1)
CellIterator addCell(const Cell &source, long id=Element::NULL_ID)
void evalElementBoundingBox(const Element &element, std::array< double, 3 > *minPoint, std::array< double, 3 > *maxPoint) const
PartitioningStatus getPartitioningStatus(bool global=false) const
void resetCellAlterationFlags(long id, AlterationFlags flags=FLAG_NONE)
void updateInterfaces(bool forcedUpdated=false)
void displayTopologyStats(std::ostream &out, unsigned int padding=0) const
virtual void rotate(const std::array< double, 3 > &n0, const std::array< double, 3 > &n1, double angle)
CellIterator cellEnd()
ConstProxyVector< std::array< double BITPIT_COMMA 3 > > getElementVertexCoordinates(const Element &element) const
long countBorderCells() const
void setTol(double tolerance)
CellConstIterator cellConstEnd() const
AlterationFlags getCellAlterationFlags(long id) const
void scale(const std::array< double, 3 > &scaling)
virtual void setDimension(int dimension)
void markCellForCoarsening(long id)
const std::unordered_map< int BITPIT_COMMA std::vector< long > > & getGhostExchangeTargets() const
void restoreVertices(std::istream &stream)
virtual void translate(const std::array< double, 3 > &translation)
VertexConstIterator internalVertexConstEnd() const
bool arePartitioningInfoDirty(bool global=true) const
CellIterator internalCell2GhostCell(long id, int owner, int haloLayer)
@ PARTITIONING_ENABLED
No partitioning can be performed.
void setInterfaceAutoIndexing(bool enabled)
CellIterator ghostCell2InternalCell(long id)
void setAdaptionStatus(AdaptionStatus status)
bool isInterfaceAutoIndexingEnabled() const
PatchKernel & operator=(PatchKernel &&other)
void setAdjacenciesBuildStrategy(AdjacenciesBuildStrategy status)
VertexConstIterator ghostVertexConstEnd() const
VertexIterator internalVertexEnd()
bool reserveInterfaces(size_t nInterfaces)
int getOwner(bool allowDirty=false) const
long countDuplicateCells() const
bool areInterfacesDirty(bool global=false) const
double evalPartitioningUnbalance() const
CellConstIterator ghostCellConstEnd() const
static const double DEFAULT_PARTITIONING_WEIGTH
void setInterfacesBuildStrategy(InterfacesBuildStrategy status)
VertexIterator ghostVertexBegin()
const MPI_Comm & getCommunicator() const
void updateBoundingBox(bool forcedUpdated=false)
void processCellFaceNeighbours(long seedId, int nLayers, Function function) const
std::size_t getHaloSize() const
void restore(std::istream &stream, bool reregister=false)
virtual void resetVertices()
virtual int findAdjoinNeighFace(const Cell &cell, int cellFace, const Cell &neigh) const
virtual void resetInterfaces()
CellConstIterator ghostConstEnd() const
long countFreeVertices() const
virtual void reset()
void setHaloSize(std::size_t haloSize)
long countFreeCells() const
virtual void _updateAdjacencies()
void initializeInterfaces(InterfacesBuildStrategy strategy=INTERFACES_AUTOMATIC)
CellIterator internalEnd()
void enableCellBalancing(long id, bool enabled)
virtual void _resetAdjacencies(bool release)
void setVTKWriteTarget(WriteTarget targetCells)
void displayCells(std::ostream &out, unsigned int padding=0) const
bool testInterfaceAlterationFlags(long id, AlterationFlags flags) const
Vertex & getVertex(long id)
virtual std::size_t _getMaxHaloSize()
std::vector< adaption::Info > partitioningPrepare(MPI_Comm communicator, const std::unordered_map< long, int > &cellRanks, bool trackPartitioning, std::size_t haloSize=1)
std::set< int > getInternalCellPIDs()
long countOrphanVertices() const
void restoreInterfaces(std::istream &stream)
virtual long getCellCount() const
VertexIterator ghostVertex2InternalVertex(long id)
void unsetInterfaceAlterationFlags(AlterationFlags flags)
VertexConstIterator ghostVertexConstBegin() const
const std::unordered_map< int, std::vector< long > > & getGhostVertexExchangeTargets() const
virtual void _adaptionCleanup()
bool isBoundingBoxFrozen() const
void removePointFromBoundingBox(const std::array< double, 3 > &point)
void resetCellAdaptionMarker(long id)
std::vector< long > findCellVertexNeighs(long id, bool complete=true) const
void setAdaptionMode(AdaptionMode mode)
bool isDistributed(bool allowDirty=false) const
std::vector< adaption::Info > partitioningAlter(bool trackPartitioning=true, bool squeezeStorage=false)
bool deleteCell(long id)
virtual void _resetTol()
int getCellHaloLayer(long id) const
const std::unordered_map< int, std::vector< long > > & getGhostVertexExchangeSources() const
bool deleteInterface(long id)
virtual void _findCellNeighs(long id, const std::vector< long > *blackList, std::vector< long > *neighs) const
const std::array< double, 3 > & getVertexCoords(long id) const
virtual long _getCellNativeIndex(long id) const
long countOrphanCells() const
bool dump(std::ostream &stream)
void write(VTKWriteMode mode=VTKWriteMode::DEFAULT)
virtual void _updateInterfaces()
virtual void evalInterfaceBoundingBox(long id, std::array< double, 3 > *minPoint, std::array< double, 3 > *maxPoint) const
void updateAdjacencies(bool forcedUpdated=false)
ConstProxyVector< std::array< double BITPIT_COMMA 3 > > getInterfaceVertexCoordinates(long id) const
bool areAdjacenciesDirty(bool global=false) const
virtual void settleAdaptionMarkers()
VertexConstIterator internalVertexConstBegin() const
void setPartitioningMode(PartitioningMode mode)
const std::unordered_map< int BITPIT_COMMA std::vector< long > > & getGhostExchangeSources() const
long countBorderVertices() const
InterfaceIterator getInterfaceIterator(long id)
AdaptionStatus getAdaptionStatus(bool global=false) const
long getInternalCount() const
bool isDirty(bool global=false) const
std::vector< long > findCellEdgeNeighs(long id, bool complete=true) const
VertexIterator internalVertexBegin()
long countFreeInterfaces() const
bool isTolCustomized() const
CellConstIterator internalConstEnd() const
void consecutiveRenumberVertices(long offset=0)
void extractEnvelope(PatchKernel &envelope) const
long getInternalCellCount() const
std::vector< long > findOrphanCells() const
void displayInterfaces(std::ostream &out, unsigned int padding=0) const
bool isInterfaceOrphan(long id) const
virtual std::vector< adaption::Info > _adaptionPrepare(bool trackAdaption)
PiercedVector< Interface > & getInterfaces()
virtual ElementType getCellType(long id) const
void setInterfaceAlterationFlags(AlterationFlags flags)
bool testCellAlterationFlags(long id, AlterationFlags flags) const
void setPartitioningStatus(PartitioningStatus status)
CellConstIterator ghostCellConstBegin() const
long getInternalVertexCount() const
long countBorderFaces() const
bool findFaceNeighCell(long cellId, long neighId, int *cellFace, int *cellAdjacencyId) const
void consecutiveRenumber(long offsetVertices, long offsetCells, long offsetInterfaces)
CellConstIterator ghostConstBegin() const
virtual std::vector< adaption::Info > _adaptionAlter(bool trackAdaption)
double getTol() const
std::array< double, 3 > evalElementCentroid(const Element &element) const
virtual std::array< double, 3 > evalCellCentroid(long id) const
bool deleteCells(const IdStorage &ids)
std::vector< long > findOrphanVertices()
bool reserveCells(size_t nCells)
AlterationFlags getInterfaceAlterationFlags(long id) const
bool isBoundingBoxDirty(bool global=false) const
void resetInterfaceAlterationFlags(long id, AlterationFlags flags=FLAG_NONE)
std::vector< int > getNeighbourRanks()
void setCellAlterationFlags(AlterationFlags flags)
bool reserveVertices(size_t nVertices)
void setBoundingBox(const std::array< double, 3 > &minPoint, const std::array< double, 3 > &maxPoint)
void processCellsFaceNeighbours(const SeedContainer &seedIds, int nLayers, Function function) const
long countFreeFaces() const
virtual void _setHaloSize(std::size_t haloSize)
std::vector< long > findCellVertexOneRing(long id, int vertex) const
CellIterator internalCellEnd()
bool isAdaptionSupported() const
CellConstIterator getCellConstIterator(long id) const
ConstProxyVector< std::array< double BITPIT_COMMA 3 > > getCellVertexCoordinates(long id) const
virtual void resetCells()
int getCellOwner(long id) const
virtual void setCommunicator(MPI_Comm communicator)
std::vector< adaption::Info > adaptionAlter(bool trackAdaption=true, bool squeezeStorage=false)
virtual void _findCellVertexNeighs(long id, int vertex, const std::vector< long > *blackList, std::vector< long > *neighs) const
bool isThreeDimensional() const
bool deleteInterfaces(const IdStorage &ids)
Metafunction for generating a pierced vector.
PiercedVectorStorage< value_t, id_t >::const_range const_range
PiercedVectorStorage< value_t, id_t >::const_iterator const_iterator
PiercedVectorStorage< value_t, id_t >::range range
PiercedVectorStorage< value_t, id_t >::iterator iterator
The base class to be used to derive VTK streamers form.
Definition VTK.hpp:209
VTK input output for Unstructured Meshes.
Definition VTK.hpp:433
The Vertex class defines the vertexs.
Definition vertex.hpp:42
VTKFormat
Definition VTK.hpp:92
VTKWriteMode
Definition VTK.hpp:49
#define BITPIT_DEPRECATED(func)
Definition compiler.hpp:87
#define BITPIT_COMMA
Definition compiler.hpp:97
The namespace 'adaption' contains the routines and the data structures for handling patch adaption.
Definition adaption.cpp:38
The namespace 'patch' contains routines for interacting with the patch manager.