Stan  2.10.0
probability, sampling & optimization
Classes | Namespaces | Functions
util.hpp File Reference
#include <stan/interface_callbacks/writer/base_writer.hpp>
#include <stan/math/fwd/scal/fun/square.hpp>
#include <stan/math/fwd/core.hpp>
#include <stan/math/prim/mat/fun/Eigen.hpp>
#include <stan/math/rev/mat/fun/grad.hpp>
#include <stan/math/rev/core.hpp>
#include <stan/math/mix/mat/functor/derivative.hpp>
#include <stan/math/mix/mat/functor/grad_hessian.hpp>
#include <stan/math/mix/mat/functor/grad_tr_mat_times_hessian.hpp>
#include <stan/math/rev/mat/functor/gradient.hpp>
#include <stan/math/fwd/mat/functor/gradient.hpp>
#include <stan/math/mix/mat/functor/gradient_dot_vector.hpp>
#include <stan/math/mix/mat/functor/hessian.hpp>
#include <stan/math/mix/mat/functor/hessian_times_vector.hpp>
#include <stan/math/fwd/mat/functor/jacobian.hpp>
#include <stan/math/rev/mat/functor/jacobian.hpp>
#include <stan/math/mix/mat/functor/partial_derivative.hpp>
#include <cmath>
#include <iomanip>
#include <iostream>
#include <string>
#include <vector>

Go to the source code of this file.

Classes

struct  stan::model::model_functional< M >
 

Namespaces

 stan
 Probability, optimization and sampling library.
 
 stan::model
 For compiling models.
 

Functions

template<bool jacobian_adjust_transform, class M >
double stan::model::log_prob_propto (const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::ostream *msgs=0)
 Helper function to calculate log probability for double scalars up to a proportion. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
double stan::model::log_prob_grad (const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::ostream *msgs=0)
 Compute the gradient using reverse-mode automatic differentiation, writing the result into the specified gradient, using the specified perturbation. More...
 
template<bool jacobian_adjust_transform, class M >
double stan::model::log_prob_propto (const M &model, Eigen::VectorXd &params_r, std::ostream *msgs=0)
 Helper function to calculate log probability for double scalars up to a proportion. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
double stan::model::log_prob_grad (const M &model, Eigen::VectorXd &params_r, Eigen::VectorXd &gradient, std::ostream *msgs=0)
 Compute the gradient using reverse-mode automatic differentiation, writing the result into the specified gradient, using the specified perturbation. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
void stan::model::finite_diff_grad (const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &grad, double epsilon=1e-6, std::ostream *msgs=0)
 Compute the gradient using finite differences for the specified parameters, writing the result into the specified gradient, using the specified perturbation. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
int stan::model::test_gradients (const M &model, std::vector< double > &params_r, std::vector< int > &params_i, double epsilon, double error, stan::interface_callbacks::writer::base_writer &writer)
 Test the log_prob_grad() function's ability to produce accurate gradients using finite differences. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
double stan::model::grad_hess_log_prob (const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::vector< double > &hessian, std::ostream *msgs=0)
 Evaluate the log-probability, its gradient, and its Hessian at params_r. More...
 
template<class M >
void stan::model::gradient (const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, std::ostream *msgs=0)
 
template<class M >
void stan::model::gradient (const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, stan::interface_callbacks::writer::base_writer &writer)
 
template<class M >
void stan::model::hessian (const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &hess_f, std::ostream *msgs=0)
 
template<class M >
void stan::model::gradient_dot_vector (const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &f, double &grad_f_dot_v, std::ostream *msgs=0)
 
template<class M >
void stan::model::hessian_times_vector (const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &hess_f_dot_v, std::ostream *msgs=0)
 
template<class M >
void stan::model::grad_tr_mat_times_hessian (const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &X, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_tr_X_hess_f, std::ostream *msgs=0)
 

     [ Stan Home Page ] © 2011–2016, Stan Development Team.