00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
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* ) 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
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 ) 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
00128
00129
00130
00131
00132
00133
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