opm-simulators
Loading...
Searching...
No Matches
Opm::BlackoilModel< TypeTag > Class Template Reference

A model implementation for three-phase black oil. More...

#include <BlackoilModel.hpp>

Public Types

using Simulator = GetPropType<TypeTag, Properties::Simulator>
using Grid = GetPropType<TypeTag, Properties::Grid>
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>
using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>
using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>
using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>
using Indices = GetPropType<TypeTag, Properties::Indices>
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>
using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>
using Scalar = GetPropType<TypeTag, Properties::Scalar>
using ModelParameters = BlackoilModelParameters<Scalar>
using VectorBlockType = Dune::FieldVector<Scalar, numEq>
using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock
using Mat = typename SparseMatrixAdapter::IstlMatrix
using BVector = Dune::BlockVector<VectorBlockType>
using ComponentName = ::Opm::ComponentName<FluidSystem,Indices>

Public Member Functions

 BlackoilModel (Simulator &simulator, const ModelParameters &param, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
 Construct the model.
bool isParallel () const
const EclipseState & eclState () const
SimulatorReportSingle prepareStep (const SimulatorTimerInterface &timer)
 Called once before each time step.
void initialLinearization (SimulatorReportSingle &report, const int iteration, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
template<class NonlinearSolverType>
SimulatorReportSingle nonlinearIteration (const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
 Called once per nonlinear iteration.
template<class NonlinearSolverType>
SimulatorReportSingle nonlinearIterationNewton (const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
SimulatorReportSingle afterStep (const SimulatorTimerInterface &)
 Called once after each time step.
SimulatorReportSingle assembleReservoir (const SimulatorTimerInterface &, const int iterationIdx)
 Assemble the residual and Jacobian of the nonlinear system.
Scalar relativeChange () const
int linearIterationsLastSolve () const
 Number of linear iterations used in last call to solveJacobianSystem().
double & linearSolveSetupTime ()
void solveJacobianSystem (BVector &x)
 Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.
void updateSolution (const BVector &dx)
 Apply an update to the primary variables.
bool terminalOutputEnabled () const
 Return true if output to cout is wanted.
std::tuple< Scalar, Scalar > convergenceReduction (Parallel::Communication comm, const Scalar pvSumLocal, const Scalar numAquiferPvSumLocal, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg)
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.
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 based on CNV point measures.
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 convergence criteria at the last time step.
void updateTUNING (const Tuning &tuning)
ConvergenceReport getReservoirConvergence (const double reportTime, const double dt, const int iteration, const int maxIter, std::vector< Scalar > &B_avg, std::vector< Scalar > &residual_norms)
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).
int numPhases () const
 The number of active fluid phases in the model.
template<class T>
std::vector< std::vector< Scalar > > computeFluidInPlace (const T &, const std::vector< int > &fipnum) const
 Wrapper required due to not following generic API.
std::vector< std::vector< Scalar > > computeFluidInPlace (const std::vector< int > &) const
 Should not be called.
const Simulator & simulator () const
Simulator & simulator ()
const SimulatorReportSinglefailureReport () const
 return the statistics if the nonlinearIteration() method failed
const SimulatorReportlocalAccumulatedReports () const
 return the statistics of local solves accumulated for this rank
const std::vector< SimulatorReport > & domainAccumulatedReports () const
 return the statistics of local solves accumulated for each domain on this rank
void writeNonlinearIterationsPerCell (const std::filesystem::path &odir) const
 Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.
const std::vector< StepReport > & stepReports () const
void writePartitions (const std::filesystem::path &odir) const
BlackoilWellModel< TypeTag > & wellModel ()
 return the StandardWells object
const BlackoilWellModel< TypeTag > & wellModel () const
void beginReportStep ()
void endReportStep ()
template<class FluidState, class Residual>
void getMaxCoeff (const unsigned cell_idx, const IntensiveQuantities &intQuants, const FluidState &fs, const Residual &modelResid, const Scalar pvValue, std::vector< Scalar > &B_avg, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< int > &maxCoeffCell)
const ModelParameters & param () const
 Returns const reference to model parameters.
const ComponentName & compNames () const
 Returns const reference to component names.
bool hasNlddSolver () const
 Returns true if an NLDD solver exists.

Static Public Attributes

static constexpr int numEq = Indices::numEq
static constexpr int contiSolventEqIdx = Indices::contiSolventEqIdx
static constexpr int contiZfracEqIdx = Indices::contiZfracEqIdx
static constexpr int contiPolymerEqIdx = Indices::contiPolymerEqIdx
static constexpr int contiEnergyEqIdx = Indices::contiEnergyEqIdx
static constexpr int contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx
static constexpr int contiFoamEqIdx = Indices::contiFoamEqIdx
static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx
static constexpr int contiMicrobialEqIdx = Indices::contiMicrobialEqIdx
static constexpr int contiOxygenEqIdx = Indices::contiOxygenEqIdx
static constexpr int contiUreaEqIdx = Indices::contiUreaEqIdx
static constexpr int contiBiofilmEqIdx = Indices::contiBiofilmEqIdx
static constexpr int contiCalciteEqIdx = Indices::contiCalciteEqIdx
static constexpr int solventSaturationIdx = Indices::solventSaturationIdx
static constexpr int zFractionIdx = Indices::zFractionIdx
static constexpr int polymerConcentrationIdx = Indices::polymerConcentrationIdx
static constexpr int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx
static constexpr int temperatureIdx = Indices::temperatureIdx
static constexpr int foamConcentrationIdx = Indices::foamConcentrationIdx
static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx
static constexpr int microbialConcentrationIdx = Indices::microbialConcentrationIdx
static constexpr int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx
static constexpr int ureaConcentrationIdx = Indices::ureaConcentrationIdx
static constexpr int biofilmVolumeFractionIdx = Indices::biofilmVolumeFractionIdx
static constexpr int calciteVolumeFractionIdx = Indices::calciteVolumeFractionIdx

Protected Attributes

Simulator & simulator_
const Grid & grid_
ModelParameters param_
SimulatorReportSingle failureReport_
BlackoilWellModel< TypeTag > & well_model_
bool terminal_output_
 Whether we print something to std::cout.
long int global_nc_
 The number of cells of the global grid.
std::vector< std::vector< Scalar > > residual_norms_history_
Scalar current_relaxation_
BVector dx_old_
std::vector< StepReportconvergence_reports_
ComponentName compNames_ {}
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
 Non-linear DD solver.
BlackoilModelConvergenceMonitor< Scalar > conv_monitor_

Static Protected Attributes

static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>()
static constexpr bool has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>()
static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>()
static constexpr bool has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>()
static constexpr bool has_energy_ = getPropValue<TypeTag, Properties::EnableEnergy>()
static constexpr bool has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>()
static constexpr bool has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>()
static constexpr bool has_bioeffects_ = getPropValue<TypeTag, Properties::EnableBioeffects>()
static constexpr bool has_micp_ = Indices::enableMICP

Detailed Description

template<class TypeTag>
class Opm::BlackoilModel< TypeTag >

A model implementation for three-phase black oil.

The simulator is capable of handling three-phase problems where gas can be dissolved in oil and vice versa. It uses an industry-standard TPFA discretization with per-phase upwind weighting of mobilities.

Constructor & Destructor Documentation

◆ BlackoilModel()

template<class TypeTag>
Opm::BlackoilModel< TypeTag >::BlackoilModel ( Simulator & simulator,
const ModelParameters & param,
BlackoilWellModel< TypeTag > & well_model,
const bool terminal_output )

Construct the model.

It will retain references to the arguments of this functions, and they are expected to remain in scope for the lifetime of the solver.

Parameters
simulatorReference to main simulator
[in]paramparameters
[in]well_modelReference to well model
[in]terminal_outputrequest output to cout/cerr

Member Function Documentation

◆ afterStep()

template<class TypeTag>
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::afterStep ( const SimulatorTimerInterface & )

Called once after each time step.

In this class, this function does nothing.

Parameters
[in]timersimulation timer

◆ getConvergence()

template<class TypeTag>
ConvergenceReport Opm::BlackoilModel< TypeTag >::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).

Parameters
[in]timersimulation timer
[in]iterationcurrent iteration number
[in]maxItermaximum number of iterations
[out]residual_normsCNV residuals by phase

◆ localConvergenceData()

template<class TypeTag>
std::pair< typename BlackoilModel< TypeTag >::Scalar, typename BlackoilModel< TypeTag >::Scalar > Opm::BlackoilModel< TypeTag >::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.

Returns
A pair of the local pore volume of interior cells and the pore volumes of the cells associated with a numerical aquifer.

◆ nonlinearIteration()

template<class TypeTag>
template<class NonlinearSolverType>
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::nonlinearIteration ( const int iteration,
const SimulatorTimerInterface & timer,
NonlinearSolverType & nonlinear_solver )

Called once per nonlinear iteration.

This model will perform a Newton-Raphson update, changing reservoir_state and well_state. It will also use the nonlinear_solver to do relaxation of updates if necessary.

Parameters
[in]iterationshould be 0 for the first call of a new timestep
[in]timersimulation timer
[in]nonlinear_solvernonlinear solver used (for oscillation/relaxation control)

◆ prepareStep()

template<class TypeTag>
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::prepareStep ( const SimulatorTimerInterface & timer)

Called once before each time step.

Parameters
[in]timersimulation timer

The documentation for this class was generated from the following files: