CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/Validation/RecoHI/plugins/HiBasicGenTest.cc

Go to the documentation of this file.
00001 #include "HepMC/GenEvent.h"
00002 #include "HepMC/GenParticle.h"
00003 #include "Validation/RecoHI/plugins/HiBasicGenTest.h"
00004 #include "DataFormats/Candidate/interface/Candidate.h"
00005 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00006 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00007 
00008 #include "HepMC/GenEvent.h"
00009 #include "HepMC/HeavyIon.h"
00010 
00011 #include <TString.h>
00012 #include <TMath.h>
00013 
00014 using namespace edm;
00015 using namespace HepMC;
00016 
00017 HiBasicGenTest::HiBasicGenTest(const edm::ParameterSet& iPSet)
00018 {
00019   dbe = 0;
00020   dbe = edm::Service<DQMStore>().operator->();
00021 }
00022 
00023 HiBasicGenTest::~HiBasicGenTest() {}
00024 
00025 void HiBasicGenTest::beginJob()
00026 {
00027   if(dbe){
00029     dbe->setCurrentFolder("Generator/Particles");
00030     
00032     for(int ibin=0; ibin<3; ibin++) {
00033       dnchdeta[ibin] = dbe->book1D(Form("dnchdeta%d",ibin), ";#eta;dN^{ch}/d#eta", 100, -6.0, 6.0);
00034       dnchdpt[ibin] = dbe->book1D(Form("dnchdpt%d",ibin), ";p_{T};dN^{ch}/dp_{T}", 200, 0.0, 100.0);
00035       b[ibin] = dbe->book1D(Form("b%d",ibin),";b[fm];events",100, 0.0, 20.0);
00036       dnchdphi[ibin] = dbe->book1D(Form("dnchdphi%d",ibin),";#phi;dN^{ch}/d#phi",100, -3.2, 3.2);
00037 
00038       dbe->tag(dnchdeta[ibin]->getFullname(),1+ibin*4);
00039       dbe->tag(dnchdpt[ibin]->getFullname(),2+ibin*4);
00040       dbe->tag(b[ibin]->getFullname(),3+ibin*4);
00041       dbe->tag(dnchdphi[ibin]->getFullname(),4+ibin*4);
00042     }
00043 
00044     rp = dbe->book1D("phi0",";#phi_{RP};events",100,-3.2,3.2);
00045     dbe->tag(rp->getFullname(),13);
00046 
00047 
00048  }
00049   return;
00050 }
00051 
00052 void HiBasicGenTest::endJob(){
00053   // normalization of histograms can be done here (or in post-processor)
00054   return;
00055 }
00056 
00057 void HiBasicGenTest::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00058 {
00059   iSetup.getData(pdt);
00060   return;
00061 }
00062 
00063 void HiBasicGenTest::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00064 
00065 void HiBasicGenTest::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00066 { 
00067 
00068   Handle<HepMCProduct> mc;
00069   iEvent.getByLabel("generator",mc);
00070   const HepMC::GenEvent *evt = mc->GetEvent();
00071   const HepMC::HeavyIon *hi = evt->heavy_ion();
00072   
00073   double ip = hi->impact_parameter();
00074   double phi0 = hi->event_plane_angle();
00075 
00076   // fill reaction plane distribution
00077   rp->Fill(phi0);
00078   
00079   // if the event is in one of the centrality bins of interest fill hists
00080   int cbin=-1;
00081   if(ip < 5.045) cbin=0;
00082   else if (ip < 7.145 && ip > 5.045) cbin=1;
00083   else if (ip < 15.202 && ip > 14.283) cbin=2;
00084   if(cbin<0) return;
00085 
00086   // fill impact parameter distributions
00087   b[cbin]->Fill(ip);
00088 
00089   // loop over particles
00090   HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
00091   HepMC::GenEvent::particle_const_iterator end = evt->particles_end();
00092   for(HepMC::GenEvent::particle_const_iterator it = begin; it != end; ++it){
00093     
00094     // only fill hists for status=1 particles
00095     if((*it)->status() != 1) continue;
00096 
00097     // only fill hists for charged particles
00098     int pdg_id = (*it)->pdg_id();
00099     const ParticleData * part = pdt->particle(pdg_id);
00100     int charge = static_cast<int>(part->charge());
00101     if(charge==0) continue;
00102     
00103     float eta = (*it)->momentum().eta();
00104     float phi = (*it)->momentum().phi();
00105     float pt = (*it)->momentum().perp();
00106     
00107     dnchdeta[cbin]->Fill(eta);
00108     dnchdpt[cbin]->Fill(pt);
00109     
00110     double pi = TMath::Pi();
00111     double p = phi-phi0;
00112     if(p > pi) p = p - 2*pi;
00113     if(p < -1*pi) p = p + 2*pi;
00114     dnchdphi[cbin]->Fill(p);
00115 
00116   } 
00117 
00118   return;
00119 
00120 }
00121 
00122 #include "FWCore/Framework/interface/MakerMacros.h"
00123 DEFINE_FWK_MODULE(HiBasicGenTest);
00124 
00125