CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/CommonTools/RecoAlgos/interface/TrackingParticleSelector.h

Go to the documentation of this file.
00001 #ifndef RecoSelectors_TrackingParticleSelector_h
00002 #define RecoSelectors_TrackingParticleSelector_h
00003 /* \class TrackingParticleSelector
00004  *
00005  * \author Giuseppe Cerati, INFN
00006  *
00007  *  $Date: 2009/03/04 13:11:28 $
00008  *  $Revision: 1.1 $
00009  *
00010  */
00011 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00012 
00013 class TrackingParticleSelector {
00014 
00015 public:
00016   TrackingParticleSelector(){}
00017   TrackingParticleSelector ( double ptMin,double minRapidity,double maxRapidity,
00018                              double tip,double lip,int minHit, bool signalOnly, bool chargedOnly, 
00019                              std::vector<int> pdgId = std::vector<int>()) :
00020     ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
00021     tip_( tip ), lip_( lip ), minHit_( minHit ), signalOnly_(signalOnly), chargedOnly_(chargedOnly), pdgId_( pdgId ) { }
00022   
00024   bool operator()( const TrackingParticle & tp ) const { 
00025 
00026     //quickly reject if it is from pile-up
00027     if (signalOnly_ && !(tp.eventId().bunchCrossing()==0 && tp.eventId().event()==0) )return false;
00028 
00029     if (chargedOnly_ && tp.charge()==0) return false;//select only if charge!=0
00030     bool testId = false;
00031     unsigned int idSize = pdgId_.size();
00032     if (idSize==0) testId = true;
00033     else for (unsigned int it=0;it!=idSize;++it){
00034       if (tp.pdgId()==pdgId_[it]) testId = true;
00035     }
00036     return (
00037             tp.matchedHit() >= minHit_ &&
00038             sqrt(tp.momentum().perp2()) >= ptMin_ && 
00039             tp.momentum().eta() >= minRapidity_ && tp.momentum().eta() <= maxRapidity_ && 
00040             sqrt(tp.vertex().perp2()) <= tip_ &&
00041             fabs(tp.vertex().z()) <= lip_ &&
00042             testId
00043             );
00044   }
00045   
00046 private:
00047   double ptMin_;
00048   double minRapidity_;
00049   double maxRapidity_;
00050   double tip_;
00051   double lip_;
00052   int    minHit_;
00053   bool signalOnly_;
00054   bool chargedOnly_;
00055   std::vector<int> pdgId_;
00056 
00057 };
00058 
00059 #include "CommonTools/UtilAlgos/interface/ParameterAdapter.h"
00060 
00061 namespace reco {
00062   namespace modules {
00063     
00064     template<>
00065     struct ParameterAdapter<TrackingParticleSelector> {
00066       static TrackingParticleSelector make( const edm::ParameterSet & cfg ) {
00067         return TrackingParticleSelector(    
00068           cfg.getParameter<double>( "ptMin" ),
00069           cfg.getParameter<double>( "minRapidity" ),
00070           cfg.getParameter<double>( "maxRapidity" ),
00071           cfg.getParameter<double>( "tip" ),
00072           cfg.getParameter<double>( "lip" ),
00073           cfg.getParameter<int>( "minHit" ), 
00074           cfg.getParameter<bool>( "signalOnly" ),
00075           cfg.getParameter<bool>( "chargedOnly" ),
00076         cfg.getParameter<std::vector<int> >( "pdgId" )); 
00077       }
00078     };
00079     
00080   }
00081 }
00082 
00083 #endif