WPILibC++ 2023.4.3
incomplete_gamma.hpp
Go to the documentation of this file.
1/*################################################################################
2 ##
3 ## Copyright (C) 2016-2022 Keith O'Hara
4 ##
5 ## This file is part of the GCE-Math C++ library.
6 ##
7 ## Licensed under the Apache License, Version 2.0 (the "License");
8 ## you may not use this file except in compliance with the License.
9 ## You may obtain a copy of the License at
10 ##
11 ## http://www.apache.org/licenses/LICENSE-2.0
12 ##
13 ## Unless required by applicable law or agreed to in writing, software
14 ## distributed under the License is distributed on an "AS IS" BASIS,
15 ## WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
16 ## See the License for the specific language governing permissions and
17 ## limitations under the License.
18 ##
19 ################################################################################*/
20
21/*
22 * compile-time (regularized) incomplete gamma function
23 */
24
25#ifndef _gcem_incomplete_gamma_HPP
26#define _gcem_incomplete_gamma_HPP
27
28namespace internal
29{
30
31// 50 point Gauss-Legendre quadrature
32
33template<typename T>
34constexpr
35T
36incomplete_gamma_quad_inp_vals(const T lb, const T ub, const int counter)
37noexcept
38{
39 return (ub-lb) * gauss_legendre_50_points[counter] / T(2) + (ub + lb) / T(2);
40}
41
42template<typename T>
43constexpr
44T
45incomplete_gamma_quad_weight_vals(const T lb, const T ub, const int counter)
46noexcept
47{
48 return (ub-lb) * gauss_legendre_50_weights[counter] / T(2);
49}
50
51template<typename T>
52constexpr
53T
54incomplete_gamma_quad_fn(const T x, const T a, const T lg_term)
55noexcept
56{
57 return exp( -x + (a-T(1))*log(x) - lg_term );
58}
59
60template<typename T>
61constexpr
62T
63incomplete_gamma_quad_recur(const T lb, const T ub, const T a, const T lg_term, const int counter)
64noexcept
65{
66 return( counter < 49 ? \
67 // if
70 + incomplete_gamma_quad_recur(lb,ub,a,lg_term,counter+1) :
71 // else
74}
75
76template<typename T>
77constexpr
78T
79incomplete_gamma_quad_lb(const T a, const T z)
80noexcept
81{
82 return( a > T(1000) ? max(T(0),min(z,a) - 11*sqrt(a)) : // break integration into ranges
83 a > T(800) ? max(T(0),min(z,a) - 11*sqrt(a)) :
84 a > T(500) ? max(T(0),min(z,a) - 10*sqrt(a)) :
85 a > T(300) ? max(T(0),min(z,a) - 10*sqrt(a)) :
86 a > T(100) ? max(T(0),min(z,a) - 9*sqrt(a)) :
87 a > T(90) ? max(T(0),min(z,a) - 9*sqrt(a)) :
88 a > T(70) ? max(T(0),min(z,a) - 8*sqrt(a)) :
89 a > T(50) ? max(T(0),min(z,a) - 7*sqrt(a)) :
90 a > T(40) ? max(T(0),min(z,a) - 6*sqrt(a)) :
91 a > T(30) ? max(T(0),min(z,a) - 5*sqrt(a)) :
92 // else
93 max(T(0),min(z,a)-4*sqrt(a)) );
94}
95
96template<typename T>
97constexpr
98T
99incomplete_gamma_quad_ub(const T a, const T z)
100noexcept
101{
102 return( a > T(1000) ? min(z, a + 10*sqrt(a)) :
103 a > T(800) ? min(z, a + 10*sqrt(a)) :
104 a > T(500) ? min(z, a + 9*sqrt(a)) :
105 a > T(300) ? min(z, a + 9*sqrt(a)) :
106 a > T(100) ? min(z, a + 8*sqrt(a)) :
107 a > T(90) ? min(z, a + 8*sqrt(a)) :
108 a > T(70) ? min(z, a + 7*sqrt(a)) :
109 a > T(50) ? min(z, a + 6*sqrt(a)) :
110 // else
111 min(z, a + 5*sqrt(a)) );
112}
113
114template<typename T>
115constexpr
116T
117incomplete_gamma_quad(const T a, const T z)
118noexcept
119{
121}
122
123// reverse cf expansion
124// see: https://functions.wolfram.com/GammaBetaErf/Gamma2/10/0003/
125
126template<typename T>
127constexpr
128T
129incomplete_gamma_cf_2_recur(const T a, const T z, const int depth)
130noexcept
131{
132 return( depth < 100 ? \
133 // if
134 (1 + (depth-1)*2 - a + z) + depth*(a - depth)/incomplete_gamma_cf_2_recur(a,z,depth+1) :
135 // else
136 (1 + (depth-1)*2 - a + z) );
137}
138
139template<typename T>
140constexpr
141T
142incomplete_gamma_cf_2(const T a, const T z)
143noexcept
144{ // lower (regularized) incomplete gamma function
145 return( T(1.0) - exp(a*log(z) - z - lgamma(a)) / incomplete_gamma_cf_2_recur(a,z,1) );
146}
147
148// cf expansion
149// see: http://functions.wolfram.com/GammaBetaErf/Gamma2/10/0009/
150
151template<typename T>
152constexpr
153T
154incomplete_gamma_cf_1_coef(const T a, const T z, const int depth)
155noexcept
156{
157 return( is_odd(depth) ? - (a - 1 + T(depth+1)/T(2)) * z : T(depth)/T(2) * z );
158}
159
160template<typename T>
161constexpr
162T
163incomplete_gamma_cf_1_recur(const T a, const T z, const int depth)
164noexcept
165{
166 return( depth < GCEM_INCML_GAMMA_MAX_ITER ? \
167 // if
168 (a + depth - 1) + incomplete_gamma_cf_1_coef(a,z,depth)/incomplete_gamma_cf_1_recur(a,z,depth+1) :
169 // else
170 (a + depth - 1) );
171}
172
173template<typename T>
174constexpr
175T
176incomplete_gamma_cf_1(const T a, const T z)
177noexcept
178{ // lower (regularized) incomplete gamma function
179 return( exp(a*log(z) - z - lgamma(a)) / incomplete_gamma_cf_1_recur(a,z,1) );
180}
181
182//
183
184template<typename T>
185constexpr
186T
187incomplete_gamma_check(const T a, const T z)
188noexcept
189{
190 return( // NaN check
191 any_nan(a, z) ? \
192 GCLIM<T>::quiet_NaN() :
193 //
194 a < T(0) ? \
195 GCLIM<T>::quiet_NaN() :
196 //
197 GCLIM<T>::min() > z ? \
198 T(0) :
199 //
200 GCLIM<T>::min() > a ? \
201 T(1) :
202 // cf or quadrature
203 (a < T(10)) && (z - a < T(10)) ?
205 (a < T(10)) || (z/a > T(3)) ?
207 // else
209}
210
211template<typename T1, typename T2, typename TC = common_return_t<T1,T2>>
212constexpr
213TC
214incomplete_gamma_type_check(const T1 a, const T2 p)
215noexcept
216{
217 return incomplete_gamma_check(static_cast<TC>(a),
218 static_cast<TC>(p));
219}
220
221}
222
223/**
224 * Compile-time regularized lower incomplete gamma function
225 *
226 * @param a a real-valued, non-negative input.
227 * @param x a real-valued, non-negative input.
228 *
229 * @return the regularized lower incomplete gamma function evaluated at (\c a, \c x),
230 * \f[ \frac{\gamma(a,x)}{\Gamma(a)} = \frac{1}{\Gamma(a)} \int_0^x t^{a-1} \exp(-t) dt \f]
231 * When \c a is not too large, the value is computed using the continued fraction representation of the upper incomplete gamma function, \f$ \Gamma(a,x) \f$, using
232 * \f[ \Gamma(a,x) = \Gamma(a) - \dfrac{x^a\exp(-x)}{a - \dfrac{ax}{a + 1 + \dfrac{x}{a + 2 - \dfrac{(a+1)x}{a + 3 + \dfrac{2x}{a + 4 - \ddots}}}}} \f]
233 * where \f$ \gamma(a,x) \f$ and \f$ \Gamma(a,x) \f$ are connected via
234 * \f[ \frac{\gamma(a,x)}{\Gamma(a)} + \frac{\Gamma(a,x)}{\Gamma(a)} = 1 \f]
235 * When \f$ a > 10 \f$, a 50-point Gauss-Legendre quadrature scheme is employed.
236 */
237
238template<typename T1, typename T2>
239constexpr
240common_return_t<T1,T2>
241incomplete_gamma(const T1 a, const T2 x)
242noexcept
243{
245}
246
247#endif
EIGEN_DEVICE_FUNC const LgammaReturnType lgamma() const
\cpp11
Definition: ArrayCwiseUnaryOps.h:620
static const long double gauss_legendre_50_points[50]
Definition: gauss_legendre_50.hpp:25
static const long double gauss_legendre_50_weights[50]
Definition: gauss_legendre_50.hpp:79
#define GCEM_INCML_GAMMA_MAX_ITER
Definition: gcem_options.hpp:173
auto sqrt(const UnitType &value) noexcept -> unit_t< square_root< typename units::traits::unit_t_traits< UnitType >::unit_type >, typename units::traits::unit_t_traits< UnitType >::underlying_type, linear_scale >
computes the square root of value
Definition: math.h:483
dimensionless::scalar_t log(const ScalarUnit x) noexcept
Compute natural logarithm.
Definition: math.h:349
constexpr common_return_t< T1, T2 > incomplete_gamma(const T1 a, const T2 x) noexcept
Compile-time regularized lower incomplete gamma function.
Definition: incomplete_gamma.hpp:241
constexpr common_t< T1, T2 > max(const T1 x, const T2 y) noexcept
Compile-time pairwise maximum function.
Definition: max.hpp:35
constexpr common_t< T1, T2 > min(const T1 x, const T2 y) noexcept
Compile-time pairwise minimum function.
Definition: min.hpp:35
Definition: Eigen_Colamd.h:50
constexpr T incomplete_gamma_quad_recur(const T lb, const T ub, const T a, const T lg_term, const int counter) noexcept
Definition: incomplete_gamma.hpp:63
constexpr T incomplete_gamma_check(const T a, const T z) noexcept
Definition: incomplete_gamma.hpp:187
constexpr T incomplete_gamma_quad_lb(const T a, const T z) noexcept
Definition: incomplete_gamma.hpp:79
constexpr T incomplete_gamma_cf_1(const T a, const T z) noexcept
Definition: incomplete_gamma.hpp:176
constexpr T incomplete_gamma_quad_ub(const T a, const T z) noexcept
Definition: incomplete_gamma.hpp:99
constexpr T incomplete_gamma_cf_1_coef(const T a, const T z, const int depth) noexcept
Definition: incomplete_gamma.hpp:154
constexpr T incomplete_gamma_quad_fn(const T x, const T a, const T lg_term) noexcept
Definition: incomplete_gamma.hpp:54
constexpr T incomplete_gamma_quad_inp_vals(const T lb, const T ub, const int counter) noexcept
Definition: incomplete_gamma.hpp:36
constexpr T incomplete_gamma_quad(const T a, const T z) noexcept
Definition: incomplete_gamma.hpp:117
constexpr T incomplete_gamma_cf_2_recur(const T a, const T z, const int depth) noexcept
Definition: incomplete_gamma.hpp:129
constexpr T incomplete_gamma_cf_2(const T a, const T z) noexcept
Definition: incomplete_gamma.hpp:142
constexpr bool is_odd(const llint_t x) noexcept
Definition: is_odd.hpp:33
constexpr T incomplete_gamma_quad_weight_vals(const T lb, const T ub, const int counter) noexcept
Definition: incomplete_gamma.hpp:45
constexpr T incomplete_gamma_cf_1_recur(const T a, const T z, const int depth) noexcept
Definition: incomplete_gamma.hpp:163
constexpr bool any_nan(const T1 x, const T2 y) noexcept
Definition: is_nan.hpp:45
constexpr TC incomplete_gamma_type_check(const T1 a, const T2 p) noexcept
Definition: incomplete_gamma.hpp:214
lb
Definition: mass.h:44