#include <RockFromDeck.hpp>

Public Member Functions

 RockFromDeck ()
 Default constructor. More...
 
 RockFromDeck (std::size_t number_of_cells)
 
void init (const Opm::EclipseState &eclState, int number_of_cells, const int *global_cell, const int *cart_dims)
 
int numDimensions () const
 
int numCells () const
 
const double * porosity () const
 
const double * permeability () const
 

Static Public Member Functions

static void extractInterleavedPermeability (const Opm::EclipseState &eclState, const int number_of_cells, const int *global_cell, const int *cart_dims, const double perm_threshold, std::vector< double > &permeability)
 

Friends

class BlackoilPropsDataHandle
 

Constructor & Destructor Documentation

◆ RockFromDeck() [1/2]

Opm::RockFromDeck::RockFromDeck ( )

Default constructor.

◆ RockFromDeck() [2/2]

Opm::RockFromDeck::RockFromDeck ( std::size_t  number_of_cells)
explicit

Creates rock properties with zero porosity and permeability

Parameters
number_of_cellsThe number of cells

Member Function Documentation

◆ extractInterleavedPermeability()

static void Opm::RockFromDeck::extractInterleavedPermeability ( const Opm::EclipseState &  eclState,
const int  number_of_cells,
const int *  global_cell,
const int *  cart_dims,
const double  perm_threshold,
std::vector< double > &  permeability 
)
static

Convert the permeabilites for the logically Cartesian grid in EclipseState to an array of size number_of_cells*dim*dim for the compressed array.

Parameters
eclStateThe EclipseState (processed deck) produced by the opm-parser code
number_of_cellsThe number of cells in the grid.
global_cellThe mapping fom local to global cell indices. global_cell[i] is the corresponding global index of i.
cart_dimsThe size of the underlying cartesian grid.
perm_thresholdThe threshold for permeability
permeabilityThe result array

◆ init()

void Opm::RockFromDeck::init ( const Opm::EclipseState &  eclState,
int  number_of_cells,
const int *  global_cell,
const int *  cart_dims 
)

Initialize from deck and cell mapping.

Parameters
eclStateThe EclipseState (processed deck) produced by the opm-parser code
number_of_cellsThe number of cells in the grid.
global_cellThe mapping fom local to global cell indices. global_cell[i] is the corresponding global index of i.
cart_dimsThe size of the underlying cartesian grid.

◆ numCells()

int Opm::RockFromDeck::numCells ( ) const
inline
Returns
N, the number of cells.

◆ numDimensions()

int Opm::RockFromDeck::numDimensions ( ) const
inline
Returns
D, the number of spatial dimensions. Always 3 for deck input.

◆ permeability()

const double * Opm::RockFromDeck::permeability ( ) const
inline
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).

◆ porosity()

const double * Opm::RockFromDeck::porosity ( ) const
inline
Returns
Array of N porosity values.

Friends And Related Function Documentation

◆ BlackoilPropsDataHandle

friend class BlackoilPropsDataHandle
friend

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