CMS 3D CMS Logo

MCvsRecoVerticesAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MCvsRecoVerticesAnalyzer
4 // Class: MCvsRecoVerticesAnalyzer
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 
49 
50 #include "TH1F.h"
51 #include "TH2F.h"
52 #include "TProfile.h"
53 
54 //
55 // class decleration
56 //
57 
58 
59 class MCvsRecoVerticesAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
60 public:
62  ~MCvsRecoVerticesAnalyzer() override;
63 
64 private:
65  void analyze(const edm::Event&, const edm::EventSetup&) override;
66 
67  // ----------member data ---------------------------
68 
69 
70 
71  const bool m_useweight;
78 
88 
92 
93 
94 };
95 
96 //
97 // constants, enums and typedefs
98 //
99 
100 //
101 // static data member definitions
102 //
103 
104 //
105 // constructors and destructor
106 //
108  : m_useweight( iConfig.getParameter< bool >( "useWeight" ) )
109  , m_useVisibleVertices( iConfig.getParameter< bool >( "useVisibleVertices" ) )
110  , m_histoParameters( iConfig.getUntrackedParameter< edm::ParameterSet >( "histoParameters", edm::ParameterSet() ) )
111  , m_doubleToken( consumes< double >( iConfig.getParameter< edm::InputTag >( "weightProduct" ) ) )
112  , m_vecPileupSummaryInfoToken( consumes< std::vector<PileupSummaryInfo> >( iConfig.getParameter< edm::InputTag >( "pileupSummaryCollection" ) ) )
113  , m_recoVertexCollectionToken( consumes< reco::VertexCollection >( iConfig.getParameter< edm::InputTag >( "pvCollection" ) ) )
114  , m_hepMCProductToken( consumes< edm::HepMCProduct >( iConfig.getParameter< edm::InputTag >( "mcTruthCollection" ) ) )
115 {
116  //now do what ever initialization is needed
117 
118  if(m_useVisibleVertices) edm::LogInfo("UseVisibleVertices") << "Only visible vertices will be used to compute Npileup";
119 
120  usesResource("TFileService");
122 
123  m_hrecovsmcnvtx2d = tfserv->make<TH2F>("recovsmcnvtx2d","Number of reco vertices vs pileup interactions",60,-0.5,59.5,60,-0.5,59.5);
124  m_hrecovsmcnvtx2d->GetXaxis()->SetTitle("Pileup Interactions"); m_hrecovsmcnvtx2d->GetYaxis()->SetTitle("Reco Vertices");
125  m_hrecovsmcnvtxprof = tfserv->make<TProfile>("recovsmcnvtxprof","Mean number of reco vs pileup vertices",60,-0.5,59.5);
126  m_hrecovsmcnvtxprof->GetXaxis()->SetTitle("Pileup Interactions"); m_hrecovsmcnvtxprof->GetYaxis()->SetTitle("Reco Vertices");
127 
128  m_hrecovsmclumi2d = tfserv->make<TH2F>("recovsmclumi2d","Number of reco vertices vs ave pileup interactions",200,0.,50.,60,-0.5,59.5);
129  m_hrecovsmclumi2d->GetXaxis()->SetTitle("Average Pileup Interactions"); m_hrecovsmclumi2d->GetYaxis()->SetTitle("Reco Vertices");
130  m_hrecovsmclumiprof = tfserv->make<TProfile>("recovsmclumiprof","Mean number of reco vs ave pileup vertices",200,0.,50.);
131  m_hrecovsmclumiprof->GetXaxis()->SetTitle("Average Pileup Interactions"); m_hrecovsmclumiprof->GetYaxis()->SetTitle("Reco Vertices");
132 
133  if(m_useweight) {
134  m_hrecovsmcnvtxweightedprof = tfserv->make<TProfile>("recovsmcnvtxweightedprof","Mean number of reco vs pileup vertices (1-w) weight",60,-0.5,59.5);
135  m_hrecovsmcnvtxweightedprof->GetXaxis()->SetTitle("Pileup Interactions"); m_hrecovsmcnvtxweightedprof->GetYaxis()->SetTitle("Reco Vertices (1-w)");
136 
137  m_hrecovsmclumiweightedprof = tfserv->make<TProfile>("recovsmclumiweightedprof","Mean number of reco vs ave pileup vertices (1-w) weight",
138  200,0.,50.);
139  m_hrecovsmclumiweightedprof->GetXaxis()->SetTitle("Average Pileup Interactions");
140  m_hrecovsmclumiweightedprof->GetYaxis()->SetTitle("Reco Vertices (1-w)");
141  }
142 
143  m_hdeltazfirst = tfserv->make<TH1F>("deltazfirst","Reco-MC vertex z position (first vertex)",
144  m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
145  m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
146  m_histoParameters.getUntrackedParameter<double>("zMax",1.));
147  m_hdeltazfirst->GetXaxis()->SetTitle("#Delta z (cm)"); m_hdeltazfirst->GetYaxis()->SetTitle("Events");
148 
149  m_hdeltazclose = tfserv->make<TH1F>("deltazclose","Reco-MC vertex z position (closest vertex)",
150  m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
151  m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
152  m_histoParameters.getUntrackedParameter<double>("zMax",1.));
153  m_hdeltazclose->GetXaxis()->SetTitle("#Delta z (cm)"); m_hdeltazclose->GetYaxis()->SetTitle("Events");
154 
155  m_hclosestvtx = tfserv->make<TH1F>("closestvtx","Closest reco vtx ID",30,-0.5,29.5);
156  m_hclosestvtx->GetXaxis()->SetTitle("Vtx ID"); m_hclosestvtx->GetYaxis()->SetTitle("Events");
157 
158  m_hdeltazfirstvsnpu = tfserv->make<TH2F>("deltazfirstvsnpu","Reco-MC vertex z position (first vertex) vs Npileup",30,-0.5,29.5,
159  m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
160  m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
161  m_histoParameters.getUntrackedParameter<double>("zMax",1.));
162  m_hdeltazfirstvsnpu->GetXaxis()->SetTitle("pileup Interactions"); m_hdeltazfirstvsnpu->GetYaxis()->SetTitle("#Delta z (cm)");
163 
164  m_hdeltazclosevsnpu = tfserv->make<TH2F>("deltazclosevsnpu","Reco-MC vertex z position (closest vertex) v Npileup",30,-0.5,29.5,
165  m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
166  m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
167  m_histoParameters.getUntrackedParameter<double>("zMax",1.));
168  m_hdeltazclosevsnpu->GetXaxis()->SetTitle("Pileup Interactions"); m_hdeltazclosevsnpu->GetYaxis()->SetTitle("#Delta z (cm)");
169 
170  m_hclosestvtxvsnpu = tfserv->make<TH2F>("closestvtxvsnpu","Closest reco vtx ID vs Npileup",30,-0.5,29.5,30,-0.5,29.5);
171  m_hclosestvtxvsnpu->GetXaxis()->SetTitle("Pileup Interactions"); m_hclosestvtxvsnpu->GetYaxis()->SetTitle("Vtx ID");
172 
173 }
174 
175 
177 {
178 
179  // do anything here that needs to be done at desctruction time
180  // (e.g. close files, deallocate resources etc.)
181 
182 }
183 
184 
185 //
186 // member functions
187 //
188 
189 // ------------ method called to for each event ------------
190 void
192 {
193 
194  double weight = 1.;
195 
196  if(m_useweight) {
197  edm::Handle< double > weightprod;
198  iEvent.getByToken( m_doubleToken, weightprod );
199 
200  weight = *weightprod;
201 
202  }
203 
205  iEvent.getByToken( m_vecPileupSummaryInfoToken, pileupinfos );
206 
207  // look for the intime PileupSummaryInfo
208 
209  std::vector<PileupSummaryInfo>::const_iterator pileupinfo;
210 
211  for(pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end() ; ++pileupinfo) {
212 
213  if(pileupinfo->getBunchCrossing()==0) break;
214 
215  }
216 
217  //
218 
220  iEvent.getByToken( m_recoVertexCollectionToken, pvcoll );
221 
222 
223  //
224 
225  if(pileupinfo->getBunchCrossing()!=0) {
226 
227  edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing();
228 
229  }
230  else {
231 
232  int npileup = pileupinfo->getPU_NumInteractions();
233 
234  if(m_useVisibleVertices) npileup = pileupinfo->getPU_zpositions().size();
235 
236  m_hrecovsmcnvtx2d->Fill(npileup,pvcoll->size(),weight);
237  m_hrecovsmcnvtxprof->Fill(npileup,pvcoll->size(),weight);
238 
239  m_hrecovsmclumi2d->Fill(pileupinfo->getTrueNumInteractions(),pvcoll->size(),weight);
240  m_hrecovsmclumiprof->Fill(pileupinfo->getTrueNumInteractions(),pvcoll->size(),weight);
241 
242  if(m_useweight) {
243  m_hrecovsmcnvtxweightedprof->Fill(npileup,pvcoll->size(),1.-weight);
244  m_hrecovsmclumiweightedprof->Fill(pileupinfo->getTrueNumInteractions(),pvcoll->size(),1.-weight);
245  }
246  //
247 
249  iEvent.getByToken( m_hepMCProductToken, EvtHandle );
250 
251  const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
252 
253  // compute the difference between the main interaction vertex z position and the first vertex of the collection
254 
255  if(!pvcoll->empty()) {
256  if(!(*pvcoll)[0].isFake()) {
257  // get the first vertex
258  if(Evt->vertices_begin() != Evt->vertices_end()) {
259  m_hdeltazfirst->Fill((*pvcoll)[0].z()-(*Evt->vertices_begin())->point3d().z()/10.,weight);
260  m_hdeltazfirstvsnpu->Fill(npileup,(*pvcoll)[0].z()-(*Evt->vertices_begin())->point3d().z()/10.,weight);
261  }
262  }
263  }
264 
265  // compute the difference between the main interaction vertex z position and the closest reco vertex
266 
267  double minabsdist = -1.;
268  double mindist = -999.;
269  int closestvtx = -1;
270 
271  for(unsigned int ivtx = 0 ; ivtx < pvcoll->size() ; ++ivtx) {
272 
273  if(closestvtx < 0 || minabsdist > std::abs((*pvcoll)[ivtx].z()-(*Evt->vertices_begin())->point3d().z()/10.)) {
274  mindist = (*pvcoll)[ivtx].z()-(*Evt->vertices_begin())->point3d().z()/10.;
275  closestvtx = ivtx;
276  minabsdist = std::abs(mindist);
277  }
278 
279  }
280  if(closestvtx >= 0) {
281  m_hdeltazclose->Fill(mindist,weight);
282  m_hdeltazclosevsnpu->Fill(npileup,mindist,weight);
283  m_hclosestvtx->Fill(closestvtx,weight);
284  m_hclosestvtxvsnpu->Fill(npileup,closestvtx,weight);
285  }
286 
287  }
288 }
289 
290 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Definition: weight.py:1
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void analyze(const edm::Event &, const edm::EventSetup &) override
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::VertexCollection > m_recoVertexCollectionToken
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< edm::HepMCProduct > m_hepMCProductToken
MCvsRecoVerticesAnalyzer(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > m_vecPileupSummaryInfoToken
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
edm::EDGetTokenT< double > m_doubleToken
fixed size matrix
HLT enums.
const edm::ParameterSet m_histoParameters