CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/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: 2012/05/17 16:12:35 $
00008  *  $Revision: 1.4 $
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, bool stableOnly,
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), stableOnly_(stableOnly), pdgId_( pdgId ) { }
00022   
00024   bool operator()( const TrackingParticle & tp ) const { 
00025     if (chargedOnly_ && tp.charge()==0) return false;//select only if charge!=0
00026     bool testId = false;
00027     unsigned int idSize = pdgId_.size();
00028     if (idSize==0) testId = true;
00029     else for (unsigned int it=0;it!=idSize;++it){
00030       if (tp.pdgId()==pdgId_[it]) testId = true;
00031     }
00032     bool signal = true;
00033     if (signalOnly_) signal = (tp.eventId().bunchCrossing()== 0 && tp.eventId().event() == 0); // signal only means no PU particles
00034     // select only stable particles
00035     bool stable = true;
00036     if (stableOnly_) {
00037       if (!signal) {
00038         stable = false; // we are not interested into PU particles among the stable ones
00039       } else {
00040         for( TrackingParticle::genp_iterator j = tp.genParticle_begin(); j != tp.genParticle_end(); ++ j ) {
00041           const HepMC::GenParticle * p = j->get();
00042           if (!p || p->status() != 1) {
00043             stable = 0; break;
00044           }
00045         }
00046        // test for remaining unstabled due to lack of genparticle pointer
00047        if(stable == 1 && tp.status() == -99 && 
00048           (fabs(tp.pdgId()) != 11 && fabs(tp.pdgId()) != 13 && fabs(tp.pdgId()) != 211 &&
00049            fabs(tp.pdgId()) != 321 && fabs(tp.pdgId()) != 2212 && fabs(tp.pdgId()) != 3112 &&
00050            fabs(tp.pdgId()) != 3222 && fabs(tp.pdgId()) != 3312 && fabs(tp.pdgId()) != 3334)) stable = 0;
00051       }
00052     }
00053     return (
00054             tp.matchedHit() >= minHit_ &&
00055             sqrt(tp.momentum().perp2()) >= ptMin_ && 
00056             tp.momentum().eta() >= minRapidity_ && tp.momentum().eta() <= maxRapidity_ && 
00057             sqrt(tp.vertex().perp2()) <= tip_ &&
00058             fabs(tp.vertex().z()) <= lip_ &&
00059             testId &&
00060             signal &&
00061             stable
00062             );
00063   }
00064   
00065 private:
00066   double ptMin_;
00067   double minRapidity_;
00068   double maxRapidity_;
00069   double tip_;
00070   double lip_;
00071   int    minHit_;
00072   bool signalOnly_;
00073   bool chargedOnly_;
00074   bool stableOnly_;
00075   std::vector<int> pdgId_;
00076 
00077 };
00078 
00079 #include "CommonTools/UtilAlgos/interface/ParameterAdapter.h"
00080 
00081 namespace reco {
00082   namespace modules {
00083     
00084     template<>
00085     struct ParameterAdapter<TrackingParticleSelector> {
00086       static TrackingParticleSelector make( const edm::ParameterSet & cfg ) {
00087         return TrackingParticleSelector(    
00088           cfg.getParameter<double>( "ptMin" ),
00089           cfg.getParameter<double>( "minRapidity" ),
00090           cfg.getParameter<double>( "maxRapidity" ),
00091           cfg.getParameter<double>( "tip" ),
00092           cfg.getParameter<double>( "lip" ),
00093           cfg.getParameter<int>( "minHit" ), 
00094           cfg.getParameter<bool>( "signalOnly" ),
00095           cfg.getParameter<bool>( "chargedOnly" ),
00096           cfg.getParameter<bool>( "stableOnly" ),
00097         cfg.getParameter<std::vector<int> >( "pdgId" )); 
00098       }
00099     };
00100     
00101   }
00102 }
00103 
00104 #endif