Logo ROOT   6.08/07
Reference Guide
RooBCPEffDecay.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 /** \class RooBCPEffDecay
18  \ingroup Roofit
19 
20 PDF describing decay time distribution of B meson including effects of standard model CP violation.
21 This function can be analytically convolved with any RooResolutionModel implementation.
22 */
23 
24 
25 #include "RooFit.h"
26 
27 #include "Riostream.h"
28 #include "Riostream.h"
29 #include "RooRealVar.h"
30 #include "RooRandom.h"
31 #include "RooBCPEffDecay.h"
32 #include "RooRealIntegral.h"
33 
34 using namespace std;
35 
37 ;
38 
39 
40 
41 ////////////////////////////////////////////////////////////////////////////////
42 /// Constructor.
43 
44 RooBCPEffDecay::RooBCPEffDecay(const char *name, const char *title,
45  RooRealVar& t, RooAbsCategory& tag,
47  RooAbsReal& avgMistag, RooAbsReal& CPeigenval,
49  RooAbsReal& effRatio, RooAbsReal& delMistag,
51  RooAbsAnaConvPdf(name,title,model,t),
52  _absLambda("absLambda","Absolute value of lambda",this,a),
53  _argLambda("argLambda","Arg(Lambda)",this,b),
54  _effRatio("effRatio","B0/B0bar efficiency ratio",this,effRatio),
55  _CPeigenval("CPeigenval","CP eigen value",this,CPeigenval),
56  _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
57  _delMistag("delMistag","Delta mistag rate",this,delMistag),
58  _t("t","time",this,t),
59  _tau("tau","decay time",this,tau),
60  _dm("dm","mixing frequency",this,dm),
61  _tag("tag","CP state",this,tag),
62  _genB0Frac(0),
63  _type(type)
64 {
65  switch(type) {
66  case SingleSided:
67  _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
68  _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
69  _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
70  break ;
71  case Flipped:
72  _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
73  _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
74  _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
75  break ;
76  case DoubleSided:
77  _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
78  _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
79  _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
80  break ;
81  }
82 }
83 
84 
85 
86 ////////////////////////////////////////////////////////////////////////////////
87 /// Copy constructor.
88 
89 RooBCPEffDecay::RooBCPEffDecay(const RooBCPEffDecay& other, const char* name) :
90  RooAbsAnaConvPdf(other,name),
91  _absLambda("absLambda",this,other._absLambda),
92  _argLambda("argLambda",this,other._argLambda),
93  _effRatio("effRatio",this,other._effRatio),
94  _CPeigenval("CPeigenval",this,other._CPeigenval),
95  _avgMistag("avgMistag",this,other._avgMistag),
96  _delMistag("delMistag",this,other._delMistag),
97  _t("t",this,other._t),
98  _tau("tau",this,other._tau),
99  _dm("dm",this,other._dm),
100  _tag("tag",this,other._tag),
101  _genB0Frac(other._genB0Frac),
102  _type(other._type),
103  _basisExp(other._basisExp),
104  _basisSin(other._basisSin),
105  _basisCos(other._basisCos)
106 {
107 }
108 
109 
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 /// Destructor.
113 
115 {
116 }
117 
118 
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 /// B0 : _tag = +1
122 ///
123 /// B0bar : _tag = -1
124 /// \param[in] basisIndex
125 
127 {
128  if (basisIndex==_basisExp) {
129  //exp term: (1 -/+ dw)(1+a^2)/2
130  return (1 - _tag*_delMistag)*(1+_absLambda*_absLambda)/2 ;
131  // = 1 + _tag*deltaDil/2
132  }
133 
134  if (basisIndex==_basisSin) {
135  //sin term: +/- (1-2w)*ImLambda(= -etaCP * |l| * sin2b)
137  // = _tag*avgDil * ...
138  }
139 
140  if (basisIndex==_basisCos) {
141  //cos term: +/- (1-2w)*(1-a^2)/2
142  return -1*_tag*(1-2*_avgMistag)*(1-_absLambda*_absLambda)/2 ;
143  // = -_tag*avgDil * ...
144  }
145 
146  return 0 ;
147 }
148 
149 
150 
151 ////////////////////////////////////////////////////////////////////////////////
152 
153 Int_t RooBCPEffDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
154 {
155  if (rangeName) return 0 ;
156 
157  if (matchArgs(allVars,analVars,_tag)) return 1 ;
158  return 0 ;
159 }
160 
161 
162 
163 ////////////////////////////////////////////////////////////////////////////////
164 
165 Double_t RooBCPEffDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
166 {
167  switch(code) {
168  // No integration
169  case 0: return coefficient(basisIndex) ;
170 
171  // Integration over 'tag'
172  case 1:
173  if (basisIndex==_basisExp) {
174  return (1+_absLambda*_absLambda) ;
175  }
176 
177  if (basisIndex==_basisSin || basisIndex==_basisCos) {
178  return 0 ;
179  }
180  break ;
181 
182  default:
183  assert(0) ;
184  }
185 
186  return 0 ;
187 }
188 
189 
190 
191 ////////////////////////////////////////////////////////////////////////////////
192 
193 Int_t RooBCPEffDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
194 {
195  if (staticInitOK) {
196  if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
197  }
198  if (matchArgs(directVars,generateVars,_t)) return 1 ;
199  return 0 ;
200 }
201 
202 
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 
207 {
208  if (code==2) {
209  // Calculate the fraction of mixed events to generate
210  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
211  _tag = 1 ;
212  Double_t b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
213  _genB0Frac = b0Int/sumInt ;
214  }
215 }
216 
217 
218 
219 ////////////////////////////////////////////////////////////////////////////////
220 /// Generates mix-state dependent.
221 /// \param[in] code
222 
224 {
225  if (code==2) {
226  Double_t rand = RooRandom::uniform() ;
227  _tag = (rand<=_genB0Frac) ? 1 : -1 ;
228  }
229 
230  // Generate delta-t dependent
231  while(1) {
232  Double_t rand = RooRandom::uniform() ;
233  Double_t tval(0) ;
234 
235  switch(_type) {
236  case SingleSided:
237  tval = -_tau*log(rand);
238  break ;
239  case Flipped:
240  tval= +_tau*log(rand);
241  break ;
242  case DoubleSided:
243  tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
244  break ;
245  }
246 
247  // Accept event if T is in generated range
248  Double_t maxDil = 1.0 ;
250  Double_t maxAcceptProb = (1+al2) + fabs(maxDil*_CPeigenval*_absLambda*_argLambda) + fabs(maxDil*(1-al2)/2);
251  Double_t acceptProb = (1+al2)/2*(1-_tag*_delMistag)
252  - (_tag*(1-2*_avgMistag))*(_CPeigenval*_absLambda*_argLambda)*sin(_dm*tval)
253  - (_tag*(1-2*_avgMistag))*(1-al2)/2*cos(_dm*tval);
254 
255  Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
256 
257  if (tval<_t.max() && tval>_t.min() && accept) {
258  _t = tval ;
259  break ;
260  }
261  }
262 
263 }
264 
virtual Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
Default implementation of function implementing advertised integrals.
RooRealProxy _t
RooRealProxy _absLambda
DecayType _type
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooRealProxy _argLambda
RooRealProxy _effRatio
RooRealProxy _dm
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
double cos(double)
PDF describing decay time distribution of B meson including effects of standard model CP violation...
Double_t _genB0Frac
void generateEvent(Int_t code)
Generates mix-state dependent.
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
const RooAbsCategory & arg() const
friend class RooArgSet
Definition: RooAbsArg.h:469
virtual ~RooBCPEffDecay()
Destructor.
double sin(double)
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooRealProxy _avgMistag
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
RooRealProxy _tau
#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
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
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 Double_t coefficient(Int_t basisIndex) const
B0 : _tag = +1.
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...
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
friend class RooRealIntegral
Definition: RooAbsPdf.h:294
RooCategoryProxy _tag
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
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 Double_t b
Definition: TRolke.cxx:630
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
const Bool_t kTRUE
Definition: Rtypes.h:91
RooRealProxy _delMistag
RooRealProxy _CPeigenval
char name[80]
Definition: TGX11.cxx:109
double log(double)