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 <opm/grid/common/GridEnums.hpp>
31#include <opm/grid/common/CartesianIndexMapper.hpp>
32#include <opm/grid/LookUpCellCentroid.hh>
33
34#include <opm/input/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquiferCell.hpp>
35#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
36
37#include <opm/models/discretization/common/fvbaseproperties.hh>
38#include <opm/models/io/basevanguard.hh>
39#include <opm/models/utils/parametersystem.hh>
40#include <opm/models/utils/propertysystem.hh>
41
44
45#include <array>
46#include <cstddef>
47#include <unordered_map>
48#include <vector>
49
50namespace Opm {
51template <class TypeTag>
52class FlowBaseVanguard;
53template<typename Grid, typename GridView> struct LookUpCellCentroid;
54}
55
56namespace Opm::Properties {
57
58namespace TTag {
60}
61
62// declare the properties required by the for the ecl simulator vanguard
63template<class TypeTag, class MyTypeTag>
64struct EquilGrid {
65 using type = UndefinedProperty;
66};
67template<class TypeTag, class MyTypeTag>
69 using type = UndefinedProperty;
70};
71template<class TypeTag, class MyTypeTag>
73 using type = UndefinedProperty;
74};
75template<class TypeTag, class MyTypeTag>
77 using type = UndefinedProperty;
78};
79template<class TypeTag, class MyTypeTag>
81 using type = UndefinedProperty;
82};
83template<class TypeTag, class MyTypeTag>
85 using type = UndefinedProperty;
86};
87template<class TypeTag, class MyTypeTag>
89 using type = UndefinedProperty;
90};
91template<class TypeTag, class MyTypeTag>
93 using type = UndefinedProperty;
94};
95
96#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
97template<class TypeTag, class MyTypeTag>
98struct NumJacobiBlocks {
99 using type = UndefinedProperty;
100};
101#endif
102
103template<class TypeTag, class MyTypeTag>
105 using type = UndefinedProperty;
106};
107
108template<class TypeTag, class MyTypeTag>
110 using type = UndefinedProperty;
111};
112
113template<class TypeTag, class MyTypeTag>
115 using type = UndefinedProperty;
116};
117
118template<class TypeTag, class MyTypeTag>
120 using type = UndefinedProperty;
121};
122
123template <class TypeTag, class MyTypeTag>
125{
126 using type = UndefinedProperty;
127};
128
129template<class TypeTag, class MyTypeTag>
131 using type = UndefinedProperty;
132};
133
134template<class TypeTag>
135struct IgnoreKeywords<TypeTag, TTag::FlowBaseVanguard> {
136 static constexpr auto value = "";
137};
138template<class TypeTag>
139struct EclDeckFileName<TypeTag, TTag::FlowBaseVanguard> {
140 static constexpr auto value = "";
141};
142template<class TypeTag>
143struct EclOutputInterval<TypeTag, TTag::FlowBaseVanguard> {
144 static constexpr int value = -1;
145};
146template<class TypeTag>
147struct EnableOpmRstFile<TypeTag, TTag::FlowBaseVanguard> {
148 static constexpr bool value = false;
149};
150template<class TypeTag>
151struct ParsingStrictness<TypeTag, TTag::FlowBaseVanguard> {
152 static constexpr auto value = "normal";
153};
154template<class TypeTag>
155struct InputSkipMode<TypeTag, TTag::FlowBaseVanguard> {
156 static constexpr auto value = "100";
157};
158template<class TypeTag>
159struct SchedRestart<TypeTag, TTag::FlowBaseVanguard> {
160 static constexpr bool value = false;
161};
162template<class TypeTag>
163struct EdgeWeightsMethod<TypeTag, TTag::FlowBaseVanguard> {
164 static constexpr int value = 1;
165};
166
167#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
168template<class TypeTag>
169struct NumJacobiBlocks<TypeTag, TTag::FlowBaseVanguard> {
170 static constexpr int value = 0;
171};
172#endif
173
174template<class TypeTag>
175struct OwnerCellsFirst<TypeTag, TTag::FlowBaseVanguard> {
176 static constexpr bool value = true;
177};
178template<class TypeTag>
179struct SerialPartitioning<TypeTag, TTag::FlowBaseVanguard> {
180 static constexpr bool value = false;
181};
182
183template<class TypeTag>
184struct ZoltanImbalanceTol<TypeTag, TTag::FlowBaseVanguard> {
185 static constexpr double value = 1.1;
186};
187
188template<class TypeTag>
189struct ZoltanParams<TypeTag,TTag::FlowBaseVanguard> {
190 static constexpr auto value = "graph";
191};
192
193template <class TypeTag>
195{
196 static constexpr auto* value = "";
197};
198
199template<class TypeTag>
201 static constexpr bool value = false;
202};
203
204template<class T1, class T2>
206
207// Same as in BlackoilModelParameters.hpp but for here.
208template<class TypeTag>
209struct UseMultisegmentWell<TypeTag, TTag::FlowBaseVanguard> {
210 static constexpr bool value = true;
211};
212} // namespace Opm::Properties
213
214namespace Opm {
215
221template <class TypeTag>
222class FlowBaseVanguard : public BaseVanguard<TypeTag>,
224{
225 using ParentType = BaseVanguard<TypeTag>;
226 using Implementation = GetPropType<TypeTag, Properties::Vanguard>;
227 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
228 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
229 using ElementMapper = GetPropType<TypeTag, Properties::ElementMapper>;
230
231 enum { enableExperiments = getPropValue<TypeTag, Properties::EnableExperiments>() };
232
233public:
234 using Grid = GetPropType<TypeTag, Properties::Grid>;
235 using GridView = GetPropType<TypeTag, Properties::GridView>;
236
237protected:
238 static const int dimension = Grid::dimension;
239 static const int dimensionworld = Grid::dimensionworld;
240 using Element = typename GridView::template Codim<0>::Entity;
242
243public:
247 static void registerParameters()
248 {
249 Parameters::registerParam<TypeTag, Properties::EclDeckFileName>
250 ("The name of the file which contains the ECL deck to be simulated");
251 Parameters::registerParam<TypeTag, Properties::EclOutputInterval>
252 ("The number of report steps that ought to be skipped between two writes of ECL results");
253 Parameters::registerParam<TypeTag, Properties::EnableOpmRstFile>
254 ("Include OPM-specific keywords in the ECL restart file to "
255 "enable restart of OPM simulators from these files");
256 Parameters::registerParam<TypeTag, Properties::IgnoreKeywords>
257 ("List of Eclipse keywords which should be ignored. As a ':' separated string.");
258 Parameters::registerParam<TypeTag, Properties::ParsingStrictness>
259 ("Set strictness of parsing process. Available options are "
260 "normal (stop for critical errors), "
261 "high (stop for all errors) and "
262 "low (as normal, except do not stop due to unsupported "
263 "keywords even if marked critical");
264 Parameters::registerParam<TypeTag, Properties::InputSkipMode>
265 ("Set compatibility mode for the SKIP100/SKIP300 keywords. Options are "
266 "100 (skip SKIP100..ENDSKIP, keep SKIP300..ENDSKIP) [default], "
267 "300 (skip SKIP300..ENDSKIP, keep SKIP100..ENDSKIP) and "
268 "all (skip both SKIP100..ENDSKIP and SKIP300..ENDSKIP) ");
269 Parameters::registerParam<TypeTag, Properties::SchedRestart>
270 ("When restarting: should we try to initialize wells and "
271 "groups from historical SCHEDULE section.");
272 Parameters::registerParam<TypeTag, Properties::EdgeWeightsMethod>
273 ("Choose edge-weighing strategy: 0=uniform, 1=trans, 2=log(trans).");
274
275#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
276 Parameters::registerParam<TypeTag, Properties::NumJacobiBlocks>
277 ("Number of blocks to be created for the Block-Jacobi preconditioner.");
278#endif
279
280 Parameters::registerParam<TypeTag, Properties::OwnerCellsFirst>
281 ("Order cells owned by rank before ghost/overlap cells.");
282#if HAVE_MPI
283 Parameters::registerParam<TypeTag, Properties::SerialPartitioning>
284 ("Perform partitioning for parallel runs on a single process.");
285 Parameters::registerParam<TypeTag, Properties::ZoltanImbalanceTol>
286 ("Tolerable imbalance of the loadbalancing provided by Zoltan (default: 1.1).");
287 Parameters::registerParam<TypeTag, Properties::ZoltanParams>
288 ("Configuration of Zoltan partitioner. "
289 "Valid options are: graph, hypergraph or scotch. "
290 "Alternatively, you can request a configuration to be read "
291 "from a JSON file by giving the filename here, ending with '.json.' "
292 "See https://sandialabs.github.io/Zoltan/ug_html/ug.html "
293 "for available Zoltan options.");
294 Parameters::hideParam<TypeTag, Properties::ZoltanParams>();
295 Parameters::registerParam<TypeTag, Properties::ExternalPartition>
296 ("Name of file from which to load an externally generated "
297 "partitioning of the model's active cells for MPI "
298 "distribution purposes. If empty, the built-in partitioning "
299 "method will be employed.");
300 Parameters::hideParam<TypeTag, Properties::ExternalPartition>();
301#endif
302 Parameters::registerParam<TypeTag, Properties::AllowDistributedWells>
303 ("Allow the perforations of a well to be distributed to interior of multiple processes");
304 // register here for the use in the tests without BlackoilModelParameters
305 Parameters::registerParam<TypeTag, Properties::UseMultisegmentWell>
306 ("Use the well model for multi-segment wells instead of the one for single-segment wells");
307 }
308
315 FlowBaseVanguard(Simulator& simulator)
316 : ParentType(simulator)
317 {
318 fileName_ = Parameters::get<TypeTag, Properties::EclDeckFileName>();
319 edgeWeightsMethod_ = Dune::EdgeWeightMethod(Parameters::get<TypeTag, Properties::EdgeWeightsMethod>());
320
321#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
322 numJacobiBlocks_ = Parameters::get<TypeTag, Properties::NumJacobiBlocks>();
323#endif
324
325 ownersFirst_ = Parameters::get<TypeTag, Properties::OwnerCellsFirst>();
326#if HAVE_MPI
327 serialPartitioning_ = Parameters::get<TypeTag, Properties::SerialPartitioning>();
328 zoltanImbalanceTol_ = Parameters::get<TypeTag, Properties::ZoltanImbalanceTol>();
329 zoltanParams_ = Parameters::get<TypeTag, Properties::ZoltanParams>();
330 externalPartitionFile_ = Parameters::get<TypeTag, Properties::ExternalPartition>();
331#endif
332 enableDistributedWells_ = Parameters::get<TypeTag, Properties::AllowDistributedWells>();
333 ignoredKeywords_ = Parameters::get<TypeTag, Properties::IgnoreKeywords>();
334 int output_param = Parameters::get<TypeTag, Properties::EclOutputInterval>();
335 if (output_param >= 0)
336 outputInterval_ = output_param;
337 useMultisegmentWell_ = Parameters::get<TypeTag, Properties::UseMultisegmentWell>();
338 enableExperiments_ = enableExperiments;
339
340 init();
341 }
342
344 { return asImp_().cartesianIndexMapper(); }
345
349 const std::array<int, dimension>& cartesianDimensions() const
350 { return asImp_().cartesianIndexMapper().cartesianDimensions(); }
351
355 int cartesianSize() const
356 { return asImp_().cartesianIndexMapper().cartesianSize(); }
357
362 { return asImp_().equilCartesianIndexMapper().cartesianSize(); }
363
367 unsigned cartesianIndex(unsigned compressedCellIdx) const
368 { return asImp_().cartesianIndexMapper().cartesianIndex(compressedCellIdx); }
369
373 unsigned cartesianIndex(const std::array<int,dimension>& coords) const
374 {
375 unsigned cartIndex = coords[0];
376 int factor = cartesianDimensions()[0];
377 for (unsigned i = 1; i < dimension; ++i) {
378 cartIndex += coords[i]*factor;
379 factor *= cartesianDimensions()[i];
380 }
381
382 return cartIndex;
383 }
384
391 int compressedIndex(int cartesianCellIdx) const
392 {
393 auto index_pair = cartesianToCompressed_.find(cartesianCellIdx);
394 if (index_pair!=cartesianToCompressed_.end())
395 return index_pair->second;
396 else
397 return -1;
398 }
399
406 int compressedIndexForInterior(int cartesianCellIdx) const
407 {
408 auto index_pair = cartesianToCompressed_.find(cartesianCellIdx);
409 if (index_pair == cartesianToCompressed_.end() ||
410 !is_interior_[index_pair->second])
411 {
412 return -1;
413 }
414 else
415 {
416 return index_pair->second;
417 }
418 }
425 void cartesianCoordinate(unsigned cellIdx, std::array<int,3>& ijk) const
426 { return asImp_().cartesianIndexMapper().cartesianCoordinate(cellIdx, ijk); }
427
431 unsigned equilCartesianIndex(unsigned compressedEquilCellIdx) const
432 { return asImp_().equilCartesianIndexMapper().cartesianIndex(compressedEquilCellIdx); }
433
440 void equilCartesianCoordinate(unsigned cellIdx, std::array<int,3>& ijk) const
441 { return asImp_().equilCartesianIndexMapper().cartesianCoordinate(cellIdx, ijk); }
442
443
450 Scalar cellCenterDepth(unsigned globalSpaceIdx) const
451 {
452 return cellCenterDepth_[globalSpaceIdx];
453 }
454
455 const std::vector<Scalar>& cellCenterDepths() const
456 {
457 return cellCenterDepth_;
458 }
459
467 Scalar cellThickness(unsigned globalSpaceIdx) const
468 {
469 assert(!cellThickness_.empty());
470 return cellThickness_[globalSpaceIdx];
471 }
472
478 std::size_t globalNumCells() const
479 {
480 const auto& grid = asImp_().grid();
481 if (grid.comm().size() == 1)
482 {
483 return grid.leafGridView().size(0);
484 }
485 const auto& gridView = grid.leafGridView();
486 constexpr int codim = 0;
487 constexpr auto Part = Dune::Interior_Partition;
488 auto local_cells = std::distance(gridView.template begin<codim, Part>(),
489 gridView.template end<codim, Part>());
490 return grid.comm().sum(local_cells);
491 }
492
495 }
496
497protected:
507 template<class CartMapper>
508 std::function<std::array<double,dimensionworld>(int)>
509 cellCentroids_(const CartMapper& cartMapper, const bool& isCpGrid) const
510 {
511 return [this, cartMapper, isCpGrid](int elemIdx) {
512 std::array<double,dimensionworld> centroid;
513 const auto rank = this->gridView().comm().rank();
514 const auto maxLevel = this->gridView().grid().maxLevel();
515 bool useEclipse = !isCpGrid || (isCpGrid && (rank == 0) && (maxLevel == 0));
516 if (useEclipse)
517 {
518 centroid = this->eclState().getInputGrid().getCellCenter(cartMapper.cartesianIndex(elemIdx));
519 }
520 else
521 {
522 LookUpCellCentroid<Grid,GridView> lookUpCellCentroid(this->gridView(), cartMapper, nullptr);
523 centroid = lookUpCellCentroid(elemIdx);
524 }
525 return centroid;
526 };
527 }
528
530 {
531 asImp_().createGrids_();
532 asImp_().filterConnections_();
533 std::string outputDir = Parameters::get<TypeTag, Properties::OutputDir>();
534 bool enableEclCompatFile = !Parameters::get<TypeTag, Properties::EnableOpmRstFile>();
535 asImp_().updateOutputDir_(outputDir, enableEclCompatFile);
536 asImp_().finalizeInit_();
537 }
538
540 {
541 std::size_t num_cells = asImp_().grid().leafGridView().size(0);
542 is_interior_.resize(num_cells);
543
544 ElementMapper elemMapper(this->gridView(), Dune::mcmgElementLayout());
545 for (const auto& element : elements(this->gridView()))
546 {
547 const auto elemIdx = elemMapper.index(element);
548 unsigned cartesianCellIdx = cartesianIndex(elemIdx);
549 cartesianToCompressed_[cartesianCellIdx] = elemIdx;
550 if (element.partitionType() == Dune::InteriorEntity)
551 {
552 is_interior_[elemIdx] = 1;
553 }
554 else
555 {
556 is_interior_[elemIdx] = 0;
557 }
558 }
559 }
560
562 {
563 int numCells = this->gridView().size(/*codim=*/0);
564 cellCenterDepth_.resize(numCells);
565
566 ElementMapper elemMapper(this->gridView(), Dune::mcmgElementLayout());
567
568 const auto num_aqu_cells = this->allAquiferCells();
569
570 for(const auto& element : elements(this->gridView())) {
571 const unsigned int elemIdx = elemMapper.index(element);
572 cellCenterDepth_[elemIdx] = cellCenterDepth(element);
573
574 if (!num_aqu_cells.empty()) {
575 const unsigned int global_index = cartesianIndex(elemIdx);
576 const auto search = num_aqu_cells.find(global_index);
577 if (search != num_aqu_cells.end()) {
578 // updating the cell depth using aquifer cell depth
579 cellCenterDepth_[elemIdx] = search->second->depth;
580 }
581 }
582 }
583 }
585 {
586 if (!this->drsdtconEnabled())
587 return;
588
589 ElementMapper elemMapper(this->gridView(), Dune::mcmgElementLayout());
590
591 int numElements = this->gridView().size(/*codim=*/0);
592 cellThickness_.resize(numElements);
593
594 for (const auto& elem : elements(this->gridView())) {
595 const unsigned int elemIdx = elemMapper.index(elem);
596 cellThickness_[elemIdx] = computeCellThickness(elem);
597 }
598 }
599
600private:
601 // computed from averaging cell corner depths
602 Scalar cellCenterDepth(const Element& element) const
603 {
604 typedef typename Element::Geometry Geometry;
605 static constexpr int zCoord = Element::dimension - 1;
606 Scalar zz = 0.0;
607
608 const Geometry& geometry = element.geometry();
609 const int corners = geometry.corners();
610 for (int i=0; i < corners; ++i)
611 zz += geometry.corner(i)[zCoord];
612
613 return zz/Scalar(corners);
614 }
615
616 Scalar computeCellThickness(const typename GridView::template Codim<0>::Entity& element) const
617 {
618 typedef typename Element::Geometry Geometry;
619 static constexpr int zCoord = Element::dimension - 1;
620 Scalar zz1 = 0.0;
621 Scalar zz2 = 0.0;
622
623 const Geometry& geometry = element.geometry();
624 // This code only works with CP-grid where the
625 // number of corners are 8 and
626 // also assumes that the first
627 // 4 corners are the top surface and
628 // the 4 next are the bottomn.
629 assert(geometry.corners() == 8);
630 for (int i=0; i < 4; ++i){
631 zz1 += geometry.corner(i)[zCoord];
632 zz2 += geometry.corner(i+4)[zCoord];
633 }
634 zz1 /=4;
635 zz2 /=4;
636 return zz2-zz1;
637 }
638
639 Implementation& asImp_()
640 { return *static_cast<Implementation*>(this); }
641
642 const Implementation& asImp_() const
643 { return *static_cast<const Implementation*>(this); }
644
645protected:
649 std::unordered_map<int,int> cartesianToCompressed_;
650
653 std::vector<Scalar> cellCenterDepth_;
654
657 std::vector<Scalar> cellThickness_;
658
661 std::vector<int> is_interior_;
662};
663
664} // namespace Opm
665
666#endif // OPM_FLOW_BASE_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:224
unsigned cartesianIndex(const std::array< int, dimension > &coords) const
Return the index of the cells in the logical Cartesian grid.
Definition: FlowBaseVanguard.hpp:373
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:539
int cartesianSize() const
Returns the overall number of cells of the logically Cartesian grid.
Definition: FlowBaseVanguard.hpp:355
int compressedIndex(int cartesianCellIdx) const
Return compressed index from cartesian index.
Definition: FlowBaseVanguard.hpp:391
void updateCellThickness_()
Definition: FlowBaseVanguard.hpp:584
std::vector< Scalar > cellCenterDepth_
Cell center depths.
Definition: FlowBaseVanguard.hpp:653
const std::array< int, dimension > & cartesianDimensions() const
Returns the number of logically Cartesian cells in each direction.
Definition: FlowBaseVanguard.hpp:349
unsigned equilCartesianIndex(unsigned compressedEquilCellIdx) const
Returns the Cartesian cell id given an element index for the grid used for equilibration.
Definition: FlowBaseVanguard.hpp:431
FlowBaseVanguard(Simulator &simulator)
Create the grid for problem data files which use the ECL file format.
Definition: FlowBaseVanguard.hpp:315
int compressedIndexForInterior(int cartesianCellIdx) const
Return compressed index from cartesian index only in interior.
Definition: FlowBaseVanguard.hpp:406
const std::vector< Scalar > & cellCenterDepths() const
Definition: FlowBaseVanguard.hpp:455
unsigned cartesianIndex(unsigned compressedCellIdx) const
Returns the Cartesian cell id for identifaction with ECL data.
Definition: FlowBaseVanguard.hpp:367
static const int dimension
Definition: FlowBaseVanguard.hpp:238
std::vector< int > is_interior_
Whether a cells is in the interior.
Definition: FlowBaseVanguard.hpp:661
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:561
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:529
const CartesianIndexMapper & cartesianMapper() const
Definition: FlowBaseVanguard.hpp:343
static const int dimensionworld
Definition: FlowBaseVanguard.hpp:239
void cartesianCoordinate(unsigned cellIdx, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: FlowBaseVanguard.hpp:425
std::unordered_map< int, int > cartesianToCompressed_
Mapping between cartesian and compressed cells. It is initialized the first time it is called.
Definition: FlowBaseVanguard.hpp:649
std::vector< Scalar > cellThickness_
Cell thickness.
Definition: FlowBaseVanguard.hpp:657
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowBaseVanguard.hpp:235
GetPropType< TypeTag, Properties::Grid > Grid
Definition: FlowBaseVanguard.hpp:234
Scalar cellCenterDepth(unsigned globalSpaceIdx) const
Returns the depth of a degree of freedom [m].
Definition: FlowBaseVanguard.hpp:450
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:478
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:440
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:509
Scalar cellThickness(unsigned globalSpaceIdx) const
Returns the thickness of a degree of freedom [m].
Definition: FlowBaseVanguard.hpp:467
typename GridView::template Codim< 0 >::Entity Element
Definition: FlowBaseVanguard.hpp:240
static void registerParameters()
Register the common run-time parameters for all ECL simulator vanguards.
Definition: FlowBaseVanguard.hpp:247
int equilCartesianSize() const
Returns the overall number of cells of the logically EquilCartesian grid.
Definition: FlowBaseVanguard.hpp:361
void setupCartesianToCompressed_()
Definition: FlowBaseVanguard.hpp:493
Definition: FlowGenericVanguard.hpp:60
bool enableDistributedWells_
Definition: FlowGenericVanguard.hpp:299
bool useMultisegmentWell_
Definition: FlowGenericVanguard.hpp:302
bool ownersFirst_
Definition: FlowGenericVanguard.hpp:292
std::string zoltanParams_
Definition: FlowGenericVanguard.hpp:296
std::unordered_map< std::size_t, const NumericalAquiferCell * > allAquiferCells() const
std::string ignoredKeywords_
Definition: FlowGenericVanguard.hpp:300
bool enableExperiments_
Definition: FlowGenericVanguard.hpp:303
bool serialPartitioning_
Definition: FlowGenericVanguard.hpp:294
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:120
bool drsdtconEnabled() const
Dune::EdgeWeightMethod edgeWeightsMethod_
Definition: FlowGenericVanguard.hpp:286
double zoltanImbalanceTol_
Definition: FlowGenericVanguard.hpp:295
std::string externalPartitionFile_
Definition: FlowGenericVanguard.hpp:297
std::optional< int > outputInterval_
Definition: FlowGenericVanguard.hpp:301
std::string fileName_
Definition: FlowGenericVanguard.hpp:285
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
Definition: FlowBaseVanguard.hpp:53
Definition: FlowBaseVanguard.hpp:130
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:131
Definition: BlackoilModelParameters.hpp:42
Definition: FlowBaseVanguard.hpp:84
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:85
Definition: FlowBaseVanguard.hpp:92
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:93
Definition: FlowBaseVanguard.hpp:68
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:69
Definition: FlowBaseVanguard.hpp:64
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:65
Definition: FlowBaseVanguard.hpp:125
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:126
Definition: FlowBaseVanguard.hpp:88
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:89
Definition: FlowBaseVanguard.hpp:76
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:77
Definition: FlowBaseVanguard.hpp:104
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:105
Definition: FlowBaseVanguard.hpp:72
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:73
Definition: FlowBaseVanguard.hpp:80
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:81
Definition: FlowBaseVanguard.hpp:109
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:110
Definition: FlowBaseVanguard.hpp:59
Definition: BlackoilModelParameters.hpp:90
Definition: FlowBaseVanguard.hpp:114
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:115
Definition: FlowBaseVanguard.hpp:119
UndefinedProperty type
Definition: FlowBaseVanguard.hpp:120