CMS 3D CMS Logo

SiPixelPhase1RecHits_cfi.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 from DQMServices.Core.DQMEDHarvester import DQMEDHarvester
4 import DQM.SiPixelPhase1Common.TriggerEventFlag_cfi as trigger
5 
6 SiPixelPhase1RecHitsNRecHits = DefaultHistoTrack.clone(
7  name = "rechits",
8  title = "RecHits",
9  range_min = 0, range_max = 30, range_nbins = 30,
10  xlabel = "rechits",
11  dimensions = 0,
12  specs = VPSet(
13 
14  StandardSpecificationTrend_Num,
15  Specification().groupBy("PXBarrel/Event")
16  .reduce("COUNT")
17  .groupBy("PXBarrel")
18  .save(nbins=100, xmin=0, xmax=5000),
19 
20  Specification().groupBy("PXForward/Event")
21  .reduce("COUNT")
22  .groupBy("PXForward")
23  .save(nbins=100, xmin=0, xmax=5000),
24 
25  Specification().groupBy("PXAll/Event")
26  .reduce("COUNT")
27  .groupBy("PXAll")
28  .save(nbins=100, xmin=0, xmax=5000)
29 
30  )
31 )
32 
33 SiPixelPhase1RecHitsClustX = DefaultHistoTrack.clone(
34  name = "clustersize_x",
35  title = "Cluster Size X (OnTrack)",
36  range_min = 0, range_max = 50, range_nbins = 50,
37  xlabel = "size[pixels]",
38  dimensions = 1,
39  specs = VPSet(
40  StandardSpecification2DProfile
41  )
42 )
43 
44 SiPixelPhase1RecHitsClustY = SiPixelPhase1RecHitsClustX.clone(
45  name = "clustersize_y",
46  title = "Cluster Size Y (OnTrack)",
47  xlabel = "size[pixels]"
48 )
49 
50 SiPixelPhase1RecHitsErrorX = DefaultHistoTrack.clone(
51  enabled=False,
52  name = "rechiterror_x",
53  title = "RecHit Error in X-direction",
54  range_min = 0, range_max = 0.02, range_nbins = 100,
55  xlabel = "X error",
56  dimensions = 1,
57  specs = VPSet(
58  StandardSpecification2DProfile
59  )
60 )
61 
62 SiPixelPhase1RecHitsErrorY = SiPixelPhase1RecHitsErrorX.clone(
63  enabled=False,
64  name = "rechiterror_y",
65  title = "RecHit Error in Y-direction",
66  xlabel = "Y error"
67 )
68 
69 SiPixelPhase1RecHitsPosition = DefaultHistoTrack.clone(
70  enabled = False,
71  name = "rechit_pos",
72  title = "Position of RecHits on Module",
73  range_min = -1, range_max = 1, range_nbins = 100,
74  range_y_min = -4, range_y_max = 4, range_y_nbins = 100,
75  xlabel = "x offset",
76  ylabel = "y offset",
77  dimensions = 2,
78  specs = VPSet(
79  Specification(PerModule).groupBy("PXBarrel/PXLayer/DetId").save(),
80  Specification(PerModule).groupBy("PXForward/PXDisk/DetId").save(),
81  )
82 )
83 
84 SiPixelPhase1RecHitsProb = DefaultHistoTrack.clone(
85  name = "clusterprob",
86  title = "Cluster Probability",
87  xlabel = "log_10(Pr)",
88  range_min = -10, range_max = 1, range_nbins = 50,
89  dimensions = 1,
90  specs = VPSet(
91 
92  Specification().groupBy("PXBarrel/PXLayer").saveAll(),
93  Specification().groupBy("PXForward/PXDisk").saveAll(),
94  StandardSpecification2DProfile,
95 
96  Specification().groupBy("PXBarrel/LumiBlock")
97  .reduce("MEAN")
98  .groupBy("PXBarrel", "EXTEND_X")
99  .save(),
100 
101  Specification().groupBy("PXForward/LumiBlock")
102  .reduce("MEAN")
103  .groupBy("PXForward", "EXTEND_X")
104  .save(),
105 
106  Specification(PerLayer1D).groupBy("PXBarrel/Shell/PXLayer").save(),
107  Specification(PerLayer1D).groupBy("PXForward/HalfCylinder/PXRing/PXDisk").save()
108  )
109 )
110 
111 SiPixelPhase1RecHitsOnEdge = DefaultHistoTrack.clone(
112  name = "onedge",
113  title = "OnEdge",
114  range_min = 0, range_max = 30, range_nbins = 30,
115  xlabel = "onedge",
116  dimensions = 0,
117  specs = VPSet(
118 
119  StandardSpecificationTrend_Num,
120  StandardSpecificationOccupancy,
121  Specification().groupBy("PXBarrel/PXLayer/Event")
122  .reduce("COUNT")
123  .groupBy("PXBarrel/PXLayer/LumiBlock")
124  .reduce("MEAN")
125  .groupBy("PXBarrel/PXLayer","EXTEND_X")
126  .save(),
127  Specification().groupBy("PXForward/PXDisk/Event")
128  .reduce("COUNT")
129  .groupBy("PXForward/PXDisk/LumiBlock")
130  .reduce("MEAN")
131  .groupBy("PXForward/PXDisk","EXTEND_X")
132  .save(),
133  Specification().groupBy("PXBarrel/PXLayer/Event")
134  .reduce("COUNT")
135  .groupBy("PXBarrel/PXLayer/")
136  .save(nbins=100, xmin=0, xmax=500),
137  Specification().groupBy("PXBarrel/Event")
138  .reduce("COUNT")
139  .groupBy("PXBarrel")
140  .save(nbins=100, xmin=0, xmax=500),
141 
142  Specification().groupBy("PXForward/PXDisk/Event")
143  .reduce("COUNT")
144  .groupBy("PXForward/PXDisk/")
145  .save(nbins=100, xmin=0, xmax=500),
146  Specification().groupBy("PXForward/PXDisk/Event")
147  .reduce("COUNT")
148  .groupBy("PXForward")
149  .save(nbins=100, xmin=0, xmax=500),
150 
151  Specification().groupBy("PXAll/Event")
152  .reduce("COUNT")
153  .groupBy("PXAll")
154  .save(nbins=100, xmin=0, xmax=500)
155 
156  )
157 )
158 
159 SiPixelPhase1RecHitsOtherBad = DefaultHistoTrack.clone(
160  name = "otherbad",
161  title = "OtherBad",
162  range_min = 0, range_max = 30, range_nbins = 30,
163  xlabel = "otherbad",
164  dimensions = 0,
165  specs = VPSet(
166 
167  StandardSpecificationTrend_Num,
168  StandardSpecificationOccupancy,
169  Specification().groupBy("PXBarrel/PXLayer/Event")
170  .reduce("COUNT")
171  .groupBy("PXBarrel/PXLayer/LumiBlock")
172  .reduce("MEAN")
173  .groupBy("PXBarrel/PXLayer","EXTEND_X")
174  .save(),
175  Specification().groupBy("PXForward/PXDisk/Event")
176  .reduce("COUNT")
177  .groupBy("PXForward/PXDisk/LumiBlock")
178  .reduce("MEAN")
179  .groupBy("PXForward/PXDisk","EXTEND_X")
180  .save(),
181  Specification().groupBy("PXBarrel/PXLayer/Event")
182  .reduce("COUNT")
183  .groupBy("PXBarrel/PXLayer/")
184  .save(nbins=100, xmin=0, xmax=500),
185  Specification().groupBy("PXBarrel/Event")
186  .reduce("COUNT")
187  .groupBy("PXBarrel")
188  .save(nbins=100, xmin=0, xmax=500),
189 
190  Specification().groupBy("PXForward/PXDisk/Event")
191  .reduce("COUNT")
192  .groupBy("PXForward/PXDisk/")
193  .save(nbins=100, xmin=0, xmax=500),
194  Specification().groupBy("PXForward/PXDisk/Event")
195  .reduce("COUNT")
196  .groupBy("PXForward")
197  .save(nbins=100, xmin=0, xmax=500),
198 
199  Specification().groupBy("PXAll/Event")
200  .reduce("COUNT")
201  .groupBy("PXAll")
202  .save(nbins=100, xmin=0, xmax=500)
203  )
204 )
205 
206 SiPixelPhase1RecHitsConf = cms.VPSet(
207  SiPixelPhase1RecHitsNRecHits,
208  SiPixelPhase1RecHitsClustX,
209  SiPixelPhase1RecHitsClustY,
210  SiPixelPhase1RecHitsErrorX,
211  SiPixelPhase1RecHitsErrorY,
212  SiPixelPhase1RecHitsPosition,
213  SiPixelPhase1RecHitsProb,
214  SiPixelPhase1RecHitsOnEdge,
215  SiPixelPhase1RecHitsOtherBad,
216 )
217 
218 from DQMServices.Core.DQMEDAnalyzer import DQMEDAnalyzer
219 SiPixelPhase1RecHitsAnalyzer = DQMEDAnalyzer('SiPixelPhase1RecHits',
220  src = cms.InputTag("generalTracks"),
221  histograms = SiPixelPhase1RecHitsConf,
222  geometry = SiPixelPhase1Geometry,
223  onlyValidHits = cms.bool(False),
224  triggerflags = trigger.SiPixelPhase1Triggers,
225  VertexCut = cms.untracked.bool(True)
226 )
227 
228 SiPixelPhase1RecHitsHarvester = DQMEDHarvester("SiPixelPhase1Harvester",
229  histograms = SiPixelPhase1RecHitsConf,
230  geometry = SiPixelPhase1Geometry
231 )
T reduce(std::vector< T > x, Op op)
Definition: conifer.h:31
save
Definition: cuy.py:1164