73enum BasisSign { Both = 0, Plus = +1, Minus = -1 };
75BasisType getBasisType(
int basisCode)
77 return static_cast<BasisType
>(basisCode == 0 ? 0 : (basisCode / 10) + 1);
110 mean(
"mean",
"Mean",this,_mean),
111 sigma(
"sigma",
"Width",this,_sigma),
112 msf(
"msf",
"Mean Scale Factor",this,_meanSF),
113 ssf(
"ssf",
"Sigma Scale Factor",this,_sigmaSF)
125 msf(
"msf",this,other.
msf),
134 std::string str =
name;
137 str.erase(remove(str.begin(),str.end(),
' '),str.end());
139 if (str ==
"exp(-@0/@1)")
return expBasisPlus ;
140 if (str ==
"exp(@0/@1)")
return expBasisMinus ;
141 if (str ==
"exp(-abs(@0)/@1)")
return expBasisSum ;
142 if (str ==
"exp(-@0/@1)*sin(@0*@2)")
return sinBasisPlus ;
143 if (str ==
"exp(@0/@1)*sin(@0*@2)")
return sinBasisMinus ;
144 if (str ==
"exp(-abs(@0)/@1)*sin(@0*@2)")
return sinBasisSum ;
145 if (str ==
"exp(-@0/@1)*cos(@0*@2)")
return cosBasisPlus ;
146 if (str ==
"exp(@0/@1)*cos(@0*@2)")
return cosBasisMinus ;
147 if (str ==
"exp(-abs(@0)/@1)*cos(@0*@2)")
return cosBasisSum ;
148 if (str ==
"(@0/@1)*exp(-@0/@1)")
return linBasisPlus ;
149 if (str ==
"(@0/@1)*(@0/@1)*exp(-@0/@1)")
return quadBasisPlus ;
150 if (str ==
"exp(-@0/@1)*cosh(@0*@2/2)")
return coshBasisPlus;
151 if (str ==
"exp(@0/@1)*cosh(@0*@2/2)")
return coshBasisMinus;
152 if (str ==
"exp(-abs(@0)/@1)*cosh(@0*@2/2)")
return coshBasisSum;
153 if (str ==
"exp(-@0/@1)*sinh(@0*@2/2)")
return sinhBasisPlus;
154 if (str ==
"exp(@0/@1)*sinh(@0*@2/2)")
return sinhBasisMinus;
155 if (str ==
"exp(-abs(@0)/@1)*sinh(@0*@2/2)")
return sinhBasisSum;
166 double param1 = arg1 ? arg1->getVal() : 0.0;
167 double param2 = arg2 ? arg2->getVal() : 0.0;
173 std::span<double> output = ctx.
output();
174 std::size_t
size = output.size();
176 auto xVals = ctx.
at(
x);
177 auto meanVals = ctx.
at(
mean);
178 auto meanSfVals = ctx.
at(
msf);
179 auto sigmaVals = ctx.
at(
sigma);
180 auto sigmaSfVals = ctx.
at(
ssf);
184 const double zeroVal = 0.0;
185 auto param1Vals = param1 ? ctx.
at(param1) : std::span<const double>{&zeroVal, 1};
186 auto param2Vals = param2 ? ctx.
at(param2) : std::span<const double>{&zeroVal, 1};
188 BasisType basisType = getBasisType(
_basisCode);
189 double basisSign =
_basisCode - 10 * (basisType - 1) - 2;
195 if (basisType == expBasis) {
196 std::array<double, 1> extraArgs{basisSign};
198 {xVals, meanVals, meanSfVals, sigmaVals, sigmaSfVals, param1Vals}, extraArgs);
203 if (xVals.size() !=
size || meanVals.size() != 1 || meanSfVals.size() != 1 || sigmaVals.size() != 1 ||
204 sigmaSfVals.size() != 1 || param1Vals.size() != 1 || param2Vals.size() != 1) {
208 for (
unsigned int i = 0; i <
size; ++i) {
209 output[i] =
evaluate(xVals[i], meanVals[0] * meanSfVals[0], sigmaVals[0] * sigmaSfVals[0], param1Vals[0],
217 static double root2(std::sqrt(2.)) ;
218 static double root2pi(std::sqrt(2.*std::atan2(0.,-1.))) ;
219 static double rootpi(std::sqrt(std::atan2(0.,-1.))) ;
221 BasisType basisType = getBasisType(
basisCode);
222 BasisSign basisSign = (BasisSign)(
basisCode - 10*(basisType-1) - 2 ) ;
224 double tau = (
basisCode!=noBasis) ? param1 : 0.0;
225 if (basisType == coshBasis &&
basisCode!=noBasis ) {
226 double dGamma = param2;
227 if (dGamma==0) basisType = expBasis;
230 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
232 double result = std::exp(-0.5*xprime*xprime)/(
sigma*root2pi) ;
233 if (
basisCode!=0 && basisSign==Both) result *= 2 ;
243 double omega = (basisType==sinBasis || basisType==cosBasis) ? param2 : 0 ;
244 double dgamma = (basisType==sinhBasis || basisType==coshBasis) ? param2 : 0 ;
245 double _x = omega *tau ;
246 double _y = tau*dgamma/2;
247 double xprime = (
x-
mean)/tau ;
248 double c =
sigma/(root2*tau) ;
249 double u = xprime/(2*
c) ;
251 if (basisType==expBasis || (basisType==cosBasis && _x==0.)) {
253 if (basisSign!=Minus) result +=
evalCerf(0,-u,
c).real();
254 if (basisSign!=Plus) result +=
evalCerf(0, u,
c).real();
259 if (basisType==sinBasis) {
261 if (_x==0.)
return result ;
262 if (basisSign!=Minus) result += -
evalCerf(-_x,-u,
c).imag();
263 if (basisSign!=Plus) result += -
evalCerf( _x, u,
c).imag();
268 if (basisType==cosBasis) {
270 if (basisSign!=Minus) result +=
evalCerf(-_x,-u,
c).real();
271 if (basisSign!=Plus) result +=
evalCerf( _x, u,
c).real();
276 if (basisType==coshBasis || basisType ==sinhBasis) {
278 int sgn = ( basisType == coshBasis ? +1 : -1 );
279 if (basisSign!=Minus) result += 0.5*(
evalCerf(0,-u,
c*(1-_y)).real()+sgn*
evalCerf(0,-u,
c*(1+_y)).real()) ;
280 if (basisSign!=Plus) result += 0.5*(sgn*
evalCerf(0, u,
c*(1-_y)).real()+
evalCerf(0, u,
c*(1+_y)).real()) ;
285 if (basisType==linBasis) {
289 double f1 = std::exp(-u*u);
290 return (xprime - 2*
c*
c)*f0 + (2*
c/rootpi)*
f1 ;
294 if (basisType==quadBasis) {
298 double f1 = std::exp(-u*u);
299 double x2c2 = xprime - 2*
c*
c;
300 return ( x2c2*x2c2*f0 + (2*
c/rootpi)*x2c2*
f1 + 2*
c*
c*f0 );
362 static const double root2 = std::sqrt(2.) ;
364 static const double rootpi = std::sqrt(std::atan2(0.0,-1.0));
369 if (code==2) ssfInt = (
ssf.max(rangeName)-
ssf.min(rangeName)) ;
372 BasisSign basisSign = (BasisSign)(
_basisCode - 10*(basisType-1) - 2 ) ;
376 if (basisType == coshBasis &&
_basisCode!=noBasis ) {
378 if (dGamma==0) basisType = expBasis;
380 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
382 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 1st form" << std::endl ;
384 double xpmin = (
x.min(rangeName)-(
mean*
msf))/xscale ;
385 double xpmax = (
x.max(rangeName)-(
mean*
msf))/xscale ;
394 if (
_basisCode!=0 && basisSign==Both) result *= 2 ;
396 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 1 " << std::endl; }
397 return result*ssfInt ;
401 double omega = ((basisType==sinBasis)||(basisType==cosBasis)) ? (
static_cast<RooAbsReal*
>(
basis().getParameter(2)))->getVal() : 0 ;
402 double dgamma =((basisType==sinhBasis)||(basisType==coshBasis)) ? (
static_cast<RooAbsReal*
>(
basis().getParameter(2)))->getVal() : 0 ;
406 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 2nd form" << std::endl ;
412 double xpmin = (
x.min(rangeName)-(
mean*
msf))/tau ;
413 double xpmax = (
x.max(rangeName)-(
mean*
msf))/tau ;
414 double umin = xpmin/(2*
c) ;
415 double umax = xpmax/(2*
c) ;
417 if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
418 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 3d form tau=" << tau << std::endl ;
420 if (basisSign!=Minus) result +=
evalCerfInt(+1,0,tau,-umin,-umax,
c).real();
421 if (basisSign!=Plus) result +=
evalCerfInt(-1,0,tau, umin, umax,
c).real();
422 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 3 " << std::endl; }
423 return result*ssfInt ;
427 double _x = omega * tau ;
428 double _y = tau*dgamma/2;
430 if (basisType==sinBasis) {
431 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 4th form omega = " << omega <<
", tau = " << tau << std::endl ;
433 if (_x==0)
return result*ssfInt ;
434 if (basisSign!=Minus) result += -1*
evalCerfInt(+1,-_x,tau,-umin,-umax,
c).imag();
435 if (basisSign!=Plus) result += -1*
evalCerfInt(-1, _x,tau, umin, umax,
c).imag();
436 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 4 " << std::endl; }
437 return result*ssfInt ;
441 if (basisType==cosBasis) {
442 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 5th form omega = " << omega <<
", tau = " << tau << std::endl ;
444 if (basisSign!=Minus) result +=
evalCerfInt(+1,-_x,tau,-umin,-umax,
c).real();
445 if (basisSign!=Plus) result +=
evalCerfInt(-1, _x,tau, umin, umax,
c).real();
446 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 5 " << std::endl; }
447 return result*ssfInt ;
452 if (basisType==coshBasis || basisType == sinhBasis) {
453 if (
verboseEval()>0) {std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 8th form tau=" << tau << std::endl ; }
455 int sgn = ( basisType == coshBasis ? +1 : -1 );
456 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());
457 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());
458 if (
TMath::IsNaN(result)) {
cxcoutE(Tracing) <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") got nan during case 6 " << std::endl; }
459 return result*ssfInt ;
463 if (basisType==linBasis) {
464 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 6th form tau=" << tau << std::endl ;
467 double f1 = std::exp(-umax*umax) - std::exp(-umin*umin);
472 double f2 = tmp1 - tmp2;
473 double f3 = xpmax*tmp1 - xpmin*tmp2;
475 double expc2 = std::exp(
c*
c);
479 (1 - 2*
c*
c)*expc2*f2 +
485 if (basisType==quadBasis) {
486 if (
verboseEval()>0) std::cout <<
"RooGaussModel::analyticalIntegral(" <<
GetName() <<
") 7th form tau=" << tau << std::endl ;
490 double tmpA1 = std::exp(-umax*umax);
491 double tmpA2 = std::exp(-umin*umin);
493 double f1 = tmpA1 - tmpA2;
494 double f2 = umax*tmpA1 - umin*tmpA2;
499 double f3 = tmpB1 - tmpB2;
500 double f4 = xpmax*tmpB1 - xpmin*tmpB2;
501 double f5 = xpmax*xpmax*tmpB1 - xpmin*xpmin*tmpB2;
503 double expc2 = std::exp(
c*
c);
506 (4*
c/rootpi)*((1-
c*
c)*
f1 +
c*f2) +
507 (2*
c*
c*(2*
c*
c-1) + 2)*expc2*f3 - (4*
c*
c-2)*expc2*f4 + expc2*f5
519 std::complex<double> diff(2., 0.);
526 diff *= std::complex<double>(1., _x);
527 diff *= tau / (1.+_x*_x);
535 return matchArgs(directVars,generateVars,
x) ? 1 : 0;
543 double xmin =
x.min();
544 double xmax =
x.max();
548 if (xgen<xmax && xgen>
xmin) {
STD::complex< double > evalCerf(double swt, double u, double c)
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Signed integer 4 bytes (int).
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
virtual void doEval(RooFit::EvalContext &) const
Base function for computing multiple values of a RooAbsReal.
bool matchArgs(const RooArgSet &allDeps, RooArgSet &analDeps, const RooArgProxy &a, const Proxies &... proxies) const
RooArgSet is a container object that can hold multiple RooAbsArg objects.
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
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
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
bool canComputeBatchWithCuda() const override
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...
double analyticalIntegral(Int_t code, const char *rangeName) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
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.
Provides static functions to create and keep track of RooRealVar constants.
Int_t _basisCode
Identifier code for selected basis function.
RooAbsRealLValue & convVar() const
Return the convolution variable of the resolution model.
RooResolutionModel()=default
const RooFormulaVar & basis() const
RooTemplateProxy< RooAbsRealLValue > x
Dependent/convolution variable.
const char * GetName() const override
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...
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
STD::complex< double > evalCerf(double swt, double u, double c)
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*...