Opm::BlackoilPropsAdFromDeck Class Reference

#include <BlackoilPropsAdFromDeck.hpp>

Inheritance diagram for Opm::BlackoilPropsAdFromDeck:
Inheritance graph

Public Types

typedef
SaturationPropsFromDeck::MaterialLawManager 
MaterialLawManager
 
typedef AutoDiffBlock< double > ADB
 
typedef ADB::V V
 
typedef std::vector< int > Cells
 
enum  PhaseIndex {
  Water = BlackoilPhases::Aqua, Oil = BlackoilPhases::Liquid, Gas = BlackoilPhases::Vapour, Aqua = BlackoilPhases::Aqua,
  Liquid = BlackoilPhases::Liquid, Vapour = BlackoilPhases::Vapour, MaxNumPhases = BlackoilPhases::MaxNumPhases
}
 Canonical named indices for each phase. More...
 
typedef ADB::M M
 

Public Member Functions

 BlackoilPropsAdFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, std::shared_ptr< MaterialLawManager > materialLawManager, const UnstructuredGrid &grid, const bool init_rock=true)
 
 BlackoilPropsAdFromDeck (Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, const UnstructuredGrid &grid, const bool init_rock=true)
 
 BlackoilPropsAdFromDeck (const BlackoilPropsAdFromDeck &props, const int number_of_cells)
 Constructor to create properties for a subgrid. More...
 
int numDimensions () const
 
int numCells () const
 
virtual const int * cellPvtRegionIndex () const
 
const double * porosity () const
 
const double * permeability () const
 
int numPhases () const
 
PhaseUsage phaseUsage () const
 
const double * surfaceDensity (const int cellIdx=0) const
 
ADB muWat (const ADB &pw, const ADB &T, const Cells &cells) const
 
ADB muOil (const ADB &po, const ADB &T, const ADB &rs, const std::vector< PhasePresence > &cond, const Cells &cells) const
 
ADB muGas (const ADB &pg, const ADB &T, const ADB &rv, const std::vector< PhasePresence > &cond, const Cells &cells) const
 
ADB bWat (const ADB &pw, const ADB &T, const Cells &cells) const
 
ADB bOil (const ADB &po, const ADB &T, const ADB &rs, const std::vector< PhasePresence > &cond, const Cells &cells) const
 
ADB bGas (const ADB &pg, const ADB &T, const ADB &rv, const std::vector< PhasePresence > &cond, const Cells &cells) const
 
V rsSat (const V &po, const Cells &cells) const
 
V rsSat (const V &po, const V &so, const Cells &cells) const
 
ADB rsSat (const ADB &po, const Cells &cells) const
 
ADB rsSat (const ADB &po, const ADB &so, const Cells &cells) const
 
ADB rvSat (const ADB &po, const Cells &cells) const
 
ADB rvSat (const ADB &po, const ADB &so, const Cells &cells) const
 
std::vector< ADBrelperm (const ADB &sw, const ADB &so, const ADB &sg, const Cells &cells) const
 
std::vector< ADBcapPress (const ADB &sw, const ADB &so, const ADB &sg, const Cells &cells) const
 
void updateSatHyst (const std::vector< double > &saturation, const std::vector< int > &cells)
 
void updateSatOilMax (const std::vector< double > &saturation)
 Update for max oil saturation. More...
 
void setSwatInitScaling (const std::vector< double > &saturation, const std::vector< double > &pc)
 

Friends

class BlackoilPropsDataHandle
 

Detailed Description

This class implements the AD-adapted fluid interface for three-phase black-oil. It requires an input deck from which it reads all relevant property data.

Most methods are available in two overloaded versions, one taking a constant vector and returning the same, and one taking an AD type and returning the same. Derivatives are not returned separately by any method, only implicitly with the AD version of the methods.

Member Typedef Documentation

typedef std::vector<int> Opm::BlackoilPropsAdFromDeck::Cells
typedef SaturationPropsFromDeck::MaterialLawManager Opm::BlackoilPropsAdFromDeck::MaterialLawManager

Member Enumeration Documentation

Canonical named indices for each phase.

Enumerator
Water 
Oil 
Gas 
Aqua 
Liquid 
Vapour 
MaxNumPhases 

Constructor & Destructor Documentation

Opm::BlackoilPropsAdFromDeck::BlackoilPropsAdFromDeck ( Opm::DeckConstPtr  deck,
Opm::EclipseStateConstPtr  eclState,
std::shared_ptr< MaterialLawManager materialLawManager,
const UnstructuredGrid &  grid,
const bool  init_rock = true 
)

Constructor to create a blackoil properties from an ECL deck.

The materialLawManager parameter represents the object from opm-material which handles the creating and updating parameter objects for the capillary pressure/relperm relations for each grid cell. This object is created internally for the constructors below, but if it is already available externally some performance can be gained by creating it only once.

Parameters
deckThe unprocessed ECL deck from opm-parser
eclStateThe processed ECL deck from opm-parser
materialLawManagerThe container for the material law parameter objects
gridThe grid upon which the simulation is run on.
init_rockIf true the rock properties (rock compressibility and reference pressure) are read from the deck
Opm::BlackoilPropsAdFromDeck::BlackoilPropsAdFromDeck ( Opm::DeckConstPtr  deck,
Opm::EclipseStateConstPtr  eclState,
const UnstructuredGrid &  grid,
const bool  init_rock = true 
)

Constructor to create a blackoil properties from an ECL deck.

Parameters
deckThe unprocessed ECL deck from opm-parser
eclStateThe processed ECL deck from opm-parser
gridThe grid upon which the simulation is run on.
init_rockIf true the rock properties (rock compressibility and reference pressure) are read from the deck
Opm::BlackoilPropsAdFromDeck::BlackoilPropsAdFromDeck ( const BlackoilPropsAdFromDeck props,
const int  number_of_cells 
)

Constructor to create properties for a subgrid.

This copies all properties that are not dependant on the grid size from an existing properties object and the number of cells. All properties that do not depend on the grid dimension will be copied. For the rest will have the correct size but the values will be undefined.

Parameters
propsThe property object to copy from. number_of_cells The number of cells of the subgrid.

Member Function Documentation

ADB Opm::BlackoilPropsAdFromDeck::bGas ( const ADB pg,
const ADB T,
const ADB rv,
const std::vector< PhasePresence > &  cond,
const Cells cells 
) const
virtual

Gas formation volume factor.

Parameters
[in]pgArray of n gas pressure values.
[in]TArray of n temperature values.
[in]rvArray of n vapor oil/gas ratio
[in]condArray of n objects, each specifying which phases are present with non-zero saturation in a cell.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n formation volume factor values.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::bOil ( const ADB po,
const ADB T,
const ADB rs,
const std::vector< PhasePresence > &  cond,
const Cells cells 
) const
virtual

Oil formation volume factor.

Parameters
[in]poArray of n oil pressure values.
[in]TArray of n temperature values.
[in]rsArray of n gas solution factor values.
[in]condArray of n objects, each specifying which phases are present with non-zero saturation in a cell.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n formation volume factor values.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::bWat ( const ADB pw,
const ADB T,
const Cells cells 
) const
virtual

Water formation volume factor.

Parameters
[in]pwArray of n water pressure values.
[in]TArray of n temperature values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n formation volume factor values.

Implements Opm::BlackoilPropsAdInterface.

std::vector<ADB> Opm::BlackoilPropsAdFromDeck::capPress ( const ADB sw,
const ADB so,
const ADB sg,
const Cells cells 
) const
virtual

Capillary pressure for all phases.

Parameters
[in]swArray of n water saturation values.
[in]soArray of n oil saturation values.
[in]sgArray of n gas saturation values.
[in]cellsArray of n cell indices to be associated with the saturation values.
Returns
An std::vector with 3 elements, each an array of n capillary pressure values, containing the offsets for each p_g, p_o, p_w. The capillary pressure between two arbitrary phases alpha and beta is then given as p_alpha - p_beta.

Implements Opm::BlackoilPropsAdInterface.

virtual const int* Opm::BlackoilPropsAdFromDeck::cellPvtRegionIndex ( ) const
inlinevirtual

Return an array containing the PVT table index for each grid cell

ADB Opm::BlackoilPropsAdFromDeck::muGas ( const ADB pg,
const ADB T,
const ADB rv,
const std::vector< PhasePresence > &  cond,
const Cells cells 
) const
virtual

Gas viscosity.

Parameters
[in]pgArray of n gas pressure values.
[in]TArray of n temperature values.
[in]rvArray of n vapor oil/gas ratios.
[in]condArray of n objects, each specifying which phases are present with non-zero saturation in a cell.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n viscosity values.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::muOil ( const ADB po,
const ADB T,
const ADB rs,
const std::vector< PhasePresence > &  cond,
const Cells cells 
) const
virtual

Oil viscosity.

Parameters
[in]poArray of n oil pressure values.
[in]TArray of n temperature values.
[in]rsArray of n gas solution factor values.
[in]condArray of n objects, each specifying which phases are present with non-zero saturation in a cell.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n viscosity values.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::muWat ( const ADB pw,
const ADB T,
const Cells cells 
) const
virtual

Water viscosity.

Parameters
[in]pwArray of n water pressure values.
[in]TArray of n temperature values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n viscosity values.

Implements Opm::BlackoilPropsAdInterface.

int Opm::BlackoilPropsAdFromDeck::numCells ( ) const
virtual
Returns
N, the number of cells.

Implements Opm::BlackoilPropsAdInterface.

int Opm::BlackoilPropsAdFromDeck::numDimensions ( ) const
virtual
Returns
D, the number of spatial dimensions.

Implements Opm::BlackoilPropsAdInterface.

int Opm::BlackoilPropsAdFromDeck::numPhases ( ) const
virtual
Returns
Number of active phases (also the number of components).

Implements Opm::BlackoilPropsAdInterface.

const double* Opm::BlackoilPropsAdFromDeck::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::BlackoilPropsAdInterface.

PhaseUsage Opm::BlackoilPropsAdFromDeck::phaseUsage ( ) const
virtual
Returns
Object describing the active phases.

Implements Opm::BlackoilPropsAdInterface.

const double* Opm::BlackoilPropsAdFromDeck::porosity ( ) const
virtual
Returns
Array of N porosity values.

Implements Opm::BlackoilPropsAdInterface.

std::vector<ADB> Opm::BlackoilPropsAdFromDeck::relperm ( const ADB sw,
const ADB so,
const ADB sg,
const Cells cells 
) const
virtual

Relative permeabilities for all phases.

Parameters
[in]swArray of n water saturation values.
[in]soArray of n oil saturation values.
[in]sgArray of n gas saturation values.
[in]cellsArray of n cell indices to be associated with the saturation values.
Returns
An std::vector with 3 elements, each an array of n relperm values, containing krw, kro, krg. Use PhaseIndex for indexing into the result.

Implements Opm::BlackoilPropsAdInterface.

V Opm::BlackoilPropsAdFromDeck::rsSat ( const V po,
const Cells cells 
) const

Bubble point curve for Rs as function of oil pressure.

Parameters
[in]poArray of n oil pressure values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n bubble point values for Rs.
V Opm::BlackoilPropsAdFromDeck::rsSat ( const V po,
const V so,
const Cells cells 
) const

Bubble point curve for Rs as function of oil pressure.

Parameters
[in]poArray of n oil pressure values.
[in]soArray of n oil saturation values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n bubble point values for Rs.
ADB Opm::BlackoilPropsAdFromDeck::rsSat ( const ADB po,
const Cells cells 
) const
virtual

Bubble point curve for Rs as function of oil pressure.

Parameters
[in]poArray of n oil pressure values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n bubble point values for Rs.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::rsSat ( const ADB po,
const ADB so,
const Cells cells 
) const
virtual

Bubble point curve for Rs as function of oil pressure.

Parameters
[in]poArray of n oil pressure values.
[in]soArray of n oil saturation values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n bubble point values for Rs.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::rvSat ( const ADB po,
const Cells cells 
) const
virtual

Condensation curve for Rv as function of oil pressure.

Parameters
[in]poArray of n oil pressure values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n condensation point values for Rv.

Implements Opm::BlackoilPropsAdInterface.

ADB Opm::BlackoilPropsAdFromDeck::rvSat ( const ADB po,
const ADB so,
const Cells cells 
) const
virtual

Condensation curve for Rv as function of oil pressure.

Parameters
[in]poArray of n oil pressure values.
[in]soArray of n oil saturation values.
[in]cellsArray of n cell indices to be associated with the pressure values.
Returns
Array of n condensation point values for Rv.

Implements Opm::BlackoilPropsAdInterface.

void Opm::BlackoilPropsAdFromDeck::setSwatInitScaling ( const std::vector< double > &  saturation,
const std::vector< double > &  pc 
)

Set capillary pressure scaling according to pressure diff. and initial water saturation.

Parameters
[in]saturationArray of n*numPhases saturation values.
[in]pcArray of n*numPhases capillary pressure values.
const double* Opm::BlackoilPropsAdFromDeck::surfaceDensity ( const int  cellIdx = 0) const
virtual

Densities of stock components at surface conditions.

Returns
Array of 3 density values.

Implements Opm::BlackoilPropsAdInterface.

void Opm::BlackoilPropsAdFromDeck::updateSatHyst ( const std::vector< double > &  saturation,
const std::vector< int > &  cells 
)
virtual

Saturation update for hysteresis behavior.

Parameters
[in]cellsArray of n cell indices to be associated with the saturation values.

Implements Opm::BlackoilPropsAdInterface.

void Opm::BlackoilPropsAdFromDeck::updateSatOilMax ( const std::vector< double > &  saturation)
virtual

Update for max oil saturation.

Implements Opm::BlackoilPropsAdInterface.

Friends And Related Function Documentation

friend class BlackoilPropsDataHandle
friend

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