Belle II Software development
DebugableSimpleBoxDivisionHoughTree< AHitPtr, AInBoxAlgorithm, divisionX, divisionY > Class Template Reference

A convenience class for adding debug information to a Simple Hough Tree. More...

#include <DebugableSimpleBoxDivisionHoughTree.h>

Inheritance diagram for DebugableSimpleBoxDivisionHoughTree< AHitPtr, AInBoxAlgorithm, divisionX, divisionY >:
SimpleBoxDivisionHoughTree< AHitPtr, AInBoxAlgorithm, divisionX, divisionY > BoxDivisionHoughTree< AHitPtr, AInBoxAlgorithm::HoughBox, divisionX, divisionY > QuadraticLegendre< AHitPointerType, AHitDecisionAlgorithm, pDivisions, qDivisions > Z0TanLambdaLegendre< AHitPointerType, AHitDecisionAlgorithm, z0Divisions, tanLambdaDivisions >

Public Types

using BoxDivision = SectoredLinearDivision< HoughBox, divisions... >
 Type of the box division strategy.
 
using HoughTree = WeightedFastHoughTree< AHitPtr, HoughBox, BoxDivision >
 Type of the fast hough tree structure.
 
using Type = typename HoughBox::template Type< I >
 Type of the coordinate I.
 
using HasType = typename HoughBox::template HasType< T >
 Predicate that the given type is indeed a coordinate of the hough space.
 
using TypeIndex = typename HoughBox::template TypeIndex< T >
 Function to get the coordinate index from its type.
 
using Node = typename HoughTree::Node
 Type of the nodes used in the tree for the search.
 

Public Member Functions

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 initialize ()
 Initialize the tree with the given values.
 
std::vector< std::pair< HoughBox, std::vector< AHitPtr > > > findSingleBest (const Weight &minWeight)
 Find only the leave with the highest weight = number of items.
 
float getMaximumX () const
 Return the maximum value in x direction.
 
float getMaximumY () const
 Return the maximum value in y direction.
 
Width< 0 > getOverlapX () const
 Return the overlap in x direction.
 
Width< 1 > getOverlapY () const
 Return the overlap in y direction.
 
size_t getDivision (size_t i) const
 Getter the number of divisions at each level for coordinate index I.
 
void constructArray (double lowerBound, double upperBound, Width< I > nBinOverlap=0, Width< I > nBinWidth=0)
 Construct the discrete value array at coordinate index I.
 
void assignArray (Array< I > array, Width< I > overlap=0)
 Provide an externally constructed array by coordinate index.
 
std::enable_if_t< HasType< T >::value, void > assignArray (Array< TypeIndex< T >::value > array, Width< TypeIndex< T >::value > overlap=0)
 Provide an externally constructed array by coordinate type.
 
void seed (const AItemPtrs &items)
 Prepare the leave finding by filling the top node with given hits.
 
void fell ()
 Terminates the processing by striping all hit information from the tree.
 
void raze ()
 Release all memory that the tree aquired during the runs.
 
HoughTreegetTree () const
 Getter for the tree used in the search in the hough plane.
 
int getMaxLevel () const
 Getter for the currently set maximal level.
 
void setMaxLevel (int maxLevel)
 Setter maximal level of the hough tree.
 
int getSectorLevelSkip () const
 Getter for number of levels to skip in first level to form a finer sectored hough space.
 
void setSectorLevelSkip (int sectorLevelSkip)
 Setter for number of levels to skip in first level to form a finer sectored hough space.
 
const Array< I > & getArray () const
 Getter for the array of discrete value for coordinate I.
 

Private Types

using Super = SimpleBoxDivisionHoughTree< AHitPtr, AInBoxAlgorithm, divisionX, divisionY >
 The Super class.
 
using HoughBox = typename AInBoxAlgorithm::HoughBox
 The HoughBox we use.
 
template<size_t I>
using Width = typename HoughBox::template Width< I >
 Type of the width in coordinate I.
 
using Array = typename Type< I >::Array
 Type of the discrete value array to coordinate index I.
 
using Arrays = TupleGenerateN< Array, sizeof...(divisions)>
 Tuple type of the discrete value arrays.
 

Private Member Functions

void fillAll ()
 Fill the tree till all nodes are touched once.
 
HoughBox constructHoughPlaneImpl (const std::index_sequence< Is... > &is)
 Construct the box of the top node of the tree. Implementation unroling the indices.
 
HoughBox constructHoughPlane ()
 Construct the box of the top node of the tree.
 

Private Attributes

float m_maximumX = 0
 The maximum value in X direction.
 
float m_maximumY = 0
 The maximum value in y direction.
 
Width< 0 > m_overlapX = 0
 The overlap in X direction.
 
Width< 1 > m_overlapY = 0
 The overlap in Y direction.
 
int m_maxLevel
 Number of the maximum tree level.
 
int m_sectorLevelSkip
 Number of levels to skip in first level to form a finer sectored hough space.
 
const std::array< size_t, sizeof ...(divisions)> m_divisions
 Array of the number of divisions at each level.
 
HoughBox::Delta m_overlaps
 An tuple of division overlaps in each coordinate.
 
Arrays m_arrays
 A tuple of value arrays providing the memory for the discrete bin bounds.
 
std::unique_ptr< HoughTreem_houghTree
 Dynamic hough tree structure traversed in the leaf search.
 

Detailed Description

template<class AHitPtr, class AInBoxAlgorithm, size_t divisionX, size_t divisionY>
class Belle2::TrackFindingCDC::DebugableSimpleBoxDivisionHoughTree< AHitPtr, AInBoxAlgorithm, divisionX, divisionY >

A convenience class for adding debug information to a Simple Hough Tree.

Definition at line 27 of file DebugableSimpleBoxDivisionHoughTree.h.

Member Typedef Documentation

◆ Array

using Array = typename Type<I>::Array
privateinherited

Type of the discrete value array to coordinate index I.

Definition at line 77 of file BoxDivisionHoughTree.h.

◆ Arrays

using Arrays = TupleGenerateN<Array, sizeof...(divisions)>
privateinherited

Tuple type of the discrete value arrays.

Definition at line 80 of file BoxDivisionHoughTree.h.

◆ BoxDivision

using BoxDivision = SectoredLinearDivision<HoughBox, divisions...>
inherited

Type of the box division strategy.

Definition at line 41 of file BoxDivisionHoughTree.h.

◆ HasType

using HasType = typename HoughBox::template HasType<T>
inherited

Predicate that the given type is indeed a coordinate of the hough space.

Definition at line 52 of file BoxDivisionHoughTree.h.

◆ HoughBox

using HoughBox = typename AInBoxAlgorithm::HoughBox
privateinherited

The HoughBox we use.

Definition at line 27 of file SimpleBoxDivisionHoughTree.h.

◆ HoughTree

using HoughTree = WeightedFastHoughTree<AHitPtr , HoughBox, BoxDivision>
inherited

Type of the fast hough tree structure.

Definition at line 44 of file BoxDivisionHoughTree.h.

◆ Node

using Node = typename HoughTree::Node
inherited

Type of the nodes used in the tree for the search.

Definition at line 63 of file BoxDivisionHoughTree.h.

◆ Super

using Super = SimpleBoxDivisionHoughTree<AHitPtr, AInBoxAlgorithm, divisionX, divisionY>
private

The Super class.

Definition at line 31 of file DebugableSimpleBoxDivisionHoughTree.h.

◆ Type

using Type = typename HoughBox::template Type<I>
inherited

Type of the coordinate I.

Definition at line 48 of file BoxDivisionHoughTree.h.

◆ TypeIndex

using TypeIndex = typename HoughBox::template TypeIndex<T>
inherited

Function to get the coordinate index from its type.

Definition at line 56 of file BoxDivisionHoughTree.h.

◆ Width

using Width = typename HoughBox::template Width<I>
privateinherited

Type of the width in coordinate I.

Definition at line 31 of file SimpleBoxDivisionHoughTree.h.

Member Function Documentation

◆ assignArray() [1/2]

void assignArray ( Array< I >  array,
Width< I >  overlap = 0 
)
inlineinherited

Provide an externally constructed array by coordinate index.

Definition at line 128 of file BoxDivisionHoughTree.h.

129 {
130 std::get<I>(m_arrays) = std::move(array);
131 std::get<I>(m_overlaps) = overlap;
132
133 // In case of a discrete axes, check whether the size of the array is sensible
134 // such that the bin width at the highest granularity level is a whole number given the overlap.
135 if (std::is_integral<Width<I> >::value) {
136 const int division = getDivision(I);
137 const int granularityLevel = m_maxLevel + m_sectorLevelSkip;
138 const long int nBins = std::pow(division, granularityLevel);
139 const long int nPositions = std::get<I>(m_arrays).size();
140 const long int nWidthTimeNBins = nPositions - 1 + (nBins - 1) * overlap;
141
142 B2ASSERT("Not enough positions in array to cover all bins.\n"
143 "Expected: positions = " << nBins - (nBins - 1) * overlap + 1 << " at least.\n"
144 "Actual: positions = " << nPositions << " (overlap = " << overlap << ", bins = " << nBins << ")\n",
145 nWidthTimeNBins >= nBins);
146
147 B2ASSERT("Number of positions in array introduces inhomogeneous width at the highest granularity level.\n"
148 "Expected: positions = 'width * bins - (bins - 1) * overlap + 1'\n"
149 "Actual: positions = " << nPositions << " (overlap = " << overlap << ", bins = " << nBins << ")\n",
150 nWidthTimeNBins % nBins == 0);
151 }
152 }
int m_sectorLevelSkip
Number of levels to skip in first level to form a finer sectored hough space.
size_t getDivision(size_t i) const
Getter the number of divisions at each level for coordinate index I.

◆ assignArray() [2/2]

std::enable_if_t< HasType< T >::value, void > assignArray ( Array< TypeIndex< T >::value >  array,
Width< TypeIndex< T >::value >  overlap = 0 
)
inlineinherited

Provide an externally constructed array by coordinate type.

Definition at line 157 of file BoxDivisionHoughTree.h.

158 {
159 assignArray<TypeIndex<T>::value>(std::move(array), overlap);
160 }

◆ constructArray()

void constructArray ( double  lowerBound,
double  upperBound,
Width< I >  nBinOverlap = 0,
Width< I >  nBinWidth = 0 
)
inlineinherited

Construct the discrete value array at coordinate index I.

This function is only applicable for discrete axes. For continuous axes assignArray should be call with an array containing only the lower and upper bound of the axes range and an optional overlap.

Parameters
lowerBoundLower bound of the value range
upperBoundUpper bound of the value range
nBinOverlapOverlap of neighboring bins. Default is no overlap. Usuallly this is counted in number of discrete values
nBinWidthWidth of the bins at lowest level. Default is width of 1. Usually this is counted in numbers of discrete values

Definition at line 104 of file BoxDivisionHoughTree.h.

108 {
109 static_assert(std::is_integral<Width<I> >::value, "Method only applicable to discrete axes");
110 const size_t division = getDivision(I);
111 const int granularityLevel = m_maxLevel + m_sectorLevelSkip;
112 const size_t nBins = std::pow(division, granularityLevel);
113
114 if (nBinWidth == 0) {
115 nBinWidth = nBinOverlap + 1;
116 }
117
118 B2ASSERT("Width " << nBinWidth << "is not bigger than overlap " << nBinOverlap,
119 nBinOverlap < nBinWidth);
120
121 const auto nPositions = (nBinWidth - nBinOverlap) * nBins + nBinOverlap + 1;
122 std::get<I>(m_arrays) = linspace<float>(lowerBound, upperBound, nPositions);
123 std::get<I>(m_overlaps) = nBinOverlap;
124 }

◆ constructHoughPlane()

HoughBox constructHoughPlane ( )
inlineprivateinherited

Construct the box of the top node of the tree.

Definition at line 239 of file BoxDivisionHoughTree.h.

240 {
241 return constructHoughPlaneImpl(std::make_index_sequence<sizeof...(divisions)>());
242 }
HoughBox constructHoughPlaneImpl(const std::index_sequence< Is... > &is)
Construct the box of the top node of the tree. Implementation unroling the indices.

◆ constructHoughPlaneImpl()

HoughBox constructHoughPlaneImpl ( const std::index_sequence< Is... > &  is)
inlineprivateinherited

Construct the box of the top node of the tree. Implementation unroling the indices.

Definition at line 233 of file BoxDivisionHoughTree.h.

234 {
235 return HoughBox(Type<Is>::getRange(std::get<Is>(m_arrays))...);
236 }

◆ drawDebugPlot()

void drawDebugPlot ( const std::vector< CDCRecoHit3D > &  allHits,
const std::vector< CDCRecoHit3D > &  foundHits,
const typename AInBoxAlgorithm::HoughBox &  node 
)
inline

Draw the results to a ROOT TCanvas.

Definition at line 98 of file DebugableSimpleBoxDivisionHoughTree.h.

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 }
double R
typedef autogenerated by FFTW
HepGeom::Vector3D< double > Vector3D
3D Vector
Definition: Cell.h:34

◆ fell()

void fell ( )
inlineinherited

Terminates the processing by striping all hit information from the tree.

Definition at line 181 of file BoxDivisionHoughTree.h.

182 {
183 m_houghTree->fell();
184 }
std::unique_ptr< HoughTree > m_houghTree
Dynamic hough tree structure traversed in the leaf search.

◆ fillAll()

void fillAll ( )
inlineprivate

Fill the tree till all nodes are touched once.

This is not for finding results, but for debug reasons.

Definition at line 171 of file DebugableSimpleBoxDivisionHoughTree.h.

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 }
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.
void fillWalk(AItemInDomainMeasure &weightItemInDomain, AIsLeafPredicate &isLeaf)
Walk through the children and fill them if necessary until isLeaf returns true.

◆ findSingleBest()

std::vector< std::pair< HoughBox, std::vector< AHitPtr > > > findSingleBest ( const Weight &  minWeight)
inlineinherited

Find only the leave with the highest weight = number of items.

Definition at line 58 of file SimpleBoxDivisionHoughTree.h.

59 {
60 AInBoxAlgorithm inBoxAlgorithm;
61 auto skipLowWeightNode = [minWeight](const typename Super::Node * node) {
62 return not(node->getWeight() >= minWeight);
63 };
64 auto found = this->getTree()->findHeaviestLeafSingle(inBoxAlgorithm, this->getMaxLevel(), skipLowWeightNode);
65
66 std::vector<std::pair<HoughBox, std::vector<AHitPtr>>> result;
67 if (found) {
68 // Move the found content over. unique_ptr still destroys the left overs.
69 result.push_back(std::move(*found));
70 }
71 return result;
72 }
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.

◆ getArray()

const Array< I > & getArray ( ) const
inlineinherited

Getter for the array of discrete value for coordinate I.

Definition at line 225 of file BoxDivisionHoughTree.h.

226 {
227 return std::get<I>(m_arrays);
228 }

◆ getDivision()

size_t getDivision ( size_t  i) const
inlineinherited

Getter the number of divisions at each level for coordinate index I.

Definition at line 84 of file BoxDivisionHoughTree.h.

85 {
86 return m_divisions[i];
87 }
const std::array< size_t, sizeof ...(divisions)> m_divisions
Array of the number of divisions at each level.

◆ getMaximumX()

float getMaximumX ( ) const
inlineinherited

Return the maximum value in x direction.

Definition at line 75 of file SimpleBoxDivisionHoughTree.h.

76 {
77 return m_maximumX;
78 }

◆ getMaximumY()

float getMaximumY ( ) const
inlineinherited

Return the maximum value in y direction.

Definition at line 81 of file SimpleBoxDivisionHoughTree.h.

82 {
83 return m_maximumY;
84 }

◆ getMaxLevel()

int getMaxLevel ( ) const
inlineinherited

Getter for the currently set maximal level.

Definition at line 200 of file BoxDivisionHoughTree.h.

201 {
202 return m_maxLevel;
203 }

◆ getOverlapX()

Width< 0 > getOverlapX ( ) const
inlineinherited

Return the overlap in x direction.

Definition at line 87 of file SimpleBoxDivisionHoughTree.h.

88 {
89 return m_overlapX;
90 }

◆ getOverlapY()

Width< 1 > getOverlapY ( ) const
inlineinherited

Return the overlap in y direction.

Definition at line 93 of file SimpleBoxDivisionHoughTree.h.

94 {
95 return m_overlapY;
96 }

◆ getSectorLevelSkip()

int getSectorLevelSkip ( ) const
inlineinherited

Getter for number of levels to skip in first level to form a finer sectored hough space.

Definition at line 212 of file BoxDivisionHoughTree.h.

213 {
214 return m_sectorLevelSkip;
215 }

◆ getTree()

HoughTree * getTree ( ) const
inlineinherited

Getter for the tree used in the search in the hough plane.

Definition at line 194 of file BoxDivisionHoughTree.h.

195 {
196 return m_houghTree.get();
197 }

◆ initialize()

void initialize ( )
inlineinherited

Initialize the tree with the given values.

Definition at line 48 of file SimpleBoxDivisionHoughTree.h.

49 {
50 Super::template constructArray<0>(-getMaximumX(), getMaximumX(), getOverlapX());
51 Super::template constructArray<1>(-getMaximumY(), getMaximumY(), getOverlapY());
52
54 }
void initialize()
Initialise the algorithm by constructing the hough tree from the parameters.
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.
Width< 1 > getOverlapY() const
Return the overlap in y direction.

◆ raze()

void raze ( )
inlineinherited

Release all memory that the tree aquired during the runs.

Definition at line 187 of file BoxDivisionHoughTree.h.

188 {
189 m_houghTree->raze();
190 }

◆ seed()

void seed ( const AItemPtrs &  items)
inlineinherited

Prepare the leave finding by filling the top node with given hits.

Definition at line 174 of file BoxDivisionHoughTree.h.

175 {
176 if (not m_houghTree) initialize();
177 m_houghTree->seed(items);
178 }

◆ setMaxLevel()

void setMaxLevel ( int  maxLevel)
inlineinherited

Setter maximal level of the hough tree.

Definition at line 206 of file BoxDivisionHoughTree.h.

207 {
208 m_maxLevel = maxLevel;
209 }

◆ setSectorLevelSkip()

void setSectorLevelSkip ( int  sectorLevelSkip)
inlineinherited

Setter for number of levels to skip in first level to form a finer sectored hough space.

Definition at line 218 of file BoxDivisionHoughTree.h.

219 {
220 m_sectorLevelSkip = sectorLevelSkip;
221 }

◆ writeDebugInfoToFile()

void writeDebugInfoToFile ( const std::string &  filename)
inline

Write out some debug information to a ROOT file with the given name.

This must be done before felling the tree. Attention: This will delete a ROOT file with the same name if already present.

Parameters
filenameThe ROOT filename.

Definition at line 42 of file DebugableSimpleBoxDivisionHoughTree.h.

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 }
void walk(AWalker &walker)
Forward walk to the top node.
Definition: DynTree.h:316

Member Data Documentation

◆ m_arrays

Arrays m_arrays
privateinherited

A tuple of value arrays providing the memory for the discrete bin bounds.

Definition at line 258 of file BoxDivisionHoughTree.h.

◆ m_divisions

const std::array<size_t, sizeof ...(divisions)> m_divisions
privateinherited

Array of the number of divisions at each level.

Definition at line 252 of file BoxDivisionHoughTree.h.

◆ m_houghTree

std::unique_ptr<HoughTree> m_houghTree
privateinherited

Dynamic hough tree structure traversed in the leaf search.

Definition at line 261 of file BoxDivisionHoughTree.h.

◆ m_maximumX

float m_maximumX = 0
privateinherited

The maximum value in X direction.

Definition at line 100 of file SimpleBoxDivisionHoughTree.h.

◆ m_maximumY

float m_maximumY = 0
privateinherited

The maximum value in y direction.

Definition at line 103 of file SimpleBoxDivisionHoughTree.h.

◆ m_maxLevel

int m_maxLevel
privateinherited

Number of the maximum tree level.

Definition at line 246 of file BoxDivisionHoughTree.h.

◆ m_overlaps

HoughBox::Delta m_overlaps
privateinherited

An tuple of division overlaps in each coordinate.

Definition at line 255 of file BoxDivisionHoughTree.h.

◆ m_overlapX

Width<0> m_overlapX = 0
privateinherited

The overlap in X direction.

Definition at line 106 of file SimpleBoxDivisionHoughTree.h.

◆ m_overlapY

Width<1> m_overlapY = 0
privateinherited

The overlap in Y direction.

Definition at line 109 of file SimpleBoxDivisionHoughTree.h.

◆ m_sectorLevelSkip

int m_sectorLevelSkip
privateinherited

Number of levels to skip in first level to form a finer sectored hough space.

Definition at line 249 of file BoxDivisionHoughTree.h.


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