Automatic Differentiation
 
Loading...
Searching...
No Matches
scaled_inv_chi_square_lpdf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_OPENCL_PRIM_SCALED_INV_CHI_SQUARE_LPDF_HPP
2#define STAN_MATH_OPENCL_PRIM_SCALED_INV_CHI_SQUARE_LPDF_HPP
3#ifdef STAN_OPENCL
4
12
13namespace stan {
14namespace math {
15
38template <
39 bool propto, typename T_y_cl, typename T_dof_cl, typename T_scale_cl,
41 T_scale_cl>* = nullptr,
42 require_any_not_stan_scalar_t<T_y_cl, T_dof_cl, T_scale_cl>* = nullptr>
44 const T_y_cl& y, const T_dof_cl& nu, const T_scale_cl& s) {
45 static constexpr const char* function = "scaled_inv_chi_square_lpdf(OpenCL)";
47 using std::isfinite;
48 using std::isnan;
49
50 check_consistent_sizes(function, "Random variable", y,
51 "Degrees of freedom parameter", nu, "Scale parameter",
52 s);
53 const size_t N = max_size(y, nu, s);
54 if (N == 0) {
55 return 0.0;
56 }
58 return 0.0;
59 }
60
61 const auto& y_col = as_column_vector_or_scalar(y);
62 const auto& nu_col = as_column_vector_or_scalar(nu);
63 const auto& s_col = as_column_vector_or_scalar(s);
64
65 const auto& y_val = value_of(y_col);
66 const auto& nu_val = value_of(nu_col);
67 const auto& s_val = value_of(s_col);
68
69 auto check_y_not_nan
70 = check_cl(function, "Random variable", y_val, "not NaN");
71 auto y_not_nan = !isnan(y_val);
72 auto check_nu_positive_finite = check_cl(
73 function, "Degrees of freedom parameter", nu_val, "positive finite");
74 auto nu_positive_finite = isfinite(nu_val) && 0 < nu_val;
75 auto check_s_positive_finite
76 = check_cl(function, "Scale parameter", s_val, "positive finite");
77 auto s_positive_finite = isfinite(s_val) && 0 < s_val;
78
79 auto any_y_nonpositive = colwise_max(cast<char>(y_val <= 0.0));
80 auto half_nu = 0.5 * nu_val;
81 auto log_y = log(y_val);
82 auto inv_y = elt_divide(1.0, y_val);
83 auto log_s = log(s_val);
84 auto log_half_nu = log(half_nu);
85 auto s_square = elt_multiply(s_val, s_val);
86
87 auto logp1 = static_select<include_summand<propto, T_dof_cl>::value>(
88 elt_multiply(half_nu, log_half_nu) - lgamma(half_nu), constant(0, N, 1));
89 auto logp2
90 = static_select<include_summand<propto, T_dof_cl, T_scale_cl>::value>(
91 logp1 + elt_multiply(nu_val, log_s), logp1);
92 auto logp3 = static_select<include_summand<propto, T_dof_cl, T_y_cl>::value>(
93 logp2 - elt_multiply(half_nu + 1.0, log_y), logp2);
94 auto logp_expr = colwise_sum(
97 logp3 - elt_multiply(elt_multiply(half_nu, s_square), inv_y), logp3));
98
99 auto y_deriv = -elt_multiply(half_nu + 1.0, inv_y)
100 + elt_multiply(elt_multiply(half_nu, s_square),
101 elt_multiply(inv_y, inv_y));
102 auto nu_deriv = 0.5
103 * (log_half_nu - digamma(half_nu) - log_y
104 - elt_multiply(s_square, inv_y))
105 + log_s + 0.5;
106 auto s_deriv = elt_divide(nu_val, s_val)
107 - elt_multiply(elt_multiply(nu_val, inv_y), s_val);
108
109 matrix_cl<char> any_y_nonpositive_cl;
110 matrix_cl<double> logp_cl;
111 matrix_cl<double> nu_deriv_cl;
112 matrix_cl<double> y_deriv_cl;
113 matrix_cl<double> s_deriv_cl;
114
115 results(check_y_not_nan, check_nu_positive_finite, check_s_positive_finite,
116 any_y_nonpositive_cl, logp_cl, y_deriv_cl, nu_deriv_cl, s_deriv_cl)
117 = expressions(y_not_nan, nu_positive_finite, s_positive_finite,
118 any_y_nonpositive, logp_expr,
122
123 if (from_matrix_cl(any_y_nonpositive_cl).any()) {
124 return LOG_ZERO;
125 }
126
127 T_partials_return logp = sum(from_matrix_cl(logp_cl));
128
129 auto ops_partials = make_partials_propagator(y_col, nu_col, s_col);
130
132 partials<0>(ops_partials) = std::move(y_deriv_cl);
133 }
135 partials<1>(ops_partials) = std::move(nu_deriv_cl);
136 }
138 partials<2>(ops_partials) = std::move(s_deriv_cl);
139 }
140 return ops_partials.build(logp);
141}
142
143} // namespace math
144} // namespace stan
145#endif
146#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)
auto constant(const T a, int rows, int cols)
Matrix of repeated values in kernel generator expressions.
Definition constant.hpp:130
auto colwise_max(T &&a)
Column wise max - reduction of a kernel generator expression.
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_scale_cl > scaled_inv_chi_square_lpdf(const T_y_cl &y, const T_dof_cl &nu, const T_scale_cl &s)
The log of a scaled inverse chi-squared density for y with the specified degrees of freedom 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.
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
static constexpr double LOG_ZERO
The natural logarithm of 0, .
Definition constants.hpp:68
size_t max_size(const T1 &x1, const Ts &... xs)
Calculate the size of the largest input.
Definition max_size.hpp:19
constexpr bool any(T x)
Return true if any values in the input are true.
Definition any.hpp:21
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.
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 > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
Definition lgamma.hpp:21
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 ...
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...