RooBMixDecay.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooBMixDecay.cxx 36230 2010-10-09 20:21:02Z 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 RooBMixDecay is a RooAbsAnaConvPdf implementation that describes
00021 // the decay of B mesons with the effects of B0/B0bar mixing. 
00022 // This function can be analytically convolved with any RooResolutionModel implementation
00023 // END_HTML
00024 //
00025 
00026 #include "RooFit.h"
00027 
00028 #include "Riostream.h"
00029 #include "TMath.h"
00030 #include "RooRealVar.h"
00031 #include "RooBMixDecay.h"
00032 #include "RooRealIntegral.h"
00033 #include "RooRandom.h"
00034 
00035 ClassImp(RooBMixDecay) 
00036 ;
00037 
00038 
00039 
00040 //_____________________________________________________________________________
00041 RooBMixDecay::RooBMixDecay(const char *name, const char *title, 
00042                            RooRealVar& t, RooAbsCategory& mixState,
00043                            RooAbsCategory& tagFlav,
00044                            RooAbsReal& tau, RooAbsReal& dm,                        
00045                            RooAbsReal& mistag, RooAbsReal& delMistag,
00046                            const RooResolutionModel& model, 
00047                            DecayType type) :
00048   RooAbsAnaConvPdf(name,title,model,t), 
00049   _type(type),
00050   _mistag("mistag","Mistag rate",this,mistag),
00051   _delMistag("delMistag","Delta mistag rate",this,delMistag),
00052   _mixState("mixState","Mixing state",this,mixState),
00053   _tagFlav("tagFlav","Flavour of tagged B0",this,tagFlav),
00054   _tau("tau","Mixing life time",this,tau),
00055   _dm("dm","Mixing frequency",this,dm),
00056   _t("_t","time",this,t), _genMixFrac(0)
00057 {
00058   // Constructor
00059   switch(type) {
00060   case SingleSided:
00061     _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
00062     _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
00063     break ;
00064   case Flipped:
00065     _basisExp = declareBasis("exp(@0/@1)",RooArgList(tau,dm)) ;
00066     _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
00067     break ;
00068   case DoubleSided:
00069     _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
00070     _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
00071     break ;
00072   }
00073 }
00074 
00075 
00076 
00077 //_____________________________________________________________________________
00078 RooBMixDecay::RooBMixDecay(const RooBMixDecay& other, const char* name) : 
00079   RooAbsAnaConvPdf(other,name), 
00080   _type(other._type),
00081   _mistag("mistag",this,other._mistag),
00082   _delMistag("delMistag",this,other._delMistag),
00083   _mixState("mixState",this,other._mixState),
00084   _tagFlav("tagFlav",this,other._tagFlav),
00085   _tau("tau",this,other._tau),
00086   _dm("dm",this,other._dm),
00087   _t("t",this,other._t),
00088   _basisExp(other._basisExp),
00089   _basisCos(other._basisCos),
00090   _genMixFrac(other._genMixFrac),
00091   _genFlavFrac(other._genFlavFrac),
00092   _genFlavFracMix(other._genFlavFracMix),
00093   _genFlavFracUnmix(other._genFlavFracUnmix)
00094 {
00095   // Copy constructor
00096 }
00097 
00098 
00099 
00100 //_____________________________________________________________________________
00101 RooBMixDecay::~RooBMixDecay()
00102 {
00103   // Destructor
00104 }
00105 
00106 
00107 
00108 //_____________________________________________________________________________
00109 Double_t RooBMixDecay::coefficient(Int_t basisIndex) const 
00110 {
00111   // Comp with tFit MC: must be (1 - tagFlav*...)
00112   if (basisIndex==_basisExp) {
00113     return (1 - _tagFlav*_delMistag) ; 
00114   }
00115 
00116   if (basisIndex==_basisCos) {
00117     return _mixState*(1-2*_mistag) ;   
00118   }
00119   
00120   return 0 ;
00121 }
00122 
00123 
00124 
00125 //_____________________________________________________________________________
00126 Int_t RooBMixDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const 
00127 {
00128 //   cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
00129   if (rangeName) {
00130     return 0 ;
00131   }
00132 
00133   if (matchArgs(allVars,analVars,_mixState,_tagFlav)) return 3 ;
00134   if (matchArgs(allVars,analVars,_mixState)) return 2 ;
00135   if (matchArgs(allVars,analVars,_tagFlav)) return 1 ;
00136   return 0 ;
00137 }
00138 
00139 
00140 
00141 //_____________________________________________________________________________
00142 Double_t RooBMixDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const 
00143 {  
00144   switch(code) {
00145     // No integration
00146   case 0: return coefficient(basisIndex) ;
00147 
00148     // Integration over 'mixState' and 'tagFlav' 
00149   case 3:
00150     if (basisIndex==_basisExp) {
00151       return 4.0 ;
00152     }    
00153     if (basisIndex==_basisCos) {
00154       return 0.0 ;
00155     }
00156     break ;
00157 
00158     // Integration over 'mixState'
00159   case 2:
00160     if (basisIndex==_basisExp) {
00161       return 2.0*coefficient(basisIndex) ;
00162     }    
00163     if (basisIndex==_basisCos) {
00164       return 0.0 ;
00165     }
00166     break ;
00167 
00168     // Integration over 'tagFlav'
00169   case 1:
00170     if (basisIndex==_basisExp) {
00171       return 2.0 ;
00172     }    
00173     if (basisIndex==_basisCos) {
00174       return 2.0*coefficient(basisIndex) ;
00175     }
00176     break ;
00177 
00178   default:
00179     assert(0) ;
00180   }
00181     
00182   return 0 ;
00183 }
00184 
00185 
00186 
00187 //_____________________________________________________________________________
00188 Int_t RooBMixDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
00189 {
00190   if (staticInitOK) {
00191     if (matchArgs(directVars,generateVars,_t,_mixState,_tagFlav)) return 4 ;  
00192     if (matchArgs(directVars,generateVars,_t,_mixState)) return 3 ;  
00193     if (matchArgs(directVars,generateVars,_t,_tagFlav)) return 2 ;  
00194   }
00195 
00196   if (matchArgs(directVars,generateVars,_t)) return 1 ;  
00197   return 0 ;
00198 }
00199 
00200 
00201 
00202 //_____________________________________________________________________________
00203 void RooBMixDecay::initGenerator(Int_t code)
00204 {
00205   switch (code) {
00206   case 2:
00207     {
00208       // Calculate the fraction of B0bar events to generate
00209       Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
00210       _tagFlav = 1 ; // B0 
00211       Double_t flavInt = RooRealIntegral("flavInt","flav integral",*this,RooArgSet(_t.arg())).getVal() ;
00212       _genFlavFrac = flavInt/sumInt ;
00213       break ;
00214     }  
00215   case 3:
00216     {
00217       // Calculate the fraction of mixed events to generate
00218       Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg())).getVal() ;
00219       _mixState = -1 ; // mixed
00220       Double_t mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
00221       _genMixFrac = mixInt/sumInt ;
00222       break ;
00223     }  
00224   case 4:
00225     {
00226       // Calculate the fraction of mixed events to generate
00227       Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg(),_tagFlav.arg())).getVal() ;
00228       _mixState = -1 ; // mixed
00229       Double_t mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
00230       _genMixFrac = mixInt/sumInt ;
00231       
00232       // Calculate the fractio of B0bar tags for mixed and unmixed
00233       RooRealIntegral dtInt("mixInt","mix integral",*this,RooArgSet(_t.arg())) ;
00234       _mixState = -1 ; // Mixed
00235       _tagFlav  =  1 ; // B0
00236       _genFlavFracMix   = dtInt.getVal() / mixInt ;
00237       _mixState =  1 ; // Unmixed
00238       _tagFlav  =  1 ; // B0
00239       _genFlavFracUnmix = dtInt.getVal() / (sumInt - mixInt) ;
00240       break ;
00241     }
00242   }
00243 }
00244 
00245 
00246 
00247 
00248 //_____________________________________________________________________________
00249 void RooBMixDecay::generateEvent(Int_t code)
00250 {
00251   // Generate mix-state dependent
00252   switch(code) {
00253   case 2:
00254     {
00255       Double_t rand = RooRandom::uniform() ;
00256       _tagFlav = (Int_t) ((rand<=_genFlavFrac) ?  1 : -1) ;
00257       break ;
00258     }
00259   case 3:
00260     {
00261       Double_t rand = RooRandom::uniform() ;
00262       _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
00263       break ;
00264     }
00265   case 4:
00266     {
00267       Double_t rand = RooRandom::uniform() ;
00268       _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
00269 
00270       rand = RooRandom::uniform() ;
00271       Double_t genFlavFrac = (_mixState==-1) ? _genFlavFracMix : _genFlavFracUnmix ;
00272       _tagFlav = (Int_t) ((rand<=genFlavFrac) ?  1 : -1) ;
00273       break ;
00274     }
00275   }
00276 
00277   // Generate delta-t dependent
00278   while(1) {
00279     Double_t rand = RooRandom::uniform() ;
00280     Double_t tval(0) ;
00281 
00282     switch(_type) {
00283     case SingleSided:
00284       tval = -_tau*log(rand);
00285       break ;
00286     case Flipped:
00287       tval= +_tau*log(rand);
00288       break ;
00289     case DoubleSided:
00290       tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
00291       break ;
00292     }
00293 
00294     // Accept event if T is in generated range
00295     Double_t dil = 1-2.*_mistag ;
00296     Double_t maxAcceptProb = 1 + TMath::Abs(_delMistag) + TMath::Abs(dil) ;
00297     Double_t acceptProb = (1-_tagFlav*_delMistag) + _mixState*dil*cos(_dm*tval);
00298     Bool_t mixAccept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
00299     
00300     if (tval<_t.max() && tval>_t.min() && mixAccept) {
00301       _t = tval ;
00302       break ;
00303     }
00304   }  
00305 }

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