Belle II Software development
AxialHitQuadTreeProcessor Class Referenceabstract

A QuadTreeProcessor for TrackHits. More...

#include <AxialHitQuadTreeProcessor.h>

Inheritance diagram for AxialHitQuadTreeProcessor:
QuadTreeProcessor< long, float, const TrackingUtilities::CDCWireHit >

Public Types

using Item
 The QuadTree will only see items of this type.
 
using QuadTree
 The used QuadTree.
 
using XSpan
 This pair describes the span in X for a node.
 
using YSpan
 This pair describes the span in Y for a node.
 
using XYSpans
 This pair of spans describes the span of a node.
 
using QuadTreeChildren
 Alias for the QuadTree Children.
 
using CandidateReceiver
 This lambda function can be used for postprocessing.
 

Public Member Functions

 AxialHitQuadTreeProcessor (int lastLevel, int seedLevel, const XYSpans &ranges, PrecisionUtil::PrecisionFunction precisionFunction)
 Constructor.
 
 AxialHitQuadTreeProcessor (const TrackingUtilities::Vector2D &localOrigin, const YSpan &curvSpan, const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > *cosSinLookupTable)
 Constructor for the quad tree processor used in the off-origin extension.
 
void drawHits (std::vector< const TrackingUtilities::CDCWireHit * > hits, unsigned int color=46) const
 Draw QuadTree node.
 
void drawNode (QuadTree *node) const
 Draw QuadTree node.
 
void clear ()
 Delete all the QuadTreeItems in the tree and clear the tree.
 
void seed (const std::vector< const TrackingUtilities::CDCWireHit * > &datas)
 Fill in the items in the given vector.
 
std::vector< const TrackingUtilities::CDCWireHit * > getAssignedItems ()
 Get items that have been assigned to the seed level The returned elements are unique even if items are assigned multiple times.
 
void fill (const CandidateReceiver &candidateReceiver, int nHitsThreshold)
 Start filling the already created tree.
 
void fill (const CandidateReceiver &candidateReceiver, int nHitsThreshold, float yLimit)
 Fill vector of QuadTree instances with hits.
 
virtual void afterFillDebugHook (QuadTreeChildren &children)
 Override that function if you want to receive debug output whenever the children of a node are filled the first time Maybe you want to make some nice plots or statistics.
 
const std::map< std::pair< long, float >, std::vector< Item * > > & getDebugInformation () const
 Return the debug information if collected.
 

Static Public Member Functions

static const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > & getCosSinLookupTable ()
 Get the standard lookup table containing equally spaces unit vectors (cos, sin)
 
static std::vector< float > createCurvBound (YSpan curvSpan, int lastLevel)
 Constructs an array with the curvature bounds as generated by the default bin divisions.
 

Protected Member Functions

bool isLeaf (QuadTree *node) const final
 lastLevel depends on curvature of the track candidate
 
XYSpans createChild (QuadTree *node, int i, int j) const final
 Return the new ranges.
 
bool isInNode (QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const final
 Check whether hit belongs to the quadtree node:
 
bool checkDerivative (QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const
 Check derivative of the sinogram.
 
bool checkExtremum (QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const
 Checks whether extreme point is located within QuadTree node's ranges.
 
virtual XYSpans createChild (QuadTree *node, int iX, int iY) const
 Implement that function if you want to provide a new processor.
 
virtual bool isInNode (QuadTree *node, const TrackingUtilities::CDCWireHit *item) const=0
 Implement that function if you want to provide a new processor.
 
virtual bool isLeaf (QuadTree *node) const
 Function which checks if given node is leaf Implemented as virtual to keep possibility of changing lastLevel values depending on region is phase-space (i.e.
 
int getLastLevel () const
 Return the parameter last level.
 

Protected Attributes

std::unique_ptr< QuadTreem_quadTree
 The quad tree we work with.
 
std::deque< Itemm_items
 Storage space for the items that are referenced by the quad tree nodes.
 
std::vector< QuadTree * > m_seededTrees
 Vector of QuadTrees QuadTree instances (which are filled in the vector) cover the whole Legendre phase-space; each instance is processes independently.
 

Private Member Functions

void fillGivenTree (QuadTree *node, const CandidateReceiver &candidateReceiver, int nItemsThreshold, float yLimit)
 Internal function to do the real quad tree search: fill the nodes, check which of the n*m bins we need to process further and go one level deeper.
 
void createChildren (QuadTree *node, QuadTreeChildren &m_children) const
 Creates the sub node of a given node.
 
void fillChildren (QuadTree *node, const std::vector< Item * > &items)
 This function is called by fillGivenTree and fills the items into the corresponding children.
 
void callResultFunction (QuadTree *node, const CandidateReceiver &candidateReceiver) const
 When a node is accepted as a result, we extract a vector with the items (back transformed to AData*) and pass it together with the result node to the candidate receiver function.
 

Private Attributes

PrecisionUtil::PrecisionFunction m_precisionFunction
 Lambda which holds resolution function for the quadtree.
 
TrackingUtilities::Vector2D m_localOrigin
 Local origin on which the phase space coordinates are centered.
 
const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > * m_cosSinLookupTable
 Pinned lookup table for precomputed cosine and sine values.
 
const double c_curlCurv = 0.02
 The curvature above which the trajectory is considered a curler.
 
bool m_twoSidedPhaseSpace
 Indicator whether the two sided phases space insertion check should be used This option should automatically split back to back tracks in the low curvature regions.
 
int m_lastLevel
 The last level to be filled.
 
int m_seedLevel
 The first level to be filled, effectively skip forward to this higher granularity level.
 
bool m_debugOutput
 A flag to control the creation of the debug output.
 
std::map< std::pair< long, float >, std::vector< Item * > > m_debugOutputMap
 The calculated debug map.
 

Detailed Description

A QuadTreeProcessor for TrackHits.

Definition at line 31 of file AxialHitQuadTreeProcessor.h.

Member Typedef Documentation

◆ CandidateReceiver

using CandidateReceiver
inherited

This lambda function can be used for postprocessing.

Definition at line 62 of file QuadTreeProcessor.h.

◆ Item

using Item
inherited

The QuadTree will only see items of this type.

Definition at line 44 of file QuadTreeProcessor.h.

◆ QuadTree

using QuadTree
inherited

The used QuadTree.

Definition at line 47 of file QuadTreeProcessor.h.

◆ QuadTreeChildren

using QuadTreeChildren
inherited

Alias for the QuadTree Children.

Definition at line 59 of file QuadTreeProcessor.h.

◆ XSpan

using XSpan
inherited

This pair describes the span in X for a node.

Definition at line 50 of file QuadTreeProcessor.h.

◆ XYSpans

using XYSpans
inherited

This pair of spans describes the span of a node.

Definition at line 56 of file QuadTreeProcessor.h.

◆ YSpan

using YSpan
inherited

This pair describes the span in Y for a node.

Definition at line 53 of file QuadTreeProcessor.h.

Constructor & Destructor Documentation

◆ AxialHitQuadTreeProcessor() [1/2]

AxialHitQuadTreeProcessor ( int lastLevel,
int seedLevel,
const XYSpans & ranges,
PrecisionUtil::PrecisionFunction precisionFunction )

Constructor.

Definition at line 97 of file AxialHitQuadTreeProcessor.cc.

101 : QuadTreeProcessor(lastLevel, seedLevel, ranges)
102 , m_precisionFunction(precisionFunction)
103 , m_localOrigin(0.0, 0.0)
105{
106 m_twoSidedPhaseSpace = m_quadTree->getYMin() * m_quadTree->getYMax() < 0;
107}
bool m_twoSidedPhaseSpace
Indicator whether the two sided phases space insertion check should be used This option should automa...
TrackingUtilities::Vector2D m_localOrigin
Local origin on which the phase space coordinates are centered.
PrecisionUtil::PrecisionFunction m_precisionFunction
Lambda which holds resolution function for the quadtree.
static const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > & getCosSinLookupTable()
Get the standard lookup table containing equally spaces unit vectors (cos, sin)
const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > * m_cosSinLookupTable
Pinned lookup table for precomputed cosine and sine values.
QuadTreeProcessor(int lastLevel, int seedLevel, const XYSpans &xySpans, bool debugOutput=false)

◆ AxialHitQuadTreeProcessor() [2/2]

AxialHitQuadTreeProcessor ( const TrackingUtilities::Vector2D & localOrigin,
const YSpan & curvSpan,
const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > * cosSinLookupTable )

Constructor for the quad tree processor used in the off-origin extension.

Currently only used in zero level mode to collect hits that are in a phase space part with respect to the given point.

Definition at line 109 of file AxialHitQuadTreeProcessor.cc.

112 : QuadTreeProcessor(0, 0, { {0, cosSinLookupTable->getNPoints() - 1}, curvSpan})
113, m_localOrigin(localOrigin)
114, m_cosSinLookupTable(cosSinLookupTable)
115{
116 // Never use two sided mode in off origin extension
117 m_twoSidedPhaseSpace = false;
118}

Member Function Documentation

◆ afterFillDebugHook()

virtual void afterFillDebugHook ( QuadTreeChildren & children)
inlinevirtualinherited

Override that function if you want to receive debug output whenever the children of a node are filled the first time Maybe you want to make some nice plots or statistics.

Definition at line 359 of file QuadTreeProcessor.h.

360 {
361 if (not m_debugOutput) return;
362 for (QuadTree& childNode : children) {
363 if (childNode.getLevel() != getLastLevel()) continue; // Only write the lowest level
364 //m_debugOutputMap[ {childNode.getXMean(), childNode.getYMean()}] = childNode.getItems();
365 }
366 }

◆ callResultFunction()

void callResultFunction ( QuadTree * node,
const CandidateReceiver & candidateReceiver ) const
inlineprivateinherited

When a node is accepted as a result, we extract a vector with the items (back transformed to AData*) and pass it together with the result node to the candidate receiver function.

Definition at line 291 of file QuadTreeProcessor.h.

292 {
293 const std::vector<Item*>& foundItems = node->getItems();
294 std::vector<AData*> candidate;
295 candidate.reserve(foundItems.size());
296
297 for (Item* item : foundItems) {
298 item->setUsedFlag();
299 candidate.push_back(item->getPointer());
300 }
301
302 candidateReceiver(candidate, node);
303 }

◆ checkDerivative()

bool checkDerivative ( QuadTree * node,
const TrackingUtilities::CDCWireHit * wireHit ) const
protected

Check derivative of the sinogram.

Parameters
nodeQuadTree node
wireHitpointer to the hit to check
Returns
returns true in cases:
- positive derivative and no extremum in the node's ranges or
- extremum located in the node's ranges
returns false in other cases (namely negative derivative

Definition at line 262 of file AxialHitQuadTreeProcessor.cc.

263{
264 const Vector2D& pos2D = wireHit->getRefPos2D() - m_localOrigin;
265
266 long thetaMin = node->getXMin();
267 long thetaMax = node->getXMax();
268
269 const Vector2D& thetaVecMin = m_cosSinLookupTable->at(thetaMin);
270 const Vector2D& thetaVecMax = m_cosSinLookupTable->at(thetaMax);
271
272 float rMinD = thetaVecMin.cross(pos2D);
273 float rMaxD = thetaVecMax.cross(pos2D);
274
275 // Does not really make sense...
276 if ((rMinD > 0) && (rMaxD * rMinD >= 0)) return true;
277 if ((rMaxD * rMinD < 0)) return true;
278 return false;
279}
double cross(const Vector2D &rhs) const
Calculated the two dimensional cross product.
Definition Vector2D.h:189

◆ checkExtremum()

bool checkExtremum ( QuadTree * node,
const TrackingUtilities::CDCWireHit * wireHit ) const
protected

Checks whether extreme point is located within QuadTree node's ranges.

Parameters
nodeQuadTree node
wireHithit to check
Returns
true or false

Definition at line 281 of file AxialHitQuadTreeProcessor.cc.

282{
283 const double& l = wireHit->getRefDriftLength();
284 const Vector2D& pos2D = wireHit->getRefPos2D() - m_localOrigin;
285 double r2 = pos2D.normSquared() - l * l;
286
287 // get left and right borders of the node
288 long thetaMin = node->getXMin();
289 long thetaMax = node->getXMax();
290
291 const Vector2D& thetaVecMin = m_cosSinLookupTable->at(thetaMin);
292 const Vector2D& thetaVecMax = m_cosSinLookupTable->at(thetaMax);
293
294 if (not pos2D.isBetween(thetaVecMin, thetaVecMax)) return false;
295
296 // compute sinograms at the position
297 double r = pos2D.norm();
298 float rRight = r - l;
299 float rLeft = r + l;
300
301 // get top and bottom borders of the node
302 float rMin = node->getYMin() * r2 / 2;
303 float rMax = node->getYMax() * r2 / 2;
304
305 bool crossesRight = (rMin - rRight) * (rMax - rRight) < 0;
306 bool crossesLeft = (rMin - rLeft) * (rMax - rLeft) < 0;
307 return crossesRight or crossesLeft;
308}
bool isBetween(const Vector2D &lower, const Vector2D &upper) const
Checks if this vector is between two other vectors Between means here that when rotating the lower ve...
Definition Vector2D.h:556
double normSquared() const
Calculates .
Definition Vector2D.h:195
double norm() const
Calculates the length of the vector.
Definition Vector2D.h:206

◆ clear()

void clear ( )
inlineinherited

Delete all the QuadTreeItems in the tree and clear the tree.

Definition at line 95 of file QuadTreeProcessor.h.

96 {
97 m_seededTrees.clear();
98 m_quadTree->clearChildren();
99 m_quadTree->clearItems();
100 m_items.clear();
101 }

◆ createChild() [1/2]

AxialHitQuadTreeProcessor::XYSpans createChild ( QuadTree * node,
int i,
int j ) const
finalprotected

Return the new ranges.

We do not use the standard ranges for the lower levels.

Parameters
nodequadtree node
itheta index of the child
jrho index of the child
Returns
returns ranges of the (i;j) child

Definition at line 136 of file AxialHitQuadTreeProcessor.cc.

137{
138 const int nodeLevel = node->getLevel();
139 const int lastLevel = getLastLevel();
140 const float meanCurv = std::fabs(node->getYMax() + node->getYMin()) / 2;
141
142 // Expand bins for all nodes 7 levels before the last level (for lastLevel = 12 starting at 6)
143 // but only in a curvature region higher than 0.005. Lower than that use always standard.
144 bool standardBinning = (nodeLevel <= lastLevel - 7) or (meanCurv <= 0.005);
145
146 if (standardBinning) {
147 float r1 = node->getYLowerBound(j);
148 float r2 = node->getYUpperBound(j);
149 long theta1 = node->getXLowerBound(i);
150 long theta2 = node->getXUpperBound(i);
151
152 // Standard bin division
153 return XYSpans({theta1, theta2}, {r1, r2});
154 }
155
156 // Non-standard binning
157 // For level 6 to 7 only expand 1 / 4, for higher levels expand 1 / 8.
158 // (assuming last level == 12)
159 if (nodeLevel < lastLevel - 5) {
160 float r1 = node->getYLowerBound(j) - node->getYBinWidth(j) / 4.;
161 float r2 = node->getYUpperBound(j) + node->getYBinWidth(j) / 4.;
162
163 // long extension = pow(2, lastLevel - nodeLevel) / 4; is same as:
164 long extension = pow(2, lastLevel - nodeLevel - 2);
165
166 long theta1 = node->getXLowerBound(i) - extension;
167 if (theta1 < 0) theta1 = 0;
168
169 long theta2 = node->getXUpperBound(i) + extension;
170 if (theta2 >= m_cosSinLookupTable->getNPoints()) {
171 theta2 = m_cosSinLookupTable->getNPoints() - 1;
172 }
173
174 return XYSpans({theta1, theta2}, {r1, r2});
175 } else {
176 float r1 = node->getYLowerBound(j) - node->getYBinWidth(j) / 8.;
177 float r2 = node->getYUpperBound(j) + node->getYBinWidth(j) / 8.;
178
179 // long extension = pow(2, lastLevel - nodeLevel) / 8; is same as
180 long extension = pow(2, lastLevel - nodeLevel - 3);
181
182 long theta1 = node->getXLowerBound(i) - extension;
183 if (theta1 < 0) theta1 = 0;
184
185 long theta2 = node->getXUpperBound(i) + extension;
186 if (theta2 >= m_cosSinLookupTable->getNPoints()) {
187 theta2 = m_cosSinLookupTable->getNPoints() - 1;
188 }
189
190 return XYSpans({theta1, theta2}, {r1, r2});
191 }
192}

◆ createChild() [2/2]

virtual XYSpans createChild ( QuadTree * node,
int iX,
int iY ) const
inlineprotectedvirtualinherited

Implement that function if you want to provide a new processor.

It decides which node-spans the n * m children of the node should have. It is called when creating the nodes. The two indices iX and iY tell you where the new node will be created (as node.children[iX][iY]). You can check some information on the level or the x- or y-values by using the methods implemented for node.

Returns
a XYSpan pair of a x- and a y-span that the new child should have. If you don nt want to provide custom spans, just return XYSpans(XSpan(node->getXBinBound(iX), node->getXBinBound(iX + 1)), YSpan(node->getYBinBound(iY), node->getYBinBound(iY + 1)));

Definition at line 314 of file QuadTreeProcessor.h.

315 {
316 AX xMin = node->getXLowerBound(iX);
317 AX xMax = node->getXUpperBound(iX);
318 AY yMin = node->getYLowerBound(iY);
319 AY yMax = node->getYUpperBound(iY);
320 return XYSpans({xMin, xMax}, {yMin, yMax});
321 }

◆ createChildren()

void createChildren ( QuadTree * node,
QuadTreeChildren & m_children ) const
inlineprivateinherited

Creates the sub node of a given node.

This function is called by fillGivenTree. To calculate the spans of the children nodes the user-defined function createChiildWithParent is used.

Definition at line 252 of file QuadTreeProcessor.h.

253 {
254 for (int i = 0; i < node->getXNbins(); ++i) {
255 for (int j = 0; j < node->getYNbins(); ++j) {
256 const XYSpans& xySpans = createChild(node, i, j);
257 const XSpan& xSpan = xySpans.first;
258 const YSpan& ySpan = xySpans.second;
259 m_children.push_back(QuadTree(xSpan, ySpan, node->getLevel() + 1, node));
260 }
261 }
262 }

◆ createCurvBound()

std::vector< float > createCurvBound ( YSpan curvSpan,
int lastLevel )
static

Constructs an array with the curvature bounds as generated by the default bin divisions.

Definition at line 65 of file AxialHitQuadTreeProcessor.cc.

66{
67 std::vector<YSpan> spans{{curvSpan}};
68
69 std::vector<YSpan> nextSpans;
70 for (int level = 1; level <= lastLevel; ++level) {
71 nextSpans.clear();
72 for (const YSpan& span : spans) {
73 nextSpans.push_back(splitCurvSpan(span, level, lastLevel, 0));
74 nextSpans.push_back(splitCurvSpan(span, level, lastLevel, 1));
75 }
76 spans.swap(nextSpans);
77 }
78
79 std::vector<float> bounds;
80 for (const YSpan& span : spans) {
81 bounds.push_back(span[0]);
82 bounds.push_back(span[1]);
83 }
84
85 assert(bounds.size() == std::pow(2, lastLevel));
86 return bounds;
87}

◆ drawHits()

void drawHits ( std::vector< const TrackingUtilities::CDCWireHit * > hits,
unsigned int color = 46 ) const

Draw QuadTree node.

Definition at line 310 of file AxialHitQuadTreeProcessor.cc.

311{
312 static int nevent(0);
313
314 TCanvas* canv = new TCanvas("canv", "legendre transform", 0, 0, 1200, 600);
315 canv->cd(1);
316 TGraph* dummyGraph = new TGraph();
317 dummyGraph->SetPoint(1, -M_PI, 0);
318 dummyGraph->SetPoint(2, M_PI, 0);
319 dummyGraph->Draw("AP");
320 dummyGraph->GetXaxis()->SetTitle("#theta");
321 dummyGraph->GetYaxis()->SetTitle("#rho");
322 dummyGraph->GetXaxis()->SetRangeUser(-M_PI, M_PI);
323 dummyGraph->GetYaxis()->SetRangeUser(-0.02, 0.15);
324
325 for (const CDCWireHit* wireHit : hits) {
326 const double& l = wireHit->getRefDriftLength();
327 const Vector2D& pos2D = wireHit->getRefPos2D() - m_localOrigin;
328 double x = pos2D.x();
329 double y = pos2D.y();
330 double r2 = pos2D.normSquared() - l * l;
331
332 TF1* concaveHitLegendre = new TF1("concaveHitLegendre", "2*([0]/[3])*cos(x) + 2*([1]/[3])*sin(x) + 2*([2]/[3])", -M_PI, M_PI);
333 TF1* convexHitLegendre = new TF1("convexHitLegendre", "2*([0]/[3])*cos(x) + 2*([1]/[3])*sin(x) - 2*([2]/[3])", -M_PI, M_PI);
334 concaveHitLegendre->SetLineWidth(1);
335 convexHitLegendre->SetLineWidth(1);
336 concaveHitLegendre->SetLineColor(color);
337 convexHitLegendre->SetLineColor(color);
338
339 concaveHitLegendre->SetParameters(x, y, l, r2);
340 convexHitLegendre->SetParameters(x, y, l, r2);
341 concaveHitLegendre->Draw("CSAME");
342 convexHitLegendre->Draw("CSAME");
343 }
344// canv->Print(Form("legendreHits_%i.root", nevent));
345// canv->Print(Form("legendreHits_%i.eps", nevent));
346 canv->Print(Form("legendreHits_%i.png", nevent));
347 delete canv;
348
349 nevent++;
350}
double x() const
Getter for the x coordinate.
Definition Vector2D.h:626
double y() const
Getter for the y coordinate.
Definition Vector2D.h:641

◆ drawNode()

void drawNode ( QuadTree * node) const

Draw QuadTree node.

Definition at line 352 of file AxialHitQuadTreeProcessor.cc.

353{
354 std::vector<const CDCWireHit*> hits;
355 for (Item* item : node->getItems()) {
356 const CDCWireHit* wireHit = item->getPointer();
357 hits.push_back(wireHit);
358 }
359 drawHits(hits);
360}
void drawHits(std::vector< const TrackingUtilities::CDCWireHit * > hits, unsigned int color=46) const
Draw QuadTree node.

◆ fill() [1/2]

void fill ( const CandidateReceiver & candidateReceiver,
int nHitsThreshold )
inlineinherited

Start filling the already created tree.

Parameters
candidateReceiverthe lambda function to call after a node was selected
nHitsThresholdthe threshold on the number of items

Definition at line 171 of file QuadTreeProcessor.h.

172 {
173 fill(candidateReceiver, nHitsThreshold, std::numeric_limits<AY>::max());
174 }

◆ fill() [2/2]

void fill ( const CandidateReceiver & candidateReceiver,
int nHitsThreshold,
float yLimit )
inlineinherited

Fill vector of QuadTree instances with hits.

Parameters
candidateReceiverthe lambda function to call after a node was selected
nHitsThresholdthe threshold on the number of items
yLimitthe threshold in the rho (curvature) variable

Definition at line 182 of file QuadTreeProcessor.h.

183 {
184 std::vector<QuadTree*> quadTrees = m_seededTrees;
185 std::sort(quadTrees.begin(), quadTrees.end(), [](QuadTree * quadTree1, QuadTree * quadTree2) {
186 return quadTree1->getNItems() > quadTree2->getNItems();
187 });
188
189 for (QuadTree* tree : quadTrees) {
190 erase_remove_if(tree->getItems(), [](Item * hit) { return hit->isUsed(); });
191 fillGivenTree(tree, candidateReceiver, nHitsThreshold, yLimit);
192 }
193 }

◆ fillChildren()

void fillChildren ( QuadTree * node,
const std::vector< Item * > & items )
inlineprivateinherited

This function is called by fillGivenTree and fills the items into the corresponding children.

For this the user-defined method isInNode is called.

Definition at line 268 of file QuadTreeProcessor.h.

269 {
270 const size_t neededSize = 2 * items.size();
271 for (QuadTree& child : node->getChildren()) {
272 child.reserveItems(neededSize);
273 }
274
275 for (Item* item : items) {
276 if (item->isUsed()) continue;
277
278 for (QuadTree& child : node->getChildren()) {
279 if (isInNode(&child, item->getPointer())) {
280 child.insertItem(item);
281 }
282 }
283 }
284 afterFillDebugHook(node->getChildren());
285 }

◆ fillGivenTree()

void fillGivenTree ( QuadTree * node,
const CandidateReceiver & candidateReceiver,
int nItemsThreshold,
float yLimit )
inlineprivateinherited

Internal function to do the real quad tree search: fill the nodes, check which of the n*m bins we need to process further and go one level deeper.

Definition at line 200 of file QuadTreeProcessor.h.

204 {
205 if (node->getNItems() < nItemsThreshold) {
206 return;
207 }
208
209 if ((node->getYMin() > yLimit) or (-node->getYMax() > yLimit)) {
210 return;
211 }
212
213 if (isLeaf(node)) {
214 callResultFunction(node, candidateReceiver);
215 return;
216 }
217
218 if (node->getChildren().empty()) {
219 this->createChildren(node, node->getChildren());
220 }
221
222 if (!node->checkFilled()) {
223 fillChildren(node, node->getItems());
224 node->setFilled();
225 }
226
227 std::vector<QuadTree*> children;
228 for (QuadTree& child : node->getChildren()) {
229 children.push_back(&child);
230 }
231 const auto compareNItems = [](const QuadTree * lhs, const QuadTree * rhs) {
232 return lhs->getNItems() < rhs->getNItems();
233 };
234
235 // Explicitly count down the children
236 const int nChildren = children.size();
237 for (int iChild = 0; iChild < nChildren; ++iChild) {
238 auto itHeaviestChild = std::max_element(children.begin(), children.end(), compareNItems);
239 QuadTree* heaviestChild = *itHeaviestChild;
240 children.erase(itHeaviestChild);
241 // After we have processed some children we need to get rid of the already used hits in all the children,
242 // because this can change the number of items drastically
243 erase_remove_if(heaviestChild->getItems(), [&](Item * hit) { return hit->isUsed(); });
244 this->fillGivenTree(heaviestChild, candidateReceiver, nItemsThreshold, yLimit);
245 }
246 }

◆ getAssignedItems()

std::vector< const TrackingUtilities::CDCWireHit * > getAssignedItems ( )
inlineinherited

Get items that have been assigned to the seed level The returned elements are unique even if items are assigned multiple times.

Definition at line 152 of file QuadTreeProcessor.h.

153 {
154 std::vector<const TrackingUtilities::CDCWireHit*> result;
155 for (QuadTree* seededTree : m_seededTrees) {
156 for (Item* item : seededTree->getItems()) {
157 result.push_back(item->getPointer());
158 }
159 }
160 std::sort(result.begin(), result.end());
161 result.erase(std::unique(result.begin(), result.end()), result.end());
162 return result;
163 }

◆ getCosSinLookupTable()

const LookupTable< Vector2D > & getCosSinLookupTable ( )
static

Get the standard lookup table containing equally spaces unit vectors (cos, sin)

It contains 2**16 + 1 sampling points between -pi and pi.

Definition at line 89 of file AxialHitQuadTreeProcessor.cc.

90{
91 static const int maxLevel = PrecisionUtil::getLookupGridLevel();
92 static const int nBins = std::pow(2, maxLevel);
93 static LookupTable<Vector2D> trigonometricLookUpTable(&Vector2D::Phi, nBins, -M_PI, M_PI);
94 return trigonometricLookUpTable;
95}
static constexpr int getLookupGridLevel()
Returns desired deepness of the trigonometrical lookup table. Used as template parameter for the Trig...
static Vector2D Phi(const double phi)
Constructs a unit vector with azimuth angle equal to phi.
Definition Vector2D.h:82

◆ getDebugInformation()

const std::map< std::pair< long, float >, std::vector< Item * > > & getDebugInformation ( ) const
inlineinherited

Return the debug information if collected.

Definition at line 371 of file QuadTreeProcessor.h.

372 {
373 return m_debugOutputMap;
374 }

◆ getLastLevel()

int getLastLevel ( ) const
inlineprotectedinherited

Return the parameter last level.

Definition at line 349 of file QuadTreeProcessor.h.

350 {
351 return m_lastLevel;
352 }

◆ isInNode() [1/2]

bool isInNode ( QuadTree * node,
const TrackingUtilities::CDCWireHit * wireHit ) const
finalprotected

Check whether hit belongs to the quadtree node:

Parameters
nodequadtree node
wireHithit being checked
Returns
returns true if sinogram of the hit crosses (geometrically) borders of the node

Definition at line 194 of file AxialHitQuadTreeProcessor.cc.

195{
196 // Check whether the hit lies in the forward direction
197 if (node->getLevel() <= 4 and m_twoSidedPhaseSpace and node->getYMin() > -c_curlCurv and
198 node->getYMax() < c_curlCurv) {
199 if (not checkDerivative(node, wireHit)) return false;
200 }
201
202 const double& l = wireHit->getRefDriftLength();
203 const Vector2D& pos2D = wireHit->getRefPos2D() - m_localOrigin;
204 double r2 = pos2D.normSquared() - l * l;
205
206 using Quadlet = std::array<std::array<float, 2>, 2>;
207 Quadlet distRight{};
208 Quadlet distLeft{};
209
210 // get top and bottom borders of the node
211 float rMin = node->getYMin() * r2 / 2;
212 float rMax = node->getYMax() * r2 / 2;
213
214 // get left and right borders of the node
215 long thetaMin = node->getXMin();
216 long thetaMax = node->getXMax();
217
218 const Vector2D& thetaVecMin = m_cosSinLookupTable->at(thetaMin);
219 const Vector2D& thetaVecMax = m_cosSinLookupTable->at(thetaMax);
220
221 float rHitMin = thetaVecMin.dot(pos2D);
222 float rHitMax = thetaVecMax.dot(pos2D);
223
224 // compute sinograms at the left and right borders of the node
225 float rHitMinRight = rHitMin - l;
226 float rHitMaxRight = rHitMax - l;
227
228 float rHitMinLeft = rHitMin + l;
229 float rHitMaxLeft = rHitMax + l;
230
231 // Compute distance from the sinograms to bottom and top borders of the node
232 distRight[0][0] = rMin - rHitMinRight;
233 distRight[0][1] = rMin - rHitMaxRight;
234 distRight[1][0] = rMax - rHitMinRight;
235 distRight[1][1] = rMax - rHitMaxRight;
236
237 distLeft[0][0] = rMin - rHitMinLeft;
238 distLeft[0][1] = rMin - rHitMaxLeft;
239 distLeft[1][0] = rMax - rHitMinLeft;
240 distLeft[1][1] = rMax - rHitMaxLeft;
241
242 // Compare distance signs from sinograms to the node
243 // Check right
244 if (not sameSign(distRight[0][0], distRight[0][1], distRight[1][0], distRight[1][1])) {
245 return true;
246 }
247
248 // Check left
249 if (not sameSign(distLeft[0][0], distLeft[0][1], distLeft[1][0], distLeft[1][1])) {
250 return true;
251 }
252
253 // Check the extremum
254 float rHitMinExtr = thetaVecMin.cross(pos2D);
255 float rHitMaxExtr = thetaVecMax.cross(pos2D);
256 if (rHitMinExtr * rHitMaxExtr < 0.) return checkExtremum(node, wireHit);
257
258 // Not contained
259 return false;
260}
const double c_curlCurv
The curvature above which the trajectory is considered a curler.
bool checkDerivative(QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const
Check derivative of the sinogram.
bool checkExtremum(QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const
Checks whether extreme point is located within QuadTree node's ranges.
double dot(const Vector2D &rhs) const
Calculates the two dimensional dot product.
Definition Vector2D.h:178
bool sameSign(double expected, double actual)
Predicate checking that two values have the same sign.

◆ isInNode() [2/2]

virtual bool isInNode ( QuadTree * node,
const TrackingUtilities::CDCWireHit * item ) const
protectedpure virtualinherited

Implement that function if you want to provide a new processor.

It is called when filling the quad tree after creation. For every item in a node and every child node this function gets called and should decide, if the item should go into this child node or not.

Parameters
nodechild node
itemitem to be filled into the child node or not
Returns
true if this item belongs into this node.

◆ isLeaf() [1/2]

bool isLeaf ( QuadTree * node) const
finalprotected

lastLevel depends on curvature of the track candidate

Definition at line 121 of file AxialHitQuadTreeProcessor.cc.

122{
123 if (node->getLevel() <= 6) return false;
124 if (node->getLevel() >= getLastLevel()) return true;
125
126 const double nodeResolution = fabs(node->getYMin() - node->getYMax());
127 const double meanCurv = (node->getYMax() + node->getYMin()) / 2;
128
129 const double resolution = m_precisionFunction(meanCurv);
130 if (resolution >= nodeResolution) return true;
131
132 return false;
133}

◆ isLeaf() [2/2]

virtual bool isLeaf ( QuadTree * node) const
inlineprotectedvirtualinherited

Function which checks if given node is leaf Implemented as virtual to keep possibility of changing lastLevel values depending on region is phase-space (i.e.

setting lastLevel as a function of Y-variable)

Definition at line 337 of file QuadTreeProcessor.h.

338 {
339 if (node->getLevel() >= m_lastLevel) {
340 return true;
341 } else {
342 return false;
343 }
344 }

◆ seed()

void seed ( const std::vector< const TrackingUtilities::CDCWireHit * > & datas)
inlineinherited

Fill in the items in the given vector.

They are transformed to QuadTreeItems internally.

Definition at line 106 of file QuadTreeProcessor.h.

107 {
108 // Create the items
109 for (AData* data : datas) {
110 m_items.emplace_back(data);
111 }
112
113 // Creating the seed level
114 long nSeedBins = pow(2, m_seedLevel);
115 m_seededTrees.reserve(nSeedBins * nSeedBins);
116
117 // Expand the first levels to the seed sectors
118 m_seededTrees.push_back(m_quadTree.get());
119 std::vector<QuadTree*> nextSeededTrees;
120
121 for (int level = 0; level < m_seedLevel; ++level) {
122 for (QuadTree* node : m_seededTrees) {
123 if (node->getChildren().empty()) {
124 this->createChildren(node, node->getChildren());
125 }
126 for (QuadTree& child : node->getChildren()) {
127 nextSeededTrees.push_back(&child);
128 }
129 }
130 std::swap(nextSeededTrees, m_seededTrees);
131 nextSeededTrees.clear();
132 }
133
134 // Fill the seed level with the items
135 for (QuadTree* seededTree : m_seededTrees) {
136 seededTree->reserveItems(m_items.size());
137
138 for (Item& item : m_items) {
139 if (item.isUsed()) continue;
140 if (isInNode(seededTree, item.getPointer())) {
141 seededTree->insertItem(&item);
142 }
143 }
144 }
145 }

Member Data Documentation

◆ c_curlCurv

const double c_curlCurv = 0.02
private

The curvature above which the trajectory is considered a curler.

Definition at line 122 of file AxialHitQuadTreeProcessor.h.

◆ m_cosSinLookupTable

const TrackingUtilities::LookupTable<TrackingUtilities::Vector2D>* m_cosSinLookupTable
private

Pinned lookup table for precomputed cosine and sine values.

Definition at line 119 of file AxialHitQuadTreeProcessor.h.

◆ m_debugOutput

bool m_debugOutput
privateinherited

A flag to control the creation of the debug output.

Definition at line 398 of file QuadTreeProcessor.h.

◆ m_debugOutputMap

std::map<std::pair<long, float>, std::vector<Item*> > m_debugOutputMap
privateinherited

The calculated debug map.

Definition at line 401 of file QuadTreeProcessor.h.

◆ m_items

std::deque<Item> m_items
protectedinherited

Storage space for the items that are referenced by the quad tree nodes.

Definition at line 381 of file QuadTreeProcessor.h.

◆ m_lastLevel

int m_lastLevel
privateinherited

The last level to be filled.

Definition at line 392 of file QuadTreeProcessor.h.

◆ m_localOrigin

TrackingUtilities::Vector2D m_localOrigin
private

Local origin on which the phase space coordinates are centered.

Definition at line 116 of file AxialHitQuadTreeProcessor.h.

◆ m_precisionFunction

PrecisionUtil::PrecisionFunction m_precisionFunction
private

Lambda which holds resolution function for the quadtree.

Definition at line 113 of file AxialHitQuadTreeProcessor.h.

◆ m_quadTree

std::unique_ptr<QuadTree> m_quadTree
protectedinherited

The quad tree we work with.

Definition at line 378 of file QuadTreeProcessor.h.

◆ m_seededTrees

std::vector<QuadTree*> m_seededTrees
protectedinherited

Vector of QuadTrees QuadTree instances (which are filled in the vector) cover the whole Legendre phase-space; each instance is processes independently.

Definition at line 388 of file QuadTreeProcessor.h.

◆ m_seedLevel

int m_seedLevel
privateinherited

The first level to be filled, effectively skip forward to this higher granularity level.

Definition at line 395 of file QuadTreeProcessor.h.

◆ m_twoSidedPhaseSpace

bool m_twoSidedPhaseSpace
private

Indicator whether the two sided phases space insertion check should be used This option should automatically split back to back tracks in the low curvature regions.

Definition at line 128 of file AxialHitQuadTreeProcessor.h.


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