CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/FastSimulation/Muons/plugins/FastTSGFromIOHit.cc

Go to the documentation of this file.
00001 #include "FastSimulation/Muons/plugins/FastTSGFromIOHit.h"
00002 
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSMatchedRecHit2DCollection.h"
00014 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
00015 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00016 #include "DataFormats/Math/interface/deltaPhi.h"
00017 
00018 
00019 FastTSGFromIOHit::FastTSGFromIOHit(const edm::ParameterSet & iConfig) : theConfig (iConfig)
00020 {
00021   theCategory = "FastSimulation|Muons||FastTSGFromIOHit";
00022 
00023   thePtCut = iConfig.getParameter<double>("PtCut");
00024 
00025   theSeedCollectionLabels = iConfig.getParameter<std::vector<edm::InputTag> >("SeedCollectionLabels");
00026 
00027   theSimTrackCollectionLabel  = iConfig.getParameter<edm::InputTag>("SimTrackCollectionLabel");
00028 
00029 }
00030 
00031 FastTSGFromIOHit::~FastTSGFromIOHit()
00032 {
00033 
00034   LogTrace(theCategory) << " FastTSGFromIOHit dtor called ";
00035 
00036 }
00037 
00038 void FastTSGFromIOHit::trackerSeeds(const TrackCand& staMuon, const TrackingRegion& region, std::vector<TrajectorySeed> & result) {
00039   
00040   // Retrieve the Monte Carlo truth (SimTracks)
00041   edm::Handle<edm::SimTrackContainer> theSimTracks;
00042   getEvent()->getByLabel(theSimTrackCollectionLabel,theSimTracks);
00043     
00044   // Retrieve Seed collection
00045   unsigned seedCollections = theSeedCollectionLabels.size();
00046   std::vector<edm::Handle<edm::View<TrajectorySeed> > > theSeeds;
00047   theSeeds.resize(seedCollections);
00048   unsigned seed_size = 0;
00049   for ( unsigned iseed=0; iseed<seedCollections; ++iseed ) { 
00050     getEvent()->getByLabel(theSeedCollectionLabels[iseed], theSeeds[iseed]);
00051     seed_size += theSeeds[iseed]->size();
00052   }
00053   
00054   // Make a ref to l2 muon
00055   reco::TrackRef muRef(staMuon.second);
00056   
00057   // Cut on muons with low momenta
00058   if ( muRef->pt() < thePtCut 
00059        || muRef->innerMomentum().Rho() < thePtCut 
00060        || muRef->innerMomentum().R() < 2.5 ) return;
00061   
00062   // Copy the collection of seeds (ahem, this is time consuming!)
00063   std::vector<TrajectorySeed> tkSeeds;
00064   std::set<unsigned> tkIds;
00065   tkSeeds.reserve(seed_size);
00066   for ( unsigned iseed=0; iseed<seedCollections; ++iseed ) { 
00067     edm::Handle<edm::View<TrajectorySeed> > aSeedCollection = theSeeds[iseed];
00068     unsigned nSeeds = aSeedCollection->size();
00069     for (unsigned seednr = 0; seednr < nSeeds; ++seednr) {
00070       
00071       // The seed
00072       const BasicTrajectorySeed* aSeed = &((*aSeedCollection)[seednr]);
00073       
00074       // Find the first hit of the Seed
00075       TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
00076       const SiTrackerGSMatchedRecHit2D * theFirstSeedingRecHit = 
00077         (const SiTrackerGSMatchedRecHit2D*) (&(*(theSeedingRecHitRange.first)));
00078       
00079       // The SimTrack id associated to that recHit
00080       int simTrackId = theFirstSeedingRecHit->simtrackId();
00081       
00082       // Track already associated to a seed
00083       std::set<unsigned>::iterator tkId = tkIds.find(simTrackId);
00084       if( tkId != tkIds.end() ) continue;       
00085       
00086       const SimTrack& theSimTrack = (*theSimTracks)[simTrackId]; 
00087        
00088       const RectangularEtaPhiTrackingRegion & regionRef = dynamic_cast<const RectangularEtaPhiTrackingRegion & > (region);
00089   
00090       if( clean(muRef,regionRef,aSeed,theSimTrack) ) tkSeeds.push_back(*aSeed);
00091       tkIds.insert(simTrackId);
00092       
00093     } // End loop on seeds
00094  
00095   } // End loop on seed collections
00096   
00097   // Now create the Muon Trajectory Seed
00098   unsigned int is=0;
00099   unsigned int isMax=tkSeeds.size();
00100   for (;is!=isMax;++is){
00101     result.push_back( L3MuonTrajectorySeed(tkSeeds[is], muRef));
00102   } // End of tk seed loop
00103   
00104 }
00105 
00106 bool
00107 FastTSGFromIOHit::clean(reco::TrackRef muRef,
00108                         const RectangularEtaPhiTrackingRegion& region,
00109                         const BasicTrajectorySeed* aSeed,
00110                         const SimTrack& theSimTrack) 
00111 {
00112   //  return true; 
00113   //}
00114 
00115   // Eta cleaner   
00116   const PixelRecoRange<float>& etaRange = region.etaRange() ;  
00117   //  return true;
00118   //}
00119 
00120   double etaSeed = theSimTrack.momentum().Eta();
00121   double etaLimit  = (fabs(fabs(etaRange.max())-fabs(etaRange.mean())) <0.05) ? 
00122     0.05 : fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) ;
00123   bool inEtaRange = 
00124     etaSeed >= (etaRange.mean() - etaLimit) && 
00125     etaSeed <= (etaRange.mean() + etaLimit) ;
00126   if  ( !inEtaRange ) return false;
00127   
00128   // Phi cleaner
00129   const TkTrackingRegionsMargin<float>& phiMargin = region.phiMargin();
00130   double phiSeed = theSimTrack.momentum().Phi(); 
00131   double phiLimit  = (phiMargin.right() < 0.05 ) ? 0.05 : phiMargin.right(); 
00132   bool inPhiRange = 
00133     (fabs(deltaPhi(phiSeed,double(region.direction().phi()))) < phiLimit );
00134   if  ( !inPhiRange ) return false;
00135   
00136   // pt cleaner
00137   double ptSeed  = std::sqrt(theSimTrack.momentum().Perp2());
00138   double ptMin   = (region.ptMin()>3.5) ? 3.5: region.ptMin();  
00139   bool inPtRange = ptSeed >= ptMin &&  ptSeed<= 2*(muRef->pt());
00140   return inPtRange;
00141   
00142 }