43 RooBDecay::RooBDecay(
const char *name,
const char* title,
44 RooRealVar& t, RooAbsReal& tau, RooAbsReal& dgamma,
45 RooAbsReal& f0, RooAbsReal& f1, RooAbsReal& f2, RooAbsReal& f3,
46 RooAbsReal& dm,
const RooResolutionModel& model, DecayType type) :
47 RooAbsAnaConvPdf(name, title, model, t),
48 _t(
"t",
"time", this, t),
49 _tau(
"tau",
"Average Decay Time", this, tau),
50 _dgamma(
"dgamma",
"Delta Gamma", this, dgamma),
51 _f0(
"f0",
"Cosh Coefficient", this, f0),
52 _f1(
"f1",
"Sinh Coefficient", this, f1),
53 _f2(
"f2",
"Cos Coefficient", this, f2),
54 _f3(
"f3",
"Sin Coefficient", this, f3),
55 _dm(
"dm",
"Delta Mass", this, dm),
63 _basisCosh = declareBasis(
"exp(-@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
64 _basisSinh = declareBasis(
"exp(-@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
65 _basisCos = declareBasis(
"exp(-@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
66 _basisSin = declareBasis(
"exp(-@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
69 _basisCosh = declareBasis(
"exp(@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
70 _basisSinh = declareBasis(
"exp(@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
71 _basisCos = declareBasis(
"exp(@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
72 _basisSin = declareBasis(
"exp(@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
75 _basisCosh = declareBasis(
"exp(-abs(@0)/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
76 _basisSinh = declareBasis(
"exp(-abs(@0)/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
77 _basisCos = declareBasis(
"exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau, dm));
78 _basisSin = declareBasis(
"exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau, dm));
86 RooBDecay::RooBDecay(
const RooBDecay& other,
const char* name) :
87 RooAbsAnaConvPdf(other, name),
88 _t(
"t", this, other._t),
89 _tau(
"tau", this, other._tau),
90 _dgamma(
"dgamma", this, other._dgamma),
91 _f0(
"f0", this, other._f0),
92 _f1(
"f1", this, other._f1),
93 _f2(
"f2", this, other._f2),
94 _f3(
"f3", this, other._f3),
95 _dm(
"dm", this, other._dm),
96 _basisCosh(other._basisCosh),
97 _basisSinh(other._basisSinh),
98 _basisCos(other._basisCos),
99 _basisSin(other._basisSin),
107 RooBDecay::~RooBDecay()
113 Double_t RooBDecay::coefficient(Int_t basisIndex)
const
115 if(basisIndex == _basisCosh)
119 if(basisIndex == _basisSinh)
123 if(basisIndex == _basisCos)
127 if(basisIndex == _basisSin)
137 RooArgSet* RooBDecay::coefVars(Int_t basisIndex)
const
139 if(basisIndex == _basisCosh)
141 return _f0.arg().getVariables();
143 if(basisIndex == _basisSinh)
145 return _f1.arg().getVariables();
147 if(basisIndex == _basisCos)
149 return _f2.arg().getVariables();
151 if(basisIndex == _basisSin)
153 return _f3.arg().getVariables();
161 Int_t RooBDecay::getCoefAnalyticalIntegral(Int_t coef, RooArgSet& allVars, RooArgSet& analVars,
const char* rangeName)
const
163 if(coef == _basisCosh)
165 return _f0.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
167 if(coef == _basisSinh)
169 return _f1.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
171 if(coef == _basisCos)
173 return _f2.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
175 if(coef == _basisSin)
177 return _f3.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
185 Double_t RooBDecay::coefAnalyticalIntegral(Int_t coef, Int_t code,
const char* rangeName)
const
187 if(coef == _basisCosh)
189 return _f0.arg().analyticalIntegral(code,rangeName) ;
191 if(coef == _basisSinh)
193 return _f1.arg().analyticalIntegral(code,rangeName) ;
195 if(coef == _basisCos)
197 return _f2.arg().analyticalIntegral(code,rangeName) ;
199 if(coef == _basisSin)
201 return _f3.arg().analyticalIntegral(code,rangeName) ;
209 Int_t RooBDecay::getGenerator(
const RooArgSet& directVars, RooArgSet &generateVars, Bool_t )
const
211 if (matchArgs(directVars, generateVars, _t))
return 1;
217 void RooBDecay::generateEvent(Int_t code)
220 Double_t gammamin = 1/_tau-TMath::Abs(_dgamma)/2;
222 Double_t t = -log(RooRandom::uniform())/gammamin;
223 if (_type == Flipped || (_type == DoubleSided && RooRandom::uniform() <0.5) ) t *= -1;
224 if ( t<_t.min() || t>_t.max() )
continue;
226 Double_t dgt = _dgamma*t/2;
227 Double_t dmt = _dm*t;
228 Double_t ft = fabs(t);
229 Double_t f = exp(-ft/_tau)*(_f0*cosh(dgt)+_f1*sinh(dgt)+_f2*cos(dmt)+_f3*sin(dmt));
231 cout <<
"RooBDecay::generateEvent(" << GetName() <<
") ERROR: PDF value less than zero" << endl;
234 Double_t w = 1.001*exp(-ft*gammamin)*(TMath::Abs(_f0)+TMath::Abs(_f1)+sqrt(_f2*_f2+_f3*_f3));
236 cout <<
"RooBDecay::generateEvent(" << GetName() <<
") ERROR: Envelope function less than p.d.f. " << endl;
243 if(w*RooRandom::uniform() > f)
continue;