46 mean(
"mean",
"Mean",this,_mean),
47 sigma(
"sigma",
"Width",this,_sigma),
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)
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)
87 _flatSFInt(other._flatSFInt),
88 _asympInt(other._asympInt),
89 mean(
"mean",this,other.mean),
91 msf(
"msf",this,other.msf),
92 ssf(
"ssf",this,other.ssf)
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.))) ;
142 if (dGamma==0) basisType =
expBasis;
147 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 1st form" << endl ;
156 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 2nd form" << endl ;
170 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 3d form tau=" << tau << endl ;
174 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::getVal(" <<
GetName() <<
") got nan during case 1 " << endl; }
180 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 4th form omega = " << omega <<
", tau = " << tau << endl ;
182 if (_x==0.)
return result ;
191 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 5th form omega = " << omega <<
", tau = " << tau << endl ;
201 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 8th form tau = " << tau << endl ;
203 int sgn = ( basisType ==
coshBasis ? +1 : -1 );
212 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 6th form tau = " << tau << endl ;
217 return (xprime - 2*
c*
c)*f0 + (2*
c/rootpi)*
f1 ;
222 if (
verboseEval()>2) cout <<
"RooGaussModel::evaluate(" <<
GetName() <<
") 7th form tau = " << tau << endl ;
228 return ( x2c2*x2c2*f0 + (2*
c/rootpi)*x2c2*
f1 + 2*
c*
c*f0 );
290 static const Double_t root2 = std::sqrt(2.) ;
292 static const Double_t rootpi = std::sqrt(std::atan2(0.0,-1.0));
297 if (code==2) ssfInt = (
ssf.
max(rangeName)-
ssf.
min(rangeName)) ;
306 if (dGamma==0) basisType =
expBasis;
310 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 1st form" << endl ;
324 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 1 " << endl; }
325 return result*ssfInt ;
334 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 2nd form" << endl ;
346 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 3d form tau=" << tau << endl ;
350 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 3 " << endl; }
351 return result*ssfInt ;
359 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 4th form omega = " << omega <<
", tau = " << tau << endl ;
361 if (_x==0)
return result*ssfInt ;
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 ;
370 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 5th form omega = " << omega <<
", tau = " << tau << endl ;
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 ;
381 if (
verboseEval()>0) {cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 8th form tau=" << tau << endl ; }
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 ;
392 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 6th form tau=" << tau << endl ;
395 Double_t f1 = std::exp(-umax*umax) - std::exp(-umin*umin);
401 Double_t f3 = xpmax*tmp1 - xpmin*tmp2;
407 (1 - 2*
c*
c)*expc2*f2 +
414 if (
verboseEval()>0) cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 7th form tau=" << tau << endl ;
418 Double_t tmpA1 = std::exp(-umax*umax);
419 Double_t tmpA2 = std::exp(-umin*umin);
422 Double_t f2 = umax*tmpA1 - umin*tmpA2;
428 Double_t f4 = xpmax*tmpB1 - xpmin*tmpB2;
429 Double_t f5 = xpmax*xpmax*tmpB1 - xpmin*xpmin*tmpB2;
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
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);
458 return 2. * (ev * (mez2zcrootpi + 1.));
465 std::complex<Double_t> diff(2., 0.);
472 diff *= std::complex<Double_t>(1., _x);
473 diff *= tau / (1.+_x*_x);
481 if (
matchArgs(directVars,generateVars,
x))
return 1 ;
495 if (xgen<xmax && xgen>
xmin) {
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...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
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.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
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)
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.
static std::complex< double > erfc(const std::complex< double > z)
complex erfc function
static std::complex< double > erf(const std::complex< double > z)
complex erf function
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
RooRealConstant provides static functions to create and keep track of RooRealVar constants.
RooRealVar represents a variable that can be changed from the outside.
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.
This is the base class for the ROOT Random number generators.
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...