CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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,
TrackPSimHitRef
SimHitTPPair
 
enum  SimToRecoDenomType { denomnone, denomsim, denomreco }
 

Public Member Functions

virtual reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Reco To Sim with Collections. More...
 
virtual 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...
 
virtual
reco::RecoToSimCollectionSeed 
associateRecoToSim (const edm::Handle< edm::View< TrajectorySeed > > &, const edm::Handle< TrackingParticleCollection > &) const override
 
virtual reco::SimToRecoCollection associateSimToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Sim To Reco with Collections. More...
 
virtual 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...
 
virtual
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::RecoToSimCollectionTCandidate 
associateRecoToSim (const edm::Handle< TrackCandidateCollection > &, const edm::Handle< TrackingParticleCollection > &) const
 
virtual
reco::SimToRecoCollectionTCandidate 
associateSimToReco (const edm::Handle< TrackCandidateCollection > &, const edm::Handle< 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
< TrackerHitAssociator
associate
 
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

Definition at line 30 of file TrackAssociatorByHitsImpl.h.

Definition at line 29 of file TrackAssociatorByHitsImpl.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

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.

45  :
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
EDProductGetter const * productGetter(std::atomic< void const * > const &iCache)
SimHitTPAssociationList const * simHitsTPAssoc
def move
Definition: eostools.py:510

Member Function Documentation

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 93 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_, HLT_25ns10e33_v2_cff::quality, lumiQTWidget::t, and ThreeHitTracksAreSpecial.

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

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

Reimplemented from reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 61 of file TrackAssociatorByHitsImpl.h.

62  {
63  return TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim(tCH,tPCH);
64  }
RecoToSimCollectionSeed TrackAssociatorByHitsImpl::associateRecoToSim ( const edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
const edm::Handle< TrackingParticleCollection > &  TPCollectionH 
) const
overridevirtual

Reimplemented from reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 324 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_, HLT_25ns10e33_v2_cff::quality, fileCollector::seed, edm::AssociationMap< Tag >::size(), lumiQTWidget::t, and ThreeHitTracksAreSpecial.

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

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

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

Reimplemented from reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 68 of file TrackAssociatorByHitsImpl.h.

69  {
70  return TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(tCH,tPCH);
71  }
SimToRecoCollectionSeed TrackAssociatorByHitsImpl::associateSimToReco ( const edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
const edm::Handle< TrackingParticleCollection > &  TPCollectionH 
) const
overridevirtual

Reimplemented from reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 388 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_, HLT_25ns10e33_v2_cff::quality, quality_SimToReco, fileCollector::seed, edm::AssociationMap< Tag >::size(), lumiQTWidget::t, and ThreeHitTracksAreSpecial.

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

Definition at line 528 of file TrackAssociatorByHitsImpl.cc.

References TrackerHitAssociator::associateHitId(), end, spr::find(), and getHitPtr().

531  {
532  int doublecount = 0 ;
533  std::vector<SimHitIdpr> SimTrackIdsDC;
534  // cout<<begin-end<<endl;
535  for (iter it = begin; it != end; it++){
536  int idcount = 0;
537  SimTrackIdsDC.clear();
538  associate->associateHitId(*getHitPtr(it), SimTrackIdsDC);
539  // cout<<SimTrackIdsDC.size()<<endl;
540  if(SimTrackIdsDC.size()>1){
541  // cout<<(t.g4Track_end()-t.g4Track_begin())<<endl;
542  for (TrackingParticle::g4t_iterator g4T = t . g4Track_begin(); g4T != t . g4Track_end(); ++g4T) {
543  if(find(SimTrackIdsDC.begin(), SimTrackIdsDC.end(),SimHitIdpr((*g4T).trackId(), SimTrackIdsDC.begin()->second)) != SimTrackIdsDC.end() ){
544  idcount++;
545  }
546  }
547  }
548  if (idcount>1) doublecount+=(idcount-1);
549  }
550  return doublecount;
551 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
#define end
Definition: vmac.h:37
std::vector< SimTrack >::const_iterator g4t_iterator
const TrackingRecHit * getHitPtr(edm::OwnVector< TrackingRecHit >::const_iterator iter) const
std::pair< uint32_t, EncodedEventId > SimHitIdpr
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
#define begin
Definition: vmac.h:30
const TrackingRecHit* TrackAssociatorByHitsImpl::getHitPtr ( edm::OwnVector< TrackingRecHit >::const_iterator  iter) const
inlineprivate

Definition at line 115 of file TrackAssociatorByHitsImpl.h.

Referenced by getDoubleCount(), and getMatchedIds().

115 {return &*iter;}
const TrackingRecHit* TrackAssociatorByHitsImpl::getHitPtr ( trackingRecHit_iterator  iter) const
inlineprivate

Definition at line 116 of file TrackAssociatorByHitsImpl.h.

116 {return &**iter;}
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 445 of file TrackAssociatorByHitsImpl.cc.

References TrackerHitAssociator::associateHitId(), end, TrackingRecHit::geographicalId(), getHitPtr(), TrackingRecHit::isValid(), j, LogTrace, and DetId::rawId().

450  {
451  matchedIds.clear();
452  ri=0;//valid rechits
453  for (iter it = begin; it != end; it++){
454  const TrackingRecHit *hit=getHitPtr(it);
455  if (hit->isValid()){
456  ri++;
457  uint32_t t_detID= hit->geographicalId().rawId();
458  SimTrackIds.clear();
459  associate->associateHitId(*hit, SimTrackIds);
460  //save all the id of matched simtracks
461  if(!SimTrackIds.empty()){
462  for(size_t j=0; j<SimTrackIds.size(); j++){
463  LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << hit->isValid()
464  << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
465  << " evt=" << SimTrackIds[j].second.event()
466  << " bc=" << SimTrackIds[j].second.bunchCrossing();
467  matchedIds.push_back(SimTrackIds[j]);
468  }
469  }
471  //****to be used when the crossing frame is in the event and with flag TrackAssociatorByHitsImplESProducer.associateRecoTracks = false
472  //std::vector<PSimHit> tmpSimHits = associate->associateHit(*getHitPtr(it));
474  //for(size_t j=0; j<tmpSimHits.size(); j++) {
475  // LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << getHitPtr(it)->isValid()
476  // << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
477  // << " evt=" << SimTrackIds[j].second.event()
478  // << " bc=" << SimTrackIds[j].second.bunchCrossing()
479  // << " process=" << tmpSimHits[j].processType()
480  // << " eloss=" << tmpSimHits[j].energyLoss()
481  // << " part type=" << tmpSimHits[j].particleType()
482  // << " pabs=" << tmpSimHits[j].pabs()
483  // << " trId=" << tmpSimHits[j].trackId();
484  //}
486  }else{
487  LogTrace("TrackAssociator") <<"\t\t Invalid Hit On "<<hit->geographicalId().rawId();
488  }
489  }//trackingRecHit loop
490 }
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
#define LogTrace(id)
const TrackingRecHit * getHitPtr(edm::OwnVector< TrackingRecHit >::const_iterator iter) const
bool isValid() const
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
#define begin
Definition: vmac.h:30
DetId geographicalId() const
int TrackAssociatorByHitsImpl::getShared ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  idcachev,
TrackingParticle const &  t 
) const
private

Definition at line 493 of file TrackAssociatorByHitsImpl.cc.

References KineDebug3::count(), TrackingParticle::eventId(), spr::find(), j, and TrackingParticle::numberOfHits().

Referenced by associateRecoToSim(), and associateSimToReco().

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

Member Data Documentation

const bool TrackAssociatorByHitsImpl::AbsoluteNumberOfHits
private

Definition at line 113 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

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

Definition at line 101 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double TrackAssociatorByHitsImpl::cut_RecoToSim
private

Definition at line 108 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim().

edm::EDProductGetter const* TrackAssociatorByHitsImpl::productGetter_
private

Definition at line 100 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double TrackAssociatorByHitsImpl::purity_SimToReco
private

Definition at line 107 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

const double TrackAssociatorByHitsImpl::quality_SimToReco
private

Definition at line 106 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

SimHitTPAssociationList const* TrackAssociatorByHitsImpl::simHitsTPAssoc
private

Definition at line 103 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

SimToRecoDenomType TrackAssociatorByHitsImpl::SimToRecoDenominator
private

Definition at line 105 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHitsImpl::ThreeHitTracksAreSpecial
private

Definition at line 112 of file TrackAssociatorByHitsImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

TrackerTopology const* TrackAssociatorByHitsImpl::tTopo
private

Definition at line 102 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHitsImpl::UseGrouped
private

Definition at line 110 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHitsImpl::UsePixels
private

Definition at line 109 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHitsImpl::UseSplitting
private

Definition at line 111 of file TrackAssociatorByHitsImpl.h.

Referenced by associateSimToReco().