Stan Math Library
4.9.0
Automatic Differentiation
|
Eigen::Matrix< var, Eigen::Dynamic, 1 > stan::math::solve_newton_tol | ( | const F & | f, |
const T & | x, | ||
const double | scaling_step_size, | ||
const double | function_tolerance, | ||
const int64_t | max_num_steps, | ||
std::ostream *const | msgs, | ||
const T_Args &... | args | ||
) |
Return the solution to the specified system of algebraic equations given an initial guess, and parameters and data, which get passed into the algebraic system.
Use the KINSOL solver from the SUNDIALS suite.
The user can also specify the scaled step size, the function tolerance, and the maximum number of steps.
This overload handles var parameters.
The Jacobian (J_{xy}) (i.e., Jacobian of unknown (x) w.r.t. the parameter (y)) is calculated given the solution as follows. Since [ f(x, y) = 0, ] we have ((J_{pq}) being the Jacobian matrix (\tfrac {dq} {dq})) [
J_{fx} J_{xy} = J_{fy}. ] Let (eta) be the adjoint with respect to (x); then to calculate [ \eta J_{xy}, ] we solve [
F | type of equation system function. |
T | type of initial guess vector. |
Args | types of additional parameters to the equation system functor |
[in] | f | Functor that evaluated the system of equations. |
[in] | x | Vector of starting values. |
[in,out] | msgs | The print stream for warning messages. |
[in] | scaling_step_size | Scaled-step stopping tolerance. If a Newton step is smaller than the scaling step tolerance, the code breaks, assuming the solver is no longer making significant progress (i.e. is stuck) |
[in] | function_tolerance | determines whether roots are acceptable. |
[in] | max_num_steps | maximum number of function evaluations. |
[in] | args | Additional parameters to the equation system functor. |
<code>std::invalid_argument</code> | if x has size zero. |
<code>std::invalid_argument</code> | if x has non-finite elements. |
<code>std::invalid_argument</code> | if scaled_step_size is strictly negative. |
<code>std::invalid_argument</code> | if function_tolerance is strictly negative. |
<code>std::invalid_argument</code> | if max_num_steps is not positive. |
<code>std::domain_error | if solver exceeds max_num_steps. |
Definition at line 138 of file solve_newton.hpp.