1#ifndef STAN_MATH_OPENCL_PRIM_SKEW_DOUBLE_EXPONENTIAL_LPDF_HPP 
    2#define STAN_MATH_OPENCL_PRIM_SKEW_DOUBLE_EXPONENTIAL_LPDF_HPP 
   33template <
bool propto, 
typename T_y_cl, 
typename T_loc_cl, 
typename T_scale_cl,
 
   34          typename T_skewness_cl,
 
   36              T_y_cl, T_loc_cl, T_scale_cl, T_skewness_cl>* = 
nullptr,
 
   38                                        T_skewness_cl>* = 
nullptr>
 
   39inline return_type_t<T_y_cl, T_loc_cl, T_scale_cl, T_skewness_cl>
 
   41                             const T_scale_cl& sigma,
 
   42                             const T_skewness_cl& tau) {
 
   43  static constexpr const char* function
 
   44      = 
"skew_double_exponential_lpdf(OpenCL)";
 
   45  using T_partials_return
 
   51                         mu, 
"Scale parameter", sigma, 
"Inv_scale paramter",
 
   53  const size_t N = 
max_size(y, mu, sigma, tau);
 
   58                       T_skewness_cl>::value) {
 
   68  const auto& mu_val = 
value_of(mu_col);
 
   69  const auto& sigma_val = 
value_of(sigma_col);
 
   70  const auto& tau_val = 
value_of(tau_col);
 
   73      = 
check_cl(function, 
"Random variable", y_val, 
"not_nan");
 
   74  auto y_not_nan_expr = !isnan(y_val);
 
   76      = 
check_cl(function, 
"Location parameter", mu_val, 
"finite");
 
   77  auto mu_finite_expr = 
isfinite(mu_val);
 
   78  auto check_sigma_positive_finite
 
   79      = 
check_cl(function, 
"Scale parameter", sigma_val, 
"positive finite");
 
   80  auto sigma_positive_finite_expr = 
isfinite(sigma_val) && sigma_val > 0;
 
   81  auto check_tau_bounded = 
check_cl(function, 
"Skewness parameter", tau_val,
 
   82                                    "in the interval [0, 1]");
 
   83  auto tau_bounded_expr = 0.0 < tau_val && tau_val < 1.0;
 
   86  auto y_m_mu = y_val - mu_val;
 
   87  auto diff_sign = 
sign(y_m_mu);
 
   88  auto diff_sign_smaller_0 = diff_sign < 0.0;
 
   89  auto abs_diff_y_mu = 
fabs(y_m_mu);
 
   90  auto abs_diff_y_mu_over_sigma = 
elt_multiply(abs_diff_y_mu, inv_sigma);
 
   91  auto tmp = diff_sign_smaller_0 + 
elt_multiply(diff_sign, tau_val);
 
   94  auto logp1 = -2.0 * expo;
 
   95  auto logp2 = static_select<include_summand<propto, T_scale_cl>::value>(
 
   96      logp1 - 
log(sigma_val), logp1);
 
   97  auto logp3 = static_select<include_summand<propto, T_skewness_cl>::value>(
 
   98      logp2 + 
log(tau_val) + 
log1m(tau_val), logp2);
 
  102  auto y_deriv = -mu_deriv;
 
  103  auto sigma_deriv = -inv_sigma + 2.0 * 
elt_multiply(expo, inv_sigma);
 
  105                   - 
elt_multiply(diff_sign * 2.0, abs_diff_y_mu_over_sigma);
 
  113  results(check_y_not_nan, check_mu_finite, check_sigma_positive_finite,
 
  114          check_tau_bounded, logp_cl, y_deriv_cl, mu_deriv_cl, sigma_deriv_cl,
 
  116      = 
expressions(y_not_nan_expr, mu_finite_expr, sigma_positive_finite_expr,
 
  117                    tau_bounded_expr, logp_expr,
 
  118                    calc_if<is_autodiff_v<T_y_cl>>(y_deriv),
 
  119                    calc_if<is_autodiff_v<T_loc_cl>>(mu_deriv),
 
  120                    calc_if<is_autodiff_v<T_scale_cl>>(sigma_deriv),
 
  121                    calc_if<is_autodiff_v<T_skewness_cl>>(tau_deriv));
 
  131  if constexpr (is_autodiff_v<T_y_cl>) {
 
  132    partials<0>(ops_partials) = std::move(y_deriv_cl);
 
  134  if constexpr (is_autodiff_v<T_loc_cl>) {
 
  135    partials<1>(ops_partials) = std::move(mu_deriv_cl);
 
  137  if constexpr (is_autodiff_v<T_scale_cl>) {
 
  138    partials<2>(ops_partials) = std::move(sigma_deriv_cl);
 
  140  if constexpr (is_autodiff_v<T_skewness_cl>) {
 
  141    partials<3>(ops_partials) = std::move(tau_deriv_cl);
 
  143  return ops_partials.build(logp);
 
Represents an arithmetic matrix on the OpenCL device.
 
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.
 
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)
 
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_loc_cl, T_scale_cl, T_skewness_cl > skew_double_exponential_lpdf(const T_y_cl &y, const T_loc_cl &mu, const T_scale_cl &sigma, const T_skewness_cl &tau)
Returns the log PMF of the skew double exponential distribution.
 
auto from_matrix_cl(const T &src)
Copies the source matrix that is stored on the OpenCL device to the destination Eigen matrix.
 
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.
 
auto sign(const T &x)
Returns signs of the arguments.
 
T value_of(const fvar< T > &v)
Return the value of the specified variable.
 
fvar< T > log(const fvar< T > &x)
 
static constexpr double LOG_TWO
The natural logarithm of 2, .
 
void check_consistent_sizes(const char *)
Trivial no input case, this function is a no-op.
 
auto sum(const std::vector< T > &m)
Return the sum of the entries of the specified standard vector.
 
int64_t max_size(const T1 &x1, const Ts &... xs)
Calculate the size of the largest input.
 
fvar< T > log1m(const fvar< T > &x)
 
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
 
fvar< T > fabs(const fvar< T > &x)
 
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.
 
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...