fsh.h File Reference
Include dependency graph for fsh.h:

Go to the source code of this file.

Classes

struct  fsh_data
 

Functions

void fsh_destroy (struct fsh_data *h)
 
struct fsh_datacfsh_construct (struct UnstructuredGrid *G, well_t *W)
 
void cfsh_assemble (struct FlowBoundaryConditions *bc, const double *src, const double *Binv, const double *Biv, const double *P, const double *gpress, well_control_t *wctrl, const double *WI, const double *BivW, const double *wdp, struct fsh_data *h)
 
struct fsh_dataifsh_construct (struct UnstructuredGrid *G, well_t *W)
 
void ifsh_assemble (struct FlowBoundaryConditions *bc, const double *src, const double *Binv, const double *gpress, well_control_t *wctrl, const double *WI, const double *wdp, struct fsh_data *h)
 
void fsh_press_flux (struct UnstructuredGrid *G, const double *Binv, const double *gpress, struct fsh_data *h, double *cpress, double *fflux, double *wpress, double *wflux)
 

Detailed Description

Routines and data structures to support the construction and formation of hybridized pressure solvers based on Schur complement reductions.

Pressure solvers based on this strategy will be structured according to the following scheme

  1. Construct FSH data object suitable for the particular problem using either of the functions cfsh_construct() or ifsh_construct() for compressible or incompressible flow, respectively
  2. Compute static discretisation quantities, for instance using functions mim_ip_simple_all() and mim_ip_compute_gpress()
  3. For each time step or non-linear iteration:
    1. Compute dynamic discretisation quantities incorporating effects of multiple phases and non-linear feedback.
    2. Assemble system of simultaneous linear equations using functions cfsh_assemble() or ifsh_assemble()
    3. Solve the resulting system of linear equations, available in the A and b objects of the FSH data object, using some linear solver software. The solution should be placed in the x object of the FSH object.
    4. Reconstruct derived quantities such as cell pressures and interface fluxes using function fsh_press_flux(). Function fsh_press_flux() relies on the solution to the linear system being stored in the x object.
  4. Release resources using function fsh_destroy() at end of simulation.

Function Documentation

void cfsh_assemble ( struct FlowBoundaryConditions bc,
const double *  src,
const double *  Binv,
const double *  Biv,
const double *  P,
const double *  gpress,
well_control_t wctrl,
const double *  WI,
const double *  BivW,
const double *  wdp,
struct fsh_data h 
)

Form Schur-complement system of simultaneous linear equations arising in compressible flow using a hybridized formulation.

Upon returning from function cfsh_assemble(), the resulting system of simultaneous linear equations is stored in h->A and h->b.

Parameters
[in]bcBoundary conditions.
[in]srcExplicit source terms.
[in]BinvInverse of block-diagonal matrix $B$ Typically computed using functions mim_ip_simple_all() and mim_ip_mobility_update().
[in]Biv$B^{-1}v$.
[in]PCompressible accumulation term.
[in]gpressGravity pressure.
[in]wctrlWell controls. Ignored.
[in]WIWell indices. Ignored.
[in]BivW$B^{-1}v$ for wells. Ignored.
[in]wdpGravity pressure along well track. Ignored.
[in,out]hData object.
struct fsh_data* cfsh_construct ( struct UnstructuredGrid G,
well_t W 
)

Construct compressible hybrid flow-solver data object for a given grid and well pattern.

Parameters
[in]GGrid.
[in]WWell topology. Ignored.
Returns
Fully formed data object suitable for use in a compressible pressure solver. NULL in case of construction failure.
void fsh_destroy ( struct fsh_data h)

Dispose of all memory associated to FSH object.

Parameters
[in,out]hFSH object. Following a call to function fsh_destroy(), the pointer is invalid.
void fsh_press_flux ( struct UnstructuredGrid G,
const double *  Binv,
const double *  gpress,
struct fsh_data h,
double *  cpress,
double *  fflux,
double *  wpress,
double *  wflux 
)

Compute cell pressures (cpress) and interface fluxes (fflux) from current solution of system of linear equations, h->x. Back substitution process, projected half-contact fluxes.

Parameters
[in]GGrid.
[in]BinvInverse of block-diagonal matrix $B$ Must coincide with the matrix used to form the system of linear equations currently stored in the data object.
[in]gpressGravity pressure. Must coincide with the array used to form the system of linear equations.
[in]hData object.
[out]cpressCell pressure.
[out]ffluxInterface fluxes.
[out]wpressWell pressure.
[out]wfluxWell perforation fluxes.
void ifsh_assemble ( struct FlowBoundaryConditions bc,
const double *  src,
const double *  Binv,
const double *  gpress,
well_control_t wctrl,
const double *  WI,
const double *  wdp,
struct fsh_data h 
)

Form Schur-complement system of simultaneous linear equations arising in compressible flow using a hybridized formulation.

Upon returning from function cfsh_assemble(), the resulting system of simultaneous linear equations is stored in h->A and h->b.

Parameters
[in]bcBoundary conditions.
[in]srcExplicit source terms.
[in]BinvInverse of block-diagonal matrix $B$ Typically computed using functions mim_ip_simple_all() and mim_ip_mobility_update().
[in]gpressGravity pressure.
[in]wctrlWell controls.
[in]WIWell indices.
[in]wdpGravity pressure along well track.
[in,out]hData object.
struct fsh_data* ifsh_construct ( struct UnstructuredGrid G,
well_t W 
)

Construct incompressible hybrid flow-solver data object for a given grid and well pattern.

Parameters
GGrid.
WWell topology.
Returns
Fully formed data object suitable for use in an incompressible pressure solver. NULL in case of construction failure.