Belle II Software  release-08-01-10
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/trackFindingCDC/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 
23 namespace Belle2 {
28  namespace TrackFindingCDC {
29 
36  template<typename AX, typename AY, class AData>
38 
39  public:
42 
45 
47  using XSpan = typename QuadTree::XSpan;
48 
50  using YSpan = typename QuadTree::YSpan;
51 
53  using XYSpans = std::pair<XSpan, YSpan>;
54 
57 
59  using CandidateReceiver = std::function<void(const std::vector<AData*>&, QuadTree*)>;
60 
61  public:
69  QuadTreeProcessor(int lastLevel,
70  int seedLevel,
71  const XYSpans& xySpans,
72  bool debugOutput = false)
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  }
80 
85  {
86  clear();
87  }
88 
92  void clear()
93  {
94  m_seededTrees.clear();
95  m_quadTree->clearChildren();
96  m_quadTree->clearItems();
97  m_items.clear();
98  }
99 
103  void seed(const std::vector<AData*>& datas)
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  }
143 
144  public:
149  std::vector<AData*> getAssignedItems()
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  }
161 
162  public:
168  void fill(const CandidateReceiver& candidateReceiver, int nHitsThreshold)
169  {
170  fill(candidateReceiver, nHitsThreshold, std::numeric_limits<AY>::max());
171  }
172 
179  void fill(const CandidateReceiver& candidateReceiver, int nHitsThreshold, float yLimit)
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  }
191 
192  private:
198  const CandidateReceiver& candidateReceiver,
199  int nItemsThreshold,
200  AY yLimit)
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  }
244 
249  void createChildren(QuadTree* node, QuadTreeChildren& m_children) const
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  }
260 
265  void fillChildren(QuadTree* node, const std::vector<Item*>& items)
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  }
282  }
283 
288  void callResultFunction(QuadTree* node, const CandidateReceiver& candidateReceiver) const
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  }
301 
302  protected: // Section of specialisable functions
311  virtual XYSpans createChild(QuadTree* node, int iX, int iY) const
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  }
319 
327  virtual bool isInNode(QuadTree* node, AData* item) const = 0;
328 
334  virtual bool isLeaf(QuadTree* node) const
335  {
336  if (node->getLevel() >= m_lastLevel) {
337  return true;
338  } else {
339  return false;
340  }
341  }
342 
346  int getLastLevel() const
347  {
348  return m_lastLevel;
349  }
350 
351  public: // debug stuff
356  virtual void afterFillDebugHook(QuadTreeChildren& children)
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  }
364 
368  const std::map<std::pair<AX, AY>, std::vector<Item*>>& getDebugInformation() const
369  {
370  return m_debugOutputMap;
371  }
372 
373  protected:
375  std::unique_ptr<QuadTree> m_quadTree;
376 
378  std::deque<Item> m_items;
379 
385  std::vector<QuadTree*> m_seededTrees;
386 
387  private:
390 
393 
396 
398  std::map<std::pair<AX, AY>, std::vector<Item*>> m_debugOutputMap;
399  };
400  }
402 }
This class serves as a wrapper around all things that should go into a QuadTree.
Definition: QuadTreeItem.h:27
Class which holds quadtree structure.
Definition: QuadTreeNode.h:29
std::vector< This > Children
Type of the child node structure for this node.
Definition: QuadTreeNode.h:48
AY getYMax() const
Get maximal "r" value of the node.
Definition: QuadTreeNode.h:198
std::array< AY, 2 > YSpan
Type for a span in the Y direction that is covered by the tree.
Definition: QuadTreeNode.h:39
void setFilled()
Set status of node to "filled" (children nodes has been filled)
Definition: QuadTreeNode.h:138
AY getYMin() const
Get minimal "r" value of the node.
Definition: QuadTreeNode.h:192
int getLevel() const
Returns level of the node in tree (i.e., how much ancestors the node has)
Definition: QuadTreeNode.h:126
AX getXLowerBound(int iBin) const
Get lower "Theta" value of given bin.
Definition: QuadTreeNode.h:174
std::vector< AItem * > & getItems()
Get items from node.
Definition: QuadTreeNode.h:91
AY getYUpperBound(int iBin) const
Get upper "r" value of given bin.
Definition: QuadTreeNode.h:215
std::array< AX, 2 > XSpan
Type for a span in the X direction that is covered by the tree.
Definition: QuadTreeNode.h:36
AY getYLowerBound(int iBin) const
Get lower "r" value of given bin.
Definition: QuadTreeNode.h:209
bool checkFilled() const
Check whether node has been processed, i.e.
Definition: QuadTreeNode.h:132
int getNItems() const
Check if the node passes threshold on number of hits.
Definition: QuadTreeNode.h:97
constexpr int getXNbins() const
Get number of bins in "Theta" direction.
Definition: QuadTreeNode.h:150
constexpr int getYNbins() const
Get number of bins in "r" direction.
Definition: QuadTreeNode.h:186
Children & getChildren()
Returns the children structure of this node.
Definition: QuadTreeNode.h:109
AX getXUpperBound(int iBin) const
Get upper "Theta" value of given bin.
Definition: QuadTreeNode.h:180
This abstract class serves as a base class for all implementations of track processors.
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.
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.
const std::map< std::pair< AX, AY >, std::vector< Item * > > & getDebugInformation() const
Return the debug information if collected.
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.
Abstract base class for different kinds of events.