CpGrid.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: CpGrid.hpp
4//
5// Created: Fri May 29 20:26:36 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// B�rd Skaflestad <bard.skaflestad@sintef.no>
9// Antonella Ritorto <antonella.ritorto@opm-op.com>
10//
11// $Date$
12//
13// $Revision$
14//
15//===========================================================================
16
17/*
18 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
19 Copyright 2009, 2010, 2014, 2022-2023 Equinor ASA.
20 Copyright 2014, 2015 Dr. Blatt - HPC-Simulartion-Software & Services
21 Copyright 2015 NTNU
22
23 This file is part of The Open Porous Media project (OPM).
24
25 OPM is free software: you can redistribute it and/or modify
26 it under the terms of the GNU General Public License as published by
27 the Free Software Foundation, either version 3 of the License, or
28 (at your option) any later version.
29
30 OPM is distributed in the hope that it will be useful,
31 but WITHOUT ANY WARRANTY; without even the implied warranty of
32 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
33 GNU General Public License for more details.
34
35 You should have received a copy of the GNU General Public License
36 along with OPM. If not, see <http://www.gnu.org/licenses/>.
37*/
38
39#ifndef OPM_CPGRID_HEADER
40#define OPM_CPGRID_HEADER
41
42// Warning suppression for Dune includes.
44
45#include <dune/grid/common/grid.hh>
46
48
50
54
56
58
59#include <set>
60
61namespace Opm
62{
63struct NNCdata;
64class EclipseGrid;
65class EclipseState;
66}
67
68namespace Dune
69{
70 class CpGrid;
71
72 namespace cpgrid
73 {
74 class CpGridData;
75 template <int> class Entity;
76 template<int,int> class Geometry;
77 class HierarchicIterator;
78 class IntersectionIterator;
79 template<int, PartitionIteratorType> class Iterator;
80 class LevelGlobalIdSet;
81 class GlobalIdSet;
82 class Intersection;
83 class IntersectionIterator;
84 class IndexSet;
85 class IdSet;
86
87 }
88}
89
90namespace Dune
91{
92
94 //
95 // CpGridTraits
96 //
98
100 {
102 typedef CpGrid Grid;
103
112
115
118 template <int cd>
119 struct Codim
120 {
123 typedef cpgrid::Geometry<3-cd, 3> Geometry;
124 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> Geometry;
127 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> LocalGeometry;
130
133
136
139
142 template <PartitionIteratorType pitype>
144 {
149 };
150 };
151
154 template <PartitionIteratorType pitype>
156 {
158 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid> > LevelGridView;
160 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> > LeafGridView;
161
162 };
163
165 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid>> LevelGridView;
167 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>> LeafGridView;
168
177
181 };
182
184 //
185 // CpGridFamily
186 //
188
190 {
192 };
193
195 //
196 // CpGrid
197 //
199
201 class CpGrid
202 : public GridDefaultImplementation<3, 3, double, CpGridFamily>
203 {
204 friend class cpgrid::CpGridData;
205 friend class cpgrid::Entity<0>;
206 friend class cpgrid::Entity<1>;
207 friend class cpgrid::Entity<2>;
208 friend class cpgrid::Entity<3>;
209 template<int dim>
210 friend cpgrid::Entity<dim> createEntity(const CpGrid&,int,bool);
211
212 public:
213
214 // --- Typedefs ---
215
216
219
220
221 // --- Methods ---
222
223
226
227 explicit CpGrid(MPIHelper::MPICommunicator comm);
228
229#if HAVE_ECL_INPUT
266 std::vector<std::size_t>
267 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
268 Opm::EclipseState* ecl_state,
269 bool periodic_extension,
270 bool turn_normals,
271 bool clip_z,
272 bool pinchActive,
273 bool edge_conformal);
274
312 std::vector<std::size_t>
313 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
314 Opm::EclipseState* ecl_state,
315 bool periodic_extension,
316 bool turn_normals = false,
317 bool clip_z = false,
318 bool edge_conformal = false);
319#endif // HAVE_ECL_INPUT
320
335 void processEclipseFormat(const grdecl& input_data,
336 bool remove_ij_boundary,
337 bool turn_normals = false,
338 bool edge_conformal = false);
339
341
347
348
355 void createCartesian(const std::array<int, 3>& dims,
356 const std::array<double, 3>& cellsize,
357 const std::array<int, 3>& shift = {0,0,0});
358
362 const std::array<int, 3>& logicalCartesianSize() const;
363
371 const std::vector<int>& globalCell() const;
372
374 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData() const;
375
377 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData();
378
381
384
392 void getIJK(const int c, std::array<int,3>& ijk) const;
394
398 bool uniqueBoundaryIds() const;
399
402 void setUniqueBoundaryIds(bool uids);
403
404
405 // --- Dune interface below ---
406
408 // \@{
413 std::string name() const;
414
416 int maxLevel() const;
417
419 template<int codim>
420 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
422 template<int codim>
423 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
424
426 template<int codim, PartitionIteratorType PiType>
427 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
429 template<int codim, PartitionIteratorType PiType>
430 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
431
433 template<int codim>
434 typename Traits::template Codim<codim>::LeafIterator leafbegin() const;
436 template<int codim>
437 typename Traits::template Codim<codim>::LeafIterator leafend() const;
438
440 template<int codim, PartitionIteratorType PiType>
441 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const;
443 template<int codim, PartitionIteratorType PiType>
444 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const;
445
447 int size (int level, int codim) const;
448
450 int size (int codim) const;
451
453 int size (int level, GeometryType type) const;
454
456 int size (GeometryType type) const;
457
459 const Traits::GlobalIdSet& globalIdSet() const;
460
462 const Traits::LocalIdSet& localIdSet() const;
463
465 const Traits::LevelIndexSet& levelIndexSet(int level) const;
466
468 const Traits::LeafIndexSet& leafIndexSet() const;
469
474 void globalRefine (int refCount, bool throwOnFailure = false);
475
476 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
477
479 template <int codim>
481
500 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
501 const std::vector<std::array<int,3>>& startIJK_vec,
502 const std::vector<std::array<int,3>>& endIJK_vec,
503 const std::vector<std::string>& lgr_name_vec,
504 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
505
512 void autoRefine(const std::array<int,3>& nxnynz);
513
514 // @brief TO BE DONE
515 const std::map<std::string,int>& getLgrNameToLevel() const;
516
517 // @breif Compute center of an entity/element/cell in the Eclipse way:
518 // - Average of the 4 corners of the bottom face.
519 // - Average of the 4 corners of the top face.
520 // Return average of the previous computations.
521 // @param [in] int Index of a cell.
522 // @return 'eclipse centroid'
523 std::array<double,3> getEclCentroid(const int& idx) const;
524
525 // @breif Compute center of an entity/element/cell in the Eclipse way:
526 // - Average of the 4 corners of the bottom face.
527 // - Average of the 4 corners of the top face.
528 // Return average of the previous computations.
529 // @param [in] Entity<0> Entity
530 // @return 'eclipse centroid'
531 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
532
533 // @brief Return parent (coarse) intersection (face) of a refined face on the leaf grid view, whose neighboring cells
534 // are two: one coarse cell (equivalent to its origin cell from level 0), and one refined cell
535 // from certain LGR.
536 // Used in Transmissibility_impl.hpp
538
559 bool mark(int refCount, const cpgrid::Entity<0>& element, bool throwOnFailure = false);
560
564 int getMark(const cpgrid::Entity<0>& element) const;
565
568 bool preAdapt();
569
572 bool adapt();
573
587 bool refineAndUpdateGrid(const std::vector<std::array<int,3>>& cells_per_dim_vec,
588 const std::vector<int>& assignRefinedLevel,
589 const std::vector<std::string>& lgr_name_vec,
590 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
591 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
592
594 void postAdapt();
596
597 private:
598 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
599 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
600 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
601 const int& corner_count,
602 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
603 const int& preAdaptMaxLevel,
604 const int& newLevels);
605
606 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
607 const std::vector<std::array<int,3>>& cells_per_dim_vec,
608 const std::vector<int>& lgr_with_at_least_one_active_cell);
609
616 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
617 public:
624 private:
625
638 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
639
643 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
644
645 private:
659 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
660 const std::vector<std::array<int,3>>& endIJK_vec,
661 std::vector<int>& assignRefinedLevel,
662 std::vector<int>& lgr_with_at_least_one_active_cell);
663
665 void populateCellIndexSetRefinedGrid(int level);
666
668 void populateCellIndexSetLeafGridView();
669
671 void populateLeafGlobalIdSet();
672
673 public:
674
679 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
680
682 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
683
685 unsigned int overlapSize(int) const;
686
687
689 unsigned int ghostSize(int) const;
690
692 unsigned int overlapSize(int, int) const;
693
695 unsigned int ghostSize(int, int) const;
696
698 unsigned int numBoundarySegments() const;
699
700 void setPartitioningParams(const std::map<std::string,std::string>& params);
701
702 // loadbalance is not part of the grid interface therefore we skip it.
703
713 bool loadBalance(int overlapLayers=1,
714 int partitionMethod = Dune::PartitionMethod::zoltan,
715 double imbalanceTol = 1.1,
716 int level =-1)
717 {
718 using std::get;
719 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
720 /* ownersFirst = */ false,
721 /* wells = */ nullptr,
722 /* possibleFutureConnections = */ {},
723 /* serialPartitioning = */ false,
724 /* transmissibilities = */ nullptr,
725 /* addCornerCells = */ true,
726 overlapLayers,
727 partitionMethod,
728 imbalanceTol,
729 /* allowDistributedWells = */ false,
730 /* input_cell_part = */ {},
731 level));
732 }
733
734 // loadbalance is not part of the grid interface therefore we skip it.
735
746 bool loadBalanceSerial(int overlapLayers=1,
747 int partitionMethod = Dune::PartitionMethod::zoltan,
748 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
749 double imbalanceTol = 1.1,
750 int level = -1)
751 {
752 using std::get;
753 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
754 /* ownersFirst = */ false,
755 /* wells = */ nullptr,
756 /* possibleFutureConnections = */ {},
757 /* serialPartitioning = */ true,
758 /* transmissibilities = */ nullptr,
759 /* addCornerCells = */ true,
760 overlapLayers,
761 partitionMethod,
762 imbalanceTol,
763 /* allowDistributedWells = */ false,
764 /* input_cell_part = */ {},
765 level));
766 }
767
768 // loadbalance is not part of the grid interface therefore we skip it.
769
798 std::pair<bool,std::vector<std::pair<std::string,bool>>>
799 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
800 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
801 const double* transmissibilities = nullptr,
802 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
803 int level = -1)
804 {
805 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
806 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
807 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
808 /* input_cell_part = */ {}, level);
809 }
810
811 // loadbalance is not part of the grid interface therefore we skip it.
812
846 std::pair<bool,std::vector<std::pair<std::string,bool>>>
847 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
848 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
849 const double* transmissibilities = nullptr, bool ownersFirst=false,
850 bool addCornerCells=false, int overlapLayers=1,
851 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
852 double imbalanceTol = 1.1,
853 int level = -1)
854 {
855 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
856 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
857 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
858 }
859
888 template<class DataHandle>
889 std::pair<bool, std::vector<std::pair<std::string,bool> > >
890 loadBalance(DataHandle& data,
891 const std::vector<cpgrid::OpmWellType> * wells,
892 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
893 const double* transmissibilities = nullptr,
894 int overlapLayers=1, int partitionMethod = 1, int level =-1)
895 {
896 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
897 using std::get;
898 if (get<0>(ret))
899 {
901 }
902 return ret;
903 }
904
939 template<class DataHandle>
940 std::pair<bool, std::vector<std::pair<std::string,bool> > >
941 loadBalance(DataHandle& data, EdgeWeightMethod method,
942 const std::vector<cpgrid::OpmWellType> * wells,
943 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
944 bool serialPartitioning,
945 const double* transmissibilities = nullptr, bool ownersFirst=false,
946 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
947 double imbalanceTol = 1.1,
948 bool allowDistributedWells = false)
949 {
950 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
951 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
952 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
953 using std::get;
954 if (get<0>(ret))
955 {
957 }
958 return ret;
959 }
960
985 template<class DataHandle>
986 std::pair<bool, std::vector<std::pair<std::string,bool> > >
987 loadBalance(DataHandle& data, const std::vector<int>& parts,
988 const std::vector<cpgrid::OpmWellType> * wells,
989 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
990 bool ownersFirst=false,
991 bool addCornerCells=false, int overlapLayers=1)
992 {
993 using std::get;
994 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
995 possibleFutureConnections,
996 /* serialPartitioning = */ false,
997 /* transmissibilities = */ {},
998 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
999 /* imbalanceTol (ignored) = */ 0.0,
1000 /* allowDistributedWells = */ true, parts, /* level = */ 0);
1001 using std::get;
1002 if (get<0>(ret))
1003 {
1005 }
1006 return ret;
1007 }
1008
1016 template<class DataHandle>
1017 bool loadBalance(DataHandle& data,
1018 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1019 {
1020 // decltype usage needed to tell the compiler not to use this function if first
1021 // argument is std::vector but rather loadbalance by parts
1022 bool ret = loadBalance(overlapLayers, partitionMethod);
1023 if (ret)
1024 {
1026 }
1027 return ret;
1028 }
1029
1041 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1042 bool addCornerCells=false, int overlapLayers=1)
1043 {
1044 using std::get;
1045 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1046 {},
1047 /* serialPartitioning = */ false,
1048 /* trabsmissibilities = */ {},
1049 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1050 /* imbalanceTol (ignored) = */ 0.0,
1051 /* allowDistributedWells = */ true, parts));
1052 }
1053
1066 template<class DataHandle>
1067 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1068 bool addCornerCells=false, int overlapLayers=1)
1069 {
1070 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1071 if (ret)
1072 {
1074 }
1075 return ret;
1076 }
1077
1091 std::vector<int>
1092 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1093 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1094 const double* transmissibilities,
1095 const int numParts,
1096 const double imbalanceTol) const;
1097
1105 template<class DataHandle>
1106 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1107 {
1108 communicate(data, iftype, dir);
1109 }
1110
1118 template<class DataHandle>
1119 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1120
1122 const typename CpGridTraits::Communication& comm () const;
1124
1125 // ------------ End of Dune interface, start of simplified interface --------------
1126
1132
1133 // enum { dimension = 3 }; // already defined
1134
1135 typedef Dune::FieldVector<double, 3> Vector;
1136
1137
1138 const std::vector<double>& zcornData() const;
1139
1140
1141 // Topology
1146 int numCells(int level = -1) const;
1147
1152 int numFaces(int level = -1) const;
1153
1155 int numVertices() const;
1156
1157
1166 int numCellFaces(int cell, int level = -1) const;
1167
1174 int cellFace(int cell, int local_index, int level = -1) const;
1175
1179
1199 int faceCell(int face, int local_index, int level = -1) const;
1200
1207 int numCellFaces() const;
1208
1209 int numFaceVertices(int face) const;
1210
1215 int faceVertex(int face, int local_index) const;
1216
1219 double cellCenterDepth(int cell_index) const;
1220
1221
1222 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1223
1224 const Vector faceAreaNormalEcl(int face) const;
1225
1226
1227 // Geometry
1231 const Vector& vertexPosition(int vertex) const;
1232
1235 double faceArea(int face) const;
1236
1239 const Vector& faceCentroid(int face) const;
1240
1244 const Vector& faceNormal(int face) const;
1245
1248 double cellVolume(int cell) const;
1249
1252 const Vector& cellCentroid(int cell) const;
1253
1256 template<int codim>
1258 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1259 FieldVector<double, 3>,
1260 const FieldVector<double, 3>&, int>
1261 {
1262 public:
1264 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1269 : iter_(iter)
1270 {}
1271
1272 const FieldVector<double,3>& dereference() const
1273 {
1274 return iter_->center();
1275 }
1277 {
1278 ++iter_;
1279 }
1280 const FieldVector<double,3>& elementAt(int n)
1281 {
1282 return iter_[n]->center();
1283 }
1284 void advance(int n){
1285 iter_+=n;
1286 }
1288 {
1289 --iter_;
1290 }
1292 {
1293 return o-iter_;
1294 }
1295 bool equals(const CentroidIterator& o) const{
1296 return o==iter_;
1297 }
1298 private:
1300 GeometryIterator iter_;
1301 };
1302
1305
1308
1309 // Extra
1310 int boundaryId(int face) const;
1311
1318 template<class Cell2FacesRowIterator>
1319 int
1320 faceTag(const Cell2FacesRowIterator& cell_face) const;
1321
1323
1324 // ------------ End of simplified interface --------------
1325
1326 //------------- methods not in the DUNE grid interface.
1327
1332
1333
1342 template<class DataHandle>
1343 void scatterData(DataHandle& handle) const;
1344
1351 template<class DataHandle>
1352 void gatherData(DataHandle& handle) const;
1353
1356
1386
1390
1393
1397
1398#if HAVE_MPI
1403
1406
1411
1413
1415
1417
1419#endif
1420
1422 const std::vector<int>& sortedNumAquiferCells() const;
1423
1424 private:
1456 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1457 scatterGrid(EdgeWeightMethod method,
1458 bool ownersFirst,
1459 const std::vector<cpgrid::OpmWellType> * wells,
1460 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1461 bool serialPartitioning,
1462 const double* transmissibilities,
1463 bool addCornerCells,
1464 int overlapLayers,
1465 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1466 double imbalanceTol = 1.1,
1467 bool allowDistributedWells = true,
1468 const std::vector<int>& input_cell_part = {},
1469 int level = -1);
1470
1475 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1477 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1479 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1481 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
1487 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1488 /*
1489 * @brief Interface for scattering and gathering point data.
1490 *
1491 * @warning Will only update owner cells
1492 */
1493 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1497 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1498
1499
1503 std::map<std::string,std::string> partitioningParams;
1504
1505 }; // end Class CpGrid
1506
1507} // end namespace Dune
1508
1512
1513
1514namespace Dune
1515{
1516
1517 namespace Capabilities
1518 {
1520 template <>
1521 struct hasEntity<CpGrid, 0>
1522 {
1523 static const bool v = true;
1524 };
1525
1527 template <>
1528 struct hasEntity<CpGrid, 3>
1529 {
1530 static const bool v = true;
1531 };
1532
1533 template<>
1534 struct canCommunicate<CpGrid,0>
1535 {
1536 static const bool v = true;
1537 };
1538
1539 template<>
1540 struct canCommunicate<CpGrid,3>
1541 {
1542 static const bool v = true;
1543 };
1544
1546 template <>
1547 struct hasBackupRestoreFacilities<CpGrid>
1548 {
1549 static const bool v = false;
1550 };
1551
1552 }
1553
1554 template<class DataHandle>
1555 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1556 {
1557 current_data_->back()->communicate(data, iftype, dir);
1558 }
1559
1560
1561 template<class DataHandle>
1562 void CpGrid::scatterData([[maybe_unused]] DataHandle& handle) const
1563 {
1564#if HAVE_MPI
1565 if (distributed_data_.empty()) {
1566 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
1567 } else {
1568 distributed_data_[0]->scatterData(handle, data_[0].get(),
1569 distributed_data_[0].get(),
1572 }
1573#endif
1574 }
1575
1576 template<class DataHandle>
1577 void CpGrid::gatherData([[maybe_unused]] DataHandle& handle) const
1578 {
1579#if HAVE_MPI
1580 if (distributed_data_.empty()) {
1581 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
1582 } else {
1583 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1584 }
1585#endif
1586 }
1587
1588
1589 template<class Cell2FacesRowIterator>
1590 int
1591 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
1592 {
1593 // Note that this relies on the following implementation detail:
1594 // The grid is always constructed such that the interior faces constructed
1595 // with orientation set to true are
1596 // oriented along the positive IJK direction. Oriented means that
1597 // the first cell attached to face has the lower index.
1598 // For faces along the boundary (only one cell, always attached at index 0)
1599 // the orientation has to be determined by the orientation of the cell.
1600 // If it is true then in UnstructuredGrid it would be stored at index 0,
1601 // otherwise at index 1.
1602 const int cell = cell_face.getCellIndex();
1603 const int face = *cell_face;
1604 assert (0 <= cell); assert (cell < numCells());
1605 assert (0 <= face); assert (face < numFaces());
1606
1608
1609 const cpgrid::EntityRep<1> f(face, true);
1610 const F2C& f2c = current_data_->back()->face_to_cell_[f];
1611 const face_tag tag = current_data_->back()->face_tag_[f];
1612
1613 assert ((f2c.size() == 1) || (f2c.size() == 2));
1614
1615 int inside_cell = 0;
1616
1617 if ( f2c.size() == 2 ) // Two cells => interior
1618 {
1619 if ( f2c[1].index() == cell )
1620 {
1621 inside_cell = 1;
1622 }
1623 }
1624 const bool normal_is_in = ! f2c[inside_cell].orientation();
1625
1626 switch (tag) {
1627 case I_FACE:
1628 // LEFT : RIGHT
1629 return normal_is_in ? 0 : 1; // min(I) : max(I)
1630 case J_FACE:
1631 // BACK : FRONT
1632 return normal_is_in ? 2 : 3; // min(J) : max(J)
1633 case K_FACE:
1634 // Note: TOP at min(K) as 'z' measures *depth*.
1635 // TOP : BOTTOM
1636 return normal_is_in ? 4 : 5; // min(K) : max(K)
1637 case NNC_FACE:
1638 // For nnc faces we return the otherwise unused value -1.
1639 return -1;
1640 default:
1641 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
1642 }
1643 }
1644
1645 template<int dim>
1647
1648} // namespace Dune
1649
1652#include "cpgrid/Intersection.hpp"
1653#include "cpgrid/Geometry.hpp"
1654#include "cpgrid/Indexsets.hpp"
1655
1656#endif // OPM_CPGRID_HEADER
DataHandle & data
Definition: CpGridData.hpp:1166
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:1261
void increment()
Definition: CpGrid.hpp:1276
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:1268
const FieldVector< double, 3 > & dereference() const
Definition: CpGrid.hpp:1272
int distanceTo(const CentroidIterator &o)
Definition: CpGrid.hpp:1291
void decrement()
Definition: CpGrid.hpp:1287
void advance(int n)
Definition: CpGrid.hpp:1284
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:1265
bool equals(const CentroidIterator &o) const
Definition: CpGrid.hpp:1295
const FieldVector< double, 3 > & elementAt(int n)
Definition: CpGrid.hpp:1280
[ 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:987
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:1577
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:941
cpgrid::CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGrid.hpp:1400
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:1591
const std::vector< double > & zcornData() const
ParallelIndexSet & getCellIndexSet()
void setUniqueBoundaryIds(bool uids)
std::array< double, 3 > getEclCentroid(const cpgrid::Entity< 0 > &elem) const
bool refineAndUpdateGrid(const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< int > &assignRefinedLevel, const std::vector< std::string > &lgr_name_vec, const std::vector< std::array< int, 3 > > &startIJK_vec=std::vector< std::array< int, 3 > >{}, const std::vector< std::array< int, 3 > > &endIJK_vec=std::vector< std::array< int, 3 > >{})
Triggers the grid refinement process, allowing to select diffrent refined level grids.
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, double imbalanceTol=1.1, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:847
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:1355
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:890
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:746
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:713
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:1402
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:1017
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:799
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:1135
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:1067
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:1041
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:1106
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:1405
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:1562
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:118
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:99
Definition: Geometry.hpp:76
The global id set for Dune.
Definition: Indexsets.hpp:487
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:118
Definition: Indexsets.hpp:56
Definition: Intersection.hpp:279
Definition: Intersection.hpp:66
Definition: Iterators.hpp:60
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
@ simple
Use simple approach based on rectangular partitioning the underlying cartesian grid.
Definition: GridEnums.hpp:46
@ zoltanGoG
use Zoltan on GraphOfGrid for partitioning
Definition: GridEnums.hpp:52
@ zoltan
Use Zoltan for partitioning.
Definition: GridEnums.hpp:48
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Me...
Definition: GridEnums.hpp:34
@ defaultTransEdgeWgt
Use the transmissibilities as edge weights.
Definition: GridEnums.hpp:38
void getFirstChildGlobalIds(const Dune::CpGrid &grid, std::vector< int > &parentToFirstChildGlobalIds)
Retrieves the global ids of the first child for each parent cell in the grid.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
face_tag
Definition: preprocess.h:66
@ K_FACE
Definition: preprocess.h:69
@ J_FACE
Definition: preprocess.h:68
@ NNC_FACE
Definition: preprocess.h:70
@ I_FACE
Definition: preprocess.h:67
Definition: CpGrid.hpp:190
CpGridTraits Traits
Definition: CpGrid.hpp:191
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:144
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition: CpGrid.hpp:146
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition: CpGrid.hpp:148
Traits associated with a specific codim.
Definition: CpGrid.hpp:120
cpgrid::Entity< cd > Entity
The type of the entity.
Definition: CpGrid.hpp:129
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity. IMPORTANT: Codim<codim>::Geometry == Geometry<di...
Definition: CpGrid.hpp:123
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition: CpGrid.hpp:126
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition: CpGrid.hpp:135
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition: CpGrid.hpp:132
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition: CpGrid.hpp:138
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:156
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:160
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:158
Definition: CpGrid.hpp:100
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition: CpGrid.hpp:170
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition: CpGrid.hpp:109
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:167
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:165
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition: CpGrid.hpp:174
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition: CpGrid.hpp:111
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition: CpGrid.hpp:172
cpgrid::CpGridDataTraits::CollectiveCommunication CollectiveCommunication
Definition: CpGrid.hpp:180
GlobalIdSet LocalIdSet
The type of the local id set.
Definition: CpGrid.hpp:176
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGrid.hpp:179
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition: CpGrid.hpp:114
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition: CpGrid.hpp:107
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition: CpGrid.hpp:105
CpGrid Grid
The type that implements the grid.
Definition: CpGrid.hpp:102
Dune::RemoteIndices< ParallelIndexSet > RemoteIndices
The type of the remote indices information.
Definition: CpGridDataTraits.hpp:83
typename CommunicationType::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGridDataTraits.hpp:80
Dune::Communication< MPICommunicator > CollectiveCommunication
Definition: CpGridDataTraits.hpp:59
Dune::OwnerOverlapCopyCommunication< int, int > CommunicationType
type of OwnerOverlap communication for cells
Definition: CpGridDataTraits.hpp:77
Dune::Communication< MPICommunicator > Communication
Definition: CpGridDataTraits.hpp:58
Communicator::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGridDataTraits.hpp:74
Definition: preprocess.h:56