CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
ElectronMomentumCorrector Class Reference

#include <ElectronMomentumCorrector.h>

Public Member Functions

void correct (reco::GsfElectron &, TrajectoryStateOnSurface &)
 
 ElectronMomentumCorrector ()
 
 ~ElectronMomentumCorrector ()
 

Private Member Functions

float energyError (float E, float *par) const
 

Private Attributes

float errorEnergy_
 
float errorTrackMomentum_
 
math::XYZTLorentzVector newMomentum_
 

Detailed Description

Definition at line 19 of file ElectronMomentumCorrector.h.

Constructor & Destructor Documentation

ElectronMomentumCorrector::ElectronMomentumCorrector ( )
inline

Definition at line 24 of file ElectronMomentumCorrector.h.

References newMomentum_.

XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
math::XYZTLorentzVector newMomentum_
ElectronMomentumCorrector::~ElectronMomentumCorrector ( )
inline

Definition at line 26 of file ElectronMomentumCorrector.h.

26 {}

Member Function Documentation

void ElectronMomentumCorrector::correct ( reco::GsfElectron electron,
TrajectoryStateOnSurface vtxTsos 
)

Definition at line 33 of file ElectronMomentumCorrector.cc.

References reco::GsfElectron::classification(), reco::GsfElectron::correctMomentum(), reco::GsfElectron::ecalEnergy(), reco::GsfElectron::ecalEnergyError(), errorEnergy_, errorTrackMomentum_, reco::GsfElectron::GAP, reco::GsfElectron::GOLDEN, reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), reco::GsfElectron::isGap(), reco::GsfElectron::isMomentumCorrected(), GaussianSumUtilities1D::mode(), MultiGaussianStateTransform::multiState1D(), newMomentum_, reco::LeafCandidate::p4(), mathSSE::return(), reco::GsfElectron::SHOWERING, mathSSE::sqrt(), reco::btau::trackMomentum, reco::GsfElectron::trackMomentumAtVtx(), reco::GsfElectron::UNKNOWN, and SingleGaussianState1D::variance().

Referenced by GsfElectronAlgo::createElectron().

33  {
34 
35  if (electron.isMomentumCorrected())
36  {
37  edm::LogWarning("ElectronMomentumCorrector::correct")<<"already done" ;
38  return ;
39  }
40 
41  newMomentum_ = electron.p4() ; // default
42  int elClass = electron.classification() ;
43 
44  // irrelevant classification
45  if ( (elClass <= reco::GsfElectron::UNKNOWN) ||
46  (elClass>reco::GsfElectron::GAP) )
47  {
48  edm::LogWarning("ElectronMomentumCorrector::correct")<<"unexpected classification" ;
49  return ;
50  }
51 
52  float scEnergy = electron.ecalEnergy() ;
53  errorEnergy_ = electron.ecalEnergyError() ;
54 
55  float trackMomentum = electron.trackMomentumAtVtx().R() ;
56  errorTrackMomentum_ = 999. ;
57 
58  // retreive momentum error
60  GaussianSumUtilities1D qpUtils(qpState);
61  errorTrackMomentum_ = trackMomentum*trackMomentum*sqrt(qpUtils.mode().variance());
62 
63  float finalMomentum = electron.p4().t(); // initial
64  float finalMomentumError = 999.;
65 
66 
67  // first check for large errors
68 
69  if (errorTrackMomentum_/trackMomentum > 0.5 && errorEnergy_/scEnergy <= 0.5) {
70  finalMomentum = scEnergy; finalMomentumError = errorEnergy_;
71  }
72  else if (errorTrackMomentum_/trackMomentum <= 0.5 && errorEnergy_/scEnergy > 0.5){
73  finalMomentum = trackMomentum; finalMomentumError = errorTrackMomentum_;
74  }
75  else if (errorTrackMomentum_/trackMomentum > 0.5 && errorEnergy_/scEnergy > 0.5){
76  if (errorTrackMomentum_/trackMomentum < errorEnergy_/scEnergy) {
77  finalMomentum = trackMomentum; finalMomentumError = errorTrackMomentum_;
78  }
79  else{
80  finalMomentum = scEnergy; finalMomentumError = errorEnergy_;
81  }
82  }
83 
84  // then apply the combination algorithm
85  else {
86 
87  // calculate E/p and corresponding error
88  float eOverP = scEnergy / trackMomentum;
89  float errorEOverP = sqrt(
90  (errorEnergy_/trackMomentum)*(errorEnergy_/trackMomentum) +
91  (scEnergy*errorTrackMomentum_/trackMomentum/trackMomentum)*
92  (scEnergy*errorTrackMomentum_/trackMomentum/trackMomentum));
93 
94  if ( eOverP > 1 + 2.5*errorEOverP )
95  {
96  finalMomentum = scEnergy; finalMomentumError = errorEnergy_;
97  if ((elClass==reco::GsfElectron::GOLDEN) && electron.isEB() && (eOverP<1.15))
98  {
99  if (scEnergy<15) {finalMomentum = trackMomentum ; finalMomentumError = errorTrackMomentum_;}
100  }
101  }
102  else if ( eOverP < 1 - 2.5*errorEOverP )
103  {
104  finalMomentum = scEnergy; finalMomentumError = errorEnergy_;
105  if (elClass==reco::GsfElectron::SHOWERING)
106  {
107  if (electron.isEB())
108  {
109  if(scEnergy<18) {finalMomentum = trackMomentum; finalMomentumError = errorTrackMomentum_;}
110  }
111  else if (electron.isEE())
112  {
113  if(scEnergy<13) {finalMomentum = trackMomentum; finalMomentumError = errorTrackMomentum_;}
114  }
115  else
116  { edm::LogWarning("ElectronMomentumCorrector::correct")<<"nor barrel neither endcap electron ?!" ; }
117  }
118  else if (electron.isGap())
119  {
120  if(scEnergy<60) {finalMomentum = trackMomentum; finalMomentumError = errorTrackMomentum_;}
121  }
122  }
123  else
124  {
125  // combination
126  finalMomentum = (scEnergy/errorEnergy_/errorEnergy_ + trackMomentum/errorTrackMomentum_/errorTrackMomentum_) /
128  float finalMomentumVariance = 1 / (1/errorEnergy_/errorEnergy_ + 1/errorTrackMomentum_/errorTrackMomentum_);
129  finalMomentumError = sqrt(finalMomentumVariance);
130  }
131 
132  }
133 
134  math::XYZTLorentzVector oldMomentum = electron.p4() ;
136  ( oldMomentum.x()*finalMomentum/oldMomentum.t(),
137  oldMomentum.y()*finalMomentum/oldMomentum.t(),
138  oldMomentum.z()*finalMomentum/oldMomentum.t(),
139  finalMomentum ) ;
140 
141 
142  // final set
143  electron.correctMomentum(newMomentum_,errorTrackMomentum_,finalMomentumError);
144 
145  }
bool isEE() const
Definition: GsfElectron.h:335
bool isEB() const
Definition: GsfElectron.h:334
return((rh^lh)&mask)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
bool isGap() const
Definition: GsfElectron.h:336
bool isMomentumCorrected() const
Definition: GsfElectron.h:589
MultiGaussianState1D multiState1D(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &, unsigned int)
math::XYZVector trackMomentumAtVtx() const
Definition: GsfElectron.h:244
T sqrt(T t)
Definition: SSEVec.h:28
void correctMomentum(const LorentzVector &momentum, float trackMomentumError, float electronMomentumError)
Definition: GsfElectron.cc:153
float ecalEnergyError() const
Definition: GsfElectron.h:588
Classification classification() const
Definition: GsfElectron.h:530
float ecalEnergy() const
Definition: GsfElectron.h:587
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
math::XYZTLorentzVector newMomentum_
float ElectronMomentumCorrector::energyError ( float  E,
float *  par 
) const
private

Member Data Documentation

float ElectronMomentumCorrector::errorEnergy_
private

Definition at line 43 of file ElectronMomentumCorrector.h.

Referenced by correct().

float ElectronMomentumCorrector::errorTrackMomentum_
private

Definition at line 44 of file ElectronMomentumCorrector.h.

Referenced by correct().

math::XYZTLorentzVector ElectronMomentumCorrector::newMomentum_
private

Definition at line 42 of file ElectronMomentumCorrector.h.

Referenced by correct(), and ElectronMomentumCorrector().