Opm::BlackoilModelBase< Grid, Implementation > Class Template Reference

#include <BlackoilModelBase.hpp>

Inheritance diagram for Opm::BlackoilModelBase< Grid, Implementation >:
Inheritance graph

Classes

struct  ReservoirResidualQuant
 
struct  WellOps
 

Public Types

typedef AutoDiffBlock< double > ADB
 
typedef ADB::V V
 
typedef ADB::M M
 
typedef ModelTraits
< Implementation >
::ReservoirState 
ReservoirState
 
typedef ModelTraits
< Implementation >::WellState 
WellState
 
typedef ModelTraits
< Implementation >
::ModelParameters 
ModelParameters
 
typedef ModelTraits
< Implementation >
::SolutionState 
SolutionState
 

Public Member Functions

 BlackoilModelBase (const ModelParameters &param, const Grid &grid, const BlackoilPropsAdInterface &fluid, const DerivedGeology &geo, const RockCompressibility *rock_comp_props, const Wells *wells, const NewtonIterationBlackoilInterface &linsolver, Opm::EclipseStateConstPtr eclState, const bool has_disgas, const bool has_vapoil, const bool terminal_output)
 
void setThresholdPressures (const std::vector< double > &threshold_pressures_by_face)
 Set threshold pressures that prevent or reduce flow. This prevents flow across faces if the potential difference is less than the threshold. If the potential difference is greater, the threshold value is subtracted before calculating flow. This is treated symmetrically, so flow is prevented or reduced in both directions equally. More...
 
void prepareStep (const double dt, ReservoirState &reservoir_state, WellState &well_state)
 
void afterStep (const double dt, ReservoirState &reservoir_state, WellState &well_state)
 
void assemble (const ReservoirState &reservoir_state, WellState &well_state, const bool initial_assembly)
 
std::vector< double > computeResidualNorms () const
 Compute the residual norms of the mass balance for each phase, the well flux, and the well equation. More...
 
int sizeNonLinear () const
 The size (number of unknowns) of the nonlinear system of equations. More...
 
int linearIterationsLastSolve () const
 Number of linear iterations used in last call to solveJacobianSystem(). More...
 
V solveJacobianSystem () const
 
void updateState (const V &dx, ReservoirState &reservoir_state, WellState &well_state)
 
bool terminalOutputEnabled () const
 Return true if output to cout is wanted. More...
 
bool getConvergence (const double dt, const int iteration)
 
int numPhases () const
 The number of active fluid phases in the model. More...
 
int numMaterials () const
 
const std::string & materialName (int material_index) const
 
void updateEquationsScaling ()
 Update the scaling factors for mass balance equations. More...
 

Protected Types

typedef Eigen::Array< double,
Eigen::Dynamic, Eigen::Dynamic,
Eigen::RowMajor > 
DataBlock
 

Protected Member Functions

Implementation & asImpl ()
 
const Implementation & asImpl () const
 
bool wellsActive () const
 
bool localWellsActive () const
 
const Wells & wells () const
 
void makeConstantState (SolutionState &state) const
 
SolutionState variableState (const ReservoirState &x, const WellState &xw) const
 
std::vector< VvariableStateInitials (const ReservoirState &x, const WellState &xw) const
 
void variableReservoirStateInitials (const ReservoirState &x, std::vector< V > &vars0) const
 
void variableWellStateInitials (const WellState &xw, std::vector< V > &vars0) const
 
std::vector< int > variableStateIndices () const
 
std::vector< int > variableWellStateIndices () const
 
SolutionState variableStateExtractVars (const ReservoirState &x, const std::vector< int > &indices, std::vector< ADB > &vars) const
 
void variableStateExtractWellsVars (const std::vector< int > &indices, std::vector< ADB > &vars, SolutionState &state) const
 
void computeAccum (const SolutionState &state, const int aix)
 
void computeWellConnectionPressures (const SolutionState &state, const WellState &xw)
 
void assembleMassBalanceEq (const SolutionState &state)
 
void solveWellEq (const std::vector< ADB > &mob_perfcells, const std::vector< ADB > &b_perfcells, SolutionState &state, WellState &well_state)
 
void computeWellFlux (const SolutionState &state, const std::vector< ADB > &mob_perfcells, const std::vector< ADB > &b_perfcells, V &aliveWells, std::vector< ADB > &cq_s)
 
void updatePerfPhaseRatesAndPressures (const std::vector< ADB > &cq_s, const SolutionState &state, WellState &xw)
 
void addWellFluxEq (const std::vector< ADB > &cq_s, const SolutionState &state)
 
void addWellContributionToMassBalanceEq (const std::vector< ADB > &cq_s, const SolutionState &state, const WellState &xw)
 
void addWellControlEq (const SolutionState &state, const WellState &xw, const V &aliveWells)
 
void updateWellControls (WellState &xw) const
 
void updateWellState (const V &dwells, WellState &well_state)
 
bool getWellConvergence (const int iteration)
 
bool isVFPActive () const
 
std::vector< ADBcomputePressures (const ADB &po, const ADB &sw, const ADB &so, const ADB &sg) const
 
V computeGasPressure (const V &po, const V &sw, const V &so, const V &sg) const
 
std::vector< ADBcomputeRelPerm (const SolutionState &state) const
 
void computeMassFlux (const int actph, const V &transi, const ADB &kr, const ADB &p, const SolutionState &state)
 
void applyThresholdPressures (ADB &dp)
 
ADB fluidViscosity (const int phase, const ADB &p, const ADB &temp, const ADB &rs, const ADB &rv, const std::vector< PhasePresence > &cond) const
 
ADB fluidReciprocFVF (const int phase, const ADB &p, const ADB &temp, const ADB &rs, const ADB &rv, const std::vector< PhasePresence > &cond) const
 
ADB fluidDensity (const int phase, const ADB &b, const ADB &rs, const ADB &rv) const
 
V fluidRsSat (const V &p, const V &so, const std::vector< int > &cells) const
 
ADB fluidRsSat (const ADB &p, const ADB &so, const std::vector< int > &cells) const
 
V fluidRvSat (const V &p, const V &so, const std::vector< int > &cells) const
 
ADB fluidRvSat (const ADB &p, const ADB &so, const std::vector< int > &cells) const
 
ADB poroMult (const ADB &p) const
 
ADB transMult (const ADB &p) const
 
const std::vector< PhasePresence > phaseCondition () const
 
void classifyCondition (const ReservoirState &state)
 
void updatePrimalVariableFromState (const ReservoirState &state)
 
void updatePhaseCondFromPrimalVariable ()
 Update the phaseCondition_ member based on the primalVariable_ member. More...
 
double convergenceReduction (const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &B, const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &tempV, const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &R, std::vector< double > &R_sum, std::vector< double > &maxCoeff, std::vector< double > &B_avg, std::vector< double > &maxNormWell, int nc, int nw) const
 Compute the reduction within the convergence check. More...
 
double dpMaxRel () const
 
double dsMax () const
 
double drMaxRel () const
 
double maxResidualAllowed () const
 

Protected Attributes

const Grid & grid_
 
const BlackoilPropsAdInterfacefluid_
 
const DerivedGeologygeo_
 
const RockCompressibility * rock_comp_props_
 
const Wells * wells_
 
VFPProperties vfp_properties_
 
const
NewtonIterationBlackoilInterface
linsolver_
 
const std::vector< bool > active_
 
const std::vector< int > canph_
 
const std::vector< int > cells_
 
HelperOps ops_
 
const WellOps wops_
 
const bool has_disgas_
 
const bool has_vapoil_
 
ModelParameters param_
 
bool use_threshold_pressure_
 
bool wells_active_
 
V threshold_pressures_by_interior_face_
 
std::vector
< ReservoirResidualQuant
rq_
 
std::vector< PhasePresence > phaseCondition_
 
V isRs_
 
V isRv_
 
V isSg_
 
V well_perforation_densities_
 
V well_perforation_pressure_diffs_
 
LinearisedBlackoilResidual residual_
 
bool terminal_output_
 Whether we print something to std::cout. More...
 
std::vector< int > primalVariable_
 
V pvdt_
 
std::vector< std::string > material_name_
 

Detailed Description

template<class Grid, class Implementation>
class Opm::BlackoilModelBase< Grid, Implementation >

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.

It uses automatic differentiation via the class AutoDiffBlock to simplify assembly of the jacobian matrix.

Template Parameters
GridUnstructuredGrid or CpGrid.
ImplementationProvides concrete state types.

Member Typedef Documentation

template<class Grid, class Implementation>
typedef AutoDiffBlock<double> Opm::BlackoilModelBase< Grid, Implementation >::ADB
template<class Grid, class Implementation>
typedef Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> Opm::BlackoilModelBase< Grid, Implementation >::DataBlock
protected
template<class Grid, class Implementation>
typedef ADB::M Opm::BlackoilModelBase< Grid, Implementation >::M
template<class Grid, class Implementation>
typedef ModelTraits<Implementation>::ModelParameters Opm::BlackoilModelBase< Grid, Implementation >::ModelParameters
template<class Grid, class Implementation>
typedef ModelTraits<Implementation>::ReservoirState Opm::BlackoilModelBase< Grid, Implementation >::ReservoirState
template<class Grid, class Implementation>
typedef ModelTraits<Implementation>::SolutionState Opm::BlackoilModelBase< Grid, Implementation >::SolutionState
template<class Grid, class Implementation>
typedef ADB::V Opm::BlackoilModelBase< Grid, Implementation >::V
template<class Grid, class Implementation>
typedef ModelTraits<Implementation>::WellState Opm::BlackoilModelBase< Grid, Implementation >::WellState

Constructor & Destructor Documentation

template<class Grid, class Implementation >
Opm::BlackoilModelBase< Grid, Implementation >::BlackoilModelBase ( const ModelParameters param,
const Grid &  grid,
const BlackoilPropsAdInterface fluid,
const DerivedGeology geo,
const RockCompressibility *  rock_comp_props,
const Wells *  wells,
const NewtonIterationBlackoilInterface linsolver,
Opm::EclipseStateConstPtr  eclState,
const bool  has_disgas,
const bool  has_vapoil,
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
[in]paramparameters
[in]gridgrid data structure
[in]fluidfluid properties
[in]georock properties
[in]rock_comp_propsif non-null, rock compressibility properties
[in]wellswell structure
[in]vfp_propertiesVertical flow performance tables
[in]linsolverlinear solver
[in]eclStateeclipse state
[in]has_disgasturn on dissolved gas
[in]has_vapoilturn on vaporized oil feature
[in]terminal_outputrequest output to cout/cerr

References Opm::AutoDiffBlock< double >::null(), and Opm::BlackoilPropsAdInterface::numPhases().

Member Function Documentation

template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::addWellContributionToMassBalanceEq ( const std::vector< ADB > &  cq_s,
const SolutionState state,
const WellState xw 
)
protected
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::afterStep ( const double  dt,
ReservoirState reservoir_state,
WellState well_state 
)

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

Parameters
[in]dttime step size
[in,out]reservoir_statereservoir state variables
[in,out]well_statewell state variables
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::applyThresholdPressures ( ADB dp)
protected
template<class Grid, class Implementation>
const Implementation& Opm::BlackoilModelBase< Grid, Implementation >::asImpl ( ) const
inlineprotected

Access the most-derived class used for static polymorphism (CRTP).

template<class Grid , class Implementation >
std::vector< ADB > Opm::BlackoilModelBase< Grid, Implementation >::computePressures ( const ADB po,
const ADB sw,
const ADB so,
const ADB sg 
) const
protected
template<class Grid , class Implementation >
std::vector< double > Opm::BlackoilModelBase< Grid, Implementation >::computeResidualNorms ( ) const

Compute the residual norms of the mass balance for each phase, the well flux, and the well equation.

Returns
a vector that contains for each phase the norm of the mass balance and afterwards the norm of the residual of the well flux and the well equation.

References Opm::detail::infinityNorm(), Opm::detail::infinityNormWell(), Opm::BlackoilModelBase< Grid, Implementation >::linsolver_, Opm::LinearisedBlackoilResidual::material_balance_eq, Opm::NewtonIterationBlackoilInterface::parallelInformation(), Opm::BlackoilModelBase< Grid, Implementation >::residual_, Opm::LinearisedBlackoilResidual::well_eq, and Opm::LinearisedBlackoilResidual::well_flux_eq.

template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::computeWellConnectionPressures ( const SolutionState state,
const WellState xw 
)
protected
template<class Grid , class Implementation >
double Opm::BlackoilModelBase< Grid, Implementation >::convergenceReduction ( const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &  B,
const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &  tempV,
const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &  R,
std::vector< double > &  R_sum,
std::vector< double > &  maxCoeff,
std::vector< double > &  B_avg,
std::vector< double > &  maxNormWell,
int  nc,
int  nw 
) const
protected

Compute the reduction within the convergence check.

Parameters
[in]BA matrix with MaxNumPhases columns and the same number rows as the number of cells of the grid. B.col(i) contains the values for phase i.
[in]tempVA matrix with MaxNumPhases columns and the same number rows as the number of cells of the grid. tempV.col(i) contains the values for phase i.
[in]RA matrix with MaxNumPhases columns and the same number rows as the number of cells of the grid. B.col(i) contains the values for phase i.
[out]R_sumAn array of size MaxNumPhases where entry i contains the sum of R for the phase i.
[out]maxCoeffAn array of size MaxNumPhases where entry i contains the maximum of tempV for the phase i.
[out]B_avgAn array of size MaxNumPhases where entry i contains the average of B for the phase i.
[out]maxNormWellThe maximum of the well equations for each phase.
[in]ncThe number of cells of the local grid.
[in]nwThe number of wells on the local grid.
Returns
The total pore volume over all cells.

References Opm::BlackoilModelBase< Grid, Implementation >::asImpl(), Opm::BlackoilModelBase< Grid, Implementation >::geo_, Opm::BlackoilModelBase< Grid, Implementation >::linsolver_, Opm::NewtonIterationBlackoilInterface::parallelInformation(), Opm::DerivedGeology::poreVolume(), Opm::BlackoilModelBase< Grid, Implementation >::residual_, Opm::AutoDiffBlock< Scalar >::value(), and Opm::LinearisedBlackoilResidual::well_flux_eq.

Referenced by Opm::BlackoilModelBase< Grid, Implementation >::getConvergence(), and Opm::BlackoilModelBase< Grid, Implementation >::getWellConvergence().

template<class Grid, class Implementation>
double Opm::BlackoilModelBase< Grid, Implementation >::dpMaxRel ( ) const
inlineprotected
template<class Grid, class Implementation>
double Opm::BlackoilModelBase< Grid, Implementation >::drMaxRel ( ) const
inlineprotected
template<class Grid, class Implementation>
double Opm::BlackoilModelBase< Grid, Implementation >::dsMax ( ) const
inlineprotected
template<class Grid , class Implementation >
ADB Opm::BlackoilModelBase< Grid, Implementation >::fluidDensity ( const int  phase,
const ADB b,
const ADB rs,
const ADB rv 
) const
protected
template<class Grid , class Implementation >
ADB Opm::BlackoilModelBase< Grid, Implementation >::fluidReciprocFVF ( const int  phase,
const ADB p,
const ADB temp,
const ADB rs,
const ADB rv,
const std::vector< PhasePresence > &  cond 
) const
protected
template<class Grid , class Implementation >
ADB Opm::BlackoilModelBase< Grid, Implementation >::fluidRsSat ( const ADB p,
const ADB so,
const std::vector< int > &  cells 
) const
protected
template<class Grid , class Implementation >
ADB Opm::BlackoilModelBase< Grid, Implementation >::fluidRvSat ( const ADB p,
const ADB so,
const std::vector< int > &  cells 
) const
protected
template<class Grid , class Implementation >
ADB Opm::BlackoilModelBase< Grid, Implementation >::fluidViscosity ( const int  phase,
const ADB p,
const ADB temp,
const ADB rs,
const ADB rv,
const std::vector< PhasePresence > &  cond 
) const
protected
template<class Grid , class Implementation >
int Opm::BlackoilModelBase< Grid, Implementation >::linearIterationsLastSolve ( ) const

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

template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::makeConstantState ( SolutionState state) const
protected
template<class Grid , class Implementation >
const std::string & Opm::BlackoilModelBase< Grid, Implementation >::materialName ( int  material_index) const

The name of an active material in the model. It is required that material_index < numMaterials().

Referenced by Opm::BlackoilModelBase< Grid, Implementation >::getConvergence(), and Opm::BlackoilModelBase< Grid, Implementation >::getWellConvergence().

template<class Grid, class Implementation>
double Opm::BlackoilModelBase< Grid, Implementation >::maxResidualAllowed ( ) const
inlineprotected
template<class Grid , class Implementation >
int Opm::BlackoilModelBase< Grid, Implementation >::numMaterials ( ) const

The number of active materials in the model. This should be equal to the number of material balance equations.

template<class Grid , class Implementation >
int Opm::BlackoilModelBase< Grid, Implementation >::numPhases ( ) const

The number of active fluid phases in the model.

template<class Grid, class Implementation>
const std::vector<PhasePresence> Opm::BlackoilModelBase< Grid, Implementation >::phaseCondition ( ) const
inlineprotected
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::prepareStep ( const double  dt,
ReservoirState reservoir_state,
WellState well_state 
)

Called once before each time step.

Parameters
[in]dttime step size
[in,out]reservoir_statereservoir state variables
[in,out]well_statewell state variables

References Opm::Gas.

template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::setThresholdPressures ( const std::vector< double > &  threshold_pressures_by_face)

Set threshold pressures that prevent or reduce flow. This prevents flow across faces if the potential difference is less than the threshold. If the potential difference is greater, the threshold value is subtracted before calculating flow. This is treated symmetrically, so flow is prevented or reduced in both directions equally.

Parameters
[in]threshold_pressures_by_facearray of size equal to the number of faces of the grid passed in the constructor.

Referenced by Opm::SimulatorFullyImplicitBlackoilSolvent< GridT >::createSolver().

template<class Grid , class Implementation >
int Opm::BlackoilModelBase< Grid, Implementation >::sizeNonLinear ( ) const

The size (number of unknowns) of the nonlinear system of equations.

template<class Grid , class Implementation >
V Opm::BlackoilModelBase< Grid, Implementation >::solveJacobianSystem ( ) const
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::solveWellEq ( const std::vector< ADB > &  mob_perfcells,
const std::vector< ADB > &  b_perfcells,
SolutionState state,
WellState well_state 
)
protected
template<class Grid , class Implementation >
bool Opm::BlackoilModelBase< Grid, Implementation >::terminalOutputEnabled ( ) const

Return true if output to cout is wanted.

template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::updatePerfPhaseRatesAndPressures ( const std::vector< ADB > &  cq_s,
const SolutionState state,
WellState xw 
)
protected
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::updatePhaseCondFromPrimalVariable ( )
protected
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::updateState ( const V dx,
ReservoirState reservoir_state,
WellState well_state 
)

Apply an update to the primary variables, chopped if appropriate.

Parameters
[in]dxupdates to apply to primary variables
[in,out]reservoir_statereservoir state variables
[in,out]well_statewell state variables

References Opm::BlackoilModelBase< Grid, Implementation >::active_, Opm::BlackoilModelBase< Grid, Implementation >::cells_, Opm::BlackoilModelBase< Grid, Implementation >::computeGasPressure(), Opm::BlackoilModelBase< Grid, Implementation >::dpMaxRel(), Opm::BlackoilModelBase< Grid, Implementation >::drMaxRel(), Opm::BlackoilModelBase< Grid, Implementation >::dsMax(), Opm::BlackoilModelBase< Grid, Implementation >::fluid_, Opm::BlackoilModelBase< Grid, Implementation >::fluidRsSat(), Opm::BlackoilModelBase< Grid, Implementation >::fluidRvSat(), Opm::Gas, Opm::BlackoilModelBase< Grid, Implementation >::grid_, Opm::BlackoilModelBase< Grid, Implementation >::has_disgas_, Opm::BlackoilModelBase< Grid, Implementation >::has_vapoil_, Opm::BlackoilModelBase< Grid, Implementation >::isRs_, Opm::BlackoilModelBase< Grid, Implementation >::isRv_, Opm::BlackoilModelBase< Grid, Implementation >::isSg_, Opm::BlackoilModelBase< Grid, Implementation >::localWellsActive(), Opm::BlackoilPropsAdInterface::numPhases(), Opm::Oil, Opm::BlackoilPropsAdInterface::phaseUsage(), Opm::BlackoilModelBase< Grid, Implementation >::primalVariable_, Opm::RS, Opm::RV, Opm::Sg, Opm::sign(), Opm::sqrt(), Opm::subset(), Opm::BlackoilModelBase< Grid, Implementation >::updatePhaseCondFromPrimalVariable(), Opm::BlackoilModelBase< Grid, Implementation >::updateWellState(), Opm::Water, and Opm::BlackoilModelBase< Grid, Implementation >::wells().

template<class Grid , class Implementation >
BlackoilModelBase< Grid, Implementation >::SolutionState Opm::BlackoilModelBase< Grid, Implementation >::variableState ( const ReservoirState x,
const WellState xw 
) const
protected
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::variableStateExtractWellsVars ( const std::vector< int > &  indices,
std::vector< ADB > &  vars,
SolutionState state 
) const
protected
template<class Grid , class Implementation >
std::vector< V > Opm::BlackoilModelBase< Grid, Implementation >::variableStateInitials ( const ReservoirState x,
const WellState xw 
) const
protected
template<class Grid , class Implementation >
std::vector< int > Opm::BlackoilModelBase< Grid, Implementation >::variableWellStateIndices ( ) const
protected
template<class Grid , class Implementation >
void Opm::BlackoilModelBase< Grid, Implementation >::variableWellStateInitials ( const WellState xw,
std::vector< V > &  vars0 
) const
protected
template<class Grid, class Implementation>
bool Opm::BlackoilModelBase< Grid, Implementation >::wellsActive ( ) const
inlineprotected

Member Data Documentation

template<class Grid, class Implementation>
const std::vector<bool> Opm::BlackoilModelBase< Grid, Implementation >::active_
protected
template<class Grid, class Implementation>
const std::vector<int> Opm::BlackoilModelBase< Grid, Implementation >::canph_
protected
template<class Grid, class Implementation>
const BlackoilPropsAdInterface& Opm::BlackoilModelBase< Grid, Implementation >::fluid_
protected
template<class Grid, class Implementation>
std::vector<std::string> Opm::BlackoilModelBase< Grid, Implementation >::material_name_
protected
template<class Grid, class Implementation>
V Opm::BlackoilModelBase< Grid, Implementation >::pvdt_
protected
template<class Grid, class Implementation>
const RockCompressibility* Opm::BlackoilModelBase< Grid, Implementation >::rock_comp_props_
protected
template<class Grid, class Implementation>
V Opm::BlackoilModelBase< Grid, Implementation >::threshold_pressures_by_interior_face_
protected
template<class Grid, class Implementation>
bool Opm::BlackoilModelBase< Grid, Implementation >::use_threshold_pressure_
protected
template<class Grid, class Implementation>
const Wells* Opm::BlackoilModelBase< Grid, Implementation >::wells_
protected
template<class Grid, class Implementation>
bool Opm::BlackoilModelBase< Grid, Implementation >::wells_active_
protected

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