Belle II Software development
AxialHitQuadTreeProcessor Class Referenceabstract

A QuadTreeProcessor for TrackHits. More...

#include <AxialHitQuadTreeProcessor.h>

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

Public Types

using Item = QuadTreeItem< const CDCWireHit >
 The QuadTree will only see items of this type.
 
using QuadTree = QuadTreeNode< long, float, Item >
 The used QuadTree.
 
using XSpan = typename QuadTree::XSpan
 This pair describes the span in X for a node.
 
using YSpan = typename QuadTree::YSpan
 This pair describes the span in Y for a node.
 
using XYSpans = std::pair< XSpan, YSpan >
 This pair of spans describes the span of a node.
 
using QuadTreeChildren = typename QuadTree::Children
 Alias for the QuadTree Children.
 
using CandidateReceiver = std::function< void(const std::vector< const CDCWireHit * > &, QuadTree *)>
 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 Vector2D &localOrigin, const YSpan &curvSpan, const LookupTable< Vector2D > *cosSinLookupTable)
 Constructor for the quad tree processor used in the off-origin extension.
 
void drawHits (std::vector< const 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 CDCWireHit * > &datas)
 Fill in the items in the given vector.
 
std::vector< const 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 LookupTable< 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 CDCWireHit *wireHit) const final
 Check whether hit belongs to the quadtree node:
 
bool checkDerivative (QuadTree *node, const CDCWireHit *wireHit) const
 Check derivative of the sinogram.
 
bool checkExtremum (QuadTree *node, const CDCWireHit *wireHit) const
 Checks whether extremum point is located whithin 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 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

const LookupTable< Vector2D > * m_cosSinLookupTable
 Pinned lookup table for precompute cosine and sine values.
 
Vector2D m_localOrigin
 Local origin on which the phase space coordinates are centered.
 
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.
 
const double c_curlCurv = 0.02
 The curvature above which the trajectory is considered a curler.
 
PrecisionUtil::PrecisionFunction m_precisionFunction
 Lambda which holds resolution function for the quadtree.
 
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 26 of file AxialHitQuadTreeProcessor.h.

Member Typedef Documentation

◆ CandidateReceiver

using CandidateReceiver = std::function<void(const std::vector<const CDCWireHit *>&, QuadTree*)>
inherited

This lambda function can be used for postprocessing.

Definition at line 59 of file QuadTreeProcessor.h.

◆ Item

using Item = QuadTreeItem<const CDCWireHit >
inherited

The QuadTree will only see items of this type.

Definition at line 41 of file QuadTreeProcessor.h.

◆ QuadTree

using QuadTree = QuadTreeNode<long , float , Item>
inherited

The used QuadTree.

Definition at line 44 of file QuadTreeProcessor.h.

◆ QuadTreeChildren

using QuadTreeChildren = typename QuadTree::Children
inherited

Alias for the QuadTree Children.

Definition at line 56 of file QuadTreeProcessor.h.

◆ XSpan

using XSpan = typename QuadTree::XSpan
inherited

This pair describes the span in X for a node.

Definition at line 47 of file QuadTreeProcessor.h.

◆ XYSpans

using XYSpans = std::pair<XSpan, YSpan>
inherited

This pair of spans describes the span of a node.

Definition at line 53 of file QuadTreeProcessor.h.

◆ YSpan

using YSpan = typename QuadTree::YSpan
inherited

This pair describes the span in Y for a node.

Definition at line 50 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 96 of file AxialHitQuadTreeProcessor.cc.

100 : QuadTreeProcessor(lastLevel, seedLevel, ranges)
102 , m_localOrigin(0.0, 0.0)
103 , m_precisionFunction(precisionFunction)
104{
105 m_twoSidedPhaseSpace = m_quadTree->getYMin() * m_quadTree->getYMax() < 0;
106}
bool m_twoSidedPhaseSpace
Indicator whether the two sided phases space insertion check should be used This option should automa...
Vector2D m_localOrigin
Local origin on which the phase space coordinates are centered.
const LookupTable< Vector2D > * m_cosSinLookupTable
Pinned lookup table for precompute cosine and sine values.
PrecisionUtil::PrecisionFunction m_precisionFunction
Lambda which holds resolution function for the quadtree.
static const LookupTable< Vector2D > & getCosSinLookupTable()
Get the standard lookup table containing equally spaces unit vectors (cos, sin)
This abstract class serves as a base class for all implementations of track processors.
std::unique_ptr< QuadTree > m_quadTree
The quad tree we work with.

◆ AxialHitQuadTreeProcessor() [2/2]

AxialHitQuadTreeProcessor ( const Vector2D localOrigin,
const YSpan curvSpan,
const LookupTable< 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 108 of file AxialHitQuadTreeProcessor.cc.

111 : QuadTreeProcessor(0, 0, { {0, cosSinLookupTable->getNPoints() - 1}, curvSpan})
112, m_cosSinLookupTable(cosSinLookupTable)
113, m_localOrigin(localOrigin)
114{
115 // Never use two sided mode in off origin extension
116 m_twoSidedPhaseSpace = false;
117}
int getNPoints() const
Return the number of finite sampling points in this lookup table.
Definition: LookupTable.h:108

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 356 of file QuadTreeProcessor.h.

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

◆ 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 288 of file QuadTreeProcessor.h.

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

◆ checkDerivative()

bool checkDerivative ( QuadTree node,
const 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 261 of file AxialHitQuadTreeProcessor.cc.

262{
263 const Vector2D& pos2D = wireHit->getRefPos2D() - m_localOrigin;
264
265 long thetaMin = node->getXMin();
266 long thetaMax = node->getXMax();
267
268 const Vector2D& thetaVecMin = m_cosSinLookupTable->at(thetaMin);
269 const Vector2D& thetaVecMax = m_cosSinLookupTable->at(thetaMax);
270
271 float rMinD = thetaVecMin.cross(pos2D);
272 float rMaxD = thetaVecMax.cross(pos2D);
273
274 // Does not really make sense...
275 if ((rMinD > 0) && (rMaxD * rMinD >= 0)) return true;
276 if ((rMaxD * rMinD < 0)) return true;
277 return false;
278}
A two dimensional vector which is equipped with functions for correct handeling of orientation relat...
Definition: Vector2D.h:32
double cross(const Vector2D &rhs) const
Calculated the two dimensional cross product.
Definition: Vector2D.h:163

◆ checkExtremum()

bool checkExtremum ( QuadTree node,
const CDCWireHit wireHit 
) const
protected

Checks whether extremum point is located whithin QuadTree node's ranges.

Parameters
nodeQuadTree node
wireHithit to check
Returns
true or false

Definition at line 280 of file AxialHitQuadTreeProcessor.cc.

281{
282 const double& l = wireHit->getRefDriftLength();
283 const Vector2D& pos2D = wireHit->getRefPos2D() - m_localOrigin;
284 double r2 = pos2D.normSquared() - l * l;
285
286 // get left and right borders of the node
287 long thetaMin = node->getXMin();
288 long thetaMax = node->getXMax();
289
290 const Vector2D& thetaVecMin = m_cosSinLookupTable->at(thetaMin);
291 const Vector2D& thetaVecMax = m_cosSinLookupTable->at(thetaMax);
292
293 if (not pos2D.isBetween(thetaVecMin, thetaVecMax)) return false;
294
295 // compute sinograms at the position
296 double r = pos2D.norm();
297 float rRight = r - l;
298 float rLeft = r + l;
299
300 // get top and bottom borders of the node
301 float rMin = node->getYMin() * r2 / 2;
302 float rMax = node->getYMax() * r2 / 2;
303
304 bool crossesRight = (rMin - rRight) * (rMax - rRight) < 0;
305 bool crossesLeft = (rMin - rLeft) * (rMax - rLeft) < 0;
306 return crossesRight or crossesLeft;
307}
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:525
double normSquared() const
Calculates .
Definition: Vector2D.h:169
double norm() const
Calculates the length of the vector.
Definition: Vector2D.h:175

◆ clear()

void clear ( )
inlineinherited

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

Definition at line 92 of file QuadTreeProcessor.h.

93 {
94 m_seededTrees.clear();
95 m_quadTree->clearChildren();
96 m_quadTree->clearItems();
97 m_items.clear();
98 }
std::vector< QuadTree * > m_seededTrees
Vector of QuadTrees QuadTree instances (which are filled in the vector) cover the whole Legendre phas...
std::deque< Item > m_items
Storage space for the items that are referenced by the quad tree nodes.

◆ 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 135 of file AxialHitQuadTreeProcessor.cc.

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

◆ 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 311 of file QuadTreeProcessor.h.

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

◆ 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 249 of file QuadTreeProcessor.h.

250 {
251 for (int i = 0; i < node->getXNbins(); ++i) {
252 for (int j = 0; j < node->getYNbins(); ++j) {
253 const XYSpans& xySpans = createChild(node, i, j);
254 const XSpan& xSpan = xySpans.first;
255 const YSpan& ySpan = xySpans.second;
256 m_children.push_back(QuadTree(xSpan, ySpan, node->getLevel() + 1, node));
257 }
258 }
259 }
QuadTreeNode< long, float, Item > QuadTree
The used QuadTree.
virtual XYSpans createChild(QuadTree *node, int iX, int iY) const
Implement that function if you want to provide a new processor.
typename QuadTree::XSpan XSpan
This pair describes the span in X for a node.
typename QuadTree::YSpan YSpan
This pair describes the span in Y for a node.

◆ 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 64 of file AxialHitQuadTreeProcessor.cc.

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

◆ drawHits()

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

Draw QuadTree node.

Definition at line 309 of file AxialHitQuadTreeProcessor.cc.

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

◆ drawNode()

void drawNode ( QuadTree node) const

Draw QuadTree node.

Definition at line 351 of file AxialHitQuadTreeProcessor.cc.

352{
353 std::vector<const CDCWireHit*> hits;
354 for (Item* item : node->getItems()) {
355 const CDCWireHit* wireHit = item->getPointer();
356 hits.push_back(wireHit);
357 }
358 drawHits(hits);
359}
void drawHits(std::vector< const CDCWireHit * > hits, unsigned int color=46) const
Draw QuadTree node.
QuadTreeItem< const CDCWireHit > Item
The QuadTree will only see items of this type.

◆ 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 168 of file QuadTreeProcessor.h.

169 {
170 fill(candidateReceiver, nHitsThreshold, std::numeric_limits<AY>::max());
171 }
void fill(const CandidateReceiver &candidateReceiver, int nHitsThreshold)
Start filling the already created tree.

◆ 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 179 of file QuadTreeProcessor.h.

180 {
181 std::vector<QuadTree*> quadTrees = m_seededTrees;
182 std::sort(quadTrees.begin(), quadTrees.end(), [](QuadTree * quadTree1, QuadTree * quadTree2) {
183 return quadTree1->getNItems() > quadTree2->getNItems();
184 });
185
186 for (QuadTree* tree : quadTrees) {
187 erase_remove_if(tree->getItems(), [](Item * hit) { return hit->isUsed(); });
188 fillGivenTree(tree, candidateReceiver, nHitsThreshold, yLimit);
189 }
190 }
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 nee...

◆ 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 265 of file QuadTreeProcessor.h.

266 {
267 const size_t neededSize = 2 * items.size();
268 for (QuadTree& child : node->getChildren()) {
269 child.reserveItems(neededSize);
270 }
271
272 for (Item* item : items) {
273 if (item->isUsed()) continue;
274
275 for (QuadTree& child : node->getChildren()) {
276 if (isInNode(&child, item->getPointer())) {
277 child.insertItem(item);
278 }
279 }
280 }
281 afterFillDebugHook(node->getChildren());
282 }
virtual bool isInNode(QuadTree *node, const CDCWireHit *item) const=0
Implement that function if you want to provide a new processor.
virtual void afterFillDebugHook(QuadTreeChildren &children)
Override that function if you want to receive debug output whenever the children of a node are filled...

◆ 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 197 of file QuadTreeProcessor.h.

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

◆ getAssignedItems()

std::vector< const 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 149 of file QuadTreeProcessor.h.

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

◆ 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 88 of file AxialHitQuadTreeProcessor.cc.

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

◆ getDebugInformation()

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

Return the debug information if collected.

Definition at line 368 of file QuadTreeProcessor.h.

369 {
370 return m_debugOutputMap;
371 }
std::map< std::pair< long, float >, std::vector< Item * > > m_debugOutputMap
The calculated debug map.

◆ getLastLevel()

int getLastLevel ( ) const
inlineprotectedinherited

Return the parameter last level.

Definition at line 346 of file QuadTreeProcessor.h.

347 {
348 return m_lastLevel;
349 }

◆ isInNode() [1/2]

virtual bool isInNode ( QuadTree node,
const 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.

◆ isInNode() [2/2]

bool isInNode ( QuadTree node,
const 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 193 of file AxialHitQuadTreeProcessor.cc.

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

◆ isLeaf() [1/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 334 of file QuadTreeProcessor.h.

335 {
336 if (node->getLevel() >= m_lastLevel) {
337 return true;
338 } else {
339 return false;
340 }
341 }

◆ isLeaf() [2/2]

bool isLeaf ( QuadTree node) const
finalprotected

lastLevel depends on curvature of the track candidate

Definition at line 120 of file AxialHitQuadTreeProcessor.cc.

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

◆ seed()

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

Fill in the items in the given vector.

They are transformed to QuadTreeItems internally.

Definition at line 103 of file QuadTreeProcessor.h.

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

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 120 of file AxialHitQuadTreeProcessor.h.

◆ m_cosSinLookupTable

const LookupTable<Vector2D>* m_cosSinLookupTable
private

Pinned lookup table for precompute cosine and sine values.

Definition at line 108 of file AxialHitQuadTreeProcessor.h.

◆ m_debugOutput

bool m_debugOutput
privateinherited

A flag to control the creation of the debug output.

Definition at line 395 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 398 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 378 of file QuadTreeProcessor.h.

◆ m_lastLevel

int m_lastLevel
privateinherited

The last level to be filled.

Definition at line 389 of file QuadTreeProcessor.h.

◆ m_localOrigin

Vector2D m_localOrigin
private

Local origin on which the phase space coordinates are centered.

Definition at line 111 of file AxialHitQuadTreeProcessor.h.

◆ m_precisionFunction

PrecisionUtil::PrecisionFunction m_precisionFunction
private

Lambda which holds resolution function for the quadtree.

Definition at line 123 of file AxialHitQuadTreeProcessor.h.

◆ m_quadTree

std::unique_ptr<QuadTree> m_quadTree
protectedinherited

The quad tree we work with.

Definition at line 375 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 385 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 392 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 117 of file AxialHitQuadTreeProcessor.h.


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