diff --git a/Documentation/ProjectFile/prj/nonlinear_solvers/nonlinear_solver/t_prefix.md b/Documentation/ProjectFile/prj/nonlinear_solvers/nonlinear_solver/t_prefix.md new file mode 100644 index 0000000000000000000000000000000000000000..a0c291ee737ad8f5ef18e4adc8ac42964dee7c4d --- /dev/null +++ b/Documentation/ProjectFile/prj/nonlinear_solvers/nonlinear_solver/t_prefix.md @@ -0,0 +1,5 @@ +Sets the options prefix for the PETSc SNES solver. + +See the PETSc +[SNESSetOptionsPrefix](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetOptionsPrefix.html) +documentation for details. diff --git a/NumLib/ODESolver/NonlinearSolver.cpp b/NumLib/ODESolver/NonlinearSolver.cpp index 87892a7155e93eb679f1d664479da903102d02f5..910fcbd3cfc4d00a1e25d9a78449c7256ae41ee8 100644 --- a/NumLib/ODESolver/NonlinearSolver.cpp +++ b/NumLib/ODESolver/NonlinearSolver.cpp @@ -439,9 +439,13 @@ createNonlinearSolver(GlobalLinearSolver& linear_solver, #ifdef USE_PETSC if (boost::iequals(type, "PETScSNES")) { + auto prefix = + //! \ogs_file_param{prj__nonlinear_solvers__nonlinear_solver__prefix} + config.getConfigParameter<std::string>("prefix", ""); auto const tag = NonlinearSolverTag::Newton; using ConcreteNLS = PETScNonlinearSolver; - return std::make_pair(std::make_unique<ConcreteNLS>(linear_solver), + return std::make_pair(std::make_unique<ConcreteNLS>( + linear_solver, max_iter, std::move(prefix)), tag); } diff --git a/NumLib/ODESolver/PETScNonlinearSolver.cpp b/NumLib/ODESolver/PETScNonlinearSolver.cpp index eda70e72af287a6abf5780172d666def061ddfbf..dbbdd3b78a8032028b396f18f180d7eac62f35b4 100644 --- a/NumLib/ODESolver/PETScNonlinearSolver.cpp +++ b/NumLib/ODESolver/PETScNonlinearSolver.cpp @@ -78,9 +78,25 @@ PetscErrorCode updateJacobian(SNES /*snes*/, Vec /*x*/, Mat J, namespace NumLib { PETScNonlinearSolver::PETScNonlinearSolver( - GlobalLinearSolver& /*linear_solver*/) + GlobalLinearSolver& /*linear_solver*/, int const maxiter, + std::string prefix) { SNESCreate(PETSC_COMM_WORLD, &_snes_solver); + if (!prefix.empty()) + { + prefix = prefix + "_"; + SNESSetOptionsPrefix(_snes_solver, prefix.c_str()); + } + // force SNESSolve() to take at least one iteration regardless of the + // initial residual norm + SNESSetForceIteration(_snes_solver, PETSC_TRUE); + + // Set the maximum iterations. + PetscReal atol, rtol, stol; + PetscInt maxf; + SNESGetTolerances(_snes_solver, &atol, &rtol, &stol, nullptr, &maxf); + SNESSetTolerances(_snes_solver, atol, rtol, stol, maxiter, maxf); + SNESSetFromOptions(_snes_solver); #ifndef NDEBUG PetscOptionsView(nullptr, PETSC_VIEWER_STDOUT_WORLD); diff --git a/NumLib/ODESolver/PETScNonlinearSolver.h b/NumLib/ODESolver/PETScNonlinearSolver.h index 85b71434d4e2a7c1c7c24e93ab69e509a009343f..7f423e62facd20e5c91f108bc18993a2b7ef6260 100644 --- a/NumLib/ODESolver/PETScNonlinearSolver.h +++ b/NumLib/ODESolver/PETScNonlinearSolver.h @@ -30,8 +30,13 @@ public: /*! Constructs a new instance. * * \param linear_solver the linear solver used by this nonlinear solver. + * \param maxiter the maximum number of iterations used to solve the + * equation. + * \param prefix used to set the SNES options. */ - explicit PETScNonlinearSolver(GlobalLinearSolver& linear_solver); + PETScNonlinearSolver(GlobalLinearSolver& linear_solver, + int maxiter, + std::string prefix); //! Set the nonlinear equation system that will be solved. void setEquationSystem(System& eq, ConvergenceCriterion& conv_crit);