#include <FML3PtSmearer.h>
Public Member Functions | |
FML3PtSmearer (const RandomEngine *engine) | |
Constructor. | |
math::XYZTLorentzVector | smear (math::XYZTLorentzVector simP4, math::XYZVector recP3) const |
smear the transverse momentum of a reco::Muon | |
~FML3PtSmearer () | |
Destructor. | |
Private Member Functions | |
double | error (double thePt, double theEta) const |
smear the transverse momentum of a SimplL3MuGMTCand | |
double | funShift (double x) const |
double | funSigma (double eta, double pt) const |
double | funSigmaEta (double x) const |
double | funSigmaPt (double x) const |
double | shift (double thePt, double theEta) const |
Private Attributes | |
const RandomEngine * | random |
Static Private Attributes | |
static double | MuonMassSquared_ = 0.10565837*0.10565837 |
Class to deal with the 'smearing' of the L3 muon transverse momentum. The output momentum is generated according to the probablility that a MC generated muon with the same pt leads to that momentum value in the GMT.
Definition at line 20 of file FML3PtSmearer.h.
FML3PtSmearer::FML3PtSmearer | ( | const RandomEngine * | engine | ) |
FML3PtSmearer::~FML3PtSmearer | ( | ) |
double FML3PtSmearer::error | ( | double | thePt, |
double | theEta | ||
) | const [private] |
smear the transverse momentum of a SimplL3MuGMTCand
Definition at line 54 of file FML3PtSmearer.cc.
References funSigma().
Referenced by smear().
{ return funSigma(fabs(theEta),thePt)/thePt; }
double FML3PtSmearer::funShift | ( | double | x | ) | const [private] |
Definition at line 62 of file FML3PtSmearer.cc.
References ExpressReco_HICollisions_FallBack::x.
Referenced by shift().
double FML3PtSmearer::funSigma | ( | double | eta, |
double | pt | ||
) | const [private] |
Definition at line 68 of file FML3PtSmearer.cc.
References funSigmaEta(), funSigmaPt(), and ExpressReco_HICollisions_FallBack::sigma.
Referenced by error().
{ double sigma = funSigmaPt(pt) * funSigmaEta(eta); return sigma; }
double FML3PtSmearer::funSigmaEta | ( | double | x | ) | const [private] |
double FML3PtSmearer::funSigmaPt | ( | double | x | ) | const [private] |
Definition at line 73 of file FML3PtSmearer.cc.
References ExpressReco_HICollisions_FallBack::x.
Referenced by funSigma().
double FML3PtSmearer::shift | ( | double | thePt, |
double | theEta | ||
) | const [private] |
Definition at line 58 of file FML3PtSmearer.cc.
References funShift().
Referenced by smear().
{ return funShift(fabs(theEta)); }
math::XYZTLorentzVector FML3PtSmearer::smear | ( | math::XYZTLorentzVector | simP4, |
math::XYZVector | recP3 | ||
) | const |
smear the transverse momentum of a reco::Muon
Definition at line 19 of file FML3PtSmearer.cc.
References error(), RandomEngine::gaussShoot(), MuonMassSquared_, random, shift(), and mathSSE::sqrt().
Referenced by ParamL3MuonProducer::produce().
{ double ptSim = std::sqrt(simP4.perp2()); double invPtSim; if (ptSim>0.) invPtSim = 1./ptSim; else { // Better if we throw an exception here... simP4.SetPx(recP3.X()); simP4.SetPy(recP3.Y()); simP4.SetPz(recP3.Z()); double muonEnergy=std::sqrt(simP4.P()*simP4.P()+MuonMassSquared_); simP4.SetE(muonEnergy); return simP4; } double etaSim = simP4.eta(); double invPtNew = random->gaussShoot(invPtSim,error(ptSim,etaSim)); invPtNew /= ( 1. + shift(ptSim,etaSim)*invPtNew); if (invPtNew>0.) { double invPtRec = std::sqrt(recP3.perp2()); if (invPtRec>0) invPtRec = 1./invPtRec; else invPtRec=invPtSim; simP4.SetPx(recP3.x()*invPtRec/invPtNew); simP4.SetPy(recP3.y()*invPtRec/invPtNew); simP4.SetPz(recP3.z()*invPtRec/invPtNew); double muonEnergy=std::sqrt(simP4.P()*simP4.P()+MuonMassSquared_); simP4.SetE(muonEnergy); } else { simP4.SetPx(recP3.x()); simP4.SetPy(recP3.y()); simP4.SetPz(recP3.z()); double muonEnergy=std::sqrt(simP4.P()*simP4.P()+MuonMassSquared_); simP4.SetE(muonEnergy); } return simP4; }
double FML3PtSmearer::MuonMassSquared_ = 0.10565837*0.10565837 [static, private] |
Definition at line 35 of file FML3PtSmearer.h.
Referenced by smear().
const RandomEngine* FML3PtSmearer::random [private] |
Definition at line 37 of file FML3PtSmearer.h.
Referenced by smear().