00001 // 00002 // Select the partons status 2 and 3 for MC Jet Flavour 00003 // Author: Attilio 00004 // Date: 10.10.2007 00005 // 00006 00007 //======================================================================= 00008 00009 // user include files 00010 #include "FWCore/Framework/interface/EDProducer.h" 00011 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00012 #include "FWCore/ParameterSet/interface/InputTag.h" 00013 00014 #include "FWCore/Framework/interface/Event.h" 00015 #include "FWCore/Framework/interface/EventSetup.h" 00016 #include "FWCore/Framework/interface/MakerMacros.h" 00017 #include "FWCore/Framework/interface/ESHandle.h" 00018 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00019 00020 #include "DataFormats/Common/interface/Handle.h" 00021 #include "DataFormats/Common/interface/Ref.h" 00022 //#include "DataFormats/Candidate/interface/CandidateFwd.h" 00023 //#include "DataFormats/HepMCCandidate/interface/GenParticleCandidate.h" 00024 #include "DataFormats/HepMCCandidate/interface/GenParticle.h" 00025 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h" 00026 00027 #include <memory> 00028 #include <string> 00029 #include <iostream> 00030 #include <vector> 00031 00032 using namespace std; 00033 using namespace reco; 00034 using namespace edm; 00035 00036 class PartonSelector : public edm::EDProducer 00037 { 00038 public: 00039 PartonSelector( const edm::ParameterSet & ); 00040 ~PartonSelector(); 00041 00042 private: 00043 00044 virtual void produce(edm::Event&, const edm::EventSetup& ); 00045 bool withLeptons; // Optionally specify leptons 00046 bool withTop; // Optionally include top quarks in the list 00047 }; 00048 //========================================================================= 00049 00050 PartonSelector::PartonSelector( const edm::ParameterSet& iConfig ) 00051 { 00052 produces<reco::GenParticleRefVector>(); 00053 withLeptons = iConfig.getParameter<bool>("withLeptons"); 00054 if ( iConfig.exists("withTop") ) { 00055 withTop = iConfig.getParameter<bool>("withTop"); 00056 } else { 00057 withTop = false; 00058 } 00059 } 00060 00061 //========================================================================= 00062 00063 PartonSelector::~PartonSelector() 00064 { 00065 } 00066 00067 // ------------ method called to produce the data ------------ 00068 00069 void PartonSelector::produce( Event& iEvent, const EventSetup& iEs ) 00070 { 00071 00072 //edm::Handle <reco::CandidateView> particles; 00073 edm::Handle <reco::GenParticleCollection> particles; 00074 iEvent.getByLabel ("genParticles", particles ); 00075 edm::LogVerbatim("PartonSelector") << "=== GenParticle size:" << particles->size(); 00076 int nPart=0; 00077 00078 auto_ptr<GenParticleRefVector> thePartons ( new GenParticleRefVector); 00079 00080 for (size_t m = 0; m < particles->size(); m++) { 00081 00082 // Don't take into account first 6 particles in generator list 00083 if (m<6) continue; 00084 00085 const GenParticle & aParticle = (*particles)[ m ]; 00086 00087 bool isAParton = false; 00088 bool isALepton = false; 00089 int flavour = abs(aParticle.pdgId()); 00090 if(flavour == 1 || 00091 flavour == 2 || 00092 flavour == 3 || 00093 flavour == 4 || 00094 flavour == 5 || 00095 (flavour == 6 && withTop) || 00096 flavour == 21 ) isAParton = true; 00097 if(flavour == 11 || 00098 flavour == 12 || 00099 flavour == 13 || 00100 flavour == 14 || 00101 flavour == 15 || 00102 flavour == 16 ) isALepton = true; 00103 00104 00105 //Add Partons status 3 00106 if( aParticle.status() == 3 && isAParton ) { 00107 thePartons->push_back( GenParticleRef( particles, m ) ); 00108 nPart++; 00109 } 00110 00111 //Add Partons stauts 2 00112 if( aParticle.numberOfDaughters() > 0 && 00113 ( aParticle.daughter(0)->pdgId() == 91 || aParticle.daughter(0)->pdgId() == 92 ) ) { 00114 thePartons->push_back( GenParticleRef( particles, m ) ); 00115 nPart++; 00116 } 00117 00118 //Add Leptons 00119 // Here you have to decide what to do with taus.... 00120 // Now all leptons, including e and mu from leptonic tau decays, are added 00121 if( withLeptons && aParticle.status() == 3 && isALepton ) { 00122 thePartons->push_back( GenParticleRef( particles, m ) ); 00123 nPart++; 00124 } 00125 } 00126 00127 edm::LogVerbatim("PartonSelector") << "=== GenParticle selected:" << nPart; 00128 iEvent.put( thePartons ); 00129 00130 } 00131 00132 DEFINE_FWK_MODULE(PartonSelector); 00133