1#ifndef STAN_MATH_PRIM_PROB_WIENER_FULL_LPDF_HPP
2#define STAN_MATH_PRIM_PROB_WIENER_FULL_LPDF_HPP
32template <
typename T_y,
typename T_a,
typename T_v,
typename T_w,
typename T_sv,
33 typename T_sw,
typename T_err>
35 const T_w& w,
const T_sv& sv,
const T_sw& sw,
37 auto low = w - sw / 2.0;
38 const auto lower_value
39 = wiener5_density<GradientCalc::ON>(y, a, v, low, sv, log_error);
40 auto high = w + sw / 2.0;
41 const auto upper_value
42 = wiener5_density<GradientCalc::ON>(y, a, v, high, sv, log_error);
43 return 0.5 * (lower_value + upper_value) / sw;
74template <
GradientCalc GradSW,
typename F,
typename T_y,
typename T_a,
75 typename T_v,
typename T_w,
typename T_sv,
typename T_sw,
76 typename T_err, std::enable_if_t<!GradSW>* =
nullptr>
78 T_w&& w, T_sv&& sv, T_sw&& sw,
80 return functor(y_diff, a, v, w, sv, log_error);
111template <
GradientCalc GradSW,
typename F,
typename T_y,
typename T_a,
112 typename T_v,
typename T_w,
typename T_sv,
typename T_sw,
113 typename T_err, std::enable_if_t<GradSW>* =
nullptr>
115 T_w&& w, T_sv&& sv, T_sw&& sw,
117 return functor(y_diff, a, v, w, sv, sw, log_error);
137 typename Wiener7FunctorT,
typename T_err,
typename... TArgs>
139 T_err&& hcubature_err, TArgs&&... args) {
140 const auto functor = [&wiener7_functor](
auto&&... integration_args) {
142 [&wiener7_functor](
auto&& x,
auto&& y,
auto&& a,
auto&& v,
auto&& w,
143 auto&& t0,
auto&& sv,
auto&& sw,
auto&& st0,
145 using ret_t =
return_type_t<
decltype(x),
decltype(a),
decltype(v),
146 decltype(w),
decltype(t0),
decltype(sv),
147 decltype(sw),
decltype(st0),
148 decltype(st0),
decltype(log_error)>;
150 const auto sw_val = GradSW ? 0 : sw;
151 const auto new_t0 = t0 + st0 * x_vec[(sw_val != 0) ? 1 : 0];
152 if (y - new_t0 <= 0) {
155 const auto new_w = w + sw_val * (x_vec[0] - 0.5);
156 return ret_t(conditionally_grad_sw<GradSW>(
157 wiener7_functor, y - new_t0, a, v, new_w, sv, sw, log_error));
160 integration_args...);
162 return estimate_with_err_check<0, 8, GradW7, GradientCalc::ON>(
163 functor, hcubature_err, args...);
317template <
bool propto =
false,
typename T_y,
typename T_a,
typename T_t0,
318 typename T_w,
typename T_v,
typename T_sv,
typename T_sw,
320inline auto wiener_lpdf(
const T_y& y,
const T_a& a,
const T_t0& t0,
321 const T_w& w,
const T_v& v,
const T_sv& sv,
322 const T_sw& sw,
const T_st0& st0,
323 const double& precision_derivatives = 1
e-4) {
339 using T_partials_return
342 static constexpr const char* function_name =
"wiener_lpdf";
344 "Boundary separation", a,
"Drift rate", v,
345 "A-priori bias", w,
"Nondecision time", t0,
346 "Inter-trial variability in drift rate", sv,
347 "Inter-trial variability in A-priori bias", sw,
348 "Inter-trial variability in Nondecision time", st0);
354 T_t0_ref t0_ref = t0;
355 T_sv_ref sv_ref = sv;
356 T_sw_ref sw_ref = sw;
357 T_st0_ref st0_ref = st0;
370 check_less(function_name,
"A-priori bias", w_val, 1);
373 check_finite(function_name,
"Nondecision time", t0_val);
376 check_finite(function_name,
"Inter-trial variability in drift rate", sv_val);
377 check_bounded(function_name,
"Inter-trial variability in A-priori bias",
380 "Inter-trial variability in Nondecision time", st0_val);
381 check_finite(function_name,
"Inter-trial variability in Nondecision time",
384 const size_t N =
max_size(y, a, v, w, t0, sv, sw, st0);
396 const size_t N_y_t0 =
max_size(y, t0, st0);
398 for (
size_t i = 0; i < N_y_t0; ++i) {
399 if (y_vec[i] <= t0_vec[i]) {
401 std::stringstream msg;
402 msg <<
", but must be greater than nondecision time = " << t0_vec[i];
403 std::string msg_str(msg.str());
410 for (
size_t i = 0; i < N_beta_sw; ++i) {
411 if (
unlikely(w_vec[i] - .5 * sw_vec[i] <= 0)) {
413 std::stringstream msg;
414 msg <<
", but must be smaller than 2*(A-priori bias) = "
416 std::string msg_str(msg.str());
418 "Inter-trial variability in A-priori bias",
419 sw_vec[i],
" = ", msg_str.c_str());
422 if (
unlikely(w_vec[i] + .5 * sw_vec[i] >= 1)) {
424 std::stringstream msg;
425 msg <<
", but must be smaller than 2*(1-A-priori bias) = "
426 << 2 * (1 - w_vec[i]);
427 std::string msg_str(msg.str());
429 "Inter-trial variability in A-priori bias",
430 sw_vec[i],
" = ", msg_str.c_str());
435 const T_partials_return log_error_density =
log(1
e-6);
437 const auto error_bound = precision_derivatives;
438 const auto log_error_derivative =
log(error_bound);
439 const T_partials_return absolute_error_hcubature = 0.0;
441 const T_partials_return relative_error_hcubature = .9 * error_bound;
442 const T_partials_return log_error_absolute =
log(1
e-12);
443 const int maximal_evaluations_hcubature = 6000;
444 T_partials_return log_density = 0.0;
446 v_ref, sv_ref, sw_ref, st0_ref);
450 for (
size_t i = 0; i < N; i++) {
451 if (sw_vec[i] == 0 && st0_vec[i] == 0) {
455 result += wiener_lpdf<propto>(y_vec[i], a_vec[i], t0_vec[i], w_vec[i],
456 v_vec[i], sv_vec[i], precision_derivatives);
459 const T_partials_return y_value = y_vec.val(i);
460 const T_partials_return a_value = a_vec.val(i);
461 const T_partials_return v_value = v_vec.val(i);
462 const T_partials_return w_value = w_vec.val(i);
463 const T_partials_return t0_value = t0_vec.val(i);
464 const T_partials_return sv_value = sv_vec.val(i);
465 const T_partials_return sw_value = sw_vec.val(i);
466 const T_partials_return st0_value = st0_vec.val(i);
467 const int dim = (sw_value != 0) + (st0_value != 0);
469 "(Inter-trial variability in A-priori bias) + "
470 "(Inter-trial variability in nondecision time)",
473 Eigen::Matrix<T_partials_return, -1, 1> xmin = Eigen::VectorXd::Zero(dim);
474 Eigen::Matrix<T_partials_return, -1, 1> xmax = Eigen::VectorXd::Ones(dim);
475 if (st0_value != 0) {
476 xmax[dim - 1] =
fmin(1.0, (y_value - t0_value) / st0_value);
479 T_partials_return hcubature_err
480 = log_error_absolute - log_error_density +
LOG_TWO + 1;
482 const auto params = std::make_tuple(y_value, a_value, v_value, w_value,
483 t0_value, sv_value, sw_value, st0_value,
485 T_partials_return density
486 = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::OFF>(
488 return internal::wiener5_density<GradientCalc::ON>(args...);
490 hcubature_err, params, dim, xmin, xmax,
491 maximal_evaluations_hcubature, absolute_error_hcubature,
492 relative_error_hcubature / 2);
493 log_density +=
log(density);
494 hcubature_err = log_error_absolute - log_error_derivative
499 const T_partials_return deriv_t_7
500 = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::OFF>(
502 return internal::wiener5_grad_t<GradientCalc::ON>(args...);
504 hcubature_err, params, dim, xmin, xmax,
505 maximal_evaluations_hcubature, absolute_error_hcubature,
506 relative_error_hcubature / 2)
512 partials<0>(ops_partials)[i] = deriv_t_7;
515 partials<1>(ops_partials)[i]
516 = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::OFF>(
518 return internal::wiener5_grad_a<GradientCalc::ON>(args...);
520 hcubature_err, params, dim, xmin, xmax,
521 maximal_evaluations_hcubature, absolute_error_hcubature,
522 relative_error_hcubature / 2)
526 partials<2>(ops_partials)[i] = -deriv_t_7;
529 partials<3>(ops_partials)[i]
530 = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::ON>(
532 return internal::wiener5_grad_w<GradientCalc::ON>(args...);
534 hcubature_err, params, dim, xmin, xmax,
535 maximal_evaluations_hcubature, absolute_error_hcubature,
536 relative_error_hcubature / 2)
540 partials<4>(ops_partials)[i]
541 = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::OFF>(
543 return internal::wiener5_grad_v<GradientCalc::ON>(args...);
545 hcubature_err, params, dim, xmin, xmax,
546 maximal_evaluations_hcubature, absolute_error_hcubature,
547 relative_error_hcubature / 2)
551 partials<5>(ops_partials)[i]
552 = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::OFF>(
554 return internal::wiener5_grad_sv<GradientCalc::ON>(args...);
556 hcubature_err, params, dim, xmin, xmax,
557 maximal_evaluations_hcubature, absolute_error_hcubature,
558 relative_error_hcubature / 2)
563 partials<6>(ops_partials)[i] = 0;
565 if (st0_value == 0) {
567 6, 0, GradientCalc::OFF, GradientCalc::ON>(
569 hcubature_err, y_value - t0_value, a_value, v_value, w_value,
570 sv_value, sw_value, log_error_absolute -
LOG_TWO);
575 hcubature_err, params, 1, xmin, xmax,
576 maximal_evaluations_hcubature, absolute_error_hcubature,
577 relative_error_hcubature / 2);
579 partials<6>(ops_partials)[i] =
derivative / density - 1.0 / sw_value;
584 if (st0_value == 0) {
585 partials<7>(ops_partials)[i] = 0;
586 }
else if (y_value - (t0_value + st0_value) <= 0) {
587 partials<7>(ops_partials)[i] = -1 / st0_value;
589 const T_partials_return t0_st0 = t0_value + st0_value;
594 return internal::wiener5_density<GradientCalc::ON>(args...);
596 log_error_derivative +
log(st0_value), y_value - t0_st0, a_value,
597 v_value, w_value, sv_value, log_error_absolute -
LOG_TWO);
599 const T_partials_return new_error = log_error_absolute -
LOG_TWO;
601 = std::make_tuple(y_value, a_value, v_value, w_value, t0_st0,
602 sv_value, sw_value, 0.0, new_error);
603 f = internal::wiener7_integrate<GradientCalc::OFF, GradientCalc::OFF>(
605 return internal::wiener5_density<GradientCalc::ON>(args...);
607 hcubature_err, params_st, 1, xmin, xmax,
608 maximal_evaluations_hcubature, absolute_error_hcubature,
609 relative_error_hcubature / 2.0);
611 partials<7>(ops_partials)[i] = -1 / st0_value + f / st0_value / density;
615 return result + ops_partials.build(log_density);
scalar_seq_view provides a uniform sequence-like wrapper around either a scalar or a sequence of scal...
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
auto wiener7_integrate(const Wiener7FunctorT &wiener7_functor, T_err &&hcubature_err, TArgs &&... args)
Implementation function for preparing arguments and functor to be passed to the hcubature() function ...
auto wiener7_grad_sw(const T_y &y, const T_a &a, const T_v &v, const T_w &w, const T_sv &sv, const T_sw &sw, T_err log_error)
Calculate the derivative of the wiener7 density w.r.t.
auto conditionally_grad_sw(F &&functor, T_y &&y_diff, T_a &&a, T_v &&v, T_w &&w, T_sv &&sv, T_sw &&sw, T_err &&log_error)
Helper function for agnostically calling wiener5 functions (to be integrated over) or directly callin...
auto estimate_with_err_check(F &&functor, T_err &&err, ArgsTupleT &&... args_tuple)
Utility function for estimating a function with a given set of arguments, checking the result against...
void derivative(const F &f, const T &x, T &fx, T &dfx_dx)
Return the derivative of the specified univariate function at the specified argument.
void check_nonnegative(const char *function, const char *name, const T_y &y)
Check if y is non-negative.
fvar< T > fmin(const fvar< T > &x1, const fvar< T > &x2)
void check_bounded(const char *function, const char *name, const T_y &y, const T_low &low, const T_high &high)
Check if the value is between the low and high values, inclusively.
static constexpr double e()
Return the base of the natural logarithm.
auto hcubature(const F &integrand, const ParsTuple &pars, const int dim, const Eigen::Matrix< T_a, Eigen::Dynamic, 1 > &a, const Eigen::Matrix< T_b, Eigen::Dynamic, 1 > &b, const int max_eval, const TAbsErr reqAbsError, const TRelErr reqRelError)
Compute the [dim]-dimensional integral of the function from to within specified relative and absol...
fvar< T > log(const fvar< T > &x)
void throw_domain_error(const char *function, const char *name, const T &y, const char *msg1, const char *msg2)
Throw a domain error with a consistently formatted message.
static constexpr double LOG_TWO
The natural logarithm of 2, .
auto as_value_column_array_or_scalar(T &&a)
Extract the value from an object and for eigen vectors and std::vectors convert to an eigen column ar...
void check_consistent_sizes(const char *)
Trivial no input case, this function is a no-op.
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
void check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
auto wiener_lpdf(const T_y &y, const T_a &a, const T_t0 &t0, const T_w &w, const T_v &v, const T_sv &sv, const double &precision_derivatives=1e-4)
Log-density function for the 5-parameter Wiener density.
void check_less(const char *function, const char *name, const T_y &y, const T_high &high, Idxs... idxs)
Throw an exception if y is not strictly less than high.
int64_t max_size(const T1 &x1, const Ts &... xs)
Calculate the size of the largest input.
void check_greater(const char *function, const char *name, const T_y &y, const T_low &low, Idxs... idxs)
Throw an exception if y is not strictly greater than low.
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
void check_positive_finite(const char *function, const char *name, const T_y &y)
Check if y is positive and finite.
fvar< T > fabs(const fvar< T > &x)
typename ref_type_if< true, T >::type ref_type_t
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 ...
Extends std::true_type when instantiated with zero or more template parameters, all of which extend t...
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...