CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
TrackAssociatorByHitsImpl Class Reference

#include <TrackAssociatorByHitsImpl.h>

Inheritance diagram for TrackAssociatorByHitsImpl:
reco::TrackToTrackingParticleAssociatorBaseImpl

Public Types

typedef std::vector< SimHitTPPairSimHitTPAssociationList
 
typedef std::pair< TrackingParticleRef, TrackPSimHitRefSimHitTPPair
 
enum  SimToRecoDenomType { denomnone, denomsim, denomreco }
 

Public Member Functions

reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Reco To Sim with Collections. More...
 
reco::RecoToSimCollection associateRecoToSim (const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const override
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
reco::RecoToSimCollectionSeed associateRecoToSim (const edm::Handle< edm::View< TrajectorySeed > > &, const edm::Handle< TrackingParticleCollection > &) const override
 
reco::SimToRecoCollection associateSimToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Sim To Reco with Collections. More...
 
reco::SimToRecoCollection associateSimToReco (const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const override
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
reco::SimToRecoCollectionSeed associateSimToReco (const edm::Handle< edm::View< TrajectorySeed > > &, const edm::Handle< TrackingParticleCollection > &) const override
 
 TrackAssociatorByHitsImpl (edm::EDProductGetter const &productGetter, std::unique_ptr< TrackerHitAssociator > iAssociate, TrackerTopology const *iTopo, SimHitTPAssociationList const *iSimHitsTPAssoc, SimToRecoDenomType iSimToRecoDenominator, double iQuality_SimToReco, double iPurity_SimToReco, double iCut_RecoToSim, bool iUsePixels, bool iUseGrouped, bool iUseSplitting, bool ThreeHitTracksAreSpecial, bool AbsoluteNumberOfHits)
 
- Public Member Functions inherited from reco::TrackToTrackingParticleAssociatorBaseImpl
virtual reco::RecoToSimCollection associateRecoToSim (const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
 
virtual reco::RecoToSimCollectionSeed associateRecoToSim (const edm::Handle< edm::View< TrajectorySeed >> &, const edm::Handle< TrackingParticleCollection > &) const
 
virtual reco::RecoToSimCollectionTCandidate associateRecoToSim (const edm::Handle< TrackCandidateCollection > &, const edm::RefVector< TrackingParticleCollection > &) const
 
virtual reco::SimToRecoCollection associateSimToReco (const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
 
virtual reco::SimToRecoCollectionSeed associateSimToReco (const edm::Handle< edm::View< TrajectorySeed >> &, const edm::Handle< TrackingParticleCollection > &) const
 
virtual reco::SimToRecoCollectionTCandidate associateSimToReco (const edm::Handle< TrackCandidateCollection > &, const edm::RefVector< TrackingParticleCollection > &) const
 
 TrackToTrackingParticleAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~TrackToTrackingParticleAssociatorBaseImpl ()
 Destructor. More...
 

Private Member Functions

template<typename iter >
int getDoubleCount (iter, iter, TrackerHitAssociator *, TrackingParticle const &) const
 
const TrackingRecHitgetHitPtr (edm::OwnVector< TrackingRecHit >::const_iterator iter) const
 
const TrackingRecHitgetHitPtr (trackingRecHit_iterator iter) const
 
template<typename iter >
void getMatchedIds (std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, int &, iter, iter, TrackerHitAssociator *) const
 
int getShared (std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticle const &) const
 

Private Attributes

const bool AbsoluteNumberOfHits
 
std::unique_ptr< TrackerHitAssociatorassociate
 
const double cut_RecoToSim
 
edm::EDProductGetter const * productGetter_
 
const double purity_SimToReco
 
const double quality_SimToReco
 
SimHitTPAssociationList const * simHitsTPAssoc
 
SimToRecoDenomType SimToRecoDenominator
 
const bool ThreeHitTracksAreSpecial
 
TrackerTopology const * tTopo
 
const bool UseGrouped
 
const bool UsePixels
 
const bool UseSplitting
 

Detailed Description

Definition at line 23 of file TrackAssociatorByHitsImpl.h.

Member Typedef Documentation

◆ SimHitTPAssociationList

Definition at line 28 of file TrackAssociatorByHitsImpl.h.

◆ SimHitTPPair

Definition at line 27 of file TrackAssociatorByHitsImpl.h.

Member Enumeration Documentation

◆ SimToRecoDenomType

Constructor & Destructor Documentation

◆ TrackAssociatorByHitsImpl()

TrackAssociatorByHitsImpl::TrackAssociatorByHitsImpl ( edm::EDProductGetter const &  productGetter,
std::unique_ptr< TrackerHitAssociator iAssociate,
TrackerTopology const *  iTopo,
SimHitTPAssociationList const *  iSimHitsTPAssoc,
SimToRecoDenomType  iSimToRecoDenominator,
double  iQuality_SimToReco,
double  iPurity_SimToReco,
double  iCut_RecoToSim,
bool  iUsePixels,
bool  iUseGrouped,
bool  iUseSplitting,
bool  ThreeHitTracksAreSpecial,
bool  AbsoluteNumberOfHits 
)

Definition at line 33 of file TrackAssociatorByHitsImpl.cc.

47  associate(std::move(iAssociate)),
48  tTopo(iTopo),
49  simHitsTPAssoc(iSimHitsTPAssoc),
50  SimToRecoDenominator(iSimToRecoDenominator),
51  quality_SimToReco(iQuality_SimToReco),
52  purity_SimToReco(iPurity_SimToReco),
53  cut_RecoToSim(iCut_RecoToSim),
54  UsePixels(iUsePixels),
55  UseGrouped(iUseGrouped),
56  UseSplitting(iUseSplitting),
57  ThreeHitTracksAreSpecial(iThreeHitTracksAreSpecial),
58  AbsoluteNumberOfHits(iAbsoluteNumberOfHits) {}
edm::EDProductGetter const * productGetter_
std::unique_ptr< TrackerHitAssociator > associate
SimHitTPAssociationList const * simHitsTPAssoc
EDProductGetter const * productGetter(std::atomic< void const *> const &iCache)
def move(src, dest)
Definition: eostools.py:511

Member Function Documentation

◆ associateRecoToSim() [1/3]

RecoToSimCollection TrackAssociatorByHitsImpl::associateRecoToSim ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  TPCollectionH 
) const
overridevirtual

Association Reco To Sim with Collections.

Implements reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 92 of file TrackAssociatorByHitsImpl.cc.

References funct::abs(), AbsoluteNumberOfHits, associate, edm::RefToBaseVector< T >::begin(), cut_RecoToSim, edm::RefToBaseVector< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), edm::AssociationMap< Tag >::post_insert(), productGetter_, quality, submitPVValidationJobs::t, ThreeHitTracksAreSpecial, and HLT_2023v12_cff::track.

94  {
95  //edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHitsImpl::associateRecoToSim - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
96  int nshared = 0;
97  float quality = 0; //fraction or absolute number of shared hits
98  std::vector<SimHitIdpr> SimTrackIds;
99  std::vector<SimHitIdpr> matchedIds;
100  RecoToSimCollection outputCollection(productGetter_);
101 
102  //dereference the edm::Refs only once
103  std::vector<TrackingParticle const*> tPC;
104  tPC.reserve(tPC.size());
105  for (auto const& ref : TPCollectionH) {
106  tPC.push_back(&(*ref));
107  }
108 
109  //get the ID of the recotrack by hits
110  int tindex = 0;
111  for (edm::RefToBaseVector<reco::Track>::const_iterator track = tC.begin(); track != tC.end(); track++, tindex++) {
112  matchedIds.clear();
113  int ri = 0; //valid rechits
114  //LogTrace("TrackAssociator") << "\nNEW TRACK - track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found();
115  getMatchedIds<trackingRecHit_iterator>(
116  matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate.get());
117 
118  //LogTrace("TrackAssociator") << "MATCHED IDS LIST BEGIN" ;
119  //for(size_t j=0; j<matchedIds.size(); j++){
120  // LogTrace("TrackAssociator") << "matchedIds[j].first=" << matchedIds[j].first;
121  //}
122  //LogTrace("TrackAssociator") << "MATCHED IDS LIST END" ;
123  //LogTrace("TrackAssociator") << "#matched ids=" << matchedIds.size() << " #tps=" << tPC.size();
124 
125  //save id for the track
126  std::vector<SimHitIdpr> idcachev;
127  if (!matchedIds.empty()) {
128  int tpindex = 0;
129  for (auto t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
130  //int nsimhit = (*t)->trackPSimHit(DetId::Tracker).size();
131  //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << (*t)->pdgId() << " with number of PSimHits: " << nsimhit;
132  idcachev.clear();
133  nshared = getShared(matchedIds, idcachev, **t);
134 
135  //if electron subtract double counting
136  if (abs((*t)->pdgId()) == 11 && ((*t)->g4Track_end() - (*t)->g4Track_begin()) > 1) {
137  nshared -= getDoubleCount<trackingRecHit_iterator>(
138  (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate.get(), **t);
139  }
140 
142  quality = static_cast<double>(nshared);
143  else if (ri != 0)
144  quality = (static_cast<double>(nshared) / static_cast<double>(ri));
145  else
146  quality = 0;
147  //cut on the fraction
148  //float purity = 1.0*nshared/ri;
149  if (quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri == 3 && nshared < 3)) {
150  //if a track has just 3 hits we require that all 3 hits are shared
151  outputCollection.insert(tC[tindex], std::make_pair(TPCollectionH[tpindex], quality));
152  // LogTrace("TrackAssociator") << "reco::Track number " << tindex << " with #hits=" << ri <<" pt=" << (*track)->pt()
153  // << " associated to TP (pdgId, nb segments, p) = "
154  // << (*t).pdgId() << " " << (*t).g4Tracks().size()
155  // << " " << (*t).momentum() << " #hits=" << nsimhit
156  // << " TP index=" << tpindex<< " with quality =" << quality;
157  } else {
158  // LogTrace("TrackAssociator") <<"reco::Track number " << tindex << " with #hits=" << ri
159  // << " NOT associated with quality =" << quality;
160  }
161  } //TP loop
162  }
163  }
164  //LogTrace("TrackAssociator") << "% of Assoc Tracks=" << ((double)outputCollection.size())/((double)tC.size());
165  outputCollection.post_insert();
166  return outputCollection;
167 }
edm::EDProductGetter const * productGetter_
std::unique_ptr< TrackerHitAssociator > associate
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticle const &) const
string quality
const_iterator begin() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const_iterator end() const

◆ associateRecoToSim() [2/3]

reco::RecoToSimCollection TrackAssociatorByHitsImpl::associateRecoToSim ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< TrackingParticleCollection > &  tPCH 
) const
inlineoverride

compare reco to sim the handle of reco::Track and TrackingParticle collections

Definition at line 57 of file TrackAssociatorByHitsImpl.h.

58  {
59  return TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim(tCH, tPCH);
60  }

◆ associateRecoToSim() [3/3]

RecoToSimCollectionSeed TrackAssociatorByHitsImpl::associateRecoToSim ( const edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
const edm::Handle< TrackingParticleCollection > &  TPCollectionH 
) const
override

Definition at line 325 of file TrackAssociatorByHitsImpl.cc.

References funct::abs(), AbsoluteNumberOfHits, associate, edm::View< T >::begin(), cut_RecoToSim, edm::View< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LogTrace, edm::AssociationMap< Tag >::post_insert(), edm::Handle< T >::product(), productGetter_, quality, fileCollector::seed, edm::AssociationMap< Tag >::size(), submitPVValidationJobs::t, and ThreeHitTracksAreSpecial.

327  {
328  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHitsImpl::associateRecoToSim - #seeds="
329  << seedCollectionH->size() << " #TPs=" << TPCollectionH->size();
330  int nshared = 0;
331  float quality = 0; //fraction or absolute number of shared hits
332  std::vector<SimHitIdpr> SimTrackIds;
333  std::vector<SimHitIdpr> matchedIds;
334  RecoToSimCollectionSeed outputCollection(productGetter_);
335 
336  const TrackingParticleCollection& tPC = *(TPCollectionH.product());
337 
338  const edm::View<TrajectorySeed> sC = *(seedCollectionH.product());
339 
340  //get the ID of the recotrack by hits
341  int tindex = 0;
342  for (edm::View<TrajectorySeed>::const_iterator seed = sC.begin(); seed != sC.end(); seed++, tindex++) {
343  matchedIds.clear();
344  int ri = 0; //valid rechits
345  int nsimhit = seed->nHits();
346  LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << nsimhit;
347  getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(
348  matchedIds, SimTrackIds, ri, seed->recHits().begin(), seed->recHits().end(), associate.get());
349 
350  //save id for the track
351  std::vector<SimHitIdpr> idcachev;
352  if (!matchedIds.empty()) {
353  int tpindex = 0;
354  for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
355  LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId()
356  << " with number of PSimHits: " << nsimhit;
357  idcachev.clear();
358  nshared = getShared(matchedIds, idcachev, *t);
359 
360  //if electron subtract double counting
361  if (abs(t->pdgId()) == 11 && (t->g4Track_end() - t->g4Track_begin()) > 1) {
362  nshared -= getDoubleCount<edm::OwnVector<TrackingRecHit>::const_iterator>(
363  seed->recHits().begin(), seed->recHits().end(), associate.get(), *t);
364  }
365 
367  quality = static_cast<double>(nshared);
368  else if (ri != 0)
369  quality = (static_cast<double>(nshared) / static_cast<double>(ri));
370  else
371  quality = 0;
372  //cut on the fraction
373  if (quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri == 3 && nshared < 3)) {
374  //if a track has just 3 hits we require that all 3 hits are shared
375  outputCollection.insert(
376  edm::RefToBase<TrajectorySeed>(seedCollectionH, tindex),
377  std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex), quality));
378  LogTrace("TrackAssociator") << "Seed number " << tindex << " with #hits=" << ri
379  << "associated to TP (pdgId, nb segments, p) = " << (*t).pdgId() << " "
380  << (*t).g4Tracks().size() << " " << (*t).momentum() << " number " << tpindex
381  << " with quality =" << quality;
382  } else {
383  LogTrace("TrackAssociator") << "Seed number " << tindex << " with #hits=" << ri
384  << " NOT associated with quality =" << quality;
385  }
386  } //TP loop
387  }
388  }
389  LogTrace("TrackAssociator") << "% of Assoc Seeds="
390  << ((double)outputCollection.size()) / ((double)seedCollectionH->size());
391 
392  outputCollection.post_insert();
393  return outputCollection;
394 }
Log< level::Info, true > LogVerbatim
T const * product() const
Definition: Handle.h:70
edm::EDProductGetter const * productGetter_
std::unique_ptr< TrackerHitAssociator > associate
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticle const &) const
#define LogTrace(id)
string quality
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
const_iterator begin() const
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator end() const

◆ associateSimToReco() [1/3]

SimToRecoCollection TrackAssociatorByHitsImpl::associateSimToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  TPCollectionH 
) const
overridevirtual

Association Sim To Reco with Collections.

Implements reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 169 of file TrackAssociatorByHitsImpl.cc.

References AbsoluteNumberOfHits, associate, edm::RefToBaseVector< T >::begin(), denomreco, denomsim, edm::RefToBaseVector< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), TrackerTopology::layer(), SiStripDetId::partnerDetId(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, edm::AssociationMap< Tag >::post_insert(), productGetter_, purity_SimToReco, quality, quality_SimToReco, FastTimerService_cff::range, DetId::rawId(), simHitsTPAssoc, SimHitTPAssociationProducer::simHitTPAssociationListGreater(), SimToRecoDenominator, DetId::subdetId(), submitPVValidationJobs::t, SiStripDetId::TEC, ThreeHitTracksAreSpecial, SiStripDetId::TIB, SiStripDetId::TID, SiStripDetId::TOB, HLT_2023v12_cff::track, tTopo, UseGrouped, UsePixels, and UseSplitting.

171  {
172  //edm::ESHandle<TrackerTopology> tTopoHand;
173  //setup->get<IdealGeometryRecord>().get(tTopoHand);
174 
175  // edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHitsImpl::associateSimToReco - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
176  float quality = 0; //fraction or absolute number of shared hits
177  int nshared = 0;
178  std::vector<SimHitIdpr> SimTrackIds;
179  std::vector<SimHitIdpr> matchedIds;
180  SimToRecoCollection outputCollection(productGetter_);
181 
182  //dereferene the edm::Refs only once
183  std::vector<TrackingParticle const*> tPC;
184  tPC.reserve(tPC.size());
185  for (auto const& ref : TPCollectionH) {
186  tPC.push_back(&(*ref));
187  }
188 
189  //for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t) {
190  // LogTrace("TrackAssociator") << "NEW TP DUMP";
191  // for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin();g4T != t -> g4Track_end(); ++g4T) {
192  // LogTrace("TrackAssociator") << "(*g4T).trackId()=" <<(*g4T).trackId() ;
193  // }
194  //}
195 
196  //cdj edm::Handle<SimHitTPAssociationProducer::SimHitTPAssociationList> simHitsTPAssoc;
197  //warning: make sure the TP collection used in the map is the same used in the associator!
198  //e->getByLabel(_simHitTpMapTag,simHitsTPAssoc);
199 
200  //get the ID of the recotrack by hits
201  int tindex = 0;
202  for (edm::RefToBaseVector<reco::Track>::const_iterator track = tC.begin(); track != tC.end(); track++, tindex++) {
203  //LogTrace("TrackAssociator") << "\nNEW TRACK - hits of track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found();
204  int ri = 0; //valid rechits
205  getMatchedIds<trackingRecHit_iterator>(
206  matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate.get());
207 
208  //save id for the track
209  std::vector<SimHitIdpr> idcachev;
210  if (!matchedIds.empty()) {
211  int tpindex = 0;
212  for (auto t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
213  idcachev.clear();
214  float totsimhit = 0;
215 
216  //int nsimhit = trackerPSimHit.size();
217  std::vector<PSimHit> tphits;
218  //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << (*t)->pdgId() << " with number of PSimHits: " << nsimhit;
219 
220  nshared = getShared(matchedIds, idcachev, **t);
221 
222  //for(std::vector<PSimHit>::const_iterator TPhit = (*t)->trackerPSimHit_begin(); TPhit != (*t)->trackerPSimHit_end(); TPhit++){
223  // unsigned int detid = TPhit->detUnitId();
224  // DetId detId = DetId(TPhit->detUnitId());
225  // LogTrace("TrackAssociator") << " hit trackId= " << TPhit->trackId() << " det ID = " << detid
226  // << " SUBDET = " << detId.subdetId() << " layer = " << LayerFromDetid(detId);
227  //}
228 
229  if (nshared != 0) { //do not waste time recounting when it is not needed!!!!
230 
231  //count the TP simhit
232  //LogTrace("TrackAssociator") << "recounting of tp hits";
233 
234  std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(
235  TPCollectionH[tpindex], TrackPSimHitRef()); //SimHit is dummy: for simHitTPAssociationListGreater
236  // sorting only the cluster is needed
237  auto range = std::equal_range(simHitsTPAssoc->begin(),
238  simHitsTPAssoc->end(),
239  clusterTPpairWithDummyTP,
241  for (auto ip = range.first; ip != range.second; ++ip) {
242  TrackPSimHitRef TPhit = ip->second;
243  DetId dId = DetId(TPhit->detUnitId());
244 
245  unsigned int subdetId = static_cast<unsigned int>(dId.subdetId());
246  if (!UsePixels && (subdetId == PixelSubdetector::PixelBarrel || subdetId == PixelSubdetector::PixelEndcap))
247  continue;
248 
249  //unsigned int dRawId = dId.rawId();
250  SiStripDetId* stripDetId = nullptr;
251  if (subdetId == SiStripDetId::TIB || subdetId == SiStripDetId::TOB || subdetId == SiStripDetId::TID ||
252  subdetId == SiStripDetId::TEC)
253  stripDetId = new SiStripDetId(dId);
254  //LogTrace("TrackAssociator") << "consider hit SUBDET = " << subdetId
255  // << " layer = " << LayerFromDetid(dId)
256  // << " id = " << dId.rawId();
257  bool newhit = true;
258  for (std::vector<PSimHit>::const_iterator TPhitOK = tphits.begin(); TPhitOK != tphits.end(); TPhitOK++) {
259  DetId dIdOK = DetId(TPhitOK->detUnitId());
260  //unsigned int dRawIdOK = dIdOK.rawId();
261  //LogTrace("TrackAssociator") << "\t\tcompare with SUBDET = " << dIdOK.subdetId()
262  // << " layer = " << LayerFromDetid(dIdOK)
263  // << " id = " << dIdOK.rawId();
264  //no grouped, no splitting
265  if (!UseGrouped && !UseSplitting)
266  if (tTopo->layer(dId) == tTopo->layer(dIdOK) && dId.subdetId() == dIdOK.subdetId())
267  newhit = false;
268  //no grouped, splitting
269  if (!UseGrouped && UseSplitting)
270  if (tTopo->layer(dId) == tTopo->layer(dIdOK) && dId.subdetId() == dIdOK.subdetId() &&
271  (stripDetId == nullptr || stripDetId->partnerDetId() != dIdOK.rawId()))
272  newhit = false;
273  //grouped, no splitting
274  if (UseGrouped && !UseSplitting)
275  if (tTopo->layer(dId) == tTopo->layer(dIdOK) && dId.subdetId() == dIdOK.subdetId() &&
276  stripDetId != nullptr && stripDetId->partnerDetId() == dIdOK.rawId())
277  newhit = false;
278  //grouped, splitting
279  if (UseGrouped && UseSplitting)
280  newhit = true;
281  }
282  if (newhit) {
283  //LogTrace("TrackAssociator") << "\t\tok";
284  tphits.push_back(*TPhit);
285  }
286  //else LogTrace("TrackAssociator") << "\t\tno";
287  delete stripDetId;
288  }
289  totsimhit = tphits.size();
290  }
291 
293  quality = static_cast<double>(nshared);
294  else if (SimToRecoDenominator == denomsim && totsimhit != 0)
295  quality = ((double)nshared) / ((double)totsimhit);
296  else if (SimToRecoDenominator == denomreco && ri != 0)
297  quality = ((double)nshared) / ((double)ri);
298  else
299  quality = 0;
300  //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit << " re-counted = " << totsimhit
301  //<< " nshared = " << nshared << " nrechit = " << ri;
302 
303  float purity = 1.0 * nshared / ri;
304  if (quality > quality_SimToReco && !(ThreeHitTracksAreSpecial && totsimhit == 3 && nshared < 3) &&
305  (AbsoluteNumberOfHits || (purity > purity_SimToReco))) {
306  //if a track has just 3 hits we require that all 3 hits are shared
307  outputCollection.insert(TPCollectionH[tpindex], std::make_pair(tC[tindex], quality));
308  // LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
309  // << " re-counted = " << totsimhit << " nshared = " << nshared
310  // << " associated to track number " << tindex << " with pt=" << (*track)->pt()
311  // << " with hit quality =" << quality ;
312  } else {
313  // LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
314  // << " re-counted = " << totsimhit << " nshared = " << nshared
315  // << " NOT associated with quality =" << quality;
316  }
317  }
318  }
319  }
320  //LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH.size());
321  outputCollection.post_insert();
322  return outputCollection;
323 }
static bool simHitTPAssociationListGreater(SimHitTPPair i, SimHitTPPair j)
static constexpr auto TID
Definition: SiStripDetId.h:38
edm::EDProductGetter const * productGetter_
std::unique_ptr< TrackerHitAssociator > associate
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticle const &) const
unsigned int layer(const DetId &id) const
SimHitTPAssociationList const * simHitsTPAssoc
string quality
const_iterator begin() const
const_iterator end() const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
static constexpr auto TOB
Definition: SiStripDetId.h:39
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:18
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
static constexpr auto TIB
Definition: SiStripDetId.h:37
uint32_t partnerDetId() const
Definition: SiStripDetId.h:170
edm::Ref< edm::PSimHitContainer > TrackPSimHitRef
static constexpr auto TEC
Definition: SiStripDetId.h:40

◆ associateSimToReco() [2/3]

reco::SimToRecoCollection TrackAssociatorByHitsImpl::associateSimToReco ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< TrackingParticleCollection > &  tPCH 
) const
inlineoverride

compare reco to sim the handle of reco::Track and TrackingParticle collections

Definition at line 64 of file TrackAssociatorByHitsImpl.h.

65  {
66  return TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(tCH, tPCH);
67  }

◆ associateSimToReco() [3/3]

SimToRecoCollectionSeed TrackAssociatorByHitsImpl::associateSimToReco ( const edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
const edm::Handle< TrackingParticleCollection > &  TPCollectionH 
) const
override

Definition at line 396 of file TrackAssociatorByHitsImpl.cc.

References AbsoluteNumberOfHits, associate, edm::View< T >::begin(), edm::View< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LogTrace, edm::AssociationMap< Tag >::post_insert(), edm::Handle< T >::product(), productGetter_, quality, quality_SimToReco, fileCollector::seed, edm::AssociationMap< Tag >::size(), submitPVValidationJobs::t, and ThreeHitTracksAreSpecial.

398  {
399  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHitsImpl::associateSimToReco - #seeds="
400  << seedCollectionH->size() << " #TPs=" << TPCollectionH->size();
401  float quality = 0; //fraction or absolute number of shared hits
402  int nshared = 0;
403  std::vector<SimHitIdpr> SimTrackIds;
404  std::vector<SimHitIdpr> matchedIds;
405  SimToRecoCollectionSeed outputCollection(productGetter_);
406 
407  const TrackingParticleCollection& tPC = *TPCollectionH.product();
408 
409  const edm::View<TrajectorySeed> sC = *(seedCollectionH.product());
410 
411  //get the ID of the recotrack by hits
412  int tindex = 0;
413  for (edm::View<TrajectorySeed>::const_iterator seed = sC.begin(); seed != sC.end(); seed++, tindex++) {
414  int ri = 0; //valid rechits
415  LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << seed->nHits();
416  getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(
417  matchedIds, SimTrackIds, ri, seed->recHits().begin(), seed->recHits().end(), associate.get());
418 
419  //save id for the track
420  std::vector<SimHitIdpr> idcachev;
421  if (!matchedIds.empty()) {
422  int tpindex = 0;
423  for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
424  idcachev.clear();
425  int nsimhit = t->numberOfTrackerHits();
426  LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId()
427  << " with number of PSimHits: " << nsimhit;
428  nshared = getShared(matchedIds, idcachev, *t);
429 
431  quality = static_cast<double>(nshared);
432  else if (ri != 0)
433  quality = ((double)nshared) / ((double)ri);
434  else
435  quality = 0;
436  //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit
437  //<< " nshared = " << nshared
438  //<< " nrechit = " << ri;
439  if (quality > quality_SimToReco && !(ThreeHitTracksAreSpecial && ri == 3 && nshared < 3)) {
440  outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),
441  std::make_pair(edm::RefToBase<TrajectorySeed>(seedCollectionH, tindex), quality));
442  LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
443  << " associated to seed number " << tindex << " with #hits=" << ri
444  << " with hit quality =" << quality;
445  } else {
446  LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
447  << " NOT associated with quality =" << quality;
448  }
449  }
450  }
451  }
452  LogTrace("TrackAssociator") << "% of Assoc TPs="
453  << ((double)outputCollection.size()) / ((double)TPCollectionH->size());
454  outputCollection.post_insert();
455  return outputCollection;
456 }
Log< level::Info, true > LogVerbatim
T const * product() const
Definition: Handle.h:70
edm::EDProductGetter const * productGetter_
std::unique_ptr< TrackerHitAssociator > associate
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticle const &) const
#define LogTrace(id)
string quality
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
const_iterator begin() const
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator end() const

◆ getDoubleCount()

template<typename iter >
int TrackAssociatorByHitsImpl::getDoubleCount ( iter  begin,
iter  end,
TrackerHitAssociator associate,
TrackingParticle const &  t 
) const
private

Definition at line 540 of file TrackAssociatorByHitsImpl.cc.

References associate, mps_fire::end, spr::find(), getHitPtr(), and submitPVValidationJobs::t.

543  {
544  int doublecount = 0;
545  std::vector<SimHitIdpr> SimTrackIdsDC;
546  // cout<<begin-end<<endl;
547  for (iter it = begin; it != end; it++) {
548  int idcount = 0;
549  SimTrackIdsDC.clear();
550  associate->associateHitId(*getHitPtr(it), SimTrackIdsDC);
551  // cout<<SimTrackIdsDC.size()<<endl;
552  if (SimTrackIdsDC.size() > 1) {
553  // cout<<(t.g4Track_end()-t.g4Track_begin())<<endl;
554  for (TrackingParticle::g4t_iterator g4T = t.g4Track_begin(); g4T != t.g4Track_end(); ++g4T) {
555  if (find(SimTrackIdsDC.begin(),
556  SimTrackIdsDC.end(),
557  SimHitIdpr((*g4T).trackId(), SimTrackIdsDC.begin()->second)) != SimTrackIdsDC.end()) {
558  idcount++;
559  }
560  }
561  }
562  if (idcount > 1)
563  doublecount += (idcount - 1);
564  }
565  return doublecount;
566 }
std::unique_ptr< TrackerHitAssociator > associate
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< SimTrack >::const_iterator g4t_iterator
const TrackingRecHit * getHitPtr(edm::OwnVector< TrackingRecHit >::const_iterator iter) const
std::pair< uint32_t, EncodedEventId > SimHitIdpr

◆ getHitPtr() [1/2]

const TrackingRecHit* TrackAssociatorByHitsImpl::getHitPtr ( edm::OwnVector< TrackingRecHit >::const_iterator  iter) const
inlineprivate

Definition at line 102 of file TrackAssociatorByHitsImpl.h.

Referenced by getDoubleCount(), and getMatchedIds().

102 { return &*iter; }

◆ getHitPtr() [2/2]

const TrackingRecHit* TrackAssociatorByHitsImpl::getHitPtr ( trackingRecHit_iterator  iter) const
inlineprivate

Definition at line 103 of file TrackAssociatorByHitsImpl.h.

103 { return &**iter; }

◆ getMatchedIds()

template<typename iter >
void TrackAssociatorByHitsImpl::getMatchedIds ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  SimTrackIds,
int &  ri,
iter  begin,
iter  end,
TrackerHitAssociator associate 
) const
private

Definition at line 459 of file TrackAssociatorByHitsImpl.cc.

References associate, mps_fire::end, getHitPtr(), dqmiolumiharvest::j, and LogTrace.

464  {
465  matchedIds.clear();
466  ri = 0; //valid rechits
467  for (iter it = begin; it != end; it++) {
468  const TrackingRecHit* hit = getHitPtr(it);
469  if (hit->isValid()) {
470  ri++;
471  uint32_t t_detID = hit->geographicalId().rawId();
472  SimTrackIds.clear();
473  associate->associateHitId(*hit, SimTrackIds);
474  //save all the id of matched simtracks
475  if (!SimTrackIds.empty()) {
476  for (size_t j = 0; j < SimTrackIds.size(); j++) {
477  LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << hit->isValid() << " det id = " << t_detID
478  << " SimId " << SimTrackIds[j].first << " evt=" << SimTrackIds[j].second.event()
479  << " bc=" << SimTrackIds[j].second.bunchCrossing();
480  matchedIds.push_back(SimTrackIds[j]);
481  }
482  }
484  //****to be used when the crossing frame is in the event and with flag TrackAssociatorByHitsImplESProducer.associateRecoTracks = false
485  //std::vector<PSimHit> tmpSimHits = associate->associateHit(*getHitPtr(it));
487  //for(size_t j=0; j<tmpSimHits.size(); j++) {
488  // LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << getHitPtr(it)->isValid()
489  // << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
490  // << " evt=" << SimTrackIds[j].second.event()
491  // << " bc=" << SimTrackIds[j].second.bunchCrossing()
492  // << " process=" << tmpSimHits[j].processType()
493  // << " eloss=" << tmpSimHits[j].energyLoss()
494  // << " part type=" << tmpSimHits[j].particleType()
495  // << " pabs=" << tmpSimHits[j].pabs()
496  // << " trId=" << tmpSimHits[j].trackId();
497  //}
499  } else {
500  LogTrace("TrackAssociator") << "\t\t Invalid Hit On " << hit->geographicalId().rawId();
501  }
502  } //trackingRecHit loop
503 }
std::unique_ptr< TrackerHitAssociator > associate
#define LogTrace(id)
const TrackingRecHit * getHitPtr(edm::OwnVector< TrackingRecHit >::const_iterator iter) const

◆ getShared()

int TrackAssociatorByHitsImpl::getShared ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  idcachev,
TrackingParticle const &  t 
) const
private

Definition at line 505 of file TrackAssociatorByHitsImpl.cc.

References submitPVResolutionJobs::count, spr::find(), dqmiolumiharvest::j, and submitPVValidationJobs::t.

Referenced by associateRecoToSim(), and associateSimToReco().

507  {
508  int nshared = 0;
509  if (t.numberOfHits() == 0)
510  return nshared; //should use trackerPSimHit but is not const
511 
512  for (size_t j = 0; j < matchedIds.size(); j++) {
513  //LogTrace("TrackAssociator") << "now matchedId=" << matchedIds[j].first;
514  if (find(idcachev.begin(), idcachev.end(), matchedIds[j]) == idcachev.end()) {
515  //only the first time we see this ID
516  idcachev.push_back(matchedIds[j]);
517 
518  for (TrackingParticle::g4t_iterator g4T = t.g4Track_begin(); g4T != t.g4Track_end(); ++g4T) {
519  // LogTrace("TrackAssociator") << " TP (ID, Ev, BC) = " << (*g4T).trackId()
520  // << ", " << t.eventId().event() << ", "<< t.eventId().bunchCrossing()
521  // << " Match(ID, Ev, BC) = " << matchedIds[j].first
522  // << ", " << matchedIds[j].second.event() << ", "
523  // << matchedIds[j].second.bunchCrossing() ;
524  //<< "\t G4 Track Momentum " << (*g4T).momentum()
525  //<< " \t reco Track Momentum " << track->momentum();
526  if ((*g4T).trackId() == matchedIds[j].first && t.eventId() == matchedIds[j].second) {
527  int countedhits = std::count(matchedIds.begin(), matchedIds.end(), matchedIds[j]);
528  nshared += countedhits;
529 
530  // LogTrace("TrackAssociator") << "hits shared by this segment : " << countedhits;
531  // LogTrace("TrackAssociator") << "hits shared so far : " << nshared;
532  }
533  } //g4Tracks loop
534  }
535  }
536  return nshared;
537 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< SimTrack >::const_iterator g4t_iterator

Member Data Documentation

◆ AbsoluteNumberOfHits

const bool TrackAssociatorByHitsImpl::AbsoluteNumberOfHits
private

Definition at line 100 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

◆ associate

std::unique_ptr<TrackerHitAssociator> TrackAssociatorByHitsImpl::associate
private

◆ cut_RecoToSim

const double TrackAssociatorByHitsImpl::cut_RecoToSim
private

Definition at line 95 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim().

◆ productGetter_

edm::EDProductGetter const* TrackAssociatorByHitsImpl::productGetter_
private

Definition at line 87 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

◆ purity_SimToReco

const double TrackAssociatorByHitsImpl::purity_SimToReco
private

Definition at line 94 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ quality_SimToReco

const double TrackAssociatorByHitsImpl::quality_SimToReco
private

Definition at line 93 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ simHitsTPAssoc

SimHitTPAssociationList const* TrackAssociatorByHitsImpl::simHitsTPAssoc
private

Definition at line 90 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ SimToRecoDenominator

SimToRecoDenomType TrackAssociatorByHitsImpl::SimToRecoDenominator
private

Definition at line 92 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ ThreeHitTracksAreSpecial

const bool TrackAssociatorByHitsImpl::ThreeHitTracksAreSpecial
private

Definition at line 99 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

◆ tTopo

TrackerTopology const* TrackAssociatorByHitsImpl::tTopo
private

Definition at line 89 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ UseGrouped

const bool TrackAssociatorByHitsImpl::UseGrouped
private

Definition at line 97 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ UsePixels

const bool TrackAssociatorByHitsImpl::UsePixels
private

Definition at line 96 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

◆ UseSplitting

const bool TrackAssociatorByHitsImpl::UseSplitting
private

Definition at line 98 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().