Belle II Software development
AxialHitQuadTreeProcessor.cc
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#include <tracking/trackFindingCDC/legendre/quadtree/AxialHitQuadTreeProcessor.h>
9
10#include <tracking/trackingUtilities/eventdata/hits/CDCWireHit.h>
11
12#include <array>
13#include <vector>
14
15#include <TF1.h>
16#include <TCanvas.h>
17#include <TGraph.h>
18#include <TAxis.h>
19
20using namespace Belle2;
21using namespace TrackFindingCDC;
22using namespace TrackingUtilities;
23
24namespace {
25 bool sameSign(double n1, double n2, double n3, double n4)
26 {
27 return ((n1 > 0 && n2 > 0 && n3 > 0 && n4 > 0) || (n1 < 0 && n2 < 0 && n3 < 0 && n4 < 0));
28 }
29
30 using YSpan = AxialHitQuadTreeProcessor::YSpan;
31 YSpan splitCurvSpan(const YSpan& curvSpan, int nodeLevel, int lastLevel, int j)
32 {
33 const float meanCurv = curvSpan[0] + (curvSpan[1] - curvSpan[0]) / 2.0;
34 const std::array<float, 3> binBounds{curvSpan[0], meanCurv, curvSpan[1]};
35 const float binWidth = binBounds[j + 1] - binBounds[j];
36
37 const bool standardBinning = (nodeLevel <= lastLevel - 7) or (std::fabs(meanCurv) <= 0.005);
38 if (standardBinning) {
39 // B2INFO("Case 1* " << meanCurv << " " << (meanCurv <= 0.005));
40 float curv1 = binBounds[j];
41 float curv2 = binBounds[j + 1];
42
43 // Standard bin division
44 return {curv1, curv2};
45 }
46
47 // Non-standard binning
48 // For level 6 to 7 only expand 1 / 4, for higher levels expand 1 / 8.
49 // (assuming last level == 12)
50 if (nodeLevel < lastLevel - 5) {
51 // B2INFO("Case 2*");
52 float curv1 = binBounds[j] - binWidth / 4.;
53 float curv2 = binBounds[j + 1] + binWidth / 4.;
54 return {curv1, curv2};
55 } else {
56 // B2INFO("Case 3*");
57 float curv1 = binBounds[j] - binWidth / 8.;
58 float curv2 = binBounds[j + 1] + binWidth / 8.;
59 return {curv1, curv2};
60 }
61 }
62
63}
64
65std::vector<float> AxialHitQuadTreeProcessor::createCurvBound(YSpan curvSpan, int lastLevel)
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}
88
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}
96
98 int seedLevel,
99 const XYSpans& ranges,
100 PrecisionUtil::PrecisionFunction precisionFunction)
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}
108
110 const YSpan& curvSpan,
111 const LookupTable<Vector2D>* cosSinLookupTable)
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}
119
120
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}
134
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}
193
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}
261
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}
280
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}
309
310void AxialHitQuadTreeProcessor::drawHits(std::vector<const CDCWireHit*> hits, unsigned int color) const
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}
351
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}
XYSpans createChild(QuadTree *node, int i, int j) const final
Return the new ranges.
const double c_curlCurv
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 automa...
void drawHits(std::vector< const TrackingUtilities::CDCWireHit * > hits, unsigned int color=46) const
Draw QuadTree node.
TrackingUtilities::Vector2D m_localOrigin
Local origin on which the phase space coordinates are centered.
bool checkDerivative(QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const
Check derivative of the sinogram.
bool isLeaf(QuadTree *node) const final
lastLevel depends on curvature of the track candidate
bool isInNode(QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const final
Check whether hit belongs to the quadtree node:
bool checkExtremum(QuadTree *node, const TrackingUtilities::CDCWireHit *wireHit) const
Checks whether extreme point is located within QuadTree node's ranges.
AxialHitQuadTreeProcessor(int lastLevel, int seedLevel, const XYSpans &ranges, PrecisionUtil::PrecisionFunction precisionFunction)
Constructor.
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)
void drawNode(QuadTree *node) const
Draw QuadTree node.
static std::vector< float > createCurvBound(YSpan curvSpan, int lastLevel)
Constructs an array with the curvature bounds as generated by the default bin divisions.
const TrackingUtilities::LookupTable< TrackingUtilities::Vector2D > * m_cosSinLookupTable
Pinned lookup table for precomputed cosine and sine values.
static constexpr int getLookupGridLevel()
Returns desired deepness of the trigonometrical lookup table. Used as template parameter for the Trig...
std::function< double(double)> PrecisionFunction
Function type which is used for resolution calculations (resolution=f(curvature)) Takes a curvature v...
std::vector< AItem * > & getItems()
Get items from node.
AY getYBinWidth(int iBin)
Getter for the width of the iBin bin in "r" direction.
AX getXMax() const
Get maximal "Theta" value of the node.
AY getYMax() const
Get maximal "r" value of the node.
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.
AX getXMin() const
Get minimal "Theta" value of the node.
AX getXUpperBound(int iBin) const
Get upper "Theta" value of given bin.
QuadTreeProcessor(int lastLevel, int seedLevel, const XYSpans &xySpans, bool debugOutput=false)
Class representing a hit wire in the central drift chamber.
Definition CDCWireHit.h:58
Class which holds precomputed values of a function.
Definition LookupTable.h:50
int getNPoints() const
Return the number of finite sampling points in this lookup table.
A two dimensional vector which is equipped with functions for correct handling of orientation relate...
Definition Vector2D.h:36
double dot(const Vector2D &rhs) const
Calculates the two dimensional dot product.
Definition Vector2D.h:178
double cross(const Vector2D &rhs) const
Calculated the two dimensional cross product.
Definition Vector2D.h:189
double x() const
Getter for the x coordinate.
Definition Vector2D.h:626
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 y() const
Getter for the y coordinate.
Definition Vector2D.h:641
double norm() const
Calculates the length of the vector.
Definition Vector2D.h:206
static Vector2D Phi(const double phi)
Constructs a unit vector with azimuth angle equal to phi.
Definition Vector2D.h:82
bool sameSign(double expected, double actual)
Predicate checking that two values have the same sign.
Abstract base class for different kinds of events.