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
|
|
|
|
2015-01-20 05:55:46 -06:00
|
|
|
#include <cassert>
|
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
#include <opm/autodiff/AutoDiffBlock.hpp>
|
|
|
|
#include <opm/autodiff/AutoDiffHelpers.hpp>
|
|
|
|
#include <opm/autodiff/BlackoilPropsAdInterface.hpp>
|
2014-04-08 08:56:08 -05:00
|
|
|
#include <opm/autodiff/LinearisedBlackoilResidual.hpp>
|
2014-04-08 09:11:08 -05:00
|
|
|
#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2015-01-20 07:21:45 -06:00
|
|
|
#include <array>
|
|
|
|
|
2013-05-24 04:00:55 -05:00
|
|
|
struct UnstructuredGrid;
|
|
|
|
struct Wells;
|
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
namespace Opm {
|
|
|
|
|
2014-05-16 11:02:55 -05:00
|
|
|
namespace parameter { class ParameterGroup; }
|
2013-05-24 04:00:55 -05:00
|
|
|
class DerivedGeology;
|
2013-06-03 07:14:48 -05:00
|
|
|
class RockCompressibility;
|
2014-04-08 09:11:08 -05:00
|
|
|
class NewtonIterationBlackoilInterface;
|
2013-05-24 04:00:55 -05:00
|
|
|
class BlackoilState;
|
2014-03-18 05:27:40 -05:00
|
|
|
class WellStateFullyImplicitBlackoil;
|
2013-05-24 04:00:55 -05:00
|
|
|
|
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.
|
2014-02-20 06:17:18 -06:00
|
|
|
template<class T>
|
2013-05-24 04:00:55 -05:00
|
|
|
class FullyImplicitBlackoilSolver
|
|
|
|
{
|
2013-05-23 11:28:50 -05:00
|
|
|
public:
|
2014-10-01 05:48:41 -05:00
|
|
|
// the Newton relaxation type
|
|
|
|
enum RelaxType { DAMPEN, SOR };
|
|
|
|
|
|
|
|
// class holding the solver parameters
|
|
|
|
struct SolverParameter
|
|
|
|
{
|
|
|
|
double dp_max_rel_;
|
|
|
|
double ds_max_;
|
2014-11-05 06:03:00 -06:00
|
|
|
double dr_max_rel_;
|
2014-10-01 05:48:41 -05:00
|
|
|
enum RelaxType relax_type_;
|
|
|
|
double relax_max_;
|
|
|
|
double relax_increment_;
|
|
|
|
double relax_rel_tol_;
|
2015-01-21 08:43:18 -06:00
|
|
|
double max_residual_allowed_;
|
2015-03-03 06:23:33 -06:00
|
|
|
double tolerance_mb_;
|
|
|
|
double tolerance_cnv_;
|
2015-03-04 05:31:04 -06:00
|
|
|
double tolerance_wells_;
|
2014-10-01 05:48:41 -05:00
|
|
|
int max_iter_;
|
|
|
|
|
|
|
|
SolverParameter( const parameter::ParameterGroup& param );
|
|
|
|
SolverParameter();
|
|
|
|
|
|
|
|
void reset();
|
|
|
|
};
|
|
|
|
|
2014-02-20 06:17:18 -06:00
|
|
|
/// \brief The type of the grid that we use.
|
|
|
|
typedef T Grid;
|
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.
|
2014-10-01 05:48:41 -05:00
|
|
|
/// \param[in] param parameters
|
2013-09-20 07:55:43 -05:00
|
|
|
/// \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
|
2014-10-01 05:48:41 -05:00
|
|
|
FullyImplicitBlackoilSolver(const SolverParameter& param,
|
2014-05-16 11:02:55 -05:00
|
|
|
const Grid& grid ,
|
2013-05-24 03:49:59 -05:00
|
|
|
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,
|
2015-01-19 07:14:18 -06:00
|
|
|
const Wells* wells,
|
2014-05-27 06:36:22 -05:00
|
|
|
const NewtonIterationBlackoilInterface& linsolver,
|
|
|
|
const bool has_disgas,
|
2015-03-04 06:42:27 -06:00
|
|
|
const bool has_vapoil,
|
|
|
|
const bool terminal_output);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2014-08-27 04:30:42 -05:00
|
|
|
/// \brief Set threshold pressures that prevent or reduce flow.
|
|
|
|
/// This prevents flow across faces if the potential
|
|
|
|
/// difference is less than the threshold. If the potential
|
|
|
|
/// difference is greater, the threshold value is subtracted
|
|
|
|
/// before calculating flow. This is treated symmetrically, so
|
|
|
|
/// flow is prevented or reduced in both directions equally.
|
2014-08-27 07:21:34 -05:00
|
|
|
/// \param[in] threshold_pressures_by_face array of size equal to the number of faces
|
2014-08-27 04:30:42 -05:00
|
|
|
/// of the grid passed in the constructor.
|
2014-08-27 07:21:34 -05:00
|
|
|
void setThresholdPressures(const std::vector<double>& threshold_pressures_by_face);
|
2014-08-27 04:30:42 -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
|
2014-10-20 07:47:45 -05:00
|
|
|
/// \return number of linear iterations used
|
2014-10-03 07:18:31 -05:00
|
|
|
int
|
2013-05-24 08:27:19 -05:00
|
|
|
step(const double dt ,
|
|
|
|
BlackoilState& state ,
|
2014-03-18 05:27:40 -05:00
|
|
|
WellStateFullyImplicitBlackoil& wstate);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2015-03-04 08:02:00 -06:00
|
|
|
unsigned int newtonIterations () const { return newtonIterations_; }
|
|
|
|
unsigned int linearIterations () const { return linearIterations_; }
|
|
|
|
|
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;
|
2014-11-20 05:31:50 -06:00
|
|
|
ADB temperature;
|
2013-05-23 11:28:50 -05:00
|
|
|
std::vector<ADB> saturation;
|
2013-06-02 01:19:21 -05:00
|
|
|
ADB rs;
|
2014-01-10 07:19:37 -06:00
|
|
|
ADB rv;
|
2013-06-02 01:58:30 -05:00
|
|
|
ADB qs;
|
2015-01-19 07:14:18 -06:00
|
|
|
ADB bhp;
|
2015-03-06 04:43:16 -06:00
|
|
|
// Below are quantities stored in the state for optimization purposes.
|
|
|
|
std::vector<ADB> canonical_phase_pressures; // Always has 3 elements, even if only 2 phases active.
|
2013-05-23 11:28:50 -05:00
|
|
|
};
|
|
|
|
|
2013-05-24 10:22:35 -05:00
|
|
|
struct WellOps {
|
2015-01-19 07:14:18 -06:00
|
|
|
WellOps(const Wells* wells);
|
2013-05-24 10:22:35 -05:00
|
|
|
M w2p; // well -> perf (scatter)
|
|
|
|
M p2w; // perf -> well (gather)
|
|
|
|
};
|
|
|
|
|
2015-01-19 13:07:31 -06:00
|
|
|
enum { Water = BlackoilPropsAdInterface::Water,
|
|
|
|
Oil = BlackoilPropsAdInterface::Oil ,
|
|
|
|
Gas = BlackoilPropsAdInterface::Gas ,
|
|
|
|
MaxNumPhases = BlackoilPropsAdInterface::MaxNumPhases
|
|
|
|
};
|
2013-05-24 04:14:05 -05:00
|
|
|
|
2014-05-27 08:09:58 -05:00
|
|
|
enum PrimalVariables { Sg = 0, RS = 1, RV = 2 };
|
2014-05-19 11:41:38 -05:00
|
|
|
|
2013-05-24 04:14:05 -05:00
|
|
|
// Member data
|
2014-02-20 06:17:18 -06:00
|
|
|
const Grid& grid_;
|
2013-05-23 11:28:50 -05:00
|
|
|
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_;
|
2015-01-19 07:14:18 -06:00
|
|
|
const Wells* wells_;
|
2014-04-08 09:11:08 -05:00
|
|
|
const NewtonIterationBlackoilInterface& linsolver_;
|
2013-05-24 03:39:10 -05:00
|
|
|
// For each canonical phase -> true if active
|
|
|
|
const std::vector<bool> active_;
|
2014-06-03 03:59:50 -05:00
|
|
|
// Size = # active phases. Maps active -> canonical phase indices.
|
2013-05-24 03:39:10 -05:00
|
|
|
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_;
|
2014-05-27 06:36:22 -05:00
|
|
|
const bool has_disgas_;
|
|
|
|
const bool has_vapoil_;
|
2014-10-01 05:48:41 -05:00
|
|
|
|
|
|
|
SolverParameter param_;
|
2014-08-27 03:27:49 -05:00
|
|
|
bool use_threshold_pressure_;
|
2014-08-27 06:46:23 -05:00
|
|
|
V threshold_pressures_by_interior_face_;
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
std::vector<ReservoirResidualQuant> rq_;
|
2014-01-10 07:19:37 -06:00
|
|
|
std::vector<PhasePresence> phaseCondition_;
|
2014-03-18 02:48:34 -05:00
|
|
|
V well_perforation_pressure_diffs_; // Diff to bhp for each well perforation.
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2014-04-08 08:56:08 -05:00
|
|
|
LinearisedBlackoilResidual residual_;
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2015-02-20 09:03:08 -06:00
|
|
|
/// \brief Whether we print something to std::cout
|
2015-02-20 09:02:06 -06:00
|
|
|
bool terminal_output_;
|
2015-03-04 08:02:00 -06:00
|
|
|
unsigned int newtonIterations_;
|
|
|
|
unsigned int linearIterations_;
|
2015-02-20 04:35:47 -06:00
|
|
|
|
2014-05-27 08:09:58 -05:00
|
|
|
std::vector<int> primalVariable_;
|
|
|
|
|
2013-05-24 04:14:05 -05:00
|
|
|
// Private methods.
|
2015-01-19 07:14:18 -06:00
|
|
|
|
|
|
|
// return true if wells are available
|
|
|
|
bool wellsActive() const { return wells_ ? wells_->number_of_wells > 0 : false ; }
|
|
|
|
// return wells object
|
2015-01-20 05:55:46 -06:00
|
|
|
const Wells& wells () const { assert( bool(wells_ != 0) ); return *wells_; }
|
2015-01-19 07:14:18 -06:00
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
SolutionState
|
2013-05-24 16:20:15 -05:00
|
|
|
constantState(const BlackoilState& x,
|
2014-03-18 05:27:40 -05:00
|
|
|
const WellStateFullyImplicitBlackoil& xw);
|
2013-05-23 11:28:50 -05:00
|
|
|
|
|
|
|
SolutionState
|
2013-05-24 16:20:15 -05:00
|
|
|
variableState(const BlackoilState& x,
|
2014-03-18 05:27:40 -05:00
|
|
|
const WellStateFullyImplicitBlackoil& 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
|
|
|
|
2014-03-18 02:48:34 -05:00
|
|
|
void computeWellConnectionPressures(const SolutionState& state,
|
2014-03-20 04:04:45 -05:00
|
|
|
const WellStateFullyImplicitBlackoil& xw);
|
2014-03-18 02:48:34 -05:00
|
|
|
|
2014-03-11 08:28:16 -05:00
|
|
|
void
|
2014-03-25 05:14:11 -05:00
|
|
|
addWellControlEq(const SolutionState& state,
|
2014-04-24 11:50:49 -05:00
|
|
|
const WellStateFullyImplicitBlackoil& xw,
|
|
|
|
const V& aliveWells);
|
2014-03-11 08:28:16 -05:00
|
|
|
|
|
|
|
void
|
2014-04-24 11:50:49 -05:00
|
|
|
addWellEq(const SolutionState& state,
|
2014-05-08 06:32:49 -05:00
|
|
|
WellStateFullyImplicitBlackoil& xw,
|
2014-04-24 11:50:49 -05:00
|
|
|
V& aliveWells);
|
2014-03-11 08:28:16 -05:00
|
|
|
|
2014-04-22 06:52:42 -05:00
|
|
|
void updateWellControls(ADB& bhp,
|
|
|
|
ADB& well_phase_flow_rate,
|
2014-03-25 08:31:06 -05:00
|
|
|
WellStateFullyImplicitBlackoil& xw) const;
|
2014-03-25 05:14:11 -05:00
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
void
|
2013-05-24 08:27:19 -05:00
|
|
|
assemble(const V& dtpv,
|
2014-03-18 05:27:40 -05:00
|
|
|
const BlackoilState& x,
|
2014-03-25 05:14:11 -05:00
|
|
|
WellStateFullyImplicitBlackoil& 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,
|
2014-03-18 05:27:40 -05:00
|
|
|
WellStateFullyImplicitBlackoil& well_state);
|
2013-05-26 04:49:44 -05:00
|
|
|
|
2013-11-26 06:51:10 -06:00
|
|
|
std::vector<ADB>
|
|
|
|
computePressures(const SolutionState& state) const;
|
|
|
|
|
2015-01-07 04:24:11 -06:00
|
|
|
std::vector<ADB>
|
|
|
|
computePressures(const ADB& po,
|
|
|
|
const ADB& sw,
|
|
|
|
const ADB& so,
|
|
|
|
const ADB& sg) const;
|
|
|
|
|
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,
|
2013-11-26 06:51:10 -06:00
|
|
|
const ADB& kr ,
|
|
|
|
const ADB& p ,
|
2013-05-24 04:14:05 -05:00
|
|
|
const SolutionState& state );
|
2013-05-23 11:28:50 -05:00
|
|
|
|
2014-08-27 03:27:49 -05:00
|
|
|
void applyThresholdPressures(ADB& dp);
|
|
|
|
|
2013-05-24 04:38:17 -05:00
|
|
|
double
|
|
|
|
residualNorm() const;
|
|
|
|
|
2015-01-27 04:57:36 -06:00
|
|
|
/// \brief Compute the residual norms of the mass balance for each phase,
|
|
|
|
/// the well flux, and the well equation.
|
|
|
|
/// \return a vector that contains for each phase the norm of the mass balance
|
2015-01-27 07:05:38 -06:00
|
|
|
/// and afterwards the norm of the residual of the well flux and the well equation.
|
2015-01-27 04:57:36 -06:00
|
|
|
std::vector<double> computeResidualNorms() const;
|
2014-05-19 03:41:23 -05:00
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
ADB
|
|
|
|
fluidViscosity(const int phase,
|
|
|
|
const ADB& p ,
|
2014-11-20 05:31:50 -06:00
|
|
|
const ADB& temp ,
|
2013-05-30 07:43:32 -05:00
|
|
|
const ADB& rs ,
|
2014-01-10 07:19:37 -06:00
|
|
|
const ADB& rv ,
|
2013-12-03 11:12:54 -06:00
|
|
|
const std::vector<PhasePresence>& cond,
|
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 ,
|
2014-11-20 05:31:50 -06:00
|
|
|
const ADB& temp ,
|
2013-05-30 07:43:32 -05:00
|
|
|
const ADB& rs ,
|
2014-01-10 07:19:37 -06:00
|
|
|
const ADB& rv ,
|
2013-12-03 11:12:54 -06:00
|
|
|
const std::vector<PhasePresence>& cond,
|
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 ,
|
2014-11-20 05:31:50 -06:00
|
|
|
const ADB& temp ,
|
2013-05-30 07:43:32 -05:00
|
|
|
const ADB& rs ,
|
2014-01-10 07:19:37 -06:00
|
|
|
const ADB& rv ,
|
2013-12-03 11:12:54 -06:00
|
|
|
const std::vector<PhasePresence>& cond,
|
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
|
2014-01-10 07:19:37 -06:00
|
|
|
fluidRsSat(const V& p,
|
2014-07-05 08:06:12 -05:00
|
|
|
const V& so,
|
2013-05-30 07:43:32 -05:00
|
|
|
const std::vector<int>& cells) const;
|
|
|
|
|
2013-05-27 03:29:04 -05:00
|
|
|
ADB
|
2014-01-10 07:19:37 -06:00
|
|
|
fluidRsSat(const ADB& p,
|
2014-07-05 08:06:12 -05:00
|
|
|
const ADB& so,
|
2014-01-10 07:19:37 -06:00
|
|
|
const std::vector<int>& cells) const;
|
|
|
|
|
|
|
|
V
|
|
|
|
fluidRvSat(const V& p,
|
2014-07-05 08:06:12 -05:00
|
|
|
const V& so,
|
2014-01-10 07:19:37 -06:00
|
|
|
const std::vector<int>& cells) const;
|
|
|
|
|
|
|
|
ADB
|
|
|
|
fluidRvSat(const ADB& p,
|
2014-07-05 08:06:12 -05:00
|
|
|
const ADB& so,
|
2013-05-27 03:29:04 -05:00
|
|
|
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
|
2013-12-03 11:12:54 -06:00
|
|
|
classifyCondition(const SolutionState& state,
|
|
|
|
std::vector<PhasePresence>& cond ) const;
|
2014-01-10 07:19:37 -06:00
|
|
|
|
|
|
|
const std::vector<PhasePresence>
|
|
|
|
phaseCondition() const {return phaseCondition_;}
|
|
|
|
|
|
|
|
void
|
|
|
|
classifyCondition(const BlackoilState& state);
|
|
|
|
|
2014-06-03 06:30:35 -05:00
|
|
|
|
|
|
|
/// update the primal variable for Sg, Rv or Rs. The Gas phase must
|
|
|
|
/// be active to call this method.
|
2014-05-27 08:09:58 -05:00
|
|
|
void
|
|
|
|
updatePrimalVariableFromState(const BlackoilState& state);
|
|
|
|
|
2014-07-17 07:34:07 -05:00
|
|
|
/// Update the phaseCondition_ member based on the primalVariable_ member.
|
|
|
|
void
|
|
|
|
updatePhaseCondFromPrimalVariable();
|
|
|
|
|
2014-05-06 08:47:08 -05:00
|
|
|
/// Compute convergence based on total mass balance (tol_mb) and maximum
|
|
|
|
/// residual mass balance (tol_cnv).
|
2014-11-10 02:29:13 -06:00
|
|
|
bool getConvergence(const double dt, const int iteration);
|
2014-05-06 08:47:08 -05:00
|
|
|
|
2015-01-27 05:49:55 -06:00
|
|
|
/// \brief Compute the reduction within the convergence check.
|
2015-01-28 09:53:16 -06:00
|
|
|
/// \param[in] B A matrix with MaxNumPhases columns and the same number rows
|
|
|
|
/// as the number of cells of the grid. B.col(i) contains the values
|
|
|
|
/// for phase i.
|
|
|
|
/// \param[in] tempV A matrix with MaxNumPhases columns and the same number rows
|
|
|
|
/// as the number of cells of the grid. tempV.col(i) contains the
|
|
|
|
/// values
|
|
|
|
/// for phase i.
|
|
|
|
/// \param[in] R A matrix with MaxNumPhases columns and the same number rows
|
|
|
|
/// as the number of cells of the grid. B.col(i) contains the values
|
|
|
|
/// for phase i.
|
|
|
|
/// \param[out] R_sum An array of size MaxNumPhases where entry i contains the sum
|
|
|
|
/// of R for the phase i.
|
2015-01-28 12:08:48 -06:00
|
|
|
/// \param[out] maxCoeff An array of size MaxNumPhases where entry i contains the
|
2015-01-28 12:17:25 -06:00
|
|
|
/// maximum of tempV for the phase i.
|
2015-01-28 12:08:48 -06:00
|
|
|
/// \param[out] B_avg An array of size MaxNumPhases where entry i contains the average
|
|
|
|
/// of B for the phase i.
|
2015-01-28 09:53:16 -06:00
|
|
|
/// \param[in] nc The number of cells of the local grid.
|
|
|
|
/// \return The total pore volume over all cells.
|
2015-01-27 05:49:55 -06:00
|
|
|
double
|
2015-01-20 07:21:45 -06:00
|
|
|
convergenceReduction(const Eigen::Array<double, Eigen::Dynamic, MaxNumPhases>& B,
|
|
|
|
const Eigen::Array<double, Eigen::Dynamic, MaxNumPhases>& tempV,
|
|
|
|
const Eigen::Array<double, Eigen::Dynamic, MaxNumPhases>& R,
|
|
|
|
std::array<double,MaxNumPhases>& R_sum,
|
2015-01-28 12:08:48 -06:00
|
|
|
std::array<double,MaxNumPhases>& maxCoeff,
|
|
|
|
std::array<double,MaxNumPhases>& B_avg,
|
2015-01-20 07:21:45 -06:00
|
|
|
int nc) const;
|
|
|
|
|
2014-05-23 04:39:26 -05:00
|
|
|
void detectNewtonOscillations(const std::vector<std::vector<double>>& residual_history,
|
2014-05-19 08:43:56 -05:00
|
|
|
const int it, const double relaxRelTol,
|
2014-05-23 04:39:26 -05:00
|
|
|
bool& oscillate, bool& stagnate) const;
|
2014-05-19 08:43:56 -05:00
|
|
|
|
2014-05-23 04:39:26 -05:00
|
|
|
void stablizeNewton(V& dx, V& dxOld, const double omega, const RelaxType relax_type) const;
|
2014-05-22 02:59:50 -05:00
|
|
|
|
2014-10-01 05:48:41 -05:00
|
|
|
double dpMaxRel() const { return param_.dp_max_rel_; }
|
|
|
|
double dsMax() const { return param_.ds_max_; }
|
2014-11-05 06:03:00 -06:00
|
|
|
double drMaxRel() const { return param_.dr_max_rel_; }
|
2014-10-01 05:48:41 -05:00
|
|
|
enum RelaxType relaxType() const { return param_.relax_type_; }
|
|
|
|
double relaxMax() const { return param_.relax_max_; };
|
|
|
|
double relaxIncrement() const { return param_.relax_increment_; };
|
|
|
|
double relaxRelTol() const { return param_.relax_rel_tol_; };
|
2015-01-21 08:43:18 -06:00
|
|
|
double maxIter() const { return param_.max_iter_; }
|
|
|
|
double maxResidualAllowed() const { return param_.max_residual_allowed_; }
|
2014-01-10 07:19:37 -06:00
|
|
|
|
2013-05-23 11:28:50 -05:00
|
|
|
};
|
|
|
|
} // namespace Opm
|
|
|
|
|
2014-02-20 06:17:18 -06:00
|
|
|
#include "FullyImplicitBlackoilSolver_impl.hpp"
|
2013-05-24 03:49:59 -05:00
|
|
|
|
|
|
|
#endif // OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED
|