2015-06-01 06:33:37 -05:00
|
|
|
/*
|
|
|
|
Copyright 2015 SINTEF ICT, Applied Mathematics.
|
|
|
|
Copyright 2015 IRIS AS
|
|
|
|
Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services
|
|
|
|
Copyright 2015 NTNU
|
|
|
|
Copyright 2015 Statoil AS
|
2015-10-07 10:48:27 -05:00
|
|
|
Copyright 2015 IRIS AS
|
2015-06-01 06:33:37 -05:00
|
|
|
|
|
|
|
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-06-16 04:35:01 -05:00
|
|
|
#ifndef OPM_NEWTONITERATIONBLACKOILINTERLEAVED_HEADER_INCLUDED
|
|
|
|
#define OPM_NEWTONITERATIONBLACKOILINTERLEAVED_HEADER_INCLUDED
|
|
|
|
|
2018-02-05 05:43:20 -06:00
|
|
|
#include <opm/autodiff/CPRPreconditioner.hpp>
|
2015-06-01 06:33:37 -05:00
|
|
|
#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
|
2018-01-30 06:25:35 -06:00
|
|
|
#include <opm/common/utility/parameters/ParameterGroup.hpp>
|
2015-06-16 04:35:01 -05:00
|
|
|
|
2015-10-08 07:04:05 -05:00
|
|
|
#include <array>
|
2015-06-01 06:33:37 -05:00
|
|
|
#include <memory>
|
|
|
|
|
|
|
|
namespace Opm
|
|
|
|
{
|
2015-10-09 05:03:58 -05:00
|
|
|
/// This class carries all parameters for the NewtonIterationBlackoilInterleaved class
|
|
|
|
struct NewtonIterationBlackoilInterleavedParameters
|
2018-02-05 05:43:20 -06:00
|
|
|
: public CPRParameter
|
2015-10-09 05:03:58 -05:00
|
|
|
{
|
|
|
|
double linear_solver_reduction_;
|
2017-05-18 10:49:50 -05:00
|
|
|
double ilu_relaxation_;
|
2015-10-09 05:03:58 -05:00
|
|
|
int linear_solver_maxiter_;
|
|
|
|
int linear_solver_restart_;
|
|
|
|
int linear_solver_verbosity_;
|
2017-06-02 07:52:49 -05:00
|
|
|
int ilu_fillin_level_;
|
2015-10-09 05:03:58 -05:00
|
|
|
bool newton_use_gmres_;
|
2016-05-11 08:13:52 -05:00
|
|
|
bool require_full_sparsity_pattern_;
|
2016-05-27 05:55:46 -05:00
|
|
|
bool ignoreConvergenceFailure_;
|
2016-11-21 10:18:24 -06:00
|
|
|
bool linear_solver_use_amg_;
|
2018-02-05 05:43:20 -06:00
|
|
|
bool use_cpr_;
|
2015-10-09 05:03:58 -05:00
|
|
|
|
|
|
|
NewtonIterationBlackoilInterleavedParameters() { reset(); }
|
|
|
|
// read values from parameter class
|
2017-04-28 08:36:25 -05:00
|
|
|
NewtonIterationBlackoilInterleavedParameters( const ParameterGroup& param )
|
2018-02-05 05:43:20 -06:00
|
|
|
: CPRParameter(param)
|
2015-10-09 05:03:58 -05:00
|
|
|
{
|
|
|
|
// set default parameters
|
|
|
|
reset();
|
|
|
|
|
|
|
|
// read parameters (using previsouly set default values)
|
|
|
|
newton_use_gmres_ = param.getDefault("newton_use_gmres", newton_use_gmres_ );
|
|
|
|
linear_solver_reduction_ = param.getDefault("linear_solver_reduction", linear_solver_reduction_ );
|
|
|
|
linear_solver_maxiter_ = param.getDefault("linear_solver_maxiter", linear_solver_maxiter_);
|
|
|
|
linear_solver_restart_ = param.getDefault("linear_solver_restart", linear_solver_restart_);
|
|
|
|
linear_solver_verbosity_ = param.getDefault("linear_solver_verbosity", linear_solver_verbosity_);
|
2016-05-11 08:13:52 -05:00
|
|
|
require_full_sparsity_pattern_ = param.getDefault("require_full_sparsity_pattern", require_full_sparsity_pattern_);
|
2016-05-27 05:55:46 -05:00
|
|
|
ignoreConvergenceFailure_ = param.getDefault("linear_solver_ignoreconvergencefailure", ignoreConvergenceFailure_);
|
2016-11-21 10:18:24 -06:00
|
|
|
linear_solver_use_amg_ = param.getDefault("linear_solver_use_amg", linear_solver_use_amg_ );
|
2017-05-18 10:49:50 -05:00
|
|
|
ilu_relaxation_ = param.getDefault("ilu_relaxation", ilu_relaxation_ );
|
2017-06-02 07:52:49 -05:00
|
|
|
ilu_fillin_level_ = param.getDefault("ilu_fillin_level", ilu_fillin_level_ );
|
2018-02-05 05:43:20 -06:00
|
|
|
|
|
|
|
// Check whether to use cpr approach
|
|
|
|
const std::string cprSolver = "cpr";
|
|
|
|
use_cpr_ = ( param.getDefault("solver_approach", std::string()) == cprSolver );
|
2015-10-09 05:03:58 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
// set default values
|
|
|
|
void reset()
|
|
|
|
{
|
2018-02-05 05:43:20 -06:00
|
|
|
use_cpr_ = false;
|
2015-10-09 05:03:58 -05:00
|
|
|
newton_use_gmres_ = false;
|
|
|
|
linear_solver_reduction_ = 1e-2;
|
2016-11-14 06:26:38 -06:00
|
|
|
linear_solver_maxiter_ = 150;
|
2015-10-09 05:03:58 -05:00
|
|
|
linear_solver_restart_ = 40;
|
|
|
|
linear_solver_verbosity_ = 0;
|
2016-05-11 08:13:52 -05:00
|
|
|
require_full_sparsity_pattern_ = false;
|
2016-05-27 05:55:46 -05:00
|
|
|
ignoreConvergenceFailure_ = false;
|
2016-11-21 10:18:24 -06:00
|
|
|
linear_solver_use_amg_ = false;
|
2017-06-02 07:52:49 -05:00
|
|
|
ilu_fillin_level_ = 0;
|
2017-05-18 10:49:50 -05:00
|
|
|
ilu_relaxation_ = 0.9;
|
2015-10-09 05:03:58 -05:00
|
|
|
}
|
|
|
|
};
|
|
|
|
|
2015-06-01 06:33:37 -05:00
|
|
|
|
|
|
|
/// This class solves the fully implicit black-oil system by
|
|
|
|
/// solving the reduced system (after eliminating well variables)
|
|
|
|
/// as a block-structured matrix (one block for all cell variables).
|
|
|
|
class NewtonIterationBlackoilInterleaved : public NewtonIterationBlackoilInterface
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
|
|
|
|
/// Construct a system solver.
|
2015-10-09 05:06:29 -05:00
|
|
|
/// \param[in] param parameters controlling the behaviour of the linear solvers
|
2015-06-01 06:33:37 -05:00
|
|
|
/// \param[in] parallelInformation In the case of a parallel run
|
2015-10-09 05:06:29 -05:00
|
|
|
/// with dune-istl the information about the parallelization.
|
2017-04-28 08:36:25 -05:00
|
|
|
NewtonIterationBlackoilInterleaved(const ParameterGroup& param,
|
2015-06-01 06:33:37 -05:00
|
|
|
const boost::any& parallelInformation=boost::any());
|
|
|
|
|
|
|
|
/// Solve the system of linear equations Ax = b, with A being the
|
|
|
|
/// combined derivative matrix of the residual and b
|
|
|
|
/// being the residual itself.
|
|
|
|
/// \param[in] residual residual object containing A and b.
|
|
|
|
/// \return the solution x
|
|
|
|
virtual SolutionVector computeNewtonIncrement(const LinearisedBlackoilResidual& residual) const;
|
|
|
|
|
|
|
|
/// \copydoc NewtonIterationBlackoilInterface::iterations
|
|
|
|
virtual int iterations () const { return iterations_; }
|
|
|
|
|
|
|
|
/// \copydoc NewtonIterationBlackoilInterface::parallelInformation
|
|
|
|
virtual const boost::any& parallelInformation() const;
|
|
|
|
|
|
|
|
private:
|
2015-10-08 07:04:05 -05:00
|
|
|
// max number of equations supported, increase if necessary
|
|
|
|
static const int maxNumberEquations_ = 6 ;
|
2015-06-01 06:33:37 -05:00
|
|
|
|
2016-02-08 08:52:32 -06:00
|
|
|
mutable std::array< std::unique_ptr< NewtonIterationBlackoilInterface >, maxNumberEquations_+1 > newtonIncrementDoublePrecision_;
|
2016-02-08 06:28:28 -06:00
|
|
|
mutable std::array< std::unique_ptr< NewtonIterationBlackoilInterface >, maxNumberEquations_+1 > newtonIncrementSinglePrecision_;
|
2015-10-09 05:03:58 -05:00
|
|
|
NewtonIterationBlackoilInterleavedParameters parameters_;
|
2015-06-01 06:33:37 -05:00
|
|
|
boost::any parallelInformation_;
|
2015-10-07 10:48:27 -05:00
|
|
|
mutable int iterations_;
|
2015-06-01 06:33:37 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
|
|
|
|
#endif // OPM_NEWTONITERATIONBLACKOILINTERLEAVED_HEADER_INCLUDED
|