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 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 #include <numeric>
23 
24 #include <vector>
25 
26 // user include files
29 
33 
35 
37 
40 
42 
46 
47 
48 #include "TH1F.h"
49 #include "TH2F.h"
50 #include "TProfile.h"
51 
52 //
53 // class decleration
54 //
55 
56 
58 public:
59  explicit MCVerticesAnalyzer(const edm::ParameterSet&);
61 
62 private:
63  virtual void beginJob() ;
64  virtual void beginRun(const edm::Run&, const edm::EventSetup&);
65  virtual void endRun(const edm::Run&, const edm::EventSetup&);
66  virtual void analyze(const edm::Event&, const edm::EventSetup&);
67  virtual void endJob() ;
68 
69  // ----------member data ---------------------------
70 
71 
72 
73  const bool m_useweight;
74 
78 
79  TH1F* m_hnvtx;
80  TH2F* m_hnvtxvsbx;
81  TH1F* m_hlumi;
84  TProfile* m_hnvtxweightprof;
85  TH1F* m_hmainvtxx;
86  TH1F* m_hmainvtxy;
87  TH1F* m_hmainvtxz;
89 
90 };
91 
92 //
93 // constants, enums and typedefs
94 //
95 
96 //
97 // static data member definitions
98 //
99 
100 //
101 // constructors and destructor
102 //
104  : m_useweight( iConfig.getParameter< bool >( "useWeight" ) )
105  , m_doubleToken( consumes< double >( iConfig.getParameter< edm::InputTag >( "weightProduct" ) ) )
106  , m_vecPileupSummaryInfoToken( consumes< std::vector<PileupSummaryInfo> >( iConfig.getParameter< edm::InputTag >( "pileupSummaryCollection" ) ) )
107  , m_hepMCProductToken( consumes< edm::HepMCProduct >( iConfig.getParameter< edm::InputTag >( "mcTruthCollection" ) ) )
108 {
109  //now do what ever initialization is needed
110 
111 
112 
114 
115  m_hnvtx = tfserv->make<TH1F>("nvtx","Number of pileup vertices",60,-0.5,59.5);
116  m_hnvtx->GetXaxis()->SetTitle("Number of Interactions");
117 
118  m_hnvtxvsbx = tfserv->make<TH2F>("nvtxvsbx","Number of pileup vertices vs BX",9,-4.5,4.5,60,-0.5,59.5);
119  m_hnvtxvsbx->GetXaxis()->SetTitle("BX number");
120  m_hnvtxvsbx->GetYaxis()->SetTitle("Number of Interactions");
121 
122  m_hlumi = tfserv->make<TH1F>("lumi","BX luminosity*xsect",200,0.,50.);
123  m_hlumi->GetXaxis()->SetTitle("Average Number of Interactions");
124 
125  m_hnvtxvslumi = tfserv->make<TH2F>("nvtxvslumi","Npileup vs BX luminosity*xsect",200,0.,50.,60,-0.5,59.5);
126  m_hnvtxvslumi->GetXaxis()->SetTitle("Average Number of Interactions"); m_hnvtxvslumi->GetYaxis()->SetTitle("Number of Interactions");
127 
128  if(m_useweight) {
129  m_hnvtxweight = tfserv->make<TH1F>("nvtxweight","Number of pileup vertices (1-w)",60,-0.5,59.5);
130  m_hnvtxweight->GetXaxis()->SetTitle("Number of Interactions");
131  m_hnvtxweightprof = tfserv->make<TProfile>("nvtxweightprof","Mean (1-w) vs Number of pileup interactions",60,-0.5,59.5);
132  m_hnvtxweightprof->GetXaxis()->SetTitle("Number of Interactions");
133  }
134 
135  m_hmainvtxx = tfserv->make<TH1F>("mainvtxx","Main vertex x position",200,-.5,.5);
136  m_hmainvtxx->GetXaxis()->SetTitle("X (cm)");
137  m_hmainvtxy = tfserv->make<TH1F>("mainvtxy","Main vertex y position",200,-.5,.5);
138  m_hmainvtxy->GetXaxis()->SetTitle("Y (cm)");
139  m_hmainvtxz = tfserv->make<TH1F>("mainvtxz","Main vertex z position",600,-30.,30.);
140  m_hmainvtxz->GetXaxis()->SetTitle("Z (cm)");
141  m_hpileupvtxz = tfserv->make<TH1F>("pileupvtxz","PileUp vertices z position",600,-30.,30.);
142  m_hpileupvtxz->GetXaxis()->SetTitle("Z (cm)");
143 
144 }
145 
146 
148 {
149 
150  // do anything here that needs to be done at desctruction time
151  // (e.g. close files, deallocate resources etc.)
152 
153 }
154 
155 
156 //
157 // member functions
158 //
159 
160 // ------------ method called to for each event ------------
161 void
163 {
164 
165  double weight = 1.;
166 
167  if(m_useweight) {
168  edm::Handle<double> weightprod;
169  iEvent.getByToken( m_doubleToken, weightprod );
170 
171  weight = *weightprod;
172 
173  }
174 
175 
177  iEvent.getByToken( m_vecPileupSummaryInfoToken, pileupinfos );
178 
179  //
180 
181  if(pileupinfos.isValid()) {
182 
183  // look for the intime PileupSummaryInfo
184 
185  std::vector<PileupSummaryInfo>::const_iterator pileupinfo;
186 
187  for(pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end() ; ++pileupinfo) {
188  m_hnvtxvsbx->Fill(pileupinfo->getBunchCrossing(),pileupinfo->getPU_NumInteractions(),weight);
189  }
190 
191 
192  for(pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end() ; ++pileupinfo) {
193  if(pileupinfo->getBunchCrossing()==0) break;
194  }
195 
196  //
197 
198  if(pileupinfo->getBunchCrossing()!=0) {
199  edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing();
200  }
201  else {
202 
203  m_hlumi->Fill(pileupinfo->getTrueNumInteractions(),weight);
204  m_hnvtx->Fill(pileupinfo->getPU_NumInteractions(),weight);
205  m_hnvtxvslumi->Fill(pileupinfo->getTrueNumInteractions(),pileupinfo->getPU_NumInteractions(),weight);
206 
207  if(m_useweight) {
208  m_hnvtxweight->Fill(pileupinfo->getPU_NumInteractions(),1.-weight);
209  m_hnvtxweightprof->Fill(pileupinfo->getPU_NumInteractions(),1.-weight);
210  }
211 
212  const std::vector<float>& zpositions = pileupinfo->getPU_zpositions();
213 
214  for(std::vector<float>::const_iterator zpos = zpositions.begin() ; zpos != zpositions.end() ; ++zpos) {
215 
216  m_hpileupvtxz->Fill(*zpos,weight);
217 
218  }
219  }
220  }
221  // main interaction part
222 
224  iEvent.getByToken( m_hepMCProductToken, EvtHandle );
225 
226  if(EvtHandle.isValid()) {
227 
228  const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
229 
230  // get the first vertex
231 
232  if(Evt->vertices_begin() != Evt->vertices_end()) {
233 
234  m_hmainvtxx->Fill((*Evt->vertices_begin())->point3d().x()/10.,weight);
235  m_hmainvtxy->Fill((*Evt->vertices_begin())->point3d().y()/10.,weight);
236  m_hmainvtxz->Fill((*Evt->vertices_begin())->point3d().z()/10.,weight);
237 
238  }
239  }
240 }
241 
242 void
244 {
245 }
246 
247 void
249 {
250 }
251 
252 
253 
254 // ------------ method called once each job just before starting event loop ------------
255 void
257 {
258 }
259 
260 // ------------ method called once each job just after ending the event loop ------------
261 void
263 {
264 }
265 
266 //define this as a plug-in
edm::EDGetTokenT< double > m_doubleToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > m_vecPileupSummaryInfoToken
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
bool isValid() const
Definition: HandleBase.h:76
MCVerticesAnalyzer(const edm::ParameterSet &)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
int weight
Definition: histoStyle.py:50
edm::EDGetTokenT< edm::HepMCProduct > m_hepMCProductToken
virtual void endRun(const edm::Run &, const edm::EventSetup &)
Definition: Run.h:41