CMS 3D CMS Logo

ShortenedTrackValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/OfflineValidation
4 // Class: ShortenedTrackValidation
5 //
6 /*
7  *\class ShortenedTrackValidation ShortenedTrackValidation.cc Alignment/OfflineValidation/plugins/ShortenedTrackValidation.cc
8 
9  Description: This module is meant to monitor the track pT resolution using the amputated tracks method, by comparing the performance using different alignments.
10 
11  Implementation: The implemenation takes advantage of the existing implementation in the DQM/TrackingMonitorSource.
12 
13 */
14 //
15 // Original Author: Marco Musich
16 // Created: Fri, 05 Jan 2023 11:41:00 GMT
17 //
18 //
19 
20 // ROOT includes files
21 #include "TMath.h"
22 #include "TFile.h"
23 #include "TH1D.h"
24 #include "TH1I.h"
25 #include "TH2D.h"
26 #include "TProfile.h"
27 #include "TLorentzVector.h"
28 
29 // standard includes
30 #include <fmt/printf.h>
31 
32 // user includes
47 #include "FWCore/Utilities/interface/transform.h" // for edm::vector_transform
48 
49 #define CREATE_HIST_1D(varname, nbins, first, last, fs) fs.make<TH1D>(#varname, #varname, nbins, first, last)
50 
51 #define CREATE_HIST_2D(varname, nbins, first, last, fs) \
52  fs.make<TH2D>(#varname, #varname, nbins, first, last, nbins, first, last)
53 
56 
57 class ShortenedTrackValidation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
58  class trackingMon {
59  public:
61  ~trackingMon() = default;
62 
63  void book(const TFileDirectory &fs) {
64  h_chi2ndof = CREATE_HIST_1D(h_chi2ndof, 100, 0.0, 10.0, fs);
67  h_trkOriAlgo =
69  h_P = CREATE_HIST_1D(h_P, 100, 0.0, 200.0, fs);
70  h_Pt = CREATE_HIST_1D(h_Pt, 100, 0.0, 100.0, fs);
71  h_nHit = CREATE_HIST_1D(h_nHit, 50, -0.5, 49.5, fs);
72  h_nHit2D = CREATE_HIST_1D(h_nHit2D, 20, -0.5, 19.5, fs);
73  h_Charge = CREATE_HIST_1D(h_Charge, 3, -1.5, 1.5, fs);
74  h_QoverP = CREATE_HIST_1D(h_QoverP, 100, -1.0, 1.0, fs);
75  h_QoverPZoom = CREATE_HIST_1D(h_QoverPZoom, 100, -0.1, 0.1, fs);
76  h_Eta = CREATE_HIST_1D(h_Eta, 100, -3., 3., fs);
77  h_Phi = CREATE_HIST_1D(h_Phi, 100, -M_PI, M_PI, fs);
78  h_vx = CREATE_HIST_1D(h_vx, 100, -0.5, 0.5, fs);
79  h_vy = CREATE_HIST_1D(h_vy, 100, -0.5, 0.5, fs);
80  h_vz = CREATE_HIST_1D(h_vz, 100, -20.0, 20.0, fs);
81  h_d0 = CREATE_HIST_1D(h_d0, 100, -0.5, 0.5, fs);
82  h_dz = CREATE_HIST_1D(h_dz, 100, -20.0, 20.0, fs);
83  h_dxy = CREATE_HIST_1D(h_dxy, 100, -0.5, 0.5, fs);
84  h_nhpxb = CREATE_HIST_1D(h_nhpxb, 10, -0.5, 9.5, fs);
85  h_nhpxe = CREATE_HIST_1D(h_nhpxe, 10, -0.5, 9.5, fs);
86  h_nhTIB = CREATE_HIST_1D(h_nhTIB, 20, -0.5, 19.5, fs);
87  h_nhTID = CREATE_HIST_1D(h_nhTID, 20, -0.5, 19.5, fs);
88  h_nhTOB = CREATE_HIST_1D(h_nhTOB, 20, -0.5, 19.5, fs);
89  h_nhTEC = CREATE_HIST_1D(h_nhTEC, 20, -0.5, 19.5, fs);
90  h_dxyBS = CREATE_HIST_1D(h_dxyBS, 100, -0.05, 0.05, fs);
91  h_d0BS = CREATE_HIST_1D(h_d0BS, 100, -0.05, 0.05, fs);
92  h_dzBS = CREATE_HIST_1D(h_dzBS, 100, -20.0, 20., fs);
93  h_dxyPV = CREATE_HIST_1D(h_dxyPV, 100, -0.05, 0.05, fs);
94  h_d0PV = CREATE_HIST_1D(h_d0PV, 100, -0.05, 0.05, fs);
95  h_dzPV = CREATE_HIST_1D(h_dzPV, 100, -0.05, 0.05, fs);
96 
97  edm::LogInfo("trackingMonitoring") << "done booking";
98  }
99 
100  //____________________________________________________________
102  int myquality = -99;
104  myquality = -1;
105  if (track.quality(reco::TrackBase::loose))
106  myquality = 0;
107  if (track.quality(reco::TrackBase::tight))
108  myquality = 1;
109  if (track.quality(reco::TrackBase::highPurity))
110  myquality = 2;
112  myquality = 3;
113 
114  return myquality;
115  }
116 
117  //____________________________________________________________
118  static bool isHit2D(const TrackingRecHit &hit) {
119  if (hit.dimension() < 2) {
120  return false; // some (muon...) stuff really has RecHit1D
121  } else {
122  const DetId detId(hit.geographicalId());
123  if (detId.det() == DetId::Tracker) {
124  if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
125  return true; // pixel is always 2D
126  } else { // should be SiStrip now
127  if (dynamic_cast<const SiStripRecHit2D *>(&hit))
128  return false; // normal hit
129  else if (dynamic_cast<const SiStripMatchedRecHit2D *>(&hit))
130  return true; // matched is 2D
131  else if (dynamic_cast<const ProjectedSiStripRecHit2D *>(&hit))
132  return false; // crazy hit...
133  else {
134  edm::LogError("UnknownType") << "@SUB=CalibrationTrackSelector::isHit2D"
135  << "Tracker hit not in pixel and neither SiStripRecHit2D nor "
136  << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
137  return false;
138  }
139  }
140  } else { // not tracker??
141  edm::LogWarning("DetectorMismatch") << "@SUB=CalibrationTrackSelector::isHit2D"
142  << "Hit not in tracker with 'official' dimension >=2.";
143  return true; // dimension() >= 2 so accept that...
144  }
145  }
146  // never reached...
147  }
148 
149  //____________________________________________________________
150  unsigned int count2DHits(const reco::Track &track) {
151  unsigned int nHit2D = 0;
152  for (auto iHit = track.recHitsBegin(); iHit != track.recHitsEnd(); ++iHit) {
153  if (isHit2D(**iHit)) {
154  ++nHit2D;
155  }
156  }
157  return nHit2D;
158  }
159 
160  //____________________________________________________________
161  void fill(const reco::Track &track, const reco::BeamSpot &beamSpot, const reco::Vertex &pvtx) {
162  h_chi2ndof->Fill(track.normalizedChi2());
163  h_trkQuality->Fill(trackQual(track));
164  h_trkAlgo->Fill(static_cast<float>(track.algo()));
165  h_trkOriAlgo->Fill(static_cast<float>(track.originalAlgo()));
166  h_P->Fill(track.p());
167  h_Pt->Fill(track.pt());
168  h_nHit->Fill(track.numberOfValidHits());
169  h_nHit2D->Fill(count2DHits(track));
170  h_Charge->Fill(track.charge());
171  h_QoverP->Fill(track.qoverp());
172  h_QoverPZoom->Fill(track.qoverp());
173  h_Eta->Fill(track.eta());
174  h_Phi->Fill(track.phi());
175  h_vx->Fill(track.vx());
176  h_vy->Fill(track.vy());
177  h_vz->Fill(track.vz());
178  h_d0->Fill(track.d0());
179  h_dz->Fill(track.dz());
180  h_dxy->Fill(track.dxy());
181  h_nhpxb->Fill(track.hitPattern().numberOfValidPixelBarrelHits());
182  h_nhpxe->Fill(track.hitPattern().numberOfValidPixelEndcapHits());
183  h_nhTIB->Fill(track.hitPattern().numberOfValidStripTIBHits());
184  h_nhTID->Fill(track.hitPattern().numberOfValidStripTIDHits());
185  h_nhTOB->Fill(track.hitPattern().numberOfValidStripTOBHits());
186  h_nhTEC->Fill(track.hitPattern().numberOfValidStripTECHits());
187 
188  math::XYZPoint BS(beamSpot.x0(), beamSpot.y0(), beamSpot.z0());
189  h_dxyBS->Fill(track.dxy(BS));
190  h_d0BS->Fill(-track.dxy(BS));
191  h_dzBS->Fill(track.dz(BS));
192 
193  math::XYZPoint PV(pvtx.x(), pvtx.y(), pvtx.z());
194  h_dxyPV->Fill(track.dxy(PV));
195  h_d0PV->Fill(-track.dxy(PV));
196  h_dzPV->Fill(track.dz(PV));
197  }
198 
199  private:
200  TH1D *h_chi2ndof;
202  TH1D *h_trkAlgo;
204  TH1D *h_P;
205  TH1D *h_Pt;
206  TH1D *h_nHit;
207  TH1D *h_nHit2D;
208  TH1D *h_Charge;
209  TH1D *h_QoverP;
211  TH1D *h_Eta;
212  TH1D *h_Phi;
213  TH1D *h_vx;
214  TH1D *h_vy;
215  TH1D *h_vz;
216  TH1D *h_d0;
217  TH1D *h_dz;
218  TH1D *h_dxy;
219  TH1D *h_nhpxb;
220  TH1D *h_nhpxe;
221  TH1D *h_nhTIB;
222  TH1D *h_nhTID;
223  TH1D *h_nhTOB;
224  TH1D *h_nhTEC;
225  TH1D *h_dxyBS;
226  TH1D *h_d0BS;
227  TH1D *h_dzBS;
228  TH1D *h_dxyPV;
229  TH1D *h_d0PV;
230  TH1D *h_dzPV;
231  };
232 
234  public:
236  ~trackComparator() = default;
237 
238  //__________________________________________________
239  void book(const TFileDirectory &fs) {
240  h2_chi2ndof = CREATE_HIST_2D(h2_chi2ndof, 100, 0.0, 10.0, fs);
242  h2_trkOriAlgo =
244  h2_P = CREATE_HIST_2D(h2_P, 100, 0.0, 200.0, fs);
245  h2_Pt = CREATE_HIST_2D(h2_Pt, 100, 0.0, 100.0, fs);
246  h2_nHit = CREATE_HIST_2D(h2_nHit, 50, -0.5, 49.5, fs);
247  h2_Charge = CREATE_HIST_2D(h2_Charge, 3, -1.5, 1.5, fs);
248  h2_QoverPZoom = CREATE_HIST_2D(h2_QoverPZoom, 100, -0.1, 0.1, fs);
249  h2_Eta = CREATE_HIST_2D(h2_Eta, 100, -3., 3., fs);
250  h2_Phi = CREATE_HIST_2D(h2_Phi, 100, -M_PI, M_PI, fs);
251  h2_vx = CREATE_HIST_2D(h2_vx, 100, -0.5, 0.5, fs);
252  h2_vy = CREATE_HIST_2D(h2_vy, 100, -0.5, 0.5, fs);
253  h2_vz = CREATE_HIST_2D(h2_vz, 100, -20.0, 20.0, fs);
254  h2_d0 = CREATE_HIST_2D(h2_d0, 100, -0.5, 0.5, fs);
255  h2_dz = CREATE_HIST_2D(h2_dz, 100, -20.0, 20.0, fs);
256  h2_nhpxb = CREATE_HIST_2D(h2_nhpxb, 10, -0.5, 9.5, fs);
257  h2_nhpxe = CREATE_HIST_2D(h2_nhpxe, 10, -0.5, 9.5, fs);
258  h2_nhTIB = CREATE_HIST_2D(h2_nhTIB, 20, -0.5, 19.5, fs);
259  h2_nhTID = CREATE_HIST_2D(h2_nhTID, 20, -0.5, 19.5, fs);
260  h2_nhTOB = CREATE_HIST_2D(h2_nhTOB, 20, -0.5, 19.5, fs);
261  h2_nhTEC = CREATE_HIST_2D(h2_nhTEC, 20, -0.5, 19.5, fs);
262  }
263 
264  //__________________________________________________
265  void fill(const reco::Track &tk1, const reco::Track &tk2) {
266  h2_chi2ndof->Fill(tk1.normalizedChi2(), tk2.normalizedChi2());
267  h2_trkAlgo->Fill(static_cast<float>(tk1.algo()), static_cast<float>(tk2.algo()));
268  h2_trkOriAlgo->Fill(static_cast<float>(tk1.originalAlgo()), static_cast<float>(tk2.originalAlgo()));
269  h2_P->Fill(tk1.p(), tk2.p());
270  h2_Pt->Fill(tk1.pt(), tk2.p());
271  h2_nHit->Fill(tk1.numberOfValidHits(), tk2.numberOfValidHits());
272  h2_Charge->Fill(tk1.charge(), tk2.charge());
273  h2_QoverPZoom->Fill(tk1.qoverp(), tk2.qoverp());
274  h2_Eta->Fill(tk1.eta(), tk2.eta());
275  h2_Phi->Fill(tk1.phi(), tk2.phi());
276  h2_vx->Fill(tk1.vx(), tk2.vx());
277  h2_vy->Fill(tk1.vy(), tk2.vy());
278  h2_vz->Fill(tk1.vz(), tk2.vz());
279  h2_d0->Fill(tk1.d0(), tk2.d0());
280  h2_dz->Fill(tk2.dz(), tk2.dz());
287  }
288 
289  private:
290  TH2D *h2_chi2ndof;
291  TH2D *h2_trkAlgo;
293  TH2D *h2_P;
294  TH2D *h2_Pt;
295  TH2D *h2_nHit;
296  TH2D *h2_Charge;
298  TH2D *h2_Eta;
299  TH2D *h2_Phi;
300  TH2D *h2_vx;
301  TH2D *h2_vy;
302  TH2D *h2_vz;
303  TH2D *h2_d0;
304  TH2D *h2_dz;
305  TH2D *h2_nhpxb;
306  TH2D *h2_nhpxe;
307  TH2D *h2_nhTIB;
308  TH2D *h2_nhTID;
309  TH2D *h2_nhTOB;
310  TH2D *h2_nhTEC;
311  };
312 
313 public:
315  ~ShortenedTrackValidation() override = default;
316  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
317 
318 private:
319  template <typename T, typename... Args>
320  T *book(const TFileDirectory &dir, const Args &...args) const;
321  void beginJob() override;
322  void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
323 
324  // ----------member data ---------------------------
327  const std::vector<std::string> hitsRemain_;
328  const double minTracksEta_;
329  const double maxTracksEta_;
330  const double minTracksPt_;
331  const double maxTracksPt_;
332 
333  const double maxDr_;
335  const std::vector<edm::InputTag> tracksRerecoTag_;
339  const std::vector<edm::EDGetTokenT<std::vector<reco::Track>>> tracksRerecoToken_;
342 
343  // monitoring histograms
344  std::vector<TH1F *> histsPtRatioAll_;
345  std::vector<TH1F *> histsPtDiffAll_;
346  std::vector<TH1F *> histsEtaDiffAll_;
347  std::vector<TH1F *> histsPhiDiffAll_;
348  std::vector<TH2F *> histsPtRatioVsDeltaRAll_;
349  std::vector<TH1F *> histsDeltaPtOverPtAll_;
350  std::vector<TH1F *> histsPtAll_;
351  std::vector<TH1F *> histsNhitsAll_;
352  std::vector<TH1F *> histsDeltaRAll_;
353 
355  std::vector<trackComparator *> comparators_;
356  static constexpr double muMass = 0.105658;
357 };
358 
359 // -----------------------------
360 // constructors and destructor
361 // -----------------------------
363  : folderName_(ps.getUntrackedParameter<std::string>("folderName", "TrackRefitting")),
364  hitsRemain_(ps.getUntrackedParameter<std::vector<std::string>>("hitsRemainInput")),
365  minTracksEta_(ps.getUntrackedParameter<double>("minTracksEtaInput", 0.0)),
366  maxTracksEta_(ps.getUntrackedParameter<double>("maxTracksEtaInput", 2.2)),
367  minTracksPt_(ps.getUntrackedParameter<double>("minTracksPtInput", 15.0)),
368  maxTracksPt_(ps.getUntrackedParameter<double>("maxTracksPtInput", 99999.9)),
369  maxDr_(ps.getUntrackedParameter<double>("maxDrInput", 0.01)),
370  tracksTag_(ps.getUntrackedParameter<edm::InputTag>("tracksInputTag", edm::InputTag("generalTracks", "", "DQM"))),
371  tracksRerecoTag_(ps.getUntrackedParameter<std::vector<edm::InputTag>>("tracksRerecoInputTag")),
372  BeamSpotTag_(ps.getUntrackedParameter<edm::InputTag>("BeamSpotTag", edm::InputTag("offlineBeamSpot"))),
373  VerticesTag_(ps.getUntrackedParameter<edm::InputTag>("VerticesTag", edm::InputTag("offlinePrimaryVertices"))),
374  tracksToken_(consumes<std::vector<reco::Track>>(tracksTag_)),
375  tracksRerecoToken_(edm::vector_transform(
376  tracksRerecoTag_, [this](edm::InputTag const &tag) { return consumes<std::vector<reco::Track>>(tag); })),
377  beamspotToken_(consumes<reco::BeamSpot>(BeamSpotTag_)),
378  vertexToken_(consumes<reco::VertexCollection>(VerticesTag_)) {
379  usesResource(TFileService::kSharedResource);
380  histsPtRatioAll_.clear();
381  histsPtDiffAll_.clear();
382  histsEtaDiffAll_.clear();
383  histsPhiDiffAll_.clear();
384  histsPtRatioVsDeltaRAll_.clear();
385  histsDeltaPtOverPtAll_.clear();
386  histsPtAll_.clear();
387  histsNhitsAll_.clear();
388  histsDeltaRAll_.clear();
389  comparators_.clear();
390 
391  comparators_.reserve(hitsRemain_.size());
392  for (unsigned int i = 0; i < hitsRemain_.size(); ++i) {
393  comparators_.push_back(new trackComparator());
394  }
395 }
396 
397 //__________________________________________________________________________________
398 template <typename T, typename... Args>
399 T *ShortenedTrackValidation::book(const TFileDirectory &dir, const Args &...args) const {
400  T *t = dir.make<T>(args...);
401  return t;
402 }
403 
404 //__________________________________________________________________________________
406  std::string currentFolder = folderName_ + "/Resolutions";
407  TFileDirectory ShortTrackResolution = fs_->mkdir(currentFolder);
408  currentFolder = folderName_ + "/Tracks";
409  TFileDirectory TrackQuals = fs_->mkdir(currentFolder);
410 
411  for (unsigned int i = 0; i < hitsRemain_.size(); ++i) {
412  histsPtRatioAll_.push_back(
413  book<TH1F>(ShortTrackResolution,
414  fmt::sprintf("trackPtRatio_%s", hitsRemain_[i]).c_str(),
415  fmt::sprintf("Short Track p_{T} / Full Track p_{T} - %s layers;p_{T}^{short}/p_{T}^{full};n. tracks",
416  hitsRemain_[i])
417  .c_str(),
418  100,
419  0.5,
420  1.5));
421 
422  histsPtDiffAll_.push_back(book<TH1F>(
423  ShortTrackResolution,
424  fmt::sprintf("trackPtDiff_%s", hitsRemain_[i]).c_str(),
425  fmt::sprintf("Short Track p_{T} - Full Track p_{T} - %s layers;p_{T}^{short} - p_{T}^{full} [GeV];n. tracks",
426  hitsRemain_[i])
427  .c_str(),
428  100,
429  -10.,
430  10.));
431 
432  histsEtaDiffAll_.push_back(
433  book<TH1F>(ShortTrackResolution,
434  fmt::sprintf("trackEtaDiff_%s", hitsRemain_[i]).c_str(),
435  fmt::sprintf("Short Track #eta - Full Track #eta - %s layers;#eta^{short} - #eta^{full};n. tracks",
436  hitsRemain_[i])
437  .c_str(),
438  100,
439  -0.001,
440  0.001));
441 
442  histsPhiDiffAll_.push_back(
443  book<TH1F>(ShortTrackResolution,
444  fmt::sprintf("trackPhiDiff_%s", hitsRemain_[i]).c_str(),
445  fmt::sprintf("Short Track #phi - Full Track #phi - %s layers;#phi^{short} - #phi^{full};n. tracks",
446  hitsRemain_[i])
447  .c_str(),
448  100,
449  -0.001,
450  0.001));
451 
452  histsPtRatioVsDeltaRAll_.push_back(
453  book<TH2F>(ShortTrackResolution,
454  fmt::sprintf("trackPtRatioVsDeltaR_%s", hitsRemain_[i]).c_str(),
455  fmt::sprintf("Short Track p_{T} / Full Track p_{T} - %s layers vs "
456  "#DeltaR;#DeltaR(short,full);p_{T}^{short}/p_{T}^{full} [GeV];n. tracks",
457  hitsRemain_[i])
458  .c_str(),
459  100,
460  0.,
461  0.01,
462  101,
463  -0.05,
464  2.05));
465 
466  histsDeltaPtOverPtAll_.push_back(
467  book<TH1F>(ShortTrackResolution,
468  fmt::sprintf("trackDeltaPtOverPt_%s", hitsRemain_[i]).c_str(),
469  fmt::sprintf("Short Track p_{T} - Full Track p_{T} / Full Track p_{T} - %s layers;p_{T}^{short} - "
470  "p_{T}^{full} / p^{full}_{T};n. tracks",
471  hitsRemain_[i])
472  .c_str(),
473  101,
474  -5.,
475  5.));
476 
477  histsPtAll_.push_back(
478  book<TH1F>(TrackQuals,
479  fmt::sprintf("trackPt_%s", hitsRemain_[i]).c_str(),
480  fmt::sprintf("Short Track p_{T} - %s layers;p_{T}^{short} [GeV];n. tracks", hitsRemain_[i]).c_str(),
481  100,
482  0.,
483  100.));
484 
485  histsNhitsAll_.push_back(
486  book<TH1F>(TrackQuals,
487  fmt::sprintf("trackNhits_%s", hitsRemain_[i]).c_str(),
488  fmt::sprintf("Short Track n. hits - %s layers; n. hits per track;n. tracks", hitsRemain_[i]).c_str(),
489  20,
490  -0.5,
491  19.5));
492 
493  histsDeltaRAll_.push_back(book<TH1F>(
494  TrackQuals,
495  fmt::sprintf("trackDeltaR_%s", hitsRemain_[i]).c_str(),
496  fmt::sprintf("Short Track / Long Track #DeltaR %s layers;#DeltaR(short,long);n. tracks", hitsRemain_[i]).c_str(),
497  100,
498  0.,
499  0.005));
500 
501  currentFolder = fmt::sprintf("%s/Compare_%sHit", folderName_, hitsRemain_[i]);
502  comparators_[i]->book(fs_->mkdir(currentFolder));
503  }
504 
505  currentFolder = folderName_ + "/OriginalTrack";
506  TFileDirectory original = fs_->mkdir(currentFolder);
508 }
509 
510 //__________________________________________________________________________________
512  const auto &tracks = iEvent.getHandle(tracksToken_);
513 
514  if (!tracks.isValid()) {
515  edm::LogError("ShortenedTrackValidation") << "Missing input track collection " << tracksTag_.encode() << std::endl;
516  return;
517  }
518 
520  edm::Handle<reco::BeamSpot> beamSpotHandle = iEvent.getHandle(beamspotToken_);
521  if (beamSpotHandle.isValid()) {
522  beamSpot = *beamSpotHandle;
523  } else {
525  }
526 
527  reco::Vertex pvtx;
528  edm::Handle<reco::VertexCollection> vertexHandle = iEvent.getHandle(vertexToken_);
529  if (vertexHandle.isValid()) {
530  pvtx = (*vertexHandle).at(0);
531  } else {
532  pvtx = reco::Vertex();
533  }
534 
535  // the original long track
536  for (const auto &track : *tracks) {
537  const reco::HitPattern &hp = track.hitPattern();
538  if (int(int(hp.numberOfValidHits()) - int(hp.numberOfAllHits(reco::HitPattern::TRACK_HITS))) != 0) {
539  break;
540  }
541 
542  // fill the original track properties monitoring
544 
545  TLorentzVector tvec;
546  tvec.SetPtEtaPhiM(track.pt(), track.eta(), track.phi(), muMass);
547 
548  int i = 0; // token index
549  // loop on the re-recoed shortened track collections
550  for (const auto &token : tracksRerecoToken_) {
551  const auto &tracks_rereco = iEvent.getHandle(token);
552 
553  for (const auto &track_rereco : *tracks_rereco) {
554  TLorentzVector trerecovec;
555  trerecovec.SetPtEtaPhiM(track_rereco.pt(), track_rereco.eta(), track_rereco.phi(), 0.0);
556  double deltaR = tvec.DeltaR(trerecovec);
557 
558  if (deltaR < maxDr_) {
559  if (track_rereco.pt() >= minTracksPt_ && track_rereco.pt() <= maxTracksPt_ &&
560  std::abs(track_rereco.eta()) >= minTracksEta_ && std::abs(track_rereco.eta()) <= maxTracksEta_) {
561  // fill the 2D comparisons per track
562  comparators_[i]->fill(track, track_rereco);
563 
564  histsPtRatioAll_[i]->Fill(1.0 * track_rereco.pt() / track.pt());
565  histsPtDiffAll_[i]->Fill(track_rereco.pt() - track.pt());
566  histsDeltaPtOverPtAll_[i]->Fill((track_rereco.pt() - track.pt()) / track.pt());
567  histsEtaDiffAll_[i]->Fill(track_rereco.eta() - track.eta());
568  histsPhiDiffAll_[i]->Fill(track_rereco.phi() - track.phi());
569  histsPtRatioVsDeltaRAll_[i]->Fill(deltaR, track_rereco.pt() / track.pt());
570  histsPtAll_[i]->Fill(track_rereco.pt());
571  histsNhitsAll_[i]->Fill(track_rereco.numberOfValidHits());
572  histsDeltaRAll_[i]->Fill(deltaR);
573  }
574  }
575  }
576  ++i;
577  }
578  }
579 }
580 
581 //__________________________________________________________________________________
584  desc.addUntracked<std::string>("folderName", "TrackRefitting");
585  desc.addUntracked<std::vector<std::string>>("hitsRemainInput", {});
586  desc.addUntracked<double>("minTracksEtaInput", 0.0);
587  desc.addUntracked<double>("maxTracksEtaInput", 2.2);
588  desc.addUntracked<double>("minTracksPtInput", 15.0);
589  desc.addUntracked<double>("maxTracksPtInput", 99999.9);
590  desc.addUntracked<double>("maxDrInput", 0.01);
591  desc.addUntracked<edm::InputTag>("tracksInputTag", edm::InputTag("generalTracks", "", "DQM"));
592  desc.addUntracked<std::vector<edm::InputTag>>("tracksRerecoInputTag", {});
593  desc.addUntracked<edm::InputTag>("BeamSpotTag", edm::InputTag("offlineBeamSpot"));
594  desc.addUntracked<edm::InputTag>("VerticesTag", edm::InputTag("offlinePrimaryVertices"));
595  descriptions.addWithDefaultLabel(desc);
596 }
597 
598 // Define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
double qoverp() const
q / p
Definition: TrackBase.h:599
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void fill(const reco::Track &tk1, const reco::Track &tk2)
double z() const
z coordinate
Definition: Vertex.h:134
ShortenedTrackValidation(const edm::ParameterSet &)
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:655
const std::vector< std::string > hitsRemain_
edm::Service< TFileService > fs_
std::string encode() const
Definition: InputTag.cc:159
double p() const
momentum vector magnitude
Definition: TrackBase.h:631
T * book(const TFileDirectory &dir, const Args &...args) const
std::vector< TH1F * > histsNhitsAll_
const int kBPIX
int numberOfValidStripTIBHits() const
Definition: HitPattern.h:847
int numberOfValidStripTECHits() const
Definition: HitPattern.h:859
Log< level::Error, false > LogError
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< TH1F * > histsDeltaPtOverPtAll_
std::vector< TH1F * > histsPhiDiffAll_
static constexpr double muMass
double pt() const
track transverse momentum
Definition: TrackBase.h:637
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:548
std::vector< TH1F * > histsPtRatioAll_
int charge() const
track electric charge
Definition: TrackBase.h:596
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
int iEvent
Definition: GenABIO.cc:224
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:661
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.h:839
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int numberOfValidStripTIDHits() const
Definition: HitPattern.h:851
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< TH2F * > histsPtRatioVsDeltaRAll_
double x() const
x coordinate
Definition: Vertex.h:130
const edm::EDGetTokenT< std::vector< reco::Track > > tracksToken_
double y() const
y coordinate
Definition: Vertex.h:132
#define M_PI
unsigned int count2DHits(const reco::Track &track)
std::vector< TH1F * > histsDeltaRAll_
Log< level::Info, false > LogInfo
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
Definition: DetId.h:17
TrackAlgorithm algo() const
Definition: TrackBase.h:547
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:611
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const std::vector< edm::InputTag > tracksRerecoTag_
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
~ShortenedTrackValidation() override=default
bool isValid() const
Definition: HandleBase.h:70
#define CREATE_HIST_1D(varname, nbins, first, last, fs)
std::vector< trackComparator * > comparators_
std::vector< TH1F * > histsEtaDiffAll_
fixed size matrix
HLT enums.
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:658
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.h:835
const edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
#define CREATE_HIST_2D(varname, nbins, first, last, fs)
Log< level::Warning, false > LogWarning
const std::vector< edm::EDGetTokenT< std::vector< reco::Track > > > tracksRerecoToken_
long double T
Definition: PV.py:1
const int kFPIX
const edm::EDGetTokenT< reco::VertexCollection > vertexToken_
std::vector< TH1F * > histsPtDiffAll_
int numberOfValidStripTOBHits() const
Definition: HitPattern.h:855
static bool isHit2D(const TrackingRecHit &hit)
void fill(const reco::Track &track, const reco::BeamSpot &beamSpot, const reco::Vertex &pvtx)