Automatic Differentiation
 
Loading...
Searching...
No Matches
neg_binomial_2_log_lpmf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_OPENCL_PRIM_NEG_BINOMIAL_2_LOG_LPMF_HPP
2#define STAN_MATH_OPENCL_PRIM_NEG_BINOMIAL_2_LOG_LPMF_HPP
3#ifdef STAN_OPENCL
4
14
15namespace stan {
16namespace math {
17
36template <bool propto, typename T_n_cl, typename T_log_location_cl,
37 typename T_precision_cl,
39 T_n_cl, T_log_location_cl, T_precision_cl>* = nullptr,
40 require_any_not_stan_scalar_t<T_n_cl, T_log_location_cl,
41 T_precision_cl>* = nullptr>
42inline return_type_t<T_n_cl, T_log_location_cl, T_precision_cl>
43neg_binomial_2_log_lpmf(const T_n_cl& n, const T_log_location_cl& eta,
44 const T_precision_cl& phi) {
45 static constexpr const char* function = "neg_binomial_2_log_lpmf(OpenCL)";
46 using T_partials_return
48 using std::isfinite;
49 using std::isnan;
50
51 check_consistent_sizes(function, "Failures variable", n,
52 "Log location parameter", eta, "Precision parameter",
53 phi);
54 const size_t N = max_size(n, eta, phi);
55 if (N == 0) {
56 return 0.0;
57 }
58 if (!include_summand<propto, T_n_cl, T_log_location_cl,
59 T_precision_cl>::value) {
60 return 0.0;
61 }
62
63 const auto& eta_col = as_column_vector_or_scalar(eta);
64 const auto& phi_col = as_column_vector_or_scalar(phi);
65
66 const auto& eta_val = value_of(eta_col);
67 const auto& phi_val = value_of(phi_col);
68
69 auto check_n_nonnegative
70 = check_cl(function, "Failures variable", n, "nonnegative");
71 auto n_nonnegative = n >= 0;
72 auto check_eta_finite
73 = check_cl(function, "Log location parameter", eta_val, "finite");
74 auto eta_finite = isfinite(eta_val);
75 auto check_phi_positive_finite
76 = check_cl(function, "Precision parameter", phi_val, "positive finite");
77 auto phi_positive_finite = 0 < phi_val && isfinite(phi_val);
78
79 auto log_phi = log(phi_val);
80 auto exp_eta = exp(eta_val);
81 auto exp_eta_over_exp_eta_phi
82 = elt_divide(1.0, elt_divide(phi_val, exp_eta) + 1.0);
83 auto log1p_exp_eta_m_logphi = log1p_exp(eta_val - log_phi);
84 auto n_plus_phi = n + phi_val;
85
86 auto logp1 = -elt_multiply(phi_val, log1p_exp_eta_m_logphi)
87 - elt_multiply(n, log_phi + log1p_exp_eta_m_logphi);
88 auto logp2 = static_select<include_summand<propto, T_precision_cl>::value>(
89 logp1 + binomial_coefficient_log(n_plus_phi - 1, n), logp1);
90 auto logp_expr = colwise_sum(
92 logp2 + elt_multiply(n, eta_val), logp2));
93
94 auto eta_deriv = n - elt_multiply(n_plus_phi, exp_eta_over_exp_eta_phi);
95 auto phi_deriv = exp_eta_over_exp_eta_phi - elt_divide(n, exp_eta + phi_val)
96 - log1p_exp_eta_m_logphi - digamma(phi_val)
97 + digamma(n_plus_phi);
98
99 matrix_cl<double> logp_cl;
100 matrix_cl<double> eta_deriv_cl;
101 matrix_cl<double> phi_deriv_cl;
102
103 results(check_n_nonnegative, check_eta_finite, check_phi_positive_finite,
104 logp_cl, eta_deriv_cl, phi_deriv_cl)
105 = expressions(n_nonnegative, eta_finite, phi_positive_finite, logp_expr,
108
109 T_partials_return logp = sum(from_matrix_cl(logp_cl));
110
111 auto ops_partials = make_partials_propagator(eta_col, phi_col);
112
114 partials<0>(ops_partials) = std::move(eta_deriv_cl);
115 }
117 partials<1>(ops_partials) = std::move(phi_deriv_cl);
118 }
119 return ops_partials.build(logp);
120}
121
122} // namespace math
123} // namespace stan
124#endif
125#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)
isfinite_< as_operation_cl_t< T > > isfinite(T &&a)
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.
binomial_coefficient_log_< as_operation_cl_t< T1 >, as_operation_cl_t< T2 > > binomial_coefficient_log(T1 &&a, T2 &&b)
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.
return_type_t< T_n_cl, T_log_location_cl, T_precision_cl > neg_binomial_2_log_lpmf(const T_n_cl &n, const T_log_location_cl &eta, const T_precision_cl &phi)
The log of the log transformed negative binomial density for the specified scalars given the specifie...
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
require_all_t< is_prim_or_rev_kernel_expression< std::decay_t< Types > >... > require_all_prim_or_rev_kernel_expression_t
Require type satisfies is_prim_or_rev_kernel_expression.
require_any_not_t< is_stan_scalar< std::decay_t< Types > >... > require_any_not_stan_scalar_t
Require at least one of the types do not satisfy is_stan_scalar.
size_t max_size(const T1 &x1, const Ts &... xs)
Calculate the size of the largest input.
Definition max_size.hpp:19
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:15
T1 static_select(T1 &&a, T2 &&b)
Returns one of the arguments that can be of different type, depending on the compile time condition.
fvar< T > log1p_exp(const fvar< T > &x)
Definition log1p_exp.hpp:13
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
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
fvar< T > digamma(const fvar< T > &x)
Return the derivative of the log gamma function at the specified argument.
Definition digamma.hpp:23
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
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...