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 Attributes
reco::tau::RecoTauEnergyRecoveryPlugin Class Reference
Inheritance diagram for reco::tau::RecoTauEnergyRecoveryPlugin:
reco::tau::RecoTauModifierPlugin reco::tau::RecoTauEventHolderPlugin reco::tau::RecoTauNamedPlugin

Public Member Functions

virtual void beginEvent ()
 
void operator() (PFTau &) const
 
 RecoTauEnergyRecoveryPlugin (const edm::ParameterSet &)
 
virtual ~RecoTauEnergyRecoveryPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauModifierPlugin
 RecoTauModifierPlugin (const edm::ParameterSet &pset)
 
virtual ~RecoTauModifierPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauEventHolderPlugin
const edm::Eventevt () const
 
const edm::EventSetupevtSetup () const
 
 RecoTauEventHolderPlugin (const edm::ParameterSet &pset)
 
void setup (const edm::Event &, const edm::EventSetup &)
 
virtual ~RecoTauEventHolderPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauNamedPlugin
const std::string & name () const
 
 RecoTauNamedPlugin (const edm::ParameterSet &pset)
 
virtual ~RecoTauNamedPlugin ()
 

Private Attributes

unsigned corrLevel_
 
double lev1EtaWindow_
 
double lev1PhiWindow_
 
RecoTauQualityCutsqcuts_
 
edm::InputTag srcVertices_
 
const reco::VertextheEventVertex_
 

Detailed Description

Definition at line 31 of file RecoTauEnergyRecoveryPlugin.cc.

Constructor & Destructor Documentation

reco::tau::RecoTauEnergyRecoveryPlugin::RecoTauEnergyRecoveryPlugin ( const edm::ParameterSet cfg)
explicit

Definition at line 53 of file RecoTauEnergyRecoveryPlugin.cc.

References edm::ParameterSet::getParameter(), qcuts_, and srcVertices_.

54  : RecoTauModifierPlugin(cfg),
55  theEventVertex_(0),
56  corrLevel_(cfg.getParameter<unsigned>("corrLevel")),
57  lev1PhiWindow_(cfg.getParameter<double>("lev1PhiWindow")),
58  lev1EtaWindow_(cfg.getParameter<double>("lev1EtaWindow"))
59 {
60  edm::ParameterSet cfgQualityCuts = cfg.getParameter<edm::ParameterSet>("qualityCuts");
61  srcVertices_ = cfgQualityCuts.getParameter<edm::InputTag>("primaryVertexSrc");
62  qcuts_ = new RecoTauQualityCuts(cfgQualityCuts.getParameter<edm::ParameterSet>("signalQualityCuts"));
63 }
T getParameter(std::string const &) const
RecoTauModifierPlugin(const edm::ParameterSet &pset)
reco::tau::RecoTauEnergyRecoveryPlugin::~RecoTauEnergyRecoveryPlugin ( )
virtual

Definition at line 65 of file RecoTauEnergyRecoveryPlugin.cc.

References qcuts_.

66 {
67  delete qcuts_;
68 }

Member Function Documentation

void reco::tau::RecoTauEnergyRecoveryPlugin::beginEvent ( )
virtual

Reimplemented from reco::tau::RecoTauModifierPlugin.

Definition at line 70 of file RecoTauEnergyRecoveryPlugin.cc.

References reco::tau::RecoTauEventHolderPlugin::evt(), edm::Event::getByLabel(), qcuts_, reco::tau::RecoTauQualityCuts::setPV(), srcVertices_, and theEventVertex_.

71 {
73  evt()->getByLabel(srcVertices_, vertices);
74 
75  if ( vertices->size() >= 1 ) {
76  qcuts_->setPV(reco::VertexRef(vertices, 0));
77  theEventVertex_ = &vertices->at(0);
78  } else {
79  theEventVertex_ = 0;
80  }
81 }
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
void reco::tau::RecoTauEnergyRecoveryPlugin::operator() ( PFTau tau) const
virtual

Implements reco::tau::RecoTauModifierPlugin.

Definition at line 108 of file RecoTauEnergyRecoveryPlugin.cc.

References corrLevel_, reco::TrackBase::dxy(), reco::TrackBase::dz(), reco::LeafCandidate::energy(), reco::LeafCandidate::eta(), reco::tau::RecoTauQualityCuts::filter(), reco::tau::getTrack(), reco::PFCandidate::h0, reco::TrackBase::hitPattern(), edm::RefToBase< T >::isNonnull(), reco::tau::isTauSignalPFCandidate(), reco::PFTau::jetRef(), lev1EtaWindow_, lev1PhiWindow_, max(), reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidPixelHits(), reco::LeafCandidate::p(), reco::TrackBase::p(), reco::LeafCandidate::phi(), reco::Vertex::position(), reco::TrackBase::pt(), reco::TrackBase::ptError(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), qcuts_, reco::BaseTau::setalternatLorentzVect(), mathSSE::sqrt(), reco::tau::square(), and theEventVertex_.

109 {
110  double tauEnergyCorr = 0.;
111 
112  if ( corrLevel_ & 1 && theEventVertex_ ) {
113 
114  bool needsCorrLevel1 = false;
115 
116  std::vector<reco::PFCandidatePtr> pfJetConstituents = tau.jetRef()->getPFConstituents();
117  for ( std::vector<reco::PFCandidatePtr>::const_iterator pfJetConstituent = pfJetConstituents.begin();
118  pfJetConstituent != pfJetConstituents.end(); ++pfJetConstituent ) {
119  const reco::TrackBaseRef track = getTrack(**pfJetConstituent);
120  if ( track.isNonnull() ) {
121  double dIP = fabs(track->dxy(theEventVertex_->position()));
122  double dZ = fabs(track->dz(theEventVertex_->position()));
123  if ( track->pt() > 2.0 && fabs(tau.eta() - (*pfJetConstituent)->eta()) < lev1EtaWindow_ &&
124  !isTauSignalPFCandidate(tau, *pfJetConstituent) && (dZ < 0.2 || dIP > 0.10) ) needsCorrLevel1 = true;
125  }
126  }
127 
128  if ( needsCorrLevel1 ) {
129  std::vector<reco::PFCandidatePtr> pfJetConstituents = tau.jetRef()->getPFConstituents();
130  for ( std::vector<reco::PFCandidatePtr>::const_iterator pfJetConstituent = pfJetConstituents.begin();
131  pfJetConstituent != pfJetConstituents.end(); ++pfJetConstituent ) {
132  if ( fabs(tau.eta() - (*pfJetConstituent)->eta()) < lev1EtaWindow_ &&
133  fabs(tau.phi() - (*pfJetConstituent)->phi()) < lev1PhiWindow_ ) {
134  if ( (*pfJetConstituent)->particleId() == reco::PFCandidate::h0 ) {
135  tauEnergyCorr += (*pfJetConstituent)->energy();
136  } else {
137  if ( !isTauSignalPFCandidate(tau, *pfJetConstituent) ) {
138  double caloEn = (*pfJetConstituent)->ecalEnergy() + (*pfJetConstituent)->hcalEnergy();
139  tauEnergyCorr += caloEn;
140  }
141  }
142  }
143  }
144  }
145  }
146 
147  if ( corrLevel_ & 2 && theEventVertex_ ) {
148 
149  double leadTrackMom = 0.;
150  double leadTrackMomErr = 0.;
151  double jetCaloEn = 0.;
152 
153  std::vector<reco::PFCandidatePtr> pfJetConstituents = tau.jetRef()->getPFConstituents();
154  for ( std::vector<reco::PFCandidatePtr>::const_iterator pfJetConstituent = pfJetConstituents.begin();
155  pfJetConstituent != pfJetConstituents.end(); ++pfJetConstituent ) {
156  const reco::TrackBaseRef track = getTrack(**pfJetConstituent);
157  if ( track.isNonnull() ) {
158  double trackPt = track->pt();
159  double trackPtErr = track->ptError();
160  if ( qcuts_->filter(**pfJetConstituent) &&
161  trackPtErr < (0.20*trackPt) && track->normalizedChi2() < 5.0 && track->hitPattern().numberOfValidPixelHits() >= 1 &&
162  (trackPt - 3.*trackPtErr) > (*pfJetConstituent)->pt() && trackPt < (3.*tau.jetRef()->pt()) ) {
163  if ( track->p() > leadTrackMom ) {
164  leadTrackMom = track->p();
165  leadTrackMomErr = leadTrackMom*(trackPtErr/trackPt);
166  }
167  }
168  }
169 
170  double caloEn = (*pfJetConstituent)->ecalEnergy() + (*pfJetConstituent)->hcalEnergy();
171  jetCaloEn += caloEn;
172  }
173 
174  if ( leadTrackMom > tau.p() ) {
175  const double chargedPionMass = 0.13957; // GeV
176  double leadTrackEn = sqrt(square(leadTrackMom) + square(chargedPionMass));
177  double jetCaloEnErr = 1.00*sqrt(std::max(jetCaloEn, leadTrackEn));
178  double combEn = ((1./square(jetCaloEnErr))*jetCaloEn + (1./square(leadTrackMomErr))*leadTrackEn)/
179  ((1./square(jetCaloEnErr)) + (1./square(leadTrackMomErr)));
180  tauEnergyCorr = std::max(tauEnergyCorr, combEn - tau.energy());
181  }
182  }
183 
184  if ( tau.energy() > 0. ) {
185  double tauEnergy_corrected = tau.energy() + tauEnergyCorr;
186  double scaleFactor = tauEnergy_corrected/tau.energy();
187  double tauPx_corrected = scaleFactor*tau.px();
188  double tauPy_corrected = scaleFactor*tau.py();
189  double tauPz_corrected = scaleFactor*tau.pz();
190  tau.setalternatLorentzVect(reco::Candidate::LorentzVector(tauPx_corrected, tauPy_corrected, tauPz_corrected, tauEnergy_corrected));
191  }
192 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:111
bool filter(const reco::PFCandidate &cand) const
Filter a single PFCandidate.
const Point & position() const
position
Definition: Vertex.h:93
const T & max(const T &a, const T &b)
bool isTauSignalPFCandidate(const reco::PFTau &tau, const reco::PFCandidatePtr &pfJetConstituent)
T sqrt(T t)
Definition: SSEVec.h:46
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:194
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:127
double square(double x)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
const reco::TrackBaseRef getTrack(const reco::PFCandidate &cand)
int numberOfValidPixelHits() const
Definition: HitPattern.h:566
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:279

Member Data Documentation

unsigned reco::tau::RecoTauEnergyRecoveryPlugin::corrLevel_
private

Definition at line 47 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by operator()().

double reco::tau::RecoTauEnergyRecoveryPlugin::lev1EtaWindow_
private

Definition at line 50 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by operator()().

double reco::tau::RecoTauEnergyRecoveryPlugin::lev1PhiWindow_
private

Definition at line 49 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by operator()().

RecoTauQualityCuts* reco::tau::RecoTauEnergyRecoveryPlugin::qcuts_
private
edm::InputTag reco::tau::RecoTauEnergyRecoveryPlugin::srcVertices_
private

Definition at line 42 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by beginEvent(), and RecoTauEnergyRecoveryPlugin().

const reco::Vertex* reco::tau::RecoTauEnergyRecoveryPlugin::theEventVertex_
private

Definition at line 43 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by beginEvent(), and operator()().