LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ClusterMergeAlg.h
Go to the documentation of this file.
1 // \file ClusterMergeAlg.h
3 //
4 // \brief ClusterMergeAlg header file
5 //
6 // \author david.kaleko@gmail.com, kazuhiro@nevis.columbia.edu
7 //
9 
10 #ifndef CLUSTERMERGEALG_H
11 #define CLUSTERMERGEALG_H
12 
13 // ART includes
14 namespace fhicl {
15  class ParameterSet;
16 }
17 
19 
20 // LArSoft
24 
25 // STL
26 #include <vector>
27 
28 class TTree;
29 
30 namespace cluster {
31 
37 
38  unsigned int cluster_index;
41 
42  float start_wire;
43  float start_time;
44  float end_wire;
45  float end_time;
46 
47  double start_wire_err;
48  double start_time_err;
49  double end_wire_err;
50  double end_time_err;
51 
52  float angle;
53 
55  cluster_merge_info() : planeID()
56  {
57 
58  cluster_index = 0xffffffff;
59  view = geo::kUnknown;
60  start_wire = start_time = end_wire = end_time = -1;
61  start_wire_err = start_time_err = end_wire_err = end_time_err = -1;
62  };
63 
65  explicit cluster_merge_info(const recob::Cluster& cl)
66  : cluster_index(cl.ID())
67  , view(cl.View())
68  , planeID(cl.Plane())
69  , start_wire(cl.StartWire())
70  , start_time(cl.StartTick())
71  , end_wire(cl.EndWire())
72  , end_time(cl.EndTick())
73  , angle(cl.StartAngle())
74  {}
75  };
76 
78 
79  public:
82 
84  void VerboseMode(bool on) { _verbose = on; }
85 
87  void ReportConfig() const;
88 
90  void SetAngleCut(double angle) { _max_allowed_2D_angle_diff = angle; }
91 
93  void SetSquaredDistanceCut(double d) { _max_2D_dist2 = d; }
94 
96  void AppendClusterInfo(const recob::Cluster& in_cluster,
97  const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
98 
100  void AppendClusterInfo(const art::Ptr<recob::Cluster> in_cluster,
101  const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
102 
104  void ClearEventInfo();
105 
111  void ProcessMergeAlg();
112 
114  const std::vector<std::vector<unsigned int>> GetClusterSets() const { return _cluster_sets_v; };
115 
117  bool CompareClusters(const cluster_merge_info& clus_info_A,
118  const cluster_merge_info& clus_info_B);
119 
125  bool Angle2DCompatibility(const cluster_merge_info& clus_info_A,
126  const cluster_merge_info& clus_info_B) const;
127 
134  bool ShortestDistanceCompatibility(const cluster_merge_info& clus_info_A,
135  const cluster_merge_info& clus_info_B) const;
136 
141  double ShortestDistanceSquared(double point_x,
142  double point_y,
143  double start_x,
144  double start_y,
145  double end_x,
146  double end_y) const;
147 
152  void PrintClusterVars(cluster_merge_info clus_info) const;
153 
159  int isInClusterSets(unsigned int index) const;
160 
161  protected:
163  void SetWire2Cm(double f) { _wire_2_cm = f; }
164 
166  void SetTime2Cm(double f) { _time_2_cm = f; }
167 
169  void ClearOutputInfo();
170 
172  void ClearInputInfo();
173 
175  void ClearTTreeInfo();
176 
178  void PrepareTTree();
179 
181  void PrepareDetParams();
182 
184  void AppendHitInfo(cluster_merge_info& ci, const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
185 
190  void BuildClusterSets(const cluster_merge_info& clus_info_A,
191  const cluster_merge_info& clus_info_B);
192 
197  void FinalizeClusterSets();
198 
200  int AppendToClusterSets(unsigned int cluster_index, int merged_index = -1);
201 
202  protected:
203  bool _verbose;
205  TTree* _merge_tree;
206 
214  std::vector<int> _cluster_merged_index;
215 
221  std::vector<std::vector<unsigned int>> _cluster_sets_v;
222 
223  std::vector<cluster::cluster_merge_info> _u_clusters;
224  std::vector<cluster::cluster_merge_info> _v_clusters;
225  std::vector<cluster::cluster_merge_info> _w_clusters;
226 
227  double _wire_2_cm;
228  double _time_2_cm;
229  double _max_allowed_2D_angle_diff; //in degrees
230  double _max_2D_dist2; //in cm^2
231  double _min_distance_unit; //in cm^2
232  }; // class ClusterMergeAlg
233 
234 } //namespace cluster
235 #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.
geo::PlaneID planeID
plane ID
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Unknown view.
Definition: geo_types.h:142
Declaration of signal hit object.
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.
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
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.
TTree * _merge_tree
Quality Control TTree pointer.
Set of hits with a 2D structure.
Definition: Cluster.h:69
std::vector< int > _cluster_merged_index
Cluster finding and building.
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.
TFile f
Definition: plotHisto.C:6
float start_wire
Vertex wire.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
void SetTime2Cm(double f)
Method to set a conversion factor from time to cm scale.
void SetAngleCut(double angle)
Method to set cut value in degrees for angle compatibility test.
float start_time
Vertex time.
parameter set interface
Float_t d
Definition: plot.C:235
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.
std::vector< cluster::cluster_merge_info > _u_clusters
Input U-plane clusters&#39; information.
float end_wire
End point wire.
geo::View_t view
Wire plane ID.
recob::tracking::Plane Plane
Definition: TrackState.h:17
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.