ROOT  6.06/09
Reference Guide
RooNonCPEigenDecay.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * AH, Andreas Hoecker, Orsay, hoecker@slac.stanford.edu *
7  * SL, Sandrine Laplace, Orsay, laplace@slac.stanford.edu *
8  * JS, Jan Stark, Paris, stark@slac.stanford.edu *
9  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
10  * *
11  * Copyright (c) 2000-2005, Regents of the University of California, *
12  * IN2P3. All rights reserved. *
13  * *
14  * History *
15  * Nov-2001 WV Created initial version *
16  * Dec-2001 SL mischarge correction, direct CPV *
17  * Jan-2002 AH built dedicated generator + code cleaning *
18  * Mar-2002 JS committed debugged version to CVS *
19  * Apr-2002 AH allow prompt (ie, non-Pdf) mischarge treatment *
20  * May-2002 JS Changed the set of CP parameters (mathematically equiv.) *
21  * *
22  * Redistribution and use in source and binary forms, *
23  * with or without modification, are permitted according to the terms *
24  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
25  *****************************************************************************/
26 
27 //////////////////////////////////////////////////////////////////////////////
28 //
29 // BEGIN_HTML
30 // Time-dependent RooAbsAnaConvPdf for CP violating decays
31 // to Non-CP eigenstates (eg, B0 -> rho+- pi-+).
32 // For a description of the physics model see the
33 // BaBar Physics Book, section 6.5.2.3 .
34 // The set of CP parameters used in this class is equivalent to
35 // the one used in the Physics Book, but it is not exactly the
36 // same. Starting from the set in the BaBar Book, in order to
37 // get the parameters used here you have to change the sign of both
38 // a_c^+ and a_c^-, and then substitute:
39 // <pre>
40 // a_s^Q = S + Q* deltaS
41 // a_c^Q = C + Q*deltaC
42 // </pre>
43 // where Q denotes the charge of the rho.
44 // END_HTML
45 //
46 
47 #include "RooFit.h"
48 
49 #include "Riostream.h"
50 #include "Riostream.h"
51 #include "RooRealVar.h"
52 #include "RooRandom.h"
53 #include "RooNonCPEigenDecay.h"
54 #include "TMath.h"
55 #include "RooRealIntegral.h"
56 
57 using namespace std;
58 
60 
61 #define Debug_RooNonCPEigenDecay 1
62 
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 
66 RooNonCPEigenDecay::RooNonCPEigenDecay( const char *name, const char *title,
67  RooRealVar& t,
68  RooAbsCategory& tag,
69  RooAbsReal& tau,
70  RooAbsReal& dm,
71  RooAbsReal& avgW,
72  RooAbsReal& delW,
73  RooAbsCategory& rhoQ,
74  RooAbsReal& correctQ,
75  RooAbsReal& wQ,
76  RooAbsReal& acp,
77  RooAbsReal& C,
78  RooAbsReal& delC,
79  RooAbsReal& S,
80  RooAbsReal& delS,
81  const RooResolutionModel& model,
82  DecayType type )
83  : RooAbsAnaConvPdf( name, title, model, t ),
84  _acp ( "acp", "acp", this, acp ),
85  _avgC ( "C", "C", this, C ),
86  _delC ( "delC", "delC", this, delC ),
87  _avgS ( "S", "S", this, S ),
88  _delS ( "delS", "delS", this, delS ),
89  _avgW ( "avgW", "Average mistag rate",this, avgW ),
90  _delW ( "delW", "Shift mistag rate", this, delW ),
91  _t ( "t", "time", this, t ),
92  _tau ( "tau", "decay time", this, tau ),
93  _dm ( "dm", "mixing frequency", this, dm ),
94  _tag ( "tag", "CP state", this, tag ),
95  _rhoQ ( "rhoQ", "Charge of the rho", this, rhoQ ),
96  _correctQ ( "correctQ", "correction of rhoQ", this, correctQ ),
97  _wQ ( "wQ", "mischarge", this, wQ ),
98  _genB0Frac ( 0 ),
99  _genRhoPlusFrac( 0 ),
100  _type ( type )
101 {
102  // Constructor
103  switch(type) {
104  case SingleSided:
105  _basisExp = declareBasis( "exp(-@0/@1)", RooArgList( tau ) );
106  _basisSin = declareBasis( "exp(-@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
107  _basisCos = declareBasis( "exp(-@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
108  break;
109  case Flipped:
110  _basisExp = declareBasis( "exp(@0)/@1)", RooArgList( tau ) );
111  _basisSin = declareBasis( "exp(@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
112  _basisCos = declareBasis( "exp(@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
113  break;
114  case DoubleSided:
115  _basisExp = declareBasis( "exp(-abs(@0)/@1)", RooArgList( tau ) );
116  _basisSin = declareBasis( "exp(-abs(@0)/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
117  _basisCos = declareBasis( "exp(-abs(@0)/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
118  break;
119  }
120 }
121 
122 
123 ////////////////////////////////////////////////////////////////////////////////
124 
125 RooNonCPEigenDecay::RooNonCPEigenDecay( const char *name, const char *title,
126  RooRealVar& t,
127  RooAbsCategory& tag,
128  RooAbsReal& tau,
129  RooAbsReal& dm,
130  RooAbsReal& avgW,
131  RooAbsReal& delW,
132  RooAbsCategory& rhoQ,
133  RooAbsReal& correctQ,
134  RooAbsReal& acp,
135  RooAbsReal& C,
136  RooAbsReal& delC,
137  RooAbsReal& S,
138  RooAbsReal& delS,
139  const RooResolutionModel& model,
140  DecayType type )
141  : RooAbsAnaConvPdf( name, title, model, t ),
142  _acp ( "acp", "acp", this, acp ),
143  _avgC ( "C", "C", this, C ),
144  _delC ( "delC", "delC", this, delC ),
145  _avgS ( "S", "S", this, S ),
146  _delS ( "delS", "delS", this, delS ),
147  _avgW ( "avgW", "Average mistag rate",this, avgW ),
148  _delW ( "delW", "Shift mistag rate", this, delW ),
149  _t ( "t", "time", this, t ),
150  _tau ( "tau", "decay time", this, tau ),
151  _dm ( "dm", "mixing frequency", this, dm ),
152  _tag ( "tag", "CP state", this, tag ),
153  _rhoQ ( "rhoQ", "Charge of the rho", this, rhoQ ),
154  _correctQ ( "correctQ", "correction of rhoQ", this, correctQ ),
155  _genB0Frac ( 0 ),
156  _genRhoPlusFrac( 0 ),
157  _type ( type )
158 {
159  // dummy mischarge (must be set to zero!)
160  _wQ = RooRealProxy( "wQ", "mischarge", this, *(new RooRealVar( "wQ", "wQ", 0 )) );
161 
162  switch(type) {
163  case SingleSided:
164  _basisExp = declareBasis( "exp(-@0/@1)", RooArgList( tau ) );
165  _basisSin = declareBasis( "exp(-@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
166  _basisCos = declareBasis( "exp(-@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
167  break;
168  case Flipped:
169  _basisExp = declareBasis( "exp(@0)/@1)", RooArgList( tau ) );
170  _basisSin = declareBasis( "exp(@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
171  _basisCos = declareBasis( "exp(@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
172  break;
173  case DoubleSided:
174  _basisExp = declareBasis( "exp(-abs(@0)/@1)", RooArgList( tau ) );
175  _basisSin = declareBasis( "exp(-abs(@0)/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
176  _basisCos = declareBasis( "exp(-abs(@0)/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
177  break;
178  }
179 }
180 
181 
182 ////////////////////////////////////////////////////////////////////////////////
183 /// Copy constructor
184 
186  : RooAbsAnaConvPdf( other, name ),
187  _acp ( "acp", this, other._acp ),
188  _avgC ( "C", this, other._avgC ),
189  _delC ( "delC", this, other._delC ),
190  _avgS ( "S", this, other._avgS ),
191  _delS ( "delS", this, other._delS ),
192  _avgW ( "avgW", this, other._avgW ),
193  _delW ( "delW", this, other._delW ),
194  _t ( "t", this, other._t ),
195  _tau ( "tau", this, other._tau ),
196  _dm ( "dm", this, other._dm ),
197  _tag ( "tag", this, other._tag ),
198  _rhoQ ( "rhoQ", this, other._rhoQ ),
199  _correctQ ( "correctQ", this, other._correctQ ),
200  _wQ ( "wQ", this, other._wQ ),
201  _genB0Frac ( other._genB0Frac ),
202  _genRhoPlusFrac( other._genRhoPlusFrac ),
203  _type ( other._type ),
204  _basisExp ( other._basisExp ),
205  _basisSin ( other._basisSin ),
206  _basisCos ( other._basisCos )
207 {
208 }
209 
210 
211 ////////////////////////////////////////////////////////////////////////////////
212 /// Destructor
213 
215 {
216 }
217 
218 
219 ////////////////////////////////////////////////////////////////////////////////
220 /// B0 : _tag == -1
221 /// B0bar : _tag == +1
222 /// rho+ : _rhoQ == +1
223 /// rho- : _rhoQ == -1
224 /// the charge corrrection factor "_correctQ" serves to implement mis-charges
225 
227 {
228  Int_t rhoQc = _rhoQ * int(_correctQ);
229  assert( rhoQc == 1 || rhoQc == -1 );
230 
231  Double_t a_sin_p = _avgS + _delS;
232  Double_t a_sin_m = _avgS - _delS;
233  Double_t a_cos_p = _avgC + _delC;
234  Double_t a_cos_m = _avgC - _delC;
235 
236  if (basisIndex == _basisExp) {
237  if (rhoQc == -1 || rhoQc == +1)
238  return (1 + rhoQc*_acp*(1 - 2*_wQ))*(1 + 0.5*_tag*(2*_delW));
239  else
240  return 1;
241  }
242 
243  if (basisIndex == _basisSin) {
244 
245  if (rhoQc == -1)
246 
247  return - ((1 - _acp)*a_sin_m*(1 - _wQ) + (1 + _acp)*a_sin_p*_wQ)*(1 - 2*_avgW)*_tag;
248 
249  else if (rhoQc == +1)
250 
251  return - ((1 + _acp)*a_sin_p*(1 - _wQ) + (1 - _acp)*a_sin_m*_wQ)*(1 - 2*_avgW)*_tag;
252 
253  else
254  return - _tag*((a_sin_p + a_sin_m)/2)*(1 - 2*_avgW);
255  }
256 
257  if (basisIndex == _basisCos) {
258 
259  if ( rhoQc == -1)
260 
261  return + ((1 - _acp)*a_cos_m*(1 - _wQ) + (1 + _acp)*a_cos_p*_wQ)*(1 - 2*_avgW)*_tag;
262 
263  else if (rhoQc == +1)
264 
265  return + ((1 + _acp)*a_cos_p*(1 - _wQ) + (1 - _acp)*a_cos_m*_wQ)*(1 - 2*_avgW)*_tag;
266 
267  else
268  return _tag*((a_cos_p + a_cos_m)/2)*(1 - 2*_avgW);
269  }
270 
271  return 0;
272 }
273 
274 // advertise analytical integration
275 
276 ////////////////////////////////////////////////////////////////////////////////
277 
279  RooArgSet& analVars, const char* rangeName ) const
280 {
281  if (rangeName) return 0 ;
282 
283  if (matchArgs( allVars, analVars, _tag, _rhoQ )) return 3;
284  if (matchArgs( allVars, analVars, _rhoQ )) return 2;
285  if (matchArgs( allVars, analVars, _tag )) return 1;
286 
287  return 0;
288 }
289 
290 
291 ////////////////////////////////////////////////////////////////////////////////
292 /// correct for the right/wrong charge...
293 
295  Int_t code, const char* /*rangeName*/ ) const
296 {
297  Int_t rhoQc = _rhoQ*int(_correctQ);
298 
299  Double_t a_sin_p = _avgS + _delS;
300  Double_t a_sin_m = _avgS - _delS;
301  Double_t a_cos_p = _avgC + _delC;
302  Double_t a_cos_m = _avgC - _delC;
303 
304  switch(code) {
305 
306  // No integration
307  case 0: return coefficient(basisIndex);
308 
309  // Integration over 'tag'
310  case 1:
311  if (basisIndex == _basisExp) return 2*(1 + rhoQc*_acp*(1 - 2*_wQ));
312  if (basisIndex == _basisSin || basisIndex==_basisCos) return 0;
313  assert( kFALSE );
314 
315  // Integration over 'rhoQ'
316  case 2:
317  if (basisIndex == _basisExp) return 2*(1 + 0.5*_tag*(2.*_delW));
318 
319  if (basisIndex == _basisSin)
320 
321  return - ( (1 - _acp)*a_sin_m + (1 + _acp)*a_sin_p )*(1 - 2*_avgW)*_tag;
322 
323  if (basisIndex == _basisCos)
324 
325  return + ( (1 - _acp)*a_cos_m + (1 + _acp)*a_cos_p )*(1 - 2*_avgW)*_tag;
326 
327  assert( kFALSE );
328 
329  // Integration over 'tag' and 'rhoQ'
330  case 3:
331  if (basisIndex == _basisExp) return 2*2; // for both: tag and charge
332  if (basisIndex == _basisSin || basisIndex==_basisCos) return 0;
333  assert( kFALSE );
334 
335  default:
336  assert( kFALSE );
337  }
338 
339  return 0;
340 }
341 
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 
346  RooArgSet& generateVars, Bool_t staticInitOK ) const
347 {
348  if (staticInitOK) {
349  if (matchArgs( directVars, generateVars, _t, _tag, _rhoQ )) return 4;
350  if (matchArgs( directVars, generateVars, _t, _rhoQ )) return 3;
351  if (matchArgs( directVars, generateVars, _t, _tag )) return 2;
352  }
353  if (matchArgs( directVars, generateVars, _t )) return 1;
354  return 0;
355 }
356 
357 
358 ////////////////////////////////////////////////////////////////////////////////
359 
361 {
362  if (code == 2 || code == 4) {
363  // Calculate the fraction of mixed events to generate
364  Double_t sumInt1 = RooRealIntegral( "sumInt1", "sum integral1", *this,
365  RooArgSet( _t.arg(), _tag.arg(), _rhoQ.arg() )
366  ).getVal();
367  _tag = -1;
368  Double_t b0Int1 = RooRealIntegral( "mixInt1", "mix integral1", *this,
369  RooArgSet( _t.arg(), _rhoQ.arg() )
370  ).getVal();
371  _genB0Frac = b0Int1/sumInt1;
372 
373  if (Debug_RooNonCPEigenDecay == 1)
374  cout << " o RooNonCPEigenDecay::initgenerator: genB0Frac : "
375  << _genB0Frac
376  << ", tag dilution: " << (1 - 2*_avgW)
377  << endl;
378  }
379 
380  if (code == 3 || code == 4) {
381  // Calculate the fraction of positive rho's to generate
382  Double_t sumInt2 = RooRealIntegral( "sumInt2", "sum integral2", *this,
383  RooArgSet( _t.arg(), _tag.arg(), _rhoQ.arg() )
384  ).getVal();
385  _rhoQ = 1;
386  Double_t b0Int2 = RooRealIntegral( "mixInt2", "mix integral2", *this,
387  RooArgSet( _t.arg(), _tag.arg() )
388  ).getVal();
389  _genRhoPlusFrac = b0Int2/sumInt2;
390 
391  if (Debug_RooNonCPEigenDecay == 1)
392  cout << " o RooNonCPEigenDecay::initgenerator: genRhoPlusFrac: "
393  << _genRhoPlusFrac << endl;
394  }
395 }
396 
397 
398 
399 ////////////////////////////////////////////////////////////////////////////////
400 
402 {
403  // Generate delta-t dependent
404  while (kTRUE) {
405 
406  // B flavor and rho charge (we do not use the integrated weights)
407  if (code != 1) {
408  if (code != 3) _tag = (RooRandom::uniform()<=0.5) ? -1 : +1;
409  if (code != 2) _rhoQ = (RooRandom::uniform()<=0.5) ? 1 : -1;
410  }
411 
412  // opposite charge?
413  // Int_t rhoQc = _rhoQ*int(_correctQ);
414 
415  Double_t a_sin_p = _avgS + _delS;
416  Double_t a_sin_m = _avgS - _delS;
417  Double_t a_cos_p = _avgC + _delC;
418  Double_t a_cos_m = _avgC - _delC;
419 
420  // maximum probability density
421  double a1 = 1 + sqrt(TMath::Power(a_cos_m, 2) + TMath::Power(a_sin_m, 2));
422  double a2 = 1 + sqrt(TMath::Power(a_cos_p, 2) + TMath::Power(a_sin_p, 2));
423 
424  Double_t maxAcceptProb = (1.10 + TMath::Abs(_acp)) * (a1 > a2 ? a1 : a2);
425  // The 1.10 in the above line is a security feature to prevent crashes close to the limit at 1.00
426 
427  Double_t rand = RooRandom::uniform();
428  Double_t tval(0);
429 
430  switch(_type) {
431 
432  case SingleSided:
433  tval = -_tau*log(rand);
434  break;
435 
436  case Flipped:
437  tval = +_tau*log(rand);
438  break;
439 
440  case DoubleSided:
441  tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5));
442  break;
443  }
444 
445  // get coefficients
446  Double_t expC = coefficient( _basisExp );
447  Double_t sinC = coefficient( _basisSin );
448  Double_t cosC = coefficient( _basisCos );
449 
450  // probability density
451  Double_t acceptProb = expC + sinC*sin(_dm*tval) + cosC*cos(_dm*tval);
452 
453  // sanity check...
454  assert( acceptProb <= maxAcceptProb );
455 
456  // hit or miss...
457  Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE;
458 
459  if (accept && tval<_t.max() && tval>_t.min()) {
460  _t = tval;
461  break;
462  }
463  }
464 }
465 
#define Debug_RooNonCPEigenDecay
virtual Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
correct for the right/wrong charge...
#define assert(cond)
Definition: unittest.h:542
virtual Double_t coefficient(Int_t basisIndex) const
B0 : _tag == -1 B0bar : _tag == +1 rho+ : _rhoQ == +1 rho- : _rhoQ == -1 the charge corrrection facto...
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
ClassImp(RooNonCPEigenDecay)
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
double cos(double)
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
double sqrt(double)
RooCategoryProxy _tag
friend class RooArgSet
Definition: RooAbsArg.h:469
double sin(double)
void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
static double C[]
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:83
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
static const float S
Definition: mandel.cpp:113
virtual ~RooNonCPEigenDecay(void)
Destructor.
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
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...
RooCategoryProxy _rhoQ
#define name(a, b)
Definition: linkTestLib0.cpp:5
friend class RooRealProxy
Definition: RooAbsReal.h:403
friend class RooRealIntegral
Definition: RooAbsPdf.h:294
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
const RooAbsCategory & arg() const
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
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
double log(double)