TpfaCompressibleAssembler Class Reference

Encapsulates the cfs_tpfa (= compressible flow solver two-point flux approximation) solver modules. More...

#include <TpfaCompressibleAssembler.hpp>

Classes

struct  LinearSystem
 Encapsulate a sparse linear system in CSR format. More...
 

Public Types

enum  FlowBCTypes { FBC_UNSET = BC_NOFLOW, FBC_PRESSURE = BC_PRESSURE, FBC_FLUX = BC_FLUX_TOTVOL }
 Boundary condition types. More...
 

Public Member Functions

 TpfaCompressibleAssembler ()
 Default constructor, does nothing. More...
 
 ~TpfaCompressibleAssembler ()
 Destructor. More...
 
template<class Grid , class Wells >
void init (const Grid &grid, const Wells &wells, const double *perm, const double *porosity, const typename Grid::Vector &gravity)
 Initialize the solver's structures for a given grid, for well setup also call initWells(). More...
 
template<class Grid >
void init (const Grid &grid, const double *perm, const double *porosity, const typename Grid::Vector &gravity)
 Initialize the solver's structures for a given grid, for well setup also call initWells(). More...
 
void assemble (const double *sources, const FlowBCTypes *bctypes, const double *bcvalues, const double dt, const double *totcompr, const double *voldiscr, const double *cellA, const double *faceA, const double *wellperfA, const double *phasemobf, const double *phasemobwellperf, const double *cell_pressure, const double *gravcapf, const double *wellperf_gpot, const double *surf_dens)
 Assemble the sparse system. You must call init() prior to calling assemble(). More...
 
void linearSystem (LinearSystem &s)
 Access the linear system assembled. You must call assemble() prior to calling linearSystem(). More...
 
void computePressuresAndFluxes (std::vector< double > &cell_pressures, std::vector< double > &face_pressures, std::vector< double > &face_fluxes, std::vector< double > &well_pressures, std::vector< double > &well_fluxes)
 Compute cell pressures and face fluxes. You must call assemble() (and solve the linear system accessed by calling linearSystem()) prior to calling computePressuresAndFluxes(). More...
 
double explicitTimestepLimit (const double *faceA, const double *phasemobf, const double *phasemobf_deriv, const double *surf_dens)
 Explicit IMPES time step limit. More...
 
void explicitTransport (const double dt, double *cell_surfvols)
 Explicit IMPES transport. More...
 
void faceFluxToCellFlux (const std::vector< double > &face_fluxes, std::vector< double > &cell_fluxes)
 Compute cell fluxes from face fluxes. You must call assemble() (and solve the linear system accessed by calling linearSystem()) prior to calling faceFluxToCellFlux(). More...
 
const std::vector< int > & numCellFaces () const
 Access the number of connections (faces) per cell. Deprecated, will be removed. More...
 
const std::vector< double > & faceTransmissibilities () const
 

Detailed Description

Encapsulates the cfs_tpfa (= compressible flow solver two-point flux approximation) solver modules.

Member Enumeration Documentation

Boundary condition types.

Enumerator
FBC_UNSET 
FBC_PRESSURE 
FBC_FLUX 

Constructor & Destructor Documentation

TpfaCompressibleAssembler::TpfaCompressibleAssembler ( )
inline

Default constructor, does nothing.

TpfaCompressibleAssembler::~TpfaCompressibleAssembler ( )
inline

Destructor.

Member Function Documentation

void TpfaCompressibleAssembler::assemble ( const double *  sources,
const FlowBCTypes bctypes,
const double *  bcvalues,
const double  dt,
const double *  totcompr,
const double *  voldiscr,
const double *  cellA,
const double *  faceA,
const double *  wellperfA,
const double *  phasemobf,
const double *  phasemobwellperf,
const double *  cell_pressure,
const double *  gravcapf,
const double *  wellperf_gpot,
const double *  surf_dens 
)
inline

Assemble the sparse system. You must call init() prior to calling assemble().

Parameters
sourcesSource terms, one per cell. Positive numbers are sources, negative are sinks.
total_mobilitiesScalar total mobilities, one per cell.
omegasGravity term, one per cell. In a multi-phase flow setting this is equal to

\[ \omega = \sum_{p} \frac{\lambda_p}{\lambda_t} \rho_p \]

where $\lambda_p$ is a phase mobility, $\rho_p$ is a phase density and $\lambda_t$ is the total mobility.
void TpfaCompressibleAssembler::computePressuresAndFluxes ( std::vector< double > &  cell_pressures,
std::vector< double > &  face_pressures,
std::vector< double > &  face_fluxes,
std::vector< double > &  well_pressures,
std::vector< double > &  well_fluxes 
)
inline

Compute cell pressures and face fluxes. You must call assemble() (and solve the linear system accessed by calling linearSystem()) prior to calling computePressuresAndFluxes().

Parameters
[out]cell_pressuresCell pressure values.
[out]face_areasFace flux values.
double TpfaCompressibleAssembler::explicitTimestepLimit ( const double *  faceA,
const double *  phasemobf,
const double *  phasemobf_deriv,
const double *  surf_dens 
)
inline
void TpfaCompressibleAssembler::explicitTransport ( const double  dt,
double *  cell_surfvols 
)
inline
void TpfaCompressibleAssembler::faceFluxToCellFlux ( const std::vector< double > &  face_fluxes,
std::vector< double > &  cell_fluxes 
)
inline

Compute cell fluxes from face fluxes. You must call assemble() (and solve the linear system accessed by calling linearSystem()) prior to calling faceFluxToCellFlux().

Parameters
face_fluxes
face_areasFace flux values (usually output from computePressuresAndFluxes()).
[out]cell_fluxesCell-wise flux values. They are given in cell order, and for each cell there is one value for each adjacent face (in the same order as the cell-face topology of the grid). Positive values represent fluxes out of the cell.
const std::vector<double>& TpfaCompressibleAssembler::faceTransmissibilities ( ) const
inline
template<class Grid , class Wells >
void TpfaCompressibleAssembler::init ( const Grid &  grid,
const Wells &  wells,
const double *  perm,
const double *  porosity,
const typename Grid::Vector &  gravity 
)
inline

Initialize the solver's structures for a given grid, for well setup also call initWells().

Template Parameters
GridThis must conform to the SimpleGrid concept.
WellsThis must conform to the SimpleWells concept.
Parameters
gridThe grid object.
wellsWell specifications.
permPermeability. It should contain dim*dim entries (a full tensor) for each cell.
permPorosity by cell.
gravityArray containing gravity acceleration vector. It should contain dim entries.

Referenced by Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::setup().

template<class Grid >
void TpfaCompressibleAssembler::init ( const Grid &  grid,
const double *  perm,
const double *  porosity,
const typename Grid::Vector &  gravity 
)
inline

Initialize the solver's structures for a given grid, for well setup also call initWells().

Template Parameters
GridThis must conform to the SimpleGrid concept.
Parameters
gridThe grid object.
permPermeability. It should contain dim*dim entries (a full tensor) for each cell.
gravityArray containing gravity acceleration vector. It should contain dim entries.
void TpfaCompressibleAssembler::linearSystem ( LinearSystem s)
inline

Access the linear system assembled. You must call assemble() prior to calling linearSystem().

Parameters
[out]sThe linear system encapsulation to modify. After this call, s will point to linear system structures that are owned and allocated internally.

References TpfaCompressibleAssembler::LinearSystem::b, TpfaCompressibleAssembler::LinearSystem::ia, TpfaCompressibleAssembler::LinearSystem::ja, TpfaCompressibleAssembler::LinearSystem::n, TpfaCompressibleAssembler::LinearSystem::nnz, TpfaCompressibleAssembler::LinearSystem::sa, and TpfaCompressibleAssembler::LinearSystem::x.

const std::vector<int>& TpfaCompressibleAssembler::numCellFaces ( ) const
inline

Access the number of connections (faces) per cell. Deprecated, will be removed.


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