CMS 3D CMS Logo

PatTrackAnalyzer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cmath>
3 #include <vector>
4 #include <string>
5 
6 #include <TH1.h>
7 #include <TProfile.h>
8 
16 
18 
24 
25 class PatTrackAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
26 public:
29  ~PatTrackAnalyzer() override;
30 
31  // virtual methods called from base class EDAnalyzer
32  void beginJob() override;
33  void analyze(const edm::Event &event, const edm::EventSetup &es) override;
34 
35 private:
36  // configuration parameters
40 
41  // the list of track quality cuts to demand from the tracking
42  std::vector<std::string> qualities_;
43 
44  // holder for the histograms, one set per quality flag
45  struct Plots {
46  TH1 *eta, *phi;
47  TH1 *pt, *ptErr;
48  TH1 *invPt, *invPtErr;
49  TH1 *d0, *d0Err;
50  TH1 *nHits;
51 
52  TProfile *pxbHitsEta, *pxeHitsEta;
53  TProfile *tibHitsEta, *tobHitsEta;
54  TProfile *tidHitsEta, *tecHitsEta;
55  };
56 
57  std::vector<Plots> plots_;
58 };
59 
61  : srcTracksToken_(mayConsume<edm::View<reco::Track> >(params.getParameter<edm::InputTag>("src"))),
62  srcMuonsToken_(mayConsume<pat::MuonCollection>(params.getParameter<edm::InputTag>("src"))),
63  beamSpotToken_(consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"))),
64  qualities_(params.getParameter<std::vector<std::string> >("qualities")) {
65  usesResource(TFileService::kSharedResource);
66 }
67 
69 
71  // retrieve handle to auxiliary service
72  // used for storing histograms into ROOT file
74 
75  // now book the histograms, for each category
76  unsigned int nQualities = qualities_.size();
77 
78  plots_.resize(nQualities);
79 
80  for (unsigned int i = 0; i < nQualities; ++i) {
81  // the name of the quality flag
82  const char *quality = qualities_[i].c_str();
83 
84  // the set of plots
85  Plots &plots = plots_[i];
86 
87  plots.eta = fs->make<TH1F>(Form("eta_%s", quality), Form("track \\eta (%s)", quality), 100, -3, 3);
88  plots.phi = fs->make<TH1F>(Form("phi_%s", quality), Form("track \\phi (%s)", quality), 100, -M_PI, +M_PI);
89  plots.pt = fs->make<TH1F>(Form("pt_%s", quality), Form("track p_{T} (%s)", quality), 100, 0, 10);
90  plots.ptErr = fs->make<TH1F>(Form("ptErr_%s", quality), Form("track p_{T} error (%s)", quality), 100, 0, 1);
91  plots.invPt = fs->make<TH1F>(Form("invPt_%s", quality), Form("track 1/p_{T} (%s)", quality), 100, -5, 5);
92  plots.invPtErr =
93  fs->make<TH1F>(Form("invPtErr_%s", quality), Form("track 1/p_{T} error (%s)", quality), 100, 0, 0.1);
94  plots.d0 = fs->make<TH1F>(Form("d0_%s", quality), Form("track d0 (%s)", quality), 100, 0, 0.1);
95  plots.d0Err = fs->make<TH1F>(Form("d0Err_%s", quality), Form("track d0 error (%s)", quality), 100, 0, 0.1);
96  plots.nHits =
97  fs->make<TH1F>(Form("nHits_%s", quality), Form("track number of total hits (%s)", quality), 60, 0, 60);
98 
99  plots.pxbHitsEta =
100  fs->make<TProfile>(Form("pxbHitsEta_%s", quality), Form("#hits in Pixel Barrel (%s)", quality), 100, 0, 3);
101  plots.pxeHitsEta =
102  fs->make<TProfile>(Form("pxeHitsEta_%s", quality), Form("#hits in Pixel Endcap (%s)", quality), 100, 0, 3);
103  plots.tibHitsEta = fs->make<TProfile>(
104  Form("tibHitsEta_%s", quality), Form("#hits in Tracker Inner Barrel (%s)", quality), 100, 0, 3);
105  plots.tobHitsEta = fs->make<TProfile>(
106  Form("tobHitsEta_%s", quality), Form("#hits in Tracker Outer Barrel (%s)", quality), 100, 0, 3);
107  plots.tidHitsEta = fs->make<TProfile>(
108  Form("tidHitsEta_%s", quality), Form("#hits in Tracker Inner Disk (%s)", quality), 100, 0, 3);
109  plots.tecHitsEta =
110  fs->make<TProfile>(Form("tecHitsEta_%s", quality), Form("#hits in Tracker Endcap (%s)", quality), 100, 0, 3);
111  }
112 }
113 
115  // handles to kinds of data we might want to read
117  edm::Handle<edm::View<reco::Track> > tracksHandle;
119 
120  // read the beam spot
121  event.getByToken(beamSpotToken_, beamSpot);
122 
123  // our internal copy of track points
124  // (we need this in order to able to simultaneously access tracks
125  // directly or embedded in PAT objects, like muons, normally you
126  // would iterate over the handle directly)
127  std::vector<const reco::Track *> tracks;
128 
129  event.getByToken(srcTracksToken_, tracksHandle);
130  if (tracksHandle.isValid()) {
131  // framework was able to read the collection as a view of
132  // tracks, no copy them to our "tracks" variable
133  for (edm::View<reco::Track>::const_iterator iter = tracksHandle->begin(); iter != tracksHandle->end(); ++iter)
134  tracks.push_back(&*iter);
135  } else {
136  // does not exist or is not a track collection
137  // let's assume it is a collection of PAT muons
138  event.getByToken(srcMuonsToken_, muonsHandle);
139 
140  // and copy them over
141  // NOTE: We are using ->globalTrack() here
142  // This means we are using the global fit over both
143  // the inner tracker and the muon stations!
144  // other alternatives are: innerTrack(), outerTrack()
145  for (pat::MuonCollection::const_iterator iter = muonsHandle->begin(); iter != muonsHandle->end(); ++iter) {
146  reco::TrackRef track = iter->globalTrack();
147  // the muon might not be a "global" muon
148  if (track.isNonnull())
149  tracks.push_back(&*track);
150  }
151  }
152 
153  // we are done filling the tracks into our "tracks" vector.
154  // now analyze them, once for each track quality category
155 
156  unsigned int nQualities = qualities_.size();
157  for (unsigned int i = 0; i < nQualities; ++i) {
158  // we convert the quality flag from its name as a string
159  // to the enumeration value used by the tracking code
160  // (which is essentially an integer number)
162 
163  // our set of plots
164  Plots &plots = plots_[i];
165 
166  // now loop over the tracks
167  for (std::vector<const reco::Track *>::const_iterator iter = tracks.begin(); iter != tracks.end(); ++iter) {
168  // this is our track
169  const reco::Track &track = **iter;
170 
171  // ignore tracks that fail the quality cut
172  if (!track.quality(quality))
173  continue;
174 
175  // and fill all the plots
176  plots.eta->Fill(track.eta());
177  plots.phi->Fill(track.phi());
178 
179  plots.pt->Fill(track.pt());
180  plots.ptErr->Fill(track.ptError());
181 
182  plots.invPt->Fill(track.qoverp());
183  plots.invPtErr->Fill(track.qoverpError());
184 
185  // the transverse IP is taken with respect to
186  // the beam spot instead of (0, 0)
187  // because the beam spot in CMS is not at (0, 0)
188  plots.d0->Fill(track.dxy(beamSpot->position()));
189  plots.d0Err->Fill(track.dxyError());
190 
191  plots.nHits->Fill(track.numberOfValidHits());
192 
193  // the hit pattern contains information about
194  // which modules of the detector have been hit
195  const reco::HitPattern &hits = track.hitPattern();
196 
197  double absEta = std::abs(track.eta());
198  // now fill the number of hits in a layer depending on eta
199  plots.pxbHitsEta->Fill(absEta, hits.numberOfValidPixelBarrelHits());
200  plots.pxeHitsEta->Fill(absEta, hits.numberOfValidPixelEndcapHits());
201  plots.tibHitsEta->Fill(absEta, hits.numberOfValidStripTIBHits());
202  plots.tobHitsEta->Fill(absEta, hits.numberOfValidStripTOBHits());
203  plots.tidHitsEta->Fill(absEta, hits.numberOfValidStripTIDHits());
204  plots.tecHitsEta->Fill(absEta, hits.numberOfValidStripTECHits());
205  }
206  }
207 }
208 
210 
static const std::string kSharedResource
Definition: TFileService.h:76
void beginJob() override
~PatTrackAnalyzer() override
TrackQuality
track quality
Definition: TrackBase.h:150
std::vector< Plots > plots_
std::vector< std::string > qualities_
edm::EDGetTokenT< pat::MuonCollection > srcMuonsToken_
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
Definition: HeavyIon.h:7
string quality
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PatTrackAnalyzer(const edm::ParameterSet &params)
constructor and destructor
#define M_PI
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
void analyze(const edm::Event &event, const edm::EventSetup &es) override
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
bool isValid() const
Definition: HandleBase.h:70
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< edm::View< reco::Track > > srcTracksToken_
Definition: event.py:1