Go to the documentation of this file.00001 #ifndef RecoSelectors_TrackingParticleSelector_h
00002 #define RecoSelectors_TrackingParticleSelector_h
00003
00004
00005
00006
00007
00008
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;
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);
00034
00035
00036
00037 bool stable = 1;
00038 if (stableOnly_) {
00039 for( TrackingParticle::genp_iterator j = tp.genParticle_begin(); j != tp.genParticle_end(); ++ j ) {
00040 const HepMC::GenParticle * p = j->get();
00041 if (p->status() != 1) {
00042 stable = 0; break;
00043 }
00044 }
00045
00046 if(stable == 1 && tp.status() == -99 &&
00047 (fabs(tp.pdgId()) != 11 && fabs(tp.pdgId()) != 13 && fabs(tp.pdgId()) != 211 &&
00048 fabs(tp.pdgId()) != 321 && fabs(tp.pdgId()) != 2212 && fabs(tp.pdgId()) != 3112 &&
00049 fabs(tp.pdgId()) != 3222 && fabs(tp.pdgId()) != 3312 && fabs(tp.pdgId()) != 3334)) stable = 0;
00050 }
00051 return (
00052 tp.matchedHit() >= minHit_ &&
00053 sqrt(tp.momentum().perp2()) >= ptMin_ &&
00054 tp.momentum().eta() >= minRapidity_ && tp.momentum().eta() <= maxRapidity_ &&
00055 sqrt(tp.vertex().perp2()) <= tip_ &&
00056 fabs(tp.vertex().z()) <= lip_ &&
00057 testId &&
00058 signal &&
00059 stable
00060 );
00061 }
00062
00063 private:
00064 double ptMin_;
00065 double minRapidity_;
00066 double maxRapidity_;
00067 double tip_;
00068 double lip_;
00069 int minHit_;
00070 bool signalOnly_;
00071 bool chargedOnly_;
00072 bool stableOnly_;
00073 std::vector<int> pdgId_;
00074
00075 };
00076
00077 #include "CommonTools/UtilAlgos/interface/ParameterAdapter.h"
00078
00079 namespace reco {
00080 namespace modules {
00081
00082 template<>
00083 struct ParameterAdapter<TrackingParticleSelector> {
00084 static TrackingParticleSelector make( const edm::ParameterSet & cfg ) {
00085 return TrackingParticleSelector(
00086 cfg.getParameter<double>( "ptMin" ),
00087 cfg.getParameter<double>( "minRapidity" ),
00088 cfg.getParameter<double>( "maxRapidity" ),
00089 cfg.getParameter<double>( "tip" ),
00090 cfg.getParameter<double>( "lip" ),
00091 cfg.getParameter<int>( "minHit" ),
00092 cfg.getParameter<bool>( "signalOnly" ),
00093 cfg.getParameter<bool>( "chargedOnly" ),
00094 cfg.getParameter<bool>( "stableOnly" ),
00095 cfg.getParameter<std::vector<int> >( "pdgId" ));
00096 }
00097 };
00098
00099 }
00100 }
00101
00102 #endif