20 #ifndef OPM_BLACKOILPROPERTIESBASIC_HEADER_INCLUDED
21 #define OPM_BLACKOILPROPERTIESBASIC_HEADER_INCLUDED
24 #include <opm/core/props/BlackoilPropertiesInterface.hpp>
25 #include <opm/core/props/rock/RockBasic.hpp>
26 #include <opm/core/props/pvt/PvtPropertiesBasic.hpp>
27 #include <opm/core/props/satfunc/SaturationPropsBasic.hpp>
28 #include <opm/core/utility/parameters/ParameterGroup.hpp>
68 virtual const double*
porosity()
const;
111 virtual void matrix(
const int n,
128 virtual void density(
const int n,
147 virtual void relperm(
const int n,
151 double* dkrds)
const;
167 double* dpcds)
const;
203 #endif // OPM_BLACKOILPROPERTIESBASIC_HEADER_INCLUDED
Class encapsulating basic saturation function behaviour, by which we mean constant, linear or quadratic relative permeability functions for a maximum of two phases, and zero capillary pressure.
Definition: SaturationPropsBasic.hpp:36
virtual int numPhases() const
Definition: BlackoilPropertiesBasic.cpp:81
virtual const double * surfaceDensity(int cellIdx=0) const
Densities of stock components at surface conditions.
Definition: BlackoilPropertiesBasic.cpp:186
virtual void swatInitScaling(const int cell, const double pcow, double &swat)
Update capillary pressure scaling according to pressure diff.
Definition: BlackoilPropertiesBasic.cpp:250
virtual void matrix(const int n, const double *p, const double *T, const double *z, const int *cells, double *A, double *dAdp) const
Definition: BlackoilPropertiesBasic.cpp:126
virtual void density(const int n, const double *A, const int *cells, double *rho) const
Densities of stock components at reservoir conditions.
Definition: BlackoilPropertiesBasic.cpp:166
Concrete class implementing the blackoil property interface, reading all necessary input from paramet...
Definition: BlackoilPropertiesBasic.hpp:35
virtual int numDimensions() const
Definition: BlackoilPropertiesBasic.cpp:52
virtual const double * permeability() const
Definition: BlackoilPropertiesBasic.cpp:72
virtual void capPress(const int n, const double *s, const int *cells, double *pc, double *dpcds) const
Definition: BlackoilPropertiesBasic.cpp:220
Class collecting simple pvt properties for 1-3 phases.
Definition: PvtPropertiesBasic.hpp:36
BlackoilPropertiesBasic(const ParameterGroup ¶m, const int dim, const int num_cells)
Construct from parameters.
Definition: BlackoilPropertiesBasic.cpp:29
virtual int numCells() const
Definition: BlackoilPropertiesBasic.cpp:58
virtual void viscosity(const int n, const double *p, const double *T, const double *z, const int *cells, double *mu, double *dmudp) const
Definition: BlackoilPropertiesBasic.cpp:100
virtual void satRange(const int n, const int *cells, double *smin, double *smax) const
Obtain the range of allowable saturation values.
Definition: BlackoilPropertiesBasic.cpp:237
Abstract base class for blackoil fluid and reservoir properties.
Definition: BlackoilPropertiesInterface.hpp:37
Definition: BlackoilPhases.hpp:36
Definition: RockBasic.hpp:30
virtual void relperm(const int n, const double *s, const int *cells, double *kr, double *dkrds) const
Definition: BlackoilPropertiesBasic.cpp:201
virtual PhaseUsage phaseUsage() const
Definition: BlackoilPropertiesBasic.cpp:87
virtual const double * porosity() const
Definition: BlackoilPropertiesBasic.cpp:64
ParameterGroup is a class that is used to provide run-time parameters.
Definition: ParameterGroup.hpp:81
virtual ~BlackoilPropertiesBasic()
Destructor.
Definition: BlackoilPropertiesBasic.cpp:46
virtual const int * cellPvtRegionIndex() const
Return an array containing the PVT table index for each grid cell.
Definition: BlackoilPropertiesBasic.hpp:64