Belle II Software development
SimpleBoxDivisionHoughTree3D.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/BoxDivisionHoughTree.h>
10#include <tracking/trackFindingCDC/eventdata/hits/CDCRecoHit3D.h>
11
12#include <TGraph.h>
13#include <TF1.h>
14#include <TCanvas.h>
15#include <TAxis.h>
16
17namespace Belle2 {
22 namespace TrackFindingCDC {
23
26 template<class AHitPtr, class AInBoxAlgorithm, size_t divisionX, size_t divisionY, size_t divisionZ>
28 BoxDivisionHoughTree<AHitPtr, typename AInBoxAlgorithm::HoughBox, divisionX, divisionY, divisionZ> {
29
30 private:
33
35 using HoughBox = typename AInBoxAlgorithm::HoughBox;
36
38 template <size_t I>
39 using Width = typename HoughBox::template Width<I>;
40
41 public:
44 float maximumY,
45 float maximumZ,
46 Width<0> overlapX = 0,
47 Width<1> overlapY = 0,
48 Width<2> overlapZ = 0)
49 : Super(0)
50 , m_maximumX(maximumX)
51 , m_maximumY(maximumY)
52 , m_maximumZ(maximumZ)
53 , m_overlapX(overlapX)
54 , m_overlapY(overlapY)
55 , m_overlapZ(overlapZ)
56 {
57 }
58
61 {
62 Super::template constructArray<0>(-getMaximumX(), getMaximumX(), getOverlapX());
63 Super::template constructArray<1>(-getMaximumY(), getMaximumY(), getOverlapY());
64 Super::template constructArray<2>(-getMaximumZ(), getMaximumZ(), getOverlapZ());
65
67 }
68
70 std::vector<std::pair<HoughBox, std::vector<AHitPtr>>>
71 findSingleBest(const Weight& minWeight)
72 {
73 AInBoxAlgorithm inBoxAlgorithm;
74 auto skipLowWeightNode = [minWeight](const typename Super::Node * node) {
75 return not(node->getWeight() >= minWeight);
76 };
77 auto found = this->getTree()->findHeaviestLeafSingle(inBoxAlgorithm, this->getMaxLevel(), skipLowWeightNode);
78
79 std::vector<std::pair<HoughBox, std::vector<AHitPtr>>> result;
80 if (found) {
81 // Move the found content over. unique_ptr still destroys the left overs.
82 result.push_back(std::move(*found));
83 }
84 return result;
85 }
86
88 void writeDebugInfoToFile(const std::string& filename __attribute__((unused)))
89 {
90 //do nothing;
91 }
92
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 centerX = (AInBoxAlgorithm::BoxAlgorithm::centerX(node));
133 const double deltaX = (AInBoxAlgorithm::BoxAlgorithm::deltaX(node));
134 const double centerY = (AInBoxAlgorithm::BoxAlgorithm::centerY(node));
135 const double centerZ = (AInBoxAlgorithm::BoxAlgorithm::centerZ(node));
136
137 TF1* candidateL = new TF1("candL", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
138 TF1* candidateH = new TF1("candH", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
139 TF1* candidateMean = new TF1("candMean", AInBoxAlgorithm::BoxAlgorithm::debugLine(), 0, 120);
140
141 candidateL->SetParameters(centerX - deltaX, centerY, centerZ - 100.0 * deltaX);
142 candidateH->SetParameters(centerX + deltaX, centerY, centerZ + 100.0 * deltaX);
143 candidateMean->SetParameters(centerX, centerY, centerZ);
144
145 candidateL->SetLineColor(9);
146 candidateH->SetLineColor(41);
147 candidateMean->SetLineColor(2);
148
149 candidateL->Draw("same");
150 candidateH->Draw("same");
151 candidateMean->Draw("same");
152 canv.SaveAs(Form("CDCRLHits_%i.png", nevent));
153 nevent++;
154 }
155
157 float getMaximumX() const
158 {
159 return m_maximumX;
160 }
161
163 float getMaximumY() const
164 {
165 return m_maximumY;
166 }
167
169 float getMaximumZ() const
170 {
171 return m_maximumZ;
172 }
173
176 {
177 return m_overlapX;
178 }
179
182 {
183 return m_overlapY;
184 }
185
188 {
189 return m_overlapZ;
190 }
191
192 private:
194 float m_maximumX = 0;
195
197 float m_maximumY = 0;
198
200 float m_maximumZ = 0;
201
204
207
210 };
211 }
213}
double R
typedef autogenerated by FFTW
A fast hough algorithm with rectangular boxes, which are split linearly by a fixed number of division...
void initialize()
Initialise the algorithm by constructing the hough tree from the parameters.
typename HoughTree::Node Node
Type of the nodes used in the tree for the search.
Class representing a three dimensional reconstructed hit.
Definition: CDCRecoHit3D.h:52
Convenience class for the typical usage-case: A box divisioned hough tree with maximum and minimum va...
SimpleBoxDivisionHoughTree3D(float maximumX, float maximumY, float maximumZ, Width< 0 > overlapX=0, Width< 1 > overlapY=0, Width< 2 > overlapZ=0)
Constructor using the given maximal level.
float getMaximumZ() const
Return the maximum value in Z direction.
void initialize()
Initialize the tree with the given values.
float getMaximumY() const
Return the maximum value in y direction.
Width< 0 > getOverlapX() const
Return the overlap in x direction.
float getMaximumX() const
Return the maximum value in x direction.
typename AInBoxAlgorithm::HoughBox HoughBox
The HoughBox we use.
void writeDebugInfoToFile(const std::string &filename)
Write debug information into a ROOT file; not implemented.
typename HoughBox::template Width< I > Width
Type of the width in coordinate I.
Width< 1 > getOverlapY() const
Return the overlap in y direction.
std::vector< std::pair< HoughBox, std::vector< AHitPtr > > > findSingleBest(const Weight &minWeight)
Find only the leaf with the highest weight (~= number of items)
void drawDebugPlot(const std::vector< CDCRecoHit3D > &allHits, const std::vector< CDCRecoHit3D > &foundHits, const typename AInBoxAlgorithm::HoughBox &node)
Draws found hits and node boundaries FIXME this is a copy-paste from DebugableSimpleBoxDivisionHoughT...
Width< 2 > getOverlapZ() const
Return the overlap in y 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
std::unique_ptr< std::pair< ADomain, std::vector< T > > > findHeaviestLeafSingle(AItemInDomainMeasure &weightItemInDomain, int maxLevel, ASkipNodePredicate &skipNode)
Go through all children until the maxLevel is reached and find the leaf with the highest weight.
Abstract base class for different kinds of events.