FlowBaseVanguard.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_FLOW_BASE_VANGUARD_HPP
28#define OPM_FLOW_BASE_VANGUARD_HPP
29
30#include <dune/grid/common/partitionset.hh>
31
32#include <opm/grid/common/GridEnums.hpp>
33#include <opm/grid/common/CartesianIndexMapper.hpp>
34#include <opm/grid/LookUpCellCentroid.hh>
35
36#include <opm/input/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquiferCell.hpp>
37#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
38
44
47
48#include <array>
49#include <cstddef>
50#include <optional>
51#include <unordered_map>
52#include <vector>
53
54namespace Opm {
55template <class TypeTag>
56class FlowBaseVanguard;
57template<typename Grid, typename GridView> struct LookUpCellCentroid;
58}
59
60namespace Opm::Properties {
61
62namespace TTag {
63
65
66}
67
68// declare the properties required by the for the ecl simulator vanguard
69template<class TypeTag, class MyTypeTag>
70struct EquilGrid { using type = UndefinedProperty; };
71
72} // namespace Opm::Properties
73
74namespace Opm {
75
81template <class TypeTag>
82class FlowBaseVanguard : public BaseVanguard<TypeTag>,
84{
86 using Implementation = GetPropType<TypeTag, Properties::Vanguard>;
90
91 enum { enableExperiments = getPropValue<TypeTag, Properties::EnableExperiments>() };
92
93public:
96
97protected:
98 static const int dimension = Grid::dimension;
99 static const int dimensionworld = Grid::dimensionworld;
100 using Element = typename GridView::template Codim<0>::Entity;
102
103public:
107 static void registerParameters()
108 {
109 FlowGenericVanguard::registerParameters_<Scalar>();
110 }
111
118 explicit FlowBaseVanguard(Simulator& simulator)
119 : ParentType(simulator)
120 {
121#if HAVE_MPI
122 imbalanceTol_ = Parameters::Get<Parameters::ImbalanceTol<Scalar>>();
123
124 zoltanImbalanceTolSet_ = Parameters::IsSet<Parameters::ZoltanImbalanceTol<Scalar>>();
125 zoltanImbalanceTol_ = Parameters::Get<Parameters::ZoltanImbalanceTol<Scalar>>();
126#endif
127 enableExperiments_ = enableExperiments;
128
129 init();
130 }
131
133 { return asImp_().cartesianIndexMapper(); }
134
138 const std::array<int, dimension>& cartesianDimensions() const
139 { return asImp_().cartesianIndexMapper().cartesianDimensions(); }
140
144 int cartesianSize() const
145 { return asImp_().cartesianIndexMapper().cartesianSize(); }
146
151 { return asImp_().equilCartesianIndexMapper().cartesianSize(); }
152
156 unsigned cartesianIndex(unsigned compressedCellIdx) const
157 { return asImp_().cartesianIndexMapper().cartesianIndex(compressedCellIdx); }
158
162 unsigned cartesianIndex(const std::array<int,dimension>& coords) const
163 {
164 unsigned cartIndex = coords[0];
165 int factor = cartesianDimensions()[0];
166 for (unsigned i = 1; i < dimension; ++i) {
167 cartIndex += coords[i]*factor;
168 factor *= cartesianDimensions()[i];
169 }
170
171 return cartIndex;
172 }
173
180 int compressedIndex(int cartesianCellIdx) const
181 {
182 auto index_pair = cartesianToCompressed_.find(cartesianCellIdx);
183 if (index_pair!=cartesianToCompressed_.end())
184 return index_pair->second;
185 else
186 return -1;
187 }
188
195 int compressedIndexForInterior(int cartesianCellIdx) const
196 {
197 auto index_pair = cartesianToCompressed_.find(cartesianCellIdx);
198 if (index_pair == cartesianToCompressed_.end() ||
199 !is_interior_[index_pair->second])
200 {
201 return -1;
202 }
203 else
204 {
205 return index_pair->second;
206 }
207 }
208
209 virtual int compressedIndexForInteriorLGR([[maybe_unused]] const std::string& lgr_tag,
210 [[maybe_unused]] const Connection& conn) const
211 {
212 throw std::runtime_error("compressedIndexForInteriorLGR not implemented");
213 }
214
221 void cartesianCoordinate(unsigned cellIdx, std::array<int,3>& ijk) const
222 { return asImp_().cartesianIndexMapper().cartesianCoordinate(cellIdx, ijk); }
223
227 unsigned equilCartesianIndex(unsigned compressedEquilCellIdx) const
228 { return asImp_().equilCartesianIndexMapper().cartesianIndex(compressedEquilCellIdx); }
229
236 void equilCartesianCoordinate(unsigned cellIdx, std::array<int,3>& ijk) const
237 { return asImp_().equilCartesianIndexMapper().cartesianCoordinate(cellIdx, ijk); }
238
239
246 Scalar cellCenterDepth(unsigned globalSpaceIdx) const
247 {
248 return cellCenterDepth_[globalSpaceIdx];
249 }
250
251 const std::vector<Scalar>& cellCenterDepths() const
252 {
253 return cellCenterDepth_;
254 }
255
263 Scalar cellThickness(unsigned globalSpaceIdx) const
264 {
265 assert(!cellThickness_.empty());
266 return cellThickness_[globalSpaceIdx];
267 }
268
274 std::size_t globalNumCells() const
275 {
276 const auto& grid = asImp_().grid();
277 if (grid.comm().size() == 1)
278 {
279 return grid.leafGridView().size(0);
280 }
281 const auto& gridView = grid.leafGridView();
282 constexpr int codim = 0;
283 constexpr auto Part = Dune::Interior_Partition;
284 auto local_cells = std::distance(gridView.template begin<codim, Part>(),
285 gridView.template end<codim, Part>());
286 return grid.comm().sum(local_cells);
287 }
288
289protected:
299 template<class CartMapper>
300 std::function<std::array<double,dimensionworld>(int)>
301 cellCentroids_(const CartMapper& cartMapper, const bool& isCpGrid) const
302 {
303 return [this, cartMapper, isCpGrid](int elemIdx) {
304 std::array<double,dimensionworld> centroid;
305 const auto rank = this->gridView().comm().rank();
306 const auto maxLevel = this->gridView().grid().maxLevel();
307 bool useEclipse = !isCpGrid || (isCpGrid && (rank == 0) && (maxLevel == 0));
308 if (useEclipse)
309 {
310 centroid = this->eclState().getInputGrid().getCellCenter(cartMapper.cartesianIndex(elemIdx));
311 }
312 else
313 {
314 LookUpCellCentroid<Grid,GridView> lookUpCellCentroid(this->gridView(), cartMapper, nullptr);
315 centroid = lookUpCellCentroid(elemIdx);
316 }
317 return centroid;
318 };
319 }
320
322 {
323 asImp_().createGrids_();
324 asImp_().filterConnections_();
325 std::string outputDir = Parameters::Get<Parameters::OutputDir>();
326 bool enableEclCompatFile = !Parameters::Get<Parameters::EnableOpmRstFile>();
327 asImp_().updateOutputDir_(outputDir, enableEclCompatFile);
328 const std::string& dryRunString = Parameters::Get<Parameters::EnableDryRun>();
329 asImp_().updateNOSIM_(dryRunString);
330 asImp_().finalizeInit_();
331 }
332
334 {
335 std::size_t num_cells = asImp_().grid().leafGridView().size(0);
336 is_interior_.resize(num_cells);
337
338 ElementMapper elemMapper(this->gridView(), Dune::mcmgElementLayout());
339 for (const auto& element : elements(this->gridView()))
340 {
341 const auto elemIdx = elemMapper.index(element);
342 unsigned cartesianCellIdx = cartesianIndex(elemIdx);
343 cartesianToCompressed_[cartesianCellIdx] = elemIdx;
344 if (element.partitionType() == Dune::InteriorEntity)
345 {
346 is_interior_[elemIdx] = 1;
347 }
348 else
349 {
350 is_interior_[elemIdx] = 0;
351 }
352 }
353 }
354
356 {
357 int numCells = this->gridView().size(/*codim=*/0);
358 cellCenterDepth_.resize(numCells);
359
360 ElementMapper elemMapper(this->gridView(), Dune::mcmgElementLayout());
361
362 const auto num_aqu_cells = this->allAquiferCells();
363
364 for(const auto& element : elements(this->gridView())) {
365 const unsigned int elemIdx = elemMapper.index(element);
366 cellCenterDepth_[elemIdx] = cellCenterDepth(element);
367
368 if (!num_aqu_cells.empty()) {
369 const unsigned int global_index = cartesianIndex(elemIdx);
370 const auto search = num_aqu_cells.find(global_index);
371 if (search != num_aqu_cells.end()) {
372 // updating the cell depth using aquifer cell depth
373 cellCenterDepth_[elemIdx] = search->second->depth;
374 }
375 }
376 }
377 }
379 {
380 if (!this->drsdtconEnabled())
381 return;
382
383 ElementMapper elemMapper(this->gridView(), Dune::mcmgElementLayout());
384
385 int numElements = this->gridView().size(/*codim=*/0);
386 cellThickness_.resize(numElements);
387
388 for (const auto& elem : elements(this->gridView())) {
389 const unsigned int elemIdx = elemMapper.index(elem);
390 cellThickness_[elemIdx] = computeCellThickness(elem);
391 }
392 }
393
394private:
395 // computed from averaging cell corner depths
396 Scalar cellCenterDepth(const Element& element) const
397 {
398 typedef typename Element::Geometry Geometry;
399 static constexpr int zCoord = Element::dimension - 1;
400 Scalar zz = 0.0;
401
402 const Geometry& geometry = element.geometry();
403 const int corners = geometry.corners();
404 for (int i=0; i < corners; ++i)
405 zz += geometry.corner(i)[zCoord];
406
407 return zz/Scalar(corners);
408 }
409
410 Scalar computeCellThickness(const typename GridView::template Codim<0>::Entity& element) const
411 {
412 typedef typename Element::Geometry Geometry;
413 static constexpr int zCoord = Element::dimension - 1;
414 Scalar zz1 = 0.0;
415 Scalar zz2 = 0.0;
416
417 const Geometry& geometry = element.geometry();
418 // This code only works with CP-grid where the
419 // number of corners are 8 and
420 // also assumes that the first
421 // 4 corners are the top surface and
422 // the 4 next are the bottomn.
423 assert(geometry.corners() == 8);
424 for (int i=0; i < 4; ++i){
425 zz1 += geometry.corner(i)[zCoord];
426 zz2 += geometry.corner(i+4)[zCoord];
427 }
428 zz1 /=4;
429 zz2 /=4;
430 return zz2-zz1;
431 }
432
433 Implementation& asImp_()
434 { return *static_cast<Implementation*>(this); }
435
436 const Implementation& asImp_() const
437 { return *static_cast<const Implementation*>(this); }
438
439protected:
443 std::unordered_map<int,int> cartesianToCompressed_;
447 mutable std::optional<std::vector<std::unordered_map<std::size_t, std::size_t>>> lgrMappers_;
450 std::vector<Scalar> cellCenterDepth_;
451
454 std::vector<Scalar> cellThickness_;
455
458 std::vector<int> is_interior_;
459};
460
461} // namespace Opm
462
463#endif // OPM_FLOW_BASE_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
Provides the base class for most (all?) simulator vanguards.
Definition: basevanguard.hh:50
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition: basevanguard.hh:70
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:84
unsigned cartesianIndex(const std::array< int, dimension > &coords) const
Return the index of the cells in the logical Cartesian grid.
Definition: FlowBaseVanguard.hpp:162
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:333
int cartesianSize() const
Returns the overall number of cells of the logically Cartesian grid.
Definition: FlowBaseVanguard.hpp:144
int compressedIndex(int cartesianCellIdx) const
Return compressed index from cartesian index.
Definition: FlowBaseVanguard.hpp:180
void updateCellThickness_()
Definition: FlowBaseVanguard.hpp:378
std::vector< Scalar > cellCenterDepth_
Cell center depths.
Definition: FlowBaseVanguard.hpp:450
const std::array< int, dimension > & cartesianDimensions() const
Returns the number of logically Cartesian cells in each direction.
Definition: FlowBaseVanguard.hpp:138
unsigned equilCartesianIndex(unsigned compressedEquilCellIdx) const
Returns the Cartesian cell id given an element index for the grid used for equilibration.
Definition: FlowBaseVanguard.hpp:227
FlowBaseVanguard(Simulator &simulator)
Create the grid for problem data files which use the ECL file format.
Definition: FlowBaseVanguard.hpp:118
int compressedIndexForInterior(int cartesianCellIdx) const
Return compressed index from cartesian index only in interior.
Definition: FlowBaseVanguard.hpp:195
const std::vector< Scalar > & cellCenterDepths() const
Definition: FlowBaseVanguard.hpp:251
unsigned cartesianIndex(unsigned compressedCellIdx) const
Returns the Cartesian cell id for identifaction with ECL data.
Definition: FlowBaseVanguard.hpp:156
static const int dimension
Definition: FlowBaseVanguard.hpp:98
std::vector< int > is_interior_
Whether a cells is in the interior.
Definition: FlowBaseVanguard.hpp:458
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:355
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:321
const CartesianIndexMapper & cartesianMapper() const
Definition: FlowBaseVanguard.hpp:132
static const int dimensionworld
Definition: FlowBaseVanguard.hpp:99
void cartesianCoordinate(unsigned cellIdx, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: FlowBaseVanguard.hpp:221
std::unordered_map< int, int > cartesianToCompressed_
Mapping between cartesian and compressed cells. It is initialized the first time it is called.
Definition: FlowBaseVanguard.hpp:443
std::vector< Scalar > cellThickness_
Cell thickness.
Definition: FlowBaseVanguard.hpp:454
std::optional< std::vector< std::unordered_map< std::size_t, std::size_t > > > lgrMappers_
Mapping between LGR cartesian and compressed cells. It is initialized as it is called.
Definition: FlowBaseVanguard.hpp:447
Scalar cellCenterDepth(unsigned globalSpaceIdx) const
Returns the depth of a degree of freedom [m].
Definition: FlowBaseVanguard.hpp:246
virtual int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const
Definition: FlowBaseVanguard.hpp:209
std::size_t globalNumCells() const
Get the number of cells in the global leaf grid view. \warn This is a collective operation that needs...
Definition: FlowBaseVanguard.hpp:274
void equilCartesianCoordinate(unsigned cellIdx, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell of the grid used for EQUIL.
Definition: FlowBaseVanguard.hpp:236
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:301
Scalar cellThickness(unsigned globalSpaceIdx) const
Returns the thickness of a degree of freedom [m].
Definition: FlowBaseVanguard.hpp:263
typename GridView::template Codim< 0 >::Entity Element
Definition: FlowBaseVanguard.hpp:100
static void registerParameters()
Register the common run-time parameters for all ECL simulator vanguards.
Definition: FlowBaseVanguard.hpp:107
int equilCartesianSize() const
Returns the overall number of cells of the logically EquilCartesian grid.
Definition: FlowBaseVanguard.hpp:150
Definition: FlowGenericVanguard.hpp:107
std::unordered_map< std::size_t, const NumericalAquiferCell * > allAquiferCells() const
bool zoltanImbalanceTolSet_
Definition: FlowGenericVanguard.hpp:383
bool enableExperiments_
Definition: FlowGenericVanguard.hpp:400
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:167
bool drsdtconEnabled() const
double zoltanImbalanceTol_
Definition: FlowGenericVanguard.hpp:384
double imbalanceTol_
Definition: FlowGenericVanguard.hpp:381
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilmodel.hh:79
Definition: blackoilboundaryratevector.hh:39
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
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: FlowBaseVanguard.hpp:57
Definition: FlowBaseVanguard.hpp:70
Definition: FlowBaseVanguard.hpp:64
a tag to mark properties as undefined
Definition: propertysystem.hh:38