ROOT  6.06/09
Reference Guide
RooBDecay.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * PL, Parker C Lund, UC Irvine *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
9  * *
10  * Copyright (c) 2000-2005, Regents of the University of California *
11  * and Stanford University. All rights reserved. *
12  * *
13  * Redistribution and use in source and binary forms, *
14  * with or without modification, are permitted according to the terms *
15  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
16  *****************************************************************************/
17 
18 
19 //////////////////////////////////////////////////////////////////////////////
20 //
21 // BEGIN_HTML
22 // Most general description of B decay time distribution with effects
23 // of CP violation, mixing and life time differences. This function can
24 // be analytically convolved with any RooResolutionModel implementation
25 // END_HTML
26 //
27 
28 
29 #include "RooFit.h"
30 
31 #include "Riostream.h"
32 #include "TMath.h"
33 #include "RooBDecay.h"
34 #include "RooRealVar.h"
35 #include "RooRandom.h"
36 
37 #include "TError.h"
38 
39 using namespace std;
40 
42 
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 
46 RooBDecay::RooBDecay(const char *name, const char* title,
47  RooRealVar& t, RooAbsReal& tau, RooAbsReal& dgamma,
49  RooAbsReal& dm, const RooResolutionModel& model, DecayType type) :
50  RooAbsAnaConvPdf(name, title, model, t),
51  _t("t", "time", this, t),
52  _tau("tau", "Average Decay Time", this, tau),
53  _dgamma("dgamma", "Delta Gamma", this, dgamma),
54  _f0("f0", "Cosh Coefficient", this, f0),
55  _f1("f1", "Sinh Coefficient", this, f1),
56  _f2("f2", "Cos Coefficient", this, f2),
57  _f3("f3", "Sin Coefficient", this, f3),
58  _dm("dm", "Delta Mass", this, dm),
59  _type(type)
60 
61 {
62  //Constructor
63  switch(type)
64  {
65  case SingleSided:
66  _basisCosh = declareBasis("exp(-@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
67  _basisSinh = declareBasis("exp(-@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
68  _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
69  _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
70  break;
71  case Flipped:
72  _basisCosh = declareBasis("exp(@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
73  _basisSinh = declareBasis("exp(@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
74  _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
75  _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
76  break;
77  case DoubleSided:
78  _basisCosh = declareBasis("exp(-abs(@0)/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
79  _basisSinh = declareBasis("exp(-abs(@0)/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
80  _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau, dm));
81  _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau, dm));
82  break;
83  }
84 }
85 
86 ////////////////////////////////////////////////////////////////////////////////
87 ///Copy constructor
88 
89 RooBDecay::RooBDecay(const RooBDecay& other, const char* name) :
90  RooAbsAnaConvPdf(other, name),
91  _t("t", this, other._t),
92  _tau("tau", this, other._tau),
93  _dgamma("dgamma", this, other._dgamma),
94  _f0("f0", this, other._f0),
95  _f1("f1", this, other._f1),
96  _f2("f2", this, other._f2),
97  _f3("f3", this, other._f3),
98  _dm("dm", this, other._dm),
99  _basisCosh(other._basisCosh),
100  _basisSinh(other._basisSinh),
101  _basisCos(other._basisCos),
102  _basisSin(other._basisSin),
103  _type(other._type)
104 {
105 }
106 
107 
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 ///Destructor
111 
113 {
114 }
115 
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 
120 {
121  if(basisIndex == _basisCosh)
122  {
123  return _f0;
124  }
125  if(basisIndex == _basisSinh)
126  {
127  return _f1;
128  }
129  if(basisIndex == _basisCos)
130  {
131  return _f2;
132  }
133  if(basisIndex == _basisSin)
134  {
135  return _f3;
136  }
137 
138  return 0 ;
139 }
140 
141 
142 ////////////////////////////////////////////////////////////////////////////////
143 
145 {
146  if(basisIndex == _basisCosh)
147  {
148  return _f0.arg().getVariables();
149  }
150  if(basisIndex == _basisSinh)
151  {
152  return _f1.arg().getVariables();
153  }
154  if(basisIndex == _basisCos)
155  {
156  return _f2.arg().getVariables();
157  }
158  if(basisIndex == _basisSin)
159  {
160  return _f3.arg().getVariables();
161  }
162 
163  return 0 ;
164 }
165 
166 
167 
168 ////////////////////////////////////////////////////////////////////////////////
169 
170 Int_t RooBDecay::getCoefAnalyticalIntegral(Int_t coef, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
171 {
172  if(coef == _basisCosh)
173  {
174  return _f0.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
175  }
176  if(coef == _basisSinh)
177  {
178  return _f1.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
179  }
180  if(coef == _basisCos)
181  {
182  return _f2.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
183  }
184  if(coef == _basisSin)
185  {
186  return _f3.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
187  }
188 
189  return 0 ;
190 }
191 
192 
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 
196 Double_t RooBDecay::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* rangeName) const
197 {
198  if(coef == _basisCosh)
199  {
200  return _f0.arg().analyticalIntegral(code,rangeName) ;
201  }
202  if(coef == _basisSinh)
203  {
204  return _f1.arg().analyticalIntegral(code,rangeName) ;
205  }
206  if(coef == _basisCos)
207  {
208  return _f2.arg().analyticalIntegral(code,rangeName) ;
209  }
210  if(coef == _basisSin)
211  {
212  return _f3.arg().analyticalIntegral(code,rangeName) ;
213  }
214 
215  return 0 ;
216 }
217 
218 
219 
220 ////////////////////////////////////////////////////////////////////////////////
221 
222 Int_t RooBDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
223 {
224  if (matchArgs(directVars, generateVars, _t)) return 1;
225  return 0;
226 }
227 
228 
229 
230 ////////////////////////////////////////////////////////////////////////////////
231 
233 {
234  R__ASSERT(code==1);
235  Double_t gammamin = 1/_tau-TMath::Abs(_dgamma)/2;
236  while(1) {
237  Double_t t = -log(RooRandom::uniform())/gammamin;
238  if (_type == Flipped || (_type == DoubleSided && RooRandom::uniform() <0.5) ) t *= -1;
239  if ( t<_t.min() || t>_t.max() ) continue;
240 
241  Double_t dgt = _dgamma*t/2;
242  Double_t dmt = _dm*t;
243  Double_t ft = fabs(t);
244  Double_t f = exp(-ft/_tau)*(_f0*cosh(dgt)+_f1*sinh(dgt)+_f2*cos(dmt)+_f3*sin(dmt));
245  if(f < 0) {
246  cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: PDF value less than zero" << endl;
247  ::abort();
248  }
249  Double_t w = 1.001*exp(-ft*gammamin)*(TMath::Abs(_f0)+TMath::Abs(_f1)+sqrt(_f2*_f2+_f3*_f3));
250  if(w < f) {
251  cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: Envelope function less than p.d.f. " << endl;
252  cout << _f0 << endl;
253  cout << _f1 << endl;
254  cout << _f2 << endl;
255  cout << _f3 << endl;
256  ::abort();
257  }
258  if(w*RooRandom::uniform() > f) continue;
259  _t = t;
260  break;
261  }
262 }
263 
264 
265 
266 
267 
268 
269 
270 
271 
272 
273 
274 
275 
276 
277 
278 
279 
280 
281 
282 
283 
284 
285 
286 
287 
288 
289 
290 
291 
292 
293 
294 
295 
296 
297 
298 
299 
300 
301 
302 
303 
304 
305 
306 
307 
308 
309 
310 
311 
312 
313 
314 
virtual ~RooBDecay()
Destructor.
Definition: RooBDecay.cxx:112
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Definition: RooAbsReal.cxx:337
#define R__ASSERT(e)
Definition: TError.h:98
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
ClassImp(RooBDecay)
Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
Definition: RooBDecay.cxx:170
STL namespace.
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
double cos(double)
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
RooRealProxy _f0
Definition: RooBDecay.h:61
double sqrt(double)
RooArgSet * coefVars(Int_t coefIdx) const
Return set of parameters with are used exclusively by the coefficient functions.
Definition: RooBDecay.cxx:144
double sinh(double)
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2081
RooRealProxy _f2
Definition: RooBDecay.h:63
DecayType _type
Definition: RooBDecay.h:71
RooRealProxy _dm
Definition: RooBDecay.h:65
double sin(double)
RooRealProxy _f1
Definition: RooBDecay.h:62
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...
Definition: RooBDecay.cxx:222
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Definition: RooAbsReal.cxx:363
double cosh(double)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
Definition: RooBDecay.cxx:232
Int_t _basisSinh
Definition: RooBDecay.h:67
Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
Default implementation of function implementing advertised integrals.
Definition: RooBDecay.cxx:196
RooRealProxy _tau
Definition: RooBDecay.h:59
double f(double x)
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
RooRealProxy _t
Definition: RooBDecay.h:58
RooRealProxy _f3
Definition: RooBDecay.h:64
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.
RooRealProxy _dgamma
Definition: RooBDecay.h:60
#define name(a, b)
Definition: linkTestLib0.cpp:5
Int_t _basisSin
Definition: RooBDecay.h:69
Int_t _basisCosh
Definition: RooBDecay.h:66
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
double f2(const double *x)
TF1 * f1
Definition: legend1.C:11
RooBDecay()
Definition: RooBDecay.h:32
Int_t _basisCos
Definition: RooBDecay.h:68
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
virtual Double_t coefficient(Int_t basisIndex) const
Definition: RooBDecay.cxx:119
double exp(double)
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
double log(double)