#include <BlackoilModel.hpp>

Public Types

using ModelParameters = BlackoilModelParameters< TypeTag >
 
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 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)
 
bool isParallel () const
 
const EclipseState & eclState () const
 
SimulatorReportSingle prepareStep (const SimulatorTimerInterface &timer)
 
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)
 
template<class NonlinearSolverType >
SimulatorReportSingle nonlinearIterationNewton (const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
 
SimulatorReportSingle afterStep (const SimulatorTimerInterface &)
 
SimulatorReportSingle assembleReservoir (const SimulatorTimerInterface &, const int iterationIdx)
 Assemble the residual and Jacobian of the nonlinear system. More...
 
double relativeChange () const
 
int linearIterationsLastSolve () const
 Number of linear iterations used in last call to solveJacobianSystem(). More...
 
double & linearSolveSetupTime ()
 
void solveJacobianSystem (BVector &x)
 
void updateSolution (const BVector &dx)
 Apply an update to the primary variables. More...
 
bool terminalOutputEnabled () const
 Return true if output to cout is wanted. More...
 
std::tuple< double, double > convergenceReduction (Parallel::Communication comm, const double pvSumLocal, const double numAquiferPvSumLocal, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg)
 
std::pair< double, double > 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. More...
 
double computeCnvErrorPv (const std::vector< Scalar > &B_avg, double dt)
 Compute the total pore volume of cells violating CNV that are not part of a numerical aquifer. More...
 
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< double > &residual_norms)
 
int numPhases () const
 The number of active fluid phases in the model. More...
 
template<class T >
std::vector< std::vector< double > > computeFluidInPlace (const T &, const std::vector< int > &fipnum) const
 Wrapper required due to not following generic API. More...
 
std::vector< std::vector< double > > computeFluidInPlace (const std::vector< int > &) const
 Should not be called. More...
 
const Simulatorsimulator () const
 
Simulatorsimulator ()
 
const SimulatorReportSinglefailureReport () const
 return the statistics if the nonlinearIteration() method failed More...
 
SimulatorReportSingle localAccumulatedReports () const
 return the statistics if the nonlinearIteration() method failed More...
 
const std::vector< StepReport > & stepReports () const
 
void writePartitions (const std::filesystem::path &odir) const
 
const std::vector< std::vector< int > > & getConvCells () const
 
BlackoilWellModel< TypeTag > & wellModel ()
 return the StandardWells object More...
 
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 ModelParametersparam () const
 Returns const reference to model parameters. More...
 
const ComponentNamecompNames () const
 Returns const reference to component names. More...
 

Public Attributes

std::vector< bool > wasSwitched_
 

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 biofilmConcentrationIdx = Indices::biofilmConcentrationIdx
 
static constexpr int calciteConcentrationIdx = Indices::calciteConcentrationIdx
 

Protected Attributes

Simulatorsimulator_
 
const Gridgrid_
 
const PhaseUsage phaseUsage_
 
ModelParameters param_
 
SimulatorReportSingle failureReport_
 
BlackoilWellModel< TypeTag > & well_model_
 
RSTConv rst_conv_
 Helper class for RPTRST CONV. More...
 
bool terminal_output_
 Whether we print something to std::cout. More...
 
long int global_nc_
 The number of cells of the global grid. More...
 
std::vector< std::vector< double > > residual_norms_history_
 
double current_relaxation_
 
BVector dx_old_
 
std::vector< StepReportconvergence_reports_
 
ComponentName compNames_ {}
 
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
 Non-linear DD solver. More...
 

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_micp_ = getPropValue<TypeTag, Properties::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.

Member Typedef Documentation

◆ BVector

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::BVector = Dune::BlockVector<VectorBlockType>

◆ ComponentName

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::ComponentName = ::Opm::ComponentName<FluidSystem,Indices>

◆ ElementContext

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::ElementContext = GetPropType<TypeTag, Properties::ElementContext>

◆ FluidSystem

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>

◆ Grid

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::Grid = GetPropType<TypeTag, Properties::Grid>

◆ Indices

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::Indices = GetPropType<TypeTag, Properties::Indices>

◆ IntensiveQuantities

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>

◆ Mat

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::Mat = typename SparseMatrixAdapter::IstlMatrix

◆ MaterialLaw

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>

◆ MaterialLawParams

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>

◆ MatrixBlockType

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock

◆ ModelParameters

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::ModelParameters = BlackoilModelParameters<TypeTag>

◆ PrimaryVariables

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>

◆ Scalar

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::Scalar = GetPropType<TypeTag, Properties::Scalar>

◆ Simulator

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::Simulator = GetPropType<TypeTag, Properties::Simulator>

◆ SolutionVector

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>

◆ SparseMatrixAdapter

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>

◆ VectorBlockType

template<class TypeTag >
using Opm::BlackoilModel< TypeTag >::VectorBlockType = Dune::FieldVector<Scalar, numEq>

Constructor & Destructor Documentation

◆ BlackoilModel()

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

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
[in]paramparameters
[in]gridgrid data structure
[in]wellswell structure
[in]vfp_propertiesVertical flow performance tables
[in]linsolverlinear solver
[in]eclStateeclipse state
[in]terminal_outputrequest output to cout/cerr

References Opm::BlackoilModel< TypeTag >::convergence_reports_, Opm::detail::countGlobalCells(), Opm::BlackoilModel< TypeTag >::global_nc_, Opm::BlackoilModel< TypeTag >::grid_, Opm::BlackoilModelParameters< TypeTag >::matrix_add_well_contributions_, Opm::BlackoilModel< TypeTag >::nlddSolver_, Opm::BlackoilModelParameters< TypeTag >::nonlinear_solver_, and Opm::BlackoilModel< TypeTag >::param_.

Member Function Documentation

◆ afterStep()

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

Called once after each time step. In this class, this function does nothing.

Parameters
[in]timersimulation timer

References Opm::RSTConv::getData(), Opm::SimulatorReportSingle::pre_post_time, Opm::BlackoilModel< TypeTag >::rst_conv_, and Opm::BlackoilModel< TypeTag >::simulator_.

◆ assembleReservoir()

template<class TypeTag >
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::assembleReservoir ( const SimulatorTimerInterface ,
const int  iterationIdx 
)
inline

Assemble the residual and Jacobian of the nonlinear system.

References Opm::BlackoilModel< TypeTag >::simulator_, and Opm::BlackoilModel< TypeTag >::wellModel().

Referenced by Opm::BlackoilModel< TypeTag >::initialLinearization().

◆ beginReportStep()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::beginReportStep ( )
inline

◆ compNames()

template<class TypeTag >
const ComponentName & Opm::BlackoilModel< TypeTag >::compNames ( ) const
inline

Returns const reference to component names.

References Opm::BlackoilModel< TypeTag >::compNames_.

◆ computeCnvErrorPv()

◆ computeFluidInPlace() [1/2]

template<class TypeTag >
std::vector< std::vector< double > > Opm::BlackoilModel< TypeTag >::computeFluidInPlace ( const std::vector< int > &  ) const
inline

Should not be called.

References Opm::BlackoilModel< TypeTag >::computeFluidInPlace().

◆ computeFluidInPlace() [2/2]

template<class TypeTag >
template<class T >
std::vector< std::vector< double > > Opm::BlackoilModel< TypeTag >::computeFluidInPlace ( const T &  ,
const std::vector< int > &  fipnum 
) const
inline

Wrapper required due to not following generic API.

References Opm::BlackoilModel< TypeTag >::computeFluidInPlace().

Referenced by Opm::BlackoilModel< TypeTag >::computeFluidInPlace().

◆ convergenceReduction()

template<class TypeTag >
std::tuple< double, double > Opm::BlackoilModel< TypeTag >::convergenceReduction ( Parallel::Communication  comm,
const double  pvSumLocal,
const double  numAquiferPvSumLocal,
std::vector< Scalar > &  R_sum,
std::vector< Scalar > &  maxCoeff,
std::vector< Scalar > &  B_avg 
)
inline

◆ eclState()

template<class TypeTag >
const EclipseState & Opm::BlackoilModel< TypeTag >::eclState ( ) const
inline

◆ endReportStep()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::endReportStep ( )
inline

◆ failureReport()

template<class TypeTag >
const SimulatorReportSingle & Opm::BlackoilModel< TypeTag >::failureReport ( ) const
inline

return the statistics if the nonlinearIteration() method failed

References Opm::BlackoilModel< TypeTag >::failureReport_.

◆ getConvCells()

template<class TypeTag >
const std::vector< std::vector< int > > & Opm::BlackoilModel< TypeTag >::getConvCells ( ) const
inline

◆ getConvergence()

template<class TypeTag >
ConvergenceReport Opm::BlackoilModel< TypeTag >::getConvergence ( const SimulatorTimerInterface timer,
const int  iteration,
const int  maxIter,
std::vector< double > &  residual_norms 
)
inline

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

References Opm::SimulatorTimerInterface::currentStepLength(), Opm::BlackoilModel< TypeTag >::getConvergence(), Opm::BlackoilModel< TypeTag >::getReservoirConvergence(), Opm::BlackoilModel< TypeTag >::numEq, Opm::SimulatorTimerInterface::simulationTimeElapsed(), and Opm::BlackoilModel< TypeTag >::wellModel().

Referenced by Opm::BlackoilModel< TypeTag >::getConvergence(), and Opm::BlackoilModel< TypeTag >::initialLinearization().

◆ getMaxCoeff()

◆ getReservoirConvergence()

◆ initialLinearization()

◆ isParallel()

template<class TypeTag >
bool Opm::BlackoilModel< TypeTag >::isParallel ( ) const
inline

◆ linearIterationsLastSolve()

template<class TypeTag >
int Opm::BlackoilModel< TypeTag >::linearIterationsLastSolve ( ) const
inline

Number of linear iterations used in last call to solveJacobianSystem().

References Opm::BlackoilModel< TypeTag >::simulator_.

Referenced by Opm::BlackoilModel< TypeTag >::nonlinearIterationNewton().

◆ linearSolveSetupTime()

template<class TypeTag >
double & Opm::BlackoilModel< TypeTag >::linearSolveSetupTime ( )
inline

◆ localAccumulatedReports()

template<class TypeTag >
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::localAccumulatedReports ( ) const
inline

return the statistics if the nonlinearIteration() method failed

References Opm::BlackoilModel< TypeTag >::nlddSolver_.

◆ localConvergenceData()

template<class TypeTag >
std::pair< double, double > Opm::BlackoilModel< TypeTag >::localConvergenceData ( std::vector< Scalar > &  R_sum,
std::vector< Scalar > &  maxCoeff,
std::vector< Scalar > &  B_avg,
std::vector< int > &  maxCoeffCell 
)
inline

◆ nonlinearIteration()

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

◆ nonlinearIterationNewton()

◆ numPhases()

template<class TypeTag >
int Opm::BlackoilModel< TypeTag >::numPhases ( ) const
inline

The number of active fluid phases in the model.

References Opm::PhaseUsage::num_phases, and Opm::BlackoilModel< TypeTag >::phaseUsage_.

◆ param()

template<class TypeTag >
const ModelParameters & Opm::BlackoilModel< TypeTag >::param ( ) const
inline

Returns const reference to model parameters.

References Opm::BlackoilModel< TypeTag >::param_.

◆ prepareStep()

◆ relativeChange()

template<class TypeTag >
double Opm::BlackoilModel< TypeTag >::relativeChange ( ) const
inline

◆ simulator() [1/2]

template<class TypeTag >
Simulator & Opm::BlackoilModel< TypeTag >::simulator ( )
inline

◆ simulator() [2/2]

◆ solveJacobianSystem()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::solveJacobianSystem ( BVector x)
inline

◆ stepReports()

template<class TypeTag >
const std::vector< StepReport > & Opm::BlackoilModel< TypeTag >::stepReports ( ) const
inline

◆ terminalOutputEnabled()

template<class TypeTag >
bool Opm::BlackoilModel< TypeTag >::terminalOutputEnabled ( ) const
inline

Return true if output to cout is wanted.

References Opm::BlackoilModel< TypeTag >::terminal_output_.

Referenced by Opm::BlackoilModel< TypeTag >::nonlinearIterationNewton().

◆ updateSolution()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::updateSolution ( const BVector dx)
inline

◆ updateTUNING()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::updateTUNING ( const Tuning &  tuning)
inline

◆ wellModel() [1/2]

◆ wellModel() [2/2]

template<class TypeTag >
const BlackoilWellModel< TypeTag > & Opm::BlackoilModel< TypeTag >::wellModel ( ) const
inline

◆ writePartitions()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::writePartitions ( const std::filesystem::path &  odir) const
inline

Member Data Documentation

◆ biofilmConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::biofilmConcentrationIdx = Indices::biofilmConcentrationIdx
staticconstexpr

◆ calciteConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::calciteConcentrationIdx = Indices::calciteConcentrationIdx
staticconstexpr

◆ compNames_

template<class TypeTag >
ComponentName Opm::BlackoilModel< TypeTag >::compNames_ {}
protected

◆ contiBiofilmEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiBiofilmEqIdx = Indices::contiBiofilmEqIdx
staticconstexpr

◆ contiBrineEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiBrineEqIdx = Indices::contiBrineEqIdx
staticconstexpr

◆ contiCalciteEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiCalciteEqIdx = Indices::contiCalciteEqIdx
staticconstexpr

◆ contiEnergyEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiEnergyEqIdx = Indices::contiEnergyEqIdx
staticconstexpr

◆ contiFoamEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiFoamEqIdx = Indices::contiFoamEqIdx
staticconstexpr

◆ contiMicrobialEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiMicrobialEqIdx = Indices::contiMicrobialEqIdx
staticconstexpr

◆ contiOxygenEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiOxygenEqIdx = Indices::contiOxygenEqIdx
staticconstexpr

◆ contiPolymerEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiPolymerEqIdx = Indices::contiPolymerEqIdx
staticconstexpr

◆ contiPolymerMWEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx
staticconstexpr

◆ contiSolventEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiSolventEqIdx = Indices::contiSolventEqIdx
staticconstexpr

◆ contiUreaEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiUreaEqIdx = Indices::contiUreaEqIdx
staticconstexpr

◆ contiZfracEqIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::contiZfracEqIdx = Indices::contiZfracEqIdx
staticconstexpr

◆ convergence_reports_

◆ current_relaxation_

template<class TypeTag >
double Opm::BlackoilModel< TypeTag >::current_relaxation_
protected

◆ dx_old_

◆ failureReport_

◆ foamConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::foamConcentrationIdx = Indices::foamConcentrationIdx
staticconstexpr

◆ global_nc_

template<class TypeTag >
long int Opm::BlackoilModel< TypeTag >::global_nc_
protected

◆ grid_

◆ has_brine_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>()
staticconstexprprotected

◆ has_energy_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_energy_ = getPropValue<TypeTag, Properties::EnableEnergy>()
staticconstexprprotected

◆ has_extbo_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>()
staticconstexprprotected

◆ has_foam_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>()
staticconstexprprotected

◆ has_micp_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_micp_ = getPropValue<TypeTag, Properties::EnableMICP>()
staticconstexprprotected

◆ has_polymer_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>()
staticconstexprprotected

◆ has_polymermw_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>()
staticconstexprprotected

◆ has_solvent_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>()
staticconstexprprotected

◆ microbialConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::microbialConcentrationIdx = Indices::microbialConcentrationIdx
staticconstexpr

◆ nlddSolver_

◆ numEq

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::numEq = Indices::numEq
staticconstexpr

◆ oxygenConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::oxygenConcentrationIdx = Indices::oxygenConcentrationIdx
staticconstexpr

◆ param_

◆ phaseUsage_

template<class TypeTag >
const PhaseUsage Opm::BlackoilModel< TypeTag >::phaseUsage_
protected

◆ polymerConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::polymerConcentrationIdx = Indices::polymerConcentrationIdx
staticconstexpr

◆ polymerMoleWeightIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::polymerMoleWeightIdx = Indices::polymerMoleWeightIdx
staticconstexpr

◆ residual_norms_history_

template<class TypeTag >
std::vector<std::vector<double> > Opm::BlackoilModel< TypeTag >::residual_norms_history_
protected

◆ rst_conv_

◆ saltConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::saltConcentrationIdx = Indices::saltConcentrationIdx
staticconstexpr

◆ simulator_

◆ solventSaturationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::solventSaturationIdx = Indices::solventSaturationIdx
staticconstexpr

◆ temperatureIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::temperatureIdx = Indices::temperatureIdx
staticconstexpr

◆ terminal_output_

◆ ureaConcentrationIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::ureaConcentrationIdx = Indices::ureaConcentrationIdx
staticconstexpr

◆ wasSwitched_

template<class TypeTag >
std::vector<bool> Opm::BlackoilModel< TypeTag >::wasSwitched_

◆ well_model_

template<class TypeTag >
BlackoilWellModel<TypeTag>& Opm::BlackoilModel< TypeTag >::well_model_
protected

◆ zFractionIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::zFractionIdx = Indices::zFractionIdx
staticconstexpr

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