opm-simulators
Loading...
Searching...
No Matches
BlackoilModel.hpp
1/*
2 Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3 Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
4 Copyright 2014, 2015 Statoil ASA.
5 Copyright 2015 NTNU
6 Copyright 2015, 2016, 2017 IRIS AS
7
8 This file is part of the Open Porous Media project (OPM).
9
10 OPM is free software: you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation, either version 3 of the License, or
13 (at your option) any later version.
14
15 OPM is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 GNU General Public License for more details.
19
20 You should have received a copy of the GNU General Public License
21 along with OPM. If not, see <http://www.gnu.org/licenses/>.
22*/
23
24#ifndef OPM_BLACKOILMODEL_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_HEADER_INCLUDED
26
27#include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
28
29#include <opm/simulators/flow/BlackoilModelConvergenceMonitor.hpp>
30#include <opm/simulators/flow/BlackoilModelNldd.hpp>
31#include <opm/simulators/flow/BlackoilModelProperties.hpp>
33#include <opm/simulators/flow/RSTConv.hpp>
34
35#include <opm/simulators/linalg/ISTLSolver.hpp>
36
37#include <opm/simulators/timestepping/ConvergenceReport.hpp>
38#include <opm/simulators/timestepping/SimulatorReport.hpp>
39#include <opm/simulators/timestepping/SimulatorTimer.hpp>
40
41#include <opm/simulators/utils/ComponentName.hpp>
42
43#include <opm/simulators/wells/BlackoilWellModel.hpp>
44
45#include <memory>
46#include <tuple>
47#include <vector>
48
49#include <fmt/format.h>
50
51namespace Opm {
52
59template <class TypeTag>
61{
62public:
63 // --------- Types and enums ---------
76 using ModelParameters = BlackoilModelParameters<Scalar>;
77
78 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
79
80 static constexpr int numEq = Indices::numEq;
81 static constexpr int contiSolventEqIdx = Indices::contiSolventEqIdx;
82 static constexpr int contiZfracEqIdx = Indices::contiZfracEqIdx;
83 static constexpr int contiPolymerEqIdx = Indices::contiPolymerEqIdx;
84 static constexpr int contiEnergyEqIdx = Indices::contiEnergyEqIdx;
85 static constexpr int contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx;
86 static constexpr int contiFoamEqIdx = Indices::contiFoamEqIdx;
87 static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx;
88 static constexpr int contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
89 static constexpr int contiOxygenEqIdx = Indices::contiOxygenEqIdx;
90 static constexpr int contiUreaEqIdx = Indices::contiUreaEqIdx;
91 static constexpr int contiBiofilmEqIdx = Indices::contiBiofilmEqIdx;
92 static constexpr int contiCalciteEqIdx = Indices::contiCalciteEqIdx;
93 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
94 static constexpr int zFractionIdx = Indices::zFractionIdx;
95 static constexpr int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
96 static constexpr int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
97 static constexpr int temperatureIdx = Indices::temperatureIdx;
98 static constexpr int foamConcentrationIdx = Indices::foamConcentrationIdx;
99 static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx;
100 static constexpr int microbialConcentrationIdx = Indices::microbialConcentrationIdx;
101 static constexpr int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
102 static constexpr int ureaConcentrationIdx = Indices::ureaConcentrationIdx;
103 static constexpr int biofilmVolumeFractionIdx = Indices::biofilmVolumeFractionIdx;
104 static constexpr int calciteVolumeFractionIdx = Indices::calciteVolumeFractionIdx;
105
106 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
107 using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock;
108 using Mat = typename SparseMatrixAdapter::IstlMatrix;
109 using BVector = Dune::BlockVector<VectorBlockType>;
110
111 using ComponentName = ::Opm::ComponentName<FluidSystem,Indices>;
112
113 // Helper structs
115 std::pair<std::vector<double>, std::vector<int>> cnvPvSplit;
116 std::vector<unsigned> ixCells;
117 };
118
120 Scalar dPMax = 0.0;
121 Scalar dSMax = 0.0;
122 Scalar dRsMax = 0.0;
123 Scalar dRvMax = 0.0;
124 };
125
126 // Output debug flags for which tolerances used
127 enum class DebugFlags {
128 STRICT = 0,
129 RELAXED = 1,
130 TUNINGDP = 2
131 };
132
133 // --------- Public methods ---------
134
142 BlackoilModel(Simulator& simulator,
143 const ModelParameters& param,
144 BlackoilWellModel<TypeTag>& well_model,
145 const bool terminal_output);
146
147 bool isParallel() const
148 { return grid_.comm().size() > 1; }
149
150 const EclipseState& eclState() const
151 { return simulator_.vanguard().eclState(); }
152
155 SimulatorReportSingle prepareStep(const SimulatorTimerInterface& timer);
156
157 void initialLinearization(SimulatorReportSingle& report,
158 const int minIter,
159 const int maxIter,
160 const SimulatorTimerInterface& timer);
161
169 template <class NonlinearSolverType>
170 SimulatorReportSingle nonlinearIteration(const SimulatorTimerInterface& timer,
171 NonlinearSolverType& nonlinear_solver);
172
173 template <class NonlinearSolverType>
174 SimulatorReportSingle nonlinearIterationNewton(const SimulatorTimerInterface& timer,
175 NonlinearSolverType& nonlinear_solver);
176
178 SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface& timer);
179
180 // compute the "relative" change of the solution between time steps
181 Scalar relativeChange() const;
182
185 { return simulator_.model().newtonMethod().linearSolver().iterations (); }
186
187 // Obtain reference to linear solver setup time
188 double& linearSolveSetupTime()
189 { return linear_solve_setup_time_; }
190
193 void solveJacobianSystem(BVector& x);
194
196 void updateSolution(const BVector& dx);
197
200 void storeSolutionUpdate(const BVector& dx);
201 MaxSolutionUpdateData getMaxSolutionUpdate(const std::vector<unsigned>& ixCells);
202
205 { return terminal_output_; }
206
207 std::tuple<Scalar,Scalar>
208 convergenceReduction(Parallel::Communication comm,
209 const Scalar pvSumLocal,
210 const Scalar numAquiferPvSumLocal,
211 std::vector<Scalar>& R_sum,
212 std::vector<Scalar>& maxCoeff,
213 std::vector<Scalar>& B_avg);
214
218 std::pair<Scalar,Scalar>
219 localConvergenceData(std::vector<Scalar>& R_sum,
220 std::vector<Scalar>& maxCoeff,
221 std::vector<Scalar>& B_avg,
222 std::vector<int>& maxCoeffCell);
223
228 CnvPvSplitData characteriseCnvPvSplit(const std::vector<Scalar>& B_avg, const double dt);
229
232 void convergencePerCell(const std::vector<Scalar>& B_avg,
233 const double dt,
234 const double tol_cnv,
235 const double tol_cnv_energy);
236
237 void updateTUNING(const Tuning& tuning);
238 void updateTUNINGDP(const TuningDp& tuning_dp);
239
241 getReservoirConvergence(const double reportTime,
242 const double dt,
243 const int maxIter,
244 std::vector<Scalar>& B_avg,
245 std::vector<Scalar>& residual_norms);
246
254 const int maxIter,
255 std::vector<Scalar>& residual_norms);
256
258 int numPhases() const
259 { return Indices::numPhases; }
260
262 template<class T>
263 std::vector<std::vector<Scalar> >
264 computeFluidInPlace(const T&, const std::vector<int>& fipnum) const
265 { return this->computeFluidInPlace(fipnum); }
266
268 std::vector<std::vector<Scalar> >
269 computeFluidInPlace(const std::vector<int>& /*fipnum*/) const;
270
271 const Simulator& simulator() const
272 { return simulator_; }
273
274 Simulator& simulator()
275 { return simulator_; }
276
279 { return failureReport_; }
280
283
285 const std::vector<SimulatorReport>& domainAccumulatedReports() const;
286
288 void writeNonlinearIterationsPerCell(const std::filesystem::path& odir) const;
289
290 const std::vector<StepReport>& stepReports() const
291 { return convergence_reports_; }
292
295 {
296 convergence_reports_.back().report.pop_back();
297 residual_norms_history_.pop_back();
298 }
299
300 void writePartitions(const std::filesystem::path& odir) const;
301
305 { return well_model_; }
306
308 wellModel() const
309 { return well_model_; }
310
311 void beginReportStep()
312 { simulator_.problem().beginEpisode(); }
313
314 void endReportStep()
315 { simulator_.problem().endEpisode(); }
316
317 template<class FluidState, class Residual>
318 void getMaxCoeff(const unsigned cell_idx,
319 const IntensiveQuantities& intQuants,
320 const FluidState& fs,
321 const Residual& modelResid,
322 const Scalar pvValue,
323 std::vector<Scalar>& B_avg,
324 std::vector<Scalar>& R_sum,
325 std::vector<Scalar>& maxCoeff,
326 std::vector<int>& maxCoeffCell);
327
329 const ModelParameters& param() const
330 { return param_; }
331
333 const ComponentName& compNames() const
334 { return compNames_; }
335
337 bool hasNlddSolver() const
338 { return nlddSolver_ != nullptr; }
339
340protected:
341 // --------- Data members ---------
342 Simulator& simulator_;
343 const Grid& grid_;
344 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
345 static constexpr bool has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>();
346 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
347 static constexpr bool has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>();
348 static constexpr bool has_energy_ = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal;
349 static constexpr bool has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>();
350 static constexpr bool has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>();
351 static constexpr bool has_bioeffects_ = getPropValue<TypeTag, Properties::EnableBioeffects>();
352 static constexpr bool has_micp_ = Indices::enableMICP;
353
354 ModelParameters param_;
355 SimulatorReportSingle failureReport_;
356
357 // Well Model
358 BlackoilWellModel<TypeTag>& well_model_;
359
363 long int global_nc_;
364
365 std::vector<std::vector<Scalar>> residual_norms_history_;
366 Scalar current_relaxation_;
367 BVector dx_old_;
368
369 SolutionVector solUpd_;
370
371 std::vector<StepReport> convergence_reports_;
372 ComponentName compNames_{};
373
374 std::unique_ptr<BlackoilModelNldd<TypeTag>> nlddSolver_;
376
377private:
378 Scalar dpMaxRel() const { return param_.dp_max_rel_; }
379 Scalar dsMax() const { return param_.ds_max_; }
380 Scalar drMaxRel() const { return param_.dr_max_rel_; }
381 Scalar maxResidualAllowed() const { return param_.max_residual_allowed_; }
382 double linear_solve_setup_time_;
383 std::vector<bool> wasSwitched_;
384};
385
386} // namespace Opm
387
388#include <opm/simulators/flow/BlackoilModel_impl.hpp>
389
390#endif // OPM_BLACKOILMODEL_HEADER_INCLUDED
Implementation of penalty cards for three-phase black oil.
Definition BlackoilModelConvergenceMonitor.hpp:35
BlackoilModel(Simulator &simulator, const ModelParameters &param, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Construct the model.
Definition BlackoilModel_impl.hpp:78
std::vector< std::vector< Scalar > > computeFluidInPlace(const T &, const std::vector< int > &fipnum) const
Wrapper required due to not following generic API.
Definition BlackoilModel.hpp:264
int linearIterationsLastSolve() const
Number of linear iterations used in last call to solveJacobianSystem().
Definition BlackoilModel.hpp:184
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir) const
Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.
Definition BlackoilModel_impl.hpp:1232
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition BlackoilModel_impl.hpp:1222
std::pair< Scalar, Scalar > localConvergenceData(std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg, std::vector< int > &maxCoeffCell)
Get reservoir quantities on this process needed for convergence calculations.
Definition BlackoilModel_impl.hpp:718
bool terminalOutputEnabled() const
Return true if output to cout is wanted.
Definition BlackoilModel.hpp:204
int numPhases() const
The number of active fluid phases in the model.
Definition BlackoilModel.hpp:258
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition BlackoilModel_impl.hpp:1211
const ComponentName & compNames() const
Returns const reference to component names.
Definition BlackoilModel.hpp:333
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition BlackoilModel_impl.hpp:531
long int global_nc_
The number of cells of the global grid.
Definition BlackoilModel.hpp:363
const ModelParameters & param() const
Returns const reference to model parameters.
Definition BlackoilModel.hpp:329
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition BlackoilModel.hpp:374
bool hasNlddSolver() const
Returns true if an NLDD solver exists.
Definition BlackoilModel.hpp:337
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int maxIter, std::vector< Scalar > &residual_norms)
Compute convergence based on total mass balance (tol_mb) and maximum residual mass balance (tol_cnv).
Definition BlackoilModel_impl.hpp:1175
const SimulatorReportSingle & failureReport() const
return the statistics if the nonlinearIteration() method failed
Definition BlackoilModel.hpp:278
CnvPvSplitData characteriseCnvPvSplit(const std::vector< Scalar > &B_avg, const double dt)
Compute pore-volume/cell count split among "converged", "relaxed converged", "unconverged" cells base...
Definition BlackoilModel_impl.hpp:769
BlackoilWellModel< TypeTag > & wellModel()
return the StandardWells object
Definition BlackoilModel.hpp:304
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &timer)
Assemble the residual and Jacobian of the nonlinear system.
Definition BlackoilModel_impl.hpp:370
void convergencePerCell(const std::vector< Scalar > &B_avg, const double dt, const double tol_cnv, const double tol_cnv_energy)
Compute the number of Newtons required by each cell in order to satisfy the solution change convergen...
Definition BlackoilModel_impl.hpp:1129
SimulatorReportSingle nonlinearIteration(const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Called once per nonlinear iteration.
Definition BlackoilModel_impl.hpp:245
void solveJacobianSystem(BVector &x)
Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.
Definition BlackoilModel_impl.hpp:469
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Called once before each time step.
Definition BlackoilModel_impl.hpp:113
void prepareStoringSolutionUpdate()
Get solution update vector as a PrimaryVariable.
Definition BlackoilModel_impl.hpp:566
void popLastConvergenceReport()
Remove the last convergence report entry and residual norms history entry.
Definition BlackoilModel.hpp:294
bool terminal_output_
Whether we print something to std::cout.
Definition BlackoilModel.hpp:361
Class for handling the blackoil well model.
Definition BlackoilWellModel.hpp:98
Definition ComponentName.hpp:34
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Interface class for SimulatorTimer objects, to be improved.
Definition SimulatorTimerInterface.hpp:34
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void writePartitions(const std::filesystem::path &odir, const std::vector< Domain > &domains, const Grid &grid, const ElementMapper &elementMapper, const CartMapper &cartMapper)
Writes the partition vector to a file in ResInsight compatible format and a partition file for each r...
Definition NlddReporting.hpp:164
Solver parameters for the BlackoilModel.
Definition BlackoilModelParameters.hpp:194
Definition BlackoilModel.hpp:114
Definition BlackoilModel.hpp:119
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34
Definition SimulatorReport.hpp:122