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