LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
Utilities.h
Go to the documentation of this file.
1 
12 #ifndef Utilities_h
13 #define Utilities_h
14 
17 
18 #include <functional>
19 
20 #include "TVector2.h"
21 #include "TVector3.h"
22 
23 #include "Math/GenVector/DisplacementVector2D.h"
24 #include "Math/GenVector/DisplacementVector3D.h"
25 
26 namespace pma
27 {
28  typedef ROOT::Math::DisplacementVector2D< ROOT::Math::Cartesian2D<double> > Vector2D;
30 
31  typedef std::map< size_t, std::vector<double> > dedx_map;
32 
33  class Hit3D;
34  class TrkCandidate;
35  class bSegmentProjLess;
36  class bDistCenterLess2D;
37  class bDistCenterLess3D;
39  struct bTrajectory3DDistLess;
40  struct bTrack3DLonger;
41 
42  double Dist2(const TVector2 & v1, const TVector2 & v2);
43  double Dist2(const Vector2D & v1, const Vector2D & v2);
44 
45  template <typename T, typename U> double Dist2(const T & v1, const U & v2) {
46  double dx = v1.X() - v2.X(), dy = v1.Y() - v2.Y(), dz = v1.Z() - v2.Z();
47  return dx * dx + dy * dy + dz * dz;
48  }
49 
50  size_t GetHitsCount(const std::vector< pma::Hit3D* >& hits, unsigned int view);
51  double GetSummedADC(const std::vector< pma::Hit3D* >& hits, unsigned int view = geo::kUnknown);
52  double GetSummedAmpl(const std::vector< pma::Hit3D* >& hits, unsigned int view = geo::kUnknown);
53 
54  double GetHitsRadius3D(const std::vector< pma::Hit3D* >& hits, bool exact = false);
55  double GetHitsRadius2D(const std::vector< pma::Hit3D* >& hits, bool exact = false);
56 
57  double GetSegmentProjVector(const TVector2& p, const TVector2& p0, const TVector2& p1);
58  double GetSegmentProjVector(const Vector2D& p, const Vector2D& p0, const Vector2D& p1);
59 
60  double GetSegmentProjVector(const TVector3& p, const TVector3& p0, const TVector3& p1);
61  double GetSegmentProjVector(const Vector3D& p, const Vector3D& p0, const Vector3D& p1);
62 
63  TVector2 GetProjectionToSegment(const TVector2& p, const TVector2& p0, const TVector2& p1);
64  TVector3 GetProjectionToSegment(const TVector3& p, const TVector3& p0, const TVector3& p1);
65 
66  double SolveLeastSquares3D(const std::vector< std::pair<TVector3, TVector3> >& lines, TVector3& result);
67 
68  TVector2 GetProjectionToPlane(const TVector3& p, unsigned int plane, unsigned int tpc, unsigned int cryo);
69  TVector2 GetVectorProjectionToPlane(const TVector3& v, unsigned int plane, unsigned int tpc, unsigned int cryo);
70  TVector2 WireDriftToCm(unsigned int wire, float drift, unsigned int plane, unsigned int tpc, unsigned int cryo);
71  TVector2 CmToWireDrift(float xw, float yd, unsigned int plane, unsigned int tpc, unsigned int cryo);
72 }
73 
74 
76  public std::binary_function<pma::Hit3D*, pma::Hit3D*, bool>
77 {
79 };
80 
82  public std::binary_function<pma::Hit3D*, pma::Hit3D*, bool>
83 {
85 };
86 
88  public std::binary_function<const pma::TrkCandidate &, const pma::TrkCandidate &, bool>
89 {
90  bool operator() (const pma::TrkCandidate & t1, const pma::TrkCandidate & t2);
91 };
92 
94  public std::binary_function<TVector3*, TVector3*, bool>
95 {
96 public:
97  bSegmentProjLess(const TVector3& s0, const TVector3& s1);
98 
99  bool operator() (TVector3* p1, TVector3* p2)
100  {
101  if (p1 && p2)
102  {
103  double b1 = pma::GetSegmentProjVector(*p1, segStart, segStop);
104  double b2 = pma::GetSegmentProjVector(*p1, segStart, segStop);
105  return b1 < b2;
106  }
107  else return false;
108  }
109 
110 private:
111  TVector3 segStart, segStop;
112 };
113 
115  public std::binary_function<TVector2, TVector2, bool>
116 {
117 public:
118  bDistCenterLess2D(const TVector2& c) : center(c) {}
119 
120  bool operator() (TVector2 p1, TVector2 p2)
121  {
122  double b1 = pma::Dist2(p1, center);
123  double b2 = pma::Dist2(p2, center);
124  return b1 < b2;
125  }
126 
127 private:
128  TVector2 center;
129 };
130 
132  public std::binary_function<TVector3, TVector3, bool>
133 {
134 public:
135  bDistCenterLess3D(const TVector3& c) : center(c) {}
136 
137  bool operator() (TVector3 p1, TVector3 p2)
138  {
139  double b1 = pma::Dist2(p1, center);
140  double b2 = pma::Dist2(p2, center);
141  return b1 < b2;
142  }
143 
144 private:
145  TVector3 center;
146 };
147 
148 #endif
149 
double GetSummedAmpl(const std::vector< pma::Hit3D * > &hits, unsigned int view=geo::kUnknown)
Definition: Utilities.cxx:47
TTree * t1
Definition: plottest35.C:26
double GetHitsRadius2D(const std::vector< pma::Hit3D * > &hits, bool exact=false)
Definition: Utilities.cxx:78
double Dist2(const TVector2 &v1, const TVector2 &v2)
Definition: Utilities.cxx:19
Unknown view.
Definition: geo_types.h:83
bDistCenterLess3D(const TVector3 &c)
Definition: Utilities.h:135
double GetHitsRadius3D(const std::vector< pma::Hit3D * > &hits, bool exact=false)
Definition: Utilities.cxx:56
recob::tracking::Vector_t Vector3D
Definition: Utilities.h:29
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< Coord_t >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space. See recob::tracking::Coord_t for more details on the ...
Definition: TrackingTypes.h:29
void hits()
Definition: readHits.C:15
size_t GetHitsCount(const std::vector< pma::Hit3D * > &hits, unsigned int view)
Definition: Utilities.cxx:30
TVector2 GetProjectionToSegment(const TVector2 &p, const TVector2 &p0, const TVector2 &p1)
Definition: Utilities.cxx:156
ROOT::Math::DisplacementVector2D< ROOT::Math::Cartesian2D< double > > Vector2D
Definition: Utilities.h:28
bDistCenterLess2D(const TVector2 &c)
Definition: Utilities.h:118
TVector2 GetVectorProjectionToPlane(const TVector3 &v, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:287
double SolveLeastSquares3D(const std::vector< std::pair< TVector3, TVector3 > > &lines, TVector3 &result)
Definition: Utilities.cxx:176
TTree * t2
Definition: plottest35.C:36
TH1F * h2
Definition: plot.C:46
bool operator()(pma::Hit3D *h1, pma::Hit3D *h2)
Definition: Utilities.cxx:318
double GetSegmentProjVector(const TVector2 &p, const TVector2 &p0, const TVector2 &p1)
Definition: Utilities.cxx:100
TVector2 WireDriftToCm(unsigned int wire, float drift, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:296
TVector2 GetProjectionToPlane(const TVector3 &p, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:280
TH1F * h1
Definition: plot.C:43
TVector2 CmToWireDrift(float xw, float yd, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:307
std::map< size_t, std::vector< double > > dedx_map
Definition: Utilities.h:31
double GetSummedADC(const std::vector< pma::Hit3D * > &hits, unsigned int view=geo::kUnknown)
Definition: Utilities.cxx:38
art framework interface to geometry description