CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/GenFilters/src/GenLeadTrackFilter.cc

Go to the documentation of this file.
00001 /***********************************************************
00002 *                 GenLeadTrackFilter                       *
00003 *                 ------------------                       *
00004 *                                                          *
00005 * Original Author: Souvik Das, Cornell University          *
00006 * Created        : 7 August 2009                           *
00007 *                                                          *
00008 *  Allows events which have at least one generator level   *
00009 *  charged particle with pT greater than X GeV within      *
00010 *  |eta| less than Y, where X and Y are specified in the   *
00011 *  cfi configuration file.                                 *
00012 ***********************************************************/
00013 
00014 #include "GeneratorInterface/GenFilters/interface/GenLeadTrackFilter.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 
00017 using namespace edm;
00018 using namespace reco;
00019 using namespace math;
00020 
00021 GenLeadTrackFilter::GenLeadTrackFilter(const edm::ParameterSet& iConfig)
00022 {
00023   hepMCProduct_label_   = iConfig.getParameter<InputTag>("HepMCProduct");
00024   genLeadTrackPt_       = iConfig.getParameter<double>("GenLeadTrackPt");
00025   genEta_               = iConfig.getParameter<double>("GenEta");
00026 }
00027 
00028 GenLeadTrackFilter::~GenLeadTrackFilter()
00029 {
00030 }
00031 
00032 bool GenLeadTrackFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00033 {
00034   bool allow=false;
00035   
00036   ESHandle<HepPDT::ParticleDataTable> pdt;
00037   iSetup.getData(pdt);
00038 
00039   Handle<HepMCProduct> genEvent;
00040   iEvent.getByLabel(hepMCProduct_label_, genEvent);
00041   if (genEvent.isValid())
00042   {
00043     float genLeadTrackPt=-100;
00044     for (HepMC::GenEvent::particle_const_iterator iter=(*(genEvent->GetEvent())).particles_begin();
00045          iter!=(*(genEvent->GetEvent())).particles_end(); 
00046          ++iter)
00047     {
00048       HepMC::GenParticle* theParticle=*iter;
00049       double pt=pow(pow(theParticle->momentum().px(),2)+pow(theParticle->momentum().py(),2), 0.5);
00050       double charge=pdt->particle(theParticle->pdg_id())->charge();
00051       if (theParticle->status()==1 &&
00052           charge!=0 &&
00053           fabs(theParticle->momentum().eta())<genEta_ &&
00054           pt>genLeadTrackPt)
00055       {
00056         genLeadTrackPt=pt;
00057       }
00058     }
00059     if (genLeadTrackPt>genLeadTrackPt_) allow=true; else allow=false;
00060   }
00061   else 
00062   {
00063     std::cout<<"genEvent in not valid!"<<std::endl;
00064   }
00065   return allow;
00066 }
00067 
00068 // ------------ method called once each job just before starting event loop  ------------
00069 
00070 // ------------ method called once each job just after ending the event loop  ------------
00071 void GenLeadTrackFilter::endJob() 
00072 {
00073 }
00074 
00075 //define this as a plug-in
00076 DEFINE_FWK_MODULE(GenLeadTrackFilter);