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