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