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.
43#include <opm/grid/utility/platform_dependent/disable_warnings.h> // Not really needed it seems, but alas.
44
45#include <dune/grid/common/grid.hh>
50#include <opm/grid/utility/platform_dependent/reenable_warnings.h> // Not really needed it seems, but alas.
51#include "common/GridEnums.hpp"
53
54#include <set>
55
56namespace Opm
57{
58struct NNCdata;
59class EclipseGrid;
60class EclipseState;
61}
62
63namespace Dune
64{
65 class CpGrid;
66
67 namespace cpgrid
68 {
69 class CpGridData;
70 template <int> class Entity;
71 template<int,int> class Geometry;
72 class HierarchicIterator;
73 class IntersectionIterator;
74 template<int, PartitionIteratorType> class Iterator;
75 class LevelGlobalIdSet;
76 class GlobalIdSet;
77 class Intersection;
78 class IntersectionIterator;
79 class IndexSet;
80 class IdSet;
81
82 }
83}
84
85namespace Dune
86{
87
89 //
90 // CpGridTraits
91 //
93
95 {
97 typedef CpGrid Grid;
98
107
110
113 template <int cd>
114 struct Codim
115 {
118 typedef cpgrid::Geometry<3-cd, 3> Geometry;
119 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> Geometry;
122 //typedef Dune::Geometry<3-cd, 3, CpGrid, cpgrid::Geometry> LocalGeometry;
125
128
131
134
137 template <PartitionIteratorType pitype>
139 {
144 };
145 };
146
149 template <PartitionIteratorType pitype>
151 {
153 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid> > LevelGridView;
155 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> > LeafGridView;
156
157 };
158
160 typedef Dune::GridView<DefaultLevelGridViewTraits<CpGrid>> LevelGridView;
162 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>> LeafGridView;
163
172
176 };
177
179 //
180 // CpGridFamily
181 //
183
185 {
187 };
188
190 //
191 // CpGrid
192 //
194
196 class CpGrid
197 : public GridDefaultImplementation<3, 3, double, CpGridFamily>
198 {
199 friend class cpgrid::CpGridData;
200 friend class cpgrid::Entity<0>;
201 friend class cpgrid::Entity<1>;
202 friend class cpgrid::Entity<2>;
203 friend class cpgrid::Entity<3>;
204 template<int dim>
205 friend cpgrid::Entity<dim> createEntity(const CpGrid&,int,bool);
206
207 public:
208
209 // --- Typedefs ---
210
211
214
215
216 // --- Methods ---
217
218
221
222 explicit CpGrid(MPIHelper::MPICommunicator comm);
223
224#if HAVE_ECL_INPUT
243 std::vector<std::size_t>
244 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
245 Opm::EclipseState* ecl_state,
246 bool periodic_extension, bool turn_normals, bool clip_z,
247 bool pinchActive);
248
268 std::vector<std::size_t>
269 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
270 Opm::EclipseState* ecl_state,
271 bool periodic_extension, bool turn_normals = false, bool clip_z = false);
272
273#endif
274
278 void processEclipseFormat(const grdecl& input_data, bool remove_ij_boundary, bool turn_normals = false);
279
281
287
288
295 void createCartesian(const std::array<int, 3>& dims,
296 const std::array<double, 3>& cellsize,
297 const std::array<int, 3>& shift = {0,0,0});
298
302 const std::array<int, 3>& logicalCartesianSize() const;
303
311 const std::vector<int>& globalCell() const;
312
314 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData() const;
315
317 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& currentData();
318
326 void getIJK(const int c, std::array<int,3>& ijk) const;
328
332 bool uniqueBoundaryIds() const;
333
336 void setUniqueBoundaryIds(bool uids);
337
338
339 // --- Dune interface below ---
340
342 // \@{
347 std::string name() const;
348
350 int maxLevel() const;
351
353 template<int codim>
354 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
356 template<int codim>
357 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
358
360 template<int codim, PartitionIteratorType PiType>
361 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
363 template<int codim, PartitionIteratorType PiType>
364 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
365
367 template<int codim>
368 typename Traits::template Codim<codim>::LeafIterator leafbegin() const;
370 template<int codim>
371 typename Traits::template Codim<codim>::LeafIterator leafend() const;
372
374 template<int codim, PartitionIteratorType PiType>
375 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const;
377 template<int codim, PartitionIteratorType PiType>
378 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const;
379
381 int size (int level, int codim) const;
382
384 int size (int codim) const;
385
387 int size (int level, GeometryType type) const;
388
390 int size (GeometryType type) const;
391
393 const Traits::GlobalIdSet& globalIdSet() const;
394
396 const Traits::LocalIdSet& localIdSet() const;
397
399 const Traits::LevelIndexSet& levelIndexSet(int level) const;
400
402 const Traits::LeafIndexSet& leafIndexSet() const;
403
407 void globalRefine (int refCount);
408
409 const std::vector<Dune::GeometryType>& geomTypes(const int) const;
410
412 template <int codim>
414
431 void addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_per_dim_vec,
432 const std::vector<std::array<int,3>>& startIJK_vec,
433 const std::vector<std::array<int,3>>& endIJK_vec,
434 const std::vector<std::string>& lgr_name_vec);
435
436 // @brief TO BE DONE
437 const std::map<std::string,int>& getLgrNameToLevel() const;
438
439 // @breif Compute center of an entity/element/cell in the Eclipse way:
440 // - Average of the 4 corners of the bottom face.
441 // - Average of the 4 corners of the top face.
442 // Return average of the previous computations.
443 // @param [in] int Index of a cell.
444 // @return 'eclipse centroid'
445 std::array<double,3> getEclCentroid(const int& idx) const;
446
447 // @breif Compute center of an entity/element/cell in the Eclipse way:
448 // - Average of the 4 corners of the bottom face.
449 // - Average of the 4 corners of the top face.
450 // Return average of the previous computations.
451 // @param [in] Entity<0> Entity
452 // @return 'eclipse centroid'
453 std::array<double,3> getEclCentroid(const cpgrid::Entity<0>& elem) const;
454
455 // @brief Return parent (coarse) intersection (face) of a refined face on the leaf grid view, whose neighboring cells
456 // are two: one coarse cell (equivalent to its origin cell from level 0), and one refined cell
457 // from certain LGR
459
479 bool mark(int refCount, const cpgrid::Entity<0>& element);
480
484 int getMark(const cpgrid::Entity<0>& element) const;
485
488 bool preAdapt();
489
492 bool adapt();
493
507 bool adapt(const std::vector<std::array<int,3>>& cells_per_dim_vec,
508 const std::vector<int>& assignRefinedLevel,
509 const std::vector<std::string>& lgr_name_vec,
510 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
511 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
512
514 void postAdapt();
516
517 private:
518
520
579 void refineAndProvideMarkedRefinedRelations(/* Marked elements parameters */
580 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
581 int& markedElem_count,
582 std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
583 std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
584 std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
585 /* Refined cells parameters */
586 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCell_to_refinedLevelAdRefinedCell,
587 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
588 std::vector<int>& refined_cell_count_vec,
589 const std::vector<int>& assignRefinedLevel,
590 std::vector<std::vector<std::tuple<int,std::vector<int>>>>& preAdapt_parent_to_children_cells_vec,
591 /* Adapted cells parameters */
592 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCell_to_adaptedCell,
593 std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
594 int& cell_count,
595 std::vector<std::vector<int>>& preAdapt_level_to_leaf_cells_vec,
596 /* Additional parameters */
597 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
598
620 std::tuple< std::vector<std::vector<std::array<int,2>>>,
621 std::vector<std::vector<int>>,
622 std::vector<std::array<int,2>>,
623 std::vector<int>> defineChildToParentAndIdxInParentCell( const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
624 const std::vector<int>& refined_cell_count_vec,
625 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
626 const int& cell_count) const;
627
654 std::pair<std::vector<std::vector<int>>, std::vector<std::array<int,2>>>
655 defineLevelToLeafAndLeafToLevelCells(const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCell_to_refinedLevelAndRefinedCell,
656 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
657 const std::vector<int>& refined_cell_count_vec,
658 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCell_to_adaptedCell,
659 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
660 const int& cell_count) const;
661
688 void identifyRefinedCornersPerLevel(std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
689 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner,
690 std::vector<int>& refined_corner_count_vec,
691 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
692 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
693 const std::vector<int>& assignRefinedLevel,
694 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
695 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
696 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
697
715 void identifyRefinedFacesPerLevel(std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
716 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
717 std::vector<int>& refined_face_count_vec,
718 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
719 const std::vector<int>& assignRefinedLevel,
720 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
721 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
722
743 void identifyLeafGridCorners(std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
744 std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
745 int& corner_count,
746 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
747 const std::vector<int>& assignRefinedLevel,
748 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
749 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
750 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
751 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
752
770 void identifyLeafGridFaces(std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
771 std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
772 int& face_count,
773 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
774 const std::vector<int>& assignRefinedLevel,
775 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
776 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
777
779 void populateRefinedCorners(std::vector<Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<0,3>>>& refined_corners_vec,
780 const std::vector<int>& refined_corner_count_vec,
781 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
782 const int& preAdaptMaxLevel,
783 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner) const;
784
786 void populateRefinedFaces(std::vector<Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<2,3>>>& refined_faces_vec,
787 std::vector<Dune::cpgrid::EntityVariableBase<enum face_tag>>& mutable_refined_face_tags_vec,
788 std::vector<Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>>& mutable_refine_face_normals_vec,
789 std::vector<Opm::SparseTable<int>>& refined_face_to_point_vec,
790 const std::vector<int>& refined_face_count_vec,
791 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
792 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
793 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
794 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
795 const int& preAdaptMaxLevel,
796 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
797 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner) const;
798
800 void populateRefinedCells(std::vector<Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<3,3>>>& refined_cells_vec,
801 std::vector<std::vector<std::array<int,8>>>& refined_cell_to_point_vec,
802 std::vector<std::vector<int>>& refined_global_cell_vec,
803 const std::vector<int>& refined_cell_count_vec,
804 std::vector<cpgrid::OrientedEntityTable<0,1>>& refined_cell_to_face_vec,
805 std::vector<cpgrid::OrientedEntityTable<1,0>>& refined_face_to_cell_vec,
806 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
807 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
808 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
809 const std::vector<Dune::cpgrid::DefaultGeometryPolicy>& refined_geometries_vec,
810 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
811 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
812 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
813 const std::vector<int>& assignRefinedLevel,
814 const int& preAdaptMaxLevel,
815 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
816 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
817 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
818
820 void setRefinedLevelGridsGeometries( /* Refined corner arguments */
821 std::vector<Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<0,3>>>& refined_corners_vec,
822 const std::vector<int>& refined_corner_count_vec,
823 /* Refined face arguments */
824 std::vector<Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<2,3>>>& refined_faces_vec,
825 std::vector<Dune::cpgrid::EntityVariableBase<enum face_tag>>& mutable_refined_face_tags_vec,
826 std::vector<Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>>& mutable_refine_face_normals_vec,
827 std::vector<Opm::SparseTable<int>>& refined_face_to_point_vec,
828 const std::vector<int>& refined_face_count_vec,
829 /* Refined cell argumets */
830 std::vector<Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<3,3>>>& refined_cells_vec,
831 std::vector<std::vector<std::array<int,8>>>& refined_cell_to_point_vec,
832 std::vector<std::vector<int>>& refined_global_cell_vec,
833 const std::vector<int>& refined_cell_count_vec,
834 std::vector<cpgrid::OrientedEntityTable<0,1>>& refined_cell_to_face_vec,
835 std::vector<cpgrid::OrientedEntityTable<1,0>>& refined_face_to_cell_vec,
836 /* Auxiliary arguments */
837 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
838 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
839 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner,
840 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
841 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
842 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
843 const std::vector<Dune::cpgrid::DefaultGeometryPolicy>& refined_geometries_vec,
844 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
845 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
846 const std::vector<int>& assignRefinedLevel,
847 const int& preAdaptMaxLevel,
848 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
849 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
850 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
851
853 void populateLeafGridCorners(Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<0,3>>& adapted_corners,
854 const int& corners_count,
855 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
856 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner) const;
857
859 void populateLeafGridFaces(Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<2,3>>& adapted_faces,
861 Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>& mutable_face_normals,
862 Opm::SparseTable<int>& adapted_face_to_point,
863 const int& face_count,
864 const std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
865 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
866 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
867 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
868 const std::vector<int>& assignRefinedLevel,
869 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
870 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
871 const std::vector<std::array<int,3>>& cells_per_dim_vec,
872 const int& preAdaptMaxLevel) const;
873
875 void populateLeafGridCells(Dune::cpgrid::EntityVariableBase<cpgrid::Geometry<3,3>>& adapted_cells,
876 std::vector<std::array<int,8>>& adapted_cell_to_point,
877 const int& cell_count,
878 cpgrid::OrientedEntityTable<0,1>& adapted_cell_to_face,
879 cpgrid::OrientedEntityTable<1,0>& adapted_face_to_cell,
880 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
881 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
882 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
883 const Dune::cpgrid::DefaultGeometryPolicy& adapted_geometries,
884 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
885 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
886 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
887 const std::vector<int>& assignRefinedLevel,
888 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
889 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
890 const std::vector<std::array<int,3>>& cells_per_dim_vec,
891 const int& preAdaptMaxLevel) const;
892
894 void updateLeafGridViewGeometries( /* Leaf grid View Corners arguments */
896 const int& corner_count,
897 /* Leaf grid View Faces arguments */
900 Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>& mutable_face_normals,
901 Opm::SparseTable<int>& adapted_face_to_point,
902 const int& face_count,
903 /* Leaf grid View Cells argumemts */
905 std::vector<std::array<int,8>>& adapted_cell_to_point,
906 const int& cell_count,
907 cpgrid::OrientedEntityTable<0,1>& adapted_cell_to_face,
908 cpgrid::OrientedEntityTable<1,0>& adapted_face_to_cell,
909 /* Auxiliary arguments */
910 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
911 const std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
912 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
913 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
914 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
915 const Dune::cpgrid::DefaultGeometryPolicy& adapted_geometries,
916 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
917 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
918 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
919 const std::vector<int>& assignRefinedLevel,
920 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
921 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
922 const std::vector<std::array<int,3>>& cells_per_dim_vec,
923 const int& preAdaptMaxLevel) const;
924
925 void updateCornerHistoryLevels(const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
926 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
927 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
928 const int& corner_count,
929 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
930 const int& preAdaptMaxLevel,
931 const int& newLevels);
932
933 void globalIdsPartitionTypesLgrAndLeafGrids(const std::vector<int>& assignRefinedLevel,
934 const std::vector<std::array<int,3>>& cells_per_dim_vec,
935 const std::vector<int>& lgr_with_at_least_one_active_cell);
936
943 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
944 public:
951 private:
952
965 void computeGlobalCellLgr(const int& level, const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
966
970 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
971
981 std::array<int,3> getRefinedCornerIJK(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr) const;
982
997 std::array<int,3> getRefinedFaceIJK(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
998 const std::shared_ptr<cpgrid::CpGridData>& elemLgr_ptr) const;
999
1004 bool isRefinedCornerInInteriorLgr(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr) const;
1005
1011 bool isRefinedFaceInInteriorLgr(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
1012 const std::shared_ptr<cpgrid::CpGridData>& elemLgr_ptr) const;
1013
1019 bool isRefinedNewBornCornerOnLgrBoundary(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr) const;
1020
1026 bool newRefinedCornerLiesOnEdge(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr) const;
1027
1033 bool isRefinedFaceOnLgrBoundary(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
1034 const std::shared_ptr<cpgrid::CpGridData>& elemLgr_ptr) const;
1035
1041 int getParentFaceWhereNewRefinedCornerLiesOn(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr, int elemLgr) const;
1042
1048 std::array<int,2> getParentFacesAssocWithNewRefinedCornLyingOnEdge(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr, int elemLgr) const;
1049
1050 protected:
1057 int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array<int,3>& cells_per_dim_lgr1, int cornerIdxLgr1, const std::array<int,3>& cells_per_dim_lgr2) const;
1058
1066 int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array<int,3>& cells_per_dim_lgr1, int cornerIdxLgr1, int elemLgr1, int parentFaceLastAppearanceIdx,
1067 const std::array<int,3>& cells_per_dim_lgr2) const;
1068
1076 int replaceLgr1FaceIdxByLgr2FaceIdx(const std::array<int,3>& cells_per_dim_lgr1, int faceIdxInLgr1,
1077 const std::shared_ptr<cpgrid::CpGridData>& elemLgr1_ptr,
1078 const std::array<int,3>& cells_per_dim_lgr2) const;
1079 private:
1080
1087 int getParentFaceWhereNewRefinedFaceLiesOn(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
1088 const std::shared_ptr<cpgrid::CpGridData>& elemLgr_ptr,
1089 int elemLgr) const;
1090
1092
1100 bool nonNNCsSelectedCellsLGR( const std::vector<std::array<int,3>>& startIJK_vec,
1101 const std::vector<std::array<int,3>>& endIJK_vec) const;
1102
1113 void detectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
1114 const std::vector<std::array<int,3>>& endIJK_vec,
1115 std::vector<int>& lgr_with_at_least_one_active_cell);
1116
1130 void markElemAssignLevelDetectActiveLgrs(const std::vector<std::array<int,3>>& startIJK_vec,
1131 const std::vector<std::array<int,3>>& endIJK_vec,
1132 std::vector<int>& assignRefinedLevel,
1133 std::vector<int>& lgr_with_at_least_one_active_cell);
1134
1140 template<class T>
1141 void computeOnLgrParents(const std::vector<std::array<int,3>>& startIJK_vec,
1142 const std::vector<std::array<int,3>>& endIJK_vec,
1143 T func);
1144
1158 void predictMinCellAndPointGlobalIdPerProcess(const std::vector<int>& assignRefinedLevel,
1159 const std::vector<std::array<int,3>>& cells_per_dim_vec,
1160 const std::vector<int>& lgr_with_at_least_one_active_cell,
1161 int& min_globalId_cell_in_proc,
1162 int& min_globalId_point_in_proc) const;
1163
1172 void assignCellIdsAndCandidatePointIds( std::vector<std::vector<int>>& localToGlobal_cells_per_level,
1173 std::vector<std::vector<int>>& localToGlobal_points_per_level,
1174 int min_globalId_cell_in_proc,
1175 int min_globalId_point_in_proc,
1176 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
1177
1191 void selectWinnerPointIds(std::vector<std::vector<int>>& localToGlobal_points_per_level,
1192 const std::vector<std::tuple<int,std::vector<int>>>& parent_to_children,
1193 const std::vector<std::array<int,3>>& cells_per_dim_vec) const;
1194
1196 void populateCellIndexSetRefinedGrid(int level);
1197
1199 void populateCellIndexSetLeafGridView();
1200
1202 void populateLeafGlobalIdSet();
1203
1204 public:
1205
1210 std::vector<std::unordered_map<std::size_t, std::size_t>> mapLocalCartesianIndexSetsToLeafIndexSet() const;
1211
1213 std::vector<std::array<int,2>> mapLeafIndexSetToLocalCartesianIndexSets() const;
1214
1216 unsigned int overlapSize(int) const;
1217
1218
1220 unsigned int ghostSize(int) const;
1221
1223 unsigned int overlapSize(int, int) const;
1224
1226 unsigned int ghostSize(int, int) const;
1227
1229 unsigned int numBoundarySegments() const;
1230
1231 void setPartitioningParams(const std::map<std::string,std::string>& params);
1232
1233 // loadbalance is not part of the grid interface therefore we skip it.
1234
1244 bool loadBalance(int overlapLayers=1,
1245 int partitionMethod = Dune::PartitionMethod::zoltan,
1246 double imbalanceTol = 1.1,
1247 int level =-1)
1248 {
1249 using std::get;
1250 return get<0>(scatterGrid(/* edgeWeightMethod = */ defaultTransEdgeWgt,
1251 /* ownersFirst = */ false,
1252 /* wells = */ nullptr,
1253 /* possibleFutureConnections = */ {},
1254 /* serialPartitioning = */ false,
1255 /* transmissibilities = */ nullptr,
1256 /* addCornerCells = */ true,
1257 overlapLayers,
1258 partitionMethod,
1259 imbalanceTol,
1260 /* allowDistributedWells = */ false,
1261 /* input_cell_part = */ {},
1262 level));
1263 }
1264
1265 // loadbalance is not part of the grid interface therefore we skip it.
1266
1277 bool loadBalanceSerial(int overlapLayers=1,
1278 int partitionMethod = Dune::PartitionMethod::zoltan,
1279 int edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt,
1280 double imbalanceTol = 1.1,
1281 int level = -1)
1282 {
1283 using std::get;
1284 return get<0>(scatterGrid(EdgeWeightMethod(edgeWeightMethod),
1285 /* ownersFirst = */ false,
1286 /* wells = */ nullptr,
1287 /* possibleFutureConnections = */ {},
1288 /* serialPartitioning = */ true,
1289 /* transmissibilities = */ nullptr,
1290 /* addCornerCells = */ true,
1291 overlapLayers,
1292 partitionMethod,
1293 imbalanceTol,
1294 /* allowDistributedWells = */ false,
1295 /* input_cell_part = */ {},
1296 level));
1297 }
1298
1299 // loadbalance is not part of the grid interface therefore we skip it.
1300
1329 std::pair<bool,std::vector<std::pair<std::string,bool>>>
1330 loadBalance(const std::vector<cpgrid::OpmWellType> * wells,
1331 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
1332 const double* transmissibilities = nullptr,
1333 int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG,
1334 int level = -1)
1335 {
1336 return scatterGrid(defaultTransEdgeWgt, /* ownersFirst = */ false, wells, possibleFutureConnections,
1337 /* serialPartitioning = */ false, transmissibilities, /* addCornerCells = */ false,
1338 overlapLayers, partitionMethod, /* imbalanceTol = */ 1.1, /* allowDistributeWells = */ false,
1339 /* input_cell_part = */ {}, level);
1340 }
1341
1342 // loadbalance is not part of the grid interface therefore we skip it.
1343
1377 std::pair<bool,std::vector<std::pair<std::string,bool>>>
1378 loadBalance(EdgeWeightMethod method, const std::vector<cpgrid::OpmWellType> * wells,
1379 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
1380 const double* transmissibilities = nullptr, bool ownersFirst=false,
1381 bool addCornerCells=false, int overlapLayers=1,
1382 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1383 double imbalanceTol = 1.1,
1384 int level = -1)
1385 {
1386 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections, /* serialPartitioning = */ false,
1387 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
1388 /* allowDistributeWells = */ false, /* input_cell_part = */ {}, level);
1389 }
1390
1419 template<class DataHandle>
1420 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1421 loadBalance(DataHandle& data,
1422 const std::vector<cpgrid::OpmWellType> * wells,
1423 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
1424 const double* transmissibilities = nullptr,
1425 int overlapLayers=1, int partitionMethod = 1, int level =-1)
1426 {
1427 auto ret = loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
1428 using std::get;
1429 if (get<0>(ret))
1430 {
1432 }
1433 return ret;
1434 }
1435
1470 template<class DataHandle>
1471 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1472 loadBalance(DataHandle& data, EdgeWeightMethod method,
1473 const std::vector<cpgrid::OpmWellType> * wells,
1474 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1475 bool serialPartitioning,
1476 const double* transmissibilities = nullptr, bool ownersFirst=false,
1477 bool addCornerCells=false, int overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1478 double imbalanceTol = 1.1,
1479 bool allowDistributedWells = false)
1480 {
1481 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
1482 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
1483 /* input_cell_parts = */ std::vector<int>{}, /* level = */ 0);
1484 using std::get;
1485 if (get<0>(ret))
1486 {
1488 }
1489 return ret;
1490 }
1491
1508 template<class DataHandle>
1509 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1510 loadBalance(DataHandle& data, const std::vector<int>& parts,
1511 const std::vector<cpgrid::OpmWellType> * wells,
1512 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
1513 bool ownersFirst=false,
1514 bool addCornerCells=false, int overlapLayers=1)
1515 {
1516 using std::get;
1517 auto ret = scatterGrid(defaultTransEdgeWgt, ownersFirst, wells,
1518 possibleFutureConnections,
1519 /* serialPartitioning = */ false,
1520 /* transmissibilities = */ {},
1521 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1522 /* imbalanceTol (ignored) = */ 0.0,
1523 /* allowDistributedWells = */ true, parts, /* level = */ 0);
1524 using std::get;
1525 if (get<0>(ret))
1526 {
1528 }
1529 return ret;
1530 }
1531
1539 template<class DataHandle>
1540 bool loadBalance(DataHandle& data,
1541 decltype(data.fixedSize(0,0)) overlapLayers=1, int partitionMethod = Dune::PartitionMethod::zoltan)
1542 {
1543 // decltype usage needed to tell the compiler not to use this function if first
1544 // argument is std::vector but rather loadbalance by parts
1545 bool ret = loadBalance(overlapLayers, partitionMethod);
1546 if (ret)
1547 {
1549 }
1550 return ret;
1551 }
1552
1564 bool loadBalance(const std::vector<int>& parts, bool ownersFirst=false,
1565 bool addCornerCells=false, int overlapLayers=1)
1566 {
1567 using std::get;
1568 return get<0>(scatterGrid(defaultTransEdgeWgt, ownersFirst, /* wells = */ {},
1569 {},
1570 /* serialPartitioning = */ false,
1571 /* trabsmissibilities = */ {},
1572 addCornerCells, overlapLayers, /* partitionMethod =*/ Dune::PartitionMethod::simple,
1573 /* imbalanceTol (ignored) = */ 0.0,
1574 /* allowDistributedWells = */ true, parts));
1575 }
1576
1589 template<class DataHandle>
1590 bool loadBalance(DataHandle& data, const std::vector<int>& parts, bool ownersFirst=false,
1591 bool addCornerCells=false, int overlapLayers=1)
1592 {
1593 bool ret = loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1594 if (ret)
1595 {
1597 }
1598 return ret;
1599 }
1600
1613 std::vector<int>
1614 zoltanPartitionWithoutScatter(const std::vector<cpgrid::OpmWellType>* wells,
1615 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1616 const double* transmissibilities,
1617 const int numParts,
1618 const double imbalanceTol) const;
1619
1627 template<class DataHandle>
1628 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int /*level*/) const
1629 {
1630 communicate(data, iftype, dir);
1631 }
1632
1640 template<class DataHandle>
1641 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const;
1642
1644 const typename CpGridTraits::Communication& comm () const;
1646
1647 // ------------ End of Dune interface, start of simplified interface --------------
1648
1654
1655 // enum { dimension = 3 }; // already defined
1656
1657 typedef Dune::FieldVector<double, 3> Vector;
1658
1659
1660 const std::vector<double>& zcornData() const;
1661
1662
1663 // Topology
1668 int numCells(int level = -1) const;
1669
1674 int numFaces(int level = -1) const;
1675
1677 int numVertices() const;
1678
1679
1688 int numCellFaces(int cell, int level = -1) const;
1689
1696 int cellFace(int cell, int local_index, int level = -1) const;
1697
1701
1714 int faceCell(int face, int local_index, int level = -1) const;
1715
1722 int numCellFaces() const;
1723
1724 int numFaceVertices(int face) const;
1725
1730 int faceVertex(int face, int local_index) const;
1731
1734 double cellCenterDepth(int cell_index) const;
1735
1736
1737 const Vector faceCenterEcl(int cell_index, int face, const Dune::cpgrid::Intersection& intersection) const;
1738
1739 const Vector faceAreaNormalEcl(int face) const;
1740
1741
1742 // Geometry
1746 const Vector& vertexPosition(int vertex) const;
1747
1750 double faceArea(int face) const;
1751
1754 const Vector& faceCentroid(int face) const;
1755
1759 const Vector& faceNormal(int face) const;
1760
1763 double cellVolume(int cell) const;
1764
1767 const Vector& cellCentroid(int cell) const;
1768
1771 template<int codim>
1773 : public RandomAccessIteratorFacade<CentroidIterator<codim>,
1774 FieldVector<double, 3>,
1775 const FieldVector<double, 3>&, int>
1776 {
1777 public:
1779 typedef typename std::vector<cpgrid::Geometry<3-codim, 3> >::const_iterator
1784 : iter_(iter)
1785 {}
1786
1787 const FieldVector<double,3>& dereference() const
1788 {
1789 return iter_->center();
1790 }
1792 {
1793 ++iter_;
1794 }
1795 const FieldVector<double,3>& elementAt(int n)
1796 {
1797 return iter_[n]->center();
1798 }
1799 void advance(int n){
1800 iter_+=n;
1801 }
1803 {
1804 --iter_;
1805 }
1807 {
1808 return o-iter_;
1809 }
1810 bool equals(const CentroidIterator& o) const{
1811 return o==iter_;
1812 }
1813 private:
1815 GeometryIterator iter_;
1816 };
1817
1820
1823
1824 // Extra
1825 int boundaryId(int face) const;
1826
1833 template<class Cell2FacesRowIterator>
1834 int
1835 faceTag(const Cell2FacesRowIterator& cell_face) const;
1836
1838
1839 // ------------ End of simplified interface --------------
1840
1841 //------------- methods not in the DUNE grid interface.
1842
1847
1848
1857 template<class DataHandle>
1858 void scatterData(DataHandle& handle) const;
1859
1866 template<class DataHandle>
1867 void gatherData(DataHandle& handle) const;
1868
1871
1901
1905
1908
1912
1913#if HAVE_MPI
1918
1921
1926
1928
1930
1932
1934#endif
1935
1937 const std::vector<int>& sortedNumAquiferCells() const;
1938
1939 private:
1971 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1972 scatterGrid(EdgeWeightMethod method,
1973 bool ownersFirst,
1974 const std::vector<cpgrid::OpmWellType> * wells,
1975 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1976 bool serialPartitioning,
1977 const double* transmissibilities,
1978 bool addCornerCells,
1979 int overlapLayers,
1980 int partitionMethod = Dune::PartitionMethod::zoltanGoG,
1981 double imbalanceTol = 1.1,
1982 bool allowDistributedWells = true,
1983 const std::vector<int>& input_cell_part = {},
1984 int level = -1);
1985
1990 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1992 cpgrid::CpGridData* current_view_data_;
1994 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1996 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1998 std::map<std::string,int> lgr_names_ = {{"GLOBAL", 0}};
2004 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
2005 /*
2006 * @brief Interface for scattering and gathering point data.
2007 *
2008 * @warning Will only update owner cells
2009 */
2010 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
2014 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
2015
2016
2020 std::map<std::string,std::string> partitioningParams;
2021
2022 }; // end Class CpGrid
2023
2024} // end namespace Dune
2025
2029
2030
2031namespace Dune
2032{
2033
2034 namespace Capabilities
2035 {
2037 template <>
2038 struct hasEntity<CpGrid, 0>
2039 {
2040 static const bool v = true;
2041 };
2042
2044 template <>
2045 struct hasEntity<CpGrid, 3>
2046 {
2047 static const bool v = true;
2048 };
2049
2050 template<>
2051 struct canCommunicate<CpGrid,0>
2052 {
2053 static const bool v = true;
2054 };
2055
2056 template<>
2057 struct canCommunicate<CpGrid,3>
2058 {
2059 static const bool v = true;
2060 };
2061
2063 template <>
2064 struct hasBackupRestoreFacilities<CpGrid>
2065 {
2066 static const bool v = false;
2067 };
2068
2069 }
2070
2071 template<class DataHandle>
2072 void CpGrid::communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
2073 {
2074 current_view_data_->communicate(data, iftype, dir);
2075 }
2076
2077
2078 template<class DataHandle>
2079 void CpGrid::scatterData(DataHandle& handle) const
2080 {
2081#if HAVE_MPI
2082 if(distributed_data_.empty())
2083 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balanced grid!");
2084 distributed_data_[0]->scatterData(handle, data_[0].get(), distributed_data_[0].get(), cellScatterGatherInterface(),
2086#else
2087 // Suppress warnings for unused argument.
2088 (void) handle;
2089#endif
2090 }
2091
2092 template<class DataHandle>
2093 void CpGrid::gatherData(DataHandle& handle) const
2094 {
2095#if HAVE_MPI
2096 if(distributed_data_.empty())
2097 OPM_THROW(std::runtime_error, "Moving Data only allowed with a load balance grid!");
2098 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
2099#else
2100 // Suppress warnings for unused argument.
2101 (void) handle;
2102#endif
2103 }
2104
2105
2106 template<class Cell2FacesRowIterator>
2107 int
2108 CpGrid::faceTag(const Cell2FacesRowIterator& cell_face) const
2109 {
2110 // Note that this relies on the following implementation detail:
2111 // The grid is always constructed such that the interior faces constructed
2112 // with orientation set to true are
2113 // oriented along the positive IJK direction. Oriented means that
2114 // the first cell attached to face has the lower index.
2115 // For faces along the boundary (only one cell, always attached at index 0)
2116 // the orientation has to be determined by the orientation of the cell.
2117 // If it is true then in UnstructuredGrid it would be stored at index 0,
2118 // otherwise at index 1.
2119 const int cell = cell_face.getCellIndex();
2120 const int face = *cell_face;
2121 assert (0 <= cell); assert (cell < numCells());
2122 assert (0 <= face); assert (face < numFaces());
2123
2125
2126 const cpgrid::EntityRep<1> f(face, true);
2127 const F2C& f2c = current_view_data_->face_to_cell_[f];
2128 const face_tag tag = current_view_data_->face_tag_[f];
2129
2130 assert ((f2c.size() == 1) || (f2c.size() == 2));
2131
2132 int inside_cell = 0;
2133
2134 if ( f2c.size() == 2 ) // Two cells => interior
2135 {
2136 if ( f2c[1].index() == cell )
2137 {
2138 inside_cell = 1;
2139 }
2140 }
2141 const bool normal_is_in = ! f2c[inside_cell].orientation();
2142
2143 switch (tag) {
2144 case I_FACE:
2145 // LEFT : RIGHT
2146 return normal_is_in ? 0 : 1; // min(I) : max(I)
2147 case J_FACE:
2148 // BACK : FRONT
2149 return normal_is_in ? 2 : 3; // min(J) : max(J)
2150 case K_FACE:
2151 // Note: TOP at min(K) as 'z' measures *depth*.
2152 // TOP : BOTTOM
2153 return normal_is_in ? 4 : 5; // min(K) : max(K)
2154 case NNC_FACE:
2155 // For nnc faces we return the otherwise unused value -1.
2156 return -1;
2157 default:
2158 OPM_THROW(std::logic_error, "Unhandled face tag. This should never happen!");
2159 }
2160 }
2161
2162 template<int dim>
2164
2165} // namespace Dune
2166
2169#include "cpgrid/Intersection.hpp"
2170#include "cpgrid/Geometry.hpp"
2171#include "cpgrid/Indexsets.hpp"
2172
2173#endif // OPM_CPGRID_HEADER
DataHandle & data
Definition: CpGridData.hpp:1250
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:1776
void increment()
Definition: CpGrid.hpp:1791
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:1783
const FieldVector< double, 3 > & dereference() const
Definition: CpGrid.hpp:1787
int distanceTo(const CentroidIterator &o)
Definition: CpGrid.hpp:1806
void decrement()
Definition: CpGrid.hpp:1802
void advance(int n)
Definition: CpGrid.hpp:1799
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:1780
bool equals(const CentroidIterator &o) const
Definition: CpGrid.hpp:1810
const FieldVector< double, 3 > & elementAt(int n)
Definition: CpGrid.hpp:1795
[ provides Dune::Grid ]
Definition: CpGrid.hpp:198
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:1510
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:213
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:2093
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:1472
cpgrid::CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGrid.hpp:1915
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:2108
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:1378
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
void processEclipseFormat(const grdecl &input_data, bool remove_ij_boundary, bool turn_normals=false)
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:1870
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:1421
int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array< int, 3 > &cells_per_dim_lgr1, int cornerIdxLgr1, const std::array< int, 3 > &cells_per_dim_lgr2) const
A refined corner appears in two single-cell-refinements. Given the corner index in the first single-c...
Traits::template Codim< codim >::LevelIterator lend(int level) const
one past the end on this level
int numCellFaces() const
Get the sum of all faces attached to all cells.
const std::vector< int > & globalCell() const
bool loadBalanceSerial(int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan, int edgeWeightMethod=Dune::EdgeWeightMethod::defaultTransEdgeWgt, double imbalanceTol=1.1, int level=-1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1277
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:1244
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
const CpGridTraits::Communication & comm() const
Get the collective communication object.
double faceArea(int face) const
Get the area of a face.
const std::vector< Dune::GeometryType > & geomTypes(const int) const
void postAdapt()
Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'.
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
int size(int level, int codim) const
Number of grid entities per level and codim.
const std::array< int, 3 > & logicalCartesianSize() const
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
unsigned int ghostSize(int, int) const
Size of the ghost cell layer on a given level.
cpgrid::CpGridDataTraits::RemoteIndices RemoteIndices
The type of the remote indices information.
Definition: CpGrid.hpp:1917
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
int numVertices() const
Get The number of vertices.
Dune::cpgrid::Intersection getParentIntersectionFromLgrBoundaryFace(const Dune::cpgrid::Intersection &intersection) const
void getIJK(const int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
void syncDistributedGlobalCellIds()
Synchronizes cell global ids across processes after load balancing.
bool uniqueBoundaryIds() const
int faceCell(int face, int local_index, int level=-1) const
Get the index identifying a cell attached to a face.
void addLgrsUpdateLeafView(const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, const std::vector< std::string > &lgr_name_vec)
Create a grid out of a coarse one and (at most) 2 refinements(LGRs) of selected block-shaped disjoint...
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
void setPartitioningParams(const std::map< std::string, std::string > &params)
int replaceLgr1FaceIdxByLgr2FaceIdx(const std::array< int, 3 > &cells_per_dim_lgr1, int faceIdxInLgr1, const std::shared_ptr< cpgrid::CpGridData > &elemLgr1_ptr, const std::array< int, 3 > &cells_per_dim_lgr2) const
A new refined face lays on the boudndary of a single-cell-refinement appears in at most two single-ce...
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
CpGrid(MPIHelper::MPICommunicator comm)
const ParallelIndexSet & getCellIndexSet() const
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData()
Returns either data_ or distributed_data_(if non empty).
int boundaryId(int face) const
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim and PartitionIteratorType.
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
Mark entity for refinement (or coarsening).
bool loadBalance(DataHandle &data, decltype(data.fixedSize(0, 0)) overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan)
Distributes this grid and data over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1540
unsigned int overlapSize(int, int) const
Size of the overlap on a given level.
int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array< int, 3 > &cells_per_dim_lgr1, int cornerIdxLgr1, int elemLgr1, int parentFaceLastAppearanceIdx, const std::array< int, 3 > &cells_per_dim_lgr2) const
A new refined corner lays on an edge and appears in at least two single-cell-refinements....
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:1330
friend cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
const Vector & faceNormal(int face) const
Get the unit normal of a face.
Dune::FieldVector< double, 3 > Vector
Definition: CpGrid.hpp:1657
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:1590
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
const cpgrid::OrientedEntityTable< 0, 1 >::row_type cellFaceRow(int cell) const
Get a list of indices identifying all faces of a cell.
bool loadBalance(const std::vector< int > &parts, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:1564
std::vector< int > zoltanPartitionWithoutScatter(const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const double *transmissibilities, const int numParts, const double imbalanceTol) const
Partitions the grid using Zoltan without decomposing and distributing it among processes.
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
std::array< double, 3 > getEclCentroid(const int &idx) const
int numCells(int level=-1) const
Get the number of cells.
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
void globalRefine(int refCount)
Refine the grid refCount times using the default refinement rule. This behaves like marking all eleme...
bool adapt(const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< int > &assignRefinedLevel, const std::vector< std::string > &lgr_name_vec, const std::vector< std::array< int, 3 > > &startIJK_vec=std::vector< std::array< int, 3 > >{}, const std::vector< std::array< int, 3 > > &endIJK_vec=std::vector< std::array< int, 3 > >{})
Triggers the grid refinement process, allowing to select diffrent refined level grids.
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
communicate objects for all codims on a given level
Definition: CpGrid.hpp:1628
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:1920
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.
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:2079
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:138
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:1028
Definition: DefaultGeometryPolicy.hpp:53
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:99
Base class for EntityVariable and SignedEntityVariable. Forwards a restricted subset of the std::vect...
Definition: EntityRep.hpp:219
Definition: Geometry.hpp:75
The global id set for Dune.
Definition: Indexsets.hpp:450
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:108
Definition: Indexsets.hpp:56
Definition: Intersection.hpp:278
Definition: Intersection.hpp:66
Definition: Iterators.hpp:60
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
cpgrid::Entity< dim > createEntity(const CpGrid &, int, bool)
@ simple
Use simple approach based on rectangular partitioning the underlying cartesian grid.
Definition: GridEnums.hpp:46
@ zoltanGoG
use Zoltan on GraphOfGrid for partitioning
Definition: GridEnums.hpp:52
@ zoltan
Use Zoltan for partitioning.
Definition: GridEnums.hpp:48
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Me...
Definition: GridEnums.hpp:34
@ defaultTransEdgeWgt
Use the transmissibilities as edge weights.
Definition: GridEnums.hpp:38
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:185
CpGridTraits Traits
Definition: CpGrid.hpp:186
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:139
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition: CpGrid.hpp:141
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition: CpGrid.hpp:143
Traits associated with a specific codim.
Definition: CpGrid.hpp:115
cpgrid::Entity< cd > Entity
The type of the entity.
Definition: CpGrid.hpp:124
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity. IMPORTANT: Codim<codim>::Geometry == Geometry<di...
Definition: CpGrid.hpp:118
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition: CpGrid.hpp:121
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition: CpGrid.hpp:130
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition: CpGrid.hpp:127
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition: CpGrid.hpp:133
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:151
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:155
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:153
Definition: CpGrid.hpp:95
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition: CpGrid.hpp:165
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition: CpGrid.hpp:104
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:162
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:160
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition: CpGrid.hpp:169
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition: CpGrid.hpp:106
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition: CpGrid.hpp:167
cpgrid::CpGridDataTraits::CollectiveCommunication CollectiveCommunication
Definition: CpGrid.hpp:175
GlobalIdSet LocalIdSet
The type of the local id set.
Definition: CpGrid.hpp:171
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGrid.hpp:174
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition: CpGrid.hpp:109
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition: CpGrid.hpp:102
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition: CpGrid.hpp:100
CpGrid Grid
The type that implements the grid.
Definition: CpGrid.hpp:97
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