Automatic Differentiation
 
Loading...
Searching...
No Matches
student_t_lpdf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_OPENCL_PRIM_STUDENT_T_LPDF_HPP
2#define STAN_MATH_OPENCL_PRIM_STUDENT_T_LPDF_HPP
3#ifdef STAN_OPENCL
4
12
13namespace stan {
14namespace math {
15
45template <bool propto, typename T_y_cl, typename T_dof_cl, typename T_loc_cl,
46 typename T_scale_cl,
48 T_y_cl, T_dof_cl, T_loc_cl, T_scale_cl>* = nullptr,
49 require_any_not_stan_scalar_t<T_y_cl, T_dof_cl, T_loc_cl,
50 T_scale_cl>* = nullptr>
52 const T_y_cl& y, const T_dof_cl& nu, const T_loc_cl& mu,
53 const T_scale_cl& sigma) {
54 static constexpr const char* function = "student_t_lpdf(OpenCL)";
55 using T_partials_return
57 using std::isfinite;
58 using std::isnan;
59
60 check_consistent_sizes(function, "Random variable", y,
61 "Degrees of freedom parameter", nu,
62 "Location parameter", mu, "Scale parameter", sigma);
63 const size_t N = max_size(y, mu, sigma);
64 if (N == 0) {
65 return 0.0;
66 }
68 return 0.0;
69 }
70
71 const auto& y_col = as_column_vector_or_scalar(y);
72 const auto& nu_col = as_column_vector_or_scalar(nu);
73 const auto& mu_col = as_column_vector_or_scalar(mu);
74 const auto& sigma_col = as_column_vector_or_scalar(sigma);
75
76 const auto& y_val = value_of(y_col);
77 const auto& nu_val = value_of(nu_col);
78 const auto& mu_val = value_of(mu_col);
79 const auto& sigma_val = value_of(sigma_col);
80
81 auto check_y_not_nan
82 = check_cl(function, "Random variable", y_val, "not NaN");
83 auto y_not_nan = !isnan(y_val);
84 auto check_nu_positive_finite = check_cl(
85 function, "Degrees of freedom parameter", nu_val, "positive finite");
86 auto nu_positive_finite = 0 < nu_val && isfinite(nu_val);
87 auto check_mu_finite
88 = check_cl(function, "Location parameter", mu_val, "finite");
89 auto mu_finite = isfinite(mu_val);
90 auto check_sigma_positive_finite
91 = check_cl(function, "Scale parameter", sigma_val, "positive finite");
92 auto sigma_positive_finite = 0 < sigma_val && isfinite(sigma_val);
93
94 auto half_nu = 0.5 * nu_val;
95 auto y_scaled = elt_divide(y_val - mu_val, sigma_val);
96 auto square_y_scaled = elt_multiply(y_scaled, y_scaled);
97 auto square_y_scaled_over_nu = elt_divide(square_y_scaled, nu_val);
98 auto log1p_val = log1p(square_y_scaled_over_nu);
99
100 auto logp1 = -elt_multiply((half_nu + 0.5), log1p_val);
101 auto logp2 = static_select<include_summand<propto, T_dof_cl>::value>(
102 logp1 + lgamma(half_nu + 0.5) - lgamma(half_nu) - 0.5 * log(nu_val),
103 logp1);
104 auto logp_expr
106 logp2 - log(sigma_val), logp2));
107
108 auto square_sigma = elt_multiply(sigma_val, sigma_val);
109 auto deriv_y_mu = elt_divide(
110 elt_multiply(nu_val + 1, y_val - mu_val),
111 elt_multiply(elt_multiply(1 + square_y_scaled_over_nu, square_sigma),
112 nu_val));
113 auto rep_deriv = elt_divide(elt_multiply(nu_val + 1, square_y_scaled_over_nu),
114 1 + square_y_scaled_over_nu)
115 - 1;
116 auto nu_deriv = 0.5
117 * (digamma(half_nu + 0.5) - digamma(half_nu) - log1p_val
118 + elt_divide(rep_deriv, nu_val));
119 auto sigma_deriv = elt_divide(rep_deriv, sigma_val);
120
121 matrix_cl<double> logp_cl;
122 matrix_cl<double> y_deriv_cl;
123 matrix_cl<double> nu_deriv_cl;
124 matrix_cl<double> mu_deriv_cl;
125 matrix_cl<double> sigma_deriv_cl;
126
127 results(check_y_not_nan, check_nu_positive_finite, check_mu_finite,
128 check_sigma_positive_finite, logp_cl, y_deriv_cl, nu_deriv_cl,
129 mu_deriv_cl, sigma_deriv_cl)
130 = expressions(y_not_nan, nu_positive_finite, mu_finite,
131 sigma_positive_finite, logp_expr,
132 calc_if<!is_constant<T_y_cl>::value>(-deriv_y_mu),
136
137 T_partials_return logp = sum(from_matrix_cl(logp_cl));
138
140 logp -= LOG_SQRT_PI * N;
141 }
142
143 auto ops_partials
144 = make_partials_propagator(y_col, nu_col, mu_col, sigma_col);
145
147 partials<0>(ops_partials) = std::move(y_deriv_cl);
148 }
150 partials<1>(ops_partials) = std::move(nu_deriv_cl);
151 }
153 partials<2>(ops_partials) = std::move(mu_deriv_cl);
154 }
156 partials<3>(ops_partials) = std::move(sigma_deriv_cl);
157 }
158 return ops_partials.build(logp);
159}
160
161} // namespace math
162} // namespace stan
163#endif
164#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.
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_y_cl, T_dof_cl, T_loc_cl, T_scale_cl > student_t_lpdf(const T_y_cl &y, const T_dof_cl &nu, const T_loc_cl &mu, const T_scale_cl &sigma)
The log of the Student-t density for the given y, nu, mean, and scale parameter.
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.
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
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.
void check_consistent_sizes(const char *)
Trivial no input case, this function is a no-op.
static constexpr double LOG_SQRT_PI
The natural logarithm of the square root of , .
fvar< T > log1p(const fvar< T > &x)
Definition log1p.hpp:12
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
Definition lgamma.hpp:21
auto sum(const std::vector< T > &m)
Return the sum of the entries of the specified standard vector.
Definition sum.hpp:23
int64_t max_size(const T1 &x1, const Ts &... xs)
Calculate the size of the largest input.
Definition max_size.hpp:20
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
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 ...
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...