CpGridVanguard.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_CPGRID_VANGUARD_HPP
28#define OPM_CPGRID_VANGUARD_HPP
29
30#include <opm/common/TimingMacros.hpp>
31
38
39#include <array>
40#include <functional>
41#include <memory>
42#include <stdexcept>
43#include <tuple>
44#include <vector>
45
46namespace Opm {
47template <class TypeTag>
48class CpGridVanguard;
49}
50
51namespace Opm::Properties {
52
53namespace TTag {
55 using InheritsFrom = std::tuple<FlowBaseVanguard>;
56};
57}
58
59// declare the properties
60template<class TypeTag>
61struct Vanguard<TypeTag, TTag::CpGridVanguard> {
63};
64template<class TypeTag>
65struct Grid<TypeTag, TTag::CpGridVanguard> {
66 using type = Dune::CpGrid;
67};
68template<class TypeTag>
69struct EquilGrid<TypeTag, TTag::CpGridVanguard> {
71};
72
73} // namespace Opm::Properties
74
75namespace Opm {
76
84template <class TypeTag>
85class CpGridVanguard : public FlowBaseVanguard<TypeTag>
86 , public GenericCpGridVanguard<GetPropType<TypeTag, Properties::ElementMapper>,
87 GetPropType<TypeTag, Properties::GridView>,
88 GetPropType<TypeTag, Properties::Scalar>>
89{
90 friend class FlowBaseVanguard<TypeTag>;
92
96
97public:
103 static constexpr int dimensionworld = Grid::dimensionworld;
105 static constexpr bool waterEnabled = Indices::waterEnabled;
106 static constexpr bool gasEnabled = Indices::gasEnabled;
107 static constexpr bool oilEnabled = Indices::oilEnabled;
108private:
109 using Element = typename GridView::template Codim<0>::Entity;
110
111public:
112 explicit CpGridVanguard(Simulator& simulator)
113 : FlowBaseVanguard<TypeTag>(simulator)
114 {
115 this->checkConsistency();
117 }
118
119 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
120 {
121 const std::array<int,3> lgr_ijk = {conn.getI(), conn.getJ(), conn.getK()};
122 const auto& lgr_level = this->grid().getLgrNameToLevel().at(lgr_tag);
123 if (ParentType::lgrMappers_.has_value() == false) {
124 ParentType::lgrMappers_.emplace(this->grid().mapLocalCartesianIndexSetsToLeafIndexSet());
125 }
126 const auto& lgr_dim = this->grid().currentData()[lgr_level]->logicalCartesianSize();
127 const auto lgr_cartesian_index = (lgr_ijk[2]*lgr_dim[0]*lgr_dim[1]) + (lgr_ijk[1]*lgr_dim[0]) + (lgr_ijk[0]);
128 return ParentType::lgrMappers_.value()[lgr_level].at(lgr_cartesian_index);
129 }
134 {
135 const auto& runspec = this->eclState().runspec();
136 const auto& config = this->eclState().getSimulationConfig();
137 const auto& phases = runspec.phases();
138
139 // check for correct module setup
140 if (config.isThermal()) {
141 if (getPropValue<TypeTag, Properties::EnergyModuleType>() != EnergyModules::FullyImplicitThermal) {
142 throw std::runtime_error("Input specifies energy while simulator has disabled it, try xxx_energy");
143 }
144 } else {
145 if (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal) {
146 throw std::runtime_error("Input specifies no energy while simulator has energy, try run without _energy");
147 }
148 }
149
150 if (config.isDiffusive()) {
151 if (getPropValue<TypeTag, Properties::EnableDiffusion>() == false) {
152 throw std::runtime_error("Input specifies diffusion while simulator has disabled it, try xxx_diffusion");
153 }
154 }
155
156 if (runspec.micp()) {
157 if (getPropValue<TypeTag, Properties::EnableBioeffects>() == false) {
158 throw std::runtime_error("Input specifies MICP while simulator has it disabled");
159 }
160 }
161
162 if (runspec.biof()) {
163 if (getPropValue<TypeTag, Properties::EnableBioeffects>() == false) {
164 throw std::runtime_error("Input specifies Biofilm while simulator has it disabled");
165 }
166 }
167
168 if (phases.active(Phase::BRINE)) {
169 if (getPropValue<TypeTag, Properties::EnableBrine>() == false) {
170 throw std::runtime_error("Input specifies Brine while simulator has it disabled");
171 }
172 }
173
174 if (phases.active(Phase::POLYMER)) {
175 if (getPropValue<TypeTag, Properties::EnablePolymer>() == false) {
176 throw std::runtime_error("Input specifies Polymer while simulator has it disabled");
177 }
178 }
179
180 // checking for correct phases is more difficult TODO!
181 if (phases.active(Phase::ZFRACTION)) {
182 if (getPropValue<TypeTag, Properties::EnableExtbo>() == false) {
183 throw std::runtime_error("Input specifies ExBo while simulator has it disabled");
184 }
185 }
186 if (phases.active(Phase::FOAM)) {
187 if (getPropValue<TypeTag, Properties::EnableFoam>() == false) {
188 throw std::runtime_error("Input specifies Foam while simulator has it disabled");
189 }
190 }
191
192 if (phases.active(Phase::SOLVENT)) {
193 if (getPropValue<TypeTag, Properties::EnableSolvent>() == false) {
194 throw std::runtime_error("Input specifies Solvent while simulator has it disabled");
195 }
196 }
197 if(phases.active(Phase::WATER)){
198 if(waterEnabled == false){
199 throw std::runtime_error("Input specifies water while simulator has it disabled");
200 }
201 }
202 if(phases.active(Phase::GAS)){
203 if(gasEnabled == false){
204 throw std::runtime_error("Input specifies gas while simulator has it disabled");
205 }
206 }
207 if(phases.active(Phase::OIL)){
208 if(oilEnabled == false){
209 throw std::runtime_error("Input specifies oil while simulator has it disabled");
210 }
211 }
212
213 }
214
221 {
222 globalTrans_.reset();
223 }
224
226 {
227 assert( globalTrans_ != nullptr );
228 return *globalTrans_;
229 }
230
237 {
238#if HAVE_MPI
239 if (const auto& extPFile = this->externalPartitionFile();
240 !extPFile.empty() && (extPFile != "none"))
241 {
243 }
244
245 this->doLoadBalance_(this->edgeWeightsMethod(), this->ownersFirst(),
246 this->addCorners(), this->numOverlap(),
247 this->partitionMethod(), this->serialPartitioning(),
250 this->imbalanceTol(),
251 this->gridView(), this->schedule(),
252 this->eclState(), this->parallelWells_,
253 this->numJacobiBlocks(), this->enableEclOutput());
254#endif
255
256 this->updateGridView_();
258 this->updateCellDepths_();
259 this->updateCellThickness_();
260
261#if HAVE_MPI
262 this->distributeFieldProps_(this->eclState());
263#endif
264 }
265
269 void addLgrs()
270 {
271 // Check if input file contains Lgrs. Add them, if any.
272 // In a parallel run, this adds the LGRs on the distributed simulation grid.
273 if (const auto& lgrs = this->eclState().getLgrs(); lgrs.size() > 0) {
274 OpmLog::info("\nAdding LGRs to the grid and updating its leaf grid view");
275 this->addLgrsUpdateLeafView(lgrs, lgrs.size(), *this->grid_);
276
277 this->updateGridView_();
278 this->updateCellDepths_();
279 this->updateCellThickness_();
280
281 if (this->grid_->comm().size()>1) {
282 // Add LGRs and update the leaf grid view in the global (undistributed) simulation grid.
283 // Purpose: To enable synchronization of cell ids in 'serial mode',
284 // we rely on the "parent-to-children" cell id mapping.
285 OpmLog::info("\nAdding LGRs to the global view and updating its leaf grid view");
286 this->grid_->switchToGlobalView();
287 this->addLgrsUpdateLeafView(lgrs, lgrs.size(), *this->grid_);
288 this->grid_->switchToDistributedView();
289 this->grid_->syncDistributedGlobalCellIds();
290 }
291 }
292 }
293
294 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const {
295 return elemIndex;
296 }
297
298 unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const {
299 return elemIndex;
300 }
308 std::function<std::array<double,dimensionworld>(int)>
310 {
311 return this->cellCentroids_(this->cartesianIndexMapper(), true);
312 }
313
314 const std::vector<int>& globalCell()
315 {
316 return this->grid().globalCell();
317 }
318
319protected:
321 {
322 this->doCreateGrids_(this->edgeConformal(), this->eclState());
323 }
324
325 void allocTrans() override
326 {
327 OPM_TIMEBLOCK(allocateTrans);
329 this->gridView(),
330 this->cartesianIndexMapper(),
331 this->grid(),
332 this->cellCentroids(),
333 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal ||
334 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal,
335 getPropValue<TypeTag, Properties::EnableDiffusion>(),
336 getPropValue<TypeTag, Properties::EnableDispersion>()));
338 }
339
340 double getTransmissibility(unsigned I, unsigned J) const override
341 {
342 return globalTrans_->transmissibility(I,J);
343 }
344
345#if HAVE_MPI
346 const std::string& zoltanParams() const override
347 {
348 return this->zoltanParams_;
349 }
350
351 double zoltanPhgEdgeSizeThreshold() const override
352 {
353 return this->zoltanPhgEdgeSizeThreshold_;
354 }
355
356 const std::string& metisParams() const override
357 {
358 return this->metisParams_;
359 }
360#endif
361
362 // \Note: this globalTrans_ is used for domain decomposition and INIT file output.
363 // It only contains trans_ due to permeability and does not contain thermalHalfTrans_,
364 // diffusivity_ abd dispersivity_. The main reason is to reduce the memory usage for rank 0
365 // during parallel running.
366 std::unique_ptr<TransmissibilityType> globalTrans_;
367};
368
369} // namespace Opm
370
371#endif // OPM_CPGRID_VANGUARD_HPP
This file ensures that flow can be compiled in the presence of dune-fem.
Declares the properties required by the black oil model.
Definition: CollectDataOnIORank.hpp:49
void updateGridView_()
Definition: basevanguard.hh:132
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: CpGridVanguard.hpp:89
GetPropType< TypeTag, Properties::Indices > Indices
Definition: CpGridVanguard.hpp:104
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition: CpGridVanguard.hpp:220
void addLgrs()
Add LGRs and update Leaf Grid View in the simulation grid.
Definition: CpGridVanguard.hpp:269
const std::string & zoltanParams() const override
Definition: CpGridVanguard.hpp:346
double zoltanPhgEdgeSizeThreshold() const override
Definition: CpGridVanguard.hpp:351
const TransmissibilityType & globalTransmissibility() const
Definition: CpGridVanguard.hpp:225
void allocTrans() override
Definition: CpGridVanguard.hpp:325
static constexpr int dimensionworld
Definition: CpGridVanguard.hpp:103
static constexpr bool waterEnabled
Definition: CpGridVanguard.hpp:105
CpGridVanguard(Simulator &simulator)
Definition: CpGridVanguard.hpp:112
void createGrids_()
Definition: CpGridVanguard.hpp:320
unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const
Definition: CpGridVanguard.hpp:298
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: CpGridVanguard.hpp:100
static constexpr bool gasEnabled
Definition: CpGridVanguard.hpp:106
int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const override
Definition: CpGridVanguard.hpp:119
void checkConsistency()
Definition: CpGridVanguard.hpp:133
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
Definition: CpGridVanguard.hpp:294
const std::string & metisParams() const override
Definition: CpGridVanguard.hpp:356
std::function< std::array< double, dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition: CpGridVanguard.hpp:309
std::unique_ptr< TransmissibilityType > globalTrans_
Definition: CpGridVanguard.hpp:366
static constexpr bool oilEnabled
Definition: CpGridVanguard.hpp:107
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition: CpGridVanguard.hpp:236
Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar > TransmissibilityType
Definition: CpGridVanguard.hpp:102
const std::vector< int > & globalCell()
Definition: CpGridVanguard.hpp:314
double getTransmissibility(unsigned I, unsigned J) const override
Definition: CpGridVanguard.hpp:340
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::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:450
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
int numJacobiBlocks() const
Number of blocks in the Block-Jacobi preconditioner.
Definition: FlowGenericVanguard.hpp:248
const Schedule & schedule() const
Return a reference to the object that managages the ECL schedule.
Definition: FlowGenericVanguard.hpp:177
std::string zoltanParams_
Definition: FlowGenericVanguard.hpp:392
bool addCorners() const
Definition: FlowGenericVanguard.hpp:267
bool serialPartitioning() const
Parameter that decides if partitioning for parallel runs should be performed on a single process only...
Definition: FlowGenericVanguard.hpp:283
ParallelWellStruct parallelWells_
Information about wells in parallel.
Definition: FlowGenericVanguard.hpp:431
bool enableDistributedWells() const
Whether perforations of a well might be distributed.
Definition: FlowGenericVanguard.hpp:311
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:168
double imbalanceTol() const
Parameter that sets the imbalance tolarance, depending on the chosen partition method.
Definition: FlowGenericVanguard.hpp:290
double zoltanPhgEdgeSizeThreshold_
Definition: FlowGenericVanguard.hpp:391
const std::string & externalPartitionFile() const
Definition: FlowGenericVanguard.hpp:302
bool enableEclOutput() const
Whether or not to emit result files that are compatible with a commercial reservoir simulator.
Definition: FlowGenericVanguard.hpp:318
bool ownersFirst() const
Parameter that decide if cells owned by rank are ordered before ghost cells.
Definition: FlowGenericVanguard.hpp:260
Dune::EdgeWeightMethod edgeWeightsMethod() const
Parameter deciding the edge-weight strategy of the load balancer.
Definition: FlowGenericVanguard.hpp:242
Dune::PartitionMethod partitionMethod() const
Parameter deciding which partition method to use.
Definition: FlowGenericVanguard.hpp:276
bool edgeConformal() const
Definition: FlowGenericVanguard.hpp:263
bool allow_splitting_inactive_wells_
Definition: FlowGenericVanguard.hpp:401
std::string metisParams_
Definition: FlowGenericVanguard.hpp:394
int numOverlap() const
Definition: FlowGenericVanguard.hpp:270
Definition: GenericCpGridVanguard.hpp:79
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
void doLoadBalance_(const Dune::EdgeWeightMethod edgeWeightsMethod, const bool ownersFirst, const bool addCorners, const int numOverlap, const Dune::PartitionMethod partitionMethod, const bool serialPartitioning, const bool enableDistributedWells, const bool allowSplittingInactiveWells, const double imbalanceTol, const GetPropType< TypeTag, Properties::GridView > &gridView, const Schedule &schedule, EclipseState &eclState, FlowGenericVanguard::ParallelWellStruct &parallelWells, const int numJacobiBlocks, const bool enableEclOutput)
Distribute the simulation grid over multiple processes.
static void setExternalLoadBalancer(const std::function< std::vector< int >(const Dune::CpGrid &)> &loadBalancer)
Sets a function that returns external load balancing information when passed the grid.
Definition: GenericCpGridVanguard.hpp:125
Definition: Transmissibility.hpp:54
Definition: GenericCpGridVanguard.hpp:57
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilmodel.hh:75
Definition: blackoilbioeffectsmodules.hh:45
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
GetPropType< TypeTag, Properties::Grid > type
Definition: CpGridVanguard.hpp:70
Definition: FlowBaseVanguard.hpp:70
Dune::CpGrid type
Definition: CpGridVanguard.hpp:66
The type of the DUNE grid.
Definition: basicproperties.hh:104
Definition: CpGridVanguard.hpp:54
std::tuple< FlowBaseVanguard > InheritsFrom
Definition: CpGridVanguard.hpp:55
Property which provides a Vanguard (manages grids)
Definition: basicproperties.hh:100