Automatic Differentiation
 
Loading...
Searching...
No Matches
bernoulli_logit_lpmf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_OPENCL_PRIM_BERNOULLI_LOGIT_LPMF_HPP
2#define STAN_MATH_OPENCL_PRIM_BERNOULLI_LOGIT_LPMF_HPP
3#ifdef STAN_OPENCL
4
9
10namespace stan {
11namespace math {
12
25template <
26 bool propto, typename T_n_cl, typename T_prob_cl,
27 require_all_prim_or_rev_kernel_expression_t<T_n_cl, T_prob_cl>* = nullptr,
28 require_any_not_stan_scalar_t<T_n_cl, T_prob_cl>* = nullptr>
30 const T_prob_cl& theta) {
31 static constexpr const char* function = "bernoulli_logit_lpmf(OpenCL)";
32 using T_partials_return = partials_return_t<T_prob_cl>;
33 using std::isnan;
34 constexpr bool is_n_vector = !is_stan_scalar<T_n_cl>::value;
35
36 check_consistent_sizes(function, "Random variable", n,
37 "Probability parameter", theta);
38 const size_t N = is_n_vector ? math::size(n) : math::size(theta);
39 if (N == 0) {
40 return 0.0;
41 }
43 return 0.0;
44 }
45
46 const auto& theta_col = as_column_vector_or_scalar(theta);
47 const auto& theta_val = value_of(theta_col);
48
49 auto check_n_bounded = check_cl(function, "n", n, "in the interval [0, 1]");
50 auto n_bounded_expr = 0 <= n && n <= 1;
51 auto check_theta_not_nan
52 = check_cl(function, "Logit transformed probability parameter", theta_val,
53 "not NaN");
54 auto theta_not_nan_expr = !isnan(theta_val);
55
56 auto signs_expr = 2 * n - 1.0; // subtracting 1.0 converts int to double
57 auto ntheta_expr = elt_multiply(signs_expr, theta_val);
58 auto exp_m_ntheta_expr = exp(-ntheta_expr);
59 static constexpr double cutoff = 20.0;
60 auto condition1_expr = ntheta_expr > cutoff;
61 auto condition2_expr = ntheta_expr < -cutoff;
62 auto logp_expr = colwise_sum(
63 select(condition1_expr, -exp_m_ntheta_expr,
64 select(condition2_expr, ntheta_expr, -log1p(exp_m_ntheta_expr))));
65 auto deriv_expr = select(
66 condition1_expr, -exp_m_ntheta_expr,
67 select(condition2_expr, signs_expr,
68 elt_multiply(signs_expr, elt_divide(exp_m_ntheta_expr,
69 (exp_m_ntheta_expr + 1)))));
70
71 matrix_cl<double> logp_cl;
72 matrix_cl<double> deriv_cl;
73
74 results(logp_cl, deriv_cl, check_n_bounded, check_theta_not_nan)
75 = expressions(logp_expr,
77 n_bounded_expr, theta_not_nan_expr);
78
79 T_partials_return logp = sum(from_matrix_cl(logp_cl));
80 auto ops_partials = make_partials_propagator(theta_col);
81
83 partials<0>(ops_partials) = deriv_cl;
84 }
85
86 return ops_partials.build(logp);
87}
88
89} // namespace math
90} // namespace stan
91#endif
92#endif
Represents an arithmetic matrix on the OpenCL device.
Definition matrix_cl.hpp:47
elt_multiply_< as_operation_cl_t< T_a >, as_operation_cl_t< T_b > > elt_multiply(T_a &&a, T_b &&b)
select_< as_operation_cl_t< T_condition >, as_operation_cl_t< T_then >, as_operation_cl_t< T_else > > select(T_condition &&condition, T_then &&then, T_else &&els)
Selection operation on kernel generator expressions.
Definition select.hpp:148
auto check_cl(const char *function, const char *var_name, T &&y, const char *must_be)
Constructs a check on opencl matrix or expression.
Definition check_cl.hpp:219
results_cl< T_results... > results(T_results &&... results)
Deduces types for constructing results_cl object.
auto as_column_vector_or_scalar(T &&a)
as_column_vector_or_scalar of a kernel generator expression.
elt_divide_< as_operation_cl_t< T_a >, as_operation_cl_t< T_b > > elt_divide(T_a &&a, T_b &&b)
calc_if_< true, as_operation_cl_t< T > > calc_if(T &&a)
Definition calc_if.hpp:121
auto colwise_sum(T &&a)
Column wise sum - reduction of a kernel generator expression.
expressions_cl< T_expressions... > expressions(T_expressions &&... expressions)
Deduces types for constructing expressions_cl object.
auto from_matrix_cl(const T &src)
Copies the source matrix that is stored on the OpenCL device to the destination Eigen matrix.
Definition copy.hpp:61
return_type_t< T_prob_cl > bernoulli_logit_lpmf(const T_n_cl &n, const T_prob_cl &theta)
Returns the log PMF of the logit-parametrized Bernoulli distribution.
size_t size(const T &m)
Returns the size (number of the elements) of a matrix_cl or var_value<matrix_cl<T>>.
Definition size.hpp:18
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
T value_of(const fvar< T > &v)
Return the value of the specified variable.
Definition value_of.hpp:18
void check_consistent_sizes(const char *)
Trivial no input case, this function is a no-op.
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Definition sum.hpp:22
fvar< T > log1p(const fvar< T > &x)
Definition log1p.hpp:12
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
fvar< T > exp(const fvar< T > &x)
Definition exp.hpp:13
typename partials_return_type< Args... >::type partials_return_t
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
Definition fvar.hpp:9
bool isnan(const stan::math::var &a)
Checks if the given number is NaN.
Definition std_isnan.hpp:18
Checks if decayed type is a var, fvar, or arithmetic.
Extends std::true_type when instantiated with zero or more template parameters, all of which extend t...
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...