Document public interface of common FSH module.

This commit is contained in:
Bård Skaflestad
2012-06-24 16:50:16 +02:00
parent 14f9faa011
commit bfaf6e0b87

View File

@@ -20,6 +20,37 @@
#ifndef OPM_FSH_HEADER_INCLUDED #ifndef OPM_FSH_HEADER_INCLUDED
#define OPM_FHS_HEADER_INCLUDED #define OPM_FHS_HEADER_INCLUDED
/**
* \file
* 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
* -# Construct @c FSH data object suitable for the particular
* problem using either of the functions cfsh_construct() or
* ifsh_construct() for compressible or incompressible flow,
* respectively
* -# Compute static discretisation quantities, for instance
* using functions mim_ip_simple_all() and mim_ip_compute_gpress()
* -# For each time step or non-linear iteration:
* -# Compute dynamic discretisation quantities incorporating
* effects of multiple phases and non-linear feedback.
* -# Assemble system of simultaneous linear equations using
* functions cfsh_assemble() or ifsh_assemble()
* -# Solve the resulting system of linear equations, available
* in the @c A and @c b objects of the @c FSH data object,
* using some linear solver software. The solution should
* be placed in the @c x object of the @c FSH object.
* -# 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 @c x object.
* -# Release resources using function fsh_destroy() at end of
* simulation.
*/
#include <opm/core/grid.h> #include <opm/core/grid.h>
#include <opm/core/well.h> #include <opm/core/well.h>
#include <opm/core/pressure/flow_bc.h> #include <opm/core/pressure/flow_bc.h>
@@ -29,54 +60,93 @@ extern "C" {
#endif #endif
/***************************************************************/
/* Data type common to compressible and incompressible solver. */
/***************************************************************/
struct CSRMatrix; struct CSRMatrix;
struct fsh_impl; struct fsh_impl;
/** Contains the linear system for assembly, as well as internal data /**
* for the assembly routines. * Main data structure of hybridized pressure solvers based on Schur
* complement reductions. Mainly intended to present a common view
* of a Schur complement system of simultaneous linear equations
* and to hold the solution of said system.
*/ */
struct fsh_data { struct fsh_data {
/* Let \f$n_i\f$ be the number of connections/faces of grid cell /**
* number \f$i\f$. Then max_ngconn = \f$\max_i n_i\f$ * Maximum number of connections in any grid cell,
* \f[
* \mathit{max\_ngconn} = \max_c \{ n_c \}
* \f]
* in which \f$n_c\f$ denotes the number connections
* (i.e., faces) of cell \f$c\f$.
*/ */
int max_ngconn; int max_ngconn;
/* With n_i as above, sum_ngconn2 = \f$\sum_i n_i^2\f$ */
/**
* Sum of squared number of connections in all grid cells,
* \f[
* \mathit{sum\_ngconn2} = \sum_c n_c^2.
* \f]
*/
size_t sum_ngconn2; size_t sum_ngconn2;
/* Linear system */ /* Linear system */
struct CSRMatrix *A; /* Coefficient matrix */ struct CSRMatrix *A; /**< Coefficient matrix */
double *b; /* System RHS */ double *b; /**< System RHS */
double *x; /* Solution */ double *x; /**< Solution */
/* Private implementational details. */ /** Private implementational details. */
struct fsh_impl *pimpl; struct fsh_impl *pimpl;
}; };
/** Destroys the fsh data object */ /**
* Dispose of all memory associated to <CODE>FSH</CODE> object.
*
* @param[in,out] h <CODE>FSH</CODE> object. Following a call
* to function fsh_destroy(), the pointer is
* invalid.
*/
void void
fsh_destroy(struct fsh_data *h); fsh_destroy(struct fsh_data *h);
/*********************************/ /**
/* Compressible solver routines. */ * Construct compressible hybrid flow-solver data object for a
/*********************************/ * given grid and well pattern.
*
* @param[in] G Grid.
/** Constructs compressible hybrid flow-solver data object for a * @param[in] W Well topology. Ignored.
* given grid and well pattern. * @return Fully formed data object suitable for use in a
* compressible pressure solver. @c NULL in case of construction
* failure.
*/ */
struct fsh_data * struct fsh_data *
cfsh_construct(struct UnstructuredGrid *G, well_t *W); cfsh_construct(struct UnstructuredGrid *G, well_t *W);
/** Assembles the hybridized linear system for face pressures. /**
* 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
* <CODE>h->A</CODE> and <CODE>h->b</CODE>.
*
* @param[in] bc Boundary conditions.
* @param[in] src Explicit source terms.
* @param[in] Binv Inverse of block-diagonal matrix \f$B\f$
* Typically computed using functions
* mim_ip_simple_all() and
* mim_ip_mobility_update().
* @param[in] Biv \f$B^{-1}v\f$.
* @param[in] P Compressible accumulation term.
* @param[in] gpress Gravity pressure.
* @param[in] wctrl Well controls. Ignored.
* @param[in] WI Well indices. Ignored.
* @param[in] BivW \f$B^{-1}v\f$ for wells. Ignored.
* @param[in] wdp Gravity pressure along well track. Ignored.
* @param[in,out] h Data object.
*/ */
void void
cfsh_assemble(struct FlowBoundaryConditions *bc, cfsh_assemble(struct FlowBoundaryConditions *bc,
@@ -92,44 +162,39 @@ cfsh_assemble(struct FlowBoundaryConditions *bc,
struct fsh_data *h); struct fsh_data *h);
/***********************************/ /**
/* Incompressible solver routines. */ * Construct incompressible hybrid flow-solver data object for a
/***********************************/ * given grid and well pattern.
/** Constructs incompressible hybrid flow-solver data object for a
* given grid and well pattern.
* *
* @param G The grid * @param G Grid.
* @param W The wells * @param W Well topology.
* @return Fully formed data object suitable for use in an
* incompressible pressure solver. @c NULL in case of construction
* failure.
*/ */
struct fsh_data * struct fsh_data *
ifsh_construct(struct UnstructuredGrid *G, well_t *W); ifsh_construct(struct UnstructuredGrid *G, well_t *W);
/**
/** Assembles the hybridized linear system for face pressures. * Form Schur-complement system of simultaneous linear equations
* arising in compressible flow using a hybridized formulation.
* *
* This routine produces no output, other than changing the linear * Upon returning from function cfsh_assemble(), the resulting
* system embedded in the ifsh_data object. * system of simultaneous linear equations is stored in
* @param bc Boundary conditions. * <CODE>h->A</CODE> and <CODE>h->b</CODE>.
* @param src Per-cell source terms (volume per second). Positive *
* values flow are sources, negative values are sinks. * @param[in] bc Boundary conditions.
* @param Binv The cell-wise effective inner products to employ in * @param[in] src Explicit source terms.
* assembly. This should be an array of length equal to * @param[in] Binv Inverse of block-diagonal matrix \f$B\f$
* sum_ngconn2 of the ifsh_data object. For each cell i, * Typically computed using functions
* there are \f$n_i^2\f$ entries, giving the inner product for * mim_ip_simple_all() and
* that cell. The inner products may for example be * mim_ip_mobility_update().
* computed by the functions of mimetic.h. * @param[in] gpress Gravity pressure.
* @param gpress Effective gravity terms. This should be an array of length * @param[in] wctrl Well controls.
* \f$\sum_i n_i\f$. For each cell, the \f$n_i\f$ elements * @param[in] WI Well indices.
* corresponding to cell \f$i\f$ should be given by * @param[in] wdp Gravity pressure along well track.
* \f$\omega g \cdot (f_c - c_c)\f$ where the symbols * @param[in,out] h Data object.
* represent the fractional-flow-weighted densities,
* the gravity vector, face centroid and cell centroid.
* @param wctrl \TODO
* @param WI \TODO
* @param wdp \TODO
* @param h The fsh_data object to use (and whose linear system will
* be modified). Must already be constructed.
*/ */
void void
ifsh_assemble(struct FlowBoundaryConditions *bc, ifsh_assemble(struct FlowBoundaryConditions *bc,
@@ -144,21 +209,24 @@ ifsh_assemble(struct FlowBoundaryConditions *bc,
/**********************************/ /**
/* Common postprocessing routine. */ * Compute cell pressures (cpress) and interface fluxes (fflux) from
/**********************************/ * current solution of system of linear equations, <CODE>h->x</CODE>.
* Back substitution process, projected half-contact fluxes.
/** Computes cell pressures, face fluxes, well pressures and well
* fluxes from face pressures.
* *
* @param G The grid. * @param[in] G Grid.
* @param h The fsh_data object. You must have called [ic]fsh_assemble() * @param[in] Binv Inverse of block-diagonal matrix \f$B\f$
* prior to this, and solved the embedded linear system of * Must coincide with the matrix used to
* this object before you call fsh_press_flux(). * form the system of linear equations
* @param cpress[out] Cell pressures. * currently stored in the data object.
* @param fflux[out] Oriented face fluxes. * @param[in] gpress Gravity pressure. Must coincide with
* @param wpress[out] \TODO * the array used to form the system of
* @param wflux[out] \TODO * linear equations.
* @param[in] h Data object.
* @param[out] cpress Cell pressure.
* @param[out] fflux Interface fluxes.
* @param[out] wpress Well pressure.
* @param[out] wflux Well perforation fluxes.
*/ */
void void
fsh_press_flux(struct UnstructuredGrid *G, fsh_press_flux(struct UnstructuredGrid *G,