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