Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "TMath.h"
26#include "Riostream.h"
27#include "RooGaussModel.h"
28#include "RooMath.h"
29#include "RooRealConstant.h"
30#include "RooRandom.h"
31#include "RooBatchCompute.h"
32
33#include "TError.h"
34
36
37using namespace std;
38
41
43
44////////////////////////////////////////////////////////////////////////////////
45
46RooGaussModel::RooGaussModel(const char *name, const char *title, RooAbsRealLValue& xIn,
47 RooAbsReal& _mean, RooAbsReal& _sigma) :
48 RooResolutionModel(name,title,xIn),
49 _flatSFInt(false),
50 _asympInt(false),
51 mean("mean","Mean",this,_mean),
52 sigma("sigma","Width",this,_sigma),
53 msf("msf","Mean Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
54 ssf("ssf","Sigma Scale Factor",this,(RooRealVar&)RooRealConstant::value(1))
55{
56}
57
58////////////////////////////////////////////////////////////////////////////////
59
60RooGaussModel::RooGaussModel(const char *name, const char *title, RooAbsRealLValue& xIn,
61 RooAbsReal& _mean, RooAbsReal& _sigma,
62 RooAbsReal& _msSF) :
63 RooResolutionModel(name,title,xIn),
64 _flatSFInt(false),
65 _asympInt(false),
66 mean("mean","Mean",this,_mean),
67 sigma("sigma","Width",this,_sigma),
68 msf("msf","Mean Scale Factor",this,_msSF),
69 ssf("ssf","Sigma Scale Factor",this,_msSF)
70{
71}
72
73////////////////////////////////////////////////////////////////////////////////
74
75RooGaussModel::RooGaussModel(const char *name, const char *title, RooAbsRealLValue& xIn,
76 RooAbsReal& _mean, RooAbsReal& _sigma,
77 RooAbsReal& _meanSF, RooAbsReal& _sigmaSF) :
78 RooResolutionModel(name,title,xIn),
79 _flatSFInt(false),
80 _asympInt(false),
81 mean("mean","Mean",this,_mean),
82 sigma("sigma","Width",this,_sigma),
83 msf("msf","Mean Scale Factor",this,_meanSF),
84 ssf("ssf","Sigma Scale Factor",this,_sigmaSF)
85{
86}
87
88////////////////////////////////////////////////////////////////////////////////
89
92 _flatSFInt(other._flatSFInt),
93 _asympInt(other._asympInt),
94 mean("mean",this,other.mean),
95 sigma("sigma",this,other.sigma),
96 msf("msf",this,other.msf),
97 ssf("ssf",this,other.ssf)
98{
99}
100
101////////////////////////////////////////////////////////////////////////////////
102/// Destructor
103
105{
106}
107
108////////////////////////////////////////////////////////////////////////////////
109
111{
112 if (!TString("exp(-@0/@1)").CompareTo(name)) return expBasisPlus ;
113 if (!TString("exp(@0/@1)").CompareTo(name)) return expBasisMinus ;
114 if (!TString("exp(-abs(@0)/@1)").CompareTo(name)) return expBasisSum ;
115 if (!TString("exp(-@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisPlus ;
116 if (!TString("exp(@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisMinus ;
117 if (!TString("exp(-abs(@0)/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisSum ;
118 if (!TString("exp(-@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisPlus ;
119 if (!TString("exp(@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisMinus ;
120 if (!TString("exp(-abs(@0)/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisSum ;
121 if (!TString("(@0/@1)*exp(-@0/@1)").CompareTo(name)) return linBasisPlus ;
122 if (!TString("(@0/@1)*(@0/@1)*exp(-@0/@1)").CompareTo(name)) return quadBasisPlus ;
123 if (!TString("exp(-@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisPlus;
124 if (!TString("exp(@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisMinus;
125 if (!TString("exp(-abs(@0)/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisSum;
126 if (!TString("exp(-@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisPlus;
127 if (!TString("exp(@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisMinus;
128 if (!TString("exp(-abs(@0)/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisSum;
129 return 0 ;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
135{
136 auto arg1 = static_cast<RooAbsReal*>(basis().getParameter(1));
137 auto arg2 = static_cast<RooAbsReal*>(basis().getParameter(2));
138 double param1 = arg1 ? arg1->getVal() : 0.0;
139 double param2 = arg2 ? arg2->getVal() : 0.0;
140 return evaluate(x, mean * msf, sigma * ssf, param1, param2, _basisCode);
141}
142
144 RooFit::Detail::DataMap const &dataMap) const
145{
146 auto xVals = dataMap.at(x);
147 auto meanVals = dataMap.at(mean);
148 auto meanSfVals = dataMap.at(msf);
149 auto sigmaVals = dataMap.at(sigma);
150 auto sigmaSfVals = dataMap.at(ssf);
151
152 auto param1 = static_cast<RooAbsReal *>(basis().getParameter(1));
153 auto param2 = static_cast<RooAbsReal *>(basis().getParameter(2));
154 const double zeroVal = 0.0;
155 auto param1Vals = param1 ? dataMap.at(param1) : std::span<const double>{&zeroVal, 1};
156 auto param2Vals = param2 ? dataMap.at(param2) : std::span<const double>{&zeroVal, 1};
157
158 BasisType basisType = getBasisType(_basisCode);
159 double basisSign = _basisCode - 10 * (basisType - 1) - 2;
160
161 // We have an implementation also for CUDA right now only for the most used
162 // basis type, which is expBasis. If the need to support other basis types
163 // arises, they can be implemented following this example. Remember to also
164 // adapt RooGaussModel::canComputeBatchWithCuda().
165 if (basisType == expBasis) {
166 RooBatchCompute::ArgVector extraArgs{basisSign};
168 {xVals, meanVals, meanSfVals, sigmaVals, sigmaSfVals, param1Vals}, extraArgs);
169 return;
170 }
171
172 // For now, if the arrays don't have the expected input shape, fall back to the scalar mode
173 if (xVals.size() != size || meanVals.size() != 1 || meanSfVals.size() != 1 || sigmaVals.size() != 1 ||
174 sigmaSfVals.size() != 1 || param1Vals.size() != 1 || param2Vals.size() != 1) {
175 return RooAbsPdf::computeBatch(output, size, dataMap);
176 }
177
178 for (unsigned int i = 0; i < size; ++i) {
179 output[i] = evaluate(xVals[i], meanVals[0] * meanSfVals[0], sigmaVals[0] * sigmaSfVals[0], param1Vals[0],
180 param2Vals[0], _basisCode);
181 }
182}
183
184double RooGaussModel::evaluate(double x, double mean, double sigma, double param1, double param2, int basisCode)
185{
186 // *** 1st form: Straight Gaussian, used for unconvoluted PDF or expBasis with 0 lifetime ***
187 static double root2(std::sqrt(2.)) ;
188 static double root2pi(std::sqrt(2.*std::atan2(0.,-1.))) ;
189 static double rootpi(std::sqrt(std::atan2(0.,-1.))) ;
190
191 BasisType basisType = getBasisType(basisCode);
192 BasisSign basisSign = (BasisSign)( basisCode - 10*(basisType-1) - 2 ) ;
193
194 double tau = (basisCode!=noBasis) ? param1 : 0.0;
195 if (basisType == coshBasis && basisCode!=noBasis ) {
196 double dGamma = param2;
197 if (dGamma==0) basisType = expBasis;
198 }
199
200 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
201 double xprime = (x-mean)/sigma ;
202 double result = std::exp(-0.5*xprime*xprime)/(sigma*root2pi) ;
203 if (basisCode!=0 && basisSign==Both) result *= 2 ;
204 return result ;
205 }
206
207 // *** 2nd form: 0, used for sinBasis, linBasis, and quadBasis with tau=0 ***
208 if (tau==0) {
209 return 0. ;
210 }
211
212 // *** 3nd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
213 double omega = (basisType==sinBasis || basisType==cosBasis) ? param2 : 0 ;
214 double dgamma = (basisType==sinhBasis || basisType==coshBasis) ? param2 : 0 ;
215 double _x = omega *tau ;
216 double _y = tau*dgamma/2;
217 double xprime = (x-mean)/tau ;
218 double c = sigma/(root2*tau) ;
219 double u = xprime/(2*c) ;
220
221 if (basisType==expBasis || (basisType==cosBasis && _x==0.)) {
222 double result(0) ;
223 if (basisSign!=Minus) result += evalCerf(0,-u,c).real();
224 if (basisSign!=Plus) result += evalCerf(0, u,c).real();
225 return result ;
226 }
227
228 // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
229 if (basisType==sinBasis) {
230 double result(0) ;
231 if (_x==0.) return result ;
232 if (basisSign!=Minus) result += -evalCerf(-_x,-u,c).imag();
233 if (basisSign!=Plus) result += -evalCerf( _x, u,c).imag();
234 return result ;
235 }
236
237 // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
238 if (basisType==cosBasis) {
239 double result(0) ;
240 if (basisSign!=Minus) result += evalCerf(-_x,-u,c).real();
241 if (basisSign!=Plus) result += evalCerf( _x, u,c).real();
242 return result ;
243 }
244
245 // ***8th form: Convolution with exp(-|t|/tau)*cosh(dgamma*t/2), used for coshBasisSum ***
246 if (basisType==coshBasis || basisType ==sinhBasis) {
247 double result(0);
248 int sgn = ( basisType == coshBasis ? +1 : -1 );
249 if (basisSign!=Minus) result += 0.5*( evalCerf(0,-u,c*(1-_y)).real()+sgn*evalCerf(0,-u,c*(1+_y)).real()) ;
250 if (basisSign!=Plus) result += 0.5*(sgn*evalCerf(0, u,c*(1-_y)).real()+ evalCerf(0, u,c*(1+_y)).real()) ;
251 return result ;
252 }
253
254 // *** 6th form: Convolution with (t/tau)*exp(-t/tau), used for linBasis ***
255 if (basisType==linBasis) {
256 R__ASSERT(basisSign==Plus); // This should only be for positive times
257
258 double f0 = std::exp(-xprime+c*c) * RooMath::erfc(-u+c);
259 double f1 = std::exp(-u*u);
260 return (xprime - 2*c*c)*f0 + (2*c/rootpi)*f1 ;
261 }
262
263 // *** 7th form: Convolution with (t/tau)^2*exp(-t/tau), used for quadBasis ***
264 if (basisType==quadBasis) {
265 R__ASSERT(basisSign==Plus); // This should only be for positive times
266
267 double f0 = std::exp(-xprime+c*c) * RooMath::erfc(-u+c);
268 double f1 = std::exp(-u*u);
269 double x2c2 = xprime - 2*c*c;
270 return ( x2c2*x2c2*f0 + (2*c/rootpi)*x2c2*f1 + 2*c*c*f0 );
271 }
272
273 R__ASSERT(0) ;
274 return 0 ;
275}
276
277////////////////////////////////////////////////////////////////////////////////
278
279Int_t RooGaussModel::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
280{
281 switch(_basisCode) {
282
283 // Analytical integration capability of raw PDF
284 case noBasis:
285
286 // Optionally advertise flat integral over sigma scale factor
287 if (_flatSFInt) {
288 if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) return 2 ;
289 }
290
291 if (matchArgs(allVars,analVars,convVar())) return 1 ;
292 break ;
293
294 // Analytical integration capability of convoluted PDF
295 case expBasisPlus:
296 case expBasisMinus:
297 case expBasisSum:
298 case sinBasisPlus:
299 case sinBasisMinus:
300 case sinBasisSum:
301 case cosBasisPlus:
302 case cosBasisMinus:
303 case cosBasisSum:
304 case linBasisPlus:
305 case quadBasisPlus:
306 case coshBasisMinus:
307 case coshBasisPlus:
308 case coshBasisSum:
309 case sinhBasisMinus:
310 case sinhBasisPlus:
311 case sinhBasisSum:
312
313 // Optionally advertise flat integral over sigma scale factor
314 if (_flatSFInt) {
315
316 if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) {
317 return 2 ;
318 }
319 }
320
321 if (matchArgs(allVars,analVars,convVar())) return 1 ;
322 break ;
323 }
324
325 return 0 ;
326}
327
328////////////////////////////////////////////////////////////////////////////////
329
330double RooGaussModel::analyticalIntegral(Int_t code, const char* rangeName) const
331{
332 static const double root2 = std::sqrt(2.) ;
333 //static double rootPiBy2 = std::sqrt(std::atan2(0.0,-1.0)/2.0);
334 static const double rootpi = std::sqrt(std::atan2(0.0,-1.0));
335 double ssfInt(1.0) ;
336
337 // Code must be 1 or 2
338 R__ASSERT(code==1||code==2) ;
339 if (code==2) ssfInt = (ssf.max(rangeName)-ssf.min(rangeName)) ;
340
341 BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
342 BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
343
344 // *** 1st form: Straight Gaussian, used for unconvoluted PDF or expBasis with 0 lifetime ***
345 double tau = (_basisCode!=noBasis)?((RooAbsReal*)basis().getParameter(1))->getVal():0 ;
346 if (basisType == coshBasis && _basisCode!=noBasis ) {
347 double dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
348 if (dGamma==0) basisType = expBasis;
349 }
350 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
351 double xscale = root2*(sigma*ssf);
352 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 1st form" << endl ;
353
354 double xpmin = (x.min(rangeName)-(mean*msf))/xscale ;
355 double xpmax = (x.max(rangeName)-(mean*msf))/xscale ;
356
357 double result ;
358 if (_asympInt) { // modified FMV, 07/24/03
359 result = 1.0 ;
360 } else {
361 result = 0.5*(RooMath::erf(xpmax)-RooMath::erf(xpmin)) ;
362 }
363
364 if (_basisCode!=0 && basisSign==Both) result *= 2 ;
365 //cout << "Integral 1st form " << " result= " << result*ssfInt << endl;
366 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 1 " << endl; }
367 return result*ssfInt ;
368 }
369
370
371 double omega = ((basisType==sinBasis)||(basisType==cosBasis)) ? ((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
372 double dgamma =((basisType==sinhBasis)||(basisType==coshBasis)) ? ((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
373
374 // *** 2nd form: unity, used for sinBasis and linBasis with tau=0 (PDF is zero) ***
375 if (tau==0) {
376 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 2nd form" << endl ;
377 return 0. ;
378 }
379
380 // *** 3rd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
381 double c = (sigma*ssf)/(root2*tau) ;
382 double xpmin = (x.min(rangeName)-(mean*msf))/tau ;
383 double xpmax = (x.max(rangeName)-(mean*msf))/tau ;
384 double umin = xpmin/(2*c) ;
385 double umax = xpmax/(2*c) ;
386
387 if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
388 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 3d form tau=" << tau << endl ;
389 double result(0) ;
390 if (basisSign!=Minus) result += evalCerfInt(+1,0,tau,-umin,-umax,c).real();
391 if (basisSign!=Plus) result += evalCerfInt(-1,0,tau, umin, umax,c).real();
392 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 3 " << endl; }
393 return result*ssfInt ;
394 }
395
396 // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
397 double _x = omega * tau ;
398 double _y = tau*dgamma/2;
399
400 if (basisType==sinBasis) {
401 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 4th form omega = " << omega << ", tau = " << tau << endl ;
402 double result(0) ;
403 if (_x==0) return result*ssfInt ;
404 if (basisSign!=Minus) result += -1*evalCerfInt(+1,-_x,tau,-umin,-umax,c).imag();
405 if (basisSign!=Plus) result += -1*evalCerfInt(-1, _x,tau, umin, umax,c).imag();
406 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 4 " << endl; }
407 return result*ssfInt ;
408 }
409
410 // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
411 if (basisType==cosBasis) {
412 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
413 double result(0) ;
414 if (basisSign!=Minus) result += evalCerfInt(+1,-_x,tau,-umin,-umax,c).real();
415 if (basisSign!=Plus) result += evalCerfInt(-1, _x,tau, umin, umax,c).real();
416 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 5 " << endl; }
417 return result*ssfInt ;
418 }
419
420 // *** 8th form: Convolution with exp(-|t|/tau)*cosh(dgamma*t/2), used for coshBasis ***
421 // *** 9th form: Convolution with exp(-|t|/tau)*sinh(dgamma*t/2), used for sinhBasis ***
422 if (basisType==coshBasis || basisType == sinhBasis) {
423 if (verboseEval()>0) {cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 8th form tau=" << tau << endl ; }
424 double result(0) ;
425 int sgn = ( basisType == coshBasis ? +1 : -1 );
426 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());
427 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());
428 if (TMath::IsNaN(result)) { cxcoutE(Tracing) << "RooGaussModel::analyticalIntegral(" << GetName() << ") got nan during case 6 " << endl; }
429 return result*ssfInt ;
430 }
431
432 // *** 6th form: Convolution with (t/tau)*exp(-t/tau), used for linBasis ***
433 if (basisType==linBasis) {
434 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 6th form tau=" << tau << endl ;
435
436 double f0 = RooMath::erf(-umax) - RooMath::erf(-umin);
437 double f1 = std::exp(-umax*umax) - std::exp(-umin*umin);
438
439 double tmp1 = std::exp(-xpmax)*RooMath::erfc(-umax + c);
440 double tmp2 = std::exp(-xpmin)*RooMath::erfc(-umin + c);
441
442 double f2 = tmp1 - tmp2;
443 double f3 = xpmax*tmp1 - xpmin*tmp2;
444
445 double expc2 = std::exp(c*c);
446
447 return -tau*( f0 +
448 (2*c/rootpi)*f1 +
449 (1 - 2*c*c)*expc2*f2 +
450 expc2*f3
451 )*ssfInt;
452 }
453
454 // *** 7th form: Convolution with (t/tau)*(t/tau)*exp(-t/tau), used for quadBasis ***
455 if (basisType==quadBasis) {
456 if (verboseEval()>0) cout << "RooGaussModel::analyticalIntegral(" << GetName() << ") 7th form tau=" << tau << endl ;
457
458 double f0 = RooMath::erf(-umax) - RooMath::erf(-umin);
459
460 double tmpA1 = std::exp(-umax*umax);
461 double tmpA2 = std::exp(-umin*umin);
462
463 double f1 = tmpA1 - tmpA2;
464 double f2 = umax*tmpA1 - umin*tmpA2;
465
466 double tmpB1 = std::exp(-xpmax)*RooMath::erfc(-umax + c);
467 double tmpB2 = std::exp(-xpmin)*RooMath::erfc(-umin + c);
468
469 double f3 = tmpB1 - tmpB2;
470 double f4 = xpmax*tmpB1 - xpmin*tmpB2;
471 double f5 = xpmax*xpmax*tmpB1 - xpmin*xpmin*tmpB2;
472
473 double expc2 = std::exp(c*c);
474
475 return -tau*( 2*f0 +
476 (4*c/rootpi)*((1-c*c)*f1 + c*f2) +
477 (2*c*c*(2*c*c-1) + 2)*expc2*f3 - (4*c*c-2)*expc2*f4 + expc2*f5
478 )*ssfInt;
479 }
480 R__ASSERT(0) ;
481 return 0 ;
482}
483
484
485////////////////////////////////////////////////////////////////////////////////
486
487std::complex<double> RooGaussModel::evalCerfInt(double sign, double _x, double tau, double umin, double umax, double c) const
488{
489 std::complex<double> diff(2., 0.);
490 if (!_asympInt) {
491 diff = evalCerf(_x,umin,c);
492 diff -= evalCerf(_x,umax,c);
493 diff += RooMath::erf(umin) - RooMath::erf(umax);
494 diff *= sign;
495 }
496 diff *= std::complex<double>(1., _x);
497 diff *= tau / (1.+_x*_x);
498 return diff;
499}
500
501////////////////////////////////////////////////////////////////////////////////
502
503Int_t RooGaussModel::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
504{
505 if (matchArgs(directVars,generateVars,x)) return 1 ;
506 return 0 ;
507}
508
509////////////////////////////////////////////////////////////////////////////////
510
512{
513 R__ASSERT(code==1) ;
514 double xmin = x.min();
515 double xmax = x.max();
516 TRandom *generator = RooRandom::randomGenerator();
517 while(true) {
518 double xgen = generator->Gaus(mean*msf,sigma*ssf);
519 if (xgen<xmax && xgen>xmin) {
520 x = xgen ;
521 return ;
522 }
523 }
524}
#define c(i)
Definition RSha256.hxx:101
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define cxcoutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:118
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition DataMap.cxx:40
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:22
RooAbsArg * getParameter(const char *name) const
Return pointer to parameter with given name.
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
RooRealProxy sigma
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
std::complex< double > evalCerfInt(double sign, double wt, double tau, double umin, double umax, double c) const
RooRealProxy msf
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
~RooGaussModel() override
Destructor.
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
Int_t basisCode(const char *name) const override
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const override
Base function for computing multiple values of a RooAbsReal.
RooRealProxy mean
double analyticalIntegral(Int_t code, const char *rangeName) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy ssf
static BasisType getBasisType(int basisCode)
static std::complex< double > erfc(const std::complex< double > z)
complex erfc function
Definition RooMath.cxx:40
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition RooMath.cxx:59
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:51
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:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
Int_t _basisCode
Identifier code for selected basis function.
RooAbsRealLValue & convVar() const
Return the convolution variable of the resolution model.
const RooFormulaVar & basis() const
RooTemplateProxy< RooAbsRealLValue > x
Dependent/convolution variable.
double max(const char *rname=nullptr) 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.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const char * GetName() const override
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:274
Basic string class.
Definition TString.h:139
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
TF1 * f1
Definition legend1.C:11
std::vector< double > ArgVector
void compute(Config cfg, Computer comp, RestrictArr output, size_t size, const VarVector &vars, ArgVector &extraArgs)
__roohost__ __roodevice__ STD::complex< double > evalCerfApprox(double _x, double u, double c)
use the approximation: erf(z) = exp(-z*z)/(STD::sqrt(pi)*z) to explicitly cancel the divergent exp(y*...
__roohost__ __roodevice__ STD::complex< double > evalCerf(double swt, double u, double c)
Bool_t IsNaN(Double_t x)
Definition TMath.h:892
static void output()