Belle II Software development
PXDValidationAlgorithm.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 <pxd/calibration/PXDValidationAlgorithm.h>
10
11#include <string>
12#include <algorithm>
13#include <set>
14
15#include <sstream>
16#include <iostream>
17
18#include <boost/format.hpp>
19
20//ROOT
21#include <TRandom.h>
22#include <TH1.h>
23#include <TF1.h>
24
25using namespace std;
26using boost::format;
27using namespace Belle2;
28
29
30// Anonymous namespace for data objects used by PXDValidationAlgorithm class
31namespace {
32
34 void getNumberOfBins(const std::shared_ptr<TH1I>& histo_ptr, unsigned short& nBinsU, unsigned short& nBinsV)
35 {
36 set<unsigned short> uBinSet;
37 set<unsigned short> vBinSet;
38
39 // Loop over all bins of input histo
40 for (auto histoBin = 1; histoBin <= histo_ptr->GetXaxis()->GetNbins(); histoBin++) {
41 // The bin label contains the vxdid, uBin and vBin
42 string label = histo_ptr->GetXaxis()->GetBinLabel(histoBin);
43
44 // Parse label string format to read sensorID, uBin and vBin
45 istringstream stream(label);
46 string token;
47 getline(stream, token, '_');
48 getline(stream, token, '_');
49 unsigned short uBin = std::stoi(token);
50
51 getline(stream, token, '_');
52 unsigned short vBin = std::stoi(token);
53
54 uBinSet.insert(uBin);
55 vBinSet.insert(vBin);
56 }
57
58 if (uBinSet.empty() || vBinSet.empty()) {
59 B2FATAL("Not able to determine the grid size. Something is wrong with collected data.");
60 } else {
61 nBinsU = *uBinSet.rbegin() + 1;
62 nBinsV = *vBinSet.rbegin() + 1;
63 }
64 }
65
67 unsigned short getNumberOfSensors(const std::shared_ptr<TH1I>& histo_ptr)
68 {
69 set<unsigned short> sensorSet;
70
71 // Loop over all bins of input histo
72 for (auto histoBin = 1; histoBin <= histo_ptr->GetXaxis()->GetNbins(); histoBin++) {
73 // The bin label contains the vxdid, uBin and vBin
74 string label = histo_ptr->GetXaxis()->GetBinLabel(histoBin);
75
76 // Parse label string format to read sensorID, uBin and vBin
77 istringstream stream(label);
78 string token;
79 getline(stream, token, '_');
80 VxdID sensorID(token);
81 sensorSet.insert(sensorID.getID());
82 }
83 return sensorSet.size();
84 }
85
86}
87
88
90 CalibrationAlgorithm("PXDCDSTCollector")
91 , minTrackPoints(1000), save2DHists(false), saveD0Z0(false), binsD0Z0(600)
92 , m_exp(-1), m_run(-1), m_hD0(nullptr), m_hZ0(nullptr)
93 , m_hTrackPointsLayer1(nullptr), m_hTrackClustersLayer1(nullptr)
94 , m_hTrackPointsLayer2(nullptr), m_hTrackClustersLayer2(nullptr)
95{
97 " -------------------------- PXDValidationAlgorithm ---------------------------------\n"
98 " \n"
99 " Algorithm for filling validation histograms. \n"
100 " ----------------------------------------------------------------------------------------\n"
101 );
102
103
104}
105
107{
108}
109
111{
112 // Turn off ROOT GUI
113 gROOT->SetBatch();
114
115 // Current Exp No. and Run No.
116 auto expRuns = getRunList();
117 auto expRunsAll = getRunListFromAllData();
118 // size should be 1
119 if (!expRuns.size())
120 return c_Failure;
121
122 m_exp = expRuns.back().first;
123 m_run = expRuns.back().second;
124 B2INFO("Current ExpRuns: [("
125 << expRuns.front().first << ", " << expRuns.front().second << "), ("
126 << m_exp << ", " << m_run << ")]");
127
128 // Get counter histograms and set pointers
129 auto cluster_counter = getObjectPtr<TH1I>("PXDTrackClusterCounter");
130 auto point_counter = getObjectPtr<TH1I>("PXDTrackPointCounter");
131 auto selected_point_counter = getObjectPtr<TH1I>("PXDSelTrackPointCounter"); // can be empty
132 auto selected_cluster_counter = getObjectPtr<TH1I>("PXDSelTrackClusterCounter"); // can be empty
133 if (!cluster_counter) return c_NotEnoughData;
134 if (!point_counter) return c_NotEnoughData;
135
136 // Extract number of sensors from counter histograms
137 auto nSensors = getNumberOfSensors(cluster_counter);
138
139 // Extract the number of grid bins from counter histograms
140 unsigned short nBinsU = 0;
141 unsigned short nBinsV = 0;
142 getNumberOfBins(cluster_counter, nBinsU, nBinsV);
143
144
145 B2INFO("Start info collection using a " << nBinsU << "x" << nBinsV << " grid per sensor.");
146 B2INFO("Number of collected track points is " << point_counter->GetEntries()
147 << " in " << nSensors << " sensors.");
148
149 if (point_counter->GetEntries() < minTrackPoints) {
150 B2WARNING("Not enough Data: Only " << point_counter->GetEntries() << " hits were collected but " << minTrackPoints
151 << " needed!");
152 return c_NotEnoughData;
153 }
154
155
156 auto hPassedHitsLayer1 = getObjectPtr<TH2F>("hPassedHitsLayer1");
157 auto hTotalHitsLayer1 = getObjectPtr<TH2F>("hTotalHitsLayer1");
158 auto hPassedHitsLayer2 = getObjectPtr<TH2F>("hPassedHitsLayer2");
159 auto hTotalHitsLayer2 = getObjectPtr<TH2F>("hTotalHitsLayer2");
160 if (!hPassedHitsLayer1) return c_NotEnoughData;
161 if (!hTotalHitsLayer1) return c_NotEnoughData;
162 if (!hPassedHitsLayer2) return c_NotEnoughData;
163 if (!hTotalHitsLayer2) return c_NotEnoughData;
164
165 // Save the current directory to change back later
166 TDirectory* currentDir = gDirectory;
167 //if (!currentDir) currentDir = gDirectory;
168
169 // Create TFile if not exist
170 if (!m_file) {
171 std::string fileName = (this->getPrefix()) + "Validation.root";
172 B2INFO("Creating file " << fileName);
173 m_file = std::make_shared<TFile>(fileName.c_str(), "RECREATE");
174 }
175
176 // Create TTree if not exist
177 if (!m_tree) {
178
179 B2INFO("Creating TTree.");
180 m_tree = std::make_shared<TTree>("tree", "PXD validation data");
181 // Define histograms out of m_file
182 currentDir->cd();
183 m_hD0 = new TH1F("hD0", "Corrected d0;#Delta d0/#sqrt{2} [cm];Counts", binsD0Z0, -0.03, 0.03);
184 m_hZ0 = new TH1F("hZ0", "Corrected z0;#Delta z0/#sqrt{2} [cm];Counts", binsD0Z0, -0.03, 0.03);
185 m_file->cd();
186
187 m_tree->Branch<int>("exp", &m_exp);
188 m_tree->Branch<int>("run", &m_run);
189 m_tree->Branch("pxdid", &m_pxdid);
190 m_tree->Branch("uBin", &m_uBin);
191 m_tree->Branch("vBin", &m_vBin);
192 m_tree->Branch("nTrackPoints", &m_nTrackPoints);
193 m_tree->Branch("nTrackClusters", &m_nTrackClusters);
194 m_tree->Branch("nSelTrackPoints", &m_nSelTrackPoints);
195 m_tree->Branch("nSelTrackClusters", &m_nSelTrackClusters);
196 m_tree->Branch<TH1F>("hD0", &m_hD0, 32000, 0);
197 m_tree->Branch<TH1F>("hZ0", &m_hZ0, 32000, 0);
198
199 if (saveD0Z0) {
200 m_tree->Branch("d0", &m_d0);
201 m_tree->Branch("z0", &m_z0);
202 }
203
204 if (save2DHists) {
205 m_hTrackPointsLayer1 = (TH2F*)hTotalHitsLayer1->Clone();
206 m_hTrackClustersLayer1 = (TH2F*)hPassedHitsLayer1->Clone();
207 m_hTrackPointsLayer2 = (TH2F*)hTotalHitsLayer2->Clone();
208 m_hTrackClustersLayer2 = (TH2F*)hPassedHitsLayer2->Clone();
209
210 // associated histograms
211 m_tree->Branch<TH2F>("hTrackPointsLayer1", &m_hTrackPointsLayer1, 32000, 0);
212 m_tree->Branch<TH2F>("hTrackClustersLayer1", &m_hTrackClustersLayer1, 32000, 0);
213 m_tree->Branch<TH2F>("hTrackPointsLayer2", &m_hTrackPointsLayer2, 32000, 0);
214 m_tree->Branch<TH2F>("hTrackClustersLayer2", &m_hTrackClustersLayer2, 32000, 0);
215 }
216 }
217
218 m_file->cd();
219
220 if (save2DHists) {
221 // Just update 2D histograms
222 *m_hTrackPointsLayer1 = *hTotalHitsLayer1;
223 *m_hTrackClustersLayer1 = *hPassedHitsLayer1;
224 *m_hTrackPointsLayer2 = *hTotalHitsLayer2;
225 *m_hTrackClustersLayer2 = *hPassedHitsLayer2;
226
227 }
228
229 // Clear
230 m_pxdid.clear();
231 m_uBin.clear();
232 m_vBin.clear();
233 m_nTrackClusters.clear();
234 m_nTrackPoints.clear();
235 m_nSelTrackPoints.clear();
236 m_nSelTrackClusters.clear();
237 m_d0.clear();
238 m_z0.clear();
239
240 // Get resolution trees and create histograms
241 auto tree_d0z0 = getObjectPtr<TTree>("tree_d0z0");
242 if (!tree_d0z0) return c_NotEnoughData;
243 float d0, z0;
244 tree_d0z0->SetBranchAddress("d0", &d0);
245 tree_d0z0->SetBranchAddress("z0", &z0);
246 m_hD0->Reset();
247 m_hZ0->Reset();
248
249 // Fill histograms from tree
250 for (int i = 0; i < tree_d0z0->GetEntries(); i++) {
251 tree_d0z0->GetEntry(i);
252 if (saveD0Z0) {
253 m_d0.emplace_back(d0);
254 m_z0.emplace_back(z0);
255 }
256 if (fabs(d0) > 0.03 || fabs(z0) > 0.03)
257 continue;
258 m_hD0->Fill(d0);
259 m_hZ0->Fill(z0);
260 }
261
262
263 // Loop over all bins of input histo
264 for (auto histoBin = 1; histoBin <= cluster_counter->GetXaxis()->GetNbins(); histoBin++) {
265 // The bin label contains the vxdid, uBin and vBin
266 string label = cluster_counter->GetXaxis()->GetBinLabel(histoBin);
267
268 // Parse label string format to read sensorID, uBin and vBin
269 istringstream stream(label);
270 string token;
271 getline(stream, token, '_');
272 VxdID sensorID(token);
273
274 getline(stream, token, '_');
275 unsigned short uBin = std::stoi(token);
276
277 getline(stream, token, '_');
278 unsigned short vBin = std::stoi(token);
279
280 // Read back the counters for number of collected clusters
281 int numberOfClusters = cluster_counter->GetBinContent(histoBin);
282 int numberOfPoints = point_counter->GetBinContent(histoBin);
283
284 m_pxdid.emplace_back(sensorID.getLayerNumber() * 1000 + sensorID.getLadderNumber() * 10 + sensorID.getSensorNumber());
285 m_uBin.emplace_back(uBin);
286 m_vBin.emplace_back(vBin);
287 m_nTrackPoints.emplace_back(numberOfPoints);
288 m_nSelTrackPoints.push_back(selected_point_counter->GetBinContent(histoBin));
289 m_nSelTrackClusters.push_back(selected_cluster_counter->GetBinContent(histoBin));
290 m_nTrackClusters.emplace_back(numberOfClusters);
291 }
292
293 m_tree->Fill();
294
295 // At the last run of the data
296 if (m_exp == expRunsAll.back().first &&
297 m_run == expRunsAll.back().second) {
298 B2INFO("Reached Final ExpRun: (" << m_exp << ", " << m_run << ")");
299 m_file->cd();
300 m_tree->Write();
301 B2INFO("Writing Successful.");
302 //if (m_file) m_file->Close();// Don't do it! We are using shared_ptr!
303 if (m_hD0) delete m_hD0;
304 if (m_hZ0) delete m_hZ0;
305 if (save2DHists) {
310 }
311 m_file->Close();
312
313 }
314 currentDir->cd();
315
316 // Always return ok to have run-by-run info
317 return c_OK;
318}
319
320
321bool PXDValidationAlgorithm::isBoundaryRequired(const Calibration::ExpRun& /*currentRun*/)
322{
323 // First run in data as we iterate, but our boundaries weren't set manually already?
324 // Just set the first run to be a boundary and we are done.
325 if (m_boundaries.empty()) {
326 B2INFO("This is the first run encountered, let's say it is a boundary.");
327 return true;
328 } else {
329 //return true;// -> simple run-by-run
330 return false;
331 }
332}
Base class for calibration algorithms.
std::vector< Calibration::ExpRun > m_boundaries
When using the boundaries functionality from isBoundaryRequired, this is used to store the boundaries...
std::vector< Calibration::ExpRun > getRunListFromAllData() const
Get the complete list of runs from inspection of collected data.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
std::string getPrefix() const
Get the prefix used for getting calibration data.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
@ c_Failure
Failed =3 in Python.
TH1F * m_hD0
Histogram of corrected d0 for each 2-track event.
TH2F * m_hTrackPointsLayer2
Histogram of intersection points for layer 2.
TH2F * m_hTrackClustersLayer1
Histogram of track matched clusters for layer 1.
TH2F * m_hTrackPointsLayer1
Histogram of intersection points for layer 1.
std::vector< unsigned long > m_nTrackClusters
Vector of number of track matched clusters.
virtual ~PXDValidationAlgorithm() override
Destructor.
TH2F * m_hTrackClustersLayer2
Histogram of track matched clusters for layer 2.
std::vector< unsigned long > m_nTrackPoints
Vector of number of track points.
std::vector< float > m_z0
Vector of delta z0 over sqrt(2)
bool saveD0Z0
Flag to save delta d0 (z0) over sqrt(2) to TTree;.
std::vector< unsigned long > m_nSelTrackPoints
Vector of number of track points outside of defective pixels.
PXDValidationAlgorithm()
Constructor set the prefix to PXDValidationAlgorithm.
int minTrackPoints
Minimum number of track points per sensor.
std::shared_ptr< TFile > m_file
Pointer for TFile.
bool save2DHists
Flag to save 2D histograms for efficiency;.
std::vector< unsigned short > m_uBin
Vector of uBin.
std::vector< float > m_d0
Vector of delta d0 over sqrt(2)
virtual EResult calibrate() override
Run algo on data.
int binsD0Z0
Bin size for d0/z0 histogram.
std::vector< unsigned long > m_nSelTrackClusters
Vector of number of track clusters outside of defective pixels.
std::vector< unsigned short > m_vBin
Vector of vBin.
std::shared_ptr< TTree > m_tree
Pointer for TTree of the validation info.
virtual bool isBoundaryRequired(const Calibration::ExpRun &) override
Decide if a run should be a payload boundary. Only used in certain Python Algorithm Starategies.
TH1F * m_hZ0
Histogram of corrected z0 for each 2-track event.
std::vector< unsigned short > m_pxdid
Vector of PXD module id (DHE id)
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
unsigned short getNumberOfSensors(const std::shared_ptr< TH1I > &histo_ptr)
Helper function to extract number of sensors from counter histogram labels.
void getNumberOfBins(const std::shared_ptr< TH1I > &histo_ptr, unsigned short &nBinsU, unsigned short &nBinsV)
Helper function to extract number of bins along u side and v side from counter histogram labels.
Abstract base class for different kinds of events.
STL namespace.