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 
8 
10 
11 using namespace reco;
12 using namespace edm;
13 
15 {
16  moduleLabel = pSet.getParameter<InputTag>("vertexLabel");
17  beamSpotLabel = pSet.getParameter<InputTag>("beamSpotLabel");
18 
19  //
20  // Book all histograms.
21  //
22 
23  // get the store
24  dqmStore_ = edm::Service<DQMStore>().operator->();
25  dqmLabel = "OfflinePV/"+moduleLabel.label();
26  dqmStore_->setCurrentFolder(dqmLabel);
27 
28 // xPos = dqmStore_->book1D ("xPos","x Coordinate" ,100, -0.1, 0.1);
29 
30  nbvtx = dqmStore_->book1D("vtxNbr","Reconstructed Vertices in Event",50,-0.5,49.5);
31 
32  nbtksinvtx[0] = dqmStore_->book1D("otherVtxTrksNbr","Reconstructed Tracks in Vertex (other Vtx)",40,-0.5,99.5);
33  trksWeight[0] = dqmStore_->book1D("otherVtxTrksWeight","Total weight of Tracks in Vertex (other Vtx)",40,0,100.);
34  vtxchi2[0] = dqmStore_->book1D("otherVtxChi2","#chi^{2} (other Vtx)",100,0.,200.);
35  vtxndf[0] = dqmStore_->book1D("otherVtxNdf","ndof (other Vtx)",100,0.,200.);
36  vtxprob[0] = dqmStore_->book1D("otherVtxProb","#chi^{2} probability (other Vtx)",100,0.,1.);
37  nans[0] = dqmStore_->book1D("otherVtxNans","Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (other Vtx)",9,0.5,9.5);
38 
39  nbtksinvtx[1] = dqmStore_->book1D("tagVtxTrksNbr","Reconstructed Tracks in Vertex (tagged Vtx)",100,-0.5,99.5);
40  trksWeight[1] = dqmStore_->book1D("tagVtxTrksWeight","Total weight of Tracks in Vertex (tagged Vtx)",100,0,100.);
41  vtxchi2[1] = dqmStore_->book1D("tagVtxChi2","#chi^{2} (tagged Vtx)",100,0.,200.);
42  vtxndf[1] = dqmStore_->book1D("tagVtxNdf","ndof (tagged Vtx)",100,0.,200.);
43  vtxprob[1] = dqmStore_->book1D("tagVtxProb","#chi^{2} probability (tagged Vtx)",100,0.,1.);
44  nans[1] = dqmStore_->book1D("tagVtxNans","Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (tagged Vtx)",9,0.5,9.5);
45 
46  xrec[0] = dqmStore_->book1D("otherPosX","Position x Coordinate (other Vtx)",100,-0.1,0.1);
47  yrec[0] = dqmStore_->book1D("otherPosY","Position y Coordinate (other Vtx)",100,-0.1,0.1);
48  zrec[0] = dqmStore_->book1D("otherPosZ","Position z Coordinate (other Vtx)",100,-20.,20.);
49  xDiff[0] = dqmStore_->book1D("otherDiffX","X distance from BeamSpot (other Vtx)",100,-500,500);
50  yDiff[0] = dqmStore_->book1D("otherDiffY","Y distance from BeamSpot (other Vtx)",100,-500,500);
51  xerr[0] = dqmStore_->book1D("otherErrX","Uncertainty x Coordinate (other Vtx)",100,-0.1,0.1);
52  yerr[0] = dqmStore_->book1D("otherErrY","Uncertainty y Coordinate (other Vtx)",100,-0.1,0.1);
53  zerr[0] = dqmStore_->book1D("otherErrZ","Uncertainty z Coordinate (other Vtx)",100,-20.,20.);
54  xerrVsTrks[0] = dqmStore_->book2D("otherErrVsWeightX","Uncertainty x Coordinate vs. track weight (other Vtx)",100,0,100.,100,-0.1,0.1);
55  yerrVsTrks[0] = dqmStore_->book2D("otherErrVsWeightY","Uncertainty y Coordinate vs. track weight (other Vtx)",100,0,100.,100,-0.1,0.1);
56  zerrVsTrks[0] = dqmStore_->book2D("otherErrVsWeightZ","Uncertainty z Coordinate vs. track weight (other Vtx)",100,0,100.,100,-0.1,0.1);
57 
58 
59  xrec[1] = dqmStore_->book1D("tagPosX","Position x Coordinate (tagged Vtx)",100,-0.1,0.1);
60  yrec[1] = dqmStore_->book1D("tagPosY","Position y Coordinate (tagged Vtx)",100,-0.1,0.1);
61  zrec[1] = dqmStore_->book1D("tagPosZ","Position z Coordinate (tagged Vtx)",100,-20.,20.);
62  xDiff[1] = dqmStore_->book1D("tagDiffX","X distance from BeamSpot (tagged Vtx)",100,-500, 500);
63  yDiff[1] = dqmStore_->book1D("tagDiffY","Y distance from BeamSpot (tagged Vtx)",100,-500, 500);
64  xerr[1] = dqmStore_->book1D("tagErrX","Uncertainty x Coordinate (tagged Vtx)",100,0.,100);
65  yerr[1] = dqmStore_->book1D("tagErrY","Uncertainty y Coordinate (tagged Vtx)",100,0.,100);
66  zerr[1] = dqmStore_->book1D("tagErrZ","Uncertainty z Coordinate (tagged Vtx)",100,0.,100);
67  xerrVsTrks[1] = dqmStore_->book2D("tagErrVsWeightX","Uncertainty x Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
68  yerrVsTrks[1] = dqmStore_->book2D("tagErrVsWeightY","Uncertainty y Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
69  zerrVsTrks[1] = dqmStore_->book2D("tagErrVsWeightZ","Uncertainty z Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
70 
71  type[0] = dqmStore_->book1D("otherType","Vertex type (other Vtx)",3,-0.5,2.5);
72  type[1] = dqmStore_->book1D("tagType","Vertex type (tagged Vtx)",3,-0.5,2.5);
73  for (int i=0;i<2;++i){
74  type[i]->getTH1F()->GetXaxis()->SetBinLabel(1,"Valid, real");
75  type[i]->getTH1F()->GetXaxis()->SetBinLabel(2,"Valid, fake");
76  type[i]->getTH1F()->GetXaxis()->SetBinLabel(3,"Invalid");
77  }
78 
79  bsX = dqmStore_->book1D("bsX", "BeamSpot x0", 100,-0.1,0.1);
80  bsY = dqmStore_->book1D("bsY", "BeamSpot y0", 100,-0.1,0.1);
81  bsZ = dqmStore_->book1D("bsZ", "BeamSpot z0", 100,-2.,2.);
82  bsSigmaZ = dqmStore_->book1D("bsSigmaZ", "BeamSpot sigmaZ", 100, 0., 10. );
83  bsDxdz = dqmStore_->book1D("bsDxdz", "BeamSpot dxdz", 100, -0.0003, 0.0003);
84  bsDydz = dqmStore_->book1D("bsDydz", "BeamSpot dydz", 100, -0.0003, 0.0003);
85  bsBeamWidthX = dqmStore_->book1D("bsBeamWidthX", "BeamSpot BeamWidthX", 100, 0., 100.);
86  bsBeamWidthY = dqmStore_->book1D("bsBeamWidthY", "BeamSpot BeamWidthY", 100, 0., 100.);
87  bsType = dqmStore_->book1D("bsType", "BeamSpot type", 4, -1.5, 2.5);
88  bsType->getTH1F()->GetXaxis()->SetBinLabel(1,"Unknown");
89  bsType->getTH1F()->GetXaxis()->SetBinLabel(2,"Fake");
90  bsType->getTH1F()->GetXaxis()->SetBinLabel(3,"LHC");
91  bsType->getTH1F()->GetXaxis()->SetBinLabel(4,"Tracker");
92 }
93 
94 
96 {}
97 
99 {
101  iEvent.getByLabel(moduleLabel, recVtxs);
102 
103  edm::Handle<reco::BeamSpot> beamSpotHandle;
104  iEvent.getByLabel(beamSpotLabel,beamSpotHandle);
105 
106  //
107  // check for absent products and simply "return" in that case
108  //
109  if (recVtxs.isValid() == false || beamSpotHandle.isValid()== false){
110  edm::LogWarning("PrimaryVertexMonitor")
111  <<" Some products not available in the event: VertexCollection "
112  <<moduleLabel<<" "
113  <<recVtxs.isValid() <<" BeamSpot "
114  <<beamSpotLabel<<" "
115  <<beamSpotHandle.isValid()<<". Skipping plots for this event";
116  return;
117  }
118 
119  BeamSpot beamSpot = *beamSpotHandle;
120 
121  nbvtx->Fill(recVtxs->size()*1.);
122 
123  vertexPlots(recVtxs->front(), beamSpot, 1);
124 
125  for(reco::VertexCollection::const_iterator v=recVtxs->begin()+1;
126  v!=recVtxs->end(); ++v){
127  vertexPlots(*v, beamSpot, 0);
128  }
129  // Beamline plots:
130  bsX->Fill(beamSpot.x0());
131  bsY->Fill(beamSpot.y0());
132  bsZ->Fill(beamSpot.z0());
133  bsSigmaZ->Fill(beamSpot.sigmaZ());
134  bsDxdz->Fill(beamSpot.dxdz());
135  bsDydz->Fill(beamSpot.dydz());
136  bsBeamWidthX->Fill(beamSpot.BeamWidthX()*10000);
137  bsBeamWidthY->Fill(beamSpot.BeamWidthY()*10000);
138  // bsType->Fill(beamSpot.type());
139 
140 }
141 
143 {
144 
145  if (i < 0 || i > 1) return;
146  if (!v.isValid()) type[i]->Fill(2.);
147  else if (v.isFake()) type[i]->Fill(1.);
148  else type[i]->Fill(0.);
149 
150  if (v.isValid() && !v.isFake()) {
151  float weight = 0;
153  t!=v.tracks_end(); t++) weight+= v.trackWeight(*t);
154  trksWeight[i]->Fill(weight);
155  nbtksinvtx[i]->Fill(v.tracksSize());
156 
157  vtxchi2[i]->Fill(v.chi2());
158  vtxndf[i]->Fill(v.ndof());
159  vtxprob[i]->Fill(ChiSquaredProbability(v.chi2() ,v.ndof()));
160 
161  xrec[i]->Fill(v.position().x());
162  yrec[i]->Fill(v.position().y());
163  zrec[i]->Fill(v.position().z());
164 
165  float xb = beamSpot.x0() + beamSpot.dxdz() * (v.position().z() - beamSpot.z0());
166  float yb = beamSpot.y0() + beamSpot.dydz() * (v.position().z() - beamSpot.z0());
167  xDiff[i]->Fill((v.position().x() - xb)*10000);
168  yDiff[i]->Fill((v.position().y() - yb)*10000);
169 
170  xerr[i]->Fill(v.xError()*10000);
171  yerr[i]->Fill(v.yError()*10000);
172  zerr[i]->Fill(v.zError()*10000);
173  xerrVsTrks[i]->Fill(weight, v.xError()*10000);
174  yerrVsTrks[i]->Fill(weight, v.yError()*10000);
175  zerrVsTrks[i]->Fill(weight, v.zError()*10000);
176 
177  nans[i]->Fill(1.,edm::isNotFinite(v.position().x())*1.);
178  nans[i]->Fill(2.,edm::isNotFinite(v.position().y())*1.);
179  nans[i]->Fill(3.,edm::isNotFinite(v.position().z())*1.);
180 
181  int index = 3;
182  for (int k = 0; k != 3; k++) {
183  for (int j = k; j != 3; j++) {
184  index++;
185  nans[i]->Fill(index*1., edm::isNotFinite(v.covariance(k, j))*1.);
186  // in addition, diagonal element must be positive
187  if (j == k && v.covariance(k, j) < 0) {
188  nans[i]->Fill(index*1., 1.);
189  }
190  }
191  }
192  }
193 }
194 
195 
197 {
198 }
199 
200 
201 //define this as a plug-in
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:69
int i
Definition: DBlmapReader.cc:9
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
void vertexPlots(const reco::Vertex &v, const reco::BeamSpot &beamSpot, int i)
double zError() const
error on z
Definition: Vertex.h:105
void zerr(int)
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:61
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:110
const Point & position() const
position
Definition: Vertex.h:93
double dydz() const
dydz slope
Definition: BeamSpot.h:85
int iEvent
Definition: GenABIO.cc:243
bool isNotFinite(T x)
Definition: isFinite.h:10
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
double chi2() const
chi-squares
Definition: Vertex.h:82
int j
Definition: DBlmapReader.cc:9
float ChiSquaredProbability(double chiSquared, double nrDOF)
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:87
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
double dxdz() const
dxdz slope
Definition: BeamSpot.h:83
double ndof() const
Definition: Vertex.h:89
int k[5][pyjets_maxn]
DEFINE_FWK_MODULE(CosmicTrackingParticleSelector)
double xError() const
error on x
Definition: Vertex.h:101
bool isFake() const
Definition: Vertex.h:65
double sigmaZ() const
sigma z
Definition: BeamSpot.h:81
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:89
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
double y0() const
y coordinate
Definition: BeamSpot.h:67
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
int weight
Definition: histoStyle.py:50
double yError() const
error on y
Definition: Vertex.h:103
virtual void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup)
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:35
PrimaryVertexMonitor(const edm::ParameterSet &pSet)
double x0() const
x coordinate
Definition: BeamSpot.h:65