CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
26  public:
28  PatTrackAnalyzer(const edm::ParameterSet &params);
30 
31  // virtual methods called from base class EDAnalyzer
32  virtual void beginJob() override;
33  virtual 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 
60 
62  srcTracksToken_(mayConsume<edm::View<reco::Track> >(params.getParameter<edm::InputTag>("src"))),
63  srcMuonsToken_(mayConsume<pat::MuonCollection>(params.getParameter<edm::InputTag>("src"))),
64  beamSpotToken_(consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"))),
65  qualities_(params.getParameter< std::vector<std::string> >("qualities"))
66 {
67 }
68 
70 {
71 }
72 
74 {
75  // retrieve handle to auxiliary service
76  // used for storing histograms into ROOT file
78 
79  // now book the histograms, for each category
80  unsigned int nQualities = qualities_.size();
81 
82  plots_.resize(nQualities);
83 
84  for(unsigned int i = 0; i < nQualities; ++i) {
85  // the name of the quality flag
86  const char *quality = qualities_[i].c_str();
87 
88  // the set of plots
89  Plots &plots = plots_[i];
90 
91  plots.eta = fs->make<TH1F>(Form("eta_%s", quality),
92  Form("track \\eta (%s)", quality),
93  100, -3, 3);
94  plots.phi = fs->make<TH1F>(Form("phi_%s", quality),
95  Form("track \\phi (%s)", quality),
96  100, -M_PI, +M_PI);
97  plots.pt = fs->make<TH1F>(Form("pt_%s", quality),
98  Form("track p_{T} (%s)", quality),
99  100, 0, 10);
100  plots.ptErr = fs->make<TH1F>(Form("ptErr_%s", quality),
101  Form("track p_{T} error (%s)", quality),
102  100, 0, 1);
103  plots.invPt = fs->make<TH1F>(Form("invPt_%s", quality),
104  Form("track 1/p_{T} (%s)", quality),
105  100, -5, 5);
106  plots.invPtErr = fs->make<TH1F>(Form("invPtErr_%s", quality),
107  Form("track 1/p_{T} error (%s)", quality),
108  100, 0, 0.1);
109  plots.d0 = fs->make<TH1F>(Form("d0_%s", quality),
110  Form("track d0 (%s)", quality),
111  100, 0, 0.1);
112  plots.d0Err = fs->make<TH1F>(Form("d0Err_%s", quality),
113  Form("track d0 error (%s)", quality),
114  100, 0, 0.1);
115  plots.nHits = fs->make<TH1F>(Form("nHits_%s", quality),
116  Form("track number of total hits (%s)", quality),
117  60, 0, 60);
118 
119  plots.pxbHitsEta = fs->make<TProfile>(Form("pxbHitsEta_%s", quality),
120  Form("#hits in Pixel Barrel (%s)", quality),
121  100, 0, 3);
122  plots.pxeHitsEta = fs->make<TProfile>(Form("pxeHitsEta_%s", quality),
123  Form("#hits in Pixel Endcap (%s)", quality),
124  100, 0, 3);
125  plots.tibHitsEta = fs->make<TProfile>(Form("tibHitsEta_%s", quality),
126  Form("#hits in Tracker Inner Barrel (%s)", quality),
127  100, 0, 3);
128  plots.tobHitsEta = fs->make<TProfile>(Form("tobHitsEta_%s", quality),
129  Form("#hits in Tracker Outer Barrel (%s)", quality),
130  100, 0, 3);
131  plots.tidHitsEta = fs->make<TProfile>(Form("tidHitsEta_%s", quality),
132  Form("#hits in Tracker Inner Disk (%s)", quality),
133  100, 0, 3);
134  plots.tecHitsEta = fs->make<TProfile>(Form("tecHitsEta_%s", quality),
135  Form("#hits in Tracker Endcap (%s)", quality),
136  100, 0, 3);
137  }
138 }
139 
141 {
142  // handles to kinds of data we might want to read
146 
147  // read the beam spot
148  event.getByToken(beamSpotToken_, beamSpot);
149 
150  // our internal copy of track points
151  // (we need this in order to able to simultaneously access tracks
152  // directly or embedded in PAT objects, like muons, normally you
153  // would iterate over the handle directly)
154  std::vector<const reco::Track*> tracks;
155 
156  event.getByToken(srcTracksToken_, tracksHandle);
157  if (tracksHandle.isValid()) {
158  // framework was able to read the collection as a view of
159  // tracks, no copy them to our "tracks" variable
160  for(edm::View<reco::Track>::const_iterator iter = tracksHandle->begin();
161  iter != tracksHandle->end(); ++iter)
162  tracks.push_back(&*iter);
163  } else {
164  // does not exist or is not a track collection
165  // let's assume it is a collection of PAT muons
166  event.getByToken(srcMuonsToken_, muonsHandle);
167 
168  // and copy them over
169  // NOTE: We are using ->globalTrack() here
170  // This means we are using the global fit over both
171  // the inner tracker and the muon stations!
172  // other alternatives are: innerTrack(), outerTrack()
173  for(pat::MuonCollection::const_iterator iter = muonsHandle->begin();
174  iter != muonsHandle->end(); ++iter) {
175  reco::TrackRef track = iter->globalTrack();
176  // the muon might not be a "global" muon
177  if (track.isNonnull())
178  tracks.push_back(&*track);
179  }
180  }
181 
182  // we are done filling the tracks into our "tracks" vector.
183  // now analyze them, once for each track quality category
184 
185  unsigned int nQualities = qualities_.size();
186  for(unsigned int i = 0; i < nQualities; ++i) {
187  // we convert the quality flag from its name as a string
188  // to the enumeration value used by the tracking code
189  // (which is essentially an integer number)
191 
192  // our set of plots
193  Plots &plots = plots_[i];
194 
195  // now loop over the tracks
196  for(std::vector<const reco::Track*>::const_iterator iter = tracks.begin();
197  iter != tracks.end(); ++iter) {
198  // this is our track
199  const reco::Track &track = **iter;
200 
201  // ignore tracks that fail the quality cut
202  if (!track.quality(quality))
203  continue;
204 
205  // and fill all the plots
206  plots.eta->Fill(track.eta());
207  plots.phi->Fill(track.phi());
208 
209  plots.pt->Fill(track.pt());
210  plots.ptErr->Fill(track.ptError());
211 
212  plots.invPt->Fill(track.qoverp());
213  plots.invPtErr->Fill(track.qoverpError());
214 
215  // the transverse IP is taken with respect to
216  // the beam spot instead of (0, 0)
217  // because the beam spot in CMS is not at (0, 0)
218  plots.d0->Fill(track.dxy(beamSpot->position()));
219  plots.d0Err->Fill(track.dxyError());
220 
221  plots.nHits->Fill(track.numberOfValidHits());
222 
223  // the hit pattern contains information about
224  // which modules of the detector have been hit
225  const reco::HitPattern &hits = track.hitPattern();
226 
227  double absEta = std::abs(track.eta());
228  // now fill the number of hits in a layer depending on eta
229  plots.pxbHitsEta->Fill(absEta, hits.numberOfValidPixelBarrelHits());
230  plots.pxeHitsEta->Fill(absEta, hits.numberOfValidPixelEndcapHits());
231  plots.tibHitsEta->Fill(absEta, hits.numberOfValidStripTIBHits());
232  plots.tobHitsEta->Fill(absEta, hits.numberOfValidStripTOBHits());
233  plots.tidHitsEta->Fill(absEta, hits.numberOfValidStripTIDHits());
234  plots.tecHitsEta->Fill(absEta, hits.numberOfValidStripTECHits());
235  }
236  }
237 }
238 
240 
double qoverp() const
q / p
Definition: TrackBase.h:526
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
virtual void beginJob() override
TrackQuality
track quality
Definition: TrackBase.h:139
double dxyError() const
error on dxy
Definition: TrackBase.h:749
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Plots > plots_
std::vector< std::string > qualities_
edm::EDGetTokenT< pat::MuonCollection > srcMuonsToken_
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:598
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
int numberOfValidStripTOBHits() const
Definition: HitPattern.h:779
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:604
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.h:754
double pt() const
track transverse momentum
Definition: TrackBase.h:574
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:716
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int numberOfValidStripTIDHits() const
Definition: HitPattern.h:774
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:773
int numberOfValidStripTECHits() const
Definition: HitPattern.h:784
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:75
PatTrackAnalyzer(const edm::ParameterSet &params)
constructor and destructor
#define M_PI
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:710
virtual void analyze(const edm::Event &event, const edm::EventSetup &es) override
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:114
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:411
tuple tracks
Definition: testEve_cfg.py:39
int numberOfValidStripTIBHits() const
Definition: HitPattern.h:769
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:463
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.h:759
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::EDGetTokenT< edm::View< reco::Track > > srcTracksToken_
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:544