CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauIsolationMasking.cc
Go to the documentation of this file.
1 #include <boost/foreach.hpp>
8 
9 namespace reco { namespace tau {
10 
11 namespace {
12 class DRSorter {
13  public:
14  DRSorter(const reco::Candidate::LorentzVector& axis):axis_(axis){}
15  template<typename T>
16  bool operator()(const T& t1, const T& t2) const {
17  return reco::deltaR(t1->p4(), axis_) < reco::deltaR(t2->p4(), axis_);
18  }
19  private:
21 };
22 // Sort by descending pt
23 class PtSorter {
24  public:
25  PtSorter(){}
26  template<typename T>
27  bool operator()(const T& t1, const T& t2) const {
28  return t1->pt() > t2->pt();
29  }
30 };
31 
32 // Check if an object is within DR of a track collection
33 class MultiTrackDRFilter {
34  public:
35  MultiTrackDRFilter(double deltaR, const std::vector<reco::PFCandidatePtr>& trks)
36  :deltaR_(deltaR),tracks_(trks){}
37  template <typename T>
38  bool operator()(const T& t) const {
39  BOOST_FOREACH(const reco::PFCandidatePtr& trk, tracks_) {
40  if (reco::deltaR(trk->p4(), t->p4()) < deltaR_)
41  return true;
42  }
43  return false;
44  }
45  private:
46  double deltaR_;
47  const std::vector<reco::PFCandidatePtr>& tracks_;
48 };
49 
50 double square(double x) { return x*x; }
51 
52 template<typename T>
53 std::vector<reco::PFCandidatePtr> convertRefCollection(const T& coll) {
54  std::vector<reco::PFCandidatePtr> output;
55  output.reserve(coll.size());
56  BOOST_FOREACH(const reco::PFCandidatePtr cand, coll) {
57  output.push_back(cand);
58  }
59  return output;
60 }
61 }
62 
64  const edm::ParameterSet& pset):resolutions_(new PFEnergyResolution) {
65  ecalCone_ = pset.getParameter<double>("ecalCone");
66  hcalCone_ = pset.getParameter<double>("hcalCone");
67  finalHcalCone_ = pset.getParameter<double>("finalHcalCone");
68  maxSigmas_ = pset.getParameter<double>("maxSigmas");
69 }
70 
71 // Need to explicitly define this in the .cc so we can use auto_ptr + forward
73 
77 
78  typedef std::list<reco::PFCandidatePtr> PFCandList;
79  // Copy original iso collections.
80  std::copy(tau.isolationPFGammaCands().begin(),
81  tau.isolationPFGammaCands().end(), std::back_inserter(output.gammas));
83  tau.isolationPFNeutrHadrCands().end(),
84  std::back_inserter(output.h0s));
85 
86  std::vector<PFCandList*> courses;
87  courses.push_back(&(output.h0s));
88  courses.push_back(&(output.gammas));
89  // Mask using each one of the tracks
90  BOOST_FOREACH(const reco::PFCandidatePtr& track,
92  double trackerEnergy = track->energy();
93  double linkedEcalEnergy = track->ecalEnergy();
94  double linkedHcalEnergy = track->hcalEnergy();
95  math::XYZPointF posAtCalo = track->positionAtECALEntrance();
96  // Get the linked calo energies & their errors
97  double linkedSumErrSquared = 0;
98  linkedSumErrSquared += square(
99  resolutions_->getEnergyResolutionEm(linkedEcalEnergy, posAtCalo.eta()));
100  linkedSumErrSquared += square(
101  resolutions_->getEnergyResolutionHad(
102  linkedHcalEnergy, posAtCalo.eta(), posAtCalo.phi()));
103 
104  // energyDelta is the difference between associated Calo and tracker energy
105  double energyDelta = linkedEcalEnergy + linkedHcalEnergy - trackerEnergy;
106 
107  // Sort the neutral hadrons by DR to the track
108  //DRSorter sorter(track->p4());
109  PtSorter sorter;
110 
111  BOOST_FOREACH(PFCandList* course, courses) {
112  // Sort by deltaR to the current track
113  course->sort(sorter);
114  PFCandList::iterator toEatIter = course->begin();
115  // While there are still candidates to eat in this course and they are
116  // within the cone.
117  while (toEatIter != course->end()) {
118  const reco::PFCandidate& toEat = **toEatIter;
119  double toEatEnergy = toEat.energy();
120  double toEatErrorSq = square(resolution(toEat));
121  // Check if we can absorb this candidate into the track.
122  if (inCone(*track, **toEatIter) &&
123  (energyDelta + toEatEnergy)/std::sqrt(
124  linkedSumErrSquared + toEatErrorSq) < maxSigmas_ ) {
125  energyDelta += toEatEnergy;
126  linkedSumErrSquared += toEatErrorSq;
127  toEatIter = course->erase(toEatIter);
128  } else {
129  // otherwise skip to the next one
130  ++toEatIter;
131  }
132  }
133  }
134  }
135  // Filter out any final HCAL objects with in cones about the tracks.
136  // This removes upward fluctuating HCAL objects
137  if (finalHcalCone_ > 0) {
138  MultiTrackDRFilter hcalFinalFilter(finalHcalCone_,
140  std::remove_if(output.h0s.begin(), output.h0s.end(), hcalFinalFilter);
141  }
142  return output;
143 }
144 
146  const reco::PFCandidate& cand) const {
147  if (cand.particleId() == reco::PFCandidate::h0) {
148  // NB for HCAL it returns relative energy
149  return cand.energy()*resolutions_->getEnergyResolutionHad(cand.energy(),
150  cand.eta(), cand.phi());
151  } else if (cand.particleId() == reco::PFCandidate::gamma) {
152  return resolutions_->getEnergyResolutionEm(cand.energy(), cand.eta());
153  } else if (cand.particleId() == reco::PFCandidate::e) {
154  // FIXME what is the electron resolution??
155  return 0.15;
156  } else {
157  edm::LogWarning("IsoMask::res (bad pf id)")
158  << "Unknown PF ID: " << cand.particleId();
159  }
160  return -1;
161 }
162 
164  const reco::PFCandidate& cand) const {
165  double openingDR = reco::deltaR(track.positionAtECALEntrance(), cand.p4());
166  if (cand.particleId() == reco::PFCandidate::h0) {
167  return (openingDR < hcalCone_);
168  } else if (cand.particleId() == reco::PFCandidate::gamma ||
169  cand.particleId() == reco::PFCandidate::e) {
170  return openingDR < ecalCone_;
171  } else {
172  edm::LogWarning("IsoMask::inCone (bad pf id)")
173  << "Unknown PF ID: " << cand.particleId()
174  << " " << reco::PFCandidate::e;
175  }
176  return -1;
177 }
178 
179 }} // end namespace reco::tau
T getParameter(std::string const &) const
const std::vector< reco::PFCandidatePtr > & tracks_
virtual double energy() const final
energy
virtual double phi() const final
momentum azimuthal angle
double deltaR_
double resolution(const reco::PFCandidate &cand) const
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
std::auto_ptr< PFEnergyResolution > resolutions_
const math::XYZPointF & positionAtECALEntrance() const
Definition: PFCandidate.h:367
const std::vector< reco::PFCandidatePtr > & isolationPFGammaCands() const
Gamma candidates in isolation region.
Definition: PFTau.cc:93
RecoTauIsolationMasking(const edm::ParameterSet &pset)
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
T sqrt(T t)
Definition: SSEVec.h:18
const reco::Candidate::LorentzVector & axis_
bool inCone(const reco::PFCandidate &track, const reco::PFCandidate &cand) const
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
const std::vector< reco::PFCandidatePtr > & isolationPFNeutrHadrCands() const
Definition: PFTau.cc:91
JetCorrectorParametersCollection coll
Definition: classes.h:10
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
double square(const double a)
virtual double eta() const final
momentum pseudorapidity
virtual ParticleType particleId() const
Definition: PFCandidate.h:373
long double T
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
const std::vector< reco::PFCandidatePtr > & signalPFChargedHadrCands() const
Charged hadrons in signal region.
Definition: PFTau.cc:80
IsoMaskResult mask(const reco::PFTau &) const
Return a new isolation collections with masking applied.