opm-grid
Loading...
Searching...
No Matches
CpGrid.hpp
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.
43#include <opm/grid/utility/platform_dependent/disable_warnings.h>
44
45#include <dune/grid/common/grid.hh>
46
47#include <opm/grid/utility/platform_dependent/reenable_warnings.h>
48
49#include <opm/grid/common/GridEnums.hpp>
50
51#include <opm/grid/cpgrid/CpGridDataTraits.hpp>
52#include <opm/grid/cpgrid/DefaultGeometryPolicy.hpp>
53#include <opm/grid/cpgrid/OrientedEntityTable.hpp>
54
56
57#include <opm/grid/utility/OpmWellType.hpp>
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;
79 template<int, PartitionIteratorType> class Iterator;
80 class LevelGlobalIdSet;
81 class GlobalIdSet;
82 class Intersection;
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>
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
179 using Communication = cpgrid::CpGridDataTraits::Communication;
180 using CollectiveCommunication = cpgrid::CpGridDataTraits::CollectiveCommunication;
181 };
182
184 //
185 // CpGridFamily
186 //
188
190 {
191 typedef CpGridTraits Traits;
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
225 CpGrid();
226
227 explicit CpGrid(MPIHelper::MPICommunicator comm);
228
229#if HAVE_ECL_INPUT
266 std::vector<std::size_t>
267 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
268 Opm::EclipseState* ecl_state,
269 bool periodic_extension,
270 bool turn_normals,
271 bool clip_z,
272 bool pinchActive,
273 bool edge_conformal);
274
312 std::vector<std::size_t>
313 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
314 Opm::EclipseState* ecl_state,
315 bool periodic_extension,
316 bool turn_normals = false,
317 bool clip_z = false,
318 bool edge_conformal = false);
319#endif // HAVE_ECL_INPUT
320
335 void processEclipseFormat(const grdecl& input_data,
336 bool remove_ij_boundary,
337 bool turn_normals = false,
338 bool edge_conformal = false);
339
341
347
348
355 void createCartesian(const std::array<int, 3>& dims,
356 const std::array<double, 3>& cellsize,
357 const std::array<int, 3>& shift = {0,0,0});
358
362 const std::array<int, 3>& logicalCartesianSize() const;
363
371 const std::vector<int>& globalCell() const;
372
374 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData() const;
375
377 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData();
378
386 void getIJK(const int c, std::array<int,3>& ijk) const;
388
392 bool uniqueBoundaryIds() const;
393
396 void setUniqueBoundaryIds(bool uids);
397
398
399 // --- Dune interface below ---
400
402 // \@{
407 std::string name() const;
408
410 int maxLevel() const;
411
413 template<int codim>
414 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
416 template<int codim>
417 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
418
420 template<int codim, PartitionIteratorType PiType>
421 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
423 template<int codim, PartitionIteratorType PiType>
424 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
425
427 template<int codim>
428 typename Traits::template Codim<codim>::LeafIterator leafbegin() const;
430 template<int codim>
431 typename Traits::template Codim<codim>::LeafIterator leafend() const;
432
434 template<int codim, PartitionIteratorType PiType>
435 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const;
437 template<int codim, PartitionIteratorType PiType>
438 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const;
439
441 int size (int level, int codim) const;
442
444 int size (int codim) const;
445
447 int size (int level, GeometryType type) const;
448
450 int size (GeometryType type) const;
451
453 const Traits::GlobalIdSet& globalIdSet() const;
454
456 const Traits::LocalIdSet& localIdSet() const;
457
459 const Traits::LevelIndexSet& levelIndexSet(int level) const;
460
462 const Traits::LeafIndexSet& leafIndexSet() const;
463
467 void globalRefine (int refCount);
468
469 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
470
472 template <int codim>
474
493 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
494 const std::vector<std::array<int,3>>& startIJK_vec,
495 const std::vector<std::array<int,3>>& endIJK_vec,
496 const std::vector<std::string>& lgr_name_vec,
497 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
498
505 void autoRefine(const std::array<int,3>& nxnynz);
506
507 // @brief TO BE DONE
508 const std::map<std::string,int>& getLgrNameToLevel() const;
509
510 // @breif Compute center of an entity/element/cell in the Eclipse way:
511 // - Average of the 4 corners of the bottom face.
512 // - Average of the 4 corners of the top face.
513 // Return average of the previous computations.
514 // @param [in] int Index of a cell.
515 // @return 'eclipse centroid'
516 std::array<double,3> getEclCentroid(const int& idx) 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] Entity<0> Entity
523 // @return 'eclipse centroid'
524 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
525
526 // @brief Return parent (coarse) intersection (face) of a refined face on the leaf grid view, whose neighboring cells
527 // are two: one coarse cell (equivalent to its origin cell from level 0), and one refined cell
528 // from certain LGR.
529 // Used in Transmissibility_impl.hpp
531
551 bool mark(int refCount, const cpgrid::Entity<0>& element);
552
556 int getMark(const cpgrid::Entity<0>& element) const;
557
560 bool preAdapt();
561
564 bool adapt();
565
579 bool refineAndUpdateGrid(const std::vector<std::array<int,3>>& cells_per_dim_vec,
580 const std::vector<int>& assignRefinedLevel,
581 const std::vector<std::string>& lgr_name_vec,
582 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
583 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
584
586 void postAdapt();
588
589 private:
590 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
591 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
592 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
593 const int& corner_count,
594 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
595 const int& preAdaptMaxLevel,
596 const int& newLevels);
597
598 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
599 const std::vector<std::array<int,3>>& cells_per_dim_vec,
600 const std::vector<int>& lgr_with_at_least_one_active_cell);
601
608 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
609 public:
616 private:
617
630 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
631
635 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
636
637 private:
651 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
652 const std::vector<std::array<int,3>>& endIJK_vec,
653 std::vector<int>& assignRefinedLevel,
654 std::vector<int>& lgr_with_at_least_one_active_cell);
655
657 void populateCellIndexSetRefinedGrid(int level);
658
660 void populateCellIndexSetLeafGridView();
661
663 void populateLeafGlobalIdSet();
664
665 public:
666
671 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
672
674 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
675
677 unsigned int overlapSize(int) const;
678
679
681 unsigned int ghostSize(int) const;
682
684 unsigned int overlapSize(int, int) const;
685
687 unsigned int ghostSize(int, int) const;
688
690 unsigned int numBoundarySegments() const;
691
692 void setPartitioningParams(const std::map<std::string,std::string>& params);
693
694 // loadbalance is not part of the grid interface therefore we skip it.
695
705 bool loadBalance(int overlapLayers=1,
706 int partitionMethod = Dune::PartitionMethod::zoltan,
707 double imbalanceTol = 1.1,
708 int level =-1)
709 {
710 using std::get;
711 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
712 /* ownersFirst = */ false,
713 /* wells = */ nullptr,
714 /* possibleFutureConnections = */ {},
715 /* serialPartitioning = */ false,
716 /* transmissibilities = */ nullptr,
717 /* addCornerCells = */ true,
718 overlapLayers,
719 partitionMethod,
720 imbalanceTol,
721 /* allowDistributedWells = */ false,
722 /* input_cell_part = */ {},
723 level));
724 }
725
726 // loadbalance is not part of the grid interface therefore we skip it.
727
738 bool loadBalanceSerial(int overlapLayers=1,
739 int partitionMethod = Dune::PartitionMethod::zoltan,
740 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
741 double imbalanceTol = 1.1,
742 int level = -1)
743 {
744 using std::get;
745 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
746 /* ownersFirst = */ false,
747 /* wells = */ nullptr,
748 /* possibleFutureConnections = */ {},
749 /* serialPartitioning = */ true,
750 /* transmissibilities = */ nullptr,
751 /* addCornerCells = */ true,
752 overlapLayers,
753 partitionMethod,
754 imbalanceTol,
755 /* allowDistributedWells = */ false,
756 /* input_cell_part = */ {},
757 level));
758 }
759
760 // loadbalance is not part of the grid interface therefore we skip it.
761
790 std::pair<bool,std::vector<std::pair<std::string,bool>>>
791 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
792 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
793 const double* transmissibilities = nullptr,
794 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
795 int level = -1)
796 {
797 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
798 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
799 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
800 /* input_cell_part = */ {}, level);
801 }
802
803 // loadbalance is not part of the grid interface therefore we skip it.
804
838 std::pair<bool,std::vector<std::pair<std::string,bool>>>
839 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
840 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
841 const double* transmissibilities = nullptr, bool ownersFirst=false,
842 bool addCornerCells=false, int overlapLayers=1,
843 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
844 double imbalanceTol = 1.1,
845 int level = -1)
846 {
847 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
848 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
849 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
850 }
851
880 template<class DataHandle>
881 std::pair<bool, std::vector<std::pair<std::string,bool> > >
882 loadBalance(DataHandle& data,
883 const std::vector<cpgrid::OpmWellType> * wells,
884 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
885 const double* transmissibilities = nullptr,
886 int overlapLayers=1, int partitionMethod = 1, int level =-1)
887 {
888 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
889 using std::get;
890 if (get<0>(ret))
891 {
892 scatterData(data);
893 }
894 return ret;
895 }
896
931 template<class DataHandle>
932 std::pair<bool, std::vector<std::pair<std::string,bool> > >
933 loadBalance(DataHandle& data, EdgeWeightMethod method,
934 const std::vector<cpgrid::OpmWellType> * wells,
935 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
936 bool serialPartitioning,
937 const double* transmissibilities = nullptr, bool ownersFirst=false,
938 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
939 double imbalanceTol = 1.1,
940 bool allowDistributedWells = false)
941 {
942 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
943 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
944 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
945 using std::get;
946 if (get<0>(ret))
947 {
948 scatterData(data);
949 }
950 return ret;
951 }
952
977 template<class DataHandle>
978 std::pair<bool, std::vector<std::pair<std::string,bool> > >
979 loadBalance(DataHandle& data, const std::vector<int>& parts,
980 const std::vector<cpgrid::OpmWellType> * wells,
981 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
982 bool ownersFirst=false,
983 bool addCornerCells=false, int overlapLayers=1)
984 {
985 using std::get;
986 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
987 possibleFutureConnections,
988 /* serialPartitioning = */ false,
989 /* transmissibilities = */ {},
990 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
991 /* imbalanceTol (ignored) = */ 0.0,
992 /* allowDistributedWells = */ true, parts, /* level = */ 0);
993 using std::get;
994 if (get<0>(ret))
995 {
996 scatterData(data);
997 }
998 return ret;
999 }
1000
1008 template<class DataHandle>
1009 bool loadBalance(DataHandle& data,
1010 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1011 {
1012 // decltype usage needed to tell the compiler not to use this function if first
1013 // argument is std::vector but rather loadbalance by parts
1014 bool ret = loadBalance(overlapLayers, partitionMethod);
1015 if (ret)
1016 {
1017 scatterData(data);
1018 }
1019 return ret;
1020 }
1021
1033 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1034 bool addCornerCells=false, int overlapLayers=1)
1035 {
1036 using std::get;
1037 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1038 {},
1039 /* serialPartitioning = */ false,
1040 /* trabsmissibilities = */ {},
1041 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1042 /* imbalanceTol (ignored) = */ 0.0,
1043 /* allowDistributedWells = */ true, parts));
1044 }
1045
1058 template<class DataHandle>
1059 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1060 bool addCornerCells=false, int overlapLayers=1)
1061 {
1062 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1063 if (ret)
1064 {
1065 scatterData(data);
1066 }
1067 return ret;
1068 }
1069
1083 std::vector<int>
1084 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1085 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1086 const double* transmissibilities,
1087 const int numParts,
1088 const double imbalanceTol) const;
1089
1097 template<class DataHandle>
1098 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1099 {
1100 communicate(data, iftype, dir);
1101 }
1102
1110 template<class DataHandle>
1111 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1112
1114 const typename CpGridTraits::Communication& comm () const;
1116
1117 // ------------ End of Dune interface, start of simplified interface --------------
1118
1124
1125 // enum { dimension = 3 }; // already defined
1126
1127 typedef Dune::FieldVector<double, 3> Vector;
1128
1129
1130 const std::vector<double>& zcornData() const;
1131
1132
1133 // Topology
1138 int numCells(int level = -1) const;
1139
1144 int numFaces(int level = -1) const;
1145
1147 int numVertices() const;
1148
1149
1158 int numCellFaces(int cell, int level = -1) const;
1159
1166 int cellFace(int cell, int local_index, int level = -1) const;
1167
1170 const cpgrid::OrientedEntityTable<0,1>::row_type cellFaceRow(int cell) const;
1171
1191 int faceCell(int face, int local_index, int level = -1) const;
1192
1199 int numCellFaces() const;
1200
1201 int numFaceVertices(int face) const;
1202
1207 int faceVertex(int face, int local_index) const;
1208
1211 double cellCenterDepth(int cell_index) const;
1212
1213
1214 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1215
1216 const Vector faceAreaNormalEcl(int face) const;
1217
1218
1219 // Geometry
1223 const Vector& vertexPosition(int vertex) const;
1224
1227 double faceArea(int face) const;
1228
1231 const Vector& faceCentroid(int face) const;
1232
1236 const Vector& faceNormal(int face) const;
1237
1240 double cellVolume(int cell) const;
1241
1244 const Vector& cellCentroid(int cell) const;
1245
1248 template<int codim>
1250 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1251 FieldVector<double, 3>,
1252 const FieldVector<double, 3>&, int>
1253 {
1254 public:
1256 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1261 : iter_(iter)
1262 {}
1263
1264 const FieldVector<double,3>& dereference() const
1265 {
1266 return iter_->center();
1267 }
1268 void increment()
1269 {
1270 ++iter_;
1271 }
1272 const FieldVector<double,3>& elementAt(int n)
1273 {
1274 return iter_[n]->center();
1275 }
1276 void advance(int n){
1277 iter_+=n;
1278 }
1279 void decrement()
1280 {
1281 --iter_;
1282 }
1283 int distanceTo(const CentroidIterator& o)
1284 {
1285 return o-iter_;
1286 }
1287 bool equals(const CentroidIterator& o) const{
1288 return o==iter_;
1289 }
1290 private:
1292 GeometryIterator iter_;
1293 };
1294
1297
1300
1301 // Extra
1302 int boundaryId(int face) const;
1303
1310 template<class Cell2FacesRowIterator>
1311 int
1312 faceTag(const Cell2FacesRowIterator& cell_face) const;
1313
1315
1316 // ------------ End of simplified interface --------------
1317
1318 //------------- methods not in the DUNE grid interface.
1319
1324
1325
1334 template<class DataHandle>
1335 void scatterData(DataHandle& handle) const;
1336
1343 template<class DataHandle>
1344 void gatherData(DataHandle& handle) const;
1345
1347 using InterfaceMap = cpgrid::CpGridDataTraits::InterfaceMap;
1348
1378
1382
1384 void switchToGlobalView();
1385
1389
1390#if HAVE_MPI
1392 using ParallelIndexSet = cpgrid::CpGridDataTraits::ParallelIndexSet;
1394 using RemoteIndices = cpgrid::CpGridDataTraits::RemoteIndices;
1395
1397 using CommunicationType = cpgrid::CpGridDataTraits::CommunicationType;
1398
1402 const CommunicationType& cellCommunication() const;
1403
1404 ParallelIndexSet& getCellIndexSet();
1405
1406 RemoteIndices& getCellRemoteIndices();
1407
1408 const ParallelIndexSet& getCellIndexSet() const;
1409
1410 const RemoteIndices& getCellRemoteIndices() const;
1411#endif
1412
1414 const std::vector<int>& sortedNumAquiferCells() const;
1415
1416 private:
1448 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1449 scatterGrid(EdgeWeightMethod method,
1450 bool ownersFirst,
1451 const std::vector<cpgrid::OpmWellType> * wells,
1452 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1453 bool serialPartitioning,
1454 const double* transmissibilities,
1455 bool addCornerCells,
1456 int overlapLayers,
1457 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1458 double imbalanceTol = 1.1,
1459 bool allowDistributedWells = true,
1460 const std::vector<int>& input_cell_part = {},
1461 int level = -1);
1462
1467 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1469 cpgrid::CpGridData* current_view_data_;
1471 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1473 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1475 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
1481 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1482 /*
1483 * @brief Interface for scattering and gathering point data.
1484 *
1485 * @warning Will only update owner cells
1486 */
1487 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1491 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1492
1493
1497 std::map<std::string,std::string> partitioningParams;
1498
1499 }; // end Class CpGrid
1500
1501} // end namespace Dune
1502
1503#include <opm/grid/cpgrid/Entity.hpp>
1504#include <opm/grid/cpgrid/Iterators.hpp>
1505#include <opm/grid/cpgrid/CpGridData.hpp>
1506
1507
1508namespace Dune
1509{
1510
1511 namespace Capabilities
1512 {
1514 template <>
1515 struct hasEntity<CpGrid, 0>
1516 {
1517 static const bool v = true;
1518 };
1519
1521 template <>
1522 struct hasEntity<CpGrid, 3>
1523 {
1524 static const bool v = true;
1525 };
1526
1527 template<>
1528 struct canCommunicate<CpGrid,0>
1529 {
1530 static const bool v = true;
1531 };
1532
1533 template<>
1534 struct canCommunicate<CpGrid,3>
1535 {
1536 static const bool v = true;
1537 };
1538
1540 template <>
1541 struct hasBackupRestoreFacilities<CpGrid>
1542 {
1543 static const bool v = false;
1544 };
1545
1546 }
1547
1548 template<class DataHandle>
1549 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1550 {
1551 current_view_data_->communicate(data, iftype, dir);
1552 }
1553
1554
1555 template<class DataHandle>
1556 void CpGrid::scatterData([[maybe_unused]] DataHandle& handle) const
1557 {
1558#if HAVE_MPI
1559 if (distributed_data_.empty()) {
1560 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
1561 } else {
1562 distributed_data_[0]->scatterData(handle, data_[0].get(),
1563 distributed_data_[0].get(),
1566 }
1567#endif
1568 }
1569
1570 template<class DataHandle>
1571 void CpGrid::gatherData([[maybe_unused]] DataHandle& handle) const
1572 {
1573#if HAVE_MPI
1574 if (distributed_data_.empty()) {
1575 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
1576 } else {
1577 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1578 }
1579#endif
1580 }
1581
1582
1583 template<class Cell2FacesRowIterator>
1584 int
1585 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
1586 {
1587 // Note that this relies on the following implementation detail:
1588 // The grid is always constructed such that the interior faces constructed
1589 // with orientation set to true are
1590 // oriented along the positive IJK direction. Oriented means that
1591 // the first cell attached to face has the lower index.
1592 // For faces along the boundary (only one cell, always attached at index 0)
1593 // the orientation has to be determined by the orientation of the cell.
1594 // If it is true then in UnstructuredGrid it would be stored at index 0,
1595 // otherwise at index 1.
1596 const int cell = cell_face.getCellIndex();
1597 const int face = *cell_face;
1598 assert (0 <= cell); assert (cell < numCells());
1599 assert (0 <= face); assert (face < numFaces());
1600
1601 typedef cpgrid::OrientedEntityTable<1,0>::row_type F2C;
1602
1603 const cpgrid::EntityRep<1> f(face, true);
1604 const F2C& f2c = current_view_data_->face_to_cell_[f];
1605 const face_tag tag = current_view_data_->face_tag_[f];
1606
1607 assert ((f2c.size() == 1) || (f2c.size() == 2));
1608
1609 int inside_cell = 0;
1610
1611 if ( f2c.size() == 2 ) // Two cells => interior
1612 {
1613 if ( f2c[1].index() == cell )
1614 {
1615 inside_cell = 1;
1616 }
1617 }
1618 const bool normal_is_in = ! f2c[inside_cell].orientation();
1619
1620 switch (tag) {
1621 case I_FACE:
1622 // LEFT : RIGHT
1623 return normal_is_in ? 0 : 1; // min(I) : max(I)
1624 case J_FACE:
1625 // BACK : FRONT
1626 return normal_is_in ? 2 : 3; // min(J) : max(J)
1627 case K_FACE:
1628 // Note: TOP at min(K) as 'z' measures *depth*.
1629 // TOP : BOTTOM
1630 return normal_is_in ? 4 : 5; // min(K) : max(K)
1631 case NNC_FACE:
1632 // For nnc faces we return the otherwise unused value -1.
1633 return -1;
1634 default:
1635 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
1636 }
1637 }
1638
1639 template<int dim>
1640 cpgrid::Entity<dim> createEntity(const CpGrid&, int, bool);
1641
1642} // namespace Dune
1643
1644#include <opm/grid/cpgrid/PersistentContainer.hpp>
1645#include <opm/grid/cpgrid/CartesianIndexMapper.hpp>
1646#include "cpgrid/Intersection.hpp"
1647#include "cpgrid/Geometry.hpp"
1648#include "cpgrid/Indexsets.hpp"
1649
1650#endif // OPM_CPGRID_HEADER
An iterator over the centroids of the geometry of the entities.
Definition CpGrid.hpp:1253
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition CpGrid.hpp:1260
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition CpGrid.hpp:1257
[ 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:979
std::string name() const
Get the grid name.
Definition CpGrid.cpp:761
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...
Definition CpGrid.cpp:726
void switchToGlobalView()
Switch to the global view.
Definition CpGrid.cpp:1660
unsigned int numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition CpGrid.cpp:1111
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:1571
int numFaces(int level=-1) const
Get the number of faces.
Definition CpGrid.cpp:1153
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition CpGrid.cpp:1018
std::vector< std::array< int, 2 > > mapLeafIndexSetToLocalCartesianIndexSets() const
Reverse map: from leaf index cell to { level, local/level Cartesian index of the cell }...
Definition CpGrid.cpp:736
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData() const
Returns either data_ or distributed_data_(if non empty).
Definition CpGrid.cpp:670
CentroidIterator< 0 > beginCellCentroids() const
Get an iterator over the cell centroids positioned at the first one.
Definition CpGrid.cpp:1594
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition CpGrid.cpp:1013
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:933
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition CpGrid.hpp:1585
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGrid.cpp:756
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.
Definition CpGrid.cpp:1846
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:839
void createCartesian(const std::array< int, 3 > &dims, const std::array< double, 3 > &cellsize, const std::array< int, 3 > &shift={0, 0, 0})
Create a cartesian grid.
Definition CpGrid.cpp:591
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.
Definition CpGrid.cpp:1089
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition CpGrid.hpp:1347
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:882
Traits::template Codim< codim >::LevelIterator lend(int level) const
one past the end on this level
const Vector & faceNormal(int face) const
Get the unit normal of a face.
Definition CpGrid.cpp:1579
int numCellFaces() const
Get the sum of all faces attached to all cells.
Definition CpGrid.cpp:1221
const std::vector< int > & globalCell() const
Retrieve mapping from internal ("compressed") active grid cells to external ("uncompressed") cells.
Definition CpGrid.cpp:680
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:738
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Definition CpGrid.cpp:1589
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:705
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
Definition CpGrid.cpp:766
const CpGridTraits::Communication & comm() const
Get the collective communication object.
Definition CpGrid.cpp:1136
double faceArea(int face) const
Get the area of a face.
Definition CpGrid.cpp:1569
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
Definition CpGrid.cpp:1574
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition CpGrid.cpp:1023
void postAdapt()
Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'.
Definition CpGrid.cpp:2385
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition CpGrid.cpp:989
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the global grid.
Definition CpGrid.cpp:663
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
Definition CpGrid.cpp:1231
int numVertices() const
Get The number of vertices.
Definition CpGrid.cpp:1159
Dune::cpgrid::Intersection getParentIntersectionFromLgrBoundaryFace(const Dune::cpgrid::Intersection &intersection) const
Definition CpGrid.cpp:1236
void getIJK(const int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition CpGrid.cpp:746
void processEclipseFormat(const grdecl &input_data, bool remove_ij_boundary, bool turn_normals=false, bool edge_conformal=false)
Read the Eclipse grid format ('grdecl').
Definition CpGrid.cpp:1743
void syncDistributedGlobalCellIds()
Synchronizes cell global ids across processes after load balancing.
Definition CpGrid.cpp:2530
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGrid.cpp:751
int faceCell(int face, int local_index, int level=-1) const
Get the index identifying a cell attached to a face.
Definition CpGrid.cpp:1183
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim and PartitionIteratorType.
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
------------— Adaptivity (begin) ------------—
Definition CpGrid.cpp:1776
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...
Definition CpGrid.cpp:2593
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:1009
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
Definition CpGrid.cpp:1095
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:791
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
Definition CpGrid.cpp:1079
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition CpGrid.cpp:1030
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
const InterfaceMap & pointScatterGatherInterface() const
Get an interface for gathering/scattering data attached to points with communication.
Definition CpGrid.cpp:1655
CpGrid()
Default constructor.
Definition CpGrid.cpp:167
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::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
Definition CpGrid.cpp:1604
bool adapt()
Triggers the grid refinement process.
Definition CpGrid.cpp:1812
void switchToDistributedView()
Switch to the distributed view.
Definition CpGrid.cpp:1666
CentroidIterator< 1 > beginFaceCentroids() const
Get an iterator over the face centroids positioned at the first one.
Definition CpGrid.cpp:1599
bool preAdapt()
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called af...
Definition CpGrid.cpp:1800
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
Definition CpGrid.cpp:1795
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:1059
const cpgrid::OrientedEntityTable< 0, 1 >::row_type cellFaceRow(int cell) const
Get a list of indices identifying all faces of a cell.
Definition CpGrid.cpp:1178
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
Definition CpGrid.cpp:1650
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:1033
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.
Definition CpGrid.cpp:194
int numCells(int level=-1) const
Get the number of cells.
Definition CpGrid.cpp:1147
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
Definition CpGrid.cpp:1564
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
Definition CpGrid.cpp:1422
void globalRefine(int refCount)
Refine the grid refCount times using the default refinement rule.
Definition CpGrid.cpp:1035
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
The new communication interface.
Definition CpGrid.hpp:1098
int cellFace(int cell, int local_index, int level=-1) const
Get a specific face of a cell.
Definition CpGrid.cpp:1171
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
double cellVolume(int cell) const
Get the volume of the cell.
Definition CpGrid.cpp:1584
void autoRefine(const std::array< int, 3 > &nxnynz)
Global refine the grid with different refinement factors in each direction.
Definition CpGrid.cpp:2750
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition CpGrid.hpp:1556
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 Entity.hpp:72
This class encapsulates geometry for vertices, intersections, and cells.
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:201
Definition Indexsets.hpp:56
Definition Intersection.hpp:279
Definition Intersection.hpp:66
Iterator intended to be used as LeafIterator and LevelIterator (no difference due to no adaptivity) f...
Definition Iterators.hpp:60
Definition Indexsets.hpp:371
Copyright 2019 Equinor AS.
Definition GridPartitioning.cpp:673
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
@ 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
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68
Low-level corner-point processing routines and supporting data structures.
face_tag
Connection taxonomy.
Definition preprocess.h:66
@ K_FACE
Connection topologically normal to I-J plane.
Definition preprocess.h:69
@ J_FACE
Connection topologically normal to I-K plane.
Definition preprocess.h:68
@ NNC_FACE
Arbitrary non-neighbouring connection.
Definition preprocess.h:70
@ I_FACE
Connection topologically normal to J-K plane.
Definition preprocess.h:67
Definition CpGrid.hpp:190
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.
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
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
Raw corner-point specification of a particular geological model.
Definition preprocess.h:56