2015-05-18 07:18:03 -05:00
|
|
|
/*
|
|
|
|
Copyright 2015 SINTEF ICT, Applied Mathematics.
|
|
|
|
Copyright 2015 Statoil ASA.
|
|
|
|
|
|
|
|
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/>.
|
|
|
|
*/
|
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
#ifndef OPM_NONLINEARSOLVER_HEADER_INCLUDED
|
|
|
|
#define OPM_NONLINEARSOLVER_HEADER_INCLUDED
|
2015-05-18 07:18:03 -05:00
|
|
|
|
2015-05-19 03:19:30 -05:00
|
|
|
#include <opm/autodiff/AutoDiffBlock.hpp>
|
2015-05-19 03:57:06 -05:00
|
|
|
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
2016-07-05 05:20:19 -05:00
|
|
|
#include <opm/core/simulator/SimulatorTimerInterface.hpp>
|
2015-05-29 07:57:49 -05:00
|
|
|
#include <memory>
|
2015-05-18 07:18:03 -05:00
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
2015-05-18 09:10:31 -05:00
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
/// A nonlinear solver class suitable for general fully-implicit models,
|
|
|
|
/// as well as pressure, transport and sequential models.
|
2015-05-18 09:10:31 -05:00
|
|
|
template <class PhysicalModel>
|
2015-10-02 02:19:07 -05:00
|
|
|
class NonlinearSolver
|
2015-05-18 07:18:03 -05:00
|
|
|
{
|
|
|
|
public:
|
2015-05-19 03:19:30 -05:00
|
|
|
// --------- Types and enums ---------
|
|
|
|
typedef AutoDiffBlock<double> ADB;
|
|
|
|
typedef ADB::V V;
|
|
|
|
typedef ADB::M M;
|
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
// Available relaxation scheme types.
|
2015-05-19 03:19:30 -05:00
|
|
|
enum RelaxType { DAMPEN, SOR };
|
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
// Solver parameters controlling nonlinear process.
|
2015-05-21 01:54:56 -05:00
|
|
|
struct SolverParameters
|
2015-05-19 03:19:30 -05:00
|
|
|
{
|
|
|
|
enum RelaxType relax_type_;
|
|
|
|
double relax_max_;
|
|
|
|
double relax_increment_;
|
|
|
|
double relax_rel_tol_;
|
2015-10-02 02:19:07 -05:00
|
|
|
int max_iter_; // max nonlinear iterations
|
|
|
|
int min_iter_; // min nonlinear iterations
|
2015-05-19 03:19:30 -05:00
|
|
|
|
2015-05-21 01:54:56 -05:00
|
|
|
explicit SolverParameters( const parameter::ParameterGroup& param );
|
|
|
|
SolverParameters();
|
2015-05-19 03:19:30 -05:00
|
|
|
|
|
|
|
void reset();
|
|
|
|
};
|
|
|
|
|
2015-05-18 09:10:31 -05:00
|
|
|
// Forwarding types from PhysicalModel.
|
|
|
|
typedef typename PhysicalModel::ReservoirState ReservoirState;
|
|
|
|
typedef typename PhysicalModel::WellState WellState;
|
|
|
|
|
2015-05-19 03:19:30 -05:00
|
|
|
// --------- Public methods ---------
|
2015-05-18 09:10:31 -05:00
|
|
|
|
2015-05-19 03:19:30 -05:00
|
|
|
/// Construct solver for a given model.
|
2015-05-28 10:26:02 -05:00
|
|
|
///
|
2015-05-29 07:57:49 -05:00
|
|
|
/// The model is a std::unique_ptr because the object to which model points to is
|
2015-10-02 02:19:07 -05:00
|
|
|
/// not allowed to be deleted as long as the NonlinearSolver object exists.
|
2015-05-28 10:26:02 -05:00
|
|
|
///
|
2015-10-02 02:19:07 -05:00
|
|
|
/// \param[in] param parameters controlling nonlinear process
|
2015-05-28 10:26:02 -05:00
|
|
|
/// \param[in, out] model physical simulation model.
|
2015-10-02 02:19:07 -05:00
|
|
|
explicit NonlinearSolver(const SolverParameters& param,
|
2015-05-29 07:57:49 -05:00
|
|
|
std::unique_ptr<PhysicalModel> model);
|
2015-05-19 03:19:30 -05:00
|
|
|
|
|
|
|
/// Take a single forward step, after which the states will be modified
|
|
|
|
/// according to the physical model.
|
2016-07-05 05:20:19 -05:00
|
|
|
/// \param[in] timer simulation timer
|
2016-01-15 08:16:40 -06:00
|
|
|
/// \param[in, out] reservoir_state reservoir state variables
|
|
|
|
/// \param[in, out] well_state well state variables
|
|
|
|
/// \return number of linear iterations used
|
2015-05-18 09:10:31 -05:00
|
|
|
int
|
2016-07-05 05:20:19 -05:00
|
|
|
step(const SimulatorTimerInterface& timer,
|
2015-05-19 03:19:30 -05:00
|
|
|
ReservoirState& reservoir_state,
|
|
|
|
WellState& well_state);
|
2015-05-18 09:10:31 -05:00
|
|
|
|
2016-01-15 08:16:40 -06:00
|
|
|
/// Take a single forward step, after which the states will be modified
|
|
|
|
/// according to the physical model. This version allows for the
|
|
|
|
/// states passed as in/out arguments to be different from the initial
|
|
|
|
/// states.
|
2016-07-05 05:20:19 -05:00
|
|
|
/// \param[in] timer simulation timer
|
2016-01-15 08:16:40 -06:00
|
|
|
/// \param[in] initial_reservoir_state reservoir state variables at start of timestep
|
|
|
|
/// \param[in] initial_well_state well state variables at start of timestep
|
|
|
|
/// \param[in, out] reservoir_state reservoir state variables
|
|
|
|
/// \param[in, out] well_state well state variables
|
|
|
|
/// \return number of linear iterations used
|
|
|
|
int
|
2016-07-05 05:20:19 -05:00
|
|
|
step(const SimulatorTimerInterface& timer,
|
2016-01-15 08:16:40 -06:00
|
|
|
const ReservoirState& initial_reservoir_state,
|
|
|
|
const WellState& initial_well_state,
|
|
|
|
ReservoirState& reservoir_state,
|
|
|
|
WellState& well_state);
|
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
/// Number of nonlinear solver iterations used in all calls to step().
|
2016-06-30 08:31:07 -05:00
|
|
|
int nonlinearIterations() const;
|
2015-05-18 09:10:31 -05:00
|
|
|
|
|
|
|
/// Number of linear solver iterations used in all calls to step().
|
2016-06-30 08:31:07 -05:00
|
|
|
int linearIterations() const;
|
2015-05-18 09:10:31 -05:00
|
|
|
|
2016-06-20 02:31:25 -05:00
|
|
|
/// Number of well iterations used in all calls to step().
|
2016-06-30 08:31:07 -05:00
|
|
|
int wellIterations() const;
|
2016-06-20 02:31:25 -05:00
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
/// Number of nonlinear solver iterations used in the last call to step().
|
2016-06-30 08:31:07 -05:00
|
|
|
int nonlinearIterationsLastStep() const;
|
2015-05-21 02:09:16 -05:00
|
|
|
|
|
|
|
/// Number of linear solver iterations used in the last call to step().
|
2016-06-30 08:31:07 -05:00
|
|
|
int linearIterationsLastStep() const;
|
2015-05-21 02:09:16 -05:00
|
|
|
|
2016-07-19 21:44:12 -05:00
|
|
|
/// Number of well iterations used in all calls to step().
|
|
|
|
int wellIterationsLastStep() const;
|
|
|
|
|
2015-10-02 06:51:40 -05:00
|
|
|
/// Reference to physical model.
|
2015-10-31 06:32:54 -05:00
|
|
|
const PhysicalModel& model() const;
|
|
|
|
|
2015-12-01 07:16:11 -06:00
|
|
|
/// Mutable reference to physical model.
|
|
|
|
PhysicalModel& model();
|
|
|
|
|
2015-10-02 06:51:40 -05:00
|
|
|
/// Detect oscillation or stagnation in a given residual history.
|
|
|
|
void detectOscillations(const std::vector<std::vector<double>>& residual_history,
|
|
|
|
const int it, bool& oscillate, bool& stagnate) const;
|
|
|
|
|
|
|
|
/// Apply a stabilization to dx, depending on dxOld and relaxation parameters.
|
|
|
|
void stabilizeNonlinearUpdate(V& dx, V& dxOld, const double omega) const;
|
|
|
|
|
|
|
|
/// The greatest relaxation factor (i.e. smallest factor) allowed.
|
|
|
|
double relaxMax() const { return param_.relax_max_; }
|
|
|
|
|
|
|
|
/// The step-change size for the relaxation factor.
|
|
|
|
double relaxIncrement() const { return param_.relax_increment_; }
|
|
|
|
|
2015-11-24 08:21:37 -06:00
|
|
|
/// The relaxation type (DAMPEN or SOR).
|
|
|
|
enum RelaxType relaxType() const { return param_.relax_type_; }
|
|
|
|
|
|
|
|
/// The relaxation relative tolerance.
|
|
|
|
double relaxRelTol() const { return param_.relax_rel_tol_; }
|
|
|
|
|
|
|
|
/// The maximum number of nonlinear iterations allowed.
|
|
|
|
double maxIter() const { return param_.max_iter_; }
|
|
|
|
|
|
|
|
/// The minimum number of nonlinear iterations allowed.
|
|
|
|
double minIter() const { return param_.min_iter_; }
|
|
|
|
|
2015-05-18 09:10:31 -05:00
|
|
|
private:
|
2015-05-19 03:19:30 -05:00
|
|
|
// --------- Data members ---------
|
2015-05-21 01:54:56 -05:00
|
|
|
SolverParameters param_;
|
2015-05-29 07:57:49 -05:00
|
|
|
std::unique_ptr<PhysicalModel> model_;
|
2016-06-30 08:31:07 -05:00
|
|
|
int nonlinearIterations_;
|
|
|
|
int linearIterations_;
|
|
|
|
int wellIterations_;
|
|
|
|
int nonlinearIterationsLast_;
|
|
|
|
int linearIterationsLast_;
|
2016-07-19 21:44:12 -05:00
|
|
|
int wellIterationsLast_;
|
2015-05-18 07:18:03 -05:00
|
|
|
};
|
|
|
|
} // namespace Opm
|
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
#include "NonlinearSolver_impl.hpp"
|
2015-05-18 07:18:03 -05:00
|
|
|
|
2015-10-02 02:19:07 -05:00
|
|
|
#endif // OPM_NONLINEARSOLVER_HEADER_INCLUDED
|