opm-upscaling
Loading...
Searching...
No Matches
Opm::ReservoirPropertyCommon< dim, RPImpl, RockType > Class Template Reference

A property class for incompressible two-phase flow. More...

#include <ReservoirPropertyCommon.hpp>

Public Types

enum  { NumberOfPhases = 2 }
 The number of phases.
typedef ImmutableCMatrix PermTensor
 Tensor type for read-only access to permeability.
typedef OwnCMatrix MutablePermTensor
 Tensor type to be used for holding copies of permeability tensors.
typedef SharedCMatrix SharedPermTensor
 Tensor type for read and write access to permeability.

Public Member Functions

 ReservoirPropertyCommon ()
 Default constructor.
void init (const Opm::Deck &, const std::vector< int > &global_cell, const double perm_threshold=0.0, const std::string *rock_list_filename=0, const bool use_jfunction_scaling=true, const double sigma=1.0, const double theta=0.0)
 Initialize from a grdecl file.
void init (const int num_cells, const double uniform_poro=0.2, const double uniform_perm=100.0 *Opm::prefix::milli *Opm::unit::darcy)
 Initialize a uniform reservoir.
void setViscosities (double v1, double v2)
 Set viscosities of both faces.
void setDensities (double d1, double d2)
 Set densitities of both faces.
double viscosityFirstPhase () const
 Viscosity of first (water) phase.
double viscositySecondPhase () const
 Viscosity of second (oil) phase.
double densityFirstPhase () const
 Density of first (water) phase.
double densitySecondPhase () const
 Density of second (oil) phase.
double porosity (int cell_index) const
 Read-access to porosity.
double ntg (int cell_index) const
 Read-access to ntg.
double swcr (int cell_index) const
 Read-access to swcr.
double sowcr (int cell_index) const
 Read-access to sowcr.
PermTensor permeability (int cell_index) const
 Read-access to permeability.
SharedPermTensor permeabilityModifiable (int cell_index)
 Read- and write-access to permeability.
template<class Vector>
void phaseDensities (int, Vector &density) const
 Densities for both phases.
double densityDifference () const
 Difference of densities.
double cflFactor () const
 A factor useful in cfl computations.
double cflFactorGravity () const
 A factor useful in gravity cfl computations.
double cflFactorCapillary () const
 A factor useful in gravity cfl computations.
double capillaryPressure (int cell_index, double saturation) const
 Capillary pressure.
double capillaryPressureDeriv (int c, double s) const
 Derivative of Capillary pressure.
double s_min (int c) const
double s_max (int c) const
double saturationFromCapillaryPressure (int cell_index, double cap_press) const
 Inverse of the capillary pressure function.
void writeSintefLegacyFormat (const std::string &grid_prefix) const
 Write permeability and porosity in the Sintef legacy format.

Protected Member Functions

void assignPorosity (const Opm::Deck &deck, const std::vector< int > &global_cell)
void assignNTG (const Opm::Deck &deck, const std::vector< int > &global_cell)
void assignSWCR (const Opm::Deck &deck, const std::vector< int > &global_cell)
void assignSOWCR (const Opm::Deck &deck, const std::vector< int > &global_cell)
void assignPermeability (const Opm::Deck &deck, const std::vector< int > &global_cell, const double perm_threshold)
void assignRockTable (const Opm::Deck &deck, const std::vector< int > &global_cell)
void readRocks (const std::string &rock_list_file)
RPImpl & asImpl ()

Protected Attributes

std::vector< double > porosity_
std::vector< double > ntg_
std::vector< double > swcr_
std::vector< double > sowcr_
std::vector< double > permeability_
std::vector< unsigned char > permfield_valid_
double density1_
double density2_
double viscosity1_
double viscosity2_
double cfl_factor_
double cfl_factor_gravity_
double cfl_factor_capillary_
std::vector< RockType > rock_
std::vector< int > cell_to_rock_
PermeabilityKind permeability_kind_

Detailed Description

template<int dim, class RPImpl, class RockType>
class Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >

A property class for incompressible two-phase flow.

Template Parameters
dimthe dimension of the space, used for giving permeability tensors the right size.

Member Function Documentation

◆ capillaryPressure()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::capillaryPressure ( int cell_index,
double saturation ) const

Capillary pressure.

Parameters
cell_indexindex of a grid cell.
saturationa saturation value.
Returns
capillary pressure at the given cell and saturation.

◆ capillaryPressureDeriv()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::capillaryPressureDeriv ( int c,
double s ) const

Derivative of Capillary pressure.

Parameters
cindex of a grid cell.
ssaturation value.
Returns
capillary pressure at the given cell and saturation.

◆ cflFactor()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::cflFactor ( ) const

A factor useful in cfl computations.

Returns
the cfl factor.

◆ cflFactorCapillary()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::cflFactorCapillary ( ) const

A factor useful in gravity cfl computations.

Returns
the capillary cfl factor.

◆ cflFactorGravity()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::cflFactorGravity ( ) const

A factor useful in gravity cfl computations.

Returns
the gravity cfl factor.

◆ densityDifference()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::densityDifference ( ) const

Difference of densities.

Returns
densityFirstPhase() - densitySecondPhase()

◆ densityFirstPhase()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::densityFirstPhase ( ) const

Density of first (water) phase.

Returns
the density value.

◆ densitySecondPhase()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::densitySecondPhase ( ) const

Density of second (oil) phase.

Returns
the density value.

◆ init() [1/2]

template<int dim, class RPImpl, class RockType>
void Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::init ( const int num_cells,
const double uniform_poro = 0.2,
const double uniform_perm = 100.0*Opm::prefix::milli*Opm::unit::darcy )

Initialize a uniform reservoir.

Parameters
num_cellsnumber of cells in the grid.
uniform_porothe uniform porosity.
uniform_permthe uniform (scalar) permeability.

◆ init() [2/2]

template<int dim, class RPImpl, class RockType>
void Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::init ( const Opm::Deck & deck,
const std::vector< int > & global_cell,
const double perm_threshold = 0.0,
const std::string * rock_list_filename = 0,
const bool use_jfunction_scaling = true,
const double sigma = 1.0,
const double theta = 0.0 )

Initialize from a grdecl file.

Parameters
deckthe deck holding the grdecl data.
global_cellthe mapping from cell indices to the logical cartesian indices of the grdecl file.
perm_thresholdlower threshold for permeability.
rock_list_filenameif non-null, the referred string gives the filename for the rock list.
use_jfunction_scalingif true, use j-function scaling of capillary pressure, if applicable.
sigmainterface tension for j-scaling, if applicable.
thetaangle for j-scaling, if applicable.

◆ ntg()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::ntg ( int cell_index) const

Read-access to ntg.

Parameters
cell_indexindex of a grid cell.
Returns
ntg value of the cell.

◆ permeability()

template<int dim, class RPImpl, class RockType>
ReservoirPropertyCommon< dim, RPImpl, RockType >::PermTensor Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::permeability ( int cell_index) const

Read-access to permeability.

Parameters
cell_indexindex of a grid cell.
Returns
permeability value of the cell.

◆ permeabilityModifiable()

template<int dim, class RPImpl, class RockType>
ReservoirPropertyCommon< dim, RPImpl, RockType >::SharedPermTensor Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::permeabilityModifiable ( int cell_index)

Read- and write-access to permeability.

Use with caution.

Parameters
cell_indexindex of a grid cell.
Returns
permeability value of the cell.

◆ phaseDensities()

template<int dim, class RPImpl, class RockType>
template<class Vector>
void Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::phaseDensities ( int ,
Vector & density ) const

Densities for both phases.

Template Parameters
Vectora class with size() and operator[].
Parameters
cell_indexindex of a grid cell (not used).
[out]densitythe phase densities. Expected to be of size 2 before (and after) the call.

◆ porosity()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::porosity ( int cell_index) const

Read-access to porosity.

Parameters
cell_indexindex of a grid cell.
Returns
porosity value of the cell.

◆ s_max()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::s_max ( int c) const
Parameters
cindex of a grid cell.
Returns
maximum saturation in given cell.

◆ s_min()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::s_min ( int c) const
Parameters
cindex of a grid cell.
Returns
minimum saturation in given cell.

◆ saturationFromCapillaryPressure()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::saturationFromCapillaryPressure ( int cell_index,
double cap_press ) const

Inverse of the capillary pressure function.

Parameters
cell_indexindex of a grid cell.
cap_pressa capillary pressure value.
Returns
the saturation at the given cell has the given capillary pressure.

◆ setDensities()

template<int dim, class RPImpl, class RockType>
void Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::setDensities ( double d1,
double d2 )

Set densitities of both faces.

Parameters
d1the densitity of the first (water) phase.
d2the densitity of the second (oil) phase.

◆ setViscosities()

template<int dim, class RPImpl, class RockType>
void Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::setViscosities ( double v1,
double v2 )

Set viscosities of both faces.

Parameters
v1the viscosity of the first (water) phase.
v2the viscosity of the second (oil) phase.

◆ sowcr()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::sowcr ( int cell_index) const

Read-access to sowcr.

Parameters
cell_indexindex of a grid cell.
Returns
sowcr value of the cell.

◆ swcr()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::swcr ( int cell_index) const

Read-access to swcr.

Parameters
cell_indexindex of a grid cell.
Returns
swcr value of the cell.

◆ viscosityFirstPhase()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::viscosityFirstPhase ( ) const

Viscosity of first (water) phase.

Returns
the viscosity value.

◆ viscositySecondPhase()

template<int dim, class RPImpl, class RockType>
double Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::viscositySecondPhase ( ) const

Viscosity of second (oil) phase.

Returns
the viscosity value.

◆ writeSintefLegacyFormat()

template<int dim, class RPImpl, class RockType>
void Opm::ReservoirPropertyCommon< dim, RPImpl, RockType >::writeSintefLegacyFormat ( const std::string & grid_prefix) const

Write permeability and porosity in the Sintef legacy format.

Parameters
grid_prefixthe prefix of all files output by this function.

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