Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooGamma.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2 * Project: RooFit *
3 * *
4 * Simple Gamma distribution
5 * authors: Stefan A. Schmitz, Gregory Schott
6 * *
7 *****************************************************************************/
8
9/** \class RooGamma
10 \ingroup Roofit
11
12Implementation of the Gamma PDF for RooFit/RooStats.
13\f[
14f(x) = \frac{(x-\mu)^{\gamma-1} \cdot \exp^{(-(x-mu) / \beta}}{\Gamma(\gamma) \cdot \beta^{\gamma}}
15\f]
16defined for \f$ x \geq 0 \f$ if \f$ \mu = 0 \f$
17
18Notes from Kyle Cranmer:
19
20Wikipedia and several sources refer to the Gamma distribution as
21
22\f[
23G(\mu,\alpha,\beta) = \beta^\alpha \mu^{(\alpha-1)} \frac{e^{(-\beta \mu)}}{\Gamma(\alpha)}
24\f]
25
26Below is the correspondence:
27
28| Wikipedia | This Implementation |
29|-----------------|--------------------------|
30| \f$ \alpha \f$ | \f$ \gamma \f$ |
31| \f$ \beta \f$ | \f$ \frac{1}{\beta} \f$ |
32| \f$ \mu \f$ | x |
33| 0 | \f$ \mu \f$ |
34
35
36Note, that for a model Pois(N|mu), a uniform prior on mu, and a measurement N
37the posterior is in the Wikipedia parameterization Gamma(mu, alpha=N+1, beta=1)
38thus for this implementation it is:
39
40`RooGamma(_x=mu,_gamma=N+1,_beta=1,_mu=0)`
41
42Also note, that in this case it is equivalent to
43RooPoison(N,mu) and treating the function as a PDF in mu.
44**/
45
46#include "RooGamma.h"
47
48#include "RooRandom.h"
49#include "RooHelpers.h"
50#include "RooBatchCompute.h"
51
52#include "TMath.h"
54
55#include <cmath>
56
57
58////////////////////////////////////////////////////////////////////////////////
59
60RooGamma::RooGamma(const char *name, const char *title,
61 RooAbsReal& _x, RooAbsReal& _gamma,
62 RooAbsReal& _beta, RooAbsReal& _mu) :
63 RooAbsPdf(name,title),
64 x("x","Observable",this,_x),
65 gamma("gamma","Mean",this,_gamma),
66 beta("beta","Width",this,_beta),
67 mu("mu","Para",this,_mu)
68{
69 RooHelpers::checkRangeOfParameters(this, {&_gamma, &_beta}, 0.);
70}
71
72////////////////////////////////////////////////////////////////////////////////
73
74RooGamma::RooGamma(const RooGamma& other, const char* name) :
75 RooAbsPdf(other,name), x("x",this,other.x), gamma("gamma",this,other.gamma),
76 beta("beta",this,other.beta), mu("mu",this,other.mu)
77{
78}
79
80////////////////////////////////////////////////////////////////////////////////
81
82double RooGamma::evaluate() const
83{
84 return TMath::GammaDist(x, gamma, mu, beta) ;
85}
86
87////////////////////////////////////////////////////////////////////////////////
88/// Compute multiple values of Gamma PDF.
90{
92 {ctx.at(x), ctx.at(gamma), ctx.at(beta), ctx.at(mu)});
93}
94
95////////////////////////////////////////////////////////////////////////////////
96
97Int_t RooGamma::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
98{
99 if (matchArgs(allVars,analVars,x)) return 1 ;
100 return 0 ;
101}
102
103////////////////////////////////////////////////////////////////////////////////
104
105double RooGamma::analyticalIntegral(Int_t /*code*/, const char *rangeName) const
106{
107 // integral of the Gamma distribution via ROOT::Math
110}
111
112namespace {
113
114inline double randomGamma(double gamma, double beta, double mu, double xmin, double xmax)
115{
116 while (true) {
117
118 double d = gamma - 1. / 3.;
119 double c = 1. / std::sqrt(9. * d);
120 double xgen = 0;
121 double v = 0;
122
123 while (v <= 0.) {
125 v = 1. + c * xgen;
126 }
127 v = v * v * v;
128 double u = RooRandom::randomGenerator()->Uniform();
129 if (u < 1. - .0331 * (xgen * xgen) * (xgen * xgen)) {
130 double x = ((d * v) * beta + mu);
131 if ((x < xmax) && (x > xmin)) {
132 return x;
133 }
134 }
135 if (std::log(u) < 0.5 * xgen * xgen + d * (1. - v + TMath::Log(v))) {
136 double x = ((d * v) * beta + mu);
137 if ((x < xmax) && (x > xmin)) {
138 return x;
139 }
140 }
141 }
142}
143
144} // namespace
145
146////////////////////////////////////////////////////////////////////////////////
147
149{
151 return 1;
152 return 0;
153}
154
155////////////////////////////////////////////////////////////////////////////////
156/// algorithm adapted from code example in:
157/// Marsaglia, G. and Tsang, W. W.
158/// A Simple Method for Generating Gamma Variables
159/// ACM Transactions on Mathematical Software, Vol. 26, No. 3, September 2000
160///
161/// The speed of this algorithm depends on the speed of generating normal variates.
162/// The algorithm is limited to \f$ \gamma \geq 0 \f$ !
163
165{
166 if (gamma >= 1) {
167 x = randomGamma(gamma, beta, mu, x.min(), x.max());
168 return;
169 }
170
171 double xVal = 0.0;
172 bool accepted = false;
173
174 while (!accepted) {
175 double u = RooRandom::randomGenerator()->Uniform();
176 double tmp = randomGamma(1 + gamma, beta, mu, 0, std::numeric_limits<double>::infinity());
177 xVal = tmp * std::pow(u, 1.0 / gamma);
178 accepted = xVal < x.max() && xVal > x.min();
179 }
180
181 x = xVal;
182}
#define d(i)
Definition RSha256.hxx:102
#define c(i)
Definition RSha256.hxx:101
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
Implementation of the Gamma PDF for RooFit/RooStats.
Definition RooGamma.h:20
RooRealProxy beta
Definition RooGamma.h:43
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Definition RooGamma.cxx:148
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Gamma PDF.
Definition RooGamma.cxx:89
RooGamma()
Definition RooGamma.h:22
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition RooGamma.cxx:105
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition RooGamma.cxx:82
RooRealProxy x
Definition RooGamma.h:41
void generateEvent(Int_t code) override
algorithm adapted from code example in: Marsaglia, G.
Definition RooGamma.cxx:164
RooRealProxy gamma
Definition RooGamma.h:42
RooRealProxy mu
Definition RooGamma.h:44
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Definition RooGamma.cxx:97
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double gamma_cdf(double x, double alpha, double theta, double x0=0)
Cumulative distribution function of the gamma distribution (lower tail).
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
void checkRangeOfParameters(const RooAbsReal *callingClass, std::initializer_list< const RooAbsReal * > pars, double min=-std::numeric_limits< double >::max(), double max=std::numeric_limits< double >::max(), bool limitsInAllowedRange=false, std::string const &extraMessage="")
Check if the parameters have a range, and warn if the range extends below / above the set limits.
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:760
Double_t GammaDist(Double_t x, Double_t gamma, Double_t mu=0, Double_t beta=1)
Computes the density function of Gamma distribution at point x.
Definition TMath.cxx:2347