1 #ifndef STAN__AGRAD__AUTO_DIFF_HPP
2 #define STAN__AGRAD__AUTO_DIFF_HPP
23 template <
typename T,
typename F>
47 template <
typename T,
typename F>
50 const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
54 Eigen::Matrix<fvar<T>,Eigen::Dynamic,1> x_fvar(x.size());
55 for (
int i = 0; i < x.size(); ++i)
94 const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
96 Eigen::Matrix<double,Eigen::Dynamic,1>& grad_fx) {
100 Eigen::Matrix<var,Eigen::Dynamic,1> x_var(x.size());
101 for (
int i = 0; i < x.size(); ++i)
103 var fx_var = f(x_var);
105 grad_fx.resize(x.size());
107 for (
int i = 0; i < x.size(); ++i)
108 grad_fx(i) = x_var(i).adj();
109 }
catch (
const std::exception& ) {
115 template <
typename T,
typename F>
118 const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
120 Eigen::Matrix<T,Eigen::Dynamic,1>& grad_fx) {
121 Eigen::Matrix<fvar<T>,Eigen::Dynamic,1> x_fvar(x.size());
122 grad_fx.resize(x.size());
123 for (
int i = 0; i < x.size(); ++i) {
124 for (
int k = 0; k < x.size(); ++k)
125 x_fvar(k) =
fvar<T>(x(k),k==i);
127 if (i == 0) fx = fx_fvar.
val_;
128 grad_fx(i) = fx_fvar.
d_;
133 template <
typename F>
136 const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
137 Eigen::Matrix<double,Eigen::Dynamic,1>& fx,
138 Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& J) {
139 using Eigen::Matrix;
using Eigen::Dynamic;
143 Matrix<var,Dynamic,1> x_var(x.size());
144 for (
int k = 0; k < x.size(); ++k)
146 Matrix<var,Dynamic,1> fx_var = f(x_var);
147 fx.resize(fx_var.size());
148 for (
int i = 0; i < fx_var.size(); ++i)
149 fx(i) = fx_var(i).val();
150 J.resize(x.size(), fx_var.size());
151 for (
int i = 0; i < fx_var.size(); ++i) {
155 for (
int k = 0; k < x.size(); ++k)
156 J(k,i) = x_var(k).adj();
158 }
catch (
const std::exception&
e) {
164 template <
typename T,
typename F>
167 const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
168 Eigen::Matrix<T,Eigen::Dynamic,1>& fx,
169 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& J) {
170 using Eigen::Matrix;
using Eigen::Dynamic;
172 Matrix<fvar<T>,Dynamic,1> x_fvar(x.size());
173 for (
int i = 0; i < x.size(); ++i) {
174 for (
int k = 0; k < x.size(); ++k)
175 x_fvar(k) =
fvar<T>(x(k), i == k);
176 Matrix<fvar<T>,Dynamic,1> fx_fvar
179 J.resize(x.size(),fx_fvar.size());
180 fx.resize(fx_fvar.size());
181 for (
int k = 0; k < fx_fvar.size(); ++k)
182 fx(k) = fx_fvar(k).
val_;
184 for (
int k = 0; k < fx_fvar.size(); ++k) {
185 J(i,k) = fx_fvar(k).
d_;
192 template <
typename F>
195 const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
197 Eigen::Matrix<double,Eigen::Dynamic,1>&
grad,
198 Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& H) {
201 H.resize(x.size(), x.size());
202 grad.resize(x.size());
203 for (
int i = 0; i < x.size(); ++i) {
204 Eigen::Matrix<fvar<var>, Eigen::Dynamic, 1> x_fvar(x.size());
205 for (
int j = 0; j < x.size(); ++j)
209 if (i == 0) fx = fx_fvar.
val_.
val();
211 for (
int j = 0; j < x.size(); ++j)
212 H(i,j) = x_fvar(j).val_.adj();
214 }
catch (
const std::exception&
e) {
221 template <
typename T,
typename F>
224 const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
226 Eigen::Matrix<T,Eigen::Dynamic,1>&
grad,
227 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& H) {
228 H.resize(x.size(), x.size());
229 grad.resize(x.size());
230 Eigen::Matrix<fvar<fvar<T> >,Eigen::Dynamic,1> x_fvar(x.size());
231 for (
int i = 0; i < x.size(); ++i) {
232 for (
int j = i; j < x.size(); ++j) {
233 for (
int k = 0; k < x.size(); ++k)
238 fx = fx_fvar.
val_.val_;
240 grad(i) = fx_fvar.
d_.val_;
241 H(i,j) = fx_fvar.
d_.d_;
250 template <
typename T1,
typename T2,
typename F>
253 const Eigen::Matrix<T1,Eigen::Dynamic,1>& x,
254 const Eigen::Matrix<T2,Eigen::Dynamic,1>& v,
260 using Eigen::Dynamic;
261 Matrix<fvar<T1>,Dynamic,1> x_fvar(x.size());
262 for (
int i = 0; i < x.size(); ++i)
266 grad_fx_dot_v = fx_fvar.
d_;
272 template <
typename F>
275 const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
276 const Eigen::Matrix<double,Eigen::Dynamic,1>& v,
278 Eigen::Matrix<double,Eigen::Dynamic,1>& Hv) {
282 using Eigen::Dynamic;
285 Matrix<var,Dynamic,1> x_var(x.size());
286 for (
int i = 0; i < x_var.size(); ++i)
289 var grad_fx_var_dot_v;
294 for (
int i = 0; i < x.size(); ++i)
295 Hv(i) = x_var(i).adj();
296 }
catch (
const std::exception&
e) {
302 template <
typename T,
typename F>
305 const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
306 const Eigen::Matrix<T,Eigen::Dynamic,1>& v,
308 Eigen::Matrix<T,Eigen::Dynamic,1>& Hv) {
310 using Eigen::Dynamic;
311 Matrix<T,Dynamic,1>
grad;
312 Matrix<T,Dynamic,Dynamic> H;
319 template <
typename F>
323 const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
324 const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& M,
325 Eigen::Matrix<double,Eigen::Dynamic,1>& grad_tr_MH) {
327 using Eigen::Dynamic;
331 grad_tr_MH.resize(x.size());
333 Matrix<var,Dynamic,1> x_var(x.size());
334 for (
int i = 0; i < x.size(); ++i)
337 Matrix<fvar<var>,Dynamic,1> x_fvar(x.size());
340 Matrix<double,Dynamic,1> M_n(x.size());
341 for (
int n = 0; n < x.size(); ++n) {
342 for (
int k = 0; k < x.size(); ++k)
344 for (
int k = 0; k < x.size(); ++k)
348 gradient_dot_vector<fvar<var>,
double>(f,x_fvar,M_n,fx,grad_fx_dot_v);
349 sum += grad_fx_dot_v.
d_;
353 for (
int i = 0; i < x.size(); ++i)
354 grad_tr_MH(i) = x_var(i).adj();
355 }
catch (
const std::exception&
e) {
void gradient(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_fx)
Calculate the value and the gradient of the specified function at the specified argument.
static void set_zero_all_adjoints()
Reset all adjoint values in the stack to zero.
static void start_nested()
Record the current position so that recover_memory_nested() can find it.
double val() const
Return the value of this variable.
static void grad(chainable *vi)
Compute the gradient for all variables starting from the specified root variable implementation.
void partial_derivative(const F &f, const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, int n, T &fx, T &dfx_dxn)
Return the partial derivative of the specified multiivariate function at the specified argument...
void hessian_times_vector(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &Hv)
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
vari * vi_
Pointer to the implementation of this variable.
void derivative(const F &f, const T &x, T &fx, T &dfx_dx)
Return the derivative of the specified univariate function at the specified argument.
void grad_tr_mat_times_hessian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &M, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_tr_MH)
Independent (input) and dependent (output) variables for gradients.
double e()
Return the base of the natural logarithm.
void jacobian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, Eigen::Matrix< double, Eigen::Dynamic, 1 > &fx, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &J)
void gradient_dot_vector(const F &f, const Eigen::Matrix< T1, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< T2, Eigen::Dynamic, 1 > &v, T1 &fx, T1 &grad_fx_dot_v)
static void recover_memory_nested()
Recover only the memory used for the top nested call.
void hessian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &H)