opm-grid
Loading...
Searching...
No Matches
CpGridData.hpp
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/EclipseGrid.hpp>
63#include <opm/input/eclipse/EclipseState/Grid/NNC.hpp>
64#endif
65
67
68#include "Entity2IndexDataHandle.hpp"
69#include "CpGridDataTraits.hpp"
70//#include "DataHandleWrappers.hpp"
71//#include "GlobalIdMapping.hpp"
72#include "Geometry.hpp"
73
74#include <array>
75#include <initializer_list>
76#include <set>
77#include <vector>
78
79namespace Opm
80{
81class EclipseState;
82}
83namespace Dune
84{
85class CpGrid;
86
87namespace cpgrid
88{
89
90class IndexSet;
91class IdSet;
92class LevelGlobalIdSet;
93class PartitionTypeIndicator;
94template<int,int> class Geometry;
95template<int> class Entity;
96template<int> class EntityRep;
97}
98}
99
100void refine_and_check(const Dune::cpgrid::Geometry<3, 3>&,
101 const std::array<int, 3>&,
102 bool);
103
104namespace Dune
105{
106namespace cpgrid
107{
108namespace mover
109{
110template<class T, int i> struct Mover;
111}
112
118{
119 template<class T, int i> friend struct mover::Mover;
120 friend class GlobalIdSet;
121 friend class HierarchicIterator;
122 friend class Dune::cpgrid::IndexSet;
123 friend class Dune::cpgrid::IdSet;
125
126 friend
127 void ::refine_and_check(const Dune::cpgrid::Geometry<3, 3>&,
128 const std::array<int, 3>&,
129 bool);
130
131private:
132 CpGridData(const CpGridData& g);
133
134public:
135 enum{
136#ifndef MAX_DATA_COMMUNICATED_PER_ENTITY
144#else
149 MAX_DATA_PER_CELL = MAX_DATA_COMMUNICATED_PER_ENTITY
150#endif
151 };
152
153 CpGridData() = delete;
154
159 explicit CpGridData(MPIHelper::MPICommunicator comm, std::vector<std::shared_ptr<CpGridData>>& data);
160
161
162
164 explicit CpGridData(std::vector<std::shared_ptr<CpGridData>>& data);
166 ~CpGridData();
167
168
169
170
172 int size(int codim) const;
173
175 int size (GeometryType type) const
176 {
177 if (type.isCube()) {
178 return size(3 - type.dim());
179 } else {
180 return 0;
181 }
182 }
183
199 void readEclipseFormat(const std::string& filename,
200 bool periodic_extension,
201 bool turn_normals = false,
202 bool edge_conformal = false);
203
204#if HAVE_ECL_INPUT
226 void processEclipseFormat(const Opm::Deck& deck,
227 bool periodic_extension,
228 bool turn_normals = false,
229 bool clip_z = false,
230 const std::vector<double>& poreVolume = std::vector<double>{},
231 bool edge_conformal = false);
232
267 std::vector<std::size_t>
268 processEclipseFormat(const Opm::EclipseGrid* ecl_grid,
269 Opm::EclipseState* ecl_state,
270 bool periodic_extension,
271 bool turn_normals = false,
272 bool clip_z = false,
273 bool pinchActive = true,
274 bool edge_conformal = false);
275#endif
276
303 void processEclipseFormat(const grdecl& input_data,
304#if HAVE_ECL_INPUT
305 Opm::EclipseState* ecl_state,
306#endif
307 std::array<std::set<std::pair<int, int>>, 2>& nnc,
308 bool remove_ij_boundary,
309 bool turn_normals,
310 bool pinchActive,
311 double tolerance_unique_points,
312 bool edge_conformal);
313
321 void getIJK(int c, std::array<int,3>& ijk) const;
322
323 int cellFace(int cell, int local_index) const
324 {
325 return cell_to_face_[cpgrid::EntityRep<0>(cell, true)][local_index].index();
326 }
327
328 auto cellToFace(int cellIdx) const
329 {
330 return cell_to_face_[cpgrid::EntityRep<0>(cellIdx, true)];
331 }
332
333 const auto& cellToPoint() const
334 {
335 return cell_to_point_;
336 }
337
338 const auto& cellToPoint(int cellIdx) const
339 {
340 return cell_to_point_[cellIdx];
341 }
342
343 int faceToCellSize(int face) const {
344 Dune::cpgrid::EntityRep<1> faceRep(face, true);
345 return face_to_cell_[faceRep].size();
346 }
347
348 auto faceTag(int faceIdx) const
349 {
350 Dune::cpgrid::EntityRep<1> faceRep(faceIdx, true);
351 return face_tag_[faceRep];
352 }
353
354 auto faceNormals(int faceIdx) const
355 {
356 Dune::cpgrid::EntityRep<1> faceRep(faceIdx, true);
357 return face_normals_[faceRep];
358 }
359
360 auto faceToPoint(int faceIdx) const
361 {
362 return face_to_point_[faceIdx];
363 }
364
365 int numFaces() const
366 {
367 return face_to_cell_.size();
368 }
369
370 auto cornerHistorySize() const
371 {
372 return corner_history_.size();
373 }
374
375 const auto& getCornerHistory(int cornerIdx) const
376 {
377 if(cornerHistorySize()) {
378 return corner_history_[cornerIdx];
379 }
380 else {
381 OPM_THROW(std::logic_error, "Vertex has no history record.\n");
382 }
383 }
384
391 const std::vector<int>& globalCell() const
392 {
393 return global_cell_;
394 }
395
398 bool hasNNCs(const std::vector<int>& cellIndices) const;
399
412 bool mark(int refCount, const cpgrid::Entity<0>& element);
413
417 int getMark(const cpgrid::Entity<0>& element) const;
418
428 bool preAdapt();
429
431 bool adapt();
432
434 void postAdapt();
435
436private:
437 std::array<Dune::FieldVector<double,3>,8> getReferenceRefinedCorners(int idx_in_parent_cell, const std::array<int,3>& cells_per_dim) const;
438
439public:
441 int getGridIdx() const {
442 // Not the nicest way of checking if "this" points at the leaf grid view of a mixed grid (with coarse and refined cells).
443 // 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.
444 // 2. Unfortunately, level_ is default initialized by 0. This implies, in particular, that if someone wants to check the value of
445 // "this->level_" when "this" points at the leaf grid view of a grid that has been refined, this value is - unfortunately - equal to 0.
446 // 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
447 // reason we check if child_to_parent_cells_.empty() [true for actual level 0 grid, false for the leaf grid view].
448 // --- TO BE IMPROVED ---
449 if ((level_data_ptr_ ->size() >1) && (level_ == 0) && (!child_to_parent_cells_.empty())) {
450 return level_data_ptr_->size() -1;
451 }
452 return level_;
453 }
455 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& levelData() const
456 {
457 if (level_data_ptr_->empty()) {
458 OPM_THROW(std::logic_error, "Level data has not been initialized\n");
459 }
460 return *level_data_ptr_;
461 }
462
470 const std::tuple<int,std::vector<int>>& getChildrenLevelAndIndexList(int elemIdx) const {
471 return parent_to_children_cells_[elemIdx];
472 }
473
474 const std::vector<std::tuple<int,std::vector<int>>>& getParentToChildren() const {
475 return parent_to_children_cells_;
476 }
477
478 const cpgrid::DefaultGeometryPolicy getGeometry() const
479 {
480 return geometry_;
481 }
482
483 int getLeafIdxFromLevelIdx(int level_cell_idx) const
484 {
485 if (level_to_leaf_cells_.empty()) {
486 OPM_THROW(std::logic_error, "Grid has no LGRs. No mapping to the leaf.\n");
487 }
488 return level_to_leaf_cells_[level_cell_idx];
489 }
490
512 std::tuple< const std::shared_ptr<CpGridData>,
513 const std::vector<std::array<int,2>>, // parent_to_refined_corners(~boundary_old_to_new_corners)
514 const std::vector<std::tuple<int,std::vector<int>>>, // parent_to_children_faces (~boundary_old_to_new_faces)
515 const std::tuple<int, std::vector<int>>, // parent_to_children_cells
516 const std::vector<std::array<int,2>>, // child_to_parent_faces
517 const std::vector<std::array<int,2>>> // child_to_parent_cells
518 refineSingleCell(const std::array<int,3>& cells_per_dim, const int& parent_idx) const;
519
520 // @breif Compute center of an entity/element/cell in the Eclipse way:
521 // - Average of the 4 corners of the bottom face.
522 // - Average of the 4 corners of the top face.
523 // Return average of the previous computations.
524 // @param [in] int Index of a cell.
525 // @return 'eclipse centroid'
526 std::array<double,3> computeEclCentroid(const int idx) const;
527
528 // @breif Compute center of an entity/element/cell in the Eclipse way:
529 // - Average of the 4 corners of the bottom face.
530 // - Average of the 4 corners of the top face.
531 // Return average of the previous computations.
532 // @param [in] Entity<0> Entity
533 // @return 'eclipse centroid'
534 std::array<double,3> computeEclCentroid(const Entity<0>& elem) const;
535
536 // Make unique boundary ids for all intersections.
537 void computeUniqueBoundaryIds();
538
542 bool uniqueBoundaryIds() const
543 {
544 return use_unique_boundary_ids_;
545 }
546
549 void setUniqueBoundaryIds(bool uids)
550 {
551 use_unique_boundary_ids_ = uids;
552 if (use_unique_boundary_ids_ && unique_boundary_ids_.empty()) {
553 computeUniqueBoundaryIds();
554 }
555 }
556
560 const std::vector<double>& zcornData() const {
561 return zcorn;
562 }
563
564
567 const IndexSet& indexSet() const
568 {
569 return *index_set_;
570 }
571
574 {
575 return *local_id_set_;
576 }
577
580 {
581 return *global_id_set_;
582 }
583
587 const std::array<int, 3>& logicalCartesianSize() const
588 {
589 return logical_cartesian_size_;
590 }
591
595 void distributeGlobalGrid(CpGrid& grid,
596 const CpGridData& view_data,
597 const std::vector<int>& cell_part);
598
604 template<class DataHandle>
605 void communicate(DataHandle& data, InterfaceType iftype, CommunicationDirection dir);
606
607 void computeCellPartitionType();
608
609 void computePointPartitionType();
610
611 void computeCommunicationInterfaces(int noexistingPoints);
612
616 using Communication = CpGridDataTraits::Communication;
617 using CollectiveCommunication = CpGridDataTraits::CollectiveCommunication;
618
621#if HAVE_MPI
623 using Communicator = CpGridDataTraits::Communicator;
624
626 using InterfaceMap = CpGridDataTraits::InterfaceMap;
627
629 using CommunicationType = CpGridDataTraits::CommunicationType;
630
632 using ParallelIndexSet = CpGridDataTraits::ParallelIndexSet;
633
635 using RemoteIndices = CpGridDataTraits::RemoteIndices;
636
640 CommunicationType& cellCommunication()
641 {
642 return cell_comm_;
643 }
644
648 const CommunicationType& cellCommunication() const
649 {
650 return cell_comm_;
651 }
652
653 ParallelIndexSet& cellIndexSet()
654 {
655 return cellCommunication().indexSet();
656 }
657
658 const ParallelIndexSet& cellIndexSet() const
659 {
660 return cellCommunication().indexSet();
661 }
662
663 RemoteIndices& cellRemoteIndices()
664 {
665 return cellCommunication().remoteIndices();
666 }
667
668 const RemoteIndices& cellRemoteIndices() const
669 {
670 return cellCommunication().remoteIndices();
671 }
672#endif
673
675 const std::vector<int>& sortedNumAquiferCells() const
676 {
677 return aquifer_cells_;
678 }
679
680private:
681
683 void populateGlobalCellIndexSet();
684
685#if HAVE_MPI
686
692 template<class DataHandle>
693 void gatherData(DataHandle& data, CpGridData* global_view,
694 CpGridData* distributed_view);
695
696
703 template<int codim, class DataHandle>
704 void gatherCodimData(DataHandle& data, CpGridData* global_data,
705 CpGridData* distributed_data);
706
713 template<class DataHandle>
714 void scatterData(DataHandle& data, const CpGridData* global_data,
715 const CpGridData* distributed_data, const InterfaceMap& cell_inf,
716 const InterfaceMap& point_inf);
717
725 template<int codim, class DataHandle>
726 void scatterCodimData(DataHandle& data, CpGridData* global_data,
727 CpGridData* distributed_data);
728
737 template<int codim, class DataHandle>
738 void communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
739 const Interface& interface);
740
749 template<int codim, class DataHandle>
750 void communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
751 const InterfaceMap& interface);
752
753#endif
754
755 void computeGeometry(const CpGrid& grid,
756 const DefaultGeometryPolicy& globalGeometry,
757 const std::vector<int>& globalAquiferCells,
758 const OrientedEntityTable<0, 1>& globalCell2Faces,
759 DefaultGeometryPolicy& geometry,
760 std::vector<int>& aquiferCells,
761 const OrientedEntityTable<0, 1>& cell2Faces,
762 const std::vector< std::array<int,8> >& cell2Points);
763
764 // Representing the topology
776 Opm::SparseTable<int> face_to_point_;
778 std::vector< std::array<int,8> > cell_to_point_;
785 std::array<int, 3> logical_cartesian_size_{};
792 std::vector<int> global_cell_;
798 typedef FieldVector<double, 3> PointType;
802 cpgrid::EntityVariable<int, 1> unique_boundary_ids_;
804 std::unique_ptr<cpgrid::IndexSet> index_set_;
806 std::shared_ptr<const cpgrid::IdSet> local_id_set_;
808 std::shared_ptr<LevelGlobalIdSet> global_id_set_;
810 std::shared_ptr<PartitionTypeIndicator> partition_type_indicator_;
812 std::vector<int> mark_;
814 int level_{0};
816 std::vector<std::shared_ptr<CpGridData>>* level_data_ptr_;
817 // SUITABLE FOR ALL LEVELS EXCEPT FOR LEAFVIEW
819 std::vector<int> level_to_leaf_cells_; // In entry 'level cell index', we store 'leafview cell index'
// {level LGR, {child0, child1, ...}}
821 std::vector<std::tuple<int,std::vector<int>>> parent_to_children_cells_;
// {# children in x-direction, ... y-, ... z-}
823 std::array<int,3> cells_per_dim_;
824 // SUITABLE ONLY FOR LEAFVIEW
// {level, cell index in that level}
826 std::vector<std::array<int,2>> leaf_to_level_cells_;
828 std::vector<std::array<int,2>> corner_history_;
829 // SUITABLE FOR ALL LEVELS INCLUDING LEAFVIEW
// {level parent cell, parent cell index}
831 std::vector<std::array<int,2>> child_to_parent_cells_;
834 std::vector<int> cell_to_idxInParentCell_;
836 int refinement_max_level_{0};
837
839 Communication ccobj_;
840
841 // Boundary information (optional).
842 bool use_unique_boundary_ids_;
843
849 std::vector<double> zcorn;
850
852 std::vector<int> aquifer_cells_;
853
854#if HAVE_MPI
855
857 CommunicationType cell_comm_;
858
860 std::tuple<Interface,Interface,Interface,Interface,Interface> cell_interfaces_;
861 /*
862 // code deactivated, because users cannot access face indices and therefore
863 // communication on faces makes no sense!
865 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
866 face_interfaces_;
867 */
869 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
870 point_interfaces_;
871
872#endif
873
874 // Return the geometry vector corresponding to the given codim.
875 template <int codim>
876 const EntityVariable<Geometry<3 - codim, 3>, codim>& geomVector() const
877 {
878 return geometry_.geomVector<codim>();
879 }
880
881 friend class Dune::CpGrid;
882 template<int> friend class Entity;
883 template<int> friend class EntityRep;
884 friend class Intersection;
885 friend class PartitionTypeIndicator;
886};
887
888
889
890#if HAVE_MPI
891
892namespace
893{
898template<class T>
899T& getInterface(InterfaceType iftype,
900 std::tuple<T,T,T,T,T>& interfaces)
901{
902 switch(iftype)
903 {
904 case 0:
905 return std::get<0>(interfaces);
906 case 1:
907 return std::get<1>(interfaces);
908 case 2:
909 return std::get<2>(interfaces);
910 case 3:
911 return std::get<3>(interfaces);
912 case 4:
913 return std::get<4>(interfaces);
914 }
915 OPM_THROW(std::runtime_error, "Invalid Interface type was used during communication");
916}
917
918} // end unnamed namespace
919
920template<int codim, class DataHandle>
921void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
922 const Interface& interface)
923{
924 this->template communicateCodim<codim>(data, dir, interface.interfaces());
925}
926
927template<int codim, class DataHandle>
928void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data_wrapper, CommunicationDirection dir,
929 const InterfaceMap& interface)
930{
931 Communicator comm(ccobj_, interface);
932
933 if(dir==ForwardCommunication)
934 comm.forward(data_wrapper);
935 else
936 comm.backward(data_wrapper);
937}
938#endif
939
940template<class DataHandle>
941void CpGridData::communicate(DataHandle& data, InterfaceType iftype,
942 CommunicationDirection dir)
943{
944#if HAVE_MPI
945 if(data.contains(3,0))
946 {
947 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*this, data);
948 communicateCodim<0>(data_wrapper, dir, getInterface(iftype, cell_interfaces_));
949 }
950 if(data.contains(3,3))
951 {
952 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*this, data);
953 communicateCodim<3>(data_wrapper, dir, getInterface(iftype, point_interfaces_));
954 }
955#else
956 // Suppress warnings for unused arguments.
957 (void) data;
958 (void) iftype;
959 (void) dir;
960#endif
961}
962}}
963
964#if HAVE_MPI
965#include <opm/grid/cpgrid/Entity.hpp>
966#include <opm/grid/cpgrid/Indexsets.hpp>
967
968namespace Dune {
969namespace cpgrid {
970
971namespace mover
972{
973template<class T>
974class MoveBuffer
975{
976 friend class Dune::cpgrid::CpGridData;
977public:
978 void read(T& data)
979 {
980 data=buffer_[index_++];
981 }
982 void write(const T& data)
983 {
984 buffer_[index_++]=data;
985 }
986 void reset()
987 {
988 index_=0;
989 }
990 void resize(std::size_t size)
991 {
992 buffer_.resize(size);
993 index_=0;
994 }
995private:
996 std::vector<T> buffer_;
997 typename std::vector<T>::size_type index_;
998};
999template<class DataHandle,int codim>
1000struct Mover
1001{
1002};
1003
1004template<class DataHandle>
1005struct BaseMover
1006{
1007 explicit BaseMover(DataHandle& data)
1008 : data_(data)
1009 {}
1010 template<class E>
1011 void moveData(const E& from, const E& to)
1012 {
1013 std::size_t size=data_.size(from);
1014 buffer.resize(size);
1015 data_.gather(buffer, from);
1016 buffer.reset();
1017 data_.scatter(buffer, to, size);
1018 }
1019 DataHandle& data_;
1020 MoveBuffer<typename DataHandle::DataType> buffer;
1021};
1022
1023
1024template<class DataHandle>
1025struct Mover<DataHandle,0> : public BaseMover<DataHandle>
1026{
1027 Mover(DataHandle& data, CpGridData* gatherView,
1028 CpGridData* scatterView)
1029 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1030 {}
1031
1032 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1033 {
1034 Entity<0> from_entity=Entity<0>(*gatherView_, from_cell_index, true);
1035 Entity<0> to_entity=Entity<0>(*scatterView_, to_cell_index, true);
1036 this->moveData(from_entity, to_entity);
1037 }
1038 CpGridData* gatherView_;
1039 CpGridData* scatterView_;
1040};
1041
1042template<class DataHandle>
1043struct Mover<DataHandle,1> : public BaseMover<DataHandle>
1044{
1045 Mover(DataHandle& data, CpGridData* gatherView,
1046 CpGridData* scatterView)
1047 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1048 {}
1049
1050 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1051 {
1052 typedef typename OrientedEntityTable<0,1>::row_type row_type;
1053 EntityRep<0> from_cell=EntityRep<0>(from_cell_index, true);
1054 EntityRep<0> to_cell=EntityRep<0>(to_cell_index, true);
1055 const OrientedEntityTable<0,1>& table = gatherView_->cell_to_face_;
1056 row_type from_faces=table.operator[](from_cell);
1057 row_type to_faces=scatterView_->cell_to_face_[to_cell];
1058
1059 for(int i=0; i<from_faces.size(); ++i)
1060 this->moveData(from_faces[i], to_faces[i]);
1061 }
1062 CpGridData *gatherView_;
1063 CpGridData *scatterView_;
1064};
1065
1066template<class DataHandle>
1067struct Mover<DataHandle,3> : public BaseMover<DataHandle>
1068{
1069 Mover(DataHandle& data, CpGridData* gatherView,
1070 CpGridData* scatterView)
1071 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1072 {}
1073 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1074 {
1075 const std::array<int,8>& from_cell_points=
1076 gatherView_->cell_to_point_[from_cell_index];
1077 const std::array<int,8>& to_cell_points=
1078 scatterView_->cell_to_point_[to_cell_index];
1079 for(std::size_t i=0; i<8; ++i)
1080 {
1081 this->moveData(Entity<3>(*gatherView_, from_cell_points[i], true),
1082 Entity<3>(*scatterView_, to_cell_points[i], true));
1083 }
1084 }
1085 CpGridData* gatherView_;
1086 CpGridData* scatterView_;
1087};
1088
1089} // end mover namespace
1090
1091template<class DataHandle>
1092void CpGridData::scatterData(DataHandle& data, const CpGridData* global_data,
1093 const CpGridData* distributed_data, const InterfaceMap& cell_inf,
1094 const InterfaceMap& point_inf)
1095{
1096#if HAVE_MPI
1097 if(data.contains(3,0))
1098 {
1099 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*global_data, *distributed_data, data);
1100 communicateCodim<0>(data_wrapper, ForwardCommunication, cell_inf);
1101 }
1102 if(data.contains(3,3))
1103 {
1104 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*global_data, *distributed_data, data);
1105 communicateCodim<3>(data_wrapper, ForwardCommunication, point_inf);
1106 }
1107#endif
1108}
1109
1110template<int codim, class DataHandle>
1111void CpGridData::scatterCodimData(DataHandle& data, CpGridData* global_data,
1112 CpGridData* distributed_data)
1113{
1114 CpGridData *gather_view, *scatter_view;
1115 gather_view=global_data;
1116 scatter_view=distributed_data;
1117
1118 mover::Mover<DataHandle,codim> mover(data, gather_view, scatter_view);
1119
1120
1121 for(auto index=distributed_data->cellIndexSet().begin(),
1122 end = distributed_data->cellIndexSet().end();
1123 index!=end; ++index)
1124 {
1125 std::size_t from=index->global();
1126 std::size_t to=index->local();
1127 mover(from,to);
1128 }
1129}
1130
1131namespace
1132{
1133
1134template<int codim, class T, class F>
1135void visitInterior(CpGridData& distributed_data, T begin, T endit, F& func)
1136{
1137 for(T it=begin; it!=endit; ++it)
1138 {
1139 Entity<codim> entity(distributed_data, it-begin, true);
1140 PartitionType pt = entity.partitionType();
1141 if(pt==Dune::InteriorEntity)
1142 {
1143 func(*it, entity);
1144 }
1145 }
1146}
1147
1148template<class DataHandle>
1149struct GlobalIndexSizeGatherer
1150{
1151 GlobalIndexSizeGatherer(DataHandle& data_,
1152 std::vector<int>& ownedGlobalIndices_,
1153 std::vector<int>& ownedSizes_)
1154 : data(data_), ownedGlobalIndices(ownedGlobalIndices_), ownedSizes(ownedSizes_)
1155 {}
1156
1157 template<class T, class E>
1158 void operator()(T& i, E& entity)
1159 {
1160 ownedGlobalIndices.push_back(i);
1161 ownedSizes.push_back(data.size(entity));
1162 }
1163 DataHandle& data;
1164 std::vector<int>& ownedGlobalIndices;
1165 std::vector<int>& ownedSizes;
1166};
1167
1168template<class DataHandle>
1169struct DataGatherer
1170{
1171 DataGatherer(mover::MoveBuffer<typename DataHandle::DataType>& buffer_,
1172 DataHandle& data_)
1173 : buffer(buffer_), data(data_)
1174 {}
1175
1176 template<class T, class E>
1177 void operator()(T& /* it */, E& entity)
1178 {
1179 data.gather(buffer, entity);
1180 }
1181 mover::MoveBuffer<typename DataHandle::DataType>& buffer;
1182 DataHandle& data;
1183};
1184
1185}
1186
1187template<class DataHandle>
1188void CpGridData::gatherData(DataHandle& data, CpGridData* global_data,
1189 CpGridData* distributed_data)
1190{
1191#if HAVE_MPI
1192 if(data.contains(3,0))
1193 gatherCodimData<0>(data, global_data, distributed_data);
1194 if(data.contains(3,3))
1195 gatherCodimData<3>(data, global_data, distributed_data);
1196#endif
1197}
1198
1199template<int codim, class DataHandle>
1200void CpGridData::gatherCodimData(DataHandle& data, CpGridData* global_data,
1201 CpGridData* distributed_data)
1202{
1203#if HAVE_MPI
1204 // Get the mapping to global index from the global id set
1205 const std::vector<int>& mapping =
1206 distributed_data->global_id_set_->getMapping<codim>();
1207
1208 // Get the global indices and data size for the entities whose data is
1209 // to be sent, i.e. the ones that we own.
1210 std::vector<int> owned_global_indices;
1211 std::vector<int> owned_sizes;
1212 owned_global_indices.reserve(mapping.size());
1213 owned_sizes.reserve(mapping.size());
1214
1215 GlobalIndexSizeGatherer<DataHandle> gisg(data, owned_global_indices, owned_sizes);
1216 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gisg);
1217
1218 // communicate the number of indices that each processor sends
1219 int no_indices=owned_sizes.size();
1220 // We will take the address of the first elemet for MPI_Allgather below.
1221 // Make sure the containers have such an element.
1222 if ( owned_global_indices.empty() )
1223 owned_global_indices.resize(1);
1224 if ( owned_sizes.empty() )
1225 owned_sizes.resize(1);
1226 std::vector<int> no_indices_to_recv(distributed_data->ccobj_.size());
1227 distributed_data->ccobj_.allgather(&no_indices, 1, &(no_indices_to_recv[0]));
1228 // compute size of the vector capable for receiving all indices
1229 // and allgather the global indices and the sizes.
1230 // calculate displacements
1231 std::vector<int> displ(distributed_data->ccobj_.size()+1, 0);
1232 std::transform(displ.begin(), displ.end()-1, no_indices_to_recv.begin(), displ.begin()+1,
1233 std::plus<int>());
1234 int global_size=displ[displ.size()-1];//+no_indices_to_recv[displ.size()-1];
1235 std::vector<int> global_indices(global_size);
1236 std::vector<int> global_sizes(global_size);
1237 MPI_Allgatherv(&(owned_global_indices[0]), no_indices, MPITraits<int>::getType(),
1238 &(global_indices[0]), &(no_indices_to_recv[0]), &(displ[0]),
1239 MPITraits<int>::getType(),
1240 distributed_data->ccobj_);
1241 MPI_Allgatherv(&(owned_sizes[0]), no_indices, MPITraits<int>::getType(),
1242 &(global_sizes[0]), &(no_indices_to_recv[0]), &(displ[0]),
1243 MPITraits<int>::getType(),
1244 distributed_data->ccobj_);
1245 std::vector<int>().swap(owned_global_indices); // free data for reuse.
1246 // Compute the number of data items to send
1247 std::vector<int> no_data_send(distributed_data->ccobj_.size());
1248 for(typename std::vector<int>::iterator begin=no_data_send.begin(),
1249 i=begin, end=no_data_send.end(); i!=end; ++i)
1250 *i = std::accumulate(global_sizes.begin()+displ[i-begin],
1251 global_sizes.begin()+displ[i-begin+1], std::size_t());
1252 // free at least some memory that can be reused.
1253 std::vector<int>().swap(owned_sizes);
1254 // compute the displacements for receiving with allgatherv
1255 displ[0]=0;
1256 std::transform(displ.begin(), displ.end()-1, no_data_send.begin(), displ.begin()+1,
1257 std::plus<std::size_t>());
1258 // Compute the number of data items we will receive
1259 int no_data_recv = displ[displ.size()-1];//+global_sizes[displ.size()-1];
1260
1261 // Collect the data to send, gather it
1262 mover::MoveBuffer<typename DataHandle::DataType> local_data_buffer, global_data_buffer;
1263 if ( no_data_send[distributed_data->ccobj_.rank()] )
1264 {
1265 local_data_buffer.resize(no_data_send[distributed_data->ccobj_.rank()]);
1266 }
1267 else
1268 {
1269 local_data_buffer.resize(1);
1270 }
1271 global_data_buffer.resize(no_data_recv);
1272
1273 DataGatherer<DataHandle> gatherer(local_data_buffer, data);
1274 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gatherer);
1275 MPI_Allgatherv(&(local_data_buffer.buffer_[0]), no_data_send[distributed_data->ccobj_.rank()],
1276 MPITraits<typename DataHandle::DataType>::getType(),
1277 &(global_data_buffer.buffer_[0]), &(no_data_send[0]), &(displ[0]),
1278 MPITraits<typename DataHandle::DataType>::getType(),
1279 distributed_data->ccobj_);
1280 Entity2IndexDataHandle<DataHandle, codim> edata(*global_data, data);
1281 int offset=0;
1282 for(int i=0; i< codim; ++i)
1283 offset+=global_data->size(i);
1284
1285 typename std::vector<int>::const_iterator s=global_sizes.begin();
1286 for(typename std::vector<int>::const_iterator i=global_indices.begin(),
1287 end=global_indices.end();
1288 i!=end; ++s, ++i)
1289 {
1290 edata.scatter(global_data_buffer, *i-offset, *s);
1291 }
1292#endif
1293}
1294
1295} // end namespace cpgrid
1296} // end namespace Dune
1297
1298#endif
1299
1300#endif
[ provides Dune::Grid ]
Definition CpGrid.hpp:203
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:118
void postAdapt()
Clean up refinement/coarsening markers - set every element to the mark 0 which represents 'doing noth...
Definition CpGridData.cpp:1950
const cpgrid::LevelGlobalIdSet & globalIdSet() const
Get the global index set.
Definition CpGridData.hpp:579
@ MAX_DATA_PER_CELL
The maximum data items allowed per cell (DUNE < 2.5.2)
Definition CpGridData.hpp:143
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition CpGridData.hpp:175
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the grid.
Definition CpGridData.hpp:587
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition CpGridData.hpp:941
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGridData.hpp:542
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:470
int size(int codim) const
number of leaf entities per codim in this process
Definition CpGridData.cpp:149
void readEclipseFormat(const std::string &filename, bool periodic_extension, bool turn_normals=false, bool edge_conformal=false)
Read the Eclipse grid format ('grdecl').
const std::vector< int > & globalCell() const
Return global_cell_ of any level grid, or the leaf grid view (in presence of refinement).
Definition CpGridData.hpp:391
CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition CpGridData.hpp:616
~CpGridData()
Destructor.
Definition CpGridData.cpp:102
void getIJK(int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition CpGridData.cpp:1713
const IndexSet & indexSet() const
Get the index set.
Definition CpGridData.hpp:567
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
Definition CpGridData.cpp:1901
CpGridDataTraits::MPICommunicator MPICommunicator
The type of the mpi communicator.
Definition CpGridData.hpp:614
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.
Definition CpGridData.cpp:1763
void distributeGlobalGrid(CpGrid &grid, const CpGridData &view_data, const std::vector< int > &cell_part)
Redistribute a global grid.
Definition CpGridData.cpp:1489
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:455
const std::vector< double > & zcornData() const
Return the internalized zcorn copy from the grid processing, if no cells were adjusted during the min...
Definition CpGridData.hpp:560
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGridData.hpp:549
bool preAdapt()
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called af...
Definition CpGridData.cpp:1906
int getGridIdx() const
Add doc/or remove method and replace it with better approach.
Definition CpGridData.hpp:441
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, bool edge_conformal)
Read the Eclipse grid format ('grdecl').
Definition processEclipseFormat.cpp:438
bool hasNNCs(const std::vector< int > &cellIndices) const
Check all cells selected for refinement have no NNCs (no neighbor connections).
Definition CpGridData.cpp:1733
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
Mark entity for refinement or coarsening.
Definition CpGridData.cpp:1884
const cpgrid::IdSet & localIdSet() const
Get the local index set.
Definition CpGridData.hpp:573
bool adapt()
TO DO: Documentation. Triggers the grid refinement process - Currently, returns preAdapt()
Definition CpGridData.cpp:1945
const std::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
Definition CpGridData.hpp:675
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
A class design to hold a variable with a value for each entity of the given codimension,...
Definition EntityRep.hpp:272
Definition Entity.hpp:72
This class encapsulates geometry for vertices, intersections, and cells.
Definition Geometry.hpp:76
The global id set for Dune.
Definition Indexsets.hpp:487
Only needs to provide interface for doing nothing.
Definition Iterators.hpp:118
Definition Indexsets.hpp:201
Definition Indexsets.hpp:56
Definition Indexsets.hpp:371
Represents the topological relationships between sets of entities, for example cells and faces.
Definition OrientedEntityTable.hpp:139
int size() const
Returns the number of rows in the table.
Definition SparseTable.hpp:121
A class design to hold a variable with a value for each entity of the given codimension,...
Definition EntityRep.hpp:307
A SparseTable stores a table with rows of varying size as efficiently as possible.
Definition SparseTable.hpp:55
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68
Low-level corner-point processing routines and supporting data structures.
MPIHelper::MPICommunicator MPICommunicator
The type of the collective communication.
Definition CpGridDataTraits.hpp:56
AttributeSet
The type of the set of the attributes.
Definition CpGridDataTraits.hpp:66
Definition CpGridData.hpp:110
Raw corner-point specification of a particular geological model.
Definition preprocess.h:56