Logo ROOT   6.07/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 \file RooBMixDecay.cxx
19 \class RooBMixDecay
20 \ingroup Roofit
21 
22 Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes
23 the decay of B mesons with the effects of B0/B0bar mixing.
24 This function can be analytically convolved with any RooResolutionModel implementation
25 **/
26 
27 #include "RooFit.h"
28 
29 #include "Riostream.h"
30 #include "TMath.h"
31 #include "RooRealVar.h"
32 #include "RooBMixDecay.h"
33 #include "RooRealIntegral.h"
34 #include "RooRandom.h"
35 
36 using namespace std;
37 
39 ;
40 
41 
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 /// Constructor
45 
46 RooBMixDecay::RooBMixDecay(const char *name, const char *title,
47  RooRealVar& t, RooAbsCategory& mixState,
48  RooAbsCategory& tagFlav,
49  RooAbsReal& tau, RooAbsReal& dm,
50  RooAbsReal& mistag, RooAbsReal& delMistag,
51  const RooResolutionModel& model,
52  DecayType type) :
53  RooAbsAnaConvPdf(name,title,model,t),
54  _type(type),
55  _mistag("mistag","Mistag rate",this,mistag),
56  _delMistag("delMistag","Delta mistag rate",this,delMistag),
57  _mixState("mixState","Mixing state",this,mixState),
58  _tagFlav("tagFlav","Flavour of tagged B0",this,tagFlav),
59  _tau("tau","Mixing life time",this,tau),
60  _dm("dm","Mixing frequency",this,dm),
61  _t("_t","time",this,t), _genMixFrac(0)
62 {
63  switch(type) {
64  case SingleSided:
65  _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
66  _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
67  break ;
68  case Flipped:
69  _basisExp = declareBasis("exp(@0/@1)",RooArgList(tau,dm)) ;
70  _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
71  break ;
72  case DoubleSided:
73  _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
74  _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
75  break ;
76  }
77 }
78 
79 
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 /// Copy constructor
83 
84 RooBMixDecay::RooBMixDecay(const RooBMixDecay& other, const char* name) :
85  RooAbsAnaConvPdf(other,name),
86  _type(other._type),
87  _mistag("mistag",this,other._mistag),
88  _delMistag("delMistag",this,other._delMistag),
89  _mixState("mixState",this,other._mixState),
90  _tagFlav("tagFlav",this,other._tagFlav),
91  _tau("tau",this,other._tau),
92  _dm("dm",this,other._dm),
93  _t("t",this,other._t),
94  _basisExp(other._basisExp),
95  _basisCos(other._basisCos),
96  _genMixFrac(other._genMixFrac),
100 {
101 }
102 
103 
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 /// Destructor
107 
109 {
110 }
111 
112 
113 
114 ////////////////////////////////////////////////////////////////////////////////
115 /// Comp with tFit MC: must be (1 - tagFlav*...)
116 
118 {
119  if (basisIndex==_basisExp) {
120  return (1 - _tagFlav*_delMistag) ;
121  }
122 
123  if (basisIndex==_basisCos) {
124  return _mixState*(1-2*_mistag) ;
125  }
126 
127  return 0 ;
128 }
129 
130 
131 
132 ////////////////////////////////////////////////////////////////////////////////
133 /// cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
134 
135 Int_t RooBMixDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
136 {
137  if (rangeName) {
138  return 0 ;
139  }
140 
141  if (matchArgs(allVars,analVars,_mixState,_tagFlav)) return 3 ;
142  if (matchArgs(allVars,analVars,_mixState)) return 2 ;
143  if (matchArgs(allVars,analVars,_tagFlav)) return 1 ;
144  return 0 ;
145 }
146 
147 
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 
151 Double_t RooBMixDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
152 {
153  switch(code) {
154  // No integration
155  case 0: return coefficient(basisIndex) ;
156 
157  // Integration over 'mixState' and 'tagFlav'
158  case 3:
159  if (basisIndex==_basisExp) {
160  return 4.0 ;
161  }
162  if (basisIndex==_basisCos) {
163  return 0.0 ;
164  }
165  break ;
166 
167  // Integration over 'mixState'
168  case 2:
169  if (basisIndex==_basisExp) {
170  return 2.0*coefficient(basisIndex) ;
171  }
172  if (basisIndex==_basisCos) {
173  return 0.0 ;
174  }
175  break ;
176 
177  // Integration over 'tagFlav'
178  case 1:
179  if (basisIndex==_basisExp) {
180  return 2.0 ;
181  }
182  if (basisIndex==_basisCos) {
183  return 2.0*coefficient(basisIndex) ;
184  }
185  break ;
186 
187  default:
188  assert(0) ;
189  }
190 
191  return 0 ;
192 }
193 
194 
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 
198 Int_t RooBMixDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
199 {
200  if (staticInitOK) {
201  if (matchArgs(directVars,generateVars,_t,_mixState,_tagFlav)) return 4 ;
202  if (matchArgs(directVars,generateVars,_t,_mixState)) return 3 ;
203  if (matchArgs(directVars,generateVars,_t,_tagFlav)) return 2 ;
204  }
205 
206  if (matchArgs(directVars,generateVars,_t)) return 1 ;
207  return 0 ;
208 }
209 
210 
211 
212 ////////////////////////////////////////////////////////////////////////////////
213 
215 {
216  switch (code) {
217  case 2:
218  {
219  // Calculate the fraction of B0bar events to generate
220  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
221  _tagFlav = 1 ; // B0
222  Double_t flavInt = RooRealIntegral("flavInt","flav integral",*this,RooArgSet(_t.arg())).getVal() ;
223  _genFlavFrac = flavInt/sumInt ;
224  break ;
225  }
226  case 3:
227  {
228  // Calculate the fraction of mixed events to generate
229  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg())).getVal() ;
230  _mixState = -1 ; // mixed
231  Double_t mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
232  _genMixFrac = mixInt/sumInt ;
233  break ;
234  }
235  case 4:
236  {
237  // Calculate the fraction of mixed events to generate
238  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg(),_tagFlav.arg())).getVal() ;
239  _mixState = -1 ; // mixed
240  Double_t mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
241  _genMixFrac = mixInt/sumInt ;
242 
243  // Calculate the fractio of B0bar tags for mixed and unmixed
244  RooRealIntegral dtInt("mixInt","mix integral",*this,RooArgSet(_t.arg())) ;
245  _mixState = -1 ; // Mixed
246  _tagFlav = 1 ; // B0
247  _genFlavFracMix = dtInt.getVal() / mixInt ;
248  _mixState = 1 ; // Unmixed
249  _tagFlav = 1 ; // B0
250  _genFlavFracUnmix = dtInt.getVal() / (sumInt - mixInt) ;
251  break ;
252  }
253  }
254 }
255 
256 
257 
258 
259 ////////////////////////////////////////////////////////////////////////////////
260 /// Generate mix-state dependent
261 
263 {
264  switch(code) {
265  case 2:
266  {
267  Double_t rand = RooRandom::uniform() ;
268  _tagFlav = (Int_t) ((rand<=_genFlavFrac) ? 1 : -1) ;
269  break ;
270  }
271  case 3:
272  {
273  Double_t rand = RooRandom::uniform() ;
274  _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
275  break ;
276  }
277  case 4:
278  {
279  Double_t rand = RooRandom::uniform() ;
280  _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
281 
282  rand = RooRandom::uniform() ;
283  Double_t genFlavFrac = (_mixState==-1) ? _genFlavFracMix : _genFlavFracUnmix ;
284  _tagFlav = (Int_t) ((rand<=genFlavFrac) ? 1 : -1) ;
285  break ;
286  }
287  }
288 
289  // Generate delta-t dependent
290  while(1) {
291  Double_t rand = RooRandom::uniform() ;
292  Double_t tval(0) ;
293 
294  switch(_type) {
295  case SingleSided:
296  tval = -_tau*log(rand);
297  break ;
298  case Flipped:
299  tval= +_tau*log(rand);
300  break ;
301  case DoubleSided:
302  tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
303  break ;
304  }
305 
306  // Accept event if T is in generated range
307  Double_t dil = 1-2.*_mistag ;
308  Double_t maxAcceptProb = 1 + TMath::Abs(_delMistag) + TMath::Abs(dil) ;
309  Double_t acceptProb = (1-_tagFlav*_delMistag) + _mixState*dil*cos(_dm*tval);
310  Bool_t mixAccept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
311 
312  if (tval<_t.max() && tval>_t.min() && mixAccept) {
313  _t = tval ;
314  break ;
315  }
316  }
317 }
RooRealProxy _t
Definition: RooBMixDecay.h:58
Double_t _genFlavFracUnmix
Definition: RooBMixDecay.h:65
Int_t _basisCos
Definition: RooBMixDecay.h:60
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.
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects The class perfor...
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
Definition: RooBMixDecay.h:23
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
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
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
#define ClassImp(name)
Definition: Rtypes.h:279
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:84
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.
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
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
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
char name[80]
Definition: TGX11.cxx:109
double log(double)