Belle II Software development
SectorMapComparer.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 <tracking/trackFindingVXD/sectorMapTools/SectorMapComparer.h>
10#include <tracking/dataobjects/FullSecID.h>
11
12#include <framework/logging/LogConfig.h>
13
14#include <TFile.h>
15#include <TLeaf.h>
16#include <TLeafI.h>
17#include <TLeafD.h>
18#include <TCanvas.h>
19#include <TObjArray.h>
20#include <TObject.h>
21#include <TLegend.h>
22
23#include <vector>
24#include <iostream>
25
26
27
28
29using namespace Belle2;
30
31SectorMapComparer::SectorMapComparer(const std::string& SMFileFirst,
32 const std::string& SMFileSecond) : m_SMFileName_first(SMFileFirst),
33 m_SMFileName_second(SMFileSecond)
34{
36};
37
38
39std::string
40SectorMapComparer::getHash(long l1, long l2, long l3)
41{
42 std::string str = std::to_string(l1) + "-" + std::to_string(l2) + "-" + std::to_string(l3);
43 return str;
44}
45
46
47
48// sets the addresses of the leafes
49// assumes each leaf name is unique! (Not sure if root takes care for that)
50void
51SectorMapComparer::setLeafAddresses(TTree* t, std::unordered_map<std::string, double>& filterVals,
52 std::unordered_map<std::string, uint>& SecIDVals)
53{
54 TObjArray* leafList = t->GetListOfLeaves();
55 for (TObject* o : *leafList) {
56 TLeaf* l = (TLeaf*)o;
57 std::string name = l->GetName();
58 // cannot use TClass pointer here as it is TLeaf and not TLeafD or TLeafI, ClassName is virtual (gives the name of the derived)
59 TString classname = l->ClassName();
60
61 // filter leafs of type TLeafI (sector IDs) from the
62 if (classname.EqualTo("TLeafI")) {
63 SecIDVals[name] = 0;
64 l->SetAddress(&(SecIDVals[name]));
65 } else if (classname.EqualTo("TLeafD")) {
66 filterVals[name] = 0.0;
67 l->SetAddress(&(filterVals[name]));
68 } else {
69 B2WARNING("Unsupported TLeaf type: " << l->ClassName() << ". Will skip this TLeaf.");
70 }
71
72 }
73}
74
75
76
77
78// returns a map of sector combinations to the position in the branch
79// WARNING: this messes up the addresses!!!!!!!
80void
81SectorMapComparer::fillSectorToTreeIndexMap(TTree* tree, std::unordered_map<std::string, long>& map)
82{
83 // the branchnames which store the sectorIDs
84 std::vector<std::string> bNames = {"innerFullSecID", "centerFullSecID", "outerFullSecID"};
85 std::vector<uint> bVals = {0, 0, 0}; // lets hope there is no inner sector id 0! But need to be 0 here as later no check for existance in the map for segments, so default value is used.
86
87 for (uint i = 0; i < bNames.size(); i++) {
88 TLeaf* l = tree->GetLeaf(bNames[i].c_str());
89 if (l) {
90 l->SetAddress(&(bVals[i]));
91 }
92 }
93
94 // now loop to create the map
95 for (long i = 0; i < tree->GetEntries(); i++) {
96 tree->GetEntry(i);
97 std::string hash = getHash(bVals[0], bVals[1], bVals[2]);
98 //std::cout << hash << " " << i << std::endl;
99 map[hash] = i;
100 }
101
102}
103
104
105
106
107
108
109// actually compares the trees (segment and triplet filters)
110void
111SectorMapComparer::compareTrees(TTree* t_first, TTree* t_second, bool unmatchedEntries)
112{
113
114 // index the second tree
115 std::unordered_map<std::string, long> indexmap_t_second;
116 // WARNING: this messes with leaf addresses of the sector ids
117 fillSectorToTreeIndexMap(t_second, indexmap_t_second);
118
119
120 std::unordered_map<std::string, double> vals_t_first;
121 std::unordered_map<std::string, uint> ids_t_first;
122 setLeafAddresses(t_first, vals_t_first, ids_t_first);
123
124 std::unordered_map<std::string, double> vals_t_second;
125 std::unordered_map<std::string, uint> ids_t_second;
126 setLeafAddresses(t_second, vals_t_second, ids_t_second);
127
128
129 // identify Three hit combinations by existens of the branch centerID
130 bool isTriplet = t_first->GetBranch("centerFullSecID") != nullptr;
131
132 // some cross checks
133 if (vals_t_second.size() != vals_t_first.size()) {
134 B2WARNING("Number of filters stored in the two SectorMaps seem to differ! This is per se not dangerous, but some cuts may be compared to zero.");
135 }
136
137 // Creating histograms, using the list of leaves for indexing the histograms
138 TObjArray* leafList = t_first->GetListOfLeaves();
139 for (TObject* o : *leafList) {
140 TLeaf* l_first = (TLeaf*)o;
141 // no histograms for the sector ids
142 TString leafName = l_first->GetName();
143 if (leafName.Contains("FullSecID")) continue;
144
145 double min = std::min(t_first->GetMinimum(leafName), t_second->GetMinimum(leafName));
146 double max = std::max(t_first->GetMaximum(leafName), t_second->GetMaximum(leafName));
147 double range = std::fabs(min - max);
148 // TODO: make this configurable
149 int Nbins = 100;
150
151 // As in the current version the maps are used for several trees, cross check that histograms are not overwritten
152 if (m_histo_map_first.find(leafName.Data()) != m_histo_map_first.end() ||
153 m_histo_map_second.find(leafName.Data()) != m_histo_map_second.end() ||
154 m_histo_map_diff.find(leafName.Data()) != m_histo_map_diff.end()) {
155
156 B2ERROR("ERROR: Histogram for key " << leafName.Data() << " already exists");
157 }
158
159
160 m_histo_map_first[ leafName.Data() ] = TH1F(leafName + "_first", leafName, Nbins, min, max);
161 m_histo_map_first[ leafName.Data() ].SetDirectory(nullptr);
162 m_histo_map_second[ leafName.Data() ] = TH1F(leafName + "_second", leafName, Nbins, min, max);
163 m_histo_map_second[ leafName.Data() ].SetDirectory(nullptr);
164 m_histo_map_diff[ leafName.Data() ] = TH1F(leafName + "_diff", "difference " + leafName + " (SM1 - SM2)", Nbins, -0.2 * range,
165 0.2 * range);
166 m_histo_map_diff[ leafName.Data() ].SetDirectory(nullptr);
167
168 // histograms containing the max minus the min
169 TString niceName = leafName;
170 niceName.ReplaceAll("_min", "");
171 niceName.ReplaceAll("_max", "");
172 m_histo_map_range_first[ leafName.Data() ] = TH1F(leafName + "_range_first", niceName + " range (max - min)", Nbins, 0, range);
173 m_histo_map_range_first[ leafName.Data() ].SetDirectory(nullptr);
174 m_histo_map_range_second[ leafName.Data() ] = TH1F(leafName + "_range_second", niceName + " range (max - min)", Nbins, 0, range);
175 m_histo_map_range_second[ leafName.Data() ].SetDirectory(nullptr);
176 }
177
178
179
180 // loop over t1, and get corresponding value for t2
181 for (long i = 0; i < t_first->GetEntries(); i++) {
182 t_first->GetEntry(i);
183
184 std::string hash = getHash(ids_t_first["innerFullSecID"], ids_t_first["centerFullSecID"], ids_t_first["outerFullSecID"]);
185
186
187 // count all connections
188 if (isTriplet) m_map_N3HitCombs[ ids_t_first["outerFullSecID"] ]++;
189 else m_map_N2HitCombs[ ids_t_first["outerFullSecID"] ]++;
190
191
192 // filter connections inside both sector maps, if unmatchedEntries is true only the unmatched entries will be selected
193 if (indexmap_t_second.find(hash) == indexmap_t_second.end()) {
194 // case this sector combination not matched, can be thrown away if unmatchedEntries is set false
195 if (!unmatchedEntries) continue;
196 } else {
197 // sector combination is matched, read second tree and skip if unmatchedEntries is set
198 t_second->GetEntry(indexmap_t_second[hash]);
199 if (unmatchedEntries) continue;
200 }
201
202 // only count the matched ones
203 if (isTriplet) m_map_N3HitCombs_matched[ ids_t_first["outerFullSecID"] ]++;
204 else m_map_N2HitCombs_matched[ ids_t_first["outerFullSecID"] ]++;
205
206 // now fill the histograms
207 // NOTE: in case unmatchedEntries==true the values in the the second tree do not make any sense, so not fill any histograms with them
208 for (TObject* o : *leafList) {
209 TString leafName = o->GetName();
210 // TODO: find a better way to destinguish the secid leaves from the others
211 if (leafName.Contains("FullSecID")) continue;
212 m_histo_map_first[ leafName.Data() ].Fill(vals_t_first[leafName.Data()]);
213 if (!unmatchedEntries) {
214 m_histo_map_second[ leafName.Data() ].Fill(vals_t_second[leafName.Data()]);
215 m_histo_map_diff[ leafName.Data() ].Fill(vals_t_first[leafName.Data()] - vals_t_second[leafName.Data()]);
216 }
217 // calculate the range only for the max
218 if (leafName.Contains("_max")) {
219 TString minLeafName = leafName;
220 minLeafName.ReplaceAll("_max", "_min");
221 m_histo_map_range_first[ leafName.Data() ].Fill(vals_t_first[leafName.Data()] - vals_t_first[minLeafName.Data()]);
222 if (!unmatchedEntries) m_histo_map_range_second[ leafName.Data() ].Fill(vals_t_second[leafName.Data()] -
223 vals_t_second[minLeafName.Data()]);
224 }
225 }
226
227 }
228
229}
230
231void
232SectorMapComparer::showSetups(TString secmapFileName)
233{
234 TFile* f = TFile::Open(secmapFileName);
235 if (!f->IsOpen()) {
236 B2WARNING("File not opened: " << secmapFileName);
237 return;
238 }
239
240 TTree* t = (TTree*)f->Get(m_setupsTreeName.c_str());
241 if (t == nullptr) {
242 B2WARNING("tree not found! tree name: " << m_setupsTreeName);
243 return;
244 }
245
246 TString* setupKeyName = nullptr;
247 t->SetBranchAddress(m_setupsBranchName.c_str(), & setupKeyName);
248 if (setupKeyName == nullptr) {
249 B2WARNING("setupKeyName not found");
250 return;
251 } else {
252 std::cout << "Following setups found for file: " << secmapFileName << std::endl;
253 std::cout << "================================ " << std::endl;
254 for (long i = 0; i < t->GetEntries(); ++i) {
255 t->GetEntry(i);
256 std::cout << (*setupKeyName) << std::endl;
257 }
258
259 delete setupKeyName;
260 }
261}
262
263// fills the listOfTrees with the names of all trees in this directory (including their full root-path)
264// loops recursively over all subdirectories
265void
266SectorMapComparer::findTrees(TDirectory* aDir, std::vector<std::string>& listOfTrees)
267{
268 TList* keys = aDir->GetListOfKeys();
269 for (TObject* akey : *keys) {
270
271 // there should be no check needed as each key should be attached to an object
272 TObject* o = aDir->Get(akey->GetName());
273
274
275 if (o->InheritsFrom(TDirectory::Class())) findTrees((TDirectory*)o, listOfTrees);
276 if (o->InheritsFrom(TTree::Class())) {
277 listOfTrees.emplace_back(std::string(aDir->GetPath()) + std::string("/") + std::string(o->GetName()));
278 }
279 }
280
281}
282
283
284void
285SectorMapComparer::plot(bool logScale, TString pdfFileName, bool drawLegend)
286{
287
288 if (!m_isCompared) {
289 B2WARNING("You need to run SectorMapComparer::CompareMaps before this function!");
290 return;
291 }
292
293 int n = m_histo_map_first.size();
294 int counter = 1;
295 for (const auto& hist : m_histo_map_first) {
296 std::string aName = hist.first;
297 TCanvas* can = new TCanvas((aName + "_can").c_str(), aName.c_str());
298 can->Divide(1, 3);
299
300 can->cd(1)->SetLogy((int)logScale);
301 // TODO: pretty sure this is a memory leak, but right now I do not care
302 // draw legend only for first
303 TLegend* l = new TLegend(0.15, 0.7, 0.3, 0.9);
304 TH1* hbuff = m_histo_map_first[ aName ].DrawCopy();
305 l->AddEntry(hbuff, "first SectorMap", "lpf");
306 m_histo_map_second[ aName ].SetLineColor(kRed);
307 hbuff = m_histo_map_second[ aName ].DrawCopy("same");
308 l->AddEntry(hbuff, "second SectorMap", "lpf");
309 if (drawLegend) l->Draw();
310
311 can->cd(2)->SetLogy((int)logScale);;
312 m_histo_map_diff[ aName ].DrawCopy();
313
314 if (TString(aName).Contains("_max")) {
315 can->cd(3)->SetLogy((int)logScale);
316 m_histo_map_range_first[ aName ].DrawCopy();
317 m_histo_map_range_second[ aName ].SetLineColor(kRed);
318 m_histo_map_range_second[ aName ].DrawCopy("same");
319 }
320
321 // if file name is specified plot directly to pdf file
322 if (pdfFileName != "") {
323 if (n == 1) can->SaveAs(pdfFileName);
324 else {
325 if (counter == 1) can->SaveAs(pdfFileName + "(");
326 else if (counter == n) can->SaveAs(pdfFileName + ")");
327 else can->SaveAs(pdfFileName);
328 }
329 }
330 counter++;
331 }
332
333}
334
335
336void
337SectorMapComparer::compareMaps(TString setupName, bool unmatchedEntries)
338{
339 B2INFO("comparing SectorMaps in the following files");
340 B2INFO("first: " << m_SMFileName_first);
341 B2INFO("second: " << m_SMFileName_second);
342
343 TFile* f_first = TFile::Open(m_SMFileName_first.c_str());
344 TFile* f_second = TFile::Open(m_SMFileName_second.c_str());
345
346 if (!f_first->IsOpen() || !f_first->IsOpen()) {
347 B2ERROR("ERROR: one of the files not open");
348
349 if (f_first) f_first->Close();
350 if (f_second) f_second->Close();
351 return;
352 }
353
354 // get the list of trees from the first file (for second file assume same trees)
355 std::vector<std::string> listOfTrees;
356 findTrees(f_first, listOfTrees);
357
358 // clear the maps, else old results will be mixed with this run
359 clearMaps();
360
361 for (const std::string& tname_first : listOfTrees) {
362
363 // only look at the trees for the specified Setup
364 if (!TString(tname_first).Contains(setupName)) continue;
365
366
367 // assume second file has
368 std::string tname_second = tname_first;
369 tname_second.replace(0, m_SMFileName_first.length(), m_SMFileName_second);
370
371
372
373 TTree* t_first = (TTree*)f_first->Get(tname_first.c_str());
374 TTree* t_second = (TTree*)f_second->Get(tname_second.c_str());
375
376 // filter out the sector id and test for consistency
377 if (TString(tname_first).Contains("CompactSecIDs")) {
378 if (t_first->GetEntries() != t_second->GetEntries())
379 B2FATAL("The number of entries in CompactSecIDs should be identical! This indicates that two sectormaps with"
380 " different sectors on sensors are compared which will lead to wrong results! " << std::endl <<
381 "# CompactSecIDs SM1: " << t_first->GetEntries() << std::endl <<
382 "# CompactSecIDs SM2: " << t_second->GetEntries() << std::endl
383 );
384 // we dont compare those
385 continue;
386 }
387
388 B2INFO("Comparing tree: " << tname_first << " with tree: " << tname_second << std::endl);
389
390 if (!t_first || !t_second) {
391 B2ERROR("One of the trees not found!");
392 continue;
393 }
394
395 compareTrees(t_first, t_second, unmatchedEntries);
396 }
397
398 f_first->Close();
399 f_second->Close();
400
401 m_isCompared = true;
402}
403
404
405uint
406SectorMapComparer::countConnections(const std::unordered_map<uint, uint>& map, int layer, int ladder, int sensor, int sector)
407{
408 uint result = 0;
409 for (const auto& entry : map) {
410 // convert to FullSecID
411 FullSecID fullSID(entry.first);
412 if (layer >= 0 && fullSID.getLayerNumber() != layer) continue;
413 if (ladder >= 0 && fullSID.getLadderNumber() != ladder) continue;
414 if (sensor >= 0 && fullSID.getVxdID().getSensorNumber() != sensor) continue;
415 if (sector >= 0 && fullSID.getSecID() != sector) continue;
416
417 result += entry.second;
418 }
419 return result;
420}
421
422
423void
424SectorMapComparer::plotSectorStats(bool logScale, TString pdfFileName, bool drawLegend)
425{
426 if (!m_isCompared) {
427 B2WARNING("You need to run SectorMapComparer::CompareMaps before this function!");
428 return;
429 }
430
431 TH1F h_2hits("h_2hits", "number of two hit connections starting at ith layer; layer; number connections", 7, -0.5, 6.5);
432 // stats box does not give any useful information, so deactivate it
433 h_2hits.SetStats(false);
434 // make some copies
435 TH1F h_2hits_matched(h_2hits);
436 h_2hits_matched.SetName("h2hits_matched");
437 h_2hits_matched.SetLineColor(kRed);
438 // make copies for the three hit connections
439 TH1F h_3hits(h_2hits);
440 h_3hits.SetName("h_3hits");
441 h_3hits.SetTitle("number of three hit connections starting at ith layer; layer; number connections");
442 TH1F h_3hits_matched(h_3hits);
443 h_3hits_matched.SetName("h_3hits_matched");
444 h_3hits_matched.SetLineColor(kRed);
445
446 for (int ilayer = 1; ilayer < 7; ilayer++) {
447
448 // the bin is the same for all
449 int bin = h_2hits.FindBin(ilayer);
450 h_2hits.SetBinContent(bin, countConnections(m_map_N2HitCombs, ilayer, -1, -1, -1));
451 h_2hits_matched.SetBinContent(bin, countConnections(m_map_N2HitCombs_matched, ilayer, -1, -1, -1));
452 h_3hits.SetBinContent(bin, countConnections(m_map_N3HitCombs, ilayer, -1, -1, -1));
453 h_3hits_matched.SetBinContent(bin, countConnections(m_map_N3HitCombs_matched, ilayer, -1, -1, -1));
454 }
455
456 TCanvas* can = new TCanvas("secStats", "secStats");
457 can->Divide(1, 2);
458 can->cd(1)->SetLogy((int)logScale);
459 // TODO: potential memory leak, but I dont care
460 TLegend* l = new TLegend(0.15, 0.7, 0.3, 0.9);
461
462 TH1* hbuff = h_2hits.DrawCopy();
463 l->AddEntry(hbuff, "all from first map", "lpf");
464 hbuff = h_2hits_matched.DrawCopy("same");
465 l->AddEntry(hbuff, "only matched ", "lpf");
466 if (drawLegend) l->Draw();
467
468 can->cd(2)->SetLogy((int)logScale);
469 h_3hits.DrawCopy();
470 h_3hits_matched.DrawCopy("same");
471
472
473
474 // a stupid way of counting the number of sensors contained without having to invoke the geometry
475 std::vector< std::unordered_map<std::string, int> > sensorsInLayer2Hits(6);
476 for (const auto& entry : m_map_N2HitCombs) {
477 VxdID id = FullSecID(entry.first).getVxdID();
478 sensorsInLayer2Hits[id.getLayerNumber() - 1][(std::string)id ] = 1;
479 }
480 std::vector< std::unordered_map<std::string, int> > sensorsInLayer3Hits(6);
481 for (const auto& entry : m_map_N3HitCombs) {
482 VxdID id = FullSecID(entry.first).getVxdID();
483 sensorsInLayer3Hits[id.getLayerNumber() - 1][(std::string)id] = 1;
484 }
485 // Reusing old histograms!
486 for (int i = 1; i <= 6; i++) {
487 // assumes binning is still the same for all hists!
488 int bin = 0;
489 bin = h_3hits.FindBin(i);
490
491 int nSensors2Hits = sensorsInLayer2Hits[i - 1].size();
492 int nSensors3Hits = sensorsInLayer3Hits[i - 1].size();
493
494 if (nSensors2Hits > 0) {
495 h_2hits.SetBinContent(bin, h_2hits.GetBinContent(bin) / nSensors2Hits);
496 h_2hits_matched.SetBinContent(bin, h_2hits_matched.GetBinContent(bin) / nSensors2Hits);
497 } else {
498 B2WARNING("No sensors found in layer " << i);
499 }
500
501 if (nSensors3Hits > 0) {
502 h_3hits.SetBinContent(bin, h_3hits.GetBinContent(bin) / nSensors3Hits);
503 h_3hits_matched.SetBinContent(bin, h_3hits_matched.GetBinContent(bin) / nSensors3Hits);
504 } else {
505 B2WARNING("No sensors found in layer " << i);
506 }
507 }
508
509 TCanvas* can2 = new TCanvas("ConPerSensor", "Connections per sensor");
510 can2->Divide(1, 2);
511 can2->cd(1)->SetLogy((int)logScale);
512 // TODO: potential memory leak, but I dont care
513 TLegend* l2 = new TLegend(0.15, 0.7, 0.3, 0.9);
514
515 h_2hits.SetTitle("average number of 2Hit combinations per sensor starting at ith layer; layer; avg. connections/sensor");
516 hbuff = h_2hits.DrawCopy();
517 l2->AddEntry(hbuff, "all from first map", "lpf");
518 hbuff = h_2hits_matched.DrawCopy("same");
519 l2->AddEntry(hbuff, "only matched ", "lpf");
520 if (drawLegend) l2->Draw();
521
522 can2->cd(2)->SetLogy((int)logScale);
523 h_3hits.SetTitle("average number of 3Hit combinations per sensor starting at ith layer; layer; avg. connections/sensor");
524 h_3hits.DrawCopy();
525 h_3hits_matched.DrawCopy("same");
526
527
528
529 if (pdfFileName != "") {
530 can->SaveAs(pdfFileName + "(");
531 can2->SaveAs(pdfFileName + ")");
532 }
533
534}
535
536
537
538
Class to identify a sector inside of the VXD.
Definition: FullSecID.h:33
int getLadderNumber() const
returns LadderID compatible with basf2 standards
Definition: FullSecID.h:130
VxdID getVxdID() const
returns VxdID of sensor.
Definition: FullSecID.h:138
short int getLayerNumber() const
returns LayerID compatible with basf2 standards.
Definition: FullSecID.h:122
short int getSecID() const
returns SecID of current FullSecID (only unique for each sensor).
Definition: FullSecID.h:146
@ c_Info
Info: for informational messages, e.g.
Definition: LogConfig.h:27
void setLogLevel(ELogLevel logLevel)
Configure the log level.
Definition: LogConfig.cc:25
LogConfig * getLogConfig()
Returns global log system configuration.
Definition: LogSystem.h:78
static LogSystem & Instance()
Static method to get a reference to the LogSystem instance.
Definition: LogSystem.cc:31
std::string m_SMFileName_first
name of the file containing the first sectormap
void findTrees(TDirectory *aDir, std::vector< std::string > &listOfTrees)
helper function that returns a list of names for all trees contained in the TDirectory including sub-...
std::unordered_map< uint, uint > m_map_N3HitCombs_matched
count the number of 3Hit connections for each sector, only if sector combination was "matched" (Note:...
std::unordered_map< uint, uint > m_map_N3HitCombs
count the number of 3Hit connections for each sector, can be used to do some statistics
void compareMaps(TString setup, bool unmatchedEntries=false)
runs the comparison for all filters (and only the filters) in the sectormap files
std::string getHash(long l1, long l2, long l3)
I was too lacy to implement a proper hash function, so I convert the three numbers into a string and ...
std::string m_setupsBranchName
name of the branch in the setups tree that holds the names of the setups
std::string m_SMFileName_second
name of the file containing the second sectormap
std::unordered_map< std::string, TH1F > m_histo_map_second
histograms for the second sector map
void fillSectorToTreeIndexMap(TTree *tree, std::unordered_map< std::string, long > &map)
fills a map that maps sector combinations in the tree to string which is a hash combining the three (...
void clearMaps()
helper that clears all the maps used
std::unordered_map< uint, uint > m_map_N2HitCombs
count the number of 2Hit connections for each sector,
void showSetups()
lists all setups of sectormaps included in the two files
std::unordered_map< std::string, TH1F > m_histo_map_range_second
histograms that store the range of that filter (max - min)
void plot(bool logScale=true, TString pdfFileName="", bool drawLegend=true)
will create lots of Canvases!!!
std::unordered_map< std::string, TH1F > m_histo_map_diff
histograms that store differences between the two sector maps
std::unordered_map< std::string, TH1F > m_histo_map_range_first
histograms that store the range of that filter (max - min)
std::unordered_map< std::string, TH1F > m_histo_map_first
histograms for the first sector map
void setLeafAddresses(TTree *t, std::unordered_map< std::string, double > &filterVals, std::unordered_map< std::string, uint > &SecIDVals)
automatically sets the TLeaf addresses for a tree.
SectorMapComparer()=delete
dont use this constructor
void compareTrees(TTree *t_first, TTree *t_second, bool unmatchedEntries=false)
makes the comparison of two trees, and fills histograms and certain maps
uint countConnections(const std::unordered_map< uint, uint > &map, int layer=-1, int ladder=-1, int sensor=-1, int sector=-1)
counts all connections for given layer, ladder, sensor, sector.
std::unordered_map< uint, uint > m_map_N2HitCombs_matched
count the number of 2Hit connections for each sector, only if sector combination was "matched" (Note:...
std::string m_setupsTreeName
name of the tree all the setups are stored
bool m_isCompared
remember if CompareMaps already has been run, to give warnings if functions are called that need Comp...
void plotSectorStats(bool logScale=true, TString pdfFileName="", bool drawLegend=true)
Will create canvases showing statistics.
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
Abstract base class for different kinds of events.