Belle II Software development
DebugableSimpleBoxDivisionHoughTree.h
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#pragma once
9#include <tracking/trackFindingCDC/hough/trees/SimpleBoxDivisionHoughTree.h>
10#include <TTree.h>
11#include <TFile.h>
12
13#include <TGraph.h>
14#include <TF1.h>
15#include <TCanvas.h>
16#include <TAxis.h>
17
18namespace Belle2 {
23 namespace TrackFindingCDC {
24
26 template<class AHitPtr, class AInBoxAlgorithm, size_t divisionX, size_t divisionY>
28 public SimpleBoxDivisionHoughTree<AHitPtr, AInBoxAlgorithm, divisionX, divisionY> {
29 private:
32 public:
34 using Super::Super;
35
42 void writeDebugInfoToFile(const std::string& filename)
43 {
44 fillAll();
45
46 TFile openedRootFile(filename.c_str(), "RECREATE");
47 TTree weightTTree("weightTree", "A tree with the weights of the box items.");
48 TTree eventTTree("eventTree", "A tree with event information.");
49
50 double lowerX, upperX, lowerY, upperY, weight, level;
51 weightTTree.Branch("lowerX", &lowerX);
52 weightTTree.Branch("upperY", &upperY);
53 weightTTree.Branch("lowerY", &lowerY);
54 weightTTree.Branch("upperX", &upperX);
55 weightTTree.Branch("weight", &weight);
56 weightTTree.Branch("level", &level);
57
58 double lowerLimX = -Super::getMaximumX();
59 double upperLimX = Super::getMaximumX();
60 double lowerLimY = -Super::getMaximumY();
61 double upperLimY = Super::getMaximumY();
62 double maxLevel = Super::getMaxLevel();
63 eventTTree.Branch("lowerLimX", &lowerLimX);
64 eventTTree.Branch("upperLimX", &upperLimX);
65 eventTTree.Branch("lowerLimY", &lowerLimY);
66 eventTTree.Branch("upperLimY", &upperLimY);
67 eventTTree.Branch("maxLevel", &maxLevel);
68 eventTTree.Fill();
69 eventTTree.Write();
70
71 auto walker = [&](const typename Super::Node * node) -> bool {
72 // cppcheck-suppress unreadVariable
73 lowerX = node->getLowerX();
74 // cppcheck-suppress unreadVariable
75 upperX = node->getUpperX();
76 // cppcheck-suppress unreadVariable
77 lowerY = node->getLowerY();
78 // cppcheck-suppress unreadVariable
79 upperY = node->getUpperY();
80 // cppcheck-suppress unreadVariable
81 weight = node->getWeight();
82 // cppcheck-suppress unreadVariable
83 level = node->getLevel();
84
85 weightTTree.Fill();
86
87 // Always return true, as we want to access every node
88 return true;
89 };
90 Super::getTree()->walk(walker);
91
92 weightTTree.Write();
93 openedRootFile.Write();
94 openedRootFile.Close();
95 }
96
98 void drawDebugPlot(const std::vector<CDCRecoHit3D>& allHits,
99 const std::vector<CDCRecoHit3D>& foundHits,
100 const typename AInBoxAlgorithm::HoughBox& node)
101 {
102 TGraph* allHitsGraph = new TGraph();
103 allHitsGraph->SetLineWidth(2);
104 allHitsGraph->SetLineColor(9);
105
106 for (const CDCRecoHit3D& recoHit3D : allHits) {
107 const Vector3D& recoPos3D = recoHit3D.getRecoPos3D();
108 const double R = std::sqrt(recoPos3D.x() * recoPos3D.x() + recoPos3D.y() * recoPos3D.y());
109 const double Z = recoPos3D.z();
110 allHitsGraph->SetPoint(allHitsGraph->GetN(), R, Z);
111 }
112
113 static int nevent(0);
114 TCanvas canv("trackCanvas", "CDC stereo hits in an event", 0, 0, 1600, 1200);
115 canv.cd();
116 allHitsGraph->Draw("APL*");
117 allHitsGraph->GetXaxis()->SetLimits(0, 120);
118 allHitsGraph->GetYaxis()->SetRangeUser(-180, 180);
119
120 TGraph* foundHitsGraph = new TGraph();
121 foundHitsGraph->SetMarkerStyle(8);
122 foundHitsGraph->SetMarkerColor(2);
123
124 for (const CDCRecoHit3D& recoHit3D : foundHits) {
125 const Vector3D& recoPos3D = recoHit3D.getRecoPos3D();
126 const double R = std::sqrt(recoPos3D.x() * recoPos3D.x() + recoPos3D.y() * recoPos3D.y());
127 const double Z = recoPos3D.z();
128 foundHitsGraph->SetPoint(foundHitsGraph->GetN(), R, Z);
129 }
130 foundHitsGraph->Draw("P");
131
132 const double xMean = (node.getLowerX() + node.getUpperX()) / 2.0; //Z0 or Z1
133 const double yMean = (node.getLowerY() + node.getUpperY()) / 2.0; //tanLambda or Z2
134 const double xLow = node.getLowerX();
135 const double yLow = node.getLowerY();
136 const double xHigh = node.getUpperX();
137 const double yHigh = node.getUpperY();
138
139 TF1* candidateLL = new TF1("candLL", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
140 TF1* candidateLH = new TF1("candLH", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
141 TF1* candidateHL = new TF1("candHL", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
142 TF1* candidateHH = new TF1("candHH", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
143 TF1* candidateMean = new TF1("candMean", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
144
145 candidateLL->SetParameters(xLow, yLow);
146 candidateLH->SetParameters(xLow, yHigh);
147 candidateHL->SetParameters(xHigh, yLow);
148 candidateHH->SetParameters(xHigh, yHigh);
149 candidateMean->SetParameters(xMean, yMean);
150
151 candidateLL->SetLineColor(9);
152 candidateLH->SetLineColor(30);
153 candidateHL->SetLineColor(46);
154 candidateHH->SetLineColor(41);
155 candidateMean->SetLineColor(2);
156
157 candidateLL->Draw("same");
158 candidateHL->Draw("same");
159 candidateLH->Draw("same");
160 candidateHH->Draw("same");
161 candidateMean->Draw("same");
162 canv.SaveAs(Form("CDCRLHits_%i.png", nevent));
163 nevent++;
164 }
165
166 private:
171 void fillAll()
172 {
173 AInBoxAlgorithm inBoxAlgorithm;
174
175 auto isLeaf = [this](typename Super::Node * node) {
176 return (node->getLevel() > this->getMaxLevel());
177 };
178
179 this->getTree()->fillWalk(inBoxAlgorithm, isLeaf);
180 }
181 };
182 }
184}
double R
typedef autogenerated by FFTW
A fast hough algorithm with rectangular boxes, which are split linearly by a fixed number of division...
typename HoughTree::Node Node
Type of the nodes used in the tree for the search.
HoughTree * getTree() const
Getter for the tree used in the search in the hough plane.
Class representing a three dimensional reconstructed hit.
Definition: CDCRecoHit3D.h:52
A convenience class for adding debug information to a Simple Hough Tree.
void writeDebugInfoToFile(const std::string &filename)
Write out some debug information to a ROOT file with the given name.
void drawDebugPlot(const std::vector< CDCRecoHit3D > &allHits, const std::vector< CDCRecoHit3D > &foundHits, const typename AInBoxAlgorithm::HoughBox &node)
Draw the results to a ROOT TCanvas.
void walk(AWalker &walker)
Forward walk to the top node.
Definition: DynTree.h:316
Convenience class for the typical usage-case: A box divisioned hough tree with maximum and minimum va...
float getMaximumY() const
Return the maximum value in y direction.
float getMaximumX() const
Return the maximum value in x direction.
A three dimensional vector.
Definition: Vector3D.h:33
double x() const
Getter for the x coordinate.
Definition: Vector3D.h:472
double y() const
Getter for the y coordinate.
Definition: Vector3D.h:484
double z() const
Getter for the z coordinate.
Definition: Vector3D.h:496
void fillWalk(AItemInDomainMeasure &weightItemInDomain, AIsLeafPredicate &isLeaf)
Walk through the children and fill them if necessary until isLeaf returns true.
Abstract base class for different kinds of events.