CMS 3D CMS Logo

GsfSeedCleaner.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    PFTracking
00004 // Class:      GsfSeedCleaner
00005 // 
00006 // Original Author:  Michele Pioppi
00007 
00008 #include "RecoParticleFlow/PFTracking/interface/GsfSeedCleaner.h"
00009 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/ParameterSet/interface/FileInPath.h"
00012 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00013 
00014 using namespace edm;
00015 using namespace std;
00016 using namespace reco;
00017 
00018 
00019 GsfSeedCleaner::GsfSeedCleaner(const ParameterSet& iConfig):
00020   conf_(iConfig)
00021 {
00022   LogInfo("GsfSeedCleaner")<<"Seed cleaning for Gsf tracks started  ";
00023   
00024   //now do what ever initialization is needed
00025  
00026   tracksContainers_ = 
00027     iConfig.getParameter< vector < InputTag > >("TkColList");
00028  
00029 
00030 
00031   preIdLabel_ =
00032     iConfig.getParameter<InputTag>("PreIdSeedLabel");
00033  
00034   produces<TrajectorySeedCollection>();
00035 
00036 }
00037 
00038 
00039 
00040 //
00041 // member functions
00042 //
00043 
00044 // ------------ method called to produce the data  ------------
00045 void
00046 GsfSeedCleaner::produce(Event& iEvent, const EventSetup& iSetup)
00047 {
00048   
00049   LogDebug("GsfSeedCleaner")<<"START event: "<<iEvent.id().event()
00050                               <<" in run "<<iEvent.id().run();
00051   
00052   //Create empty output collections
00053   auto_ptr<TrajectorySeedCollection> output(new TrajectorySeedCollection);
00054 
00055   Handle<TrajectorySeedCollection> SeedCollection;
00056   iEvent.getByLabel(preIdLabel_,SeedCollection);
00057 
00058   TrajectorySeedCollection::const_iterator isc=SeedCollection->begin();
00059   TrajectorySeedCollection::const_iterator isc_end=SeedCollection->end();
00060 
00061   for (;isc!=isc_end;++isc){
00062     bool seed_not_used =true;
00063     //Vector of track collections
00064     for (uint istr=0; istr<tracksContainers_.size();istr++){
00065       
00066       //Track collection
00067       Handle<GsfElectronCollection> ElecCollection;
00068       iEvent.getByLabel(tracksContainers_[istr], ElecCollection);
00069 
00070 
00071       GsfElectronCollection::const_iterator itk = ElecCollection->begin();
00072       GsfElectronCollection::const_iterator itk_end = ElecCollection->end();
00073       for(;itk!=itk_end;++itk){
00074 
00075         if (seed_not_used){
00076           seed_not_used=CompareHits(*itk->gsfTrack(),*isc);
00077         }
00078       }    
00079       
00080     } //end loop on the vector of Gsf track collections
00081     if (seed_not_used)
00082       output->push_back(*isc);
00083   }//end loop on the preid seeds
00084 
00085   iEvent.put(output);
00086   
00087 }
00088 
00089 bool GsfSeedCleaner::CompareHits(const GsfTrack tk,const TrajectorySeed s){
00090 
00091  
00092   TrajectorySeed::const_iterator sh = s.recHits().first;
00093   TrajectorySeed::const_iterator sh_end = s.recHits().second;
00094   int hitinseed=0;
00095   int hitshared=0;
00096   for (;sh!=sh_end;++sh){
00097     if(!(sh->isValid())) continue;
00098     hitinseed++;
00099     trackingRecHit_iterator  ghit= tk.recHitsBegin();
00100     trackingRecHit_iterator  ghit_end= tk.recHitsEnd();
00101     for (;ghit!=ghit_end;++ghit){      
00102       if (!((*ghit)->isValid())) continue;
00103         if(((*sh).geographicalId()==(*ghit)->geographicalId())&&
00104            (((*sh).localPosition()-(*ghit)->localPosition()).mag()<0.01)) hitshared++;
00105     }
00106   }
00107 
00108   if ((hitinseed<3) &&(hitinseed==hitshared)) return false;
00109   if ((hitinseed>2) &&((hitinseed-hitshared)<2)) return false; 
00110   return true;
00111 }

Generated on Tue Jun 9 17:44:50 2009 for CMSSW by  doxygen 1.5.4