Loading [MathJax]/jax/input/TeX/config.js
Automatic Differentiation
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
laplace_marginal_bernoulli_logit_lpmf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_MIX_PROB_LAPLACE_MARGINAL_BERNOULLI_LOGIT_LPMF_HPP
2#define STAN_MATH_MIX_PROB_LAPLACE_MARGINAL_BERNOULLI_LOGIT_LPMF_HPP
3
18
19namespace stan {
20namespace math {
21
23 template <typename ThetaVec, typename YVec, typename Mean>
24 inline auto operator()(const ThetaVec& theta, const YVec& y,
25 const std::vector<int>& delta_int, Mean&& mean,
26 std::ostream* pstream) const {
27 auto theta_offset = to_ref(add(theta, mean));
28 return sum(
29 elt_multiply(theta_offset, y)
30 - elt_multiply(to_vector(delta_int), log(add(1.0, exp(theta_offset)))));
31 }
32};
33
54template <bool propto = false, typename ThetaVec, typename Mean,
55 typename CovarFun, typename CovarArgs,
58 const std::vector<int>& y, const std::vector<int>& n_samples, Mean&& mean,
59 CovarFun&& covariance_function, CovarArgs&& covar_args,
60 const ThetaVec& theta_0, double tolerance, int max_num_steps,
61 const int hessian_block_size, const int solver,
62 const int max_steps_line_search, std::ostream* msgs) {
63 laplace_options_user_supplied ops{hessian_block_size, solver,
64 max_steps_line_search, tolerance,
65 max_num_steps, value_of(theta_0)};
68 std::forward_as_tuple(to_vector(y), n_samples, std::forward<Mean>(mean)),
69 std::forward<CovarFun>(covariance_function),
70 std::forward<CovarArgs>(covar_args), ops, msgs);
71}
72
90template <bool propto = false, typename Mean, typename CovarFun,
91 typename CovarArgs>
93 const std::vector<int>& y, const std::vector<int>& n_samples, Mean&& mean,
94 CovarFun&& covariance_function, CovarArgs&& covar_args,
95 std::ostream* msgs) {
98 std::forward_as_tuple(to_vector(y), n_samples, std::forward<Mean>(mean)),
99 std::forward<CovarFun>(covariance_function),
100 std::forward<CovarArgs>(covar_args), laplace_options_default{}, msgs);
101}
102
103} // namespace math
104} // namespace stan
105
106#endif
require_t< is_eigen_vector< std::decay_t< T > > > require_eigen_vector_t
Require type satisfies is_eigen_vector.
elt_multiply_< as_operation_cl_t< T_a >, as_operation_cl_t< T_b > > elt_multiply(T_a &&a, T_b &&b)
addition_< as_operation_cl_t< T_a >, as_operation_cl_t< T_b > > add(T_a &&a, T_b &&b)
auto to_vector(T_x &&x)
Returns input matrix reshaped into a vector.
Definition to_vector.hpp:21
Reference for calculations of marginal and its gradients: Margossian et al (2020),...
scalar_type_t< T > mean(const T &m)
Returns the sample mean (i.e., average) of the coefficients in the specified std vector,...
Definition mean.hpp:20
T value_of(const fvar< T > &v)
Return the value of the specified variable.
Definition value_of.hpp:18
fvar< T > log(const fvar< T > &x)
Definition log.hpp:18
auto laplace_marginal_tol_bernoulli_logit_lpmf(const std::vector< int > &y, const std::vector< int > &n_samples, Mean &&mean, CovarFun &&covariance_function, CovarArgs &&covar_args, const ThetaVec &theta_0, double tolerance, int max_num_steps, const int hessian_block_size, const int solver, const int max_steps_line_search, std::ostream *msgs)
Wrapper function around the laplace_marginal function for a logistic Bernoulli likelihood.
auto sum(const std::vector< T > &m)
Return the sum of the entries of the specified standard vector.
Definition sum.hpp:23
auto laplace_marginal_bernoulli_logit_lpmf(const std::vector< int > &y, const std::vector< int > &n_samples, Mean &&mean, CovarFun &&covariance_function, CovarArgs &&covar_args, std::ostream *msgs)
Wrapper function around the laplace_marginal function for a logistic Bernoulli likelihood.
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:18
fvar< T > exp(const fvar< T > &x)
Definition exp.hpp:15
double laplace_marginal_density(LLFun &&ll_fun, LLTupleArgs &&ll_args, CovarFun &&covariance_function, CovarArgs &&covar_args, const laplace_options< InitTheta > &options, std::ostream *msgs)
For a latent Gaussian model with global parameters phi, latent variables theta, and observations y,...
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
auto operator()(const ThetaVec &theta, const YVec &y, const std::vector< int > &delta_int, Mean &&mean, std::ostream *pstream) const