Belle II Software development
AlignDQMModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <alignment/modules/AlignmentDQM/AlignDQMModule.h>
10#include <alignment/modules/AlignmentDQM/AlignDQMEventProcessor.h>
11#include <tracking/dqmUtils/HistogramFactory.h>
12
13#include <TDirectory.h>
14
15using namespace Belle2;
16using namespace Belle2::HistogramFactory;
17using namespace std;
18using boost::format;
19
20//-----------------------------------------------------------------
21// Register the Module
22//-----------------------------------------------------------------
23
24REG_MODULE(AlignDQM);
25
26//-----------------------------------------------------------------
27// Implementation
28//-----------------------------------------------------------------
29
31{
32 setDescription("DQM of Alignment for off line "
33 "residuals per sensor, layer, "
34 "keep also On-Line DQM from tracking: "
35 "their momentum, "
36 "Number of hits in tracks, "
37 "Number of tracks. "
38 );
39}
40
41//------------------------------------------------------------------
42// Function to define histograms
43//-----------------------------------------------------------------
44
46{
48
49 if (VXD::GeoCache::getInstance().getGeoTools()->getNumberOfLayers() == 0)
50 B2WARNING("Missing geometry for VXD.");
51
52 TDirectory* originalDirectory = gDirectory;
53 TDirectory* alignmentDirectory = originalDirectory->mkdir("AlignmentDQM", "", true);
54 TDirectory* sensorsDirectory = originalDirectory->mkdir("AlignmentDQMSensors", "", true);
55 TDirectory* layersDirectory = originalDirectory->mkdir("AlignmentDQMLayers", "", true);
56
57 alignmentDirectory->cd();
59 DefineHits();
67
68 sensorsDirectory->cd();
70
71 layersDirectory->cd();
73
74 originalDirectory->cd();
75
76 for (auto change : m_histogramParameterChanges)
77 ProcessHistogramParameterChange(get<0>(change), get<1>(change), get<2>(change));
78}
79
81{
84 return;
85
87
88 eventProcessor.Run();
89}
90
92{
94 auto gTools = geo.getGeoTools();
95
96 for (VxdID layer : geo.getLayers()) {
97 int layerIndex = gTools->getLayerIndex(layer.getLayerNumber());
98 m_ResMeanUPhiThetaLayer[layerIndex]->Divide(m_ResMeanPhiThetaLayerCounts[layerIndex]);
99 m_ResMeanVPhiThetaLayer[layerIndex]->Divide(m_ResMeanPhiThetaLayerCounts[layerIndex]);
100
101 ComputeMean(m_ResMeanUPhiLayer[layerIndex], m_ResUPhiLayer[layerIndex]);
102 ComputeMean(m_ResMeanVPhiLayer[layerIndex], m_ResVPhiLayer[layerIndex]);
103
104 ComputeMean(m_ResMeanUThetaLayer[layerIndex], m_ResUThetaLayer[layerIndex]);
105 ComputeMean(m_ResMeanVThetaLayer[layerIndex], m_ResVThetaLayer[layerIndex]);
106 }
107
108 for (int sensorIndex = 0; sensorIndex < gTools->getNumberOfSensors(); sensorIndex++) {
109 m_ResMeanUPosUVSens[sensorIndex]->Divide(m_ResMeanPosUVSensCounts[sensorIndex]);
110 m_ResMeanVPosUVSens[sensorIndex]->Divide(m_ResMeanPosUVSensCounts[sensorIndex]);
111
112 ComputeMean(m_ResMeanUPosUSens[sensorIndex], m_ResUPosUSens[sensorIndex]);
113 ComputeMean(m_ResMeanVPosUSens[sensorIndex], m_ResVPosUSens[sensorIndex]);
114
115 ComputeMean(m_ResMeanUPosVSens[sensorIndex], m_ResUPosVSens[sensorIndex]);
116 ComputeMean(m_ResMeanVPosVSens[sensorIndex], m_ResVPosVSens[sensorIndex]);
117 }
118}
119
121{
122 TDirectory* originalDirectory = gDirectory;
123
124 TDirectory* helixParameters = originalDirectory->mkdir("HelixPars", "", true);
125 TDirectory* helixCorrelations = originalDirectory->mkdir("HelixCorrelations", "", true);
126
127 int iZ0Range = 100;
128 double fZ0Range = 10.0; // Half range in cm
129 int iD0Range = 100;
130 double fD0Range = 1.0; // Half range in cm
131 int iMomRangeBig = 600;
132 int iMomRangeSmall = 60;
133 double fMomRange = 6.0;
134 int iPhiRange = 180;
135 double fPhiRange = 180.0; // Half range in deg
136 int iLambdaRange = 100;
137 double fLambdaRange = 4.0;
138 int iOmegaRange = 100;
139 double fOmegaRange = 0.1;
140
141 auto phi = Axis(iPhiRange, -fPhiRange, fPhiRange, "#phi [deg]");
142 auto D0 = Axis(iD0Range, -fD0Range, fD0Range, "d0 [cm]");
143 auto Z0 = Axis(iZ0Range, -fZ0Range, fZ0Range, "z0 [cm]");
144 auto tanLambda = Axis(iLambdaRange, -fLambdaRange, fLambdaRange, "Tan Lambda");
145 auto omega = Axis(iOmegaRange, -fOmegaRange, fOmegaRange, "Omega");
146 auto momentumBig = Axis(2 * iMomRangeBig, 0.0, fMomRange, "Momentum");
147 auto momentumSmall = Axis(2 * iMomRangeSmall, 0.0, fMomRange, "Momentum");
148
149 auto factory = Factory(this);
150
151 helixParameters->cd();
152
153 factory.yTitleDefault("Arb. Units");
154
155 m_Z0 = factory.xAxis(Z0).CreateTH1F("Z0", "z0 - the z coordinate of the perigee (beam spot position)");
156 m_D0 = factory.xAxis(D0).CreateTH1F("D0", "d0 - the signed distance to the IP in the r-phi plane");
157 m_Phi = factory.xAxis(phi).CreateTH1F("Phi",
158 "Phi - angle of the transverse momentum in the r-phi plane, with CDF naming convention");
159 m_Omega = factory.xAxis(omega).CreateTH1F("Omega",
160 "Omega - the curvature of the track. It's sign is defined by the charge of the particle");
161 m_TanLambda = factory.xAxis(tanLambda).CreateTH1F("TanLambda", "TanLambda - the slope of the track in the r-z plane");
162 m_MomPt = factory.xAxis(momentumBig).yTitle("counts").CreateTH1F("TrackMomentumPt", "Track Momentum pT");
163
164 helixCorrelations->cd();
165
166 factory.zTitleDefault("Arb. Units");
167
168 m_PhiD0 = factory.xAxis(phi).yAxis(D0).CreateTH2F("PhiD0",
169 "Phi - angle of the transverse momentum in the r-phi plane vs. d0 - signed distance to the IP in r-phi");
170 m_PhiZ0 = factory.xAxis(phi).yAxis(Z0).CreateTH2F("PhiZ0",
171 "Phi - angle of the transverse momentum in the r-phi plane vs. z0 of the perigee (to see primary vertex shifts along R or z)");
172 m_PhiMomPt = factory.xAxis(phi).yAxis(momentumSmall).CreateTH2F("PhiMomPt",
173 "Phi - angle of the transverse momentum in the r-phi plane vs. Track momentum Pt");
174 m_PhiOmega = factory.xAxis(phi).yAxis(omega).CreateTH2F("PhiOmega",
175 "Phi - angle of the transverse momentum in the r-phi plane vs. Omega - the curvature of the track");
176 m_PhiTanLambda = factory.xAxis(phi).yAxis(tanLambda).CreateTH2F("PhiTanLambda",
177 "dPhi - angle of the transverse momentum in the r-phi plane vs. TanLambda - the slope of the track in the r-z plane");
178 m_D0Z0 = factory.xAxis(D0).yAxis(Z0).CreateTH2F("D0Z0",
179 "d0 - signed distance to the IP in r-phi vs. z0 of the perigee (to see primary vertex shifts along R or z)");
180 m_D0MomPt = factory.xAxis(D0).yAxis(momentumSmall).CreateTH2F("D0MomPt",
181 "d0 - signed distance to the IP in r-phi vs. Track momentum Pt");
182 m_D0Omega = factory.xAxis(D0).yAxis(omega).CreateTH2F("D0Omega",
183 "d0 - signed distance to the IP in r-phi vs. Omega - the curvature of the track");
184 m_D0TanLambda = factory.xAxis(D0).yAxis(tanLambda).CreateTH2F("D0TanLambda",
185 "d0 - signed distance to the IP in r-phi vs. TanLambda - the slope of the track in the r-z plane");
186 m_Z0MomPt = factory.xAxis(Z0).yAxis(momentumSmall).CreateTH2F("Z0MomPt",
187 "z0 - the z0 coordinate of the perigee vs. Track momentum Pt");
188 m_Z0Omega = factory.xAxis(Z0).yAxis(omega).CreateTH2F("Z0Omega",
189 "z0 - the z0 coordinate of the perigee vs. Omega - the curvature of the track");
190 m_Z0TanLambda = factory.xAxis(Z0).yAxis(tanLambda).CreateTH2F("Z0TanLambda",
191 "z0 - the z0 coordinate of the perigee vs. TanLambda - the slope of the track in the r-z plane");
192 m_MomPtOmega = factory.xAxis(momentumSmall).yAxis(omega).CreateTH2F("MomPtOmega",
193 "Track momentum Pt vs. Omega - the curvature of the track");
194 m_MomPtTanLambda = factory.xAxis(momentumSmall).yAxis(tanLambda).CreateTH2F("MomPtTanLambda",
195 "Track momentum Pt vs. TanLambda - the slope of the track in the r-z plane");
196 m_OmegaTanLambda = factory.xAxis(omega).yAxis(tanLambda).CreateTH2F("OmegaTanLambda",
197 "Omega - the curvature of the track vs. TanLambda - the slope of the track in the r-z plane");
198
199 originalDirectory->cd();
200}
201
203{
204 TDirectory* originalDirectory = gDirectory;
205
206 TDirectory* resMeanUPosUV = originalDirectory->mkdir("ResidMeanUPositUV", "", true);
207 TDirectory* resMeanVPosUV = originalDirectory->mkdir("ResidMeanVPositUV", "", true);
208 TDirectory* resMeanPosUVCounts = originalDirectory->mkdir("ResidMeanPositUVCounts", "", true);
209 TDirectory* resMeanUPosU = originalDirectory->mkdir("ResidMeanUPositU", "", true);
210 TDirectory* resMeanVPosU = originalDirectory->mkdir("ResidMeanVPositU", "", true);
211 TDirectory* resMeanUPosV = originalDirectory->mkdir("ResidMeanUPositV", "", true);
212 TDirectory* resMeanVPosV = originalDirectory->mkdir("ResidMeanVPositV", "", true);
213 TDirectory* resUPosU = originalDirectory->mkdir("ResidUPositU", "", true);
214 TDirectory* resVPosU = originalDirectory->mkdir("ResidVPositU", "", true);
215 TDirectory* resUPosV = originalDirectory->mkdir("ResidUPositV", "", true);
216 TDirectory* resVPosV = originalDirectory->mkdir("ResidVPositV", "", true);
217 TDirectory* resids2D = originalDirectory->mkdir("Residuals2D", "", true);
218 TDirectory* resids1D = originalDirectory->mkdir("Residuals1D", "", true);
219
220 int iSizeBins = 20;
221 double fSizeMin = -50; // in mm
222 double fSizeMax = -fSizeMin;
223 double residualRange = 400; // in um
224
225 auto positionU = Axis(iSizeBins, fSizeMin, fSizeMax, "position U [mm]");
226 auto positionV = Axis(positionU).title("position V [mm]");
227 auto residualU = Axis(200, -residualRange, residualRange, "residual U [#mum]");
228 auto residualV = Axis(residualU).title("residual V [#mum]");
229
230 auto factory = Factory(this);
231
232 factory.xAxisDefault(positionU).yAxisDefault(positionV);
233
234 resMeanUPosUV->cd();
235 m_ResMeanUPosUVSens = factory.zTitle("residual U [#mum]").CreateSensorsTH2F(format("ResMeanUPosUVSens_%1%"),
236 format("Residual Mean U in Position UV, %1%"));
237 resMeanVPosUV->cd();
238 m_ResMeanVPosUVSens = factory.zTitle("residual V [#mum]").CreateSensorsTH2F(format("ResMeanVPosUVSens_%1%"),
239 format("Residual Mean V in Position UV, %1%"));
240 resMeanPosUVCounts->cd();
241 m_ResMeanPosUVSensCounts = factory.zTitle("counts").CreateSensorsTH2F(format("ResMeanPosUVCountsSens_%1%"),
242 format("Residual Mean Counts in Position UV, %1%"));
243 resMeanUPosU->cd();
244 m_ResMeanUPosUSens = factory.yTitle("residual mean U [#mum]").CreateSensorsTH1F(format("ResMeanUPosUSens_%1%"),
245 format("Residual Mean U in Position U, %1%"));
246 resMeanVPosU->cd();
247 m_ResMeanVPosUSens = factory.yTitle("residual mean V [#mum]").CreateSensorsTH1F(format("ResMeanVPosUSens_%1%"),
248 format("Residual Mean V in Position U, %1%"));
249
250 factory.xAxisDefault(positionV);
251
252 resMeanUPosV->cd();
253 m_ResMeanUPosVSens = factory.yTitle("residual mean U [#mum]").CreateSensorsTH1F(format("ResMeanUPosVSens_%1%"),
254 format("Residual Mean U in Position V, %1%"));
255 resMeanVPosV->cd();
256 m_ResMeanVPosVSens = factory.yTitle("residual mean V [#mum]").CreateSensorsTH1F(format("ResMeanVPosVSens_%1%"),
257 format("Residual Mean V in Position V, %1%"));
258
259 factory.xAxisDefault(positionU).zTitleDefault("counts");
260
261 resUPosU->cd();
262 m_ResUPosUSens = factory.yAxis(residualU).CreateSensorsTH2F(format("ResUPosUSensor_%1%"),
263 format("Residual U in Position U, %1%"));
264 resVPosU->cd();
265 m_ResVPosUSens = factory.yAxis(residualV).CreateSensorsTH2F(format("ResVPosUSensor_%1%"),
266 format("Residual V in Position U, %1%"));
267
268 factory.xAxisDefault(positionV);
269
270 resUPosV->cd();
271 m_ResUPosVSens = factory.yAxis(residualU).CreateSensorsTH2F(format("ResUPosVSensor_%1%"),
272 format("Residual U in Position V, %1%"));
273 resVPosV->cd();
274 m_ResVPosVSens = factory.yAxis(residualV).CreateSensorsTH2F(format("ResVPosVSensor_%1%"),
275 format("Residual V in Position V, %1%"));
276
277 resids2D->cd();
278 m_UBResidualsSensor = factory.xAxis(residualU).yAxis(residualV).CreateSensorsTH2F(format("UBResiduals_%1%"),
279 format("PXD Unbiased residuals for sensor %1%"));
280
281 factory.yTitleDefault("counts");
282
283 resids1D->cd();
284 m_UBResidualsSensorU = factory.xAxis(residualU).CreateSensorsTH1F(format("UBResidualsU_%1%"),
285 format("PXD Unbiased U residuals for sensor %1%"));
286 m_UBResidualsSensorV = factory.xAxis(residualV).CreateSensorsTH1F(format("UBResidualsV_%1%"),
287 format("PXD Unbiased V residuals for sensor %1%"));
288
289 originalDirectory->cd();
290}
291
293{
294 TDirectory* originalDirectory = gDirectory;
295
296 TDirectory* resMeanUPosUV = originalDirectory->mkdir("ResidLayerMeanUPositPhiTheta", "", true);
297 TDirectory* resMeanVPosUV = originalDirectory->mkdir("ResidLayerMeanVPositPhiTheta", "", true);
298 TDirectory* resMeanPosUVCounts = originalDirectory->mkdir("ResidLayerMeanPositPhiThetaCounts", "", true);
299 TDirectory* resMeanUPosU = originalDirectory->mkdir("ResidLayerMeanUPositPhi", "", true);
300 TDirectory* resMeanVPosU = originalDirectory->mkdir("ResidLayerMeanVPositPhi", "", true);
301 TDirectory* resMeanUPosV = originalDirectory->mkdir("ResidLayerMeanUPositTheta", "", true);
302 TDirectory* resMeanVPosV = originalDirectory->mkdir("ResidLayerMeanVPositTheta", "", true);
303 TDirectory* resUPosU = originalDirectory->mkdir("ResidLayerUPositPhi", "", true);
304 TDirectory* resVPosU = originalDirectory->mkdir("ResidLayerVPositPhi", "", true);
305 TDirectory* resUPosV = originalDirectory->mkdir("ResidLayerUPositTheta", "", true);
306 TDirectory* resVPosV = originalDirectory->mkdir("ResidLayerVPositTheta", "", true);
307
308 int iPhiGran = 90;
309 int iThetGran = iPhiGran / 2;
310 int iYResGran = 200;
311 double residualRange = 400; // in um
312
313 auto phi = Axis(iPhiGran, -180, 180, "Phi [deg]");
314 auto theta = Axis(iThetGran, 0, 180, "Theta [deg]");
315 auto residual = Axis(iYResGran, -residualRange, residualRange, "residual [#mum]");
316
317 auto factory = Factory(this);
318 factory.xAxisDefault(phi).yAxisDefault(theta).zTitleDefault("counts");
319
320 resMeanUPosUV->cd();
321 m_ResMeanUPhiThetaLayer = factory.zTitle("residual [#mum]").CreateLayersTH2F(format("ResMeanUPhiThetaLayer_%1%"),
322 format("Residuals Mean U in Phi Theta, Layer %1%"));
323 resMeanVPosUV->cd();
324 m_ResMeanVPhiThetaLayer = factory.zTitle("residual [#mum]").CreateLayersTH2F(format("ResMeanVPhiThetaLayer_%1%"),
325 format("Residuals Mean V in Phi Theta, Layer %1%"));
326 resMeanPosUVCounts->cd();
327 m_ResMeanPhiThetaLayerCounts = factory.CreateLayersTH2F(format("ResCounterPhiThetaLayer_%1%"),
328 format("Residuals counter in Phi Theta, Layer %1%"));
329
330 factory.yAxisDefault(residual);
331
332 resUPosU->cd();
333 m_ResUPhiLayer = factory.CreateLayersTH2F(format("ResUPhiLayer_%1%"), format("Residuals U in Phi, Layer %1%"));
334 resVPosU->cd();
335 m_ResVPhiLayer = factory.CreateLayersTH2F(format("ResVPhiLayer_%1%"), format("Residuals V in Phi, Layer %1%"));
336
337 factory.xAxisDefault(theta);
338
339 resUPosV->cd();
340 m_ResUThetaLayer = factory.CreateLayersTH2F(format("ResUThetaLayer_%1%"), format("Residuals U in Theta, Layer %1%"));
341 resVPosV->cd();
342 m_ResVThetaLayer = factory.CreateLayersTH2F(format("ResVThetaLayer_%1%"), format("Residuals V in Theta, Layer %1%"));
343 resMeanUPosV->cd();
344 m_ResMeanUThetaLayer = factory.CreateLayersTH1F(format("ResMeanUThetaLayer_%1%"),
345 format("Residuals Mean U in Theta, Layer %1%"));
346 resMeanVPosV->cd();
347 m_ResMeanVThetaLayer = factory.CreateLayersTH1F(format("ResMeanVThetaLayer_%1%"),
348 format("Residuals Mean V in Theta, Layer %1%"));
349
350 factory.xAxisDefault(phi).yTitleDefault("residual [#mum]");
351
352 resMeanUPosU->cd();
353 m_ResMeanUPhiLayer = factory.CreateLayersTH1F(format("ResMeanUPhiLayer_%1%"),
354 format("Residuals Mean U in Phi, Layer %1%"));
355 resMeanVPosU->cd();
356 m_ResMeanVPhiLayer = factory.CreateLayersTH1F(format("ResMeanVPhiLayer_%1%"),
357 format("Residuals Mean V in Phi, Layer %1%"));
358
359 originalDirectory->cd();
360}
361
363{
365
366 m_PhiZ0->Fill(tfr->getPhi0() / Unit::deg, tfr->getZ0());
367 m_PhiMomPt->Fill(tfr->getPhi0() / Unit::deg, tfr->getMomentum().Rho());
368 m_PhiOmega->Fill(tfr->getPhi0() / Unit::deg, tfr->getOmega());
369 m_PhiTanLambda->Fill(tfr->getPhi0() / Unit::deg, tfr->getTanLambda());
370 m_D0MomPt->Fill(tfr->getD0(), tfr->getMomentum().Rho());
371 m_D0Omega->Fill(tfr->getD0(), tfr->getOmega());
372 m_D0TanLambda->Fill(tfr->getD0(), tfr->getTanLambda());
373 m_Z0MomPt->Fill(tfr->getZ0(), tfr->getMomentum().Rho());
374 m_Z0Omega->Fill(tfr->getZ0(), tfr->getOmega());
375 m_Z0TanLambda->Fill(tfr->getZ0(), tfr->getTanLambda());
376 m_MomPtOmega->Fill(tfr->getMomentum().Rho(), tfr->getOmega());
377 m_MomPtTanLambda->Fill(tfr->getMomentum().Rho(), tfr->getTanLambda());
378 m_OmegaTanLambda->Fill(tfr->getOmega(), tfr->getTanLambda());
379}
380
381void AlignDQMModule::FillPositionSensors(ROOT::Math::XYZVector residual_um, ROOT::Math::XYZVector position, int sensorIndex)
382{
383 float positionU_mm = position.X() / Unit::mm;
384 float positionV_mm = position.Y() / Unit::mm;
385
386 m_ResMeanPosUVSensCounts[sensorIndex]->Fill(positionU_mm, positionV_mm);
387 m_ResMeanUPosUVSens[sensorIndex]->Fill(positionU_mm, positionV_mm, residual_um.X());
388 m_ResMeanVPosUVSens[sensorIndex]->Fill(positionU_mm, positionV_mm, residual_um.Y());
389 m_ResUPosUSens[sensorIndex]->Fill(positionU_mm, residual_um.X());
390 m_ResUPosVSens[sensorIndex]->Fill(positionV_mm, residual_um.X());
391 m_ResVPosUSens[sensorIndex]->Fill(positionU_mm, residual_um.Y());
392 m_ResVPosVSens[sensorIndex]->Fill(positionV_mm, residual_um.Y());
393}
394
395void AlignDQMModule::FillLayers(ROOT::Math::XYZVector residual_um, float phi_deg, float theta_deg, int layerIndex)
396{
397 m_ResMeanPhiThetaLayerCounts[layerIndex]->Fill(phi_deg, theta_deg);
398 m_ResMeanUPhiThetaLayer[layerIndex]->Fill(phi_deg, theta_deg, residual_um.X());
399 m_ResMeanVPhiThetaLayer[layerIndex]->Fill(phi_deg, theta_deg, residual_um.Y());
400 m_ResUPhiLayer[layerIndex]->Fill(phi_deg, residual_um.X());
401 m_ResVPhiLayer[layerIndex]->Fill(phi_deg, residual_um.Y());
402 m_ResUThetaLayer[layerIndex]->Fill(theta_deg, residual_um.X());
403 m_ResVThetaLayer[layerIndex]->Fill(theta_deg, residual_um.Y());
404}
405
406TH1F* AlignDQMModule::Create(string name, string title, int nbinsx, double xlow, double xup, string xTitle, string yTitle)
407{
408 return DQMHistoModuleBase::Create("Alig_" + name, title, nbinsx, xlow, xup, xTitle, yTitle);
409}
410
411TH2F* AlignDQMModule::Create(string name, string title, int nbinsx, double xlow, double xup, int nbinsy, double ylow, double yup,
412 string xTitle, string yTitle, string zTitle)
413{
414 return DQMHistoModuleBase::Create("Alig_" + name, title, nbinsx, xlow, xup, nbinsy, ylow, yup, xTitle, yTitle, zTitle);
415}
The purpose of this class is to process one event() in AlignDQMModule.
TH1F ** m_ResMeanVPosVSens
ResidaulMeanV vs V for sensor.
void DefineSensors()
Define histograms which depend on position for individual sensors.
TH1F ** m_ResMeanUPosVSens
ResidaulMeanU vs V for sensor.
TH2F ** m_ResVPosUSens
ResidaulV vs U for sensor.
virtual void FillHelixParametersAndCorrelations(const TrackFitResult *tfr) override
Fill histograms with helix parameters and their correlations.
TH2F * m_MomPtTanLambda
Track momentum Pt vs.
TH2F * m_PhiMomPt
Phi - the angle of the transverse momentum in the r-phi plane vs.
TH2F * m_PhiTanLambda
Phi - the angle of the transverse momentum in the r-phi plane vs.
TH2F ** m_ResMeanVPhiThetaLayer
ResidaulMeanU vs Phi vs Theta for Layer.
virtual void event() override
Module function event.
TH2F ** m_ResUPhiLayer
ResidaulU vs Phi for Layer.
virtual TH1F * Create(std::string name, std::string title, int nbinsx, double xlow, double xup, std::string xTitle, std::string yTitle) override
Function to create TH1F and add it to the vector of histograms (m_histograms).
virtual void endRun() override
Module function endRun.
TH2F ** m_ResMeanPhiThetaLayerCounts
Special Alignment related: Layer level.
TH1F ** m_ResMeanUPhiLayer
ResidaulMeanU vs Phi for Layer.
TH2F ** m_ResMeanVPosUVSens
ResidaulMeanU vs U vs V for sensor.
TH2F * m_MomPtOmega
Track momentum Pt vs.
TH2F ** m_ResMeanPosUVSensCounts
Special Alignment related: Sensor level.
TH2F * m_Z0TanLambda
z0 - the z0 coordinate of the perigee vs.
TH2F ** m_ResMeanUPhiThetaLayer
ResidaulMeanU vs Phi vs Theta for Layer.
TH1F ** m_ResMeanVPhiLayer
ResidaulMeanV vs Phi for Layer.
virtual void FillPositionSensors(ROOT::Math::XYZVector residual_um, ROOT::Math::XYZVector position, int sensorIndex)
Fill histograms which depend on position for individual sensors.
AlignDQMModule()
Constructor.
TH2F ** m_ResVThetaLayer
ResidaulV vs Theta for Layer.
TH2F * m_Z0Omega
z0 - the z0 coordinate of the perigee vs.
TH2F ** m_ResMeanUPosUVSens
ResidaulMeanU vs U vs V for sensor.
TH2F * m_OmegaTanLambda
Omega - the curvature of the track vs.
virtual void DefineHelixParametersAndCorrelations() override
All the following Define- functions should be used in the defineHisto() function to define histograms...
virtual void FillLayers(ROOT::Math::XYZVector residual_um, float phi_deg, float theta_deg, int layerIndex)
Fill histograms which depend on layerIndex.
TH1F ** m_ResMeanUThetaLayer
ResidaulMeanU vs Theta for Layer.
TH2F * m_PhiOmega
Phi - the angle of the transverse momentum in the r-phi plane vs.
TH2F * m_PhiZ0
helix parameters and their corellations
TH2F ** m_ResVPosVSens
ResidaulV vs V for sensor.
TH1F ** m_ResMeanUPosUSens
ResidaulMeanU vs U for sensor.
TH2F ** m_ResUPosVSens
ResidaulU vs V for sensor.
TH2F * m_Z0MomPt
z0 - the z0 coordinate of the perigee vs.
TH2F ** m_ResUPosUSens
ResidaulU vs U for sensor.
TH2F * m_D0MomPt
d0 - signed distance to the IP in r-phi vs.
TH1F ** m_ResMeanVThetaLayer
ResidaulMeanV vs Theta for Layer.
TH2F ** m_ResVPhiLayer
ResidaulV vs Phi for Layer.
TH2F * m_D0Omega
d0 - signed distance to the IP in r-phi vs.
TH1F ** m_ResMeanVPosUSens
ResidaulMeanV vs U for sensor.
TH2F * m_D0TanLambda
d0 - signed distance to the IP in r-phi vs.
virtual void DefineLayers()
Define histograms which depend on layerIndex.
TH2F ** m_ResUThetaLayer
ResidaulU vs Theta for Layer.
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
virtual void Run()
Call this to start processing the event data and filling histograms.
This class serves as a base for the TrackDQMModule and AlignDQMModule (and possibly other DQM histogr...
bool histogramsDefined
True if the defineHisto() was called.
TH1F * m_MomPt
Track momentum Pt.
TH1F * m_TanLambda
TanLambda - the slope of the track in the r-z plane.
TH1F ** m_UBResidualsSensorU
Unbiased residuals for PXD and SVD u per sensor.
std::vector< std::tuple< std::string, std::string, std::string > > m_histogramParameterChanges
Used for changing parameters of histograms via the ProcessHistogramParameterChange function.
virtual void DefineHalfShellsVXD()
Define histograms with unbiased residuals for half-shells for PXD and SVD sensors.
virtual void DefineUBResidualsVXD()
Define histograms with unbiased residuals in PXD and SVD sensors.
virtual void DefineTRClusters()
Define histograms with correlations between neighbor layers and cluster hitmap in IP angle range.
void ProcessHistogramParameterChange(const std::string &name, const std::string &parameter, const std::string &value)
Process one change in histogram parameters.
TH1F * m_Omega
Omega - the curvature of the track.
virtual void event() override
This method is called for each event.
virtual void DefineTrackFitStatus()
Define histograms which require FitStatus.
std::string m_tracksStoreArrayName
StoreArray name where Tracks are written.
TH2F * m_PhiD0
d0 vs Phi - the signed distance to the IP in the r-phi plane
virtual TH1F * Create(std::string name, std::string title, int nbinsx, double xlow, double xup, std::string xTitle, std::string yTitle)
Function to create TH1F and add it to the vector of histograms (m_histograms).
TH1F * m_Phi
Phi - the angle of the transverse momentum in the r-phi plane, with CDF naming convention.
virtual void DefineHits()
Define histograms with numbers of hits.
virtual void FillHelixParametersAndCorrelations(const TrackFitResult *tfr)
Fill histograms with helix parameters and their correlations.
TH2F ** m_UBResidualsSensor
Unbiased residuals for PXD and SVD u vs v per sensor.
virtual void DefineMomentumCoordinates()
Define histograms with track momentum Pt.
TH1F ** m_UBResidualsSensorV
Unbiased residuals for PXD and SVD v per sensor.
static void ComputeMean(TH1F *output, TH2F *input, bool onX=true)
Creates a graph of means by given axis from given TH2F histogram.
virtual void DefineTracks()
All the following Define- functions should be used in the defineHisto() function to define histograms...
TH2F * m_D0Z0
z0 vs d0 - signed distance to the IP in r-phi vs.
virtual void DefineMomentumAngles()
Define histograms with track momentum Pt.
std::string m_recoTracksStoreArrayName
StoreArray name where RecoTracks are written.
TH1F * m_D0
d0 - the signed distance to the IP in the r-phi plane
TH1F * m_Z0
z0 - the z0 coordinate of the perigee (beam spot position)
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
This class unites some parameters for Factory which describe one axis of histogram.
Axis & title(std::string title)
Set value of title.
This class is used for creating TH1F and TH2F objects.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Values of the result of a track fit with a given particle hypothesis.
double getOmega() const
Getter for omega.
double getD0() const
Getter for d0.
double getTanLambda() const
Getter for tanLambda.
double getZ0() const
Getter for z0.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
double getPhi0() const
Getter for phi0.
static const double mm
[millimeters]
Definition: Unit.h:70
static const double deg
degree to radians
Definition: Unit.h:109
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:176
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:142
int getLayerIndex(unsigned short layer) const
Return index of layer in plots.
Definition: GeoTools.h:426
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
STL namespace.