LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
TrackShowerSeparationAlg.h
Go to the documentation of this file.
1 // Class: TrackShowerSeparationAlg
3 // File: TrackShowerSeparationAlg.h
4 // Author: Mike Wallbank (m.wallbank@sheffield.ac.uk), November 2015
5 //
6 // Track/shower separation class.
7 // Provides methods for removing hits associated with track-like
8 // objects.
9 // To be run after track reconstruction, before shower reconstruction.
11 
12 #ifndef TrackShowerSeparationAlg_hxx
13 #define TrackShowerSeparationAlg_hxx
14 
15 // Framework
18 namespace fhicl {
19  class ParameterSet;
20 }
21 
22 // larsoft
26 
27 // ROOT
28 #include "TMath.h"
29 #include "TVector3.h"
30 
31 namespace shower {
32  class TrackShowerSeparationAlg;
33  class ReconTrack;
34 }
35 
36 #include <algorithm>
37 #include <iterator>
38 #include <map>
39 #include <vector>
40 
42 public:
43  ReconTrack(int id)
44  {
45  fID = id;
46  fTrack = false;
47  fShower = false;
48  fShowerTrack = false;
49  fShowerCone = false;
50  }
51 
52  // Setters
53  void SetVertex(TVector3 vertex) { fVertex = vertex; }
54  void SetEnd(TVector3 end) { fEnd = end; }
55  void SetLength(double length) { fLength = length; }
56  void SetVertexDir(TVector3 vertexDir) { fVertexDir = vertexDir; }
57  void SetDirection(TVector3 direction) { fDirection = direction; }
60  {
61  fSpacePoints = spacePoints;
62  }
63 
65  {
66  if (std::find(fForwardConeTracks.begin(), fForwardConeTracks.end(), track) ==
67  fForwardConeTracks.end())
68  fForwardConeTracks.push_back(track);
69  }
71  {
72  if (std::find(fBackwardConeTracks.begin(), fBackwardConeTracks.end(), track) ==
73  fBackwardConeTracks.end())
74  fBackwardConeTracks.push_back(track);
75  }
76  void AddShowerTrack(int track) { fShowerTracks.push_back(track); }
77 
78  void AddForwardSpacePoint(int spacePoint) { fForwardSpacePoints.push_back(spacePoint); }
79  void AddBackwardSpacePoint(int spacePoint) { fBackwardSpacePoints.push_back(spacePoint); }
80  void AddCylinderSpacePoint(int spacePoint) { fCylinderSpacePoints.push_back(spacePoint); }
81  void AddSphereSpacePoint(int spacePoint) { fSphereSpacePoints.push_back(spacePoint); }
82  void AddIsolationSpacePoint(int spacePoint, double distance)
83  {
84  fIsolationSpacePoints[spacePoint] = distance;
85  }
86 
87  // Getters
88  int ID() const { return fID; }
89  TVector3 Vertex() const { return fVertex; }
90  TVector3 End() const { return fEnd; }
91  double Length() const { return fLength; }
92  TVector3 VertexDirection() const { return fVertexDir; }
93  TVector3 Direction() const { return fDirection; }
94  const std::vector<art::Ptr<recob::Hit>>& Hits() const { return fHits; }
95  const std::vector<art::Ptr<recob::SpacePoint>>& SpacePoints() const { return fSpacePoints; }
96 
97  void FlipTrack()
98  {
99  TVector3 tmp = fEnd;
100  fEnd = fVertex;
101  fVertex = tmp;
102  fDirection *= -1;
103  }
104 
105  void MakeShower()
106  {
107  if (fTrack)
108  this->MakeShowerTrack();
109  else
110  this->MakeShowerCone();
111  }
113  {
114  fShower = true;
115  fShowerTrack = true;
116  fShowerCone = false;
117  fTrack = false;
118  }
120  {
121  fShower = true;
122  fShowerCone = true;
123  fShowerTrack = false;
124  fTrack = false;
125  }
126  void MakeTrack()
127  {
128  fTrack = true;
129  fShower = false;
130  fShowerTrack = false;
131  fShowerCone = false;
132  }
133 
134  bool IsShower() const { return fShower; }
135  bool IsShowerTrack() const { return fShowerTrack; }
136  bool IsShowerCone() const { return fShowerCone; }
137  bool IsTrack() const { return fTrack; }
138  bool IsUndetermined() const { return !fTrack and !fShower; }
139 
140  int TrackConeSize() const
141  {
142  return (int)fForwardConeTracks.size() - (int)fBackwardConeTracks.size();
143  }
144  bool ShowerTrackCandidate() const { return TrackConeSize() > 5; }
145  const std::vector<int>& ShowerTracks() const { return fShowerTracks; }
146  const std::vector<int>& ForwardConeTracks() const { return fForwardConeTracks; }
147 
148  int ConeSize() const
149  {
150  return (int)fForwardSpacePoints.size() - (int)fBackwardSpacePoints.size();
151  }
152  int ForwardSpacePoints() const { return fForwardSpacePoints.size(); }
153  int NumCylinderSpacePoints() const { return fCylinderSpacePoints.size(); }
154  double CylinderSpacePointRatio() const
155  {
156  return (double)fCylinderSpacePoints.size() / (double)fSpacePoints.size();
157  }
158  int NumSphereSpacePoints() const { return fSphereSpacePoints.size(); }
159  double SphereSpacePointDensity(double scale) const
160  {
161  return (double)fSphereSpacePoints.size() /
162  (4 * TMath::Pi() * TMath::Power((scale * fLength / 2.), 3) / 3.);
163  }
165  {
166  std::vector<double> distances;
167  std::transform(fIsolationSpacePoints.begin(),
168  fIsolationSpacePoints.end(),
169  std::back_inserter(distances),
170  [](const std::pair<int, double>& p) { return p.second; });
171  return TMath::Mean(distances.begin(), distances.end());
172  }
173 
174 private:
175  int fID;
176  TVector3 fVertex;
177  TVector3 fEnd;
178  double fLength;
179  TVector3 fVertexDir;
180  TVector3 fDirection;
181  std::vector<art::Ptr<recob::Hit>> fHits;
182  std::vector<art::Ptr<recob::SpacePoint>> fSpacePoints;
183 
184  std::vector<int> fForwardConeTracks;
185  std::vector<int> fBackwardConeTracks;
186  std::vector<int> fShowerTracks;
187 
188  std::vector<int> fForwardSpacePoints;
189  std::vector<int> fBackwardSpacePoints;
190  std::vector<int> fCylinderSpacePoints;
191  std::vector<int> fSphereSpacePoints;
192  std::map<int, double> fIsolationSpacePoints;
193 
194  bool fShower;
197  bool fTrack;
198 };
199 
201 public:
203 
205  void reconfigure(fhicl::ParameterSet const& pset);
206 
207  std::vector<art::Ptr<recob::Hit>> SelectShowerHits(
208  int event,
210  const std::vector<art::Ptr<recob::Track>>& tracks,
211  const std::vector<art::Ptr<recob::SpacePoint>>& spacePoints,
212  const art::FindManyP<recob::Hit>& fmht,
213  const art::FindManyP<recob::Track>& fmth,
215  const art::FindManyP<recob::Track>& fmtsp) const;
216 
217 private:
219  std::vector<int> InitialTrackLikeSegment(
220  std::map<int, std::unique_ptr<ReconTrack>>& reconTracks) const;
221 
223  TVector3 Gradient(const std::vector<TVector3>& points,
224  const std::unique_ptr<TVector3>& dir) const;
225 
227  TVector3 Gradient(const art::Ptr<recob::Track>& track) const;
228 
230  TVector3 Gradient(const std::vector<art::Ptr<recob::SpacePoint>>& spacePoints) const;
231 
234  TVector3 ProjPoint(const TVector3& point,
235  const TVector3& direction,
236  const TVector3& origin = TVector3(0, 0, 0)) const;
237 
239  TVector3 SpacePointPos(const art::Ptr<recob::SpacePoint>& spacePoint) const;
240 
242  double SpacePointsRMS(const std::vector<art::Ptr<recob::SpacePoint>>& spacePoints) const;
243 
244  // Parameters
245  int fDebug;
246 
247  // Properties
248  double fConeAngle;
250 
251  // Cuts
253  double fCylinderCut;
255 };
256 
257 #endif
void SetDirection(TVector3 direction)
const std::vector< art::Ptr< recob::SpacePoint > > & SpacePoints() const
void SetEnd(TVector3 end)
const std::vector< int > & ForwardConeTracks() const
double CylinderSpacePointRatio() const
const std::vector< int > & ShowerTracks() const
void AddForwardSpacePoint(int spacePoint)
Declaration of signal hit object.
void SetLength(double length)
void AddBackwardTrack(int track)
double IsolationSpacePointDistance() const
Float_t tmp
Definition: plot.C:35
void AddForwardTrack(int track)
std::vector< int > fBackwardSpacePoints
std::map< int, double > fIsolationSpacePoints
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
void SetSpacePoints(std::vector< art::Ptr< recob::SpacePoint >> spacePoints)
Double_t scale
Definition: plot.C:24
void SetVertex(TVector3 vertex)
void hits()
Definition: readHits.C:15
std::vector< int > fForwardConeTracks
parameter set interface
std::vector< art::Ptr< recob::SpacePoint > > fSpacePoints
const std::vector< art::Ptr< recob::Hit > > & Hits() const
Provides recob::Track data product.
std::vector< int > fSphereSpacePoints
void SetVertexDir(TVector3 vertexDir)
void AddCylinderSpacePoint(int spacePoint)
std::vector< int > fShowerTracks
std::vector< int > fForwardSpacePoints
TVector3 VertexDirection() const
TDirectory * dir
Definition: macro.C:5
double SphereSpacePointDensity(double scale) const
std::vector< art::Ptr< recob::Hit > > fHits
void AddSphereSpacePoint(int spacePoint)
void AddIsolationSpacePoint(int spacePoint, double distance)
TVector3 Direction() const
void AddShowerTrack(int track)
Float_t track
Definition: plot.C:35
void SetHits(std::vector< art::Ptr< recob::Hit >> hits)
void AddBackwardSpacePoint(int spacePoint)
std::vector< int > fCylinderSpacePoints
constexpr Point origin()
Returns a origin position with a point of the specified type.
Definition: geo_vectors.h:229
Event finding and building.
std::vector< int > fBackwardConeTracks
vertex reconstruction