CMS 3D CMS Logo

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

Generated on Tue Jun 9 17:47:55 2009 for CMSSW by  doxygen 1.5.4