21#ifdef DEAL_II_WITH_SUNDIALS
28# ifdef DEAL_II_WITH_TRILINOS
32# ifdef DEAL_II_WITH_PETSC
40# include <sundials/sundials_config.h>
42# include <kinsol/kinsol_direct.h>
43# include <sunlinsol/sunlinsol_dense.h>
44# include <sunmatrix/sunmatrix_dense.h>
56 template <
typename VectorType>
59 const unsigned int maximum_non_linear_iterations,
60 const double function_tolerance,
61 const double step_tolerance,
62 const bool no_init_setup,
63 const unsigned int maximum_setup_calls,
64 const double maximum_newton_step,
65 const double dq_relative_error,
66 const unsigned int maximum_beta_failures,
67 const unsigned int anderson_subspace_size)
69 , maximum_non_linear_iterations(maximum_non_linear_iterations)
70 , function_tolerance(function_tolerance)
71 , step_tolerance(step_tolerance)
72 , no_init_setup(no_init_setup)
73 , maximum_setup_calls(maximum_setup_calls)
74 , maximum_newton_step(maximum_newton_step)
75 , dq_relative_error(dq_relative_error)
76 , maximum_beta_failures(maximum_beta_failures)
77 , anderson_subspace_size(anderson_subspace_size)
82 template <
typename VectorType>
86 static std::string strategy_str(
"newton");
89 "Choose among newton|linesearch|fixed_point|picard",
91 "newton|linesearch|fixed_point|picard"));
92 prm.
add_action(
"Solution strategy", [&](
const std::string &value) {
93 if (value ==
"newton")
95 else if (value ==
"linesearch")
96 strategy = linesearch;
97 else if (value ==
"fixed_point")
98 strategy = fixed_point;
99 else if (value ==
"picard")
105 maximum_non_linear_iterations);
106 prm.
add_parameter(
"Function norm stopping tolerance", function_tolerance);
107 prm.
add_parameter(
"Scaled step stopping tolerance", step_tolerance);
112 maximum_setup_calls);
113 prm.
add_parameter(
"Maximum allowable scaled length of the Newton step",
114 maximum_newton_step);
115 prm.
add_parameter(
"Relative error for different quotient computation",
120 prm.
add_parameter(
"Maximum number of beta-condition failures",
121 maximum_beta_failures);
127 anderson_subspace_size);
134 template <
typename VectorType>
136 residual_or_iteration_callback(N_Vector yy, N_Vector FF,
void *user_data)
152 err = solver.
residual(*src_yy, *dst_FF);
165# if DEAL_II_SUNDIALS_VERSION_LT(4, 1, 0)
166 template <
typename VectorType>
189 template <
typename VectorType>
191 solve_with_jacobian_callback(KINMem
kinsol_mem,
197 KINSOL<VectorType> &solver =
198 *
static_cast<KINSOL<VectorType> *
>(
kinsol_mem->kin_user_data);
202 solver.reinit_vector(*src_ycur);
205 solver.reinit_vector(*src_fcur);
211 solver.reinit_vector(*src);
214 solver.reinit_vector(*dst);
218 int err = solver.solve_jacobian_system(*src_ycur, *src_fcur, *src, *dst);
231 template <
typename VectorType>
233 setup_jacobian_callback(N_Vector u,
243 const KINSOL<VectorType> &solver =
244 *
static_cast<const KINSOL<VectorType> *
>(user_data);
251 solver.reinit_vector(*ycur);
252 solver.reinit_vector(*fcur);
258 solver.setup_jacobian(*ycur, *fcur);
265 template <
typename VectorType>
276 const KINSOL<VectorType> &solver =
277 *
static_cast<const KINSOL<VectorType> *
>(LS->content);
281 if (solver.solve_with_jacobian)
289 solver.reinit_vector(*src_b);
290 solver.reinit_vector(*dst_x);
294 const int err = solver.solve_with_jacobian(*src_b, *dst_x, tol);
315 solver.reinit_vector(*src_b);
316 solver.reinit_vector(*dst_x);
326 solver.solve_jacobian_system(*src_ycur, *src_fcur, *src_b, *dst_x);
339 template <
typename VectorType>
356 template <
typename VectorType>
360 KINFree(&kinsol_mem);
362# ifdef DEAL_II_WITH_MPI
365 const int ierr = MPI_Comm_free(&communicator);
374 template <
typename VectorType>
378 unsigned int system_size = initial_guess_and_solution.size();
383# ifdef DEAL_II_WITH_MPI
386 const IndexSet is = initial_guess_and_solution.locally_owned_elements();
387 const unsigned int local_system_size = is.
n_elements();
390 N_VNew_Parallel(communicator, local_system_size, system_size);
392 u_scale = N_VNew_Parallel(communicator, local_system_size, system_size);
393 N_VConst_Parallel(1.e0, u_scale);
395 f_scale = N_VNew_Parallel(communicator, local_system_size, system_size);
396 N_VConst_Parallel(1.e0, f_scale);
403 "Trying to use a serial code with a parallel vector."));
404 solution = N_VNew_Serial(system_size);
405 u_scale = N_VNew_Serial(system_size);
406 N_VConst_Serial(1.e0, u_scale);
407 f_scale = N_VNew_Serial(system_size);
408 N_VConst_Serial(1.e0, f_scale);
411 if (get_solution_scaling)
412 copy(u_scale, get_solution_scaling());
414 if (get_function_scaling)
415 copy(f_scale, get_function_scaling());
417 copy(solution, initial_guess_and_solution);
420 KINFree(&kinsol_mem);
422 kinsol_mem = KINCreate();
425 KINInit(kinsol_mem, residual_or_iteration_callback<VectorType>, solution);
429 status = KINSetUserData(kinsol_mem,
static_cast<void *
>(
this));
432 status = KINSetNumMaxIters(kinsol_mem, data.maximum_non_linear_iterations);
435 status = KINSetFuncNormTol(kinsol_mem, data.function_tolerance);
438 status = KINSetScaledStepTol(kinsol_mem, data.step_tolerance);
441 status = KINSetMaxSetupCalls(kinsol_mem, data.maximum_setup_calls);
444 status = KINSetNoInitSetup(kinsol_mem, data.no_init_setup);
447 status = KINSetMaxNewtonStep(kinsol_mem, data.maximum_newton_step);
450 status = KINSetMaxBetaFails(kinsol_mem, data.maximum_beta_failures);
453 status = KINSetMAA(kinsol_mem, data.anderson_subspace_size);
456 status = KINSetRelErrFunc(kinsol_mem, data.dq_relative_error);
459 SUNMatrix J =
nullptr;
462 if (solve_jacobian_system ||
467# if DEAL_II_SUNDIALS_VERSION_LT(4, 1, 0)
468 auto KIN_mem =
static_cast<KINMem
>(kinsol_mem);
469 KIN_mem->kin_lsolve = solve_with_jacobian_callback<VectorType>;
472 KIN_mem->kin_lsetup = setup_jacobian_callback<VectorType>;
475# elif DEAL_II_SUNDIALS_VERSION_LT(5, 0, 0)
495 return SUNLINEARSOLVER_MATRIX_ITERATIVE;
501 LS->content =
nullptr;
513 LS->ops->solve = solve_with_jacobian_callback<VectorType>;
520 J = SUNMatNewEmpty();
523 J->ops->getid = [](SUNMatrix ) -> SUNMatrix_ID {
524 return SUNMATRIX_CUSTOM;
527 J->ops->destroy = [](SUNMatrix
A) {
530 A->content =
nullptr;
542 status = KINSetLinearSolver(kinsol_mem, LS, J);
551 KINSetJacFn(kinsol_mem, &setup_jacobian_callback<VectorType>);
558 J = SUNDenseMatrix(system_size, system_size);
559 LS = SUNDenseLinearSolver(u_scale, J);
560 status = KINDlsSetLinearSolver(kinsol_mem, LS, J);
564 if (data.strategy == AdditionalData::newton ||
565 data.strategy == AdditionalData::linesearch)
566 Assert(residual, ExcFunctionNotProvided(
"residual"));
568 if (data.strategy == AdditionalData::fixed_point ||
569 data.strategy == AdditionalData::picard)
570 Assert(iteration_function, ExcFunctionNotProvided(
"iteration_function"));
573 status = KINSol(kinsol_mem, solution, data.strategy, u_scale, f_scale);
576 copy(initial_guess_and_solution, solution);
579# ifdef DEAL_II_WITH_MPI
582 N_VDestroy_Parallel(solution);
583 N_VDestroy_Parallel(u_scale);
584 N_VDestroy_Parallel(f_scale);
589 N_VDestroy_Serial(solution);
590 N_VDestroy_Serial(u_scale);
591 N_VDestroy_Serial(f_scale);
595 status = KINGetNumNonlinSolvIters(kinsol_mem, &nniters);
600 KINFree(&kinsol_mem);
602 return static_cast<unsigned int>(nniters);
605 template <
typename VectorType>
609 reinit_vector = [](VectorType &) {
610 AssertThrow(
false, ExcFunctionNotProvided(
"reinit_vector"));
620# ifdef DEAL_II_WITH_MPI
622# ifdef DEAL_II_WITH_TRILINOS
627# ifdef DEAL_II_WITH_PETSC
628# ifndef PETSC_USE_COMPLEX
size_type n_elements() const
void add_parameter(const std::string &entry, ParameterType ¶meter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern(), const bool has_to_be_set=false)
void add_action(const std::string &entry, const std::function< void(const std::string &value)> &action)
void enter_subsection(const std::string &subsection)
AdditionalData(const SolutionStrategy &strategy=linesearch, const unsigned int maximum_non_linear_iterations=200, const double function_tolerance=0.0, const double step_tolerance=0.0, const bool no_init_setup=false, const unsigned int maximum_setup_calls=0, const double maximum_newton_step=0.0, const double dq_relative_error=0.0, const unsigned int maximum_beta_failures=0, const unsigned int anderson_subspace_size=0)
void add_parameters(ParameterHandler &prm)
std::function< int(const VectorType ¤t_u, const VectorType ¤t_f)> setup_jacobian
void set_functions_to_trigger_an_assert()
KINSOL(const AdditionalData &data=AdditionalData(), const MPI_Comm &mpi_comm=MPI_COMM_WORLD)
unsigned int solve(VectorType &initial_guess_and_solution)
GrowingVectorMemory< VectorType > mem
std::function< int(const VectorType &src, VectorType &dst)> residual
std::function< void(VectorType &)> reinit_vector
std::function< int(const VectorType &src, VectorType &dst)> iteration_function
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertNothrow(cond, exc)
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
#define AssertKINSOL(code)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SUNLinearSolver SUNLinSolNewEmpty()
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
void copy(const T *begin, const T *end, U *dest)