Automatic Differentiation
 
Loading...
Searching...
No Matches
algebra_system.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUNCTOR_ALGEBRA_SYSTEM_HPP
2#define STAN_MATH_REV_FUNCTOR_ALGEBRA_SYSTEM_HPP
3
8#include <iostream>
9#include <string>
10#include <vector>
11
12namespace stan {
13namespace math {
14
23template <typename T, int NX = Eigen::Dynamic, int NY = Eigen::Dynamic>
25 const int m_inputs, m_values;
26
28
30
31 int inputs() const { return m_inputs; }
32 int values() const { return m_values; }
33};
34
44template <typename S>
47 S fs_;
48
50 Eigen::MatrixXd J_;
51
52 explicit hybrj_functor_solver(const S& fs) : fs_(fs) {}
53
61 int operator()(const Eigen::VectorXd& iv, Eigen::VectorXd& fvec) {
62 jacobian(fs_, iv, fvec, J_);
63 return 0;
64 }
65
72 int df(const Eigen::VectorXd& iv, Eigen::MatrixXd& fjac) const {
73 fjac = J_;
74 return 0;
75 }
76
83 Eigen::MatrixXd get_jacobian(const Eigen::VectorXd& iv) {
84 Eigen::VectorXd fvec;
85 jacobian(fs_, iv, fvec, J_);
86 return J_;
87 }
88
95 Eigen::VectorXd get_value(const Eigen::VectorXd& iv) const { return fs_(iv); }
96};
97
98// TODO(jgaeb): Remove this when the chain method of the fixed point solver is
99// updated.
100template <typename T1, typename T2>
101void algebra_solver_check(const Eigen::Matrix<T1, Eigen::Dynamic, 1>& x,
102 const Eigen::Matrix<T2, Eigen::Dynamic, 1> y,
103 const std::vector<double>& dat,
104 const std::vector<int>& dat_int,
105 double function_tolerance,
106 long int max_num_steps) { // NOLINT(runtime/int)
107 check_nonzero_size("algebra_solver", "initial guess", x);
108 check_finite("algebra_solver", "initial guess", x);
109 check_finite("algebra_solver", "parameter vector", y);
110 check_finite("algebra_solver", "continuous data", dat);
111 check_finite("algebra_solver", "integer data", dat_int);
112
113 check_nonnegative("algebra_solver", "function_tolerance", function_tolerance);
114 check_positive("algebra_solver", "max_num_steps", max_num_steps);
115}
116
117} // namespace math
118} // namespace stan
119
120#endif
void jacobian(const F &f, const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, Eigen::Matrix< T, Eigen::Dynamic, 1 > &fx, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &J)
Definition jacobian.hpp:11
void check_nonnegative(const char *function, const char *name, const T_y &y)
Check if y is non-negative.
void algebra_solver_check(const Eigen::Matrix< T1, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< T2, Eigen::Dynamic, 1 > y, const std::vector< double > &dat, const std::vector< int > &dat_int, double function_tolerance, long int max_num_steps)
void check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
void check_nonzero_size(const char *function, const char *name, const T_y &y)
Check if the specified matrix/vector is of non-zero size.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
Definition fvar.hpp:9
int operator()(const Eigen::VectorXd &iv, Eigen::VectorXd &fvec)
Computes the value the algebraic function, f, when pluging in the independent variables,...
Eigen::MatrixXd J_
Jacobian of algebraic function wrt unknowns.
Eigen::MatrixXd get_jacobian(const Eigen::VectorXd &iv)
Performs the same task as the operator(), but returns the Jacobian, instead of saving it inside an ar...
Eigen::VectorXd get_value(const Eigen::VectorXd &iv) const
Performs the same task as df(), but returns the value of algebraic function, instead of saving it ins...
S fs_
Wrapper around algebraic system.
int df(const Eigen::VectorXd &iv, Eigen::MatrixXd &fjac) const
Assign the Jacobian to fjac.
A functor with the required operators to call Eigen's algebraic solver.
nlo_functor(int inputs, int values)
A structure which gets passed to Eigen's dogleg algebraic solver.