CMS 3D CMS Logo

MultiTrackValidator.cc
Go to the documentation of this file.
3 
7 
26 
32 #include<type_traits>
33 
34 
35 #include "TMath.h"
36 #include <TF1.h>
39 //#include <iostream>
40 
41 using namespace std;
42 using namespace edm;
43 
45 namespace {
46  bool trackSelected(unsigned char mask, unsigned char qual) {
47  return mask & 1<<qual;
48  }
49 
50 }
51 
53  associators(pset.getUntrackedParameter< std::vector<edm::InputTag> >("associators")),
54  label(pset.getParameter< std::vector<edm::InputTag> >("label")),
55  parametersDefiner(pset.getParameter<std::string>("parametersDefiner")),
56  parametersDefinerIsCosmic_(parametersDefiner == "CosmicParametersDefinerForTP"),
57  ignoremissingtkcollection_(pset.getUntrackedParameter<bool>("ignoremissingtrackcollection",false)),
58  useAssociators_(pset.getParameter< bool >("UseAssociators")),
59  calculateDrSingleCollection_(pset.getUntrackedParameter<bool>("calculateDrSingleCollection")),
60  doPlotsOnlyForTruePV_(pset.getUntrackedParameter<bool>("doPlotsOnlyForTruePV")),
61  doSummaryPlots_(pset.getUntrackedParameter<bool>("doSummaryPlots")),
62  doSimPlots_(pset.getUntrackedParameter<bool>("doSimPlots")),
63  doSimTrackPlots_(pset.getUntrackedParameter<bool>("doSimTrackPlots")),
64  doRecoTrackPlots_(pset.getUntrackedParameter<bool>("doRecoTrackPlots")),
65  dodEdxPlots_(pset.getUntrackedParameter<bool>("dodEdxPlots")),
66  doPVAssociationPlots_(pset.getUntrackedParameter<bool>("doPVAssociationPlots")),
67  doSeedPlots_(pset.getUntrackedParameter<bool>("doSeedPlots")),
68  doMVAPlots_(pset.getUntrackedParameter<bool>("doMVAPlots")),
69  simPVMaxZ_(pset.getUntrackedParameter<double>("simPVMaxZ"))
70 {
71  const edm::InputTag& label_tp_effic_tag = pset.getParameter< edm::InputTag >("label_tp_effic");
72  const edm::InputTag& label_tp_fake_tag = pset.getParameter< edm::InputTag >("label_tp_fake");
73 
74  if(pset.getParameter<bool>("label_tp_effic_refvector")) {
75  label_tp_effic_refvector = consumes<TrackingParticleRefVector>(label_tp_effic_tag);
76  }
77  else {
78  label_tp_effic = consumes<TrackingParticleCollection>(label_tp_effic_tag);
79  }
80  if(pset.getParameter<bool>("label_tp_fake_refvector")) {
81  label_tp_fake_refvector = consumes<TrackingParticleRefVector>(label_tp_fake_tag);
82  }
83  else {
84  label_tp_fake = consumes<TrackingParticleCollection>(label_tp_fake_tag);
85  }
86  label_pileupinfo = consumes<std::vector<PileupSummaryInfo> >(pset.getParameter< edm::InputTag >("label_pileupinfo"));
87  for(const auto& tag: pset.getParameter<std::vector<edm::InputTag>>("sim")) {
88  simHitTokens_.push_back(consumes<std::vector<PSimHit>>(tag));
89  }
90 
91  std::vector<edm::InputTag> doResolutionPlotsForLabels = pset.getParameter<std::vector<edm::InputTag> >("doResolutionPlotsForLabels");
92  doResolutionPlots_.reserve(label.size());
93  for (auto& itag : label) {
94  labelToken.push_back(consumes<edm::View<reco::Track> >(itag));
95  const bool doResol = doResolutionPlotsForLabels.empty() || (std::find(cbegin(doResolutionPlotsForLabels), cend(doResolutionPlotsForLabels), itag) != cend(doResolutionPlotsForLabels));
96  doResolutionPlots_.push_back(doResol);
97  }
98  { // check for duplicates
99  auto labelTmp = edm::vector_transform(label, [&](const edm::InputTag& tag) { return tag.label(); });
100  std::sort(begin(labelTmp), end(labelTmp));
102  const std::string* prev = &empty;
103  for(const std::string& l: labelTmp) {
104  if(l == *prev) {
105  throw cms::Exception("Configuration") << "Duplicate InputTag in labels: " << l;
106  }
107  prev = &l;
108  }
109  }
110 
112  bsSrc = consumes<reco::BeamSpot>(beamSpotTag);
113 
114  ParameterSet psetForHistoProducerAlgo = pset.getParameter<ParameterSet>("histoProducerAlgoBlock");
115  histoProducerAlgo_ = std::make_unique<MTVHistoProducerAlgoForTracker>(psetForHistoProducerAlgo, doSeedPlots_);
116 
117  dirName_ = pset.getParameter<std::string>("dirName");
118 
119  tpNLayersToken_ = consumes<edm::ValueMap<unsigned int> >(pset.getParameter<edm::InputTag>("label_tp_nlayers"));
120  tpNPixelLayersToken_ = consumes<edm::ValueMap<unsigned int> >(pset.getParameter<edm::InputTag>("label_tp_npixellayers"));
121  tpNStripStereoLayersToken_ = consumes<edm::ValueMap<unsigned int> >(pset.getParameter<edm::InputTag>("label_tp_nstripstereolayers"));
122 
123  if(dodEdxPlots_) {
124  m_dEdx1Tag = consumes<edm::ValueMap<reco::DeDxData> >(pset.getParameter< edm::InputTag >("dEdx1Tag"));
125  m_dEdx2Tag = consumes<edm::ValueMap<reco::DeDxData> >(pset.getParameter< edm::InputTag >("dEdx2Tag"));
126  }
127 
128  label_tv = consumes<TrackingVertexCollection>(pset.getParameter< edm::InputTag >("label_tv"));
130  recoVertexToken_ = consumes<edm::View<reco::Vertex> >(pset.getUntrackedParameter<edm::InputTag>("label_vertex"));
131  vertexAssociatorToken_ = consumes<reco::VertexToTrackingVertexAssociator>(pset.getUntrackedParameter<edm::InputTag>("vertexAssociator"));
132  }
133 
134  if(doMVAPlots_) {
136  auto mvaPSet = pset.getUntrackedParameter<edm::ParameterSet>("mvaLabels");
137  for(size_t iIter=0; iIter<labelToken.size(); ++iIter) {
139  labelsForToken(labelToken[iIter], labels);
140  if(mvaPSet.exists(labels.module)) {
141  mvaQualityCollectionTokens_[iIter] = edm::vector_transform(mvaPSet.getUntrackedParameter<std::vector<std::string> >(labels.module),
142  [&](const std::string& tag) {
143  return std::make_tuple(consumes<MVACollection>(edm::InputTag(tag, "MVAValues")),
144  consumes<QualityMaskCollection>(edm::InputTag(tag, "QualityMasks")));
145  });
146  }
147  }
148  }
149 
150  tpSelector = TrackingParticleSelector(pset.getParameter<double>("ptMinTP"),
151  pset.getParameter<double>("ptMaxTP"),
152  pset.getParameter<double>("minRapidityTP"),
153  pset.getParameter<double>("maxRapidityTP"),
154  pset.getParameter<double>("tipTP"),
155  pset.getParameter<double>("lipTP"),
156  pset.getParameter<int>("minHitTP"),
157  pset.getParameter<bool>("signalOnlyTP"),
158  pset.getParameter<bool>("intimeOnlyTP"),
159  pset.getParameter<bool>("chargedOnlyTP"),
160  pset.getParameter<bool>("stableOnlyTP"),
161  pset.getParameter<std::vector<int> >("pdgIdTP"));
162 
164  pset.getParameter<double>("minRapidityTP"),
165  pset.getParameter<double>("maxRapidityTP"),
166  pset.getParameter<double>("tipTP"),
167  pset.getParameter<double>("lipTP"),
168  pset.getParameter<int>("minHitTP"),
169  pset.getParameter<bool>("chargedOnlyTP"),
170  pset.getParameter<std::vector<int> >("pdgIdTP"));
171 
172 
173  ParameterSet psetVsPhi = psetForHistoProducerAlgo.getParameter<ParameterSet>("TpSelectorForEfficiencyVsPhi");
174  dRtpSelector = TrackingParticleSelector(psetVsPhi.getParameter<double>("ptMin"),
175  psetVsPhi.getParameter<double>("ptMax"),
176  psetVsPhi.getParameter<double>("minRapidity"),
177  psetVsPhi.getParameter<double>("maxRapidity"),
178  psetVsPhi.getParameter<double>("tip"),
179  psetVsPhi.getParameter<double>("lip"),
180  psetVsPhi.getParameter<int>("minHit"),
181  psetVsPhi.getParameter<bool>("signalOnly"),
182  psetVsPhi.getParameter<bool>("intimeOnly"),
183  psetVsPhi.getParameter<bool>("chargedOnly"),
184  psetVsPhi.getParameter<bool>("stableOnly"),
185  psetVsPhi.getParameter<std::vector<int> >("pdgId"));
186 
188 
189  useGsf = pset.getParameter<bool>("useGsf");
190 
191  _simHitTpMapTag = mayConsume<SimHitTPAssociationProducer::SimHitTPAssociationList>(pset.getParameter<edm::InputTag>("simHitTpMapTag"));
192 
194  labelTokenForDrCalculation = consumes<edm::View<reco::Track> >(pset.getParameter<edm::InputTag>("trackCollectionForDrCalculation"));
195  }
196 
197  if(useAssociators_) {
198  for (auto const& src: associators) {
199  associatorTokens.push_back(consumes<reco::TrackToTrackingParticleAssociator>(src));
200  }
201  } else {
202  for (auto const& src: associators) {
203  associatormapStRs.push_back(consumes<reco::SimToRecoCollection>(src));
204  associatormapRtSs.push_back(consumes<reco::RecoToSimCollection>(src));
205  }
206  }
207 }
208 
209 
211 
212 
214 
215  const auto minColl = -0.5;
216  const auto maxColl = label.size()-0.5;
217  const auto nintColl = label.size();
218 
219  auto binLabels = [&](ConcurrentMonitorElement me) {
220  for(size_t i=0; i<label.size(); ++i) {
221  me.setBinLabel(i+1, label[i].label());
222  }
223  me.disableAlphanumeric();
224  return me;
225  };
226 
227  //Booking histograms concerning with simulated tracks
228  if(doSimPlots_) {
229  ibook.cd();
230  ibook.setCurrentFolder(dirName_ + "simulation");
231 
232  histoProducerAlgo_->bookSimHistos(ibook, histograms.histoProducerAlgo);
233 
234  ibook.cd();
235  ibook.setCurrentFolder(dirName_);
236  }
237 
238  for (unsigned int ww=0;ww<associators.size();ww++){
239  ibook.cd();
240  // FIXME: these need to be moved to a subdirectory whose name depends on the associator
241  ibook.setCurrentFolder(dirName_);
242 
243  if(doSummaryPlots_) {
244  if(doSimTrackPlots_) {
245  histograms.h_assoc_coll.push_back(binLabels( ibook.book1D("num_assoc(simToReco)_coll", "N of associated (simToReco) tracks vs track collection", nintColl, minColl, maxColl) ));
246  histograms.h_simul_coll.push_back(binLabels( ibook.book1D("num_simul_coll", "N of simulated tracks vs track collection", nintColl, minColl, maxColl) ));
247  }
248  if(doRecoTrackPlots_) {
249  histograms.h_reco_coll.push_back(binLabels( ibook.book1D("num_reco_coll", "N of reco track vs track collection", nintColl, minColl, maxColl) ));
250  histograms.h_assoc2_coll.push_back(binLabels( ibook.book1D("num_assoc(recoToSim)_coll", "N of associated (recoToSim) tracks vs track collection", nintColl, minColl, maxColl) ));
251  histograms.h_looper_coll.push_back(binLabels( ibook.book1D("num_duplicate_coll", "N of associated (recoToSim) looper tracks vs track collection", nintColl, minColl, maxColl) ));
252  histograms.h_pileup_coll.push_back(binLabels( ibook.book1D("num_pileup_coll", "N of associated (recoToSim) pileup tracks vs track collection", nintColl, minColl, maxColl) ));
253  }
254  }
255 
256  for (unsigned int www=0;www<label.size();www++){
257  ibook.cd();
258  InputTag algo = label[www];
259  string dirName=dirName_;
260  if (!algo.process().empty())
261  dirName+=algo.process()+"_";
262  if(!algo.label().empty())
263  dirName+=algo.label()+"_";
264  if(!algo.instance().empty())
265  dirName+=algo.instance()+"_";
266  if (dirName.find("Tracks")<dirName.length()){
267  dirName.replace(dirName.find("Tracks"),6,"");
268  }
269  string assoc= associators[ww].label();
270  if (assoc.find("Track")<assoc.length()){
271  assoc.replace(assoc.find("Track"),5,"");
272  }
273  dirName+=assoc;
274  std::replace(dirName.begin(), dirName.end(), ':', '_');
275 
276  ibook.setCurrentFolder(dirName);
277 
278  const bool doResolutionPlots = doResolutionPlots_[www];
279 
280  if(doSimTrackPlots_) {
281  histoProducerAlgo_->bookSimTrackHistos(ibook, histograms.histoProducerAlgo, doResolutionPlots);
282  if(doPVAssociationPlots_) histoProducerAlgo_->bookSimTrackPVAssociationHistos(ibook, histograms.histoProducerAlgo);
283  }
284 
285  //Booking histograms concerning with reconstructed tracks
286  if(doRecoTrackPlots_) {
287  histoProducerAlgo_->bookRecoHistos(ibook, histograms.histoProducerAlgo, doResolutionPlots);
288  if (dodEdxPlots_) histoProducerAlgo_->bookRecodEdxHistos(ibook, histograms.histoProducerAlgo);
289  if (doPVAssociationPlots_) histoProducerAlgo_->bookRecoPVAssociationHistos(ibook, histograms.histoProducerAlgo);
290  if (doMVAPlots_) histoProducerAlgo_->bookMVAHistos(ibook, histograms.histoProducerAlgo, mvaQualityCollectionTokens_[www].size());
291  }
292 
293  if(doSeedPlots_) {
294  histoProducerAlgo_->bookSeedHistos(ibook, histograms.histoProducerAlgo);
295  }
296  }//end loop www
297  }// end loop ww
298 }
299 
300 namespace {
301  void ensureEffIsSubsetOfFake(const TrackingParticleRefVector& eff, const TrackingParticleRefVector& fake) {
302  // If efficiency RefVector is empty, don't check the product ids
303  // as it will be 0:0 for empty. This covers also the case where
304  // both are empty. The case of fake being empty and eff not is an
305  // error.
306  if(eff.empty())
307  return;
308 
309  // First ensure product ids
310  if(eff.id() != fake.id()) {
311  throw cms::Exception("Configuration") << "Efficiency and fake TrackingParticle (refs) point to different collections (eff " << eff.id() << " fake " << fake.id() << "). This is not supported. Efficiency TP set must be the same or a subset of the fake TP set.";
312  }
313 
314  // Same technique as in associationMapFilterValues
315  edm::IndexSet fakeKeys;
316  fakeKeys.reserve(fake.size());
317  for(const auto& ref: fake) {
318  fakeKeys.insert(ref.key());
319  }
320 
321  for(const auto& ref: eff) {
322  if(!fakeKeys.has(ref.key())) {
323  throw cms::Exception("Configuration") << "Efficiency TrackingParticle " << ref.key() << " is not found from the set of fake TPs. This is not supported. The efficiency TP set must be the same or a subset of the fake TP set.";
324  }
325  }
326  }
327 }
328 
330  for(const auto& simV: *htv) {
331  if(simV.eventId().bunchCrossing() != 0) continue; // remove OOTPU
332  if(simV.eventId().event() != 0) continue; // pick the PV of hard scatter
333  return &(simV.position());
334  }
335  return nullptr;
336 }
337 
340  event.getByToken(recoVertexToken_, hvertex);
341 
343  event.getByToken(vertexAssociatorToken_, hvassociator);
344 
345  auto v_r2s = hvassociator->associateRecoToSim(hvertex, htv);
346  auto pvPtr = hvertex->refAt(0);
347  if(pvPtr->isFake() || pvPtr->ndof() < 0) // skip junk vertices
348  return nullptr;
349 
350  auto pvFound = v_r2s.find(pvPtr);
351  if(pvFound == v_r2s.end())
352  return nullptr;
353 
354  for(const auto& vertexRefQuality: pvFound->val) {
355  const TrackingVertex& tv = *(vertexRefQuality.first);
356  if(tv.eventId().event() == 0 && tv.eventId().bunchCrossing() == 0) {
357  return &(pvPtr->position());
358  }
359  }
360 
361  return nullptr;
362 }
363 
365  const TrackingParticleRefVector& tPCeff,
366  const ParametersDefinerForTP& parametersDefinerTP,
367  const edm::Event& event, const edm::EventSetup& setup,
368  const reco::BeamSpot& bs,
369  std::vector<std::tuple<TrackingParticle::Vector, TrackingParticle::Point> >& momVert_tPCeff,
370  std::vector<size_t>& selected_tPCeff) const {
371  selected_tPCeff.reserve(tPCeff.size());
372  momVert_tPCeff.reserve(tPCeff.size());
373  int nIntimeTPs = 0;
375  for(size_t j=0; j<tPCeff.size(); ++j) {
376  const TrackingParticleRef& tpr = tPCeff[j];
377 
378  TrackingParticle::Vector momentum = parametersDefinerTP.momentum(event,setup,tpr);
379  TrackingParticle::Point vertex = parametersDefinerTP.vertex(event,setup,tpr);
380  if(doSimPlots_) {
381  histoProducerAlgo_->fill_generic_simTrack_histos(histograms.histoProducerAlgo, momentum, vertex, tpr->eventId().bunchCrossing());
382  }
383  if(tpr->eventId().bunchCrossing() == 0)
384  ++nIntimeTPs;
385 
386  if(cosmictpSelector(tpr,&bs,event,setup)) {
387  selected_tPCeff.push_back(j);
388  momVert_tPCeff.emplace_back(momentum, vertex);
389  }
390  }
391  }
392  else {
393  size_t j=0;
394  for(auto const& tpr: tPCeff) {
395  const TrackingParticle& tp = *tpr;
396 
397  // TODO: do we want to fill these from all TPs that include IT
398  // and OOT (as below), or limit to IT+OOT TPs passing tpSelector
399  // (as it was before)? The latter would require another instance
400  // of tpSelector with intimeOnly=False.
401  if(doSimPlots_) {
402  histoProducerAlgo_->fill_generic_simTrack_histos(histograms.histoProducerAlgo, tp.momentum(), tp.vertex(), tp.eventId().bunchCrossing());
403  }
404  if(tp.eventId().bunchCrossing() == 0)
405  ++nIntimeTPs;
406 
407  if(tpSelector(tp)) {
408  selected_tPCeff.push_back(j);
409  TrackingParticle::Vector momentum = parametersDefinerTP.momentum(event,setup,tpr);
410  TrackingParticle::Point vertex = parametersDefinerTP.vertex(event,setup,tpr);
411  momVert_tPCeff.emplace_back(momentum, vertex);
412  }
413  ++j;
414  }
415  }
416  if(doSimPlots_) {
417  histoProducerAlgo_->fill_simTrackBased_histos(histograms.histoProducerAlgo, nIntimeTPs);
418  }
419 }
420 
421 
423  const std::vector<size_t>& selected_tPCeff,
424  DynArray<float>& dR_tPCeff) const {
425  float etaL[tPCeff.size()], phiL[tPCeff.size()];
426  size_t n_selTP_dr = 0;
427  for(size_t iTP: selected_tPCeff) {
428  //calculare dR wrt inclusive collection (also with PU, low pT, displaced)
429  auto const& tp2 = *(tPCeff[iTP]);
430  auto && p = tp2.momentum();
431  etaL[iTP] = etaFromXYZ(p.x(),p.y(),p.z());
432  phiL[iTP] = atan2f(p.y(),p.x());
433  }
434  for(size_t iTP1: selected_tPCeff) {
435  auto const& tp = *(tPCeff[iTP1]);
437  if(dRtpSelector(tp)) {//only for those needed for efficiency!
438  ++n_selTP_dr;
439  float eta = etaL[iTP1];
440  float phi = phiL[iTP1];
441  for(size_t iTP2: selected_tPCeff) {
442  //calculare dR wrt inclusive collection (also with PU, low pT, displaced)
443  if (iTP1==iTP2) {continue;}
444  auto dR_tmp = reco::deltaR2(eta, phi, etaL[iTP2], phiL[iTP2]);
445  if (dR_tmp<dR) dR=dR_tmp;
446  } // ttp2 (iTP)
447  }
448  dR_tPCeff[iTP1] = std::sqrt(dR);
449  } // tp
450  return n_selTP_dr;
451 }
452 
454  int i=0;
455  float etaL[trackCollectionDr.size()];
456  float phiL[trackCollectionDr.size()];
457  bool validL[trackCollectionDr.size()];
458  for (auto const & track2 : trackCollectionDr) {
459  auto && p = track2.momentum();
460  etaL[i] = etaFromXYZ(p.x(),p.y(),p.z());
461  phiL[i] = atan2f(p.y(),p.x());
462  validL[i] = !trackFromSeedFitFailed(track2);
463  ++i;
464  }
465  for(View<reco::Track>::size_type i=0; i<trackCollection.size(); ++i){
466  auto const & track = trackCollection[i];
469  auto && p = track.momentum();
470  float eta = etaFromXYZ(p.x(),p.y(),p.z());
471  float phi = atan2f(p.y(),p.x());
472  for(View<reco::Track>::size_type j=0; j<trackCollectionDr.size(); ++j){
473  if(!validL[j]) continue;
474  auto dR_tmp = reco::deltaR2(eta, phi, etaL[j], phiL[j]);
475  if ( (dR_tmp<dR) & (dR_tmp>std::numeric_limits<float>::min())) dR=dR_tmp;
476  }
477  }
478  dR_trk[i] = std::sqrt(dR);
479  }
480 }
481 
482 
484  using namespace reco;
485 
486  LogDebug("TrackValidator") << "\n====================================================" << "\n"
487  << "Analyzing new event" << "\n"
488  << "====================================================\n" << "\n";
489 
490 
491  edm::ESHandle<ParametersDefinerForTP> parametersDefinerTPHandle;
492  setup.get<TrackAssociatorRecord>().get(parametersDefiner,parametersDefinerTPHandle);
493  //Since we modify the object, we must clone it
494  auto parametersDefinerTP = parametersDefinerTPHandle->clone();
495 
497  setup.get<TrackerTopologyRcd>().get(httopo);
498  const TrackerTopology& ttopo = *httopo;
499 
500  // FIXME: we really need to move to edm::View for reading the
501  // TrackingParticles... Unfortunately it has non-trivial
502  // consequences on the associator/association interfaces etc.
503  TrackingParticleRefVector tmpTPeff;
504  TrackingParticleRefVector tmpTPfake;
505  const TrackingParticleRefVector *tmpTPeffPtr = nullptr;
506  const TrackingParticleRefVector *tmpTPfakePtr = nullptr;
507 
509  edm::Handle<TrackingParticleRefVector> TPCollectionHeffRefVector;
510 
511  const bool tp_effic_refvector = label_tp_effic.isUninitialized();
512  if(!tp_effic_refvector) {
513  event.getByToken(label_tp_effic, TPCollectionHeff);
514  for(size_t i=0, size=TPCollectionHeff->size(); i<size; ++i) {
515  tmpTPeff.push_back(TrackingParticleRef(TPCollectionHeff, i));
516  }
517  tmpTPeffPtr = &tmpTPeff;
518  }
519  else {
520  event.getByToken(label_tp_effic_refvector, TPCollectionHeffRefVector);
521  tmpTPeffPtr = TPCollectionHeffRefVector.product();
522  }
524  edm::Handle<TrackingParticleCollection> TPCollectionHfake ;
525  event.getByToken(label_tp_fake,TPCollectionHfake);
526  for(size_t i=0, size=TPCollectionHfake->size(); i<size; ++i) {
527  tmpTPfake.push_back(TrackingParticleRef(TPCollectionHfake, i));
528  }
529  tmpTPfakePtr = &tmpTPfake;
530  }
531  else {
532  edm::Handle<TrackingParticleRefVector> TPCollectionHfakeRefVector;
533  event.getByToken(label_tp_fake_refvector, TPCollectionHfakeRefVector);
534  tmpTPfakePtr = TPCollectionHfakeRefVector.product();
535  }
536 
537  TrackingParticleRefVector const & tPCeff = *tmpTPeffPtr;
538  TrackingParticleRefVector const & tPCfake = *tmpTPfakePtr;
539 
540  ensureEffIsSubsetOfFake(tPCeff, tPCfake);
541 
544  //warning: make sure the TP collection used in the map is the same used in the MTV!
545  event.getByToken(_simHitTpMapTag,simHitsTPAssoc);
546  parametersDefinerTP->initEvent(simHitsTPAssoc);
547  cosmictpSelector.initEvent(simHitsTPAssoc);
548  }
549 
550  // Find the sim PV and tak its position
552  event.getByToken(label_tv, htv);
553  const TrackingVertex::LorentzVector *theSimPVPosition = getSimPVPosition(htv);
554  if(simPVMaxZ_ >= 0) {
555  if(!theSimPVPosition) return;
556  if(std::abs(theSimPVPosition->z()) > simPVMaxZ_) return;
557  }
558 
559  // Check, when necessary, if reco PV matches to sim PV
560  const reco::Vertex::Point *thePVposition = nullptr;
562  thePVposition = getRecoPVPosition(event, htv);
563  if(doPlotsOnlyForTruePV_ && !thePVposition)
564  return;
565 
566  // Rest of the code assumes that if thePVposition is non-null, the
567  // PV-association histograms get filled. In above, the "nullness"
568  // is used to deliver the information if the reco PV is matched to
569  // the sim PV.
571  thePVposition = nullptr;
572  }
573 
574  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
575  event.getByToken(bsSrc,recoBeamSpotHandle);
576  reco::BeamSpot const & bs = *recoBeamSpotHandle;
577 
579  event.getByToken(label_pileupinfo,puinfoH);
580  PileupSummaryInfo puinfo;
581 
582  for (unsigned int puinfo_ite=0;puinfo_ite<(*puinfoH).size();++puinfo_ite){
583  if ((*puinfoH)[puinfo_ite].getBunchCrossing()==0){
584  puinfo=(*puinfoH)[puinfo_ite];
585  break;
586  }
587  }
588 
589  // Number of 3D layers for TPs
591  event.getByToken(tpNLayersToken_, tpNLayersH);
592  const auto& nLayers_tPCeff = *tpNLayersH;
593 
594  event.getByToken(tpNPixelLayersToken_, tpNLayersH);
595  const auto& nPixelLayers_tPCeff = *tpNLayersH;
596 
597  event.getByToken(tpNStripStereoLayersToken_, tpNLayersH);
598  const auto& nStripMonoAndStereoLayers_tPCeff = *tpNLayersH;
599 
600  // Precalculate TP selection (for efficiency), and momentum and vertex wrt PCA
601  //
602  // TODO: ParametersDefinerForTP ESProduct needs to be changed to
603  // EDProduct because of consumes.
604  //
605  // In principle, we could just precalculate the momentum and vertex
606  // wrt PCA for all TPs for once and put that to the event. To avoid
607  // repetitive calculations those should be calculated only once for
608  // each TP. That would imply that we should access TPs via Refs
609  // (i.e. View) in here, since, in general, the eff and fake TP
610  // collections can be different (and at least HI seems to use that
611  // feature). This would further imply that the
612  // RecoToSimCollection/SimToRecoCollection should be changed to use
613  // View<TP> instead of vector<TP>, and migrate everything.
614  //
615  // Or we could take only one input TP collection, and do another
616  // TP-selection to obtain the "fake" collection like we already do
617  // for "efficiency" TPs.
618  std::vector<size_t> selected_tPCeff;
619  std::vector<std::tuple<TrackingParticle::Vector, TrackingParticle::Point>> momVert_tPCeff;
620  tpParametersAndSelection(histograms, tPCeff, *parametersDefinerTP, event, setup, bs, momVert_tPCeff, selected_tPCeff);
621 
622  //calculate dR for TPs
623  declareDynArray(float, tPCeff.size(), dR_tPCeff);
624  size_t n_selTP_dr = tpDR(tPCeff, selected_tPCeff, dR_tPCeff);
625 
628  event.getByToken(labelTokenForDrCalculation, trackCollectionForDrCalculation);
629  }
630 
631  // dE/dx
632  // at some point this could be generalized, with a vector of tags and a corresponding vector of Handles
633  // I'm writing the interface such to take vectors of ValueMaps
634  std::vector<const edm::ValueMap<reco::DeDxData> *> v_dEdx;
635  if(dodEdxPlots_) {
638  event.getByToken(m_dEdx1Tag, dEdx1Handle);
639  event.getByToken(m_dEdx2Tag, dEdx2Handle);
640  v_dEdx.push_back(dEdx1Handle.product());
641  v_dEdx.push_back(dEdx2Handle.product());
642  }
643 
644  std::vector<const MVACollection *> mvaCollections;
645  std::vector<const QualityMaskCollection *> qualityMaskCollections;
646  std::vector<float> mvaValues;
647 
648  int w=0; //counter counting the number of sets of histograms
649  for (unsigned int ww=0;ww<associators.size();ww++){
650  // run value filtering of recoToSim map already here as it depends only on the association, not track collection
651  reco::SimToRecoCollection const * simRecCollPFull=nullptr;
652  reco::RecoToSimCollection const * recSimCollP=nullptr;
653  reco::RecoToSimCollection recSimCollL;
654  if(!useAssociators_) {
655  Handle<reco::SimToRecoCollection > simtorecoCollectionH;
656  event.getByToken(associatormapStRs[ww], simtorecoCollectionH);
657  simRecCollPFull = simtorecoCollectionH.product();
658 
659  Handle<reco::RecoToSimCollection > recotosimCollectionH;
660  event.getByToken(associatormapRtSs[ww],recotosimCollectionH);
661  recSimCollP = recotosimCollectionH.product();
662 
663  // We need to filter the associations of the fake-TrackingParticle
664  // collection only from RecoToSim collection, otherwise the
665  // RecoToSim histograms get false entries
666  recSimCollL = associationMapFilterValues(*recSimCollP, tPCfake);
667  recSimCollP = &recSimCollL;
668  }
669 
670  for (unsigned int www=0;www<label.size();www++, w++){ // need to increment w here, since there are many continues in the loop body
671  //
672  //get collections from the event
673  //
674  edm::Handle<View<Track> > trackCollectionHandle;
675  if(!event.getByToken(labelToken[www], trackCollectionHandle)&&ignoremissingtkcollection_)continue;
676  const edm::View<Track>& trackCollection = *trackCollectionHandle;
677 
678  reco::SimToRecoCollection const * simRecCollP=nullptr;
679  reco::SimToRecoCollection simRecCollL;
680 
681  //associate tracks
682  LogTrace("TrackValidator") << "Analyzing "
683  << label[www] << " with "
684  << associators[ww] <<"\n";
685  if(useAssociators_){
687  event.getByToken(associatorTokens[ww], theAssociator);
688 
689  // The associator interfaces really need to be fixed...
691  for(edm::View<Track>::size_type i=0; i<trackCollection.size(); ++i) {
692  trackRefs.push_back(trackCollection.refAt(i));
693  }
694 
695 
696  LogTrace("TrackValidator") << "Calling associateRecoToSim method" << "\n";
697  recSimCollL = theAssociator->associateRecoToSim(trackRefs, tPCfake);
698  recSimCollP = &recSimCollL;
699  LogTrace("TrackValidator") << "Calling associateSimToReco method" << "\n";
700  // It is necessary to do the association wrt. fake TPs,
701  // because this SimToReco association is used also for
702  // duplicates. Since the set of efficiency TPs are required to
703  // be a subset of the set of fake TPs, for efficiency
704  // histograms it doesn't matter if the association contains
705  // associations of TPs not in the set of efficiency TPs.
706  simRecCollL = theAssociator->associateSimToReco(trackRefs, tPCfake);
707  simRecCollP = &simRecCollL;
708  }
709  else{
710  // We need to filter the associations of the current track
711  // collection only from SimToReco collection, otherwise the
712  // SimToReco histograms get false entries. The filtering must
713  // be done separately for each track collection.
714  simRecCollL = associationMapFilterValues(*simRecCollPFull, trackCollection);
715  simRecCollP = &simRecCollL;
716  }
717 
718  reco::RecoToSimCollection const & recSimColl = *recSimCollP;
719  reco::SimToRecoCollection const & simRecColl = *simRecCollP;
720 
721  // read MVA collections
725  for(const auto& tokenTpl: mvaQualityCollectionTokens_[www]) {
726  event.getByToken(std::get<0>(tokenTpl), hmva);
727  event.getByToken(std::get<1>(tokenTpl), hqual);
728 
729  mvaCollections.push_back(hmva.product());
730  qualityMaskCollections.push_back(hqual.product());
731  if(mvaCollections.back()->size() != trackCollection.size()) {
732  throw cms::Exception("Configuration") << "Inconsistency in track collection and MVA sizes. Track collection " << www << " has " << trackCollection.size() << " tracks, whereas the MVA " << (mvaCollections.size()-1) << " for it has " << mvaCollections.back()->size() << " entries. Double-check your configuration.";
733  }
734  if(qualityMaskCollections.back()->size() != trackCollection.size()) {
735  throw cms::Exception("Configuration") << "Inconsistency in track collection and quality mask sizes. Track collection " << www << " has " << trackCollection.size() << " tracks, whereas the quality mask " << (qualityMaskCollections.size()-1) << " for it has " << qualityMaskCollections.back()->size() << " entries. Double-check your configuration.";
736  }
737  }
738  }
739 
740  // ########################################################
741  // fill simulation histograms (LOOP OVER TRACKINGPARTICLES)
742  // ########################################################
743 
744  //compute number of tracks per eta interval
745  //
746  LogTrace("TrackValidator") << "\n# of TrackingParticles: " << tPCeff.size() << "\n";
747  int ats(0); //This counter counts the number of simTracks that are "associated" to recoTracks
748  int st(0); //This counter counts the number of simulated tracks passing the MTV selection (i.e. tpSelector(tp) )
749 
750  //loop over already-selected TPs for tracking efficiency
751  for(size_t i=0; i<selected_tPCeff.size(); ++i) {
752  size_t iTP = selected_tPCeff[i];
753  const TrackingParticleRef& tpr = tPCeff[iTP];
754  const TrackingParticle& tp = *tpr;
755 
756  auto const& momVert = momVert_tPCeff[i];
757  TrackingParticle::Vector momentumTP;
758  TrackingParticle::Point vertexTP;
759 
760  double dxySim(0);
761  double dzSim(0);
762  double dxyPVSim = 0;
763  double dzPVSim = 0;
764  double dR=dR_tPCeff[iTP];
765 
766  //---------- THIS PART HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
767  //If the TrackingParticle is collison like, get the momentum and vertex at production state
769  {
770  momentumTP = tp.momentum();
771  vertexTP = tp.vertex();
772  //Calcualte the impact parameters w.r.t. PCA
773  const TrackingParticle::Vector& momentum = std::get<TrackingParticle::Vector>(momVert);
774  const TrackingParticle::Point& vertex = std::get<TrackingParticle::Point>(momVert);
775  dxySim = TrackingParticleIP::dxy(vertex, momentum, bs.position());
776  dzSim = TrackingParticleIP::dz(vertex, momentum, bs.position());
777 
778  if(theSimPVPosition) {
779  dxyPVSim = TrackingParticleIP::dxy(vertex, momentum, *theSimPVPosition);
780  dzPVSim = TrackingParticleIP::dz(vertex, momentum, *theSimPVPosition);
781  }
782  }
783  //If the TrackingParticle is comics, get the momentum and vertex at PCA
784  else
785  {
786  momentumTP = std::get<TrackingParticle::Vector>(momVert);
787  vertexTP = std::get<TrackingParticle::Point>(momVert);
788  dxySim = TrackingParticleIP::dxy(vertexTP, momentumTP, bs.position());
789  dzSim = TrackingParticleIP::dz(vertexTP, momentumTP, bs.position());
790 
791  // Do dxy and dz vs. PV make any sense for cosmics? I guess not
792  }
793  //---------- THE PART ABOVE HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
794 
795  // in the coming lines, histos are filled using as input
796  // - momentumTP
797  // - vertexTP
798  // - dxySim
799  // - dzSim
800  if(!doSimTrackPlots_)
801  continue;
802 
803  // ##############################################
804  // fill RecoAssociated SimTracks' histograms
805  // ##############################################
806  const reco::Track *matchedTrackPointer = nullptr;
807  const reco::Track *matchedSecondTrackPointer = nullptr;
808  unsigned int selectsLoose = mvaCollections.size();
809  unsigned int selectsHP = mvaCollections.size();
810  if(simRecColl.find(tpr) != simRecColl.end()){
811  auto const & rt = simRecColl[tpr];
812  if (!rt.empty()) {
813  ats++; //This counter counts the number of simTracks that have a recoTrack associated
814  // isRecoMatched = true; // UNUSED
815  matchedTrackPointer = rt.begin()->first.get();
816  if(rt.size() >= 2) {
817  matchedSecondTrackPointer = (rt.begin()+1)->first.get();
818  }
819  LogTrace("TrackValidator") << "TrackingParticle #" << st
820  << " with pt=" << sqrt(momentumTP.perp2())
821  << " associated with quality:" << rt.begin()->second <<"\n";
822 
823  if(doMVAPlots_) {
824  // for each MVA we need to take the value of the track
825  // with largest MVA value (for the cumulative histograms)
826  //
827  // also identify the first MVA that possibly selects any
828  // track matched to this TrackingParticle, separately
829  // for loose and highPurity qualities
830  for(size_t imva=0; imva<mvaCollections.size(); ++imva) {
831  const auto& mva = *(mvaCollections[imva]);
832  const auto& qual = *(qualityMaskCollections[imva]);
833 
834  auto iMatch = rt.begin();
835  float maxMva = mva[iMatch->first.key()];
836  for(; iMatch!=rt.end(); ++iMatch) {
837  auto itrk = iMatch->first.key();
838  maxMva = std::max(maxMva, mva[itrk]);
839 
840  if(selectsLoose >= imva && trackSelected(qual[itrk], reco::TrackBase::loose))
841  selectsLoose = imva;
842  if(selectsHP >= imva && trackSelected(qual[itrk], reco::TrackBase::highPurity))
843  selectsHP = imva;
844  }
845  mvaValues.push_back(maxMva);
846  }
847  }
848  }
849  }else{
850  LogTrace("TrackValidator")
851  << "TrackingParticle #" << st
852  << " with pt,eta,phi: "
853  << sqrt(momentumTP.perp2()) << " , "
854  << momentumTP.eta() << " , "
855  << momentumTP.phi() << " , "
856  << " NOT associated to any reco::Track" << "\n";
857  }
858 
859 
860 
861 
862  int nSimHits = tp.numberOfTrackerHits();
863  int nSimLayers = nLayers_tPCeff[tpr];
864  int nSimPixelLayers = nPixelLayers_tPCeff[tpr];
865  int nSimStripMonoAndStereoLayers = nStripMonoAndStereoLayers_tPCeff[tpr];
866  histoProducerAlgo_->fill_recoAssociated_simTrack_histos(histograms.histoProducerAlgo,w,tp,momentumTP,vertexTP,dxySim,dzSim,dxyPVSim,dzPVSim,nSimHits,nSimLayers,nSimPixelLayers,nSimStripMonoAndStereoLayers,matchedTrackPointer,puinfo.getPU_NumInteractions(), dR, thePVposition, theSimPVPosition, bs.position(), mvaValues, selectsLoose, selectsHP);
867  mvaValues.clear();
868 
869  if(matchedTrackPointer && matchedSecondTrackPointer) {
870  histoProducerAlgo_->fill_duplicate_histos(histograms.histoProducerAlgo,w, *matchedTrackPointer, *matchedSecondTrackPointer);
871  }
872 
873  if(doSummaryPlots_) {
874  if(dRtpSelector(tp)) {
875  histograms.h_simul_coll[ww].fill(www);
876  if (matchedTrackPointer) {
877  histograms.h_assoc_coll[ww].fill(www);
878  }
879  }
880  }
881 
882 
883 
884 
885  } // End for (TrackingParticleCollection::size_type i=0; i<tPCeff.size(); i++){
886 
887  // ##############################################
888  // fill recoTracks histograms (LOOP OVER TRACKS)
889  // ##############################################
890  if(!doRecoTrackPlots_)
891  continue;
892  LogTrace("TrackValidator") << "\n# of reco::Tracks with "
893  << label[www].process()<<":"
894  << label[www].label()<<":"
895  << label[www].instance()
896  << ": " << trackCollection.size() << "\n";
897 
898  int sat(0); //This counter counts the number of recoTracks that are associated to SimTracks from Signal only
899  int at(0); //This counter counts the number of recoTracks that are associated to SimTracks
900  int rT(0); //This counter counts the number of recoTracks in general
901  int seed_fit_failed = 0;
902  size_t n_selTrack_dr = 0;
903 
904  //calculate dR for tracks
905  const edm::View<Track> *trackCollectionDr = &trackCollection;
907  trackCollectionDr = trackCollectionForDrCalculation.product();
908  }
909  declareDynArray(float, trackCollection.size(), dR_trk);
910  trackDR(trackCollection, *trackCollectionDr, dR_trk);
911 
912  for(View<Track>::size_type i=0; i<trackCollection.size(); ++i){
913  auto track = trackCollection.refAt(i);
914  rT++;
915  if(trackFromSeedFitFailed(*track)) ++seed_fit_failed;
916  if((*dRTrackSelector)(*track, bs.position())) ++n_selTrack_dr;
917 
918  bool isSigSimMatched(false);
919  bool isSimMatched(false);
920  bool isChargeMatched(true);
921  int numAssocRecoTracks = 0;
922  int nSimHits = 0;
923  double sharedFraction = 0.;
924 
925  auto tpFound = recSimColl.find(track);
926  isSimMatched = tpFound != recSimColl.end();
927  if (isSimMatched) {
928  const auto& tp = tpFound->val;
929  nSimHits = tp[0].first->numberOfTrackerHits();
930  sharedFraction = tp[0].second;
931  if (tp[0].first->charge() != track->charge()) isChargeMatched = false;
932  if(simRecColl.find(tp[0].first) != simRecColl.end()) numAssocRecoTracks = simRecColl[tp[0].first].size();
933  at++;
934  for (unsigned int tp_ite=0;tp_ite<tp.size();++tp_ite){
935  TrackingParticle trackpart = *(tp[tp_ite].first);
936  if ((trackpart.eventId().event() == 0) && (trackpart.eventId().bunchCrossing() == 0)){
937  isSigSimMatched = true;
938  sat++;
939  break;
940  }
941  }
942  LogTrace("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
943  << " associated with quality:" << tp.begin()->second <<"\n";
944  } else {
945  LogTrace("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
946  << " NOT associated to any TrackingParticle" << "\n";
947  }
948 
949  // set MVA values for this track
950  // take also the indices of first MVAs to select by loose and
951  // HP quality
952  unsigned int selectsLoose = mvaCollections.size();
953  unsigned int selectsHP = mvaCollections.size();
954  if(doMVAPlots_) {
955  for(size_t imva=0; imva<mvaCollections.size(); ++imva) {
956  const auto& mva = *(mvaCollections[imva]);
957  const auto& qual = *(qualityMaskCollections[imva]);
958  mvaValues.push_back(mva[i]);
959 
960  if(selectsLoose >= imva && trackSelected(qual[i], reco::TrackBase::loose))
961  selectsLoose = imva;
962  if(selectsHP >= imva && trackSelected(qual[i], reco::TrackBase::highPurity))
963  selectsHP = imva;
964  }
965  }
966 
967  double dR=dR_trk[i];
968  histoProducerAlgo_->fill_generic_recoTrack_histos(histograms.histoProducerAlgo,w,*track, ttopo, bs.position(), thePVposition, theSimPVPosition, isSimMatched,isSigSimMatched, isChargeMatched, numAssocRecoTracks, puinfo.getPU_NumInteractions(), nSimHits, sharedFraction, dR, mvaValues, selectsLoose, selectsHP);
969  mvaValues.clear();
970 
971  if(doSummaryPlots_) {
972  histograms.h_reco_coll[ww].fill(www);
973  if(isSimMatched) {
974  histograms.h_assoc2_coll[ww].fill(www);
975  if(numAssocRecoTracks>1) {
976  histograms.h_looper_coll[ww].fill(www);
977  }
978  if(!isSigSimMatched) {
979  histograms.h_pileup_coll[ww].fill(www);
980  }
981  }
982  }
983 
984  // dE/dx
985  if (dodEdxPlots_) histoProducerAlgo_->fill_dedx_recoTrack_histos(histograms.histoProducerAlgo,w,track, v_dEdx);
986 
987 
988  //Fill other histos
989  if (!isSimMatched) continue;
990 
991  histoProducerAlgo_->fill_simAssociated_recoTrack_histos(histograms.histoProducerAlgo,w,*track);
992 
993  /* TO BE FIXED LATER
994  if (associators[ww]=="trackAssociatorByChi2"){
995  //association chi2
996  double assocChi2 = -tp.begin()->second;//in association map is stored -chi2
997  h_assochi2[www]->Fill(assocChi2);
998  h_assochi2_prob[www]->Fill(TMath::Prob((assocChi2)*5,5));
999  }
1000  else if (associators[ww]=="quickTrackAssociatorByHits"){
1001  double fraction = tp.begin()->second;
1002  h_assocFraction[www]->Fill(fraction);
1003  h_assocSharedHit[www]->Fill(fraction*track->numberOfValidHits());
1004  }
1005  */
1006 
1007 
1008  if(doResolutionPlots_[www]) {
1009  //Get tracking particle parameters at point of closest approach to the beamline
1010  TrackingParticleRef tpr = tpFound->val.begin()->first;
1011  TrackingParticle::Vector momentumTP = parametersDefinerTP->momentum(event,setup,tpr);
1012  TrackingParticle::Point vertexTP = parametersDefinerTP->vertex(event,setup,tpr);
1013  int chargeTP = tpr->charge();
1014 
1015  histoProducerAlgo_->fill_ResoAndPull_recoTrack_histos(histograms.histoProducerAlgo,w,momentumTP,vertexTP,chargeTP,
1016  *track,bs.position());
1017  }
1018 
1019 
1020  //TO BE FIXED
1021  //std::vector<PSimHit> simhits=tpr.get()->trackPSimHit(DetId::Tracker);
1022  //nrecHit_vs_nsimHit_rec2sim[w]->Fill(track->numberOfValidHits(), (int)(simhits.end()-simhits.begin() ));
1023 
1024  } // End of for(View<Track>::size_type i=0; i<trackCollection.size(); ++i){
1025  mvaCollections.clear();
1026  qualityMaskCollections.clear();
1027 
1028  histoProducerAlgo_->fill_trackBased_histos(histograms.histoProducerAlgo,w,at,rT, n_selTrack_dr, n_selTP_dr);
1029  // Fill seed-specific histograms
1030  if(doSeedPlots_) {
1031  histoProducerAlgo_->fill_seed_histos(histograms.histoProducerAlgo,www, seed_fit_failed, trackCollection.size());
1032  }
1033 
1034 
1035  LogTrace("TrackValidator") << "Collection " << www << "\n"
1036  << "Total Simulated (selected): " << n_selTP_dr << "\n"
1037  << "Total Reconstructed (selected): " << n_selTrack_dr << "\n"
1038  << "Total Reconstructed: " << rT << "\n"
1039  << "Total Associated (recoToSim): " << at << "\n"
1040  << "Total Fakes: " << rT-at << "\n";
1041  } // End of for (unsigned int www=0;www<label.size();www++){
1042  } //END of for (unsigned int ww=0;ww<associators.size();ww++){
1043 
1044 }
#define LogDebug(id)
size
Write out results.
void trackDR(const edm::View< reco::Track > &trackCollection, const edm::View< reco::Track > &trackCollectionDr, DynArray< float > &dR_trk) const
T getParameter(std::string const &) const
unsigned int size_type
Definition: View.h:90
T getUntrackedParameter(std::string const &, T const &) const
std::vector< ConcurrentMonitorElement > h_pileup_coll
edm::EDGetTokenT< SimHitTPAssociationProducer::SimHitTPAssociationList > _simHitTpMapTag
int event() const
get the contents of the subdetector field (should be protected?)
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx2Tag
bool trackFromSeedFitFailed(const reco::Track &track)
std::vector< edm::EDGetTokenT< reco::SimToRecoCollection > > associatormapStRs
const double w
Definition: UKUtility.cc:23
CosmicTrackingParticleSelector cosmictpSelector
Vector momentum() const
spatial momentum vector
const_iterator end() const
last iterator over the map (read only)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
const_iterator find(const key_type &k) const
find element with specified reference key
std::vector< std::vector< std::tuple< edm::EDGetTokenT< MVACollection >, edm::EDGetTokenT< QualityMaskCollection > > > > mvaQualityCollectionTokens_
std::vector< ConcurrentMonitorElement > h_simul_coll
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNStripStereoLayersToken_
def replace(string, replacements)
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNLayersToken_
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
size_type size() const
T_AssociationMap associationMapFilterValues(const T_AssociationMap &map, const T_RefVector &valueRefs)
TrackingParticleSelector dRtpSelector
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
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > > associatorTokens
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
const bool doPlotsOnlyForTruePV_
edm::EDGetTokenT< reco::BeamSpot > bsSrc
const reco::Vertex::Point * getRecoPVPosition(const edm::Event &event, const edm::Handle< TrackingVertexCollection > &htv) const
std::vector< ConcurrentMonitorElement > h_looper_coll
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void tpParametersAndSelection(const Histograms &histograms, const TrackingParticleRefVector &tPCeff, const ParametersDefinerForTP &parametersDefinerTP, const edm::Event &event, const edm::EventSetup &setup, const reco::BeamSpot &bs, std::vector< std::tuple< TrackingParticle::Vector, TrackingParticle::Point > > &momVert_tPCeff, std::vector< size_t > &selected_tPCeff) const
std::vector< ConcurrentMonitorElement > h_reco_coll
~MultiTrackValidator() override
Destructor.
RefToBase< value_type > refAt(size_type i) const
math::XYZPointD Point
point in the space
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNPixelLayersToken_
math::XYZTLorentzVectorD LorentzVector
std::vector< edm::InputTag > associators
SingleObjectSelector< TrackingParticleCollection,::TrackingParticleSelector > TrackingParticleSelector
std::vector< edm::EDGetTokenT< edm::View< reco::Track > > > labelToken
bool has(unsigned int index) const
Check if an element (=index) is in the set.
Definition: IndexSet.h:55
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx1Tag
ProductID id() const
Accessor for product ID.
Definition: RefVector.h:122
std::vector< edm::EDGetTokenT< reco::RecoToSimCollection > > associatormapRtSs
auto dz(const T_Vertex &vertex, const T_Momentum &momentum, const T_Point &point)
T sqrt(T t)
Definition: SSEVec.h:18
size_t tpDR(const TrackingParticleRefVector &tPCeff, const std::vector< size_t > &selected_tPCeff, DynArray< float > &dR_tPCeff) const
int bunchCrossing() const
get the detector field from this detid
TrackingParticleSelector tpSelector
std::vector< bool > doResolutionPlots_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< ConcurrentMonitorElement > h_assoc_coll
ConcurrentMonitorElement book1D(Args &&...args)
Definition: DQMStore.h:223
virtual TrackingParticle::Vector momentum(const edm::Event &iEvent, const edm::EventSetup &iSetup, const Charge ch, const Point &vtx, const LorentzVector &lv) const
edm::EDGetTokenT< edm::View< reco::Track > > labelTokenForDrCalculation
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > label_pileupinfo
std::vector< edm::InputTag > label
void initEvent(edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList > simHitsTPAssocToSet) const
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
char const * module
Definition: ProductLabels.h:5
virtual std::unique_ptr< ParametersDefinerForTP > clone() const
MultiTrackValidator(const edm::ParameterSet &pset)
Constructor.
#define LogTrace(id)
ObjectSelector< CosmicTrackingParticleSelector > CosmicTrackingParticleSelector
edm::EDGetTokenT< TrackingParticleRefVector > label_tp_fake_refvector
std::unique_ptr< MTVHistoProducerAlgoForTracker > histoProducerAlgo_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< TrackingVertexCollection > label_tv
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
edm::EDGetTokenT< reco::VertexToTrackingVertexAssociator > vertexAssociatorToken_
std::unique_ptr< RecoTrackSelectorBase > dRTrackSelector
const bool ignoremissingtkcollection_
MTVHistoProducerAlgoForTrackerHistograms histoProducerAlgo
const int getPU_NumInteractions() const
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
const bool doPVAssociationPlots_
std::string const & label() const
Definition: InputTag.h:36
EncodedEventId eventId() const
Signal source, crossing number.
Point vertex() const
Parent vertex position.
edm::EDGetTokenT< edm::View< reco::Vertex > > recoVertexToken_
std::vector< edm::EDGetTokenT< std::vector< PSimHit > > > simHitTokens_
std::string const & process() const
Definition: InputTag.h:40
void dqmAnalyze(const edm::Event &, const edm::EventSetup &, const Histograms &) const override
Method called once per event.
const EncodedEventId & eventId() const
const bool calculateDrSingleCollection_
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
auto dxy(const T_Vertex &vertex, const T_Momentum &momentum, const T_Point &point)
void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, Histograms &) const override
Method called to book the DQM histograms.
T get() const
Definition: EventSetup.h:63
void push_back(const RefToBase< T > &)
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
edm::EDGetTokenT< TrackingParticleCollection > label_tp_fake
Monte Carlo truth information used for tracking validation.
void insert(unsigned int index)
Insert an element (=index) to the set.
Definition: IndexSet.h:48
bool isUninitialized() const
Definition: EDGetToken.h:73
const Point & position() const
position
Definition: BeamSpot.h:62
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
reco::SimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
math::XYZVectorD Vector
point in the space
edm::EDGetTokenT< TrackingParticleRefVector > label_tp_effic_refvector
reco::VertexRecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Vertex > > &vCH, const edm::Handle< TrackingVertexCollection > &tVCH) const
compare reco to sim the handle of reco::Vertex and TrackingVertex collections
doResolutionPlotsForLabels
do resolution plots only for these labels (or all if empty)
static std::unique_ptr< RecoTrackSelectorBase > makeRecoTrackSelectorFromTPSelectorParameters(const edm::ParameterSet &pset)
void reserve(unsigned int size)
Reserve memory for the set.
Definition: IndexSet.h:35
std::vector< ConcurrentMonitorElement > h_assoc2_coll
edm::EDGetTokenT< TrackingParticleCollection > label_tp_effic
std::string const & instance() const
Definition: InputTag.h:37
edm::Ref< TrackingParticleCollection > TrackingParticleRef
Definition: event.py:1
Definition: Run.h:44
virtual TrackingParticle::Point vertex(const edm::Event &iEvent, const edm::EventSetup &iSetup, const Charge ch, const Point &vtx, const LorentzVector &lv) const
const TrackingVertex::LorentzVector * getSimPVPosition(const edm::Handle< TrackingVertexCollection > &htv) const
const bool parametersDefinerIsCosmic_