All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TrackShowerSeparationAlg.h
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // 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.
10 ////////////////////////////////////////////////////////////////////////
11 
12 #ifndef TrackShowerSeparationAlg_hxx
13 #define TrackShowerSeparationAlg_hxx
14 
15 // Framework
16 #include "canvas/Persistency/Common/Ptr.h"
17 #include "canvas/Persistency/Common/FindManyP.h"
18 namespace fhicl { class ParameterSet; }
19 
20 // larsoft
24 
25 // ROOT
26 #include "TMath.h"
27 #include "TVector3.h"
28 
29 namespace shower {
30  class TrackShowerSeparationAlg;
31  class ReconTrack;
32 }
33 
34 #include <algorithm>
35 #include <iterator>
36 #include <map>
37 #include <vector>
38 
40  public:
41 
42  ReconTrack(int id) {
43  fID = id;
44  fTrack = false;
45  fShower = false;
46  fShowerTrack = false;
47  fShowerCone = false;
48  }
49 
50  // Setters
51  void SetVertex(TVector3 vertex) { fVertex = vertex; }
52  void SetEnd(TVector3 end) { fEnd = end; }
53  void SetLength(double length) { fLength = length; }
54  void SetVertexDir(TVector3 vertexDir) { fVertexDir = vertexDir; }
55  void SetDirection(TVector3 direction) { fDirection = direction; }
56  void SetHits(std::vector<art::Ptr<recob::Hit> > hits) { fHits = hits; }
57  void SetSpacePoints(std::vector<art::Ptr<recob::SpacePoint> > spacePoints) { fSpacePoints = spacePoints; }
58 
59  void AddForwardTrack(int track) { if (std::find(fForwardConeTracks.begin(), fForwardConeTracks.end(), track) == fForwardConeTracks.end()) fForwardConeTracks.push_back(track); }
60  void AddBackwardTrack(int track) { if (std::find(fBackwardConeTracks.begin(), fBackwardConeTracks.end(), track) == fBackwardConeTracks.end()) fBackwardConeTracks.push_back(track); }
61  void AddShowerTrack(int track) { fShowerTracks.push_back(track); }
62 
63  void AddForwardSpacePoint(int spacePoint) { fForwardSpacePoints.push_back(spacePoint); }
64  void AddBackwardSpacePoint(int spacePoint) { fBackwardSpacePoints.push_back(spacePoint); }
65  void AddCylinderSpacePoint(int spacePoint) { fCylinderSpacePoints.push_back(spacePoint); }
66  void AddSphereSpacePoint(int spacePoint) { fSphereSpacePoints.push_back(spacePoint); }
67  void AddIsolationSpacePoint(int spacePoint, double distance) { fIsolationSpacePoints[spacePoint] = distance; }
68 
69  // Getters
70  int ID() const { return fID; }
71  TVector3 Vertex() const { return fVertex; }
72  TVector3 End() const { return fEnd; }
73  double Length() const { return fLength; }
74  TVector3 VertexDirection() const { return fVertexDir; }
75  TVector3 Direction() const { return fDirection; }
76  const std::vector<art::Ptr<recob::Hit> >& Hits() const { return fHits; }
77  const std::vector<art::Ptr<recob::SpacePoint> >& SpacePoints() const { return fSpacePoints; }
78 
79  void FlipTrack() {
80  TVector3 tmp = fEnd;
81  fEnd = fVertex;
82  fVertex = tmp;
83  fDirection *= -1;
84  }
85 
86  void MakeShower() {
87  if (fTrack)
88  this->MakeShowerTrack();
89  else
90  this->MakeShowerCone();
91  }
92  void MakeShowerTrack() {
93  fShower = true;
94  fShowerTrack = true;
95  fShowerCone = false;
96  fTrack = false;
97  }
98  void MakeShowerCone() {
99  fShower = true;
100  fShowerCone = true;
101  fShowerTrack = false;
102  fTrack = false;
103  }
104  void MakeTrack() {
105  fTrack = true;
106  fShower = false;
107  fShowerTrack = false;
108  fShowerCone = false;
109  }
110 
111  bool IsShower() const { return fShower; }
112  bool IsShowerTrack() const { return fShowerTrack; }
113  bool IsShowerCone() const { return fShowerCone; }
114  bool IsTrack() const { return fTrack; }
115  bool IsUndetermined() const { return !fTrack and !fShower; }
116 
117  int TrackConeSize() const { return (int)fForwardConeTracks.size() - (int)fBackwardConeTracks.size(); }
118  bool ShowerTrackCandidate() const { return TrackConeSize() > 5; }
119  const std::vector<int>& ShowerTracks() const { return fShowerTracks; }
120  const std::vector<int>& ForwardConeTracks() const { return fForwardConeTracks; }
121 
122  int ConeSize() const { return (int)fForwardSpacePoints.size() - (int)fBackwardSpacePoints.size(); }
123  int ForwardSpacePoints() const { return fForwardSpacePoints.size(); }
124  int NumCylinderSpacePoints() const { return fCylinderSpacePoints.size(); }
125  double CylinderSpacePointRatio() const { return (double)fCylinderSpacePoints.size()/(double)fSpacePoints.size(); }
126  int NumSphereSpacePoints() const { return fSphereSpacePoints.size(); }
127  double SphereSpacePointDensity(double scale) const { return (double)fSphereSpacePoints.size()/(4*TMath::Pi()*TMath::Power((scale*fLength/2.),3)/3.); }
128  double IsolationSpacePointDistance() const { std::vector<double> distances;
129  std::transform(fIsolationSpacePoints.begin(), fIsolationSpacePoints.end(), std::back_inserter(distances), [](const std::pair<int,double>& p){return p.second;});
130  return TMath::Mean(distances.begin(), distances.end()); }
131 
132  private:
133 
134  int fID;
135  TVector3 fVertex;
136  TVector3 fEnd;
137  double fLength;
138  TVector3 fVertexDir;
139  TVector3 fDirection;
140  std::vector<art::Ptr<recob::Hit> > fHits;
141  std::vector<art::Ptr<recob::SpacePoint> > fSpacePoints;
142 
143  std::vector<int> fForwardConeTracks;
144  std::vector<int> fBackwardConeTracks;
145  std::vector<int> fShowerTracks;
146 
147  std::vector<int> fForwardSpacePoints;
148  std::vector<int> fBackwardSpacePoints;
149  std::vector<int> fCylinderSpacePoints;
150  std::vector<int> fSphereSpacePoints;
151  std::map<int,double> fIsolationSpacePoints;
152 
153  bool fShower;
156  bool fTrack;
157 
158 };
159 
161  public:
162 
163  TrackShowerSeparationAlg(fhicl::ParameterSet const& pset);
164 
165  /// Read in configurable parameters from provided parameter set
166  void reconfigure(fhicl::ParameterSet const& pset);
167 
168  std::vector<art::Ptr<recob::Hit> > SelectShowerHits(int event,
169  const std::vector<art::Ptr<recob::Hit> >& hits,
170  const std::vector<art::Ptr<recob::Track> >& tracks,
171  const std::vector<art::Ptr<recob::SpacePoint> >& spacePoints,
172  const art::FindManyP<recob::Hit>& fmht,
173  const art::FindManyP<recob::Track>& fmth,
174  const art::FindManyP<recob::SpacePoint>& fmspt,
175  const art::FindManyP<recob::Track>& fmtsp) const;
176 
177  private:
178 
179  ///
180  std::vector<int> InitialTrackLikeSegment(std::map<int,std::unique_ptr<ReconTrack> >& reconTracks) const;
181 
182  ///
183  TVector3 Gradient(const std::vector<TVector3>& points, const std::unique_ptr<TVector3>& dir) const;
184 
185  ///
186  TVector3 Gradient(const art::Ptr<recob::Track>& track) const;
187 
188  ///
189  TVector3 Gradient(const std::vector<art::Ptr<recob::SpacePoint> >& spacePoints) const;
190 
191  /// Projects the 3D point given along the line defined by the specified direction
192  /// Coordinate system is assumed to be centred at the origin unless a difference point is specified
193  TVector3 ProjPoint(const TVector3& point, const TVector3& direction, const TVector3& origin = TVector3(0,0,0)) const;
194 
195  /// Return 3D point of this space point
196  TVector3 SpacePointPos(const art::Ptr<recob::SpacePoint>& spacePoint) const;
197 
198  ///
199  double SpacePointsRMS(const std::vector<art::Ptr<recob::SpacePoint> >& spacePoints) const;
200 
201  // Parameters
202  int fDebug;
203 
204  // Properties
205  double fConeAngle;
207 
208  // Cuts
210  double fCylinderCut;
212 
213 };
214 
215 #endif
process_name vertex
Definition: cheaterreco.fcl:51
void SetDirection(TVector3 direction)
const std::vector< art::Ptr< recob::SpacePoint > > & SpacePoints() const
std::vector< int > InitialTrackLikeSegment(std::map< int, std::unique_ptr< ReconTrack > > &reconTracks) const
void SetEnd(TVector3 end)
const std::vector< int > & ForwardConeTracks() const
ClusterModuleLabel join with tracks
void reconfigure(fhicl::ParameterSet const &pset)
Read in configurable parameters from provided parameter set.
static constexpr Sample_t transform(Sample_t sample)
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
pdgs p
Definition: selectors.fcl:22
std::map< int, double > fIsolationSpacePoints
TVector3 ProjPoint(const TVector3 &point, const TVector3 &direction, const TVector3 &origin=TVector3(0, 0, 0)) const
process_name use argoneut_mc_hitfinder track
void SetSpacePoints(std::vector< art::Ptr< recob::SpacePoint > > spacePoints)
void AddForwardTrack(int track)
process_name shower
Definition: cheaterreco.fcl:51
std::vector< int > fBackwardSpacePoints
TVector3 SpacePointPos(const art::Ptr< recob::SpacePoint > &spacePoint) const
Return 3D point of this space point.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
void SetVertex(TVector3 vertex)
double SpacePointsRMS(const std::vector< art::Ptr< recob::SpacePoint > > &spacePoints) const
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
std::vector< int > fForwardConeTracks
TrackShowerSeparationAlg(fhicl::ParameterSet const &pset)
std::vector< art::Ptr< recob::Hit > > SelectShowerHits(int event, const std::vector< art::Ptr< recob::Hit > > &hits, const std::vector< art::Ptr< recob::Track > > &tracks, const std::vector< art::Ptr< recob::SpacePoint > > &spacePoints, const art::FindManyP< recob::Hit > &fmht, const art::FindManyP< recob::Track > &fmth, const art::FindManyP< recob::SpacePoint > &fmspt, const art::FindManyP< recob::Track > &fmtsp) const
std::vector< art::Ptr< recob::SpacePoint > > fSpacePoints
const std::vector< art::Ptr< recob::Hit > > & Hits() const
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
return match has_match and(match.match_pdg==11 or match.match_pdg==-11)
std::vector< int > fSphereSpacePoints
void SetVertexDir(TVector3 vertexDir)
Provides recob::Track data product.
tuple dir
Definition: dropbox.py:28
void AddCylinderSpacePoint(int spacePoint)
std::vector< int > fShowerTracks
std::vector< int > fForwardSpacePoints
TVector3 VertexDirection() const
double SphereSpacePointDensity(double scale) const
std::vector< art::Ptr< recob::Hit > > fHits
TVector3 Gradient(const std::vector< TVector3 > &points, const std::unique_ptr< TVector3 > &dir) const
void AddSphereSpacePoint(int spacePoint)
void AddIsolationSpacePoint(int spacePoint, double distance)
void SetHits(std::vector< art::Ptr< recob::Hit > > hits)
TVector3 Direction() const
void AddShowerTrack(int track)
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:227
std::vector< int > fBackwardConeTracks