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
386 void getIJK(const int c, std::array<int,3>& ijk) const;
388
392 bool uniqueBoundaryIds() const;
393
396 void setUniqueBoundaryIds(bool uids);
397
398
399 // --- Dune interface below ---
400
402 // \@{
407 std::string name() const;
408
410 int maxLevel() const;
411
413 template<int codim>
414 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
416 template<int codim>
417 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
418
420 template<int codim, PartitionIteratorType PiType>
421 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
423 template<int codim, PartitionIteratorType PiType>
424 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
425
427 template<int codim>
428 typename Traits::template Codim<codim>::LeafIterator leafbegin() const;
430 template<int codim>
431 typename Traits::template Codim<codim>::LeafIterator leafend() const;
432
434 template<int codim, PartitionIteratorType PiType>
435 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const;
437 template<int codim, PartitionIteratorType PiType>
438 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const;
439
441 int size (int level, int codim) const;
442
444 int size (int codim) const;
445
447 int size (int level, GeometryType type) const;
448
450 int size (GeometryType type) const;
451
453 const Traits::GlobalIdSet& globalIdSet() const;
454
456 const Traits::LocalIdSet& localIdSet() const;
457
459 const Traits::LevelIndexSet& levelIndexSet(int level) const;
460
462 const Traits::LeafIndexSet& leafIndexSet() const;
463
467 void globalRefine (int refCount);
468
469 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
470
472 template <int codim>
474
491 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
492 const std::vector<std::array<int,3>>& startIJK_vec,
493 const std::vector<std::array<int,3>>& endIJK_vec,
494 const std::vector<std::string>& lgr_name_vec);
495
502 void autoRefine(const std::array<int,3>& nxnynz);
503
504 // @brief TO BE DONE
505 const std::map<std::string,int>& getLgrNameToLevel() const;
506
507 // @breif Compute center of an entity/element/cell in the Eclipse way:
508 // - Average of the 4 corners of the bottom face.
509 // - Average of the 4 corners of the top face.
510 // Return average of the previous computations.
511 // @param [in] int Index of a cell.
512 // @return 'eclipse centroid'
513 std::array<double,3> getEclCentroid(const int& idx) const;
514
515 // @breif Compute center of an entity/element/cell in the Eclipse way:
516 // - Average of the 4 corners of the bottom face.
517 // - Average of the 4 corners of the top face.
518 // Return average of the previous computations.
519 // @param [in] Entity<0> Entity
520 // @return 'eclipse centroid'
521 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
522
523 // @brief Return parent (coarse) intersection (face) of a refined face on the leaf grid view, whose neighboring cells
524 // are two: one coarse cell (equivalent to its origin cell from level 0), and one refined cell
525 // from certain LGR.
526 // Used in Transmissibility_impl.hpp
528
548 bool mark(int refCount, const cpgrid::Entity<0>& element);
549
553 int getMark(const cpgrid::Entity<0>& element) const;
554
557 bool preAdapt();
558
561 bool adapt();
562
576 bool adapt(const std::vector<std::array<int,3>>& cells_per_dim_vec,
577 const std::vector<int>& assignRefinedLevel,
578 const std::vector<std::string>& lgr_name_vec,
579 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
580 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
581
583 void postAdapt();
585
586 private:
587 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
588 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
589 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
590 const int& corner_count,
591 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
592 const int& preAdaptMaxLevel,
593 const int& newLevels);
594
595 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
596 const std::vector<std::array<int,3>>& cells_per_dim_vec,
597 const std::vector<int>& lgr_with_at_least_one_active_cell);
598
605 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
606 public:
613 private:
614
627 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
628
632 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
633
634 private:
642 bool nonNNCsSelectedCellsLGR( const std::vector<std::array<int,3>>& startIJK_vec,
643 const std::vector<std::array<int,3>>& endIJK_vec) const;
644
658 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
659 const std::vector<std::array<int,3>>& endIJK_vec,
660 std::vector<int>& assignRefinedLevel,
661 std::vector<int>& lgr_with_at_least_one_active_cell);
662
664 void populateCellIndexSetRefinedGrid(int level);
665
667 void populateCellIndexSetLeafGridView();
668
670 void populateLeafGlobalIdSet();
671
672 public:
673
678 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
679
681 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
682
684 unsigned int overlapSize(int) const;
685
686
688 unsigned int ghostSize(int) const;
689
691 unsigned int overlapSize(int, int) const;
692
694 unsigned int ghostSize(int, int) const;
695
697 unsigned int numBoundarySegments() const;
698
699 void setPartitioningParams(const std::map<std::string,std::string>& params);
700
701 // loadbalance is not part of the grid interface therefore we skip it.
702
712 bool loadBalance(int overlapLayers=1,
713 int partitionMethod = Dune::PartitionMethod::zoltan,
714 double imbalanceTol = 1.1,
715 int level =-1)
716 {
717 using std::get;
718 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
719 /* ownersFirst = */ false,
720 /* wells = */ nullptr,
721 /* possibleFutureConnections = */ {},
722 /* serialPartitioning = */ false,
723 /* transmissibilities = */ nullptr,
724 /* addCornerCells = */ true,
725 overlapLayers,
726 partitionMethod,
727 imbalanceTol,
728 /* allowDistributedWells = */ false,
729 /* input_cell_part = */ {},
730 level));
731 }
732
733 // loadbalance is not part of the grid interface therefore we skip it.
734
745 bool loadBalanceSerial(int overlapLayers=1,
746 int partitionMethod = Dune::PartitionMethod::zoltan,
747 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
748 double imbalanceTol = 1.1,
749 int level = -1)
750 {
751 using std::get;
752 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
753 /* ownersFirst = */ false,
754 /* wells = */ nullptr,
755 /* possibleFutureConnections = */ {},
756 /* serialPartitioning = */ true,
757 /* transmissibilities = */ nullptr,
758 /* addCornerCells = */ true,
759 overlapLayers,
760 partitionMethod,
761 imbalanceTol,
762 /* allowDistributedWells = */ false,
763 /* input_cell_part = */ {},
764 level));
765 }
766
767 // loadbalance is not part of the grid interface therefore we skip it.
768
797 std::pair<bool,std::vector<std::pair<std::string,bool>>>
798 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
799 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
800 const double* transmissibilities = nullptr,
801 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
802 int level = -1)
803 {
804 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
805 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
806 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
807 /* input_cell_part = */ {}, level);
808 }
809
810 // loadbalance is not part of the grid interface therefore we skip it.
811
845 std::pair<bool,std::vector<std::pair<std::string,bool>>>
846 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
847 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
848 const double* transmissibilities = nullptr, bool ownersFirst=false,
849 bool addCornerCells=false, int overlapLayers=1,
850 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
851 double imbalanceTol = 1.1,
852 int level = -1)
853 {
854 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
855 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
856 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
857 }
858
887 template<class DataHandle>
888 std::pair<bool, std::vector<std::pair<std::string,bool> > >
889 loadBalance(DataHandle& data,
890 const std::vector<cpgrid::OpmWellType> * wells,
891 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
892 const double* transmissibilities = nullptr,
893 int overlapLayers=1, int partitionMethod = 1, int level =-1)
894 {
895 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
896 using std::get;
897 if (get<0>(ret))
898 {
900 }
901 return ret;
902 }
903
938 template<class DataHandle>
939 std::pair<bool, std::vector<std::pair<std::string,bool> > >
940 loadBalance(DataHandle& data, EdgeWeightMethod method,
941 const std::vector<cpgrid::OpmWellType> * wells,
942 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
943 bool serialPartitioning,
944 const double* transmissibilities = nullptr, bool ownersFirst=false,
945 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
946 double imbalanceTol = 1.1,
947 bool allowDistributedWells = false)
948 {
949 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
950 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
951 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
952 using std::get;
953 if (get<0>(ret))
954 {
956 }
957 return ret;
958 }
959
984 template<class DataHandle>
985 std::pair<bool, std::vector<std::pair<std::string,bool> > >
986 loadBalance(DataHandle& data, const std::vector<int>& parts,
987 const std::vector<cpgrid::OpmWellType> * wells,
988 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
989 bool ownersFirst=false,
990 bool addCornerCells=false, int overlapLayers=1)
991 {
992 using std::get;
993 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
994 possibleFutureConnections,
995 /* serialPartitioning = */ false,
996 /* transmissibilities = */ {},
997 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
998 /* imbalanceTol (ignored) = */ 0.0,
999 /* allowDistributedWells = */ true, parts, /* level = */ 0);
1000 using std::get;
1001 if (get<0>(ret))
1002 {
1004 }
1005 return ret;
1006 }
1007
1015 template<class DataHandle>
1016 bool loadBalance(DataHandle& data,
1017 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1018 {
1019 // decltype usage needed to tell the compiler not to use this function if first
1020 // argument is std::vector but rather loadbalance by parts
1021 bool ret = loadBalance(overlapLayers, partitionMethod);
1022 if (ret)
1023 {
1025 }
1026 return ret;
1027 }
1028
1040 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1041 bool addCornerCells=false, int overlapLayers=1)
1042 {
1043 using std::get;
1044 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1045 {},
1046 /* serialPartitioning = */ false,
1047 /* trabsmissibilities = */ {},
1048 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1049 /* imbalanceTol (ignored) = */ 0.0,
1050 /* allowDistributedWells = */ true, parts));
1051 }
1052
1065 template<class DataHandle>
1066 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1067 bool addCornerCells=false, int overlapLayers=1)
1068 {
1069 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1070 if (ret)
1071 {
1073 }
1074 return ret;
1075 }
1076
1090 std::vector<int>
1091 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1092 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1093 const double* transmissibilities,
1094 const int numParts,
1095 const double imbalanceTol) const;
1096
1104 template<class DataHandle>
1105 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1106 {
1107 communicate(data, iftype, dir);
1108 }
1109
1117 template<class DataHandle>
1118 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1119
1121 const typename CpGridTraits::Communication& comm () const;
1123
1124 // ------------ End of Dune interface, start of simplified interface --------------
1125
1131
1132 // enum { dimension = 3 }; // already defined
1133
1134 typedef Dune::FieldVector<double, 3> Vector;
1135
1136
1137 const std::vector<double>& zcornData() const;
1138
1139
1140 // Topology
1145 int numCells(int level = -1) const;
1146
1151 int numFaces(int level = -1) const;
1152
1154 int numVertices() const;
1155
1156
1165 int numCellFaces(int cell, int level = -1) const;
1166
1173 int cellFace(int cell, int local_index, int level = -1) const;
1174
1178
1198 int faceCell(int face, int local_index, int level = -1) const;
1199
1206 int numCellFaces() const;
1207
1208 int numFaceVertices(int face) const;
1209
1214 int faceVertex(int face, int local_index) const;
1215
1218 double cellCenterDepth(int cell_index) const;
1219
1220
1221 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1222
1223 const Vector faceAreaNormalEcl(int face) const;
1224
1225
1226 // Geometry
1230 const Vector& vertexPosition(int vertex) const;
1231
1234 double faceArea(int face) const;
1235
1238 const Vector& faceCentroid(int face) const;
1239
1243 const Vector& faceNormal(int face) const;
1244
1247 double cellVolume(int cell) const;
1248
1251 const Vector& cellCentroid(int cell) const;
1252
1255 template<int codim>
1257 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1258 FieldVector<double, 3>,
1259 const FieldVector<double, 3>&, int>
1260 {
1261 public:
1263 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1268 : iter_(iter)
1269 {}
1270
1271 const FieldVector<double,3>& dereference() const
1272 {
1273 return iter_->center();
1274 }
1276 {
1277 ++iter_;
1278 }
1279 const FieldVector<double,3>& elementAt(int n)
1280 {
1281 return iter_[n]->center();
1282 }
1283 void advance(int n){
1284 iter_+=n;
1285 }
1287 {
1288 --iter_;
1289 }
1291 {
1292 return o-iter_;
1293 }
1294 bool equals(const CentroidIterator& o) const{
1295 return o==iter_;
1296 }
1297 private:
1299 GeometryIterator iter_;
1300 };
1301
1304
1307
1308 // Extra
1309 int boundaryId(int face) const;
1310
1317 template<class Cell2FacesRowIterator>
1318 int
1319 faceTag(const Cell2FacesRowIterator& cell_face) const;
1320
1322
1323 // ------------ End of simplified interface --------------
1324
1325 //------------- methods not in the DUNE grid interface.
1326
1331
1332
1341 template<class DataHandle>
1342 void scatterData(DataHandle& handle) const;
1343
1350 template<class DataHandle>
1351 void gatherData(DataHandle& handle) const;
1352
1355
1385
1389
1392
1396
1397#if HAVE_MPI
1402
1405
1410
1412
1414
1416
1418#endif
1419
1421 const std::vector<int>& sortedNumAquiferCells() const;
1422
1423 private:
1455 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1456 scatterGrid(EdgeWeightMethod method,
1457 bool ownersFirst,
1458 const std::vector<cpgrid::OpmWellType> * wells,
1459 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1460 bool serialPartitioning,
1461 const double* transmissibilities,
1462 bool addCornerCells,
1463 int overlapLayers,
1464 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1465 double imbalanceTol = 1.1,
1466 bool allowDistributedWells = true,
1467 const std::vector<int>& input_cell_part = {},
1468 int level = -1);
1469
1474 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1476 cpgrid::CpGridData* current_view_data_;
1478 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1480 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1482 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
1488 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1489 /*
1490 * @brief Interface for scattering and gathering point data.
1491 *
1492 * @warning Will only update owner cells
1493 */
1494 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1498 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1499
1500
1504 std::map<std::string,std::string> partitioningParams;
1505
1506 }; // end Class CpGrid
1507
1508} // end namespace Dune
1509
1513
1514
1515namespace Dune
1516{
1517
1518 namespace Capabilities
1519 {
1521 template <>
1522 struct hasEntity<CpGrid, 0>
1523 {
1524 static const bool v = true;
1525 };
1526
1528 template <>
1529 struct hasEntity<CpGrid, 3>
1530 {
1531 static const bool v = true;
1532 };
1533
1534 template<>
1535 struct canCommunicate<CpGrid,0>
1536 {
1537 static const bool v = true;
1538 };
1539
1540 template<>
1541 struct canCommunicate<CpGrid,3>
1542 {
1543 static const bool v = true;
1544 };
1545
1547 template <>
1548 struct hasBackupRestoreFacilities<CpGrid>
1549 {
1550 static const bool v = false;
1551 };
1552
1553 }
1554
1555 template<class DataHandle>
1556 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1557 {
1558 current_view_data_->communicate(data, iftype, dir);
1559 }
1560
1561
1562 template<class DataHandle>
1563 void CpGrid::scatterData([[maybe_unused]] DataHandle& handle) const
1564 {
1565#if HAVE_MPI
1566 if (distributed_data_.empty()) {
1567 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
1568 } else {
1569 distributed_data_[0]->scatterData(handle, data_[0].get(),
1570 distributed_data_[0].get(),
1573 }
1574#endif
1575 }
1576
1577 template<class DataHandle>
1578 void CpGrid::gatherData([[maybe_unused]] DataHandle& handle) const
1579 {
1580#if HAVE_MPI
1581 if (distributed_data_.empty()) {
1582 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
1583 } else {
1584 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1585 }
1586#endif
1587 }
1588
1589
1590 template<class Cell2FacesRowIterator>
1591 int
1592 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
1593 {
1594 // Note that this relies on the following implementation detail:
1595 // The grid is always constructed such that the interior faces constructed
1596 // with orientation set to true are
1597 // oriented along the positive IJK direction. Oriented means that
1598 // the first cell attached to face has the lower index.
1599 // For faces along the boundary (only one cell, always attached at index 0)
1600 // the orientation has to be determined by the orientation of the cell.
1601 // If it is true then in UnstructuredGrid it would be stored at index 0,
1602 // otherwise at index 1.
1603 const int cell = cell_face.getCellIndex();
1604 const int face = *cell_face;
1605 assert (0 <= cell); assert (cell < numCells());
1606 assert (0 <= face); assert (face < numFaces());
1607
1609
1610 const cpgrid::EntityRep<1> f(face, true);
1611 const F2C& f2c = current_view_data_->face_to_cell_[f];
1612 const face_tag tag = current_view_data_->face_tag_[f];
1613
1614 assert ((f2c.size() == 1) || (f2c.size() == 2));
1615
1616 int inside_cell = 0;
1617
1618 if ( f2c.size() == 2 ) // Two cells => interior
1619 {
1620 if ( f2c[1].index() == cell )
1621 {
1622 inside_cell = 1;
1623 }
1624 }
1625 const bool normal_is_in = ! f2c[inside_cell].orientation();
1626
1627 switch (tag) {
1628 case I_FACE:
1629 // LEFT : RIGHT
1630 return normal_is_in ? 0 : 1; // min(I) : max(I)
1631 case J_FACE:
1632 // BACK : FRONT
1633 return normal_is_in ? 2 : 3; // min(J) : max(J)
1634 case K_FACE:
1635 // Note: TOP at min(K) as 'z' measures *depth*.
1636 // TOP : BOTTOM
1637 return normal_is_in ? 4 : 5; // min(K) : max(K)
1638 case NNC_FACE:
1639 // For nnc faces we return the otherwise unused value -1.
1640 return -1;
1641 default:
1642 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
1643 }
1644 }
1645
1646 template<int dim>
1648
1649} // namespace Dune
1650
1653#include "cpgrid/Intersection.hpp"
1654#include "cpgrid/Geometry.hpp"
1655#include "cpgrid/Indexsets.hpp"
1656
1657#endif // OPM_CPGRID_HEADER
DataHandle & data
Definition: CpGridData.hpp:1183
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:1260
void increment()
Definition: CpGrid.hpp:1275
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:1267
const FieldVector< double, 3 > & dereference() const
Definition: CpGrid.hpp:1271
int distanceTo(const CentroidIterator &o)
Definition: CpGrid.hpp:1290
void decrement()
Definition: CpGrid.hpp:1286
void advance(int n)
Definition: CpGrid.hpp:1283
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:1264
bool equals(const CentroidIterator &o) const
Definition: CpGrid.hpp:1294
const FieldVector< double, 3 > & elementAt(int n)
Definition: CpGrid.hpp:1279
[ 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:986
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:1578
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:940
cpgrid::CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGrid.hpp:1399
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:1592
const std::vector< double > & zcornData() const
ParallelIndexSet & getCellIndexSet()
void setUniqueBoundaryIds(bool uids)
std::array< double, 3 > getEclCentroid(const cpgrid::Entity< 0 > &elem) const
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:846
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:1354
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:889
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 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:745
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:712
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:1401
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.
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)
Create a grid out of a coarse one and (at most) 2 refinements(LGRs) of selected block-shaped disjoint...
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).
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:1016
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:798
friend cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
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:1134
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:1066
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:1040
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...
bool adapt(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.
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
communicate objects for all codims on a given level
Definition: CpGrid.hpp:1105
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:1404
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:1563
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:118
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:961
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