RooGExpModel.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooGExpModel.cxx 24286 2008-06-16 15:47:04Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 //
00019 // BEGIN_HTML
00020 // Class RooGExpModel is a RooResolutionModel implementation that models
00021 // a resolution function that is the convolution of a Gaussian with
00022 // a one-sided exponential. Object of class RooGExpModel can be used
00023 // for analytical convolutions with classes inheriting from RooAbsAnaConvPdf
00024 // END_HTML
00025 //
00026 
00027 #include "RooFit.h"
00028 
00029 #include "Riostream.h"
00030 #include "Riostream.h"
00031 #include "RooGExpModel.h"
00032 #include "RooMath.h"
00033 #include "RooRealConstant.h"
00034 #include "RooRandom.h"
00035 #include "RooMath.h"
00036 #include "TMath.h"
00037 
00038 ClassImp(RooGExpModel) 
00039 ;
00040 
00041 
00042 
00043 //_____________________________________________________________________________
00044 RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn, 
00045                            RooAbsReal& _sigma, RooAbsReal& _rlife, 
00046                            Bool_t nlo, Type type) : 
00047   RooResolutionModel(name,title,xIn), 
00048   sigma("sigma","Width",this,_sigma),
00049   rlife("rlife","Life time",this,_rlife),
00050   ssf("ssf","Sigma Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
00051   rsf("rsf","RLife Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
00052   _flip(type==Flipped),_nlo(nlo), _flatSFInt(kFALSE), _asympInt(kFALSE)
00053 {  
00054 }
00055 
00056 
00057 
00058 //_____________________________________________________________________________
00059 RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn, 
00060                            RooAbsReal& _sigma, RooAbsReal& _rlife, 
00061                            RooAbsReal& _rsSF,
00062                            Bool_t nlo, Type type) : 
00063   RooResolutionModel(name,title,xIn), 
00064   sigma("sigma","Width",this,_sigma),
00065   rlife("rlife","Life time",this,_rlife),
00066   ssf("ssf","Sigma Scale Factor",this,_rsSF),
00067   rsf("rsf","RLife Scale Factor",this,_rsSF),
00068   _flip(type==Flipped),
00069   _nlo(nlo), 
00070   _flatSFInt(kFALSE),
00071   _asympInt(kFALSE)
00072 {  
00073 }
00074 
00075 
00076 
00077 //_____________________________________________________________________________
00078 RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn, 
00079                            RooAbsReal& _sigma, RooAbsReal& _rlife, 
00080                            RooAbsReal& _sigmaSF, RooAbsReal& _rlifeSF,
00081                            Bool_t nlo, Type type) : 
00082   RooResolutionModel(name,title,xIn), 
00083   sigma("sigma","Width",this,_sigma),
00084   rlife("rlife","Life time",this,_rlife),
00085   ssf("ssf","Sigma Scale Factor",this,_sigmaSF),
00086   rsf("rsf","RLife Scale Factor",this,_rlifeSF),
00087   _flip(type==Flipped),
00088   _nlo(nlo), 
00089   _flatSFInt(kFALSE),
00090   _asympInt(kFALSE)
00091 {  
00092 }
00093 
00094 
00095 
00096 //_____________________________________________________________________________
00097 RooGExpModel::RooGExpModel(const RooGExpModel& other, const char* name) : 
00098   RooResolutionModel(other,name),
00099   sigma("sigma",this,other.sigma),
00100   rlife("rlife",this,other.rlife),
00101   ssf("ssf",this,other.ssf),
00102   rsf("rsf",this,other.rsf),
00103   _flip(other._flip),
00104   _nlo(other._nlo),
00105   _flatSFInt(other._flatSFInt),
00106   _asympInt(other._asympInt)
00107 {
00108 }
00109 
00110 
00111 
00112 //_____________________________________________________________________________
00113 RooGExpModel::~RooGExpModel()
00114 {
00115   // Destructor
00116 }
00117 
00118 
00119 
00120 //_____________________________________________________________________________
00121 Int_t RooGExpModel::basisCode(const char* name) const 
00122 {
00123   if (!TString("exp(-@0/@1)").CompareTo(name)) return expBasisPlus ;
00124   if (!TString("exp(@0/@1)").CompareTo(name)) return expBasisMinus ;
00125   if (!TString("exp(-abs(@0)/@1)").CompareTo(name)) return expBasisSum ;
00126   if (!TString("exp(-@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisPlus ;
00127   if (!TString("exp(@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisMinus ;
00128   if (!TString("exp(-abs(@0)/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisSum ;
00129   if (!TString("exp(-@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisPlus ;
00130   if (!TString("exp(@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisMinus ;
00131   if (!TString("exp(-abs(@0)/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisSum ;
00132   if (!TString("exp(-@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisPlus;
00133   if (!TString("exp(@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisMinus;
00134   if (!TString("exp(-abs(@0)/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisSum;
00135   if (!TString("exp(-@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisPlus;
00136   if (!TString("exp(@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisMinus;
00137   if (!TString("exp(-abs(@0)/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisSum;
00138   return 0 ;
00139 } 
00140 
00141 
00142 
00143 //_____________________________________________________________________________
00144 Double_t RooGExpModel::evaluate() const 
00145 {  
00146   static Double_t root2(sqrt(2.)) ;
00147 //   static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
00148 //   static Double_t rootpi(sqrt(atan2(0.,-1.)));
00149 
00150   BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
00151   BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
00152 
00153   Double_t fsign = _flip?-1:1 ;
00154 
00155   Double_t sig = sigma*ssf ;
00156   Double_t rtau = rlife*rsf ;
00157  
00158   Double_t tau = (_basisCode!=noBasis)?((RooAbsReal*)basis().getParameter(1))->getVal():0. ;
00159   // added, FMV 07/27/03
00160   if (basisType == coshBasis && _basisCode!=noBasis ) {
00161      Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
00162      if (dGamma==0) basisType = expBasis;
00163   }
00164 
00165   // *** 1st form: Straight GExp, used for unconvoluted PDF or expBasis with 0 lifetime ***
00166   if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
00167     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 1st form" << endl ;    
00168 
00169     Double_t expArg = sig*sig/(2*rtau*rtau) + fsign*x/rtau ;
00170 
00171     Double_t result ;
00172     if (expArg<300) {
00173       result = 1/(2*rtau) * exp(expArg) * RooMath::erfc(sig/(root2*rtau) + fsign*x/(root2*sig));
00174     } else {
00175       // If exponent argument is very large, bring canceling RooMath::erfc() term inside exponent
00176       // to avoid floating point over/underflows of intermediate calculations
00177       result = 1/(2*rtau) * exp(expArg + logErfC(sig/(root2*rtau) + fsign*x/(root2*sig))) ;
00178     }
00179 
00180 //     Double_t result = 1/(2*rtau)
00181 //                     * exp(sig*sig/(2*rtau*rtau) + fsign*x/rtau)
00182 //                     * RooMath::erfc(sig/(root2*rtau) + fsign*x/(root2*sig));
00183 
00184     // equivalent form, added FMV, 07/24/03
00185     //Double_t xprime = x/rtau ;
00186     //Double_t c = sig/(root2*rtau) ;
00187     //Double_t u = xprime/(2*c) ;
00188     //Double_t result = 0.5*evalCerfRe(fsign*u,c) ;  // sign=-1 ! 
00189 
00190     if (_basisCode!=0 && basisSign==Both) result *= 2 ;
00191     //cout << "1st form " << "x= " << x << " result= " << result << endl;
00192     return result ;    
00193   }
00194   
00195   // *** 2nd form: 0, used for sinBasis and cosBasis with tau=0 ***
00196   if (tau==0) {
00197     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 2nd form" << endl ;
00198     return 0. ;
00199   }
00200 
00201   Double_t omega = (basisType!=expBasis)?((RooAbsReal*)basis().getParameter(2))->getVal():0. ;
00202 
00203   // *** 3nd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
00204   if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
00205     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 3d form tau=" << tau << endl ;
00206     Double_t result(0) ;
00207     if (basisSign!=Minus) result += calcDecayConv(+1,tau,sig,rtau,fsign) ;  // modified FMV,08/13/03
00208     if (basisSign!=Plus)  result += calcDecayConv(-1,tau,sig,rtau,fsign) ;  // modified FMV,08/13/03
00209     //cout << "3rd form " << "x= " << x << " result= " << result << endl;
00210     return result ;
00211   }
00212   
00213   // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
00214   Double_t wt = omega *tau ;
00215   if (basisType==sinBasis) {
00216     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 4th form omega = " 
00217                              << omega << ", tau = " << tau << endl ;
00218     Double_t result(0) ;
00219     if (wt==0.) return result ;
00220     if (basisSign!=Minus) result += -1*calcSinConv(+1,sig,tau,omega,rtau,fsign).im() ;
00221     if (basisSign!=Plus) result += -1*calcSinConv(-1,sig,tau,omega,rtau,fsign).im() ;
00222     //cout << "4th form " << "x= " << x << " result= " << result << endl;
00223     return result ;
00224   }
00225 
00226   // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
00227   if (basisType==cosBasis) {
00228     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() 
00229                              << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
00230     Double_t result(0) ;
00231     if (basisSign!=Minus) result += calcSinConv(+1,sig,tau,omega,rtau,fsign).re() ;
00232     if (basisSign!=Plus) result += calcSinConv(-1,sig,tau,omega,rtau,fsign).re() ;
00233     //cout << "5th form " << "x= " << x << " result= " << result << endl;
00234     return result ;  
00235   }
00236 
00237 
00238   // *** 6th form: Convolution with exp(-t/tau)*sinh(dgamma*t/2), used for sinhBasis ***
00239   if (basisType==sinhBasis) {
00240     Double_t dgamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
00241    
00242     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
00243                              << ") 6th form = " << dgamma << ", tau = " << tau << endl;
00244     Double_t result(0);
00245     //if (basisSign!=Minus) result += calcSinhConv(+1,+1,-1,tau,dgamma,sig,rtau,fsign);
00246     //if (basisSign!=Plus) result += calcSinhConv(-1,-1,+1,tau,dgamma,sig,rtau,fsign);
00247     // better form, since it also accounts for the numerical divergence region, added FMV, 07/24/03
00248     Double_t tau1 = 1/(1/tau-dgamma/2) ; 
00249     Double_t tau2 = 1/(1/tau+dgamma/2) ;
00250     if (basisSign!=Minus) result += 0.5*(calcDecayConv(+1,tau1,sig,rtau,fsign)-calcDecayConv(+1,tau2,sig,rtau,fsign));  
00251           // modified FMV,08/13/03
00252     if (basisSign!=Plus) result += 0.5*(calcDecayConv(-1,tau2,sig,rtau,fsign)-calcDecayConv(-1,tau1,sig,rtau,fsign));
00253           // modified FMV,08/13/03
00254     //cout << "6th form " << "x= " << x << " result= " << result << endl;
00255     return result;
00256   }
00257 
00258   // *** 7th form: Convolution with exp(-t/tau)*cosh(dgamma*t/2), used for coshBasis ***
00259   if (basisType==coshBasis) {
00260     Double_t dgamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
00261     
00262     if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
00263                          << ") 7th form = " << dgamma << ", tau = " << tau << endl;
00264     Double_t result(0);
00265     //if (basisSign!=Minus) result += calcCoshConv(+1,tau,dgamma,sig,rtau,fsign);
00266     //if (basisSign!=Plus) result += calcCoshConv(-1,tau,dgamma,sig,rtau,fsign);
00267     // better form, since it also accounts for the numerical divergence region, added FMV, 07/24/03
00268     Double_t tau1 = 1/(1/tau-dgamma/2) ; 
00269     Double_t tau2 = 1/(1/tau+dgamma/2) ;
00270     if (basisSign!=Minus) result += 0.5*(calcDecayConv(+1,tau1,sig,rtau,fsign)+calcDecayConv(+1,tau2,sig,rtau,fsign));
00271           // modified FMV,08/13/03
00272     if (basisSign!=Plus) result += 0.5*(calcDecayConv(-1,tau1,sig,rtau,fsign)+calcDecayConv(-1,tau2,sig,rtau,fsign));
00273           // modified FMV,08/13/03
00274     //cout << "7th form " << "x= " << x << " result= " << result << endl;
00275     return result;
00276   }
00277   assert(0) ;
00278   return 0 ;
00279   }
00280 
00281 
00282 
00283 //_____________________________________________________________________________
00284 Double_t RooGExpModel::logErfC(Double_t xx) const
00285 {
00286   // Approximation of the log of the complex error function
00287   Double_t t,z,ans;
00288   z=fabs(xx);
00289   t=1.0/(1.0+0.5*z);
00290   
00291   if(xx >= 0.0)
00292     ans=log(t)+(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+t*(-0.18628806+
00293         t*(0.27886807+t*(-1.13520398+t*(1.48851587+t*(-0.82215223+t*0.17087277)))))))));
00294   else
00295     ans=log(2.0-t*exp(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+t*(-0.18628806+
00296         t*(0.27886807+t*(-1.13520398+t*(1.48851587+t*(-0.82215223+t*0.17087277))))))))));
00297   
00298   return ans;
00299 }
00300 
00301 
00302 
00303 //_____________________________________________________________________________
00304 RooComplex RooGExpModel::calcSinConv(Double_t sign, Double_t sig, Double_t tau, Double_t omega, Double_t rtau, Double_t fsign) const
00305 {
00306   static Double_t root2(sqrt(2.)) ;
00307 
00308   Double_t s1= -sign*x/tau;
00309   //Double_t s1= x/tau;
00310   Double_t c1= sig/(root2*tau);
00311   Double_t u1= s1/(2*c1);  
00312   Double_t s2= x/rtau;
00313   Double_t c2= sig/(root2*rtau);
00314   Double_t u2= fsign*s2/(2*c2) ;
00315   //Double_t u2= s2/(2*c2) ;
00316 
00317   RooComplex eins(1,0);
00318   RooComplex k(1/tau,sign*omega);  
00319   //return (evalCerf(-sign*omega*tau,u1,c1)+evalCerf(0,u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
00320 
00321   return (evalCerf(-sign*omega*tau,u1,c1)+RooComplex(evalCerfRe(u2,c2),0)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
00322   // equivalent form, added FMV, 07/24/03
00323   //return (evalCerf(-sign*omega*tau,-sign*u1,c1)+evalCerf(0,fsign*u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
00324 }
00325 
00326 
00327 // added FMV,08/18/03
00328 
00329 //_____________________________________________________________________________
00330 Double_t RooGExpModel::calcSinConv(Double_t sign, Double_t sig, Double_t tau, Double_t rtau, Double_t fsign) const
00331 {
00332   static Double_t root2(sqrt(2.)) ;
00333 
00334   Double_t s1= -sign*x/tau;
00335   //Double_t s1= x/tau;
00336   Double_t c1= sig/(root2*tau);
00337   Double_t u1= s1/(2*c1);  
00338   Double_t s2= x/rtau;
00339   Double_t c2= sig/(root2*rtau);
00340   Double_t u2= fsign*s2/(2*c2) ;
00341   //Double_t u2= s2/(2*c2) ;
00342 
00343   Double_t eins(1);
00344   Double_t k(1/tau);  
00345   return (evalCerfRe(u1,c1)+evalCerfRe(u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
00346   // equivalent form, added FMV, 07/24/03
00347   //return (evalCerfRe(-sign*u1,c1)+evalCerfRe(fsign*u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
00348 }
00349 
00350 
00351 
00352 //_____________________________________________________________________________
00353 Double_t RooGExpModel::calcDecayConv(Double_t sign, Double_t tau, Double_t sig, Double_t rtau, Double_t fsign) const
00354 // modified FMV,08/13/03
00355 {
00356   static Double_t root2(sqrt(2.)) ;
00357   static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
00358   static Double_t rootpi(sqrt(atan2(0.,-1.)));
00359 
00360   // Process flip status
00361   Double_t xp(x) ;
00362   //if (_flip) {
00363   //  xp   *= -1 ;
00364   //  sign *= -1 ;
00365   //}
00366   xp *= fsign ;    // modified FMV,08/13/03
00367   sign *= fsign ;  // modified FMV,08/13/03
00368 
00369   Double_t cFly;
00370   if ((sign<0)&&(fabs(tau-rtau)<tau/260)) {
00371 
00372     Double_t MeanTau=0.5*(tau+rtau);
00373     if (fabs(xp/MeanTau)>300) {
00374       return 0 ;
00375     }
00376 
00377     cFly=1./(MeanTau*MeanTau*root2pi) *
00378       exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
00379       *(sig*exp(-1/(2*sig*sig)*TMath::Power((sig*sig/MeanTau+xp),2)) 
00380         -(sig*sig/MeanTau+xp)*(rootpi/root2)*RooMath::erfc(sig/(root2*MeanTau)+xp/(root2*sig)));
00381     
00382     if(_nlo) {
00383       Double_t epsilon=0.5*(tau-rtau);
00384       Double_t a=sig/(root2*MeanTau)+xp/(root2*sig);
00385       cFly += 1./(MeanTau*MeanTau)
00386         *exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
00387         *0.5/MeanTau*epsilon*epsilon*
00388         (exp(-a*a)*(sig/MeanTau*root2/rootpi
00389                     -(4*a*sig*sig)/(2*rootpi*MeanTau*MeanTau)
00390                     +(-4/rootpi+8*a*a/rootpi)/6
00391                     *TMath::Power(sig/(root2*MeanTau),3)
00392                     +2/rootpi*(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
00393                     (sig/(root2*MeanTau)-a*(sig*sig)/(2*MeanTau*MeanTau))
00394                     +2/rootpi*((3*sig*sig)/(2*MeanTau*MeanTau)+xp/MeanTau+
00395                                0.5*TMath::Power(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,2))*sig/(root2*MeanTau))
00396          -(2*sig*sig/(MeanTau*MeanTau)+xp/MeanTau+(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
00397            (3*sig*sig/(2*MeanTau*MeanTau)+xp/MeanTau)
00398            +TMath::Power(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,3)/6)*RooMath::erfc(a));
00399     }
00400     
00401   } else {
00402 
00403     Double_t expArg1 = sig*sig/(2*tau*tau)-sign*xp/tau ;
00404     Double_t expArg2 = sig*sig/(2*rtau*rtau)+xp/rtau ;
00405 
00406     Double_t term1, term2 ;
00407     if (expArg1<300) {
00408       term1 = exp(expArg1) *RooMath::erfc(sig/(root2*tau)-sign*xp/(root2*sig)) ;
00409     } else {
00410       term1 = exp(expArg1+logErfC(sig/(root2*tau)-sign*xp/(root2*sig))) ; ;
00411     }
00412     if (expArg2<300) {
00413       term2 = exp(expArg2) *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)) ;
00414     } else {
00415       term2 = exp(expArg2+logErfC(sig/(root2*rtau)+xp/(root2*sig))) ; ;
00416     }
00417 
00418     cFly=(term1+sign*term2)/(2*(tau+sign*rtau));
00419     
00420     // WVE prevent numeric underflows 
00421     if (cFly<1e-100) {
00422       cFly = 0 ;
00423     }
00424 
00425     // equivalent form, added FMV, 07/24/03 
00426     //cFly = calcSinConv(sign, sig, tau, rtau, fsign)/(2*tau) ;
00427   }
00428 
00429   return cFly*2*tau ;    
00430 }
00431 
00432 /* commented FMV, 07/24/03
00433 
00434 //_____________________________________________________________________________
00435 Double_t RooGExpModel::calcCoshConv(Double_t sign, Double_t tau, Double_t dgamma, Double_t sig, Double_t rtau, Double_t fsign) const
00436 {
00437   
00438   
00439   
00440   static Double_t root2(sqrt(2.)) ;
00441   static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
00442   static Double_t rootpi(sqrt(atan2(0.,-1.)));
00443   Double_t tau1 = 1/(1/tau-dgamma/2);
00444   Double_t tau2 = 1/(1/tau+dgamma/2);
00445   Double_t cFly;
00446   Double_t xp(x);
00447 
00448   //if (_flip) {
00449   //  xp   *= -1 ;
00450   //  sign *= -1 ;
00451   //}
00452   xp *= fsign ;    // modified FMV,08/13/03
00453   sign *= fsign ;  // modified FMV,08/13/03
00454 
00455   cFly=tau1*(exp(sig*sig/(2*tau1*tau1)-sign*xp/tau1)
00456           *RooMath::erfc(sig/(root2*tau1)-sign*xp/(root2*sig))
00457           +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
00458           *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau1+sign*rtau))
00459     +tau2*(exp(sig*sig/(2*tau2*tau2)-sign*xp/tau2)
00460           *RooMath::erfc(sig/(root2*tau2)-sign*xp/(root2*sig))
00461           +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
00462           *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau2+sign*rtau));;
00463   return cFly;
00464 }
00465 */
00466 
00467 /* commented FMV, 07/24/03
00468 
00469 //_____________________________________________________________________________
00470 Double_t RooGExpModel::calcSinhConv(Double_t sign, Double_t sign1, Double_t sign2, Double_t tau, Double_t dgamma, Double_t sig, Double_t rtau, Double_t fsign) const
00471 {
00472   static Double_t root2(sqrt(2.)) ;
00473   static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
00474   static Double_t rootpi(sqrt(atan2(0.,-1.)));
00475   Double_t tau1 = 1/(1/tau-dgamma/2);
00476   Double_t tau2 = 1/(1/tau+dgamma/2);
00477   Double_t cFly;
00478   Double_t xp(x);
00479   
00480   //if (_flip) {
00481   //  xp   *= -1 ;
00482   //  sign1 *= -1 ;
00483   //  sign2 *= -1 ;
00484   //}
00485   xp *= fsign ;    // modified FMV,08/13/03
00486   sign1 *= fsign ;  // modified FMV,08/13/03
00487   sign2 *= fsign ;  // modified FMV,08/13/03
00488 
00489   cFly=sign1*tau1*(exp(sig*sig/(2*tau1*tau1)-sign*xp/tau1)
00490           *RooMath::erfc(sig/(root2*tau1)-sign*xp/(root2*sig))
00491           +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
00492           *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau1+sign*rtau))
00493     +sign2*tau2*(exp(sig*sig/(2*tau2*tau2)-sign*xp/tau2)
00494           *RooMath::erfc(sig/(root2*tau2)-sign*xp/(root2*sig))
00495           +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
00496           *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau2+sign*rtau));;
00497   return cFly;
00498 }
00499 */
00500 
00501 
00502 //_____________________________________________________________________________
00503 Int_t RooGExpModel::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const 
00504 {
00505   switch(_basisCode) {
00506 
00507   // Analytical integration capability of raw PDF
00508   case noBasis:
00509     if (matchArgs(allVars,analVars,convVar())) return 1 ;
00510     break ;
00511 
00512   // Analytical integration capability of convoluted PDF
00513   case expBasisPlus:
00514   case expBasisMinus:
00515   case expBasisSum:
00516   case sinBasisPlus:
00517   case sinBasisMinus:
00518   case sinBasisSum:
00519   case cosBasisPlus:
00520   case cosBasisMinus:
00521   case cosBasisSum:
00522   case sinhBasisPlus:
00523   case sinhBasisMinus:
00524   case sinhBasisSum:
00525   case coshBasisPlus:
00526   case coshBasisMinus:
00527   case coshBasisSum:
00528     
00529     // Optionally advertise flat integral over sigma scale factor
00530     if (_flatSFInt) {
00531       if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) {
00532         return 2 ;
00533       }
00534     }
00535 
00536     if (matchArgs(allVars,analVars,convVar())) return 1 ;
00537     break ;
00538   }
00539   
00540   return 0 ;
00541 }
00542 
00543 
00544 
00545 //_____________________________________________________________________________
00546 Double_t RooGExpModel::analyticalIntegral(Int_t code, const char* rangeName) const 
00547 {
00548   static Double_t root2 = sqrt(2.) ;
00549 //   static Double_t rootPiBy2 = sqrt(atan2(0.0,-1.0)/2.0);
00550   Double_t ssfInt(1.0) ;
00551 
00552   // Code must be 1 or 2
00553   assert(code==1||code==2) ;
00554   if (code==2) {
00555     ssfInt = (ssf.max(rangeName)-ssf.min(rangeName)) ;
00556   }
00557 
00558   BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
00559   BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
00560    
00561   Double_t tau = (_basisCode!=noBasis)?((RooAbsReal*)basis().getParameter(1))->getVal():0 ;
00562 
00563   // added FMV, 07/24/03
00564   if (basisType == coshBasis && _basisCode!=noBasis ) {
00565      Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
00566      if (dGamma==0) basisType = expBasis;
00567   }
00568   Double_t fsign = _flip?-1:1 ;
00569   Double_t sig = sigma*ssf ;
00570   Double_t rtau = rlife*rsf ;
00571 
00572   // *** 1st form????
00573   if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
00574     if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 1st form" << endl ;
00575 
00576     //Double_t result = 1.0 ; // WVE inferred from limit(tau->0) of cosBasisNorm
00577     // finite+asymtotic normalization, added FMV, 07/24/03
00578     Double_t xpmin = x.min(rangeName)/rtau ;
00579     Double_t xpmax = x.max(rangeName)/rtau ;
00580     Double_t c = sig/(root2*rtau) ;
00581     Double_t umin = xpmin/(2*c) ;
00582     Double_t umax = xpmax/(2*c) ;
00583     Double_t result ;
00584     if (_asympInt) {
00585       result = 1.0 ;
00586     } else {
00587       result = 0.5*evalCerfInt(-fsign,rtau,fsign*umin,fsign*umax,c)/rtau ;  //WVEFIX add 1/rtau
00588     }
00589 
00590     if (_basisCode!=0 && basisSign==Both) result *= 2 ;    
00591     //cout << "Integral 1st form " << " result= " << result*ssfInt << endl;
00592     return result*ssfInt ;    
00593   }
00594 
00595   Double_t omega = (basisType!=expBasis) ?((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
00596 
00597   // *** 2nd form: unity, used for sinBasis and cosBasis with tau=0 (PDF is zero) ***
00598   //if (tau==0&&omega!=0) {
00599   if (tau==0) {  // modified, FMV 07/24/03
00600     if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 2nd form" << endl ;
00601     return 0. ;
00602   }
00603 
00604   // *** 3rd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
00605   if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
00606     //Double_t result = 2*tau ;
00607     //if (basisSign==Both) result *= 2 ;
00608     // finite+asymtotic normalization, added FMV, 07/24/03
00609     Double_t result(0.);
00610     if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,sig,rtau,fsign,rangeName); 
00611     if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,sig,rtau,fsign,rangeName);  
00612     //cout << "Integral 3rd form " << " result= " << result*ssfInt << endl;
00613     return result*ssfInt ;
00614   }
00615   
00616   // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
00617   Double_t wt = omega * tau ;    
00618   if (basisType==sinBasis) {    
00619     if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 4th form omega = " 
00620                              << omega << ", tau = " << tau << endl ;
00621     //cout << "sin integral" << endl;
00622     Double_t result(0) ;
00623     if (wt==0) return result ;
00624     //if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega).im() ;
00625     //if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega).im() ;
00626     // finite+asymtotic normalization, added FMV, 07/24/03
00627     if (basisSign!=Minus) result += -1*calcSinConvNorm(+1,tau,omega,sig,rtau,fsign,rangeName).im();
00628     if (basisSign!=Plus) result += -1*calcSinConvNorm(-1,tau,omega,sig,rtau,fsign,rangeName).im();
00629     //cout << "Integral 4th form " << " result= " << result*ssfInt << endl;
00630     return result*ssfInt ;
00631   }
00632  
00633   // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
00634   if (basisType==cosBasis) {
00635     if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() 
00636                              << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
00637     //cout << "cos integral" << endl;
00638     Double_t result(0) ;
00639     //if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega).re() ;
00640     //if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega).re() ;
00641     // finite+asymtotic normalization, added FMV, 07/24/03
00642     if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega,sig,rtau,fsign,rangeName).re();
00643     if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega,sig,rtau,fsign,rangeName).re();
00644     //cout << "Integral 5th form " << " result= " << result*ssfInt << endl;
00645     return result*ssfInt ;
00646   }
00647   
00648   Double_t dgamma = ((basisType==coshBasis)||(basisType==sinhBasis))?((RooAbsReal*)basis().getParameter(2))->getVal():0 ;  
00649  
00650   // *** 6th form: Convolution with exp(-t/tau)*sinh(dgamma*t/2), used for sinhBasis ***
00651   if (basisType==sinhBasis) {
00652     if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() 
00653                              << ") 6th form dgamma = " << dgamma << ", tau = " << tau << endl ;
00654     Double_t tau1 = 1/(1/tau-dgamma/2);
00655     Double_t tau2 = 1/(1/tau+dgamma/2);
00656     //cout << "sinh integral" << endl;
00657     Double_t result(0) ;
00658     //if (basisSign!=Minus) result += tau1-tau2 ;
00659     //if (basisSign!=Plus) result += tau2-tau1 ;
00660     // finite+asymtotic normalization, added FMV, 07/24/03
00661     if (basisSign!=Minus) result += 0.5*(calcSinConvNorm(+1,tau1,sig,rtau,fsign,rangeName)-
00662                                          calcSinConvNorm(+1,tau2,sig,rtau,fsign,rangeName));
00663     if (basisSign!=Plus) result += 0.5*(calcSinConvNorm(-1,tau2,sig,rtau,fsign,rangeName)-
00664                                         calcSinConvNorm(-1,tau1,sig,rtau,fsign,rangeName));
00665     //cout << "Integral 6th form " << " result= " << result*ssfInt << endl;
00666     return result;
00667     }
00668 
00669   // ** 7th form: Convolution with exp(-t/tau)*cosh(dgamma*t/2), used for coshBasis ***
00670   if (basisType==coshBasis) {
00671     if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() 
00672                              << ") 6th form dgamma = " << dgamma << ", tau = " << tau << endl ;
00673     //cout << "cosh integral" << endl;
00674     Double_t tau1 = 1/(1/tau-dgamma/2);
00675     Double_t tau2 = 1/(1/tau+dgamma/2);
00676     //Double_t result = (tau1+tau2) ;
00677     //if (basisSign==Both) result *= 2 ;
00678     // finite+asymtotic normalization, added FMV, 07/24/03
00679     Double_t result(0);
00680     if (basisSign!=Minus) result += 0.5*(calcSinConvNorm(+1,tau1,sig,rtau,fsign,rangeName)+
00681                                          calcSinConvNorm(+1,tau2,sig,rtau,fsign,rangeName));
00682     if (basisSign!=Plus) result += 0.5*(calcSinConvNorm(-1,tau1,sig,rtau,fsign,rangeName)+
00683                                         calcSinConvNorm(-1,tau2,sig,rtau,fsign,rangeName));
00684     //cout << "Integral 7th form " << " result= " << result*ssfInt << endl;
00685     return result;
00686   
00687     }
00688 
00689   assert(0) ;
00690   return 1 ;
00691 }
00692 
00693 
00694 // modified FMV, 07/24/03. Finite+asymtotic normalization
00695 
00696 //_____________________________________________________________________________
00697 RooComplex RooGExpModel::calcSinConvNorm(Double_t sign, Double_t tau, Double_t omega, 
00698                                          Double_t sig, Double_t rtau, Double_t fsign, const char* rangeName) const
00699 {
00700   //  old code (asymptotic normalization only)
00701   //  RooComplex z(1/tau,sign*omega);
00702   //  return z*2/(omega*omega+1/(tau*tau));
00703 
00704   static Double_t root2(sqrt(2.)) ;
00705 
00706   Double_t smin1= x.min(rangeName)/tau;
00707   Double_t smax1= x.max(rangeName)/tau;
00708   Double_t c1= sig/(root2*tau);
00709   Double_t umin1= smin1/(2*c1);  
00710   Double_t umax1= smax1/(2*c1);  
00711   Double_t smin2= x.min(rangeName)/rtau;
00712   Double_t smax2= x.max(rangeName)/rtau;
00713   Double_t c2= sig/(root2*rtau);
00714   Double_t umin2= smin2/(2*c2) ;
00715   Double_t umax2= smax2/(2*c2) ;
00716 
00717   RooComplex eins(1,0);
00718   RooComplex k(1/tau,sign*omega);
00719   RooComplex term1 = evalCerfInt(sign,-sign*omega*tau, tau, -sign*umin1, -sign*umax1, c1);
00720   //RooComplex term2 = evalCerfInt(-fsign,0., rtau, fsign*umin2, fsign*umax2, c2)*RooComplex(fsign*sign,0);
00721   RooComplex term2 = RooComplex(evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2, c2)*fsign*sign,0);
00722   return (term1+term2)/(eins + k*fsign*sign*rtau) ;
00723 }
00724 
00725 
00726 // added FMV, 08/17/03
00727 
00728 //_____________________________________________________________________________
00729 Double_t RooGExpModel::calcSinConvNorm(Double_t sign, Double_t tau, Double_t sig, Double_t rtau, Double_t fsign, const char* rangeName) const
00730 {
00731   static Double_t root2(sqrt(2.)) ;
00732 
00733   Double_t smin1= x.min(rangeName)/tau;
00734   Double_t smax1= x.max(rangeName)/tau;
00735   Double_t c1= sig/(root2*tau);
00736   Double_t umin1= smin1/(2*c1);  
00737   Double_t umax1= smax1/(2*c1);  
00738   Double_t smin2= x.min(rangeName)/rtau;
00739   Double_t smax2= x.max(rangeName)/rtau;
00740   Double_t c2= sig/(root2*rtau);
00741   Double_t umin2= smin2/(2*c2) ;
00742   Double_t umax2= smax2/(2*c2) ;
00743 
00744   Double_t eins(1);
00745   Double_t k(1/tau);
00746   Double_t term1 = evalCerfInt(sign, tau, -sign*umin1, -sign*umax1, c1);
00747   Double_t term2 = evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2, c2)*fsign*sign;
00748 
00749   // WVE Handle 0/0 numeric divergence 
00750   if (fabs(tau-rtau)<1e-10 && fabs(term1+term2)<1e-10) {
00751     cout << "epsilon method" << endl ;
00752     static Double_t epsilon = 1e-4 ;
00753     return calcSinConvNorm(sign,tau+epsilon,sig,rtau-epsilon,fsign,rangeName) ;
00754   }
00755   return (term1+term2)/(eins + k*fsign*sign*rtau) ;
00756 }
00757 
00758 
00759 
00760 // added FMV, 07/24/03
00761 //_____________________________________________________________________________
00762 RooComplex RooGExpModel::evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
00763 {
00764   RooComplex diff;
00765   if (_asympInt) {
00766     diff = RooComplex(2,0) ;
00767   } else {
00768     diff = RooComplex(sign,0.)*(evalCerf(wt,umin,c) - evalCerf(wt,umax,c) + RooMath::erf(umin) - RooMath::erf(umax));
00769   }
00770   return RooComplex(tau/(1.+wt*wt),0)*RooComplex(1,wt)*diff;
00771 }
00772 // added FMV, 08/17/03. Modified FMV, 08/30/03
00773 
00774 //_____________________________________________________________________________
00775 Double_t RooGExpModel::evalCerfInt(Double_t sign, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
00776 {
00777   Double_t diff;
00778   if (_asympInt) {
00779     diff = 2. ;
00780   } else {
00781     if ((umin<-8 && umax>8)||(umax<-8 && umin>8)) {
00782       // If integral is over >8 sigma, approximate with full integral
00783       diff = 2. ;
00784     } else {
00785       diff = sign*(evalCerfRe(umin,c) - evalCerfRe(umax,c) + RooMath::erf(umin) - RooMath::erf(umax));
00786     }
00787   }
00788   return tau*diff;
00789 }
00790 
00791 
00792 
00793 //_____________________________________________________________________________
00794 RooComplex RooGExpModel::evalCerfApprox(Double_t swt, Double_t u, Double_t c) const
00795 {
00796   // use the approximation: erf(z) = exp(-z*z)/(sqrt(pi)*z)
00797   // to explicitly cancel the divergent exp(y*y) behaviour of
00798   // CWERF for z = x + i y with large negative y
00799 
00800   static Double_t rootpi= sqrt(atan2(0.,-1.));
00801   RooComplex z(swt*c,u+c);  
00802   RooComplex zc(u+c,-swt*c);
00803   RooComplex zsq= z*z;
00804   RooComplex v= -zsq - u*u;
00805 
00806   return v.exp()*(-zsq.exp()/(zc*rootpi) + 1)*2 ;
00807 }
00808 
00809 
00810 
00811 //_____________________________________________________________________________
00812 Int_t RooGExpModel::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
00813 {
00814   if (matchArgs(directVars,generateVars,x)) return 1 ; 
00815   return 0 ;
00816 }
00817 
00818 
00819 
00820 //_____________________________________________________________________________
00821 void RooGExpModel::generateEvent(Int_t code)
00822 {
00823   assert(code==1) ;
00824   Double_t xgen ;
00825   while(1) {
00826     Double_t xgau = RooRandom::randomGenerator()->Gaus(0,(sigma*ssf));
00827     Double_t xexp = RooRandom::uniform();
00828     if (!_flip) xgen= xgau + (rlife*rsf)*log(xexp);  // modified, FMV 08/13/03
00829     else xgen= xgau - (rlife*rsf)*log(xexp);
00830     if (xgen<x.max() && xgen>x.min()) {
00831       x = xgen ;
00832       return ;
00833     }
00834   }
00835 }
00836 
00837 
00838 
00839 
00840 
00841 
00842 
00843 
00844 
00845 
00846 
00847 
00848 
00849 
00850 
00851 
00852 
00853 
00854 
00855 
00856 
00857 
00858 
00859 
00860 
00861 

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