LArSoft  v06_85_00
Liquid Argon Software toolkit - http://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 #include "fhiclcpp/ParameterSet.h"
24 
25 // LArSoft
34 // STL
35 #include <set>
36 #include <vector>
37 
38 // ROOT
39 #include <TString.h>
40 #include <TTree.h>
41 
42 namespace cluster
43 {
45 
46  public:
47 
49  enum MatchMethod_t { kRoughZ = 0,
54  };
55 
66 
67  unsigned short cluster_index;
69  unsigned int nhits;
70  unsigned short wire_max;
71  unsigned short wire_min;
72  double start_time_max;
73  double peak_time_max;
74  double end_time_max;
75  double start_time_min;
76  double peak_time_min;
77  double end_time_min;
78  double sum_charge;
79 
81  cluster_match_info(unsigned short index) {
82  cluster_index = index;
83  view = geo::kUnknown;
84  wire_max = 0;
85  wire_min = 0xffff;
86  start_time_max = peak_time_max = end_time_max = -1.;
87  start_time_min = peak_time_min = end_time_min = 1.e9;
88  sum_charge = -1.;
89  };
90 
93  cluster_match_info(0xffff);
94  };
95  };
96 
97  public:
98 
101  //ClusterMatchAlg(fhicl::ParameterSet const& pset, art::ActivityRegistry& reg);
102 
104  virtual ~ClusterMatchAlg(){};
105 
107  void ReportConfig() const;
108 
110  void SetMCTruthModName(std::string name)
111  {_ModName_MCTruth = name;}
112 
114  void FillMCInfo(const art::Event &evt);
115 
117  void AppendClusterInfo(const recob::Cluster &in_cluster,
118  const std::vector<art::Ptr<recob::Hit> > &in_hit_v);
119 
121  void AppendClusterInfo(const art::Ptr<recob::Cluster> in_cluster,
122  const std::vector<art::Ptr<recob::Hit> > &in_hit_v);
130  void MatchThreePlanes();
131 
133  void MatchTwoPlanes();
134 
136  std::vector<std::vector<unsigned int> > GetMatchedClusters() const;
137 
139  const std::vector<std::vector<recob::SpacePoint> >& GetMatchedSpacePoints() const { return _matched_sps_v; };
140 
142  bool StoreSpacePoints() const { return _store_sps; }
143 
145  void ClearEventInfo();
146 
147  protected:
148 
150  void ClearMatchInputInfo();
151 
153  void ClearMatchOutputInfo();
154 
156  void ClearTTreeInfo();
157 
159  void PrepareDetParams();
160 
162  void PrepareTTree();
163 
165  art::PtrVector<recob::Hit> &out_hit_v,
166  const std::vector<art::Ptr<recob::Hit> > &in_hit_v);
167 
170 
175  bool Match_RoughZ(const cluster_match_info &ci1, const cluster_match_info &ci2,
176  const geo::View_t v1, const geo::View_t v2 ) const;
177 
179  bool Match_RoughTime(const cluster_match_info &ci1, const cluster_match_info &ci2);
180 
186  //bool Match_RoughTime(const cluster_match_info &ci1, const cluster_match_info &ci2, const cluster_match_info &ci3);
187 
192  bool Match_SumCharge(const cluster_match_info &uc, const cluster_match_info &vc);
193 
200  bool Match_SpacePoint(const size_t uindex, const size_t vindex, const size_t windex, std::vector<recob::SpacePoint> &sps_v);
201 
202  //
203  // Cut parameter values
204  //
205  size_t _num_sps_cut;
207  double _qratio_cut;
208 
209  //
210  // Resulting data holder for matched cluster indexes
211  //
212  std::vector<unsigned int> _matched_uclusters_v;
213  std::vector<unsigned int> _matched_vclusters_v;
214  std::vector<unsigned int> _matched_wclusters_v;
215  std::vector<std::vector<recob::SpacePoint> > _matched_sps_v;
216 
217  //
218  // Run control variables
219  //
222  bool _debug_mode;
223  bool _store_sps;
224  unsigned int _tot_planes;
228 
229  std::string _ModName_MCTruth;
230 
231  std::vector<art::PtrVector<recob::Hit> > _uhits_v;
232  std::vector<art::PtrVector<recob::Hit> > _vhits_v;
233  std::vector<art::PtrVector<recob::Hit> > _whits_v;
234 
235  std::vector<cluster_match_info> _ucluster_v;
236  std::vector<cluster_match_info> _vcluster_v;
237  std::vector<cluster_match_info> _wcluster_v;
238 
240 
241  //
242  // Quality control parameters to be saved in the TTree
243  //
244  TTree* _match_tree;
245  // Event-wise variables
246  double _mc_E;
247  double _mc_Px;
248  double _mc_Py;
249  double _mc_Pz;
250  double _mc_Vx;
251  double _mc_Vy;
252  double _mc_Vz;
253  int _pdgid;
254  unsigned short _tot_u;
255  unsigned short _tot_v;
256  unsigned short _tot_w;
257  unsigned short _tot_pass_qsum;
258  unsigned short _tot_pass_t;
259  unsigned short _tot_pass_z;
260  unsigned short _tot_pass_sps;
261 
262  // Cluster-combination-wise variable
263  std::vector<uint16_t> _u_nhits_v;
264  std::vector<uint16_t> _v_nhits_v;
265  std::vector<uint16_t> _w_nhits_v;
266  std::vector<uint16_t> _nsps;
267  std::vector<double> _qratio_v;
268  std::vector<double> _uv_tratio_v;
269  std::vector<double> _vw_tratio_v;
270  std::vector<double> _wu_tratio_v;
271 
272  // Cluster-wise variable
275  std::vector<uint16_t> _view_v;
276  std::vector<double> _charge_v;
277  std::vector<uint16_t> _nhits_v;
278  std::vector<double> _tstart_min_v;
279  std::vector<double> _tstart_max_v;
280  std::vector<double> _tpeak_min_v;
281  std::vector<double> _tpeak_max_v;
282  std::vector<double> _tend_min_v;
283  std::vector<double> _tend_max_v;
284 
285  }; // class ClusterMatchAlg
286 
287 } //namespace cluster
288 #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 "peak time" 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 "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:83
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
void ClearEventInfo()
Method to clear event-wise information.
void FillHitInfo(cluster_match_info &ci, art::PtrVector< recob::Hit > &out_hit_v, const std::vector< art::Ptr< recob::Hit > > &in_hit_v)
double _qratio_cut
Maximum difference among clusters&#39; charge sum used in Match_SumCharge method.
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.
void MatchTwoPlanes()
Two plane version of cluster matching method.
void AppendClusterInfo(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.
std::vector< cluster_match_info > _wcluster_v
Local cluster data container... W-plane.
Cluster finding and building.
Particle class.
unsigned short _tot_pass_sps
bool StoreSpacePoints() const
Method to check if it is configured to store SpacePoint.
std::vector< double > _tpeak_max_v
double _overlay_tratio_cut
Minimum overlayed time fraction among two clusters used in Match_RoughTime method.
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 "peak time" 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.
bool Match_SpacePoint(const size_t uindex, const size_t vindex, const size_t windex, std::vector< recob::SpacePoint > &sps_v)
std::vector< art::PtrVector< recob::Hit > > _uhits_v
Local Hit pointer vector container ... U-plane.
std::vector< uint16_t > _u_nhits_v
Declaration of cluster object.
void AppendClusterTreeVariables(const cluster_match_info &ci)
Internal method to fill cluster-info tree.
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)
size_t _num_sps_cut
Number of SpacePoint used to cut in Match_SpacePoint method.
std::vector< uint16_t > _view_v
unsigned short _tot_pass_qsum
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.
double start_time_max
Maximum "start time" among all hits in this cluster.
bool _match_methods[kMATCH_METHOD_MAX]
Boolean list for enabled algorithms.
void PrepareDetParams()
Internal method, called only once, to fill detector-wise information.
void ReportConfig() const
Method to report the current configuration.
Algorithm for generating space points from hits.
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.
Rough-Z comparison method ... see Match_RoughZ() description.
art framework interface to geometry description
std::vector< double > _tpeak_min_v
std::vector< cluster_match_info > _ucluster_v
Local cluster data container... U-plane.
virtual ~ClusterMatchAlg()
Default destructor.
double start_time_min
Minimum "start time" among all hits in this cluster.