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>
2018-07-02 05:29:55 -05:00
# include <opm/autodiff/ParallelOverlappingILU0.hpp>
2015-06-16 04:35:01 -05:00
2018-06-21 05:14:17 -05:00
# include <ewoms/common/parametersystem.hh>
2015-10-08 07:04:05 -05:00
# include <array>
2015-06-01 06:33:37 -05:00
# include <memory>
2018-06-21 05:14:17 -05:00
BEGIN_PROPERTIES
NEW_TYPE_TAG ( FlowIstlSolverParams ) ;
NEW_PROP_TAG ( Scalar ) ;
2018-08-06 08:59:35 -05:00
NEW_PROP_TAG ( LinearSolverReduction ) ;
NEW_PROP_TAG ( IluRelaxation ) ;
NEW_PROP_TAG ( LinearSolverMaxIter ) ;
NEW_PROP_TAG ( LinearSolverRestart ) ;
2018-06-21 05:14:17 -05:00
NEW_PROP_TAG ( FlowLinearSolverVerbosity ) ;
2018-08-06 08:59:35 -05:00
NEW_PROP_TAG ( IluFillinLevel ) ;
NEW_PROP_TAG ( MIluVariant ) ;
NEW_PROP_TAG ( IluRedblack ) ;
NEW_PROP_TAG ( IluReorderSpheres ) ;
NEW_PROP_TAG ( UseGmres ) ;
NEW_PROP_TAG ( LinearSolverRequireFullSparsityPattern ) ;
NEW_PROP_TAG ( LinearSolverIgnoreConvergenceFailure ) ;
NEW_PROP_TAG ( UseAmg ) ;
NEW_PROP_TAG ( UseCpr ) ;
SET_SCALAR_PROP ( FlowIstlSolverParams , LinearSolverReduction , 1e-2 ) ;
SET_SCALAR_PROP ( FlowIstlSolverParams , IluRelaxation , 0.9 ) ;
SET_INT_PROP ( FlowIstlSolverParams , LinearSolverMaxIter , 150 ) ;
SET_INT_PROP ( FlowIstlSolverParams , LinearSolverRestart , 40 ) ;
2018-06-21 05:14:17 -05:00
SET_INT_PROP ( FlowIstlSolverParams , FlowLinearSolverVerbosity , 0 ) ;
2018-08-06 08:59:35 -05:00
SET_INT_PROP ( FlowIstlSolverParams , IluFillinLevel , 0 ) ;
SET_STRING_PROP ( FlowIstlSolverParams , MIluVariant , " ILU " ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , IluRedblack , false ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , IluReorderSpheres , false ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , UseGmres , false ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , LinearSolverRequireFullSparsityPattern , false ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , LinearSolverIgnoreConvergenceFailure , false ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , UseAmg , false ) ;
SET_BOOL_PROP ( FlowIstlSolverParams , UseCpr , false ) ;
2018-06-21 05:14:17 -05:00
END_PROPERTIES
2015-06-01 06:33:37 -05:00
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_ ;
2018-07-02 05:29:55 -05:00
Opm : : MILU_VARIANT ilu_milu_ ;
2018-06-12 12:38:55 -05:00
bool ilu_redblack_ ;
bool ilu_reorder_sphere_ ;
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
2018-06-21 05:14:17 -05:00
template < class TypeTag >
void init ( )
{
// TODO: these parameters have undocumented non-trivial dependencies
2018-08-06 08:59:35 -05:00
linear_solver_reduction_ = EWOMS_GET_PARAM ( TypeTag , double , LinearSolverReduction ) ;
ilu_relaxation_ = EWOMS_GET_PARAM ( TypeTag , double , IluRelaxation ) ;
linear_solver_maxiter_ = EWOMS_GET_PARAM ( TypeTag , int , LinearSolverMaxIter ) ;
linear_solver_restart_ = EWOMS_GET_PARAM ( TypeTag , int , LinearSolverRestart ) ;
2018-06-21 05:14:17 -05:00
linear_solver_verbosity_ = EWOMS_GET_PARAM ( TypeTag , int , FlowLinearSolverVerbosity ) ;
2018-08-06 08:59:35 -05:00
ilu_fillin_level_ = EWOMS_GET_PARAM ( TypeTag , int , IluFillinLevel ) ;
ilu_milu_ = convertString2Milu ( EWOMS_GET_PARAM ( TypeTag , std : : string , MIluVariant ) ) ;
ilu_redblack_ = EWOMS_GET_PARAM ( TypeTag , bool , IluRedblack ) ;
ilu_reorder_sphere_ = EWOMS_GET_PARAM ( TypeTag , bool , IluReorderSpheres ) ;
newton_use_gmres_ = EWOMS_GET_PARAM ( TypeTag , bool , UseGmres ) ;
require_full_sparsity_pattern_ = EWOMS_GET_PARAM ( TypeTag , bool , LinearSolverRequireFullSparsityPattern ) ;
ignoreConvergenceFailure_ = EWOMS_GET_PARAM ( TypeTag , bool , LinearSolverIgnoreConvergenceFailure ) ;
linear_solver_use_amg_ = EWOMS_GET_PARAM ( TypeTag , bool , UseAmg ) ;
use_cpr_ = EWOMS_GET_PARAM ( TypeTag , bool , UseCpr ) ;
2018-06-21 05:14:17 -05:00
}
template < class TypeTag >
static void registerParameters ( )
{
2018-08-06 08:59:35 -05:00
EWOMS_REGISTER_PARAM ( TypeTag , double , LinearSolverReduction , " The minimum reduction of the residual which the linear solver must achieve " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , double , IluRelaxation , " The relaxation factor of the linear solver's ILU preconditioner " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , int , LinearSolverMaxIter , " The maximum number of iterations of the linear solver " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , int , LinearSolverRestart , " The number of iterations after which GMRES is restarted " ) ;
2018-06-21 05:14:17 -05:00
EWOMS_REGISTER_PARAM ( TypeTag , int , FlowLinearSolverVerbosity , " The verbosity level of the linear solver (0: off, 2: all) " ) ;
2018-08-06 08:59:35 -05:00
EWOMS_REGISTER_PARAM ( TypeTag , int , IluFillinLevel , " The fill-in level of the linear solver's ILU preconditioner " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , std : : string , MIluVariant , " Specify which variant of the modified-ILU preconditioner ought to be used " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , IluRedblack , " Use red-black partioning for the ILU preconditioner " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , IluReorderSpheres , " Reorder the entries of the matrix in the ILU preconditioner " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , UseGmres , " Use GMRES as the linear solver " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , LinearSolverRequireFullSparsityPattern , " Produce the full sparsity pattern for the linear solver " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , LinearSolverIgnoreConvergenceFailure , " Continue with the simulation like nothing happened after the linear solver did not converge " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , UseAmg , " Use AMG as the linear solver's preconditioner " ) ;
EWOMS_REGISTER_PARAM ( TypeTag , bool , UseCpr , " Use CPR as the linear solver's preconditioner " ) ;
2018-06-21 05:14:17 -05:00
}
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-06-12 12:38:55 -05:00
ilu_redblack_ = param . getDefault ( " ilu_redblack " , cpr_ilu_redblack_ ) ;
ilu_reorder_sphere_ = param . getDefault ( " ilu_reorder_sphere " , cpr_ilu_reorder_sphere_ ) ;
2018-07-02 05:29:55 -05:00
std : : string milu ( " ILU " ) ;
ilu_milu_ = convertString2Milu ( param . getDefault ( " ilu_milu " , milu ) ) ;
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 ;
2018-07-02 05:29:55 -05:00
ilu_milu_ = MILU_VARIANT : : ILU ;
2018-06-12 12:38:55 -05:00
ilu_redblack_ = false ;
ilu_reorder_sphere_ = true ;
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