All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ClusterMatchAlg.h
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // \file ClusterMatchAlg.h
3 //
4 // \brief ClusterMatchAlg header file
5 //
6 // \author kazuhiro@nevis.columbia.edu
7 //
8 ////////////////////////////////////////////////////////////////////////
9 
10 #ifndef CLUSTERMATCHALG_H
11 #define CLUSTERMATCHALG_H
12 
13 // ART includes
14 namespace art {
15  class Event;
16 }
17 
18 #include "canvas/Persistency/Common/Ptr.h"
19 #include "canvas/Persistency/Common/PtrVector.h"
20 
21 namespace fhicl {
22  class ParameterSet;
23 }
24 
25 // LArSoft
26 namespace recob {
27  class Cluster;
28 }
32 
33 namespace detinfo {
34  class DetectorClocksData;
35  class DetectorPropertiesData;
36 }
37 namespace trkf {
38  class SpacePointAlg;
39 }
40 
41 // STL
42 #include <string>
43 #include <vector>
44 
45 class TTree;
46 
47 namespace cluster {
49 
50  public:
51  /// Enum switch for various matching methods
53  kRoughZ = 0, ///< Rough-Z comparison method ... see Match_RoughZ() description
54  kRoughT, ///< Rough-Time comparison method ... see Match_RoughTime() description
55  kSpacePoint, ///< Use SpacePoint finder algorithm ... see Match_SpacePoint() description
56  kSumCharge, ///< Use summed charge comparison ... see Match_SumCharge() description
58  };
59 
60  /**
61  Local struct data container to store cluster's basic information
62  based on hits that consist the cluster. Looping over hit pointer
63  occurs when we create art::PtrVector<recob::Hit> from input file.
64  All information that is based on hits and is used for cluster-matching
65  should be extracted from there to maximize I/O efficiency as looping
66  over hits takes time.
67  In other words... all hits related variables should be stored here!
68  */
70 
71  unsigned short cluster_index; ///< Cluster's index position in the input cluster vector array
72  geo::View_t view; ///< Wire plane ID
73  unsigned int nhits; ///< Number of hits
74  unsigned short wire_max; ///< Maximum wire number in this cluster
75  unsigned short wire_min; ///< Minimum wire number in this cluster
76  double start_time_max; ///< Maximum "start time" among all hits in this cluster
77  double peak_time_max; ///< Maximum "peak time" among all hits in this cluster
78  double end_time_max; ///< Maximum "end time" among all hits in this cluster
79  double start_time_min; ///< Minimum "start time" among all hits in this cluster
80  double peak_time_min; ///< Minimum "peak time" among all hits in this cluster
81  double end_time_min; ///< Minimum "end time" among all hits in this cluster
82  double sum_charge; ///< Summed charge among all hits in this cluster
83 
84  /// Constructor with cluster's index ID
85  cluster_match_info(unsigned short index)
86  {
87  cluster_index = index;
89  wire_max = 0;
90  wire_min = 0xffff;
93  sum_charge = -1.;
94  };
95 
96  /// Default constructor
98  };
99 
100  public:
101  /// Default constructor with fhicl parameters
102  ClusterMatchAlg(fhicl::ParameterSet const& pset);
103 
104  /// Method to report the current configuration
105  void ReportConfig() const;
106 
107  /// Method to specify input MCTruth's module name (optional)
108  void
109  SetMCTruthModName(std::string name)
110  {
112  }
113 
114  /// Internal method to fill MCTruth information when available
115  void FillMCInfo(const art::Event& evt);
116 
117  /// Method to fill cluster information to be used for matching
119  const recob::Cluster& in_cluster,
120  const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
121 
122  /**
123  Method to run matching algorithms for three planes.
124  Event info must be provided prior to this function call through FillEventInfo()
125  function call. If the function is called more than once w/o supplying the new
126  art::Event object, it does not perform any new matching unless a user explicitly
127  calls ClearEventInfo() and then fill event info again though FillEventInfo().
128  */
129  void MatchThreePlanes(detinfo::DetectorClocksData const& clock_data,
130  detinfo::DetectorPropertiesData const& det_prop);
131 
132  /// Two plane version of cluster matching method.
133  void MatchTwoPlanes(detinfo::DetectorClocksData const& clock_data,
134  detinfo::DetectorPropertiesData const& det_prop);
135 
136  /// Method to retrieve matched cluster combinations. The format is [wire_plane][cluster_index]
137  std::vector<std::vector<unsigned int>> GetMatchedClusters() const;
138 
139  /// Method to retrieve matched SpacePoint for each combinations.
140  const std::vector<std::vector<recob::SpacePoint>>&
142  {
143  return _matched_sps_v;
144  };
145 
146  /// Method to check if it is configured to store SpacePoint
147  bool
149  {
150  return _store_sps;
151  }
152 
153  /// Method to clear event-wise information
154  void ClearEventInfo();
155 
156  protected:
157  /// Method to clear input cluster information
158  void ClearMatchInputInfo();
159 
160  /// Method to clear output matched cluster information
161  void ClearMatchOutputInfo();
162 
163  /// Method to clear TTree variables
164  void ClearTTreeInfo();
165 
166  /// Internal method, called only once, to fill detector-wise information
167  void PrepareDetParams(detinfo::DetectorPropertiesData const& clockData);
168 
169  /// Internal method to create output TTree for quality checking of the algorithm
170  void PrepareTTree();
171 
172  void FillHitInfo(cluster_match_info& ci,
173  art::PtrVector<recob::Hit>& out_hit_v,
174  const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
175 
176  /// Internal method to fill cluster-info tree
177  void AppendClusterTreeVariables(const cluster_match_info& ci);
178 
179  /**
180  Match clusters based on min/max Z boundary information. It checks clusters' overlap
181  along Z spatial coordinate based on 2 input cluster information.
182  */
183  bool Match_RoughZ(const cluster_match_info& ci1,
184  const cluster_match_info& ci2,
185  const geo::View_t v1,
186  const geo::View_t v2) const;
187 
188  /// Checks min/max hit timing among two clusters and make sure there is an overlap
189  bool Match_RoughTime(const cluster_match_info& ci1, const cluster_match_info& ci2);
190 
191  /**
192  Checks min/max hit timing among three clusters and make sure there is an overlap.
193  If _overlay_tratio_cut is set, then overlapped-time / cluster-timespan fraction
194  is compared to the set cut value to claim a match.
195  */
196  //bool Match_RoughTime(const cluster_match_info &ci1, const cluster_match_info &ci2, const cluster_match_info &ci3);
197 
198  /**
199  Checks the ratio of two clusters' summed charge. If the ratio is within 1 +/- set cut value,
200  two clusters are considered to match.
201  */
202  bool Match_SumCharge(const cluster_match_info& uc, const cluster_match_info& vc);
203 
204  /**
205  Cluster matching using space points. This can be slow as we run SpacePointFinder
206  algorithm per cluster pair. Three clusters (U, V, W) are considerd to "match" if
207  there found N spacepoints using hits in them and N > min_nsps where "min_nsps" is
208  the cut value you can set in SetNSpacePointCut() method.
209  */
210  bool Match_SpacePoint(detinfo::DetectorClocksData const& clockData,
212  const size_t uindex,
213  const size_t vindex,
214  const size_t windex,
215  std::vector<recob::SpacePoint>& sps_v);
216 
217  //
218  // Cut parameter values
219  //
220  size_t _num_sps_cut; ///< Number of SpacePoint used to cut in
221  ///< Match_SpacePoint method
222  double _overlay_tratio_cut; ///< Minimum overlayed time fraction among two
223  ///< clusters used in Match_RoughTime method
224  double _qratio_cut; ///< Maximum difference among clusters' charge sum used
225  ///< in Match_SumCharge method
226 
227  //
228  // Resulting data holder for matched cluster indexes
229  //
230  std::vector<unsigned int> _matched_uclusters_v; ///< U plane matched clusters' index
231  std::vector<unsigned int> _matched_vclusters_v; ///< V plane matched clusters' index
232  std::vector<unsigned int> _matched_wclusters_v; ///< W plane matched clusters' index
233  std::vector<std::vector<recob::SpacePoint>>
234  _matched_sps_v; ///< Local SpacePoint vector container
235 
236  //
237  // Run control variables
238  //
239  bool _match_methods[kMATCH_METHOD_MAX]; ///< Boolean list for enabled algorithms
241  bool _debug_mode; ///< Boolean to enable debug mode (call all enabled matching methods)
242  bool _store_sps; ///< Boolean to enable storage of SpacePoint vector
243  unsigned int _tot_planes;
247 
248  std::string _ModName_MCTruth; ///< MCTruth producer's module name
249 
250  std::vector<art::PtrVector<recob::Hit>>
251  _uhits_v; ///< Local Hit pointer vector container ... U-plane
252  std::vector<art::PtrVector<recob::Hit>>
253  _vhits_v; ///< Local Hit pointer vector container ... V-plane
254  std::vector<art::PtrVector<recob::Hit>>
255  _whits_v; ///< Local Hit pointer vector container ... W-plane
256 
257  std::vector<cluster_match_info> _ucluster_v; ///< Local cluster data container... U-plane
258  std::vector<cluster_match_info> _vcluster_v; ///< Local cluster data container... V-plane
259  std::vector<cluster_match_info> _wcluster_v; ///< Local cluster data container... W-plane
260 
261  trkf::SpacePointAlg* _sps_algo; ///< SpacePointFinder algorithm pointer
262 
263  //
264  // Quality control parameters to be saved in the TTree
265  //
266  TTree* _match_tree;
267  // Event-wise variables
268  double _mc_E;
269  double _mc_Px;
270  double _mc_Py;
271  double _mc_Pz;
272  double _mc_Vx;
273  double _mc_Vy;
274  double _mc_Vz;
275  int _pdgid;
276  unsigned short _tot_u;
277  unsigned short _tot_v;
278  unsigned short _tot_w;
279  unsigned short _tot_pass_qsum;
280  unsigned short _tot_pass_t;
281  unsigned short _tot_pass_z;
282  unsigned short _tot_pass_sps;
283 
284  // Cluster-combination-wise variable
285  std::vector<uint16_t> _u_nhits_v;
286  std::vector<uint16_t> _v_nhits_v;
287  std::vector<uint16_t> _w_nhits_v;
288  std::vector<uint16_t> _nsps;
289  std::vector<double> _qratio_v;
290  std::vector<double> _uv_tratio_v;
291  std::vector<double> _vw_tratio_v;
292  std::vector<double> _wu_tratio_v;
293 
294  // Cluster-wise variable
297  std::vector<uint16_t> _view_v;
298  std::vector<double> _charge_v;
299  std::vector<uint16_t> _nhits_v;
300  std::vector<double> _tstart_min_v;
301  std::vector<double> _tstart_max_v;
302  std::vector<double> _tpeak_min_v;
303  std::vector<double> _tpeak_max_v;
304  std::vector<double> _tend_min_v;
305  std::vector<double> _tend_max_v;
306 
307  }; // class ClusterMatchAlg
308 
309 } //namespace cluster
310 #endif
std::vector< uint16_t > _w_nhits_v
Use summed charge comparison ... see Match_SumCharge() description.
void ClearMatchInputInfo()
Method to clear input cluster information.
std::vector< std::vector< unsigned int > > GetMatchedClusters() const
Method to retrieve matched cluster combinations. The format is [wire_plane][cluster_index].
trkf::SpacePointAlg * _sps_algo
SpacePointFinder algorithm pointer.
double peak_time_max
Maximum &quot;peak time&quot; among all hits in this cluster.
std::vector< art::PtrVector< recob::Hit > > _vhits_v
Local Hit pointer vector container ... V-plane.
void FillMCInfo(const art::Event &evt)
Internal method to fill MCTruth information when available.
void SetMCTruthModName(std::string name)
Method to specify input MCTruth&#39;s module name (optional)
bool Match_RoughTime(const cluster_match_info &ci1, const cluster_match_info &ci2)
Checks min/max hit timing among two clusters and make sure there is an overlap.
double end_time_min
Minimum &quot;end time&quot; among all hits in this cluster.
std::vector< double > _vw_tratio_v
void AppendClusterInfo(detinfo::DetectorPropertiesData const &det_prop, const recob::Cluster &in_cluster, const std::vector< art::Ptr< recob::Hit >> &in_hit_v)
Method to fill cluster information to be used for matching.
process_name cluster
Definition: cheaterreco.fcl:51
unsigned short wire_min
Minimum wire number in this cluster.
std::vector< art::PtrVector< recob::Hit > > _whits_v
Local Hit pointer vector container ... W-plane.
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Unknown view.
Definition: geo_types.h:136
double end_time_max
Maximum &quot;end time&quot; among all hits in this cluster.
bool Match_SpacePoint(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const size_t uindex, const size_t vindex, const size_t windex, std::vector< recob::SpacePoint > &sps_v)
Declaration of signal hit object.
std::vector< double > _wu_tratio_v
std::vector< double > _qratio_v
void ClearEventInfo()
Method to clear event-wise information.
void ClearMatchOutputInfo()
Method to clear output matched cluster information.
bool _store_sps
Boolean to enable storage of SpacePoint vector.
bool Match_RoughZ(const cluster_match_info &ci1, const cluster_match_info &ci2, const geo::View_t v1, const geo::View_t v2) const
MatchMethod_t
Enum switch for various matching methods.
Set of hits with a 2D structure.
Definition: Cluster.h:71
Use SpacePoint finder algorithm ... see Match_SpacePoint() description.
std::vector< cluster_match_info > _wcluster_v
Local cluster data container... W-plane.
unsigned short _tot_pass_sps
bool StoreSpacePoints() const
Method to check if it is configured to store SpacePoint.
std::vector< double > _tpeak_max_v
std::vector< double > _charge_v
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
std::vector< uint16_t > _nhits_v
unsigned short cluster_index
Cluster&#39;s index position in the input cluster vector array.
void PrepareTTree()
Internal method to create output TTree for quality checking of the algorithm.
std::vector< double > _tend_max_v
std::vector< double > _tend_min_v
std::vector< uint16_t > _nsps
double peak_time_min
Minimum &quot;peak time&quot; among all hits in this cluster.
double sum_charge
Summed charge among all hits in this cluster.
Rough-Time comparison method ... see Match_RoughTime() description.
std::string _ModName_MCTruth
MCTruth producer&#39;s module name.
std::vector< art::PtrVector< recob::Hit > > _uhits_v
Local Hit pointer vector container ... U-plane.
std::vector< uint16_t > _u_nhits_v
void AppendClusterTreeVariables(const cluster_match_info &ci)
Internal method to fill cluster-info tree.
Definition of data types for geometry description.
std::vector< double > _tstart_max_v
std::vector< cluster_match_info > _vcluster_v
Local cluster data container... V-plane.
std::vector< double > _tstart_min_v
std::vector< double > _uv_tratio_v
bool _debug_mode
Boolean to enable debug mode (call all enabled matching methods)
ClusterMatchAlg(fhicl::ParameterSet const &pset)
Default constructor with fhicl parameters.
const std::vector< std::vector< recob::SpacePoint > > & GetMatchedSpacePoints() const
Method to retrieve matched SpacePoint for each combinations.
void ClearTTreeInfo()
Method to clear TTree variables.
bool Match_SumCharge(const cluster_match_info &uc, const cluster_match_info &vc)
std::vector< uint16_t > _view_v
unsigned short _tot_pass_qsum
Contains all timing reference information for the detector.
unsigned short wire_max
Maximum wire number in this cluster.
std::vector< uint16_t > _v_nhits_v
std::vector< unsigned int > _matched_uclusters_v
U plane matched clusters&#39; index.
std::vector< std::vector< recob::SpacePoint > > _matched_sps_v
Local SpacePoint vector container.
then echo fcl name
double start_time_max
Maximum &quot;start time&quot; among all hits in this cluster.
bool _match_methods[kMATCH_METHOD_MAX]
Boolean list for enabled algorithms.
TCEvent evt
Definition: DataStructs.cxx:8
void MatchTwoPlanes(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop)
Two plane version of cluster matching method.
void ReportConfig() const
Method to report the current configuration.
cluster_match_info(unsigned short index)
Constructor with cluster&#39;s index ID.
void PrepareDetParams(detinfo::DetectorPropertiesData const &clockData)
Internal method, called only once, to fill detector-wise information.
std::vector< unsigned int > _matched_wclusters_v
W plane matched clusters&#39; index.
std::vector< unsigned int > _matched_vclusters_v
V plane matched clusters&#39; index.
void MatchThreePlanes(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop)
void FillHitInfo(cluster_match_info &ci, art::PtrVector< recob::Hit > &out_hit_v, const std::vector< art::Ptr< recob::Hit >> &in_hit_v)
Rough-Z comparison method ... see Match_RoughZ() description.
auto const detProp
std::vector< double > _tpeak_min_v
std::vector< cluster_match_info > _ucluster_v
Local cluster data container... U-plane.
double start_time_min
Minimum &quot;start time&quot; among all hits in this cluster.