2013-05-23 11:28:50 -05:00
|
|
|
/*
|
|
|
|
Copyright 2013 SINTEF ICT, Applied Mathematics.
|
|
|
|
|
2013-05-24 03:49:59 -05:00
|
|
|
This file is part of the Open Porous Media project (OPM).
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
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/>.
|
|
|
|
*/
|
|
|
|
|
2013-05-24 03:49:59 -05:00
|
|
|
#ifndef OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED
|
|
|
|
#define OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
#include <opm/autodiff/AutoDiffBlock.hpp>
|
|
|
|
#include <opm/autodiff/AutoDiffHelpers.hpp>
|
|
|
|
#include <opm/autodiff/BlackoilPropsAdInterface.hpp>
|
|
|
|
|
2013-05-24 04:00:55 -05:00
|
|
|
struct UnstructuredGrid;
|
|
|
|
struct Wells;
|
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
namespace Opm {
|
|
|
|
|
2013-05-24 04:00:55 -05:00
|
|
|
class DerivedGeology;
|
2013-06-03 07:14:48 -05:00
|
|
|
class RockCompressibility;
|
2013-05-24 04:00:55 -05:00
|
|
|
class LinearSolverInterface;
|
|
|
|
class BlackoilState;
|
|
|
|
class WellState;
|
|
|
|
|
2013-05-24 04:14:05 -05:00
|
|
|
|
2013-09-20 07:55:43 -05:00
|
|
|
/// A fully implicit solver for the black-oil problem.
|
|
|
|
///
|
|
|
|
/// The simulator is capable of handling three-phase problems
|
|
|
|
/// where gas can be dissolved in oil (but not vice versa). It
|
|
|
|
/// uses an industry-standard TPFA discretization with per-phase
|
|
|
|
/// upwind weighting of mobilities.
|
|
|
|
///
|
|
|
|
/// It uses automatic differentiation via the class AutoDiffBlock
|
|
|
|
/// to simplify assembly of the jacobian matrix.
|
2013-05-24 04:00:55 -05:00
|
|
|
class FullyImplicitBlackoilSolver
|
|
|
|
{
|
2013-05-23 11:28:50 -05:00
|
|
|
public:
|
2013-09-20 07:55:43 -05:00
|
|
|
/// Construct a solver. It will retain references to the
|
|
|
|
/// arguments of this functions, and they are expected to
|
|
|
|
/// remain in scope for the lifetime of the solver.
|
|
|
|
/// \param[in] grid grid data structure
|
|
|
|
/// \param[in] fluid fluid properties
|
|
|
|
/// \param[in] geo rock properties
|
|
|
|
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
|
|
|
/// \param[in] wells well structure
|
|
|
|
/// \param[in] linsolver linear solver
|
2013-05-24 03:49:59 -05:00
|
|
|
FullyImplicitBlackoilSolver(const UnstructuredGrid& grid ,
|
|
|
|
const BlackoilPropsAdInterface& fluid,
|
2013-05-24 08:27:19 -05:00
|
|
|
const DerivedGeology& geo ,
|
2013-06-03 07:14:48 -05:00
|
|
|
const RockCompressibility* rock_comp_props,
|
2013-05-26 04:49:44 -05:00
|
|
|
const Wells& wells,
|
|
|
|
const LinearSolverInterface& linsolver);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2013-05-24 04:14:05 -05:00
|
|
|
/// Take a single forward step, modifiying
|
|
|
|
/// state.pressure()
|
|
|
|
/// state.faceflux()
|
|
|
|
/// state.saturation()
|
2013-05-27 04:32:35 -05:00
|
|
|
/// state.gasoilratio()
|
|
|
|
/// wstate.bhp()
|
2013-09-20 07:55:43 -05:00
|
|
|
/// \param[in] dt time step size
|
|
|
|
/// \param[in] state reservoir state
|
|
|
|
/// \param[in] wstate well state
|
2013-05-23 11:28:50 -05:00
|
|
|
void
|
2013-05-24 08:27:19 -05:00
|
|
|
step(const double dt ,
|
|
|
|
BlackoilState& state ,
|
|
|
|
WellState& wstate);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
private:
|
2013-05-24 04:14:05 -05:00
|
|
|
// Types and enums
|
2013-09-19 05:53:28 -05:00
|
|
|
typedef AutoDiffBlock<double> ADB;
|
2013-05-23 11:28:50 -05:00
|
|
|
typedef ADB::V V;
|
|
|
|
typedef ADB::M M;
|
|
|
|
typedef Eigen::Array<double,
|
|
|
|
Eigen::Dynamic,
|
|
|
|
Eigen::Dynamic,
|
|
|
|
Eigen::RowMajor> DataBlock;
|
|
|
|
|
|
|
|
struct ReservoirResidualQuant {
|
2013-05-24 04:14:05 -05:00
|
|
|
ReservoirResidualQuant();
|
2013-05-23 11:28:50 -05:00
|
|
|
std::vector<ADB> accum; // Accumulations
|
|
|
|
ADB mflux; // Mass flux (surface conditions)
|
|
|
|
ADB b; // Reciprocal FVF
|
|
|
|
ADB head; // Pressure drop across int. interfaces
|
|
|
|
ADB mob; // Phase mobility (per cell)
|
|
|
|
};
|
|
|
|
|
|
|
|
struct SolutionState {
|
2013-05-24 04:14:05 -05:00
|
|
|
SolutionState(const int np);
|
2013-05-23 11:28:50 -05:00
|
|
|
ADB pressure;
|
|
|
|
std::vector<ADB> saturation;
|
2013-06-02 01:19:21 -05:00
|
|
|
ADB rs;
|
2013-06-02 01:58:30 -05:00
|
|
|
ADB qs;
|
2013-05-24 16:20:15 -05:00
|
|
|
ADB bhp;
|
2013-05-23 11:28:50 -05:00
|
|
|
};
|
|
|
|
|
2013-05-24 10:22:35 -05:00
|
|
|
struct WellOps {
|
|
|
|
WellOps(const Wells& wells);
|
|
|
|
M w2p; // well -> perf (scatter)
|
|
|
|
M p2w; // perf -> well (gather)
|
|
|
|
};
|
|
|
|
|
2013-05-24 04:14:05 -05:00
|
|
|
enum { Water = BlackoilPropsAdInterface::Water,
|
|
|
|
Oil = BlackoilPropsAdInterface::Oil ,
|
|
|
|
Gas = BlackoilPropsAdInterface::Gas };
|
|
|
|
|
|
|
|
// Member data
|
2013-05-23 11:28:50 -05:00
|
|
|
const UnstructuredGrid& grid_;
|
|
|
|
const BlackoilPropsAdInterface& fluid_;
|
2013-05-24 04:00:55 -05:00
|
|
|
const DerivedGeology& geo_;
|
2013-06-03 07:14:48 -05:00
|
|
|
const RockCompressibility* rock_comp_props_;
|
2013-05-24 08:27:19 -05:00
|
|
|
const Wells& wells_;
|
2013-05-26 04:49:44 -05:00
|
|
|
const LinearSolverInterface& linsolver_;
|
2013-05-24 03:39:10 -05:00
|
|
|
// For each canonical phase -> true if active
|
|
|
|
const std::vector<bool> active_;
|
|
|
|
// Size = # active faces. Maps active -> canonical phase indices.
|
|
|
|
const std::vector<int> canph_;
|
2013-05-23 11:28:50 -05:00
|
|
|
const std::vector<int> cells_; // All grid cells
|
|
|
|
HelperOps ops_;
|
2013-05-24 10:22:35 -05:00
|
|
|
const WellOps wops_;
|
2013-05-23 11:28:50 -05:00
|
|
|
const M grav_;
|
|
|
|
|
|
|
|
std::vector<ReservoirResidualQuant> rq_;
|
|
|
|
|
2013-05-24 15:16:06 -05:00
|
|
|
// The mass_balance vector has one element for each active phase,
|
|
|
|
// each of which has size equal to the number of cells.
|
|
|
|
// The well_eq has size equal to the number of wells.
|
2013-05-23 11:28:50 -05:00
|
|
|
struct {
|
2013-05-24 15:16:06 -05:00
|
|
|
std::vector<ADB> mass_balance;
|
2013-05-27 04:32:35 -05:00
|
|
|
ADB rs_or_sg_eq; // Only used if both gas and oil present
|
2013-06-02 01:58:30 -05:00
|
|
|
ADB well_flux_eq;
|
2013-05-24 15:16:06 -05:00
|
|
|
ADB well_eq;
|
2013-05-23 11:28:50 -05:00
|
|
|
} residual_;
|
|
|
|
|
2013-05-24 04:14:05 -05:00
|
|
|
// Private methods.
|
2013-05-23 11:28:50 -05:00
|
|
|
SolutionState
|
2013-05-24 16:20:15 -05:00
|
|
|
constantState(const BlackoilState& x,
|
|
|
|
const WellState& xw);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
SolutionState
|
2013-05-24 16:20:15 -05:00
|
|
|
variableState(const BlackoilState& x,
|
|
|
|
const WellState& xw);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
void
|
|
|
|
computeAccum(const SolutionState& state,
|
2013-05-24 04:14:05 -05:00
|
|
|
const int aix );
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
void
|
2013-05-24 08:27:19 -05:00
|
|
|
assemble(const V& dtpv,
|
|
|
|
const BlackoilState& x ,
|
|
|
|
const WellState& xw );
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2013-05-30 07:43:32 -05:00
|
|
|
V solveJacobianSystem() const;
|
|
|
|
|
|
|
|
void updateState(const V& dx,
|
|
|
|
BlackoilState& state,
|
|
|
|
WellState& well_state) const;
|
2013-05-26 04:49:44 -05:00
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
std::vector<ADB>
|
2013-05-25 03:47:22 -05:00
|
|
|
computeRelPerm(const SolutionState& state) const;
|
|
|
|
|
|
|
|
std::vector<ADB>
|
|
|
|
computeRelPermWells(const SolutionState& state,
|
|
|
|
const DataBlock& well_s,
|
|
|
|
const std::vector<int>& well_cells) const;
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
void
|
|
|
|
computeMassFlux(const int actph ,
|
|
|
|
const V& transi,
|
|
|
|
const std::vector<ADB>& kr ,
|
2013-05-24 04:14:05 -05:00
|
|
|
const SolutionState& state );
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2013-05-24 04:38:17 -05:00
|
|
|
double
|
|
|
|
residualNorm() const;
|
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
ADB
|
|
|
|
fluidViscosity(const int phase,
|
|
|
|
const ADB& p ,
|
2013-05-30 07:43:32 -05:00
|
|
|
const ADB& rs ,
|
2013-11-28 04:27:25 -06:00
|
|
|
const bool* isSat,
|
2013-05-24 04:14:05 -05:00
|
|
|
const std::vector<int>& cells) const;
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
ADB
|
|
|
|
fluidReciprocFVF(const int phase,
|
|
|
|
const ADB& p ,
|
2013-05-30 07:43:32 -05:00
|
|
|
const ADB& rs ,
|
2013-11-28 04:27:25 -06:00
|
|
|
const bool* isSat,
|
2013-05-24 04:14:05 -05:00
|
|
|
const std::vector<int>& cells) const;
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
ADB
|
|
|
|
fluidDensity(const int phase,
|
|
|
|
const ADB& p ,
|
2013-05-30 07:43:32 -05:00
|
|
|
const ADB& rs ,
|
2013-11-28 04:27:25 -06:00
|
|
|
const bool* isSat,
|
2013-05-24 04:14:05 -05:00
|
|
|
const std::vector<int>& cells) const;
|
2013-05-27 03:29:04 -05:00
|
|
|
|
2013-05-30 07:43:32 -05:00
|
|
|
V
|
|
|
|
fluidRsMax(const V& p,
|
|
|
|
const std::vector<int>& cells) const;
|
|
|
|
|
2013-05-27 03:29:04 -05:00
|
|
|
ADB
|
|
|
|
fluidRsMax(const ADB& p,
|
|
|
|
const std::vector<int>& cells) const;
|
2013-06-03 07:14:48 -05:00
|
|
|
|
|
|
|
ADB
|
|
|
|
poroMult(const ADB& p) const;
|
|
|
|
|
|
|
|
ADB
|
|
|
|
transMult(const ADB& p) const;
|
2013-11-28 04:27:25 -06:00
|
|
|
|
|
|
|
void
|
|
|
|
getSaturatedCells(const SolutionState& state, bool* isSat) const;
|
2013-05-23 11:28:50 -05:00
|
|
|
};
|
|
|
|
} // namespace Opm
|
|
|
|
|
2013-05-24 03:49:59 -05:00
|
|
|
|
|
|
|
#endif // OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED
|