LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ClusterMatchAlg.h
Go to the documentation of this file.
1 // \file ClusterMatchAlg.h
3 //
4 // \brief ClusterMatchAlg header file
5 //
6 // \author kazuhiro@nevis.columbia.edu
7 //
9 
10 #ifndef CLUSTERMATCHALG_H
11 #define CLUSTERMATCHALG_H
12 
13 // ART includes
14 namespace art {
15  class Event;
16 }
17 
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:
53  kRoughZ = 0,
57  kMATCH_METHOD_MAX
58  };
59 
70 
71  unsigned short cluster_index;
73  unsigned int nhits;
74  unsigned short wire_max;
75  unsigned short wire_min;
76  double start_time_max;
77  double peak_time_max;
78  double end_time_max;
79  double start_time_min;
80  double peak_time_min;
81  double end_time_min;
82  double sum_charge;
83 
85  cluster_match_info(unsigned short index)
86  {
87  cluster_index = index;
88  view = geo::kUnknown;
89  wire_max = 0;
90  wire_min = 0xffff;
91  start_time_max = peak_time_max = end_time_max = -1.;
92  start_time_min = peak_time_min = end_time_min = 1.e9;
93  sum_charge = -1.;
94  };
95 
98  };
99 
100  public:
103 
105  void ReportConfig() const;
106 
108  void SetMCTruthModName(std::string name) { _ModName_MCTruth = name; }
109 
111  void FillMCInfo(const art::Event& evt);
112 
114  void AppendClusterInfo(detinfo::DetectorPropertiesData const& det_prop,
115  const recob::Cluster& in_cluster,
116  const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
117 
125  void MatchThreePlanes(detinfo::DetectorClocksData const& clock_data,
126  detinfo::DetectorPropertiesData const& det_prop);
127 
129  void MatchTwoPlanes(detinfo::DetectorClocksData const& clock_data,
130  detinfo::DetectorPropertiesData const& det_prop);
131 
133  std::vector<std::vector<unsigned int>> GetMatchedClusters() const;
134 
136  const std::vector<std::vector<recob::SpacePoint>>& GetMatchedSpacePoints() const
137  {
138  return _matched_sps_v;
139  };
140 
142  bool StoreSpacePoints() const { return _store_sps; }
143 
145  void ClearEventInfo();
146 
147  protected:
149  void ClearMatchInputInfo();
150 
152  void ClearMatchOutputInfo();
153 
155  void ClearTTreeInfo();
156 
158  void PrepareDetParams(detinfo::DetectorPropertiesData const& clockData);
159 
161  void PrepareTTree();
162 
163  void FillHitInfo(cluster_match_info& ci,
164  art::PtrVector<recob::Hit>& out_hit_v,
165  const std::vector<art::Ptr<recob::Hit>>& in_hit_v);
166 
168  void AppendClusterTreeVariables(const cluster_match_info& ci);
169 
174  bool Match_RoughZ(const cluster_match_info& ci1,
175  const cluster_match_info& ci2,
176  const geo::View_t v1,
177  const geo::View_t v2) const;
178 
180  bool Match_RoughTime(const cluster_match_info& ci1, const cluster_match_info& ci2);
181 
187  //bool Match_RoughTime(const cluster_match_info &ci1, const cluster_match_info &ci2, const cluster_match_info &ci3);
188 
193  bool Match_SumCharge(const cluster_match_info& uc, const cluster_match_info& vc);
194 
201  bool Match_SpacePoint(detinfo::DetectorClocksData const& clockData,
202  detinfo::DetectorPropertiesData const& detProp,
203  const size_t uindex,
204  const size_t vindex,
205  const size_t windex,
206  std::vector<recob::SpacePoint>& sps_v);
207 
208  //
209  // Cut parameter values
210  //
211  size_t _num_sps_cut;
212  double _overlay_tratio_cut;
214  double _qratio_cut;
216 
218  //
219  // Resulting data holder for matched cluster indexes
220  //
221  std::vector<unsigned int> _matched_uclusters_v;
222  std::vector<unsigned int> _matched_vclusters_v;
223  std::vector<unsigned int> _matched_wclusters_v;
224  std::vector<std::vector<recob::SpacePoint>>
226 
227  //
228  // Run control variables
229  //
230  bool _match_methods[kMATCH_METHOD_MAX];
232  bool _debug_mode;
233  bool _store_sps;
234  unsigned int _tot_planes;
238 
239  std::string _ModName_MCTruth;
240 
241  std::vector<art::PtrVector<recob::Hit>>
243  std::vector<art::PtrVector<recob::Hit>>
245  std::vector<art::PtrVector<recob::Hit>>
247 
248  std::vector<cluster_match_info> _ucluster_v;
249  std::vector<cluster_match_info> _vcluster_v;
250  std::vector<cluster_match_info> _wcluster_v;
251 
253 
254  //
255  // Quality control parameters to be saved in the TTree
256  //
257  TTree* _match_tree;
258  // Event-wise variables
259  double _mc_E;
260  double _mc_Px;
261  double _mc_Py;
262  double _mc_Pz;
263  double _mc_Vx;
264  double _mc_Vy;
265  double _mc_Vz;
266  int _pdgid;
267  unsigned short _tot_u;
268  unsigned short _tot_v;
269  unsigned short _tot_w;
270  unsigned short _tot_pass_qsum;
271  unsigned short _tot_pass_t;
272  unsigned short _tot_pass_z;
273  unsigned short _tot_pass_sps;
274 
275  // Cluster-combination-wise variable
276  std::vector<uint16_t> _u_nhits_v;
277  std::vector<uint16_t> _v_nhits_v;
278  std::vector<uint16_t> _w_nhits_v;
279  std::vector<uint16_t> _nsps;
280  std::vector<double> _qratio_v;
281  std::vector<double> _uv_tratio_v;
282  std::vector<double> _vw_tratio_v;
283  std::vector<double> _wu_tratio_v;
284 
285  // Cluster-wise variable
288  std::vector<uint16_t> _view_v;
289  std::vector<double> _charge_v;
290  std::vector<uint16_t> _nhits_v;
291  std::vector<double> _tstart_min_v;
292  std::vector<double> _tstart_max_v;
293  std::vector<double> _tpeak_min_v;
294  std::vector<double> _tpeak_max_v;
295  std::vector<double> _tend_min_v;
296  std::vector<double> _tend_max_v;
297 
298  }; // class ClusterMatchAlg
299 
300 } //namespace cluster
301 #endif
std::vector< uint16_t > _w_nhits_v
Use summed charge comparison ... see Match_SumCharge() description.
trkf::SpacePointAlg * _sps_algo
SpacePointFinder algorithm pointer.
double peak_time_max
Maximum "peak time" among all hits in this cluster.
std::vector< art::PtrVector< recob::Hit > > _vhits_v
Local Hit pointer vector container ... V-plane.
void SetMCTruthModName(std::string name)
Method to specify input MCTruth&#39;s module name (optional)
Reconstruction base classes.
double end_time_min
Minimum "end time" among all hits in this cluster.
std::vector< double > _vw_tratio_v
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:142
double end_time_max
Maximum "end time" among all hits in this cluster.
Declaration of signal hit object.
std::vector< double > _wu_tratio_v
std::vector< double > _qratio_v
bool _store_sps
Boolean to enable storage of SpacePoint vector.
MatchMethod_t
Enum switch for various matching methods.
Set of hits with a 2D structure.
Definition: Cluster.h:69
Use SpacePoint finder algorithm ... see Match_SpacePoint() description.
std::vector< cluster_match_info > _wcluster_v
Local cluster data container... W-plane.
Cluster finding and building.
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:289
std::vector< uint16_t > _nhits_v
unsigned short cluster_index
Cluster&#39;s index position in the input cluster vector array.
std::vector< double > _tend_max_v
std::vector< double > _tend_min_v
std::vector< uint16_t > _nsps
double peak_time_min
Minimum "peak time" among all hits in this cluster.
parameter set interface
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
General LArSoft Utilities.
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)
const std::vector< std::vector< recob::SpacePoint > > & GetMatchedSpacePoints() const
Method to retrieve matched SpacePoint for each combinations.
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< std::vector< recob::SpacePoint > > _matched_sps_v
Local SpacePoint vector container.
std::vector< unsigned int > _matched_uclusters_v
U plane matched clusters&#39; index.
Definition: MVAAlg.h:12
double start_time_max
Maximum "start time" among all hits in this cluster.
TCEvent evt
Definition: DataStructs.cxx:8
cluster_match_info(unsigned short index)
Constructor with cluster&#39;s index ID.
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.
std::vector< double > _tpeak_min_v
std::vector< cluster_match_info > _ucluster_v
Local cluster data container... U-plane.
double start_time_min
Minimum "start time" among all hits in this cluster.