Automatic Differentiation
 
Loading...
Searching...
No Matches
kinsol_solve.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUNCTOR_KINSOL_SOLVE_HPP
2#define STAN_MATH_REV_FUNCTOR_KINSOL_SOLVE_HPP
3
9#include <sundials/sundials_context.h>
10#include <kinsol/kinsol.h>
11#include <sunmatrix/sunmatrix_dense.h>
12#include <sunlinsol/sunlinsol_dense.h>
13#include <nvector/nvector_serial.h>
14#include <vector>
15
16namespace stan {
17namespace math {
18
60template <typename F1, typename... Args>
61Eigen::VectorXd kinsol_solve(const F1& f, const Eigen::VectorXd& x,
62 const double scaling_step_tol, // = 1e-3
63 const double function_tolerance, // = 1e-6
64 const int64_t max_num_steps, // = 200
65 const bool custom_jacobian, // = 1
66 const int steps_eval_jacobian, // = 10
67 const int global_line_search, // = KIN_LINESEARCH
68 std::ostream* const msgs, const Args&... args) {
69 int N = x.size();
70 typedef kinsol_system_data<F1, Args...> system_data;
71 system_data kinsol_data(f, x, msgs, args...);
72
73 CHECK_KINSOL_CALL(KINInit(kinsol_data.kinsol_memory_,
74 &system_data::kinsol_f_system, kinsol_data.nv_x_));
75
76 N_Vector scaling = N_VNew_Serial(N, kinsol_data.sundials_context_);
77 N_Vector nv_x = N_VNew_Serial(N, kinsol_data.sundials_context_);
78 Eigen::VectorXd x_solution(N);
79
80 try {
81 N_VConst_Serial(1.0, scaling); // no scaling
82
84 KINSetNumMaxIters(kinsol_data.kinsol_memory_, max_num_steps));
86 KINSetFuncNormTol(kinsol_data.kinsol_memory_, function_tolerance));
88 KINSetScaledStepTol(kinsol_data.kinsol_memory_, scaling_step_tol));
90 KINSetMaxSetupCalls(kinsol_data.kinsol_memory_, steps_eval_jacobian));
91
92 // CHECK
93 // The default value is 1000 * ||u_0||_D where ||u_0|| is the initial guess.
94 // So we run into issues if ||u_0|| = 0.
95 // If the norm is non-zero, use kinsol's default (accessed with 0),
96 // else use the dimension of x -- CHECK - find optimal length.
97 double max_newton_step = (x.norm() == 0) ? x.size() : 0;
99 KINSetMaxNewtonStep(kinsol_data.kinsol_memory_, max_newton_step));
100 CHECK_KINSOL_CALL(KINSetUserData(kinsol_data.kinsol_memory_,
101 static_cast<void*>(&kinsol_data)));
102
103 // construct Linear solver
104 CHECK_KINSOL_CALL(KINSetLinearSolver(kinsol_data.kinsol_memory_,
105 kinsol_data.LS_, kinsol_data.J_));
106
107 if (custom_jacobian)
108 CHECK_KINSOL_CALL(KINSetJacFn(kinsol_data.kinsol_memory_,
109 &system_data::kinsol_jacobian));
110
111 for (int i = 0; i < N; i++)
112 NV_Ith_S(nv_x, i) = x(i);
113
114 kinsol_check(KINSol(kinsol_data.kinsol_memory_, nv_x, global_line_search,
115 scaling, scaling),
116 "KINSol", max_num_steps);
117
118 for (int i = 0; i < N; i++)
119 x_solution(i) = NV_Ith_S(nv_x, i);
120 } catch (const std::exception& e) {
121 N_VDestroy(nv_x);
122 N_VDestroy(scaling);
123 throw;
124 }
125
126 N_VDestroy(nv_x);
127 N_VDestroy(scaling);
128
129 return x_solution;
130}
131
132} // namespace math
133} // namespace stan
134#endif
#define CHECK_KINSOL_CALL(call)
KINSOL algebraic system data holder.
static constexpr double e()
Return the base of the natural logarithm.
Definition constants.hpp:20
void kinsol_check(int flag, const char *func_name)
Throws an exception message when the functions in KINSOL fails.
Eigen::VectorXd kinsol_solve(const F1 &f, const Eigen::VectorXd &x, const double scaling_step_tol, const double function_tolerance, const int64_t max_num_steps, const bool custom_jacobian, const int steps_eval_jacobian, const int global_line_search, std::ostream *const msgs, const Args &... args)
Return the solution to the specified algebraic system, given an initial guess.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...