2018-11-15 04:07:47 -06:00
|
|
|
/*
|
2020-10-12 09:41:09 -05:00
|
|
|
Copyright 2015, 2020 SINTEF Digital, Mathematics and Cybernetics.
|
2018-11-15 04:07:47 -06:00
|
|
|
Copyright 2015 IRIS AS
|
|
|
|
Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services
|
|
|
|
Copyright 2015 NTNU
|
|
|
|
Copyright 2015 Statoil AS
|
|
|
|
|
|
|
|
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/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_FLOWLINEARSOLVERPARAMETERS_HEADER_INCLUDED
|
|
|
|
#define OPM_FLOWLINEARSOLVERPARAMETERS_HEADER_INCLUDED
|
|
|
|
|
2022-08-16 03:45:32 -05:00
|
|
|
#include <opm/simulators/linalg/MILU.hpp>
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
#include <opm/simulators/linalg/linalgparameters.hh>
|
2020-05-19 03:04:35 -05:00
|
|
|
#include <opm/simulators/linalg/linalgproperties.hh>
|
2019-09-16 03:58:20 -05:00
|
|
|
#include <opm/models/utils/parametersystem.hh>
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2018-11-12 04:03:54 -06:00
|
|
|
namespace Opm {
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2018-11-12 04:03:54 -06:00
|
|
|
template <class TypeTag>
|
2024-01-31 07:14:50 -06:00
|
|
|
class ISTLSolverBda;
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2023-08-09 07:54:10 -05:00
|
|
|
template <class TypeTag>
|
2024-01-31 07:14:50 -06:00
|
|
|
class ISTLSolver;
|
2018-11-12 04:03:54 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
}
|
2023-08-09 07:54:10 -05:00
|
|
|
|
2020-08-21 06:42:08 -05:00
|
|
|
namespace Opm::Properties {
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2020-08-27 03:30:29 -05:00
|
|
|
namespace TTag {
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 03:30:29 -05:00
|
|
|
struct FlowIstlSolverParams {};
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 03:30:29 -05:00
|
|
|
}
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
// Set the backend to be used.
|
|
|
|
template<class TypeTag>
|
|
|
|
struct LinearSolverBackend<TypeTag, TTag::FlowIstlSolverParams>
|
|
|
|
{
|
|
|
|
#if COMPILE_BDA_BRIDGE
|
|
|
|
using type = ISTLSolverBda<TypeTag>;
|
|
|
|
#else
|
|
|
|
using type = ISTLSolver<TypeTag>;
|
|
|
|
#endif
|
|
|
|
};
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
namespace Opm::Parameters {
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverReduction { using type = Properties::UndefinedProperty; };
|
2023-03-28 04:43:00 -05:00
|
|
|
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct RelaxedLinearSolverReduction { using type = Properties::UndefinedProperty; };
|
2023-08-09 07:54:10 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverMaxIter { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverRestart { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluRelaxation { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2023-06-09 06:46:45 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluFillinLevel { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct MiluVariant { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluRedblack { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluReorderSpheres { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct UseGmres { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverIgnoreConvergenceFailure { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct ScaleLinearSystem { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolver { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverPrintJsonDefinition { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2022-05-10 03:16:18 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct CprReuseSetup { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct CprReuseInterval { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct AcceleratorMode { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct BdaDeviceId { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct OpenclPlatformId { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-10-15 03:52:10 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct OpenclIluParallel { using type = Properties::UndefinedProperty; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverReduction<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{
|
2024-06-28 05:17:13 -05:00
|
|
|
using type = GetPropType<TypeTag, Properties::Scalar>;
|
2020-08-27 04:38:38 -05:00
|
|
|
static constexpr type value = 1e-2;
|
|
|
|
};
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct RelaxedLinearSolverReduction<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{
|
2024-06-28 05:17:13 -05:00
|
|
|
using type = GetPropType<TypeTag, Properties::Scalar>;
|
2023-03-28 04:43:00 -05:00
|
|
|
static constexpr type value = 1e-2;
|
2023-08-09 07:54:10 -05:00
|
|
|
};
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2023-03-28 04:43:00 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverMaxIter<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 200; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverRestart<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 40; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluRelaxation<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{
|
2024-06-28 05:17:13 -05:00
|
|
|
using type = GetPropType<TypeTag, Properties::Scalar>;
|
2023-06-09 06:46:45 -05:00
|
|
|
static constexpr type value = 0.9;
|
|
|
|
};
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2023-06-09 06:46:45 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluFillinLevel<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 0; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct MiluVariant<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr auto value = "ILU"; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluRedblack<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr bool value = false; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct IluReorderSpheres<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr bool value = false; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct UseGmres<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr bool value = false; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverIgnoreConvergenceFailure<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr bool value = false; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct ScaleLinearSystem<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr bool value = false; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolver<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-02-02 09:40:18 -06:00
|
|
|
{ static constexpr auto value = "cprw"; };
|
2024-06-28 05:17:13 -05:00
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct LinearSolverPrintJsonDefinition<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr auto value = true; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct CprReuseSetup<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 4; };
|
|
|
|
|
2022-05-10 03:16:18 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct CprReuseInterval<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 30; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct AcceleratorMode<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr auto value = "none"; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct BdaDeviceId<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 0; };
|
|
|
|
|
2020-08-27 04:38:38 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct OpenclPlatformId<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr int value = 0; };
|
|
|
|
|
2020-10-15 03:52:10 -05:00
|
|
|
template<class TypeTag>
|
2024-06-28 05:17:13 -05:00
|
|
|
struct OpenclIluParallel<TypeTag, Properties::TTag::FlowIstlSolverParams>
|
2024-06-28 05:17:13 -05:00
|
|
|
{ static constexpr bool value = true; }; // note: false should only be used in debug
|
2018-11-12 04:03:54 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
} // namespace Opm::Parameters
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
/// This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
|
|
|
|
struct FlowLinearSolverParameters
|
|
|
|
{
|
|
|
|
double linear_solver_reduction_;
|
|
|
|
double relaxed_linear_solver_reduction_;
|
|
|
|
int linear_solver_maxiter_;
|
|
|
|
int linear_solver_restart_;
|
|
|
|
int linear_solver_verbosity_;
|
|
|
|
double ilu_relaxation_;
|
|
|
|
int ilu_fillin_level_;
|
|
|
|
MILU_VARIANT ilu_milu_;
|
|
|
|
bool ilu_redblack_;
|
|
|
|
bool ilu_reorder_sphere_;
|
|
|
|
bool newton_use_gmres_;
|
|
|
|
bool ignoreConvergenceFailure_;
|
|
|
|
bool scale_linear_system_;
|
|
|
|
std::string linsolver_;
|
|
|
|
bool linear_solver_print_json_definition_;
|
|
|
|
int cpr_reuse_setup_;
|
|
|
|
int cpr_reuse_interval_;
|
|
|
|
std::string accelerator_mode_;
|
|
|
|
int bda_device_id_;
|
|
|
|
int opencl_platform_id_;
|
|
|
|
bool opencl_ilu_parallel_;
|
|
|
|
|
|
|
|
template <class TypeTag>
|
|
|
|
void init(bool cprRequestedInDataFile)
|
|
|
|
{
|
|
|
|
// TODO: these parameters have undocumented non-trivial dependencies
|
2024-06-28 05:17:13 -05:00
|
|
|
linear_solver_reduction_ = Parameters::get<TypeTag, Parameters::LinearSolverReduction>();
|
|
|
|
relaxed_linear_solver_reduction_ = Parameters::get<TypeTag, Parameters::RelaxedLinearSolverReduction>();
|
|
|
|
linear_solver_maxiter_ = Parameters::get<TypeTag, Parameters::LinearSolverMaxIter>();
|
|
|
|
linear_solver_restart_ = Parameters::get<TypeTag, Parameters::LinearSolverRestart>();
|
2024-07-06 03:22:47 -05:00
|
|
|
linear_solver_verbosity_ = Parameters::Get<Parameters::LinearSolverVerbosity>();
|
2024-06-28 05:17:13 -05:00
|
|
|
ilu_relaxation_ = Parameters::get<TypeTag, Parameters::IluRelaxation>();
|
|
|
|
ilu_fillin_level_ = Parameters::get<TypeTag, Parameters::IluFillinLevel>();
|
|
|
|
ilu_milu_ = convertString2Milu(Parameters::get<TypeTag, Parameters::MiluVariant>());
|
|
|
|
ilu_redblack_ = Parameters::get<TypeTag, Parameters::IluRedblack>();
|
|
|
|
ilu_reorder_sphere_ = Parameters::get<TypeTag, Parameters::IluReorderSpheres>();
|
|
|
|
newton_use_gmres_ = Parameters::get<TypeTag, Parameters::UseGmres>();
|
|
|
|
ignoreConvergenceFailure_ = Parameters::get<TypeTag, Parameters::LinearSolverIgnoreConvergenceFailure>();
|
|
|
|
scale_linear_system_ = Parameters::get<TypeTag, Parameters::ScaleLinearSystem>();
|
|
|
|
linsolver_ = Parameters::get<TypeTag, Parameters::LinearSolver>();
|
|
|
|
linear_solver_print_json_definition_ = Parameters::get<TypeTag, Parameters::LinearSolverPrintJsonDefinition>();
|
|
|
|
cpr_reuse_setup_ = Parameters::get<TypeTag, Parameters::CprReuseSetup>();
|
|
|
|
cpr_reuse_interval_ = Parameters::get<TypeTag, Parameters::CprReuseInterval>();
|
|
|
|
|
|
|
|
if (!Parameters::isSet<TypeTag, Parameters::LinearSolver>() && cprRequestedInDataFile) {
|
2024-06-28 05:17:13 -05:00
|
|
|
linsolver_ = "cpr";
|
|
|
|
} else {
|
2024-06-28 05:17:13 -05:00
|
|
|
linsolver_ = Parameters::get<TypeTag, Parameters::LinearSolver>();
|
2018-11-15 04:07:47 -06:00
|
|
|
}
|
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
accelerator_mode_ = Parameters::get<TypeTag, Parameters::AcceleratorMode>();
|
|
|
|
bda_device_id_ = Parameters::get<TypeTag, Parameters::BdaDeviceId>();
|
|
|
|
opencl_platform_id_ = Parameters::get<TypeTag, Parameters::OpenclPlatformId>();
|
|
|
|
opencl_ilu_parallel_ = Parameters::get<TypeTag, Parameters::OpenclIluParallel>();
|
2024-06-28 05:17:13 -05:00
|
|
|
}
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
template <class TypeTag>
|
|
|
|
static void registerParameters()
|
|
|
|
{
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::LinearSolverReduction>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The minimum reduction of the residual which the linear solver must achieve");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::RelaxedLinearSolverReduction>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The minimum reduction of the residual which the linear solver need to "
|
|
|
|
"achieve for the solution to be accepted");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::LinearSolverMaxIter>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The maximum number of iterations of the linear solver");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::LinearSolverRestart>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The number of iterations after which GMRES is restarted");
|
2024-07-06 03:22:47 -05:00
|
|
|
Parameters::Register<Parameters::LinearSolverVerbosity>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The verbosity level of the linear solver (0: off, 2: all)");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::IluRelaxation>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The relaxation factor of the linear solver's ILU preconditioner");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::IluFillinLevel>
|
2024-06-28 05:17:13 -05:00
|
|
|
("The fill-in level of the linear solver's ILU preconditioner");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::MiluVariant>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Specify which variant of the modified-ILU preconditioner ought to be used. "
|
|
|
|
"Possible variants are: ILU (default, plain ILU), "
|
|
|
|
"MILU_1 (lump diagonal with dropped row entries), "
|
|
|
|
"MILU_2 (lump diagonal with the sum of the absolute values of the dropped row entries), "
|
|
|
|
"MILU_3 (if diagonal is positive add sum of dropped row entries, otherwise subtract them), "
|
|
|
|
"MILU_4 (if diagonal is positive add sum of dropped row entries, otherwise do nothing");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::IluRedblack>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Use red-black partitioning for the ILU preconditioner");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::IluReorderSpheres>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Whether to reorder the entries of the matrix in the red-black "
|
|
|
|
"ILU preconditioner in spheres starting at an edge. "
|
|
|
|
"If false the original ordering is preserved in each color. "
|
|
|
|
"Otherwise why try to ensure D4 ordering (in a 2D structured grid, "
|
|
|
|
"the diagonal elements are consecutive).");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::UseGmres>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Use GMRES as the linear solver");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::LinearSolverIgnoreConvergenceFailure>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Continue with the simulation like nothing happened "
|
|
|
|
"after the linear solver did not converge");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::ScaleLinearSystem>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Scale linear system according to equation scale and primary variable types");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::LinearSolver>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Configuration of solver. Valid options are: ilu0 (default), "
|
|
|
|
"dilu, cprw, cpr (an alias for cprw), cpr_quasiimpes, "
|
|
|
|
"cpr_trueimpes, cpr_trueimpesanalytic, amg or hybrid (experimental). "
|
|
|
|
"Alternatively, you can request a configuration to be read from a "
|
|
|
|
"JSON file by giving the filename here, ending with '.json.'");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::LinearSolverPrintJsonDefinition>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Write the JSON definition of the linear solver setup to the DBG file.");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::CprReuseSetup>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Reuse preconditioner setup. Valid options are "
|
|
|
|
"0: recreate the preconditioner for every linear solve, "
|
|
|
|
"1: recreate once every timestep, "
|
|
|
|
"2: recreate if last linear solve took more than 10 iterations, "
|
|
|
|
"3: never recreate, "
|
|
|
|
"4: recreated every CprReuseInterval");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::CprReuseInterval>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Reuse preconditioner interval. Used when CprReuseSetup is set to 4, "
|
|
|
|
"then the preconditioner will be fully recreated instead of reused "
|
|
|
|
"every N linear solve, where N is this parameter.");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::AcceleratorMode>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Choose a linear solver, usage: "
|
|
|
|
"'--accelerator-mode=[none|cusparse|opencl|amgcl|rocalution|rocsparse]'");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::BdaDeviceId>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Choose device ID for cusparseSolver or openclSolver, "
|
|
|
|
"use 'nvidia-smi' or 'clinfo' to determine valid IDs");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::OpenclPlatformId>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Choose platform ID for openclSolver, use 'clinfo' "
|
|
|
|
"to determine valid platform IDs");
|
2024-06-28 05:17:13 -05:00
|
|
|
Parameters::registerParam<TypeTag, Parameters::OpenclIluParallel>
|
2024-06-28 05:17:13 -05:00
|
|
|
("Parallelize ILU decomposition and application on GPU");
|
2024-07-06 03:22:47 -05:00
|
|
|
|
|
|
|
Parameters::SetDefault<Parameters::LinearSolverVerbosity>(0);
|
2024-06-28 05:17:13 -05:00
|
|
|
}
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
FlowLinearSolverParameters() { reset(); }
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
// set default values
|
|
|
|
void reset()
|
|
|
|
{
|
|
|
|
relaxed_linear_solver_reduction_ = 1e-2;
|
|
|
|
linear_solver_reduction_ = 1e-2;
|
|
|
|
linear_solver_maxiter_ = 200;
|
|
|
|
linear_solver_restart_ = 40;
|
|
|
|
linear_solver_verbosity_ = 0;
|
|
|
|
ilu_relaxation_ = 0.9;
|
|
|
|
ilu_fillin_level_ = 0;
|
|
|
|
ilu_milu_ = MILU_VARIANT::ILU;
|
|
|
|
ilu_redblack_ = false;
|
|
|
|
ilu_reorder_sphere_ = false;
|
|
|
|
newton_use_gmres_ = false;
|
|
|
|
ignoreConvergenceFailure_ = false;
|
|
|
|
scale_linear_system_ = false;
|
2024-05-15 06:05:26 -05:00
|
|
|
linsolver_ = "cprw";
|
2024-06-28 05:17:13 -05:00
|
|
|
linear_solver_print_json_definition_ = true;
|
|
|
|
cpr_reuse_setup_ = 4;
|
|
|
|
cpr_reuse_interval_ = 30;
|
|
|
|
accelerator_mode_ = "none";
|
|
|
|
bda_device_id_ = 0;
|
|
|
|
opencl_platform_id_ = 0;
|
|
|
|
opencl_ilu_parallel_ = true;
|
|
|
|
}
|
|
|
|
};
|
2018-11-15 04:07:47 -06:00
|
|
|
|
2024-06-28 05:17:13 -05:00
|
|
|
} // namespace Opm
|
2018-11-15 04:07:47 -06:00
|
|
|
|
|
|
|
#endif // OPM_FLOWLINEARSOLVERPARAMETERS_HEADER_INCLUDED
|