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 <set>
75#include <vector>
76
77namespace Opm
78{
79class EclipseState;
80}
81namespace Dune
82{
83class CpGrid;
84
85namespace cpgrid
86{
87
88class IndexSet;
89class IdSet;
90class LevelGlobalIdSet;
91class PartitionTypeIndicator;
92template<int,int> class Geometry;
93template<int> class Entity;
94template<int> class EntityRep;
95}
96}
97
99 const std::array<int, 3>&,
100 bool);
101
102void refinePatch_and_check(const std::array<int,3>&,
103 const std::array<int,3>&,
104 const std::array<int,3>&);
105
107 const std::vector<std::array<int,3>>&,
108 const std::vector<std::array<int,3>>&,
109 const std::vector<std::array<int,3>>&,
110 const std::vector<std::string>&);
111
113 const Dune::CpGrid&);
114
115void fieldProp_check(const Dune::CpGrid& grid, Opm::EclipseGrid eclGrid, std::string deck_string);
116
117namespace Dune
118{
119namespace cpgrid
120{
121namespace mover
122{
123template<class T, int i> struct Mover;
124}
125
131{
132 template<class T, int i> friend struct mover::Mover;
133 friend class GlobalIdSet;
134 friend class HierarchicIterator;
136
137 friend
139 const std::array<int, 3>&,
140 bool);
141 friend
142 void ::refinePatch_and_check(const std::array<int,3>&,
143 const std::array<int,3>&,
144 const std::array<int,3>&);
145
146 friend
148 const std::vector<std::array<int,3>>&,
149 const std::vector<std::array<int,3>>&,
150 const std::vector<std::array<int,3>>&,
151 const std::vector<std::string>&);
152
153 friend
155 const Dune::CpGrid&);
156
157 friend
158 void ::fieldProp_check(const Dune::CpGrid& grid, Opm::EclipseGrid eclGrid, std::string deck_string);
159
160private:
161 CpGridData(const CpGridData& g);
162
163public:
164 enum{
165#ifndef MAX_DATA_COMMUNICATED_PER_ENTITY
173#else
178 MAX_DATA_PER_CELL = MAX_DATA_COMMUNICATED_PER_ENTITY
179#endif
180 };
181
185 explicit CpGridData(MPIHelper::MPICommunicator comm, std::vector<std::shared_ptr<CpGridData>>& data);
186
187
188
190 explicit CpGridData(std::vector<std::shared_ptr<CpGridData>>& data);
193
194
195
196
198 int size(int codim) const;
199
201 int size (GeometryType type) const
202 {
203 if (type.isCube()) {
204 return size(3 - type.dim());
205 } else {
206 return 0;
207 }
208 }
212 void readSintefLegacyFormat(const std::string& grid_prefix);
213
217 void writeSintefLegacyFormat(const std::string& grid_prefix) const;
218
224 void readEclipseFormat(const std::string& filename, bool periodic_extension, bool turn_normals = false);
225
226#if HAVE_ECL_INPUT
235 void processEclipseFormat(const Opm::Deck& deck, bool periodic_extension, bool turn_normals = false, bool clip_z = false,
236 const std::vector<double>& poreVolume = std::vector<double>());
237
250 std::vector<std::size_t> processEclipseFormat(const Opm::EclipseGrid* ecl_grid, Opm::EclipseState* ecl_state,
251 bool periodic_extension, bool turn_normals = false, bool clip_z = false,
252 bool pinchActive = true);
253#endif
254
266 void processEclipseFormat(const grdecl& input_data,
267#if HAVE_ECL_INPUT
268 Opm::EclipseState* ecl_state,
269#endif
270 std::array<std::set<std::pair<int, int>>, 2>& nnc,
271 bool remove_ij_boundary, bool turn_normals, bool pinchActive,
272 double tolerance_unique_points);
273
281 void getIJK(int c, std::array<int,3>& ijk) const
282 {
283 int gc = global_cell_[c];
284 ijk[0] = gc % logical_cartesian_size_[0]; gc /= logical_cartesian_size_[0];
285 ijk[1] = gc % logical_cartesian_size_[1];
286 ijk[2] = gc / logical_cartesian_size_[1];
287 }
288
294 bool disjointPatches(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
295
303 std::vector<int>
304 getPatchesCells(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
305
308 bool hasNNCs(const std::vector<int>& cellIndices) const;
309
316 void validStartEndIJKs(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
317
319 void checkCuboidShape(const std::vector<int>& cellIdx_vec) const;
320
326 bool patchesShareFace(const std::vector<std::array<int,3>>& startIJK_vec, const std::vector<std::array<int,3>>& endIJK_vec) const;
327
328private:
336 std::array<int,3> getPatchDim(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
337
345 std::vector<int> getPatchCorners(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
346
354 std::vector<int> getPatchFaces(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
355
363 std::vector<int> getPatchCells(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
364
372 std::vector<int> getPatchBoundaryCorners(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
373
381 std::array<std::vector<int>,6> getBoundaryPatchFaces(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
382
384 std::array<std::vector<double>,3> getWidthsLengthsHeights(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
385
400 Geometry<3,3> cellifyPatch(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK,
401 const std::vector<int>& patch_cells, DefaultGeometryPolicy& cellifiedPatch_geometry,
402 std::array<int,8>& cellifiedPatch_to_point,
403 std::array<int,8>& allcorners_cellifiedPatch) const;
404
405 // @brief Compute the average of array<double,3>.
406 //
407 // @param [in] vector of array<double,3>
408 // @return array<double,3> (average of the entries of the given vector).
409 std::array<double,3> getAverageArr(const std::vector<std::array<double,3>>& vec) const;
410
411public:
433 std::tuple< const std::shared_ptr<CpGridData>,
434 const std::vector<std::array<int,2>>, // parent_to_refined_corners(~boundary_old_to_new_corners)
435 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_faces (~boundary_old_to_new_faces)
436 const std::tuple<int, std::vector<int>>, // parent_to_children_cells
437 const std::vector<std::array<int,2>>, // child_to_parent_faces
438 const std::vector<std::array<int,2>>> // child_to_parent_cells
439 refineSingleCell(const std::array<int,3>& cells_per_dim, const int& parent_idx) const;
440
453 std::tuple< std::shared_ptr<CpGridData>,
454 const std::vector<std::array<int,2>>, // boundary_old_to_new_corners
455 const std::vector<std::tuple<int,std::vector<int>>>, // boundary_old_to_new_faces
456 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_faces
457 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_cell
458 const std::vector<std::array<int,2>>, // child_to_parent_faces
459 const std::vector<std::array<int,2>>> // child_to_parent_cells
460 refinePatch(const std::array<int,3>& cells_per_dim, const std::array<int,3>& startIJK, const std::array<int,3>& endIJK) const;
461
462 // @breif Compute center of an entity/element/cell in the Eclipse way:
463 // - Average of the 4 corners of the bottom face.
464 // - Average of the 4 corners of the top face.
465 // Return average of the previous computations.
466 // @param [in] int Index of a cell.
467 // @return 'eclipse centroid'
468 std::array<double,3> computeEclCentroid(const int idx) const;
469
470 // @breif Compute center of an entity/element/cell in the Eclipse way:
471 // - Average of the 4 corners of the bottom face.
472 // - Average of the 4 corners of the top face.
473 // Return average of the previous computations.
474 // @param [in] Entity<0> Entity
475 // @return 'eclipse centroid'
476 std::array<double,3> computeEclCentroid(const Entity<0>& elem) const;
477
478 // Make unique boundary ids for all intersections.
480
484 bool uniqueBoundaryIds() const
485 {
486 return use_unique_boundary_ids_;
487 }
488
491 void setUniqueBoundaryIds(bool uids)
492 {
493 use_unique_boundary_ids_ = uids;
494 if (use_unique_boundary_ids_ && unique_boundary_ids_.empty()) {
496 }
497 }
498
502 const std::vector<double>& zcornData() const {
503 return zcorn;
504 }
505
506
509 const IndexSet& indexSet() const
510 {
511 return *index_set_;
512 }
513
516 {
517 return *local_id_set_;
518 }
519
523 const std::array<int, 3>& logicalCartesianSize() const
524 {
525 return logical_cartesian_size_;
526 }
527
532 const CpGridData& view_data,
533 const std::vector<int>& cell_part);
534
540 template<class DataHandle>
541 void communicate(DataHandle& data, InterfaceType iftype, CommunicationDirection dir);
542
548#if HAVE_MPI
551
554
557
560
563
566
571 {
572 return cell_comm_;
573 }
574
579 {
580 return cell_comm_;
581 }
582
584 {
585 return cellCommunication().indexSet();
586 }
587
589 {
590 return cellCommunication().indexSet();
591 }
592
594 {
595 return cellCommunication().remoteIndices();
596 }
597
599 {
600 return cellCommunication().remoteIndices();
601 }
602#endif
603
605 const std::vector<int>& sortedNumAquiferCells() const
606 {
607 return aquifer_cells_;
608 }
609
610private:
611
613 void populateGlobalCellIndexSet();
614
615#if HAVE_MPI
616
622 template<class DataHandle>
623 void gatherData(DataHandle& data, CpGridData* global_view,
624 CpGridData* distributed_view);
625
626
633 template<int codim, class DataHandle>
634 void gatherCodimData(DataHandle& data, CpGridData* global_data,
635 CpGridData* distributed_data);
636
643 template<class DataHandle>
644 void scatterData(DataHandle& data, CpGridData* global_data,
645 CpGridData* distributed_data, const InterfaceMap& cell_inf,
646 const InterfaceMap& point_inf);
647
655 template<int codim, class DataHandle>
656 void scatterCodimData(DataHandle& data, CpGridData* global_data,
657 CpGridData* distributed_data);
658
667 template<int codim, class DataHandle>
668 void communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
669 const Interface& interface);
670
679 template<int codim, class DataHandle>
680 void communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
681 const InterfaceMap& interface);
682
683#endif
684
685 void computeGeometry(CpGrid& grid,
686 const DefaultGeometryPolicy& globalGeometry,
687 const std::vector<int>& globalAquiferCells,
688 const OrientedEntityTable<0, 1>& globalCell2Faces,
689 DefaultGeometryPolicy& geometry,
690 std::vector<int>& aquiferCells,
692 const std::vector< std::array<int,8> >& cell2Points);
693
694 // Representing the topology
706 Opm::SparseTable<int> face_to_point_;
708 std::vector< std::array<int,8> > cell_to_point_;
715 std::array<int, 3> logical_cartesian_size_;
722 std::vector<int> global_cell_;
728 typedef FieldVector<double, 3> PointType;
732 cpgrid::EntityVariable<int, 1> unique_boundary_ids_;
734 std::unique_ptr<cpgrid::IndexSet> index_set_;
736 std::shared_ptr<const cpgrid::IdSet> local_id_set_;
738 std::shared_ptr<LevelGlobalIdSet> global_id_set_;
740 std::shared_ptr<PartitionTypeIndicator> partition_type_indicator_;
742 int level_{0};
744 std::vector<std::shared_ptr<CpGridData>>* level_data_ptr_;
745 // SUITABLE FOR ALL LEVELS EXCEPT FOR LEAFVIEW
747 std::vector<int> level_to_leaf_cells_; // In entry 'level cell index', we store 'leafview cell index' // {level LGR, {child0, child1, ...}}
749 std::vector<std::tuple<int,std::vector<int>>> parent_to_children_cells_; // {# children in x-direction, ... y-, ... z-}
751 std::array<int,3> cells_per_dim_;
752 // SUITABLE ONLY FOR LEAFVIEW // {level, cell index in that level}
754 std::vector<std::array<int,2>> leaf_to_level_cells_;
755 // SUITABLE FOR ALL LEVELS INCLUDING LEAFVIEW // {level parent cell, parent cell index}
757 std::vector<std::array<int,2>> child_to_parent_cells_;
758
759
761 Communication ccobj_;
762
763 // Boundary information (optional).
764 bool use_unique_boundary_ids_;
765
771 std::vector<double> zcorn;
772
774 std::vector<int> aquifer_cells_;
775
776#if HAVE_MPI
777
779 CommunicationType cell_comm_;
780
782 std::tuple<Interface,Interface,Interface,Interface,Interface> cell_interfaces_;
783 /*
784 // code deactivated, because users cannot access face indices and therefore
785 // communication on faces makes no sense!
787 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
788 face_interfaces_;
789 */
791 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
792 point_interfaces_;
793
794#endif
795
796 // Return the geometry vector corresponding to the given codim.
797 template <int codim>
798 const EntityVariable<Geometry<3 - codim, 3>, codim>& geomVector() const
799 {
800 return geometry_.geomVector<codim>();
801 }
802
803 friend class Dune::CpGrid;
804 template<int> friend class Entity;
805 template<int> friend class EntityRep;
806 friend class Intersection;
808};
809
810
811
812#if HAVE_MPI
813
814namespace
815{
820template<class T>
821T& getInterface(InterfaceType iftype,
822 std::tuple<T,T,T,T,T>& interfaces)
823{
824 switch(iftype)
825 {
826 case 0:
827 return std::get<0>(interfaces);
828 case 1:
829 return std::get<1>(interfaces);
830 case 2:
831 return std::get<2>(interfaces);
832 case 3:
833 return std::get<3>(interfaces);
834 case 4:
835 return std::get<4>(interfaces);
836 }
837 OPM_THROW(std::runtime_error, "Invalid Interface type was used during communication");
838}
839
840} // end unnamed namespace
841
842template<int codim, class DataHandle>
843void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
844 const Interface& interface)
845{
846 this->template communicateCodim<codim>(data, dir, interface.interfaces());
847}
848
849template<int codim, class DataHandle>
850void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data_wrapper, CommunicationDirection dir,
851 const InterfaceMap& interface)
852{
853 Communicator comm(ccobj_, interface);
854
855 if(dir==ForwardCommunication)
856 comm.forward(data_wrapper);
857 else
858 comm.backward(data_wrapper);
859}
860#endif
861
862template<class DataHandle>
863void CpGridData::communicate(DataHandle& data, InterfaceType iftype,
864 CommunicationDirection dir)
865{
866#if HAVE_MPI
867 if(data.contains(3,0))
868 {
869 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*this, data);
870 communicateCodim<0>(data_wrapper, dir, getInterface(iftype, cell_interfaces_));
871 }
872 if(data.contains(3,3))
873 {
874 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*this, data);
875 communicateCodim<3>(data_wrapper, dir, getInterface(iftype, point_interfaces_));
876 }
877#else
878 // Suppress warnings for unused arguments.
879 (void) data;
880 (void) iftype;
881 (void) dir;
882#endif
883}
884}}
885
886#if HAVE_MPI
889
890namespace Dune {
891namespace cpgrid {
892
893namespace mover
894{
895template<class T>
897{
899public:
900 void read(T& data)
901 {
902 data=buffer_[index_++];
903 }
904 void write(const T& data)
905 {
906 buffer_[index_++]=data;
907 }
908 void reset()
909 {
910 index_=0;
911 }
912 void resize(std::size_t size)
913 {
914 buffer_.resize(size);
915 index_=0;
916 }
917private:
918 std::vector<T> buffer_;
919 typename std::vector<T>::size_type index_;
920};
921template<class DataHandle,int codim>
922struct Mover
923{
924};
925
926template<class DataHandle>
928{
929 BaseMover(DataHandle& data)
930 : data_(data)
931 {}
932 template<class E>
933 void moveData(const E& from, const E& to)
934 {
935 std::size_t size=data_.size(from);
936 buffer.resize(size);
937 data_.gather(buffer, from);
938 buffer.reset();
939 data_.scatter(buffer, to, size);
940 }
941 DataHandle& data_;
943};
944
945
946template<class DataHandle>
948{
949 Mover<DataHandle,0>(DataHandle& data, CpGridData* gatherView,
950 CpGridData* scatterView)
951 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
952 {}
953
954 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
955 {
956 Entity<0> from_entity=Entity<0>(*gatherView_, from_cell_index, true);
957 Entity<0> to_entity=Entity<0>(*scatterView_, to_cell_index, true);
958 this->moveData(from_entity, to_entity);
959 }
961 CpGridData* scatterView_;
962};
963
964template<class DataHandle>
966{
967 Mover<DataHandle,1>(DataHandle& data, CpGridData* gatherView,
968 CpGridData* scatterView)
969 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
970 {}
971
972 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
973 {
974 typedef typename OrientedEntityTable<0,1>::row_type row_type;
975 EntityRep<0> from_cell=EntityRep<0>(from_cell_index, true);
976 EntityRep<0> to_cell=EntityRep<0>(to_cell_index, true);
977 OrientedEntityTable<0,1>& table = gatherView_->cell_to_face_;
978 row_type from_faces=table.operator[](from_cell);
979 row_type to_faces=scatterView_->cell_to_face_[to_cell];
980
981 for(int i=0; i<from_faces.size(); ++i)
982 this->moveData(from_faces[i], to_faces[i]);
983 }
985 CpGridData *scatterView_;
986};
987
988template<class DataHandle>
990{
991 Mover<DataHandle,3>(DataHandle& data, CpGridData* gatherView,
992 CpGridData* scatterView)
993 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
994 {}
995 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
996 {
997 const std::array<int,8>& from_cell_points=
998 gatherView_->cell_to_point_[from_cell_index];
999 const std::array<int,8>& to_cell_points=
1000 scatterView_->cell_to_point_[to_cell_index];
1001 for(std::size_t i=0; i<8; ++i)
1002 {
1003 this->moveData(Entity<3>(*gatherView_, from_cell_points[i], true),
1004 Entity<3>(*scatterView_, to_cell_points[i], true));
1005 }
1006 }
1009};
1010
1011} // end mover namespace
1012
1013template<class DataHandle>
1014void CpGridData::scatterData(DataHandle& data, CpGridData* global_data,
1015 CpGridData* distributed_data, const InterfaceMap& cell_inf,
1016 const InterfaceMap& point_inf)
1017{
1018#if HAVE_MPI
1019 if(data.contains(3,0))
1020 {
1021 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*global_data, *distributed_data, data);
1022 communicateCodim<0>(data_wrapper, ForwardCommunication, cell_inf);
1023 }
1024 if(data.contains(3,3))
1025 {
1026 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*global_data, *distributed_data, data);
1027 communicateCodim<3>(data_wrapper, ForwardCommunication, point_inf);
1028 }
1029#endif
1030}
1031
1032template<int codim, class DataHandle>
1033void CpGridData::scatterCodimData(DataHandle& data, CpGridData* global_data,
1034 CpGridData* distributed_data)
1035{
1036 CpGridData *gather_view, *scatter_view;
1037 gather_view=global_data;
1038 scatter_view=distributed_data;
1039
1040 mover::Mover<DataHandle,codim> mover(data, gather_view, scatter_view);
1041
1042
1043 for(auto index=distributed_data->cellIndexSet().begin(),
1044 end = distributed_data->cellIndexSet().end();
1045 index!=end; ++index)
1046 {
1047 std::size_t from=index->global();
1048 std::size_t to=index->local();
1049 mover(from,to);
1050 }
1051}
1052
1053namespace
1054{
1055
1056template<int codim, class T, class F>
1057void visitInterior(CpGridData& distributed_data, T begin, T endit, F& func)
1058{
1059 for(T it=begin; it!=endit; ++it)
1060 {
1061 Entity<codim> entity(distributed_data, it-begin, true);
1062 PartitionType pt = entity.partitionType();
1063 if(pt==Dune::InteriorEntity)
1064 {
1065 func(*it, entity);
1066 }
1067 }
1068}
1069
1070template<class DataHandle>
1071struct GlobalIndexSizeGatherer
1072{
1073 GlobalIndexSizeGatherer(DataHandle& data_,
1074 std::vector<int>& ownedGlobalIndices_,
1075 std::vector<int>& ownedSizes_)
1076 : data(data_), ownedGlobalIndices(ownedGlobalIndices_), ownedSizes(ownedSizes_)
1077 {}
1078
1079 template<class T, class E>
1080 void operator()(T& i, E& entity)
1081 {
1082 ownedGlobalIndices.push_back(i);
1083 ownedSizes.push_back(data.size(entity));
1084 }
1085 DataHandle& data;
1086 std::vector<int>& ownedGlobalIndices;
1087 std::vector<int>& ownedSizes;
1088};
1089
1090template<class DataHandle>
1091struct DataGatherer
1092{
1093 DataGatherer(mover::MoveBuffer<typename DataHandle::DataType>& buffer_,
1094 DataHandle& data_)
1095 : buffer(buffer_), data(data_)
1096 {}
1097
1098 template<class T, class E>
1099 void operator()(T& /* it */, E& entity)
1100 {
1101 data.gather(buffer, entity);
1102 }
1103 mover::MoveBuffer<typename DataHandle::DataType>& buffer;
1104 DataHandle& data;
1105};
1106
1107}
1108
1109template<class DataHandle>
1110void CpGridData::gatherData(DataHandle& data, CpGridData* global_data,
1111 CpGridData* distributed_data)
1112{
1113#if HAVE_MPI
1114 if(data.contains(3,0))
1115 gatherCodimData<0>(data, global_data, distributed_data);
1116 if(data.contains(3,3))
1117 gatherCodimData<3>(data, global_data, distributed_data);
1118#endif
1119}
1120
1121template<int codim, class DataHandle>
1122void CpGridData::gatherCodimData(DataHandle& data, CpGridData* global_data,
1123 CpGridData* distributed_data)
1124{
1125#if HAVE_MPI
1126 // Get the mapping to global index from the global id set
1127 const std::vector<int>& mapping =
1128 distributed_data->global_id_set_->getMapping<codim>();
1129
1130 // Get the global indices and data size for the entities whose data is
1131 // to be sent, i.e. the ones that we own.
1132 std::vector<int> owned_global_indices;
1133 std::vector<int> owned_sizes;
1134 owned_global_indices.reserve(mapping.size());
1135 owned_sizes.reserve(mapping.size());
1136
1137 GlobalIndexSizeGatherer<DataHandle> gisg(data, owned_global_indices, owned_sizes);
1138 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gisg);
1139
1140 // communicate the number of indices that each processor sends
1141 int no_indices=owned_sizes.size();
1142 // We will take the address of the first elemet for MPI_Allgather below.
1143 // Make sure the containers have such an element.
1144 if ( owned_global_indices.empty() )
1145 owned_global_indices.resize(1);
1146 if ( owned_sizes.empty() )
1147 owned_sizes.resize(1);
1148 std::vector<int> no_indices_to_recv(distributed_data->ccobj_.size());
1149 distributed_data->ccobj_.allgather(&no_indices, 1, &(no_indices_to_recv[0]));
1150 // compute size of the vector capable for receiving all indices
1151 // and allgather the global indices and the sizes.
1152 // calculate displacements
1153 std::vector<int> displ(distributed_data->ccobj_.size()+1, 0);
1154 std::transform(displ.begin(), displ.end()-1, no_indices_to_recv.begin(), displ.begin()+1,
1155 std::plus<int>());
1156 int global_size=displ[displ.size()-1];//+no_indices_to_recv[displ.size()-1];
1157 std::vector<int> global_indices(global_size);
1158 std::vector<int> global_sizes(global_size);
1159 MPI_Allgatherv(&(owned_global_indices[0]), no_indices, MPITraits<int>::getType(),
1160 &(global_indices[0]), &(no_indices_to_recv[0]), &(displ[0]),
1161 MPITraits<int>::getType(),
1162 distributed_data->ccobj_);
1163 MPI_Allgatherv(&(owned_sizes[0]), no_indices, MPITraits<int>::getType(),
1164 &(global_sizes[0]), &(no_indices_to_recv[0]), &(displ[0]),
1165 MPITraits<int>::getType(),
1166 distributed_data->ccobj_);
1167 std::vector<int>().swap(owned_global_indices); // free data for reuse.
1168 // Compute the number of data items to send
1169 std::vector<int> no_data_send(distributed_data->ccobj_.size());
1170 for(typename std::vector<int>::iterator begin=no_data_send.begin(),
1171 i=begin, end=no_data_send.end(); i!=end; ++i)
1172 *i = std::accumulate(global_sizes.begin()+displ[i-begin],
1173 global_sizes.begin()+displ[i-begin+1], std::size_t());
1174 // free at least some memory that can be reused.
1175 std::vector<int>().swap(owned_sizes);
1176 // compute the displacements for receiving with allgatherv
1177 displ[0]=0;
1178 std::transform(displ.begin(), displ.end()-1, no_data_send.begin(), displ.begin()+1,
1179 std::plus<std::size_t>());
1180 // Compute the number of data items we will receive
1181 int no_data_recv = displ[displ.size()-1];//+global_sizes[displ.size()-1];
1182
1183 // Collect the data to send, gather it
1184 mover::MoveBuffer<typename DataHandle::DataType> local_data_buffer, global_data_buffer;
1185 if ( no_data_send[distributed_data->ccobj_.rank()] )
1186 {
1187 local_data_buffer.resize(no_data_send[distributed_data->ccobj_.rank()]);
1188 }
1189 else
1190 {
1191 local_data_buffer.resize(1);
1192 }
1193 global_data_buffer.resize(no_data_recv);
1194
1195 DataGatherer<DataHandle> gatherer(local_data_buffer, data);
1196 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gatherer);
1197 MPI_Allgatherv(&(local_data_buffer.buffer_[0]), no_data_send[distributed_data->ccobj_.rank()],
1198 MPITraits<typename DataHandle::DataType>::getType(),
1199 &(global_data_buffer.buffer_[0]), &(no_data_send[0]), &(displ[0]),
1200 MPITraits<typename DataHandle::DataType>::getType(),
1201 distributed_data->ccobj_);
1202 Entity2IndexDataHandle<DataHandle, codim> edata(*global_data, data);
1203 int offset=0;
1204 for(int i=0; i< codim; ++i)
1205 offset+=global_data->size(i);
1206
1207 typename std::vector<int>::const_iterator s=global_sizes.begin();
1208 for(typename std::vector<int>::const_iterator i=global_indices.begin(),
1209 end=global_indices.end();
1210 i!=end; ++s, ++i)
1211 {
1212 edata.scatter(global_data_buffer, *i-offset, *s);
1213 }
1214#endif
1215}
1216
1217} // end namespace cpgrid
1218} // end namespace Dune
1219
1220#endif
1221
1222#endif
DataHandle & data
Definition: CpGridData.hpp:1085
void refinePatch_and_check(const std::array< int, 3 > &, const std::array< int, 3 > &, const std::array< int, 3 > &)
mover::MoveBuffer< typename DataHandle::DataType > & buffer
Definition: CpGridData.hpp:1103
void check_global_refine(const Dune::CpGrid &, const Dune::CpGrid &)
std::vector< int > & ownedGlobalIndices
Definition: CpGridData.hpp:1086
void refine_and_check(const Dune::cpgrid::Geometry< 3, 3 > &, const std::array< int, 3 > &, bool)
std::vector< int > & ownedSizes
Definition: CpGridData.hpp:1087
void fieldProp_check(const Dune::CpGrid &grid, Opm::EclipseGrid eclGrid, std::string deck_string)
[ provides Dune::Grid ]
Definition: CpGrid.hpp:238
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:131
CpGridDataTraits::CommunicationType CommunicationType
type of OwnerOverlap communication for cells
Definition: CpGridData.hpp:559
@ MAX_DATA_PER_CELL
The maximum data items allowed per cell (DUNE < 2.5.2)
Definition: CpGridData.hpp:172
CpGridDataTraits::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGridData.hpp:562
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 size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: CpGridData.hpp:201
const std::array< int, 3 > & logicalCartesianSize() const
Definition: CpGridData.hpp:523
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:863
bool uniqueBoundaryIds() const
Definition: CpGridData.hpp:484
std::array< double, 3 > computeEclCentroid(const Entity< 0 > &elem) const
RemoteIndices & cellRemoteIndices()
Definition: CpGridData.hpp:593
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.
void readSintefLegacyFormat(const std::string &grid_prefix)
int size(int codim) const
number of leaf entities per codim in this process
CpGridDataTraits::CollectiveCommunication CollectiveCommunication
Definition: CpGridData.hpp:547
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:556
CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGridData.hpp:546
void getIJK(int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: CpGridData.hpp:281
const IndexSet & indexSet() const
Definition: CpGridData.hpp:509
CommunicationType & cellCommunication()
Get the owner-overlap-copy communication for cells.
Definition: CpGridData.hpp:570
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:565
ParallelIndexSet & cellIndexSet()
Definition: CpGridData.hpp:583
const ParallelIndexSet & cellIndexSet() const
Definition: CpGridData.hpp:588
CpGridDataTraits::MPICommunicator MPICommunicator
The type of the mpi communicator.
Definition: CpGridData.hpp:544
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< double > & zcornData() const
Definition: CpGridData.hpp:502
void setUniqueBoundaryIds(bool uids)
Definition: CpGridData.hpp:491
void writeSintefLegacyFormat(const std::string &grid_prefix) const
bool hasNNCs(const std::vector< int > &cellIndices) const
Check all cells selected for refinement have no NNCs (no neighbor connections). Assumption: all grid ...
const CommunicationType & cellCommunication() const
Get the owner-overlap-copy communication for cells.
Definition: CpGridData.hpp:578
const cpgrid::IdSet & localIdSet() const
Get the local index set.
Definition: CpGridData.hpp:515
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:553
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:605
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:598
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:304
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:108
Definition: Indexsets.hpp:199
Definition: Indexsets.hpp:55
Definition: Intersection.hpp:66
Definition: PartitionTypeIndicator.hpp:50
Definition: CpGridData.hpp:897
void write(const T &data)
Definition: CpGridData.hpp:904
void reset()
Definition: CpGridData.hpp:908
void read(T &data)
Definition: CpGridData.hpp:900
void resize(std::size_t size)
Definition: CpGridData.hpp:912
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:55
Dune::VariableSizeCommunicator<> Communicator
The type of the Communicator.
Definition: CpGridDataTraits.hpp:70
Dune::RemoteIndices< ParallelIndexSet > RemoteIndices
The type of the remote indices information.
Definition: CpGridDataTraits.hpp:82
typename CommunicationType::ParallelIndexSet ParallelIndexSet
The type of the parallel index set.
Definition: CpGridDataTraits.hpp:79
Dune::Communication< MPICommunicator > CollectiveCommunication
Definition: CpGridDataTraits.hpp:58
Dune::OwnerOverlapCopyCommunication< int, int > CommunicationType
type of OwnerOverlap communication for cells
Definition: CpGridDataTraits.hpp:76
Dune::Communication< MPICommunicator > Communication
Definition: CpGridDataTraits.hpp:57
AttributeSet
The type of the set of the attributes.
Definition: CpGridDataTraits.hpp:65
Communicator::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGridDataTraits.hpp:73
Definition: CpGridData.hpp:928
BaseMover(DataHandle &data)
Definition: CpGridData.hpp:929
void moveData(const E &from, const E &to)
Definition: CpGridData.hpp:933
MoveBuffer< typename DataHandle::DataType > buffer
Definition: CpGridData.hpp:942
DataHandle & data_
Definition: CpGridData.hpp:941
Definition: CpGridData.hpp:948
void operator()(std::size_t from_cell_index, std::size_t to_cell_index)
Definition: CpGridData.hpp:954
CpGridData * gatherView_
Definition: CpGridData.hpp:960
Definition: CpGridData.hpp:966
CpGridData * gatherView_
Definition: CpGridData.hpp:984
void operator()(std::size_t from_cell_index, std::size_t to_cell_index)
Definition: CpGridData.hpp:972
Definition: CpGridData.hpp:990
CpGridData * scatterView_
Definition: CpGridData.hpp:1008
CpGridData * gatherView_
Definition: CpGridData.hpp:1007
void operator()(std::size_t from_cell_index, std::size_t to_cell_index)
Definition: CpGridData.hpp:995
Definition: CpGridData.hpp:923
Definition: preprocess.h:56