mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
gammaDistribution.hpp
Go to the documentation of this file.
1 /** \file gammaDistribution.hpp
2  * \brief The Gamma Distribution.
3  * \ingroup gen_math_files
4  * \author Jared R. Males (jaredmales@gmail.com)
5  *
6  */
7 
8 //***********************************************************************//
9 // Copyright 2023 Jared R. Males (jaredmales@gmail.com)
10 //
11 // This file is part of mxlib.
12 //
13 // mxlib is free software: you can redistribute it and/or modify
14 // it under the terms of the GNU General Public License as published by
15 // the Free Software Foundation, either version 3 of the License, or
16 // (at your option) any later version.
17 //
18 // mxlib is distributed in the hope that it will be useful,
19 // but WITHOUT ANY WARRANTY; without even the implied warranty of
20 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 // GNU General Public License for more details.
22 //
23 // You should have received a copy of the GNU General Public License
24 // along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25 //***********************************************************************//
26 
27 #ifndef gammaDistribution_hpp
28 #define gammaDistribution_hpp
29 
30 #include <cmath>
31 
32 #include "gamma.hpp"
33 
34 namespace mx
35 {
36 namespace math
37 {
38 namespace func
39 {
40 
41 /// The denominator of the Gamma Distribution
42 /** Can be used to avoid repeated calculations when the parameters are constant
43  *
44  * \tparam realT a real floating point type
45  *
46  * \returns the denominator of the Gamma Distribution.
47  *
48  * \ingroup gen_math_gammaDist
49  */
50 template<typename realT>
51 realT gammaDistributionDenom( realT k, ///< [in] shape parameter
52  realT q ///< [in] the scale parameter
53  )
54 {
55  try
56  {
57  realT d = pow(q, k) * tgamma<realT>(k);
58 
59  if(!std::isnormal(d)) return std::numeric_limits<realT>::max();
60 
61  return d;
62  }
63  catch(...)
64  {
65  return std::numeric_limits<realT>::max();
66  }
67 }
68 
69 ///The general shifted Gamma Distribution at a point using an arbitrary peak scaling parameter
70 /** Calculates the value of the Gamma Distribution at a location specified by x.
71  *
72  * \tparam realT a real floating point type
73  *
74  * \returns the value of the Gamma distribution at x.
75  *
76  * \ingroup gen_math_gammaDist
77  */
78 template<typename realT>
79 realT gammaDistribution( realT x, ///< [in] the location at which to calculate the distribution
80  realT x0, ///< [in] the location parameter
81  realT k, ///< [in] shape parameter
82  realT q, ///< [in] the scale parameter
83  realT denom ///< [in] the denominator, or 1/peak-scale.
84  )
85 {
86  if(x - x0 < 0) return 0;
87 
88  realT v = pow(x-x0, k-1) * exp( -(x-x0)/q) / denom;
89 
90  if(!std::isnormal(v)) return 0;
91 
92  return v;
93 
94 }
95 
96 ///The general shifted Gamma Distribution at a point.
97 /** Calculates the value of the Gamma Distribution at a location specified by x.
98  *
99  * \tparam realT a real floating point type
100  *
101  * \returns the value of the Gamma distribution at x.
102  *
103  * \ingroup gen_math_gammaDist
104  */
105 template<typename realT>
106 realT gammaDistribution( realT x, ///< [in] the location at which to calculate the distribution
107  realT x0, ///< [in] the location parameter
108  realT k, ///< [in] shape parameter
109  realT q ///< [in] the scale parameter
110  )
111 {
112 
113  return gammaDistribution<realT>(x, x0, k, q, gammaDistributionDenom<realT>(k,q));
114 }
115 
116 /// The mean of the Gamma Distribution
117 /** Calculates the mean of the Gamma Distribution for the given parameters.
118  *
119  * \tparam realT a real floating point type
120  *
121  * \returns the mean of the Gamma Distribution.
122  *
123  * \ingroup gen_math_gammaDist
124  */
125 template<typename realT>
126 realT gammaDistributionMean( realT x0, ///< [in] the location parameter
127  realT k, ///< [in] shape parameter
128  realT theta ///< [in] the scale parameter
129  )
130 {
131  return x0 + k*theta;
132 }
133 
134 /// The mode of the Gamma Distribution
135 /** Calculates the mode of the Gamma Distribution for the given parameters.
136  *
137  * \tparam realT a real floating point type
138  *
139  * \returns the mode of the Gamma Distribution.
140  *
141  * \ingroup gen_math_gammaDist
142  */
143 template<typename realT>
144 realT gammaDistributionMode( realT x0, ///< [in] the location parameter
145  realT k, ///< [in] shape parameter
146  realT theta ///< [in] the scale parameter
147  )
148 {
149  if(k >= 1) return x0 + (k-1)*theta;
150 
151  return 0;
152 }
153 
154 /// The variance of the Gamma Distribution
155 /** Calculates the variance of the Gamma Distribution for the given parameters.
156  *
157  * \tparam realT a real floating point type
158  *
159  * \returns the variance of the Gamma Distribution.
160  *
161  * \ingroup gen_math_gammaDist
162  */
163 template<typename realT>
164 realT gammaDistributionVariance( realT k, ///< [in] shape parameter
165  realT theta ///< [in] the scale parameter
166  )
167 {
168  return k*theta*theta;
169 }
170 
171 } //namespace func
172 } //namespace math
173 } //namespace mx
174 
175 #endif //gammaDistribution_hpp
Declares and defines the gamma function.
constexpr units::realT k()
Boltzmann Constant.
Definition: constants.hpp:71
realT gammaDistributionMode(realT x0, realT k, realT theta)
The mode of the Gamma Distribution.
realT gammaDistributionDenom(realT k, realT q)
The denominator of the Gamma Distribution.
realT gammaDistributionVariance(realT k, realT theta)
The variance of the Gamma Distribution.
realT gammaDistribution(realT x, realT x0, realT k, realT q)
The general shifted Gamma Distribution at a point.
realT gammaDistributionMean(realT x0, realT k, realT theta)
The mean of the Gamma Distribution.
The mxlib c++ namespace.
Definition: mxError.hpp:107