Logo ROOT  
Reference Guide
RooGaussModel.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 RooGaussModel
18 \ingroup Roofit
19
20Class RooGaussModel implements a RooResolutionModel that models a Gaussian
21distribution. Object of class RooGaussModel can be used
22for analytical convolutions with classes inheriting from RooAbsAnaConvPdf
23**/
24
25#include "RooFit.h"
26
27#include "TMath.h"
28#include "Riostream.h"
29#include "RooGaussModel.h"
30#include "RooRealConstant.h"
31#include "RooRandom.h"
32
33#include "TError.h"
34
35using namespace std;
36
38
39////////////////////////////////////////////////////////////////////////////////
40
41RooGaussModel::RooGaussModel(const char *name, const char *title, RooAbsRealLValue& xIn,
42 RooAbsReal& _mean, RooAbsReal& _sigma) :
43 RooResolutionModel(name,title,xIn),
44 _flatSFInt(kFALSE),
45 _asympInt(kFALSE),
46 mean("mean","Mean",this,_mean),
47 sigma("sigma","Width",this,_sigma),
48 msf("msf","Mean Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
49 ssf("ssf","Sigma Scale Factor",this,(RooRealVar&)RooRealConstant::value(1))
50{
51}
52
53////////////////////////////////////////////////////////////////////////////////
54
55RooGaussModel::RooGaussModel(const char *name, const char *title, RooAbsRealLValue& xIn,
56 RooAbsReal& _mean, RooAbsReal& _sigma,
57 RooAbsReal& _msSF) :
58 RooResolutionModel(name,title,xIn),
59 _flatSFInt(kFALSE),
60 _asympInt(kFALSE),
61 mean("mean","Mean",this,_mean),
62 sigma("sigma","Width",this,_sigma),
63 msf("msf","Mean Scale Factor",this,_msSF),
64 ssf("ssf","Sigma Scale Factor",this,_msSF)
65{
66}
67
68////////////////////////////////////////////////////////////////////////////////
69
70RooGaussModel::RooGaussModel(const char *name, const char *title, RooAbsRealLValue& xIn,
71 RooAbsReal& _mean, RooAbsReal& _sigma,
72 RooAbsReal& _meanSF, RooAbsReal& _sigmaSF) :
73 RooResolutionModel(name,title,xIn),
74 _flatSFInt(kFALSE),
75 _asympInt(kFALSE),
76 mean("mean","Mean",this,_mean),
77 sigma("sigma","Width",this,_sigma),
78 msf("msf","Mean Scale Factor",this,_meanSF),
79 ssf("ssf","Sigma Scale Factor",this,_sigmaSF)
80{
81}
82
83////////////////////////////////////////////////////////////////////////////////
84
87 _flatSFInt(other._flatSFInt),
88 _asympInt(other._asympInt),
89 mean("mean",this,other.mean),
90 sigma("sigma",this,other.sigma),
91 msf("msf",this,other.msf),
92 ssf("ssf",this,other.ssf)
93{
94}
95
96////////////////////////////////////////////////////////////////////////////////
97/// Destructor
98
100{
101}
102
103////////////////////////////////////////////////////////////////////////////////
104
106{
107 if (!TString("exp(-@0/@1)").CompareTo(name)) return expBasisPlus ;
108 if (!TString("exp(@0/@1)").CompareTo(name)) return expBasisMinus ;
109 if (!TString("exp(-abs(@0)/@1)").CompareTo(name)) return expBasisSum ;
110 if (!TString("exp(-@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisPlus ;
111 if (!TString("exp(@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisMinus ;
112 if (!TString("exp(-abs(@0)/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisSum ;
113 if (!TString("exp(-@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisPlus ;
114 if (!TString("exp(@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisMinus ;
115 if (!TString("exp(-abs(@0)/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisSum ;
116 if (!TString("(@0/@1)*exp(-@0/@1)").CompareTo(name)) return linBasisPlus ;
117 if (!TString("(@0/@1)*(@0/@1)*exp(-@0/@1)").CompareTo(name)) return quadBasisPlus ;
118 if (!TString("exp(-@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisPlus;
119 if (!TString("exp(@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisMinus;
120 if (!TString("exp(-abs(@0)/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisSum;
121 if (!TString("exp(-@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisPlus;
122 if (!TString("exp(@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisMinus;
123 if (!TString("exp(-abs(@0)/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisSum;
124 return 0 ;
125}
126
127////////////////////////////////////////////////////////////////////////////////
128
130{
131 // *** 1st form: Straight Gaussian, used for unconvoluted PDF or expBasis with 0 lifetime ***
132 static Double_t root2(std::sqrt(2.)) ;
133 static Double_t root2pi(std::sqrt(2.*std::atan2(0.,-1.))) ;
134 static Double_t rootpi(std::sqrt(std::atan2(0.,-1.))) ;
135
136 BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
137 BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
138
140 if (basisType == coshBasis && _basisCode!=noBasis ) {
141 Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
142 if (dGamma==0) basisType = expBasis;
143 }
144
145 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
146 Double_t xprime = (x-(mean*msf))/(sigma*ssf) ;
147 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 1st form" << endl ;
148
149 Double_t result = std::exp(-0.5*xprime*xprime)/(sigma*ssf*root2pi) ;
150 if (_basisCode!=0 && basisSign==Both) result *= 2 ;
151 return result ;
152 }
153
154 // *** 2nd form: 0, used for sinBasis, linBasis, and quadBasis with tau=0 ***
155 if (tau==0) {
156 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 2nd form" << endl ;
157 return 0. ;
158 }
159
160 // *** 3nd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
161 Double_t omega = (basisType==sinBasis || basisType==cosBasis) ? ((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
162 Double_t dgamma = (basisType==sinhBasis || basisType==coshBasis) ? ((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
163 Double_t _x = omega *tau ;
164 Double_t _y = tau*dgamma/2;
165 Double_t xprime = (x-(mean*msf))/tau ;
166 Double_t c = (sigma*ssf)/(root2*tau) ;
167 Double_t u = xprime/(2*c) ;
168
169 if (basisType==expBasis || (basisType==cosBasis && _x==0.)) {
170 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 3d form tau=" << tau << endl ;
171 Double_t result(0) ;
172 if (basisSign!=Minus) result += evalCerf(0,-u,c).real();
173 if (basisSign!=Plus) result += evalCerf(0, u,c).real();
174 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::getVal(" << GetName() << ") got nan during case 1 " << endl; }
175 return result ;
176 }
177
178 // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
179 if (basisType==sinBasis) {
180 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 4th form omega = " << omega << ", tau = " << tau << endl ;
181 Double_t result(0) ;
182 if (_x==0.) return result ;
183 if (basisSign!=Minus) result += -evalCerf(-_x,-u,c).imag();
184 if (basisSign!=Plus) result += -evalCerf( _x, u,c).imag();
185 if (TMath::IsNaN(result)) cxcoutE(Tracing) << "RooGaussModel::getVal(" << GetName() << ") got nan during case 3 " << endl;
186 return result ;
187 }
188
189 // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
190 if (basisType==cosBasis) {
191 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
192 Double_t result(0) ;
193 if (basisSign!=Minus) result += evalCerf(-_x,-u,c).real();
194 if (basisSign!=Plus) result += evalCerf( _x, u,c).real();
195 if (TMath::IsNaN(result)) cxcoutE(Tracing) << "RooGaussModel::getVal(" << GetName() << ") got nan during case 4 " << endl;
196 return result ;
197 }
198
199 // ***8th form: Convolution with exp(-|t|/tau)*cosh(dgamma*t/2), used for coshBasisSum ***
200 if (basisType==coshBasis || basisType ==sinhBasis) {
201 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 8th form tau = " << tau << endl ;
202 Double_t result(0);
203 int sgn = ( basisType == coshBasis ? +1 : -1 );
204 if (basisSign!=Minus) result += 0.5*( evalCerf(0,-u,c*(1-_y)).real()+sgn*evalCerf(0,-u,c*(1+_y)).real()) ;
205 if (basisSign!=Plus) result += 0.5*(sgn*evalCerf(0, u,c*(1-_y)).real()+ evalCerf(0, u,c*(1+_y)).real()) ;
206 if (TMath::IsNaN(result)) cxcoutE(Tracing) << "RooGaussModel::getVal(" << GetName() << ") got nan during case 8 " << endl;
207 return result ;
208 }
209
210 // *** 6th form: Convolution with (t/tau)*exp(-t/tau), used for linBasis ***
211 if (basisType==linBasis) {
212 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 6th form tau = " << tau << endl ;
213 R__ASSERT(basisSign==Plus); // This should only be for positive times
214
215 Double_t f0 = std::exp(-xprime+c*c) * RooMath::erfc(-u+c);
216 Double_t f1 = std::exp(-u*u);
217 return (xprime - 2*c*c)*f0 + (2*c/rootpi)*f1 ;
218 }
219
220 // *** 7th form: Convolution with (t/tau)^2*exp(-t/tau), used for quadBasis ***
221 if (basisType==quadBasis) {
222 if (verboseEval()>2) cout << "RooGaussModel::evaluate(" << GetName() << ") 7th form tau = " << tau << endl ;
223 R__ASSERT(basisSign==Plus); // This should only be for positive times
224
225 Double_t f0 = std::exp(-xprime+c*c) * RooMath::erfc(-u+c);
226 Double_t f1 = std::exp(-u*u);
227 Double_t x2c2 = xprime - 2*c*c;
228 return ( x2c2*x2c2*f0 + (2*c/rootpi)*x2c2*f1 + 2*c*c*f0 );
229 }
230
231 R__ASSERT(0) ;
232 return 0 ;
233}
234
235////////////////////////////////////////////////////////////////////////////////
236
237Int_t RooGaussModel::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
238{
239 switch(_basisCode) {
240
241 // Analytical integration capability of raw PDF
242 case noBasis:
243
244 // Optionally advertise flat integral over sigma scale factor
245 if (_flatSFInt) {
246 if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) return 2 ;
247 }
248
249 if (matchArgs(allVars,analVars,convVar())) return 1 ;
250 break ;
251
252 // Analytical integration capability of convoluted PDF
253 case expBasisPlus:
254 case expBasisMinus:
255 case expBasisSum:
256 case sinBasisPlus:
257 case sinBasisMinus:
258 case sinBasisSum:
259 case cosBasisPlus:
260 case cosBasisMinus:
261 case cosBasisSum:
262 case linBasisPlus:
263 case quadBasisPlus:
264 case coshBasisMinus:
265 case coshBasisPlus:
266 case coshBasisSum:
267 case sinhBasisMinus:
268 case sinhBasisPlus:
269 case sinhBasisSum:
270
271 // Optionally advertise flat integral over sigma scale factor
272 if (_flatSFInt) {
273
274 if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) {
275 return 2 ;
276 }
277 }
278
279 if (matchArgs(allVars,analVars,convVar())) return 1 ;
280 break ;
281 }
282
283 return 0 ;
284}
285
286////////////////////////////////////////////////////////////////////////////////
287
288Double_t RooGaussModel::analyticalIntegral(Int_t code, const char* rangeName) const
289{
290 static const Double_t root2 = std::sqrt(2.) ;
291 //static Double_t rootPiBy2 = std::sqrt(std::atan2(0.0,-1.0)/2.0);
292 static const Double_t rootpi = std::sqrt(std::atan2(0.0,-1.0));
293 Double_t ssfInt(1.0) ;
294
295 // Code must be 1 or 2
296 R__ASSERT(code==1||code==2) ;
297 if (code==2) ssfInt = (ssf.max(rangeName)-ssf.min(rangeName)) ;
298
299 BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
300 BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
301
302 // *** 1st form: Straight Gaussian, used for unconvoluted PDF or expBasis with 0 lifetime ***
304 if (basisType == coshBasis && _basisCode!=noBasis ) {
305 Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
306 if (dGamma==0) basisType = expBasis;
307 }
308 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
309 Double_t xscale = root2*(sigma*ssf);
310 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 1st form" << endl ;
311
312 Double_t xpmin = (x.min(rangeName)-(mean*msf))/xscale ;
313 Double_t xpmax = (x.max(rangeName)-(mean*msf))/xscale ;
314
315 Double_t result ;
316 if (_asympInt) { // modified FMV, 07/24/03
317 result = 1.0 ;
318 } else {
319 result = 0.5*(RooMath::erf(xpmax)-RooMath::erf(xpmin)) ;
320 }
321
322 if (_basisCode!=0 && basisSign==Both) result *= 2 ;
323 //cout << "Integral 1st form " << " result= " << result*ssfInt << endl;
324 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 1 " << endl; }
325 return result*ssfInt ;
326 }
327
328
329 Double_t omega = ((basisType==sinBasis)||(basisType==cosBasis)) ? ((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
330 Double_t dgamma =((basisType==sinhBasis)||(basisType==coshBasis)) ? ((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
331
332 // *** 2nd form: unity, used for sinBasis and linBasis with tau=0 (PDF is zero) ***
333 if (tau==0) {
334 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 2nd form" << endl ;
335 return 0. ;
336 }
337
338 // *** 3rd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
339 Double_t c = (sigma*ssf)/(root2*tau) ;
340 Double_t xpmin = (x.min(rangeName)-(mean*msf))/tau ;
341 Double_t xpmax = (x.max(rangeName)-(mean*msf))/tau ;
342 Double_t umin = xpmin/(2*c) ;
343 Double_t umax = xpmax/(2*c) ;
344
345 if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
346 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 3d form tau=" << tau << endl ;
347 Double_t result(0) ;
348 if (basisSign!=Minus) result += evalCerfInt(+1,0,tau,-umin,-umax,c).real();
349 if (basisSign!=Plus) result += evalCerfInt(-1,0,tau, umin, umax,c).real();
350 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 3 " << endl; }
351 return result*ssfInt ;
352 }
353
354 // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
355 Double_t _x = omega * tau ;
356 Double_t _y = tau*dgamma/2;
357
358 if (basisType==sinBasis) {
359 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 4th form omega = " << omega << ", tau = " << tau << endl ;
360 Double_t result(0) ;
361 if (_x==0) return result*ssfInt ;
362 if (basisSign!=Minus) result += -1*evalCerfInt(+1,-_x,tau,-umin,-umax,c).imag();
363 if (basisSign!=Plus) result += -1*evalCerfInt(-1, _x,tau, umin, umax,c).imag();
364 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 4 " << endl; }
365 return result*ssfInt ;
366 }
367
368 // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
369 if (basisType==cosBasis) {
370 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
371 Double_t result(0) ;
372 if (basisSign!=Minus) result += evalCerfInt(+1,-_x,tau,-umin,-umax,c).real();
373 if (basisSign!=Plus) result += evalCerfInt(-1, _x,tau, umin, umax,c).real();
374 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 5 " << endl; }
375 return result*ssfInt ;
376 }
377
378 // *** 8th form: Convolution with exp(-|t|/tau)*cosh(dgamma*t/2), used for coshBasis ***
379 // *** 9th form: Convolution with exp(-|t|/tau)*sinh(dgamma*t/2), used for sinhBasis ***
380 if (basisType==coshBasis || basisType == sinhBasis) {
381 if (verboseEval()>0) {cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 8th form tau=" << tau << endl ; }
382 Double_t result(0) ;
383 int sgn = ( basisType == coshBasis ? +1 : -1 );
384 if (basisSign!=Minus) result += 0.5*( evalCerfInt(+1,0,tau/(1-_y),-umin,-umax,c*(1-_y)).real()+ sgn*evalCerfInt(+1,0,tau/(1+_y),-umin,-umax,c*(1+_y)).real());
385 if (basisSign!=Plus) result += 0.5*(sgn*evalCerfInt(-1,0,tau/(1-_y), umin, umax,c*(1-_y)).real()+ evalCerfInt(-1,0,tau/(1+_y), umin, umax,c*(1+_y)).real());
386 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 6 " << endl; }
387 return result*ssfInt ;
388 }
389
390 // *** 6th form: Convolution with (t/tau)*exp(-t/tau), used for linBasis ***
391 if (basisType==linBasis) {
392 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 6th form tau=" << tau << endl ;
393
394 Double_t f0 = RooMath::erf(-umax) - RooMath::erf(-umin);
395 Double_t f1 = std::exp(-umax*umax) - std::exp(-umin*umin);
396
397 Double_t tmp1 = std::exp(-xpmax)*RooMath::erfc(-umax + c);
398 Double_t tmp2 = std::exp(-xpmin)*RooMath::erfc(-umin + c);
399
400 Double_t f2 = tmp1 - tmp2;
401 Double_t f3 = xpmax*tmp1 - xpmin*tmp2;
402
403 Double_t expc2 = std::exp(c*c);
404
405 return -tau*( f0 +
406 (2*c/rootpi)*f1 +
407 (1 - 2*c*c)*expc2*f2 +
408 expc2*f3
409 )*ssfInt;
410 }
411
412 // *** 7th form: Convolution with (t/tau)*(t/tau)*exp(-t/tau), used for quadBasis ***
413 if (basisType==quadBasis) {
414 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 7th form tau=" << tau << endl ;
415
416 Double_t f0 = RooMath::erf(-umax) - RooMath::erf(-umin);
417
418 Double_t tmpA1 = std::exp(-umax*umax);
419 Double_t tmpA2 = std::exp(-umin*umin);
420
421 Double_t f1 = tmpA1 - tmpA2;
422 Double_t f2 = umax*tmpA1 - umin*tmpA2;
423
424 Double_t tmpB1 = std::exp(-xpmax)*RooMath::erfc(-umax + c);
425 Double_t tmpB2 = std::exp(-xpmin)*RooMath::erfc(-umin + c);
426
427 Double_t f3 = tmpB1 - tmpB2;
428 Double_t f4 = xpmax*tmpB1 - xpmin*tmpB2;
429 Double_t f5 = xpmax*xpmax*tmpB1 - xpmin*xpmin*tmpB2;
430
431 Double_t expc2 = std::exp(c*c);
432
433 return -tau*( 2*f0 +
434 (4*c/rootpi)*((1-c*c)*f1 + c*f2) +
435 (2*c*c*(2*c*c-1) + 2)*expc2*f3 - (4*c*c-2)*expc2*f4 + expc2*f5
436 )*ssfInt;
437 }
438 R__ASSERT(0) ;
439 return 0 ;
440}
441
442////////////////////////////////////////////////////////////////////////////////
443/// use the approximation: erf(z) = exp(-z*z)/(std::sqrt(pi)*z)
444/// to explicitly cancel the divergent exp(y*y) behaviour of
445/// CWERF for z = x + i y with large negative y
446
448{
449 static const Double_t rootpi= std::sqrt(std::atan2(0.,-1.));
450 const std::complex<Double_t> z(_x * c, u + c);
451 const std::complex<Double_t> zc(u + c, - _x * c);
452 const std::complex<Double_t> zsq((z.real() + z.imag()) * (z.real() - z.imag()),
453 2. * z.real() * z.imag());
454 const std::complex<Double_t> v(-zsq.real() - u*u, -zsq.imag());
455 const std::complex<Double_t> ev = std::exp(v);
456 const std::complex<Double_t> mez2zcrootpi = -std::exp(zsq)/(zc*rootpi);
457
458 return 2. * (ev * (mez2zcrootpi + 1.));
459}
460
461////////////////////////////////////////////////////////////////////////////////
462
463std::complex<Double_t> RooGaussModel::evalCerfInt(Double_t sign, Double_t _x, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
464{
465 std::complex<Double_t> diff(2., 0.);
466 if (!_asympInt) {
467 diff = evalCerf(_x,umin,c);
468 diff -= evalCerf(_x,umax,c);
469 diff += RooMath::erf(umin) - RooMath::erf(umax);
470 diff *= sign;
471 }
472 diff *= std::complex<Double_t>(1., _x);
473 diff *= tau / (1.+_x*_x);
474 return diff;
475}
476
477////////////////////////////////////////////////////////////////////////////////
478
479Int_t RooGaussModel::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
480{
481 if (matchArgs(directVars,generateVars,x)) return 1 ;
482 return 0 ;
483}
484
485////////////////////////////////////////////////////////////////////////////////
486
488{
489 R__ASSERT(code==1) ;
490 Double_t xmin = x.min();
491 Double_t xmax = x.max();
492 TRandom *generator = RooRandom::randomGenerator();
493 while(true) {
494 Double_t xgen = generator->Gaus(mean*msf,sigma*ssf);
495 if (xgen<xmax && xgen>xmin) {
496 x = xgen ;
497 return ;
498 }
499 }
500}
#define c(i)
Definition: RSha256.hxx:101
#define cxcoutE(a)
Definition: RooMsgService.h:97
const Bool_t kFALSE
Definition: RtypesCore.h:90
#define ClassImp(name)
Definition: Rtypes.h:361
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
double atan2(double, double)
double sqrt(double)
double exp(double)
friend class RooArgSet
Definition: RooAbsArg.h:572
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Definition: RooAbsPdf.cxx:3303
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooAbsArg * getParameter(const char *name) const
Definition: RooFormulaVar.h:40
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26
RooRealProxy sigma
Definition: RooGaussModel.h:83
Bool_t _asympInt
Definition: RooGaussModel.h:80
Bool_t _flatSFInt
Definition: RooGaussModel.h:78
RooRealProxy msf
Definition: RooGaussModel.h:84
virtual ~RooGaussModel()
Destructor.
std::complex< Double_t > evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
static std::complex< Double_t > evalCerf(Double_t swt, Double_t u, Double_t c)
Definition: RooGaussModel.h:69
static std::complex< Double_t > evalCerfApprox(Double_t swt, Double_t u, Double_t c)
use the approximation: erf(z) = exp(-z*z)/(std::sqrt(pi)*z) to explicitly cancel the divergent exp(y*...
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
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...
virtual Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
virtual Int_t basisCode(const char *name) const
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooRealProxy mean
Definition: RooGaussModel.h:82
RooRealProxy ssf
Definition: RooGaussModel.h:85
static std::complex< double > erfc(const std::complex< double > z)
complex erfc function
Definition: RooMath.cxx:556
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:580
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:53
RooRealConstant provides static functions to create and keep track of RooRealVar constants.
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
RooAbsRealLValue & convVar() const
Return the convolution variable of the resolution model.
const RooFormulaVar & basis() const
RooTemplateProxy< RooAbsRealLValue > x
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
const T & arg() const
Return reference to object held in proxy.
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:27
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:263
Basic string class.
Definition: TString.h:131
const Double_t sigma
TF1 * f1
Definition: legend1.C:11
Bool_t IsNaN(Double_t x)
Definition: TMath.h:882