24#ifndef OPM_BLACKOILMODEL_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_HEADER_INCLUDED
27#include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
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>
35#include <opm/simulators/linalg/ISTLSolver.hpp>
37#include <opm/simulators/timestepping/ConvergenceReport.hpp>
38#include <opm/simulators/timestepping/SimulatorReport.hpp>
39#include <opm/simulators/timestepping/SimulatorTimer.hpp>
41#include <opm/simulators/utils/ComponentName.hpp>
43#include <opm/simulators/wells/BlackoilWellModel.hpp>
49#include <fmt/format.h>
59template <
class TypeTag>
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;
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>;
115 std::pair<std::vector<double>, std::vector<int>> cnvPvSplit;
116 std::vector<unsigned> ixCells;
127 enum class DebugFlags {
143 const ModelParameters&
param,
145 const bool terminal_output);
147 bool isParallel()
const
148 {
return grid_.comm().size() > 1; }
150 const EclipseState& eclState()
const
151 {
return simulator_.vanguard().eclState(); }
155 SimulatorReportSingle
prepareStep(
const SimulatorTimerInterface& timer);
157 void initialLinearization(SimulatorReportSingle& report,
160 const SimulatorTimerInterface& timer);
169 template <
class NonlinearSolverType>
171 NonlinearSolverType& nonlinear_solver);
173 template <
class NonlinearSolverType>
174 SimulatorReportSingle nonlinearIterationNewton(
const SimulatorTimerInterface& timer,
175 NonlinearSolverType& nonlinear_solver);
181 Scalar relativeChange()
const;
185 {
return simulator_.model().newtonMethod().linearSolver().iterations (); }
188 double& linearSolveSetupTime()
189 {
return linear_solve_setup_time_; }
200 void storeSolutionUpdate(
const BVector& dx);
201 MaxSolutionUpdateData getMaxSolutionUpdate(
const std::vector<unsigned>& ixCells);
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);
218 std::pair<Scalar,Scalar>
220 std::vector<Scalar>& maxCoeff,
221 std::vector<Scalar>& B_avg,
222 std::vector<int>& maxCoeffCell);
234 const double tol_cnv,
235 const double tol_cnv_energy);
237 void updateTUNING(
const Tuning& tuning);
238 void updateTUNINGDP(
const TuningDp& tuning_dp);
241 getReservoirConvergence(
const double reportTime,
244 std::vector<Scalar>& B_avg,
245 std::vector<Scalar>& residual_norms);
255 std::vector<Scalar>& residual_norms);
259 {
return Indices::numPhases; }
263 std::vector<std::vector<Scalar> >
268 std::vector<std::vector<Scalar> >
272 {
return simulator_; }
274 Simulator& simulator()
275 {
return simulator_; }
279 {
return failureReport_; }
290 const std::vector<StepReport>& stepReports()
const
291 {
return convergence_reports_; }
296 convergence_reports_.back().report.pop_back();
297 residual_norms_history_.pop_back();
305 {
return well_model_; }
309 {
return well_model_; }
311 void beginReportStep()
312 { simulator_.problem().beginEpisode(); }
315 { simulator_.problem().endEpisode(); }
317 template<
class Flu
idState,
class Res
idual>
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);
329 const ModelParameters&
param()
const
334 {
return compNames_; }
352 static constexpr bool has_micp_ = Indices::enableMICP;
354 ModelParameters param_;
365 std::vector<std::vector<Scalar>> residual_norms_history_;
366 Scalar current_relaxation_;
369 SolutionVector solUpd_;
371 std::vector<StepReport> convergence_reports_;
372 ComponentName compNames_{};
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_;
388#include <opm/simulators/flow/BlackoilModel_impl.hpp>
Implementation of penalty cards for three-phase black oil.
Definition BlackoilModelConvergenceMonitor.hpp:35
BlackoilModel(Simulator &simulator, const ModelParameters ¶m, 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