CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MCVerticesAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MCVerticesAnalyzer
4 // Class: MCVerticesAnalyzer
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Thu Dec 16 16:32:56 CEST 2010
16 // $Id: MCVerticesAnalyzer.cc,v 1.6 2011/11/12 16:49:19 venturia Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 #include <numeric>
24 
25 #include <vector>
26 
27 // user include files
30 
34 
36 
38 
41 
43 
47 
48 
49 #include "TH1F.h"
50 #include "TH2F.h"
51 #include "TProfile.h"
52 
53 //
54 // class decleration
55 //
56 
57 
59 public:
60  explicit MCVerticesAnalyzer(const edm::ParameterSet&);
62 
63 private:
64  virtual void beginJob() ;
65  virtual void beginRun(const edm::Run&, const edm::EventSetup&);
66  virtual void endRun(const edm::Run&, const edm::EventSetup&);
67  virtual void analyze(const edm::Event&, const edm::EventSetup&);
68  virtual void endJob() ;
69 
70  // ----------member data ---------------------------
71 
72 
73 
76  const bool m_useweight;
78 
79  TH1F* m_hnvtx;
80  TH1F* m_hlumi;
83  TProfile* m_hnvtxweightprof;
84  TH1F* m_hmainvtxx;
85  TH1F* m_hmainvtxy;
86  TH1F* m_hmainvtxz;
88 
89 };
90 
91 //
92 // constants, enums and typedefs
93 //
94 
95 //
96 // static data member definitions
97 //
98 
99 //
100 // constructors and destructor
101 //
103  m_pileupcollection(iConfig.getParameter<edm::InputTag>("pileupSummaryCollection")),
104  m_mctruthcollection(iConfig.getParameter<edm::InputTag>("mcTruthCollection")),
105  m_useweight(iConfig.getParameter<bool>("useWeight")),
106  m_weight(iConfig.getParameter<edm::InputTag>("weightProduct"))
107 
108 
109 {
110  //now do what ever initialization is needed
111 
112 
113 
115 
116  m_hnvtx = tfserv->make<TH1F>("nvtx","Number of pileup vertices",60,-0.5,59.5);
117  m_hnvtx->GetXaxis()->SetTitle("Number of Interactions");
118 
119  m_hlumi = tfserv->make<TH1F>("lumi","BX luminosity*xsect",200,0.,50.);
120  m_hlumi->GetXaxis()->SetTitle("Average Number of Interactions");
121 
122  m_hnvtxvslumi = tfserv->make<TH2F>("nvtxvslumi","Npileup vs BX luminosity*xsect",200,0.,50.,60,-0.5,59.5);
123  m_hnvtxvslumi->GetXaxis()->SetTitle("Average Number of Interactions"); m_hnvtxvslumi->GetYaxis()->SetTitle("Number of Interactions");
124 
125  if(m_useweight) {
126  m_hnvtxweight = tfserv->make<TH1F>("nvtxweight","Number of pileup vertices (1-w)",60,-0.5,59.5);
127  m_hnvtxweight->GetXaxis()->SetTitle("Number of Interactions");
128  m_hnvtxweightprof = tfserv->make<TProfile>("nvtxweightprof","Mean (1-w) vs Number of pileup interactions",60,-0.5,59.5);
129  m_hnvtxweightprof->GetXaxis()->SetTitle("Number of Interactions");
130  }
131 
132  m_hmainvtxx = tfserv->make<TH1F>("mainvtxx","Main vertex x position",200,-.5,.5);
133  m_hmainvtxx->GetXaxis()->SetTitle("X (cm)");
134  m_hmainvtxy = tfserv->make<TH1F>("mainvtxy","Main vertex y position",200,-.5,.5);
135  m_hmainvtxy->GetXaxis()->SetTitle("Y (cm)");
136  m_hmainvtxz = tfserv->make<TH1F>("mainvtxz","Main vertex z position",600,-30.,30.);
137  m_hmainvtxz->GetXaxis()->SetTitle("Z (cm)");
138  m_hpileupvtxz = tfserv->make<TH1F>("pileupvtxz","PileUp vertices z position",600,-30.,30.);
139  m_hpileupvtxz->GetXaxis()->SetTitle("Z (cm)");
140 
141 }
142 
143 
145 {
146 
147  // do anything here that needs to be done at desctruction time
148  // (e.g. close files, deallocate resources etc.)
149 
150 }
151 
152 
153 //
154 // member functions
155 //
156 
157 // ------------ method called to for each event ------------
158 void
160 {
161  using namespace edm;
162 
163  double weight = 1.;
164 
165  if(m_useweight) {
166  Handle<double> weightprod;
167  iEvent.getByLabel(m_weight,weightprod);
168 
169  weight = *weightprod;
170 
171  }
172 
173 
175  iEvent.getByLabel(m_pileupcollection,pileupinfos);
176 
177  //
178 
179  if(pileupinfos.isValid()) {
180 
181  // look for the intime PileupSummaryInfo
182 
183  std::vector<PileupSummaryInfo>::const_iterator pileupinfo;
184  for(pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end() ; ++pileupinfo) {
185  if(pileupinfo->getBunchCrossing()==0) break;
186  }
187 
188  //
189 
190  if(pileupinfo->getBunchCrossing()!=0) {
191  edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing();
192  }
193  else {
194 
195  m_hlumi->Fill(pileupinfo->getTrueNumInteractions(),weight);
196  m_hnvtx->Fill(pileupinfo->getPU_NumInteractions(),weight);
197  m_hnvtxvslumi->Fill(pileupinfo->getTrueNumInteractions(),pileupinfo->getPU_NumInteractions(),weight);
198 
199  if(m_useweight) {
200  m_hnvtxweight->Fill(pileupinfo->getPU_NumInteractions(),1.-weight);
201  m_hnvtxweightprof->Fill(pileupinfo->getPU_NumInteractions(),1.-weight);
202  }
203 
204  const std::vector<float>& zpositions = pileupinfo->getPU_zpositions();
205 
206  for(std::vector<float>::const_iterator zpos = zpositions.begin() ; zpos != zpositions.end() ; ++zpos) {
207 
208  m_hpileupvtxz->Fill(*zpos,weight);
209 
210  }
211  }
212  }
213  // main interaction part
214 
215  Handle< HepMCProduct > EvtHandle ;
216  iEvent.getByLabel(m_mctruthcollection, EvtHandle ) ;
217 
218  if(EvtHandle.isValid()) {
219 
220  const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
221 
222  // get the first vertex
223 
224  if(Evt->vertices_begin() != Evt->vertices_end()) {
225 
226  m_hmainvtxx->Fill((*Evt->vertices_begin())->point3d().x()/10.,weight);
227  m_hmainvtxy->Fill((*Evt->vertices_begin())->point3d().y()/10.,weight);
228  m_hmainvtxz->Fill((*Evt->vertices_begin())->point3d().z()/10.,weight);
229 
230  }
231  }
232 }
233 
234 void
236 {
237 }
238 
239 void
241 {
242 }
243 
244 
245 
246 // ------------ method called once each job just before starting event loop ------------
247 void
249 {
250 }
251 
252 // ------------ method called once each job just after ending the event loop ------------
253 void
255 {
256 }
257 
258 //define this as a plug-in
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int iEvent
Definition: GenABIO.cc:243
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
edm::InputTag m_mctruthcollection
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
MCVerticesAnalyzer(const edm::ParameterSet &)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
edm::InputTag m_pileupcollection
T * make() const
make new ROOT object
int weight
Definition: histoStyle.py:50
virtual void endRun(const edm::Run &, const edm::EventSetup &)
Definition: Run.h:36