From: Wolfgang Bangerth Date: Mon, 12 Jun 2023 16:33:14 +0000 (-0600) Subject: Enforce callback return value standard for the new NonlinearSolverSelector class. X-Git-Tag: v9.5.0-rc1~125^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=d957d865a4aa18dc0aacf90a5c8bdcd8564ab4a5;p=dealii.git Enforce callback return value standard for the new NonlinearSolverSelector class. --- diff --git a/include/deal.II/numerics/nonlinear.h b/include/deal.II/numerics/nonlinear.h index 77d007599d..cf4dc0c441 100644 --- a/include/deal.II/numerics/nonlinear.h +++ b/include/deal.II/numerics/nonlinear.h @@ -288,9 +288,16 @@ public: * A function object that users should supply and that is intended to * compute the residual `dst = F(src)`. * - * This function should return an int for either failure or success. + * @note This variable represents a + * @ref GlossUserProvidedCallBack "user provided callback". + * See there for a description of how to deal with errors and other + * requirements and conventions. Some of the underlying packages + * used by this class can deal with recoverable exceptions, whereas + * others cannot. As a consequence, if a callback + * throws an exception of type RecoverableUserCallbackError, then this + * exception may or may not be treated like any other exception. */ - std::function residual; + std::function residual; /** * A function object that users may supply and that is intended to @@ -308,8 +315,17 @@ public: * approximate Jacobian matrix $L$. * * @param current_u Current value of $u$ + * + * @note This variable represents a + * @ref GlossUserProvidedCallBack "user provided callback". + * See there for a description of how to deal with errors and other + * requirements and conventions. Some of the underlying packages + * used by this class can deal with recoverable exceptions, whereas + * others cannot. As a consequence, if a callback + * throws an exception of type RecoverableUserCallbackError, then this + * exception may or may not be treated like any other exception. */ - std::function setup_jacobian; + std::function setup_jacobian; /** * A function object that users may supply and that is intended to solve @@ -322,9 +338,18 @@ public: * @param[out] dst The solution of $J^{-1} * \texttt{src}$. * @param[in] tolerance The tolerance with which to solve the linear system * of equations. + * + * @note This variable represents a + * @ref GlossUserProvidedCallBack "user provided callback". + * See there for a description of how to deal with errors and other + * requirements and conventions. Some of the underlying packages + * used by this class can deal with recoverable exceptions, whereas + * others cannot. As a consequence, if a callback + * throws an exception of type RecoverableUserCallbackError, then this + * exception may or may not be treated like any other exception. */ std::function< - int(const VectorType &rhs, VectorType &dst, const double tolerance)> + void(const VectorType &rhs, VectorType &dst, const double tolerance)> solve_with_jacobian; private: @@ -567,12 +592,12 @@ NonlinearSolverSelector::solve_with_petsc( nonlinear_solver.solve_with_jacobian = [&](const PETScWrappers::MPI::Vector &src, - PETScWrappers::MPI::Vector & dst) -> int { - // PETSc does not gives a tolerance, so we have to choose something - // reasonable to provide to the user: - const double tolerance = 1e-6; - return solve_with_jacobian(src, dst, tolerance); - }; + PETScWrappers::MPI::Vector & dst) { + // PETSc does not gives a tolerance, so we have to choose something + // reasonable to provide to the user: + const double tolerance = 1e-6; + solve_with_jacobian(src, dst, tolerance); + }; nonlinear_solver.solve(initial_guess_and_solution); } diff --git a/tests/numerics/nonlinear_solver_selector_01.cc b/tests/numerics/nonlinear_solver_selector_01.cc index e3126ba98b..958cf1aa38 100644 --- a/tests/numerics/nonlinear_solver_selector_01.cc +++ b/tests/numerics/nonlinear_solver_selector_01.cc @@ -369,22 +369,16 @@ namespace nonlinear_solver_selector_test nonlinear_solver.residual = [&](const Vector &evaluation_point, Vector & residual) { compute_residual(evaluation_point, residual); - - return 0; }; nonlinear_solver.setup_jacobian = [&](const Vector ¤t_u) { compute_and_factorize_jacobian(current_u); - - return 0; }; nonlinear_solver.solve_with_jacobian = [&](const Vector &rhs, Vector & dst, const double tolerance) { this->solve(rhs, dst, tolerance); - - return 0; }; nonlinear_solver.solve(current_solution); @@ -402,6 +396,4 @@ main() MinimalSurfaceProblem<2> laplace_problem_2d; laplace_problem_2d.run(); - - return 0; } diff --git a/tests/numerics/nonlinear_solver_selector_03.cc b/tests/numerics/nonlinear_solver_selector_03.cc index 9026f57231..1968ca6171 100644 --- a/tests/numerics/nonlinear_solver_selector_03.cc +++ b/tests/numerics/nonlinear_solver_selector_03.cc @@ -451,22 +451,16 @@ namespace MPI_nonlinear_solver_selector_test nonlinear_solver.residual = [&](const LA::MPI::Vector &evaluation_point, LA::MPI::Vector & residual) { compute_residual(evaluation_point, residual); - - return 0; }; nonlinear_solver.setup_jacobian = [&](const LA::MPI::Vector ¤t_u) { compute_and_factorize_jacobian(current_u); - - return 0; }; nonlinear_solver.solve_with_jacobian = [&](const LA::MPI::Vector &rhs, LA::MPI::Vector & dst, const double tolerance) { this->solve(rhs, dst, tolerance); - - return 0; }; nonlinear_solver.solve(current_solution); @@ -487,6 +481,4 @@ main(int argc, char *argv[]) MinimalSurfaceProblem<2> laplace_problem_2d; laplace_problem_2d.run(); - - return 0; }