/* Copyright 2013 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 . */ #ifndef OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED #define OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED #include #include #include struct UnstructuredGrid; struct Wells; namespace Opm { class DerivedGeology; class RockCompressibility; class LinearSolverInterface; class BlackoilState; class WellState; /// 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. class FullyImplicitBlackoilSolver { public: /// 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 FullyImplicitBlackoilSolver(const UnstructuredGrid& grid , const BlackoilPropsAdInterface& fluid, const DerivedGeology& geo , const RockCompressibility* rock_comp_props, const Wells& wells, const LinearSolverInterface& linsolver); /// Take a single forward step, modifiying /// state.pressure() /// state.faceflux() /// state.saturation() /// state.gasoilratio() /// wstate.bhp() /// \param[in] dt time step size /// \param[in] state reservoir state /// \param[in] wstate well state void step(const double dt , BlackoilState& state , WellState& wstate); private: // Types and enums typedef AutoDiffBlock ADB; typedef ADB::V V; typedef ADB::M M; typedef Eigen::Array DataBlock; struct ReservoirResidualQuant { ReservoirResidualQuant(); std::vector 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 { SolutionState(const int np); ADB pressure; std::vector saturation; ADB rs; ADB rv; ADB qs; ADB bhp; }; struct WellOps { WellOps(const Wells& wells); M w2p; // well -> perf (scatter) M p2w; // perf -> well (gather) }; enum { Water = BlackoilPropsAdInterface::Water, Oil = BlackoilPropsAdInterface::Oil , Gas = BlackoilPropsAdInterface::Gas }; // Member data const UnstructuredGrid& grid_; const BlackoilPropsAdInterface& fluid_; const DerivedGeology& geo_; const RockCompressibility* rock_comp_props_; const Wells& wells_; const LinearSolverInterface& linsolver_; // For each canonical phase -> true if active const std::vector active_; // Size = # active faces. Maps active -> canonical phase indices. const std::vector canph_; const std::vector cells_; // All grid cells HelperOps ops_; const WellOps wops_; const M grav_; std::vector rq_; std::vector phaseCondition_; // 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. struct { std::vector mass_balance; ADB rs_or_sg_eq; // Only used if both gas and oil present ADB well_flux_eq; ADB well_eq; } residual_; // Private methods. SolutionState constantState(const BlackoilState& x, const WellState& xw); SolutionState variableState(const BlackoilState& x, const WellState& xw); void computeAccum(const SolutionState& state, const int aix ); void assemble(const V& dtpv, const BlackoilState& x , const WellState& xw ); V solveJacobianSystem() const; void updateState(const V& dx, BlackoilState& state, WellState& well_state); std::vector computePressures(const SolutionState& state) const; std::vector computeRelPerm(const SolutionState& state) const; std::vector computeRelPermWells(const SolutionState& state, const DataBlock& well_s, const std::vector& well_cells) const; void computeMassFlux(const int actph , const V& transi, const ADB& kr , const ADB& p , const SolutionState& state ); double residualNorm() const; ADB fluidViscosity(const int phase, const ADB& p , const ADB& rs , const ADB& rv , const std::vector& cond, const std::vector& cells) const; ADB fluidReciprocFVF(const int phase, const ADB& p , const ADB& rs , const ADB& rv , const std::vector& cond, const std::vector& cells) const; ADB fluidDensity(const int phase, const ADB& p , const ADB& rs , const ADB& rv , const std::vector& cond, const std::vector& cells) const; V fluidRsSat(const V& p, const std::vector& cells) const; ADB fluidRsSat(const ADB& p, const std::vector& cells) const; V fluidRvSat(const V& p, const std::vector& cells) const; ADB fluidRvSat(const ADB& p, const std::vector& cells) const; ADB poroMult(const ADB& p) const; ADB transMult(const ADB& p) const; void classifyCondition(const SolutionState& state, std::vector& cond ) const; const std::vector phaseCondition() const {return phaseCondition_;} void classifyCondition(const BlackoilState& state); }; } // namespace Opm #endif // OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED