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#include <opm/grid/cpgrid/LevelCartesianIndexMapper.hpp>
38
40
46
47#include <array>
48#include <cstddef>
49#include <memory>
50#include <tuple>
51#include <vector>
52
53namespace Opm {
54template <class TypeTag>
56
57} // namespace Opm
58
59namespace Opm::Properties {
60
61namespace TTag {
63 using InheritsFrom = std::tuple<FlowBaseVanguard>;
64};
65}
66
67// declare the properties
68template<class TypeTag>
69struct Vanguard<TypeTag, TTag::AluGridVanguard> {
71};
72template<class TypeTag>
73struct Grid<TypeTag, TTag::AluGridVanguard> {
74#if HAVE_MPI
75 using type = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm>;
76#else
77 using type = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridNoComm>;
78#endif //HAVE_MPI
79};
80template<class TypeTag>
81struct EquilGrid<TypeTag, TTag::AluGridVanguard> {
82 using type = Dune::CpGrid;
83};
84
85} // namespace Opm::Properties
86
87namespace Opm {
88
96template <class TypeTag>
97class AluGridVanguard : public FlowBaseVanguard<TypeTag>
98{
99 friend class FlowBaseVanguard<TypeTag>;
101
105
106public:
114 using Factory = Dune::FromToGridFactory<Grid>;
115
116 static constexpr int dimension = Grid::dimension;
117 static constexpr int dimensionworld = Grid::dimensionworld;
118
119 explicit AluGridVanguard(Simulator& simulator)
120 : FlowBaseVanguard<TypeTag>(simulator)
121 {
122 this->mpiRank = FlowGenericVanguard::comm().rank();
124 }
125
130 { return *grid_; }
131
135 const Grid& grid() const
136 { return *grid_; }
137
147 const EquilGrid& equilGrid() const
148 { return *equilGrid_; }
149
158 {
161
162 delete equilGrid_;
163 equilGrid_ = nullptr;
164 }
165
172 {
173 auto gridView = grid().leafGridView();
174 auto dataHandle = cartesianIndexMapper_->dataHandle(gridView);
175 grid().loadBalance(*dataHandle);
176
177 // communicate non-interior cells values
178 // TODO: This causes an infinite loop within fixedSize()
179 // on dune-alugrid 2.11. Since this application has never really
180 // worked in parallel, simply workaround by not calling
181 // in a sequential run for now.
182 if (grid().comm().size() > 1) {
183 grid().communicate(*dataHandle,
184 Dune::InteriorBorder_All_Interface,
185 Dune::ForwardCommunication );
186 }
187
188 if (grid().size(0))
189 {
190 globalTrans_ = std::make_unique<TransmissibilityType>(this->eclState(),
191 this->gridView(),
192 this->cartesianIndexMapper(),
193 this->grid(),
194 this->cellCentroids(),
195 getPropValue<TypeTag, Properties::EnergyModuleType>()
196 == EnergyModules::FullyImplicitThermal,
197 getPropValue<TypeTag,
199 getPropValue<TypeTag,
201 // Re-ordering for ALUGrid
203 [&](unsigned int i) { return gridEquilIdxToGridIdx(i);});
204 }
205
206 }
207
208 void addLgrs()
209 {
210 // do nothing: AluGrid with LGRs not supported yet!
211 }
212
213 template<class DataHandle>
214 void scatterData(DataHandle& /*handle*/) const
215 {
216 // not existing for this type of grid yet
217 }
218
219 template<class DataHandle>
220 void gatherData(DataHandle& /*handle*/) const
221 {
222 // not existing for this type of grid yet
223 }
224
225 template<class DataHandle, class InterfaceType, class CommunicationDirection>
226 void communicate (DataHandle& /*data*/, InterfaceType /*iftype*/,
227 CommunicationDirection /*dir*/) const
228 {
229 // not existing for this type of grid yet
230 }
231
238 {
239 globalTrans_.reset();
240 }
241
247 { return *cartesianIndexMapper_; }
248
256
261 { return *equilCartesianIndexMapper_; }
262
270 std::function<std::array<double,dimensionworld>(int)>
272 {
273 return this->cellCentroids_(this->cartesianIndexMapper(), false);
274 }
275
277 {
278 assert( globalTrans_ != nullptr );
279 return *globalTrans_;
280 }
281
282 const std::vector<int>& globalCell()
283 {
284 return cartesianCellId_;
285 }
286
287 std::vector<int> cellPartition() const
288 {
289 // not required for this type of grid yet
290 return {};
291 }
292
293 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const {
294 return equilGridToGrid_[elemIndex];
295 }
296
297 unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const {
298 return ordering_[elemIndex];
299 }
300
301protected:
303 {
304 // we use separate grid objects: one for the calculation of the initial condition
305 // via EQUIL and one for the actual simulation. The reason is that the EQUIL code
306 // cannot cope with arbitrary Dune grids and is also allergic to distributed
307 // grids.
308
310 // create the EQUIL grid
312 const EclipseGrid* input_grid = nullptr;
313 std::vector<double> global_porv;
314 // At this stage the ParallelEclipseState instance is still in global
315 // view; on rank 0 we have undistributed data for the entire grid, on
316 // the other ranks the EclipseState is empty.
317 if (mpiRank == 0) {
318 // Processing grid
319 input_grid = &this->eclState().getInputGrid();
320 global_porv = this->eclState().fieldProps().porv(true);
321 OpmLog::info("\nProcessing grid");
322 }
323
324#if HAVE_MPI
325 this->equilGrid_ = std::make_unique<Dune::CpGrid>(FlowGenericVanguard::comm());
326#else
327 this->equilGrid_ = std::make_unique<Dune::CpGrid>();
328#endif
329 // Note: removed_cells is guaranteed to be empty on ranks other than 0.
330 auto removed_cells =
331 this->equilGrid_->processEclipseFormat(input_grid,
332 &this->eclState(),
333 /*isPeriodic=*/false,
334 /*flipNormals=*/false,
335 /*clipZ=*/false);
336
337 cartesianCellId_ = this->equilGrid_->globalCell();
338
339 for (unsigned i = 0; i < dimension; ++i)
340 cartesianDimension_[i] = this->equilGrid_->logicalCartesianSize()[i];
341
342 equilCartesianIndexMapper_ = std::make_unique<EquilCartesianIndexMapper>(*equilGrid_);
343
345 // create the simulation grid
347
348 factory_ = std::make_unique<Factory>();
350 OpmLog::warning("Space Filling Curve (SFC) ordering is enabled: see flow_blackoil_alugrid for more informations on disabling/enabling SFC reordering");
351 equilGridToGrid_.resize(ordering_.size());
352 for (std::size_t index = 0; index < ordering_.size(); ++index) {
353 equilGridToGrid_[ordering_[index]] = index;
354 }
355
356 cartesianIndexMapper_ = std::make_unique<CartesianIndexMapper>(*grid_, cartesianDimension_, cartesianCellId_);
357 this->updateGridView_();
359 this->updateCellDepths_();
360 this->updateCellThickness_();
361 }
362
364 {
365 // not handling the removal of completions for this type of grid yet.
366 }
367
368 std::unique_ptr<Grid> grid_;
369 std::unique_ptr<EquilGrid> equilGrid_;
370 std::vector<int> cartesianCellId_;
371 std::vector<unsigned int> ordering_;
372 std::vector<unsigned int> equilGridToGrid_;
373 std::array<int,dimension> cartesianDimension_;
374 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
375 std::unique_ptr<EquilCartesianIndexMapper> equilCartesianIndexMapper_;
376 std::unique_ptr<Factory> factory_;
377 // \Note: this globalTrans_ is used for domain decomposition and INIT file output.
378 // It only contains trans_ due to permeability and does not contain thermalHalfTrans_,
379 // diffusivity_ abd dispersivity_. The main reason is to reduce the memory usage for rank 0
380 // during parallel running.
381 std::unique_ptr<TransmissibilityType> globalTrans_;
383};
384
385} // namespace Opm
386
387#endif // OPM_ALUGRID_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
Helper class for grid instantiation of ECL file-format using problems.
Definition: AluGridVanguard.hpp:98
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition: AluGridVanguard.hpp:171
std::vector< unsigned int > ordering_
Definition: AluGridVanguard.hpp:371
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: AluGridVanguard.hpp:246
Dune::FromToGridFactory< Grid > Factory
Definition: AluGridVanguard.hpp:114
std::unique_ptr< Grid > grid_
Definition: AluGridVanguard.hpp:368
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: AluGridVanguard.hpp:108
const LevelCartesianIndexMapper levelCartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: AluGridVanguard.hpp:254
const TransmissibilityType & globalTransmissibility() const
Definition: AluGridVanguard.hpp:276
std::function< std::array< double, dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition: AluGridVanguard.hpp:271
std::array< int, dimension > cartesianDimension_
Definition: AluGridVanguard.hpp:373
static constexpr int dimensionworld
Definition: AluGridVanguard.hpp:117
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
Definition: AluGridVanguard.hpp:293
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition: AluGridVanguard.hpp:237
std::unique_ptr< CartesianIndexMapper > cartesianIndexMapper_
Definition: AluGridVanguard.hpp:374
void communicate(DataHandle &, InterfaceType, CommunicationDirection) const
Definition: AluGridVanguard.hpp:226
AluGridVanguard(Simulator &simulator)
Definition: AluGridVanguard.hpp:119
void createGrids_()
Definition: AluGridVanguard.hpp:302
GetPropType< TypeTag, Properties::Grid > Grid
Definition: AluGridVanguard.hpp:107
const Grid & grid() const
Return a reference to the simulation grid.
Definition: AluGridVanguard.hpp:135
void filterConnections_()
Definition: AluGridVanguard.hpp:363
std::unique_ptr< TransmissibilityType > globalTrans_
Definition: AluGridVanguard.hpp:381
const EquilCartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition: AluGridVanguard.hpp:260
const std::vector< int > & globalCell()
Definition: AluGridVanguard.hpp:282
std::unique_ptr< EquilCartesianIndexMapper > equilCartesianIndexMapper_
Definition: AluGridVanguard.hpp:375
std::vector< int > cartesianCellId_
Definition: AluGridVanguard.hpp:370
std::vector< int > cellPartition() const
Definition: AluGridVanguard.hpp:287
std::unique_ptr< EquilGrid > equilGrid_
Definition: AluGridVanguard.hpp:369
unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const
Definition: AluGridVanguard.hpp:297
static constexpr int dimension
Definition: AluGridVanguard.hpp:116
void gatherData(DataHandle &) const
Definition: AluGridVanguard.hpp:220
Grid & grid()
Return a reference to the simulation grid.
Definition: AluGridVanguard.hpp:129
const EquilGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition: AluGridVanguard.hpp:147
Opm::LevelCartesianIndexMapper< Grid > LevelCartesianIndexMapper
Definition: AluGridVanguard.hpp:111
std::vector< unsigned int > equilGridToGrid_
Definition: AluGridVanguard.hpp:372
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition: AluGridVanguard.hpp:157
std::unique_ptr< Factory > factory_
Definition: AluGridVanguard.hpp:376
int mpiRank
Definition: AluGridVanguard.hpp:382
void addLgrs()
Definition: AluGridVanguard.hpp:208
void scatterData(DataHandle &) const
Definition: AluGridVanguard.hpp:214
void updateGridView_()
Definition: basevanguard.hh:132
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition: basevanguard.hh:70
Definition: AluGridVanguard.hpp:55
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:84
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:336
void updateCellThickness_()
Definition: FlowBaseVanguard.hpp:381
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:358
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:325
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:305
static Parallel::Communication & comm()
Obtain global communicator.
Definition: FlowGenericVanguard.hpp:336
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:168
Definition: EclGenericWriter.hpp:50
Definition: Transmissibility.hpp:54
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilmodel.hh:80
Definition: blackoilbioeffectsmodules.hh:45
constexpr auto getPropValue()
get the value data member of a property
Definition: propertysystem.hh:240
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:233
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:91
Enable dispersive fluxes?
Definition: multiphasebaseproperties.hh:95
Dune::CpGrid type
Definition: AluGridVanguard.hpp:82
Definition: FlowBaseVanguard.hpp:70
Dune::ALUGrid< 3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm > type
Definition: AluGridVanguard.hpp:75
The type of the DUNE grid.
Definition: basicproperties.hh:104
Definition: AluGridVanguard.hpp:62
std::tuple< FlowBaseVanguard > InheritsFrom
Definition: AluGridVanguard.hpp:63
Property which provides a Vanguard (manages grids)
Definition: basicproperties.hh:100