All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ClusterMergeAlg.h
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // \file ClusterMergeAlg.h
3 //
4 // \brief ClusterMergeAlg header file
5 //
6 // \author david.kaleko@gmail.com, kazuhiro@nevis.columbia.edu
7 //
8 ////////////////////////////////////////////////////////////////////////
9 
10 #ifndef CLUSTERMERGEALG_H
11 #define CLUSTERMERGEALG_H
12 
13 // ART includes
14 namespace fhicl {
15  class ParameterSet;
16 }
17 
18 #include "canvas/Persistency/Common/Ptr.h"
19 
20 // LArSoft
24 
25 // STL
26 #include <vector>
27 
28 class TTree;
29 
30 namespace cluster
31 {
32 
33  /**
34  \struct cluster_merge_info
35  A utility struct for cluster-wise analysis information for merging
36  */
38 
39  unsigned int cluster_index; ///< Input cluster ID
40  geo::View_t view; ///< Wire plane ID
41  geo::PlaneID planeID; ///< plane ID
42 
43  float start_wire; ///< Vertex wire
44  float start_time; ///< Vertex time
45  float end_wire; ///< End point wire
46  float end_time; ///< End point time
47 
48  double start_wire_err; ///< Vertex wire error
49  double start_time_err; ///< Vertex time error
50  double end_wire_err; ///< End point wire error
51  double end_time_err; ///< End point time error
52 
53  float angle; ///< 2D starting angle (in radians)
54 
55  /// Default constructor
57 
58  cluster_index = 0xffffffff;
62 
63  };
64 
65  /// Initialization from a recob::Cluster
66  explicit cluster_merge_info(const recob::Cluster& cl)
67  :cluster_index(cl.ID())
68  ,view(cl.View())
69  ,planeID(cl.Plane())
70  ,start_wire(cl.StartWire())
71  ,start_time(cl.StartTick())
72  ,end_wire(cl.EndWire())
73  ,end_time(cl.EndTick())
74  ,angle(cl.StartAngle())
75  {}
76  };
77 
79 
80  public:
81 
82  /// Default constructor with fhicl parameters
83  ClusterMergeAlg(fhicl::ParameterSet const& pset);
84 
85  /// Method to set verbose mode
86  void VerboseMode(bool on) { _verbose = on; }
87 
88  /// Method to report the current configuration
89  void ReportConfig() const;
90 
91  /// Method to set cut value in degrees for angle compatibility test
93 
94  /// Method to set cut value in cm^2 for distance compatibility test
95  void SetSquaredDistanceCut(double d) { _max_2D_dist2 = d; }
96 
97  /// Method to add a cluster information for processing
98  void AppendClusterInfo(const recob::Cluster &in_cluster,
99  const std::vector<art::Ptr<recob::Hit> > &in_hit_v);
100 
101  /// Method to add a cluster information for processing
102  void AppendClusterInfo(const art::Ptr<recob::Cluster> in_cluster,
103  const std::vector<art::Ptr<recob::Hit> > &in_hit_v);
104 
105  /// Method to clear event-wise information (both input cluster info & output merged cluster sets)
106  void ClearEventInfo();
107 
108  /**
109  Method to execute the algorithm. All parameter configuration + adding input cluster information
110  should be done before calling this function. This function generate a resulting set of cluster IDs
111  for merging, which can be accessed through GetClusterSets() function.
112  */
113  void ProcessMergeAlg();
114 
115  /// Method to extract resulting set of cluster IDs for merging computed by ProcessMergeAlg() function.
116  const std::vector<std::vector<unsigned int> > GetClusterSets () const {return _cluster_sets_v;};
117 
118  /// Method to compare a compatibility between two clusters
119  bool CompareClusters(const cluster_merge_info &clus_info_A,
120  const cluster_merge_info &clus_info_B);
121 
122  /**
123  Function to compare the 2D angles of two clusters and return true if they are
124  within the maximum allowed parameter. Includes shifting by 180 for backwards clusters.
125  This is called within CompareClusters().
126  */
127  bool Angle2DCompatibility(const cluster_merge_info &clus_info_A,
128  const cluster_merge_info &clus_info_B) const;
129 
130  /**
131  Function to compare the 2D distance of two clusters and return true if they are
132  within the maximum allowed distance.The distance-squared is computed by another
133  function, ShortestDistanceSquared().
134  This is called within CompareClusters().
135  */
136  bool ShortestDistanceCompatibility(const cluster_merge_info &clus_info_A,
137  const cluster_merge_info &clus_info_B) const;
138 
139  /**
140  Function to compute a distance between a 2D point (point_x, point_y) to a 2D finite line segment
141  (start_x, start_y) => (end_x, end_y).
142  */
143  double ShortestDistanceSquared(double point_x, double point_y,
144  double start_x, double start_y,
145  double end_x, double end_y ) const;
146 
147  /**
148  Function to print to screen a specific cluser's info
149  from ClusterPrepAna module. Used for debugging.
150  */
151  void PrintClusterVars(cluster_merge_info clus_info) const;
152 
153  /**
154  Utility function to check if an index is already somewhere inside of _cluster_sets_v vector
155  returns the location of the element vector in _cluster_sets_v that contains the index
156  and returns -1 if the index is not in _cluster_sets_v anywhere
157  */
158  int isInClusterSets(unsigned int index) const;
159 
160  protected:
161 
162  /// Method to set a conversion factor from wire to cm scale
163  void SetWire2Cm(double f) { _wire_2_cm = f; }
164 
165  /// Method to set a conversion factor from time to cm scale
166  void SetTime2Cm(double f) { _time_2_cm = f; }
167 
168  /// Method to clear output merged cluster sets (_cluster_sets_v)
169  void ClearOutputInfo();
170 
171  /// Method to clear input cluster information
172  void ClearInputInfo();
173 
174  /// Method to clear quality control TTree variables
175  void ClearTTreeInfo();
176 
177  /// Method to prepare quality control TTree
178  void PrepareTTree();
179 
180  /// Method to prepare detector parameters
181  void PrepareDetParams();
182 
183  /// Method to fill hit-array-related information
185  const std::vector<art::Ptr<recob::Hit> > &in_hit_v);
186 
187  /**
188  For a given pair of clusters, this function calls CompareClusters() and append to the resulting
189  merged cluster sets (_cluster_sets_v) by calling AppendToClusterSets() when they are compatible.
190  */
191  void BuildClusterSets(const cluster_merge_info &clus_info_A,
192  const cluster_merge_info &clus_info_B);
193 
194  /**
195  Function to loop through _cluster_sets_v and add in the un-mergable clusters
196  individually, because BuildClusterSets wouldn't have included them anywhere
197  */
198  void FinalizeClusterSets();
199 
200  /// A function to add a cluster to a merged sets (_cluster_sets_v)
201  int AppendToClusterSets(unsigned int cluster_index, int merged_index=-1);
202 
203  protected:
204 
205  bool _verbose; ///< Verbose mode boolean
206  bool _det_params_prepared; ///< Boolean to keep track of detector parameter preparation
207  TTree* _merge_tree; ///< Quality Control TTree pointer
208 
209  /**
210  Book-keeping vector which length is same as input cluster array's length.
211  The stored value is the merged cluster sets' index (_cluster_sets_v).
212  For instance, given 5 clusters (0, 1, 2, 3, 4) as an input among which
213  1,2,3 are to be merged. _cluster_merged_index may hold contents like
214  [1,0,0,0,2] when _cluster_sets_v contents are [[1,2,3],[0],[4]].
215  */
216  std::vector<int> _cluster_merged_index;
217 
218  /**
219  The result container of ProcessMergeAlg() function.
220  The structure is such that the inner vector holds the cluster IDs to be merged into one cluster.
221  Naturally we expect multiple merged clusters, hence it's a vector of vector.
222  */
223  std::vector<std::vector<unsigned int> > _cluster_sets_v;
224 
225  std::vector<cluster::cluster_merge_info> _u_clusters; ///< Input U-plane clusters' information
226  std::vector<cluster::cluster_merge_info> _v_clusters; ///< Input V-plane clusters' information
227  std::vector<cluster::cluster_merge_info> _w_clusters; ///< Input W-plane clusters' information
228 
229  double _wire_2_cm; ///< Conversion factor from wire number to cm scale
230  double _time_2_cm; ///< Conversion factor from time to cm scale
231  double _max_allowed_2D_angle_diff; //in degrees
232  double _max_2D_dist2; //in cm^2
233  double _min_distance_unit; //in cm^2
234  }; // class ClusterMergeAlg
235 
236 } //namespace cluster
237 #endif
void VerboseMode(bool on)
Method to set verbose mode.
void SetSquaredDistanceCut(double d)
Method to set cut value in cm^2 for distance compatibility test.
BEGIN_PROLOG true icarus_rawdigitfilter FilterTools FilterPlane1 Plane
double ShortestDistanceSquared(double point_x, double point_y, double start_x, double start_y, double end_x, double end_y) const
process_name cluster
Definition: cheaterreco.fcl:51
geo::PlaneID planeID
plane ID
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Unknown view.
Definition: geo_types.h:136
Declaration of signal hit object.
void ReportConfig() const
Method to report the current configuration.
std::vector< cluster::cluster_merge_info > _v_clusters
Input V-plane clusters&#39; information.
double _time_2_cm
Conversion factor from time to cm scale.
const std::vector< std::vector< unsigned int > > GetClusterSets() const
Method to extract resulting set of cluster IDs for merging computed by ProcessMergeAlg() function...
cluster_merge_info()
Default constructor.
int AppendToClusterSets(unsigned int cluster_index, int merged_index=-1)
A function to add a cluster to a merged sets (_cluster_sets_v)
The data type to uniquely identify a Plane.
Definition: geo_types.h:472
double end_wire_err
End point wire error.
cluster_merge_info(const recob::Cluster &cl)
Initialization from a recob::Cluster.
double _wire_2_cm
Conversion factor from wire number to cm scale.
BEGIN_PROLOG StartTick
bool ShortestDistanceCompatibility(const cluster_merge_info &clus_info_A, const cluster_merge_info &clus_info_B) const
TTree * _merge_tree
Quality Control TTree pointer.
void ClearInputInfo()
Method to clear input cluster information.
Set of hits with a 2D structure.
Definition: Cluster.h:71
std::vector< int > _cluster_merged_index
ClusterMergeAlg(fhicl::ParameterSet const &pset)
Default constructor with fhicl parameters.
std::vector< cluster::cluster_merge_info > _w_clusters
Input W-plane clusters&#39; information.
float end_time
End point time.
double start_wire_err
Vertex wire error.
bool _verbose
Verbose mode boolean.
float start_wire
Vertex wire.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
void SetTime2Cm(double f)
Method to set a conversion factor from time to cm scale.
void ClearTTreeInfo()
Method to clear quality control TTree variables.
void PrepareDetParams()
Method to prepare detector parameters.
void SetAngleCut(double angle)
Method to set cut value in degrees for angle compatibility test.
void PrepareTTree()
Method to prepare quality control TTree.
float start_time
Vertex time.
bool Angle2DCompatibility(const cluster_merge_info &clus_info_A, const cluster_merge_info &clus_info_B) const
void BuildClusterSets(const cluster_merge_info &clus_info_A, const cluster_merge_info &clus_info_B)
void AppendHitInfo(cluster_merge_info &ci, const std::vector< art::Ptr< recob::Hit > > &in_hit_v)
Method to fill hit-array-related information.
Declaration of cluster object.
float angle
2D starting angle (in radians)
Definition of data types for geometry description.
std::vector< std::vector< unsigned int > > _cluster_sets_v
bool _det_params_prepared
Boolean to keep track of detector parameter preparation.
unsigned int cluster_index
Input cluster ID.
void PrintClusterVars(cluster_merge_info clus_info) const
std::vector< cluster::cluster_merge_info > _u_clusters
Input U-plane clusters&#39; information.
void ClearEventInfo()
Method to clear event-wise information (both input cluster info &amp; output merged cluster sets) ...
bool CompareClusters(const cluster_merge_info &clus_info_A, const cluster_merge_info &clus_info_B)
Method to compare a compatibility between two clusters.
float end_wire
End point wire.
void AppendClusterInfo(const recob::Cluster &in_cluster, const std::vector< art::Ptr< recob::Hit > > &in_hit_v)
Method to add a cluster information for processing.
int isInClusterSets(unsigned int index) const
geo::View_t view
Wire plane ID.
finds tracks best matching by angle
void ClearOutputInfo()
Method to clear output merged cluster sets (_cluster_sets_v)
double start_time_err
Vertex time error.
double end_time_err
End point time error.
void SetWire2Cm(double f)
Method to set a conversion factor from wire to cm scale.