Logo ROOT   6.07/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 \file RooBDecay.cxx
21 \class RooBDecay
22 \ingroup Roofit
23 
24 Most general description of B decay time distribution with effects
25 of CP violation, mixing and life time differences. This function can
26 be analytically convolved with any RooResolutionModel implementation
27 **/
28 
29 
30 #include "RooFit.h"
31 
32 #include "Riostream.h"
33 #include "TMath.h"
34 #include "RooBDecay.h"
35 #include "RooRealVar.h"
36 #include "RooRandom.h"
37 
38 #include "TError.h"
39 
40 using namespace std;
41 
43 
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 
47 RooBDecay::RooBDecay(const char *name, const char* title,
48  RooRealVar& t, RooAbsReal& tau, RooAbsReal& dgamma,
50  RooAbsReal& dm, const RooResolutionModel& model, DecayType type) :
51  RooAbsAnaConvPdf(name, title, model, t),
52  _t("t", "time", this, t),
53  _tau("tau", "Average Decay Time", this, tau),
54  _dgamma("dgamma", "Delta Gamma", this, dgamma),
55  _f0("f0", "Cosh Coefficient", this, f0),
56  _f1("f1", "Sinh Coefficient", this, f1),
57  _f2("f2", "Cos Coefficient", this, f2),
58  _f3("f3", "Sin Coefficient", this, f3),
59  _dm("dm", "Delta Mass", this, dm),
60  _type(type)
61 
62 {
63  //Constructor
64  switch(type)
65  {
66  case SingleSided:
67  _basisCosh = declareBasis("exp(-@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
68  _basisSinh = declareBasis("exp(-@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
69  _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
70  _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
71  break;
72  case Flipped:
73  _basisCosh = declareBasis("exp(@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
74  _basisSinh = declareBasis("exp(@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
75  _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
76  _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
77  break;
78  case DoubleSided:
79  _basisCosh = declareBasis("exp(-abs(@0)/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
80  _basisSinh = declareBasis("exp(-abs(@0)/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
81  _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau, dm));
82  _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau, dm));
83  break;
84  }
85 }
86 
87 ////////////////////////////////////////////////////////////////////////////////
88 ///Copy constructor
89 
90 RooBDecay::RooBDecay(const RooBDecay& other, const char* name) :
91  RooAbsAnaConvPdf(other, name),
92  _t("t", this, other._t),
93  _tau("tau", this, other._tau),
94  _dgamma("dgamma", this, other._dgamma),
95  _f0("f0", this, other._f0),
96  _f1("f1", this, other._f1),
97  _f2("f2", this, other._f2),
98  _f3("f3", this, other._f3),
99  _dm("dm", this, other._dm),
100  _basisCosh(other._basisCosh),
101  _basisSinh(other._basisSinh),
102  _basisCos(other._basisCos),
103  _basisSin(other._basisSin),
104  _type(other._type)
105 {
106 }
107 
108 
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 ///Destructor
112 
114 {
115 }
116 
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 
121 {
122  if(basisIndex == _basisCosh)
123  {
124  return _f0;
125  }
126  if(basisIndex == _basisSinh)
127  {
128  return _f1;
129  }
130  if(basisIndex == _basisCos)
131  {
132  return _f2;
133  }
134  if(basisIndex == _basisSin)
135  {
136  return _f3;
137  }
138 
139  return 0 ;
140 }
141 
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 
146 {
147  if(basisIndex == _basisCosh)
148  {
149  return _f0.arg().getVariables();
150  }
151  if(basisIndex == _basisSinh)
152  {
153  return _f1.arg().getVariables();
154  }
155  if(basisIndex == _basisCos)
156  {
157  return _f2.arg().getVariables();
158  }
159  if(basisIndex == _basisSin)
160  {
161  return _f3.arg().getVariables();
162  }
163 
164  return 0 ;
165 }
166 
167 
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 
171 Int_t RooBDecay::getCoefAnalyticalIntegral(Int_t coef, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
172 {
173  if(coef == _basisCosh)
174  {
175  return _f0.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
176  }
177  if(coef == _basisSinh)
178  {
179  return _f1.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
180  }
181  if(coef == _basisCos)
182  {
183  return _f2.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
184  }
185  if(coef == _basisSin)
186  {
187  return _f3.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
188  }
189 
190  return 0 ;
191 }
192 
193 
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 
197 Double_t RooBDecay::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* rangeName) const
198 {
199  if(coef == _basisCosh)
200  {
201  return _f0.arg().analyticalIntegral(code,rangeName) ;
202  }
203  if(coef == _basisSinh)
204  {
205  return _f1.arg().analyticalIntegral(code,rangeName) ;
206  }
207  if(coef == _basisCos)
208  {
209  return _f2.arg().analyticalIntegral(code,rangeName) ;
210  }
211  if(coef == _basisSin)
212  {
213  return _f3.arg().analyticalIntegral(code,rangeName) ;
214  }
215 
216  return 0 ;
217 }
218 
219 
220 
221 ////////////////////////////////////////////////////////////////////////////////
222 
223 Int_t RooBDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
224 {
225  if (matchArgs(directVars, generateVars, _t)) return 1;
226  return 0;
227 }
228 
229 
230 
231 ////////////////////////////////////////////////////////////////////////////////
232 
234 {
235  R__ASSERT(code==1);
236  Double_t gammamin = 1/_tau-TMath::Abs(_dgamma)/2;
237  while(1) {
238  Double_t t = -log(RooRandom::uniform())/gammamin;
239  if (_type == Flipped || (_type == DoubleSided && RooRandom::uniform() <0.5) ) t *= -1;
240  if ( t<_t.min() || t>_t.max() ) continue;
241 
242  Double_t dgt = _dgamma*t/2;
243  Double_t dmt = _dm*t;
244  Double_t ft = fabs(t);
245  Double_t f = exp(-ft/_tau)*(_f0*cosh(dgt)+_f1*sinh(dgt)+_f2*cos(dmt)+_f3*sin(dmt));
246  if(f < 0) {
247  cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: PDF value less than zero" << endl;
248  ::abort();
249  }
250  Double_t w = 1.001*exp(-ft*gammamin)*(TMath::Abs(_f0)+TMath::Abs(_f1)+sqrt(_f2*_f2+_f3*_f3));
251  if(w < f) {
252  cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: Envelope function less than p.d.f. " << endl;
253  cout << _f0 << endl;
254  cout << _f1 << endl;
255  cout << _f2 << endl;
256  cout << _f3 << endl;
257  ::abort();
258  }
259  if(w*RooRandom::uniform() > f) continue;
260  _t = t;
261  break;
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 
315 
virtual ~RooBDecay()
Destructor.
Definition: RooBDecay.cxx:113
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
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:171
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:145
Most general description of B decay time distribution with effects of CP violation, mixing and life time differences.
Definition: RooBDecay.h:24
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
double sinh(double)
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2082
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
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
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:223
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:233
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:197
RooRealProxy _tau
Definition: RooBDecay.h:59
#define ClassImp(name)
Definition: Rtypes.h:279
double f(double x)
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:84
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
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:120
double exp(double)
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
char name[80]
Definition: TGX11.cxx:109
double log(double)