RooGamma.cxx

Go to the documentation of this file.
00001  /***************************************************************************** 
00002   * Project: RooFit                                                           * 
00003   *                                                                           * 
00004   * Simple Gamma distribution 
00005   * authors: Stefan A. Schmitz, Gregory Schott
00006   *                                                                           * 
00007   *****************************************************************************/ 
00008 
00009 //
00010 // implementation of the Gamma PDF for RooFit/RooStats
00011 // $f(x) = \frac{(x-\mu)^{\gamma-1} \cdot \exp^{(-(x-mu) / \beta}}{\Gamma(\gamma) \cdot \beta^{\gamma}}$
00012 // defined for $x \geq 0$ if $\mu = 0$
00013 //
00014 
00015 // Notes from Kyle Cranmer
00016 // Wikipedia and several sources refer to the Gamma distribution as
00017 // G(mu|alpha,beta) = beta^alpha mu^(alpha-1) e^(-beta mu) / Gamma(alpha)
00018 // Below is the correspondance
00019 // Wikipedia | This Implementation
00020 //---------------------------------
00021 // alpha     | gamma
00022 // beta      | 1/beta
00023 // mu        | x
00024 // 0         | mu
00025 //
00026 // Note, that for a model Pois(N|mu), a uniform prior on mu, and a measurement N
00027 // the posterior is in the Wikipedia parametrization Gamma(mu, alpha=N+1, beta=1)
00028 // thus for this implementation it is:
00029 // RooGamma(_x=mu,_gamma=N+1,_beta=1,_mu=0)
00030 // Also note, that in this case it is equivalent to
00031 // RooPoison(N,mu) and treating the function as a PDF in mu.
00032 
00033 
00034 #include "RooFit.h"
00035 
00036 #include "Riostream.h"
00037 #include "Riostream.h"
00038 #include <math.h>
00039 
00040 #include "RooGamma.h"
00041 #include "RooAbsReal.h"
00042 #include "RooRealVar.h"
00043 #include "RooRandom.h"
00044 #include "RooMath.h"
00045 
00046 #include <iostream> 
00047 #include "TMath.h"
00048 
00049 #include <Math/SpecFuncMathCore.h>
00050 #include <Math/PdfFuncMathCore.h>
00051 #include <Math/ProbFuncMathCore.h>
00052 
00053 ClassImp(RooGamma)
00054 
00055 
00056 //_____________________________________________________________________________
00057 RooGamma::RooGamma(const char *name, const char *title,
00058                          RooAbsReal& _x, RooAbsReal& _gamma,
00059                          RooAbsReal& _beta, RooAbsReal& _mu) :
00060   RooAbsPdf(name,title),
00061   x("x","Observable",this,_x),
00062   gamma("gamma","Mean",this,_gamma),
00063   beta("beta","Width",this,_beta),
00064   mu("mu","Para",this,_mu)
00065 {
00066 }
00067 
00068 
00069 
00070 //_____________________________________________________________________________
00071 RooGamma::RooGamma(const RooGamma& other, const char* name) : 
00072   RooAbsPdf(other,name), x("x",this,other.x), gamma("mean",this,other.gamma),
00073   beta("beta",this,other.beta), mu("mu",this,other.mu)
00074 {
00075 }
00076 
00077 
00078 
00079 //_____________________________________________________________________________
00080 Double_t RooGamma::evaluate() const
00081 {
00082 
00083   Double_t arg= x ;  
00084   Double_t ret = TMath::GammaDist(arg, gamma, mu, beta) ;
00085   return ret ;
00086 }
00087 
00088 
00089 
00090 //_____________________________________________________________________________
00091 Int_t RooGamma::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const 
00092 {
00093   if (matchArgs(allVars,analVars,x)) return 1 ;
00094   return 0 ;
00095 }
00096 
00097 
00098 
00099 //_____________________________________________________________________________
00100 Double_t RooGamma::analyticalIntegral(Int_t code, const char* rangeName) const 
00101 {
00102   assert(code==1) ;
00103 
00104  //integral of the Gamma distribution via ROOT::Math
00105   Double_t integral = ROOT::Math::gamma_cdf(x.max(rangeName), gamma, beta, mu) - ROOT::Math::gamma_cdf(x.min(rangeName), gamma, beta, mu);
00106   return integral ;
00107 
00108 
00109 }
00110 
00111 
00112 
00113 
00114 //_____________________________________________________________________________
00115 Int_t RooGamma::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
00116 {
00117   if (matchArgs(directVars,generateVars,x)) return 1 ;  
00118   return 0 ;
00119 }
00120 
00121 
00122 
00123 //_____________________________________________________________________________
00124 void RooGamma::generateEvent(Int_t code)
00125 {
00126   assert(code==1) ;
00127 //algorithm adapted from code example in:
00128 //Marsaglia, G. and Tsang, W. W.
00129 //A Simple Method for Generating Gamma Variables
00130 //ACM Transactions on Mathematical Software, Vol. 26, No. 3, September 2000
00131 //
00132 //The speed of this algorithm depends on the speed of generating normal variates.
00133 //The algorithm is limited to $\gamma \geq 0$ !
00134 
00135   while(1) {
00136 
00137   double d = 0;
00138   double c = 0;
00139   double xgen = 0;
00140   double v = 0;
00141   double u = 0; 
00142   d = gamma -1./3.; c = 1./TMath::Sqrt(9.*d);
00143 
00144   while(v <= 0.){
00145     xgen = RooRandom::randomGenerator()->Gaus(); v = 1. + c*xgen;
00146   }
00147   v = v*v*v; u = RooRandom::randomGenerator()->Uniform();
00148   if( u < 1.-.0331*(xgen*xgen)*(xgen*xgen) ) {
00149     if ( (((d*v)* beta + mu ) < x.max()) && (((d*v)* beta + mu) > x.min()) ) {
00150       x = ((d*v)* beta + mu) ;
00151       break;
00152     } 
00153   } 
00154   if( TMath::Log(u) < 0.5*xgen*xgen + d*(1.-v + TMath::Log(v)) ) { 
00155     if ( (((d*v)* beta + mu ) < x.max()) && (((d*v)* beta + mu) > x.min()) ) {
00156       x = ((d*v)* beta + mu) ;
00157       break;
00158     } 
00159   } 
00160 
00161   }
00162 
00163 
00164   return;
00165 }
00166 
00167 

Generated on Tue Jul 5 14:55:12 2011 for ROOT_528-00b_version by  doxygen 1.5.1