59enum BasisType { none = 0, expBasis = 1, sinBasis = 2, cosBasis = 3, sinhBasis = 4, coshBasis = 5 };
61enum BasisSign { Both = 0, Plus = +1, Minus = -1 };
89 _mean(
"mean",
"Mean of Gaussian component", this, meanIn),
90 sigma(
"sigma",
"Width", this, sigmaIn),
91 rlife(
"rlife",
"Life time", this, rlifeIn),
92 _meanSF(
"meanSF",
"Scale factor for mean", this, meanSF),
93 ssf(
"ssf",
"Sigma Scale Factor", this, sigmaSF),
94 rsf(
"rsf",
"RLife Scale Factor", this, rlifeSF),
118 sigma(
"sigma",
"Width",this,_sigma),
119 rlife(
"rlife",
"Life time",this,_rlife),
123 _flip(
type==Flipped),_nlo(nlo), _flatSFInt(false), _asympInt(false)
144 sigma(
"sigma",
"Width",this,_sigma),
145 rlife(
"rlife",
"Life time",this,_rlife),
147 ssf(
"ssf",
"Sigma Scale Factor",this,_rsSF),
148 rsf(
"rsf",
"RLife Scale Factor",this,_rsSF),
149 _flip(
type==Flipped),
174 sigma(
"sigma",
"Width",this,_sigma),
175 rlife(
"rlife",
"Life time",this,_rlife),
177 ssf(
"ssf",
"Sigma Scale Factor",this,_sigmaSF),
178 rsf(
"rsf",
"RLife Scale Factor",this,_rlifeSF),
179 _flip(
type==Flipped),
190 _mean(
"mean", this, other._mean),
192 rlife(
"rlife",this,other.rlife),
193 _meanSF(
"meanSf", this, other._meanSF),
194 ssf(
"ssf",this,other.ssf),
195 rsf(
"rsf",this,other.rsf),
198 _flatSFInt(other._flatSFInt),
199 _asympInt(other._asympInt)
207 std::string str =
name;
210 str.erase(remove(str.begin(),str.end(),
' '),str.end());
212 if (str ==
"exp(-@0/@1)")
return expBasisPlus ;
213 if (str ==
"exp(@0/@1)")
return expBasisMinus ;
214 if (str ==
"exp(-abs(@0)/@1)")
return expBasisSum ;
215 if (str ==
"exp(-@0/@1)*sin(@0*@2)")
return sinBasisPlus ;
216 if (str ==
"exp(@0/@1)*sin(@0*@2)")
return sinBasisMinus ;
217 if (str ==
"exp(-abs(@0)/@1)*sin(@0*@2)")
return sinBasisSum ;
218 if (str ==
"exp(-@0/@1)*cos(@0*@2)")
return cosBasisPlus ;
219 if (str ==
"exp(@0/@1)*cos(@0*@2)")
return cosBasisMinus ;
220 if (str ==
"exp(-abs(@0)/@1)*cos(@0*@2)")
return cosBasisSum ;
221 if (str ==
"exp(-@0/@1)*sinh(@0*@2/2)")
return sinhBasisPlus;
222 if (str ==
"exp(@0/@1)*sinh(@0*@2/2)")
return sinhBasisMinus;
223 if (str ==
"exp(-abs(@0)/@1)*sinh(@0*@2/2)")
return sinhBasisSum;
224 if (str ==
"exp(-@0/@1)*cosh(@0*@2/2)")
return coshBasisPlus;
225 if (str ==
"exp(@0/@1)*cosh(@0*@2/2)")
return coshBasisMinus;
226 if (str ==
"exp(-abs(@0)/@1)*cosh(@0*@2/2)")
return coshBasisSum;
235double logErfC(
double xx)
245 (-z * z - 1.26551223 +
251 t * (-1.13520398 + t * (1.48851587 + t * (-0.82215223 + t * 0.17087277)))))))));
256 exp(-z * z - 1.26551223 +
261 t * (0.27886807 + t * (-1.13520398 +
262 t * (1.48851587 + t * (-0.82215223 + t * 0.17087277))))))))));
275 static double rootpi=
sqrt(
atan2(0.,-1.));
276 std::complex<double> z(swt*
c,u+
c);
277 std::complex<double> zc(u+
c,-swt*
c);
278 std::complex<double> zsq= z*z;
279 std::complex<double>
v= -zsq - u*u;
281 return std::exp(
v)*(-std::exp(zsq)/(zc*rootpi) + 1.)*2.;
286std::complex<double>
evalCerf(
double swt,
double u,
double c)
288 std::complex<double> z(swt*
c,u+
c);
295inline double evalCerfRe(
double u,
double c) {
296 double expArg = u*2*
c+
c*
c ;
300 return exp(expArg+logErfC(u+
c));
312 static double root2(sqrt(2.)) ;
315 BasisSign basisSign = (BasisSign)(
_basisCode - 10*(basisType-1) - 2 ) ;
317 double fsign =
_flip?-1:1 ;
324 if (basisType == coshBasis &&
_basisCode!=noBasis ) {
326 if (dGamma==0) basisType = expBasis;
330 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
331 if (
verboseEval()>2) std::cout <<
"RooGExpModel::evaluate(" <<
GetName() <<
") 1st form" << std::endl ;
333 double expArg = sig*sig/(2*rtau*rtau) + fsign*(
x -
_mean*
_meanSF)/rtau ;
341 result = 1/(2*rtau) * exp(expArg + logErfC(sig/(root2*rtau) + fsign*(
x -
_mean*
_meanSF)/(root2*sig))) ;
361 if (
verboseEval()>2) std::cout <<
"RooGExpModel::evaluate(" <<
GetName() <<
") 2nd form" << std::endl ;
368 if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
369 if (
verboseEval()>2) std::cout <<
"RooGExpModel::evaluate(" <<
GetName() <<
") 3d form tau=" << tau << std::endl ;
378 double wt = omega *tau ;
379 if (basisType==sinBasis) {
380 if (
verboseEval()>2) std::cout <<
"RooGExpModel::evaluate(" <<
GetName() <<
") 4th form omega = "
381 << omega <<
", tau = " << tau << std::endl ;
383 if (wt==0.)
return result ;
384 if (basisSign!=Minus)
result += -1*
calcSinConv(+1,sig,tau,omega,rtau,fsign).imag() ;
385 if (basisSign!=Plus)
result += -1*
calcSinConv(-1,sig,tau,omega,rtau,fsign).imag() ;
391 if (basisType==cosBasis) {
393 <<
") 5th form omega = " << omega <<
", tau = " << tau << std::endl ;
395 if (basisSign!=Minus)
result +=
calcSinConv(+1,sig,tau,omega,rtau,fsign).real() ;
403 if (basisType==sinhBasis) {
407 <<
") 6th form = " << dgamma <<
", tau = " << tau << std::endl;
412 double tau1 = 1/(1/tau-dgamma/2) ;
413 double tau2 = 1/(1/tau+dgamma/2) ;
423 if (basisType==coshBasis) {
427 <<
") 7th form = " << dgamma <<
", tau = " << tau << std::endl;
432 double tau1 = 1/(1/tau-dgamma/2) ;
433 double tau2 = 1/(1/tau+dgamma/2) ;
450 static double root2(sqrt(2.)) ;
454 double c1= sig/(root2*tau);
455 double u1=
s1/(2*
c1);
457 double c2= sig/(root2*rtau);
458 double u2= fsign*s2/(2*
c2) ;
461 std::complex<double> eins(1,0);
462 std::complex<double> k(1/tau,sign*omega);
465 return (evalCerf(-sign*omega*tau,u1,
c1)+std::complex<double>(evalCerfRe(u2,
c2),0)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
476 static double root2(sqrt(2.)) ;
480 double c1= sig/(root2*tau);
481 double u1=
s1/(2*
c1);
483 double c2= sig/(root2*rtau);
484 double u2= fsign*s2/(2*
c2) ;
489 return (evalCerfRe(u1,
c1)+evalCerfRe(u2,
c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
499 static double root2(sqrt(2.)) ;
500 static double root2pi(sqrt(2*atan2(0.,-1.))) ;
501 static double rootpi(sqrt(atan2(0.,-1.)));
504 double xp(
x - _mean*_meanSF) ;
513 if ((sign<0)&&(std::abs(tau-rtau)<tau/260)) {
515 double MeanTau=0.5*(tau+rtau);
516 if (std::abs(xp/MeanTau)>300) {
520 cFly=1./(MeanTau*MeanTau*root2pi) *
521 exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
522 *(sig*exp(-1/(2*sig*sig)*std::pow((sig*sig/MeanTau+xp),2))
523 -(sig*sig/MeanTau+xp)*(rootpi/root2)*
RooMath::erfc(sig/(root2*MeanTau)+xp/(root2*sig)));
526 double epsilon=0.5*(tau-rtau);
527 double a=sig/(root2*MeanTau)+xp/(root2*sig);
528 cFly += 1./(MeanTau*MeanTau)
529 *exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
530 *0.5/MeanTau*epsilon*epsilon*
531 (exp(-
a*
a)*(sig/MeanTau*root2/rootpi
532 -(4*
a*sig*sig)/(2*rootpi*MeanTau*MeanTau)
533 +(-4/rootpi+8*
a*
a/rootpi)/6
534 *std::pow(sig/(root2*MeanTau),3)
535 +2/rootpi*(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
536 (sig/(root2*MeanTau)-
a*(sig*sig)/(2*MeanTau*MeanTau))
537 +2/rootpi*((3*sig*sig)/(2*MeanTau*MeanTau)+xp/MeanTau+
538 0.5*std::pow(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,2))*sig/(root2*MeanTau))
539 -(2*sig*sig/(MeanTau*MeanTau)+xp/MeanTau+(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
540 (3*sig*sig/(2*MeanTau*MeanTau)+xp/MeanTau)
541 +std::pow(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,3)/6)*
RooMath::erfc(
a));
546 double expArg1 = sig*sig/(2*tau*tau)-sign*xp/tau ;
547 double expArg2 = sig*sig/(2*rtau*rtau)+xp/rtau ;
552 term1 = exp(expArg1) *
RooMath::erfc(sig/(root2*tau)-sign*xp/(root2*sig)) ;
554 term1 = exp(expArg1+logErfC(sig/(root2*tau)-sign*xp/(root2*sig))) ;
557 term2 = exp(expArg2) *
RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)) ;
559 term2 = exp(expArg2+logErfC(sig/(root2*rtau)+xp/(root2*sig))) ;
562 cFly=(term1+sign*term2)/(2*(tau+sign*rtau));
580double RooGExpModel::calcCoshConv(double sign, double tau, double dgamma, double sig, double rtau, double fsign) const
584 static double root2(sqrt(2.)) ;
585 static double root2pi(sqrt(2*atan2(0.,-1.))) ;
586 static double rootpi(sqrt(atan2(0.,-1.)));
587 double tau1 = 1/(1/tau-dgamma/2);
588 double tau2 = 1/(1/tau+dgamma/2);
596 xp *= fsign ; // modified FMV,08/13/03
597 sign *= fsign ; // modified FMV,08/13/03
599 cFly=tau1*(exp(sig*sig/(2*tau1*tau1)-sign*xp/tau1)
600 *RooMath::erfc(sig/(root2*tau1)-sign*xp/(root2*sig))
601 +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
602 *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau1+sign*rtau))
603 +tau2*(exp(sig*sig/(2*tau2*tau2)-sign*xp/tau2)
604 *RooMath::erfc(sig/(root2*tau2)-sign*xp/(root2*sig))
605 +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
606 *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau2+sign*rtau));
615double RooGExpModel::calcSinhConv(double sign, double sign1, double sign2, double tau, double dgamma, double sig, double rtau, double fsign) const
617 static double root2(sqrt(2.)) ;
618 static double root2pi(sqrt(2*atan2(0.,-1.))) ;
619 static double rootpi(sqrt(atan2(0.,-1.)));
620 double tau1 = 1/(1/tau-dgamma/2);
621 double tau2 = 1/(1/tau+dgamma/2);
630 xp *= fsign ; // modified FMV,08/13/03
631 sign1 *= fsign ; // modified FMV,08/13/03
632 sign2 *= fsign ; // modified FMV,08/13/03
634 cFly=sign1*tau1*(exp(sig*sig/(2*tau1*tau1)-sign*xp/tau1)
635 *RooMath::erfc(sig/(root2*tau1)-sign*xp/(root2*sig))
636 +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
637 *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau1+sign*rtau))
638 +sign2*tau2*(exp(sig*sig/(2*tau2*tau2)-sign*xp/tau2)
639 *RooMath::erfc(sig/(root2*tau2)-sign*xp/(root2*sig))
640 +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
641 *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau2+sign*rtau));
692 static double root2 = sqrt(2.) ;
703 BasisSign basisSign = (BasisSign)(
_basisCode - 10*(basisType-1) - 2 ) ;
708 if (basisType == coshBasis &&
_basisCode!=noBasis ) {
710 if (dGamma==0) basisType = expBasis;
712 double fsign =
_flip?-1:1 ;
717 if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
718 if (
verboseEval()>0) std::cout <<
"RooGExpModel::analyticalIntegral(" <<
GetName() <<
") 1st form" << std::endl ;
724 double c = sig/(root2*rtau) ;
725 double umin = xpmin/(2*
c) ;
726 double umax = xpmax/(2*
c) ;
739 double omega = (basisType!=expBasis) ?(
static_cast<RooAbsReal*
>(
basis().getParameter(2)))->
getVal() : 0 ;
744 if (
verboseEval()>0) std::cout <<
"RooGExpModel::analyticalIntegral(" <<
GetName() <<
") 2nd form" << std::endl ;
749 if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
761 double wt = omega * tau ;
762 if (basisType==sinBasis) {
763 if (
verboseEval()>0) std::cout <<
"RooGExpModel::analyticalIntegral(" <<
GetName() <<
") 4th form omega = "
764 << omega <<
", tau = " << tau << std::endl ;
767 if (wt==0)
return result ;
778 if (basisType==cosBasis) {
780 <<
") 5th form omega = " << omega <<
", tau = " << tau << std::endl ;
795 if (basisType==sinhBasis) {
797 <<
") 6th form dgamma = " << dgamma <<
", tau = " << tau << std::endl ;
798 double tau1 = 1/(1/tau-dgamma/2);
799 double tau2 = 1/(1/tau+dgamma/2);
814 if (basisType==coshBasis) {
816 <<
") 6th form dgamma = " << dgamma <<
", tau = " << tau << std::endl ;
818 double tau1 = 1/(1/tau-dgamma/2);
819 double tau2 = 1/(1/tau+dgamma/2);
845 double sig,
double rtau,
double fsign,
const char* rangeName)
const
847 static double root2(sqrt(2.)) ;
851 double c1= sig/(root2*tau);
852 double umin1= smin1/(2*
c1);
853 double umax1= smax1/(2*
c1);
856 double c2= sig/(root2*rtau);
857 double umin2= smin2/(2*
c2) ;
858 double umax2= smax2/(2*
c2) ;
860 std::complex<double> eins(1,0);
861 std::complex<double> k(1/tau,sign*omega);
862 std::complex<double> term1 =
evalCerfInt(sign,-sign*omega*tau, tau, -sign*umin1, -sign*umax1,
c1);
864 std::complex<double> term2 = std::complex<double>(
evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2,
c2)*fsign*sign,0);
865 return (term1+term2)/(eins + k*fsign*sign*rtau) ;
874 static double root2(sqrt(2.)) ;
878 double c1= sig/(root2*tau);
879 double umin1= smin1/(2*
c1);
880 double umax1= smax1/(2*
c1);
883 double c2= sig/(root2*rtau);
884 double umin2= smin2/(2*
c2) ;
885 double umax2= smax2/(2*
c2) ;
889 double term1 =
evalCerfInt(sign, tau, -sign*umin1, -sign*umax1,
c1);
890 double term2 =
evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2,
c2)*fsign*sign;
893 if (std::abs(tau-rtau)<1
e-10 && std::abs(term1+term2)<1
e-10) {
894 std::cout <<
"epsilon method" << std::endl ;
895 static double epsilon = 1
e-4 ;
896 return calcSinConvNorm(sign,tau+epsilon,sig,rtau-epsilon,fsign,rangeName) ;
898 return (term1+term2)/(eins + k*fsign*sign*rtau) ;
906 std::complex<double> diff;
908 diff = std::complex<double>(2,0) ;
912 return std::complex<double>(tau/(1.+wt*wt),0)*std::complex<double>(1,wt)*diff;
924 if ((umin<-8 && umax>8)||(umax<-8 && umin>8)) {
938 if (
matchArgs(directVars,generateVars,
x))
return 1 ;
952 xgen = xgau + (
rlife*
rsf)*log(xexp);
954 xgen = xgau - (
rlife *
rsf) * log(xexp);
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
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
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 Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
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 ...
Abstract base class for objects that represent a real value and implements functionality common to al...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
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.
The RooGExpModel is a RooResolutionModel implementation that models a resolution function that is the...
std::complex< double > evalCerfInt(double sign, double wt, double tau, double umin, double umax, double c) const
std::complex< double > calcSinConv(double sign, double sig, double tau, double omega, double rtau, double fsign) const
double calcDecayConv(double sign, double tau, double sig, double rtau, double fsign) const
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
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...
std::complex< double > calcSinConvNorm(double sign, double tau, double omega, double sig, double rtau, double fsign, const char *rangeName) const
old code (asymptotic normalization only) std::complex<double> z(1/tau,sign*omega); return z*2/(omega*...
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Int_t basisCode(const char *name) const override
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 std::complex< double > faddeeva_fast(std::complex< double > z)
evaluate Faddeeva function for complex argument (fast version)
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Provides static functions to create and keep track of RooRealVar constants.
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.
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...
RVec< PromoteType< T > > log(const RVec< T > &v)
RVec< PromoteTypes< T0, T1 > > atan2(const T0 &x, const RVec< T1 > &v)
RVec< PromoteType< T > > exp(const RVec< T > &v)
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
__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)