Automatic Differentiation
 
Loading...
Searching...
No Matches
std_normal_lcdf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_OPENCL_KERNELS_DEVICE_FUNCTIONS_STD_NORMAL_LCDF_HPP
2#define STAN_MATH_OPENCL_KERNELS_DEVICE_FUNCTIONS_STD_NORMAL_LCDF_HPP
3#ifdef STAN_OPENCL
4
6#include <string>
7
8namespace stan {
9namespace math {
10namespace opencl_kernels {
11// \cond
12static constexpr const char* std_normal_lcdf_device_function
13 = "\n"
14 "#ifndef STAN_MATH_OPENCL_KERNELS_DEVICE_FUNCTIONS_STD_NORMAL_LCDF\n"
15 "#define "
16 "STAN_MATH_OPENCL_KERNELS_DEVICE_FUNCTIONS_STD_NORMAL_LCDF\n" STRINGIFY(
24 inline double std_normal_lcdf_scaled_impl(double scaled_y) {
25 double lcdf_n;
26 if (scaled_y > 0.0) {
27 // CDF(x) = 1/2 + 1/2 erf(x) = 1 - 1/2 erfc(x)
28 lcdf_n = log1p(-0.5 * erfc(scaled_y));
29 if (isnan(lcdf_n)) {
30 lcdf_n = 0;
31 }
32 } else if (scaled_y > -20.0) {
33 // CDF(x) = 1/2 - 1/2 erf(-x) = 1/2 erfc(-x)
34 lcdf_n = log(erfc(-scaled_y)) - M_LN2;
35 } else if (10.0 * log(fabs(scaled_y)) < log(DBL_MAX)) {
36 // Need direct approximation once erfc(-x) underflows.
37 const double x2 = scaled_y * scaled_y;
38 const double x4 = pow(scaled_y, 4);
39 const double x6 = pow(scaled_y, 6);
40 const double x8 = pow(scaled_y, 8);
41 const double x10 = pow(scaled_y, 10);
42 const double temp_p = 0.000658749161529837803157
43 + 0.0160837851487422766278 / x2
44 + 0.125781726111229246204 / x4
45 + 0.360344899949804439429 / x6
46 + 0.305326634961232344035 / x8
47 + 0.0163153871373020978498 / x10;
48 const double temp_q
49 = -0.00233520497626869185443 - 0.0605183413124413191178 / x2
50 - 0.527905102951428412248 / x4 - 1.87295284992346047209 / x6
51 - 2.56852019228982242072 / x8 - 1.0 / x10;
52 lcdf_n = log(0.5 * M_2_SQRTPI + (temp_p / temp_q) / x2) - M_LN2
53 - log(-scaled_y) - x2;
54 } else {
55 lcdf_n = -INFINITY;
56 }
57 return lcdf_n;
58 }
59
68 inline double std_normal_lcdf_dscaled_impl(double scaled_y) {
69 double dnlcdf = 0.0;
70 double t = 0.0;
71 double t2 = 0.0;
72 double t4 = 0.0;
73 const double x2 = scaled_y * scaled_y;
74
75 if (scaled_y > 2.9) {
76 t = 1.0 / (1.0 + 0.3275911 * scaled_y);
77 t2 = t * t;
78 t4 = pow(t, 4);
79 dnlcdf = 0.5 * M_2_SQRTPI
80 / (exp(x2) - 0.254829592 + 0.284496736 * t
81 - 1.421413741 * t2 + 1.453152027 * t2 * t
82 - 1.061405429 * t4);
83 } else if (scaled_y > 2.5) {
84 t = scaled_y - 2.7;
85 t2 = t * t;
86 t4 = pow(t, 4);
87 dnlcdf = 0.0003849882382 - 0.002079084702 * t
88 + 0.005229340880 * t2 - 0.008029540137 * t2 * t
89 + 0.008232190507 * t4 - 0.005692364250 * t4 * t
90 + 0.002399496363 * pow(t, 6);
91 } else if (scaled_y > 2.1) {
92 t = scaled_y - 2.3;
93 t2 = t * t;
94 t4 = pow(t, 4);
95 dnlcdf = 0.002846135439 - 0.01310032351 * t + 0.02732189391 * t2
96 - 0.03326906904 * t2 * t + 0.02482478940 * t4
97 - 0.009883071924 * t4 * t - 0.0002771362254 * pow(t, 6);
98 } else if (scaled_y > 1.5) {
99 t = scaled_y - 1.85;
100 t2 = t * t;
101 t4 = pow(t, 4);
102 dnlcdf = 0.01849212058 - 0.06876280470 * t + 0.1099906382 * t2
103 - 0.09274533184 * t2 * t + 0.03543327418 * t4
104 + 0.005644855518 * t4 * t - 0.01111434424 * pow(t, 6);
105 } else if (scaled_y > 0.8) {
106 t = scaled_y - 1.15;
107 t2 = t * t;
108 t4 = pow(t, 4);
109 dnlcdf = 0.1585747034 - 0.3898677543 * t + 0.3515963775 * t2
110 - 0.09748053605 * t2 * t - 0.04347986191 * t4
111 + 0.02182506378 * t4 * t + 0.01074751427 * pow(t, 6);
112 } else if (scaled_y > 0.1) {
113 t = scaled_y - 0.45;
114 t2 = t * t;
115 t4 = pow(t, 4);
116 dnlcdf = 0.6245634904 - 0.9521866949 * t + 0.3986215682 * t2
117 + 0.04700850676 * t2 * t - 0.03478651979 * t4
118 - 0.01772675404 * t4 * t + 0.0006577254811 * pow(t, 6);
119 } else if (10.0 * log(fabs(scaled_y)) < log(DBL_MAX)) {
120 t = 1.0 / (1.0 - 0.3275911 * scaled_y);
121 t2 = t * t;
122 t4 = pow(t, 4);
123 dnlcdf
124 = M_2_SQRTPI
125 / (0.254829592 * t - 0.284496736 * t2 + 1.421413741 * t2 * t
126 - 1.453152027 * t4 + 1.061405429 * t4 * t);
127 if (scaled_y < -29.0) {
128 dnlcdf += 0.0015065154280332 * x2
129 - 0.3993154819705530 * scaled_y - 4.2919418242931700;
130 } else if (scaled_y < -17.0) {
131 dnlcdf += 0.0001263257217272 * x2 * scaled_y
132 + 0.0123586859488623 * x2
133 - 0.0860505264736028 * scaled_y - 1.252783383752970;
134 } else if (scaled_y < -7.0) {
135 dnlcdf += 0.000471585349920831 * x2 * scaled_y
136 + 0.0296839305424034 * x2
137 + 0.207402143352332 * scaled_y + 0.425316974683324;
138 } else if (scaled_y < -3.9) {
139 dnlcdf += -0.0006972280656443 * x2 * scaled_y
140 + 0.0068218494628567 * x2
141 + 0.0585761964460277 * scaled_y + 0.1034397670201370;
142 } else if (scaled_y < -2.1) {
143 dnlcdf += -0.0018742199480885 * x2 * scaled_y
144 - 0.0097119598291202 * x2
145 - 0.0170137970924080 * scaled_y - 0.0100428567412041;
146 }
147 } else {
148 dnlcdf = INFINITY;
149 }
150
151 return dnlcdf;
152 }) "\n#endif\n"; // NOLINT
153// \endcond
154
155} // namespace opencl_kernels
156} // namespace math
157} // namespace stan
158
159#endif
160#endif
std_normal_lcdf_dscaled_impl_< as_operation_cl_t< T > > std_normal_lcdf_dscaled_impl(T &&a)
std_normal_lcdf_scaled_impl_< as_operation_cl_t< T > > std_normal_lcdf_scaled_impl(T &&a)
bool isnan(double_d a)
Definition double_d.hpp:327
auto pow(const T1 &x1, const T2 &x2)
Definition pow.hpp:32
fvar< T > log(const fvar< T > &x)
Definition log.hpp:18
fvar< T > erfc(const fvar< T > &x)
Definition erfc.hpp:16
fvar< T > log1p(const fvar< T > &x)
Definition log1p.hpp:12
fvar< T > fabs(const fvar< T > &x)
Definition fabs.hpp:16
fvar< T > exp(const fvar< T > &x)
Definition exp.hpp:15
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
#define STRINGIFY(...)
Definition stringify.hpp:9