CpGrid.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: CpGrid.hpp
4//
5// Created: Fri May 29 20:26:36 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// B�rd Skaflestad <bard.skaflestad@sintef.no>
9// Antonella Ritorto <antonella.ritorto@opm-op.com>
10//
11// $Date$
12//
13// $Revision$
14//
15//===========================================================================
16
17/*
18 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
19 Copyright 2009, 2010, 2014, 2022-2023 Equinor ASA.
20 Copyright 2014, 2015 Dr. Blatt - HPC-Simulartion-Software & Services
21 Copyright 2015 NTNU
22
23 This file is part of The Open Porous Media project (OPM).
24
25 OPM is free software: you can redistribute it and/or modify
26 it under the terms of the GNU General Public License as published by
27 the Free Software Foundation, either version 3 of the License, or
28 (at your option) any later version.
29
30 OPM is distributed in the hope that it will be useful,
31 but WITHOUT ANY WARRANTY; without even the implied warranty of
32 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
33 GNU General Public License for more details.
34
35 You should have received a copy of the GNU General Public License
36 along with OPM. If not, see <http://www.gnu.org/licenses/>.
37*/
38
39#ifndef OPM_CPGRID_HEADER
40#define OPM_CPGRID_HEADER
41
42// Warning suppression for Dune includes.
44
45#include <dune/grid/common/grid.hh>
46
48
50
54
56
58
59#include <set>
60
61namespace Opm
62{
63struct NNCdata;
64class EclipseGrid;
65class EclipseState;
66}
67
68namespace Dune
69{
70 class CpGrid;
71
72 namespace cpgrid
73 {
74 class CpGridData;
75 template <int> class Entity;
76 template<int,int> class Geometry;
77 class HierarchicIterator;
78 class IntersectionIterator;
79 template<int, PartitionIteratorType> class Iterator;
80 class LevelGlobalIdSet;
81 class GlobalIdSet;
82 class Intersection;
83 class IntersectionIterator;
84 class IndexSet;
85 class IdSet;
86
87 }
88}
89
90namespace Dune
91{
92
94 //
95 // CpGridTraits
96 //
98
100 {
102 typedef CpGrid Grid;
103
112
115
118 template <int cd>
119 struct Codim
120 {
123 typedef cpgrid::Geometry<3-cd, 3> Geometry;
124 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> Geometry;
127 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> LocalGeometry;
130
133
136
139
142 template <PartitionIteratorType pitype>
144 {
149 };
150 };
151
154 template <PartitionIteratorType pitype>
156 {
158 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid> > LevelGridView;
160 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> > LeafGridView;
161
162 };
163
165 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid>> LevelGridView;
167 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>> LeafGridView;
168
177
181 };
182
184 //
185 // CpGridFamily
186 //
188
190 {
192 };
193
195 //
196 // CpGrid
197 //
199
201 class CpGrid
202 : public GridDefaultImplementation<3, 3, double, CpGridFamily>
203 {
204 friend class cpgrid::CpGridData;
205 friend class cpgrid::Entity<0>;
206 friend class cpgrid::Entity<1>;
207 friend class cpgrid::Entity<2>;
208 friend class cpgrid::Entity<3>;
209 template<int dim>
210 friend cpgrid::Entity<dim> createEntity(const CpGrid&,int,bool);
211
212 public:
213
214 // --- Typedefs ---
215
216
219
220
221 // --- Methods ---
222
223
226
227 explicit CpGrid(MPIHelper::MPICommunicator comm);
228
229#if HAVE_OPM_COMMON
266 std::vector<std::size_t>
267 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
268 Opm::EclipseState* ecl_state,
269 bool periodic_extension,
270 bool turn_normals,
271 bool clip_z,
272 bool pinchActive,
273 bool edge_conformal);
274
312 std::vector<std::size_t>
313 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
314 Opm::EclipseState* ecl_state,
315 bool periodic_extension,
316 bool turn_normals = false,
317 bool clip_z = false,
318 bool edge_conformal = false);
319#endif // HAVE_OPM_COMMON
320
335 void processEclipseFormat(const grdecl& input_data,
336 bool remove_ij_boundary,
337 bool turn_normals = false,
338 bool edge_conformal = false);
339
341
347
348
355 void createCartesian(const std::array<int, 3>& dims,
356 const std::array<double, 3>& cellsize,
357 const std::array<int, 3>& shift = {0,0,0});
358
362 const std::array<int, 3>& logicalCartesianSize() const;
363
371 const std::vector<int>& globalCell() const;
372
374 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData() const;
375
377 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData();
378
381
384
392 void getIJK(const int c, std::array<int,3>& ijk) const;
394
398 bool uniqueBoundaryIds() const;
399
402 void setUniqueBoundaryIds(bool uids);
403
404
405 // --- Dune interface below ---
406
408 // \@{
413 std::string name() const;
414
416 int maxLevel() const;
417
419 template<int codim>
420 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
422 template<int codim>
423 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
424
426 template<int codim, PartitionIteratorType PiType>
427 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
429 template<int codim, PartitionIteratorType PiType>
430 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
431
433 template<int codim>
434 typename Traits::template Codim<codim>::LeafIterator leafbegin() const;
436 template<int codim>
437 typename Traits::template Codim<codim>::LeafIterator leafend() const;
438
440 template<int codim, PartitionIteratorType PiType>
441 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const;
443 template<int codim, PartitionIteratorType PiType>
444 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const;
445
447 int size (int level, int codim) const;
448
450 int size (int codim) const;
451
453 int size (int level, GeometryType type) const;
454
456 int size (GeometryType type) const;
457
459 const Traits::GlobalIdSet& globalIdSet() const;
460
462 const Traits::LocalIdSet& localIdSet() const;
463
465 const Traits::LevelIndexSet& levelIndexSet(int level) const;
466
468 const Traits::LeafIndexSet& leafIndexSet() const;
469
475 void globalRefine (int refCount, bool throwOnFailure = false);
476
477 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
478
480 template <int codim>
482
501 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
502 const std::vector<std::array<int,3>>& startIJK_vec,
503 const std::vector<std::array<int,3>>& endIJK_vec,
504 const std::vector<std::string>& lgr_name_vec,
505 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
506
513 void autoRefine(const std::array<int,3>& nxnynz);
514
515 // @brief TO BE DONE
516 const std::map<std::string,int>& getLgrNameToLevel() const;
517
518 // @breif Compute center of an entity/element/cell in the Eclipse way:
519 // - Average of the 4 corners of the bottom face.
520 // - Average of the 4 corners of the top face.
521 // Return average of the previous computations.
522 // @param [in] int Index of a cell.
523 // @return 'eclipse centroid'
524 std::array<double,3> getEclCentroid(const int& idx) const;
525
526 // @breif Compute center of an entity/element/cell in the Eclipse way:
527 // - Average of the 4 corners of the bottom face.
528 // - Average of the 4 corners of the top face.
529 // Return average of the previous computations.
530 // @param [in] Entity<0> Entity
531 // @return 'eclipse centroid'
532 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
533
544
565 bool mark(int refCount, const cpgrid::Entity<0>& element, bool throwOnFailure = false);
566
570 int getMark(const cpgrid::Entity<0>& element) const;
571
574 bool preAdapt();
575
578 bool adapt();
579
598 bool refineAndUpdateGrid(bool throwOnFailure,
599 const std::vector<std::array<int,3>>& cells_per_dim_vec,
600 const std::vector<int>& assignRefinedLevel,
601 const std::vector<std::string>& lgr_name_vec,
602 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
603 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
604
606 void postAdapt();
608
609 private:
610 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
611 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
612 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
613 const int& corner_count,
614 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
615 const int& preAdaptMaxLevel,
616 const int& newLevels);
617
618 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
619 const std::vector<std::array<int,3>>& cells_per_dim_vec,
620 const std::vector<int>& lgr_with_at_least_one_active_cell);
621
628 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
629 public:
636 private:
637
650 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
651
655 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
656
657 private:
671 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
672 const std::vector<std::array<int,3>>& endIJK_vec,
673 std::vector<int>& assignRefinedLevel,
674 std::vector<int>& lgr_with_at_least_one_active_cell);
675
677 void populateCellIndexSetRefinedGrid(int level);
678
680 void populateCellIndexSetLeafGridView();
681
683 void populateLeafGlobalIdSet();
684
685 public:
686
691 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
692
694 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
695
697 unsigned int overlapSize(int) const;
698
699
701 unsigned int ghostSize(int) const;
702
704 unsigned int overlapSize(int, int) const;
705
707 unsigned int ghostSize(int, int) const;
708
710 unsigned int numBoundarySegments() const;
711
712 void setPartitioningParams(const std::map<std::string,std::string>& params);
713
714 // loadbalance is not part of the grid interface therefore we skip it.
715
725 bool loadBalance(int overlapLayers=1,
726 int partitionMethod = Dune::PartitionMethod::zoltan,
727 double imbalanceTol = 1.1,
728 int level =-1)
729 {
730 using std::get;
731 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
732 /* ownersFirst = */ false,
733 /* wells = */ nullptr,
734 /* possibleFutureConnections = */ {},
735 /* serialPartitioning = */ false,
736 /* transmissibilities = */ nullptr,
737 /* addCornerCells = */ true,
738 overlapLayers,
739 partitionMethod,
740 imbalanceTol,
741 /* allowDistributedWells = */ false,
742 /* input_cell_part = */ {},
743 level));
744 }
745
746 // loadbalance is not part of the grid interface therefore we skip it.
747
758 bool loadBalanceSerial(int overlapLayers=1,
759 int partitionMethod = Dune::PartitionMethod::zoltan,
760 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
761 double imbalanceTol = 1.1,
762 int level = -1)
763 {
764 using std::get;
765 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
766 /* ownersFirst = */ false,
767 /* wells = */ nullptr,
768 /* possibleFutureConnections = */ {},
769 /* serialPartitioning = */ true,
770 /* transmissibilities = */ nullptr,
771 /* addCornerCells = */ true,
772 overlapLayers,
773 partitionMethod,
774 imbalanceTol,
775 /* allowDistributedWells = */ false,
776 /* input_cell_part = */ {},
777 level));
778 }
779
780 // loadbalance is not part of the grid interface therefore we skip it.
781
810 std::pair<bool,std::vector<std::pair<std::string,bool>>>
811 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
812 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
813 const double* transmissibilities = nullptr,
814 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
815 int level = -1)
816 {
817 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
818 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
819 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
820 /* input_cell_part = */ {}, level);
821 }
822
823 // loadbalance is not part of the grid interface therefore we skip it.
824
858 std::pair<bool,std::vector<std::pair<std::string,bool>>>
859 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
860 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
861 const double* transmissibilities = nullptr, bool ownersFirst=false,
862 bool addCornerCells=false, int overlapLayers=1,
863 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
864 double imbalanceTol = 1.1,
865 int level = -1)
866 {
867 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
868 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
869 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
870 }
871
900 template<class DataHandle>
901 std::pair<bool, std::vector<std::pair<std::string,bool> > >
902 loadBalance(DataHandle& data,
903 const std::vector<cpgrid::OpmWellType> * wells,
904 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
905 const double* transmissibilities = nullptr,
906 int overlapLayers=1, int partitionMethod = 1, int level =-1)
907 {
908 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
909 using std::get;
910 if (get<0>(ret))
911 {
913 }
914 return ret;
915 }
916
951 template<class DataHandle>
952 std::pair<bool, std::vector<std::pair<std::string,bool> > >
953 loadBalance(DataHandle& data, EdgeWeightMethod method,
954 const std::vector<cpgrid::OpmWellType> * wells,
955 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
956 bool serialPartitioning,
957 const double* transmissibilities = nullptr, bool ownersFirst=false,
958 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
959 double imbalanceTol = 1.1,
960 bool allowDistributedWells = false)
961 {
962 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
963 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
964 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
965 using std::get;
966 if (get<0>(ret))
967 {
969 }
970 return ret;
971 }
972
997 template<class DataHandle>
998 std::pair<bool, std::vector<std::pair<std::string,bool> > >
999 loadBalance(DataHandle& data, const std::vector<int>& parts,
1000 const std::vector<cpgrid::OpmWellType> * wells,
1001 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
1002 bool ownersFirst=false,
1003 bool addCornerCells=false, int overlapLayers=1)
1004 {
1005 using std::get;
1006 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
1007 possibleFutureConnections,
1008 /* serialPartitioning = */ false,
1009 /* transmissibilities = */ {},
1010 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1011 /* imbalanceTol (ignored) = */ 0.0,
1012 /* allowDistributedWells = */ true, parts, /* level = */ 0);
1013 using std::get;
1014 if (get<0>(ret))
1015 {
1017 }
1018 return ret;
1019 }
1020
1028 template<class DataHandle>
1029 bool loadBalance(DataHandle& data,
1030 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1031 {
1032 // decltype usage needed to tell the compiler not to use this function if first
1033 // argument is std::vector but rather loadbalance by parts
1034 bool ret = loadBalance(overlapLayers, partitionMethod);
1035 if (ret)
1036 {
1038 }
1039 return ret;
1040 }
1041
1053 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1054 bool addCornerCells=false, int overlapLayers=1)
1055 {
1056 using std::get;
1057 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1058 {},
1059 /* serialPartitioning = */ false,
1060 /* trabsmissibilities = */ {},
1061 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1062 /* imbalanceTol (ignored) = */ 0.0,
1063 /* allowDistributedWells = */ true, parts));
1064 }
1065
1078 template<class DataHandle>
1079 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1080 bool addCornerCells=false, int overlapLayers=1)
1081 {
1082 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1083 if (ret)
1084 {
1086 }
1087 return ret;
1088 }
1089
1103 std::vector<int>
1104 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1105 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1106 const double* transmissibilities,
1107 const int numParts,
1108 const double imbalanceTol) const;
1109
1117 template<class DataHandle>
1118 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1119 {
1120 communicate(data, iftype, dir);
1121 }
1122
1130 template<class DataHandle>
1131 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1132
1134 const typename CpGridTraits::Communication& comm () const;
1136
1137 // ------------ End of Dune interface, start of simplified interface --------------
1138
1144
1145 // enum { dimension = 3 }; // already defined
1146
1147 typedef Dune::FieldVector<double, 3> Vector;
1148
1149
1150 const std::vector<double>& zcornData() const;
1151
1152
1153 // Topology
1158 int numCells(int level = -1) const;
1159
1164 int numFaces(int level = -1) const;
1165
1167 int numVertices() const;
1168
1169
1178 int numCellFaces(int cell, int level = -1) const;
1179
1186 int cellFace(int cell, int local_index, int level = -1) const;
1187
1191
1211 int faceCell(int face, int local_index, int level = -1) const;
1212
1219 int numCellFaces() const;
1220
1221 int numFaceVertices(int face) const;
1222
1227 int faceVertex(int face, int local_index) const;
1228
1231 double cellCenterDepth(int cell_index) const;
1232
1233
1234 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1235
1236 const Vector faceAreaNormalEcl(int face) const;
1237
1238
1239 // Geometry
1243 const Vector& vertexPosition(int vertex) const;
1244
1247 double faceArea(int face) const;
1248
1251 const Vector& faceCentroid(int face) const;
1252
1256 const Vector& faceNormal(int face) const;
1257
1260 double cellVolume(int cell) const;
1261
1264 const Vector& cellCentroid(int cell) const;
1265
1268 template<int codim>
1270 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1271 FieldVector<double, 3>,
1272 const FieldVector<double, 3>&, int>
1273 {
1274 public:
1276 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1281 : iter_(iter)
1282 {}
1283
1284 const FieldVector<double,3>& dereference() const
1285 {
1286 return iter_->center();
1287 }
1289 {
1290 ++iter_;
1291 }
1292 const FieldVector<double,3>& elementAt(int n)
1293 {
1294 return iter_[n]->center();
1295 }
1296 void advance(int n){
1297 iter_+=n;
1298 }
1300 {
1301 --iter_;
1302 }
1304 {
1305 return o-iter_;
1306 }
1307 bool equals(const CentroidIterator& o) const{
1308 return o==iter_;
1309 }
1310 private:
1312 GeometryIterator iter_;
1313 };
1314
1317
1320
1321 // Extra
1322 int boundaryId(int face) const;
1323
1330 template<class Cell2FacesRowIterator>
1331 int
1332 faceTag(const Cell2FacesRowIterator& cell_face) const;
1333
1335
1336 // ------------ End of simplified interface --------------
1337
1338 //------------- methods not in the DUNE grid interface.
1339
1344
1345
1354 template<class DataHandle>
1355 void scatterData(DataHandle& handle) const;
1356
1363 template<class DataHandle>
1364 void gatherData(DataHandle& handle) const;
1365
1368
1398
1402
1405
1409
1410#if HAVE_MPI
1415
1418
1423
1425
1427
1429
1431#endif
1432
1434 const std::vector<int>& sortedNumAquiferCells() const;
1435
1436 private:
1468 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1469 scatterGrid(EdgeWeightMethod method,
1470 bool ownersFirst,
1471 const std::vector<cpgrid::OpmWellType> * wells,
1472 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1473 bool serialPartitioning,
1474 const double* transmissibilities,
1475 bool addCornerCells,
1476 int overlapLayers,
1477 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1478 double imbalanceTol = 1.1,
1479 bool allowDistributedWells = true,
1480 const std::vector<int>& input_cell_part = {},
1481 int level = -1);
1482
1487 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1489 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1491 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1493 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
1499 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1500 /*
1501 * @brief Interface for scattering and gathering point data.
1502 *
1503 * @warning Will only update owner cells
1504 */
1505 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1509 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1510
1511
1515 std::map<std::string,std::string> partitioningParams;
1516
1517 }; // end Class CpGrid
1518
1519} // end namespace Dune
1520
1524
1525
1526namespace Dune
1527{
1528
1529 namespace Capabilities
1530 {
1532 template <>
1533 struct hasEntity<CpGrid, 0>
1534 {
1535 static const bool v = true;
1536 };
1537
1539 template <>
1540 struct hasEntity<CpGrid, 3>
1541 {
1542 static const bool v = true;
1543 };
1544
1545 template<>
1546 struct canCommunicate<CpGrid,0>
1547 {
1548 static const bool v = true;
1549 };
1550
1551 template<>
1552 struct canCommunicate<CpGrid,3>
1553 {
1554 static const bool v = true;
1555 };
1556
1558 template <>
1559 struct hasBackupRestoreFacilities<CpGrid>
1560 {
1561 static const bool v = false;
1562 };
1563
1564 }
1565
1566 template<class DataHandle>
1567 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1568 {
1569 current_data_->back()->communicate(data, iftype, dir);
1570 }
1571
1572
1573 template<class DataHandle>
1574 void CpGrid::scatterData([[maybe_unused]] DataHandle& handle) const
1575 {
1576#if HAVE_MPI
1577 if (distributed_data_.empty()) {
1578 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
1579 } else {
1580 distributed_data_[0]->scatterData(handle, data_[0].get(),
1581 distributed_data_[0].get(),
1584 }
1585#endif
1586 }
1587
1588 template<class DataHandle>
1589 void CpGrid::gatherData([[maybe_unused]] DataHandle& handle) const
1590 {
1591#if HAVE_MPI
1592 if (distributed_data_.empty()) {
1593 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
1594 } else {
1595 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1596 }
1597#endif
1598 }
1599
1600
1601 template<class Cell2FacesRowIterator>
1602 int
1603 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
1604 {
1605 // Note that this relies on the following implementation detail:
1606 // The grid is always constructed such that the interior faces constructed
1607 // with orientation set to true are
1608 // oriented along the positive IJK direction. Oriented means that
1609 // the first cell attached to face has the lower index.
1610 // For faces along the boundary (only one cell, always attached at index 0)
1611 // the orientation has to be determined by the orientation of the cell.
1612 // If it is true then in UnstructuredGrid it would be stored at index 0,
1613 // otherwise at index 1.
1614 const int cell = cell_face.getCellIndex();
1615 const int face = *cell_face;
1616 assert (0 <= cell); assert (cell < numCells());
1617 assert (0 <= face); assert (face < numFaces());
1618
1620
1621 const cpgrid::EntityRep<1> f(face, true);
1622 const F2C& f2c = current_data_->back()->face_to_cell_[f];
1623 const face_tag tag = current_data_->back()->face_tag_[f];
1624
1625 assert ((f2c.size() == 1) || (f2c.size() == 2));
1626
1627 int inside_cell = 0;
1628
1629 if ( f2c.size() == 2 ) // Two cells => interior
1630 {
1631 if ( f2c[1].index() == cell )
1632 {
1633 inside_cell = 1;
1634 }
1635 }
1636 const bool normal_is_in = ! f2c[inside_cell].orientation();
1637
1638 switch (tag) {
1639 case I_FACE:
1640 // LEFT : RIGHT
1641 return normal_is_in ? 0 : 1; // min(I) : max(I)
1642 case J_FACE:
1643 // BACK : FRONT
1644 return normal_is_in ? 2 : 3; // min(J) : max(J)
1645 case K_FACE:
1646 // Note: TOP at min(K) as 'z' measures *depth*.
1647 // TOP : BOTTOM
1648 return normal_is_in ? 4 : 5; // min(K) : max(K)
1649 case NNC_FACE:
1650 // For nnc faces we return the otherwise unused value -1.
1651 return -1;
1652 default:
1653 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
1654 }
1655 }
1656
1657 template<int dim>
1659
1660} // namespace Dune
1661
1664#include "cpgrid/Intersection.hpp"
1665#include "cpgrid/Geometry.hpp"
1666#include "cpgrid/Indexsets.hpp"
1667
1668#endif // OPM_CPGRID_HEADER
DataHandle & data
Definition: CpGridData.hpp:1164
#define OPM_THROW(Exception, message)
Definition: ErrorMacros.hpp:29
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:1273
void increment()
Definition: CpGrid.hpp:1288
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:1280
const FieldVector< double, 3 > & dereference() const
Definition: CpGrid.hpp:1284
int distanceTo(const CentroidIterator &o)
Definition: CpGrid.hpp:1303
void decrement()
Definition: CpGrid.hpp:1299
void advance(int n)
Definition: CpGrid.hpp:1296
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:1277
bool equals(const CentroidIterator &o) const
Definition: CpGrid.hpp:1307
const FieldVector< double, 3 > & elementAt(int n)
Definition: CpGrid.hpp:1292
[ 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:999
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:1589
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:953
cpgrid::CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGrid.hpp:1412
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:1603
const std::vector< double > & zcornData() const
ParallelIndexSet & getCellIndexSet()
void setUniqueBoundaryIds(bool uids)
std::array< double, 3 > getEclCentroid(const cpgrid::Entity< 0 > &elem) const
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, double imbalanceTol=1.1, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:859
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:1367
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:902
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:758
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:725
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:1414
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
Returns the coarse-level-zero ("parent/origin") intersection of a refined intersection (face) in the ...
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:1029
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:811
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:1147
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:1079
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
const cpgrid::OrientedEntityTable< 0, 1 >::row_type cellFaceRow(int cell) const
Get a list of indices identifying all faces of a cell.
bool refineAndUpdateGrid(bool throwOnFailure, const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< int > &assignRefinedLevel, const std::vector< std::string > &lgr_name_vec, const std::vector< std::array< int, 3 > > &startIJK_vec=std::vector< std::array< int, 3 > >{}, const std::vector< std::array< int, 3 > > &endIJK_vec=std::vector< std::array< int, 3 > >{})
Triggers the grid refinement process, allowing to select diffrent refined level grids.
bool loadBalance(const std::vector< int > &parts, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1053
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:1118
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:1417
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:1574
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:118
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:98
Definition: Geometry.hpp:76
The global id set for Dune.
Definition: Indexsets.hpp:483
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:118
Definition: Indexsets.hpp:57
Definition: Intersection.hpp:276
Definition: Intersection.hpp:63
Definition: Iterators.hpp:60
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
@ simple
Use simple approach based on rectangular partitioning the underlying cartesian grid.
Definition: GridEnums.hpp:46
@ zoltanGoG
use Zoltan on GraphOfGrid for partitioning
Definition: GridEnums.hpp:52
@ zoltan
Use Zoltan for partitioning.
Definition: GridEnums.hpp:48
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Me...
Definition: GridEnums.hpp:34
@ defaultTransEdgeWgt
Use the transmissibilities as edge weights.
Definition: GridEnums.hpp:38
void getFirstChildGlobalIds(const Dune::CpGrid &grid, std::vector< int > &parentToFirstChildGlobalIds)
Retrieves the global ids of the first child for each parent cell in the grid.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
face_tag
Definition: preprocess.h:66
@ K_FACE
Definition: preprocess.h:69
@ J_FACE
Definition: preprocess.h:68
@ NNC_FACE
Definition: preprocess.h:70
@ I_FACE
Definition: preprocess.h:67
Definition: CpGrid.hpp:190
CpGridTraits Traits
Definition: CpGrid.hpp:191
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:144
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition: CpGrid.hpp:146
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition: CpGrid.hpp:148
Traits associated with a specific codim.
Definition: CpGrid.hpp:120
cpgrid::Entity< cd > Entity
The type of the entity.
Definition: CpGrid.hpp:129
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity. IMPORTANT: Codim<codim>::Geometry == Geometry<di...
Definition: CpGrid.hpp:123
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition: CpGrid.hpp:126
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition: CpGrid.hpp:135
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition: CpGrid.hpp:132
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition: CpGrid.hpp:138
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:156
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:160
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:158
Definition: CpGrid.hpp:100
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition: CpGrid.hpp:170
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition: CpGrid.hpp:109
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:167
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:165
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition: CpGrid.hpp:174
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition: CpGrid.hpp:111
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition: CpGrid.hpp:172
cpgrid::CpGridDataTraits::CollectiveCommunication CollectiveCommunication
Definition: CpGrid.hpp:180
GlobalIdSet LocalIdSet
The type of the local id set.
Definition: CpGrid.hpp:176
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGrid.hpp:179
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition: CpGrid.hpp:114
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition: CpGrid.hpp:107
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition: CpGrid.hpp:105
CpGrid Grid
The type that implements the grid.
Definition: CpGrid.hpp:102
Dune::RemoteIndices< ParallelIndexSet > RemoteIndices
The type of the remote indices information.
Definition: CpGridDataTraits.hpp:83
typename CommunicationType::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGridDataTraits.hpp:80
Dune::Communication< MPICommunicator > CollectiveCommunication
Definition: CpGridDataTraits.hpp:59
Dune::OwnerOverlapCopyCommunication< int, int > CommunicationType
type of OwnerOverlap communication for cells
Definition: CpGridDataTraits.hpp:77
Dune::Communication< MPICommunicator > Communication
Definition: CpGridDataTraits.hpp:58
Communicator::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGridDataTraits.hpp:74
Definition: preprocess.h:56