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>
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;
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>;
121 const ModelParameters&
param,
123 const bool terminal_output);
125 bool isParallel()
const
126 {
return grid_.comm().size() > 1; }
128 const EclipseState& eclState()
const
129 {
return simulator_.vanguard().eclState(); }
148 template <
class NonlinearSolverType>
151 NonlinearSolverType& nonlinear_solver);
153 template <
class NonlinearSolverType>
156 NonlinearSolverType& nonlinear_solver);
165 const int iterationIdx);
168 Scalar relativeChange()
const;
172 {
return simulator_.model().newtonMethod().linearSolver().iterations (); }
175 double& linearSolveSetupTime()
176 {
return linear_solve_setup_time_; }
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);
200 std::pair<Scalar,Scalar>
202 std::vector<Scalar>& maxCoeff,
203 std::vector<Scalar>& B_avg,
204 std::vector<int>& maxCoeffCell);
209 std::pair<std::vector<double>, std::vector<int>>
216 const double tol_cnv,
217 const double tol_cnv_energy,
218 const int iteration);
220 void updateTUNING(
const Tuning& tuning);
223 getReservoirConvergence(
const double reportTime,
227 std::vector<Scalar>& B_avg,
228 std::vector<Scalar>& residual_norms);
240 std::vector<Scalar>& residual_norms);
244 {
return Indices::numPhases; }
248 std::vector<std::vector<Scalar> >
253 std::vector<std::vector<Scalar> >
257 {
return simulator_; }
259 Simulator& simulator()
260 {
return simulator_; }
264 {
return failureReport_; }
275 const std::vector<StepReport>& stepReports()
const
276 {
return convergence_reports_; }
278 void writePartitions(
const std::filesystem::path& odir)
const;
281 BlackoilWellModel<TypeTag>&
283 {
return well_model_; }
287 {
return well_model_; }
289 void beginReportStep()
290 { simulator_.problem().beginEpisode(); }
293 { simulator_.problem().endEpisode(); }
295 template<
class Flu
idState,
class Res
idual>
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);
307 const ModelParameters&
param()
const
312 {
return compNames_; }
330 static constexpr bool has_micp_ = Indices::enableMICP;
332 ModelParameters param_;
343 std::vector<std::vector<Scalar>> residual_norms_history_;
344 Scalar current_relaxation_;
347 std::vector<StepReport> convergence_reports_;
348 ComponentName compNames_{};
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_;
364#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: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