/* Copyright 2015, 2020 SINTEF Digital, Mathematics and Cybernetics. 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 . */ #include #include #include #include namespace Opm { void FlowLinearSolverParameters::init(bool cprRequestedInDataFile) { // TODO: these parameters have undocumented non-trivial dependencies relaxed_linear_solver_reduction_ = Parameters::Get(); linear_solver_restart_ = Parameters::Get(); linear_solver_verbosity_ = Parameters::Get(); ilu_relaxation_ = Parameters::Get(); ilu_fillin_level_ = Parameters::Get(); ilu_milu_ = convertString2Milu(Parameters::Get()); ilu_redblack_ = Parameters::Get(); ilu_reorder_sphere_ = Parameters::Get(); newton_use_gmres_ = Parameters::Get(); ignoreConvergenceFailure_ = Parameters::Get(); scale_linear_system_ = Parameters::Get(); linsolver_ = Parameters::Get(); linear_solver_print_json_definition_ = Parameters::Get(); cpr_reuse_setup_ = Parameters::Get(); cpr_reuse_interval_ = Parameters::Get(); if (!Parameters::IsSet() && cprRequestedInDataFile) { linsolver_ = "cpr"; } else { linsolver_ = Parameters::Get(); } // if local solver from nldd, use nldd local linear solver parameters if (is_nldd_local_solver_) { linsolver_ = Parameters::Get(); linear_solver_reduction_ = Parameters::Get(); linear_solver_maxiter_ = Parameters::Get(); } else { linsolver_ = Parameters::Get(); linear_solver_reduction_ = Parameters::Get(); linear_solver_maxiter_ = Parameters::Get(); } accelerator_mode_ = Parameters::Get(); gpu_device_id_ = Parameters::Get(); opencl_platform_id_ = Parameters::Get(); opencl_ilu_parallel_ = Parameters::Get(); } void FlowLinearSolverParameters::registerParameters() { Parameters::Register ("The minimum reduction of the residual which the linear solver must achieve"); Parameters::Register ("The minimum reduction of the residual which the NLDD local linear solver must achieve"); Parameters::Register ("The minimum reduction of the residual which the linear solver need to " "achieve for the solution to be accepted"); Parameters::Register ("The maximum number of iterations of the linear solver"); Parameters::Register ("The maximum number of iterations of the NLDD local linear solver"); Parameters::Register ("The number of iterations after which GMRES is restarted"); Parameters::Register ("The verbosity level of the linear solver (0: off, 2: all)"); Parameters::Register ("The relaxation factor of the linear solver's ILU preconditioner"); Parameters::Register ("The fill-in level of the linear solver's ILU preconditioner"); Parameters::Register ("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"); Parameters::Register ("Use red-black partitioning for the ILU preconditioner"); Parameters::Register ("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)."); Parameters::Register ("Use GMRES as the linear solver"); Parameters::Register ("Continue with the simulation like nothing happened " "after the linear solver did not converge"); Parameters::Register ("Scale linear system according to equation scale and primary variable types"); Parameters::Register ("Configuration of solver. Valid options are: cprw (default), " "ilu0, dilu, 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.'"); Parameters::Register ("Configuration of NLDD local linear solver. Valid options are: ilu0 (default), " "dilu, cpr_quasiimpes and amg. " "Alternatively, you can request a configuration to be read from a " "JSON file by giving the filename here, ending with '.json.'"); Parameters::Register ("Write the JSON definition of the linear solver setup to the DBG file."); Parameters::Register ("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"); Parameters::Register ("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."); Parameters::Register ("Choose a linear solver, usage: " "'--accelerator-mode=[none|cusparse|opencl|amgcl|rocalution|rocsparse]'"); Parameters::Register ("Choose device ID for cusparseSolver or openclSolver, " "use 'nvidia-smi' or 'clinfo' to determine valid IDs"); Parameters::Register ("Choose platform ID for openclSolver, use 'clinfo' " "to determine valid platform IDs"); Parameters::Register ("Parallelize ILU decomposition and application on GPU"); Parameters::SetDefault(0); } void FlowLinearSolverParameters::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; is_nldd_local_solver_ = false; linsolver_ = "cprw"; linear_solver_print_json_definition_ = true; cpr_reuse_setup_ = 4; cpr_reuse_interval_ = 30; accelerator_mode_ = "none"; gpu_device_id_ = 0; opencl_platform_id_ = 0; opencl_ilu_parallel_ = true; } } // namespace Opm