![]() |
![]() |
#include <HZZ4LRooPdfs.h>
Public Member Functions | |
virtual TObject * | clone (const char *newname) const |
Roo4lMasses2D_BkgGGZZ () | |
Roo4lMasses2D_BkgGGZZ (const char *name, const char *title, RooAbsReal &_mZstar, RooAbsReal &_mZZ, RooAbsReal &_channelVal) | |
Roo4lMasses2D_BkgGGZZ (const Roo4lMasses2D_BkgGGZZ &other, const char *name=0) | |
virtual | ~Roo4lMasses2D_BkgGGZZ () |
Protected Member Functions | |
Double_t | evaluate () const |
Double_t | UnitStep (double arg) const |
Protected Attributes | |
RooRealProxy | channelVal |
RooRealProxy | mZstar |
RooRealProxy | mZZ |
Definition at line 277 of file HZZ4LRooPdfs.h.
Roo4lMasses2D_BkgGGZZ::Roo4lMasses2D_BkgGGZZ | ( | ) |
Referenced by clone().
Roo4lMasses2D_BkgGGZZ::Roo4lMasses2D_BkgGGZZ | ( | const char * | name, |
const char * | title, | ||
RooAbsReal & | _mZstar, | ||
RooAbsReal & | _mZZ, | ||
RooAbsReal & | _channelVal | ||
) |
Definition at line 2512 of file HZZ4LRooPdfs.cc.
: RooAbsPdf(name,title), mZstar("mZstar","mZstar",this,_mZstar), mZZ("mZZ","mZZ",this,_mZZ), channelVal("channelVal","channelVal",this,_channelVal) { }
Roo4lMasses2D_BkgGGZZ::Roo4lMasses2D_BkgGGZZ | ( | const Roo4lMasses2D_BkgGGZZ & | other, |
const char * | name = 0 |
||
) |
Definition at line 2525 of file HZZ4LRooPdfs.cc.
: RooAbsPdf(other,name), mZstar("mZstar",this,other.mZstar), mZZ("mZZ",this,other.mZZ), channelVal("channelVal",this,other.channelVal) { }
virtual Roo4lMasses2D_BkgGGZZ::~Roo4lMasses2D_BkgGGZZ | ( | ) | [inline, virtual] |
Definition at line 287 of file HZZ4LRooPdfs.h.
{ }
virtual TObject* Roo4lMasses2D_BkgGGZZ::clone | ( | const char * | newname | ) | const [inline, virtual] |
Definition at line 286 of file HZZ4LRooPdfs.h.
References Roo4lMasses2D_BkgGGZZ().
{ return new Roo4lMasses2D_BkgGGZZ(*this,newname); }
double Roo4lMasses2D_BkgGGZZ::evaluate | ( | ) | const [protected] |
*
*
Definition at line 2545 of file HZZ4LRooPdfs.cc.
References beta, channelVal, gather_cfg::cout, alignCSCRings::e, create_public_lumi_plots::exp, f, Gamma, funct::integral(), mZstar, mZZ, Pi, and mathSSE::sqrt().
{ double Gamma=0, mZ=0, a0=0, a1=0, a2=0, a3=0, Gamma0=0, m0=0, f=0, f0=0; double a0_bkgd=0, a1_bkgd=0, a2_bkgd=0, a3_bkgd=0, a4_bkgd=0, a5_bkgd=0, a6_bkgd=0, a7_bkgd=0; double a8_bkgd=0, a9_bkgd=0;//, a10_bkgd=0, a11_bkgd=0, a12_bkgd=0, a13_bkgd=0; a2 = 0.; a3 = 0.; if (channelVal == 1.){ mZ = 91.2; Gamma0 = 24.3063 + -0.169814*mZZ + 0.00274393*(mZZ - 108.86)*(mZZ - 108.86); m0 = -273.264 + 6.81815*mZZ + -0.0520652*TMath::Power(mZZ,2) + 0.000129716*TMath::Power(mZZ,3); a0 = 0.00205023 + 1.49685e-05*mZZ + 2.45743e-07*(mZZ - 110.553)*(mZZ - 110.553); a1 = -34.7502 + 1.15173*mZZ +-0.0152287*TMath::Power(mZZ,2) + 0.000100538*TMath::Power(mZZ,3)+ -3.317e-07*TMath::Power(mZZ,4)+ 4.37829e-10*TMath::Power(mZZ,5); a2 = 0.0710702 + -0.00232026*mZZ + 3.00042e-05*TMath::Power(mZZ,2) + -1.92831e-07*TMath::Power(mZZ,3)+ 6.18009e-10*TMath::Power(mZZ,4)+ -7.92751e-13*TMath::Power(mZZ,5); f = 4.81681e-01*(TMath::Erf( (mZZ-1.47620e+02)/8.12397e-01 ) + 1.); Gamma = 252.406 + -5.11437*mZZ + 0.0345238*TMath::Power(mZZ,2) + -7.44284e-05*TMath::Power(mZZ,3); f0 = 3.68735e-01*(TMath::Erf( (mZZ-1.17687e+02)/2.61381e-01 ) + 1.); a0_bkgd = 138.962; a1_bkgd = 30.6644; a2_bkgd = 102.09; a3_bkgd = 0.0416272; a4_bkgd = 178.464; a5_bkgd = 12.6718; a6_bkgd = 38.1437; a7_bkgd = 0.13019; a8_bkgd = 38.1146; a9_bkgd = 0.0277349; } else if (channelVal == 2.){ mZ = 91.2; Gamma0 = 64.5791 + -0.465027*mZZ + 0.00680112*(mZZ - 109.083)*(mZZ - 109.083); m0 = -4528.09 + 98.4925*mZZ + -0.70193*TMath::Power(mZZ,2) + 0.00164952*TMath::Power(mZZ,3); a0 = 0.00205023 + 1.49685e-05*mZZ + 2.45743e-07*(mZZ - 110.553)*(mZZ - 110.553); a1 = 0.0549555 + -0.000367875*mZZ + 6.16188e-10*(mZZ - 106.028)*(mZZ - 106.028)*(mZZ - 106.028)*(mZZ - 106.028); a2 = -0.000222052 + 1.51606e-06*mZZ + -9.57278e-08*(mZZ - 144.667)*(mZZ - 144.667); f = 4.82988e-01*(TMath::Erf( (mZZ-1.47620e+02)/8.12397e-01 ) + 1.); Gamma = 252.406 + -5.11437*mZZ + 0.0345238*TMath::Power(mZZ,2) + -7.44284e-05*TMath::Power(mZZ,3); f0 = 3.68735e-01*(TMath::Erf( (mZZ-1.17687e+02)/2.61381e-01 ) + 1.); a0_bkgd = 138.962; a1_bkgd = 30.6644; a2_bkgd = 102.09; a3_bkgd = 0.0416272; a4_bkgd = 178.464; a5_bkgd = 12.6718; a6_bkgd = 38.1437; a7_bkgd = 0.13019; a8_bkgd = 38.1146; a9_bkgd = 0.0277349; } else if (channelVal == 3.){ mZ = 91.2; Gamma0 = 36.9986 + -0.250136*mZZ + 0.00474031*(mZZ - 108.732)*(mZZ - 108.732); m0 = -3266.9 + 73.5969*mZZ + -0.544448*TMath::Power(mZZ,2) + 0.00133127*TMath::Power(mZZ,3); a0 = 0.00205023 + 1.49685e-05*mZZ + 2.45743e-07*(mZZ - 110.553)*(mZZ - 110.553); a1 = 0.0206819 + -0.000168968*mZZ + 4.26934e-11*(mZZ - 40.4456)*(mZZ - 40.4456)*(mZZ - 40.4456)*(mZZ - 40.4456); a2 = -1.271e-05 + 8.742e-08*mZZ + -6.34276e-08*(mZZ - 144.032)*(mZZ - 144.032); f = 4.82988e-01*(TMath::Erf( (mZZ-1.47620e+02)/8.12397e-01 ) + 1.); Gamma = 252.406 + -5.11437*mZZ + 0.0345238*TMath::Power(mZZ,2) + -7.44284e-05*TMath::Power(mZZ,3); f0 = 3.68735e-01*(TMath::Erf( (mZZ-1.17687e+02)/2.61381e-01 ) + 1.); a0_bkgd = 138.962; a1_bkgd = 30.6644; a2_bkgd = 102.09; a3_bkgd = 0.0416272; a4_bkgd = 178.464; a5_bkgd = 12.6718; a6_bkgd = 38.1437; a7_bkgd = 0.13019; a8_bkgd = 38.1146; a9_bkgd = 0.0277349; } else{ std::cout << "Invalid paramters for this PDF!" << std::endl; } // ZZ shape double ZZshape = (.5+.5*TMath::Erf((mZZ-a0_bkgd)/a1_bkgd))*(a3_bkgd/(1+exp((mZZ-a0_bkgd)/a2_bkgd)))+(.5+.5*TMath::Erf((mZZ-a4_bkgd)/a5_bkgd))*(a7_bkgd/(1+exp((mZZ-a4_bkgd)/a6_bkgd))+a9_bkgd/(1+exp((mZZ-a4_bkgd)/a8_bkgd))); double mZDistribution, acceptance; //================= mZ Distribution ================= double beta= (1-(mZstar-mZ)*(mZstar-mZ)/(mZZ*mZZ))*(1-(mZstar+mZ)*(mZstar+mZ)/(mZZ*mZZ)); if(beta<0) return 0.0000001; else mZDistribution = beta; //================= acceptance (chebychev polynomials) acceptance = a0+a1*mZstar+a2*(2*mZstar*mZstar-1)+a3*(4*mZstar*mZstar*mZstar-3*mZstar); //================= on-shell Z contamination ======== double onShellZ = exp(-(mZstar-mZ)*(mZstar-mZ)/(2*Gamma*Gamma))/sqrt(2*3.1415*Gamma*Gamma); //double errorFunc = TMath::Erf((mZstar - m0)/Gamma); double errorFunc = exp(-(mZstar-m0)*(mZstar-m0)/(2*Gamma0*Gamma0))/sqrt(2*3.1415*Gamma0*Gamma0); //std::cout << mZZDistribution << " : " << mZDistribution << " : " << acceptance << std::endl; double final = mZDistribution*((1-f)*(f0*acceptance+(1-f0)*errorFunc)+f*onShellZ); if (final <= 0) final = 1e-12; double integral; double E = 2.71828183; // ======================================================== // integral double mZstarFix=12.; double integralT1Lo=((1 - f)*f0*((a1*Power(mZstarFix,6))/6. + (2*a2*Power(mZstarFix,7))/7. + (a0 - a2)*mZstarFix*Power(Power(mZ,2) - Power(mZZ,2),2) + (a1*Power(mZstarFix,2)*Power(Power(mZ,2) - Power(mZZ,2),2))/2. - (a1*Power(mZstarFix,4)*(Power(mZ,2) + Power(mZZ,2)))/2. + (Power(mZstarFix,5)*(a0 - a2*(1 + 4*Power(mZ,2) + 4*Power(mZZ,2))))/5. + (2*Power(mZstarFix,3)*(-(a0*(Power(mZ,2) + Power(mZZ,2))) + a2*(Power(mZ,4) + Power(mZZ,2) + Power(mZZ,4) + Power(mZ,2)*(1 - 2*Power(mZZ,2)))))/3.))/Power(mZZ,4); mZstarFix=mZZ-mZ; double integralT1Hi=((1 - f)*f0*((a1*Power(mZstarFix,6))/6. + (2*a2*Power(mZstarFix,7))/7. + (a0 - a2)*mZstarFix*Power(Power(mZ,2) - Power(mZZ,2),2) + (a1*Power(mZstarFix,2)*Power(Power(mZ,2) - Power(mZZ,2),2))/2. - (a1*Power(mZstarFix,4)*(Power(mZ,2) + Power(mZZ,2)))/2. + (Power(mZstarFix,5)*(a0 - a2*(1 + 4*Power(mZ,2) + 4*Power(mZZ,2))))/5. + (2*Power(mZstarFix,3)*(-(a0*(Power(mZ,2) + Power(mZZ,2))) + a2*(Power(mZ,4) + Power(mZZ,2) + Power(mZZ,4) + Power(mZ,2)*(1 - 2*Power(mZZ,2)))))/3.))/Power(mZZ,4); //std::cout << "IntegralHI: " << integralT1Hi << ", IntegralLO: " << integralT1Lo << std::endl; mZstarFix=12.; double t2loA = (3*Power(Gamma0,4) + Power(m0,4) + Power(mZ,4) - 2*Power(Gamma0,2)*Power(mZZ,2) + Power(mZZ,4) + Power(m0,2)*(6*Power(Gamma0,2) - 2*Power(mZ,2) - 2*Power(mZZ,2)) - 2*Power(mZ,2)*(Power(Gamma0,2) + Power(mZZ,2))); double integralT2Lo=((-1 + f)*(-1 + f0)*(-1.*((Power(Gamma0,2)*(Power(m0,3) + Power(m0,2)*mZstarFix + mZstarFix*(3*Power(Gamma0,2) - 2*Power(mZ,2) + Power(mZstarFix,2) - 2*Power(mZZ,2)) + m0*(5*Power(Gamma0,2) - 2*Power(mZ,2) + Power(mZstarFix,2) - 2*Power(mZZ,2))))/Power(E,Power(m0 - mZstarFix,2)/(2.*Power(Gamma0,2)))) - Gamma0*t2loA*sqrt(TMath::Pi()/2.)*TMath::Erf((m0 - mZstarFix)/(sqrt(2)*Gamma0))))/Power(mZZ,4); integralT2Lo/=sqrt(2*TMath::Pi()*Gamma0*Gamma0); mZstarFix=mZZ-mZ; double t2hiA = (3*Power(Gamma0,4) + Power(m0,4) + Power(mZ,4) - 2*Power(Gamma0,2)*Power(mZZ,2) + Power(mZZ,4) + Power(m0,2)*(6*Power(Gamma0,2) - 2*Power(mZ,2) - 2*Power(mZZ,2)) - 2*Power(mZ,2)*(Power(Gamma0,2) + Power(mZZ,2))); double integralT2Hi=((-1 + f)*(-1 + f0)*(-1.*((Power(Gamma0,2)*(Power(m0,3) + Power(m0,2)*mZstarFix + mZstarFix*(3*Power(Gamma0,2) - 2*Power(mZ,2) + Power(mZstarFix,2) - 2*Power(mZZ,2)) + m0*(5*Power(Gamma0,2) - 2*Power(mZ,2) + Power(mZstarFix,2) - 2*Power(mZZ,2))))/Power(E,Power(m0 - mZstarFix,2)/(2.*Power(Gamma0,2)))) - Gamma0*t2hiA*sqrt(TMath::Pi()/2.)*TMath::Erf((m0 - mZstarFix)/(sqrt(2)*Gamma0))))/Power(mZZ,4); integralT2Hi/=sqrt(2*TMath::Pi()*Gamma0*Gamma0); mZstarFix=12.; double t3loA = (3*Power(Gamma,4) + 4*Power(Gamma,2)*Power(mZ,2) - 2*(Power(Gamma,2) + 2*Power(mZ,2))*Power(mZZ,2) + Power(mZZ,4)); double integralT3Lo=(f*((Power(Gamma,2)*(Power(mZ,3) + Power(mZ,2)*mZstarFix - mZstarFix*(3*Power(Gamma,2) + Power(mZstarFix,2) - 2*Power(mZZ,2)) - mZ*(5*Power(Gamma,2) + Power(mZstarFix,2) - 2*Power(mZZ,2))))/Power(E,Power(mZ - mZstarFix,2)/(2.*Power(Gamma,2))) - Gamma*t3loA*sqrt(TMath::Pi()/2.)*TMath::Erf((mZ - mZstarFix)/(sqrt(2)*Gamma))))/Power(mZZ,4); integralT3Lo/=sqrt(2*TMath::Pi()*Gamma*Gamma); mZstarFix=mZZ-mZ; double t3hiA = (3*Power(Gamma,4) + 4*Power(Gamma,2)*Power(mZ,2) - 2*(Power(Gamma,2) + 2*Power(mZ,2))*Power(mZZ,2) + Power(mZZ,4)); double integralT3Hi=(f*((Power(Gamma,2)*(Power(mZ,3) + Power(mZ,2)*mZstarFix - mZstarFix*(3*Power(Gamma,2) + Power(mZstarFix,2) - 2*Power(mZZ,2)) - mZ*(5*Power(Gamma,2) + Power(mZstarFix,2) - 2*Power(mZZ,2))))/Power(E,Power(mZ - mZstarFix,2)/(2.*Power(Gamma,2))) - Gamma*t3hiA*sqrt(TMath::Pi()/2.)*TMath::Erf((mZ - mZstarFix)/(sqrt(2)*Gamma))))/Power(mZZ,4); integralT3Hi/=sqrt(2*TMath::Pi()*Gamma*Gamma); integral = integralT1Hi + integralT2Hi + integralT3Hi - integralT1Lo - integralT2Lo - integralT3Lo; if (integral <= 0.) { std::cout << "integral is less than zero...." << std::endl; //integral = 100000000.; } return (final/integral)*ZZshape; //*/ //return final; }
double Roo4lMasses2D_BkgGGZZ::UnitStep | ( | double | arg | ) | const [protected] |
Definition at line 2534 of file HZZ4LRooPdfs.cc.
{ if(arg<0.0){ return 0.0; } else{ return 1.0; } }
RooRealProxy Roo4lMasses2D_BkgGGZZ::channelVal [protected] |
Definition at line 293 of file HZZ4LRooPdfs.h.
Referenced by evaluate().
RooRealProxy Roo4lMasses2D_BkgGGZZ::mZstar [protected] |
Definition at line 291 of file HZZ4LRooPdfs.h.
Referenced by evaluate().
RooRealProxy Roo4lMasses2D_BkgGGZZ::mZZ [protected] |
Definition at line 292 of file HZZ4LRooPdfs.h.
Referenced by evaluate().