Belle II Software development
QuadTreeProcessor.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8#pragma once
9
10#include <tracking/trackFindingCDC/legendre/quadtree/QuadTreeNode.h>
11#include <tracking/trackFindingCDC/legendre/quadtree/QuadTreeItem.h>
12
13#include <tracking/trackingUtilities/utilities/Algorithms.h>
14
15#include <algorithm>
16#include <functional>
17#include <memory>
18#include <map>
19#include <vector>
20#include <deque>
21#include <utility>
22
23namespace Belle2 {
28 namespace TrackingUtilities {
29 class CDCWireHit;
30 }
31 namespace TrackFindingCDC {
32
39 template<typename AX, typename AY, class AData>
41
42 public:
45
48
50 using XSpan = typename QuadTree::XSpan;
51
53 using YSpan = typename QuadTree::YSpan;
54
56 using XYSpans = std::pair<XSpan, YSpan>;
57
60
62 using CandidateReceiver = std::function<void(const std::vector<AData*>&, QuadTree*)>;
63
64 public:
72 QuadTreeProcessor(int lastLevel,
73 int seedLevel,
74 const XYSpans& xySpans,
75 bool debugOutput = false)
76 : m_quadTree{std::make_unique<QuadTree>(xySpans.first, xySpans.second, 0, nullptr)}
77 , m_lastLevel(lastLevel)
78 , m_seedLevel(seedLevel)
79 , m_debugOutput(debugOutput)
81 {
82 }
83
88 {
89 clear();
90 }
91
95 void clear()
96 {
97 m_seededTrees.clear();
98 m_quadTree->clearChildren();
99 m_quadTree->clearItems();
100 m_items.clear();
101 }
102
106 void seed(const std::vector<AData*>& datas)
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 }
146
147 public:
152 std::vector<AData*> getAssignedItems()
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 }
164
165 public:
171 void fill(const CandidateReceiver& candidateReceiver, int nHitsThreshold)
172 {
173 fill(candidateReceiver, nHitsThreshold, std::numeric_limits<AY>::max());
174 }
175
182 void fill(const CandidateReceiver& candidateReceiver, int nHitsThreshold, float yLimit)
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 }
194
195 private:
201 const CandidateReceiver& candidateReceiver,
202 int nItemsThreshold,
203 AY yLimit)
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 }
247
252 void createChildren(QuadTree* node, QuadTreeChildren& m_children) const
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 }
263
268 void fillChildren(QuadTree* node, const std::vector<Item*>& items)
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 }
285 }
286
291 void callResultFunction(QuadTree* node, const CandidateReceiver& candidateReceiver) const
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 }
304
305 protected: // Section of specialisable functions
314 virtual XYSpans createChild(QuadTree* node, int iX, int iY) const
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 }
322
330 virtual bool isInNode(QuadTree* node, AData* item) const = 0;
331
337 virtual bool isLeaf(QuadTree* node) const
338 {
339 if (node->getLevel() >= m_lastLevel) {
340 return true;
341 } else {
342 return false;
343 }
344 }
345
349 int getLastLevel() const
350 {
351 return m_lastLevel;
352 }
353
354 public: // debug stuff
359 virtual void afterFillDebugHook(QuadTreeChildren& children)
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 }
367
371 const std::map<std::pair<AX, AY>, std::vector<Item*>>& getDebugInformation() const
372 {
373 return m_debugOutputMap;
374 }
375
376 protected:
378 std::unique_ptr<QuadTree> m_quadTree;
379
381 std::deque<Item> m_items;
382
388 std::vector<QuadTree*> m_seededTrees;
389
390 private:
393
396
399
401 std::map<std::pair<AX, AY>, std::vector<Item*>> m_debugOutputMap;
402 };
403 }
405}
This class serves as a wrapper around all things that should go into a QuadTree.
Class which holds quadtree structure.
Children & getChildren()
Returns the children structure of this node.
std::vector< AItem * > & getItems()
Get items from node.
AY getYMax() const
Get maximal "r" value of the node.
void setFilled()
Set status of node to "filled" (children nodes has been filled)
AY getYMin() const
Get minimal "r" value of the node.
int getLevel() const
Returns level of the node in tree (i.e., how much ancestors the node has)
AX getXLowerBound(int iBin) const
Get lower "Theta" value of given bin.
AY getYUpperBound(int iBin) const
Get upper "r" value of given bin.
AY getYLowerBound(int iBin) const
Get lower "r" value of given bin.
bool checkFilled() const
Check whether node has been processed, i.e.
int getNItems() const
Check if the node passes threshold on number of hits.
constexpr int getXNbins() const
Get number of bins in "Theta" direction.
constexpr int getYNbins() const
Get number of bins in "r" direction.
AX getXUpperBound(int iBin) const
Get upper "Theta" value of given bin.
void fill(const CandidateReceiver &candidateReceiver, int nHitsThreshold, float yLimit)
Fill vector of QuadTree instances with hits.
QuadTreeNode< AX, AY, Item > QuadTree
The used QuadTree.
typename QuadTree::Children QuadTreeChildren
Alias for the QuadTree Children.
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.
int m_lastLevel
The last level to be filled.
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...
virtual ~QuadTreeProcessor()
Destructor deletes the quad tree.
const std::map< std::pair< AX, AY >, std::vector< Item * > > & getDebugInformation() const
Return the debug information if collected.
std::pair< XSpan, YSpan > XYSpans
This pair of spans describes the span of a node.
QuadTreeProcessor(int lastLevel, int seedLevel, const XYSpans &xySpans, bool debugOutput=false)
Constructor is very simple.
std::vector< QuadTree * > m_seededTrees
Vector of QuadTrees QuadTree instances (which are filled in the vector) cover the whole Legendre phas...
bool m_debugOutput
A flag to control the creation of the debug output.
int getLastLevel() const
Return the parameter last level.
virtual XYSpans createChild(QuadTree *node, int iX, int iY) const
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 la...
typename QuadTree::XSpan XSpan
This pair describes the span in X for a node.
std::deque< Item > m_items
Storage space for the items that are referenced by the quad tree nodes.
std::vector< AData * > getAssignedItems()
Get items that have been assigned to the seed level The returned elements are unique even if items ar...
typename QuadTree::YSpan YSpan
This pair describes the span in Y for a node.
void fill(const CandidateReceiver &candidateReceiver, int nHitsThreshold)
Start filling the already created tree.
void clear()
Delete all the QuadTreeItems in the tree and clear the tree.
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*) ...
int m_seedLevel
The first level to be filled, effectively skip forward to this higher granularity level.
void seed(const std::vector< AData * > &datas)
Fill in the items in the given vector.
std::function< void(const std::vector< AData * > &, QuadTree *)> CandidateReceiver
This lambda function can be used for postprocessing.
QuadTreeItem< AData > Item
The QuadTree will only see items of this type.
virtual bool isInNode(QuadTree *node, AData *item) const =0
Implement that function if you want to provide a new processor.
std::map< std::pair< AX, AY >, std::vector< Item * > > m_debugOutputMap
The calculated debug map.
virtual void afterFillDebugHook(QuadTreeChildren &children)
Override that function if you want to receive debug output whenever the children of a node are filled...
std::unique_ptr< QuadTree > m_quadTree
The quad tree we work with.
Class representing a hit wire in the central drift chamber.
Definition CDCWireHit.h:58
Abstract base class for different kinds of events.
STL namespace.