Black oil model for coupling Flow simulations with TPSA geomechanics. More...

#include <BlackoilModelTPSA.hpp>

Inheritance diagram for Opm::BlackoilModelTPSA< TypeTag >:
Inheritance graph

Public Types

using ModelParameters = typename ParentType::ModelParameters
 
enum class  DebugFlags { STRICT = 0 , RELAXED = 1 , TUNINGDP = 2 }
 
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 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

 BlackoilModelTPSA (Simulator &simulator, const ModelParameters &param, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
 Constructor. More...
 
template<class NonlinearSolverType >
SimulatorReportSingle nonlinearIteration (const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
 Perform a nonlinear iteration updating Flow and TPSA geomechanics. More...
 
template<class NonlinearSolverType >
SimulatorReportSingle nonlinearIterationFixedStressTPSA (const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
 Perform a nonlinear iteration updating Flow and TPSA geomechanics in a fixed-stress, iterative loop. More...
 
template<class NonlinearSolverType >
SimulatorReportSingle nonlinearIterationLaggedTPSA (const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
 Perform a nonlinear iteration updating Flow and TPSA geomechanics in a lagged scheme. More...
 
bool solveTpsaEquations ()
 Solve TPSA geomechanics equations. More...
 
bool isParallel () const
 
const EclipseState & eclState () const
 
SimulatorReportSingle prepareStep (const SimulatorTimerInterface &timer)
 
void initialLinearization (SimulatorReportSingle &report, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
 
template<class NonlinearSolverType >
SimulatorReportSingle nonlinearIterationNewton (const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
 
SimulatorReportSingle assembleReservoir (const SimulatorTimerInterface &timer)
 Assemble the residual and Jacobian of the nonlinear system. More...
 
Scalar 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...
 
void prepareStoringSolutionUpdate ()
 Get solution update vector as a PrimaryVariable. More...
 
void storeSolutionUpdate (const BVector &dx)
 
MaxSolutionUpdateData getMaxSolutionUpdate (const std::vector< unsigned > &ixCells)
 
bool terminalOutputEnabled () const
 Return true if output to cout is wanted. More...
 
std::tuple< Scalar, ScalarconvergenceReduction (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, ScalarlocalConvergenceData (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...
 
CnvPvSplitData 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. Also returns list of cells where CNV is greater than its strict tolerance. More...
 
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 convergence criteria at the last time step. More...
 
void updateTUNING (const Tuning &tuning)
 
void updateTUNINGDP (const TuningDp &tuning_dp)
 
ConvergenceReport getReservoirConvergence (const double reportTime, const double dt, const int maxIter, std::vector< Scalar > &B_avg, std::vector< Scalar > &residual_norms)
 
ConvergenceReport getConvergence (const SimulatorTimerInterface &timer, const int maxIter, std::vector< Scalar > &residual_norms)
 
int numPhases () const
 The number of active fluid phases in the model. More...
 
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. More...
 
std::vector< std::vector< Scalar > > 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...
 
const SimulatorReportlocalAccumulatedReports () const
 return the statistics of local solves accumulated for this rank More...
 
const std::vector< SimulatorReport > & domainAccumulatedReports () const
 return the statistics of local solves accumulated for each domain on this rank More...
 
void writeNonlinearIterationsPerCell (const std::filesystem::path &odir) const
 Write the number of nonlinear iterations per cell to a file in ResInsight compatible format. More...
 
const std::vector< StepReport > & stepReports () const
 
void popLastConvergenceReport ()
 Remove the last convergence report entry and residual norms history entry. More...
 
void writePartitions (const std::filesystem::path &odir) 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...
 
bool hasNlddSolver () const
 Returns true if an NLDD solver exists. More...
 

Static Public Attributes

static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>()
 
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

Simulatorsimulator_
 
const Gridgrid_
 
ModelParameters param_
 
SimulatorReportSingle failureReport_
 
BlackoilWellModel< TypeTag > & well_model_
 
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< Scalar > > residual_norms_history_
 
Scalar current_relaxation_
 
BVector dx_old_
 
SolutionVector solUpd_
 
std::vector< StepReportconvergence_reports_
 
ComponentName compNames_ {}
 
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
 Non-linear DD solver. More...
 
BlackoilModelConvergenceMonitor< Scalarconv_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::EnergyModuleType>() == EnergyModules::FullyImplicitThermal
 
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::BlackoilModelTPSA< TypeTag >

Black oil model for coupling Flow simulations with TPSA geomechanics.

Member Typedef Documentation

◆ BVector

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

◆ ComponentName

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

◆ ElementContext

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

◆ FluidSystem

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

◆ Grid

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

◆ Indices

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

◆ IntensiveQuantities

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

◆ Mat

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

◆ MaterialLaw

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

◆ MaterialLawParams

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

◆ MatrixBlockType

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

◆ ModelParameters

template<class TypeTag >
using Opm::BlackoilModelTPSA< TypeTag >::ModelParameters = typename ParentType::ModelParameters

◆ PrimaryVariables

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

◆ SolutionVector

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

◆ SparseMatrixAdapter

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

◆ VectorBlockType

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

Member Enumeration Documentation

◆ DebugFlags

template<class TypeTag >
enum class Opm::BlackoilModel::DebugFlags
stronginherited
Enumerator
STRICT 
RELAXED 
TUNINGDP 

Constructor & Destructor Documentation

◆ BlackoilModelTPSA()

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

Constructor.

Parameters
simulatorReference to simulator object
paramReference to parameters for model
well_modelRefenerence to well model
terminal_outputBool for terminal output

Member Function Documentation

◆ assembleReservoir()

template<class TypeTag >
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::assembleReservoir ( const SimulatorTimerInterface timer)
inherited

Assemble the residual and Jacobian of the nonlinear system.

◆ beginReportStep()

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

◆ characteriseCnvPvSplit()

template<class TypeTag >
BlackoilModel< TypeTag >::CnvPvSplitData Opm::BlackoilModel< TypeTag >::characteriseCnvPvSplit ( const std::vector< Scalar > &  B_avg,
const double  dt 
)
inherited

Compute pore-volume/cell count split among "converged", "relaxed converged", "unconverged" cells based on CNV point measures. Also returns list of cells where CNV is greater than its strict tolerance.

References OPM_BEGIN_PARALLEL_TRY_CATCH, and OPM_END_PARALLEL_TRY_CATCH.

◆ compNames()

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

Returns const reference to component names.

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

◆ computeFluidInPlace() [1/2]

template<class TypeTag >
std::vector< std::vector< typename BlackoilModel< TypeTag >::Scalar > > Opm::BlackoilModel< TypeTag >::computeFluidInPlace ( const std::vector< int > &  ) const
inherited

Should not be called.

◆ computeFluidInPlace() [2/2]

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

Wrapper required due to not following generic API.

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

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

◆ convergencePerCell()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::convergencePerCell ( const std::vector< Scalar > &  B_avg,
const double  dt,
const double  tol_cnv,
const double  tol_cnv_energy 
)
inherited

Compute the number of Newtons required by each cell in order to satisfy the solution change convergence criteria at the last time step.

References OPM_BEGIN_PARALLEL_TRY_CATCH, and OPM_END_PARALLEL_TRY_CATCH.

◆ convergenceReduction()

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

◆ domainAccumulatedReports()

template<class TypeTag >
const std::vector< SimulatorReport > & Opm::BlackoilModel< TypeTag >::domainAccumulatedReports
inherited

return the statistics of local solves accumulated for each domain on this rank

◆ eclState()

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

◆ endReportStep()

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

◆ failureReport()

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

return the statistics if the nonlinearIteration() method failed

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

◆ getConvergence()

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

Compute convergence based on total mass balance (tol_mb) and maximum residual mass balance (tol_cnv).

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

References Opm::SimulatorTimerInterface::currentStepLength(), and Opm::SimulatorTimerInterface::simulationTimeElapsed().

◆ getMaxCoeff()

template<class TypeTag >
template<class FluidState , class Residual >
void Opm::BlackoilModel< TypeTag >::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 
)
inherited

◆ getMaxSolutionUpdate()

template<class TypeTag >
BlackoilModel< TypeTag >::MaxSolutionUpdateData Opm::BlackoilModel< TypeTag >::getMaxSolutionUpdate ( const std::vector< unsigned > &  ixCells)
inherited

◆ getReservoirConvergence()

template<class TypeTag >
ConvergenceReport Opm::BlackoilModel< TypeTag >::getReservoirConvergence ( const double  reportTime,
const double  dt,
const int  maxIter,
std::vector< Scalar > &  B_avg,
std::vector< Scalar > &  residual_norms 
)
inherited

◆ hasNlddSolver()

template<class TypeTag >
bool Opm::BlackoilModel< TypeTag >::hasNlddSolver ( ) const
inlineinherited

Returns true if an NLDD solver exists.

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

◆ initialLinearization()

◆ isParallel()

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

◆ linearIterationsLastSolve()

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

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

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

◆ linearSolveSetupTime()

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

◆ localAccumulatedReports()

template<class TypeTag >
const SimulatorReport & Opm::BlackoilModel< TypeTag >::localAccumulatedReports
inherited

return the statistics of local solves accumulated for this rank

◆ 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 
)
inherited

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.

References OPM_BEGIN_PARALLEL_TRY_CATCH, and OPM_END_PARALLEL_TRY_CATCH.

◆ nonlinearIteration()

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

Perform a nonlinear iteration updating Flow and TPSA geomechanics.

Parameters
timerSimulation timer
nonlinear_solverNonlinear solver type
Returns
Report for simulator performance
Note
Strategies of coupling Flow and TPSA currently implemented:
  • fixed-stress: fixed-stress algorithm, i.e. iteratively solving Flow and TPSA equations in sequence
  • lagged: one-way coupling where Flow is solved with TPSA info from previous time step

References Opm::BlackoilModelTPSA< TypeTag >::nonlinearIterationFixedStressTPSA(), Opm::BlackoilModelTPSA< TypeTag >::nonlinearIterationLaggedTPSA(), and Opm::BlackoilModel< TypeTag >::simulator_.

◆ nonlinearIterationFixedStressTPSA()

template<class TypeTag >
template<class NonlinearSolverType >
SimulatorReportSingle Opm::BlackoilModelTPSA< TypeTag >::nonlinearIterationFixedStressTPSA ( const SimulatorTimerInterface timer,
NonlinearSolverType &  nonlinear_solver 
)
inline

Perform a nonlinear iteration updating Flow and TPSA geomechanics in a fixed-stress, iterative loop.

Parameters
timerSimulation timer
nonlinear_solverNonlinear solver type
Returns
Report for simulator performance

References Opm::SimulatorReportSingle::converged, Opm::BlackoilModel< TypeTag >::nonlinearIteration(), Opm::BlackoilModel< TypeTag >::simulator_, and Opm::BlackoilModelTPSA< TypeTag >::solveTpsaEquations().

Referenced by Opm::BlackoilModelTPSA< TypeTag >::nonlinearIteration().

◆ nonlinearIterationLaggedTPSA()

template<class TypeTag >
template<class NonlinearSolverType >
SimulatorReportSingle Opm::BlackoilModelTPSA< TypeTag >::nonlinearIterationLaggedTPSA ( const SimulatorTimerInterface timer,
NonlinearSolverType &  nonlinear_solver 
)
inline

Perform a nonlinear iteration updating Flow and TPSA geomechanics in a lagged scheme.

Parameters
iterationFlow nonlinear iteration
timerSimulation timer
nonlinear_solverNonlinear solver type
Returns
Report for simulator performance

References Opm::BlackoilModel< TypeTag >::nonlinearIteration(), Opm::BlackoilModel< TypeTag >::simulator_, and Opm::BlackoilModelTPSA< TypeTag >::solveTpsaEquations().

Referenced by Opm::BlackoilModelTPSA< TypeTag >::nonlinearIteration().

◆ nonlinearIterationNewton()

◆ numPhases()

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

The number of active fluid phases in the model.

◆ param()

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

Returns const reference to model parameters.

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

◆ popLastConvergenceReport()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::popLastConvergenceReport ( )
inlineinherited

Remove the last convergence report entry and residual norms history entry.

References Opm::BlackoilModel< TypeTag >::convergence_reports_, and Opm::BlackoilModel< TypeTag >::residual_norms_history_.

◆ prepareStep()

◆ prepareStoringSolutionUpdate()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::prepareStoringSolutionUpdate
inherited

Get solution update vector as a PrimaryVariable.

◆ relativeChange()

template<class TypeTag >
BlackoilModel< TypeTag >::Scalar Opm::BlackoilModel< TypeTag >::relativeChange
inherited

◆ simulator() [1/2]

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

◆ simulator() [2/2]

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

◆ solveJacobianSystem()

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

Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.

◆ solveTpsaEquations()

template<class TypeTag >
bool Opm::BlackoilModelTPSA< TypeTag >::solveTpsaEquations ( )
inline

Solve TPSA geomechanics equations.

Returns
Bool indicating TPSA convergence
Note
Calls Newton method for TPSA

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

Referenced by Opm::BlackoilModelTPSA< TypeTag >::nonlinearIterationFixedStressTPSA(), and Opm::BlackoilModelTPSA< TypeTag >::nonlinearIterationLaggedTPSA().

◆ stepReports()

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

◆ storeSolutionUpdate()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::storeSolutionUpdate ( const BVector dx)
inherited

◆ terminalOutputEnabled()

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

Return true if output to cout is wanted.

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

◆ updateSolution()

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

Apply an update to the primary variables.

◆ updateTUNING()

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

◆ updateTUNINGDP()

template<class TypeTag >
void Opm::BlackoilModel< TypeTag >::updateTUNINGDP ( const TuningDp &  tuning_dp)
inherited

◆ wellModel() [1/2]

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

◆ wellModel() [2/2]

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

◆ writeNonlinearIterationsPerCell()

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

Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.

◆ writePartitions()

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

Member Data Documentation

◆ biofilmVolumeFractionIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::biofilmVolumeFractionIdx = Indices::biofilmVolumeFractionIdx
staticconstexprinherited

◆ calciteVolumeFractionIdx

template<class TypeTag >
constexpr int Opm::BlackoilModel< TypeTag >::calciteVolumeFractionIdx = Indices::calciteVolumeFractionIdx
staticconstexprinherited

◆ compNames_

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

◆ contiBiofilmEqIdx

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

◆ contiBrineEqIdx

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

◆ contiCalciteEqIdx

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

◆ contiEnergyEqIdx

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

◆ contiFoamEqIdx

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

◆ contiMicrobialEqIdx

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

◆ contiOxygenEqIdx

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

◆ contiPolymerEqIdx

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

◆ contiPolymerMWEqIdx

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

◆ contiSolventEqIdx

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

◆ contiUreaEqIdx

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

◆ contiZfracEqIdx

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

◆ conv_monitor_

template<class TypeTag >
BlackoilModelConvergenceMonitor<Scalar> Opm::BlackoilModel< TypeTag >::conv_monitor_
protectedinherited

◆ convergence_reports_

template<class TypeTag >
std::vector<StepReport> Opm::BlackoilModel< TypeTag >::convergence_reports_
protectedinherited

◆ current_relaxation_

template<class TypeTag >
Scalar Opm::BlackoilModel< TypeTag >::current_relaxation_
protectedinherited

◆ dx_old_

template<class TypeTag >
BVector Opm::BlackoilModel< TypeTag >::dx_old_
protectedinherited

◆ enableSaltPrecipitation

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>()
staticconstexprinherited

◆ failureReport_

template<class TypeTag >
SimulatorReportSingle Opm::BlackoilModel< TypeTag >::failureReport_
protectedinherited

◆ foamConcentrationIdx

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

◆ global_nc_

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

The number of cells of the global grid.

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

◆ grid_

template<class TypeTag >
const Grid& Opm::BlackoilModel< TypeTag >::grid_
protectedinherited

◆ has_bioeffects_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_bioeffects_ = getPropValue<TypeTag, Properties::EnableBioeffects>()
staticconstexprprotectedinherited

◆ has_brine_

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

◆ has_energy_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_energy_ = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal
staticconstexprprotectedinherited

◆ has_extbo_

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

◆ has_foam_

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

◆ has_micp_

template<class TypeTag >
constexpr bool Opm::BlackoilModel< TypeTag >::has_micp_ = Indices::enableMICP
staticconstexprprotectedinherited

◆ has_polymer_

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

◆ has_polymermw_

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

◆ has_solvent_

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

◆ microbialConcentrationIdx

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

◆ nlddSolver_

template<class TypeTag >
std::unique_ptr<BlackoilModelNldd<TypeTag> > Opm::BlackoilModel< TypeTag >::nlddSolver_
protectedinherited

◆ numEq

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

◆ oxygenConcentrationIdx

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

◆ param_

template<class TypeTag >
ModelParameters Opm::BlackoilModel< TypeTag >::param_
protectedinherited

◆ polymerConcentrationIdx

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

◆ polymerMoleWeightIdx

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

◆ residual_norms_history_

template<class TypeTag >
std::vector<std::vector<Scalar> > Opm::BlackoilModel< TypeTag >::residual_norms_history_
protectedinherited

◆ saltConcentrationIdx

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

◆ simulator_

◆ solUpd_

template<class TypeTag >
SolutionVector Opm::BlackoilModel< TypeTag >::solUpd_
protectedinherited

◆ solventSaturationIdx

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

◆ temperatureIdx

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

◆ terminal_output_

template<class TypeTag >
bool Opm::BlackoilModel< TypeTag >::terminal_output_
protectedinherited

Whether we print something to std::cout.

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

◆ ureaConcentrationIdx

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

◆ well_model_

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

◆ zFractionIdx

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

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