CMS 3D CMS Logo

TrackIsoCalculator.cc
Go to the documentation of this file.
2 
4 
5 using namespace edm;
6 using namespace reco;
7 using namespace std;
8 
10  : recCollection_{trackCollection}, trackQuality_{trackQuality} {}
11 
13  const double x,
14  const double threshold,
15  const double innerDR) {
16  double totalPt = 0;
17 
18  for (auto const& recTrack : recCollection_) {
19  bool goodtrack = recTrack.quality(reco::TrackBase::qualityByName(trackQuality_));
20  if (!goodtrack)
21  continue;
22 
23  double pt = recTrack.pt();
24  double dR2 = reco::deltaR2(cluster, recTrack);
25  if (dR2 >= (0.01 * x * x))
26  continue;
27  if (dR2 < innerDR * innerDR)
28  continue;
29  if (pt > threshold)
30  totalPt = totalPt + pt;
31  }
32 
33  return totalPt;
34 }
35 
37  const double x,
38  const double threshold,
39  const double innerDR) {
40  double SClusterEta = cluster.eta();
41  double totalPt = 0.0;
42 
43  for (auto const& recTrack : recCollection_) {
44  bool goodtrack = recTrack.quality(reco::TrackBase::qualityByName(trackQuality_));
45  if (!goodtrack)
46  continue;
47 
48  double pt = recTrack.pt();
49  if (std::abs(recTrack.eta() - SClusterEta) >= 0.1 * x)
50  continue;
51  if (reco::deltaR2(cluster, recTrack) < innerDR * innerDR)
52  continue;
53 
54  if (pt > threshold)
55  totalPt = totalPt + pt;
56  }
57 
58  double Tx = getTrackIso(cluster, x, threshold, innerDR);
59  double CTx = (Tx - totalPt / 40.0 * x) * (1 / (1 - x / 40.));
60 
61  return CTx;
62 }
TrackIsoCalculator(reco::TrackCollection const &trackCollection, std::string const &trackQuality)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
reco::TrackCollection const & recCollection_
std::string const & trackQuality_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
double getBkgSubTrackIso(reco::Photon const &clus, const double i, const double threshold, const double innerDR=0)
Return the background-subtracted tracker energy in a cone around the photon.
fixed size matrix
HLT enums.
float x
double getTrackIso(reco::Photon const &clus, const double i, const double threshold, const double innerDR=0)
Return the tracker energy in a cone around the photon.
double eta() const final
momentum pseudorapidity