CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EleTkIsolFromCands.h
Go to the documentation of this file.
1 #ifndef RECOEGAMMA_EGAMMAISOLATIONALGOS_ELETKISOLFROMCANDS_H
2 #define RECOEGAMMA_EGAMMAISOLATIONALGOS_ELETKISOLFROMCANDS_H
3 
11 
12 //author S. Harper (RAL)
13 //this class does a simple calculation of the track isolation for a track with eta,
14 //phi and z vtx (typically the GsfTrack of the electron). It uses
15 //PFPackedCandidates as proxies for the tracks. It now has been upgraded to
16 //take general tracks as input also
17 //
18 //Note: due to rounding and space saving issues, the isolation calculated on tracks and
19 //PFPackedCandidates will differ slightly even if the same cuts are used. These differences
20 //are thought to be inconsquencial but as such its important to remake the PFPackedCandidates
21 //in RECO/AOD if you want to be consistant with miniAOD
22 //
23 //The track version is more for variables that are stored in the electron rather than
24 //objects which are recomputed on the fly for AOD/miniAOD
25 //
26 //Note, the tracks in miniAOD have additional cuts on and are missing algo information so this
27 //has to be taken into account when using them
28 
29 //new for 9X:
30 // 1) its now possible to use generalTracks as input
31 // 2) adapted to 9X PF packed candidate improvements
32 // 2a) the PF packed candidates now store the pt of the track in addition to the pt of the
33 // candidate. This means there is no need to pass in the electron collection any more
34 // to try and undo the electron e/p combination when the candidate is an electron
35 // 2b) we now not only store the GsfTrack of the electron but also the general track of the electron
36 // there are three input collections now:
37 // packedPFCandidates : all PF candidates (with the track for ele candidates being the gsftrack)
38 // lostTracks : all tracks which were not made into a PFCandidate but passed some preselection
39 // lostTracks::eleTracks : KF electron tracks
40 // as such to avoid double counting the GSF and KF versions of an electron track, we now need to
41 // specify if the electron PF candidates are to be rejected in the sum over that collection or not.
42 // Note in all this, I'm not concerned about the electron in questions track, that will be rejected,
43 // I'm concerned about near by fake electrons which have been recoed by PF
44 // This is handled by the PIDVeto, which obviously is only used/required when using PFCandidates
45 
47 public:
48  struct TrkCuts {
49  float minPt;
50  float minDR2;
51  float maxDR2;
52  float minDEta;
53  float maxDZ;
54  float minHits;
55  float minPixelHits;
56  float maxDPtPt;
57  std::vector<reco::TrackBase::TrackQuality> allowedQualities;
58  std::vector<reco::TrackBase::TrackAlgorithm> algosToReject;
59  explicit TrkCuts(const edm::ParameterSet& para);
61  };
62 
63  struct Configuration {
64  explicit Configuration(const edm::ParameterSet& para)
65  : barrelCuts(para.getParameter<edm::ParameterSet>("barrelCuts")),
66  endcapCuts(para.getParameter<edm::ParameterSet>("endcapCuts")) {}
69  };
70 
71  enum class PIDVeto {
72  NONE = 0,
73  ELES,
74  NONELES,
75  };
76 
78  : cfg_{cfg}, tracks_{&tracks} {}
81  PIDVeto pidVeto = PIDVeto::NONE)
82  : cfg_{cfg}, cands_{&cands}, pidVeto_{pidVeto} {}
83 
85 
86  static PIDVeto pidVetoFromStr(const std::string& vetoStr);
87 
88  struct Output {
89  const int nTracks;
90  const float ptSum;
91  };
92 
93  Output operator()(const reco::TrackBase& electronTrack);
94 
95 private:
96  // For each electron, we want to try out which tracks are in a cone around
97  // it. However, this will get expensive if there are many electrons and
98  // tracks (Phase II conditions). In particular, calling
99  // reco::TrackBase::eta() many times is costy because eta is not precomputed.
100  // To solve this, we first cache the tracks in a simpler data structure in
101  // which eta is already computed (TrackTable). Furthermore, the tracks are
102  // preselected by the cuts that can already be applied without considering
103  // the electron. Note that this has to be done twice, because the required
104  // preselection is different for barrel and endcap electrons.
105 
107 
108  static bool passPIDVeto(const int pdgId, const EleTkIsolFromCands::PIDVeto pidVeto);
109 
112  TrkCuts const& cuts,
114 
115  static bool passTrackPreselection(const reco::TrackBase& trk, float trkPt, const TrkCuts& cuts);
116 
117  //no qualities specified, accept all, ORed
118  static bool passQual(const reco::TrackBase& trk, const std::vector<reco::TrackBase::TrackQuality>& quals);
119  static bool passAlgo(const reco::TrackBase& trk, const std::vector<reco::TrackBase::TrackAlgorithm>& algosToRej);
120 
122 
124 
125  // All of these member variables are related to the caching of preselected tracks
126  reco::TrackCollection const* tracks_ = nullptr;
133 };
134 
135 #endif
TrkCuts(const edm::ParameterSet &para)
tuple cfg
Definition: looper.py:296
reco::TrackCollection const * tracks_
static bool passTrackPreselection(const reco::TrackBase &trk, float trkPt, const TrkCuts &cuts)
static edm::ParameterSetDescription pSetDescript()
std::vector< pat::PackedCandidate > PackedCandidateCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
auto const & tracks
cannot be loose
static TrackTable preselectTracks(reco::TrackCollection const &tracks, TrkCuts const &cuts)
static bool passAlgo(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackAlgorithm > &algosToRej)
std::vector< reco::TrackBase::TrackQuality > allowedQualities
TkSoA const *__restrict__ CAHitNtupletGeneratorKernelsGPU::QualityCuts cuts
Output operator()(const reco::TrackBase &electronTrack)
std::vector< reco::TrackBase::TrackAlgorithm > algosToReject
static edm::ParameterSetDescription pSetDescript()
Configuration(const edm::ParameterSet &para)
static bool passQual(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackQuality > &quals)
TrackTable preselectedTracksWithBarrelCuts_
static PIDVeto pidVetoFromStr(const std::string &vetoStr)
static TrackTable preselectTracksFromCands(pat::PackedCandidateCollection const &cands, TrkCuts const &cuts, PIDVeto=PIDVeto::NONE)
Configuration const & cfg_
EleTkIsolFromCands(Configuration const &cfg, reco::TrackCollection const &tracks)
pat::PackedCandidateCollection const * cands_
EleTkIsolFromCands(Configuration const &cfg, pat::PackedCandidateCollection const &cands, PIDVeto pidVeto=PIDVeto::NONE)
TrackTable const & getPreselectedTracks(bool isBarrel)
TrackTable preselectedTracksWithEndcapCuts_
static bool passPIDVeto(const int pdgId, const EleTkIsolFromCands::PIDVeto pidVeto)