CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackToTrackComparisonHists.cc
Go to the documentation of this file.
3 
5 
9 
12 
13 //
14 // constructors and destructor
15 //
17  : monitoredTrackInputTag_(iConfig.getParameter<edm::InputTag>("monitoredTrack")),
18  referenceTrackInputTag_(iConfig.getParameter<edm::InputTag>("referenceTrack")),
19  topDirName_(iConfig.getParameter<std::string>("topDirName")),
20  dRmin_(iConfig.getParameter<double>("dRmin")),
21  pTCutForPlateau_(iConfig.getParameter<double>("pTCutForPlateau")),
22  dxyCutForPlateau_(iConfig.getParameter<double>("dxyCutForPlateau")),
23  dzWRTPvCut_(iConfig.getParameter<double>("dzWRTPvCut")),
24  requireValidHLTPaths_(iConfig.getParameter<bool>("requireValidHLTPaths")),
25  genTriggerEventFlag_(new GenericTriggerEventFlag(
26  iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this))
27 
28 {
29  initialize_parameter(iConfig);
30 
31  //now do what ever initialization is needed
32  monitoredTrackToken_ = consumes<reco::TrackCollection>(monitoredTrackInputTag_);
33  referenceTrackToken_ = consumes<reco::TrackCollection>(referenceTrackInputTag_);
34  monitoredBSToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("monitoredBeamSpot"));
35  referenceBSToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("referenceBeamSpot"));
36  monitoredPVToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("monitoredPrimaryVertices"));
37  referencePVToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("referencePrimaryVertices"));
38 
41 
44 
45  matchTracksMEs_.label = "matches";
46 }
47 
50  genTriggerEventFlag_.reset();
51 }
52 
54 
56  LogDebug("TrackToTrackComparisonHists")
57  << " requireValidHLTPaths_ " << requireValidHLTPaths_ << " hltPathsAreValid_ " << hltPathsAreValid_ << "\n";
58  // if valid HLT paths are required,
59  // analyze event only if paths are valid
60  if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
61  return;
62  }
63 
64  LogDebug("TrackToTrackComparisonHists") << " genTriggerEventFlag_->on() " << genTriggerEventFlag_->on()
65  << " accept: " << genTriggerEventFlag_->accept(iEvent, iSetup) << "\n";
66  // Filter out events if Trigger Filtering is requested
67  if (genTriggerEventFlag_->on() && !genTriggerEventFlag_->accept(iEvent, iSetup)) {
68  return;
69  }
70 
71  //
72  // Get Reference Track Info
73  //
74  edm::Handle<reco::TrackCollection> referenceTracksHandle;
75  iEvent.getByToken(referenceTrackToken_, referenceTracksHandle);
76  if (!referenceTracksHandle.isValid()) {
77  edm::LogError("TrackToTrackComparisonHists") << "referenceTracksHandle not found, skipping event";
78  return;
79  }
80  reco::TrackCollection referenceTracks = *referenceTracksHandle;
81 
82  edm::Handle<reco::BeamSpot> referenceBSHandle;
83  iEvent.getByToken(referenceBSToken_, referenceBSHandle);
84  if (!referenceBSHandle.isValid()) {
85  edm::LogError("TrackToTrackComparisonHists") << "referenceBSHandle not found, skipping event";
86  return;
87  }
88  reco::BeamSpot referenceBS = *referenceBSHandle;
89 
90  edm::Handle<reco::VertexCollection> referencePVHandle;
91  iEvent.getByToken(referencePVToken_, referencePVHandle);
92  if (!referencePVHandle.isValid()) {
93  edm::LogError("TrackToTrackComparisonHists") << "referencePVHandle not found, skipping event";
94  return;
95  }
96  if (referencePVHandle->empty()) {
97  edm::LogInfo("TrackToTrackComparisonHists") << "referencePVHandle->size is 0 ";
98  return;
99  }
100  reco::Vertex referencePV = referencePVHandle->at(0);
101 
102  //
103  // Get Monitored Track Info
104  //
105  edm::Handle<reco::TrackCollection> monitoredTracksHandle;
106  iEvent.getByToken(monitoredTrackToken_, monitoredTracksHandle);
107  if (!monitoredTracksHandle.isValid()) {
108  edm::LogError("TrackToTrackComparisonHists") << "monitoredTracksHandle not found, skipping event";
109  return;
110  }
111  reco::TrackCollection monitoredTracks = *monitoredTracksHandle;
112 
113  edm::Handle<reco::BeamSpot> monitoredBSHandle;
114  iEvent.getByToken(monitoredBSToken_, monitoredBSHandle);
115  if (!monitoredTracksHandle.isValid()) {
116  edm::LogError("TrackToTrackComparisonHists") << "monitoredBSHandle not found, skipping event";
117  return;
118  }
119  reco::BeamSpot monitoredBS = *monitoredBSHandle;
120 
121  edm::Handle<reco::VertexCollection> monitoredPVHandle;
122  iEvent.getByToken(monitoredPVToken_, monitoredPVHandle);
123  if (!monitoredPVHandle.isValid()) {
124  edm::LogError("TrackToTrackComparisonHists") << "monitoredPVHandle not found, skipping event";
125  return;
126  }
127  if (monitoredPVHandle->empty()) {
128  edm::LogInfo("TrackToTrackComparisonHists") << "monitoredPVHandle->size is 0 ";
129  return;
130  }
131  reco::Vertex monitoredPV = monitoredPVHandle->at(0);
132 
133  edm::LogInfo("TrackToTrackComparisonHists")
134  << "analyzing " << monitoredTrackInputTag_.process() << ":" << monitoredTrackInputTag_.label() << ":"
135  << monitoredTrackInputTag_.instance() << " w.r.t. " << referenceTrackInputTag_.process() << ":"
137 
138  //
139  // Build the dR maps
140  //
141  idx2idxByDoubleColl monitored2referenceColl;
142  fillMap(monitoredTracks, referenceTracks, monitored2referenceColl, dRmin_);
143 
144  idx2idxByDoubleColl reference2monitoredColl;
145  fillMap(referenceTracks, monitoredTracks, reference2monitoredColl, dRmin_);
146 
147  unsigned int nReferenceTracks(0); // Counts the number of refernce tracks
148  unsigned int nMatchedReferenceTracks(0); // Counts the number of matched refernce tracks
149  unsigned int nMonitoredTracks(0); // Counts the number of monitored tracks
150  unsigned int nUnmatchedMonitoredTracks(0); // Counts the number of unmatched monitored tracks
151 
152  //
153  // loop over reference tracks
154  //
155  LogDebug("TrackToTrackComparisonHists") << "\n# of tracks (reference): " << referenceTracks.size() << "\n";
156  for (idx2idxByDoubleColl::const_iterator pItr = reference2monitoredColl.begin(), eItr = reference2monitoredColl.end();
157  pItr != eItr;
158  ++pItr) {
159  nReferenceTracks++;
160  int trackIdx = pItr->first;
161  reco::Track track = referenceTracks.at(trackIdx);
162 
163  float dzWRTpv = track.dz(referencePV.position());
164  if (fabs(dzWRTpv) > dzWRTPvCut_)
165  continue;
166 
167  fill_generic_tracks_histos(*&referenceTracksMEs_, &track, &referenceBS, &referencePV);
168 
169  std::map<double, int> trackDRmap = pItr->second;
170  if (trackDRmap.empty()) {
173  continue;
174  }
175 
176  double dRmin = trackDRmap.begin()->first;
177  (referenceTracksMEs_.h_dRmin)->Fill(dRmin);
179 
180  bool matched = false;
181  if (dRmin < dRmin_)
182  matched = true;
183 
184  if (matched) {
185  nMatchedReferenceTracks++;
186  fill_generic_tracks_histos(*&matchedReferenceTracksMEs_, &track, &referenceBS, &referencePV);
189 
190  int matchedTrackIndex = trackDRmap[dRmin];
191  reco::Track matchedTrack = monitoredTracks.at(matchedTrackIndex);
192  fill_matching_tracks_histos(*&matchTracksMEs_, &track, &matchedTrack, &referenceBS, &referencePV);
193  }
194 
195  } // Over reference tracks
196 
197  //
198  // loop over monitoed tracks
199  //
200  LogDebug("TrackToTrackComparisonHists") << "\n# of tracks (monitored): " << monitoredTracks.size() << "\n";
201  for (idx2idxByDoubleColl::const_iterator pItr = monitored2referenceColl.begin(), eItr = monitored2referenceColl.end();
202  pItr != eItr;
203  ++pItr) {
204  nMonitoredTracks++;
205  int trackIdx = pItr->first;
206  reco::Track track = monitoredTracks.at(trackIdx);
207 
208  float dzWRTpv = track.dz(monitoredPV.position());
209  if (fabs(dzWRTpv) > dzWRTPvCut_)
210  continue;
211 
212  fill_generic_tracks_histos(*&monitoredTracksMEs_, &track, &monitoredBS, &monitoredPV);
213 
214  std::map<double, int> trackDRmap = pItr->second;
215  if (trackDRmap.empty()) {
218  continue;
219  }
220 
221  double dRmin = trackDRmap.begin()->first;
222  (monitoredTracksMEs_.h_dRmin)->Fill(dRmin);
224 
225  bool matched = false;
226  if (dRmin < dRmin_)
227  matched = true;
228 
229  if (!matched) {
230  nUnmatchedMonitoredTracks++;
231  fill_generic_tracks_histos(*&unMatchedMonitoredTracksMEs_, &track, &monitoredBS, &monitoredPV);
234  }
235 
236  } // over monitoed tracks
237 
238  edm::LogInfo("TrackToTrackComparisonHists")
239  << "Total reference tracks: " << nReferenceTracks << "\n"
240  << "Total matched reference tracks: " << nMatchedReferenceTracks << "\n"
241  << "Total monitored tracks: " << nMonitoredTracks << "\n"
242  << "Total unMatched monitored tracks: " << nUnmatchedMonitoredTracks << "\n";
243 }
244 
246  edm::Run const& iRun,
247  edm::EventSetup const& iSetup) {
249  genTriggerEventFlag_->initRun(iRun, iSetup);
250 
251  // check if every HLT path specified has a valid match in the HLT Menu
253  (genTriggerEventFlag_ && genTriggerEventFlag_->on() && genTriggerEventFlag_->allHLTPathsAreValid());
254 
255  // if valid HLT paths are required,
256  // create DQM outputs only if all paths are valid
258  return;
259  }
260 
262 
263  bookHistos(ibooker, referenceTracksMEs_, "ref", dir);
264  bookHistos(ibooker, matchedReferenceTracksMEs_, "ref_matched", dir);
265 
266  bookHistos(ibooker, monitoredTracksMEs_, "mon", dir);
267  bookHistos(ibooker, unMatchedMonitoredTracksMEs_, "mon_unMatched", dir);
268 
269  book_matching_tracks_histos(ibooker, matchTracksMEs_, "matches", dir);
270 }
271 
272 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
275 
276  desc.add<bool>("requireValidHLTPaths", true);
277 
278  desc.add<edm::InputTag>("monitoredTrack", edm::InputTag("hltMergedTracks"));
279  desc.add<edm::InputTag>("monitoredBeamSpot", edm::InputTag("hltOnlineBeamSpot"));
280  desc.add<edm::InputTag>("monitoredPrimaryVertices", edm::InputTag("hltVerticesPFSelector"));
281 
282  desc.add<edm::InputTag>("referenceTrack", edm::InputTag("generalTracks"));
283  desc.add<edm::InputTag>("referenceBeamSpot", edm::InputTag("offlineBeamSpot"));
284  desc.add<edm::InputTag>("referencePrimaryVertices", edm::InputTag("offlinePrimaryVertices"));
285 
286  desc.add<std::string>("topDirName", "HLT/Tracking/ValidationWRTOffline");
287  desc.add<double>("dRmin", 0.002);
288 
289  desc.add<double>("pTCutForPlateau", 0.9);
290  desc.add<double>("dxyCutForPlateau", 2.5);
291  desc.add<double>("dzWRTPvCut", 1e6);
292 
293  edm::ParameterSetDescription genericTriggerEventPSet;
294  genericTriggerEventPSet.add<bool>("andOr", false);
295  genericTriggerEventPSet.add<edm::InputTag>("dcsInputTag", edm::InputTag("scalersRawToDigi"));
296  genericTriggerEventPSet.add<std::vector<int> >("dcsPartitions", {24, 25, 26, 27, 28, 29}); // 24-27: strip, 28-29
297  genericTriggerEventPSet.add<bool>("andOrDcs", false);
298  genericTriggerEventPSet.add<bool>("errorReplyDcs", true);
299  genericTriggerEventPSet.add<std::string>("dbLabel", "");
300  genericTriggerEventPSet.add<bool>("andOrHlt", true); // True:=OR; False:=AND
301  genericTriggerEventPSet.add<edm::InputTag>("hltInputTag", edm::InputTag("TriggerResults::HLT"));
302  genericTriggerEventPSet.add<std::vector<std::string> >("hltPaths", {});
303  genericTriggerEventPSet.add<std::string>("hltDBKey", "");
304  genericTriggerEventPSet.add<bool>("errorReplyHlt", false);
305  genericTriggerEventPSet.add<unsigned int>("verbosityLevel", 1);
306  desc.add<edm::ParameterSetDescription>("genericTriggerEventPSet", genericTriggerEventPSet);
307 
309  fillHistoPSetDescription(histoPSet);
310  desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);
311 
312  descriptions.add("trackToTrackComparisonHists", desc);
313 }
314 
316  reco::TrackCollection tracks2,
317  idx2idxByDoubleColl& map,
318  float dRMin) {
319  //
320  // loop on tracks1
321  //
322  int i = 0;
323  for (const auto& track1 : tracks1) {
324  std::map<double, int> tmp;
325  int j = 0;
326  float smallest_dR = 1e9;
327  int smallest_dR_j = -1;
328 
329  //
330  // loop on tracks2
331  //
332  for (const auto& track2 : tracks2) {
333  double dR = reco::deltaR(track1.eta(), track1.phi(), track2.eta(), track2.phi());
334 
335  if (dR < smallest_dR) {
336  smallest_dR = dR;
337  smallest_dR_j = j;
338  }
339 
340  if (dR < dRMin) {
341  tmp[dR] = j;
342  }
343 
344  j++;
345  }
346 
347  //
348  // If there are no tracks that pass the dR store the smallest (for debugging/validating matching)
349  //
350  if (tmp.empty())
351  tmp[smallest_dR] = smallest_dR_j;
352 
353  map.push_back(std::make_pair(i, tmp));
354  i++;
355  }
356 }
357 
359  generalME& mes,
360  TString label,
361  std::string& dir) {
362  book_generic_tracks_histos(ibooker, mes, label, dir);
363 }
364 
366  generalME& mes,
367  TString label,
368  std::string& dir) {
369  ibooker.cd();
370  ibooker.setCurrentFolder(dir);
371  (mes.h_pt) = ibooker.book1D(label + "_pt", "track p_{T}", Pt_nbin, Pt_rangeMin, Pt_rangeMax);
372  (mes.h_eta) = ibooker.book1D(label + "_eta", "track pseudorapidity", Eta_nbin, Eta_rangeMin, Eta_rangeMax);
373  (mes.h_phi) = ibooker.book1D(label + "_phi", "track #phi", Phi_nbin, Phi_rangeMin, Phi_rangeMax);
374  (mes.h_dxy) =
375  ibooker.book1D(label + "_dxy", "track transverse dca to beam spot", Dxy_nbin, Dxy_rangeMin, Dxy_rangeMax);
376  (mes.h_dz) = ibooker.book1D(label + "_dz", "track longitudinal dca to beam spot", Dz_nbin, Dz_rangeMin, Dz_rangeMax);
377  (mes.h_dxyWRTpv) = ibooker.book1D(
378  label + "_dxyWRTpv", "track transverse dca to primary vertex", Dxy_nbin, Dxy_rangeMin, Dxy_rangeMax);
379  (mes.h_dzWRTpv) = ibooker.book1D(
380  label + "_dzWRTpv", "track longitudinal dca to primary vertex", Dz_nbin, 0.1 * Dz_rangeMin, 0.1 * Dz_rangeMax);
381  (mes.h_charge) = ibooker.book1D(label + "_charge", "track charge", 5, -2, 2);
382  (mes.h_hits) = ibooker.book1D(label + "_hits", "track number of hits", 35, -0.5, 34.5);
383  (mes.h_dRmin) = ibooker.book1D(label + "_dRmin", "track min dR", 100, 0., 0.01);
384  (mes.h_dRmin_l) = ibooker.book1D(label + "_dRmin_l", "track min dR", 100, 0., 0.4);
385 
386  (mes.h_pt_vs_eta) = ibooker.book2D(label + "_ptVSeta",
387  "track p_{T} vs #eta",
388  Eta_nbin,
389  Eta_rangeMin,
390  Eta_rangeMax,
391  Pt_nbin,
392  Pt_rangeMin,
393  Pt_rangeMax);
394 }
395 
397  matchingME& mes,
398  TString label,
399  std::string& dir) {
400  ibooker.cd();
401  ibooker.setCurrentFolder(dir);
402 
403  (mes.h_hits_vs_hits) = ibooker.book2D(
404  label + "_hits_vs_hits", "monitored track # hits vs reference track # hits", 35, -0.5, 34.5, 35, -0.5, 34.5);
405  (mes.h_pt_vs_pt) = ibooker.book2D(label + "_pt_vs_pt",
406  "monitored track p_{T} vs reference track p_{T}",
407  Pt_nbin,
408  Pt_rangeMin,
409  Pt_rangeMax,
410  Pt_nbin,
411  Pt_rangeMin,
412  Pt_rangeMax);
413  (mes.h_eta_vs_eta) = ibooker.book2D(label + "_eta_vs_eta",
414  "monitored track #eta vs reference track #eta",
415  Eta_nbin,
416  Eta_rangeMin,
417  Eta_rangeMax,
418  Eta_nbin,
419  Eta_rangeMin,
420  Eta_rangeMax);
421  (mes.h_phi_vs_phi) = ibooker.book2D(label + "_phi_vs_phi",
422  "monitored track #phi vs reference track #phi",
423  Phi_nbin,
424  Phi_rangeMin,
425  Phi_rangeMax,
426  Phi_nbin,
427  Phi_rangeMin,
428  Phi_rangeMax);
429 
430  (mes.h_dPt) = ibooker.book1D(label + "_dPt", "#Delta track #P_T", ptRes_nbin, ptRes_rangeMin, ptRes_rangeMax);
431  (mes.h_dEta) = ibooker.book1D(label + "_dEta", "#Delta track #eta", etaRes_nbin, etaRes_rangeMin, etaRes_rangeMax);
432  (mes.h_dPhi) = ibooker.book1D(label + "_dPhi", "#Delta track #phi", phiRes_nbin, phiRes_rangeMin, phiRes_rangeMax);
433  (mes.h_dDxy) = ibooker.book1D(
434  label + "_dDxy", "#Delta track transverse dca to beam spot", dxyRes_nbin, dxyRes_rangeMin, dxyRes_rangeMax);
435  (mes.h_dDz) = ibooker.book1D(
436  label + "_dDz", "#Delta track longitudinal dca to beam spot", dzRes_nbin, dzRes_rangeMin, dzRes_rangeMax);
437  (mes.h_dDxyWRTpv) = ibooker.book1D(label + "_dDxyWRTpv",
438  "#Delta track transverse dca to primary vertex ",
439  dxyRes_nbin,
442  (mes.h_dDzWRTpv) = ibooker.book1D(label + "_dDzWRTpv",
443  "#Delta track longitudinal dca to primary vertex",
444  dzRes_nbin,
447  (mes.h_dCharge) = ibooker.book1D(label + "_dCharge", "#Delta track charge", 5, -2.5, 2.5);
448  (mes.h_dHits) = ibooker.book1D(label + "_dHits", "#Delta track number of hits", 39, -19.5, 19.5);
449 }
450 
452  generalME& mes, reco::Track* trk, reco::BeamSpot* bs, reco::Vertex* pv, bool requirePlateau) {
453  float pt = trk->pt();
454  float eta = trk->eta();
455  float phi = trk->phi();
456  float dxy = trk->dxy(bs->position());
457  float dz = trk->dz(bs->position());
458  float dxyWRTpv = trk->dxy(pv->position());
459  float dzWRTpv = trk->dz(pv->position());
460  float charge = trk->charge();
461  float nhits = trk->hitPattern().numberOfValidHits();
462 
463  bool dxyOnPlateau = (fabs(dxyWRTpv) < dxyCutForPlateau_);
464  bool pTOnPlateau = (pt > pTCutForPlateau_);
465 
466  if (dxyOnPlateau || !requirePlateau) {
467  (mes.h_pt)->Fill(pt);
468  }
469 
470  if ((pTOnPlateau && dxyOnPlateau) || !requirePlateau) {
471  (mes.h_eta)->Fill(eta);
472  (mes.h_phi)->Fill(phi);
473  (mes.h_dz)->Fill(dz);
474  (mes.h_dzWRTpv)->Fill(dzWRTpv);
475  (mes.h_charge)->Fill(charge);
476  (mes.h_hits)->Fill(nhits);
477  }
478 
479  if (pTOnPlateau || !requirePlateau) {
480  (mes.h_dxy)->Fill(dxy);
481  (mes.h_dxyWRTpv)->Fill(dxyWRTpv);
482  }
483 
484  (mes.h_pt_vs_eta)->Fill(eta, pt);
485 }
486 
489  float mon_pt = mon->pt();
490  float mon_eta = mon->eta();
491  float mon_phi = mon->phi();
492  float mon_dxy = mon->dxy(bs->position());
493  float mon_dz = mon->dz(bs->position());
494  float mon_dxyWRTpv = mon->dxy(pv->position());
495  float mon_dzWRTpv = mon->dz(pv->position());
496  float mon_charge = mon->charge();
497  float mon_nhits = mon->hitPattern().numberOfValidHits();
498 
499  float ref_pt = ref->pt();
500  float ref_eta = ref->eta();
501  float ref_phi = ref->phi();
502  float ref_dxy = ref->dxy(bs->position());
503  float ref_dz = ref->dz(bs->position());
504  float ref_dxyWRTpv = ref->dxy(pv->position());
505  float ref_dzWRTpv = ref->dz(pv->position());
506  float ref_charge = ref->charge();
507  float ref_nhits = ref->hitPattern().numberOfValidHits();
508 
509  (mes.h_hits_vs_hits)->Fill(ref_nhits, mon_nhits);
510  (mes.h_pt_vs_pt)->Fill(ref_pt, mon_pt);
511  (mes.h_eta_vs_eta)->Fill(ref_eta, mon_eta);
512  (mes.h_phi_vs_phi)->Fill(ref_phi, mon_phi);
513 
514  (mes.h_dPt)->Fill(ref_pt - mon_pt);
515  (mes.h_dEta)->Fill(ref_eta - mon_eta);
516  (mes.h_dPhi)->Fill(ref_phi - mon_phi);
517  (mes.h_dDxy)->Fill(ref_dxy - mon_dxy);
518  (mes.h_dDz)->Fill(ref_dz - mon_dz);
519  (mes.h_dDxyWRTpv)->Fill(ref_dxyWRTpv - mon_dxyWRTpv);
520  (mes.h_dDzWRTpv)->Fill(ref_dzWRTpv - mon_dzWRTpv);
521  (mes.h_dCharge)->Fill(ref_charge - mon_charge);
522  (mes.h_dHits)->Fill(ref_nhits - mon_nhits);
523 }
524 
526  const edm::ParameterSet& pset = iConfig.getParameter<edm::ParameterSet>("histoPSet");
527 
528  Eta_rangeMin = pset.getParameter<double>("Eta_rangeMin");
529  Eta_rangeMax = pset.getParameter<double>("Eta_rangeMax");
530  Eta_nbin = pset.getParameter<unsigned int>("Eta_nbin");
531 
532  Pt_rangeMin = pset.getParameter<double>("Pt_rangeMin");
533  Pt_rangeMax = pset.getParameter<double>("Pt_rangeMax");
534  Pt_nbin = pset.getParameter<unsigned int>("Pt_nbin");
535 
536  Phi_rangeMin = pset.getParameter<double>("Phi_rangeMin");
537  Phi_rangeMax = pset.getParameter<double>("Phi_rangeMax");
538  Phi_nbin = pset.getParameter<unsigned int>("Phi_nbin");
539 
540  Dxy_rangeMin = pset.getParameter<double>("Dxy_rangeMin");
541  Dxy_rangeMax = pset.getParameter<double>("Dxy_rangeMax");
542  Dxy_nbin = pset.getParameter<unsigned int>("Dxy_nbin");
543 
544  Dz_rangeMin = pset.getParameter<double>("Dz_rangeMin");
545  Dz_rangeMax = pset.getParameter<double>("Dz_rangeMax");
546  Dz_nbin = pset.getParameter<unsigned int>("Dz_nbin");
547 
548  ptRes_rangeMin = pset.getParameter<double>("ptRes_rangeMin");
549  ptRes_rangeMax = pset.getParameter<double>("ptRes_rangeMax");
550  ptRes_nbin = pset.getParameter<unsigned int>("ptRes_nbin");
551 
552  phiRes_rangeMin = pset.getParameter<double>("phiRes_rangeMin");
553  phiRes_rangeMax = pset.getParameter<double>("phiRes_rangeMax");
554  phiRes_nbin = pset.getParameter<unsigned int>("phiRes_nbin");
555 
556  etaRes_rangeMin = pset.getParameter<double>("etaRes_rangeMin");
557  etaRes_rangeMax = pset.getParameter<double>("etaRes_rangeMax");
558  etaRes_nbin = pset.getParameter<unsigned int>("etaRes_nbin");
559 
560  dxyRes_rangeMin = pset.getParameter<double>("dxyRes_rangeMin");
561  dxyRes_rangeMax = pset.getParameter<double>("dxyRes_rangeMax");
562  dxyRes_nbin = pset.getParameter<unsigned int>("dxyRes_nbin");
563 
564  dzRes_rangeMin = pset.getParameter<double>("dzRes_rangeMin");
565  dzRes_rangeMax = pset.getParameter<double>("dzRes_rangeMax");
566  dzRes_nbin = pset.getParameter<unsigned int>("dzRes_nbin");
567 }
568 
570  pset.add<double>("Eta_rangeMin", -2.5);
571  pset.add<double>("Eta_rangeMax", 2.5);
572  pset.add<unsigned int>("Eta_nbin", 50);
573 
574  pset.add<double>("Pt_rangeMin", 0.1);
575  pset.add<double>("Pt_rangeMax", 100.0);
576  pset.add<unsigned int>("Pt_nbin", 1000);
577 
578  pset.add<double>("Phi_rangeMin", -3.1416);
579  pset.add<double>("Phi_rangeMax", 3.1416);
580  pset.add<unsigned int>("Phi_nbin", 36);
581 
582  pset.add<double>("Dxy_rangeMin", -1.0);
583  pset.add<double>("Dxy_rangeMax", 1.0);
584  pset.add<unsigned int>("Dxy_nbin", 300);
585 
586  pset.add<double>("Dz_rangeMin", -30.0);
587  pset.add<double>("Dz_rangeMax", 30.0);
588  pset.add<unsigned int>("Dz_nbin", 60);
589 
590  pset.add<double>("ptRes_rangeMin", -0.1);
591  pset.add<double>("ptRes_rangeMax", 0.1);
592  pset.add<unsigned int>("ptRes_nbin", 100);
593 
594  pset.add<double>("phiRes_rangeMin", -0.01);
595  pset.add<double>("phiRes_rangeMax", 0.01);
596  pset.add<unsigned int>("phiRes_nbin", 300);
597 
598  pset.add<double>("etaRes_rangeMin", -0.01);
599  pset.add<double>("etaRes_rangeMax", 0.01);
600  pset.add<unsigned int>("etaRes_nbin", 300);
601 
602  pset.add<double>("dxyRes_rangeMin", -0.05);
603  pset.add<double>("dxyRes_rangeMax", 0.05);
604  pset.add<unsigned int>("dxyRes_nbin", 500);
605 
606  pset.add<double>("dzRes_rangeMin", -0.05);
607  pset.add<double>("dzRes_rangeMax", 0.05);
608  pset.add<unsigned int>("dzRes_nbin", 150);
609 }
edm::EDGetTokenT< reco::TrackCollection > monitoredTrackToken_
void beginJob(const edm::EventSetup &iSetup)
edm::EDGetTokenT< reco::BeamSpot > referenceBSToken_
edm::EDGetTokenT< reco::VertexCollection > monitoredPVToken_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
void bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
int numberOfValidHits() const
Definition: HitPattern.h:811
edm::EDGetTokenT< reco::TrackCollection > referenceTrackToken_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
Log< level::Error, false > LogError
const Point & position() const
position
Definition: Vertex.h:127
char const * label
int iEvent
Definition: GenABIO.cc:224
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
std::unique_ptr< GenericTriggerEventFlag > genTriggerEventFlag_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void fill_generic_tracks_histos(generalME &mes, reco::Track *trk, reco::BeamSpot *bs, reco::Vertex *pv, bool requirePlateau=true)
double pt() const
track transverse momentum
Definition: TrackBase.h:637
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
static void fillHistoPSetDescription(edm::ParameterSetDescription &pset)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
std::vector< std::pair< int, std::map< double, int > > > idx2idxByDoubleColl
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
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
Log< level::Info, false > LogInfo
void initialize_parameter(const edm::ParameterSet &iConfig)
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void fill_matching_tracks_histos(matchingME &mes, reco::Track *mon, reco::Track *ref, reco::BeamSpot *bs, reco::Vertex *pv)
std::string const & label() const
Definition: InputTag.h:36
edm::EDGetTokenT< reco::BeamSpot > monitoredBSToken_
std::string const & process() const
Definition: InputTag.h:40
edm::EDGetTokenT< reco::VertexCollection > referencePVToken_
void fillMap(reco::TrackCollection tracks1, reco::TrackCollection tracks2, idx2idxByDoubleColl &map, float dRMin)
void book_matching_tracks_histos(DQMStore::IBooker &ibooker, matchingME &mes, TString label, std::string &dir)
void book_generic_tracks_histos(DQMStore::IBooker &ibooker, generalME &mes, TString label, std::string &dir)
int charge() const
track electric charge
Definition: TrackBase.h:596
const Point & position() const
position
Definition: BeamSpot.h:59
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
tmp
align.sh
Definition: createJobs.py:716
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
std::string const & instance() const
Definition: InputTag.h:37
void bookHistos(DQMStore::IBooker &ibooker, generalME &mes, TString label, std::string &dir)
Definition: Run.h:45
#define LogDebug(id)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
TrackToTrackComparisonHists(const edm::ParameterSet &)