CpGridData.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: CpGridData.hpp
4//
5// Created: Sep 17 21:11:41 2013
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// Bård Skaflestad <bard.skaflestad@sintef.no>
9// Markus Blatt <markus@dr-blatt.de>
10// Antonella Ritorto <antonella.ritorto@opm-op.com>
11//
12// Comment: Major parts of this file originated in dune/grid/CpGrid.hpp
13// and got transfered here during refactoring for the parallelization.
14//
15// $Date$
16//
17// $Revision$
18//
19//===========================================================================
20
21/*
22 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
23 Copyright 2009, 2010, 2013, 2022-2023 Equinor ASA.
24 Copyright 2013 Dr. Blatt - HPC-Simulation-Software & Services
25
26 This file is part of The Open Porous Media project (OPM).
27
28 OPM is free software: you can redistribute it and/or modify
29 it under the terms of the GNU General Public License as published by
30 the Free Software Foundation, either version 3 of the License, or
31 (at your option) any later version.
32
33 OPM is distributed in the hope that it will be useful,
34 but WITHOUT ANY WARRANTY; without even the implied warranty of
35 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
36 GNU General Public License for more details.
37
38 You should have received a copy of the GNU General Public License
39 along with OPM. If not, see <http://www.gnu.org/licenses/>.
40*/
48#ifndef OPM_CPGRIDDATA_HEADER
49#define OPM_CPGRIDDATA_HEADER
50
51
52#include <dune/common/parallel/mpihelper.hh>
53#ifdef HAVE_DUNE_ISTL
54#include <dune/istl/owneroverlapcopy.hh>
55#endif
56
57#include <dune/common/parallel/communication.hh>
58#include <dune/common/parallel/variablesizecommunicator.hh>
59#include <dune/grid/common/gridenums.hh>
60
61#if HAVE_ECL_INPUT
62#include <opm/input/eclipse/EclipseState/Grid/NNC.hpp>
63#endif
64
66
68#include "CpGridDataTraits.hpp"
69//#include "DataHandleWrappers.hpp"
70//#include "GlobalIdMapping.hpp"
71#include "Geometry.hpp"
72
73#include <array>
74#include <initializer_list>
75#include <set>
76#include <vector>
77
78namespace Opm
79{
80class EclipseState;
81}
82namespace Dune
83{
84class CpGrid;
85
86namespace cpgrid
87{
88
89class IndexSet;
90class IdSet;
91class LevelGlobalIdSet;
92class PartitionTypeIndicator;
93template<int,int> class Geometry;
94template<int> class Entity;
95template<int> class EntityRep;
96}
97}
98
100 const std::array<int, 3>&,
101 bool);
102
103void refinePatch_and_check(const std::array<int,3>&,
104 const std::array<int,3>&,
105 const std::array<int,3>&);
106
108 const std::vector<std::array<int,3>>&,
109 const std::vector<std::array<int,3>>&,
110 const std::vector<std::array<int,3>>&,
111 const std::vector<std::string>&);
112
114 const Dune::CpGrid&);
115
116void fieldProp_check(const Dune::CpGrid& grid, Opm::EclipseGrid eclGrid, const std::string& deck_string);
117
118void testInactiveCellsLgrs(const std::string&,
119 const std::vector<std::array<int,3>>&,
120 const std::vector<std::array<int,3>>&,
121 const std::vector<std::array<int,3>>&,
122 const std::vector<std::string>&);
123
124namespace Dune
125{
126namespace cpgrid
127{
128namespace mover
129{
130template<class T, int i> struct Mover;
131}
132
138{
139 template<class T, int i> friend struct mover::Mover;
140 friend class GlobalIdSet;
141 friend class HierarchicIterator;
145
146 friend
148 const std::array<int, 3>&,
149 bool);
150 friend
151 void ::refinePatch_and_check(const std::array<int,3>&,
152 const std::array<int,3>&,
153 const std::array<int,3>&);
154
155 friend
157 const std::vector<std::array<int,3>>&,
158 const std::vector<std::array<int,3>>&,
159 const std::vector<std::array<int,3>>&,
160 const std::vector<std::string>&);
161
162 friend
164 const Dune::CpGrid&);
165 friend
166 void ::fieldProp_check(const Dune::CpGrid& grid, Opm::EclipseGrid eclGrid, const std::string& deck_string);
167 friend
168 void ::testInactiveCellsLgrs(const std::string&,
169 const std::vector<std::array<int,3>>&,
170 const std::vector<std::array<int,3>>&,
171 const std::vector<std::array<int,3>>&,
172 const std::vector<std::string>&);
173
174private:
175 CpGridData(const CpGridData& g);
176
177public:
178 enum{
179#ifndef MAX_DATA_COMMUNICATED_PER_ENTITY
187#else
192 MAX_DATA_PER_CELL = MAX_DATA_COMMUNICATED_PER_ENTITY
193#endif
194 };
195
196 CpGridData() = delete;
197
201 explicit CpGridData(MPIHelper::MPICommunicator comm, std::vector<std::shared_ptr<CpGridData>>& data);
202
203
204
206 explicit CpGridData(std::vector<std::shared_ptr<CpGridData>>& data);
209
210
211
212
214 int size(int codim) const;
215
217 int size (GeometryType type) const
218 {
219 if (type.isCube()) {
220 return size(3 - type.dim());
221 } else {
222 return 0;
223 }
224 }
225
231 void readEclipseFormat(const std::string& filename, bool periodic_extension, bool turn_normals = false);
232
233#if HAVE_ECL_INPUT
242 void processEclipseFormat(const Opm::Deck& deck, bool periodic_extension, bool turn_normals = false, bool clip_z = false,
243 const std::vector<double>& poreVolume = std::vector<double>());
244
258 std::vector<std::size_t>
259 processEclipseFormat(const Opm::EclipseGrid* ecl_grid, Opm::EclipseState* ecl_state,
260 bool periodic_extension, bool turn_normals = false, bool clip_z = false,
261 bool pinchActive = true);
262#endif
263
275 void processEclipseFormat(const grdecl& input_data,
276#if HAVE_ECL_INPUT
277 Opm::EclipseState* ecl_state,
278#endif
279 std::array<std::set<std::pair<int, int>>, 2>& nnc,
280 bool remove_ij_boundary, bool turn_normals, bool pinchActive,
281 double tolerance_unique_points);
282
290 void getIJK(int c, std::array<int,3>& ijk) const
291 {
292 ijk = getIJK(global_cell_[c], logical_cartesian_size_);
293 }
294
295 int cellFace(int cell, int local_index) const
296 {
297 return cell_to_face_[cpgrid::EntityRep<0>(cell, true)][local_index].index();
298 }
299
300 int faceToCellSize(int face) const {
301 Dune::cpgrid::EntityRep<1> faceEntity(face, true);
302 return face_to_cell_[faceEntity].size();
303 }
304
311 const std::vector<int>& globalCell() const
312 {
313 return global_cell_;
314 }
315
322 std::array<int,3> getIJK(int idx_in_parent_cell, const std::array<int,3>& cells_per_dim) const
323 {
324 // idx = k*cells_per_dim_[0]*cells_per_dim_[1] + j*cells_per_dim_[0] + i
325 // with 0<= i < cells_per_dim_[0], 0<= j < cells_per_dim_[1], 0<= k <cells_per_dim_[2].
326 assert(cells_per_dim[0]);
327 assert(cells_per_dim[1]);
328 assert(cells_per_dim[2]);
329
330 std::array<int,3> ijk = {0,0,0};
331 ijk[0] = idx_in_parent_cell % cells_per_dim[0]; idx_in_parent_cell /= cells_per_dim[0];
332 ijk[1] = idx_in_parent_cell % cells_per_dim[1];
333 ijk[2] = idx_in_parent_cell /cells_per_dim[1];
334 return ijk;
335 }
336
342 bool disjointPatches(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
343
351 std::vector<int>
352 getPatchesCells(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
353
356 bool hasNNCs(const std::vector<int>& cellIndices) const;
357
364 void validStartEndIJKs(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
365
367 void checkCuboidShape(const std::vector<int>& cellIdx_vec) const;
368
374 bool patchesShareFace(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
375
376 int sharedFaceTag(const std::vector<std::array<int,3>>& startIJK_2Patches, const std::vector<std::array<int,3>>& endIJK_2Patches) const;
377
378
391 bool mark(int refCount, const cpgrid::Entity<0>& element);
392
396 int getMark(const cpgrid::Entity<0>& element) const;
397
407 bool preAdapt();
408
410 bool adapt();
411
413 void postAdapt();
414
415private:
431 bool compatibleSubdivisions(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) const;
434
435 std::array<Dune::FieldVector<double,3>,8> getReferenceRefinedCorners(int idx_in_parent_cell, const std::array<int,3>& cells_per_dim) const;
436
444 std::array<int,3> getPatchDim(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
445
453 std::vector<int> getPatchCorners(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
454
462 std::vector<int> getPatchFaces(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
463
471 std::vector<int> getPatchCells(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
472
480 std::vector<int> getPatchBoundaryCorners(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
481
489 std::array<std::vector<int>,6> getBoundaryPatchFaces(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
490
492 std::array<std::vector<double>,3> getWidthsLengthsHeights(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
493
508 Geometry<3,3> cellifyPatch(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK,
509 const std::vector<int>& patch_cells, DefaultGeometryPolicy& cellifiedPatch_geometry,
510 std::array<int,8>& cellifiedPatch_to_point,
511 std::array<int,8>& allcorners_cellifiedPatch) const;
512
513 // @brief Compute the average of array<double,3>.
514 //
515 // @param [in] vector of array<double,3>
516 // @return array<double,3> (average of the entries of the given vector).
517 std::array<double,3> getAverageArr(const std::vector<std::array<double,3>>& vec) const;
518
519public:
521 int getGridIdx() const {
522 // Not the nicest way of checking if "this" points at the leaf grid view of a mixed grid (with coarse and refined cells).
523 // 1. When the grid has been refined at least onece, level_data_ptr_ ->size() >1. Therefore, there is a chance of "this" pointing at the leaf grid view.
524 // 2. Unfortunately, level_ is default initialized by 0. This implies, in particular, that if someone wants to check the value of
525 // "this->level_" when "this" points at the leaf grid view of a grid that has been refined, this value is - unfortunately - equal to 0.
526 // 3. Due to 2. we need an extra bool value to distinguish between the actual level 0 grid and such a leaf grid view (with incorrect level_ == 0). For this
527 // reason we check if child_to_parent_cells_.empty() [true for actual level 0 grid, false for the leaf grid view].
528 // --- TO BE IMPROVED ---
529 if ((level_data_ptr_ ->size() >1) && (level_ == 0) && (!child_to_parent_cells_.empty())) {
530 return level_data_ptr_->size() -1;
531 }
532 return level_;
533 }
535 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& levelData() const
536 {
537 if (level_data_ptr_->empty()) {
538 OPM_THROW(std::logic_error, "Level data has not been initialized\n");
539 }
540 return *level_data_ptr_;
541 }
542
550 const std::tuple<int,std::vector<int>>& getChildrenLevelAndIndexList(int elemIdx) const {
551 return parent_to_children_cells_[elemIdx];
552 }
553
554 const std::vector<std::tuple<int,std::vector<int>>>& getParentToChildren() const {
555 return parent_to_children_cells_;
556 }
578 std::tuple< const std::shared_ptr<CpGridData>,
579 const std::vector<std::array<int,2>>, // parent_to_refined_corners(~boundary_old_to_new_corners)
580 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_faces (~boundary_old_to_new_faces)
581 const std::tuple<int, std::vector<int>>, // parent_to_children_cells
582 const std::vector<std::array<int,2>>, // child_to_parent_faces
583 const std::vector<std::array<int,2>>> // child_to_parent_cells
584 refineSingleCell(const std::array<int,3>& cells_per_dim, const int& parent_idx) const;
585
598 std::tuple< std::shared_ptr<CpGridData>,
599 const std::vector<std::array<int,2>>, // boundary_old_to_new_corners
600 const std::vector<std::tuple<int,std::vector<int>>>, // boundary_old_to_new_faces
601 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_faces
602 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_cell
603 const std::vector<std::array<int,2>>, // child_to_parent_faces
604 const std::vector<std::array<int,2>>> // child_to_parent_cells
605 refinePatch(const std::array<int,3>& cells_per_dim, const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
606
607 // @breif Compute center of an entity/element/cell in the Eclipse way:
608 // - Average of the 4 corners of the bottom face.
609 // - Average of the 4 corners of the top face.
610 // Return average of the previous computations.
611 // @param [in] int Index of a cell.
612 // @return 'eclipse centroid'
613 std::array<double,3> computeEclCentroid(const int idx) const;
614
615 // @breif Compute center of an entity/element/cell in the Eclipse way:
616 // - Average of the 4 corners of the bottom face.
617 // - Average of the 4 corners of the top face.
618 // Return average of the previous computations.
619 // @param [in] Entity<0> Entity
620 // @return 'eclipse centroid'
621 std::array<double,3> computeEclCentroid(const Entity<0>& elem) const;
622
623 // Make unique boundary ids for all intersections.
625
629 bool uniqueBoundaryIds() const
630 {
631 return use_unique_boundary_ids_;
632 }
633
636 void setUniqueBoundaryIds(bool uids)
637 {
638 use_unique_boundary_ids_ = uids;
639 if (use_unique_boundary_ids_ && unique_boundary_ids_.empty()) {
641 }
642 }
643
647 const std::vector<double>& zcornData() const {
648 return zcorn;
649 }
650
651
654 const IndexSet& indexSet() const
655 {
656 return *index_set_;
657 }
658
661 {
662 return *local_id_set_;
663 }
664
667 {
668 return *global_id_set_;
669 }
670
674 const std::array<int, 3>& logicalCartesianSize() const
675 {
676 return logical_cartesian_size_;
677 }
678
683 const CpGridData& view_data,
684 const std::vector<int>& cell_part);
685
691 template<class DataHandle>
692 void communicate(DataHandle& data, InterfaceType iftype, CommunicationDirection dir);
693
695
697
698 void computeCommunicationInterfaces(int noexistingPoints);
699
705
708#if HAVE_MPI
711
714
717
720
723
728 {
729 return cell_comm_;
730 }
731
736 {
737 return cell_comm_;
738 }
739
741 {
742 return cellCommunication().indexSet();
743 }
744
746 {
747 return cellCommunication().indexSet();
748 }
749
751 {
752 return cellCommunication().remoteIndices();
753 }
754
756 {
757 return cellCommunication().remoteIndices();
758 }
759#endif
760
762 const std::vector<int>& sortedNumAquiferCells() const
763 {
764 return aquifer_cells_;
765 }
766
767private:
768
770 void populateGlobalCellIndexSet();
771
772#if HAVE_MPI
773
779 template<class DataHandle>
780 void gatherData(DataHandle& data, CpGridData* global_view,
781 CpGridData* distributed_view);
782
783
790 template<int codim, class DataHandle>
791 void gatherCodimData(DataHandle& data, CpGridData* global_data,
792 CpGridData* distributed_data);
793
800 template<class DataHandle>
801 void scatterData(DataHandle& data, const CpGridData* global_data,
802 const CpGridData* distributed_data, const InterfaceMap& cell_inf,
803 const InterfaceMap& point_inf);
804
812 template<int codim, class DataHandle>
813 void scatterCodimData(DataHandle& data, CpGridData* global_data,
814 CpGridData* distributed_data);
815
824 template<int codim, class DataHandle>
825 void communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
826 const Interface& interface);
827
836 template<int codim, class DataHandle>
837 void communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
838 const InterfaceMap& interface);
839
840#endif
841
842 void computeGeometry(const CpGrid& grid,
843 const DefaultGeometryPolicy& globalGeometry,
844 const std::vector<int>& globalAquiferCells,
845 const OrientedEntityTable<0, 1>& globalCell2Faces,
846 DefaultGeometryPolicy& geometry,
847 std::vector<int>& aquiferCells,
849 const std::vector< std::array<int,8> >& cell2Points);
850
851 // Representing the topology
863 Opm::SparseTable<int> face_to_point_;
865 std::vector< std::array<int,8> > cell_to_point_;
872 std::array<int, 3> logical_cartesian_size_{};
879 std::vector<int> global_cell_;
885 typedef FieldVector<double, 3> PointType;
889 cpgrid::EntityVariable<int, 1> unique_boundary_ids_;
891 std::unique_ptr<cpgrid::IndexSet> index_set_;
893 std::shared_ptr<const cpgrid::IdSet> local_id_set_;
895 std::shared_ptr<LevelGlobalIdSet> global_id_set_;
897 std::shared_ptr<PartitionTypeIndicator> partition_type_indicator_;
899 std::vector<int> mark_;
901 int level_{0};
903 std::vector<std::shared_ptr<CpGridData>>* level_data_ptr_;
904 // SUITABLE FOR ALL LEVELS EXCEPT FOR LEAFVIEW
906 std::vector<int> level_to_leaf_cells_; // In entry 'level cell index', we store 'leafview cell index' // {level LGR, {child0, child1, ...}}
908 std::vector<std::tuple<int,std::vector<int>>> parent_to_children_cells_; // {# children in x-direction, ... y-, ... z-}
910 std::array<int,3> cells_per_dim_;
911 // SUITABLE ONLY FOR LEAFVIEW // {level, cell index in that level}
913 std::vector<std::array<int,2>> leaf_to_level_cells_;
915 std::vector<std::array<int,2>> corner_history_;
916 // SUITABLE FOR ALL LEVELS INCLUDING LEAFVIEW // {level parent cell, parent cell index}
918 std::vector<std::array<int,2>> child_to_parent_cells_;
921 std::vector<int> cell_to_idxInParentCell_;
922
923
924
926 Communication ccobj_;
927
928 // Boundary information (optional).
929 bool use_unique_boundary_ids_;
930
936 std::vector<double> zcorn;
937
939 std::vector<int> aquifer_cells_;
940
941#if HAVE_MPI
942
944 CommunicationType cell_comm_;
945
947 std::tuple<Interface,Interface,Interface,Interface,Interface> cell_interfaces_;
948 /*
949 // code deactivated, because users cannot access face indices and therefore
950 // communication on faces makes no sense!
952 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
953 face_interfaces_;
954 */
956 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
957 point_interfaces_;
958
959#endif
960
961 // Return the geometry vector corresponding to the given codim.
962 template <int codim>
963 const EntityVariable<Geometry<3 - codim, 3>, codim>& geomVector() const
964 {
965 return geometry_.geomVector<codim>();
966 }
967
968 friend class Dune::CpGrid;
969 template<int> friend class Entity;
970 template<int> friend class EntityRep;
971 friend class Intersection;
973};
974
975
976
977#if HAVE_MPI
978
979namespace
980{
985template<class T>
986T& getInterface(InterfaceType iftype,
987 std::tuple<T,T,T,T,T>& interfaces)
988{
989 switch(iftype)
990 {
991 case 0:
992 return std::get<0>(interfaces);
993 case 1:
994 return std::get<1>(interfaces);
995 case 2:
996 return std::get<2>(interfaces);
997 case 3:
998 return std::get<3>(interfaces);
999 case 4:
1000 return std::get<4>(interfaces);
1001 }
1002 OPM_THROW(std::runtime_error, "Invalid Interface type was used during communication");
1003}
1004
1005} // end unnamed namespace
1006
1007template<int codim, class DataHandle>
1008void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
1009 const Interface& interface)
1010{
1011 this->template communicateCodim<codim>(data, dir, interface.interfaces());
1012}
1013
1014template<int codim, class DataHandle>
1015void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data_wrapper, CommunicationDirection dir,
1016 const InterfaceMap& interface)
1017{
1018 Communicator comm(ccobj_, interface);
1019
1020 if(dir==ForwardCommunication)
1021 comm.forward(data_wrapper);
1022 else
1023 comm.backward(data_wrapper);
1024}
1025#endif
1026
1027template<class DataHandle>
1028void CpGridData::communicate(DataHandle& data, InterfaceType iftype,
1029 CommunicationDirection dir)
1030{
1031#if HAVE_MPI
1032 if(data.contains(3,0))
1033 {
1034 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*this, data);
1035 communicateCodim<0>(data_wrapper, dir, getInterface(iftype, cell_interfaces_));
1036 }
1037 if(data.contains(3,3))
1038 {
1039 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*this, data);
1040 communicateCodim<3>(data_wrapper, dir, getInterface(iftype, point_interfaces_));
1041 }
1042#else
1043 // Suppress warnings for unused arguments.
1044 (void) data;
1045 (void) iftype;
1046 (void) dir;
1047#endif
1048}
1049}}
1050
1051#if HAVE_MPI
1054
1055namespace Dune {
1056namespace cpgrid {
1057
1058namespace mover
1059{
1060template<class T>
1062{
1064public:
1065 void read(T& data)
1066 {
1067 data=buffer_[index_++];
1068 }
1069 void write(const T& data)
1070 {
1071 buffer_[index_++]=data;
1072 }
1073 void reset()
1074 {
1075 index_=0;
1076 }
1077 void resize(std::size_t size)
1078 {
1079 buffer_.resize(size);
1080 index_=0;
1081 }
1082private:
1083 std::vector<T> buffer_;
1084 typename std::vector<T>::size_type index_;
1085};
1086template<class DataHandle,int codim>
1087struct Mover
1088{
1089};
1090
1091template<class DataHandle>
1093{
1094 explicit BaseMover(DataHandle& data)
1095 : data_(data)
1096 {}
1097 template<class E>
1098 void moveData(const E& from, const E& to)
1099 {
1100 std::size_t size=data_.size(from);
1101 buffer.resize(size);
1102 data_.gather(buffer, from);
1103 buffer.reset();
1104 data_.scatter(buffer, to, size);
1105 }
1106 DataHandle& data_;
1108};
1109
1110
1111template<class DataHandle>
1113{
1114 Mover(DataHandle& data, CpGridData* gatherView,
1115 CpGridData* scatterView)
1116 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1117 {}
1118
1119 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1120 {
1121 Entity<0> from_entity=Entity<0>(*gatherView_, from_cell_index, true);
1122 Entity<0> to_entity=Entity<0>(*scatterView_, to_cell_index, true);
1123 this->moveData(from_entity, to_entity);
1124 }
1127};
1128
1129template<class DataHandle>
1131{
1132 Mover(DataHandle& data, CpGridData* gatherView,
1133 CpGridData* scatterView)
1134 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1135 {}
1136
1137 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1138 {
1139 typedef typename OrientedEntityTable<0,1>::row_type row_type;
1140 EntityRep<0> from_cell=EntityRep<0>(from_cell_index, true);
1141 EntityRep<0> to_cell=EntityRep<0>(to_cell_index, true);
1142 const OrientedEntityTable<0,1>& table = gatherView_->cell_to_face_;
1143 row_type from_faces=table.operator[](from_cell);
1144 row_type to_faces=scatterView_->cell_to_face_[to_cell];
1145
1146 for(int i=0; i<from_faces.size(); ++i)
1147 this->moveData(from_faces[i], to_faces[i]);
1148 }
1151};
1152
1153template<class DataHandle>
1155{
1156 Mover(DataHandle& data, CpGridData* gatherView,
1157 CpGridData* scatterView)
1158 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1159 {}
1160 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1161 {
1162 const std::array<int,8>& from_cell_points=
1163 gatherView_->cell_to_point_[from_cell_index];
1164 const std::array<int,8>& to_cell_points=
1165 scatterView_->cell_to_point_[to_cell_index];
1166 for(std::size_t i=0; i<8; ++i)
1167 {
1168 this->moveData(Entity<3>(*gatherView_, from_cell_points[i], true),
1169 Entity<3>(*scatterView_, to_cell_points[i], true));
1170 }
1171 }
1174};
1175
1176} // end mover namespace
1177
1178template<class DataHandle>
1179void CpGridData::scatterData(DataHandle& data, const CpGridData* global_data,
1180 const CpGridData* distributed_data, const InterfaceMap& cell_inf,
1181 const InterfaceMap& point_inf)
1182{
1183#if HAVE_MPI
1184 if(data.contains(3,0))
1185 {
1186 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*global_data, *distributed_data, data);
1187 communicateCodim<0>(data_wrapper, ForwardCommunication, cell_inf);
1188 }
1189 if(data.contains(3,3))
1190 {
1191 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*global_data, *distributed_data, data);
1192 communicateCodim<3>(data_wrapper, ForwardCommunication, point_inf);
1193 }
1194#endif
1195}
1196
1197template<int codim, class DataHandle>
1198void CpGridData::scatterCodimData(DataHandle& data, CpGridData* global_data,
1199 CpGridData* distributed_data)
1200{
1201 CpGridData *gather_view, *scatter_view;
1202 gather_view=global_data;
1203 scatter_view=distributed_data;
1204
1205 mover::Mover<DataHandle,codim> mover(data, gather_view, scatter_view);
1206
1207
1208 for(auto index=distributed_data->cellIndexSet().begin(),
1209 end = distributed_data->cellIndexSet().end();
1210 index!=end; ++index)
1211 {
1212 std::size_t from=index->global();
1213 std::size_t to=index->local();
1214 mover(from,to);
1215 }
1216}
1217
1218namespace
1219{
1220
1221template<int codim, class T, class F>
1222void visitInterior(CpGridData& distributed_data, T begin, T endit, F& func)
1223{
1224 for(T it=begin; it!=endit; ++it)
1225 {
1226 Entity<codim> entity(distributed_data, it-begin, true);
1227 PartitionType pt = entity.partitionType();
1228 if(pt==Dune::InteriorEntity)
1229 {
1230 func(*it, entity);
1231 }
1232 }
1233}
1234
1235template<class DataHandle>
1236struct GlobalIndexSizeGatherer
1237{
1238 GlobalIndexSizeGatherer(DataHandle& data_,
1239 std::vector<int>& ownedGlobalIndices_,
1240 std::vector<int>& ownedSizes_)
1241 : data(data_), ownedGlobalIndices(ownedGlobalIndices_), ownedSizes(ownedSizes_)
1242 {}
1243
1244 template<class T, class E>
1245 void operator()(T& i, E& entity)
1246 {
1247 ownedGlobalIndices.push_back(i);
1248 ownedSizes.push_back(data.size(entity));
1249 }
1250 DataHandle& data;
1251 std::vector<int>& ownedGlobalIndices;
1252 std::vector<int>& ownedSizes;
1253};
1254
1255template<class DataHandle>
1256struct DataGatherer
1257{
1258 DataGatherer(mover::MoveBuffer<typename DataHandle::DataType>& buffer_,
1259 DataHandle& data_)
1260 : buffer(buffer_), data(data_)
1261 {}
1262
1263 template<class T, class E>
1264 void operator()(T& /* it */, E& entity)
1265 {
1266 data.gather(buffer, entity);
1267 }
1268 mover::MoveBuffer<typename DataHandle::DataType>& buffer;
1269 DataHandle& data;
1270};
1271
1272}
1273
1274template<class DataHandle>
1275void CpGridData::gatherData(DataHandle& data, CpGridData* global_data,
1276 CpGridData* distributed_data)
1277{
1278#if HAVE_MPI
1279 if(data.contains(3,0))
1280 gatherCodimData<0>(data, global_data, distributed_data);
1281 if(data.contains(3,3))
1282 gatherCodimData<3>(data, global_data, distributed_data);
1283#endif
1284}
1285
1286template<int codim, class DataHandle>
1287void CpGridData::gatherCodimData(DataHandle& data, CpGridData* global_data,
1288 CpGridData* distributed_data)
1289{
1290#if HAVE_MPI
1291 // Get the mapping to global index from the global id set
1292 const std::vector<int>& mapping =
1293 distributed_data->global_id_set_->getMapping<codim>();
1294
1295 // Get the global indices and data size for the entities whose data is
1296 // to be sent, i.e. the ones that we own.
1297 std::vector<int> owned_global_indices;
1298 std::vector<int> owned_sizes;
1299 owned_global_indices.reserve(mapping.size());
1300 owned_sizes.reserve(mapping.size());
1301
1302 GlobalIndexSizeGatherer<DataHandle> gisg(data, owned_global_indices, owned_sizes);
1303 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gisg);
1304
1305 // communicate the number of indices that each processor sends
1306 int no_indices=owned_sizes.size();
1307 // We will take the address of the first elemet for MPI_Allgather below.
1308 // Make sure the containers have such an element.
1309 if ( owned_global_indices.empty() )
1310 owned_global_indices.resize(1);
1311 if ( owned_sizes.empty() )
1312 owned_sizes.resize(1);
1313 std::vector<int> no_indices_to_recv(distributed_data->ccobj_.size());
1314 distributed_data->ccobj_.allgather(&no_indices, 1, &(no_indices_to_recv[0]));
1315 // compute size of the vector capable for receiving all indices
1316 // and allgather the global indices and the sizes.
1317 // calculate displacements
1318 std::vector<int> displ(distributed_data->ccobj_.size()+1, 0);
1319 std::transform(displ.begin(), displ.end()-1, no_indices_to_recv.begin(), displ.begin()+1,
1320 std::plus<int>());
1321 int global_size=displ[displ.size()-1];//+no_indices_to_recv[displ.size()-1];
1322 std::vector<int> global_indices(global_size);
1323 std::vector<int> global_sizes(global_size);
1324 MPI_Allgatherv(&(owned_global_indices[0]), no_indices, MPITraits<int>::getType(),
1325 &(global_indices[0]), &(no_indices_to_recv[0]), &(displ[0]),
1326 MPITraits<int>::getType(),
1327 distributed_data->ccobj_);
1328 MPI_Allgatherv(&(owned_sizes[0]), no_indices, MPITraits<int>::getType(),
1329 &(global_sizes[0]), &(no_indices_to_recv[0]), &(displ[0]),
1330 MPITraits<int>::getType(),
1331 distributed_data->ccobj_);
1332 std::vector<int>().swap(owned_global_indices); // free data for reuse.
1333 // Compute the number of data items to send
1334 std::vector<int> no_data_send(distributed_data->ccobj_.size());
1335 for(typename std::vector<int>::iterator begin=no_data_send.begin(),
1336 i=begin, end=no_data_send.end(); i!=end; ++i)
1337 *i = std::accumulate(global_sizes.begin()+displ[i-begin],
1338 global_sizes.begin()+displ[i-begin+1], std::size_t());
1339 // free at least some memory that can be reused.
1340 std::vector<int>().swap(owned_sizes);
1341 // compute the displacements for receiving with allgatherv
1342 displ[0]=0;
1343 std::transform(displ.begin(), displ.end()-1, no_data_send.begin(), displ.begin()+1,
1344 std::plus<std::size_t>());
1345 // Compute the number of data items we will receive
1346 int no_data_recv = displ[displ.size()-1];//+global_sizes[displ.size()-1];
1347
1348 // Collect the data to send, gather it
1349 mover::MoveBuffer<typename DataHandle::DataType> local_data_buffer, global_data_buffer;
1350 if ( no_data_send[distributed_data->ccobj_.rank()] )
1351 {
1352 local_data_buffer.resize(no_data_send[distributed_data->ccobj_.rank()]);
1353 }
1354 else
1355 {
1356 local_data_buffer.resize(1);
1357 }
1358 global_data_buffer.resize(no_data_recv);
1359
1360 DataGatherer<DataHandle> gatherer(local_data_buffer, data);
1361 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gatherer);
1362 MPI_Allgatherv(&(local_data_buffer.buffer_[0]), no_data_send[distributed_data->ccobj_.rank()],
1363 MPITraits<typename DataHandle::DataType>::getType(),
1364 &(global_data_buffer.buffer_[0]), &(no_data_send[0]), &(displ[0]),
1365 MPITraits<typename DataHandle::DataType>::getType(),
1366 distributed_data->ccobj_);
1367 Entity2IndexDataHandle<DataHandle, codim> edata(*global_data, data);
1368 int offset=0;
1369 for(int i=0; i< codim; ++i)
1370 offset+=global_data->size(i);
1371
1372 typename std::vector<int>::const_iterator s=global_sizes.begin();
1373 for(typename std::vector<int>::const_iterator i=global_indices.begin(),
1374 end=global_indices.end();
1375 i!=end; ++s, ++i)
1376 {
1377 edata.scatter(global_data_buffer, *i-offset, *s);
1378 }
1379#endif
1380}
1381
1382} // end namespace cpgrid
1383} // end namespace Dune
1384
1385#endif
1386
1387#endif
DataHandle & data
Definition: CpGridData.hpp:1250
void refinePatch_and_check(const std::array< int, 3 > &, const std::array< int, 3 > &, const std::array< int, 3 > &)
void testInactiveCellsLgrs(const std::string &, const std::vector< std::array< int, 3 > > &, const std::vector< std::array< int, 3 > > &, const std::vector< std::array< int, 3 > > &, const std::vector< std::string > &)
mover::MoveBuffer< typename DataHandle::DataType > & buffer
Definition: CpGridData.hpp:1268
void fieldProp_check(const Dune::CpGrid &grid, Opm::EclipseGrid eclGrid, const std::string &deck_string)
void check_global_refine(const Dune::CpGrid &, const Dune::CpGrid &)
std::vector< int > & ownedGlobalIndices
Definition: CpGridData.hpp:1251
void refine_and_check(const Dune::cpgrid::Geometry< 3, 3 > &, const std::array< int, 3 > &, bool)
std::vector< int > & ownedSizes
Definition: CpGridData.hpp:1252
[ provides Dune::Grid ]
Definition: CpGrid.hpp:198
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:138
void postAdapt()
Clean up refinement/coarsening markers - set every element to the mark 0 which represents 'doing noth...
const cpgrid::LevelGlobalIdSet & globalIdSet() const
Get the global index set.
Definition: CpGridData.hpp:666
CpGridDataTraits::CommunicationType CommunicationType
type of OwnerOverlap communication for cells
Definition: CpGridData.hpp:716
@ MAX_DATA_PER_CELL
The maximum data items allowed per cell (DUNE < 2.5.2)
Definition: CpGridData.hpp:186
CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGridData.hpp:719
std::vector< int > getPatchesCells(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Compute cell indices of selected patches of cells (Cartesian grid required).
int sharedFaceTag(const std::vector< std::array< int, 3 > > &startIJK_2Patches, const std::vector< std::array< int, 3 > > &endIJK_2Patches) const
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: CpGridData.hpp:217
const std::array< int, 3 > & logicalCartesianSize() const
Definition: CpGridData.hpp:674
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:1028
bool uniqueBoundaryIds() const
Definition: CpGridData.hpp:629
const std::tuple< int, std::vector< int > > & getChildrenLevelAndIndexList(int elemIdx) const
Retrieves the level and child indices of a given parent cell.
Definition: CpGridData.hpp:550
std::array< double, 3 > computeEclCentroid(const Entity< 0 > &elem) const
std::array< int, 3 > getIJK(int idx_in_parent_cell, const std::array< int, 3 > &cells_per_dim) const
Extract Cartesian index triplet (i,j,k) given an index between 0 and NXxNYxNZ -1 where NX,...
Definition: CpGridData.hpp:322
void computeCommunicationInterfaces(int noexistingPoints)
RemoteIndices & cellRemoteIndices()
Definition: CpGridData.hpp:750
bool patchesShareFace(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Determine if a finite amount of patches (of cells) share a face.
int size(int codim) const
number of leaf entities per codim in this process
CpGridDataTraits::CollectiveCommunication CollectiveCommunication
Definition: CpGridData.hpp:704
const std::vector< int > & globalCell() const
Definition: CpGridData.hpp:311
void readEclipseFormat(const std::string &filename, bool periodic_extension, bool turn_normals=false)
CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGridData.hpp:713
CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGridData.hpp:703
const std::vector< std::tuple< int, std::vector< int > > > & getParentToChildren() const
Definition: CpGridData.hpp:554
void getIJK(int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: CpGridData.hpp:290
const IndexSet & indexSet() const
Definition: CpGridData.hpp:654
CommunicationType & cellCommunication()
Get the owner-overlap-copy communication for cells.
Definition: CpGridData.hpp:727
std::tuple< std::shared_ptr< CpGridData >, const std::vector< std::array< int, 2 > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::vector< std::array< int, 2 > >, const std::vector< std::array< int, 2 > > > refinePatch(const std::array< int, 3 > &cells_per_dim, const std::array< int, 3 > &startIJK, const std::array< int, 3 > &endIJK) const
Refine a (connected block-shaped) patch of cells. Based on the patch, a Geometry<3,...
CpGridDataTraits::RemoteIndices RemoteIndices
The type of the remote indices information.
Definition: CpGridData.hpp:722
ParallelIndexSet & cellIndexSet()
Definition: CpGridData.hpp:740
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
const ParallelIndexSet & cellIndexSet() const
Definition: CpGridData.hpp:745
CpGridDataTraits::MPICommunicator MPICommunicator
The type of the mpi communicator.
Definition: CpGridData.hpp:701
void processEclipseFormat(const grdecl &input_data, std::array< std::set< std::pair< int, int > >, 2 > &nnc, bool remove_ij_boundary, bool turn_normals, bool pinchActive, double tolerance_unique_points)
CpGridData(std::vector< std::shared_ptr< CpGridData > > &data)
Constructor.
std::tuple< const std::shared_ptr< CpGridData >, const std::vector< std::array< int, 2 > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::tuple< int, std::vector< int > >, const std::vector< std::array< int, 2 > >, const std::vector< std::array< int, 2 > > > refineSingleCell(const std::array< int, 3 > &cells_per_dim, const int &parent_idx) const
Refine a single cell and return a shared pointer of CpGridData type.
void distributeGlobalGrid(CpGrid &grid, const CpGridData &view_data, const std::vector< int > &cell_part)
Redistribute a global grid.
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & levelData() const
Add doc/or remove method and replace it with better approach.
Definition: CpGridData.hpp:535
const std::vector< double > & zcornData() const
Definition: CpGridData.hpp:647
void setUniqueBoundaryIds(bool uids)
Definition: CpGridData.hpp:636
bool preAdapt()
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called af...
int getGridIdx() const
Add doc/or remove method and replace it with better approach.
Definition: CpGridData.hpp:521
bool hasNNCs(const std::vector< int > &cellIndices) const
Check all cells selected for refinement have no NNCs (no neighbor connections). Assumption: all grid ...
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
Mark entity for refinement or coarsening.
const CommunicationType & cellCommunication() const
Get the owner-overlap-copy communication for cells.
Definition: CpGridData.hpp:735
int cellFace(int cell, int local_index) const
Definition: CpGridData.hpp:295
const cpgrid::IdSet & localIdSet() const
Get the local index set.
Definition: CpGridData.hpp:660
bool adapt()
TO DO: Documentation. Triggers the grid refinement process - Currently, returns preAdapt()
int faceToCellSize(int face) const
Definition: CpGridData.hpp:300
void checkCuboidShape(const std::vector< int > &cellIdx_vec) const
Check that every cell to be refined has cuboid shape.
CpGridDataTraits::Communicator Communicator
The type of the Communicator.
Definition: CpGridData.hpp:710
bool disjointPatches(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Determine if a finite amount of patches (of cells) are disjoint, namely, they do not share any corner...
std::array< double, 3 > computeEclCentroid(const int idx) const
const std::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
Definition: CpGridData.hpp:762
void validStartEndIJKs(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Check startIJK and endIJK of each patch of cells to be refined are valid, i.e. startIJK and endIJK ve...
const RemoteIndices & cellRemoteIndices() const
Definition: CpGridData.hpp:755
CpGridData(MPIHelper::MPICommunicator comm, std::vector< std::shared_ptr< CpGridData > > &data)
Definition: DefaultGeometryPolicy.hpp:53
const EntityVariable< cpgrid::Geometry< 3 - codim, 3 >, codim > & geomVector() const
Definition: DefaultGeometryPolicy.hpp:86
Wrapper that turns a data handle suitable for dune-grid into one based on integers instead of entitie...
Definition: Entity2IndexDataHandle.hpp:56
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:99
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:199
Definition: Indexsets.hpp:56
Definition: Intersection.hpp:66
Definition: Indexsets.hpp:350
int size() const
Returns the number of rows in the table.
Definition: SparseTable.hpp:121
Definition: PartitionTypeIndicator.hpp:50
Definition: CpGridData.hpp:1062
void write(const T &data)
Definition: CpGridData.hpp:1069
void reset()
Definition: CpGridData.hpp:1073
void read(T &data)
Definition: CpGridData.hpp:1065
void resize(std::size_t size)
Definition: CpGridData.hpp:1077
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
Dune::cpgrid::Cell2FacesContainer cell2Faces(const Dune::CpGrid &grid)
Get the cell to faces mapping of a grid.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
MPIHelper::MPICommunicator MPICommunicator
The type of the collective communication.
Definition: CpGridDataTraits.hpp:56
Dune::VariableSizeCommunicator<> Communicator
The type of the Communicator.
Definition: CpGridDataTraits.hpp:71
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
AttributeSet
The type of the set of the attributes.
Definition: CpGridDataTraits.hpp:66
Communicator::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGridDataTraits.hpp:74
Definition: CpGridData.hpp:1093
BaseMover(DataHandle &data)
Definition: CpGridData.hpp:1094
void moveData(const E &from, const E &to)
Definition: CpGridData.hpp:1098
MoveBuffer< typename DataHandle::DataType > buffer
Definition: CpGridData.hpp:1107
DataHandle & data_
Definition: CpGridData.hpp:1106
Definition: CpGridData.hpp:1113
void operator()(std::size_t from_cell_index, std::size_t to_cell_index)
Definition: CpGridData.hpp:1119
CpGridData * scatterView_
Definition: CpGridData.hpp:1126
CpGridData * gatherView_
Definition: CpGridData.hpp:1125
Mover(DataHandle &data, CpGridData *gatherView, CpGridData *scatterView)
Definition: CpGridData.hpp:1114
Definition: CpGridData.hpp:1131
Mover(DataHandle &data, CpGridData *gatherView, CpGridData *scatterView)
Definition: CpGridData.hpp:1132
CpGridData * gatherView_
Definition: CpGridData.hpp:1149
CpGridData * scatterView_
Definition: CpGridData.hpp:1150
void operator()(std::size_t from_cell_index, std::size_t to_cell_index)
Definition: CpGridData.hpp:1137
Definition: CpGridData.hpp:1155
CpGridData * scatterView_
Definition: CpGridData.hpp:1173
CpGridData * gatherView_
Definition: CpGridData.hpp:1172
Mover(DataHandle &data, CpGridData *gatherView, CpGridData *scatterView)
Definition: CpGridData.hpp:1156
void operator()(std::size_t from_cell_index, std::size_t to_cell_index)
Definition: CpGridData.hpp:1160
Definition: CpGridData.hpp:1088
Definition: preprocess.h:56