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
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
103
104public:
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,
188 getPropValue<TypeTag,
190 getPropValue<TypeTag,
192 // Re-ordering for ALUGrid
194 [&](unsigned int i) { return gridEquilIdxToGridIdx(i);});
195 }
196
197 }
198
199 template<class DataHandle>
200 void scatterData(DataHandle& /*handle*/) const
201 {
202 // not existing for this type of grid yet
203 }
204
205 template<class DataHandle>
206 void gatherData(DataHandle& /*handle*/) const
207 {
208 // not existing for this type of grid yet
209 }
210
211 template<class DataHandle, class InterfaceType, class CommunicationDirection>
212 void communicate (DataHandle& /*data*/, InterfaceType /*iftype*/,
213 CommunicationDirection /*dir*/) const
214 {
215 // not existing for this type of grid yet
216 }
217
224 {
225 globalTrans_.reset();
226 }
227
233 { return *cartesianIndexMapper_; }
234
239 { return *equilCartesianIndexMapper_; }
240
248 std::function<std::array<double,dimensionworld>(int)>
250 {
251 return this->cellCentroids_(this->cartesianIndexMapper(), false);
252 }
253
255 {
256 assert( globalTrans_ != nullptr );
257 return *globalTrans_;
258 }
259
260 const std::vector<int>& globalCell()
261 {
262 return cartesianCellId_;
263 }
264
265 std::vector<int> cellPartition() const
266 {
267 // not required for this type of grid yet
268 return {};
269 }
270
271 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const {
272 return equilGridToGrid_[elemIndex];
273 }
274
275 unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const {
276 return ordering_[elemIndex];
277 }
278
279protected:
281 {
282 // we use separate grid objects: one for the calculation of the initial condition
283 // via EQUIL and one for the actual simulation. The reason is that the EQUIL code
284 // cannot cope with arbitrary Dune grids and is also allergic to distributed
285 // grids.
286
288 // create the EQUIL grid
290 const EclipseGrid* input_grid = nullptr;
291 std::vector<double> global_porv;
292 // At this stage the ParallelEclipseState instance is still in global
293 // view; on rank 0 we have undistributed data for the entire grid, on
294 // the other ranks the EclipseState is empty.
295 if (mpiRank == 0) {
296 // Processing grid
297 input_grid = &this->eclState().getInputGrid();
298 global_porv = this->eclState().fieldProps().porv(true);
299 OpmLog::info("\nProcessing grid");
300 }
301
302#if HAVE_MPI
303 this->equilGrid_ = std::make_unique<Dune::CpGrid>(FlowGenericVanguard::comm());
304#else
305 this->equilGrid_ = std::make_unique<Dune::CpGrid>();
306#endif
307 // Note: removed_cells is guaranteed to be empty on ranks other than 0.
308 auto removed_cells =
309 this->equilGrid_->processEclipseFormat(input_grid,
310 &this->eclState(),
311 /*isPeriodic=*/false,
312 /*flipNormals=*/false,
313 /*clipZ=*/false);
314
315 cartesianCellId_ = this->equilGrid_->globalCell();
316
317 for (unsigned i = 0; i < dimension; ++i)
318 cartesianDimension_[i] = this->equilGrid_->logicalCartesianSize()[i];
319
320 equilCartesianIndexMapper_ = std::make_unique<EquilCartesianIndexMapper>(*equilGrid_);
321
323 // create the simulation grid
325
326 factory_ = std::make_unique<Factory>();
328 OpmLog::warning("Space Filling Curve Ordering is not yet supported: DISABLE_ALUGRID_SFC_ORDERING is enabled");
329 equilGridToGrid_.resize(ordering_.size());
330 for (std::size_t index = 0; index < ordering_.size(); ++index) {
331 equilGridToGrid_[ordering_[index]] = index;
332 }
333
334 cartesianIndexMapper_ = std::make_unique<CartesianIndexMapper>(*grid_, cartesianDimension_, cartesianCellId_);
335 this->updateGridView_();
337 this->updateCellDepths_();
338 this->updateCellThickness_();
339 }
340
342 {
343 // not handling the removal of completions for this type of grid yet.
344 }
345
346 std::unique_ptr<Grid> grid_;
347 std::unique_ptr<EquilGrid> equilGrid_;
348 std::vector<int> cartesianCellId_;
349 std::vector<unsigned int> ordering_;
350 std::vector<unsigned int> equilGridToGrid_;
351 std::array<int,dimension> cartesianDimension_;
352 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
353 std::unique_ptr<EquilCartesianIndexMapper> equilCartesianIndexMapper_;
354 std::unique_ptr<Factory> factory_;
355 // \Note: this globalTrans_ is used for domain decomposition and INIT file output.
356 // It only contains trans_ due to permeability and does not contain thermalHalfTrans_,
357 // diffusivity_ abd dispersivity_. The main reason is to reduce the memory usage for rank 0
358 // during parallel running.
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:349
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: AluGridVanguard.hpp:232
Dune::FromToGridFactory< Grid > Factory
Definition: AluGridVanguard.hpp:111
std::unique_ptr< Grid > grid_
Definition: AluGridVanguard.hpp:346
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: AluGridVanguard.hpp:106
const TransmissibilityType & globalTransmissibility() const
Definition: AluGridVanguard.hpp:254
std::function< std::array< double, dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition: AluGridVanguard.hpp:249
std::array< int, dimension > cartesianDimension_
Definition: AluGridVanguard.hpp:351
static constexpr int dimensionworld
Definition: AluGridVanguard.hpp:114
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
Definition: AluGridVanguard.hpp:271
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition: AluGridVanguard.hpp:223
std::unique_ptr< CartesianIndexMapper > cartesianIndexMapper_
Definition: AluGridVanguard.hpp:352
void communicate(DataHandle &, InterfaceType, CommunicationDirection) const
Definition: AluGridVanguard.hpp:212
AluGridVanguard(Simulator &simulator)
Definition: AluGridVanguard.hpp:116
void createGrids_()
Definition: AluGridVanguard.hpp:280
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:341
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:238
const std::vector< int > & globalCell()
Definition: AluGridVanguard.hpp:260
std::unique_ptr< EquilCartesianIndexMapper > equilCartesianIndexMapper_
Definition: AluGridVanguard.hpp:353
std::vector< int > cartesianCellId_
Definition: AluGridVanguard.hpp:348
std::vector< int > cellPartition() const
Definition: AluGridVanguard.hpp:265
std::unique_ptr< EquilGrid > equilGrid_
Definition: AluGridVanguard.hpp:347
unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const
Definition: AluGridVanguard.hpp:275
static constexpr int dimension
Definition: AluGridVanguard.hpp:113
void gatherData(DataHandle &) const
Definition: AluGridVanguard.hpp:206
Grid & grid()
Return a reference to the simulation grid.
Definition: AluGridVanguard.hpp:126
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:350
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:354
int mpiRank
Definition: AluGridVanguard.hpp:360
void scatterData(DataHandle &) const
Definition: AluGridVanguard.hpp:200
void updateGridView_()
Definition: basevanguard.hh:120
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition: basevanguard.hh:69
Definition: AluGridVanguard.hpp:53
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:83
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:350
void updateCellThickness_()
Definition: FlowBaseVanguard.hpp:395
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:372
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:338
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:318
static Parallel::Communication & comm()
Obtain global communicator.
Definition: FlowGenericVanguard.hpp:306
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:159
Definition: Transmissibility.hpp:54
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition: propertysystem.hh:242
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:79
Enable dispersive fluxes?
Definition: multiphasebaseproperties.hh:82
Specify whether energy should be considered as a conservation quantity or not.
Definition: multiphasebaseproperties.hh:76
Dune::CpGrid type
Definition: AluGridVanguard.hpp:80
Definition: FlowBaseVanguard.hpp:69
Dune::ALUGrid< 3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm > type
Definition: AluGridVanguard.hpp:73
The type of the DUNE grid.
Definition: basicproperties.hh:100
Definition: AluGridVanguard.hpp:60
std::tuple< FlowBaseVanguard > InheritsFrom
Definition: AluGridVanguard.hpp:61
Property which provides a Vanguard (manages grids)
Definition: basicproperties.hh:96