CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/SimTracker/TrackAssociation/src/TrackAssociatorByHits.cc

Go to the documentation of this file.
00001 //
00002 //
00003 #include "FWCore/Framework/interface/Frameworkfwd.h"
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Framework/interface/EventSetup.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "DataFormats/Common/interface/Ref.h"
00009 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByHits.h"
00010 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 
00014 //reco track
00015 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00016 #include "DataFormats/TrackReco/interface/Track.h"
00017 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00018 //TrackingParticle
00019 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00020 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
00021 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00022 //##---new stuff
00023 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00024 #include "DataFormats/DetId/interface/DetId.h"
00025 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00026 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00027 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00028 #include "DataFormats/SiStripDetId/interface/TECDetId.h" 
00029 #include "DataFormats/SiStripDetId/interface/TIBDetId.h" 
00030 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00031 #include "DataFormats/SiStripDetId/interface/TOBDetId.h" 
00032 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00033 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00034 using namespace reco;
00035 using namespace std;
00036 
00037 /* Constructor */
00038 TrackAssociatorByHits::TrackAssociatorByHits (const edm::ParameterSet& conf) :  
00039   conf_(conf),
00040   AbsoluteNumberOfHits(conf_.getParameter<bool>("AbsoluteNumberOfHits")),
00041   SimToRecoDenominator(denomnone),
00042   quality_SimToReco(conf_.getParameter<double>("Quality_SimToReco")),
00043   purity_SimToReco(conf_.getParameter<double>("Purity_SimToReco")),
00044   cut_RecoToSim(conf_.getParameter<double>("Cut_RecoToSim")),
00045   UsePixels(conf_.getParameter<bool>("UsePixels")),
00046   UseGrouped(conf_.getParameter<bool>("UseGrouped")),
00047   UseSplitting(conf_.getParameter<bool>("UseSplitting")),
00048   ThreeHitTracksAreSpecial(conf_.getParameter<bool>("ThreeHitTracksAreSpecial"))
00049 {
00050   std::string tmp = conf_.getParameter<string>("SimToRecoDenominator");
00051   if (tmp=="sim") {
00052     SimToRecoDenominator = denomsim;
00053   } else if (tmp == "reco") {
00054     SimToRecoDenominator = denomreco;
00055   } 
00056 
00057   if (SimToRecoDenominator == denomnone) {
00058     throw cms::Exception("TrackAssociatorByHits") << "SimToRecoDenominator not specified as sim or reco";
00059   }
00060 
00061 }
00062 
00063 
00064 /* Destructor */
00065 TrackAssociatorByHits::~TrackAssociatorByHits()
00066 {
00067   //do cleanup here
00068 }
00069 
00070 //
00071 //---member functions
00072 //
00073 
00074 RecoToSimCollection  
00075 TrackAssociatorByHits::associateRecoToSim(const edm::RefToBaseVector<reco::Track>& tC, 
00076                                           const edm::RefVector<TrackingParticleCollection>& TPCollectionH,
00077                                           const edm::Event * e,
00078                                           const edm::EventSetup *setup ) const{
00079 
00080   //edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateRecoToSim - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
00081   int nshared = 0;
00082   float quality=0;//fraction or absolute number of shared hits
00083   std::vector< SimHitIdpr> SimTrackIds;
00084   std::vector< SimHitIdpr> matchedIds; 
00085   RecoToSimCollection  outputCollection;
00086   
00087   TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
00088   
00089   TrackingParticleCollection tPC;
00090   if (TPCollectionH.size()!=0)  tPC = *const_cast<TrackingParticleCollection*>(TPCollectionH.product());
00091 
00092   //get the ID of the recotrack  by hits 
00093   int tindex=0;
00094   for (edm::RefToBaseVector<reco::Track>::const_iterator track=tC.begin(); track!=tC.end(); track++, tindex++){
00095     matchedIds.clear();
00096     int ri=0;//valid rechits
00097     //LogTrace("TrackAssociator") << "\nNEW TRACK - track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found(); 
00098     getMatchedIds<trackingRecHit_iterator>(matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate);
00099 
00100     //LogTrace("TrackAssociator") << "MATCHED IDS LIST BEGIN" ;
00101     //for(size_t j=0; j<matchedIds.size(); j++){
00102     //  LogTrace("TrackAssociator") << "matchedIds[j].first=" << matchedIds[j].first;
00103     //}
00104     //LogTrace("TrackAssociator") << "MATCHED IDS LIST END" ;
00105     //LogTrace("TrackAssociator") << "#matched ids=" << matchedIds.size() << " #tps=" << tPC.size();
00106 
00107     //save id for the track
00108     std::vector<SimHitIdpr> idcachev;
00109     if(!matchedIds.empty()){
00110 
00111       int tpindex =0;
00112       for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
00113         //int nsimhit = t->trackPSimHit(DetId::Tracker).size(); 
00114         //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
00115         idcachev.clear();
00116         nshared = getShared(matchedIds, idcachev, t);
00117 
00118         //if electron subtract double counting
00119         if (abs(t->pdgId())==11&&(t->g4Track_end()-t->g4Track_begin())>1){
00120           nshared-=getDoubleCount<trackingRecHit_iterator>((*track)->recHitsBegin(), (*track)->recHitsEnd(), associate, t);
00121         }
00122 
00123         if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
00124         else if(ri!=0) quality = (static_cast<double>(nshared)/static_cast<double>(ri));
00125         else quality = 0;
00126         //cut on the fraction
00127         //float purity = 1.0*nshared/ri;
00128         if(quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3)){
00129           //if a track has just 3 hits we require that all 3 hits are shared
00130           outputCollection.insert(tC[tindex],
00131                                   std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),
00132                                                  quality));
00133 //          LogTrace("TrackAssociator") << "reco::Track number " << tindex  << " with #hits=" << ri <<" pt=" << (*track)->pt() 
00134 //                                      << " associated to TP (pdgId, nb segments, p) = " 
00135 //                                      << (*t).pdgId() << " " << (*t).g4Tracks().size() 
00136 //                                      << " " << (*t).momentum() << " #hits=" << nsimhit
00137 //                                      << " TP index=" << tpindex<< " with quality =" << quality;
00138         } else {
00139 //          LogTrace("TrackAssociator") <<"reco::Track number " << tindex << " with #hits=" << ri 
00140 //                                      << " NOT associated with quality =" << quality;
00141         }
00142       }//TP loop
00143     }
00144   }
00145   //LogTrace("TrackAssociator") << "% of Assoc Tracks=" << ((double)outputCollection.size())/((double)tC.size());
00146   delete associate;
00147   outputCollection.post_insert();
00148   return outputCollection;
00149 }
00150 
00151 
00152 SimToRecoCollection  
00153 TrackAssociatorByHits::associateSimToReco(const edm::RefToBaseVector<reco::Track>& tC, 
00154                                           const edm::RefVector<TrackingParticleCollection>& TPCollectionH,
00155                                           const edm::Event * e,
00156                                           const edm::EventSetup *setup ) const{
00157 //  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateSimToReco - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
00158   float quality=0;//fraction or absolute number of shared hits
00159   int nshared = 0;
00160   std::vector< SimHitIdpr> SimTrackIds;
00161   std::vector< SimHitIdpr> matchedIds; 
00162   SimToRecoCollection  outputCollection;
00163 
00164   TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
00165   
00166   TrackingParticleCollection tPC;
00167   if (TPCollectionH.size()!=0)  tPC = *const_cast<TrackingParticleCollection*>(TPCollectionH.product());
00168 
00169   //for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t) {
00170   //  LogTrace("TrackAssociator") << "NEW TP DUMP";
00171   //  for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin();g4T !=  t -> g4Track_end(); ++g4T) {
00172   //    LogTrace("TrackAssociator") << "(*g4T).trackId()=" <<(*g4T).trackId() ;
00173   //  }
00174   //}
00175   
00176   //get the ID of the recotrack  by hits 
00177   int tindex=0;
00178   for (edm::RefToBaseVector<reco::Track>::const_iterator track=tC.begin(); track!=tC.end(); track++, tindex++){
00179     //LogTrace("TrackAssociator") << "\nNEW TRACK - hits of track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found(); 
00180     int ri=0;//valid rechits
00181     getMatchedIds<trackingRecHit_iterator>(matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate);
00182 
00183     //save id for the track
00184     std::vector<SimHitIdpr> idcachev;
00185     if(!matchedIds.empty()){
00186         
00187       int tpindex =0;
00188       for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
00189         idcachev.clear();
00190         std::vector<PSimHit> trackerPSimHit( t->trackPSimHit(DetId::Tracker) );
00191         //int nsimhit = trackerPSimHit.size();
00192         float totsimhit = 0; 
00193         std::vector<PSimHit> tphits;
00194         //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
00195 
00196         nshared = getShared(matchedIds, idcachev, t);
00197 
00198         //for(std::vector<PSimHit>::const_iterator TPhit = t->trackerPSimHit_begin(); TPhit != t->trackerPSimHit_end(); TPhit++){
00199         //  unsigned int detid = TPhit->detUnitId();
00200         //  DetId detId = DetId(TPhit->detUnitId());
00201         //  LogTrace("TrackAssociator") <<  " hit trackId= " << TPhit->trackId() << " det ID = " << detid 
00202         //                                    << " SUBDET = " << detId.subdetId() << " layer = " << LayerFromDetid(detId); 
00203         //}
00204 
00205         if (nshared!=0) {//do not waste time recounting when it is not needed!!!!
00206 
00207           //count the TP simhit
00208           //LogTrace("TrackAssociator") << "recounting of tp hits";
00209           for(std::vector<PSimHit>::const_iterator TPhit = trackerPSimHit.begin(); TPhit != trackerPSimHit.end(); TPhit++){
00210             DetId dId = DetId(TPhit->detUnitId());
00211           
00212             unsigned int subdetId = static_cast<unsigned int>(dId.subdetId());
00213             if (!UsePixels && (subdetId==PixelSubdetector::PixelBarrel || subdetId==PixelSubdetector::PixelEndcap) )
00214               continue;
00215 
00216             //unsigned int dRawId = dId.rawId();
00217             SiStripDetId* stripDetId = 0;
00218             if (subdetId==SiStripDetId::TIB||subdetId==SiStripDetId::TOB||
00219                 subdetId==SiStripDetId::TID||subdetId==SiStripDetId::TEC)
00220               stripDetId= new SiStripDetId(dId);
00221             //LogTrace("TrackAssociator") << "consider hit SUBDET = " << subdetId
00222             //                            << " layer = " << LayerFromDetid(dId) 
00223             //                            << " id = " << dId.rawId();
00224             bool newhit = true;
00225             for(std::vector<PSimHit>::const_iterator TPhitOK = tphits.begin(); TPhitOK != tphits.end(); TPhitOK++){
00226               DetId dIdOK = DetId(TPhitOK->detUnitId());
00227               //unsigned int dRawIdOK = dIdOK.rawId();
00228               //LogTrace("TrackAssociator") << "\t\tcompare with SUBDET = " << dIdOK.subdetId()
00229               //                            << " layer = " << LayerFromDetid(dIdOK)
00230               //                            << " id = " << dIdOK.rawId();
00231               //no grouped, no splitting
00232               if (!UseGrouped && !UseSplitting)
00233                 if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
00234                     dId.subdetId()==dIdOK.subdetId()) newhit = false;
00235               //no grouped, splitting
00236               if (!UseGrouped && UseSplitting)
00237                 if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
00238                     dId.subdetId()==dIdOK.subdetId() &&
00239                     (stripDetId==0 || stripDetId->partnerDetId()!=dIdOK.rawId()))
00240                   newhit = false;
00241               //grouped, no splitting
00242               if (UseGrouped && !UseSplitting)
00243                 if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
00244                     dId.subdetId()==dIdOK.subdetId() &&
00245                     stripDetId!=0 && stripDetId->partnerDetId()==dIdOK.rawId())
00246                   newhit = false;
00247               //grouped, splitting
00248               if (UseGrouped && UseSplitting)
00249                 newhit = true;
00250             }
00251             if (newhit) {
00252               //LogTrace("TrackAssociator") << "\t\tok";
00253               tphits.push_back(*TPhit);
00254             }
00255             //else LogTrace("TrackAssociator") << "\t\tno";
00256             delete stripDetId;
00257           }
00258           totsimhit = tphits.size();
00259         }
00260 
00261         if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
00262         else if(SimToRecoDenominator == denomsim && totsimhit!=0) quality = ((double) nshared)/((double)totsimhit);
00263         else if(SimToRecoDenominator == denomreco && ri!=0) quality = ((double) nshared)/((double)ri);
00264         else quality = 0;
00265         //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit << " re-counted = " << totsimhit 
00266         //<< " nshared = " << nshared << " nrechit = " << ri;
00267         
00268         float purity = 1.0*nshared/ri;
00269         if (quality>quality_SimToReco && !(ThreeHitTracksAreSpecial && totsimhit==3 && nshared<3) && (AbsoluteNumberOfHits||(purity>purity_SimToReco))) {
00270           //if a track has just 3 hits we require that all 3 hits are shared
00271           outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex), 
00272                                   std::make_pair(tC[tindex],quality));
00273 //          LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit 
00274 //                                      << " re-counted = "  << totsimhit << " nshared = " << nshared 
00275 //                                      << " associated to track number " << tindex << " with pt=" << (*track)->pt() 
00276 //                                      << " with hit quality =" << quality ;
00277         } else {
00278 //          LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit 
00279 //                                      << " re-counted = "  << totsimhit << " nshared = " << nshared 
00280 //                                      << " NOT associated with quality =" << quality;
00281         }
00282       }
00283     }
00284   }
00285   //LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH.size());
00286   delete associate;
00287   outputCollection.post_insert();
00288   return outputCollection;
00289 }
00290 
00291 int TrackAssociatorByHits::LayerFromDetid(const DetId& detId ) const
00292 {
00293   int layerNumber=0;
00294   unsigned int subdetId = static_cast<unsigned int>(detId.subdetId()); 
00295   if ( subdetId == StripSubdetector::TIB) 
00296     { 
00297       TIBDetId tibid(detId.rawId()); 
00298       layerNumber = tibid.layer();
00299     }
00300   else if ( subdetId ==  StripSubdetector::TOB )
00301     { 
00302       TOBDetId tobid(detId.rawId()); 
00303       layerNumber = tobid.layer();
00304     }
00305   else if ( subdetId ==  StripSubdetector::TID) 
00306     { 
00307       TIDDetId tidid(detId.rawId());
00308       layerNumber = tidid.wheel();
00309     }
00310   else if ( subdetId ==  StripSubdetector::TEC )
00311     { 
00312       TECDetId tecid(detId.rawId()); 
00313       layerNumber = tecid.wheel(); 
00314     }
00315   else if ( subdetId ==  PixelSubdetector::PixelBarrel ) 
00316     { 
00317       PXBDetId pxbid(detId.rawId()); 
00318       layerNumber = pxbid.layer();  
00319     }
00320   else if ( subdetId ==  PixelSubdetector::PixelEndcap ) 
00321     { 
00322       PXFDetId pxfid(detId.rawId()); 
00323       layerNumber = pxfid.disk();  
00324     }
00325   else
00326     LogTrace("TrackAssociator") << "Unknown subdetid: " <<  subdetId;
00327   
00328   return layerNumber;
00329 } 
00330 
00331 RecoToSimCollectionSeed  
00332 TrackAssociatorByHits::associateRecoToSim(edm::Handle<edm::View<TrajectorySeed> >& seedCollectionH,
00333                                           edm::Handle<TrackingParticleCollection>&  TPCollectionH,     
00334                                           const edm::Event * e,
00335                                           const edm::EventSetup *setup ) const{
00336 
00337   edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateRecoToSim - #seeds="
00338                                       <<seedCollectionH->size()<<" #TPs="<<TPCollectionH->size();
00339   int nshared = 0;
00340   float quality=0;//fraction or absolute number of shared hits
00341   std::vector< SimHitIdpr> SimTrackIds;
00342   std::vector< SimHitIdpr> matchedIds; 
00343   RecoToSimCollectionSeed  outputCollection;
00344   
00345   TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
00346   
00347   const TrackingParticleCollection tPC   = *(TPCollectionH.product());
00348 
00349   const edm::View<TrajectorySeed> sC = *(seedCollectionH.product()); 
00350   
00351   //get the ID of the recotrack  by hits 
00352   int tindex=0;
00353   for (edm::View<TrajectorySeed>::const_iterator seed=sC.begin(); seed!=sC.end(); seed++, tindex++) {
00354     matchedIds.clear();
00355     int ri=0;//valid rechits
00356     int nsimhit = seed->recHits().second-seed->recHits().first;
00357     LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << nsimhit;
00358     getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(matchedIds, SimTrackIds, ri, seed->recHits().first, seed->recHits().second, associate );
00359 
00360     //save id for the track
00361     std::vector<SimHitIdpr> idcachev;
00362     if(!matchedIds.empty()){
00363 
00364       int tpindex =0;
00365       for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
00366         LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
00367         idcachev.clear();
00368         nshared = getShared(matchedIds, idcachev, t);
00369 
00370         //if electron subtract double counting
00371         if (abs(t->pdgId())==11&&(t->g4Track_end()-t->g4Track_begin())>1){
00372           nshared-=getDoubleCount<edm::OwnVector<TrackingRecHit>::const_iterator>(seed->recHits().first, seed->recHits().second, associate, t);
00373         }
00374         
00375         if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
00376         else if(ri!=0) quality = (static_cast<double>(nshared)/static_cast<double>(ri));
00377         else quality = 0;
00378         //cut on the fraction
00379         if(quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3) ){
00380           //if a track has just 3 hits we require that all 3 hits are shared
00381           outputCollection.insert(edm::RefToBase<TrajectorySeed>(seedCollectionH,tindex), 
00382                                   std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),quality));
00383           LogTrace("TrackAssociator") << "Seed number " << tindex << " with #hits=" << ri
00384                                       << "associated to TP (pdgId, nb segments, p) = " 
00385                                       << (*t).pdgId() << " " << (*t).g4Tracks().size() 
00386                                       << " " << (*t).momentum() <<" number " << tpindex << " with quality =" << quality;
00387         } else {
00388           LogTrace("TrackAssociator") <<"Seed number " << tindex << " with #hits=" << ri << " NOT associated with quality =" << quality;
00389         }
00390       }//TP loop
00391     }
00392   }
00393   LogTrace("TrackAssociator") << "% of Assoc Seeds=" << ((double)outputCollection.size())/((double)seedCollectionH->size());
00394   delete associate;
00395   outputCollection.post_insert();
00396   return outputCollection;
00397 }
00398 
00399 
00400 SimToRecoCollectionSeed
00401 TrackAssociatorByHits::associateSimToReco(edm::Handle<edm::View<TrajectorySeed> >& seedCollectionH,
00402                                           edm::Handle<TrackingParticleCollection>& TPCollectionH, 
00403                                           const edm::Event * e,
00404                                           const edm::EventSetup *setup ) const{
00405 
00406   edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateSimToReco - #seeds="
00407                                       <<seedCollectionH->size()<<" #TPs="<<TPCollectionH->size();
00408   float quality=0;//fraction or absolute number of shared hits
00409   int nshared = 0;
00410   std::vector< SimHitIdpr> SimTrackIds;
00411   std::vector< SimHitIdpr> matchedIds; 
00412   SimToRecoCollectionSeed  outputCollection;
00413 
00414   TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
00415   
00416   TrackingParticleCollection tPC =*const_cast<TrackingParticleCollection*>(TPCollectionH.product());
00417 
00418   const edm::View<TrajectorySeed> sC = *(seedCollectionH.product()); 
00419 
00420   //get the ID of the recotrack  by hits 
00421   int tindex=0;
00422   for (edm::View<TrajectorySeed>::const_iterator seed=sC.begin(); seed!=sC.end(); seed++, tindex++) {
00423     int ri=0;//valid rechits
00424     LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << seed->recHits().second-seed->recHits().first;
00425     getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(matchedIds, SimTrackIds, ri, seed->recHits().first, seed->recHits().second, associate );
00426 
00427     //save id for the track
00428     std::vector<SimHitIdpr> idcachev;
00429     if(!matchedIds.empty()){
00430       int tpindex =0;
00431       for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
00432         idcachev.clear();
00433         int nsimhit = t->trackPSimHit(DetId::Tracker).size(); 
00434         LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
00435         nshared = getShared(matchedIds, idcachev, t);
00436         
00437         if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
00438         else if(ri!=0) quality = ((double) nshared)/((double)ri);
00439         else quality = 0;
00440         //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit 
00441         //<< " nshared = " << nshared 
00442         //<< " nrechit = " << ri;
00443         if(quality > quality_SimToReco && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3) ){
00444           outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex), 
00445                                   std::make_pair(edm::RefToBase<TrajectorySeed>(seedCollectionH,tindex), quality));
00446           LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
00447                                       << " associated to seed number " << tindex << " with #hits=" << ri 
00448                                       << " with hit quality =" << quality ;
00449         }
00450         else {
00451           LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit << " NOT associated with quality =" << quality;
00452         }
00453       }
00454     }
00455   }
00456   LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH->size());
00457   delete associate;
00458   outputCollection.post_insert();
00459   return outputCollection;
00460 }
00461 
00462 template<typename iter>
00463 void TrackAssociatorByHits::getMatchedIds(std::vector<SimHitIdpr>& matchedIds, 
00464                                           std::vector<SimHitIdpr>& SimTrackIds, 
00465                                           int& ri, 
00466                                           iter begin,
00467                                           iter end,
00468                                           TrackerHitAssociator* associate ) const {
00469     matchedIds.clear();
00470     ri=0;//valid rechits
00471     for (iter it = begin;  it != end; it++){
00472       const TrackingRecHit *hit=getHitPtr(it);
00473       if (hit->isValid()){
00474         ri++;
00475         uint32_t t_detID=  hit->geographicalId().rawId();
00476         SimTrackIds.clear();      
00477         associate->associateHitId(*hit, SimTrackIds);
00478         //save all the id of matched simtracks
00479         if(!SimTrackIds.empty()){
00480          for(size_t j=0; j<SimTrackIds.size(); j++){
00481            LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << hit->isValid() 
00482                                        << " det id = " << t_detID << " SimId " << SimTrackIds[j].first 
00483                                        << " evt=" << SimTrackIds[j].second.event() 
00484                                        << " bc=" << SimTrackIds[j].second.bunchCrossing();  
00485            matchedIds.push_back(SimTrackIds[j]);                        
00486          }
00487         }
00489         //****to be used when the crossing frame is in the event and with flag TrackAssociatorByHitsESProducer.associateRecoTracks = false
00490         //std::vector<PSimHit> tmpSimHits = associate->associateHit(*getHitPtr(it));
00492         //for(size_t j=0; j<tmpSimHits.size(); j++) {
00493         //  LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << getHitPtr(it)->isValid()
00494         //                              << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
00495         //                              << " evt=" << SimTrackIds[j].second.event()
00496         //                              << " bc=" << SimTrackIds[j].second.bunchCrossing()
00497         //                              << " process=" << tmpSimHits[j].processType()
00498         //                              << " eloss=" << tmpSimHits[j].energyLoss()
00499         //                              << " part type=" << tmpSimHits[j].particleType()
00500         //                              << " pabs=" << tmpSimHits[j].pabs()
00501         //                              << " trId=" << tmpSimHits[j].trackId();
00502         //}
00504       }else{
00505         LogTrace("TrackAssociator") <<"\t\t Invalid Hit On "<<hit->geographicalId().rawId();
00506       }
00507     }//trackingRecHit loop
00508 }
00509 
00510 
00511 int TrackAssociatorByHits::getShared(std::vector<SimHitIdpr>& matchedIds, 
00512                                      std::vector<SimHitIdpr>& idcachev,
00513                                      TrackingParticleCollection::const_iterator t) const {
00514   int nshared = 0;
00515   if (t->trackPSimHit().size()==0) return nshared;//should use trackerPSimHit but is not const
00516 
00517   for(size_t j=0; j<matchedIds.size(); j++){
00518     //LogTrace("TrackAssociator") << "now matchedId=" << matchedIds[j].first;
00519     if(find(idcachev.begin(), idcachev.end(),matchedIds[j]) == idcachev.end() ){
00520       //only the first time we see this ID 
00521       idcachev.push_back(matchedIds[j]);
00522       
00523       for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin(); g4T !=  t -> g4Track_end(); ++g4T) {
00524 //      LogTrace("TrackAssociator") << " TP   (ID, Ev, BC) = " << (*g4T).trackId() 
00525 //                                  << ", " << t->eventId().event() << ", "<< t->eventId().bunchCrossing()
00526 //                                  << " Match(ID, Ev, BC) = " <<  matchedIds[j].first
00527 //                                  << ", " << matchedIds[j].second.event() << ", "
00528 //                                  << matchedIds[j].second.bunchCrossing() ;
00529                                     //<< "\t G4  Track Momentum " << (*g4T).momentum() 
00530                                     //<< " \t reco Track Momentum " << track->momentum();             
00531         if((*g4T).trackId() == matchedIds[j].first && t->eventId() == matchedIds[j].second){
00532                 int countedhits = std::count(matchedIds.begin(), matchedIds.end(), matchedIds[j]);
00533                 nshared += countedhits;
00534 
00535 //              LogTrace("TrackAssociator") << "hits shared by this segment : " << countedhits;
00536 //              LogTrace("TrackAssociator") << "hits shared so far : " << nshared;
00537         }
00538       }//g4Tracks loop
00539     }
00540   }
00541   return nshared;
00542 }
00543 
00544 
00545 template<typename iter>
00546 int TrackAssociatorByHits::getDoubleCount(iter begin,
00547                                           iter end,
00548                                           TrackerHitAssociator* associate,
00549                                           TrackingParticleCollection::const_iterator t) const {
00550   int doublecount = 0 ;
00551   std::vector<SimHitIdpr> SimTrackIdsDC;
00552   //  cout<<begin-end<<endl;
00553   for (iter it = begin;  it != end; it++){
00554     int idcount = 0;
00555     SimTrackIdsDC.clear();
00556     associate->associateHitId(*getHitPtr(it), SimTrackIdsDC);
00557     //    cout<<SimTrackIdsDC.size()<<endl;
00558     if(SimTrackIdsDC.size()>1){
00559       //     cout<<(t->g4Track_end()-t->g4Track_begin())<<endl;
00560       for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin(); g4T !=  t -> g4Track_end(); ++g4T) {
00561         if(find(SimTrackIdsDC.begin(), SimTrackIdsDC.end(),SimHitIdpr((*g4T).trackId(), SimTrackIdsDC.begin()->second)) != SimTrackIdsDC.end() ){
00562           idcount++;
00563         }
00564       }
00565     }
00566     if (idcount>1) doublecount+=(idcount-1);
00567   }
00568   return doublecount;
00569 }