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 int numEq = Indices::numEq;
79 static constexpr int contiSolventEqIdx = Indices::contiSolventEqIdx;
80 static constexpr int contiZfracEqIdx = Indices::contiZfracEqIdx;
81 static constexpr int contiPolymerEqIdx = Indices::contiPolymerEqIdx;
82 static constexpr int contiEnergyEqIdx = Indices::contiEnergyEqIdx;
83 static constexpr int contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx;
84 static constexpr int contiFoamEqIdx = Indices::contiFoamEqIdx;
85 static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx;
86 static constexpr int contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
87 static constexpr int contiOxygenEqIdx = Indices::contiOxygenEqIdx;
88 static constexpr int contiUreaEqIdx = Indices::contiUreaEqIdx;
89 static constexpr int contiBiofilmEqIdx = Indices::contiBiofilmEqIdx;
90 static constexpr int contiCalciteEqIdx = Indices::contiCalciteEqIdx;
91 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
92 static constexpr int zFractionIdx = Indices::zFractionIdx;
93 static constexpr int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
94 static constexpr int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
95 static constexpr int temperatureIdx = Indices::temperatureIdx;
96 static constexpr int foamConcentrationIdx = Indices::foamConcentrationIdx;
97 static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx;
98 static constexpr int microbialConcentrationIdx = Indices::microbialConcentrationIdx;
99 static constexpr int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
100 static constexpr int ureaConcentrationIdx = Indices::ureaConcentrationIdx;
101 static constexpr int biofilmVolumeFractionIdx = Indices::biofilmVolumeFractionIdx;
102 static constexpr int calciteVolumeFractionIdx = Indices::calciteVolumeFractionIdx;
103
104 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
105 using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock;
106 using Mat = typename SparseMatrixAdapter::IstlMatrix;
107 using BVector = Dune::BlockVector<VectorBlockType>;
108
109 using ComponentName = ::Opm::ComponentName<FluidSystem,Indices>;
110
111 // --------- Public methods ---------
112
120 BlackoilModel(Simulator& simulator,
121 const ModelParameters& param,
122 BlackoilWellModel<TypeTag>& well_model,
123 const bool terminal_output);
124
125 bool isParallel() const
126 { return grid_.comm().size() > 1; }
127
128 const EclipseState& eclState() const
129 { return simulator_.vanguard().eclState(); }
130
134
135 void initialLinearization(SimulatorReportSingle& report,
136 const int iteration,
137 const int minIter,
138 const int maxIter,
139 const SimulatorTimerInterface& timer);
140
148 template <class NonlinearSolverType>
149 SimulatorReportSingle nonlinearIteration(const int iteration,
150 const SimulatorTimerInterface& timer,
151 NonlinearSolverType& nonlinear_solver);
152
153 template <class NonlinearSolverType>
154 SimulatorReportSingle nonlinearIterationNewton(const int iteration,
155 const SimulatorTimerInterface& timer,
156 NonlinearSolverType& nonlinear_solver);
157
162
165 const int iterationIdx);
166
167 // compute the "relative" change of the solution between time steps
168 Scalar relativeChange() const;
169
172 { return simulator_.model().newtonMethod().linearSolver().iterations (); }
173
174 // Obtain reference to linear solver setup time
175 double& linearSolveSetupTime()
176 { return linear_solve_setup_time_; }
177
180 void solveJacobianSystem(BVector& x);
181
183 void updateSolution(const BVector& dx);
184
187 { return terminal_output_; }
188
189 std::tuple<Scalar,Scalar>
190 convergenceReduction(Parallel::Communication comm,
191 const Scalar pvSumLocal,
192 const Scalar numAquiferPvSumLocal,
193 std::vector<Scalar>& R_sum,
194 std::vector<Scalar>& maxCoeff,
195 std::vector<Scalar>& B_avg);
196
200 std::pair<Scalar,Scalar>
201 localConvergenceData(std::vector<Scalar>& R_sum,
202 std::vector<Scalar>& maxCoeff,
203 std::vector<Scalar>& B_avg,
204 std::vector<int>& maxCoeffCell);
205
209 std::pair<std::vector<double>, std::vector<int>>
210 characteriseCnvPvSplit(const std::vector<Scalar>& B_avg, const double dt);
211
214 void convergencePerCell(const std::vector<Scalar>& B_avg,
215 const double dt,
216 const double tol_cnv,
217 const double tol_cnv_energy,
218 const int iteration);
219
220 void updateTUNING(const Tuning& tuning);
221
223 getReservoirConvergence(const double reportTime,
224 const double dt,
225 const int iteration,
226 const int maxIter,
227 std::vector<Scalar>& B_avg,
228 std::vector<Scalar>& residual_norms);
229
238 const int iteration,
239 const int maxIter,
240 std::vector<Scalar>& residual_norms);
241
243 int numPhases() const
244 { return Indices::numPhases; }
245
247 template<class T>
248 std::vector<std::vector<Scalar> >
249 computeFluidInPlace(const T&, const std::vector<int>& fipnum) const
250 { return this->computeFluidInPlace(fipnum); }
251
253 std::vector<std::vector<Scalar> >
254 computeFluidInPlace(const std::vector<int>& /*fipnum*/) const;
255
256 const Simulator& simulator() const
257 { return simulator_; }
258
259 Simulator& simulator()
260 { return simulator_; }
261
264 { return failureReport_; }
265
268
270 const std::vector<SimulatorReport>& domainAccumulatedReports() const;
271
273 void writeNonlinearIterationsPerCell(const std::filesystem::path& odir) const;
274
275 const std::vector<StepReport>& stepReports() const
276 { return convergence_reports_; }
277
278 void writePartitions(const std::filesystem::path& odir) const;
279
281 BlackoilWellModel<TypeTag>&
283 { return well_model_; }
284
286 wellModel() const
287 { return well_model_; }
288
289 void beginReportStep()
290 { simulator_.problem().beginEpisode(); }
291
292 void endReportStep()
293 { simulator_.problem().endEpisode(); }
294
295 template<class FluidState, class Residual>
296 void getMaxCoeff(const unsigned cell_idx,
297 const IntensiveQuantities& intQuants,
298 const FluidState& fs,
299 const Residual& modelResid,
300 const Scalar pvValue,
301 std::vector<Scalar>& B_avg,
302 std::vector<Scalar>& R_sum,
303 std::vector<Scalar>& maxCoeff,
304 std::vector<int>& maxCoeffCell);
305
307 const ModelParameters& param() const
308 { return param_; }
309
311 const ComponentName& compNames() const
312 { return compNames_; }
313
315 bool hasNlddSolver() const
316 { return nlddSolver_ != nullptr; }
317
318protected:
319 // --------- Data members ---------
320 Simulator& simulator_;
321 const Grid& grid_;
322 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
323 static constexpr bool has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>();
324 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
325 static constexpr bool has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>();
326 static constexpr bool has_energy_ = getPropValue<TypeTag, Properties::EnableEnergy>();
327 static constexpr bool has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>();
328 static constexpr bool has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>();
329 static constexpr bool has_bioeffects_ = getPropValue<TypeTag, Properties::EnableBioeffects>();
330 static constexpr bool has_micp_ = Indices::enableMICP;
331
332 ModelParameters param_;
333 SimulatorReportSingle failureReport_;
334
335 // Well Model
336 BlackoilWellModel<TypeTag>& well_model_;
337
341 long int global_nc_;
342
343 std::vector<std::vector<Scalar>> residual_norms_history_;
344 Scalar current_relaxation_;
345 BVector dx_old_;
346
347 std::vector<StepReport> convergence_reports_;
348 ComponentName compNames_{};
349
350 std::unique_ptr<BlackoilModelNldd<TypeTag>> nlddSolver_;
352
353private:
354 Scalar dpMaxRel() const { return param_.dp_max_rel_; }
355 Scalar dsMax() const { return param_.ds_max_; }
356 Scalar drMaxRel() const { return param_.dr_max_rel_; }
357 Scalar maxResidualAllowed() const { return param_.max_residual_allowed_; }
358 double linear_solve_setup_time_;
359 std::vector<bool> wasSwitched_;
360};
361
362} // namespace Opm
363
364#include <opm/simulators/flow/BlackoilModel_impl.hpp>
365
366#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:51
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:249
int linearIterationsLastSolve() const
Number of linear iterations used in last call to solveJacobianSystem().
Definition BlackoilModel.hpp:171
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Called once per nonlinear iteration.
Definition BlackoilModel_impl.hpp:214
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:1046
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &, const int iterationIdx)
Assemble the residual and Jacobian of the nonlinear system.
Definition BlackoilModel_impl.hpp:357
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int iteration, 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:988
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition BlackoilModel_impl.hpp:1036
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:604
bool terminalOutputEnabled() const
Return true if output to cout is wanted.
Definition BlackoilModel.hpp:186
int numPhases() const
The number of active fluid phases in the model.
Definition BlackoilModel.hpp:243
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition BlackoilModel_impl.hpp:1025
SimulatorReportSingle afterStep(const SimulatorTimerInterface &)
Called once after each time step.
Definition BlackoilModel_impl.hpp:344
const ComponentName & compNames() const
Returns const reference to component names.
Definition BlackoilModel.hpp:311
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition BlackoilModel_impl.hpp:520
long int global_nc_
The number of cells of the global grid.
Definition BlackoilModel.hpp:341
const ModelParameters & param() const
Returns const reference to model parameters.
Definition BlackoilModel.hpp:307
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition BlackoilModel.hpp:350
bool hasNlddSolver() const
Returns true if an NLDD solver exists.
Definition BlackoilModel.hpp:315
std::pair< std::vector< double >, std::vector< int > > 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:655
const SimulatorReportSingle & failureReport() const
return the statistics if the nonlinearIteration() method failed
Definition BlackoilModel.hpp:263
BlackoilWellModel< TypeTag > & wellModel()
return the StandardWells object
Definition BlackoilModel.hpp:282
void solveJacobianSystem(BVector &x)
Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.
Definition BlackoilModel_impl.hpp:458
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Called once before each time step.
Definition BlackoilModel_impl.hpp:86
void convergencePerCell(const std::vector< Scalar > &B_avg, const double dt, const double tol_cnv, const double tol_cnv_energy, const int iteration)
Compute the number of Newtons required by each cell in order to satisfy the solution change convergen...
Definition BlackoilModel_impl.hpp:941
bool terminal_output_
Whether we print something to std::cout.
Definition BlackoilModel.hpp:339
Class for handling the blackoil well model.
Definition BlackoilWellModel.hpp:102
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:43
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
Solver parameters for the BlackoilModel.
Definition BlackoilModelParameters.hpp:180
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34
Definition SimulatorReport.hpp:122