CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/HiGenCommon/src/HighMultiplicityGenFilter.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    HighMultiplicityGenFilter
00004 // Class:      HighMultiplicityGenFilter
00005 // 
00013 //
00014 // Original Author:  Wei Li
00015 //         Created:  Tue Dec  8 23:51:37 EST 2009
00016 // $Id: HighMultiplicityGenFilter.cc,v 1.1 2010/08/31 21:39:08 edwenger Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDFilter.h"
00027 #include "FWCore/Framework/interface/EventSetup.h"
00028 #include "FWCore/Framework/interface/ESHandle.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00033 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00034 //
00035 // class declaration
00036 //
00037 
00038 class HighMultiplicityGenFilter : public edm::EDFilter {
00039    public:
00040       explicit HighMultiplicityGenFilter(const edm::ParameterSet&);
00041       ~HighMultiplicityGenFilter();
00042 
00043    private:
00044       virtual void beginJob();
00045       virtual bool filter(edm::Event&, const edm::EventSetup&);
00046       virtual void endJob();
00047       
00048       // ----------member data ---------------------------
00049       edm::ESHandle <ParticleDataTable> pdt; 
00050       double etaMax;
00051       double ptMin;
00052       int nMin; 
00053       int nAccepted;
00054 };
00055 
00056 //
00057 // constants, enums and typedefs
00058 //
00059 
00060 //
00061 // static data member definitions
00062 //
00063 
00064 //
00065 // constructors and destructor
00066 //
00067 HighMultiplicityGenFilter::HighMultiplicityGenFilter(const edm::ParameterSet& iConfig) :
00068 etaMax(iConfig.getUntrackedParameter<double>("etaMax")),
00069 ptMin(iConfig.getUntrackedParameter<double>("ptMin")),
00070 nMin(iConfig.getUntrackedParameter<int>("nMin"))
00071 {
00072   //now do what ever initialization is needed
00073   nAccepted = 0; 
00074 }
00075 
00076 
00077 HighMultiplicityGenFilter::~HighMultiplicityGenFilter()
00078 {
00079  
00080    // do anything here that needs to be done at desctruction time
00081    // (e.g. close files, deallocate resources etc.)
00082 
00083 }
00084 
00085 
00086 //
00087 // member functions
00088 //
00089 
00090 // ------------ method called on each new Event  ------------
00091 bool
00092 HighMultiplicityGenFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00093 {
00094 
00095   bool accepted = false;
00096   edm::Handle<edm::HepMCProduct> evt;
00097   iEvent.getByLabel("generator", evt);
00098 
00099   iSetup.getData(pdt);
00100 
00101   const HepMC::GenEvent * myGenEvent = evt->GetEvent();
00102 
00103   int nMult=0;
00104   for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();   p != myGenEvent->particles_end(); ++p ) {
00105 
00106     if((*p)->status()!=1) continue;
00107 
00108     double charge = 0;
00109     int pid = (*p)->pdg_id();
00110     if(abs(pid) > 100000) { std::cout<<"pid="<<pid<<" status="<<(*p)->status()<<std::endl; continue; }
00111     const ParticleData* part = pdt->particle(pid);
00112     if(part) charge = part->charge();
00113     if(charge == 0) continue;
00114 
00115     if ( 
00116          (*p)->momentum().perp() > ptMin 
00117          && fabs((*p)->momentum().eta()) < etaMax  ) nMult++;
00118   }
00119   if(nMult>=nMin) { nAccepted++; accepted = true; }
00120   return accepted;
00121 }
00122 
00123 // ------------ method called once each job just before starting event loop  ------------
00124 void 
00125 HighMultiplicityGenFilter::beginJob()
00126 {}
00127 
00128 // ------------ method called once each job just after ending the event loop  ------------
00129 void 
00130 HighMultiplicityGenFilter::endJob() {
00131   std::cout<<"There are "<<nAccepted<<" events with multiplicity greater than "<<nMin<<std::endl;
00132 }
00133 
00134 //define this as a plug-in
00135 DEFINE_FWK_MODULE(HighMultiplicityGenFilter);