CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelRecHitModule.cc
Go to the documentation of this file.
6 // STL
7 #include <vector>
8 #include <memory>
9 #include <string>
10 #include <iostream>
11 #include <stdlib.h>
12 
13 // Data Formats
18 //
19 // Constructors
20 //
24  id_(id)
25 {
26 }
27 
28 //
29 // Destructor
30 //
32 //
33 // Book histograms
34 //
36  bool twoD, bool reducedSet) {
37 
38  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
39  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
40  bool isHalfModule = false;
41  if(barrel){
42  isHalfModule = PixelBarrelName(DetId(id_)).isHalfModule();
43  }
44 
45  std::string hid;
46  // Get collection name and instantiate Histo Id builder
47  edm::InputTag src = iConfig.getParameter<edm::InputTag>( "src" );
48  // Get DQM interface
49  DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
50 
51 
52  if(type==0){
53  SiPixelHistogramId* theHistogramId = new SiPixelHistogramId( src.label() );
54  if(!reducedSet)
55  {
56  if(twoD){
57  // XYPosition
58  hid = theHistogramId->setHistoId("xypos",id_);
59  meXYPos_ = theDMBE->book2D(hid,"XY Position",100,-1.,1,100,-4,4);
60  meXYPos_->setAxisTitle("X Position",1);
61  meXYPos_->setAxisTitle("Y Position",2);
62  }
63  else{
64  // projections of XYPosition
65  hid = theHistogramId->setHistoId("xypos",id_);
66  meXYPos_px_ = theDMBE->book1D(hid+"_px","X Position",100,-1.,1);
67  meXYPos_px_->setAxisTitle("X Position",1);
68  meXYPos_py_ = theDMBE->book1D(hid+"_py","Y Position",100,-4,4);
69  meXYPos_py_->setAxisTitle("Y Position",1);
70  }
71  }
72  hid = theHistogramId->setHistoId("ClustX",id_);
73  meClustX_ = theDMBE->book1D(hid, "RecHit X size", 10, 0., 10.);
74  meClustX_->setAxisTitle("RecHit size X dimension", 1);
75  hid = theHistogramId->setHistoId("ClustY",id_);
76  meClustY_ = theDMBE->book1D(hid, "RecHit Y size", 15, 0., 15.);
77  meClustY_->setAxisTitle("RecHit size Y dimension", 1);
78 
79  hid = theHistogramId->setHistoId("ErrorX",id_);
80  meErrorX_ = theDMBE->book1D(hid, "RecHit error X", 100,0.,0.02);
81  meErrorX_->setAxisTitle("RecHit error X", 1);
82  hid = theHistogramId->setHistoId("ErrorY",id_);
83  meErrorY_ = theDMBE->book1D(hid, "RecHit error Y", 100,0.,0.02);
84  meErrorY_->setAxisTitle("RecHit error Y", 1);
85 
86  hid = theHistogramId->setHistoId("nRecHits",id_);
87  menRecHits_ = theDMBE->book1D(hid, "# of rechits in this module", 8, 0, 8);
88  menRecHits_->setAxisTitle("number of rechits",1);
89  delete theHistogramId;
90  }
91 
92  if(type==1 && barrel){
93  uint32_t DBladder = PixelBarrelName(DetId(id_)).ladderName();
94  char sladder[80]; sprintf(sladder,"Ladder_%02i",DBladder);
95  hid = src.label() + "_" + sladder;
96  if(isHalfModule) hid += "H";
97  else hid += "F";
98  if(!reducedSet)
99  {
100  if(twoD){
101  meXYPosLad_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
102  meXYPosLad_->setAxisTitle("X Position",1);
103  meXYPosLad_->setAxisTitle("Y Position",2);
104  }
105  else{
106  // projections of XYPosition
107  meXYPosLad_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
108  meXYPosLad_px_->setAxisTitle("X Position",1);
109  meXYPosLad_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
110  meXYPosLad_py_->setAxisTitle("Y Position",1);
111  }
112  }
113  meClustXLad_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
114  meClustXLad_->setAxisTitle("RecHit size X dimension", 1);
115  meClustYLad_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
116  meClustYLad_->setAxisTitle("RecHit size Y dimension", 1);
117  meErrorXLad_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
118  meErrorXLad_->setAxisTitle("RecHit error X", 1);
119  meErrorYLad_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
120  meErrorYLad_->setAxisTitle("RecHit error Y", 1);
121  menRecHitsLad_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
122  menRecHitsLad_->setAxisTitle("number of rechits",1);
123 
124  }
125 
126  if(type==2 && barrel){
127 
128  uint32_t DBlayer = PixelBarrelName(DetId(id_)).layerName();
129  char slayer[80]; sprintf(slayer,"Layer_%i",DBlayer);
130  hid = src.label() + "_" + slayer;
131 
132  if(!reducedSet)
133  {
134  if(twoD){
135  meXYPosLay_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
136  meXYPosLay_->setAxisTitle("X Position",1);
137  meXYPosLay_->setAxisTitle("Y Position",2);
138  }
139  else{
140  // projections of XYPosition
141  meXYPosLay_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
142  meXYPosLay_px_->setAxisTitle("X Position",1);
143  meXYPosLay_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
144  meXYPosLay_py_->setAxisTitle("Y Position",1);
145  }
146  }
147 
148  meClustXLay_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
149  meClustXLay_->setAxisTitle("RecHit size X dimension", 1);
150  meClustYLay_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
151  meClustYLay_->setAxisTitle("RecHit size Y dimension", 1);
152  meErrorXLay_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
153  meErrorXLay_->setAxisTitle("RecHit error X", 1);
154  meErrorYLay_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
155  meErrorYLay_->setAxisTitle("RecHit error Y", 1);
156  menRecHitsLay_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
157  menRecHitsLay_->setAxisTitle("number of rechits",1);
158 
159  }
160 
161  if(type==3 && barrel){
162  uint32_t DBmodule = PixelBarrelName(DetId(id_)).moduleName();
163  char smodule[80]; sprintf(smodule,"Ring_%i",DBmodule);
164  hid = src.label() + "_" + smodule;
165 
166  if(!reducedSet)
167  {
168  if(twoD){
169  meXYPosPhi_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
170  meXYPosPhi_->setAxisTitle("X Position",1);
171  meXYPosPhi_->setAxisTitle("Y Position",2);
172  }
173  else{
174  // projections of XYPosition
175  meXYPosPhi_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
176  meXYPosPhi_px_->setAxisTitle("X Position",1);
177  meXYPosPhi_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
178  meXYPosPhi_py_->setAxisTitle("Y Position",1);
179  }
180  }
181  meClustXPhi_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
182  meClustXPhi_->setAxisTitle("RecHit size X dimension", 1);
183  meClustYPhi_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
184  meClustYPhi_->setAxisTitle("RecHit size Y dimension", 1);
185  meErrorXPhi_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
186  meErrorXPhi_->setAxisTitle("RecHit error X", 1);
187  meErrorYPhi_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
188  meErrorYPhi_->setAxisTitle("RecHit error Y", 1);
189  menRecHitsPhi_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
190  menRecHitsPhi_->setAxisTitle("number of rechits",1);
191 
192  }
193 
194  if(type==4 && endcap){
195  uint32_t blade= PixelEndcapName(DetId(id_)).bladeName();
196 
197  char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
198  hid = src.label() + "_" + sblade;
199 // meXYPosBlade_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
200 // meXYPosBlade_->setAxisTitle("X Position",1);
201 // meXYPosBlade_->setAxisTitle("Y Position",2);
202 
203  meClustXBlade_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
204  meClustXBlade_->setAxisTitle("RecHit size X dimension", 1);
205  meClustYBlade_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
206  meClustYBlade_->setAxisTitle("RecHit size Y dimension", 1);
207  meErrorXBlade_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
208  meErrorXBlade_->setAxisTitle("RecHit error X", 1);
209  meErrorYBlade_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
210  meErrorYBlade_->setAxisTitle("RecHit error Y", 1);
211  menRecHitsBlade_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
212  menRecHitsBlade_->setAxisTitle("number of rechits",1);
213 
214  }
215  if(type==5 && endcap){
216  uint32_t disk = PixelEndcapName(DetId(id_)).diskName();
217 
218  char sdisk[80]; sprintf(sdisk, "Disk_%i",disk);
219  hid = src.label() + "_" + sdisk;
220 // meXYPosDisk_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
221 // meXYPosDisk_->setAxisTitle("X Position",1);
222 // meXYPosDisk_->setAxisTitle("Y Position",2);
223 
224  meClustXDisk_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
225  meClustXDisk_->setAxisTitle("RecHit size X dimension", 1);
226  meClustYDisk_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
227  meClustYDisk_->setAxisTitle("RecHit size Y dimension", 1);
228  meErrorXDisk_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
229  meErrorXDisk_->setAxisTitle("RecHit error X", 1);
230  meErrorYDisk_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
231  meErrorYDisk_->setAxisTitle("RecHit error Y", 1);
232  menRecHitsDisk_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
233  menRecHitsDisk_->setAxisTitle("number of rechits",1);
234 
235  }
236 
237  if(type==6 && endcap){
238  uint32_t panel= PixelEndcapName(DetId(id_)).pannelName();
240  char slab[80]; sprintf(slab, "Panel_%i_Ring_%i",panel, module);
241  hid = src.label() + "_" + slab;
242 
243  if(!reducedSet)
244  {
245  if(twoD){
246  meXYPosRing_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
247  meXYPosRing_->setAxisTitle("X Position",1);
248  meXYPosRing_->setAxisTitle("Y Position",2);
249  }
250  else{
251  // projections of XYPosition
252  meXYPosRing_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
253  meXYPosRing_px_->setAxisTitle("X Position",1);
254  meXYPosRing_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
255  meXYPosRing_py_->setAxisTitle("Y Position",1);
256  }
257  }
258  meClustXRing_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
259  meClustXRing_->setAxisTitle("RecHit size X dimension", 1);
260  meClustYRing_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
261  meClustYRing_->setAxisTitle("RecHit size Y dimension", 1);
262  meErrorXRing_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
263  meErrorXRing_->setAxisTitle("RecHit error X", 1);
264  meErrorYRing_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
265  meErrorYRing_->setAxisTitle("RecHit error Y", 1);
266  menRecHitsRing_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
267  menRecHitsRing_->setAxisTitle("number of rechits",1);
268 
269  }
270 
271 }
272 //
273 // Fill histograms
274 //
275 void SiPixelRecHitModule::fill(const float& rechit_x, const float& rechit_y,
276  const int& sizeX, const int& sizeY,
277  const float& lerr_x, const float& lerr_y,
278  bool modon, bool ladon, bool layon, bool phion,
279  bool bladeon, bool diskon, bool ringon,
280  bool twoD, bool reducedSet) {
281 
282  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
283  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
284 
285  //std::cout << rechit_x << " " << rechit_y << " " << sizeX << " " << sizeY << std::endl;
286  if(modon){
287  meClustX_->Fill(sizeX);
288  meClustY_->Fill(sizeY);
289  meErrorX_->Fill(lerr_x);
290  meErrorY_->Fill(lerr_y);
291  if(!reducedSet)
292  {
293  if(twoD) meXYPos_->Fill(rechit_x, rechit_y);
294  else {
295  meXYPos_px_->Fill(rechit_x);
296  meXYPos_py_->Fill(rechit_y);
297  }
298  }
299  }
300  //std::cout<<"number of detector units="<<numberOfDetUnits<<std::endl;
301 
302  if(ladon && barrel){
303  meClustXLad_->Fill(sizeX);
304  meClustYLad_->Fill(sizeY);
305  meErrorXLad_->Fill(lerr_x);
306  meErrorYLad_->Fill(lerr_y);
307  if(!reducedSet)
308  {
309  if(twoD) meXYPosLad_->Fill(rechit_x, rechit_y);
310  else{
311  meXYPosLad_px_->Fill(rechit_x);
312  meXYPosLad_py_->Fill(rechit_y);
313  }
314  }
315  }
316 
317  if(layon && barrel){
318  meClustXLay_->Fill(sizeX);
319  meClustYLay_->Fill(sizeY);
320  meErrorXLay_->Fill(lerr_x);
321  meErrorYLay_->Fill(lerr_y);
322  if(!reducedSet)
323  {
324  if(twoD) meXYPosLay_->Fill(rechit_x, rechit_y);
325  else{
326  meXYPosLay_px_->Fill(rechit_x);
327  meXYPosLay_py_->Fill(rechit_y);
328  }
329  }
330  }
331 
332  if(phion && barrel){
333  meClustXPhi_->Fill(sizeX);
334  meClustYPhi_->Fill(sizeY);
335  meErrorXPhi_->Fill(lerr_x);
336  meErrorYPhi_->Fill(lerr_y);
337  if(!reducedSet)
338  {
339  if(twoD) meXYPosPhi_->Fill(rechit_x, rechit_y);
340  else{
341  meXYPosPhi_px_->Fill(rechit_x);
342  meXYPosPhi_py_->Fill(rechit_y);
343  }
344  }
345  }
346 
347  if(bladeon && endcap){
348  //meXYPosBlade_->Fill(rechit_x, rechit_y);
349  meClustXBlade_->Fill(sizeX);
350  meClustYBlade_->Fill(sizeY);
351  meErrorXBlade_->Fill(lerr_x);
352  meErrorYBlade_->Fill(lerr_y);
353  }
354 
355  if(diskon && endcap){
356  //meXYPosDisk_->Fill(rechit_x, rechit_y);
357  meClustXDisk_->Fill(sizeX);
358  meClustYDisk_->Fill(sizeY);
359  meErrorXDisk_->Fill(lerr_x);
360  meErrorYDisk_->Fill(lerr_y);
361  }
362 
363  if(ringon && endcap){
364  meClustXRing_->Fill(sizeX);
365  meClustYRing_->Fill(sizeY);
366  meErrorXRing_->Fill(lerr_x);
367  meErrorYRing_->Fill(lerr_y);
368  if(!reducedSet)
369  {
370  if(twoD) meXYPosRing_->Fill(rechit_x, rechit_y);
371  else{
372  meXYPosRing_px_->Fill(rechit_x);
373  meXYPosRing_py_->Fill(rechit_y);
374  }
375  }
376  }
377 }
378 
379 void SiPixelRecHitModule::nfill(const int& nrec, bool modon, bool ladon, bool layon, bool phion, bool bladeon, bool diskon, bool ringon) {
380 
381  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
382  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
383 
384  if(modon) menRecHits_->Fill(nrec);
385  //barrel
386  if(ladon && barrel) menRecHitsLad_->Fill(nrec);
387  if(layon && barrel) menRecHitsLay_->Fill(nrec);
388  if(phion && barrel) menRecHitsPhi_->Fill(nrec);
389  //endcap
390  if(bladeon && endcap) menRecHitsBlade_->Fill(nrec);
391  if(diskon && endcap) menRecHitsDisk_->Fill(nrec);
392  if(ringon && endcap) menRecHitsRing_->Fill(nrec);
393 }
int plaquetteName() const
plaquetteId (in pannel)
MonitorElement * menRecHitsLay_
type
Definition: HCALResponse.h:22
T getParameter(std::string const &) const
MonitorElement * meXYPosPhi_px_
MonitorElement * meClustYRing_
MonitorElement * meClustXLay_
MonitorElement * meXYPosPhi_
void nfill(const int &nrec, bool modon=true, bool ladon=false, bool layon=false, bool phion=false, bool bladeon=false, bool diskon=false, bool ringon=false)
MonitorElement * meXYPosPhi_py_
MonitorElement * meErrorXDisk_
MonitorElement * meClustYPhi_
SiPixelRecHitModule()
Default constructor.
int moduleName() const
module id (index in z)
MonitorElement * meXYPosLay_py_
MonitorElement * meClustXDisk_
MonitorElement * meClustYDisk_
MonitorElement * meErrorY_
MonitorElement * meErrorYRing_
~SiPixelRecHitModule()
Destructor.
MonitorElement * meClustXLad_
MonitorElement * meErrorYPhi_
MonitorElement * meXYPosRing_px_
MonitorElement * meXYPosLay_
MonitorElement * meErrorXLad_
MonitorElement * meClustX_
void Fill(long long x)
MonitorElement * meClustXBlade_
MonitorElement * meClustXPhi_
MonitorElement * meXYPosLad_px_
int bladeName() const
blade id
MonitorElement * meErrorXBlade_
bool isHalfModule() const
full or half module
MonitorElement * meXYPosRing_
MonitorElement * meErrorYLad_
MonitorElement * menRecHitsRing_
MonitorElement * meXYPosLay_px_
MonitorElement * meClustXRing_
MonitorElement * meXYPos_py_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
MonitorElement * meErrorYLay_
MonitorElement * meXYPos_
MonitorElement * menRecHits_
std::string setHistoId(std::string variable, uint32_t &rawId)
Set Histogram Id.
Definition: DetId.h:20
int ladderName() const
ladder id (index in phi)
MonitorElement * meXYPosLad_
int layerName() const
layer id
MonitorElement * meErrorYDisk_
MonitorElement * menRecHitsDisk_
std::string const & label() const
Definition: InputTag.h:25
MonitorElement * meErrorX_
MonitorElement * menRecHitsPhi_
void book(const edm::ParameterSet &iConfig, int type=0, bool twoD=true, bool reducedSet=false)
Book histograms.
int pannelName() const
pannel id
MonitorElement * menRecHitsLad_
MonitorElement * meClustYLad_
MonitorElement * meClustY_
MonitorElement * meXYPosRing_py_
int diskName() const
disk id
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Definition: vlib.h:209
MonitorElement * meXYPosLad_py_
MonitorElement * meErrorXLay_
MonitorElement * menRecHitsBlade_
MonitorElement * meClustYBlade_
void fill(const float &rechit_x, const float &rechit_y, const int &sizeX, const int &sizeY, const float &lerr_x, const float &lerr_y, bool modon=true, bool ladon=false, bool layon=false, bool phion=false, bool bladeon=false, bool diskon=false, bool ringon=false, bool twoD=true, bool reducedSet=false)
Fill histograms.
MonitorElement * meClustYLay_
MonitorElement * meXYPos_px_
MonitorElement * meErrorXRing_
MonitorElement * meErrorYBlade_
tuple src
Definition: align_tpl.py:87
MonitorElement * meErrorXPhi_