opm-simulators/opm/core/pressure/IncompTpfa.hpp

122 lines
5.3 KiB
C++
Raw Normal View History

2012-02-20 06:23:01 -06:00
/*
Copyright 2012 SINTEF ICT, Applied Mathematics.
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef OPM_INCOMPTPFA_HEADER_INCLUDED
#define OPM_INCOMPTPFA_HEADER_INCLUDED
#include <vector>
struct UnstructuredGrid;
struct ifs_tpfa_data;
struct FlowBoundaryConditions;
2012-02-20 06:23:01 -06:00
namespace Opm
{
class LinearSolverInterface;
2012-02-20 06:23:01 -06:00
/// Encapsulating a tpfa pressure solver for the incompressible case.
/// Supports gravity and simple sources as driving forces.
/// Below we use the shortcuts D for the number of dimensions, N
/// for the number of cells and F for the number of faces.
/// Note: we intend to add wells in the future.
class IncompTpfa
{
public:
/// Construct solver.
/// \param[in] g A 2d or 3d grid.
/// \param[in] permeability Array of permeability tensors, the array
/// should have size N*D^2, if D == g.dimensions
/// and N == g.number_of_cells.
/// \param[in] gravity Gravity vector. If nonzero, the array should
/// have D elements.
IncompTpfa(const UnstructuredGrid& g,
const double* permeability,
const double* gravity,
const LinearSolverInterface& linsolver);
2012-02-20 06:23:01 -06:00
/// Destructor.
~IncompTpfa();
/// Assemble and solve incompressible pressure system.
2012-02-20 06:23:01 -06:00
/// \param[in] totmob Must contain N total mobility values (one per cell).
/// totmob = \sum_{p} kr_p/mu_p.
/// \param[in] omega Must be empty if constructor gravity argument was null.
/// Otherwise must contain N mobility-weighted density values (one per cell).
2012-02-20 06:23:01 -06:00
/// omega = \frac{\sum_{p} mob_p rho_p}{\sum_p rho_p}.
/// \param[in] src Must contain N source rates (one per cell).
2012-02-20 06:23:01 -06:00
/// Positive values represent total inflow rates,
/// negative values represent total outflow rates.
/// \param[in] bcs If non-null, specifies boundary conditions.
/// If null, noflow conditions are assumed.
2012-02-20 06:23:01 -06:00
/// \param[out] pressure Will contain N cell-pressure values.
/// \param[out] faceflux Will contain F signed face flux values.
void solve(const std::vector<double>& totmob,
const std::vector<double>& omega,
const std::vector<double>& src,
const FlowBoundaryConditions* bcs,
2012-02-20 06:23:01 -06:00
std::vector<double>& pressure,
std::vector<double>& faceflux);
/// Assemble and solve pressure system with rock compressibility (assumed constant per cell).
/// \param[in] totmob Must contain N total mobility values (one per cell).
/// totmob = \sum_{p} kr_p/mu_p.
/// \param[in] omega Must be empty if constructor gravity argument was null.
/// Otherwise must contain N mobility-weighted density values (one per cell).
/// omega = \frac{\sum_{p} mob_p rho_p}{\sum_p rho_p}.
/// \param[in] src Must contain N source rates (one per cell).
/// Positive values represent total inflow rates,
/// negative values represent total outflow rates.
/// \param[in] bcs If non-null, specifies boundary conditions.
/// If null, noflow conditions are assumed.
/// \param[in] porevol Must contain N pore volumes.
/// \param[in] rock_comp Must contain N rock compressibilities.
/// rock_comp = (d poro / d p)*(1/poro).
/// \param[in] dt Timestep.
/// \param[out] pressure Will contain N cell-pressure values.
/// \param[out] faceflux Will contain F signed face flux values.
void solve(const std::vector<double>& totmob,
const std::vector<double>& omega,
const std::vector<double>& src,
const FlowBoundaryConditions* bcs,
const std::vector<double>& porevol,
const std::vector<double>& rock_comp,
const double dt,
std::vector<double>& pressure,
std::vector<double>& faceflux);
/// Expose read-only reference to internal half-transmissibility.
const ::std::vector<double>& getHalfTrans() const { return htrans_; }
2012-02-20 06:23:01 -06:00
private:
const UnstructuredGrid& grid_;
const LinearSolverInterface& linsolver_;
2012-02-20 06:23:01 -06:00
::std::vector<double> htrans_;
::std::vector<double> trans_ ;
::std::vector<double> gpress_;
::std::vector<double> gpress_omegaweighted_;
struct ifs_tpfa_data* h_;
};
} // namespace Opm
#endif // OPM_INCOMPTPFA_HEADER_INCLUDED