AluGridVanguard.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_ALUGRID_VANGUARD_HPP
28#define OPM_ALUGRID_VANGUARD_HPP
29
30#include <dune/alugrid/common/fromtogridfactory.hh>
31#include <dune/alugrid/dgf.hh>
32#include <dune/alugrid/grid.hh>
33
34#include <opm/common/OpmLog/OpmLog.hpp>
35
36#include <opm/grid/CpGrid.hpp>
37
38#include <opm/models/common/multiphasebaseproperties.hh>
39
44
45#include <array>
46#include <cstddef>
47#include <memory>
48#include <tuple>
49#include <vector>
50
51namespace Opm {
52template <class TypeTag>
54
55} // namespace Opm
56
57namespace Opm::Properties {
58
59namespace TTag {
61 using InheritsFrom = std::tuple<FlowBaseVanguard>;
62};
63}
64
65// declare the properties
66template<class TypeTag>
67struct Vanguard<TypeTag, TTag::AluGridVanguard> {
69};
70template<class TypeTag>
71struct Grid<TypeTag, TTag::AluGridVanguard> {
72#if HAVE_MPI
73 using type = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm>;
74#else
75 using type = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridNoComm>;
76#endif //HAVE_MPI
77};
78template<class TypeTag>
79struct EquilGrid<TypeTag, TTag::AluGridVanguard> {
80 using type = Dune::CpGrid;
81};
82
83} // namespace Opm::Properties
84
85namespace Opm {
86
94template <class TypeTag>
95class AluGridVanguard : public FlowBaseVanguard<TypeTag>
96{
97 friend class FlowBaseVanguard<TypeTag>;
99
100 using ElementMapper = GetPropType<TypeTag, Properties::ElementMapper>;
101 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
102 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
103
104public:
105 using Grid = GetPropType<TypeTag, Properties::Grid>;
106 using EquilGrid = GetPropType<TypeTag, Properties::EquilGrid>;
107 using GridView = GetPropType<TypeTag, Properties::GridView>;
111 using Factory = Dune::FromToGridFactory<Grid>;
112
113 static constexpr int dimension = Grid::dimension;
114 static constexpr int dimensionworld = Grid::dimensionworld;
115
116 AluGridVanguard(Simulator& simulator)
117 : FlowBaseVanguard<TypeTag>(simulator)
118 {
119 this->mpiRank = FlowGenericVanguard::comm().rank();
121 }
122
127 { return *grid_; }
128
132 const Grid& grid() const
133 { return *grid_; }
134
144 const EquilGrid& equilGrid() const
145 { return *equilGrid_; }
146
155 {
158
159 delete equilGrid_;
160 equilGrid_ = nullptr;
161 }
162
169 {
170 auto gridView = grid().leafGridView();
171 auto dataHandle = cartesianIndexMapper_->dataHandle(gridView);
172 grid().loadBalance(*dataHandle);
173
174 // communicate non-interior cells values
175 grid().communicate(*dataHandle,
176 Dune::InteriorBorder_All_Interface,
177 Dune::ForwardCommunication );
178
179 if (grid().size(0))
180 {
181 globalTrans_ = std::make_unique<TransmissibilityType>(this->eclState(),
182 this->gridView(),
183 this->cartesianIndexMapper(),
184 this->grid(),
185 this->cellCentroids(),
186 getPropValue<TypeTag,
187 Properties::EnableEnergy>(),
188 getPropValue<TypeTag,
189 Properties::EnableDiffusion>(),
190 getPropValue<TypeTag,
191 Properties::EnableDispersion>());
192 // Re-ordering for ALUGrid
193 globalTrans_->update(false, [&](unsigned int i) { return gridEquilIdxToGridIdx(i);});
194 }
195
196 }
197
198 template<class DataHandle>
199 void scatterData(DataHandle& /*handle*/) const
200 {
201 // not existing for this type of grid yet
202 }
203
204 template<class DataHandle>
205 void gatherData(DataHandle& /*handle*/) const
206 {
207 // not existing for this type of grid yet
208 }
209
210 template<class DataHandle, class InterfaceType, class CommunicationDirection>
211 void communicate (DataHandle& /*data*/, InterfaceType /*iftype*/,
212 CommunicationDirection /*dir*/) const
213 {
214 // not existing for this type of grid yet
215 }
216
223 {
224 globalTrans_.reset();
225 }
226
232 { return *cartesianIndexMapper_; }
233
238 { return *equilCartesianIndexMapper_; }
239
247 std::function<std::array<double,dimensionworld>(int)>
249 {
250 return this->cellCentroids_(this->cartesianIndexMapper(), false);
251 }
252
254 {
255 assert( globalTrans_ != nullptr );
256 return *globalTrans_;
257 }
258
260 {
261 globalTrans_.reset();
262 }
263
264 const std::vector<int>& globalCell()
265 {
266 return cartesianCellId_;
267 }
268
269 std::vector<int> cellPartition() const
270 {
271 // not required for this type of grid yet
272 return {};
273 }
274
275 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const {
276 return equilGridToGrid_[elemIndex];
277 }
278
279 unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const {
280 return ordering_[elemIndex];
281 }
282
283protected:
285 {
286 // we use separate grid objects: one for the calculation of the initial condition
287 // via EQUIL and one for the actual simulation. The reason is that the EQUIL code
288 // cannot cope with arbitrary Dune grids and is also allergic to distributed
289 // grids.
290
292 // create the EQUIL grid
294 const EclipseGrid* input_grid = nullptr;
295 std::vector<double> global_porv;
296 // At this stage the ParallelEclipseState instance is still in global
297 // view; on rank 0 we have undistributed data for the entire grid, on
298 // the other ranks the EclipseState is empty.
299 if (mpiRank == 0) {
300 // Processing grid
301 input_grid = &this->eclState().getInputGrid();
302 global_porv = this->eclState().fieldProps().porv(true);
303 OpmLog::info("\nProcessing grid");
304 }
305
306#if HAVE_MPI
307 this->equilGrid_ = std::make_unique<Dune::CpGrid>(FlowGenericVanguard::comm());
308#else
309 this->equilGrid_ = std::make_unique<Dune::CpGrid>();
310#endif
311 // Note: removed_cells is guaranteed to be empty on ranks other than 0.
312 auto removed_cells =
313 this->equilGrid_->processEclipseFormat(input_grid,
314 &this->eclState(),
315 /*isPeriodic=*/false,
316 /*flipNormals=*/false,
317 /*clipZ=*/false);
318
319 cartesianCellId_ = this->equilGrid_->globalCell();
320
321 for (unsigned i = 0; i < dimension; ++i)
322 cartesianDimension_[i] = this->equilGrid_->logicalCartesianSize()[i];
323
324 equilCartesianIndexMapper_ = std::make_unique<EquilCartesianIndexMapper>(*equilGrid_);
325
327 // create the simulation grid
329
330 factory_ = std::make_unique<Factory>();
332 OpmLog::warning("Space Filling Curve Ordering is not yet supported: DISABLE_ALUGRID_SFC_ORDERING is enabled");
333 equilGridToGrid_.resize(ordering_.size());
334 for (std::size_t index = 0; index < ordering_.size(); ++index) {
335 equilGridToGrid_[ordering_[index]] = index;
336 }
337
338 cartesianIndexMapper_ = std::make_unique<CartesianIndexMapper>(*grid_, cartesianDimension_, cartesianCellId_);
339 this->updateGridView_();
341 this->updateCellDepths_();
342 this->updateCellThickness_();
343 }
344
346 {
347 // not handling the removal of completions for this type of grid yet.
348 }
349
350 std::unique_ptr<Grid> grid_;
351 std::unique_ptr<EquilGrid> equilGrid_;
352 std::vector<int> cartesianCellId_;
353 std::vector<unsigned int> ordering_;
354 std::vector<unsigned int> equilGridToGrid_;
355 std::array<int,dimension> cartesianDimension_;
356 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
357 std::unique_ptr<EquilCartesianIndexMapper> equilCartesianIndexMapper_;
358 std::unique_ptr<Factory> factory_;
359 std::unique_ptr<TransmissibilityType> globalTrans_;
361};
362
363} // namespace Opm
364
365#endif // OPM_ALUGRID_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
Helper class for grid instantiation of ECL file-format using problems.
Definition: AluGridVanguard.hpp:96
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition: AluGridVanguard.hpp:168
std::vector< unsigned int > ordering_
Definition: AluGridVanguard.hpp:353
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: AluGridVanguard.hpp:231
Dune::FromToGridFactory< Grid > Factory
Definition: AluGridVanguard.hpp:111
std::unique_ptr< Grid > grid_
Definition: AluGridVanguard.hpp:350
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: AluGridVanguard.hpp:106
const TransmissibilityType & globalTransmissibility() const
Definition: AluGridVanguard.hpp:253
std::function< std::array< double, dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition: AluGridVanguard.hpp:248
std::array< int, dimension > cartesianDimension_
Definition: AluGridVanguard.hpp:355
static constexpr int dimensionworld
Definition: AluGridVanguard.hpp:114
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
Definition: AluGridVanguard.hpp:275
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition: AluGridVanguard.hpp:222
std::unique_ptr< CartesianIndexMapper > cartesianIndexMapper_
Definition: AluGridVanguard.hpp:356
void communicate(DataHandle &, InterfaceType, CommunicationDirection) const
Definition: AluGridVanguard.hpp:211
AluGridVanguard(Simulator &simulator)
Definition: AluGridVanguard.hpp:116
void createGrids_()
Definition: AluGridVanguard.hpp:284
GetPropType< TypeTag, Properties::Grid > Grid
Definition: AluGridVanguard.hpp:105
const Grid & grid() const
Return a reference to the simulation grid.
Definition: AluGridVanguard.hpp:132
void filterConnections_()
Definition: AluGridVanguard.hpp:345
std::unique_ptr< TransmissibilityType > globalTrans_
Definition: AluGridVanguard.hpp:359
const EquilCartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition: AluGridVanguard.hpp:237
const std::vector< int > & globalCell()
Definition: AluGridVanguard.hpp:264
GetPropType< TypeTag, Properties::GridView > GridView
Definition: AluGridVanguard.hpp:107
std::unique_ptr< EquilCartesianIndexMapper > equilCartesianIndexMapper_
Definition: AluGridVanguard.hpp:357
std::vector< int > cartesianCellId_
Definition: AluGridVanguard.hpp:352
std::vector< int > cellPartition() const
Definition: AluGridVanguard.hpp:269
std::unique_ptr< EquilGrid > equilGrid_
Definition: AluGridVanguard.hpp:351
unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const
Definition: AluGridVanguard.hpp:279
static constexpr int dimension
Definition: AluGridVanguard.hpp:113
void gatherData(DataHandle &) const
Definition: AluGridVanguard.hpp:205
Grid & grid()
Return a reference to the simulation grid.
Definition: AluGridVanguard.hpp:126
void releaseGlobalTransmissibility()
Definition: AluGridVanguard.hpp:259
const EquilGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition: AluGridVanguard.hpp:144
std::vector< unsigned int > equilGridToGrid_
Definition: AluGridVanguard.hpp:354
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition: AluGridVanguard.hpp:154
std::unique_ptr< Factory > factory_
Definition: AluGridVanguard.hpp:358
int mpiRank
Definition: AluGridVanguard.hpp:360
void scatterData(DataHandle &) const
Definition: AluGridVanguard.hpp:199
Definition: AluGridVanguard.hpp:53
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:216
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:526
void updateCellThickness_()
Definition: FlowBaseVanguard.hpp:571
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:548
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:516
std::function< std::array< double, dimensionworld >(int)> cellCentroids_(const CartMapper &cartMapper, const bool &isCpGrid) const
Get function to query cell centroids for a distributed grid.
Definition: FlowBaseVanguard.hpp:496
static Parallel::Communication & comm()
Obtain global communicator.
Definition: FlowGenericVanguard.hpp:255
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:120
Definition: Transmissibility.hpp:54
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
Dune::CpGrid type
Definition: AluGridVanguard.hpp:80
Definition: FlowBaseVanguard.hpp:64
Dune::ALUGrid< 3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm > type
Definition: AluGridVanguard.hpp:73
Definition: AluGridVanguard.hpp:60
std::tuple< FlowBaseVanguard > InheritsFrom
Definition: AluGridVanguard.hpp:61