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  return me;
224  };
225 
226  //Booking histograms concerning with simulated tracks
227  if(doSimPlots_) {
228  ibook.cd();
229  ibook.setCurrentFolder(dirName_ + "simulation");
230 
231  histoProducerAlgo_->bookSimHistos(ibook, histograms.histoProducerAlgo);
232 
233  ibook.cd();
234  ibook.setCurrentFolder(dirName_);
235  }
236 
237  for (unsigned int ww=0;ww<associators.size();ww++){
238  ibook.cd();
239  // FIXME: these need to be moved to a subdirectory whose name depends on the associator
240  ibook.setCurrentFolder(dirName_);
241 
242  if(doSummaryPlots_) {
243  if(doSimTrackPlots_) {
244  histograms.h_assoc_coll.push_back(binLabels( ibook.book1D("num_assoc(simToReco)_coll", "N of associated (simToReco) tracks vs track collection", nintColl, minColl, maxColl) ));
245  histograms.h_simul_coll.push_back(binLabels( ibook.book1D("num_simul_coll", "N of simulated tracks vs track collection", nintColl, minColl, maxColl) ));
246  }
247  if(doRecoTrackPlots_) {
248  histograms.h_reco_coll.push_back(binLabels( ibook.book1D("num_reco_coll", "N of reco track vs track collection", nintColl, minColl, maxColl) ));
249  histograms.h_assoc2_coll.push_back(binLabels( ibook.book1D("num_assoc(recoToSim)_coll", "N of associated (recoToSim) tracks vs track collection", nintColl, minColl, maxColl) ));
250  histograms.h_looper_coll.push_back(binLabels( ibook.book1D("num_duplicate_coll", "N of associated (recoToSim) looper tracks vs track collection", nintColl, minColl, maxColl) ));
251  histograms.h_pileup_coll.push_back(binLabels( ibook.book1D("num_pileup_coll", "N of associated (recoToSim) pileup tracks vs track collection", nintColl, minColl, maxColl) ));
252  }
253  }
254 
255  for (unsigned int www=0;www<label.size();www++){
256  ibook.cd();
257  InputTag algo = label[www];
258  string dirName=dirName_;
259  if (algo.process()!="")
260  dirName+=algo.process()+"_";
261  if(algo.label()!="")
262  dirName+=algo.label()+"_";
263  if(algo.instance()!="")
264  dirName+=algo.instance()+"_";
265  if (dirName.find("Tracks")<dirName.length()){
266  dirName.replace(dirName.find("Tracks"),6,"");
267  }
268  string assoc= associators[ww].label();
269  if (assoc.find("Track")<assoc.length()){
270  assoc.replace(assoc.find("Track"),5,"");
271  }
272  dirName+=assoc;
273  std::replace(dirName.begin(), dirName.end(), ':', '_');
274 
275  ibook.setCurrentFolder(dirName);
276 
277  const bool doResolutionPlots = doResolutionPlots_[www];
278 
279  if(doSimTrackPlots_) {
280  histoProducerAlgo_->bookSimTrackHistos(ibook, histograms.histoProducerAlgo, doResolutionPlots);
281  if(doPVAssociationPlots_) histoProducerAlgo_->bookSimTrackPVAssociationHistos(ibook, histograms.histoProducerAlgo);
282  }
283 
284  //Booking histograms concerning with reconstructed tracks
285  if(doRecoTrackPlots_) {
286  histoProducerAlgo_->bookRecoHistos(ibook, histograms.histoProducerAlgo, doResolutionPlots);
287  if (dodEdxPlots_) histoProducerAlgo_->bookRecodEdxHistos(ibook, histograms.histoProducerAlgo);
288  if (doPVAssociationPlots_) histoProducerAlgo_->bookRecoPVAssociationHistos(ibook, histograms.histoProducerAlgo);
289  if (doMVAPlots_) histoProducerAlgo_->bookMVAHistos(ibook, histograms.histoProducerAlgo, mvaQualityCollectionTokens_[www].size());
290  }
291 
292  if(doSeedPlots_) {
293  histoProducerAlgo_->bookSeedHistos(ibook, histograms.histoProducerAlgo);
294  }
295  }//end loop www
296  }// end loop ww
297 }
298 
299 namespace {
300  void ensureEffIsSubsetOfFake(const TrackingParticleRefVector& eff, const TrackingParticleRefVector& fake) {
301  // If efficiency RefVector is empty, don't check the product ids
302  // as it will be 0:0 for empty. This covers also the case where
303  // both are empty. The case of fake being empty and eff not is an
304  // error.
305  if(eff.empty())
306  return;
307 
308  // First ensure product ids
309  if(eff.id() != fake.id()) {
310  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.";
311  }
312 
313  // Same technique as in associationMapFilterValues
314  edm::IndexSet fakeKeys;
315  fakeKeys.reserve(fake.size());
316  for(const auto& ref: fake) {
317  fakeKeys.insert(ref.key());
318  }
319 
320  for(const auto& ref: eff) {
321  if(!fakeKeys.has(ref.key())) {
322  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.";
323  }
324  }
325  }
326 }
327 
329  for(const auto& simV: *htv) {
330  if(simV.eventId().bunchCrossing() != 0) continue; // remove OOTPU
331  if(simV.eventId().event() != 0) continue; // pick the PV of hard scatter
332  return &(simV.position());
333  }
334  return nullptr;
335 }
336 
339  event.getByToken(recoVertexToken_, hvertex);
340 
342  event.getByToken(vertexAssociatorToken_, hvassociator);
343 
344  auto v_r2s = hvassociator->associateRecoToSim(hvertex, htv);
345  auto pvPtr = hvertex->refAt(0);
346  if(pvPtr->isFake() || pvPtr->ndof() < 0) // skip junk vertices
347  return nullptr;
348 
349  auto pvFound = v_r2s.find(pvPtr);
350  if(pvFound == v_r2s.end())
351  return nullptr;
352 
353  for(const auto& vertexRefQuality: pvFound->val) {
354  const TrackingVertex& tv = *(vertexRefQuality.first);
355  if(tv.eventId().event() == 0 && tv.eventId().bunchCrossing() == 0) {
356  return &(pvPtr->position());
357  }
358  }
359 
360  return nullptr;
361 }
362 
364  const TrackingParticleRefVector& tPCeff,
365  const ParametersDefinerForTP& parametersDefinerTP,
366  const edm::Event& event, const edm::EventSetup& setup,
367  const reco::BeamSpot& bs,
368  std::vector<std::tuple<TrackingParticle::Vector, TrackingParticle::Point> >& momVert_tPCeff,
369  std::vector<size_t>& selected_tPCeff) const {
370  selected_tPCeff.reserve(tPCeff.size());
371  momVert_tPCeff.reserve(tPCeff.size());
372  int nIntimeTPs = 0;
374  for(size_t j=0; j<tPCeff.size(); ++j) {
375  const TrackingParticleRef& tpr = tPCeff[j];
376 
377  TrackingParticle::Vector momentum = parametersDefinerTP.momentum(event,setup,tpr);
378  TrackingParticle::Point vertex = parametersDefinerTP.vertex(event,setup,tpr);
379  if(doSimPlots_) {
380  histoProducerAlgo_->fill_generic_simTrack_histos(histograms.histoProducerAlgo, momentum, vertex, tpr->eventId().bunchCrossing());
381  }
382  if(tpr->eventId().bunchCrossing() == 0)
383  ++nIntimeTPs;
384 
385  if(cosmictpSelector(tpr,&bs,event,setup)) {
386  selected_tPCeff.push_back(j);
387  momVert_tPCeff.emplace_back(momentum, vertex);
388  }
389  }
390  }
391  else {
392  size_t j=0;
393  for(auto const& tpr: tPCeff) {
394  const TrackingParticle& tp = *tpr;
395 
396  // TODO: do we want to fill these from all TPs that include IT
397  // and OOT (as below), or limit to IT+OOT TPs passing tpSelector
398  // (as it was before)? The latter would require another instance
399  // of tpSelector with intimeOnly=False.
400  if(doSimPlots_) {
401  histoProducerAlgo_->fill_generic_simTrack_histos(histograms.histoProducerAlgo, tp.momentum(), tp.vertex(), tp.eventId().bunchCrossing());
402  }
403  if(tp.eventId().bunchCrossing() == 0)
404  ++nIntimeTPs;
405 
406  if(tpSelector(tp)) {
407  selected_tPCeff.push_back(j);
408  TrackingParticle::Vector momentum = parametersDefinerTP.momentum(event,setup,tpr);
409  TrackingParticle::Point vertex = parametersDefinerTP.vertex(event,setup,tpr);
410  momVert_tPCeff.emplace_back(momentum, vertex);
411  }
412  ++j;
413  }
414  }
415  if(doSimPlots_) {
416  histoProducerAlgo_->fill_simTrackBased_histos(histograms.histoProducerAlgo, nIntimeTPs);
417  }
418 }
419 
420 
422  const std::vector<size_t>& selected_tPCeff,
423  DynArray<float>& dR_tPCeff) const {
424  float etaL[tPCeff.size()], phiL[tPCeff.size()];
425  size_t n_selTP_dr = 0;
426  for(size_t iTP: selected_tPCeff) {
427  //calculare dR wrt inclusive collection (also with PU, low pT, displaced)
428  auto const& tp2 = *(tPCeff[iTP]);
429  auto && p = tp2.momentum();
430  etaL[iTP] = etaFromXYZ(p.x(),p.y(),p.z());
431  phiL[iTP] = atan2f(p.y(),p.x());
432  }
433  for(size_t iTP1: selected_tPCeff) {
434  auto const& tp = *(tPCeff[iTP1]);
436  if(dRtpSelector(tp)) {//only for those needed for efficiency!
437  ++n_selTP_dr;
438  float eta = etaL[iTP1];
439  float phi = phiL[iTP1];
440  for(size_t iTP2: selected_tPCeff) {
441  //calculare dR wrt inclusive collection (also with PU, low pT, displaced)
442  if (iTP1==iTP2) {continue;}
443  auto dR_tmp = reco::deltaR2(eta, phi, etaL[iTP2], phiL[iTP2]);
444  if (dR_tmp<dR) dR=dR_tmp;
445  } // ttp2 (iTP)
446  }
447  dR_tPCeff[iTP1] = std::sqrt(dR);
448  } // tp
449  return n_selTP_dr;
450 }
451 
453  int i=0;
454  float etaL[trackCollectionDr.size()];
455  float phiL[trackCollectionDr.size()];
456  bool validL[trackCollectionDr.size()];
457  for (auto const & track2 : trackCollectionDr) {
458  auto && p = track2.momentum();
459  etaL[i] = etaFromXYZ(p.x(),p.y(),p.z());
460  phiL[i] = atan2f(p.y(),p.x());
461  validL[i] = !trackFromSeedFitFailed(track2);
462  ++i;
463  }
464  for(View<reco::Track>::size_type i=0; i<trackCollection.size(); ++i){
465  auto const & track = trackCollection[i];
468  auto && p = track.momentum();
469  float eta = etaFromXYZ(p.x(),p.y(),p.z());
470  float phi = atan2f(p.y(),p.x());
471  for(View<reco::Track>::size_type j=0; j<trackCollectionDr.size(); ++j){
472  if(!validL[j]) continue;
473  auto dR_tmp = reco::deltaR2(eta, phi, etaL[j], phiL[j]);
474  if ( (dR_tmp<dR) & (dR_tmp>std::numeric_limits<float>::min())) dR=dR_tmp;
475  }
476  }
477  dR_trk[i] = std::sqrt(dR);
478  }
479 }
480 
481 
483  using namespace reco;
484 
485  LogDebug("TrackValidator") << "\n====================================================" << "\n"
486  << "Analyzing new event" << "\n"
487  << "====================================================\n" << "\n";
488 
489 
490  edm::ESHandle<ParametersDefinerForTP> parametersDefinerTPHandle;
491  setup.get<TrackAssociatorRecord>().get(parametersDefiner,parametersDefinerTPHandle);
492  //Since we modify the object, we must clone it
493  auto parametersDefinerTP = parametersDefinerTPHandle->clone();
494 
496  setup.get<TrackerTopologyRcd>().get(httopo);
497  const TrackerTopology& ttopo = *httopo;
498 
499  // FIXME: we really need to move to edm::View for reading the
500  // TrackingParticles... Unfortunately it has non-trivial
501  // consequences on the associator/association interfaces etc.
502  TrackingParticleRefVector tmpTPeff;
503  TrackingParticleRefVector tmpTPfake;
504  const TrackingParticleRefVector *tmpTPeffPtr = nullptr;
505  const TrackingParticleRefVector *tmpTPfakePtr = nullptr;
506 
508  edm::Handle<TrackingParticleRefVector> TPCollectionHeffRefVector;
509 
510  const bool tp_effic_refvector = label_tp_effic.isUninitialized();
511  if(!tp_effic_refvector) {
512  event.getByToken(label_tp_effic, TPCollectionHeff);
513  for(size_t i=0, size=TPCollectionHeff->size(); i<size; ++i) {
514  tmpTPeff.push_back(TrackingParticleRef(TPCollectionHeff, i));
515  }
516  tmpTPeffPtr = &tmpTPeff;
517  }
518  else {
519  event.getByToken(label_tp_effic_refvector, TPCollectionHeffRefVector);
520  tmpTPeffPtr = TPCollectionHeffRefVector.product();
521  }
523  edm::Handle<TrackingParticleCollection> TPCollectionHfake ;
524  event.getByToken(label_tp_fake,TPCollectionHfake);
525  for(size_t i=0, size=TPCollectionHfake->size(); i<size; ++i) {
526  tmpTPfake.push_back(TrackingParticleRef(TPCollectionHfake, i));
527  }
528  tmpTPfakePtr = &tmpTPfake;
529  }
530  else {
531  edm::Handle<TrackingParticleRefVector> TPCollectionHfakeRefVector;
532  event.getByToken(label_tp_fake_refvector, TPCollectionHfakeRefVector);
533  tmpTPfakePtr = TPCollectionHfakeRefVector.product();
534  }
535 
536  TrackingParticleRefVector const & tPCeff = *tmpTPeffPtr;
537  TrackingParticleRefVector const & tPCfake = *tmpTPfakePtr;
538 
539  ensureEffIsSubsetOfFake(tPCeff, tPCfake);
540 
543  //warning: make sure the TP collection used in the map is the same used in the MTV!
544  event.getByToken(_simHitTpMapTag,simHitsTPAssoc);
545  parametersDefinerTP->initEvent(simHitsTPAssoc);
546  cosmictpSelector.initEvent(simHitsTPAssoc);
547  }
548 
549  // Find the sim PV and tak its position
551  event.getByToken(label_tv, htv);
552  const TrackingVertex::LorentzVector *theSimPVPosition = getSimPVPosition(htv);
553  if(simPVMaxZ_ >= 0) {
554  if(!theSimPVPosition) return;
555  if(std::abs(theSimPVPosition->z()) > simPVMaxZ_) return;
556  }
557 
558  // Check, when necessary, if reco PV matches to sim PV
559  const reco::Vertex::Point *thePVposition = nullptr;
561  thePVposition = getRecoPVPosition(event, htv);
562  if(doPlotsOnlyForTruePV_ && !thePVposition)
563  return;
564 
565  // Rest of the code assumes that if thePVposition is non-null, the
566  // PV-association histograms get filled. In above, the "nullness"
567  // is used to deliver the information if the reco PV is matched to
568  // the sim PV.
570  thePVposition = nullptr;
571  }
572 
573  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
574  event.getByToken(bsSrc,recoBeamSpotHandle);
575  reco::BeamSpot const & bs = *recoBeamSpotHandle;
576 
578  event.getByToken(label_pileupinfo,puinfoH);
579  PileupSummaryInfo puinfo;
580 
581  for (unsigned int puinfo_ite=0;puinfo_ite<(*puinfoH).size();++puinfo_ite){
582  if ((*puinfoH)[puinfo_ite].getBunchCrossing()==0){
583  puinfo=(*puinfoH)[puinfo_ite];
584  break;
585  }
586  }
587 
588  // Number of 3D layers for TPs
590  event.getByToken(tpNLayersToken_, tpNLayersH);
591  const auto& nLayers_tPCeff = *tpNLayersH;
592 
593  event.getByToken(tpNPixelLayersToken_, tpNLayersH);
594  const auto& nPixelLayers_tPCeff = *tpNLayersH;
595 
596  event.getByToken(tpNStripStereoLayersToken_, tpNLayersH);
597  const auto& nStripMonoAndStereoLayers_tPCeff = *tpNLayersH;
598 
599  // Precalculate TP selection (for efficiency), and momentum and vertex wrt PCA
600  //
601  // TODO: ParametersDefinerForTP ESProduct needs to be changed to
602  // EDProduct because of consumes.
603  //
604  // In principle, we could just precalculate the momentum and vertex
605  // wrt PCA for all TPs for once and put that to the event. To avoid
606  // repetitive calculations those should be calculated only once for
607  // each TP. That would imply that we should access TPs via Refs
608  // (i.e. View) in here, since, in general, the eff and fake TP
609  // collections can be different (and at least HI seems to use that
610  // feature). This would further imply that the
611  // RecoToSimCollection/SimToRecoCollection should be changed to use
612  // View<TP> instead of vector<TP>, and migrate everything.
613  //
614  // Or we could take only one input TP collection, and do another
615  // TP-selection to obtain the "fake" collection like we already do
616  // for "efficiency" TPs.
617  std::vector<size_t> selected_tPCeff;
618  std::vector<std::tuple<TrackingParticle::Vector, TrackingParticle::Point>> momVert_tPCeff;
619  tpParametersAndSelection(histograms, tPCeff, *parametersDefinerTP, event, setup, bs, momVert_tPCeff, selected_tPCeff);
620 
621  //calculate dR for TPs
622  declareDynArray(float, tPCeff.size(), dR_tPCeff);
623  size_t n_selTP_dr = tpDR(tPCeff, selected_tPCeff, dR_tPCeff);
624 
627  event.getByToken(labelTokenForDrCalculation, trackCollectionForDrCalculation);
628  }
629 
630  // dE/dx
631  // at some point this could be generalized, with a vector of tags and a corresponding vector of Handles
632  // I'm writing the interface such to take vectors of ValueMaps
633  std::vector<const edm::ValueMap<reco::DeDxData> *> v_dEdx;
634  if(dodEdxPlots_) {
637  event.getByToken(m_dEdx1Tag, dEdx1Handle);
638  event.getByToken(m_dEdx2Tag, dEdx2Handle);
639  v_dEdx.push_back(dEdx1Handle.product());
640  v_dEdx.push_back(dEdx2Handle.product());
641  }
642 
643  std::vector<const MVACollection *> mvaCollections;
644  std::vector<const QualityMaskCollection *> qualityMaskCollections;
645  std::vector<float> mvaValues;
646 
647  int w=0; //counter counting the number of sets of histograms
648  for (unsigned int ww=0;ww<associators.size();ww++){
649  // run value filtering of recoToSim map already here as it depends only on the association, not track collection
650  reco::SimToRecoCollection const * simRecCollPFull=nullptr;
651  reco::RecoToSimCollection const * recSimCollP=nullptr;
652  reco::RecoToSimCollection recSimCollL;
653  if(!useAssociators_) {
654  Handle<reco::SimToRecoCollection > simtorecoCollectionH;
655  event.getByToken(associatormapStRs[ww], simtorecoCollectionH);
656  simRecCollPFull = simtorecoCollectionH.product();
657 
658  Handle<reco::RecoToSimCollection > recotosimCollectionH;
659  event.getByToken(associatormapRtSs[ww],recotosimCollectionH);
660  recSimCollP = recotosimCollectionH.product();
661 
662  // We need to filter the associations of the fake-TrackingParticle
663  // collection only from RecoToSim collection, otherwise the
664  // RecoToSim histograms get false entries
665  recSimCollL = associationMapFilterValues(*recSimCollP, tPCfake);
666  recSimCollP = &recSimCollL;
667  }
668 
669  for (unsigned int www=0;www<label.size();www++, w++){ // need to increment w here, since there are many continues in the loop body
670  //
671  //get collections from the event
672  //
673  edm::Handle<View<Track> > trackCollectionHandle;
674  if(!event.getByToken(labelToken[www], trackCollectionHandle)&&ignoremissingtkcollection_)continue;
675  const edm::View<Track>& trackCollection = *trackCollectionHandle;
676 
677  reco::SimToRecoCollection const * simRecCollP=nullptr;
678  reco::SimToRecoCollection simRecCollL;
679 
680  //associate tracks
681  LogTrace("TrackValidator") << "Analyzing "
682  << label[www] << " with "
683  << associators[ww] <<"\n";
684  if(useAssociators_){
686  event.getByToken(associatorTokens[ww], theAssociator);
687 
688  // The associator interfaces really need to be fixed...
690  for(edm::View<Track>::size_type i=0; i<trackCollection.size(); ++i) {
691  trackRefs.push_back(trackCollection.refAt(i));
692  }
693 
694 
695  LogTrace("TrackValidator") << "Calling associateRecoToSim method" << "\n";
696  recSimCollL = theAssociator->associateRecoToSim(trackRefs, tPCfake);
697  recSimCollP = &recSimCollL;
698  LogTrace("TrackValidator") << "Calling associateSimToReco method" << "\n";
699  // It is necessary to do the association wrt. fake TPs,
700  // because this SimToReco association is used also for
701  // duplicates. Since the set of efficiency TPs are required to
702  // be a subset of the set of fake TPs, for efficiency
703  // histograms it doesn't matter if the association contains
704  // associations of TPs not in the set of efficiency TPs.
705  simRecCollL = theAssociator->associateSimToReco(trackRefs, tPCfake);
706  simRecCollP = &simRecCollL;
707  }
708  else{
709  // We need to filter the associations of the current track
710  // collection only from SimToReco collection, otherwise the
711  // SimToReco histograms get false entries. The filtering must
712  // be done separately for each track collection.
713  simRecCollL = associationMapFilterValues(*simRecCollPFull, trackCollection);
714  simRecCollP = &simRecCollL;
715  }
716 
717  reco::RecoToSimCollection const & recSimColl = *recSimCollP;
718  reco::SimToRecoCollection const & simRecColl = *simRecCollP;
719 
720  // read MVA collections
724  for(const auto& tokenTpl: mvaQualityCollectionTokens_[www]) {
725  event.getByToken(std::get<0>(tokenTpl), hmva);
726  event.getByToken(std::get<1>(tokenTpl), hqual);
727 
728  mvaCollections.push_back(hmva.product());
729  qualityMaskCollections.push_back(hqual.product());
730  if(mvaCollections.back()->size() != trackCollection.size()) {
731  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.";
732  }
733  if(qualityMaskCollections.back()->size() != trackCollection.size()) {
734  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.";
735  }
736  }
737  }
738 
739  // ########################################################
740  // fill simulation histograms (LOOP OVER TRACKINGPARTICLES)
741  // ########################################################
742 
743  //compute number of tracks per eta interval
744  //
745  LogTrace("TrackValidator") << "\n# of TrackingParticles: " << tPCeff.size() << "\n";
746  int ats(0); //This counter counts the number of simTracks that are "associated" to recoTracks
747  int st(0); //This counter counts the number of simulated tracks passing the MTV selection (i.e. tpSelector(tp) )
748 
749  //loop over already-selected TPs for tracking efficiency
750  for(size_t i=0; i<selected_tPCeff.size(); ++i) {
751  size_t iTP = selected_tPCeff[i];
752  const TrackingParticleRef& tpr = tPCeff[iTP];
753  const TrackingParticle& tp = *tpr;
754 
755  auto const& momVert = momVert_tPCeff[i];
756  TrackingParticle::Vector momentumTP;
757  TrackingParticle::Point vertexTP;
758 
759  double dxySim(0);
760  double dzSim(0);
761  double dxyPVSim = 0;
762  double dzPVSim = 0;
763  double dR=dR_tPCeff[iTP];
764 
765  //---------- THIS PART HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
766  //If the TrackingParticle is collison like, get the momentum and vertex at production state
768  {
769  momentumTP = tp.momentum();
770  vertexTP = tp.vertex();
771  //Calcualte the impact parameters w.r.t. PCA
772  const TrackingParticle::Vector& momentum = std::get<TrackingParticle::Vector>(momVert);
773  const TrackingParticle::Point& vertex = std::get<TrackingParticle::Point>(momVert);
774  dxySim = TrackingParticleIP::dxy(vertex, momentum, bs.position());
775  dzSim = TrackingParticleIP::dz(vertex, momentum, bs.position());
776 
777  if(theSimPVPosition) {
778  dxyPVSim = TrackingParticleIP::dxy(vertex, momentum, *theSimPVPosition);
779  dzPVSim = TrackingParticleIP::dz(vertex, momentum, *theSimPVPosition);
780  }
781  }
782  //If the TrackingParticle is comics, get the momentum and vertex at PCA
783  else
784  {
785  momentumTP = std::get<TrackingParticle::Vector>(momVert);
786  vertexTP = std::get<TrackingParticle::Point>(momVert);
787  dxySim = TrackingParticleIP::dxy(vertexTP, momentumTP, bs.position());
788  dzSim = TrackingParticleIP::dz(vertexTP, momentumTP, bs.position());
789 
790  // Do dxy and dz vs. PV make any sense for cosmics? I guess not
791  }
792  //---------- THE PART ABOVE HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
793 
794  // in the coming lines, histos are filled using as input
795  // - momentumTP
796  // - vertexTP
797  // - dxySim
798  // - dzSim
799  if(!doSimTrackPlots_)
800  continue;
801 
802  // ##############################################
803  // fill RecoAssociated SimTracks' histograms
804  // ##############################################
805  const reco::Track *matchedTrackPointer = nullptr;
806  const reco::Track *matchedSecondTrackPointer = nullptr;
807  unsigned int selectsLoose = mvaCollections.size();
808  unsigned int selectsHP = mvaCollections.size();
809  if(simRecColl.find(tpr) != simRecColl.end()){
810  auto const & rt = simRecColl[tpr];
811  if (!rt.empty()) {
812  ats++; //This counter counts the number of simTracks that have a recoTrack associated
813  // isRecoMatched = true; // UNUSED
814  matchedTrackPointer = rt.begin()->first.get();
815  if(rt.size() >= 2) {
816  matchedSecondTrackPointer = (rt.begin()+1)->first.get();
817  }
818  LogTrace("TrackValidator") << "TrackingParticle #" << st
819  << " with pt=" << sqrt(momentumTP.perp2())
820  << " associated with quality:" << rt.begin()->second <<"\n";
821 
822  if(doMVAPlots_) {
823  // for each MVA we need to take the value of the track
824  // with largest MVA value (for the cumulative histograms)
825  //
826  // also identify the first MVA that possibly selects any
827  // track matched to this TrackingParticle, separately
828  // for loose and highPurity qualities
829  for(size_t imva=0; imva<mvaCollections.size(); ++imva) {
830  const auto& mva = *(mvaCollections[imva]);
831  const auto& qual = *(qualityMaskCollections[imva]);
832 
833  auto iMatch = rt.begin();
834  float maxMva = mva[iMatch->first.key()];
835  for(; iMatch!=rt.end(); ++iMatch) {
836  auto itrk = iMatch->first.key();
837  maxMva = std::max(maxMva, mva[itrk]);
838 
839  if(selectsLoose >= imva && trackSelected(qual[itrk], reco::TrackBase::loose))
840  selectsLoose = imva;
841  if(selectsHP >= imva && trackSelected(qual[itrk], reco::TrackBase::highPurity))
842  selectsHP = imva;
843  }
844  mvaValues.push_back(maxMva);
845  }
846  }
847  }
848  }else{
849  LogTrace("TrackValidator")
850  << "TrackingParticle #" << st
851  << " with pt,eta,phi: "
852  << sqrt(momentumTP.perp2()) << " , "
853  << momentumTP.eta() << " , "
854  << momentumTP.phi() << " , "
855  << " NOT associated to any reco::Track" << "\n";
856  }
857 
858 
859 
860 
861  int nSimHits = tp.numberOfTrackerHits();
862  int nSimLayers = nLayers_tPCeff[tpr];
863  int nSimPixelLayers = nPixelLayers_tPCeff[tpr];
864  int nSimStripMonoAndStereoLayers = nStripMonoAndStereoLayers_tPCeff[tpr];
865  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);
866  mvaValues.clear();
867 
868  if(matchedTrackPointer && matchedSecondTrackPointer) {
869  histoProducerAlgo_->fill_duplicate_histos(histograms.histoProducerAlgo,w, *matchedTrackPointer, *matchedSecondTrackPointer);
870  }
871 
872  if(doSummaryPlots_) {
873  if(dRtpSelector(tp)) {
874  histograms.h_simul_coll[ww].fill(www);
875  if (matchedTrackPointer) {
876  histograms.h_assoc_coll[ww].fill(www);
877  }
878  }
879  }
880 
881 
882 
883 
884  } // End for (TrackingParticleCollection::size_type i=0; i<tPCeff.size(); i++){
885 
886  // ##############################################
887  // fill recoTracks histograms (LOOP OVER TRACKS)
888  // ##############################################
889  if(!doRecoTrackPlots_)
890  continue;
891  LogTrace("TrackValidator") << "\n# of reco::Tracks with "
892  << label[www].process()<<":"
893  << label[www].label()<<":"
894  << label[www].instance()
895  << ": " << trackCollection.size() << "\n";
896 
897  int sat(0); //This counter counts the number of recoTracks that are associated to SimTracks from Signal only
898  int at(0); //This counter counts the number of recoTracks that are associated to SimTracks
899  int rT(0); //This counter counts the number of recoTracks in general
900  int seed_fit_failed = 0;
901  size_t n_selTrack_dr = 0;
902 
903  //calculate dR for tracks
904  const edm::View<Track> *trackCollectionDr = &trackCollection;
906  trackCollectionDr = trackCollectionForDrCalculation.product();
907  }
908  declareDynArray(float, trackCollection.size(), dR_trk);
909  trackDR(trackCollection, *trackCollectionDr, dR_trk);
910 
911  for(View<Track>::size_type i=0; i<trackCollection.size(); ++i){
912  auto track = trackCollection.refAt(i);
913  rT++;
914  if(trackFromSeedFitFailed(*track)) ++seed_fit_failed;
915  if((*dRTrackSelector)(*track, bs.position())) ++n_selTrack_dr;
916 
917  bool isSigSimMatched(false);
918  bool isSimMatched(false);
919  bool isChargeMatched(true);
920  int numAssocRecoTracks = 0;
921  int nSimHits = 0;
922  double sharedFraction = 0.;
923 
924  auto tpFound = recSimColl.find(track);
925  isSimMatched = tpFound != recSimColl.end();
926  if (isSimMatched) {
927  const auto& tp = tpFound->val;
928  nSimHits = tp[0].first->numberOfTrackerHits();
929  sharedFraction = tp[0].second;
930  if (tp[0].first->charge() != track->charge()) isChargeMatched = false;
931  if(simRecColl.find(tp[0].first) != simRecColl.end()) numAssocRecoTracks = simRecColl[tp[0].first].size();
932  at++;
933  for (unsigned int tp_ite=0;tp_ite<tp.size();++tp_ite){
934  TrackingParticle trackpart = *(tp[tp_ite].first);
935  if ((trackpart.eventId().event() == 0) && (trackpart.eventId().bunchCrossing() == 0)){
936  isSigSimMatched = true;
937  sat++;
938  break;
939  }
940  }
941  LogTrace("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
942  << " associated with quality:" << tp.begin()->second <<"\n";
943  } else {
944  LogTrace("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
945  << " NOT associated to any TrackingParticle" << "\n";
946  }
947 
948  // set MVA values for this track
949  // take also the indices of first MVAs to select by loose and
950  // HP quality
951  unsigned int selectsLoose = mvaCollections.size();
952  unsigned int selectsHP = mvaCollections.size();
953  if(doMVAPlots_) {
954  for(size_t imva=0; imva<mvaCollections.size(); ++imva) {
955  const auto& mva = *(mvaCollections[imva]);
956  const auto& qual = *(qualityMaskCollections[imva]);
957  mvaValues.push_back(mva[i]);
958 
959  if(selectsLoose >= imva && trackSelected(qual[i], reco::TrackBase::loose))
960  selectsLoose = imva;
961  if(selectsHP >= imva && trackSelected(qual[i], reco::TrackBase::highPurity))
962  selectsHP = imva;
963  }
964  }
965 
966  double dR=dR_trk[i];
967  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);
968  mvaValues.clear();
969 
970  if(doSummaryPlots_) {
971  histograms.h_reco_coll[ww].fill(www);
972  if(isSimMatched) {
973  histograms.h_assoc2_coll[ww].fill(www);
974  if(numAssocRecoTracks>1) {
975  histograms.h_looper_coll[ww].fill(www);
976  }
977  if(!isSigSimMatched) {
978  histograms.h_pileup_coll[ww].fill(www);
979  }
980  }
981  }
982 
983  // dE/dx
984  if (dodEdxPlots_) histoProducerAlgo_->fill_dedx_recoTrack_histos(histograms.histoProducerAlgo,w,track, v_dEdx);
985 
986 
987  //Fill other histos
988  if (!isSimMatched) continue;
989 
990  histoProducerAlgo_->fill_simAssociated_recoTrack_histos(histograms.histoProducerAlgo,w,*track);
991 
992  /* TO BE FIXED LATER
993  if (associators[ww]=="trackAssociatorByChi2"){
994  //association chi2
995  double assocChi2 = -tp.begin()->second;//in association map is stored -chi2
996  h_assochi2[www]->Fill(assocChi2);
997  h_assochi2_prob[www]->Fill(TMath::Prob((assocChi2)*5,5));
998  }
999  else if (associators[ww]=="quickTrackAssociatorByHits"){
1000  double fraction = tp.begin()->second;
1001  h_assocFraction[www]->Fill(fraction);
1002  h_assocSharedHit[www]->Fill(fraction*track->numberOfValidHits());
1003  }
1004  */
1005 
1006 
1007  if(doResolutionPlots_[www]) {
1008  //Get tracking particle parameters at point of closest approach to the beamline
1009  TrackingParticleRef tpr = tpFound->val.begin()->first;
1010  TrackingParticle::Vector momentumTP = parametersDefinerTP->momentum(event,setup,tpr);
1011  TrackingParticle::Point vertexTP = parametersDefinerTP->vertex(event,setup,tpr);
1012  int chargeTP = tpr->charge();
1013 
1014  histoProducerAlgo_->fill_ResoAndPull_recoTrack_histos(histograms.histoProducerAlgo,w,momentumTP,vertexTP,chargeTP,
1015  *track,bs.position());
1016  }
1017 
1018 
1019  //TO BE FIXED
1020  //std::vector<PSimHit> simhits=tpr.get()->trackPSimHit(DetId::Tracker);
1021  //nrecHit_vs_nsimHit_rec2sim[w]->Fill(track->numberOfValidHits(), (int)(simhits.end()-simhits.begin() ));
1022 
1023  } // End of for(View<Track>::size_type i=0; i<trackCollection.size(); ++i){
1024  mvaCollections.clear();
1025  qualityMaskCollections.clear();
1026 
1027  histoProducerAlgo_->fill_trackBased_histos(histograms.histoProducerAlgo,w,at,rT, n_selTrack_dr, n_selTP_dr);
1028  // Fill seed-specific histograms
1029  if(doSeedPlots_) {
1030  histoProducerAlgo_->fill_seed_histos(histograms.histoProducerAlgo,www, seed_fit_failed, trackCollection.size());
1031  }
1032 
1033 
1034  LogTrace("TrackValidator") << "Collection " << www << "\n"
1035  << "Total Simulated (selected): " << n_selTP_dr << "\n"
1036  << "Total Reconstructed (selected): " << n_selTrack_dr << "\n"
1037  << "Total Reconstructed: " << rT << "\n"
1038  << "Total Associated (recoToSim): " << at << "\n"
1039  << "Total Fakes: " << rT-at << "\n";
1040  } // End of for (unsigned int www=0;www<label.size();www++){
1041  } //END of for (unsigned int ww=0;ww<associators.size();ww++){
1042 
1043 }
#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:519
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_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
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 T & get() const
Definition: EventSetup.h:59
const bool ignoremissingtkcollection_
MTVHistoProducerAlgoForTrackerHistograms histoProducerAlgo
const int getPU_NumInteractions() const
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
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.
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:43
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_