Encapsulates the cfs_tpfa (= compressible flow solver two-point flux approximation) solver modules.
More...
#include <TpfaCompressibleAssembler.hpp>
|
| 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 |
|
Encapsulates the cfs_tpfa (= compressible flow solver two-point flux approximation) solver modules.
◆ FlowBCTypes
Boundary condition types.
Enumerator |
---|
FBC_UNSET | |
FBC_PRESSURE | |
FBC_FLUX | |
◆ TpfaCompressibleAssembler()
TpfaCompressibleAssembler::TpfaCompressibleAssembler |
( |
| ) |
|
|
inline |
Default constructor, does nothing.
◆ ~TpfaCompressibleAssembler()
TpfaCompressibleAssembler::~TpfaCompressibleAssembler |
( |
| ) |
|
|
inline |
◆ assemble()
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
-
◆ computePressuresAndFluxes()
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_pressures | Cell pressure values. |
[out] | face_areas | Face flux values. |
◆ explicitTimestepLimit()
double TpfaCompressibleAssembler::explicitTimestepLimit |
( |
const double * |
faceA, |
|
|
const double * |
phasemobf, |
|
|
const double * |
phasemobf_deriv, |
|
|
const double * |
surf_dens |
|
) |
| |
|
inline |
◆ explicitTransport()
void TpfaCompressibleAssembler::explicitTransport |
( |
const double |
dt, |
|
|
double * |
cell_surfvols |
|
) |
| |
|
inline |
◆ faceFluxToCellFlux()
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_areas | Face flux values (usually output from computePressuresAndFluxes()). |
[out] | cell_fluxes | Cell-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. |
◆ faceTransmissibilities()
const std::vector< double > & TpfaCompressibleAssembler::faceTransmissibilities |
( |
| ) |
const |
|
inline |
◆ init() [1/2]
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
-
Grid | This must conform to the SimpleGrid concept. |
- Parameters
-
grid | The grid object. |
perm | Permeability. It should contain dim*dim entries (a full tensor) for each cell. |
gravity | Array containing gravity acceleration vector. It should contain dim entries. |
◆ init() [2/2]
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
-
Grid | This must conform to the SimpleGrid concept. |
Wells | This must conform to the SimpleWells concept. |
- Parameters
-
grid | The grid object. |
wells | Well specifications. |
perm | Permeability. It should contain dim*dim entries (a full tensor) for each cell. |
perm | Porosity by cell. |
gravity | Array containing gravity acceleration vector. It should contain dim entries. |
References init().
Referenced by init(), and Opm::TpfaCompressible< GridInterface, RockInterface, FluidInterface, WellsInterface, BCInterface >::setup().
◆ linearSystem()
void TpfaCompressibleAssembler::linearSystem |
( |
LinearSystem & |
s | ) |
|
|
inline |
◆ numCellFaces()
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:
|