![]() |
![]() |
#include <PhysicsTools/JetMCAlgos/plugins/TauGenJetProducer.h>
Public Member Functions | |
virtual void | beginJob (const edm::EventSetup &c) |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
TauGenJetProducer (const edm::ParameterSet &) | |
~TauGenJetProducer () | |
Private Attributes | |
bool | includeNeutrinos_ |
if yes, neutrinos will be included, for debug purposes | |
edm::InputTag | inputTagGenParticles_ |
Input PFCandidates. | |
bool | verbose_ |
verbose ? |
Definition at line 23 of file TauGenJetProducer.h.
TauGenJetProducer::TauGenJetProducer | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 22 of file TauGenJetProducer.cc.
References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), includeNeutrinos_, inputTagGenParticles_, and verbose_.
00022 { 00023 00024 00025 inputTagGenParticles_ 00026 = iConfig.getParameter<InputTag>("GenParticles"); 00027 00028 includeNeutrinos_ 00029 = iConfig.getParameter<bool>("includeNeutrinos"); 00030 00031 verbose_ = 00032 iConfig.getUntrackedParameter<bool>("verbose",false); 00033 00034 00035 00036 produces<GenJetCollection>(); 00037 }
TauGenJetProducer::~TauGenJetProducer | ( | ) |
void TauGenJetProducer::beginJob | ( | const edm::EventSetup & | c | ) | [virtual] |
void TauGenJetProducer::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 48 of file TauGenJetProducer.cc.
References funct::abs(), edm::RefVector< C, T, F >::begin(), GenMuonPlsPt100GeV_cfg::cout, edm::RefVector< C, T, F >::end(), lat::endl(), err, Exception, GenParticlesHelper::findDescendents(), GenParticlesHelper::findParticles(), genParticles_cfi::genParticles, edm::Event::getByLabel(), includeNeutrinos_, inputTagGenParticles_, metsig::jet, edm::Event::put(), edm::refToPtr(), and verbose_.
00049 { 00050 00051 Handle<GenParticleCollection> genParticles; 00052 00053 bool found = iEvent.getByLabel( inputTagGenParticles_, genParticles); 00054 00055 if(!found ) { 00056 std::ostringstream err; 00057 err<<" cannot get collection: " 00058 <<inputTagGenParticles_<<std::endl; 00059 edm::LogError("TauGenJetProducer")<<err.str(); 00060 throw cms::Exception( "MissingProduct", err.str()); 00061 } 00062 00063 std::auto_ptr<GenJetCollection> 00064 pOutVisTaus(new GenJetCollection() ); 00065 00066 00067 using namespace GenParticlesHelper; 00068 00069 GenParticleRefVector allStatus2Taus; 00070 findParticles( *genParticles, 00071 allStatus2Taus, 15, 2); 00072 00073 for( IGR iTau=allStatus2Taus.begin(); iTau!=allStatus2Taus.end(); ++iTau) { 00074 00075 // look for all status 1 (stable) descendents 00076 GenParticleRefVector descendents; 00077 findDescendents( *iTau, descendents, 1); 00078 00079 // loop on descendents, and take all except neutrinos 00080 math::XYZTLorentzVector sumVisMom; 00081 Particle::Charge charge = 0; 00082 Jet::Constituents constituents; 00083 00084 if(verbose_) 00085 cout<<"tau "<<(*iTau)<<endl; 00086 00087 00088 for(IGR igr = descendents.begin(); 00089 igr!= descendents.end(); ++igr ) { 00090 00091 int absPdgId = abs( (*igr)->pdgId()); 00092 00093 // neutrinos 00094 if(!includeNeutrinos_ ) { 00095 if( absPdgId == 12 || 00096 absPdgId == 14 || 00097 absPdgId == 16 ) 00098 continue; 00099 } 00100 00101 if(verbose_) 00102 cout<<"\t"<<(*igr)<<endl; 00103 00104 charge += (*igr)->charge(); 00105 sumVisMom += (*igr)->p4(); 00106 00107 // need to convert the vector of reference to the constituents 00108 // to a vector of pointers to build the genjet 00109 constituents.push_back( refToPtr( *igr) ); 00110 } 00111 00112 math::XYZPoint vertex; 00113 GenJet::Specific specific; 00114 00115 GenJet jet( sumVisMom, vertex, specific, constituents); 00116 pOutVisTaus->push_back( jet ); 00117 00118 } 00119 iEvent.put( pOutVisTaus ); 00120 }
bool TauGenJetProducer::includeNeutrinos_ [private] |
if yes, neutrinos will be included, for debug purposes
Definition at line 40 of file TauGenJetProducer.h.
Referenced by produce(), and TauGenJetProducer().
Input PFCandidates.
Definition at line 37 of file TauGenJetProducer.h.
Referenced by produce(), and TauGenJetProducer().
bool TauGenJetProducer::verbose_ [private] |
verbose ?
Definition at line 43 of file TauGenJetProducer.h.
Referenced by produce(), and TauGenJetProducer().