8#include <tracking/trackFindingCDC/legendre/quadtree/AxialHitQuadTreeProcessor.h>
10#include <tracking/trackFindingCDC/eventdata/hits/CDCWireHit.h>
21using namespace TrackFindingCDC;
24 bool sameSign(
double n1,
double n2,
double n3,
double n4)
26 return ((n1 > 0 && n2 > 0 && n3 > 0 && n4 > 0) || (n1 < 0 && n2 < 0 && n3 < 0 && n4 < 0));
29 using YSpan = AxialHitQuadTreeProcessor::YSpan;
30 YSpan splitCurvSpan(
const YSpan& curvSpan,
int nodeLevel,
int lastLevel,
int j)
32 const float meanCurv = curvSpan[0] + (curvSpan[1] - curvSpan[0]) / 2.0;
33 const std::array<float, 3> binBounds{curvSpan[0], meanCurv, curvSpan[1]};
34 const float binWidth = binBounds[j + 1] - binBounds[j];
36 const bool standardBinning = (nodeLevel <= lastLevel - 7) or (std::fabs(meanCurv) <= 0.005);
37 if (standardBinning) {
39 float curv1 = binBounds[j];
40 float curv2 = binBounds[j + 1];
43 return {curv1, curv2};
49 if (nodeLevel < lastLevel - 5) {
51 float curv1 = binBounds[j] - binWidth / 4.;
52 float curv2 = binBounds[j + 1] + binWidth / 4.;
53 return {curv1, curv2};
56 float curv1 = binBounds[j] - binWidth / 8.;
57 float curv2 = binBounds[j + 1] + binWidth / 8.;
58 return {curv1, curv2};
66 std::vector<YSpan> spans{{curvSpan}};
68 std::vector<YSpan> nextSpans;
69 for (
int level = 1; level <= lastLevel; ++level) {
71 for (
const YSpan& span : spans) {
72 nextSpans.push_back(splitCurvSpan(span, level, lastLevel, 0));
73 nextSpans.push_back(splitCurvSpan(span, level, lastLevel, 1));
75 spans.swap(nextSpans);
78 std::vector<float> bounds;
79 for (
const YSpan& span : spans) {
80 bounds.push_back(span[0]);
81 bounds.push_back(span[1]);
84 assert(bounds.size() == std::pow(2, lastLevel));
91 static const int nBins = std::pow(2, maxLevel);
93 return trigonometricLookUpTable;
101 , m_cosSinLookupTable(&getCosSinLookupTable())
102 , m_localOrigin(0.0, 0.0)
103 , m_precisionFunction(precisionFunction)
109 const YSpan& curvSpan,
112, m_cosSinLookupTable(cosSinLookupTable)
113, m_localOrigin(localOrigin)
116 m_twoSidedPhaseSpace =
false;
122 if (node->
getLevel() <= 6)
return false;
125 const double nodeResolution = fabs(node->
getYMin() - node->
getYMax());
129 if (resolution >= nodeResolution)
return true;
137 const int nodeLevel = node->
getLevel();
139 const float meanCurv = std::fabs(node->
getYMax() + node->
getYMin()) / 2;
143 bool standardBinning = (nodeLevel <= lastLevel - 7) or (meanCurv <= 0.005);
145 if (standardBinning) {
152 return XYSpans({theta1, theta2}, {r1, r2});
158 if (nodeLevel < lastLevel - 5) {
163 long extension = pow(2, lastLevel - nodeLevel - 2);
166 if (theta1 < 0) theta1 = 0;
173 return XYSpans({theta1, theta2}, {r1, r2});
179 long extension = pow(2, lastLevel - nodeLevel - 3);
182 if (theta1 < 0) theta1 = 0;
189 return XYSpans({theta1, theta2}, {r1, r2});
201 const double& l = wireHit->getRefDriftLength();
205 using Quadlet = std::array<std::array<float, 2>, 2>;
210 float rMin = node->
getYMin() * r2 / 2;
211 float rMax = node->
getYMax() * r2 / 2;
214 long thetaMin = node->
getXMin();
215 long thetaMax = node->
getXMax();
220 float rHitMin = thetaVecMin.
dot(pos2D);
221 float rHitMax = thetaVecMax.
dot(pos2D);
224 float rHitMinRight = rHitMin - l;
225 float rHitMaxRight = rHitMax - l;
227 float rHitMinLeft = rHitMin + l;
228 float rHitMaxLeft = rHitMax + l;
231 distRight[0][0] = rMin - rHitMinRight;
232 distRight[0][1] = rMin - rHitMaxRight;
233 distRight[1][0] = rMax - rHitMinRight;
234 distRight[1][1] = rMax - rHitMaxRight;
236 distLeft[0][0] = rMin - rHitMinLeft;
237 distLeft[0][1] = rMin - rHitMaxLeft;
238 distLeft[1][0] = rMax - rHitMinLeft;
239 distLeft[1][1] = rMax - rHitMaxLeft;
243 if (not sameSign(distRight[0][0], distRight[0][1], distRight[1][0], distRight[1][1])) {
248 if (not sameSign(distLeft[0][0], distLeft[0][1], distLeft[1][0], distLeft[1][1])) {
253 float rHitMinExtr = thetaVecMin.
cross(pos2D);
254 float rHitMaxExtr = thetaVecMax.
cross(pos2D);
255 if (rHitMinExtr * rHitMaxExtr < 0.)
return checkExtremum(node, wireHit);
265 long thetaMin = node->
getXMin();
266 long thetaMax = node->
getXMax();
271 float rMinD = thetaVecMin.
cross(pos2D);
272 float rMaxD = thetaVecMax.
cross(pos2D);
275 if ((rMinD > 0) && (rMaxD * rMinD >= 0))
return true;
276 if ((rMaxD * rMinD < 0))
return true;
282 const double& l = wireHit->getRefDriftLength();
287 long thetaMin = node->
getXMin();
288 long thetaMax = node->
getXMax();
293 if (not pos2D.
isBetween(thetaVecMin, thetaVecMax))
return false;
296 double r = pos2D.
norm();
297 float rRight = r - l;
301 float rMin = node->
getYMin() * r2 / 2;
302 float rMax = node->
getYMax() * r2 / 2;
304 bool crossesRight = (rMin - rRight) * (rMax - rRight) < 0;
305 bool crossesLeft = (rMin - rLeft) * (rMax - rLeft) < 0;
306 return crossesRight or crossesLeft;
311 static int nevent(0);
313 TCanvas* canv =
new TCanvas(
"canv",
"legendre transform", 0, 0, 1200, 600);
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);
325 const double& l = wireHit->getRefDriftLength();
327 double x = pos2D.
x();
328 double y = pos2D.
y();
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);
338 concaveHitLegendre->SetParameters(x, y, l, r2);
339 convexHitLegendre->SetParameters(x, y, l, r2);
340 concaveHitLegendre->Draw(
"CSAME");
341 convexHitLegendre->Draw(
"CSAME");
345 canv->Print(Form(
"legendreHits_%i.png", nevent));
353 std::vector<const CDCWireHit*> hits;
355 const CDCWireHit* wireHit = item->getPointer();
356 hits.push_back(wireHit);
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 checkDerivative(QuadTree *node, const CDCWireHit *wireHit) const
Check derivative of the sinogram.
void drawHits(std::vector< const CDCWireHit * > hits, unsigned int color=46) const
Draw QuadTree node.
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 precomputed cosine and sine values.
bool isLeaf(QuadTree *node) const final
lastLevel depends on curvature of the track candidate
bool checkExtremum(QuadTree *node, const 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.
bool isInNode(QuadTree *node, const CDCWireHit *wireHit) const final
Check whether hit belongs to the quadtree node:
static const LookupTable< 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.
Class representing a hit wire in the central drift chamber.
Class which holds precomputed values of a function.
int getNPoints() const
Return the number of finite sampling points in this lookup table.
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...
This class serves as a wrapper around all things that should go into a QuadTree.
Class which holds quadtree structure.
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.
This abstract class serves as a base class for all implementations of track processors.
std::pair< XSpan, YSpan > XYSpans
This pair of spans describes the span of a node.
int getLastLevel() const
Return the parameter last level.
typename QuadTree::YSpan YSpan
This pair describes the span in Y for a node.
std::unique_ptr< QuadTree > m_quadTree
The quad tree we work with.
A two dimensional vector which is equipped with functions for correct handling of orientation relate...
double dot(const Vector2D &rhs) const
Calculates the two dimensional dot product.
double cross(const Vector2D &rhs) const
Calculated the two dimensional cross product.
double x() const
Getter for the x coordinate.
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...
double normSquared() const
Calculates .
double y() const
Getter for the y coordinate.
double norm() const
Calculates the length of the vector.
static Vector2D Phi(const double phi)
Constructs a unit vector with azimuth angle equal to phi.
bool sameSign(double expected, double actual)
Predicate checking that two values have the same sign.
Abstract base class for different kinds of events.