CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PrimaryVertexMonitor.cc
Go to the documentation of this file.
5 
6 
8 
10 
11 #include "TMath.h"
12 
13 using namespace reco;
14 using namespace edm;
15 
17  : conf_ ( pSet )
18  , TopFolderName_ ( pSet.getParameter<std::string>("TopFolderName") )
19  , AlignmentLabel_( pSet.getParameter<std::string>("AlignmentLabel"))
20  , nbvtx(NULL)
21  , bsX(NULL)
22  , bsY(NULL)
23  , bsZ(NULL)
24  , bsSigmaZ(NULL)
25  , bsDxdz(NULL)
26  , bsDydz(NULL)
27  , bsBeamWidthX(NULL)
28  , bsBeamWidthY(NULL)
29  , bsType(NULL)
30  , sumpt(NULL)
31  , ntracks(NULL)
32  , weight(NULL)
33  , chi2ndf(NULL)
34  , chi2prob(NULL)
35  , dxy(NULL)
36  , dxy2(NULL)
37  , dz(NULL)
38  , dxyErr(NULL)
39  , dzErr(NULL)
40  , dxyVsPhi_pt1(NULL)
41  , dzVsPhi_pt1(NULL)
42  , dxyVsEta_pt1(NULL)
43  , dzVsEta_pt1(NULL)
44  , dxyVsPhi_pt10(NULL)
45  , dzVsPhi_pt10(NULL)
46  , dxyVsEta_pt10(NULL)
47  , dzVsEta_pt10(NULL)
48 {
49  // dqmStore_ = edm::Service<DQMStore>().operator->();
50 
51 
52  vertexInputTag_ = pSet.getParameter<InputTag>("vertexLabel");
53  beamSpotInputTag_ = pSet.getParameter<InputTag>("beamSpotLabel");
54  vertexToken_ = consumes<reco::VertexCollection>(vertexInputTag_);
55  beamspotToken_ = consumes<reco::BeamSpot> (beamSpotInputTag_);
56 
57 }
58 
59 // -- BeginRun
60 //---------------------------------------------------------------------------------//
61 void
63  edm::Run const &, edm::EventSetup const &) {
64 
65  std::string dqmLabel = "";
66 
67  //
68  // Book all histograms.
69  //
70 
71  // get the store
72  dqmLabel = TopFolderName_+"/"+vertexInputTag_.label();
73  iBooker.setCurrentFolder(dqmLabel);
74 
75 // xPos = iBooker.book1D ("xPos","x Coordinate" ,100, -0.1, 0.1);
76 
77  nbvtx = iBooker.book1D("vtxNbr","Reconstructed Vertices in Event",50,-0.5,49.5);
78  nbgvtx = iBooker.book1D("goodvtxNbr","Reconstructed Good Vertices in Event",50,-0.5,49.5);
79 
80  // to be configured each year...
81  auto vposx = conf_.getParameter<double>("Xpos");
82  auto vposy = conf_.getParameter<double>("Ypos");
83 
84  nbtksinvtx[0] = iBooker.book1D("otherVtxTrksNbr","Reconstructed Tracks in Vertex (other Vtx)",40,-0.5,99.5);
85  trksWeight[0] = iBooker.book1D("otherVtxTrksWeight","Total weight of Tracks in Vertex (other Vtx)",40,0,100.);
86  vtxchi2[0] = iBooker.book1D("otherVtxChi2","#chi^{2} (other Vtx)",100,0.,200.);
87  vtxndf[0] = iBooker.book1D("otherVtxNdf","ndof (other Vtx)",100,0.,200.);
88  vtxprob[0] = iBooker.book1D("otherVtxProb","#chi^{2} probability (other Vtx)",100,0.,1.);
89  nans[0] = iBooker.book1D("otherVtxNans","Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (other Vtx)",9,0.5,9.5);
90 
91  nbtksinvtx[1] = iBooker.book1D("tagVtxTrksNbr","Reconstructed Tracks in Vertex (tagged Vtx)",100,-0.5,99.5);
92  trksWeight[1] = iBooker.book1D("tagVtxTrksWeight","Total weight of Tracks in Vertex (tagged Vtx)",100,0,100.);
93  vtxchi2[1] = iBooker.book1D("tagVtxChi2","#chi^{2} (tagged Vtx)",100,0.,200.);
94  vtxndf[1] = iBooker.book1D("tagVtxNdf","ndof (tagged Vtx)",100,0.,200.);
95  vtxprob[1] = iBooker.book1D("tagVtxProb","#chi^{2} probability (tagged Vtx)",100,0.,1.);
96  nans[1] = iBooker.book1D("tagVtxNans","Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (tagged Vtx)",9,0.5,9.5);
97 
98  xrec[0] = iBooker.book1D("otherPosX","Position x Coordinate (other Vtx)",100,vposx-0.1,vposx+0.1);
99  yrec[0] = iBooker.book1D("otherPosY","Position y Coordinate (other Vtx)",100,vposy-0.1,vposy+0.1);
100  zrec[0] = iBooker.book1D("otherPosZ","Position z Coordinate (other Vtx)",100,-20.,20.);
101  xDiff[0] = iBooker.book1D("otherDiffX","X distance from BeamSpot (other Vtx)",100,-500,500);
102  yDiff[0] = iBooker.book1D("otherDiffY","Y distance from BeamSpot (other Vtx)",100,-500,500);
103  xerr[0] = iBooker.book1D("otherErrX","Uncertainty x Coordinate (other Vtx)",100,0.,100);
104  yerr[0] = iBooker.book1D("otherErrY","Uncertainty y Coordinate (other Vtx)",100,0.,100);
105  zerr[0] = iBooker.book1D("otherErrZ","Uncertainty z Coordinate (other Vtx)",100,0.,100);
106  xerrVsTrks[0] = iBooker.book2D("otherErrVsWeightX","Uncertainty x Coordinate vs. track weight (other Vtx)",100,0,100.,100,0.,100);
107  yerrVsTrks[0] = iBooker.book2D("otherErrVsWeightY","Uncertainty y Coordinate vs. track weight (other Vtx)",100,0,100.,100,0.,100);
108  zerrVsTrks[0] = iBooker.book2D("otherErrVsWeightZ","Uncertainty z Coordinate vs. track weight (other Vtx)",100,0,100.,100,0.,100);
109 
110 
111  xrec[1] = iBooker.book1D("tagPosX","Position x Coordinate (tagged Vtx)",100,vposx-0.1,vposx+0.1);
112  yrec[1] = iBooker.book1D("tagPosY","Position y Coordinate (tagged Vtx)",100,vposx-0.1,vposy+0.1);
113  zrec[1] = iBooker.book1D("tagPosZ","Position z Coordinate (tagged Vtx)",100,-20.,20.);
114  xDiff[1] = iBooker.book1D("tagDiffX","X distance from BeamSpot (tagged Vtx)",100,-500, 500);
115  yDiff[1] = iBooker.book1D("tagDiffY","Y distance from BeamSpot (tagged Vtx)",100,-500, 500);
116  xerr[1] = iBooker.book1D("tagErrX","Uncertainty x Coordinate (tagged Vtx)",100,0.,100);
117  yerr[1] = iBooker.book1D("tagErrY","Uncertainty y Coordinate (tagged Vtx)",100,0.,100);
118  zerr[1] = iBooker.book1D("tagErrZ","Uncertainty z Coordinate (tagged Vtx)",100,0.,100);
119  xerrVsTrks[1] = iBooker.book2D("tagErrVsWeightX","Uncertainty x Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
120  yerrVsTrks[1] = iBooker.book2D("tagErrVsWeightY","Uncertainty y Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
121  zerrVsTrks[1] = iBooker.book2D("tagErrVsWeightZ","Uncertainty z Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
122 
123  type[0] = iBooker.book1D("otherType","Vertex type (other Vtx)",3,-0.5,2.5);
124  type[1] = iBooker.book1D("tagType","Vertex type (tagged Vtx)",3,-0.5,2.5);
125  for (int i=0;i<2;++i){
126  type[i]->getTH1F()->GetXaxis()->SetBinLabel(1,"Valid, real");
127  type[i]->getTH1F()->GetXaxis()->SetBinLabel(2,"Valid, fake");
128  type[i]->getTH1F()->GetXaxis()->SetBinLabel(3,"Invalid");
129  }
130 
131 
132  // get the store
133  dqmLabel = TopFolderName_+"/"+beamSpotInputTag_.label();
134  iBooker.setCurrentFolder(dqmLabel);
135 
136  bsX = iBooker.book1D("bsX", "BeamSpot x0", 100,-0.1,0.1);
137  bsY = iBooker.book1D("bsY", "BeamSpot y0", 100,-0.1,0.1);
138  bsZ = iBooker.book1D("bsZ", "BeamSpot z0", 100,-2.,2.);
139  bsSigmaZ = iBooker.book1D("bsSigmaZ", "BeamSpot sigmaZ", 100, 0., 10. );
140  bsDxdz = iBooker.book1D("bsDxdz", "BeamSpot dxdz", 100, -0.0003, 0.0003);
141  bsDydz = iBooker.book1D("bsDydz", "BeamSpot dydz", 100, -0.0003, 0.0003);
142  bsBeamWidthX = iBooker.book1D("bsBeamWidthX", "BeamSpot BeamWidthX", 100, 0., 100.);
143  bsBeamWidthY = iBooker.book1D("bsBeamWidthY", "BeamSpot BeamWidthY", 100, 0., 100.);
144  bsType = iBooker.book1D("bsType", "BeamSpot type", 4, -1.5, 2.5);
145  bsType->getTH1F()->GetXaxis()->SetBinLabel(1,"Unknown");
146  bsType->getTH1F()->GetXaxis()->SetBinLabel(2,"Fake");
147  bsType->getTH1F()->GetXaxis()->SetBinLabel(3,"LHC");
148  bsType->getTH1F()->GetXaxis()->SetBinLabel(4,"Tracker");
149 
150 
151  // get the store
152  dqmLabel = TopFolderName_+"/"+AlignmentLabel_;
153  iBooker.setCurrentFolder(dqmLabel);
154 
155  int TKNoBin = conf_.getParameter<int>( "TkSizeBin");
156  double TKNoMin = conf_.getParameter<double>("TkSizeMin");
157  double TKNoMax = conf_.getParameter<double>("TkSizeMax");
158 
159  int DxyBin = conf_.getParameter<int>( "DxyBin");
160  double DxyMin = conf_.getParameter<double>("DxyMin");
161  double DxyMax = conf_.getParameter<double>("DxyMax");
162 
163  int DzBin = conf_.getParameter<int>( "DzBin");
164  double DzMin = conf_.getParameter<double>("DzMin");
165  double DzMax = conf_.getParameter<double>("DzMax");
166 
167  int PhiBin = conf_.getParameter<int>( "PhiBin");
168  double PhiMin = conf_.getParameter<double>("PhiMin");
169  double PhiMax = conf_.getParameter<double>("PhiMax");
170 
171  int EtaBin = conf_.getParameter<int>( "EtaBin");
172  double EtaMin = conf_.getParameter<double>("EtaMin");
173  double EtaMax = conf_.getParameter<double>("EtaMax");
174 
175 
176  ntracks = iBooker.book1D("ntracks","number of PV tracks (p_{T} > 1 GeV)", TKNoBin, TKNoMin, TKNoMax);
177  ntracks->setAxisTitle("Number of PV Tracks (p_{T} > 1 GeV) per Event", 1);
178  ntracks->setAxisTitle("Number of Event", 2);
179 
180  weight = iBooker.book1D("weight","weight of PV tracks (p_{T} > 1 GeV)", 100, 0., 1.);
181  weight->setAxisTitle("weight of PV Tracks (p_{T} > 1 GeV) per Event", 1);
182  weight->setAxisTitle("Number of Event", 2);
183 
184  sumpt = iBooker.book1D("sumpt", "#Sum p_{T} of PV tracks (p_{T} > 1 GeV)", 100,-0.5,249.5);
185  chi2ndf = iBooker.book1D("chi2ndf", "PV tracks (p_{T} > 1 GeV) #chi^{2}/ndof", 100, 0., 20. );
186  chi2prob = iBooker.book1D("chi2prob","PV tracks (p_{T} > 1 GeV) #chi^{2} probability",100, 0., 1. );
187  dxy = iBooker.book1D("dxy", "PV tracks (p_{T} > 1 GeV) d_{xy} (cm)", DxyBin, DxyMin, DxyMax);
188  dxy2 = iBooker.book1D("dxyzoom", "PV tracks (p_{T} > 1 GeV) d_{xy} (cm)", DxyBin, DxyMin/5., DxyMax/5.);
189  dz = iBooker.book1D("dz", "PV tracks (p_{T} > 1 GeV) d_{z} (cm)", DzBin, DzMin, DzMax );
190  dxyErr = iBooker.book1D("dxyErr", "PV tracks (p_{T} > 1 GeV) d_{xy} error (cm)", 100, 0., 0.2 );
191  dzErr = iBooker.book1D("dzErr", "PV tracks (p_{T} > 1 GeV) d_{z} error(cm)", 100, 0., 1. );
192 
193  dxyVsPhi_pt1 = iBooker.bookProfile("dxyVsPhi_pt1", "PV tracks (p_{T} > 1 GeV) d_{xy} (cm) VS track #phi",PhiBin, PhiMin, PhiMax, DxyBin, DxyMin, DxyMax,"");
194  dxyVsPhi_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) #phi", 1);
195  dxyVsPhi_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) d_{xy}",2);
196  dzVsPhi_pt1 = iBooker.bookProfile("dzVsPhi_pt1", "PV tracks (p_{T} > 1 GeV) d_{z} (cm) VS track #phi", PhiBin, PhiMin, PhiMax, DzBin, DzMin, DzMax, "");
197  dzVsPhi_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) #phi", 1);
198  dzVsPhi_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) d_{z}",2);
199 
200  dxyVsEta_pt1 = iBooker.bookProfile("dxyVsEta_pt1", "PV tracks (p_{T} > 1 GeV) d_{xy} (cm) VS track #eta",EtaBin, EtaMin, EtaMax, DxyBin, DxyMin, DxyMax,"");
201  dxyVsEta_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) #eta", 1);
202  dxyVsEta_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) d_{xy}",2);
203  dzVsEta_pt1 = iBooker.bookProfile("dzVsEta_pt1", "PV tracks (p_{T} > 1 GeV) d_{z} (cm) VS track #eta", EtaBin, EtaMin, EtaMax, DzBin, DzMin, DzMax, "");
204  dzVsEta_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) #eta", 1);
205  dzVsEta_pt1->setAxisTitle("PV track (p_{T} > 1 GeV) d_{z}",2);
206 
207  dxyVsPhi_pt10 = iBooker.bookProfile("dxyVsPhi_pt10", "PV tracks (p_{T} > 1 GeV) d_{xy} (cm) VS track #phi",PhiBin, PhiMin, PhiMax, DxyBin, DxyMin, DxyMax,"");
208  dxyVsPhi_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) #phi", 1);
209  dxyVsPhi_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) d_{xy}",2);
210  dzVsPhi_pt10 = iBooker.bookProfile("dzVsPhi_pt10", "PV tracks (p_{T} > 10 GeV) d_{z} (cm) VS track #phi", PhiBin, PhiMin, PhiMax, DzBin, DzMin, DzMax, "");
211  dzVsPhi_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) #phi", 1);
212  dzVsPhi_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) d_{z}",2);
213 
214  dxyVsEta_pt10 = iBooker.bookProfile("dxyVsEta_pt10", "PV tracks (p_{T} > 10 GeV) d_{xy} (cm) VS track #eta",EtaBin, EtaMin, EtaMax, DxyBin, DxyMin, DxyMax,"");
215  dxyVsEta_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) #eta", 1);
216  dxyVsEta_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) d_{xy}",2);
217  dzVsEta_pt10 = iBooker.bookProfile("dzVsEta_pt10", "PV tracks (p_{T} > 10 GeV) d_{z} (cm) VS track #eta", EtaBin, EtaMin, EtaMax, DzBin, DzMin, DzMax, "");
218  dzVsEta_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) #eta", 1);
219  dzVsEta_pt10->setAxisTitle("PV track (p_{T} > 10 GeV) d_{z}",2);
220 
221 }
222 
223 
225 {}
226 
228 {
230  iEvent.getByToken(vertexToken_, recVtxs);
231 
232  edm::Handle<reco::BeamSpot> beamSpotHandle;
233  iEvent.getByToken(beamspotToken_,beamSpotHandle);
234 
235  //
236  // check for absent products and simply "return" in that case
237  //
238  if (recVtxs.isValid() == false || beamSpotHandle.isValid()== false){
239  edm::LogWarning("PrimaryVertexMonitor")
240  <<" Some products not available in the event: VertexCollection "
241  <<vertexInputTag_<<" "
242  <<recVtxs.isValid() <<" BeamSpot "
243  <<beamSpotInputTag_<<" "
244  <<beamSpotHandle.isValid()<<". Skipping plots for this event";
245  return;
246  }
247 
248  BeamSpot beamSpot = *beamSpotHandle;
249 
250  nbvtx->Fill(recVtxs->size()*1.);
251  int ng=0;
252  for (auto const & vx : (*recVtxs) )
253  if (vx.isValid() && !vx.isFake() && vx.ndof()>=4.) ++ng;
254  nbgvtx->Fill(ng*1.);
255 
256  vertexPlots(recVtxs->front(), beamSpot, 1);
257 
258  // fill PV tracks MEs (as now, for alignment)
259  pvTracksPlots(recVtxs->front());
260 
261  for(reco::VertexCollection::const_iterator v=recVtxs->begin()+1;
262  v!=recVtxs->end(); ++v){
263  vertexPlots(*v, beamSpot, 0);
264  }
265  // Beamline plots:
266  bsX->Fill(beamSpot.x0());
267  bsY->Fill(beamSpot.y0());
268  bsZ->Fill(beamSpot.z0());
269  bsSigmaZ->Fill(beamSpot.sigmaZ());
270  bsDxdz->Fill(beamSpot.dxdz());
271  bsDydz->Fill(beamSpot.dydz());
272  bsBeamWidthX->Fill(beamSpot.BeamWidthX()*10000);
273  bsBeamWidthY->Fill(beamSpot.BeamWidthY()*10000);
274  // bsType->Fill(beamSpot.type());
275 
276 }
277 
278 void
280 {
281 
282  const math::XYZPoint myVertex(v.position().x(),v.position().y(),v.position().z());
283 
284  if ( !v.isValid() ) return;
285  if ( v.isFake() ) return;
286 
287  if ( v.tracksSize() == 0 ) {
288  ntracks -> Fill ( 0 );
289  return;
290  }
291 
292  size_t nTracks = 0;
293  float sumPT = 0.;
295 
296  bool isHighPurity = (**t).quality(reco::TrackBase::highPurity);
297  if ( !isHighPurity ) continue;
298 
299  float pt = (**t).pt();
300  if ( pt < 1. ) continue;
301 
302  nTracks++;
303 
304  float eta = (**t).eta();
305  float phi = (**t).phi();
306 
307  float w = v.trackWeight(*t);
308  float chi2NDF = (**t).normalizedChi2();
309  float chi2Prob = TMath::Prob((**t).chi2(),(int)(**t).ndof());
310  float Dxy = (**t).dxy(myVertex);
311  float Dz = (**t).dz(myVertex);
312  float DxyErr = (**t).dxyError();
313  float DzErr = (**t).dzError();
314 
315  sumPT += pt*pt;
316 
317 
318  // fill MEs
319  weight -> Fill (w);
320  chi2ndf -> Fill (chi2NDF);
321  chi2prob -> Fill (chi2Prob);
322  dxy -> Fill (Dxy);
323  dxy2 -> Fill (Dxy);
324  dz -> Fill (Dz);
325  dxyErr -> Fill (DxyErr);
326  dzErr -> Fill (DzErr);
327 
328  dxyVsPhi_pt1 -> Fill (phi,Dxy);
329  dzVsPhi_pt1 -> Fill (phi,Dz);
330  dxyVsEta_pt1 -> Fill (eta,Dxy);
331  dzVsEta_pt1 -> Fill (eta,Dz);
332 
333  if ( pt < 10. ) continue;
334  dxyVsPhi_pt10 -> Fill (phi,Dxy);
335  dzVsPhi_pt10 -> Fill (phi,Dz);
336  dxyVsEta_pt10 -> Fill (eta,Dxy);
337  dzVsEta_pt10 -> Fill (eta,Dz);
338  }
339  ntracks -> Fill (float(nTracks));
340  sumpt -> Fill (sumPT);
341 
342 }
343 
345 {
346 
347  if (i < 0 || i > 1) return;
348  if (!v.isValid()) type[i]->Fill(2.);
349  else if (v.isFake()) type[i]->Fill(1.);
350  else type[i]->Fill(0.);
351 
352  if (v.isValid() && !v.isFake()) {
353  float weight = 0;
355  t!=v.tracks_end(); t++) weight+= v.trackWeight(*t);
356  trksWeight[i]->Fill(weight);
357  nbtksinvtx[i]->Fill(v.tracksSize());
358 
359  vtxchi2[i]->Fill(v.chi2());
360  vtxndf[i]->Fill(v.ndof());
362 
363  xrec[i]->Fill(v.position().x());
364  yrec[i]->Fill(v.position().y());
365  zrec[i]->Fill(v.position().z());
366 
367  float xb = beamSpot.x0() + beamSpot.dxdz() * (v.position().z() - beamSpot.z0());
368  float yb = beamSpot.y0() + beamSpot.dydz() * (v.position().z() - beamSpot.z0());
369  xDiff[i]->Fill((v.position().x() - xb)*10000);
370  yDiff[i]->Fill((v.position().y() - yb)*10000);
371 
372  xerr[i]->Fill(v.xError()*10000);
373  yerr[i]->Fill(v.yError()*10000);
374  zerr[i]->Fill(v.zError()*10000);
375  xerrVsTrks[i]->Fill(weight, v.xError()*10000);
376  yerrVsTrks[i]->Fill(weight, v.yError()*10000);
377  zerrVsTrks[i]->Fill(weight, v.zError()*10000);
378 
379  nans[i]->Fill(1.,edm::isNotFinite(v.position().x())*1.);
380  nans[i]->Fill(2.,edm::isNotFinite(v.position().y())*1.);
381  nans[i]->Fill(3.,edm::isNotFinite(v.position().z())*1.);
382 
383  int index = 3;
384  for (int k = 0; k != 3; k++) {
385  for (int j = k; j != 3; j++) {
386  index++;
387  nans[i]->Fill(index*1., edm::isNotFinite(v.covariance(k, j))*1.);
388  // in addition, diagonal element must be positive
389  if (j == k && v.covariance(k, j) < 0) {
390  nans[i]->Fill(index*1., 1.);
391  }
392  }
393  }
394  }
395 }
396 
397 //define this as a plug-in
type
Definition: HCALResponse.h:21
MonitorElement * xerrVsTrks[2]
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
MonitorElement * xerr[2]
tuple t
Definition: tree.py:139
int i
Definition: DBlmapReader.cc:9
MonitorElement * vtxndf[2]
MonitorElement * dxyVsEta_pt1
MonitorElement * dzVsPhi_pt10
MonitorElement * dzVsPhi_pt1
const double w
Definition: UKUtility.cc:23
void pvTracksPlots(const reco::Vertex &v)
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
void vertexPlots(const reco::Vertex &v, const reco::BeamSpot &beamSpot, int i)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
double zError() const
error on z
Definition: Vertex.h:118
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * bsSigmaZ
edm::InputTag beamSpotInputTag_
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:60
#define NULL
Definition: scimark2.h:8
MonitorElement * chi2ndf
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:123
MonitorElement * dxyVsPhi_pt10
MonitorElement * xDiff[2]
const Point & position() const
position
Definition: Vertex.h:106
MonitorElement * dzVsEta_pt1
edm::ParameterSet conf_
void Fill(long long x)
MonitorElement * zerrVsTrks[2]
MonitorElement * nans[2]
double dydz() const
dydz slope
Definition: BeamSpot.h:84
int iEvent
Definition: GenABIO.cc:230
bool isNotFinite(T x)
Definition: isFinite.h:10
MonitorElement * yerrVsTrks[2]
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * dxyVsEta_pt10
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * zerr[2]
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
MonitorElement * yrec[2]
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
double chi2() const
chi-squares
Definition: Vertex.h:95
int j
Definition: DBlmapReader.cc:9
float ChiSquaredProbability(double chiSquared, double nrDOF)
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
MonitorElement * zrec[2]
const double EtaMin[kNumberCalorimeter]
MonitorElement * vtxchi2[2]
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
bool isValid() const
Definition: HandleBase.h:75
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
double ndof() const
Definition: Vertex.h:102
MonitorElement * dzVsEta_pt10
double xError() const
error on x
Definition: Vertex.h:114
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
bool isFake() const
Definition: Vertex.h:64
MonitorElement * chi2prob
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
MonitorElement * nbtksinvtx[2]
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
TH1F * getTH1F(void) const
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
std::string const & label() const
Definition: InputTag.h:43
MonitorElement * yerr[2]
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
double y0() const
y coordinate
Definition: BeamSpot.h:66
MonitorElement * xrec[2]
MonitorElement * vtxprob[2]
MonitorElement * ntracks
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
MonitorElement * bsBeamWidthY
int weight
Definition: histoStyle.py:50
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * yDiff[2]
double yError() const
error on y
Definition: Vertex.h:116
MonitorElement * trksWeight[2]
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:34
PrimaryVertexMonitor(const edm::ParameterSet &pSet)
Definition: Run.h:41
MonitorElement * bsBeamWidthX
MonitorElement * dxyVsPhi_pt1
double x0() const
x coordinate
Definition: BeamSpot.h:64