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
TrackAssociatorByHits Class Reference

#include <TrackAssociatorByHits.h>

Inheritance diagram for TrackAssociatorByHits:
TrackAssociatorBase

Public Types

enum  SimToRecoDenomType { denomnone, denomsim, denomreco }
 

Public Member Functions

reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 Association Reco To Sim with Collections. More...
 
reco::RecoToSimCollection associateRecoToSim (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
reco::RecoToSimCollectionSeed associateRecoToSim (edm::Handle< edm::View< TrajectorySeed > > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 
reco::SimToRecoCollection associateSimToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 Association Sim To Reco with Collections. More...
 
reco::SimToRecoCollection associateSimToReco (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
reco::SimToRecoCollectionSeed associateSimToReco (edm::Handle< edm::View< TrajectorySeed > > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 
template<typename iter >
int getDoubleCount (iter, iter, TrackerHitAssociator *, TrackingParticleCollection::const_iterator) 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 > &, TrackingParticleCollection::const_iterator) const
 
 TrackAssociatorByHits (const edm::ParameterSet &)
 
 ~TrackAssociatorByHits ()
 
- Public Member Functions inherited from TrackAssociatorBase
virtual
reco::RecoToSimCollectionSeed 
associateRecoToSim (const edm::Handle< edm::View< TrajectorySeed > > &, const edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 
virtual
reco::RecoToSimCollectionTCandidate 
associateRecoToSim (const edm::Handle< TrackCandidateCollection > &, const edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 
virtual
reco::SimToRecoCollectionSeed 
associateSimToReco (const edm::Handle< edm::View< TrajectorySeed > > &, const edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 
virtual
reco::SimToRecoCollectionTCandidate 
associateSimToReco (const edm::Handle< TrackCandidateCollection > &, const edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 
 TrackAssociatorBase ()
 Constructor. More...
 
virtual ~TrackAssociatorBase ()
 Destructor. More...
 

Private Member Functions

const TrackingRecHitgetHitPtr (edm::OwnVector< TrackingRecHit >::const_iterator iter) const
 
const TrackingRecHitgetHitPtr (trackingRecHit_iterator iter) const
 
int LayerFromDetid (const DetId &) const
 

Private Attributes

const bool AbsoluteNumberOfHits
 
const edm::ParameterSetconf_
 
const double cut_RecoToSim
 
const double purity_SimToReco
 
const double quality_SimToReco
 
SimToRecoDenomType SimToRecoDenominator
 
const bool ThreeHitTracksAreSpecial
 
const bool UseGrouped
 
const bool UsePixels
 
const bool UseSplitting
 

Detailed Description

Definition at line 21 of file TrackAssociatorByHits.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

TrackAssociatorByHits::TrackAssociatorByHits ( const edm::ParameterSet conf)
explicit

Definition at line 38 of file TrackAssociatorByHits.cc.

References conf_, denomnone, denomreco, denomsim, edm::hlt::Exception, edm::ParameterSet::getParameter(), SimToRecoDenominator, and tmp.

38  :
39  conf_(conf),
40  AbsoluteNumberOfHits(conf_.getParameter<bool>("AbsoluteNumberOfHits")),
42  quality_SimToReco(conf_.getParameter<double>("Quality_SimToReco")),
43  purity_SimToReco(conf_.getParameter<double>("Purity_SimToReco")),
44  cut_RecoToSim(conf_.getParameter<double>("Cut_RecoToSim")),
45  UsePixels(conf_.getParameter<bool>("UsePixels")),
46  UseGrouped(conf_.getParameter<bool>("UseGrouped")),
47  UseSplitting(conf_.getParameter<bool>("UseSplitting")),
48  ThreeHitTracksAreSpecial(conf_.getParameter<bool>("ThreeHitTracksAreSpecial"))
49 {
50  std::string tmp = conf_.getParameter<string>("SimToRecoDenominator");
51  if (tmp=="sim") {
53  } else if (tmp == "reco") {
55  }
56 
58  throw cms::Exception("TrackAssociatorByHits") << "SimToRecoDenominator not specified as sim or reco";
59  }
60 
61 }
T getParameter(std::string const &) const
const edm::ParameterSet & conf_
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
SimToRecoDenomType SimToRecoDenominator
TrackAssociatorByHits::~TrackAssociatorByHits ( )

Definition at line 65 of file TrackAssociatorByHits.cc.

66 {
67  //do cleanup here
68 }

Member Function Documentation

RecoToSimCollection TrackAssociatorByHits::associateRecoToSim ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  TPCollectionH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
virtual

Association Reco To Sim with Collections.

Implements TrackAssociatorBase.

Definition at line 75 of file TrackAssociatorByHits.cc.

References abs, AbsoluteNumberOfHits, edm::RefToBaseVector< T >::begin(), conf_, cut_RecoToSim, edm::RefToBaseVector< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), edm::RefVector< C, T, F >::size(), lumiQTWidget::t, and ThreeHitTracksAreSpecial.

78  {
79 
80  //edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateRecoToSim - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
81  int nshared = 0;
82  float quality=0;//fraction or absolute number of shared hits
83  std::vector< SimHitIdpr> SimTrackIds;
84  std::vector< SimHitIdpr> matchedIds;
85  RecoToSimCollection outputCollection;
86 
87  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
88 
90  if (TPCollectionH.size()!=0) tPC = *const_cast<TrackingParticleCollection*>(TPCollectionH.product());
91 
92  //get the ID of the recotrack by hits
93  int tindex=0;
94  for (edm::RefToBaseVector<reco::Track>::const_iterator track=tC.begin(); track!=tC.end(); track++, tindex++){
95  matchedIds.clear();
96  int ri=0;//valid rechits
97  //LogTrace("TrackAssociator") << "\nNEW TRACK - track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found();
98  getMatchedIds<trackingRecHit_iterator>(matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate);
99 
100  //LogTrace("TrackAssociator") << "MATCHED IDS LIST BEGIN" ;
101  //for(size_t j=0; j<matchedIds.size(); j++){
102  // LogTrace("TrackAssociator") << "matchedIds[j].first=" << matchedIds[j].first;
103  //}
104  //LogTrace("TrackAssociator") << "MATCHED IDS LIST END" ;
105  //LogTrace("TrackAssociator") << "#matched ids=" << matchedIds.size() << " #tps=" << tPC.size();
106 
107  //save id for the track
108  std::vector<SimHitIdpr> idcachev;
109  if(!matchedIds.empty()){
110 
111  int tpindex =0;
112  for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
113  //int nsimhit = t->trackPSimHit(DetId::Tracker).size();
114  //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: " << nsimhit;
115  idcachev.clear();
116  nshared = getShared(matchedIds, idcachev, t);
117 
118  //if electron subtract double counting
119  if (abs(t->pdgId())==11&&(t->g4Track_end()-t->g4Track_begin())>1){
120  nshared-=getDoubleCount<trackingRecHit_iterator>((*track)->recHitsBegin(), (*track)->recHitsEnd(), associate, t);
121  }
122 
123  if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
124  else if(ri!=0) quality = (static_cast<double>(nshared)/static_cast<double>(ri));
125  else quality = 0;
126  //cut on the fraction
127  //float purity = 1.0*nshared/ri;
128  if(quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3)){
129  //if a track has just 3 hits we require that all 3 hits are shared
130  outputCollection.insert(tC[tindex],
131  std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),
132  quality));
133 // LogTrace("TrackAssociator") << "reco::Track number " << tindex << " with #hits=" << ri <<" pt=" << (*track)->pt()
134 // << " associated to TP (pdgId, nb segments, p) = "
135 // << (*t).pdgId() << " " << (*t).g4Tracks().size()
136 // << " " << (*t).momentum() << " #hits=" << nsimhit
137 // << " TP index=" << tpindex<< " with quality =" << quality;
138  } else {
139 // LogTrace("TrackAssociator") <<"reco::Track number " << tindex << " with #hits=" << ri
140 // << " NOT associated with quality =" << quality;
141  }
142  }//TP loop
143  }
144  }
145  //LogTrace("TrackAssociator") << "% of Assoc Tracks=" << ((double)outputCollection.size())/((double)tC.size());
146  delete associate;
147  outputCollection.post_insert();
148  return outputCollection;
149 }
std::vector< TrackingParticle > TrackingParticleCollection
C const * product() const
Accessor for product collection.
Definition: RefVector.h:272
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticleCollection::const_iterator) const
const_iterator end() const
#define abs(x)
Definition: mlp_lapack.h:159
void post_insert()
post insert action
void insert(const key_type &k, const data_type &v)
insert an association
const edm::ParameterSet & conf_
const_iterator begin() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
reco::RecoToSimCollection TrackAssociatorByHits::associateRecoToSim ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
inlinevirtual

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

Reimplemented from TrackAssociatorBase.

Definition at line 44 of file TrackAssociatorByHits.h.

References TrackAssociatorBase::associateRecoToSim(), event(), and HcalObjRepresent::setup().

47  {
48  return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event,setup);
49  }
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
RecoToSimCollectionSeed TrackAssociatorByHits::associateRecoToSim ( edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
edm::Handle< TrackingParticleCollection > &  TPCollectionH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const

Definition at line 332 of file TrackAssociatorByHits.cc.

References abs, AbsoluteNumberOfHits, edm::View< T >::begin(), conf_, cut_RecoToSim, edm::View< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LogTrace, edm::AssociationMap< Tag >::post_insert(), edm::Handle< T >::product(), edm::AssociationMap< Tag >::size(), lumiQTWidget::t, and ThreeHitTracksAreSpecial.

335  {
336 
337  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateRecoToSim - #seeds="
338  <<seedCollectionH->size()<<" #TPs="<<TPCollectionH->size();
339  int nshared = 0;
340  float quality=0;//fraction or absolute number of shared hits
341  std::vector< SimHitIdpr> SimTrackIds;
342  std::vector< SimHitIdpr> matchedIds;
343  RecoToSimCollectionSeed outputCollection;
344 
345  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
346 
347  const TrackingParticleCollection tPC = *(TPCollectionH.product());
348 
349  const edm::View<TrajectorySeed> sC = *(seedCollectionH.product());
350 
351  //get the ID of the recotrack by hits
352  int tindex=0;
353  for (edm::View<TrajectorySeed>::const_iterator seed=sC.begin(); seed!=sC.end(); seed++, tindex++) {
354  matchedIds.clear();
355  int ri=0;//valid rechits
356  int nsimhit = seed->recHits().second-seed->recHits().first;
357  LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << nsimhit;
358  getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(matchedIds, SimTrackIds, ri, seed->recHits().first, seed->recHits().second, associate );
359 
360  //save id for the track
361  std::vector<SimHitIdpr> idcachev;
362  if(!matchedIds.empty()){
363 
364  int tpindex =0;
365  for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
366  LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: " << nsimhit;
367  idcachev.clear();
368  nshared = getShared(matchedIds, idcachev, t);
369 
370  //if electron subtract double counting
371  if (abs(t->pdgId())==11&&(t->g4Track_end()-t->g4Track_begin())>1){
372  nshared-=getDoubleCount<edm::OwnVector<TrackingRecHit>::const_iterator>(seed->recHits().first, seed->recHits().second, associate, t);
373  }
374 
375  if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
376  else if(ri!=0) quality = (static_cast<double>(nshared)/static_cast<double>(ri));
377  else quality = 0;
378  //cut on the fraction
379  if(quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3) ){
380  //if a track has just 3 hits we require that all 3 hits are shared
381  outputCollection.insert(edm::RefToBase<TrajectorySeed>(seedCollectionH,tindex),
382  std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),quality));
383  LogTrace("TrackAssociator") << "Seed number " << tindex << " with #hits=" << ri
384  << "associated to TP (pdgId, nb segments, p) = "
385  << (*t).pdgId() << " " << (*t).g4Tracks().size()
386  << " " << (*t).momentum() <<" number " << tpindex << " with quality =" << quality;
387  } else {
388  LogTrace("TrackAssociator") <<"Seed number " << tindex << " with #hits=" << ri << " NOT associated with quality =" << quality;
389  }
390  }//TP loop
391  }
392  }
393  LogTrace("TrackAssociator") << "% of Assoc Seeds=" << ((double)outputCollection.size())/((double)seedCollectionH->size());
394  delete associate;
395  outputCollection.post_insert();
396  return outputCollection;
397 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::vector< TrackingParticle > TrackingParticleCollection
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticleCollection::const_iterator) const
#define abs(x)
Definition: mlp_lapack.h:159
void post_insert()
post insert action
#define LogTrace(id)
size_type size() const
map size
void insert(const key_type &k, const data_type &v)
insert an association
const edm::ParameterSet & conf_
T const * product() const
Definition: Handle.h:74
const_iterator begin() const
const_iterator end() const
SimToRecoCollection TrackAssociatorByHits::associateSimToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  TPCollectionH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
virtual

Association Sim To Reco with Collections.

Implements TrackAssociatorBase.

Definition at line 153 of file TrackAssociatorByHits.cc.

References AbsoluteNumberOfHits, edm::RefToBaseVector< T >::begin(), conf_, denomreco, denomsim, edm::RefToBaseVector< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LayerFromDetid(), SiStripDetId::partnerDetId(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), purity_SimToReco, quality_SimToReco, DetId::rawId(), SimToRecoDenominator, edm::RefVector< C, T, F >::size(), DetId::subdetId(), lumiQTWidget::t, SiStripDetId::TEC, ThreeHitTracksAreSpecial, SiStripDetId::TIB, SiStripDetId::TID, SiStripDetId::TOB, DetId::Tracker, UseGrouped, UsePixels, and UseSplitting.

156  {
157 // edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateSimToReco - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
158  float quality=0;//fraction or absolute number of shared hits
159  int nshared = 0;
160  std::vector< SimHitIdpr> SimTrackIds;
161  std::vector< SimHitIdpr> matchedIds;
162  SimToRecoCollection outputCollection;
163 
164  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
165 
167  if (TPCollectionH.size()!=0) tPC = *const_cast<TrackingParticleCollection*>(TPCollectionH.product());
168 
169  //for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t) {
170  // LogTrace("TrackAssociator") << "NEW TP DUMP";
171  // for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin();g4T != t -> g4Track_end(); ++g4T) {
172  // LogTrace("TrackAssociator") << "(*g4T).trackId()=" <<(*g4T).trackId() ;
173  // }
174  //}
175 
176  //get the ID of the recotrack by hits
177  int tindex=0;
178  for (edm::RefToBaseVector<reco::Track>::const_iterator track=tC.begin(); track!=tC.end(); track++, tindex++){
179  //LogTrace("TrackAssociator") << "\nNEW TRACK - hits of track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found();
180  int ri=0;//valid rechits
181  getMatchedIds<trackingRecHit_iterator>(matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate);
182 
183  //save id for the track
184  std::vector<SimHitIdpr> idcachev;
185  if(!matchedIds.empty()){
186 
187  int tpindex =0;
188  for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
189  idcachev.clear();
190  std::vector<PSimHit> trackerPSimHit( t->trackPSimHit(DetId::Tracker) );
191  //int nsimhit = trackerPSimHit.size();
192  float totsimhit = 0;
193  std::vector<PSimHit> tphits;
194  //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: " << nsimhit;
195 
196  nshared = getShared(matchedIds, idcachev, t);
197 
198  //for(std::vector<PSimHit>::const_iterator TPhit = t->trackerPSimHit_begin(); TPhit != t->trackerPSimHit_end(); TPhit++){
199  // unsigned int detid = TPhit->detUnitId();
200  // DetId detId = DetId(TPhit->detUnitId());
201  // LogTrace("TrackAssociator") << " hit trackId= " << TPhit->trackId() << " det ID = " << detid
202  // << " SUBDET = " << detId.subdetId() << " layer = " << LayerFromDetid(detId);
203  //}
204 
205  if (nshared!=0) {//do not waste time recounting when it is not needed!!!!
206 
207  //count the TP simhit
208  //LogTrace("TrackAssociator") << "recounting of tp hits";
209  for(std::vector<PSimHit>::const_iterator TPhit = trackerPSimHit.begin(); TPhit != trackerPSimHit.end(); TPhit++){
210  DetId dId = DetId(TPhit->detUnitId());
211 
212  unsigned int subdetId = static_cast<unsigned int>(dId.subdetId());
214  continue;
215 
216  //unsigned int dRawId = dId.rawId();
217  SiStripDetId* stripDetId = 0;
218  if (subdetId==SiStripDetId::TIB||subdetId==SiStripDetId::TOB||
219  subdetId==SiStripDetId::TID||subdetId==SiStripDetId::TEC)
220  stripDetId= new SiStripDetId(dId);
221  //LogTrace("TrackAssociator") << "consider hit SUBDET = " << subdetId
222  // << " layer = " << LayerFromDetid(dId)
223  // << " id = " << dId.rawId();
224  bool newhit = true;
225  for(std::vector<PSimHit>::const_iterator TPhitOK = tphits.begin(); TPhitOK != tphits.end(); TPhitOK++){
226  DetId dIdOK = DetId(TPhitOK->detUnitId());
227  //unsigned int dRawIdOK = dIdOK.rawId();
228  //LogTrace("TrackAssociator") << "\t\tcompare with SUBDET = " << dIdOK.subdetId()
229  // << " layer = " << LayerFromDetid(dIdOK)
230  // << " id = " << dIdOK.rawId();
231  //no grouped, no splitting
232  if (!UseGrouped && !UseSplitting)
233  if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
234  dId.subdetId()==dIdOK.subdetId()) newhit = false;
235  //no grouped, splitting
236  if (!UseGrouped && UseSplitting)
237  if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
238  dId.subdetId()==dIdOK.subdetId() &&
239  (stripDetId==0 || stripDetId->partnerDetId()!=dIdOK.rawId()))
240  newhit = false;
241  //grouped, no splitting
242  if (UseGrouped && !UseSplitting)
243  if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
244  dId.subdetId()==dIdOK.subdetId() &&
245  stripDetId!=0 && stripDetId->partnerDetId()==dIdOK.rawId())
246  newhit = false;
247  //grouped, splitting
248  if (UseGrouped && UseSplitting)
249  newhit = true;
250  }
251  if (newhit) {
252  //LogTrace("TrackAssociator") << "\t\tok";
253  tphits.push_back(*TPhit);
254  }
255  //else LogTrace("TrackAssociator") << "\t\tno";
256  delete stripDetId;
257  }
258  totsimhit = tphits.size();
259  }
260 
261  if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
262  else if(SimToRecoDenominator == denomsim && totsimhit!=0) quality = ((double) nshared)/((double)totsimhit);
263  else if(SimToRecoDenominator == denomreco && ri!=0) quality = ((double) nshared)/((double)ri);
264  else quality = 0;
265  //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit << " re-counted = " << totsimhit
266  //<< " nshared = " << nshared << " nrechit = " << ri;
267 
268  float purity = 1.0*nshared/ri;
269  if (quality>quality_SimToReco && !(ThreeHitTracksAreSpecial && totsimhit==3 && nshared<3) && (AbsoluteNumberOfHits||(purity>purity_SimToReco))) {
270  //if a track has just 3 hits we require that all 3 hits are shared
271  outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),
272  std::make_pair(tC[tindex],quality));
273 // LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
274 // << " re-counted = " << totsimhit << " nshared = " << nshared
275 // << " associated to track number " << tindex << " with pt=" << (*track)->pt()
276 // << " with hit quality =" << quality ;
277  } else {
278 // LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
279 // << " re-counted = " << totsimhit << " nshared = " << nshared
280 // << " NOT associated with quality =" << quality;
281  }
282  }
283  }
284  }
285  //LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH.size());
286  delete associate;
287  outputCollection.post_insert();
288  return outputCollection;
289 }
std::vector< TrackingParticle > TrackingParticleCollection
C const * product() const
Accessor for product collection.
Definition: RefVector.h:272
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticleCollection::const_iterator) const
const_iterator end() const
int LayerFromDetid(const DetId &) const
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
void post_insert()
post insert action
uint32_t partnerDetId() const
Definition: SiStripDetId.h:168
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
Definition: DetId.h:20
void insert(const key_type &k, const data_type &v)
insert an association
const edm::ParameterSet & conf_
const_iterator begin() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
SimToRecoDenomType SimToRecoDenominator
reco::SimToRecoCollection TrackAssociatorByHits::associateSimToReco ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
inlinevirtual

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

Reimplemented from TrackAssociatorBase.

Definition at line 52 of file TrackAssociatorByHits.h.

References TrackAssociatorBase::associateSimToReco(), event(), and HcalObjRepresent::setup().

55  {
56  return TrackAssociatorBase::associateSimToReco(tCH,tPCH,event,setup);
57  }
virtual reco::SimToRecoCollection associateSimToReco(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
SimToRecoCollectionSeed TrackAssociatorByHits::associateSimToReco ( edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
edm::Handle< TrackingParticleCollection > &  TPCollectionH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const

Definition at line 401 of file TrackAssociatorByHits.cc.

References AbsoluteNumberOfHits, edm::View< T >::begin(), conf_, edm::View< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LogTrace, edm::AssociationMap< Tag >::post_insert(), edm::Handle< T >::product(), quality_SimToReco, edm::AssociationMap< Tag >::size(), lumiQTWidget::t, ThreeHitTracksAreSpecial, and DetId::Tracker.

404  {
405 
406  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateSimToReco - #seeds="
407  <<seedCollectionH->size()<<" #TPs="<<TPCollectionH->size();
408  float quality=0;//fraction or absolute number of shared hits
409  int nshared = 0;
410  std::vector< SimHitIdpr> SimTrackIds;
411  std::vector< SimHitIdpr> matchedIds;
412  SimToRecoCollectionSeed outputCollection;
413 
414  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
415 
416  TrackingParticleCollection tPC =*const_cast<TrackingParticleCollection*>(TPCollectionH.product());
417 
418  const edm::View<TrajectorySeed> sC = *(seedCollectionH.product());
419 
420  //get the ID of the recotrack by hits
421  int tindex=0;
422  for (edm::View<TrajectorySeed>::const_iterator seed=sC.begin(); seed!=sC.end(); seed++, tindex++) {
423  int ri=0;//valid rechits
424  LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << seed->recHits().second-seed->recHits().first;
425  getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(matchedIds, SimTrackIds, ri, seed->recHits().first, seed->recHits().second, associate );
426 
427  //save id for the track
428  std::vector<SimHitIdpr> idcachev;
429  if(!matchedIds.empty()){
430  int tpindex =0;
431  for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
432  idcachev.clear();
433  int nsimhit = t->trackPSimHit(DetId::Tracker).size();
434  LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: " << nsimhit;
435  nshared = getShared(matchedIds, idcachev, t);
436 
437  if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
438  else if(ri!=0) quality = ((double) nshared)/((double)ri);
439  else quality = 0;
440  //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit
441  //<< " nshared = " << nshared
442  //<< " nrechit = " << ri;
443  if(quality > quality_SimToReco && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3) ){
444  outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),
445  std::make_pair(edm::RefToBase<TrajectorySeed>(seedCollectionH,tindex), quality));
446  LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
447  << " associated to seed number " << tindex << " with #hits=" << ri
448  << " with hit quality =" << quality ;
449  }
450  else {
451  LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit << " NOT associated with quality =" << quality;
452  }
453  }
454  }
455  }
456  LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH->size());
457  delete associate;
458  outputCollection.post_insert();
459  return outputCollection;
460 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::vector< TrackingParticle > TrackingParticleCollection
int getShared(std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticleCollection::const_iterator) const
void post_insert()
post insert action
#define LogTrace(id)
size_type size() const
map size
void insert(const key_type &k, const data_type &v)
insert an association
const edm::ParameterSet & conf_
T const * product() const
Definition: Handle.h:74
const_iterator begin() const
const_iterator end() const
template<typename iter >
int TrackAssociatorByHits::getDoubleCount ( iter  begin,
iter  end,
TrackerHitAssociator associate,
TrackingParticleCollection::const_iterator  t 
) const

Definition at line 546 of file TrackAssociatorByHits.cc.

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

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

Definition at line 98 of file TrackAssociatorByHits.h.

Referenced by getDoubleCount(), and getMatchedIds().

98 {return &*iter;}
const TrackingRecHit* TrackAssociatorByHits::getHitPtr ( trackingRecHit_iterator  iter) const
inlineprivate

Definition at line 99 of file TrackAssociatorByHits.h.

99 {return &**iter;}
template<typename iter >
void TrackAssociatorByHits::getMatchedIds ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  SimTrackIds,
int &  ri,
iter  begin,
iter  end,
TrackerHitAssociator associate 
) const

Definition at line 463 of file TrackAssociatorByHits.cc.

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

468  {
469  matchedIds.clear();
470  ri=0;//valid rechits
471  for (iter it = begin; it != end; it++){
472  const TrackingRecHit *hit=getHitPtr(it);
473  if (hit->isValid()){
474  ri++;
475  uint32_t t_detID= hit->geographicalId().rawId();
476  SimTrackIds.clear();
477  associate->associateHitId(*hit, SimTrackIds);
478  //save all the id of matched simtracks
479  if(!SimTrackIds.empty()){
480  for(size_t j=0; j<SimTrackIds.size(); j++){
481  LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << hit->isValid()
482  << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
483  << " evt=" << SimTrackIds[j].second.event()
484  << " bc=" << SimTrackIds[j].second.bunchCrossing();
485  matchedIds.push_back(SimTrackIds[j]);
486  }
487  }
489  //****to be used when the crossing frame is in the event and with flag TrackAssociatorByHitsESProducer.associateRecoTracks = false
490  //std::vector<PSimHit> tmpSimHits = associate->associateHit(*getHitPtr(it));
492  //for(size_t j=0; j<tmpSimHits.size(); j++) {
493  // LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << getHitPtr(it)->isValid()
494  // << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
495  // << " evt=" << SimTrackIds[j].second.event()
496  // << " bc=" << SimTrackIds[j].second.bunchCrossing()
497  // << " process=" << tmpSimHits[j].processType()
498  // << " eloss=" << tmpSimHits[j].energyLoss()
499  // << " part type=" << tmpSimHits[j].particleType()
500  // << " pabs=" << tmpSimHits[j].pabs()
501  // << " trId=" << tmpSimHits[j].trackId();
502  //}
504  }else{
505  LogTrace("TrackAssociator") <<"\t\t Invalid Hit On "<<hit->geographicalId().rawId();
506  }
507  }//trackingRecHit loop
508 }
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit)
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
const TrackingRecHit * getHitPtr(edm::OwnVector< TrackingRecHit >::const_iterator iter) const
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:38
#define LogTrace(id)
bool isValid() const
#define begin
Definition: vmac.h:31
DetId geographicalId() const
int TrackAssociatorByHits::getShared ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  idcachev,
TrackingParticleCollection::const_iterator  t 
) const

Definition at line 511 of file TrackAssociatorByHits.cc.

References prof2calltree::count, spr::find(), and j.

Referenced by associateRecoToSim(), and associateSimToReco().

513  {
514  int nshared = 0;
515  if (t->trackPSimHit().size()==0) return nshared;//should use trackerPSimHit but is not const
516 
517  for(size_t j=0; j<matchedIds.size(); j++){
518  //LogTrace("TrackAssociator") << "now matchedId=" << matchedIds[j].first;
519  if(find(idcachev.begin(), idcachev.end(),matchedIds[j]) == idcachev.end() ){
520  //only the first time we see this ID
521  idcachev.push_back(matchedIds[j]);
522 
523  for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin(); g4T != t -> g4Track_end(); ++g4T) {
524 // LogTrace("TrackAssociator") << " TP (ID, Ev, BC) = " << (*g4T).trackId()
525 // << ", " << t->eventId().event() << ", "<< t->eventId().bunchCrossing()
526 // << " Match(ID, Ev, BC) = " << matchedIds[j].first
527 // << ", " << matchedIds[j].second.event() << ", "
528 // << matchedIds[j].second.bunchCrossing() ;
529  //<< "\t G4 Track Momentum " << (*g4T).momentum()
530  //<< " \t reco Track Momentum " << track->momentum();
531  if((*g4T).trackId() == matchedIds[j].first && t->eventId() == matchedIds[j].second){
532  int countedhits = std::count(matchedIds.begin(), matchedIds.end(), matchedIds[j]);
533  nshared += countedhits;
534 
535 // LogTrace("TrackAssociator") << "hits shared by this segment : " << countedhits;
536 // LogTrace("TrackAssociator") << "hits shared so far : " << nshared;
537  }
538  }//g4Tracks loop
539  }
540  }
541  return nshared;
542 }
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
int TrackAssociatorByHits::LayerFromDetid ( const DetId detId) const
private

Definition at line 291 of file TrackAssociatorByHits.cc.

References PXFDetId::disk(), PXBDetId::layer(), TOBDetId::layer(), TIBDetId::layer(), align::tib::layerNumber(), LogTrace, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, DetId::rawId(), DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, StripSubdetector::TOB, TIDDetId::wheel(), and TECDetId::wheel().

Referenced by associateSimToReco().

292 {
293  int layerNumber=0;
294  unsigned int subdetId = static_cast<unsigned int>(detId.subdetId());
295  if ( subdetId == StripSubdetector::TIB)
296  {
297  TIBDetId tibid(detId.rawId());
298  layerNumber = tibid.layer();
299  }
300  else if ( subdetId == StripSubdetector::TOB )
301  {
302  TOBDetId tobid(detId.rawId());
303  layerNumber = tobid.layer();
304  }
305  else if ( subdetId == StripSubdetector::TID)
306  {
307  TIDDetId tidid(detId.rawId());
308  layerNumber = tidid.wheel();
309  }
310  else if ( subdetId == StripSubdetector::TEC )
311  {
312  TECDetId tecid(detId.rawId());
313  layerNumber = tecid.wheel();
314  }
315  else if ( subdetId == PixelSubdetector::PixelBarrel )
316  {
317  PXBDetId pxbid(detId.rawId());
318  layerNumber = pxbid.layer();
319  }
320  else if ( subdetId == PixelSubdetector::PixelEndcap )
321  {
322  PXFDetId pxfid(detId.rawId());
323  layerNumber = pxfid.disk();
324  }
325  else
326  LogTrace("TrackAssociator") << "Unknown subdetid: " << subdetId;
327 
328  return layerNumber;
329 }
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
unsigned int layerNumber(align::ID)
Layer number increases with rho from 1 to 8.
Definition: TIBNameSpace.h:85
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
#define LogTrace(id)
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50

Member Data Documentation

const bool TrackAssociatorByHits::AbsoluteNumberOfHits
private

Definition at line 87 of file TrackAssociatorByHits.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const edm::ParameterSet& TrackAssociatorByHits::conf_
private
const double TrackAssociatorByHits::cut_RecoToSim
private

Definition at line 91 of file TrackAssociatorByHits.h.

Referenced by associateRecoToSim().

const double TrackAssociatorByHits::purity_SimToReco
private

Definition at line 90 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

const double TrackAssociatorByHits::quality_SimToReco
private

Definition at line 89 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

SimToRecoDenomType TrackAssociatorByHits::SimToRecoDenominator
private

Definition at line 88 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco(), and TrackAssociatorByHits().

const bool TrackAssociatorByHits::ThreeHitTracksAreSpecial
private

Definition at line 95 of file TrackAssociatorByHits.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const bool TrackAssociatorByHits::UseGrouped
private

Definition at line 93 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHits::UsePixels
private

Definition at line 92 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHits::UseSplitting
private

Definition at line 94 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().