1#ifndef STAN_MATH_PRIM_FUN_BINOMIAL_COEFFICIENT_LOG_HPP 
    2#define STAN_MATH_PRIM_FUN_BINOMIAL_COEFFICIENT_LOG_HPP 
   79template <
typename T_n, 
typename T_k,
 
   80          require_all_stan_scalar_t<T_n, T_k>* = 
nullptr>
 
   94  const T_partials_return n_dbl = 
value_of(n);
 
   95  const T_partials_return k_dbl = 
value_of(k);
 
   96  const T_partials_return n_plus_1 = n_dbl + 1;
 
   97  const T_partials_return n_plus_1_mk = n_plus_1 - k_dbl;
 
   99  static constexpr const char* function = 
"binomial_coefficient_log";
 
  107  T_partials_return value;
 
  113    value = -
lbeta(n_plus_1_mk, k_dbl + 1) - 
log1p(n_dbl);
 
  116  if constexpr (is_any_autodiff_v<T_n, T_k>) {
 
  124    T_partials_return digamma_n_plus_1_mk = 
digamma(n_plus_1_mk);
 
  126    if constexpr (is_autodiff_v<T_n>) {
 
  129          partials<0>(ops_partials)[0] = 0;
 
  134        partials<0>(ops_partials)[0]
 
  135            = (
digamma(n_plus_1) - digamma_n_plus_1_mk);
 
  138    if constexpr (is_autodiff_v<T_k>) {
 
  139      if (k_dbl == 0 && n_dbl == -1.0) {
 
  141      } 
else if (k_dbl == -1) {
 
  142        partials<1>(ops_partials)[0] = 
INFTY;
 
  144        partials<1>(ops_partials)[0]
 
  145            = (digamma_n_plus_1_mk - 
digamma(k_dbl + 1));
 
  150  return ops_partials.build(value);
 
  163template <
typename T1, 
typename T2, require_any_container_t<T1, T2>* = 
nullptr>
 
  166      [](
auto&& c, 
auto&& d) {
 
  168                                        std::forward<
decltype(d)>(d));
 
  170      std::forward<T1>(a), std::forward<T2>(b));
 
binomial_coefficient_log_< as_operation_cl_t< T1 >, as_operation_cl_t< T2 > > binomial_coefficient_log(T1 &&a, T2 &&b)
 
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
 
double value_of_rec(const fvar< T > &v)
Return the value of the specified variable.
 
static constexpr double NOT_A_NUMBER
(Quiet) not-a-number value.
 
static constexpr double e()
Return the base of the natural logarithm.
 
fvar< T > lbeta(const fvar< T > &x1, const fvar< T > &x2)
 
auto apply_scalar_binary(F &&f, T1 &&x, T2 &&y)
Base template function for vectorization of binary scalar functions defined by applying a functor to ...
 
T value_of(const fvar< T > &v)
Return the value of the specified variable.
 
static constexpr double NEGATIVE_INFTY
Negative infinity.
 
void check_greater_or_equal(const char *function, const char *name, const T_y &y, const T_low &low, Idxs... idxs)
Throw an exception if y is not greater or equal than low.
 
fvar< T > log1p(const fvar< T > &x)
 
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
 
constexpr double lgamma_stirling_diff_useful
 
bool is_any_nan(const T &x)
Returns true if the input is NaN and false otherwise.
 
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
 
static constexpr double INFTY
Positive infinity.
 
fvar< T > digamma(const fvar< T > &x)
Return the derivative of the log gamma function at the specified argument.
 
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 ...