9#include <tracking/trackFindingVXD/sectorMapTools/SectorMapComparer.h>
10#include <tracking/dataobjects/FullSecID.h>
12#include <framework/logging/LogConfig.h>
32 const std::string& SMFileSecond) : m_SMFileName_first(SMFileFirst),
33 m_SMFileName_second(SMFileSecond)
42 std::string str = std::to_string(l1) +
"-" + std::to_string(l2) +
"-" + std::to_string(l3);
52 std::unordered_map<std::string, uint>& SecIDVals)
54 TObjArray* leafList = t->GetListOfLeaves();
55 for (TObject* o : *leafList) {
57 std::string name = l->GetName();
59 TString classname = l->ClassName();
62 if (classname.EqualTo(
"TLeafI")) {
64 l->SetAddress(&(SecIDVals[name]));
65 }
else if (classname.EqualTo(
"TLeafD")) {
66 filterVals[name] = 0.0;
67 l->SetAddress(&(filterVals[name]));
69 B2WARNING(
"Unsupported TLeaf type: " << l->ClassName() <<
". Will skip this TLeaf.");
84 std::vector<std::string> bNames = {
"innerFullSecID",
"centerFullSecID",
"outerFullSecID"};
85 std::vector<uint> bVals = {0, 0, 0};
87 for (uint i = 0; i < bNames.size(); i++) {
88 TLeaf* l = tree->GetLeaf(bNames[i].c_str());
90 l->SetAddress(&(bVals[i]));
95 for (
long i = 0; i < tree->GetEntries(); i++) {
97 std::string hash =
getHash(bVals[0], bVals[1], bVals[2]);
115 std::unordered_map<std::string, long> indexmap_t_second;
120 std::unordered_map<std::string, double> vals_t_first;
121 std::unordered_map<std::string, uint> ids_t_first;
124 std::unordered_map<std::string, double> vals_t_second;
125 std::unordered_map<std::string, uint> ids_t_second;
130 bool isTriplet = t_first->GetBranch(
"centerFullSecID") !=
nullptr;
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.");
138 TObjArray* leafList = t_first->GetListOfLeaves();
139 for (TObject* o : *leafList) {
140 TLeaf* l_first = (TLeaf*)o;
142 TString leafName = l_first->GetName();
143 if (leafName.Contains(
"FullSecID"))
continue;
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);
156 B2ERROR(
"ERROR: Histogram for key " << leafName.Data() <<
" already exists");
160 m_histo_map_first[ leafName.Data() ] = TH1F(leafName +
"_first", leafName, Nbins, min, max);
162 m_histo_map_second[ leafName.Data() ] = TH1F(leafName +
"_second", leafName, Nbins, min, max);
164 m_histo_map_diff[ leafName.Data() ] = TH1F(leafName +
"_diff",
"difference " + leafName +
" (SM1 - SM2)", Nbins, -0.2 * range,
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);
174 m_histo_map_range_second[ leafName.Data() ] = TH1F(leafName +
"_range_second", niceName +
" range (max - min)", Nbins, 0, range);
181 for (
long i = 0; i < t_first->GetEntries(); i++) {
182 t_first->GetEntry(i);
184 std::string hash =
getHash(ids_t_first[
"innerFullSecID"], ids_t_first[
"centerFullSecID"], ids_t_first[
"outerFullSecID"]);
193 if (indexmap_t_second.find(hash) == indexmap_t_second.end()) {
195 if (!unmatchedEntries)
continue;
198 t_second->GetEntry(indexmap_t_second[hash]);
199 if (unmatchedEntries)
continue;
208 for (TObject* o : *leafList) {
209 TString leafName = o->GetName();
211 if (leafName.Contains(
"FullSecID"))
continue;
213 if (!unmatchedEntries) {
215 m_histo_map_diff[ leafName.Data() ].Fill(vals_t_first[leafName.Data()] - vals_t_second[leafName.Data()]);
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()]);
223 vals_t_second[minLeafName.Data()]);
234 TFile* f = TFile::Open(secmapFileName);
236 B2WARNING(
"File not opened: " << secmapFileName);
246 TString* setupKeyName =
nullptr;
248 if (setupKeyName ==
nullptr) {
249 B2WARNING(
"setupKeyName not found");
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) {
256 std::cout << (*setupKeyName) << std::endl;
268 TList* keys = aDir->GetListOfKeys();
269 for (TObject* akey : *keys) {
272 TObject* o = aDir->Get(akey->GetName());
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()));
289 B2WARNING(
"You need to run SectorMapComparer::CompareMaps before this function!");
296 std::string aName = hist.first;
297 TCanvas* can =
new TCanvas((aName +
"_can").c_str(), aName.c_str());
300 can->cd(1)->SetLogy((
int)logScale);
303 TLegend* l =
new TLegend(0.15, 0.7, 0.3, 0.9);
305 l->AddEntry(hbuff,
"first SectorMap",
"lpf");
308 l->AddEntry(hbuff,
"second SectorMap",
"lpf");
309 if (drawLegend) l->Draw();
311 can->cd(2)->SetLogy((
int)logScale);;
314 if (TString(aName).Contains(
"_max")) {
315 can->cd(3)->SetLogy((
int)logScale);
322 if (pdfFileName !=
"") {
323 if (n == 1) can->SaveAs(pdfFileName);
325 if (counter == 1) can->SaveAs(pdfFileName +
"(");
326 else if (counter == n) can->SaveAs(pdfFileName +
")");
327 else can->SaveAs(pdfFileName);
339 B2INFO(
"comparing SectorMaps in the following files");
346 if (!f_first->IsOpen() || !f_first->IsOpen()) {
347 B2ERROR(
"ERROR: one of the files not open");
349 if (f_first) f_first->Close();
350 if (f_second) f_second->Close();
355 std::vector<std::string> listOfTrees;
361 for (
const std::string& tname_first : listOfTrees) {
364 if (!TString(tname_first).Contains(setupName))
continue;
368 std::string tname_second = tname_first;
373 TTree* t_first = (TTree*)f_first->Get(tname_first.c_str());
374 TTree* t_second = (TTree*)f_second->Get(tname_second.c_str());
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
388 B2INFO(
"Comparing tree: " << tname_first <<
" with tree: " << tname_second << std::endl);
390 if (!t_first || !t_second) {
391 B2ERROR(
"One of the trees not found!");
409 for (
const auto& entry : map) {
415 if (sector >= 0 && fullSID.
getSecID() != sector)
continue;
417 result += entry.second;
427 B2WARNING(
"You need to run SectorMapComparer::CompareMaps before this function!");
431 TH1F h_2hits(
"h_2hits",
"number of two hit connections starting at ith layer; layer; number connections", 7, -0.5, 6.5);
433 h_2hits.SetStats(
false);
435 TH1F h_2hits_matched(h_2hits);
436 h_2hits_matched.SetName(
"h2hits_matched");
437 h_2hits_matched.SetLineColor(kRed);
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);
446 for (
int ilayer = 1; ilayer < 7; ilayer++) {
449 int bin = h_2hits.FindBin(ilayer);
456 TCanvas* can =
new TCanvas(
"secStats",
"secStats");
458 can->cd(1)->SetLogy((
int)logScale);
460 TLegend* l =
new TLegend(0.15, 0.7, 0.3, 0.9);
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();
468 can->cd(2)->SetLogy((
int)logScale);
470 h_3hits_matched.DrawCopy(
"same");
475 std::vector< std::unordered_map<std::string, int> > sensorsInLayer2Hits(6);
478 sensorsInLayer2Hits[
id.getLayerNumber() - 1][(std::string)
id ] = 1;
480 std::vector< std::unordered_map<std::string, int> > sensorsInLayer3Hits(6);
483 sensorsInLayer3Hits[
id.getLayerNumber() - 1][(std::string)
id] = 1;
486 for (
int i = 1; i <= 6; i++) {
489 bin = h_3hits.FindBin(i);
491 int nSensors2Hits = sensorsInLayer2Hits[i - 1].size();
492 int nSensors3Hits = sensorsInLayer3Hits[i - 1].size();
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);
498 B2WARNING(
"No sensors found in layer " << i);
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);
505 B2WARNING(
"No sensors found in layer " << i);
509 TCanvas* can2 =
new TCanvas(
"ConPerSensor",
"Connections per sensor");
511 can2->cd(1)->SetLogy((
int)logScale);
513 TLegend* l2 =
new TLegend(0.15, 0.7, 0.3, 0.9);
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();
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");
525 h_3hits_matched.DrawCopy(
"same");
529 if (pdfFileName !=
"") {
530 can->SaveAs(pdfFileName +
"(");
531 can2->SaveAs(pdfFileName +
")");
Class to identify a sector inside of the VXD.
int getLadderNumber() const
returns LadderID compatible with basf2 standards
VxdID getVxdID() const
returns VxdID of sensor.
short int getLayerNumber() const
returns LayerID compatible with basf2 standards.
short int getSecID() const
returns SecID of current FullSecID (only unique for each sensor).
@ c_Info
Info: for informational messages, e.g.
void setLogLevel(ELogLevel logLevel)
Configure the log level.
LogConfig * getLogConfig()
Returns global log system configuration.
static LogSystem & Instance()
Static method to get a reference to the LogSystem instance.
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.
baseType getSensorNumber() const
Get the sensor id.
Abstract base class for different kinds of events.