// Class RooGExpModel is a RooResolutionModel implementation that models
// a resolution function that is the convolution of a Gaussian with
// a one-sided exponential. Object of class RooGExpModel can be used
// for analytical convolutions with classes inheriting from RooAbsAnaConvPdf
// END_HTML
#include "RooFit.h"
#include "Riostream.h"
#include "Riostream.h"
#include "RooGExpModel.h"
#include "RooMath.h"
#include "RooRealConstant.h"
#include "RooRandom.h"
#include "RooMath.h"
#include "TMath.h"
#include "TError.h"
using namespace std;
ClassImp(RooGExpModel)
;
RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn,
RooAbsReal& _sigma, RooAbsReal& _rlife,
Bool_t nlo, Type type) :
RooResolutionModel(name,title,xIn),
sigma("sigma","Width",this,_sigma),
rlife("rlife","Life time",this,_rlife),
ssf("ssf","Sigma Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
rsf("rsf","RLife Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
_flip(type==Flipped),_nlo(nlo), _flatSFInt(kFALSE), _asympInt(kFALSE)
{
}
RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn,
RooAbsReal& _sigma, RooAbsReal& _rlife,
RooAbsReal& _rsSF,
Bool_t nlo, Type type) :
RooResolutionModel(name,title,xIn),
sigma("sigma","Width",this,_sigma),
rlife("rlife","Life time",this,_rlife),
ssf("ssf","Sigma Scale Factor",this,_rsSF),
rsf("rsf","RLife Scale Factor",this,_rsSF),
_flip(type==Flipped),
_nlo(nlo),
_flatSFInt(kFALSE),
_asympInt(kFALSE)
{
}
RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn,
RooAbsReal& _sigma, RooAbsReal& _rlife,
RooAbsReal& _sigmaSF, RooAbsReal& _rlifeSF,
Bool_t nlo, Type type) :
RooResolutionModel(name,title,xIn),
sigma("sigma","Width",this,_sigma),
rlife("rlife","Life time",this,_rlife),
ssf("ssf","Sigma Scale Factor",this,_sigmaSF),
rsf("rsf","RLife Scale Factor",this,_rlifeSF),
_flip(type==Flipped),
_nlo(nlo),
_flatSFInt(kFALSE),
_asympInt(kFALSE)
{
}
RooGExpModel::RooGExpModel(const RooGExpModel& other, const char* name) :
RooResolutionModel(other,name),
sigma("sigma",this,other.sigma),
rlife("rlife",this,other.rlife),
ssf("ssf",this,other.ssf),
rsf("rsf",this,other.rsf),
_flip(other._flip),
_nlo(other._nlo),
_flatSFInt(other._flatSFInt),
_asympInt(other._asympInt)
{
}
RooGExpModel::~RooGExpModel()
{
}
Int_t RooGExpModel::basisCode(const char* name) const
{
if (!TString("exp(-@0/@1)").CompareTo(name)) return expBasisPlus ;
if (!TString("exp(@0/@1)").CompareTo(name)) return expBasisMinus ;
if (!TString("exp(-abs(@0)/@1)").CompareTo(name)) return expBasisSum ;
if (!TString("exp(-@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisPlus ;
if (!TString("exp(@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisMinus ;
if (!TString("exp(-abs(@0)/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisSum ;
if (!TString("exp(-@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisPlus ;
if (!TString("exp(@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisMinus ;
if (!TString("exp(-abs(@0)/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisSum ;
if (!TString("exp(-@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisPlus;
if (!TString("exp(@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisMinus;
if (!TString("exp(-abs(@0)/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisSum;
if (!TString("exp(-@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisPlus;
if (!TString("exp(@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisMinus;
if (!TString("exp(-abs(@0)/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisSum;
return 0 ;
}
Double_t RooGExpModel::evaluate() const
{
static Double_t root2(sqrt(2.)) ;
BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
Double_t fsign = _flip?-1:1 ;
Double_t sig = sigma*ssf ;
Double_t rtau = rlife*rsf ;
Double_t tau = (_basisCode!=noBasis)?((RooAbsReal*)basis().getParameter(1))->getVal():0. ;
if (basisType == coshBasis && _basisCode!=noBasis ) {
Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
if (dGamma==0) basisType = expBasis;
}
if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 1st form" << endl ;
Double_t expArg = sig*sig/(2*rtau*rtau) + fsign*x/rtau ;
Double_t result ;
if (expArg<300) {
result = 1/(2*rtau) * exp(expArg) * RooMath::erfc(sig/(root2*rtau) + fsign*x/(root2*sig));
} else {
result = 1/(2*rtau) * exp(expArg + logErfC(sig/(root2*rtau) + fsign*x/(root2*sig))) ;
}
if (_basisCode!=0 && basisSign==Both) result *= 2 ;
return result ;
}
if (tau==0) {
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 2nd form" << endl ;
return 0. ;
}
Double_t omega = (basisType!=expBasis)?((RooAbsReal*)basis().getParameter(2))->getVal():0. ;
if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 3d form tau=" << tau << endl ;
Double_t result(0) ;
if (basisSign!=Minus) result += calcDecayConv(+1,tau,sig,rtau,fsign) ;
if (basisSign!=Plus) result += calcDecayConv(-1,tau,sig,rtau,fsign) ;
return result ;
}
Double_t wt = omega *tau ;
if (basisType==sinBasis) {
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 4th form omega = "
<< omega << ", tau = " << tau << endl ;
Double_t result(0) ;
if (wt==0.) return result ;
if (basisSign!=Minus) result += -1*calcSinConv(+1,sig,tau,omega,rtau,fsign).imag() ;
if (basisSign!=Plus) result += -1*calcSinConv(-1,sig,tau,omega,rtau,fsign).imag() ;
return result ;
}
if (basisType==cosBasis) {
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
<< ") 5th form omega = " << omega << ", tau = " << tau << endl ;
Double_t result(0) ;
if (basisSign!=Minus) result += calcSinConv(+1,sig,tau,omega,rtau,fsign).real() ;
if (basisSign!=Plus) result += calcSinConv(-1,sig,tau,omega,rtau,fsign).real() ;
return result ;
}
if (basisType==sinhBasis) {
Double_t dgamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
<< ") 6th form = " << dgamma << ", tau = " << tau << endl;
Double_t result(0);
Double_t tau1 = 1/(1/tau-dgamma/2) ;
Double_t tau2 = 1/(1/tau+dgamma/2) ;
if (basisSign!=Minus) result += 0.5*(calcDecayConv(+1,tau1,sig,rtau,fsign)-calcDecayConv(+1,tau2,sig,rtau,fsign));
if (basisSign!=Plus) result += 0.5*(calcDecayConv(-1,tau2,sig,rtau,fsign)-calcDecayConv(-1,tau1,sig,rtau,fsign));
return result;
}
if (basisType==coshBasis) {
Double_t dgamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
<< ") 7th form = " << dgamma << ", tau = " << tau << endl;
Double_t result(0);
Double_t tau1 = 1/(1/tau-dgamma/2) ;
Double_t tau2 = 1/(1/tau+dgamma/2) ;
if (basisSign!=Minus) result += 0.5*(calcDecayConv(+1,tau1,sig,rtau,fsign)+calcDecayConv(+1,tau2,sig,rtau,fsign));
if (basisSign!=Plus) result += 0.5*(calcDecayConv(-1,tau1,sig,rtau,fsign)+calcDecayConv(-1,tau2,sig,rtau,fsign));
return result;
}
R__ASSERT(0) ;
return 0 ;
}
Double_t RooGExpModel::logErfC(Double_t xx) const
{
Double_t t,z,ans;
z=fabs(xx);
t=1.0/(1.0+0.5*z);
if(xx >= 0.0)
ans=log(t)+(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+t*(-0.18628806+
t*(0.27886807+t*(-1.13520398+t*(1.48851587+t*(-0.82215223+t*0.17087277)))))))));
else
ans=log(2.0-t*exp(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+t*(-0.18628806+
t*(0.27886807+t*(-1.13520398+t*(1.48851587+t*(-0.82215223+t*0.17087277))))))))));
return ans;
}
std::complex<Double_t> RooGExpModel::calcSinConv(Double_t sign, Double_t sig, Double_t tau, Double_t omega, Double_t rtau, Double_t fsign) const
{
static Double_t root2(sqrt(2.)) ;
Double_t s1= -sign*x/tau;
Double_t c1= sig/(root2*tau);
Double_t u1= s1/(2*c1);
Double_t s2= x/rtau;
Double_t c2= sig/(root2*rtau);
Double_t u2= fsign*s2/(2*c2) ;
std::complex<Double_t> eins(1,0);
std::complex<Double_t> k(1/tau,sign*omega);
return (evalCerf(-sign*omega*tau,u1,c1)+std::complex<Double_t>(evalCerfRe(u2,c2),0)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
}
Double_t RooGExpModel::calcSinConv(Double_t sign, Double_t sig, Double_t tau, Double_t rtau, Double_t fsign) const
{
static Double_t root2(sqrt(2.)) ;
Double_t s1= -sign*x/tau;
Double_t c1= sig/(root2*tau);
Double_t u1= s1/(2*c1);
Double_t s2= x/rtau;
Double_t c2= sig/(root2*rtau);
Double_t u2= fsign*s2/(2*c2) ;
Double_t eins(1);
Double_t k(1/tau);
return (evalCerfRe(u1,c1)+evalCerfRe(u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
}
Double_t RooGExpModel::calcDecayConv(Double_t sign, Double_t tau, Double_t sig, Double_t rtau, Double_t fsign) const
{
static Double_t root2(sqrt(2.)) ;
static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
static Double_t rootpi(sqrt(atan2(0.,-1.)));
Double_t xp(x) ;
xp *= fsign ;
sign *= fsign ;
Double_t cFly;
if ((sign<0)&&(fabs(tau-rtau)<tau/260)) {
Double_t MeanTau=0.5*(tau+rtau);
if (fabs(xp/MeanTau)>300) {
return 0 ;
}
cFly=1./(MeanTau*MeanTau*root2pi) *
exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
*(sig*exp(-1/(2*sig*sig)*TMath::Power((sig*sig/MeanTau+xp),2))
-(sig*sig/MeanTau+xp)*(rootpi/root2)*RooMath::erfc(sig/(root2*MeanTau)+xp/(root2*sig)));
if(_nlo) {
Double_t epsilon=0.5*(tau-rtau);
Double_t a=sig/(root2*MeanTau)+xp/(root2*sig);
cFly += 1./(MeanTau*MeanTau)
*exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
*0.5/MeanTau*epsilon*epsilon*
(exp(-a*a)*(sig/MeanTau*root2/rootpi
-(4*a*sig*sig)/(2*rootpi*MeanTau*MeanTau)
+(-4/rootpi+8*a*a/rootpi)/6
*TMath::Power(sig/(root2*MeanTau),3)
+2/rootpi*(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
(sig/(root2*MeanTau)-a*(sig*sig)/(2*MeanTau*MeanTau))
+2/rootpi*((3*sig*sig)/(2*MeanTau*MeanTau)+xp/MeanTau+
0.5*TMath::Power(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,2))*sig/(root2*MeanTau))
-(2*sig*sig/(MeanTau*MeanTau)+xp/MeanTau+(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
(3*sig*sig/(2*MeanTau*MeanTau)+xp/MeanTau)
+TMath::Power(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,3)/6)*RooMath::erfc(a));
}
} else {
Double_t expArg1 = sig*sig/(2*tau*tau)-sign*xp/tau ;
Double_t expArg2 = sig*sig/(2*rtau*rtau)+xp/rtau ;
Double_t term1, term2 ;
if (expArg1<300) {
term1 = exp(expArg1) *RooMath::erfc(sig/(root2*tau)-sign*xp/(root2*sig)) ;
} else {
term1 = exp(expArg1+logErfC(sig/(root2*tau)-sign*xp/(root2*sig))) ; ;
}
if (expArg2<300) {
term2 = exp(expArg2) *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)) ;
} else {
term2 = exp(expArg2+logErfC(sig/(root2*rtau)+xp/(root2*sig))) ; ;
}
cFly=(term1+sign*term2)/(2*(tau+sign*rtau));
if (cFly<1e-100) {
cFly = 0 ;
}
}
return cFly*2*tau ;
}
Int_t RooGExpModel::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* ) const
{
switch(_basisCode) {
case noBasis:
if (matchArgs(allVars,analVars,convVar())) return 1 ;
break ;
case expBasisPlus:
case expBasisMinus:
case expBasisSum:
case sinBasisPlus:
case sinBasisMinus:
case sinBasisSum:
case cosBasisPlus:
case cosBasisMinus:
case cosBasisSum:
case sinhBasisPlus:
case sinhBasisMinus:
case sinhBasisSum:
case coshBasisPlus:
case coshBasisMinus:
case coshBasisSum:
if (_flatSFInt) {
if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) {
return 2 ;
}
}
if (matchArgs(allVars,analVars,convVar())) return 1 ;
break ;
}
return 0 ;
}
Double_t RooGExpModel::analyticalIntegral(Int_t code, const char* rangeName) const
{
static Double_t root2 = sqrt(2.) ;
Double_t ssfInt(1.0) ;
R__ASSERT(code==1||code==2) ;
if (code==2) {
ssfInt = (ssf.max(rangeName)-ssf.min(rangeName)) ;
}
BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
Double_t tau = (_basisCode!=noBasis)?((RooAbsReal*)basis().getParameter(1))->getVal():0 ;
if (basisType == coshBasis && _basisCode!=noBasis ) {
Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
if (dGamma==0) basisType = expBasis;
}
Double_t fsign = _flip?-1:1 ;
Double_t sig = sigma*ssf ;
Double_t rtau = rlife*rsf ;
if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 1st form" << endl ;
Double_t xpmin = x.min(rangeName)/rtau ;
Double_t xpmax = x.max(rangeName)/rtau ;
Double_t c = sig/(root2*rtau) ;
Double_t umin = xpmin/(2*c) ;
Double_t umax = xpmax/(2*c) ;
Double_t result ;
if (_asympInt) {
result = 1.0 ;
} else {
result = 0.5*evalCerfInt(-fsign,rtau,fsign*umin,fsign*umax,c)/rtau ;
}
if (_basisCode!=0 && basisSign==Both) result *= 2 ;
return result*ssfInt ;
}
Double_t omega = (basisType!=expBasis) ?((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
if (tau==0) {
if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 2nd form" << endl ;
return 0. ;
}
if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
Double_t result(0.);
if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,sig,rtau,fsign,rangeName);
if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,sig,rtau,fsign,rangeName);
return result*ssfInt ;
}
Double_t wt = omega * tau ;
if (basisType==sinBasis) {
if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 4th form omega = "
<< omega << ", tau = " << tau << endl ;
Double_t result(0) ;
if (wt==0) return result ;
if (basisSign!=Minus) result += -1*calcSinConvNorm(+1,tau,omega,sig,rtau,fsign,rangeName).imag();
if (basisSign!=Plus) result += -1*calcSinConvNorm(-1,tau,omega,sig,rtau,fsign,rangeName).imag();
return result*ssfInt ;
}
if (basisType==cosBasis) {
if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName()
<< ") 5th form omega = " << omega << ", tau = " << tau << endl ;
Double_t result(0) ;
if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega,sig,rtau,fsign,rangeName).real();
if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega,sig,rtau,fsign,rangeName).real();
return result*ssfInt ;
}
Double_t dgamma = ((basisType==coshBasis)||(basisType==sinhBasis))?((RooAbsReal*)basis().getParameter(2))->getVal():0 ;
if (basisType==sinhBasis) {
if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName()
<< ") 6th form dgamma = " << dgamma << ", tau = " << tau << endl ;
Double_t tau1 = 1/(1/tau-dgamma/2);
Double_t tau2 = 1/(1/tau+dgamma/2);
Double_t result(0) ;
if (basisSign!=Minus) result += 0.5*(calcSinConvNorm(+1,tau1,sig,rtau,fsign,rangeName)-
calcSinConvNorm(+1,tau2,sig,rtau,fsign,rangeName));
if (basisSign!=Plus) result += 0.5*(calcSinConvNorm(-1,tau2,sig,rtau,fsign,rangeName)-
calcSinConvNorm(-1,tau1,sig,rtau,fsign,rangeName));
return result;
}
if (basisType==coshBasis) {
if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName()
<< ") 6th form dgamma = " << dgamma << ", tau = " << tau << endl ;
Double_t tau1 = 1/(1/tau-dgamma/2);
Double_t tau2 = 1/(1/tau+dgamma/2);
Double_t result(0);
if (basisSign!=Minus) result += 0.5*(calcSinConvNorm(+1,tau1,sig,rtau,fsign,rangeName)+
calcSinConvNorm(+1,tau2,sig,rtau,fsign,rangeName));
if (basisSign!=Plus) result += 0.5*(calcSinConvNorm(-1,tau1,sig,rtau,fsign,rangeName)+
calcSinConvNorm(-1,tau2,sig,rtau,fsign,rangeName));
return result;
}
R__ASSERT(0) ;
return 1 ;
}
std::complex<Double_t> RooGExpModel::calcSinConvNorm(Double_t sign, Double_t tau, Double_t omega,
Double_t sig, Double_t rtau, Double_t fsign, const char* rangeName) const
{
static Double_t root2(sqrt(2.)) ;
Double_t smin1= x.min(rangeName)/tau;
Double_t smax1= x.max(rangeName)/tau;
Double_t c1= sig/(root2*tau);
Double_t umin1= smin1/(2*c1);
Double_t umax1= smax1/(2*c1);
Double_t smin2= x.min(rangeName)/rtau;
Double_t smax2= x.max(rangeName)/rtau;
Double_t c2= sig/(root2*rtau);
Double_t umin2= smin2/(2*c2) ;
Double_t umax2= smax2/(2*c2) ;
std::complex<Double_t> eins(1,0);
std::complex<Double_t> k(1/tau,sign*omega);
std::complex<Double_t> term1 = evalCerfInt(sign,-sign*omega*tau, tau, -sign*umin1, -sign*umax1, c1);
std::complex<Double_t> term2 = std::complex<Double_t>(evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2, c2)*fsign*sign,0);
return (term1+term2)/(eins + k*fsign*sign*rtau) ;
}
Double_t RooGExpModel::calcSinConvNorm(Double_t sign, Double_t tau, Double_t sig, Double_t rtau, Double_t fsign, const char* rangeName) const
{
static Double_t root2(sqrt(2.)) ;
Double_t smin1= x.min(rangeName)/tau;
Double_t smax1= x.max(rangeName)/tau;
Double_t c1= sig/(root2*tau);
Double_t umin1= smin1/(2*c1);
Double_t umax1= smax1/(2*c1);
Double_t smin2= x.min(rangeName)/rtau;
Double_t smax2= x.max(rangeName)/rtau;
Double_t c2= sig/(root2*rtau);
Double_t umin2= smin2/(2*c2) ;
Double_t umax2= smax2/(2*c2) ;
Double_t eins(1);
Double_t k(1/tau);
Double_t term1 = evalCerfInt(sign, tau, -sign*umin1, -sign*umax1, c1);
Double_t term2 = evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2, c2)*fsign*sign;
if (fabs(tau-rtau)<1e-10 && fabs(term1+term2)<1e-10) {
cout << "epsilon method" << endl ;
static Double_t epsilon = 1e-4 ;
return calcSinConvNorm(sign,tau+epsilon,sig,rtau-epsilon,fsign,rangeName) ;
}
return (term1+term2)/(eins + k*fsign*sign*rtau) ;
}
std::complex<Double_t> RooGExpModel::evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
{
std::complex<Double_t> diff;
if (_asympInt) {
diff = std::complex<Double_t>(2,0) ;
} else {
diff = std::complex<Double_t>(sign,0.)*(evalCerf(wt,umin,c) - evalCerf(wt,umax,c) + RooMath::erf(umin) - RooMath::erf(umax));
}
return std::complex<Double_t>(tau/(1.+wt*wt),0)*std::complex<Double_t>(1,wt)*diff;
}
Double_t RooGExpModel::evalCerfInt(Double_t sign, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
{
Double_t diff;
if (_asympInt) {
diff = 2. ;
} else {
if ((umin<-8 && umax>8)||(umax<-8 && umin>8)) {
diff = 2. ;
} else {
diff = sign*(evalCerfRe(umin,c) - evalCerfRe(umax,c) + RooMath::erf(umin) - RooMath::erf(umax));
}
}
return tau*diff;
}
std::complex<Double_t> RooGExpModel::evalCerfApprox(Double_t swt, Double_t u, Double_t c)
{
static Double_t rootpi= sqrt(atan2(0.,-1.));
std::complex<Double_t> z(swt*c,u+c);
std::complex<Double_t> zc(u+c,-swt*c);
std::complex<Double_t> zsq= z*z;
std::complex<Double_t> v= -zsq - u*u;
return std::exp(v)*(-std::exp(zsq)/(zc*rootpi) + 1.)*2.;
}
Int_t RooGExpModel::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t ) const
{
if (matchArgs(directVars,generateVars,x)) return 1 ;
return 0 ;
}
void RooGExpModel::generateEvent(Int_t code)
{
R__ASSERT(code==1) ;
Double_t xgen ;
while(1) {
Double_t xgau = RooRandom::randomGenerator()->Gaus(0,(sigma*ssf));
Double_t xexp = RooRandom::uniform();
if (!_flip) xgen= xgau + (rlife*rsf)*log(xexp);
else xgen= xgau - (rlife*rsf)*log(xexp);
if (xgen<x.max() && xgen>x.min()) {
x = xgen ;
return ;
}
}
}