CMS 3D CMS Logo

EleTkIsolFromCands.h
Go to the documentation of this file.
1 #ifndef RECOEGAMMA_EGAMMAISOLATIONALGOS_ELETKISOLFROMCANDS_H
2 #define RECOEGAMMA_EGAMMAISOLATIONALGOS_ELETKISOLFROMCANDS_H
3 
9 
10 //author S. Harper (RAL)
11 //this class does a simple calculation of the track isolation for a track with eta,
12 //phi and z vtx (typically the GsfTrack of the electron). It uses
13 //PFPackedCandidates as proxies for the tracks. It now has been upgraded to
14 //take general tracks as input also
15 //
16 //Note: due to rounding and space saving issues, the isolation calculated on tracks and
17 //PFPackedCandidates will differ slightly even if the same cuts are used. These differences
18 //are thought to be inconsquencial but as such its important to remake the PFPackedCandidates
19 //in RECO/AOD if you want to be consistant with miniAOD
20 //
21 //The track version is more for variables that are stored in the electron rather than
22 //objects which are recomputed on the fly for AOD/miniAOD
23 //
24 //Note, the tracks in miniAOD have additional cuts on and are missing algo information so this
25 //has to be taken into account when using them
26 
27 //new for 9X:
28 // 1) its now possible to use generalTracks as input
29 // 2) adapted to 9X PF packed candidate improvements
30 // 2a) the PF packed candidates now store the pt of the track in addition to the pt of the
31 // candidate. This means there is no need to pass in the electron collection any more
32 // to try and undo the electron e/p combination when the candidate is an electron
33 // 2b) we now not only store the GsfTrack of the electron but also the general track of the electron
34 // there are three input collections now:
35 // packedPFCandidates : all PF candidates (with the track for ele candidates being the gsftrack)
36 // lostTracks : all tracks which were not made into a PFCandidate but passed some preselection
37 // lostTracks::eleTracks : KF electron tracks
38 // as such to avoid double counting the GSF and KF versions of an electron track, we now need to
39 // specify if the electron PF candidates are to be rejected in the sum over that collection or not.
40 // Note in all this, I'm not concerned about the electron in questions track, that will be rejected,
41 // I'm concerned about near by fake electrons which have been recoed by PF
42 // This is handled by the PIDVeto, which obviously is only used/required when using PFCandidates
43 
44 
46 public:
47  enum class PIDVeto{
48  NONE=0,
49  ELES,
50  NONELES,
51  };
52 
53 private:
54  struct TrkCuts {
55  float minPt;
56  float minDR2;
57  float maxDR2;
58  float minDEta;
59  float maxDZ;
60  float minHits;
61  float minPixelHits;
62  float maxDPtPt;
63  std::vector<reco::TrackBase::TrackQuality> allowedQualities;
64  std::vector<reco::TrackBase::TrackAlgorithm> algosToReject;
65  explicit TrkCuts(const edm::ParameterSet& para);
67  };
68 
70 
71 public:
72  explicit EleTkIsolFromCands(const edm::ParameterSet& para);
73  EleTkIsolFromCands(const EleTkIsolFromCands&)=default;
74  ~EleTkIsolFromCands()=default;
76 
78 
79  std::pair<int,double> calIsol(const reco::TrackBase& trk,const pat::PackedCandidateCollection& cands,const PIDVeto=PIDVeto::NONE)const;
80  std::pair<int,double> calIsol(const double eleEta,const double elePhi,const double eleVZ,
82 
83 
84  std::pair<int,double> calIsol(const reco::TrackBase& trk,const reco::TrackCollection& tracks)const;
85  std::pair<int,double> calIsol(const double eleEta,const double elePhi,const double eleVZ,
86  const reco::TrackCollection& tracks)const;
87 
88  //little helper function for the four calIsol functions for it to directly return the pt
89  template<typename ...Args>
90  double calIsolPt(Args && ...args)const{return calIsol(std::forward<Args>(args)...).second;}
91 
92  static PIDVeto pidVetoFromStr(const std::string& vetoStr);
93  static bool passPIDVeto(const int pdgId,const EleTkIsolFromCands::PIDVeto pidVeto);
94 
95 private:
96  static bool passTrkSel(const reco::TrackBase& trk,
97  const double trkPt,
98  const TrkCuts& cuts,
99  const double eleEta,const double elePhi,
100  const double eleVZ);
101  //no qualities specified, accept all, ORed
102  static bool passQual(const reco::TrackBase& trk,
103  const std::vector<reco::TrackBase::TrackQuality>& quals);
104  static bool passAlgo(const reco::TrackBase& trk,
105  const std::vector<reco::TrackBase::TrackAlgorithm>& algosToRej);
106 };
107 
108 
109 #endif
~EleTkIsolFromCands()=default
std::vector< pat::PackedCandidate > PackedCandidateCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
static bool passAlgo(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackAlgorithm > &algosToRej)
std::vector< reco::TrackBase::TrackQuality > allowedQualities
std::vector< reco::TrackBase::TrackAlgorithm > algosToReject
static edm::ParameterSetDescription pSetDescript()
double calIsolPt(Args &&...args) const
std::pair< int, double > calIsol(const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const PIDVeto=PIDVeto::NONE) const
static bool passQual(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackQuality > &quals)
static PIDVeto pidVetoFromStr(const std::string &vetoStr)
static bool passTrkSel(const reco::TrackBase &trk, const double trkPt, const TrkCuts &cuts, const double eleEta, const double elePhi, const double eleVZ)
EleTkIsolFromCands & operator=(const EleTkIsolFromCands &)=default
EleTkIsolFromCands(const edm::ParameterSet &para)
static bool passPIDVeto(const int pdgId, const EleTkIsolFromCands::PIDVeto pidVeto)