Belle II Software development
SectorMapComparer Class Reference

A root tool that compares two Sectormaps (local root files) and produces some statistics output. More...

#include <SectorMapComparer.h>

Inheritance diagram for SectorMapComparer:

Public Member Functions

 SectorMapComparer ()=delete
 dont use this constructor
 
 SectorMapComparer (const std::string &SMFileFirst, const std::string &SMFileSecond)
 constructor
 
 ~SectorMapComparer ()=default
 default destructor
 
void plot (bool logScale=true, TString pdfFileName="", bool drawLegend=true)
 will create lots of Canvases!!!
 
void plotSectorStats (bool logScale=true, TString pdfFileName="", bool drawLegend=true)
 Will create canvases showing statistics.
 
void compareMaps (TString setup, bool unmatchedEntries=false)
 runs the comparison for all filters (and only the filters) in the sectormap files
 
void showSetups ()
 lists all setups of sectormaps included in the two files
 
void showSetups (TString sectorMapFileName)
 lists all sector map setups in the specified file
 
void showFiles ()
 lists the names of the files which are used
 

Private Member Functions

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::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 use the string as hash for the map.
 
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-directories
 
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.
 
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 (two) sector ids
 
void compareTrees (TTree *t_first, TTree *t_second, bool unmatchedEntries=false)
 makes the comparison of two trees, and fills histograms and certain maps
 
void clearMaps ()
 helper that clears all the maps used
 

Private Attributes

std::unordered_map< std::string, TH1F > m_histo_map_first
 histograms for the first sector map
 
std::unordered_map< std::string, TH1F > m_histo_map_second
 histograms for the second sector map
 
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_range_second
 histograms that store the range of that filter (max - min)
 
std::unordered_map< uint, uint > m_map_N2HitCombs
 count the number of 2Hit connections for each 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: meaning of "matched" depends on the setting)
 
std::unordered_map< uint, uint > m_map_N3HitCombs
 count the number of 3Hit connections for each sector, can be used to do some statistics
 
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: meaning of "matched" depends on the setting)
 
std::string m_SMFileName_first
 name of the file containing the first sectormap
 
std::string m_SMFileName_second
 name of the file containing the second sectormap
 
std::string m_setupsTreeName = "Setups"
 name of the tree all the setups are stored
 
std::string m_setupsBranchName = "name"
 name of the branch in the setups tree that holds the names of the setups
 
bool m_isCompared = false
 remember if CompareMaps already has been run, to give warnings if functions are called that need CompareMaps to have run
 

Detailed Description

A root tool that compares two Sectormaps (local root files) and produces some statistics output.

Yes we old people do use root. So get over it and open a root session to create an Object to let its magic unfold. This tool assumes that both sectormaps use the same setup (i.e. number of sectors on sensor, and filter configuration). It still may work with slightly different filters as long as the naming did not change. But if the number of sectors on sensor is different the matching of sectors will not work and the results may become misleading!

The output will be several plots comparing the filter cuts (min, max) for each of the contained filters. In addition a simple counting of sector connections is provided.

Per default only sector combinations which are matched between both sectormaps are drawn. If the option unmatchedEntries==true is set when calling CompareMaps only the sector combinations which are not matched from the first sector map are plotted.

Definition at line 38 of file SectorMapComparer.h.

Constructor & Destructor Documentation

◆ SectorMapComparer()

SectorMapComparer ( const std::string &  SMFileFirst,
const std::string &  SMFileSecond 
)

constructor

Parameters
SMFileFirstroot file name the first sector map is contained in
SMFileSecondroot file name the second sector map is contained in

Definition at line 31 of file SectorMapComparer.cc.

32 : m_SMFileName_first(SMFileFirst),
33 m_SMFileName_second(SMFileSecond)
34{
36};
@ 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
std::string m_SMFileName_second
name of the file containing the second sectormap

Member Function Documentation

◆ clearMaps()

void clearMaps ( )
inlineprivate

helper that clears all the maps used

Definition at line 140 of file SectorMapComparer.h.

141 {
142 m_histo_map_first.clear();
143 m_histo_map_second.clear();
144 m_histo_map_diff.clear();
147 m_map_N2HitCombs.clear();
149 m_map_N3HitCombs.clear();
151 }
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
std::unordered_map< std::string, TH1F > m_histo_map_second
histograms for the second sector map
std::unordered_map< uint, uint > m_map_N2HitCombs
count the number of 2Hit connections for each sector,
std::unordered_map< std::string, TH1F > m_histo_map_range_second
histograms that store the range of that filter (max - min)
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
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:...

◆ compareMaps()

void compareMaps ( TString  setup,
bool  unmatchedEntries = false 
)

runs the comparison for all filters (and only the filters) in the sectormap files

Definition at line 337 of file SectorMapComparer.cc.

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}
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-...
void clearMaps()
helper that clears all the maps used
void compareTrees(TTree *t_first, TTree *t_second, bool unmatchedEntries=false)
makes the comparison of two trees, and fills histograms and certain maps
bool m_isCompared
remember if CompareMaps already has been run, to give warnings if functions are called that need Comp...

◆ compareTrees()

void compareTrees ( TTree *  t_first,
TTree *  t_second,
bool  unmatchedEntries = false 
)
private

makes the comparison of two trees, and fills histograms and certain maps

Parameters
t_firstfirst tree
t_secondsecdon tree
unmatchedEntriesif true only the unmatched entries will be filled (only for the first tree)

Definition at line 111 of file SectorMapComparer.cc.

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 distinguish 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}
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 ...
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 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.

◆ countConnections()

uint countConnections ( const std::unordered_map< uint, uint > &  map,
int  layer = -1,
int  ladder = -1,
int  sensor = -1,
int  sector = -1 
)
private

counts all connections for given layer, ladder, sensor, sector.

If -1 is given for any of those it is summed over the corresponding entries This has to be called after the CompareMaps function, else the result will be zero

Parameters
mapreference to map mapping from fullsecid to number of times
layer,ladder,sensor,sectorlayer, ladder, sensor, sector number, if -1 it will be summed over

Definition at line 406 of file SectorMapComparer.cc.

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}
Class to identify a sector inside of the VXD.
Definition: FullSecID.h:33

◆ fillSectorToTreeIndexMap()

void fillSectorToTreeIndexMap ( TTree *  tree,
std::unordered_map< std::string, long > &  map 
)
private

fills a map that maps sector combinations in the tree to string which is a hash combining the three (two) sector ids

Parameters
treethe tree to be mapped
mapthe map which will be filled in the function

Definition at line 81 of file SectorMapComparer.cc.

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 existence 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}

◆ findTrees()

void findTrees ( TDirectory *  aDir,
std::vector< std::string > &  listOfTrees 
)
private

helper function that returns a list of names for all trees contained in the TDirectory including sub-directories

Parameters
aDira root directory pointer (e.g. TFile *)
listOfTreesthis vector will be filled with the names of the tree (NOTE: the name contains the full path (including file name))

Definition at line 266 of file SectorMapComparer.cc.

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}

◆ getHash()

std::string getHash ( long  l1,
long  l2,
long  l3 
)
private

I was too lacy to implement a proper hash function, so I convert the three numbers into a string and use the string as hash for the map.

Parameters
l1,l2,l3just three numbers of type long
Returns
: a string which combines the three numbers which can be used as hash as the hash function for std::string already exists

Definition at line 40 of file SectorMapComparer.cc.

41{
42 std::string str = std::to_string(l1) + "-" + std::to_string(l2) + "-" + std::to_string(l3);
43 return str;
44}

◆ plot()

void plot ( bool  logScale = true,
TString  pdfFileName = "",
bool  drawLegend = true 
)

will create lots of Canvases!!!

Parameters
logScaleif true the y-axis will be log scale
pdfFileNamename of the output file, has to be a file format which is supported by TCanvas::SaveAs (e.g. .pdf) if empty "" no output file will be created
drawLegendif true a legend is drawn, may hide interisting parts so one can deactivate it

Definition at line 285 of file SectorMapComparer.cc.

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}

◆ plotSectorStats()

void plotSectorStats ( bool  logScale = true,
TString  pdfFileName = "",
bool  drawLegend = true 
)

Will create canvases showing statistics.

Parameters
logScaleif true the y-axis will be log scale
pdfFileNamename of the output file, has to a file format which is supported by TCanvas::SaveAs (e.g. .pdf) if empty "" no output file will be created
drawLegendif true a legend is drawn, may hide interisting parts so one can deactivate it

Definition at line 424 of file SectorMapComparer.cc.

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}
VxdID getVxdID() const
returns VxdID of sensor.
Definition: FullSecID.h:138
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.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33

◆ setLeafAddresses()

void setLeafAddresses ( TTree *  t,
std::unordered_map< std::string, double > &  filterVals,
std::unordered_map< std::string, uint > &  SecIDVals 
)
private

automatically sets the TLeaf addresses for a tree.

Currently only TLeafI and TLeafD are supported

Parameters
tthe tree for setting the addresses
filterValsthe targets for the TLeafD addresses, will be filled in the function
SecIDValsthe targets for the TLeafI addresses, will be filled in the function

Definition at line 51 of file SectorMapComparer.cc.

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}

◆ showFiles()

void showFiles ( )
inline

lists the names of the files which are used

Definition at line 85 of file SectorMapComparer.h.

86 {
87 B2INFO("Following files will be compared:");
88 B2INFO("=================================");
89 B2INFO("First : " << m_SMFileName_first);
90 B2INFO("Second: " << m_SMFileName_second);
91 }

◆ showSetups() [1/2]

void showSetups ( )
inline

lists all setups of sectormaps included in the two files

Definition at line 73 of file SectorMapComparer.h.

74 {
75 std::cout << std::endl;
77 std::cout << std::endl;
79 }
void showSetups()
lists all setups of sectormaps included in the two files

◆ showSetups() [2/2]

void showSetups ( TString  sectorMapFileName)

lists all sector map setups in the specified file

Definition at line 232 of file SectorMapComparer.cc.

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}
std::string m_setupsBranchName
name of the branch in the setups tree that holds the names of the setups
std::string m_setupsTreeName
name of the tree all the setups are stored

Member Data Documentation

◆ m_histo_map_diff

std::unordered_map<std::string, TH1F> m_histo_map_diff
private

histograms that store differences between the two sector maps

Definition at line 157 of file SectorMapComparer.h.

◆ m_histo_map_first

std::unordered_map<std::string, TH1F> m_histo_map_first
private

histograms for the first sector map

Definition at line 155 of file SectorMapComparer.h.

◆ m_histo_map_range_first

std::unordered_map<std::string, TH1F> m_histo_map_range_first
private

histograms that store the range of that filter (max - min)

Definition at line 158 of file SectorMapComparer.h.

◆ m_histo_map_range_second

std::unordered_map<std::string, TH1F> m_histo_map_range_second
private

histograms that store the range of that filter (max - min)

Definition at line 159 of file SectorMapComparer.h.

◆ m_histo_map_second

std::unordered_map<std::string, TH1F> m_histo_map_second
private

histograms for the second sector map

Definition at line 156 of file SectorMapComparer.h.

◆ m_isCompared

bool m_isCompared = false
private

remember if CompareMaps already has been run, to give warnings if functions are called that need CompareMaps to have run

Definition at line 180 of file SectorMapComparer.h.

◆ m_map_N2HitCombs

std::unordered_map<uint, uint> m_map_N2HitCombs
private

count the number of 2Hit connections for each sector,

Definition at line 162 of file SectorMapComparer.h.

◆ m_map_N2HitCombs_matched

std::unordered_map<uint, uint> m_map_N2HitCombs_matched
private

count the number of 2Hit connections for each sector, only if sector combination was "matched" (Note: meaning of "matched" depends on the setting)

Definition at line 164 of file SectorMapComparer.h.

◆ m_map_N3HitCombs

std::unordered_map<uint, uint> m_map_N3HitCombs
private

count the number of 3Hit connections for each sector, can be used to do some statistics

Definition at line 167 of file SectorMapComparer.h.

◆ m_map_N3HitCombs_matched

std::unordered_map<uint, uint> m_map_N3HitCombs_matched
private

count the number of 3Hit connections for each sector, only if sector combination was "matched" (Note: meaning of "matched" depends on the setting)

Definition at line 169 of file SectorMapComparer.h.

◆ m_setupsBranchName

std::string m_setupsBranchName = "name"
private

name of the branch in the setups tree that holds the names of the setups

Definition at line 177 of file SectorMapComparer.h.

◆ m_setupsTreeName

std::string m_setupsTreeName = "Setups"
private

name of the tree all the setups are stored

Definition at line 176 of file SectorMapComparer.h.

◆ m_SMFileName_first

std::string m_SMFileName_first
private

name of the file containing the first sectormap

Definition at line 173 of file SectorMapComparer.h.

◆ m_SMFileName_second

std::string m_SMFileName_second
private

name of the file containing the second sectormap

Definition at line 174 of file SectorMapComparer.h.


The documentation for this class was generated from the following files: