CMS 3D CMS Logo

EgammaL1TkIsolation.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_EgammaIsolationAlgos_EgammaL1TkIsolation_h
2 #define RecoEgamma_EgammaIsolationAlgos_EgammaL1TkIsolation_h
3 
9 
10 //author S. Harper (RAL/CERN)
11 //based on the work of Swagata Mukherjee and Giulia Sorrentino
12 
14 public:
15  explicit EgammaL1TkIsolation(const edm::ParameterSet& para);
16 
21  return desc;
22  }
23 
24  std::pair<int, double> calIsol(const reco::TrackBase& trk, const L1TrackCollection& l1Tks) const;
25 
26  std::pair<int, double> calIsol(const double objEta,
27  const double objPhi,
28  const double objZ,
29  const L1TrackCollection& l1Tks) const;
30 
31  //little helper function for the two calIsol functions for it to directly return the pt
32  template <typename... Args>
33  double calIsolPt(Args&&... args) const {
34  return calIsol(std::forward<Args>(args)...).second;
35  }
36 
37 private:
38  struct TrkCuts {
39  float minPt;
40  float minDR2;
41  float maxDR2;
42  float minDEta;
43  float maxDZ;
44  explicit TrkCuts(const edm::ParameterSet& para);
46  };
47 
48  size_t etaBinNr(double eta) const;
49  static bool passTrkSel(const L1Track& trk,
50  const double trkPt,
51  const TrkCuts& cuts,
52  const double objEta,
53  const double objPhi,
54  const double objZ);
55 
56  bool useAbsEta_;
57  std::vector<double> etaBoundaries_;
58  std::vector<TrkCuts> trkCuts_;
59 };
60 
61 #endif
double calIsolPt(Args &&... args) const
std::vector< double > etaBoundaries_
static bool passTrkSel(const L1Track &trk, const double trkPt, const TrkCuts &cuts, const double objEta, const double objPhi, const double objZ)
static edm::ParameterSetDescription makePSetDescription()
EgammaL1TkIsolation(const edm::ParameterSet &para)
TrkCuts(const edm::ParameterSet &para)
std::pair< int, double > calIsol(const reco::TrackBase &trk, const L1TrackCollection &l1Tks) const
static edm::ParameterSetDescription makePSetDescription()
std::vector< TrkCuts > trkCuts_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
size_t etaBinNr(double eta) const
std::vector< L1Track > L1TrackCollection
Definition: L1Track.h:9