Belle II Software development
QuadTreeProcessor< AX, AY, AData > Class Template Referenceabstract

This abstract class serves as a base class for all implementations of track processors. More...

#include <QuadTreeProcessor.h>

Public Types

using Item = QuadTreeItem< AData >
 The QuadTree will only see items of this type.
 
using QuadTree = QuadTreeNode< AX, AY, 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< AData * > &, QuadTree *)>
 This lambda function can be used for postprocessing.
 

Public Member Functions

 QuadTreeProcessor (int lastLevel, int seedLevel, const XYSpans &xySpans, bool debugOutput=false)
 Constructor is very simple.
 
virtual ~QuadTreeProcessor ()
 Destructor deletes the quad tree.
 
void clear ()
 Delete all the QuadTreeItems in the tree and clear the tree.
 
void seed (const std::vector< AData * > &datas)
 Fill in the items in the given vector.
 
std::vector< AData * > 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< AX, AY >, std::vector< Item * > > & getDebugInformation () const
 Return the debug information if collected.
 

Protected Member Functions

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, AData *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, AY 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

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< AX, AY >, std::vector< Item * > > m_debugOutputMap
 The calculated debug map.
 

Detailed Description

template<typename AX, typename AY, class AData>
class Belle2::TrackFindingCDC::QuadTreeProcessor< AX, AY, AData >

This abstract class serves as a base class for all implementations of track processors.

It provides some functions to create, fill, clear and postprocess a quad tree. If you want to use your own class as a quad tree item, you have to overload this processor. You have provide only the two functions isInNode and createChild.

Definition at line 37 of file QuadTreeProcessor.h.

Member Typedef Documentation

◆ CandidateReceiver

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

This lambda function can be used for postprocessing.

Definition at line 59 of file QuadTreeProcessor.h.

◆ Item

using Item = QuadTreeItem<AData>

The QuadTree will only see items of this type.

Definition at line 41 of file QuadTreeProcessor.h.

◆ QuadTree

using QuadTree = QuadTreeNode<AX, AY, Item>

The used QuadTree.

Definition at line 44 of file QuadTreeProcessor.h.

◆ QuadTreeChildren

Alias for the QuadTree Children.

Definition at line 56 of file QuadTreeProcessor.h.

◆ XSpan

using XSpan = typename QuadTree::XSpan

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>

This pair of spans describes the span of a node.

Definition at line 53 of file QuadTreeProcessor.h.

◆ YSpan

using YSpan = typename QuadTree::YSpan

This pair describes the span in Y for a node.

Definition at line 50 of file QuadTreeProcessor.h.

Constructor & Destructor Documentation

◆ QuadTreeProcessor()

QuadTreeProcessor ( int  lastLevel,
int  seedLevel,
const XYSpans xySpans,
bool  debugOutput = false 
)
inline

Constructor is very simple.

The QuadTree has to be constructed elsewhere.

Parameters
lastLeveldescribing the last search level for the quad tree creation
seedLevelfirst level to be filled, effectively skip forward to this higher granularity level
xySpanspair of spans describing the span of a node
debugOutputenable debug output

Definition at line 69 of file QuadTreeProcessor.h.

73 : m_quadTree{std::make_unique<QuadTree>(xySpans.first, xySpans.second, 0, nullptr)}
74 , m_lastLevel(lastLevel)
75 , m_seedLevel(seedLevel)
76 , m_debugOutput(debugOutput)
78 {
79 }
int m_lastLevel
The last level to be filled.
bool m_debugOutput
A flag to control the creation of the debug output.
int m_seedLevel
The first level to be filled, effectively skip forward to this higher granularity level.
std::map< std::pair< AX, AY >, std::vector< Item * > > m_debugOutputMap
The calculated debug map.
std::unique_ptr< QuadTree > m_quadTree
The quad tree we work with.

◆ ~QuadTreeProcessor()

virtual ~QuadTreeProcessor ( )
inlinevirtual

Destructor deletes the quad tree.

Definition at line 84 of file QuadTreeProcessor.h.

85 {
86 clear();
87 }
void clear()
Delete all the QuadTreeItems in the tree and clear the tree.

Member Function Documentation

◆ afterFillDebugHook()

virtual void afterFillDebugHook ( QuadTreeChildren children)
inlinevirtual

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 }
QuadTreeNode< AX, AY, Item > QuadTree
The used QuadTree.
int getLastLevel() const
Return the parameter last level.

◆ callResultFunction()

void callResultFunction ( QuadTree node,
const CandidateReceiver candidateReceiver 
) const
inlineprivate

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 }
QuadTreeItem< AData > Item
The QuadTree will only see items of this type.

◆ clear()

void clear ( )
inline

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()

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

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 }
std::pair< XSpan, YSpan > XYSpans
This pair of spans describes the span of a node.

◆ createChildren()

void createChildren ( QuadTree node,
QuadTreeChildren m_children 
) const
inlineprivate

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 }
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.

◆ fill() [1/2]

void fill ( const CandidateReceiver candidateReceiver,
int  nHitsThreshold 
)
inline

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 
)
inline

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, AY 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 
)
inlineprivate

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, AData *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,
AY  yLimit 
)
inlineprivate

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< AData * > getAssignedItems ( )
inline

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 }

◆ getDebugInformation()

const std::map< std::pair< AX, AY >, std::vector< Item * > > & getDebugInformation ( ) const
inline

Return the debug information if collected.

Definition at line 368 of file QuadTreeProcessor.h.

369 {
370 return m_debugOutputMap;
371 }

◆ getLastLevel()

int getLastLevel ( ) const
inlineprotected

Return the parameter last level.

Definition at line 346 of file QuadTreeProcessor.h.

347 {
348 return m_lastLevel;
349 }

◆ isInNode()

virtual bool isInNode ( QuadTree node,
AData *  item 
) const
protectedpure virtual

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()

virtual bool isLeaf ( QuadTree node) const
inlineprotectedvirtual

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 }

◆ seed()

void seed ( const std::vector< AData * > &  datas)
inline

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 }

Member Data Documentation

◆ m_debugOutput

bool m_debugOutput
private

A flag to control the creation of the debug output.

Definition at line 395 of file QuadTreeProcessor.h.

◆ m_debugOutputMap

std::map<std::pair<AX, AY>, std::vector<Item*> > m_debugOutputMap
private

The calculated debug map.

Definition at line 398 of file QuadTreeProcessor.h.

◆ m_items

std::deque<Item> m_items
protected

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
private

The last level to be filled.

Definition at line 389 of file QuadTreeProcessor.h.

◆ m_quadTree

std::unique_ptr<QuadTree> m_quadTree
protected

The quad tree we work with.

Definition at line 375 of file QuadTreeProcessor.h.

◆ m_seededTrees

std::vector<QuadTree*> m_seededTrees
protected

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
private

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

Definition at line 392 of file QuadTreeProcessor.h.


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