RooBDecay.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooBDecay.cxx 24286 2008-06-16 15:47:04Z wouter $
00005  * Authors:                                                                  *
00006  *   PL, Parker C Lund,   UC Irvine                                          *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00009  *                                                                           *
00010  * Copyright (c) 2000-2005, Regents of the University of California          *
00011  *                          and Stanford University. All rights reserved.    *
00012  *                                                                           *
00013  * Redistribution and use in source and binary forms,                        *
00014  * with or without modification, are permitted according to the terms        *
00015  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00016  *****************************************************************************/
00017 
00018 
00019 //////////////////////////////////////////////////////////////////////////////
00020 //
00021 // BEGIN_HTML
00022 // Most general description of B decay time distribution with effects
00023 // of CP violation, mixing and life time differences. This function can 
00024 // be analytically convolved with any RooResolutionModel implementation
00025 // END_HTML
00026 //
00027 
00028 
00029 #include "RooFit.h"
00030 
00031 #include "Riostream.h"
00032 #include "TMath.h"
00033 #include "RooBDecay.h"
00034 #include "RooRealVar.h"
00035 #include "RooRandom.h"
00036 
00037 ClassImp(RooBDecay);
00038 
00039 
00040 //_____________________________________________________________________________
00041 RooBDecay::RooBDecay(const char *name, const char* title, 
00042                RooRealVar& t, RooAbsReal& tau, RooAbsReal& dgamma,
00043                RooAbsReal& f0, RooAbsReal& f1, RooAbsReal& f2, RooAbsReal& f3, 
00044                RooAbsReal& dm, const RooResolutionModel& model, DecayType type) :
00045   RooAbsAnaConvPdf(name, title, model, t),
00046   _t("t", "time", this, t),
00047   _tau("tau", "Average Decay Time", this, tau),
00048   _dgamma("dgamma", "Delta Gamma", this, dgamma),
00049   _f0("f0", "Cosh Coefficient", this, f0),
00050   _f1("f1", "Sinh Coefficient", this, f1),
00051   _f2("f2", "Cos Coefficient", this, f2),
00052   _f3("f3", "Sin Coefficient", this, f3),
00053   _dm("dm", "Delta Mass", this, dm),
00054   _type(type)
00055 
00056 {
00057   //Constructor
00058   switch(type)
00059     {
00060     case SingleSided:
00061       _basisCosh = declareBasis("exp(-@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
00062       _basisSinh = declareBasis("exp(-@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
00063       _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
00064       _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
00065       break;
00066     case Flipped:
00067       _basisCosh = declareBasis("exp(@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
00068       _basisSinh = declareBasis("exp(@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
00069       _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
00070       _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
00071       break;
00072     case DoubleSided:
00073       _basisCosh = declareBasis("exp(-abs(@0)/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
00074       _basisSinh = declareBasis("exp(-abs(@0)/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
00075       _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau, dm));
00076       _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau, dm));
00077       break;
00078     }
00079 }
00080 
00081 //_____________________________________________________________________________
00082 RooBDecay::RooBDecay(const RooBDecay& other, const char* name) :
00083   RooAbsAnaConvPdf(other, name),
00084   _t("t", this, other._t),
00085   _tau("tau", this, other._tau),
00086   _dgamma("dgamma", this, other._dgamma),
00087   _f0("f0", this, other._f0),
00088   _f1("f1", this, other._f1),
00089   _f2("f2", this, other._f2),
00090   _f3("f3", this, other._f3),
00091   _dm("dm", this, other._dm),
00092   _basisCosh(other._basisCosh),
00093   _basisSinh(other._basisSinh),
00094   _basisCos(other._basisCos),
00095   _basisSin(other._basisSin),
00096   _type(other._type)
00097 {
00098   //Copy constructor
00099 }
00100 
00101 
00102 
00103 //_____________________________________________________________________________
00104 RooBDecay::~RooBDecay()
00105 {
00106   //Destructor
00107 }
00108 
00109 
00110 //_____________________________________________________________________________
00111 Double_t RooBDecay::coefficient(Int_t basisIndex) const
00112 {
00113   if(basisIndex == _basisCosh)
00114     {  
00115       return _f0;
00116     }
00117   if(basisIndex == _basisSinh)
00118     {
00119       return _f1;
00120     }
00121   if(basisIndex == _basisCos)
00122     {
00123       return _f2;
00124     }
00125   if(basisIndex == _basisSin)
00126     {
00127       return _f3;
00128     }
00129 
00130   return 0 ;
00131 }
00132 
00133 
00134 //_____________________________________________________________________________
00135 RooArgSet* RooBDecay::coefVars(Int_t basisIndex) const 
00136 {
00137   if(basisIndex == _basisCosh)
00138     {  
00139       return _f0.arg().getVariables();
00140     }
00141   if(basisIndex == _basisSinh)
00142     {
00143       return _f1.arg().getVariables();
00144     }
00145   if(basisIndex == _basisCos)
00146     {
00147       return _f2.arg().getVariables();
00148     }
00149   if(basisIndex == _basisSin)
00150     {
00151       return _f3.arg().getVariables();
00152     }
00153 
00154   return 0 ;  
00155 }
00156 
00157 
00158 
00159 //_____________________________________________________________________________
00160 Int_t RooBDecay::getCoefAnalyticalIntegral(Int_t coef, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const 
00161 {
00162   if(coef == _basisCosh)
00163     {  
00164       return _f0.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
00165     }
00166   if(coef == _basisSinh)
00167     {
00168       return _f1.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
00169     }
00170   if(coef == _basisCos)
00171     {
00172       return _f2.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
00173     }
00174   if(coef == _basisSin)
00175     {
00176       return _f3.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
00177     }
00178 
00179   return 0 ;
00180 }
00181 
00182 
00183 
00184 //_____________________________________________________________________________
00185 Double_t RooBDecay::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* rangeName) const 
00186 {
00187   if(coef == _basisCosh)
00188     {  
00189       return _f0.arg().analyticalIntegral(code,rangeName) ;
00190     }
00191   if(coef == _basisSinh)
00192     {
00193       return _f1.arg().analyticalIntegral(code,rangeName) ;
00194     }
00195   if(coef == _basisCos)
00196     {
00197       return _f2.arg().analyticalIntegral(code,rangeName) ;
00198     }
00199   if(coef == _basisSin)
00200     {
00201       return _f3.arg().analyticalIntegral(code,rangeName) ;
00202     }
00203 
00204   return 0 ;
00205 }
00206 
00207 
00208 
00209 //_____________________________________________________________________________
00210 Int_t RooBDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
00211 {
00212   if (matchArgs(directVars, generateVars, _t)) return 1;
00213   return 0;
00214 }
00215 
00216 
00217 
00218 //_____________________________________________________________________________
00219 void RooBDecay::generateEvent(Int_t code)
00220 {
00221   assert(code==1);
00222   Double_t gammamin = 1/_tau-TMath::Abs(_dgamma)/2;
00223   while(1) {
00224     Double_t t = -log(RooRandom::uniform())/gammamin;
00225     if (_type == Flipped || (_type == DoubleSided && RooRandom::uniform() <0.5) ) t *= -1;
00226     if ( t<_t.min() || t>_t.max() ) continue;
00227 
00228     Double_t dgt = _dgamma*t/2;
00229     Double_t dmt = _dm*t;
00230     Double_t ft = fabs(t);
00231     Double_t f = exp(-ft/_tau)*(_f0*cosh(dgt)+_f1*sinh(dgt)+_f2*cos(dmt)+_f3*sin(dmt));
00232     if(f < 0) {
00233       cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: PDF value less than zero" << endl;
00234       ::abort();
00235     }
00236     Double_t w = 1.001*exp(-ft*gammamin)*(TMath::Abs(_f0)+TMath::Abs(_f1)+sqrt(_f2*_f2+_f3*_f3));
00237     if(w < f) {
00238       cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: Envelope function less than p.d.f. " << endl;
00239       cout << _f0 << endl;
00240       cout << _f1 << endl;
00241       cout << _f2 << endl;
00242       cout << _f3 << endl;
00243       ::abort();
00244     }
00245     if(w*RooRandom::uniform() > f) continue;
00246     _t = t;
00247     break;
00248   }
00249 }
00250 
00251 
00252 
00253 
00254 
00255 
00256 
00257 
00258 
00259 
00260 
00261 
00262 
00263 
00264 
00265 
00266 
00267 
00268 
00269 
00270 
00271 
00272 
00273 
00274 
00275 
00276 
00277 
00278 
00279 
00280 
00281 
00282 
00283 
00284 
00285 
00286 
00287 
00288 
00289 
00290 
00291 
00292 
00293 
00294 
00295 
00296 
00297 
00298 
00299 
00300 
00301 

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