CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/Calibration/IsolatedParticles/src/MatchingSimTrack.cc

Go to the documentation of this file.
00001 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
00002 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00003 
00004 #include "Calibration/IsolatedParticles/interface/MatchingSimTrack.h"
00005 
00006 #include <iostream>
00007 
00008 namespace spr{
00009 
00010   edm::SimTrackContainer::const_iterator matchedSimTrack(const edm::Event& iEvent, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, const reco::Track* pTrack, TrackerHitAssociator& associate, bool debug){
00011   
00012     edm::SimTrackContainer::const_iterator itr = SimTk->end();;
00013     
00014     edm::SimTrackContainer::const_iterator simTrkItr;
00015     edm::SimVertexContainer::const_iterator simVtxItr;
00016  
00017     //Get the vector of PsimHits associated to TrackerRecHits and select the 
00018     //matching SimTrack on the basis of maximum occurance of trackIds
00019     std::vector<unsigned int> trkId, trkOcc;
00020     int i=0;
00021     for (trackingRecHit_iterator iTrkHit = pTrack->recHitsBegin(); iTrkHit != pTrack->recHitsEnd(); ++iTrkHit, ++i) {
00022     
00023       std::vector<PSimHit> matchedSimIds = associate.associateHit((**iTrkHit));
00024       for (unsigned int isim=0; isim<matchedSimIds.size(); isim++) {
00025         unsigned tkId = matchedSimIds[isim].trackId();
00026         bool found = false;
00027         for (unsigned int j=0; j<trkId.size(); j++) {
00028           if ( tkId == trkId[j] ) {
00029             trkOcc[j]++;
00030             found = true;
00031             break;
00032           }
00033         }
00034         if (!found) { trkId.push_back(tkId); trkOcc.push_back(1); }
00035       }
00036     }
00037 
00038     if (debug) {
00039       std::cout << "Reconstructed Track with " << i << " recHits.";
00040       for (unsigned int isim=0; isim<trkId.size(); isim++){
00041         std::cout << "\n trkId " << trkId[isim] << "  Occurance " << trkOcc[isim] << ", ";
00042       } 
00043       std::cout << std::endl;
00044     }
00045   
00046     int matchedId=0; 
00047     unsigned int matchSimTrk=0;
00048     if (trkOcc.size() > 0) {
00049       unsigned int maxTrkOcc=0, idxMax=0;
00050       for(unsigned int j=0; j<trkOcc.size(); j++) {
00051         if(trkOcc[j] > maxTrkOcc ) { maxTrkOcc = trkOcc[j]; idxMax = j; }
00052       }
00053       matchSimTrk = trkId[idxMax];
00054       for (simTrkItr = SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++) {
00055         if ( simTrkItr->trackId() == matchSimTrk ) {
00056           matchedId = simTrkItr->type();
00057           if (debug) std::cout << "matched trackId (maximum occurance) " << matchSimTrk << " type " << matchedId << std::endl;
00058           itr = simTrkItr;
00059           break;
00060         }
00061       }
00062     }
00063   
00064     if (matchedId==0 && debug) {
00065       std::cout << "Could not find matched SimTrk and track history now " << std::endl;
00066     }
00067 
00068     return itr;
00069   }
00070 
00071   //Returns a vector of TrackId originating from the matching SimTrack
00072   std::vector<int> matchedSimTrackId(const edm::Event& iEvent, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, const reco::Track* pTrack, TrackerHitAssociator& associate, bool debug) {
00073 
00074     // get the matching SimTrack
00075     edm::SimTrackContainer::const_iterator trkInfo = spr::matchedSimTrack(iEvent, SimTk, SimVtx, pTrack,associate, debug);
00076     unsigned int matchSimTrk = trkInfo->trackId();
00077     if (debug) std::cout << "matchedSimTrackId finds the SimTrk ID of the current track to be " << matchSimTrk << std::endl;
00078   
00079     std::vector<int> matchTkid;
00080     if( trkInfo->type() != 0) {
00081       edm::SimTrackContainer::const_iterator simTrkItr;
00082       for(simTrkItr = SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++){
00083         if (validSimTrack(matchSimTrk, simTrkItr, SimTk, SimVtx, false))
00084           matchTkid.push_back((int)simTrkItr->trackId());
00085       }
00086     }
00087     return matchTkid;
00088   }
00089 
00090   spr::simTkInfo matchedSimTrackInfo(unsigned int simTkId, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, bool debug) {
00091     spr::simTkInfo info;
00092     for (edm::SimTrackContainer::const_iterator simTrkItr=SimTk->begin(); 
00093          simTrkItr!= SimTk->end(); simTrkItr++) {
00094       if (simTkId == simTrkItr->trackId()) {
00095         if (spr::validSimTrack(simTkId, simTrkItr, SimTk, SimVtx, debug)) {
00096           info.found  = true;
00097           info.pdgId  = simTrkItr->type();
00098           info.charge = simTrkItr->charge();
00099         } else {
00100           edm::SimTrackContainer::const_iterator parentItr =  spr::parentSimTrack(simTrkItr, SimTk, SimVtx, debug);
00101           if (debug) {
00102             if (parentItr != SimTk->end() ) std::cout << "original parent of " << simTrkItr->trackId() << " " << parentItr->trackId() << ", " << parentItr->type() << std::endl;
00103             else                            std::cout << "original parent of " << simTrkItr->trackId() << " not found" << std::endl;
00104           }
00105           if (parentItr != SimTk->end()) {
00106             info.found  = true;
00107             info.pdgId  = parentItr->type();
00108             info.charge = parentItr->charge();
00109           }
00110         }
00111         break;
00112       }
00113     }
00114     return info;
00115   }
00116 
00117   // Checks if this SimTrack=thisTrkItr originates from the one with trackId=simTkId
00118   bool validSimTrack(unsigned int simTkId, edm::SimTrackContainer::const_iterator thisTrkItr, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, bool debug){
00119     
00120     if (debug) std::cout << "Inside validSimTrack: trackId " << thisTrkItr->trackId() << " vtxIndex " <<  thisTrkItr->vertIndex() << " to be matched to " << simTkId << std::endl;
00121     
00122     //This track originates from simTkId
00123     if (thisTrkItr->trackId() == simTkId) return true;
00124 
00125     //Otherwise trace back the history using SimTracks and SimVertices
00126     int vertIndex = thisTrkItr->vertIndex();
00127     if (vertIndex == -1 || vertIndex >= (int)SimVtx->size()) return false;
00128     
00129     edm::SimVertexContainer::const_iterator simVtxItr= SimVtx->begin();
00130     for (int iv=0; iv<vertIndex; iv++) simVtxItr++;
00131     int parent = simVtxItr->parentIndex();
00132     if (debug) std::cout << "validSimTrack:: parent index " << parent <<" ";
00133     if (parent < 0 && simVtxItr != SimVtx->begin()) {
00134       const math::XYZTLorentzVectorD pos1 = simVtxItr->position();
00135       for (simVtxItr=SimVtx->begin(); simVtxItr!=SimVtx->end(); ++simVtxItr) {
00136         if (simVtxItr->parentIndex() > 0) {
00137           const math::XYZTLorentzVectorD pos2 = pos1 - simVtxItr->position();
00138           double dist = pos2.P();
00139           if (dist < 0.001) {
00140             parent = simVtxItr->parentIndex();
00141             break;
00142           }
00143         }
00144       }
00145     }
00146     if (debug) std::cout << "final index " << parent << std::endl;;
00147     for(edm::SimTrackContainer::const_iterator simTrkItr= SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++){
00148       if ((int)simTrkItr->trackId() == parent && simTrkItr != thisTrkItr) return  validSimTrack(simTkId, simTrkItr, SimTk, SimVtx, debug) ; 
00149     }
00150   
00151     return false;
00152   
00153   }
00154 
00155   //Returns the parent of a SimTrack
00156   edm::SimTrackContainer::const_iterator parentSimTrack(edm::SimTrackContainer::const_iterator thisTrkItr, edm::Handle<edm::SimTrackContainer>& SimTk, edm::Handle<edm::SimVertexContainer>& SimVtx, bool debug){
00157   
00158     edm::SimTrackContainer::const_iterator itr = SimTk->end();
00159 
00160     int vertIndex = thisTrkItr->vertIndex();
00161     int type = thisTrkItr->type(); int charge = (int)thisTrkItr->charge();
00162     if (debug) std::cout << "SimTrackParent " << thisTrkItr->trackId() << " Vertex " <<  vertIndex << " Type " << type << " Charge " << charge << std::endl;
00163 
00164     if( vertIndex == -1 )                      return thisTrkItr;  
00165     else if (vertIndex >= (int)SimVtx->size()) return itr; 
00166 
00167     edm::SimVertexContainer::const_iterator simVtxItr= SimVtx->begin();
00168     for (int iv=0; iv<vertIndex; iv++) simVtxItr++;
00169     int parent = simVtxItr->parentIndex();
00170 
00171     if (parent < 0 && simVtxItr != SimVtx->begin()) {
00172       const math::XYZTLorentzVectorD pos1 = simVtxItr->position();
00173       for (simVtxItr=SimVtx->begin(); simVtxItr!=SimVtx->end(); ++simVtxItr) {
00174         if (simVtxItr->parentIndex() > 0) {
00175           const math::XYZTLorentzVectorD pos2 = pos1 - simVtxItr->position();
00176           double dist = pos2.P();
00177           if (dist < 0.001) {
00178             parent = simVtxItr->parentIndex();
00179             break;
00180           }
00181         }
00182       }
00183     }
00184     for (edm::SimTrackContainer::const_iterator simTrkItr= SimTk->begin(); simTrkItr!= SimTk->end(); simTrkItr++){
00185       if ((int)simTrkItr->trackId() == parent && simTrkItr != thisTrkItr) return  parentSimTrack(simTrkItr, SimTk, SimVtx, debug); 
00186     }
00187 
00188     return thisTrkItr; 
00189   }
00190 
00191 
00192 }