Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface > Class Template Reference

#include <TpfaCompressible.hpp>

Public Types

enum  ReturnCode { SolveOk, VolumeDiscrepancyTooLarge, FailedToConverge }
 
typedef TpfaCompressibleAssembler PressureAssembler
 

Public Member Functions

 TpfaCompressible ()
 Default constructor. Does nothing. More...
 
void init (const Opm::parameter::ParameterGroup &param)
 Initializes run-time parameters of the solver. More...
 
FluidInterface::CompVec inflowMixture () const
 Accessor for the inflow mixture. More...
 
void setup (const GridInterface &grid, const RockInterface &rock, const FluidInterface &fluid, const WellsInterface &wells, const typename GridInterface::Vector &grav, const BCInterface &bc, const std::vector< typename FluidInterface::PhaseVec > *face_pressure=0)
 Setup routine, does grid/rock-dependent initialization. More...
 
double volumeDiscrepancyLimit () const
 
const std::vector< double > & faceTransmissibilities ()
 
bool volumeDiscrepancyAcceptable (const std::vector< typename FluidInterface::PhaseVec > &cell_pressure, const std::vector< typename FluidInterface::PhaseVec > &face_pressure, const std::vector< double > &well_perf_pressure, const std::vector< typename FluidInterface::CompVec > &cell_z, const double dt)
 
ReturnCode solve (std::vector< typename FluidInterface::PhaseVec > &cell_pressure, std::vector< typename FluidInterface::PhaseVec > &face_pressure, const std::vector< typename FluidInterface::CompVec > &cell_z, std::vector< double > &face_flux, std::vector< double > &well_bhp_pressures, std::vector< double > &well_perf_pressures, std::vector< double > &well_perf_fluxes, const std::vector< double > &src, const double dt)
 Construct and solve system of linear equations for the phase pressure values on cells and faces, also compute total face fluxes. More...
 
double stableStepIMPES ()
 Call this function after solve(). More...
 
void doStepIMPES (std::vector< typename FluidInterface::CompVec > &cell_z, const double dt)
 Do an IMPES step using the facilities of the pressure solver. More...
 

Member Typedef Documentation

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
typedef TpfaCompressibleAssembler Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::PressureAssembler

Member Enumeration Documentation

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
enum Opm::TpfaCompressible::ReturnCode
Enumerator
SolveOk 
VolumeDiscrepancyTooLarge 
FailedToConverge 

Constructor & Destructor Documentation

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::TpfaCompressible ( )
inline

Default constructor. Does nothing.

Member Function Documentation

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
void Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::doStepIMPES ( std::vector< typename FluidInterface::CompVec > &  cell_z,
const double  dt 
)
inline

Do an IMPES step using the facilities of the pressure solver.

References TpfaCompressibleAssembler::explicitTransport().

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
const std::vector<double>& Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::faceTransmissibilities ( )
inline
template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
FluidInterface::CompVec Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::inflowMixture ( ) const
inline

Accessor for the inflow mixture.

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
void Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::init ( const Opm::parameter::ParameterGroup &  param)
inline

Initializes run-time parameters of the solver.

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
void Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::setup ( const GridInterface &  grid,
const RockInterface &  rock,
const FluidInterface &  fluid,
const WellsInterface &  wells,
const typename GridInterface::Vector &  grav,
const BCInterface &  bc,
const std::vector< typename FluidInterface::PhaseVec > *  face_pressure = 0 
)
inline

Setup routine, does grid/rock-dependent initialization.

Parameters
[in]gridThe grid.
[in]rockThe cell-wise permeabilities and porosities.
[in]fluidFluid properties.
[in]wellsWell specifications.
[in]gravGravity vector. Its Euclidian two-norm value represents the strength of the gravity field (in units of m/s^2) while its direction is the direction of gravity in the current model.
[in]bcBoundary conditions.
[in]face_pressureFace pressures. Allow optional micromanagement of pressure boundary conditions.

References TpfaCompressibleAssembler::FBC_FLUX, TpfaCompressibleAssembler::FBC_PRESSURE, TpfaCompressibleAssembler::FBC_UNSET, TpfaCompressibleAssembler::init(), Opm::BCBase< T >::isDirichlet(), Opm::BCBase< T >::isNeumann(), Opm::FlowBC::outflux(), and Opm::FlowBC::pressure().

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
ReturnCode Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::solve ( std::vector< typename FluidInterface::PhaseVec > &  cell_pressure,
std::vector< typename FluidInterface::PhaseVec > &  face_pressure,
const std::vector< typename FluidInterface::CompVec > &  cell_z,
std::vector< double > &  face_flux,
std::vector< double > &  well_bhp_pressures,
std::vector< double > &  well_perf_pressures,
std::vector< double > &  well_perf_fluxes,
const std::vector< double > &  src,
const double  dt 
)
inline

Construct and solve system of linear equations for the phase pressure values on cells and faces, also compute total face fluxes.

Parameters
[in,out]cell_pressurePhase pressures per cell.
[in,out]face_pressurePhase pressures per face.
[in,out]cell_zSurface volume per cell. Only changed if the
transport
argument is true.
[out]face_fluxTotal (summed over all phases) volume flux (signed) across each face.
[out]well_perf_pressuresPressure in each well perforation.
[out]well_perf_fluxesTotal (summed over all phases) volume flux (signed, positive meaning injection) from each well perforation.
[in]srcExplicit source terms. One scalar value for each grid cell representing the rate (in units of m^3/s) of fluid being injected into (>0) or extracted from (<0) a given grid cell.
[in]dtTimestep for pressure solver.
[in]transportIf true, modify
cell_z
by IMPES scheme.

References Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::SolveOk.

template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
double Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::stableStepIMPES ( )
inline
template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
bool Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::volumeDiscrepancyAcceptable ( const std::vector< typename FluidInterface::PhaseVec > &  cell_pressure,
const std::vector< typename FluidInterface::PhaseVec > &  face_pressure,
const std::vector< double > &  well_perf_pressure,
const std::vector< typename FluidInterface::CompVec > &  cell_z,
const double  dt 
)
inline
template<class GridInterface , class RockInterface , class FluidInterface , class WellsInterface , class BCInterface >
double Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::volumeDiscrepancyLimit ( ) const
inline

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