TrdMCClusterR Fit Classifier
AMSTracker.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <limits>
4 #include <vector>
5 
6 #include <TH1F.h>
7 #include <TH2D.h>
8 
9 #include "Line.h"
10 #include "Point.h"
11 
12 const double TrkL1_ZPos = 158.920;
13 const double TrkL1_ZPosCutoff = 2;
14 const double closest_distance_cutoff = 10;
15 
16 
17 namespace TrackerL1
18 {
19  //___________________________________________________________________________
20  bool HasL1XYHit(const TH2D * hTrackerL1_XY, const double& x, const double& y, const bool& debug = false)
21  {
23  int binX = hTrackerL1_XY->GetXaxis()->FindBin(x);
24  int binY = hTrackerL1_XY->GetYaxis()->FindBin(y);
25 
26  if (debug == true)
27  printf("\n\nFor (%f, %f): (binX, binY) = (%d, %d). Bin content = %f\n\n", x, y, binX, binY, hTrackerL1_XY->GetBinContent(binX, binY));
28 
29  return (hTrackerL1_XY->GetBinContent(binX, binY) > 0);
30  }
31 
32 
33  //___________________________________________________________________________
34  std::pair<Point::Point2DwID, bool> FindClosestTrkL1ClusterFrom3DLine(const TH2D * hTrackerL1_XY,
35  const Line::Line3DParametricForm& line,
36  const std::vector<Point::Point2DwID>& gMCClustersXY,
37  const bool& debug = false)
38  {
40  if (!hTrackerL1_XY)
41  {
42  std::cerr<<"\n\nWARNING! Bad pointer for hTrackerL1_XY\n\n";
43  return std::make_pair(Point::Point2DwID(), false);
44  }
45  if (line.direction[2] == 0)
46  return std::make_pair(Point::Point2DwID(), false);
47 
49  double t = (TrkL1_ZPos - line.point[2]) / line.direction[2];
50 
52  double x_extrapolated = line.point[0] + line.direction[0] * t;
53  double y_extrapolated = line.point[1] + line.direction[1] * t;
54 
55  if (debug == true)
56  std::cout<<"\nThe extrapolated X and Y at L1 for the 3D line are: "<<x_extrapolated<<", "<<y_extrapolated;
57 
59  bool has_L1_hit = HasL1XYHit(hTrackerL1_XY, x_extrapolated, y_extrapolated);
60  if (!has_L1_hit)
61  return std::make_pair(Point::Point2DwID(), false);
62 
64  Point::Point2DwID closestPoint;
65  double minDistance = std::numeric_limits<double>::max();
66 
67  for (const auto& MCCluster : gMCClustersXY)
68  {
69  double distance = Point::GetDistance2D(x_extrapolated, y_extrapolated, MCCluster.x, MCCluster.y);
70 
71  if (debug == true)
72  printf("\nDistance from (%7.2f, %7.2f) = %7.2f cm.", MCCluster.x, MCCluster.y, distance);
73 
74  if (distance < minDistance)
75  {
76  minDistance = distance;
77  closestPoint = MCCluster;
78  }
79  }
80 
81 
87 
88  if (minDistance <= closest_distance_cutoff)
89  {
90  closestPoint.closedis = minDistance;
91  return std::make_pair(closestPoint, true);
92  }
93 
94  else
95  return std::make_pair(Point::Point2DwID(), false);
96  }
97 }
const double TrkL1_ZPos
Definition: AMSTracker.h:12
const double closest_distance_cutoff
cm
Definition: AMSTracker.h:14
const double TrkL1_ZPosCutoff
cm
Definition: AMSTracker.h:13
double GetDistance2D(double x1, double y1, double x2, double y2)
Definition: Point.h:43
cm
Definition: AMSTracker.h:18
bool HasL1XYHit(const TH2D *hTrackerL1_XY, const double &x, const double &y, const bool &debug=false)
Definition: AMSTracker.h:20
std::pair< Point::Point2DwID, bool > FindClosestTrkL1ClusterFrom3DLine(const TH2D *hTrackerL1_XY, const Line::Line3DParametricForm &line, const std::vector< Point::Point2DwID > &gMCClustersXY, const bool &debug=false)
Definition: AMSTracker.h:34
Definition: Line.h:60
std::vector< double > point
Definition: Line.h:63
std::vector< double > direction
Vector containing the point on the line (x0, y0, z0)
Definition: Line.h:64
Definition: Point.h:18
double closedis
Definition: Point.h:25