Belle II Software development
VXDDQMExpressRecoModule.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 "vxd/modules/vxdDQM/VXDDQMExpressRecoModule.h"
10
11#include <framework/datastore/StoreArray.h>
12#include <framework/datastore/RelationArray.h>
13
14#include <svd/dataobjects/SVDShaperDigit.h>
15#include <svd/dataobjects/SVDCluster.h>
16#include <svd/geometry/SensorInfo.h>
17
18#include <pxd/dataobjects/PXDDigit.h>
19#include <pxd/dataobjects/PXDCluster.h>
20#include <pxd/geometry/SensorInfo.h>
21
22#include <vxd/geometry/SensorInfoBase.h>
23#include <vxd/geometry/GeoTools.h>
24
25#include <TMath.h>
26#include <Math/Vector3D.h>
27
28#include <boost/format.hpp>
29
30#include "TDirectory.h"
31
32using namespace std;
33using boost::format;
34using namespace Belle2;
35
36//-----------------------------------------------------------------
37// Register the Module
38//-----------------------------------------------------------------
39REG_MODULE(VXDDQMExpressReco);
40
41
42//-----------------------------------------------------------------
43// Implementation
44//-----------------------------------------------------------------
45
47{
48 //Set module properties
49 setDescription("VXD DQM module for Express Reco "
50 "Recommended Number of events for monitor is 40 kEvents or more to fill all histograms "
51 );
52
53 setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
54 addParam("UseDigits", m_UseDigits,
55 "flag <0,1> for using digits only, no cluster information will be required, default = 0 ", m_UseDigits);
56 addParam("CorrelationGranulation", m_CorrelationGranulation,
57 "Set granulation of histogram plots, default is 1 degree, min = 0.02 degree, max = 1 degree ", m_CorrelationGranulation);
58 addParam("SwapPXD", m_SwapPXD, "flag <0,1> very special case for swap of phi-theta coordinates, default = 0 ", m_SwapPXD);
59 addParam("CutCorrelationSigPXD", m_CutCorrelationSigPXD,
60 "Cut threshold of PXD signal for accepting to correlations, default = 0 ADU ", m_CutCorrelationSigPXD);
61 addParam("CutCorrelationSigUSVD", m_CutCorrelationSigUSVD,
62 "Cut threshold of SVD signal for accepting to correlations in u, default = 0 ADU ", m_CutCorrelationSigUSVD);
63 addParam("CutCorrelationSigVSVD", m_CutCorrelationSigVSVD,
64 "Cut threshold of SVD signal for accepting to correlations in v, default = 0 ADU ", m_CutCorrelationSigVSVD);
65 addParam("CutCorrelationTimeSVD", m_CutCorrelationTimeSVD,
66 "Cut threshold of SVD time window for accepting to correlations, default = 70 ns ", m_CutCorrelationTimeSVD);
67
68 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
69 std::string("VXDExpReco"));
70}
71
72
73VXDDQMExpressRecoModule::~VXDDQMExpressRecoModule()
74{
75}
76
77//------------------------------------------------------------------
78// Function to define histograms
79//-----------------------------------------------------------------
80
82{
84 if (gTools->getNumberOfLayers() == 0) {
85 B2WARNING("Missing geometry for VXD, check steering file.");
86 return;
87 }
88 if (gTools->getNumberOfPXDLayers() == 0) {
89 B2WARNING("Missing geometry for PXD.");
90 }
91 if (gTools->getNumberOfSVDLayers() == 0) {
92 B2WARNING("Missing geometry for SVD.");
93 }
94
95 // Create a separate histogram directories and cd into it.
96 TDirectory* oldDir = gDirectory;
97 if (m_histogramDirectoryName != "") {
98 oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
99 oldDir->cd(m_histogramDirectoryName.c_str());
100 }
101
102 if (m_CorrelationGranulation > 5.0) m_CorrelationGranulation = 5.0; // set maximum of gramularity to 1 degree
103 if (m_CorrelationGranulation < 0.02) m_CorrelationGranulation = 0.02; // set minimum of gramularity to 0.02 degree
104
105 // if binning go over h_MaxBins it decrease preset of range
106 int h_MaxBins = 2000; //maximal size of histogram binning:
107
109
110 // basic constants presets:
111 int nVXDLayers = gTools->getNumberOfLayers();
112
113 // Create basic histograms:
114 m_correlationsSP = (TH2F**) new TH2F*[nVXDLayers * nVXDLayers];
115 m_correlationsSP1DPhi = (TH1F**) new TH1F*[nVXDLayers * nVXDLayers];
116 m_correlationsSP1DTheta = (TH1F**) new TH1F*[nVXDLayers * nVXDLayers];
117
118
119 for (int i = 0; i < nVXDLayers; i++) {
120 for (int j = 0; j < nVXDLayers; j++) {
121 m_correlationsSP[nVXDLayers * j + i] = NULL;
122 m_correlationsSP1DPhi[nVXDLayers * j + i] = NULL;
123 m_correlationsSP1DTheta[nVXDLayers * j + i] = NULL;
124 }
125 }
126
127 string Diru = str(format("Phi"));
128 string Dirv = str(format("Theta"));
129 string Unit = str(format("degree"));
130 string AxisTitPhi = str(format("%1% position [%2%]") % Diru.c_str() % Unit.c_str()); // 0..360 degree, from x+ to y+...
131 string AxisTitTheta = str(format("%1% position [%2%]") % Dirv.c_str() % Unit.c_str()); // 0..180 degree, z+ to z-
132 if (m_UseDigits) {
133 AxisTitPhi = str(format("From digits: %1%") % AxisTitPhi);
134 AxisTitTheta = str(format("From digits: %1%") % AxisTitTheta);
135 }
136 for (VxdID layer1 : geo.getLayers()) {
137 int i = layer1.getLayerNumber() - gTools->getFirstLayer();
138 float uSize1s, uSize1e;
139 float vSize1s, vSize1e;
140 int nStripsU1, nStripsV1;
141 uSize1s = -180.0; // degree
142 uSize1e = 180.0; // degree
143 vSize1s = 0.0; // degree
144 vSize1e = 180.0; // degree
145 nStripsU1 = (uSize1e - uSize1s) / m_CorrelationGranulation;
146 if (nStripsU1 > h_MaxBins) nStripsU1 = h_MaxBins;
147 nStripsV1 = nStripsU1;
148 for (VxdID layer2 : geo.getLayers()) {
149 int j = layer2.getLayerNumber() - gTools->getFirstLayer();
150 float uSize2s, uSize2e;
151 float vSize2s, vSize2e;
152 int nStripsU2, nStripsV2;
153 uSize2s = -180.0; // degree
154 uSize2e = 180.0; // degree
155 vSize2s = 0.0; // degree
156 vSize2e = 180.0; // degree
157 nStripsU2 = nStripsU1;
158 nStripsV2 = nStripsV1;
159
160 if (i == j) { // hit maps
161 string nameSP = str(format("DQMER_VXD_Hitmap_L%1%") % layer2.getLayerNumber());
162 string titleSP = str(format("DQM ER VXD Hitmap, layer %1%") % layer2.getLayerNumber());
163 m_correlationsSP[nVXDLayers * j + i] = new TH2F(nameSP.c_str(), titleSP.c_str(),
164 nStripsU2, uSize2s, uSize2e,
165 nStripsV2, vSize2s, vSize2e);
166 m_correlationsSP[nVXDLayers * j + i]->GetXaxis()->SetTitle(AxisTitPhi.c_str());
167 m_correlationsSP[nVXDLayers * j + i]->GetYaxis()->SetTitle(AxisTitTheta.c_str());
168 m_correlationsSP[nVXDLayers * j + i]->GetZaxis()->SetTitle("hits");
169
170 nameSP = str(format("DQMER_VXD_Hitmap_%1%_L%2%") % Diru.c_str() % layer2.getLayerNumber());
171 titleSP = str(format("DQM ER VXD Hitmap in %1%, layer %2%") % Diru.c_str() % layer2.getLayerNumber());
172 m_correlationsSP1DPhi[nVXDLayers * j + i] = new TH1F(nameSP.c_str(), titleSP.c_str(),
173 nStripsU2, uSize2s, uSize2e);
174 m_correlationsSP1DPhi[nVXDLayers * j + i]->GetXaxis()->SetTitle(AxisTitPhi.c_str());
175 m_correlationsSP1DPhi[nVXDLayers * j + i]->GetYaxis()->SetTitle("hits");
176
177 nameSP = str(format("DQMER_VXD_Hitmap_%1%_L%2%") % Dirv.c_str() % layer2.getLayerNumber());
178 titleSP = str(format("DQM ER VXD Hitmap in %1%, layer %2%") % Dirv.c_str() % layer2.getLayerNumber());
179 m_correlationsSP1DTheta[nVXDLayers * j + i] = new TH1F(nameSP.c_str(), titleSP.c_str(),
180 nStripsV2, vSize2s, vSize2e);
181 m_correlationsSP1DTheta[nVXDLayers * j + i]->GetXaxis()->SetTitle(AxisTitTheta.c_str());
182 m_correlationsSP1DTheta[nVXDLayers * j + i]->GetYaxis()->SetTitle("hits");
183
184 } else if (i < j) { // correlations for Phi
185 string nameSP = str(format("DQMER_VXD_Correlations_%1%_L%2%_L%3%") % Diru.c_str() % layer1.getLayerNumber() %
186 layer2.getLayerNumber());
187 string titleSP = str(format("DQM ER VXD Correlations in %1%, layers %2% %3%") % Diru.c_str() % layer1.getLayerNumber() %
188 layer2.getLayerNumber());
189 m_correlationsSP[nVXDLayers * j + i] = new TH2F(nameSP.c_str(), titleSP.c_str(),
190 nStripsU1, uSize1s, uSize1e,
191 nStripsU2, uSize2s, uSize2e);
192 string axisxtitle = str(format("%1%, layer %2%") % AxisTitPhi.c_str() % layer1.getLayerNumber());
193 string axisytitle = str(format("%1%, layer %2%") % AxisTitPhi.c_str() % layer2.getLayerNumber());
194 m_correlationsSP[nVXDLayers * j + i]->GetXaxis()->SetTitle(axisxtitle.c_str());
195 m_correlationsSP[nVXDLayers * j + i]->GetYaxis()->SetTitle(axisytitle.c_str());
196 m_correlationsSP[nVXDLayers * j + i]->GetZaxis()->SetTitle("hits");
197
198 nameSP = str(format("DQMER_VXD_1D_Correlations_%1%_L%2%_L%3%") % Diru.c_str() % layer1.getLayerNumber() % layer2.getLayerNumber());
199 titleSP = str(format("DQM ER VXD 1D Correlations in %1%, layers %2% %3%") % Diru.c_str() % layer1.getLayerNumber() %
200 layer2.getLayerNumber());
201 m_correlationsSP1DPhi[nVXDLayers * j + i] = new TH1F(nameSP.c_str(), titleSP.c_str(),
202 nStripsU1, uSize1s, uSize1e);
203 axisxtitle = str(format("%1%, layer %2% - %3%") % AxisTitPhi.c_str() % layer1.getLayerNumber() % layer2.getLayerNumber());
204 m_correlationsSP1DPhi[nVXDLayers * j + i]->GetXaxis()->SetTitle(axisxtitle.c_str());
205 m_correlationsSP1DPhi[nVXDLayers * j + i]->GetYaxis()->SetTitle("hits");
206 } else { // correlations for Theta
207 string nameSP = str(format("DQMER_VXD_Correlations_%1%_L%2%_L%3%") % Dirv.c_str() % layer1.getLayerNumber() %
208 layer2.getLayerNumber());
209 string titleSP = str(format("DQM ER VXD Correlations in %1%, layers %2% %3%") % Dirv.c_str() % layer1.getLayerNumber() %
210 layer2.getLayerNumber());
211 m_correlationsSP[nVXDLayers * j + i] = new TH2F(nameSP.c_str(), titleSP.c_str(),
212 nStripsV1, vSize1s, vSize1e,
213 nStripsV2, vSize2s, vSize2e);
214 string axisxtitle = str(format("%1%, layer %2%") % AxisTitTheta.c_str() % layer1.getLayerNumber());
215 string axisytitle = str(format("%1%, layer %2%") % AxisTitTheta.c_str() % layer2.getLayerNumber());
216 m_correlationsSP[nVXDLayers * j + i]->GetXaxis()->SetTitle(axisxtitle.c_str());
217 m_correlationsSP[nVXDLayers * j + i]->GetYaxis()->SetTitle(axisytitle.c_str());
218 m_correlationsSP[nVXDLayers * j + i]->GetZaxis()->SetTitle("hits");
219
220 nameSP = str(format("DQMER_VXD_1D_Correlations_%1%_L%2%_L%3%") % Dirv.c_str() % layer1.getLayerNumber() % layer2.getLayerNumber());
221 titleSP = str(format("DQM ER VXD 1D Correlations in %1%, layers %2% %3%") % Dirv.c_str() % layer1.getLayerNumber() %
222 layer2.getLayerNumber());
223 m_correlationsSP1DTheta[nVXDLayers * j + i] = new TH1F(nameSP.c_str(), titleSP.c_str(),
224 nStripsV1, -vSize1e, vSize1e);
225 axisxtitle = str(format("%1%, layer %2% - %3%") % AxisTitTheta.c_str() % layer1.getLayerNumber() % layer2.getLayerNumber());
226 m_correlationsSP1DTheta[nVXDLayers * j + i]->GetXaxis()->SetTitle(axisxtitle.c_str());
227 m_correlationsSP1DTheta[nVXDLayers * j + i]->GetYaxis()->SetTitle("hits");
228 }
229 }
230 }
231
232 oldDir->cd();
233}
234
235
237{
238 // Register histograms (calls back defineHisto)
239 REG_HISTOGRAM
240
241 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
242 if (gTools->getNumberOfLayers() != 0) {
243 //Register collections
248 RelationArray relPXDClusterDigits(storePXDClusters, storePXDDigits);
249 RelationArray relSVDClusterDigits(storeSVDClusters, storeSVDShaperDigits);
250 m_storePXDClustersName = storePXDClusters.getName();
251 m_relPXDClusterDigitName = relPXDClusterDigits.getName();
252 m_storeSVDClustersName = storeSVDClusters.getName();
253 m_relSVDClusterDigitName = relSVDClusterDigits.getName();
254
255 storePXDDigits.isRequired();
256 storeSVDShaperDigits.isRequired();
257
258 //Store names to speed up creation later
259 m_storePXDDigitsName = storePXDDigits.getName();
260 m_storeSVDShaperDigitsName = storeSVDShaperDigits.getName();
261 }
262}
263
265{
266 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
267 if (gTools->getNumberOfLayers() == 0) return;
268
269 // Just to make sure, reset all the histograms.
270 int nVXDLayers = gTools->getNumberOfLayers();
271 for (int i = 0; i < nVXDLayers; i++) {
272 for (int j = 0; j < nVXDLayers; j++) {
273 if (m_correlationsSP[nVXDLayers * j + i] != NULL) m_correlationsSP[nVXDLayers * j + i]->Reset();
274 if (m_correlationsSP1DPhi[nVXDLayers * j + i] != NULL) m_correlationsSP1DPhi[nVXDLayers * j + i]->Reset();
275 if (m_correlationsSP1DTheta[nVXDLayers * j + i] != NULL) m_correlationsSP1DTheta[nVXDLayers * j + i]->Reset();
276 }
277 }
278}
279
280
282{
283 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
284 if (gTools->getNumberOfLayers() == 0) return;
285
286 const StoreArray<PXDDigit> storePXDDigits(m_storePXDDigitsName);
287 const StoreArray<SVDShaperDigit> storeSVDShaperDigits(m_storeSVDShaperDigitsName);
288
289 const StoreArray<SVDCluster> storeSVDClusters(m_storeSVDClustersName);
290 const StoreArray<PXDCluster> storePXDClusters(m_storePXDClustersName);
291
292 const RelationArray relPXDClusterDigits(storePXDClusters, storePXDDigits, m_relPXDClusterDigitName);
293 const RelationArray relSVDClusterDigits(storeSVDClusters, storeSVDShaperDigits, m_relSVDClusterDigitName);
294
295 // If there are no digits, leave
296 if (!storePXDDigits && !storeSVDShaperDigits) return;
297
298 // basic constants presets:
299 int nVXDLayers = gTools->getNumberOfLayers();
300 int firstVXDLayer = gTools->getFirstLayer();
301 int firstPXDLayer = gTools->getFirstPXDLayer();
302 int lastPXDLayer = gTools->getLastPXDLayer();
303 int firstSVDLayer = gTools->getFirstSVDLayer();
304 int lastSVDLayer = gTools->getLastSVDLayer();
305
306 int MaxHits = 0;
307 if (m_UseDigits) {
308 MaxHits = storeSVDShaperDigits.getEntries() + storePXDDigits.getEntries();
309 } else {
310 MaxHits = storeSVDClusters.getEntries() + storePXDClusters.getEntries();
311 }
312
313 // If there are no hits, leave
314 if (MaxHits == 0) return;
315
316 for (int i1 = 0; i1 < MaxHits; i1++) {
317 // preparing of first value for correlation plots with postfix "1":
318 float fTime1 = 0.0;
319 float fPosSPU1 = 0.0;
320 float fPosSPV1 = 0.0;
321 int iIsPXD1 = 0;
322 int index1 = 0;
323 int iIsU1 = 0;
324 int iIsV1 = 0;
325 int iLayer1 = 0;
326 if (m_UseDigits) { // DIGITS:
327 if (i1 < storePXDDigits.getEntries()) { // PXD digits/clusters:
328 const PXDDigit& digitPXD1 = *storePXDDigits[i1];
329 iLayer1 = digitPXD1.getSensorID().getLayerNumber();
330 if ((iLayer1 < firstPXDLayer) || (iLayer1 > lastPXDLayer)) continue;
331 index1 = iLayer1 - firstVXDLayer;
332 float fCharge1 = digitPXD1.getCharge();
333 if (fCharge1 < m_CutCorrelationSigPXD) continue;
334 VxdID sensorID1 = digitPXD1.getSensorID();
335 auto info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID1));
336 ROOT::Math::XYZVector rLocal1(info.getUCellPosition(digitPXD1.getUCellID()), info.getVCellPosition(digitPXD1.getVCellID()), 0);
337 ROOT::Math::XYZVector ral1 = info.pointToGlobal(rLocal1);
338 iIsPXD1 = 1;
339 iIsU1 = 1;
340 iIsV1 = 1;
341 fPosSPU1 = ral1.Phi() / TMath::Pi() * 180;
342 fPosSPV1 = ral1.Theta() / TMath::Pi() * 180;
343 if (m_SwapPXD) {
344 fPosSPV1 = ral1.Phi() / TMath::Pi() * 180;
345 fPosSPU1 = ral1.Theta() / TMath::Pi() * 180;
346 }
347 } else { // SVD digits/clusters:
348 const SVDShaperDigit& digitSVD1 = *storeSVDShaperDigits[i1 - storePXDDigits.getEntries()];
349 iLayer1 = digitSVD1.getSensorID().getLayerNumber();
350 if ((iLayer1 < firstSVDLayer) || (iLayer1 > lastSVDLayer)) continue;
351 index1 = iLayer1 - firstVXDLayer;
352 fTime1 = digitSVD1.getFADCTime();
353 VxdID sensorID1 = digitSVD1.getSensorID();
354 auto info = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID1));
355 SVDShaperDigit::APVFloatSamples samples = digitSVD1.getSamples();
356 if (digitSVD1.isUStrip()) {
357 int iCont = 0;
358 for (size_t i = 0; i < SVDShaperDigit::c_nAPVSamples; ++i) {
359 float fCharge1 = samples[i];
360 if (fCharge1 > m_CutCorrelationSigUSVD) iCont = 1;
361 }
362 if (iCont == 0) continue;
363 float possi = info.getUCellPosition(digitSVD1.getCellID());
364 ROOT::Math::XYZVector rLocal1(possi, 0, 0);
365 ROOT::Math::XYZVector ral1 = info.pointToGlobal(rLocal1);
366 iIsU1 = 1;
367 fPosSPU1 = ral1.Phi() / TMath::Pi() * 180;
368 } else {
369 int iCont = 0;
370 for (size_t i = 0; i < SVDShaperDigit::c_nAPVSamples; ++i) {
371 float fCharge1 = samples[i];
372 if (fCharge1 > m_CutCorrelationSigVSVD) iCont = 1;
373 }
374 if (iCont == 0) continue;
375
376 //float possi = digitSVD1.getCellPosition(); // is not work anymore
377 float possi = info.getVCellPosition(digitSVD1.getCellID());
378
379 ROOT::Math::XYZVector rLocal1(0, possi, 0);
380 ROOT::Math::XYZVector ral1 = info.pointToGlobal(rLocal1);
381 iIsV1 = 1;
382 fPosSPV1 = ral1.Theta() / TMath::Pi() * 180;
383 }
384 }
385 } else { // CLUSTERS:
386 if (i1 < storePXDClusters.getEntries()) { // PXD digits/clusters:
387 const PXDCluster& clusterPXD1 = *storePXDClusters[i1];
388 iLayer1 = clusterPXD1.getSensorID().getLayerNumber();
389 if ((iLayer1 < firstPXDLayer) || (iLayer1 > lastPXDLayer)) continue;
390 index1 = iLayer1 - firstVXDLayer;
391 float fCharge1 = clusterPXD1.getCharge();
392 if (fCharge1 < m_CutCorrelationSigPXD) continue;
393 VxdID sensorID1 = clusterPXD1.getSensorID();
394 auto info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID1));
395 ROOT::Math::XYZVector rLocal1(clusterPXD1.getU(), clusterPXD1.getV(), 0);
396 ROOT::Math::XYZVector ral1 = info.pointToGlobal(rLocal1);
397 iIsPXD1 = 1;
398 iIsU1 = 1;
399 iIsV1 = 1;
400 fPosSPU1 = ral1.Phi() / TMath::Pi() * 180;
401 fPosSPV1 = ral1.Theta() / TMath::Pi() * 180;
402 if (m_SwapPXD) {
403 fPosSPV1 = ral1.Phi() / TMath::Pi() * 180;
404 fPosSPU1 = ral1.Theta() / TMath::Pi() * 180;
405 }
406 } else { // SVD digits/clusters:
407 const SVDCluster& clusterSVD1 = *storeSVDClusters[i1 - storePXDClusters.getEntries()];
408 iLayer1 = clusterSVD1.getSensorID().getLayerNumber();
409 if ((iLayer1 < firstSVDLayer) || (iLayer1 > lastSVDLayer)) continue;
410 index1 = iLayer1 - firstVXDLayer;
411 float fCharge1 = clusterSVD1.getCharge();
412 fTime1 = clusterSVD1.getClsTime();
413 VxdID sensorID1 = clusterSVD1.getSensorID();
414 auto info = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID1));
415 if (clusterSVD1.isUCluster()) {
416 if (fCharge1 < m_CutCorrelationSigUSVD * 200) continue; // in electrons
417 ROOT::Math::XYZVector rLocal1(clusterSVD1.getPosition(), 0, 0);
418 ROOT::Math::XYZVector ral1 = info.pointToGlobal(rLocal1);
419 iIsU1 = 1;
420 fPosSPU1 = ral1.Phi() / TMath::Pi() * 180;
421 } else {
422 if (fCharge1 < m_CutCorrelationSigVSVD * 200) continue; // in electrons
423 ROOT::Math::XYZVector rLocal1(0, clusterSVD1.getPosition(), 0);
424 ROOT::Math::XYZVector ral1 = info.pointToGlobal(rLocal1);
425 iIsV1 = 1;
426 fPosSPV1 = ral1.Theta() / TMath::Pi() * 180;
427 }
428 }
429 }
430 // hit maps for PXD:
431 if ((iIsU1 == 1) && (iIsV1 == 1)) {
432 if (m_correlationsSP[nVXDLayers * index1 + index1] != NULL)
433 m_correlationsSP[nVXDLayers * index1 + index1]->Fill(fPosSPU1, fPosSPV1);
434 if (m_correlationsSP1DPhi[nVXDLayers * index1 + index1] != NULL)
435 m_correlationsSP1DPhi[nVXDLayers * index1 + index1]->Fill(fPosSPU1);
436 if (m_correlationsSP1DTheta[nVXDLayers * index1 + index1] != NULL)
437 m_correlationsSP1DTheta[nVXDLayers * index1 + index1]->Fill(fPosSPV1);
438 }
439 for (int i2 = 0; i2 < MaxHits; i2++) {
440 // preparing of second value for correlation plots with postfix "2":
441 float fTime2 = 0.0;
442 float fPosSPU2 = 0.0;
443 float fPosSPV2 = 0.0;
444 int iIsPXD2 = 0;
445 int index2 = 0;
446 int iIsU2 = 0;
447 int iIsV2 = 0;
448 int iLayer2 = 0;
449 if (m_UseDigits) { // DIGITS:
450 if (i2 < storePXDDigits.getEntries()) { // PXD digits/clusters:
451 const PXDDigit& digitPXD2 = *storePXDDigits[i2];
452 iLayer2 = digitPXD2.getSensorID().getLayerNumber();
453 if ((iLayer2 < firstPXDLayer) || (iLayer2 > lastPXDLayer)) continue;
454 index2 = iLayer2 - firstVXDLayer;
455 float fCharge2 = digitPXD2.getCharge();
456 if (fCharge2 < m_CutCorrelationSigPXD) continue;
457 VxdID sensorID2 = digitPXD2.getSensorID();
458 auto info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID2));
459 ROOT::Math::XYZVector rLocal2(info.getUCellPosition(digitPXD2.getUCellID()), info.getVCellPosition(digitPXD2.getVCellID()), 0);
460 ROOT::Math::XYZVector ral2 = info.pointToGlobal(rLocal2);
461 iIsPXD2 = 1;
462 iIsU2 = 1;
463 iIsV2 = 1;
464 fPosSPU2 = ral2.Phi() / TMath::Pi() * 180;
465 fPosSPV2 = ral2.Theta() / TMath::Pi() * 180;
466 if (m_SwapPXD) {
467 fPosSPV2 = ral2.Phi() / TMath::Pi() * 180;
468 fPosSPU2 = ral2.Theta() / TMath::Pi() * 180;
469 }
470 } else { // SVD digits/clusters:
471 const SVDShaperDigit& digitSVD2 = *storeSVDShaperDigits[i2 - storePXDDigits.getEntries()];
472 iLayer2 = digitSVD2.getSensorID().getLayerNumber();
473 if ((iLayer2 < firstSVDLayer) || (iLayer2 > lastSVDLayer)) continue;
474 index2 = iLayer2 - firstVXDLayer;
475 fTime2 = digitSVD2.getFADCTime();
476 VxdID sensorID2 = digitSVD2.getSensorID();
477 auto info = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID2));
478
479 SVDShaperDigit::APVFloatSamples samples = digitSVD2.getSamples();
480 if (digitSVD2.isUStrip()) {
481 int iCont = 0;
482 for (size_t i = 0; i < SVDShaperDigit::c_nAPVSamples; ++i) {
483 float fCharge2 = samples[i];
484 if (fCharge2 > m_CutCorrelationSigUSVD) iCont = 1;
485 }
486 if (iCont == 0) continue;
487 float possi = info.getUCellPosition(digitSVD2.getCellID());
488 ROOT::Math::XYZVector rLocal2(possi, 0, 0);
489 ROOT::Math::XYZVector ral2 = info.pointToGlobal(rLocal2);
490 iIsU2 = 1;
491 fPosSPU2 = ral2.Phi() / TMath::Pi() * 180;
492 } else {
493 int iCont = 0;
494 for (size_t i = 0; i < SVDShaperDigit::c_nAPVSamples; ++i) {
495 float fCharge2 = samples[i];
496 if (fCharge2 > m_CutCorrelationSigVSVD) iCont = 1;
497 }
498 if (iCont == 0) continue;
499 float possi = info.getVCellPosition(digitSVD2.getCellID());
500 ROOT::Math::XYZVector rLocal2(0, possi, 0);
501 ROOT::Math::XYZVector ral2 = info.pointToGlobal(rLocal2);
502 iIsV2 = 1;
503 fPosSPV2 = ral2.Theta() / TMath::Pi() * 180;
504 }
505 }
506 } else { // CLUSTERS:
507 if (i2 < storePXDClusters.getEntries()) { // PXD digits/clusters:
508 const PXDCluster& clusterPXD2 = *storePXDClusters[i2];
509 iLayer2 = clusterPXD2.getSensorID().getLayerNumber();
510 if ((iLayer2 < firstPXDLayer) || (iLayer2 > lastPXDLayer)) continue;
511 index2 = iLayer2 - firstVXDLayer;
512 float fCharge2 = clusterPXD2.getCharge();
513 if (fCharge2 < m_CutCorrelationSigPXD) continue;
514 VxdID sensorID2 = clusterPXD2.getSensorID();
515 auto info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID2));
516 ROOT::Math::XYZVector rLocal2(clusterPXD2.getU(), clusterPXD2.getV(), 0);
517 ROOT::Math::XYZVector ral2 = info.pointToGlobal(rLocal2);
518 iIsPXD2 = 1;
519 iIsU2 = 1;
520 iIsV2 = 1;
521 fPosSPU2 = ral2.Phi() / TMath::Pi() * 180;
522 fPosSPV2 = ral2.Theta() / TMath::Pi() * 180;
523 if (m_SwapPXD) {
524 fPosSPV2 = ral2.Phi() / TMath::Pi() * 180;
525 fPosSPU2 = ral2.Theta() / TMath::Pi() * 180;
526 }
527 } else { // SVD digits/clusters:
528 const SVDCluster& clusterSVD2 = *storeSVDClusters[i2 - storePXDClusters.getEntries()];
529 iLayer2 = clusterSVD2.getSensorID().getLayerNumber();
530 if ((iLayer2 < firstSVDLayer) || (iLayer2 > lastSVDLayer)) continue;
531 index2 = iLayer2 - firstVXDLayer;
532 float fCharge2 = clusterSVD2.getCharge();
533 fTime2 = clusterSVD2.getClsTime();
534 VxdID sensorID2 = clusterSVD2.getSensorID();
535 auto info = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID2));
536 if (clusterSVD2.isUCluster()) {
537 if (fCharge2 < m_CutCorrelationSigUSVD * 200) continue; // in electrons
538 ROOT::Math::XYZVector rLocal2(clusterSVD2.getPosition(), 0, 0);
539 ROOT::Math::XYZVector ral2 = info.pointToGlobal(rLocal2);
540 iIsU2 = 1;
541 fPosSPU2 = ral2.Phi() / TMath::Pi() * 180;
542 } else {
543 if (fCharge2 < m_CutCorrelationSigVSVD * 200) continue; // in electrons
544 ROOT::Math::XYZVector rLocal2(0, clusterSVD2.getPosition(), 0);
545 ROOT::Math::XYZVector ral2 = info.pointToGlobal(rLocal2);
546 iIsV2 = 1;
547 fPosSPV2 = ral2.Theta() / TMath::Pi() * 180;
548 }
549 }
550 }
551 if ((iIsPXD1 == 0) && (iIsPXD2 == 0))
552 if ((fabs(fTime1 - fTime2)) > m_CutCorrelationTimeSVD) continue;
553 // ready to fill correlation histograms and hit maps:
554 if ((index1 == index2) && (iIsU1 == 1) && (iIsV2 == 1) && (iIsPXD1 == 0) && (iIsPXD2 == 0)) {
555 // hit maps for SVD:
556 if (m_correlationsSP[nVXDLayers * index2 + index1] != NULL)
557 m_correlationsSP[nVXDLayers * index2 + index1]->Fill(fPosSPU1, fPosSPV2);
558 if (m_correlationsSP1DPhi[nVXDLayers * index2 + index1] != NULL)
559 m_correlationsSP1DPhi[nVXDLayers * index2 + index1]->Fill(fPosSPU1);
560 if (m_correlationsSP1DTheta[nVXDLayers * index2 + index1] != NULL)
561 m_correlationsSP1DTheta[nVXDLayers * index2 + index1]->Fill(fPosSPV2);
562 } else if ((index1 < index2) && (iIsU1 == iIsU2) && (iIsU1 == 1)) {
563 // correlations for u
564 if (m_correlationsSP[nVXDLayers * index2 + index1] != NULL)
565 m_correlationsSP[nVXDLayers * index2 + index1]->Fill(fPosSPU1, fPosSPU2);
566 if (m_correlationsSP1DPhi[nVXDLayers * index2 + index1] != NULL)
567 m_correlationsSP1DPhi[nVXDLayers * index2 + index1]->Fill(fPosSPU2 - fPosSPU1);
568 } else if ((index1 > index2) && (iIsV1 == iIsV2) && (iIsV1 == 1)) {
569 // correlations for v
570 if (m_correlationsSP[nVXDLayers * index2 + index1] != NULL)
571 m_correlationsSP[nVXDLayers * index2 + index1]->Fill(fPosSPV2, fPosSPV1);
572 if (m_correlationsSP1DTheta[nVXDLayers * index2 + index1] != NULL)
573 m_correlationsSP1DTheta[nVXDLayers * index2 + index1]->Fill(fPosSPV2 - fPosSPV1);
574 }
575 }
576 }
577
578}
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
float getV() const
Get v coordinate of hit position.
Definition: PXDCluster.h:136
unsigned short getCharge() const
Get collected charge.
Definition: PXDCluster.h:156
VxdID getSensorID() const
Get the sensor ID.
Definition: PXDCluster.h:126
float getU() const
Get u coordinate of hit position.
Definition: PXDCluster.h:131
The PXD digit class.
Definition: PXDDigit.h:27
unsigned short getVCellID() const
Get cell ID in v.
Definition: PXDDigit.h:74
unsigned short getUCellID() const
Get cell ID in u.
Definition: PXDDigit.h:69
unsigned short getCharge() const
Get collected charge.
Definition: PXDDigit.h:79
VxdID getSensorID() const
Get the sensor ID.
Definition: PXDDigit.h:64
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
float getClsTime() const
Get average of waveform maximum times of cluster strip signals.
Definition: SVDCluster.h:134
float getCharge() const
Get collected charge.
Definition: SVDCluster.h:144
VxdID getSensorID() const
Get the sensor ID.
Definition: SVDCluster.h:102
bool isUCluster() const
Get the direction of strips.
Definition: SVDCluster.h:110
float getPosition(double v=0) const
Get the coordinate of reconstructed hit.
Definition: SVDCluster.h:117
The SVD ShaperDigit class.
static const std::size_t c_nAPVSamples
Number of APV samples stored.
std::array< APVFloatSampleType, c_nAPVSamples > APVFloatSamples
array of APVFloatSampleType objects
VxdID getSensorID() const
Get the sensor ID.
APVFloatSamples getSamples() const
Get array of samples.
short int getCellID() const
Get strip ID.
bool isUStrip() const
Get strip direction.
float getFADCTime() const
Get digit FADCTime estimate.
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
Definition: SensorInfo.h:25
const std::string & getName() const
Return name under which the object is saved in the DataStore.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
The Unit class.
Definition: Unit.h:40
TH2F ** m_correlationsSP
Correlations and hit maps from space points.
void initialize() override final
Initialize.
std::string m_relPXDClusterDigitName
PXDClustersToPXDDigits RelationArray name.
int m_SwapPXD
flag <0,1> very special case for swap of u-v coordinates
float m_CorrelationGranulation
set granulation of histogram plots, default is 1 deg (1 mm), min = 0.02, max = 5.0
std::string m_storePXDClustersName
PXDClusters StoreArray name.
void defineHisto() override final
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
std::string m_relSVDClusterDigitName
SVDClustersToSVDDigits RelationArray name.
float m_CutCorrelationTimeSVD
Cut threshold of SVD time window for accepting to correlations, default = 70 ns.
float m_CutCorrelationSigVSVD
Cut threshold of SVD signal for accepting to correlations in v, default = 0 ADU.
std::string m_storeSVDShaperDigitsName
SVDShaperDigits StoreArray name.
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
TH1F ** m_correlationsSP1DPhi
Correlations and hit maps from space points - differencies in Phi.
float m_CutCorrelationSigPXD
Cut threshold of PXD signal for accepting to correlations, default = 0 ADU.
int m_UseDigits
flag <0,1> for using digits only, no clusters will be required, default = 0
float m_CutCorrelationSigUSVD
Cut threshold of SVD signal for accepting to correlations in u, default = 0 ADU.
void beginRun() override final
Begin run.
TH1F ** m_correlationsSP1DTheta
Correlations and hit maps from space points - differencies in Theta.
std::string m_storeSVDClustersName
SVDClusters StoreArray name.
std::string m_storePXDDigitsName
PXDDigits StoreArray name.
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
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
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
unsigned short getNumberOfLayers() const
Get number of VXD layers.
Definition: GeoTools.h:41
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.