ROOT  6.06/09
Reference Guide
RooBMixDecay.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes
21 // the decay of B mesons with the effects of B0/B0bar mixing.
22 // This function can be analytically convolved with any RooResolutionModel implementation
23 // END_HTML
24 //
25 
26 #include "RooFit.h"
27 
28 #include "Riostream.h"
29 #include "TMath.h"
30 #include "RooRealVar.h"
31 #include "RooBMixDecay.h"
32 #include "RooRealIntegral.h"
33 #include "RooRandom.h"
34 
35 using namespace std;
36 
38 ;
39 
40 
41 
42 ////////////////////////////////////////////////////////////////////////////////
43 /// Constructor
44 
45 RooBMixDecay::RooBMixDecay(const char *name, const char *title,
46  RooRealVar& t, RooAbsCategory& mixState,
47  RooAbsCategory& tagFlav,
48  RooAbsReal& tau, RooAbsReal& dm,
49  RooAbsReal& mistag, RooAbsReal& delMistag,
50  const RooResolutionModel& model,
51  DecayType type) :
52  RooAbsAnaConvPdf(name,title,model,t),
53  _type(type),
54  _mistag("mistag","Mistag rate",this,mistag),
55  _delMistag("delMistag","Delta mistag rate",this,delMistag),
56  _mixState("mixState","Mixing state",this,mixState),
57  _tagFlav("tagFlav","Flavour of tagged B0",this,tagFlav),
58  _tau("tau","Mixing life time",this,tau),
59  _dm("dm","Mixing frequency",this,dm),
60  _t("_t","time",this,t), _genMixFrac(0)
61 {
62  switch(type) {
63  case SingleSided:
64  _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
65  _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
66  break ;
67  case Flipped:
68  _basisExp = declareBasis("exp(@0/@1)",RooArgList(tau,dm)) ;
69  _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
70  break ;
71  case DoubleSided:
72  _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
73  _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
74  break ;
75  }
76 }
77 
78 
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// Copy constructor
82 
83 RooBMixDecay::RooBMixDecay(const RooBMixDecay& other, const char* name) :
84  RooAbsAnaConvPdf(other,name),
85  _type(other._type),
86  _mistag("mistag",this,other._mistag),
87  _delMistag("delMistag",this,other._delMistag),
88  _mixState("mixState",this,other._mixState),
89  _tagFlav("tagFlav",this,other._tagFlav),
90  _tau("tau",this,other._tau),
91  _dm("dm",this,other._dm),
92  _t("t",this,other._t),
93  _basisExp(other._basisExp),
94  _basisCos(other._basisCos),
95  _genMixFrac(other._genMixFrac),
96  _genFlavFrac(other._genFlavFrac),
97  _genFlavFracMix(other._genFlavFracMix),
98  _genFlavFracUnmix(other._genFlavFracUnmix)
99 {
100 }
101 
102 
103 
104 ////////////////////////////////////////////////////////////////////////////////
105 /// Destructor
106 
108 {
109 }
110 
111 
112 
113 ////////////////////////////////////////////////////////////////////////////////
114 /// Comp with tFit MC: must be (1 - tagFlav*...)
115 
117 {
118  if (basisIndex==_basisExp) {
119  return (1 - _tagFlav*_delMistag) ;
120  }
121 
122  if (basisIndex==_basisCos) {
123  return _mixState*(1-2*_mistag) ;
124  }
125 
126  return 0 ;
127 }
128 
129 
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
133 
134 Int_t RooBMixDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
135 {
136  if (rangeName) {
137  return 0 ;
138  }
139 
140  if (matchArgs(allVars,analVars,_mixState,_tagFlav)) return 3 ;
141  if (matchArgs(allVars,analVars,_mixState)) return 2 ;
142  if (matchArgs(allVars,analVars,_tagFlav)) return 1 ;
143  return 0 ;
144 }
145 
146 
147 
148 ////////////////////////////////////////////////////////////////////////////////
149 
150 Double_t RooBMixDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
151 {
152  switch(code) {
153  // No integration
154  case 0: return coefficient(basisIndex) ;
155 
156  // Integration over 'mixState' and 'tagFlav'
157  case 3:
158  if (basisIndex==_basisExp) {
159  return 4.0 ;
160  }
161  if (basisIndex==_basisCos) {
162  return 0.0 ;
163  }
164  break ;
165 
166  // Integration over 'mixState'
167  case 2:
168  if (basisIndex==_basisExp) {
169  return 2.0*coefficient(basisIndex) ;
170  }
171  if (basisIndex==_basisCos) {
172  return 0.0 ;
173  }
174  break ;
175 
176  // Integration over 'tagFlav'
177  case 1:
178  if (basisIndex==_basisExp) {
179  return 2.0 ;
180  }
181  if (basisIndex==_basisCos) {
182  return 2.0*coefficient(basisIndex) ;
183  }
184  break ;
185 
186  default:
187  assert(0) ;
188  }
189 
190  return 0 ;
191 }
192 
193 
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 
197 Int_t RooBMixDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
198 {
199  if (staticInitOK) {
200  if (matchArgs(directVars,generateVars,_t,_mixState,_tagFlav)) return 4 ;
201  if (matchArgs(directVars,generateVars,_t,_mixState)) return 3 ;
202  if (matchArgs(directVars,generateVars,_t,_tagFlav)) return 2 ;
203  }
204 
205  if (matchArgs(directVars,generateVars,_t)) return 1 ;
206  return 0 ;
207 }
208 
209 
210 
211 ////////////////////////////////////////////////////////////////////////////////
212 
214 {
215  switch (code) {
216  case 2:
217  {
218  // Calculate the fraction of B0bar events to generate
219  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
220  _tagFlav = 1 ; // B0
221  Double_t flavInt = RooRealIntegral("flavInt","flav integral",*this,RooArgSet(_t.arg())).getVal() ;
222  _genFlavFrac = flavInt/sumInt ;
223  break ;
224  }
225  case 3:
226  {
227  // Calculate the fraction of mixed events to generate
228  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg())).getVal() ;
229  _mixState = -1 ; // mixed
230  Double_t mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
231  _genMixFrac = mixInt/sumInt ;
232  break ;
233  }
234  case 4:
235  {
236  // Calculate the fraction of mixed events to generate
237  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg(),_tagFlav.arg())).getVal() ;
238  _mixState = -1 ; // mixed
239  Double_t mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
240  _genMixFrac = mixInt/sumInt ;
241 
242  // Calculate the fractio of B0bar tags for mixed and unmixed
243  RooRealIntegral dtInt("mixInt","mix integral",*this,RooArgSet(_t.arg())) ;
244  _mixState = -1 ; // Mixed
245  _tagFlav = 1 ; // B0
246  _genFlavFracMix = dtInt.getVal() / mixInt ;
247  _mixState = 1 ; // Unmixed
248  _tagFlav = 1 ; // B0
249  _genFlavFracUnmix = dtInt.getVal() / (sumInt - mixInt) ;
250  break ;
251  }
252  }
253 }
254 
255 
256 
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 /// Generate mix-state dependent
260 
262 {
263  switch(code) {
264  case 2:
265  {
266  Double_t rand = RooRandom::uniform() ;
267  _tagFlav = (Int_t) ((rand<=_genFlavFrac) ? 1 : -1) ;
268  break ;
269  }
270  case 3:
271  {
272  Double_t rand = RooRandom::uniform() ;
273  _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
274  break ;
275  }
276  case 4:
277  {
278  Double_t rand = RooRandom::uniform() ;
279  _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
280 
281  rand = RooRandom::uniform() ;
282  Double_t genFlavFrac = (_mixState==-1) ? _genFlavFracMix : _genFlavFracUnmix ;
283  _tagFlav = (Int_t) ((rand<=genFlavFrac) ? 1 : -1) ;
284  break ;
285  }
286  }
287 
288  // Generate delta-t dependent
289  while(1) {
290  Double_t rand = RooRandom::uniform() ;
291  Double_t tval(0) ;
292 
293  switch(_type) {
294  case SingleSided:
295  tval = -_tau*log(rand);
296  break ;
297  case Flipped:
298  tval= +_tau*log(rand);
299  break ;
300  case DoubleSided:
301  tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
302  break ;
303  }
304 
305  // Accept event if T is in generated range
306  Double_t dil = 1-2.*_mistag ;
307  Double_t maxAcceptProb = 1 + TMath::Abs(_delMistag) + TMath::Abs(dil) ;
308  Double_t acceptProb = (1-_tagFlav*_delMistag) + _mixState*dil*cos(_dm*tval);
309  Bool_t mixAccept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
310 
311  if (tval<_t.max() && tval>_t.min() && mixAccept) {
312  _t = tval ;
313  break ;
314  }
315  }
316 }
RooRealProxy _t
Definition: RooBMixDecay.h:58
Double_t _genFlavFracUnmix
Definition: RooBMixDecay.h:65
Int_t _basisCos
Definition: RooBMixDecay.h:60
#define assert(cond)
Definition: unittest.h:542
RooRealProxy _delMistag
Definition: RooBMixDecay.h:53
void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
RooRealProxy _mistag
Definition: RooBMixDecay.h:52
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
STL namespace.
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t _genMixFrac
Definition: RooBMixDecay.h:62
double cos(double)
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
void generateEvent(Int_t code)
Generate mix-state dependent.
friend class RooArgSet
Definition: RooAbsArg.h:469
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
Double_t _genFlavFrac
do not persist
Definition: RooBMixDecay.h:63
RooCategoryProxy _tagFlav
Definition: RooBMixDecay.h:55
ClassImp(RooBMixDecay)
RooCategoryProxy _mixState
Definition: RooBMixDecay.h:54
Double_t _genFlavFracMix
Definition: RooBMixDecay.h:64
virtual Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
Default implementation of function implementing advertised integrals.
virtual Double_t coefficient(Int_t basisIndex) const
Comp with tFit MC: must be (1 - tagFlav*...)
DecayType _type
Definition: RooBMixDecay.h:51
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:83
RooRealProxy _dm
Definition: RooBMixDecay.h:57
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
int type
Definition: TGX11.cxx:120
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
virtual ~RooBMixDecay()
Destructor.
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooRealProxy _tau
Definition: RooBMixDecay.h:56
friend class RooRealIntegral
Definition: RooAbsPdf.h:294
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
const RooAbsCategory & arg() const
Int_t _basisExp
Definition: RooBMixDecay.h:59
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
double log(double)