CMS 3D CMS Logo

GeneralPurposeTrackAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/OfflineValidation
4 // Class: GeneralPurposeTrackAnalyzer
5 //
9 //
10 // Original Author: Marco Musich
11 // Created: Mon, 13 Jun 2016 15:07:11 GMT
12 //
13 //
14 
15 // ROOT includes files
16 
17 #include "TMath.h"
18 #include "TFile.h"
19 #include "TH1D.h"
20 #include "TH1I.h"
21 #include "TH2D.h"
22 #include "TProfile.h"
23 
24 // system includes files
25 
26 #include <cmath>
27 #include <fstream>
28 #include <iostream>
29 #include <map>
30 #include <string>
31 #include <vector>
32 #include <boost/range/adaptor/indexed.hpp>
33 
34 // user include files
82 
83 // toggle to enable debugging
84 #define DEBUG 0
85 
88 
89 class GeneralPurposeTrackAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
90 public:
93  magFieldToken_(esConsumes<edm::Transition::BeginRun>()),
94  geomTokenBR_(esConsumes<edm::Transition::BeginRun>()),
97  doLatencyAnalysis_ = pset.getParameter<bool>("doLatencyAnalysis");
98  if (doLatencyAnalysis_) {
99  latencyToken_ = esConsumes<edm::Transition::BeginRun>();
100  }
101 
102  usesResource(TFileService::kSharedResource);
103 
104  TkTag_ = pset.getParameter<edm::InputTag>("TkTag");
105  theTrackCollectionToken_ = consumes<reco::TrackCollection>(TkTag_);
106 
107  TriggerResultsTag_ = pset.getParameter<edm::InputTag>("TriggerResultsTag");
108  hltresultsToken_ = consumes<edm::TriggerResults>(TriggerResultsTag_);
109 
110  BeamSpotTag_ = pset.getParameter<edm::InputTag>("BeamSpotTag");
111  beamspotToken_ = consumes<reco::BeamSpot>(BeamSpotTag_);
112 
113  VerticesTag_ = pset.getParameter<edm::InputTag>("VerticesTag");
114  vertexToken_ = consumes<reco::VertexCollection>(VerticesTag_);
115 
116  isCosmics_ = pset.getParameter<bool>("isCosmics");
117 
118  pmap = std::make_unique<TrackerMap>("Pixel");
119  pmap->onlyPixel(true);
120  pmap->setTitle("Pixel Hit entries");
121  pmap->setPalette(1);
122 
123  tmap = std::make_unique<TrackerMap>("Strip");
124  tmap->setTitle("Strip Hit entries");
125  tmap->setPalette(1);
126 
127  pixelmap = std::make_unique<Phase1PixelMaps>("COLZ0 L");
128  pixelmap->bookBarrelHistograms("entriesBarrel", "# hits", "# pixel hits");
129  pixelmap->bookForwardHistograms("entriesForward", "# hits", "# pixel hits");
130 
131  pixelrocsmap_ = std::make_unique<Phase1PixelROCMaps>("");
132  }
133 
134  ~GeneralPurposeTrackAnalyzer() override = default;
135 
137 
138  template <class OBJECT_TYPE>
139  int index(const std::vector<OBJECT_TYPE *> &vec, const TString &name) {
140  for (const auto &iter : vec | boost::adaptors::indexed(0)) {
141  if (iter.value() && iter.value()->GetName() == name) {
142  return iter.index();
143  }
144  }
145  edm::LogError("GeneralPurposeTrackAnalyzer") << "@SUB=GeneralPurposeTrackAnalyzer::index"
146  << " could not find " << name;
147  return -1;
148  }
149 
150  template <typename T, typename... Args>
151  T *book(const Args &...args) const {
152  T *t = fs->make<T>(args...);
153  return t;
154  }
155 
156 private:
157  // tokens for the event setup
161 
165 
167 
169 
171 
172  std::unique_ptr<TrackerMap> tmap;
173  std::unique_ptr<TrackerMap> pmap;
174 
175  std::unique_ptr<Phase1PixelMaps> pixelmap;
176  std::unique_ptr<Phase1PixelROCMaps> pixelrocsmap_;
177 
178  TH1D *hchi2ndof;
179  TH1D *hNtrk;
180  TH1D *hNtrkZoom;
181  TH1I *htrkQuality;
182  TH1I *htrkAlgo;
183  TH1I *htrkOriAlgo;
185  TH1D *hP;
186  TH1D *hPPlus;
187  TH1D *hPMinus;
188  TH1D *hPt;
189  TH1D *hMinPt;
190  TH1D *hPtPlus;
191  TH1D *hPtMinus;
192  TH1D *hHit;
193  TH1D *hHit2D;
194 
201 
204 
207 
208  TH1D *hHitPlus;
209  TH1D *hHitMinus;
210 
211  TH1D *hPhp;
212  TH1D *hPthp;
213  TH1D *hHithp;
214  TH1D *hEtahp;
215  TH1D *hPhihp;
216  TH1D *hchi2ndofhp;
217  TH1D *hchi2Probhp;
218 
219  TH1D *hCharge;
220  TH1D *hQoverP;
221  TH1D *hQoverPZoom;
222  TH1D *hEta;
223  TH1D *hEtaPlus;
224  TH1D *hEtaMinus;
225  TH1D *hPhi;
226  TH1D *hPhiBarrel;
231  TH1D *hPhiPlus;
232  TH1D *hPhiMinus;
233 
234  TH1D *hDeltaPhi;
235  TH1D *hDeltaEta;
236  TH1D *hDeltaR;
237 
238  TH1D *hvx;
239  TH1D *hvy;
240  TH1D *hvz;
241  TH1D *hd0;
242  TH1D *hdz;
243  TH1D *hdxy;
244 
245  TH2D *hd0PVvsphi;
246  TH2D *hd0PVvseta;
247  TH2D *hd0PVvspt;
248 
249  TH2D *hd0vsphi;
250  TH2D *hd0vseta;
251  TH2D *hd0vspt;
252 
253  TH1D *hnhpxb;
254  TH1D *hnhpxe;
255  TH1D *hnhTIB;
256  TH1D *hnhTID;
257  TH1D *hnhTOB;
258  TH1D *hnhTEC;
259 
261 
262  TProfile *pNBpixHitsVsVx;
263  TProfile *pNBpixHitsVsVy;
264  TProfile *pNBpixHitsVsVz;
265 
266  TH1D *hMultCand;
267 
268  TH1D *hdxyBS;
269  TH1D *hd0BS;
270  TH1D *hdzBS;
271  TH1D *hdxyPV;
272  TH1D *hd0PV;
273  TH1D *hdzPV;
274  TH1D *hrun;
275  TH1D *hlumi;
276 
277  TH1F *h_BSx0;
278  TH1F *h_BSy0;
279  TH1F *h_BSz0;
283  TH1F *h_BSdxdz;
284  TH1F *h_BSdydz;
285 
286  std::vector<TH1 *> vTrackHistos_;
287  std::vector<TH1 *> vTrackProfiles_;
288  std::vector<TH1 *> vTrack2DHistos_;
289 
292 
293  TH1D *modeByRun_;
294  TH1D *fieldByRun_;
295 
296  int ievt;
297  int itrks;
298  int mode;
299  float etaMax_;
301 
303 
308 
311 
316 
317  std::map<std::string, std::pair<int, int> > triggerMap_;
318  std::map<int, std::pair<int, float> > conditionsMap_;
319  std::map<int, std::pair<int, int> > runInfoMap_;
320 
321  //*************************************************************
322  void analyze(const edm::Event &event, const edm::EventSetup &setup) override
323  //*************************************************************
324  {
325  ievt++;
326 
327  // geometry setup
328  const TrackerGeometry *theGeometry = &setup.getData(geomToken_);
329 
331  const reco::TrackCollection tC = *(trackCollection.product());
332  itrks += tC.size();
333 
334  runInfoMap_[event.run()].first += 1;
335  runInfoMap_[event.run()].second += tC.size();
336 
337  if (DEBUG) {
338  edm::LogInfo("GeneralPurposeTrackAnalyzer") << "Reconstructed " << tC.size() << " tracks" << std::endl;
339  }
340  //int iCounter=0;
341 
343  if (hltresults.isValid()) {
344  const edm::TriggerNames &triggerNames_ = event.triggerNames(*hltresults);
345  int ntrigs = hltresults->size();
346  //const vector<std::string> &triggernames = triggerNames_.triggerNames();
347 
348  for (int itrig = 0; itrig != ntrigs; ++itrig) {
349  const std::string &trigName = triggerNames_.triggerName(itrig);
350  bool accept = hltresults->accept(itrig);
351  if (accept == 1) {
352  if (DEBUG) {
353  edm::LogInfo("GeneralPurposeTrackAnalyzer")
354  << trigName << " " << accept << " ,track size: " << tC.size() << std::endl;
355  }
356  triggerMap_[trigName].first += 1;
357  triggerMap_[trigName].second += tC.size();
358  // triggerInfo.push_back(pair <std::string, int> (trigName, accept));
359  }
360  }
361  }
362 
363  hrun->Fill(event.run());
364  hlumi->Fill(event.luminosityBlock());
365 
366  int nHighPurityTracks = 0;
367 
368  for (auto track = tC.cbegin(); track != tC.cend(); track++) {
369  unsigned int nHit2D = 0;
370  std::bitset<16> rocsToMask;
371  for (auto iHit = track->recHitsBegin(); iHit != track->recHitsEnd(); ++iHit) {
372  if (this->isHit2D(**iHit)) {
373  ++nHit2D;
374  }
375  // rest the ROCs for the map
376  rocsToMask.reset();
377  const DetId &detId = (*iHit)->geographicalId();
378  const GeomDet *geomDet(theGeometry->idToDet(detId));
379 
380  const SiPixelRecHit *pixhit = dynamic_cast<const SiPixelRecHit *>(*iHit);
381 
382  if (pixhit) {
383  if (pixhit->isValid()) {
384  unsigned int subid = detId.subdetId();
385  int detid_db = detId.rawId();
386 
387  // get the cluster
388  auto clustp = pixhit->cluster();
389 
390  if (clustp.isNull())
391  continue;
392  auto const &cluster = *clustp;
393  int row = cluster.x() - 0.5, col = cluster.y() - 0.5;
394  int rocId = coord_.roc(detId, std::make_pair(row, col));
395 
397  pmap->fill(detid_db, 1);
398  } else if (phase_ == SiPixelPI::phase::one) {
399  rocsToMask.set(rocId);
400  pixelrocsmap_->fillSelectedRocs(detid_db, rocsToMask, 1);
401 
402  if (subid == PixelSubdetector::PixelBarrel) {
403  pixelmap->fillBarrelBin("entriesBarrel", detid_db, 1);
404  } else {
405  pixelmap->fillForwardBin("entriesForward", detid_db, 1);
406  }
407  }
408 
409  LocalPoint lp = (*iHit)->localPosition();
410  //LocalError le = (*iHit)->localPositionError();
411 
412  GlobalPoint GP = geomDet->surface().toGlobal(lp);
413 
414  if ((subid == PixelSubdetector::PixelBarrel) || (subid == PixelSubdetector::PixelEndcap)) {
415  // 1 = PXB, 2 = PXF
416  if (subid == PixelSubdetector::PixelBarrel) {
417  hHitCountVsThetaBPix->Fill(GP.theta());
418  hHitCountVsPhiBPix->Fill(GP.phi());
419 
420  hHitCountVsZBPix->Fill(GP.z());
421  hHitCountVsXBPix->Fill(GP.x());
422  hHitCountVsYBPix->Fill(GP.y());
423 
424  } else if (subid == PixelSubdetector::PixelEndcap) {
425  hHitCountVsThetaFPix->Fill(GP.theta());
426  hHitCountVsPhiFPix->Fill(GP.phi());
427 
428  hHitCountVsZFPix->Fill(GP.z());
429  hHitCountVsXFPix->Fill(GP.x());
430  hHitCountVsYFPix->Fill(GP.y());
431  }
432  }
433  }
434  } else {
435  if ((*iHit)->isValid() && phase_ != SiPixelPI::phase::two) {
436  tmap->fill(detId.rawId(), 1);
437  }
438  }
439  }
440 
441  hHit2D->Fill(nHit2D);
442  hHit->Fill(track->numberOfValidHits());
443  hnhpxb->Fill(track->hitPattern().numberOfValidPixelBarrelHits());
444  hnhpxe->Fill(track->hitPattern().numberOfValidPixelEndcapHits());
445  hnhTIB->Fill(track->hitPattern().numberOfValidStripTIBHits());
446  hnhTID->Fill(track->hitPattern().numberOfValidStripTIDHits());
447  hnhTOB->Fill(track->hitPattern().numberOfValidStripTOBHits());
448  hnhTEC->Fill(track->hitPattern().numberOfValidStripTECHits());
449 
450  // fill hit composition histogram
451  if (track->hitPattern().numberOfValidPixelBarrelHits() != 0) {
452  hHitComposition->Fill(0., track->hitPattern().numberOfValidPixelBarrelHits());
453 
454  pNBpixHitsVsVx->Fill(track->vx(), track->hitPattern().numberOfValidPixelBarrelHits());
455  pNBpixHitsVsVy->Fill(track->vy(), track->hitPattern().numberOfValidPixelBarrelHits());
456  pNBpixHitsVsVz->Fill(track->vz(), track->hitPattern().numberOfValidPixelBarrelHits());
457  }
458  if (track->hitPattern().numberOfValidPixelEndcapHits() != 0) {
459  hHitComposition->Fill(1., track->hitPattern().numberOfValidPixelEndcapHits());
460  }
461  if (track->hitPattern().numberOfValidStripTIBHits() != 0) {
462  hHitComposition->Fill(2., track->hitPattern().numberOfValidStripTIBHits());
463  }
464  if (track->hitPattern().numberOfValidStripTIDHits() != 0) {
465  hHitComposition->Fill(3., track->hitPattern().numberOfValidStripTIDHits());
466  }
467  if (track->hitPattern().numberOfValidStripTOBHits() != 0) {
468  hHitComposition->Fill(4., track->hitPattern().numberOfValidStripTOBHits());
469  }
470  if (track->hitPattern().numberOfValidStripTECHits() != 0) {
471  hHitComposition->Fill(5., track->hitPattern().numberOfValidStripTECHits());
472  }
473 
474  hCharge->Fill(track->charge());
475  hQoverP->Fill(track->qoverp());
476  hQoverPZoom->Fill(track->qoverp());
477  hPt->Fill(track->pt());
478  hP->Fill(track->p());
479  hchi2ndof->Fill(track->normalizedChi2());
480  hEta->Fill(track->eta());
481  hPhi->Fill(track->phi());
482 
483  if (fabs(track->eta()) < 0.8) {
484  hPhiBarrel->Fill(track->phi());
485  }
486  if (track->eta() > 0.8 && track->eta() < 1.4) {
487  hPhiOverlapPlus->Fill(track->phi());
488  }
489  if (track->eta() < -0.8 && track->eta() > -1.4) {
490  hPhiOverlapMinus->Fill(track->phi());
491  }
492  if (track->eta() > 1.4) {
493  hPhiEndcapPlus->Fill(track->phi());
494  }
495  if (track->eta() < -1.4) {
496  hPhiEndcapMinus->Fill(track->phi());
497  }
498 
499  hd0->Fill(track->d0());
500  hdz->Fill(track->dz());
501  hdxy->Fill(track->dxy());
502  hvx->Fill(track->vx());
503  hvy->Fill(track->vy());
504  hvz->Fill(track->vz());
505 
506  htrkAlgo->Fill(static_cast<int>(track->algo()));
507  htrkOriAlgo->Fill(static_cast<int>(track->originalAlgo()));
508 
509  int myquality = -99;
510  if (track->quality(reco::TrackBase::undefQuality)) {
511  myquality = -1;
512  htrkQuality->Fill(myquality);
513  }
514  if (track->quality(reco::TrackBase::loose)) {
515  myquality = 0;
516  htrkQuality->Fill(myquality);
517  }
518  if (track->quality(reco::TrackBase::tight)) {
519  myquality = 1;
520  htrkQuality->Fill(myquality);
521  }
522  if (track->quality(reco::TrackBase::highPurity) && (!isCosmics_)) {
523  myquality = 2;
524  htrkQuality->Fill(myquality);
525  hPhp->Fill(track->p());
526  hPthp->Fill(track->pt());
527  hHithp->Fill(track->numberOfValidHits());
528  hEtahp->Fill(track->eta());
529  hPhihp->Fill(track->phi());
530  hchi2ndofhp->Fill(track->normalizedChi2());
531  hchi2Probhp->Fill(TMath::Prob(track->chi2(), track->ndof()));
532  nHighPurityTracks++;
533  }
534  if (track->quality(reco::TrackBase::confirmed)) {
535  myquality = 3;
536  htrkQuality->Fill(myquality);
537  }
538  if (track->quality(reco::TrackBase::goodIterative)) {
539  myquality = 4;
540  htrkQuality->Fill(myquality);
541  }
542 
543  // Fill 1D track histos
544  static const int etaindex = this->index(vTrackHistos_, "h_tracketa");
545  vTrackHistos_[etaindex]->Fill(track->eta());
546  static const int phiindex = this->index(vTrackHistos_, "h_trackphi");
547  vTrackHistos_[phiindex]->Fill(track->phi());
548  static const int numOfValidHitsindex = this->index(vTrackHistos_, "h_trackNumberOfValidHits");
549  vTrackHistos_[numOfValidHitsindex]->Fill(track->numberOfValidHits());
550  static const int numOfLostHitsindex = this->index(vTrackHistos_, "h_trackNumberOfLostHits");
551  vTrackHistos_[numOfLostHitsindex]->Fill(track->numberOfLostHits());
552 
553  GlobalPoint gPoint(track->vx(), track->vy(), track->vz());
554  double theLocalMagFieldInInverseGeV = magneticField_->inInverseGeV(gPoint).z();
555  double kappa = -track->charge() * theLocalMagFieldInInverseGeV / track->pt();
556 
557  static const int kappaindex = this->index(vTrackHistos_, "h_curvature");
558  vTrackHistos_[kappaindex]->Fill(kappa);
559  static const int kappaposindex = this->index(vTrackHistos_, "h_curvature_pos");
560  if (track->charge() > 0)
561  vTrackHistos_[kappaposindex]->Fill(fabs(kappa));
562  static const int kappanegindex = this->index(vTrackHistos_, "h_curvature_neg");
563  if (track->charge() < 0)
564  vTrackHistos_[kappanegindex]->Fill(fabs(kappa));
565 
566  double chi2Prob = TMath::Prob(track->chi2(), track->ndof());
567  double normchi2 = track->normalizedChi2();
568 
569  static const int normchi2index = this->index(vTrackHistos_, "h_normchi2");
570  vTrackHistos_[normchi2index]->Fill(normchi2);
571  static const int chi2index = this->index(vTrackHistos_, "h_chi2");
572  vTrackHistos_[chi2index]->Fill(track->chi2());
573  static const int chi2Probindex = this->index(vTrackHistos_, "h_chi2Prob");
574  vTrackHistos_[chi2Probindex]->Fill(chi2Prob);
575  static const int ptindex = this->index(vTrackHistos_, "h_pt");
576  static const int pt2index = this->index(vTrackHistos_, "h_ptrebin");
577  vTrackHistos_[ptindex]->Fill(track->pt());
578  vTrackHistos_[pt2index]->Fill(track->pt());
579  if (track->ptError() != 0.) {
580  static const int ptResolutionindex = this->index(vTrackHistos_, "h_ptResolution");
581  vTrackHistos_[ptResolutionindex]->Fill(track->ptError() / track->pt());
582  }
583  // Fill track profiles
584  static const int d0phiindex = this->index(vTrackProfiles_, "p_d0_vs_phi");
585  vTrackProfiles_[d0phiindex]->Fill(track->phi(), track->d0());
586  static const int dzphiindex = this->index(vTrackProfiles_, "p_dz_vs_phi");
587  vTrackProfiles_[dzphiindex]->Fill(track->phi(), track->dz());
588  static const int d0etaindex = this->index(vTrackProfiles_, "p_d0_vs_eta");
589  vTrackProfiles_[d0etaindex]->Fill(track->eta(), track->d0());
590  static const int dzetaindex = this->index(vTrackProfiles_, "p_dz_vs_eta");
591  vTrackProfiles_[dzetaindex]->Fill(track->eta(), track->dz());
592  static const int chiProbphiindex = this->index(vTrackProfiles_, "p_chi2Prob_vs_phi");
593  vTrackProfiles_[chiProbphiindex]->Fill(track->phi(), chi2Prob);
594  static const int chiProbabsd0index = this->index(vTrackProfiles_, "p_chi2Prob_vs_d0");
595  vTrackProfiles_[chiProbabsd0index]->Fill(fabs(track->d0()), chi2Prob);
596  static const int chiProbabsdzindex = this->index(vTrackProfiles_, "p_chi2Prob_vs_dz");
597  vTrackProfiles_[chiProbabsdzindex]->Fill(track->dz(), chi2Prob);
598  static const int chiphiindex = this->index(vTrackProfiles_, "p_chi2_vs_phi");
599  vTrackProfiles_[chiphiindex]->Fill(track->phi(), track->chi2());
600  static const int normchiphiindex = this->index(vTrackProfiles_, "p_normchi2_vs_phi");
601  vTrackProfiles_[normchiphiindex]->Fill(track->phi(), normchi2);
602  static const int chietaindex = this->index(vTrackProfiles_, "p_chi2_vs_eta");
603  vTrackProfiles_[chietaindex]->Fill(track->eta(), track->chi2());
604  static const int normchiptindex = this->index(vTrackProfiles_, "p_normchi2_vs_pt");
605  vTrackProfiles_[normchiptindex]->Fill(track->pt(), normchi2);
606  static const int normchipindex = this->index(vTrackProfiles_, "p_normchi2_vs_p");
607  vTrackProfiles_[normchipindex]->Fill(track->p(), normchi2);
608  static const int chiProbetaindex = this->index(vTrackProfiles_, "p_chi2Prob_vs_eta");
609  vTrackProfiles_[chiProbetaindex]->Fill(track->eta(), chi2Prob);
610  static const int normchietaindex = this->index(vTrackProfiles_, "p_normchi2_vs_eta");
611  vTrackProfiles_[normchietaindex]->Fill(track->eta(), normchi2);
612  static const int kappaphiindex = this->index(vTrackProfiles_, "p_kappa_vs_phi");
613  vTrackProfiles_[kappaphiindex]->Fill(track->phi(), kappa);
614  static const int kappaetaindex = this->index(vTrackProfiles_, "p_kappa_vs_eta");
615  vTrackProfiles_[kappaetaindex]->Fill(track->eta(), kappa);
616  static const int ptResphiindex = this->index(vTrackProfiles_, "p_ptResolution_vs_phi");
617  vTrackProfiles_[ptResphiindex]->Fill(track->phi(), track->ptError() / track->pt());
618  static const int ptResetaindex = this->index(vTrackProfiles_, "p_ptResolution_vs_eta");
619  vTrackProfiles_[ptResetaindex]->Fill(track->eta(), track->ptError() / track->pt());
620 
621  // Fill 2D track histos
622  static const int etaphiindex_2d = this->index(vTrack2DHistos_, "h2_phi_vs_eta");
623  vTrack2DHistos_[etaphiindex_2d]->Fill(track->eta(), track->phi());
624  static const int d0phiindex_2d = this->index(vTrack2DHistos_, "h2_d0_vs_phi");
625  vTrack2DHistos_[d0phiindex_2d]->Fill(track->phi(), track->d0());
626  static const int dzphiindex_2d = this->index(vTrack2DHistos_, "h2_dz_vs_phi");
627  vTrack2DHistos_[dzphiindex_2d]->Fill(track->phi(), track->dz());
628  static const int d0etaindex_2d = this->index(vTrack2DHistos_, "h2_d0_vs_eta");
629  vTrack2DHistos_[d0etaindex_2d]->Fill(track->eta(), track->d0());
630  static const int dzetaindex_2d = this->index(vTrack2DHistos_, "h2_dz_vs_eta");
631  vTrack2DHistos_[dzetaindex_2d]->Fill(track->eta(), track->dz());
632  static const int chiphiindex_2d = this->index(vTrack2DHistos_, "h2_chi2_vs_phi");
633  vTrack2DHistos_[chiphiindex_2d]->Fill(track->phi(), track->chi2());
634  static const int chiProbphiindex_2d = this->index(vTrack2DHistos_, "h2_chi2Prob_vs_phi");
635  vTrack2DHistos_[chiProbphiindex_2d]->Fill(track->phi(), chi2Prob);
636  static const int chiProbabsd0index_2d = this->index(vTrack2DHistos_, "h2_chi2Prob_vs_d0");
637  vTrack2DHistos_[chiProbabsd0index_2d]->Fill(fabs(track->d0()), chi2Prob);
638  static const int normchiphiindex_2d = this->index(vTrack2DHistos_, "h2_normchi2_vs_phi");
639  vTrack2DHistos_[normchiphiindex_2d]->Fill(track->phi(), normchi2);
640  static const int chietaindex_2d = this->index(vTrack2DHistos_, "h2_chi2_vs_eta");
641  vTrack2DHistos_[chietaindex_2d]->Fill(track->eta(), track->chi2());
642  static const int chiProbetaindex_2d = this->index(vTrack2DHistos_, "h2_chi2Prob_vs_eta");
643  vTrack2DHistos_[chiProbetaindex_2d]->Fill(track->eta(), chi2Prob);
644  static const int normchietaindex_2d = this->index(vTrack2DHistos_, "h2_normchi2_vs_eta");
645  vTrack2DHistos_[normchietaindex_2d]->Fill(track->eta(), normchi2);
646  static const int kappaphiindex_2d = this->index(vTrack2DHistos_, "h2_kappa_vs_phi");
647  vTrack2DHistos_[kappaphiindex_2d]->Fill(track->phi(), kappa);
648  static const int kappaetaindex_2d = this->index(vTrack2DHistos_, "h2_kappa_vs_eta");
649  vTrack2DHistos_[kappaetaindex_2d]->Fill(track->eta(), kappa);
650  static const int normchi2kappa_2d = this->index(vTrack2DHistos_, "h2_normchi2_vs_kappa");
651  vTrack2DHistos_[normchi2kappa_2d]->Fill(normchi2, kappa);
652 
653  //dxy with respect to the beamspot
655  edm::Handle<reco::BeamSpot> beamSpotHandle = event.getHandle(beamspotToken_);
656  if (beamSpotHandle.isValid()) {
657  beamSpot = *beamSpotHandle;
658 
659  double BSx0 = beamSpot.x0();
660  double BSy0 = beamSpot.y0();
661  double BSz0 = beamSpot.z0();
662  double Beamsigmaz = beamSpot.sigmaZ();
663  double Beamdxdz = beamSpot.dxdz();
664  double Beamdydz = beamSpot.dydz();
665  double BeamWidthX = beamSpot.BeamWidthX();
666  double BeamWidthY = beamSpot.BeamWidthY();
667 
668  math::XYZPoint point(BSx0, BSy0, BSz0);
669  double dxy = track->dxy(point);
670  double dz = track->dz(point);
671  hdxyBS->Fill(dxy);
672  hd0BS->Fill(-dxy);
673  hdzBS->Fill(dz);
674 
675  h_BSx0->Fill(BSx0);
676  h_BSy0->Fill(BSy0);
677  h_BSz0->Fill(BSz0);
678  h_Beamsigmaz->Fill(Beamsigmaz);
679  h_BeamWidthX->Fill(BeamWidthX);
680  h_BeamWidthY->Fill(BeamWidthY);
681  h_BSdxdz->Fill(Beamdxdz);
682  h_BSdydz->Fill(Beamdydz);
683  }
684 
685  //dxy with respect to the primary vertex
686  reco::Vertex pvtx;
687  edm::Handle<reco::VertexCollection> vertexHandle = event.getHandle(vertexToken_);
688  double mindxy = 100.;
689  double dz = 100;
690  if (vertexHandle.isValid() && !isCosmics_) {
691  for (auto pvtx = vertexHandle->cbegin(); pvtx != vertexHandle->cend(); ++pvtx) {
692  math::XYZPoint mypoint(pvtx->x(), pvtx->y(), pvtx->z());
693  if (abs(mindxy) > abs(track->dxy(mypoint))) {
694  mindxy = track->dxy(mypoint);
695  dz = track->dz(mypoint);
696  }
697  }
698 
699  hdxyPV->Fill(mindxy);
700  hd0PV->Fill(-mindxy);
701  hdzPV->Fill(dz);
702 
703  hd0PVvsphi->Fill(track->phi(), -mindxy);
704  hd0PVvseta->Fill(track->eta(), -mindxy);
705  hd0PVvspt->Fill(track->pt(), -mindxy);
706 
707  } else {
708  hdxyPV->Fill(100);
709  hd0PV->Fill(100);
710  hdzPV->Fill(100);
711 
712  hd0vsphi->Fill(track->phi(), -track->dxy());
713  hd0vseta->Fill(track->eta(), -track->dxy());
714  hd0vspt->Fill(track->pt(), -track->dxy());
715  }
716 
717  if (DEBUG) {
718  edm::LogInfo("GeneralPurposeTrackAnalyzer") << "end of track loop" << std::endl;
719  }
720  }
721 
722  hNtrk->Fill(tC.size());
723  hNtrkZoom->Fill(tC.size());
724  hNhighPurity->Fill(nHighPurityTracks);
725 
726  if (DEBUG) {
727  edm::LogInfo("GeneralPurposeTrackAnalyzer") << "end of analysis" << std::endl;
728  }
729  }
730 
731  //*************************************************************
732  void beginRun(edm::Run const &run, edm::EventSetup const &setup) override
733  //*************************************************************
734  {
735  // initialize runInfoMap_
736  runInfoMap_[run.run()].first = 0;
737  runInfoMap_[run.run()].second = 0;
738 
739  // Magnetic Field setup
740  magneticField_ = setup.getHandle(magFieldToken_);
741  float B_ = magneticField_.product()->inTesla(GlobalPoint(0, 0, 0)).mag();
742 
743  if (DEBUG) {
744  edm::LogInfo("GeneralPurposeTrackAnalyzer")
745  << "run number:" << run.run() << " magnetic field: " << B_ << " [T]" << std::endl;
746  }
747 
748  const TrackerGeometry *trackerGeometry = &setup.getData(geomTokenBR_);
749  if (trackerGeometry->isThere(GeomDetEnumerators::P2PXB) || trackerGeometry->isThere(GeomDetEnumerators::P2PXEC)) {
751  } else if (trackerGeometry->isThere(GeomDetEnumerators::P1PXB) ||
752  trackerGeometry->isThere(GeomDetEnumerators::P1PXEC)) {
754  } else {
756  }
757 
758  // if it's a phase-2 geometry there are no phase-1 conditions
759  if (phase_ == SiPixelPI::phase::two) {
760  mode = 0;
761  } else {
762  if (doLatencyAnalysis_) {
763  //SiStrip Latency
764  const SiStripLatency *apvlat = &setup.getData(latencyToken_);
765  if (apvlat->singleReadOutMode() == 1) {
766  mode = 1; // peak mode
767  } else if (apvlat->singleReadOutMode() == 0) {
768  mode = -1; // deco mode
769  }
770  } else {
771  mode = 0;
772  }
773  }
774 
775  conditionsMap_[run.run()].first = mode;
776  conditionsMap_[run.run()].second = B_;
777 
778  // init the sipixel coordinates
779  const TrackerTopology *trackerTopology = &setup.getData(trackerTopologyTokenBR_);
780  const SiPixelFedCablingMap *siPixelFedCablingMap = &setup.getData(siPixelFedCablingMapTokenBR_);
781 
782  // Pixel Phase-1 helper class
783  coord_.init(trackerTopology, trackerGeometry, siPixelFedCablingMap);
784  }
785 
786  //*************************************************************
787  void endRun(edm::Run const &, edm::EventSetup const &) override {}
788  //*************************************************************
789 
790  //*************************************************************
791  void beginJob() override
792  //*************************************************************
793  {
794  if (DEBUG) {
795  edm::LogInfo("GeneralPurposeTrackAnalyzer") << __LINE__ << std::endl;
796  }
797 
798  TH1D::SetDefaultSumw2(kTRUE);
799  etaMax_ = 3.0;
800 
801  ievt = 0;
802  itrks = 0;
803 
804  hrun = book<TH1D>("h_run", "run", 100000, 230000, 240000);
805  hlumi = book<TH1D>("h_lumi", "lumi", 1000, 0, 1000);
806 
807  hchi2ndof = book<TH1D>("h_chi2ndof", "chi2/ndf;#chi^{2}/ndf;tracks", 100, 0, 5.);
808  hCharge = book<TH1D>("h_charge", "charge;Charge of the track;tracks", 5, -2.5, 2.5);
809  hNtrk = book<TH1D>("h_Ntrk", "ntracks;Number of Tracks;events", 200, 0., 200.);
810  hNtrkZoom = book<TH1D>("h_NtrkZoom", "Number of tracks; number of tracks;events", 10, 0., 10.);
811  hNhighPurity =
812  book<TH1D>("h_NhighPurity", "n. high purity tracks;Number of high purity tracks;events", 200, 0., 200.);
813 
814  htrkAlgo = book<TH1I>("h_trkAlgo",
815  "tracking step;iterative tracking step;tracks",
817  0.,
818  double(reco::TrackBase::algoSize));
819 
820  htrkOriAlgo = book<TH1I>("h_trkOriAlgo",
821  "original tracking step;original iterative tracking step;tracks",
823  0.,
824  double(reco::TrackBase::algoSize));
825 
826  for (size_t ibin = 0; ibin < reco::TrackBase::algoSize - 1; ibin++) {
827  htrkAlgo->GetXaxis()->SetBinLabel(ibin + 1, (reco::TrackBase::algoNames[ibin]).c_str());
828  htrkOriAlgo->GetXaxis()->SetBinLabel(ibin + 1, (reco::TrackBase::algoNames[ibin]).c_str());
829  }
830 
831  htrkQuality = book<TH1I>("h_trkQuality", "track quality;track quality;tracks", 6, -1, 5);
832  std::string qualities[7] = {"undef", "loose", "tight", "highPurity", "confirmed", "goodIterative"};
833  for (int nbin = 1; nbin <= htrkQuality->GetNbinsX(); nbin++) {
834  htrkQuality->GetXaxis()->SetBinLabel(nbin, (qualities[nbin - 1]).c_str());
835  }
836 
837  hP = book<TH1D>("h_P", "Momentum;track momentum [GeV];tracks", 100, 0., 100.);
838  hQoverP = book<TH1D>("h_qoverp", "Track q/p; track q/p [GeV^{-1}];tracks", 100, -1., 1.);
839  hQoverPZoom = book<TH1D>("h_qoverpZoom", "Track q/p; track q/p [GeV^{-1}];tracks", 100, -0.1, 0.1);
840  hPt = book<TH1D>("h_Pt", "Transverse Momentum;track p_{T} [GeV];tracks", 100, 0., 100.);
841  hHit = book<TH1D>("h_nHits", "Number of hits;track n. hits;tracks", 50, -0.5, 49.5);
842  hHit2D = book<TH1D>("h_nHit2D", "Number of 2D hits; number of 2D hits;tracks", 20, 0, 20);
843 
844  hHitCountVsZBPix = book<TH1D>("h_HitCountVsZBpix", "Number of BPix hits vs z;hit global z;hits", 60, -30, 30);
845  hHitCountVsZFPix = book<TH1D>("h_HitCountVsZFpix", "Number of FPix hits vs z;hit global z;hits", 100, -100, 100);
846 
847  hHitCountVsXBPix = book<TH1D>("h_HitCountVsXBpix", "Number of BPix hits vs x;hit global x;hits", 20, -20, 20);
848  hHitCountVsXFPix = book<TH1D>("h_HitCountVsXFpix", "Number of FPix hits vs x;hit global x;hits", 20, -20, 20);
849 
850  hHitCountVsYBPix = book<TH1D>("h_HitCountVsYBpix", "Number of BPix hits vs y;hit global y;hits", 20, -20, 20);
851  hHitCountVsYFPix = book<TH1D>("h_HitCountVsYFpix", "Number of FPix hits vs y;hit global y;hits", 20, -20, 20);
852 
854  book<TH1D>("h_HitCountVsThetaBpix", "Number of BPix hits vs #theta;hit global #theta;hits", 20, 0., M_PI);
856  book<TH1D>("h_HitCountVsPhiBpix", "Number of BPix hits vs #phi;hit global #phi;hits", 20, -M_PI, M_PI);
857 
859  book<TH1D>("h_HitCountVsThetaFpix", "Number of FPix hits vs #theta;hit global #theta;hits", 40, 0., M_PI);
861  book<TH1D>("h_HitCountVsPhiFpix", "Number of FPix hits vs #phi;hit global #phi;hits", 20, -M_PI, M_PI);
862 
863  hEta = book<TH1D>("h_Eta", "Track pseudorapidity; track #eta;tracks", 100, -etaMax_, etaMax_);
864  hPhi = book<TH1D>("h_Phi", "Track azimuth; track #phi;tracks", 100, -M_PI, M_PI);
865 
866  hPhiBarrel = book<TH1D>("h_PhiBarrel", "hPhiBarrel (0<|#eta|<0.8);track #Phi;tracks", 100, -M_PI, M_PI);
868  book<TH1D>("h_PhiOverlapPlus", "hPhiOverlapPlus (0.8<#eta<1.4);track #phi;tracks", 100, -M_PI, M_PI);
870  book<TH1D>("h_PhiOverlapMinus", "hPhiOverlapMinus (-1.4<#eta<-0.8);track #phi;tracks", 100, -M_PI, M_PI);
871  hPhiEndcapPlus = book<TH1D>("h_PhiEndcapPlus", "hPhiEndcapPlus (#eta>1.4);track #phi;track", 100, -M_PI, M_PI);
872  hPhiEndcapMinus = book<TH1D>("h_PhiEndcapMinus", "hPhiEndcapMinus (#eta<1.4);track #phi;tracks", 100, -M_PI, M_PI);
873 
874  h_BSx0 = book<TH1F>("h_BSx0", "x-coordinate of reco beamspot;x^{BS}_{0};n_{events}", 100, -0.1, 0.1);
875  h_BSy0 = book<TH1F>("h_BSy0", "y-coordinate of reco beamspot;y^{BS}_{0};n_{events}", 100, -0.1, 0.1);
876  h_BSz0 = book<TH1F>("h_BSz0", "z-coordinate of reco beamspot;z^{BS}_{0};n_{events}", 100, -1., 1.);
877  h_Beamsigmaz = book<TH1F>("h_Beamsigmaz", "z-coordinate beam width;#sigma_{Z}^{beam};n_{events}", 100, 0., 1.);
878  h_BeamWidthX = book<TH1F>("h_BeamWidthX", "x-coordinate beam width;#sigma_{X}^{beam};n_{events}", 100, 0., 0.01);
879  h_BeamWidthY = book<TH1F>("h_BeamWidthY", "y-coordinate beam width;#sigma_{Y}^{beam};n_{events}", 100, 0., 0.01);
880  h_BSdxdz = book<TH1F>("h_BSdxdz", "BeamSpot dxdz;beamspot dx/dz;n_{events}", 100, -0.0003, 0.0003);
881  h_BSdydz = book<TH1F>("h_BSdydz", "BeamSpot dydz;beamspot dy/dz;n_{events}", 100, -0.0003, 0.0003);
882 
883  if (!isCosmics_) {
884  hPhp = book<TH1D>("h_P_hp", "Momentum (high purity);track momentum [GeV];tracks", 100, 0., 100.);
885  hPthp = book<TH1D>("h_Pt_hp", "Transverse Momentum (high purity);track p_{T} [GeV];tracks", 100, 0., 100.);
886  hHithp = book<TH1D>("h_nHit_hp", "Number of hits (high purity);track n. hits;tracks", 30, 0, 30);
887  hEtahp = book<TH1D>("h_Eta_hp", "Track pseudorapidity (high purity); track #eta;tracks", 100, -etaMax_, etaMax_);
888  hPhihp = book<TH1D>("h_Phi_hp", "Track azimuth (high purity); track #phi;tracks", 100, -M_PI, M_PI);
889  hchi2ndofhp = book<TH1D>("h_chi2ndof_hp", "chi2/ndf (high purity);#chi^{2}/ndf;tracks", 100, 0, 5.);
890  hchi2Probhp = book<TH1D>(
891  "h_chi2_Prob_hp", "#chi^{2} probability (high purity);#chi^{2}prob_{Track};Number of Tracks", 100, 0.0, 1.);
892 
893  hvx = book<TH1D>("h_vx", "Track v_{x} ; track v_{x} [cm];tracks", 100, -1.5, 1.5);
894  hvy = book<TH1D>("h_vy", "Track v_{y} ; track v_{y} [cm];tracks", 100, -1.5, 1.5);
895  hvz = book<TH1D>("h_vz", "Track v_{z} ; track v_{z} [cm];tracks", 100, -20., 20.);
896  hd0 = book<TH1D>("h_d0", "Track d_{0} ; track d_{0} [cm];tracks", 100, -1., 1.);
897  hdxy = book<TH1D>("h_dxy", "Track d_{xy}; track d_{xy} [cm]; tracks", 100, -0.5, 0.5);
898  hdz = book<TH1D>("h_dz", "Track d_{z} ; track d_{z} [cm]; tracks", 100, -20, 20);
899 
900  hd0PVvsphi =
901  book<TH2D>("h2_d0PVvsphi", "hd0PVvsphi;track #phi;track d_{0}(PV) [cm]", 160, -M_PI, M_PI, 100, -1., 1.);
902  hd0PVvseta =
903  book<TH2D>("h2_d0PVvseta", "hdPV0vseta;track #eta;track d_{0}(PV) [cm]", 160, -2.5, 2.5, 100, -1., 1.);
904  hd0PVvspt = book<TH2D>("h2_d0PVvspt", "hdPV0vspt;track p_{T};d_{0}(PV) [cm]", 50, 0., 100., 100, -1, 1.);
905 
906  hdxyBS = book<TH1D>("h_dxyBS", "hdxyBS; track d_{xy}(BS) [cm];tracks", 100, -0.1, 0.1);
907  hd0BS = book<TH1D>("h_d0BS", "hd0BS ; track d_{0}(BS) [cm];tracks", 100, -0.1, 0.1);
908  hdzBS = book<TH1D>("h_dzBS", "hdzBS ; track d_{z}(BS) [cm];tracks", 100, -12, 12);
909  hdxyPV = book<TH1D>("h_dxyPV", "hdxyPV; track d_{xy}(PV) [cm];tracks", 100, -0.1, 0.1);
910  hd0PV = book<TH1D>("h_d0PV", "hd0PV ; track d_{0}(PV) [cm];tracks", 100, -0.15, 0.15);
911  hdzPV = book<TH1D>("h_dzPV", "hdzPV ; track d_{z}(PV) [cm];tracks", 100, -0.1, 0.1);
912 
913  hnhTIB = book<TH1D>("h_nHitTIB", "nhTIB;# hits in TIB; tracks", 20, 0., 20.);
914  hnhTID = book<TH1D>("h_nHitTID", "nhTID;# hits in TID; tracks", 20, 0., 20.);
915  hnhTOB = book<TH1D>("h_nHitTOB", "nhTOB;# hits in TOB; tracks", 20, 0., 20.);
916  hnhTEC = book<TH1D>("h_nHitTEC", "nhTEC;# hits in TEC; tracks", 20, 0., 20.);
917 
918  } else {
919  hvx = book<TH1D>("h_vx", "Track v_{x};track v_{x} [cm];tracks", 100, -100., 100.);
920  hvy = book<TH1D>("h_vy", "Track v_{y};track v_{y} [cm];tracks", 100, -100., 100.);
921  hvz = book<TH1D>("h_vz", "Track v_{z};track v_{z} [cm];track", 100, -100., 100.);
922  hd0 = book<TH1D>("h_d0", "Track d_{0};track d_{0} [cm];track", 100, -100., 100.);
923  hdxy = book<TH1D>("h_dxy", "Track d_{xy};track d_{xy} [cm];tracks", 100, -100, 100);
924  hdz = book<TH1D>("h_dz", "Track d_{z};track d_{z} [cm];tracks", 100, -200, 200);
925 
926  hd0vsphi = book<TH2D>(
927  "h2_d0vsphi", "Track d_{0} vs #phi; track #phi;track d_{0} [cm]", 160, -3.20, 3.20, 100, -100., 100.);
928  hd0vseta = book<TH2D>(
929  "h2_d0vseta", "Track d_{0} vs #eta; track #eta;track d_{0} [cm]", 160, -3.20, 3.20, 100, -100., 100.);
930  hd0vspt =
931  book<TH2D>("h2_d0vspt", "Track d_{0} vs p_{T};track p_{T};track d_{0} [cm]", 50, 0., 100., 100, -100, 100);
932 
933  hdxyBS = book<TH1D>("h_dxyBS", "Track d_{xy}(BS);d_{xy}(BS) [cm];tracks", 100, -100., 100.);
934  hd0BS = book<TH1D>("h_d0BS", "Track d_{0}(BS);d_{0}(BS) [cm];tracks", 100, -100., 100.);
935  hdzBS = book<TH1D>("h_dzBS", "Track d_{z}(BS);d_{z}(BS) [cm];tracks", 100, -100., 100.);
936  hdxyPV = book<TH1D>("h_dxyPV", "Track d_{xy}(PV); d_{xy}(PV) [cm];tracks", 100, -100., 100.);
937  hd0PV = book<TH1D>("h_d0PV", "Track d_{0}(PV); d_{0}(PV) [cm];tracks", 100, -100., 100.);
938  hdzPV = book<TH1D>("h_dzPV", "Track d_{z}(PV); d_{z}(PV) [cm];tracks", 100, -100., 100.);
939 
940  hnhTIB = book<TH1D>("h_nHitTIB", "nhTIB;# hits in TIB; tracks", 30, 0., 30.);
941  hnhTID = book<TH1D>("h_nHitTID", "nhTID;# hits in TID; tracks", 30, 0., 30.);
942  hnhTOB = book<TH1D>("h_nHitTOB", "nhTOB;# hits in TOB; tracks", 30, 0., 30.);
943  hnhTEC = book<TH1D>("h_nHitTEC", "nhTEC;# hits in TEC; tracks", 30, 0., 30.);
944  }
945 
946  hnhpxb = book<TH1D>("h_nHitPXB", "nhpxb;# hits in Pixel Barrel; tracks", 10, 0., 10.);
947  hnhpxe = book<TH1D>("h_nHitPXE", "nhpxe;# hits in Pixel Endcap; tracks", 10, 0., 10.);
948 
949  hHitComposition = book<TH1D>("h_hitcomposition", "track hit composition;;# hits", 6, -0.5, 5.5);
950 
952  book<TProfile>("p_NpixHits_vs_Vx", "n. Barrel Pixel hits vs. v_{x};v_{x} (cm);n. BPix hits", 20, -20, 20);
953 
955  book<TProfile>("p_NpixHits_vs_Vy", "n. Barrel Pixel hits vs. v_{y};v_{y} (cm);n. BPix hits", 20, -20, 20);
956 
958  book<TProfile>("p_NpixHits_vs_Vz", "n. Barrel Pixel hits vs. v_{z};v_{z} (cm);n. BPix hits", 20, -100, 100);
959 
960  std::string dets[6] = {"PXB", "PXF", "TIB", "TID", "TOB", "TEC"};
961 
962  for (int i = 1; i <= hHitComposition->GetNbinsX(); i++) {
963  hHitComposition->GetXaxis()->SetBinLabel(i, (dets[i - 1]).c_str());
964  }
965 
966  // vector of track histograms
967 
968  // clang-format off
969 
970  vTrackHistos_.push_back(book<TH1F>("h_tracketa", "Track #eta;#eta_{Track};Number of Tracks", 90, -etaMax_, etaMax_));
971  vTrackHistos_.push_back(book<TH1F>("h_trackphi", "Track #phi;#phi_{Track};Number of Tracks", 90, -M_PI, M_PI));
972  vTrackHistos_.push_back(book<TH1F>("h_trackNumberOfValidHits", "Track # of valid hits;# of valid hits _{Track};Number of Tracks", 40, 0., 40.));
973  vTrackHistos_.push_back(book<TH1F>("h_trackNumberOfLostHits", "Track # of lost hits;# of lost hits _{Track};Number of Tracks", 10, 0., 10.));
974  vTrackHistos_.push_back(book<TH1F>("h_curvature", "Curvature #kappa;#kappa_{Track};Number of Tracks", 100, -.05, .05));
975  vTrackHistos_.push_back(book<TH1F>("h_curvature_pos", "Curvature |#kappa| Positive Tracks;|#kappa_{pos Track}|;Number of Tracks", 100, .0, .05));
976  vTrackHistos_.push_back(book<TH1F>("h_curvature_neg", "Curvature |#kappa| Negative Tracks;|#kappa_{neg Track}|;Number of Tracks", 100, .0, .05));
977  vTrackHistos_.push_back(book<TH1F>("h_diff_curvature","Curvature |#kappa| Tracks Difference;|#kappa_{Track}|;# Pos Tracks - # Neg Tracks",100, .0, .05));
978  vTrackHistos_.push_back(book<TH1F>("h_chi2", "Track #chi^{2};#chi^{2}_{Track};Number of Tracks", 500, -0.01, 500.));
979  vTrackHistos_.push_back(book<TH1F>("h_chi2Prob", "#chi^{2} probability;Track Prob(#chi^{2},ndof);Number of Tracks", 100, 0.0, 1.));
980  vTrackHistos_.push_back(book<TH1F>("h_normchi2", "#chi^{2}/ndof;#chi^{2}/ndof;Number of Tracks", 100, -0.01, 10.));
981  //variable binning for chi2/ndof vs. pT
982  double xBins[19] = {0., 0.15, 0.5, 1., 1.5, 2., 2.5, 3., 3.5, 4., 4.5, 5., 7., 10., 15., 25., 40., 100., 200.};
983  vTrackHistos_.push_back(book<TH1F>("h_pt", "Track p_{T};p_{T}^{track} [GeV];Number of Tracks", 250, 0., 250));
984  vTrackHistos_.push_back(book<TH1F>("h_ptrebin", "Track p_{T};p_{T}^{track} [GeV];Number of Tracks", 18, xBins));
985  vTrackHistos_.push_back(book<TH1F>("h_ptResolution", "#delta_{p_{T}}/p_{T}^{track};#delta_{p_{T}}/p_{T}^{track};Number of Tracks", 100, 0., 0.5));
986 
987  vTrackProfiles_.push_back(book<TProfile>("p_d0_vs_phi", "Transverse Impact Parameter vs. #phi;#phi_{Track};#LT d_{0} #GT [cm]", 100, -M_PI, M_PI));
988  vTrackProfiles_.push_back(book<TProfile>("p_dz_vs_phi", "Longitudinal Impact Parameter vs. #phi;#phi_{Track};#LT d_{z} #GT [cm]", 100, -M_PI, M_PI));
989  vTrackProfiles_.push_back(book<TProfile>("p_d0_vs_eta", "Transverse Impact Parameter vs. #eta;#eta_{Track};#LT d_{0} #GT [cm]", 100, -etaMax_, etaMax_));
990  vTrackProfiles_.push_back(book<TProfile>("p_dz_vs_eta", "Longitudinal Impact Parameter vs. #eta;#eta_{Track};#LT d_{z} #GT [cm]", 100, -etaMax_, etaMax_));
991  vTrackProfiles_.push_back(book<TProfile>("p_chi2_vs_phi", "#chi^{2} vs. #phi;#phi_{Track};#LT #chi^{2} #GT", 100, -M_PI, M_PI));
992  vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_phi", "#chi^{2} probablility vs. #phi;#phi_{Track};#LT #chi^{2} probability#GT", 100, -M_PI, M_PI));
993  vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_d0", "#chi^{2} probablility vs. |d_{0}|;|d_{0}|[cm];#LT #chi^{2} probability#GT", 100, 0, 80));
994  vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_dz", "#chi^{2} probablility vs. dz;d_{z} [cm];#LT #chi^{2} probability#GT", 100, -30, 30));
995  vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_phi", "#chi^{2}/ndof vs. #phi;#phi_{Track};#LT #chi^{2}/ndof #GT", 100, -M_PI, M_PI));
996  vTrackProfiles_.push_back(book<TProfile>("p_chi2_vs_eta", "#chi^{2} vs. #eta;#eta_{Track};#LT #chi^{2} #GT", 100, -etaMax_, etaMax_));
997  vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_pt", "norm #chi^{2} vs. p_{T}_{Track}; p_{T}_{Track};#LT #chi^{2}/ndof #GT", 18, xBins));
998  vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_p", "#chi^{2}/ndof vs. p_{Track};p_{Track};#LT #chi^{2}/ndof #GT", 18, xBins));
999  vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_eta","#chi^{2} probability vs. #eta;#eta_{Track};#LT #chi^{2} probability #GT",100,-etaMax_,etaMax_));
1000  vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_eta", "#chi^{2}/ndof vs. #eta;#eta_{Track};#LT #chi^{2}/ndof #GT", 100, -etaMax_, etaMax_));
1001  vTrackProfiles_.push_back(book<TProfile>("p_kappa_vs_phi", "#kappa vs. #phi;#phi_{Track};#kappa", 100, -M_PI, M_PI));
1002  vTrackProfiles_.push_back(book<TProfile>("p_kappa_vs_eta", "#kappa vs. #eta;#eta_{Track};#kappa", 100, -etaMax_, etaMax_));
1003  vTrackProfiles_.push_back(book<TProfile>("p_ptResolution_vs_phi","#delta_{p_{T}}/p_{T}^{track};#phi^{track};#delta_{p_{T}}/p_{T}^{track}",100,-M_PI,M_PI));
1004  vTrackProfiles_.push_back(book<TProfile>("p_ptResolution_vs_eta","#delta_{p_{T}}/p_{T}^{track};#eta^{track};#delta_{p_{T}}/p_{T}^{track}",100,-etaMax_,etaMax_));
1005 
1006  vTrack2DHistos_.push_back(book<TH2F>("h2_phi_vs_eta", "Track #phi vs. #eta;#eta_{Track};#phi_{Track}",50, -etaMax_, etaMax_, 50, -M_PI, M_PI));
1007  vTrack2DHistos_.push_back(book<TH2F>("h2_d0_vs_phi", "Transverse Impact Parameter vs. #phi;#phi_{Track};d_{0} [cm]", 100, -M_PI, M_PI, 100, -1., 1.));
1008  vTrack2DHistos_.push_back(book<TH2F>("h2_dz_vs_phi", "Longitudinal Impact Parameter vs. #phi;#phi_{Track};d_{z} [cm]",100, -M_PI, M_PI, 100, -100., 100.));
1009  vTrack2DHistos_.push_back(book<TH2F>("h2_d0_vs_eta", "Transverse Impact Parameter vs. #eta;#eta_{Track};d_{0} [cm]", 100, -etaMax_, etaMax_, 100, -1., 1.));
1010  vTrack2DHistos_.push_back(book<TH2F>("h2_dz_vs_eta", "Longitudinal Impact Parameter vs. #eta;#eta_{Track};d_{z} [cm]", 100, -etaMax_, etaMax_, 100, -100., 100.));
1011  vTrack2DHistos_.push_back(book<TH2F>("h2_chi2_vs_phi", "#chi^{2} vs. #phi;#phi_{Track};#chi^{2}", 100, -M_PI, M_PI, 500, 0., 500.));
1012  vTrack2DHistos_.push_back(book<TH2F>("h2_chi2Prob_vs_phi", "#chi^{2} probability vs. #phi;#phi_{Track};#chi^{2} probability", 100, -M_PI, M_PI, 100, 0., 1.));
1013  vTrack2DHistos_.push_back(book<TH2F>("h2_chi2Prob_vs_d0", "#chi^{2} probability vs. |d_{0}|;|d_{0}| [cm];#chi^{2} probability", 100, 0, 80, 100, 0., 1.));
1014  vTrack2DHistos_.push_back(book<TH2F>("h2_normchi2_vs_phi", "#chi^{2}/ndof vs. #phi;#phi_{Track};#chi^{2}/ndof", 100, -M_PI, M_PI, 100, 0., 10.));
1015  vTrack2DHistos_.push_back(book<TH2F>("h2_chi2_vs_eta", "#chi^{2} vs. #eta;#eta_{Track};#chi^{2}", 100, -etaMax_, etaMax_, 500, 0., 500.));
1016  vTrack2DHistos_.push_back(book<TH2F>("h2_chi2Prob_vs_eta", "#chi^{2} probaility vs. #eta;#eta_{Track};#chi^{2} probability", 100, -etaMax_, etaMax_, 100, 0., 1.));
1017  vTrack2DHistos_.push_back(book<TH2F>("h2_normchi2_vs_eta", "#chi^{2}/ndof vs. #eta;#eta_{Track};#chi^{2}/ndof", 100, -etaMax_, etaMax_, 100, 0., 10.));
1018  vTrack2DHistos_.push_back(book<TH2F>("h2_kappa_vs_phi", "#kappa vs. #phi;#phi_{Track};#kappa", 100, -M_PI, M_PI, 100, .0, .05));
1019  vTrack2DHistos_.push_back(book<TH2F>("h2_kappa_vs_eta", "#kappa vs. #eta;#eta_{Track};#kappa", 100, -etaMax_, etaMax_, 100, .0, .05));
1020  vTrack2DHistos_.push_back(book<TH2F>("h2_normchi2_vs_kappa", "#kappa vs. #chi^{2}/ndof;#chi^{2}/ndof;#kappa", 100, 0., 10, 100, -.03, .03));
1021 
1022  // clang-format on
1023 
1024  } //beginJob
1025 
1026  //*************************************************************
1027  void endJob() override
1028  //*************************************************************
1029  {
1030  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1031  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "Events run in total: " << ievt << std::endl;
1032  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "n. tracks: " << itrks << std::endl;
1033  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1034 
1035  int nFiringTriggers = triggerMap_.size();
1036  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "firing triggers: " << nFiringTriggers << std::endl;
1037  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1038 
1039  tksByTrigger_ =
1040  book<TH1D>("tksByTrigger", "tracks by HLT path;;% of # traks", nFiringTriggers, -0.5, nFiringTriggers - 0.5);
1041  evtsByTrigger_ =
1042  book<TH1D>("evtsByTrigger", "events by HLT path;;% of # events", nFiringTriggers, -0.5, nFiringTriggers - 0.5);
1043 
1044  int i = 0;
1045  for (const auto &it : triggerMap_) {
1046  i++;
1047 
1048  double trkpercent = ((it.second).second) * 100. / double(itrks);
1049  double evtpercent = ((it.second).first) * 100. / double(ievt);
1050 
1051  std::cout.precision(4);
1052 
1053  edm::LogPrint("GeneralPurposeTrackAnalyzer")
1054  << "HLT path: " << std::setw(60) << std::left << it.first << " | events firing: " << std::right
1055  << std::setw(8) << (it.second).first << " (" << std::setw(8) << std::fixed << std::setprecision(4)
1056  << evtpercent << "%)"
1057  << " | tracks collected: " << std::setw(10) << (it.second).second << " (" << std::setw(8) << std::fixed
1058  << std::setprecision(4) << trkpercent << "%)";
1059 
1060  tksByTrigger_->SetBinContent(i, trkpercent);
1061  tksByTrigger_->GetXaxis()->SetBinLabel(i, (it.first).c_str());
1062 
1063  evtsByTrigger_->SetBinContent(i, evtpercent);
1064  evtsByTrigger_->GetXaxis()->SetBinLabel(i, (it.first).c_str());
1065  }
1066 
1067  int nRuns = conditionsMap_.size();
1068 
1069  std::vector<int> theRuns_;
1070  for (const auto &it : conditionsMap_) {
1071  theRuns_.push_back(it.first);
1072  }
1073 
1074  sort(theRuns_.begin(), theRuns_.end());
1075  int runRange = theRuns_.back() - theRuns_.front() + 1;
1076 
1077  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1078  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "first run: " << theRuns_.front() << std::endl;
1079  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "last run: " << theRuns_.back() << std::endl;
1080  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "considered runs: " << nRuns << std::endl;
1081  edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1082 
1083  modeByRun_ = book<TH1D>("modeByRun",
1084  "Strip APV mode by run number;;APV mode (-1=deco,+1=peak)",
1085  runRange,
1086  theRuns_.front() - 0.5,
1087  theRuns_.back() + 0.5);
1088 
1089  fieldByRun_ = book<TH1D>("fieldByRun",
1090  "CMS B-field intensity by run number;;B-field intensity [T]",
1091  runRange,
1092  theRuns_.front() - 0.5,
1093  theRuns_.back() + 0.5);
1094 
1095  for (const auto &the_r : theRuns_) {
1096  if (conditionsMap_.find(the_r)->second.first != 0) {
1097  edm::LogPrint("GeneralPurposeTrackAnalyzer")
1098  << "run:" << the_r << " | isPeak: " << std::setw(4) << conditionsMap_.find(the_r)->second.first
1099  << "| B-field: " << conditionsMap_.find(the_r)->second.second << " [T]"
1100  << "| events: " << std::setw(10) << runInfoMap_.find(the_r)->second.first << ", tracks " << std::setw(10)
1101  << runInfoMap_.find(the_r)->second.second << std::endl;
1102  }
1103 
1104  modeByRun_->SetBinContent((the_r - theRuns_.front()) + 1, conditionsMap_.find(the_r)->second.first);
1105  fieldByRun_->SetBinContent((the_r - theRuns_.front()) + 1, conditionsMap_.find(the_r)->second.second);
1106  modeByRun_->GetXaxis()->SetBinLabel((the_r - theRuns_.front()) + 1, std::to_string(the_r).c_str());
1107  fieldByRun_->GetXaxis()->SetBinLabel((the_r - theRuns_.front()) + 1, std::to_string(the_r).c_str());
1108  }
1109 
1110  if (phase_ < SiPixelPI::phase::two) {
1111  if (phase_ == SiPixelPI::phase::zero) {
1112  pmap->save(true, 0, 0, "PixelHitMap.pdf", 600, 800);
1113  pmap->save(true, 0, 0, "PixelHitMap.png", 500, 750);
1114  }
1115 
1116  tmap->save(true, 0, 0, "StripHitMap.pdf");
1117  tmap->save(true, 0, 0, "StripHitMap.png");
1118 
1119  gStyle->SetPalette(kRainBow);
1120  pixelmap->beautifyAllHistograms();
1121 
1122  TCanvas cB("CanvBarrel", "CanvBarrel", 1200, 1000);
1123  pixelmap->drawBarrelMaps("entriesBarrel", cB);
1124  cB.SaveAs("pixelBarrelEntries.png");
1125 
1126  TCanvas cF("CanvForward", "CanvForward", 1600, 1000);
1127  pixelmap->drawForwardMaps("entriesForward", cF);
1128  cF.SaveAs("pixelForwardEntries.png");
1129 
1130  TCanvas cRocs = TCanvas("cRocs", "cRocs", 1200, 1600);
1131  pixelrocsmap_->drawMaps(cRocs, "Pixel on-track clusters occupancy");
1132  cRocs.SaveAs("Phase1PixelROCMaps_fullROCs.png");
1133  }
1134  }
1135 
1136  //*************************************************************
1138  //*************************************************************
1139  {
1140  bool countStereoHitAs2D_ = true;
1141  // we count SiStrip stereo modules as 2D if selected via countStereoHitAs2D_
1142  // (since they provide theta information)
1143  if (!hit.isValid() ||
1144  (hit.dimension() < 2 && !countStereoHitAs2D_ && !dynamic_cast<const SiStripRecHit1D *>(&hit))) {
1145  return false; // real RecHit1D - but SiStripRecHit1D depends on countStereoHitAs2D_
1146  } else {
1147  const DetId detId(hit.geographicalId());
1148  if (detId.det() == DetId::Tracker) {
1149  if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
1150  return true; // pixel is always 2D
1151  } else { // should be SiStrip now
1152  const SiStripDetId stripId(detId);
1153  if (stripId.stereo())
1154  return countStereoHitAs2D_; // stereo modules
1155  else if (dynamic_cast<const SiStripRecHit1D *>(&hit) || dynamic_cast<const SiStripRecHit2D *>(&hit))
1156  return false; // rphi modules hit
1157  //the following two are not used any more since ages...
1158  else if (dynamic_cast<const SiStripMatchedRecHit2D *>(&hit))
1159  return true; // matched is 2D
1160  else if (dynamic_cast<const ProjectedSiStripRecHit2D *>(&hit)) {
1161  const ProjectedSiStripRecHit2D *pH = static_cast<const ProjectedSiStripRecHit2D *>(&hit);
1162  return (countStereoHitAs2D_ && this->isHit2D(pH->originalHit())); // depends on original...
1163  } else {
1164  edm::LogError("UnknownType") << "@SUB=GeneralPurposeTrackAnalyzer::isHit2D"
1165  << "Tracker hit not in pixel, neither SiStripRecHit[12]D nor "
1166  << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
1167  return false;
1168  }
1169  }
1170  } else { // not tracker??
1171  edm::LogWarning("DetectorMismatch") << "@SUB=GeneralPurposeTrackAnalyzer::isHit2D"
1172  << "Hit not in tracker with 'official' dimension >=2.";
1173  return true; // dimension() >= 2 so accept that...
1174  }
1175  }
1176  // never reached...
1177  }
1178 };
1179 
1180 //*************************************************************
1182 //*************************************************************
1183 {
1185  desc.setComment("Generic track analyzer to check ALCARECO sample quantities");
1186  desc.add<edm::InputTag>("TkTag", edm::InputTag("generalTracks"));
1187  desc.add<edm::InputTag>("TriggerResultsTag", edm::InputTag("TriggerResults", "", "HLT"));
1188  desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("offlineBeamSpot"));
1189  desc.add<edm::InputTag>("VerticesTag", edm::InputTag("offlinePrimaryVertices"));
1190  desc.add<bool>("isCosmics", false);
1191  desc.add<bool>("doLatencyAnalysis", true);
1192  descriptions.addWithDefaultLabel(desc);
1193 }
1194 
static const std::string kSharedResource
Definition: TFileService.h:76
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
std::map< int, std::pair< int, float > > conditionsMap_
double z() const
z coordinate
Definition: Vertex.h:133
edm::EDGetTokenT< edm::TriggerResults > hltresultsToken_
edm::EDGetTokenT< reco::TrackCollection > theTrackCollectionToken_
T z() const
Definition: PV3DBase.h:61
Strings const & triggerNames() const
Definition: TriggerNames.cc:48
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
std::unique_ptr< TrackerMap > tmap
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:36
static void fillDescriptions(edm::ConfigurationDescriptions &)
std::unique_ptr< TrackerMap > pmap
std::string to_string(const V &value)
Definition: OMSAccess.h:71
Log< level::Error, false > LogError
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
void beginRun(edm::Run const &run, edm::EventSetup const &setup) override
const int kFPIX
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomTokenBR_
U second(std::pair< T, U > const &p)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T const * product() const
Definition: ESHandle.h:86
edm::ESHandle< MagneticField > magneticField_
std::map< std::string, std::pair< int, int > > triggerMap_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
bool isThere(GeomDetEnumerators::SubDetector subdet) const
GeneralPurposeTrackAnalyzer(const edm::ParameterSet &pset)
uint32_t stereo() const
Definition: SiStripDetId.h:168
T mag() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Transition
Definition: Transition.h:12
void init(const TrackerTopology *, const TrackerGeometry *, const SiPixelFedCablingMap *)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
int roc(const DetId &, const std::pair< int, int > &)
std::map< int, std::pair< int, int > > runInfoMap_
double x() const
x coordinate
Definition: Vertex.h:129
trackCollection
Definition: JetHT_cfg.py:51
Log< level::Warning, true > LogPrint
const TrackerGeomDet * idToDet(DetId) const override
double y() const
y coordinate
Definition: Vertex.h:131
const char * qualities[3]
#define M_PI
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
Log< level::Info, false > LogInfo
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:18
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
Definition: DetId.h:17
int16_t singleReadOutMode() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyTokenBR_
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:50
edm::ESGetToken< SiStripLatency, SiStripLatencyRcd > latencyToken_
static const std::string algoNames[]
Definition: TrackBase.h:147
void endRun(edm::Run const &, edm::EventSetup const &) override
bool isValid() const
Definition: HandleBase.h:70
std::unique_ptr< Phase1PixelROCMaps > pixelrocsmap_
HLT enums.
SiStripRecHit2D originalHit() const
col
Definition: cuy.py:1009
T * book(const Args &...args) const
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
Log< level::Warning, false > LogWarning
std::unique_ptr< Phase1PixelMaps > pixelmap
const int kBPIX
long double T
bool isHit2D(const TrackingRecHit &hit)
~GeneralPurposeTrackAnalyzer() override=default
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
const edm::ESGetToken< SiPixelFedCablingMap, SiPixelFedCablingMapRcd > siPixelFedCablingMapTokenBR_
Definition: event.py:1
Definition: Run.h:45
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
Our base class.
Definition: SiPixelRecHit.h:23
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
int index(const std::vector< OBJECT_TYPE *> &vec, const TString &name)