ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/Validation/plugins/HiBasicGenTest.cc
Revision: 1.1
Committed: Tue Mar 23 10:10:23 2010 UTC (15 years, 1 month ago) by edwenger
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Log Message:
migrated HI validation code to new directory

File Contents

# Content
1 #include "HepMC/GenEvent.h"
2 #include "HepMC/GenParticle.h"
3 #include "CmsHi/Validation/plugins/HiBasicGenTest.h"
4 #include "DataFormats/Candidate/interface/Candidate.h"
5 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
6 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
7
8 #include "HepMC/GenEvent.h"
9 #include "HepMC/HeavyIon.h"
10
11 #include <TString.h>
12 #include <TMath.h>
13
14 using namespace edm;
15 using namespace HepMC;
16
17 HiBasicGenTest::HiBasicGenTest(const edm::ParameterSet& iPSet)
18 {
19 dbe = 0;
20 dbe = edm::Service<DQMStore>().operator->();
21 }
22
23 HiBasicGenTest::~HiBasicGenTest() {}
24
25 void HiBasicGenTest::beginJob()
26 {
27 if(dbe){
28 ///Setting the DQM top directories
29 dbe->setCurrentFolder("Generator/Particles");
30
31 ///Booking the ME's
32 for(int ibin=0; ibin<3; ibin++) {
33 dnchdeta[ibin] = dbe->book1D(Form("dnchdeta%d",ibin), ";#eta;dN^{ch}/d#eta", 100, -6.0, 6.0);
34 dnchdpt[ibin] = dbe->book1D(Form("dnchdpt%d",ibin), ";p_{T};dN^{ch}/dp_{T}", 200, 0.0, 100.0);
35 b[ibin] = dbe->book1D(Form("b%d",ibin),";b[fm];events",100, 0.0, 20.0);
36 dnchdphi[ibin] = dbe->book1D(Form("dnchdphi%d",ibin),";#phi;dN^{ch}/d#phi",100, -3.2, 3.2);
37
38 dbe->tag(dnchdeta[ibin]->getFullname(),1+ibin*4);
39 dbe->tag(dnchdpt[ibin]->getFullname(),2+ibin*4);
40 dbe->tag(b[ibin]->getFullname(),3+ibin*4);
41 dbe->tag(dnchdphi[ibin]->getFullname(),4+ibin*4);
42 }
43
44 rp = dbe->book1D("phi0",";#phi_{RP};events",100,-3.2,3.2);
45 dbe->tag(rp->getFullname(),13);
46
47
48 }
49 return;
50 }
51
52 void HiBasicGenTest::endJob(){
53 // normalization of histograms can be done here (or in post-processor)
54 return;
55 }
56
57 void HiBasicGenTest::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
58 {
59 iSetup.getData(pdt);
60 return;
61 }
62
63 void HiBasicGenTest::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
64
65 void HiBasicGenTest::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
66 {
67
68 Handle<HepMCProduct> mc;
69 iEvent.getByLabel("generator",mc);
70 const HepMC::GenEvent *evt = mc->GetEvent();
71 const HepMC::HeavyIon *hi = evt->heavy_ion();
72
73 double ip = hi->impact_parameter();
74 double phi0 = hi->event_plane_angle();
75
76 // fill reaction plane distribution
77 rp->Fill(phi0);
78
79 // if the event is in one of the centrality bins of interest fill hists
80 int cbin=-1;
81 if(ip < 5.045) cbin=0;
82 else if (ip < 7.145 && ip > 5.045) cbin=1;
83 else if (ip < 15.202 && ip > 14.283) cbin=2;
84 if(cbin<0) return;
85
86 // fill impact parameter distributions
87 b[cbin]->Fill(ip);
88
89 // loop over particles
90 HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
91 HepMC::GenEvent::particle_const_iterator end = evt->particles_end();
92 for(HepMC::GenEvent::particle_const_iterator it = begin; it != end; ++it){
93
94 // only fill hists for status=1 particles
95 if((*it)->status() != 1) continue;
96
97 // only fill hists for charged particles
98 int pdg_id = (*it)->pdg_id();
99 const ParticleData * part = pdt->particle(pdg_id);
100 int charge = static_cast<int>(part->charge());
101 if(charge==0) continue;
102
103 float eta = (*it)->momentum().eta();
104 float phi = (*it)->momentum().phi();
105 float pt = (*it)->momentum().perp();
106
107 dnchdeta[cbin]->Fill(eta);
108 dnchdpt[cbin]->Fill(pt);
109
110 double pi = TMath::Pi();
111 double p = phi-phi0;
112 if(p > pi) p = p - 2*pi;
113 if(p < -1*pi) p = p + 2*pi;
114 dnchdphi[cbin]->Fill(p);
115
116 }
117
118 return;
119
120 }
121
122 #include "FWCore/Framework/interface/MakerMacros.h"
123 DEFINE_FWK_MODULE(HiBasicGenTest);
124
125