CpGrid.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: CpGrid.hpp
4//
5// Created: Fri May 29 20:26:36 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// B�rd Skaflestad <bard.skaflestad@sintef.no>
9// Antonella Ritorto <antonella.ritorto@opm-op.com>
10//
11// $Date$
12//
13// $Revision$
14//
15//===========================================================================
16
17/*
18 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
19 Copyright 2009, 2010, 2014, 2022-2023 Equinor ASA.
20 Copyright 2014, 2015 Dr. Blatt - HPC-Simulartion-Software & Services
21 Copyright 2015 NTNU
22
23 This file is part of The Open Porous Media project (OPM).
24
25 OPM is free software: you can redistribute it and/or modify
26 it under the terms of the GNU General Public License as published by
27 the Free Software Foundation, either version 3 of the License, or
28 (at your option) any later version.
29
30 OPM is distributed in the hope that it will be useful,
31 but WITHOUT ANY WARRANTY; without even the implied warranty of
32 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
33 GNU General Public License for more details.
34
35 You should have received a copy of the GNU General Public License
36 along with OPM. If not, see <http://www.gnu.org/licenses/>.
37*/
38
39#ifndef OPM_CPGRID_HEADER
40#define OPM_CPGRID_HEADER
41
42// Warning suppression for Dune includes.
44
45#include <dune/grid/common/grid.hh>
46
48
50
54
56
58
59#include <set>
60
61namespace Opm
62{
63struct NNCdata;
64class EclipseGrid;
65class EclipseState;
66}
67
68namespace Dune
69{
70 class CpGrid;
71
72 namespace cpgrid
73 {
74 class CpGridData;
75 template <int> class Entity;
76 template<int,int> class Geometry;
77 class HierarchicIterator;
78 class IntersectionIterator;
79 template<int, PartitionIteratorType> class Iterator;
80 class LevelGlobalIdSet;
81 class GlobalIdSet;
82 class Intersection;
83 class IntersectionIterator;
84 class IndexSet;
85 class IdSet;
86
87 }
88}
89
90namespace Dune
91{
92
94 //
95 // CpGridTraits
96 //
98
100 {
102 typedef CpGrid Grid;
103
112
115
118 template <int cd>
119 struct Codim
120 {
123 typedef cpgrid::Geometry<3-cd, 3> Geometry;
124 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> Geometry;
127 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> LocalGeometry;
130
133
136
139
142 template <PartitionIteratorType pitype>
144 {
149 };
150 };
151
154 template <PartitionIteratorType pitype>
156 {
158 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid> > LevelGridView;
160 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> > LeafGridView;
161
162 };
163
165 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid>> LevelGridView;
167 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>> LeafGridView;
168
177
181 };
182
184 //
185 // CpGridFamily
186 //
188
190 {
192 };
193
195 //
196 // CpGrid
197 //
199
201 class CpGrid
202 : public GridDefaultImplementation<3, 3, double, CpGridFamily>
203 {
204 friend class cpgrid::CpGridData;
205 friend class cpgrid::Entity<0>;
206 friend class cpgrid::Entity<1>;
207 friend class cpgrid::Entity<2>;
208 friend class cpgrid::Entity<3>;
209 template<int dim>
210 friend cpgrid::Entity<dim> createEntity(const CpGrid&,int,bool);
211
212 public:
213
214 // --- Typedefs ---
215
216
219
220
221 // --- Methods ---
222
223
226
227 explicit CpGrid(MPIHelper::MPICommunicator comm);
228
229#if HAVE_ECL_INPUT
266 std::vector<std::size_t>
267 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
268 Opm::EclipseState* ecl_state,
269 bool periodic_extension,
270 bool turn_normals,
271 bool clip_z,
272 bool pinchActive,
273 bool edge_conformal);
274
312 std::vector<std::size_t>
313 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
314 Opm::EclipseState* ecl_state,
315 bool periodic_extension,
316 bool turn_normals = false,
317 bool clip_z = false,
318 bool edge_conformal = false);
319#endif // HAVE_ECL_INPUT
320
335 void processEclipseFormat(const grdecl& input_data,
336 bool remove_ij_boundary,
337 bool turn_normals = false,
338 bool edge_conformal = false);
339
341
347
348
355 void createCartesian(const std::array<int, 3>& dims,
356 const std::array<double, 3>& cellsize,
357 const std::array<int, 3>& shift = {0,0,0});
358
362 const std::array<int, 3>& logicalCartesianSize() const;
363
371 const std::vector<int>& globalCell() const;
372
374 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData() const;
375
377 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData();
378
381
384
392 void getIJK(const int c, std::array<int,3>& ijk) const;
394
398 bool uniqueBoundaryIds() const;
399
402 void setUniqueBoundaryIds(bool uids);
403
404
405 // --- Dune interface below ---
406
408 // \@{
413 std::string name() const;
414
416 int maxLevel() const;
417
419 template<int codim>
420 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
422 template<int codim>
423 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
424
426 template<int codim, PartitionIteratorType PiType>
427 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
429 template<int codim, PartitionIteratorType PiType>
430 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
431
433 template<int codim>
434 typename Traits::template Codim<codim>::LeafIterator leafbegin() const;
436 template<int codim>
437 typename Traits::template Codim<codim>::LeafIterator leafend() const;
438
440 template<int codim, PartitionIteratorType PiType>
441 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const;
443 template<int codim, PartitionIteratorType PiType>
444 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const;
445
447 int size (int level, int codim) const;
448
450 int size (int codim) const;
451
453 int size (int level, GeometryType type) const;
454
456 int size (GeometryType type) const;
457
459 const Traits::GlobalIdSet& globalIdSet() const;
460
462 const Traits::LocalIdSet& localIdSet() const;
463
465 const Traits::LevelIndexSet& levelIndexSet(int level) const;
466
468 const Traits::LeafIndexSet& leafIndexSet() const;
469
473 void globalRefine (int refCount);
474
475 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
476
478 template <int codim>
480
499 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
500 const std::vector<std::array<int,3>>& startIJK_vec,
501 const std::vector<std::array<int,3>>& endIJK_vec,
502 const std::vector<std::string>& lgr_name_vec,
503 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
504
511 void autoRefine(const std::array<int,3>& nxnynz);
512
513 // @brief TO BE DONE
514 const std::map<std::string,int>& getLgrNameToLevel() const;
515
516 // @breif Compute center of an entity/element/cell in the Eclipse way:
517 // - Average of the 4 corners of the bottom face.
518 // - Average of the 4 corners of the top face.
519 // Return average of the previous computations.
520 // @param [in] int Index of a cell.
521 // @return 'eclipse centroid'
522 std::array<double,3> getEclCentroid(const int& idx) const;
523
524 // @breif Compute center of an entity/element/cell in the Eclipse way:
525 // - Average of the 4 corners of the bottom face.
526 // - Average of the 4 corners of the top face.
527 // Return average of the previous computations.
528 // @param [in] Entity<0> Entity
529 // @return 'eclipse centroid'
530 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
531
532 // @brief Return parent (coarse) intersection (face) of a refined face on the leaf grid view, whose neighboring cells
533 // are two: one coarse cell (equivalent to its origin cell from level 0), and one refined cell
534 // from certain LGR.
535 // Used in Transmissibility_impl.hpp
537
557 bool mark(int refCount, const cpgrid::Entity<0>& element);
558
562 int getMark(const cpgrid::Entity<0>& element) const;
563
566 bool preAdapt();
567
570 bool adapt();
571
585 bool refineAndUpdateGrid(const std::vector<std::array<int,3>>& cells_per_dim_vec,
586 const std::vector<int>& assignRefinedLevel,
587 const std::vector<std::string>& lgr_name_vec,
588 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
589 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
590
592 void postAdapt();
594
595 private:
596 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
597 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
598 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
599 const int& corner_count,
600 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
601 const int& preAdaptMaxLevel,
602 const int& newLevels);
603
604 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
605 const std::vector<std::array<int,3>>& cells_per_dim_vec,
606 const std::vector<int>& lgr_with_at_least_one_active_cell);
607
614 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
615 public:
622 private:
623
636 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
637
641 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
642
643 private:
657 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
658 const std::vector<std::array<int,3>>& endIJK_vec,
659 std::vector<int>& assignRefinedLevel,
660 std::vector<int>& lgr_with_at_least_one_active_cell);
661
663 void populateCellIndexSetRefinedGrid(int level);
664
666 void populateCellIndexSetLeafGridView();
667
669 void populateLeafGlobalIdSet();
670
671 public:
672
677 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
678
680 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
681
683 unsigned int overlapSize(int) const;
684
685
687 unsigned int ghostSize(int) const;
688
690 unsigned int overlapSize(int, int) const;
691
693 unsigned int ghostSize(int, int) const;
694
696 unsigned int numBoundarySegments() const;
697
698 void setPartitioningParams(const std::map<std::string,std::string>& params);
699
700 // loadbalance is not part of the grid interface therefore we skip it.
701
711 bool loadBalance(int overlapLayers=1,
712 int partitionMethod = Dune::PartitionMethod::zoltan,
713 double imbalanceTol = 1.1,
714 int level =-1)
715 {
716 using std::get;
717 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
718 /* ownersFirst = */ false,
719 /* wells = */ nullptr,
720 /* possibleFutureConnections = */ {},
721 /* serialPartitioning = */ false,
722 /* transmissibilities = */ nullptr,
723 /* addCornerCells = */ true,
724 overlapLayers,
725 partitionMethod,
726 imbalanceTol,
727 /* allowDistributedWells = */ false,
728 /* input_cell_part = */ {},
729 level));
730 }
731
732 // loadbalance is not part of the grid interface therefore we skip it.
733
744 bool loadBalanceSerial(int overlapLayers=1,
745 int partitionMethod = Dune::PartitionMethod::zoltan,
746 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
747 double imbalanceTol = 1.1,
748 int level = -1)
749 {
750 using std::get;
751 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
752 /* ownersFirst = */ false,
753 /* wells = */ nullptr,
754 /* possibleFutureConnections = */ {},
755 /* serialPartitioning = */ true,
756 /* transmissibilities = */ nullptr,
757 /* addCornerCells = */ true,
758 overlapLayers,
759 partitionMethod,
760 imbalanceTol,
761 /* allowDistributedWells = */ false,
762 /* input_cell_part = */ {},
763 level));
764 }
765
766 // loadbalance is not part of the grid interface therefore we skip it.
767
796 std::pair<bool,std::vector<std::pair<std::string,bool>>>
797 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
798 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
799 const double* transmissibilities = nullptr,
800 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
801 int level = -1)
802 {
803 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
804 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
805 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
806 /* input_cell_part = */ {}, level);
807 }
808
809 // loadbalance is not part of the grid interface therefore we skip it.
810
844 std::pair<bool,std::vector<std::pair<std::string,bool>>>
845 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
846 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
847 const double* transmissibilities = nullptr, bool ownersFirst=false,
848 bool addCornerCells=false, int overlapLayers=1,
849 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
850 double imbalanceTol = 1.1,
851 int level = -1)
852 {
853 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
854 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
855 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
856 }
857
886 template<class DataHandle>
887 std::pair<bool, std::vector<std::pair<std::string,bool> > >
888 loadBalance(DataHandle& data,
889 const std::vector<cpgrid::OpmWellType> * wells,
890 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
891 const double* transmissibilities = nullptr,
892 int overlapLayers=1, int partitionMethod = 1, int level =-1)
893 {
894 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
895 using std::get;
896 if (get<0>(ret))
897 {
899 }
900 return ret;
901 }
902
937 template<class DataHandle>
938 std::pair<bool, std::vector<std::pair<std::string,bool> > >
939 loadBalance(DataHandle& data, EdgeWeightMethod method,
940 const std::vector<cpgrid::OpmWellType> * wells,
941 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
942 bool serialPartitioning,
943 const double* transmissibilities = nullptr, bool ownersFirst=false,
944 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
945 double imbalanceTol = 1.1,
946 bool allowDistributedWells = false)
947 {
948 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
949 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
950 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
951 using std::get;
952 if (get<0>(ret))
953 {
955 }
956 return ret;
957 }
958
983 template<class DataHandle>
984 std::pair<bool, std::vector<std::pair<std::string,bool> > >
985 loadBalance(DataHandle& data, const std::vector<int>& parts,
986 const std::vector<cpgrid::OpmWellType> * wells,
987 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
988 bool ownersFirst=false,
989 bool addCornerCells=false, int overlapLayers=1)
990 {
991 using std::get;
992 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
993 possibleFutureConnections,
994 /* serialPartitioning = */ false,
995 /* transmissibilities = */ {},
996 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
997 /* imbalanceTol (ignored) = */ 0.0,
998 /* allowDistributedWells = */ true, parts, /* level = */ 0);
999 using std::get;
1000 if (get<0>(ret))
1001 {
1003 }
1004 return ret;
1005 }
1006
1014 template<class DataHandle>
1015 bool loadBalance(DataHandle& data,
1016 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1017 {
1018 // decltype usage needed to tell the compiler not to use this function if first
1019 // argument is std::vector but rather loadbalance by parts
1020 bool ret = loadBalance(overlapLayers, partitionMethod);
1021 if (ret)
1022 {
1024 }
1025 return ret;
1026 }
1027
1039 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1040 bool addCornerCells=false, int overlapLayers=1)
1041 {
1042 using std::get;
1043 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1044 {},
1045 /* serialPartitioning = */ false,
1046 /* trabsmissibilities = */ {},
1047 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1048 /* imbalanceTol (ignored) = */ 0.0,
1049 /* allowDistributedWells = */ true, parts));
1050 }
1051
1064 template<class DataHandle>
1065 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1066 bool addCornerCells=false, int overlapLayers=1)
1067 {
1068 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1069 if (ret)
1070 {
1072 }
1073 return ret;
1074 }
1075
1089 std::vector<int>
1090 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1091 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1092 const double* transmissibilities,
1093 const int numParts,
1094 const double imbalanceTol) const;
1095
1103 template<class DataHandle>
1104 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1105 {
1106 communicate(data, iftype, dir);
1107 }
1108
1116 template<class DataHandle>
1117 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1118
1120 const typename CpGridTraits::Communication& comm () const;
1122
1123 // ------------ End of Dune interface, start of simplified interface --------------
1124
1130
1131 // enum { dimension = 3 }; // already defined
1132
1133 typedef Dune::FieldVector<double, 3> Vector;
1134
1135
1136 const std::vector<double>& zcornData() const;
1137
1138
1139 // Topology
1144 int numCells(int level = -1) const;
1145
1150 int numFaces(int level = -1) const;
1151
1153 int numVertices() const;
1154
1155
1164 int numCellFaces(int cell, int level = -1) const;
1165
1172 int cellFace(int cell, int local_index, int level = -1) const;
1173
1177
1197 int faceCell(int face, int local_index, int level = -1) const;
1198
1205 int numCellFaces() const;
1206
1207 int numFaceVertices(int face) const;
1208
1213 int faceVertex(int face, int local_index) const;
1214
1217 double cellCenterDepth(int cell_index) const;
1218
1219
1220 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1221
1222 const Vector faceAreaNormalEcl(int face) const;
1223
1224
1225 // Geometry
1229 const Vector& vertexPosition(int vertex) const;
1230
1233 double faceArea(int face) const;
1234
1237 const Vector& faceCentroid(int face) const;
1238
1242 const Vector& faceNormal(int face) const;
1243
1246 double cellVolume(int cell) const;
1247
1250 const Vector& cellCentroid(int cell) const;
1251
1254 template<int codim>
1256 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1257 FieldVector<double, 3>,
1258 const FieldVector<double, 3>&, int>
1259 {
1260 public:
1262 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1267 : iter_(iter)
1268 {}
1269
1270 const FieldVector<double,3>& dereference() const
1271 {
1272 return iter_->center();
1273 }
1275 {
1276 ++iter_;
1277 }
1278 const FieldVector<double,3>& elementAt(int n)
1279 {
1280 return iter_[n]->center();
1281 }
1282 void advance(int n){
1283 iter_+=n;
1284 }
1286 {
1287 --iter_;
1288 }
1290 {
1291 return o-iter_;
1292 }
1293 bool equals(const CentroidIterator& o) const{
1294 return o==iter_;
1295 }
1296 private:
1298 GeometryIterator iter_;
1299 };
1300
1303
1306
1307 // Extra
1308 int boundaryId(int face) const;
1309
1316 template<class Cell2FacesRowIterator>
1317 int
1318 faceTag(const Cell2FacesRowIterator& cell_face) const;
1319
1321
1322 // ------------ End of simplified interface --------------
1323
1324 //------------- methods not in the DUNE grid interface.
1325
1330
1331
1340 template<class DataHandle>
1341 void scatterData(DataHandle& handle) const;
1342
1349 template<class DataHandle>
1350 void gatherData(DataHandle& handle) const;
1351
1354
1384
1388
1391
1395
1396#if HAVE_MPI
1401
1404
1409
1411
1413
1415
1417#endif
1418
1420 const std::vector<int>& sortedNumAquiferCells() const;
1421
1422 private:
1454 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1455 scatterGrid(EdgeWeightMethod method,
1456 bool ownersFirst,
1457 const std::vector<cpgrid::OpmWellType> * wells,
1458 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1459 bool serialPartitioning,
1460 const double* transmissibilities,
1461 bool addCornerCells,
1462 int overlapLayers,
1463 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1464 double imbalanceTol = 1.1,
1465 bool allowDistributedWells = true,
1466 const std::vector<int>& input_cell_part = {},
1467 int level = -1);
1468
1473 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1475 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1477 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1479 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
1485 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1486 /*
1487 * @brief Interface for scattering and gathering point data.
1488 *
1489 * @warning Will only update owner cells
1490 */
1491 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1495 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1496
1497
1501 std::map<std::string,std::string> partitioningParams;
1502
1503 }; // end Class CpGrid
1504
1505} // end namespace Dune
1506
1510
1511
1512namespace Dune
1513{
1514
1515 namespace Capabilities
1516 {
1518 template <>
1519 struct hasEntity<CpGrid, 0>
1520 {
1521 static const bool v = true;
1522 };
1523
1525 template <>
1526 struct hasEntity<CpGrid, 3>
1527 {
1528 static const bool v = true;
1529 };
1530
1531 template<>
1532 struct canCommunicate<CpGrid,0>
1533 {
1534 static const bool v = true;
1535 };
1536
1537 template<>
1538 struct canCommunicate<CpGrid,3>
1539 {
1540 static const bool v = true;
1541 };
1542
1544 template <>
1545 struct hasBackupRestoreFacilities<CpGrid>
1546 {
1547 static const bool v = false;
1548 };
1549
1550 }
1551
1552 template<class DataHandle>
1553 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1554 {
1555 current_data_->back()->communicate(data, iftype, dir);
1556 }
1557
1558
1559 template<class DataHandle>
1560 void CpGrid::scatterData([[maybe_unused]] DataHandle& handle) const
1561 {
1562#if HAVE_MPI
1563 if (distributed_data_.empty()) {
1564 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
1565 } else {
1566 distributed_data_[0]->scatterData(handle, data_[0].get(),
1567 distributed_data_[0].get(),
1570 }
1571#endif
1572 }
1573
1574 template<class DataHandle>
1575 void CpGrid::gatherData([[maybe_unused]] DataHandle& handle) const
1576 {
1577#if HAVE_MPI
1578 if (distributed_data_.empty()) {
1579 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
1580 } else {
1581 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1582 }
1583#endif
1584 }
1585
1586
1587 template<class Cell2FacesRowIterator>
1588 int
1589 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
1590 {
1591 // Note that this relies on the following implementation detail:
1592 // The grid is always constructed such that the interior faces constructed
1593 // with orientation set to true are
1594 // oriented along the positive IJK direction. Oriented means that
1595 // the first cell attached to face has the lower index.
1596 // For faces along the boundary (only one cell, always attached at index 0)
1597 // the orientation has to be determined by the orientation of the cell.
1598 // If it is true then in UnstructuredGrid it would be stored at index 0,
1599 // otherwise at index 1.
1600 const int cell = cell_face.getCellIndex();
1601 const int face = *cell_face;
1602 assert (0 <= cell); assert (cell < numCells());
1603 assert (0 <= face); assert (face < numFaces());
1604
1606
1607 const cpgrid::EntityRep<1> f(face, true);
1608 const F2C& f2c = current_data_->back()->face_to_cell_[f];
1609 const face_tag tag = current_data_->back()->face_tag_[f];
1610
1611 assert ((f2c.size() == 1) || (f2c.size() == 2));
1612
1613 int inside_cell = 0;
1614
1615 if ( f2c.size() == 2 ) // Two cells => interior
1616 {
1617 if ( f2c[1].index() == cell )
1618 {
1619 inside_cell = 1;
1620 }
1621 }
1622 const bool normal_is_in = ! f2c[inside_cell].orientation();
1623
1624 switch (tag) {
1625 case I_FACE:
1626 // LEFT : RIGHT
1627 return normal_is_in ? 0 : 1; // min(I) : max(I)
1628 case J_FACE:
1629 // BACK : FRONT
1630 return normal_is_in ? 2 : 3; // min(J) : max(J)
1631 case K_FACE:
1632 // Note: TOP at min(K) as 'z' measures *depth*.
1633 // TOP : BOTTOM
1634 return normal_is_in ? 4 : 5; // min(K) : max(K)
1635 case NNC_FACE:
1636 // For nnc faces we return the otherwise unused value -1.
1637 return -1;
1638 default:
1639 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
1640 }
1641 }
1642
1643 template<int dim>
1645
1646} // namespace Dune
1647
1650#include "cpgrid/Intersection.hpp"
1651#include "cpgrid/Geometry.hpp"
1652#include "cpgrid/Indexsets.hpp"
1653
1654#endif // OPM_CPGRID_HEADER
DataHandle & data
Definition: CpGridData.hpp:1165
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:1259
void increment()
Definition: CpGrid.hpp:1274
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:1266
const FieldVector< double, 3 > & dereference() const
Definition: CpGrid.hpp:1270
int distanceTo(const CentroidIterator &o)
Definition: CpGrid.hpp:1289
void decrement()
Definition: CpGrid.hpp:1285
void advance(int n)
Definition: CpGrid.hpp:1282
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:1263
bool equals(const CentroidIterator &o) const
Definition: CpGrid.hpp:1293
const FieldVector< double, 3 > & elementAt(int n)
Definition: CpGrid.hpp:1278
[ provides Dune::Grid ]
Definition: CpGrid.hpp:203
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(DataHandle &data, const std::vector< int > &parts, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:985
std::string name() const
Get the grid name.
std::vector< std::unordered_map< std::size_t, std::size_t > > mapLocalCartesianIndexSetsToLeafIndexSet() const
Compute for each level grid, a map from the global_cell_[ cell index in level grid ] to the leaf inde...
CentroidIterator< 0 > beginCellCentroids() const
Get an iterator over the cell centroids positioned at the first one.
void switchToGlobalView()
Switch to the global view.
const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection &intersection) const
int numFaceVertices(int face) const
unsigned int numBoundarySegments() const
returns the number of boundary segments within the macro grid
CpGridFamily GridFamily
Family typedef, why is this not defined by Grid<>?
Definition: CpGrid.hpp:218
Traits::template Codim< codim >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
void gatherData(DataHandle &handle) const
Moves data from the distributed view to the global (all data on process) view.
Definition: CpGrid.hpp:1575
CentroidIterator< 1 > beginFaceCentroids() const
Get an iterator over the face centroids positioned at the first one.
int numFaces(int level=-1) const
Get the number of faces.
std::vector< std::array< int, 2 > > mapLeafIndexSetToLocalCartesianIndexSets() const
Reverse map: from leaf index cell to { level, local/level Cartesian index of the cell }.
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData() const
Returns either data_ or distributed_data_(if non empty).
int size(GeometryType type) const
number of leaf entities per geometry type in this process
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(DataHandle &data, EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, bool serialPartitioning, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, double imbalanceTol=1.1, bool allowDistributedWells=false)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:939
cpgrid::CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGrid.hpp:1398
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:1589
const std::vector< double > & zcornData() const
ParallelIndexSet & getCellIndexSet()
void setUniqueBoundaryIds(bool uids)
std::array< double, 3 > getEclCentroid(const cpgrid::Entity< 0 > &elem) const
bool refineAndUpdateGrid(const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< int > &assignRefinedLevel, const std::vector< std::string > &lgr_name_vec, const std::vector< std::array< int, 3 > > &startIJK_vec=std::vector< std::array< int, 3 > >{}, const std::vector< std::array< int, 3 > > &endIJK_vec=std::vector< std::array< int, 3 > >{})
Triggers the grid refinement process, allowing to select diffrent refined level grids.
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, double imbalanceTol=1.1, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:845
void createCartesian(const std::array< int, 3 > &dims, const std::array< double, 3 > &cellsize, const std::array< int, 3 > &shift={0, 0, 0})
Traits::template Codim< codim >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
unsigned int overlapSize(int) const
Size of the overlap on the leaf level.
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGrid.hpp:1353
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(DataHandle &data, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, int overlapLayers=1, int partitionMethod=1, int level=-1)
Distributes this grid and data over the available nodes in a distributed machine.
Definition: CpGrid.hpp:888
Traits::template Codim< codim >::LevelIterator lend(int level) const
one past the end on this level
int numCellFaces() const
Get the sum of all faces attached to all cells.
const Dune::cpgrid::CpGridData & currentLeafData() const
Returns current view data (the leaf grid)
const std::vector< int > & globalCell() const
bool loadBalanceSerial(int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan, int edgeWeightMethod=Dune::EdgeWeightMethod::defaultTransEdgeWgt, double imbalanceTol=1.1, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:744
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
bool loadBalance(int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan, double imbalanceTol=1.1, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:711
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
const CpGridTraits::Communication & comm() const
Get the collective communication object.
double faceArea(int face) const
Get the area of a face.
const std::vector< Dune::GeometryType > & geomTypes(const int) const
void postAdapt()
Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'.
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
int size(int level, int codim) const
Number of grid entities per level and codim.
const std::array< int, 3 > & logicalCartesianSize() const
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
unsigned int ghostSize(int, int) const
Size of the ghost cell layer on a given level.
cpgrid::CpGridDataTraits::RemoteIndices RemoteIndices
The type of the remote indices information.
Definition: CpGrid.hpp:1400
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
int numVertices() const
Get The number of vertices.
Dune::cpgrid::Intersection getParentIntersectionFromLgrBoundaryFace(const Dune::cpgrid::Intersection &intersection) const
void getIJK(const int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
void processEclipseFormat(const grdecl &input_data, bool remove_ij_boundary, bool turn_normals=false, bool edge_conformal=false)
void syncDistributedGlobalCellIds()
Synchronizes cell global ids across processes after load balancing.
bool uniqueBoundaryIds() const
int faceCell(int face, int local_index, int level=-1) const
Get the index identifying a cell attached to a face.
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
void setPartitioningParams(const std::map< std::string, std::string > &params)
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
CpGrid(MPIHelper::MPICommunicator comm)
const ParallelIndexSet & getCellIndexSet() const
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData()
Returns either data_ or distributed_data_(if non empty).
int boundaryId(int face) const
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim and PartitionIteratorType.
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
Mark entity for refinement (or coarsening).
void addLgrsUpdateLeafView(const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, const std::vector< std::string > &lgr_name_vec, const std::vector< std::string > &lgr_parent_grid_name_vec=std::vector< std::string >{})
Create a grid out of a coarse one and (at most) 2 refinements(LGRs) of selected block-shaped disjoint...
bool loadBalance(DataHandle &data, decltype(data.fixedSize(0, 0)) overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan)
Distributes this grid and data over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1015
unsigned int overlapSize(int, int) const
Size of the overlap on a given level.
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:797
friend cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
Dune::cpgrid::CpGridData & currentLeafData()
Returns current view data (the leaf grid)
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
const Vector & faceNormal(int face) const
Get the unit normal of a face.
Dune::FieldVector< double, 3 > Vector
Definition: CpGrid.hpp:1133
CpGrid()
Default constructor.
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level and PartitionIteratorType.
const std::map< std::string, int > & getLgrNameToLevel() const
const std::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
bool adapt()
Triggers the grid refinement process. Returns true if the grid has changed, false otherwise.
void switchToDistributedView()
Switch to the distributed view.
const RemoteIndices & getCellRemoteIndices() const
const Vector faceAreaNormalEcl(int face) const
bool preAdapt()
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called af...
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
bool loadBalance(DataHandle &data, const std::vector< int > &parts, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1)
Distributes this grid and data over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1065
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
const cpgrid::OrientedEntityTable< 0, 1 >::row_type cellFaceRow(int cell) const
Get a list of indices identifying all faces of a cell.
bool loadBalance(const std::vector< int > &parts, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1039
std::vector< int > zoltanPartitionWithoutScatter(const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const double *transmissibilities, const int numParts, const double imbalanceTol) const
Partitions the grid using Zoltan without decomposing and distributing it among processes.
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
std::array< double, 3 > getEclCentroid(const int &idx) const
int numCells(int level=-1) const
Get the number of cells.
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
void globalRefine(int refCount)
Refine the grid refCount times using the default refinement rule. This behaves like marking all eleme...
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
communicate objects for all codims on a given level
Definition: CpGrid.hpp:1104
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
RemoteIndices & getCellRemoteIndices()
cpgrid::CpGridDataTraits::CommunicationType CommunicationType
The type of the owner-overlap-copy communication.
Definition: CpGrid.hpp:1403
int size(int codim) const
number of leaf entities per codim in this process
const InterfaceMap & pointScatterGatherInterface() const
Get an interface for gathering/scattering data attached to points with communication.
int cellFace(int cell, int local_index, int level=-1) const
Get a specific face of a cell.
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
int numCellFaces(int cell, int level=-1) const
Get the number of faces of a cell.
double cellVolume(int cell) const
Get the volume of the cell.
void autoRefine(const std::array< int, 3 > &nxnynz)
Global refine the grid with different refinement factors in each direction.
const CommunicationType & cellCommunication() const
Get the owner-overlap-copy communication for cells.
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition: CpGrid.hpp:1560
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:118
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:99
Definition: Geometry.hpp:76
The global id set for Dune.
Definition: Indexsets.hpp:487
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:118
Definition: Indexsets.hpp:56
Definition: Intersection.hpp:279
Definition: Intersection.hpp:66
Definition: Iterators.hpp:60
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
@ simple
Use simple approach based on rectangular partitioning the underlying cartesian grid.
Definition: GridEnums.hpp:46
@ zoltanGoG
use Zoltan on GraphOfGrid for partitioning
Definition: GridEnums.hpp:52
@ zoltan
Use Zoltan for partitioning.
Definition: GridEnums.hpp:48
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Me...
Definition: GridEnums.hpp:34
@ defaultTransEdgeWgt
Use the transmissibilities as edge weights.
Definition: GridEnums.hpp:38
void getFirstChildGlobalIds(const Dune::CpGrid &grid, std::vector< int > &parentToFirstChildGlobalIds)
Retrieves the global ids of the first child for each parent cell in the grid.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
face_tag
Definition: preprocess.h:66
@ K_FACE
Definition: preprocess.h:69
@ J_FACE
Definition: preprocess.h:68
@ NNC_FACE
Definition: preprocess.h:70
@ I_FACE
Definition: preprocess.h:67
Definition: CpGrid.hpp:190
CpGridTraits Traits
Definition: CpGrid.hpp:191
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:144
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition: CpGrid.hpp:146
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition: CpGrid.hpp:148
Traits associated with a specific codim.
Definition: CpGrid.hpp:120
cpgrid::Entity< cd > Entity
The type of the entity.
Definition: CpGrid.hpp:129
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity. IMPORTANT: Codim<codim>::Geometry == Geometry<di...
Definition: CpGrid.hpp:123
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition: CpGrid.hpp:126
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition: CpGrid.hpp:135
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition: CpGrid.hpp:132
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition: CpGrid.hpp:138
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:156
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:160
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:158
Definition: CpGrid.hpp:100
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition: CpGrid.hpp:170
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition: CpGrid.hpp:109
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:167
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:165
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition: CpGrid.hpp:174
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition: CpGrid.hpp:111
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition: CpGrid.hpp:172
cpgrid::CpGridDataTraits::CollectiveCommunication CollectiveCommunication
Definition: CpGrid.hpp:180
GlobalIdSet LocalIdSet
The type of the local id set.
Definition: CpGrid.hpp:176
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGrid.hpp:179
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition: CpGrid.hpp:114
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition: CpGrid.hpp:107
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition: CpGrid.hpp:105
CpGrid Grid
The type that implements the grid.
Definition: CpGrid.hpp:102
Dune::RemoteIndices< ParallelIndexSet > RemoteIndices
The type of the remote indices information.
Definition: CpGridDataTraits.hpp:83
typename CommunicationType::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGridDataTraits.hpp:80
Dune::Communication< MPICommunicator > CollectiveCommunication
Definition: CpGridDataTraits.hpp:59
Dune::OwnerOverlapCopyCommunication< int, int > CommunicationType
type of OwnerOverlap communication for cells
Definition: CpGridDataTraits.hpp:77
Dune::Communication< MPICommunicator > Communication
Definition: CpGridDataTraits.hpp:58
Communicator::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGridDataTraits.hpp:74
Definition: preprocess.h:56