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_OPM_COMMON
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_OPM_COMMON
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
475 void globalRefine (int refCount, bool throwOnFailure = false);
476
477 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
478
480 template <int codim>
482
501 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
502 const std::vector<std::array<int,3>>& startIJK_vec,
503 const std::vector<std::array<int,3>>& endIJK_vec,
504 const std::vector<std::string>& lgr_name_vec,
505 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
506
513 void autoRefine(const std::array<int,3>& nxnynz);
514
515 // @brief TO BE DONE
516 const std::map<std::string,int>& getLgrNameToLevel() const;
517
518 // @breif Compute center of an entity/element/cell in the Eclipse way:
519 // - Average of the 4 corners of the bottom face.
520 // - Average of the 4 corners of the top face.
521 // Return average of the previous computations.
522 // @param [in] int Index of a cell.
523 // @return 'eclipse centroid'
524 std::array<double,3> getEclCentroid(const int& idx) const;
525
526 // @breif Compute center of an entity/element/cell in the Eclipse way:
527 // - Average of the 4 corners of the bottom face.
528 // - Average of the 4 corners of the top face.
529 // Return average of the previous computations.
530 // @param [in] Entity<0> Entity
531 // @return 'eclipse centroid'
532 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
533
534 // @brief Return parent (coarse) intersection (face) of a refined face on the leaf grid view, whose neighboring cells
535 // are two: one coarse cell (equivalent to its origin cell from level 0), and one refined cell
536 // from certain LGR.
537 // Used in Transmissibility_impl.hpp
539
560 bool mark(int refCount, const cpgrid::Entity<0>& element, bool throwOnFailure = false);
561
565 int getMark(const cpgrid::Entity<0>& element) const;
566
569 bool preAdapt();
570
573 bool adapt();
574
593 bool refineAndUpdateGrid(bool throwOnFailure,
594 const std::vector<std::array<int,3>>& cells_per_dim_vec,
595 const std::vector<int>& assignRefinedLevel,
596 const std::vector<std::string>& lgr_name_vec,
597 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
598 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
599
601 void postAdapt();
603
604 private:
605 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
606 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
607 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
608 const int& corner_count,
609 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
610 const int& preAdaptMaxLevel,
611 const int& newLevels);
612
613 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
614 const std::vector<std::array<int,3>>& cells_per_dim_vec,
615 const std::vector<int>& lgr_with_at_least_one_active_cell);
616
623 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
624 public:
631 private:
632
645 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
646
650 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
651
652 private:
666 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
667 const std::vector<std::array<int,3>>& endIJK_vec,
668 std::vector<int>& assignRefinedLevel,
669 std::vector<int>& lgr_with_at_least_one_active_cell);
670
672 void populateCellIndexSetRefinedGrid(int level);
673
675 void populateCellIndexSetLeafGridView();
676
678 void populateLeafGlobalIdSet();
679
680 public:
681
686 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
687
689 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
690
692 unsigned int overlapSize(int) const;
693
694
696 unsigned int ghostSize(int) const;
697
699 unsigned int overlapSize(int, int) const;
700
702 unsigned int ghostSize(int, int) const;
703
705 unsigned int numBoundarySegments() const;
706
707 void setPartitioningParams(const std::map<std::string,std::string>& params);
708
709 // loadbalance is not part of the grid interface therefore we skip it.
710
720 bool loadBalance(int overlapLayers=1,
721 int partitionMethod = Dune::PartitionMethod::zoltan,
722 double imbalanceTol = 1.1,
723 int level =-1)
724 {
725 using std::get;
726 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
727 /* ownersFirst = */ false,
728 /* wells = */ nullptr,
729 /* possibleFutureConnections = */ {},
730 /* serialPartitioning = */ false,
731 /* transmissibilities = */ nullptr,
732 /* addCornerCells = */ true,
733 overlapLayers,
734 partitionMethod,
735 imbalanceTol,
736 /* allowDistributedWells = */ false,
737 /* input_cell_part = */ {},
738 level));
739 }
740
741 // loadbalance is not part of the grid interface therefore we skip it.
742
753 bool loadBalanceSerial(int overlapLayers=1,
754 int partitionMethod = Dune::PartitionMethod::zoltan,
755 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
756 double imbalanceTol = 1.1,
757 int level = -1)
758 {
759 using std::get;
760 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
761 /* ownersFirst = */ false,
762 /* wells = */ nullptr,
763 /* possibleFutureConnections = */ {},
764 /* serialPartitioning = */ true,
765 /* transmissibilities = */ nullptr,
766 /* addCornerCells = */ true,
767 overlapLayers,
768 partitionMethod,
769 imbalanceTol,
770 /* allowDistributedWells = */ false,
771 /* input_cell_part = */ {},
772 level));
773 }
774
775 // loadbalance is not part of the grid interface therefore we skip it.
776
805 std::pair<bool,std::vector<std::pair<std::string,bool>>>
806 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
807 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
808 const double* transmissibilities = nullptr,
809 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
810 int level = -1)
811 {
812 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
813 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
814 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
815 /* input_cell_part = */ {}, level);
816 }
817
818 // loadbalance is not part of the grid interface therefore we skip it.
819
853 std::pair<bool,std::vector<std::pair<std::string,bool>>>
854 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
855 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
856 const double* transmissibilities = nullptr, bool ownersFirst=false,
857 bool addCornerCells=false, int overlapLayers=1,
858 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
859 double imbalanceTol = 1.1,
860 int level = -1)
861 {
862 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
863 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
864 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
865 }
866
895 template<class DataHandle>
896 std::pair<bool, std::vector<std::pair<std::string,bool> > >
897 loadBalance(DataHandle& data,
898 const std::vector<cpgrid::OpmWellType> * wells,
899 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
900 const double* transmissibilities = nullptr,
901 int overlapLayers=1, int partitionMethod = 1, int level =-1)
902 {
903 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
904 using std::get;
905 if (get<0>(ret))
906 {
908 }
909 return ret;
910 }
911
946 template<class DataHandle>
947 std::pair<bool, std::vector<std::pair<std::string,bool> > >
948 loadBalance(DataHandle& data, EdgeWeightMethod method,
949 const std::vector<cpgrid::OpmWellType> * wells,
950 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
951 bool serialPartitioning,
952 const double* transmissibilities = nullptr, bool ownersFirst=false,
953 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
954 double imbalanceTol = 1.1,
955 bool allowDistributedWells = false)
956 {
957 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
958 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
959 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
960 using std::get;
961 if (get<0>(ret))
962 {
964 }
965 return ret;
966 }
967
992 template<class DataHandle>
993 std::pair<bool, std::vector<std::pair<std::string,bool> > >
994 loadBalance(DataHandle& data, const std::vector<int>& parts,
995 const std::vector<cpgrid::OpmWellType> * wells,
996 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
997 bool ownersFirst=false,
998 bool addCornerCells=false, int overlapLayers=1)
999 {
1000 using std::get;
1001 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
1002 possibleFutureConnections,
1003 /* serialPartitioning = */ false,
1004 /* transmissibilities = */ {},
1005 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1006 /* imbalanceTol (ignored) = */ 0.0,
1007 /* allowDistributedWells = */ true, parts, /* level = */ 0);
1008 using std::get;
1009 if (get<0>(ret))
1010 {
1012 }
1013 return ret;
1014 }
1015
1023 template<class DataHandle>
1024 bool loadBalance(DataHandle& data,
1025 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1026 {
1027 // decltype usage needed to tell the compiler not to use this function if first
1028 // argument is std::vector but rather loadbalance by parts
1029 bool ret = loadBalance(overlapLayers, partitionMethod);
1030 if (ret)
1031 {
1033 }
1034 return ret;
1035 }
1036
1048 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1049 bool addCornerCells=false, int overlapLayers=1)
1050 {
1051 using std::get;
1052 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1053 {},
1054 /* serialPartitioning = */ false,
1055 /* trabsmissibilities = */ {},
1056 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1057 /* imbalanceTol (ignored) = */ 0.0,
1058 /* allowDistributedWells = */ true, parts));
1059 }
1060
1073 template<class DataHandle>
1074 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1075 bool addCornerCells=false, int overlapLayers=1)
1076 {
1077 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1078 if (ret)
1079 {
1081 }
1082 return ret;
1083 }
1084
1098 std::vector<int>
1099 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1100 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1101 const double* transmissibilities,
1102 const int numParts,
1103 const double imbalanceTol) const;
1104
1112 template<class DataHandle>
1113 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1114 {
1115 communicate(data, iftype, dir);
1116 }
1117
1125 template<class DataHandle>
1126 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1127
1129 const typename CpGridTraits::Communication& comm () const;
1131
1132 // ------------ End of Dune interface, start of simplified interface --------------
1133
1139
1140 // enum { dimension = 3 }; // already defined
1141
1142 typedef Dune::FieldVector<double, 3> Vector;
1143
1144
1145 const std::vector<double>& zcornData() const;
1146
1147
1148 // Topology
1153 int numCells(int level = -1) const;
1154
1159 int numFaces(int level = -1) const;
1160
1162 int numVertices() const;
1163
1164
1173 int numCellFaces(int cell, int level = -1) const;
1174
1181 int cellFace(int cell, int local_index, int level = -1) const;
1182
1186
1206 int faceCell(int face, int local_index, int level = -1) const;
1207
1214 int numCellFaces() const;
1215
1216 int numFaceVertices(int face) const;
1217
1222 int faceVertex(int face, int local_index) const;
1223
1226 double cellCenterDepth(int cell_index) const;
1227
1228
1229 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1230
1231 const Vector faceAreaNormalEcl(int face) const;
1232
1233
1234 // Geometry
1238 const Vector& vertexPosition(int vertex) const;
1239
1242 double faceArea(int face) const;
1243
1246 const Vector& faceCentroid(int face) const;
1247
1251 const Vector& faceNormal(int face) const;
1252
1255 double cellVolume(int cell) const;
1256
1259 const Vector& cellCentroid(int cell) const;
1260
1263 template<int codim>
1265 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1266 FieldVector<double, 3>,
1267 const FieldVector<double, 3>&, int>
1268 {
1269 public:
1271 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1276 : iter_(iter)
1277 {}
1278
1279 const FieldVector<double,3>& dereference() const
1280 {
1281 return iter_->center();
1282 }
1284 {
1285 ++iter_;
1286 }
1287 const FieldVector<double,3>& elementAt(int n)
1288 {
1289 return iter_[n]->center();
1290 }
1291 void advance(int n){
1292 iter_+=n;
1293 }
1295 {
1296 --iter_;
1297 }
1299 {
1300 return o-iter_;
1301 }
1302 bool equals(const CentroidIterator& o) const{
1303 return o==iter_;
1304 }
1305 private:
1307 GeometryIterator iter_;
1308 };
1309
1312
1315
1316 // Extra
1317 int boundaryId(int face) const;
1318
1325 template<class Cell2FacesRowIterator>
1326 int
1327 faceTag(const Cell2FacesRowIterator& cell_face) const;
1328
1330
1331 // ------------ End of simplified interface --------------
1332
1333 //------------- methods not in the DUNE grid interface.
1334
1339
1340
1349 template<class DataHandle>
1350 void scatterData(DataHandle& handle) const;
1351
1358 template<class DataHandle>
1359 void gatherData(DataHandle& handle) const;
1360
1363
1393
1397
1400
1404
1405#if HAVE_MPI
1410
1413
1418
1420
1422
1424
1426#endif
1427
1429 const std::vector<int>& sortedNumAquiferCells() const;
1430
1431 private:
1463 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1464 scatterGrid(EdgeWeightMethod method,
1465 bool ownersFirst,
1466 const std::vector<cpgrid::OpmWellType> * wells,
1467 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1468 bool serialPartitioning,
1469 const double* transmissibilities,
1470 bool addCornerCells,
1471 int overlapLayers,
1472 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1473 double imbalanceTol = 1.1,
1474 bool allowDistributedWells = true,
1475 const std::vector<int>& input_cell_part = {},
1476 int level = -1);
1477
1482 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1484 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1486 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1488 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
1494 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1495 /*
1496 * @brief Interface for scattering and gathering point data.
1497 *
1498 * @warning Will only update owner cells
1499 */
1500 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1504 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1505
1506
1510 std::map<std::string,std::string> partitioningParams;
1511
1512 }; // end Class CpGrid
1513
1514} // end namespace Dune
1515
1519
1520
1521namespace Dune
1522{
1523
1524 namespace Capabilities
1525 {
1527 template <>
1528 struct hasEntity<CpGrid, 0>
1529 {
1530 static const bool v = true;
1531 };
1532
1534 template <>
1535 struct hasEntity<CpGrid, 3>
1536 {
1537 static const bool v = true;
1538 };
1539
1540 template<>
1541 struct canCommunicate<CpGrid,0>
1542 {
1543 static const bool v = true;
1544 };
1545
1546 template<>
1547 struct canCommunicate<CpGrid,3>
1548 {
1549 static const bool v = true;
1550 };
1551
1553 template <>
1554 struct hasBackupRestoreFacilities<CpGrid>
1555 {
1556 static const bool v = false;
1557 };
1558
1559 }
1560
1561 template<class DataHandle>
1562 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1563 {
1564 current_data_->back()->communicate(data, iftype, dir);
1565 }
1566
1567
1568 template<class DataHandle>
1569 void CpGrid::scatterData([[maybe_unused]] DataHandle& handle) const
1570 {
1571#if HAVE_MPI
1572 if (distributed_data_.empty()) {
1573 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
1574 } else {
1575 distributed_data_[0]->scatterData(handle, data_[0].get(),
1576 distributed_data_[0].get(),
1579 }
1580#endif
1581 }
1582
1583 template<class DataHandle>
1584 void CpGrid::gatherData([[maybe_unused]] DataHandle& handle) const
1585 {
1586#if HAVE_MPI
1587 if (distributed_data_.empty()) {
1588 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
1589 } else {
1590 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1591 }
1592#endif
1593 }
1594
1595
1596 template<class Cell2FacesRowIterator>
1597 int
1598 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
1599 {
1600 // Note that this relies on the following implementation detail:
1601 // The grid is always constructed such that the interior faces constructed
1602 // with orientation set to true are
1603 // oriented along the positive IJK direction. Oriented means that
1604 // the first cell attached to face has the lower index.
1605 // For faces along the boundary (only one cell, always attached at index 0)
1606 // the orientation has to be determined by the orientation of the cell.
1607 // If it is true then in UnstructuredGrid it would be stored at index 0,
1608 // otherwise at index 1.
1609 const int cell = cell_face.getCellIndex();
1610 const int face = *cell_face;
1611 assert (0 <= cell); assert (cell < numCells());
1612 assert (0 <= face); assert (face < numFaces());
1613
1615
1616 const cpgrid::EntityRep<1> f(face, true);
1617 const F2C& f2c = current_data_->back()->face_to_cell_[f];
1618 const face_tag tag = current_data_->back()->face_tag_[f];
1619
1620 assert ((f2c.size() == 1) || (f2c.size() == 2));
1621
1622 int inside_cell = 0;
1623
1624 if ( f2c.size() == 2 ) // Two cells => interior
1625 {
1626 if ( f2c[1].index() == cell )
1627 {
1628 inside_cell = 1;
1629 }
1630 }
1631 const bool normal_is_in = ! f2c[inside_cell].orientation();
1632
1633 switch (tag) {
1634 case I_FACE:
1635 // LEFT : RIGHT
1636 return normal_is_in ? 0 : 1; // min(I) : max(I)
1637 case J_FACE:
1638 // BACK : FRONT
1639 return normal_is_in ? 2 : 3; // min(J) : max(J)
1640 case K_FACE:
1641 // Note: TOP at min(K) as 'z' measures *depth*.
1642 // TOP : BOTTOM
1643 return normal_is_in ? 4 : 5; // min(K) : max(K)
1644 case NNC_FACE:
1645 // For nnc faces we return the otherwise unused value -1.
1646 return -1;
1647 default:
1648 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
1649 }
1650 }
1651
1652 template<int dim>
1654
1655} // namespace Dune
1656
1659#include "cpgrid/Intersection.hpp"
1660#include "cpgrid/Geometry.hpp"
1661#include "cpgrid/Indexsets.hpp"
1662
1663#endif // OPM_CPGRID_HEADER
DataHandle & data
Definition: CpGridData.hpp:1166
#define OPM_THROW(Exception, message)
Definition: ErrorMacros.hpp:29
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:1268
void increment()
Definition: CpGrid.hpp:1283
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:1275
const FieldVector< double, 3 > & dereference() const
Definition: CpGrid.hpp:1279
int distanceTo(const CentroidIterator &o)
Definition: CpGrid.hpp:1298
void decrement()
Definition: CpGrid.hpp:1294
void advance(int n)
Definition: CpGrid.hpp:1291
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:1272
bool equals(const CentroidIterator &o) const
Definition: CpGrid.hpp:1302
const FieldVector< double, 3 > & elementAt(int n)
Definition: CpGrid.hpp:1287
[ 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:994
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:1584
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:948
cpgrid::CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGrid.hpp:1407
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:1598
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:854
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:1362
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:897
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:753
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:720
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.
bool mark(int refCount, const cpgrid::Entity< 0 > &element, bool throwOnFailure=false)
Mark entity for refinement (or coarsening).
cpgrid::CpGridDataTraits::RemoteIndices RemoteIndices
The type of the remote indices information.
Definition: CpGrid.hpp:1409
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.
void globalRefine(int refCount, bool throwOnFailure=false)
Refine the grid refCount times using the default refinement rule. This behaves like marking all eleme...
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.
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:1024
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:806
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:1142
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:1074
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 refineAndUpdateGrid(bool throwOnFailure, 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.
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:1048
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 communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
communicate objects for all codims on a given level
Definition: CpGrid.hpp:1113
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:1412
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:1569
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:98
Definition: Geometry.hpp:76
The global id set for Dune.
Definition: Indexsets.hpp:483
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:118
Definition: Indexsets.hpp:57
Definition: Intersection.hpp:276
Definition: Intersection.hpp:63
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