#include <BlackoilPropertiesFromDeck.hpp>
|
| BlackoilPropertiesFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, const UnstructuredGrid &grid, bool init_rock=true) |
|
| BlackoilPropertiesFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, const UnstructuredGrid &grid, const parameter::ParameterGroup ¶m, bool init_rock=true) |
|
| BlackoilPropertiesFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, int number_of_cells, const int *global_cell, const int *cart_dims, bool init_rock=true) |
|
| BlackoilPropertiesFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, int number_of_cells, const int *global_cell, const int *cart_dims, const parameter::ParameterGroup ¶m, bool init_rock=true) |
|
| BlackoilPropertiesFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, std::shared_ptr< MaterialLawManager > materialLawManager, int number_of_cells, const int *global_cell, const int *cart_dims, const parameter::ParameterGroup ¶m, bool init_rock=true) |
|
virtual | ~BlackoilPropertiesFromDeck () |
| Destructor. More...
|
|
virtual int | numDimensions () const |
|
virtual int | numCells () const |
|
virtual const int * | cellPvtRegionIndex () const |
|
virtual const double * | porosity () const |
|
virtual const double * | permeability () const |
|
virtual int | numPhases () const |
|
virtual PhaseUsage | phaseUsage () const |
|
virtual void | viscosity (const int n, const double *p, const double *T, const double *z, const int *cells, double *mu, double *dmudp) const |
|
virtual void | matrix (const int n, const double *p, const double *T, const double *z, const int *cells, double *A, double *dAdp) const |
|
virtual void | density (const int n, const double *A, const int *cells, double *rho) const |
|
virtual const double * | surfaceDensity (int cellIdx=0) const |
|
virtual void | relperm (const int n, const double *s, const int *cells, double *kr, double *dkrds) const |
|
virtual void | capPress (const int n, const double *s, const int *cells, double *pc, double *dpcds) const |
|
virtual void | satRange (const int n, const int *cells, double *smin, double *smax) const |
|
virtual void | swatInitScaling (const int cell, const double pcow, double &swat) |
|
Concrete class implementing the blackoil property interface, reading all data and properties from eclipse deck input.
Opm::BlackoilPropertiesFromDeck::BlackoilPropertiesFromDeck |
( |
Opm::DeckConstPtr |
deck, |
|
|
Opm::EclipseStateConstPtr |
eclState, |
|
|
const UnstructuredGrid & |
grid, |
|
|
bool |
init_rock = true |
|
) |
| |
Initialize from deck and grid. - Parameters
-
[in] | deck | Deck input parser |
[in] | grid | Grid to which property object applies, needed for the mapping from cell indices (typically from a processed grid) to logical cartesian indices consistent with the deck. |
Opm::BlackoilPropertiesFromDeck::BlackoilPropertiesFromDeck |
( |
Opm::DeckConstPtr |
deck, |
|
|
Opm::EclipseStateConstPtr |
eclState, |
|
|
const UnstructuredGrid & |
grid, |
|
|
const parameter::ParameterGroup & |
param, |
|
|
bool |
init_rock = true |
|
) |
| |
Initialize from deck, grid and parameters. - Parameters
-
[in] | deck | Deck input parser |
[in] | grid | Grid to which property object applies, needed for the mapping from cell indices (typically from a processed grid) to logical cartesian indices consistent with the deck. |
[in] | param | Parameters. Accepted parameters include: pvt_tab_size (200) number of uniform sample points for dead-oil pvt tables. sat_tab_size (200) number of uniform sample points for saturation tables. threephase_model("simple") three-phase relperm model (accepts "simple" and "stone2"). For both size parameters, a 0 or negative value indicates that no spline fitting is to be done, and the input fluid data used directly for linear interpolation. |
Opm::BlackoilPropertiesFromDeck::BlackoilPropertiesFromDeck |
( |
Opm::DeckConstPtr |
deck, |
|
|
Opm::EclipseStateConstPtr |
eclState, |
|
|
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
const int * |
cart_dims, |
|
|
bool |
init_rock = true |
|
) |
| |
Opm::BlackoilPropertiesFromDeck::BlackoilPropertiesFromDeck |
( |
Opm::DeckConstPtr |
deck, |
|
|
Opm::EclipseStateConstPtr |
eclState, |
|
|
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
const int * |
cart_dims, |
|
|
const parameter::ParameterGroup & |
param, |
|
|
bool |
init_rock = true |
|
) |
| |
Opm::BlackoilPropertiesFromDeck::BlackoilPropertiesFromDeck |
( |
Opm::DeckConstPtr |
deck, |
|
|
Opm::EclipseStateConstPtr |
eclState, |
|
|
std::shared_ptr< MaterialLawManager > |
materialLawManager, |
|
|
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
const int * |
cart_dims, |
|
|
const parameter::ParameterGroup & |
param, |
|
|
bool |
init_rock = true |
|
) |
| |
virtual Opm::BlackoilPropertiesFromDeck::~BlackoilPropertiesFromDeck |
( |
| ) |
|
|
virtual |
virtual void Opm::BlackoilPropertiesFromDeck::capPress |
( |
const int |
n, |
|
|
const double * |
s, |
|
|
const int * |
cells, |
|
|
double * |
pc, |
|
|
double * |
dpcds |
|
) |
| const |
|
virtual |
- Parameters
-
[in] | n | Number of data points. |
[in] | s | Array of nP saturation values. |
[in] | cells | Array of n cell indices to be associated with the s values. |
[out] | pc | Array of nP capillary pressure values, array must be valid before calling. |
[out] | dpcds | If non-null: array of nP^2 derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = {dpc_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m_01 ...) |
Implements Opm::BlackoilPropertiesInterface.
virtual const int* Opm::BlackoilPropertiesFromDeck::cellPvtRegionIndex |
( |
| ) |
const |
|
inlinevirtual |
virtual void Opm::BlackoilPropertiesFromDeck::density |
( |
const int |
n, |
|
|
const double * |
A, |
|
|
const int * |
cells, |
|
|
double * |
rho |
|
) |
| const |
|
virtual |
Densities of stock components at reservoir conditions. - Parameters
-
[in] | n | Number of data points. |
[in] | A | Array of nP^2 values, where the P^2 values for a cell give the matrix A = RB^{-1} which relates z to u by z = Au. The matrices are assumed to be in Fortran order, and are typically the result of a call to the method matrix(). |
[in] | cells | The index of the grid cell of each data point. |
[out] | rho | Array of nP density values, array must be valid before calling. |
Implements Opm::BlackoilPropertiesInterface.
virtual void Opm::BlackoilPropertiesFromDeck::matrix |
( |
const int |
n, |
|
|
const double * |
p, |
|
|
const double * |
T, |
|
|
const double * |
z, |
|
|
const int * |
cells, |
|
|
double * |
A, |
|
|
double * |
dAdp |
|
) |
| const |
|
virtual |
- Parameters
-
[in] | n | Number of data points. |
[in] | p | Array of n pressure values. |
[in] | T | Array of n temperature values. |
[in] | z | Array of nP surface volume values. |
[in] | cells | Array of n cell indices to be associated with the p and z values. |
[out] | A | Array of nP^2 values, array must be valid before calling. The P^2 values for a cell give the matrix A = RB^{-1} which relates z to u by z = Au. The matrices are output in Fortran order. |
[out] | dAdp | If non-null: array of nP^2 matrix derivative values, array must be valid before calling. The matrices are output in Fortran order. |
Implements Opm::BlackoilPropertiesInterface.
virtual int Opm::BlackoilPropertiesFromDeck::numCells |
( |
| ) |
const |
|
virtual |
virtual int Opm::BlackoilPropertiesFromDeck::numDimensions |
( |
| ) |
const |
|
virtual |
virtual int Opm::BlackoilPropertiesFromDeck::numPhases |
( |
| ) |
const |
|
virtual |
virtual const double* Opm::BlackoilPropertiesFromDeck::permeability |
( |
| ) |
const |
|
virtual |
- Returns
- Array of ND^2 permeability values. The D^2 permeability values for a cell are organized as a matrix, which is symmetric (so ordering does not matter).
Implements Opm::BlackoilPropertiesInterface.
virtual PhaseUsage Opm::BlackoilPropertiesFromDeck::phaseUsage |
( |
| ) |
const |
|
virtual |
virtual const double* Opm::BlackoilPropertiesFromDeck::porosity |
( |
| ) |
const |
|
virtual |
virtual void Opm::BlackoilPropertiesFromDeck::relperm |
( |
const int |
n, |
|
|
const double * |
s, |
|
|
const int * |
cells, |
|
|
double * |
kr, |
|
|
double * |
dkrds |
|
) |
| const |
|
virtual |
- Parameters
-
[in] | n | Number of data points. |
[in] | s | Array of nP saturation values. |
[in] | cells | Array of n cell indices to be associated with the s values. |
[out] | kr | Array of nP relperm values, array must be valid before calling. |
[out] | dkrds | If non-null: array of nP^2 relperm derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = {dkr_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m_01 ...) |
Implements Opm::BlackoilPropertiesInterface.
virtual void Opm::BlackoilPropertiesFromDeck::satRange |
( |
const int |
n, |
|
|
const int * |
cells, |
|
|
double * |
smin, |
|
|
double * |
smax |
|
) |
| const |
|
virtual |
Obtain the range of allowable saturation values. In cell cells[i], saturation of phase p is allowed to be in the interval [smin[i*P + p], smax[i*P + p]]. - Parameters
-
[in] | n | Number of data points. |
[in] | cells | Array of n cell indices. |
[out] | smin | Array of nP minimum s values, array must be valid before calling. |
[out] | smax | Array of nP maximum s values, array must be valid before calling. |
Implements Opm::BlackoilPropertiesInterface.
virtual const double* Opm::BlackoilPropertiesFromDeck::surfaceDensity |
( |
int |
cellIdx = 0 | ) |
const |
|
virtual |
virtual void Opm::BlackoilPropertiesFromDeck::swatInitScaling |
( |
const int |
cell, |
|
|
const double |
pcow, |
|
|
double & |
swat |
|
) |
| |
|
virtual |
Update capillary pressure scaling according to pressure diff. and initial water saturation. - Parameters
-
[in] | cell | Cell index. |
[in] | pcow | P_oil - P_water. |
| in/out] | swat Water saturation. / Possibly modified Water saturation. |
Implements Opm::BlackoilPropertiesInterface.
virtual void Opm::BlackoilPropertiesFromDeck::viscosity |
( |
const int |
n, |
|
|
const double * |
p, |
|
|
const double * |
T, |
|
|
const double * |
z, |
|
|
const int * |
cells, |
|
|
double * |
mu, |
|
|
double * |
dmudp |
|
) |
| const |
|
virtual |
- Parameters
-
[in] | n | Number of data points. |
[in] | p | Array of n pressure values. |
[in] | T | Array of n temperature values. |
[in] | z | Array of nP surface volume values. |
[in] | cells | Array of n cell indices to be associated with the p and z values. |
[out] | mu | Array of nP viscosity values, array must be valid before calling. |
[out] | dmudp | If non-null: array of nP viscosity derivative values, array must be valid before calling. |
Implements Opm::BlackoilPropertiesInterface.
The documentation for this class was generated from the following file:
|