ClusterMatchAlg.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // ClusterMatchAlg source
4 //
5 ////////////////////////////////////////////////////////////////////////
6 
7 #include "ClusterMatchAlg.h"
8 
9 #include "TString.h"
10 #include "TTree.h"
11 
14 #include "art_root_io/TFileService.h"
17 #include "fhiclcpp/ParameterSet.h"
19 
28 
29 namespace cluster {
30 
31  ClusterMatchAlg::ClusterMatchAlg(fhicl::ParameterSet const& pset) : _ModName_MCTruth("")
32  {
33  _debug_mode = pset.get<bool>("DebugMode");
34  _store_sps = pset.get<bool>("StoreSpacePoint");
35  _num_sps_cut = pset.get<size_t>("CutParam_NumSpacePoint");
36  _overlay_tratio_cut = pset.get<double>("CutParam_OverlayTimeFraction");
37  _qratio_cut = pset.get<double>("CutParam_SumChargeRatio");
38  std::vector<size_t> algo_list = pset.get<std::vector<size_t>>("MatchAlgoList");
39 
40  _sps_algo = new trkf::SpacePointAlg(pset.get<fhicl::ParameterSet>("SpacePointAlg"));
41  _match_tree = 0;
42  _cluster_tree = 0;
43  _save_cluster_info = true;
44  _det_params_prepared = false;
45 
46  for (size_t i = 0; i < (size_t)(kMATCH_METHOD_MAX); ++i)
47 
48  _match_methods[i] = false;
49 
50  for (auto const v : algo_list) {
51 
52  if (v >= (size_t)(kMATCH_METHOD_MAX))
53 
54  mf::LogError("ClusterMatchAlg") << Form("Invalid algorithm enum: %zu", v);
55 
56  else
57  _match_methods[v] = true;
58  }
59 
60  ReportConfig();
61 
63  }
64 
65  void
67  {
68  std::ostringstream msg;
69  msg << std::endl
70  << " ClusterMatchAlg Configuration: " << std::endl
71  << "---------------------------------------------" << std::endl;
72  msg << " Debug Mode ... " << (_debug_mode ? "enabled!" : "disabled!") << std::endl;
73  msg << " RoughZ ....... " << (_match_methods[kRoughZ] ? "enabled!" : "disabled!") << std::endl;
74  msg << " RoughT ....... " << (_match_methods[kRoughT] ? "enabled!" : "disabled!") << std::endl;
75  msg << " SpacePoint ... " << (_match_methods[kSpacePoint] ? "enabled!" : "disabled!")
76  << std::endl;
77  msg << " SumCharge .... " << (_match_methods[kSumCharge] ? "enabled!" : "disabled!")
78  << std::endl;
79  msg << std::endl;
80  msg << " Overlay-Time Fraction Cut : " << _overlay_tratio_cut << std::endl
81  << " Charge-Ratio Diff. Cut : " << _qratio_cut << std::endl
82  << " Minimum # of SpacePoint : " << _num_sps_cut << std::endl
83  << std::endl;
84  msg << "---------------------------------------------" << std::endl;
85 
86  mf::LogWarning("ClusterMatchAlg") << msg.str();
87  }
88 
89  void
91  {
92  _ucluster_v.clear();
93  _vcluster_v.clear();
94  _wcluster_v.clear();
95 
96  _uhits_v.clear();
97  _vhits_v.clear();
98  _whits_v.clear();
99  }
100 
101  void
103  {
104  _matched_uclusters_v.clear();
105  _matched_vclusters_v.clear();
106  _matched_wclusters_v.clear();
107  _matched_sps_v.clear();
108  }
109 
110  void
112  {
113 
114  _mc_E = 0;
115  _mc_Px = 0;
116  _mc_Py = 0;
117  _mc_Pz = 0;
118  _mc_Vx = 0;
119  _mc_Vy = 0;
120  _mc_Vz = 0;
121  _pdgid = 0;
122  _tot_u = 0;
123  _tot_v = 0;
124  _tot_w = 0;
125  _tot_pass_z = 0;
126  _tot_pass_t = 0;
127  _tot_pass_sps = 0;
128  _tot_pass_qsum = 0;
129  _qratio_v.clear();
130  _uv_tratio_v.clear();
131  _vw_tratio_v.clear();
132  _wu_tratio_v.clear();
133  _u_nhits_v.clear();
134  _v_nhits_v.clear();
135  _w_nhits_v.clear();
136  _nsps.clear();
137 
138  _view_v.clear();
139  _charge_v.clear();
140  _nhits_v.clear();
141  _tstart_min_v.clear();
142  _tstart_max_v.clear();
143  _tpeak_min_v.clear();
144  _tpeak_max_v.clear();
145  _tend_min_v.clear();
146  _tend_max_v.clear();
147  }
148 
149  void
151  {
154  ClearTTreeInfo();
155  }
156 
157  void
159  {
160  if (!_match_tree) {
162  _match_tree = fileService->make<TTree>("match_tree", "");
163  _match_tree->Branch("mc_E", &_mc_E, "mc_E/D");
164  _match_tree->Branch("mc_Px", &_mc_Px, "mc_Px/D");
165  _match_tree->Branch("mc_Py", &_mc_Py, "mc_Py/D");
166  _match_tree->Branch("mc_Pz", &_mc_Pz, "mc_Pz/D");
167  _match_tree->Branch("mc_Vx", &_mc_Vx, "mc_Vx/D");
168  _match_tree->Branch("mc_Vy", &_mc_Vy, "mc_Vy/D");
169  _match_tree->Branch("mc_Vz", &_mc_Vz, "mc_Vz/D");
170 
171  _match_tree->Branch("pdgid", &_pdgid, "pdgid/I");
172  _match_tree->Branch("tot_u", &_tot_u, "tot_u/s");
173  _match_tree->Branch("tot_v", &_tot_v, "tot_v/s");
174  _match_tree->Branch("tot_w", &_tot_w, "tot_w/s");
175  _match_tree->Branch("tot_pass_t", &_tot_pass_t, "tot_pass_t/s");
176  _match_tree->Branch("tot_pass_z", &_tot_pass_z, "tot_pass_z/s");
177  _match_tree->Branch("tot_pass_sps", &_tot_pass_sps, "tot_pass_sps/s");
178  _match_tree->Branch("tot_pass_qsum", &_tot_pass_qsum, "tot_pass_qsum/s");
179 
180  _match_tree->Branch("uv_tratio_v", "std::vector<double>", &_uv_tratio_v);
181  _match_tree->Branch("vw_tratio_v", "std::vector<double>", &_vw_tratio_v);
182  _match_tree->Branch("wu_tratio_v", "std::vector<double>", &_wu_tratio_v);
183 
184  _match_tree->Branch("qratio_v", "std::vector<double>", &_qratio_v);
185  _match_tree->Branch("u_nhits_v", "std::vector<UShort_t>", &_u_nhits_v);
186  _match_tree->Branch("v_nhits_v", "std::vector<UShort_t>", &_v_nhits_v);
187  _match_tree->Branch("w_nhits_v", "std::vector<UShort_t>", &_w_nhits_v);
188  _match_tree->Branch("nsps", "std::vector<UShort_t>", &_nsps);
189  }
192  _cluster_tree = fileService->make<TTree>("cluster_tree", "");
193  _cluster_tree->Branch("view_v", "std::vector<uint16_t>", &_view_v);
194  _cluster_tree->Branch("charge_v", "std::vector<double>", &_charge_v);
195  _cluster_tree->Branch("nhits_v", "std::vector<uint16_t>", &_nhits_v);
196  _cluster_tree->Branch("tstart_min_v", "std::vector<double>", &_tstart_min_v);
197  _cluster_tree->Branch("tstart_max_v", "std::vector<double>", &_tstart_max_v);
198  _cluster_tree->Branch("tpeak_min_v", "std::vector<double>", &_tpeak_min_v);
199  _cluster_tree->Branch("tpeak_max_v", "std::vector<double>", &_tpeak_max_v);
200  _cluster_tree->Branch("tend_min_v", "std::vector<double>", &_tend_min_v);
201  _cluster_tree->Branch("tend_max_v", "std::vector<double>", &_tend_max_v);
202  }
203  }
204 
205  void
207  {
208  if (!_ModName_MCTruth.size()) return;
209 
210  std::vector<const simb::MCTruth*> mciArray;
211 
212  try {
213 
214  evt.getView(_ModName_MCTruth, mciArray);
215  }
216  catch (art::Exception const& e) {
217 
218  if (e.categoryCode() != art::errors::ProductNotFound) throw;
219  }
220 
221  for (size_t i = 0; i < mciArray.size(); ++i) {
222 
223  if (i == 1) {
224  mf::LogWarning("ClusterMatchAlg") << " Ignoring > 2nd MCTruth in MC generator...";
225  break;
226  }
227  const simb::MCTruth* mci_ptr(mciArray.at(i));
228 
229  for (size_t j = 0; j < (size_t)(mci_ptr->NParticles()); ++j) {
230 
231  if (j == 1) {
232  mf::LogWarning("ClusterMatchAlg") << " Ignoring > 2nd MCParticle in MC generator...";
233  break;
234  }
235 
236  const simb::MCParticle part(mci_ptr->GetParticle(j));
237 
238  _pdgid = part.PdgCode();
239  _mc_E = part.E();
240  _mc_Px = part.Px();
241  _mc_Py = part.Py();
242  _mc_Pz = part.Pz();
243  _mc_Vx = part.Vx();
244  _mc_Vy = part.Vy();
245  _mc_Vz = part.Vz();
246  }
247  }
248  }
249 
250  void
252  {
253  if (!_det_params_prepared) {
254  // Total number of planes
256  _tot_planes = geo->Nplanes();
257 
258  // Ask DetectorPrperties about time-offset among different wire planes ... used to correct timing
259  // difference among different wire planes in the following loop.
260  _time_offset_uplane = det_prop.GetXTicksOffset(geo::kU, 0, 0);
261  _time_offset_vplane = det_prop.GetXTicksOffset(geo::kV, 0, 0);
263  if (_tot_planes > 2) _time_offset_wplane = det_prop.GetXTicksOffset(geo::kW, 0, 0);
264  _det_params_prepared = true;
265  }
266  }
267 
268  void
270  const recob::Cluster& in_cluster,
271  const std::vector<art::Ptr<recob::Hit>>& in_hit_v)
272  {
273  PrepareDetParams(det_prop);
274  cluster_match_info ci(in_cluster.ID());
275  ci.view = in_cluster.View();
276 
278  FillHitInfo(ci, hit_ptrv, in_hit_v);
279 
280  // Save created art::PtrVector & cluster_match_info struct object
281  switch (ci.view) {
282  case geo::kU:
283  _uhits_v.push_back(hit_ptrv);
284  _ucluster_v.push_back(ci);
286  break;
287  case geo::kV:
288  _vhits_v.push_back(hit_ptrv);
289  _vcluster_v.push_back(ci);
291  break;
292  case geo::kW:
293  _whits_v.push_back(hit_ptrv);
294  _wcluster_v.push_back(ci);
296  break;
297  default:
298  mf::LogError("ClusterMatchAlg") << Form("Found an invalid plane ID: %d", in_cluster.View());
299  }
300  }
301 
302  void
304  art::PtrVector<recob::Hit>& out_hit_v,
305  const std::vector<art::Ptr<recob::Hit>>& in_hit_v)
306  {
307 
308  out_hit_v.reserve(in_hit_v.size());
309 
310  double time_offset = 0;
311  if (ci.view == geo::kU)
312  time_offset = _time_offset_uplane;
313  else if (ci.view == geo::kV)
314  time_offset = _time_offset_vplane;
315  else if (ci.view == geo::kW)
316  time_offset = _time_offset_wplane;
317 
318  // Loop over hits in this cluster
319  for (auto const hit : in_hit_v) {
320 
321  unsigned int wire = hit->WireID().Wire;
322  double tstart = hit->PeakTimePlusRMS(-1.) - time_offset;
323  double tpeak = hit->PeakTime() - time_offset;
324  double tend = hit->PeakTimePlusRMS(+1.) - time_offset;
325 
326  ci.sum_charge += hit->Integral();
327 
328  ci.wire_max = (ci.wire_max < wire) ? wire : ci.wire_max;
329  ci.wire_min = (ci.wire_min > wire) ? wire : ci.wire_min;
330 
331  ci.start_time_max = (ci.start_time_max < tstart) ? tstart : ci.start_time_max;
332  ci.peak_time_max = (ci.peak_time_max < tpeak) ? tpeak : ci.peak_time_max;
333  ci.end_time_max = (ci.end_time_max < tend) ? tend : ci.end_time_max;
334 
335  ci.start_time_min = (ci.start_time_min > tstart) ? tstart : ci.start_time_min;
336  ci.peak_time_min = (ci.peak_time_min > tpeak) ? tpeak : ci.peak_time_min;
337  ci.end_time_min = (ci.end_time_min > tend) ? tend : ci.end_time_min;
338 
339  out_hit_v.push_back(hit);
340  }
341 
342  ci.nhits = in_hit_v.size();
343  }
344 
345  void
347  {
348  if (_cluster_tree) {
349  _view_v.push_back(ci.view);
350  _charge_v.push_back(ci.sum_charge);
351  _nhits_v.push_back(ci.nhits);
352  _tstart_min_v.push_back(ci.start_time_min);
353  _tstart_max_v.push_back(ci.start_time_max);
354  _tpeak_min_v.push_back(ci.peak_time_min);
355  _tpeak_max_v.push_back(ci.peak_time_max);
356  _tend_min_v.push_back(ci.end_time_min);
357  _tend_max_v.push_back(ci.end_time_max);
358  }
359  }
360 
361  //########################################################################################
362  bool
364  const cluster_match_info& ci2,
365  const geo::View_t v1,
366  const geo::View_t v2) const
367  {
369  double y, z_min, z_max;
370  y = z_min = z_max = -1;
371  geo_h->IntersectionPoint(ci1.wire_min, ci2.wire_min, v1, v2, 0, 0, y, z_min);
372  geo_h->IntersectionPoint(ci1.wire_max, ci2.wire_max, v1, v2, 0, 0, y, z_max);
373  return (z_max > z_min);
374  }
375 
376  bool
378  {
379  double time_overlay = std::min(ci1.end_time_max, ci2.end_time_max) -
381 
382  double overlay_tratio =
383  time_overlay /
384  (ci1.end_time_max - ci1.start_time_min + ci2.end_time_max - ci2.start_time_min) * 2.;
385 
386  if ((ci1.view == geo::kU && ci2.view == geo::kV) ||
387  (ci1.view == geo::kV && ci2.view == geo::kU))
388  _uv_tratio_v.push_back(overlay_tratio);
389  else if ((ci1.view == geo::kV && ci2.view == geo::kW) ||
390  (ci1.view == geo::kW && ci2.view == geo::kV))
391  _vw_tratio_v.push_back(overlay_tratio);
392  else if ((ci1.view == geo::kW && ci2.view == geo::kU) ||
393  (ci1.view == geo::kU && ci2.view == geo::kW))
394  _wu_tratio_v.push_back(overlay_tratio);
395 
396  return (overlay_tratio > _overlay_tratio_cut);
397  }
398 
399  bool
401  {
402  double qratio = (uc.sum_charge) / (vc.sum_charge);
403 
404  // For quality check log
405  _qratio_v.push_back(qratio);
406 
407  return ((1 - _qratio_cut) < qratio && (qratio) < (1 + _qratio_cut));
408  }
409 
410  bool
412  detinfo::DetectorPropertiesData const& det_prop,
413  const size_t uindex,
414  const size_t vindex,
415  const size_t windex,
416  std::vector<recob::SpacePoint>& sps_v)
417  {
418  bool use_wplane = _tot_planes > 2;
419 
420  if (uindex >= _ucluster_v.size() || vindex >= _vcluster_v.size() ||
421  (use_wplane && (windex >= _wcluster_v.size()))) {
422 
423  mf::LogError("ClusterMatchAlg")
424  << std::endl
425  << Form(
426  "Requested to cluster-index (U,V,W) = (%zu,%zu,%zu) where max-length is (%zu,%zu,%zu)",
427  uindex,
428  vindex,
429  windex,
430  _ucluster_v.size(),
431  _vcluster_v.size(),
432  _wcluster_v.size())
433  << std::endl;
434  return false;
435  }
436 
437  // Define a time range in which hits are used for spacepoint finding ... here "peak time" is the relevant one
438  double trange_min =
439  std::min(_ucluster_v.at(uindex).peak_time_min, _vcluster_v.at(vindex).peak_time_min);
440  if (use_wplane) trange_min = std::min(trange_min, _wcluster_v.at(windex).peak_time_min);
441 
442  double trange_max =
443  std::max(_ucluster_v.at(uindex).peak_time_max, _vcluster_v.at(vindex).peak_time_max);
444  if (use_wplane) trange_max = std::max(trange_max, _wcluster_v.at(windex).peak_time_max);
445 
446  // Space-point algorithm applies additional dT
447  trange_min -= _sps_algo->maxDT();
448  trange_max += _sps_algo->maxDT();
449 
450  // Make PtrVector<recob::Hit> for relevant Hits
451  art::PtrVector<recob::Hit> hit_group;
452  size_t max_size = _uhits_v.at(uindex).size() + _vhits_v.at(vindex).size();
453  if (use_wplane) max_size += _whits_v.at(windex).size();
454  hit_group.reserve(max_size);
455  // Loop over hits in U-plane
456  for (auto const hit : _uhits_v.at(uindex)) {
457  if (hit->PeakTime() < trange_min) continue;
458  if (hit->PeakTime() > trange_max) continue;
459  hit_group.push_back(hit);
460  }
461  // Check if any hit found in this plane
462  size_t u_nhits = hit_group.size();
463  if (!u_nhits && !_debug_mode) return false;
464  // Loop over hits in V-plane
465  for (auto const hit : _vhits_v.at(vindex)) {
466  if (hit->PeakTime() < trange_min) continue;
467  if (hit->PeakTime() > trange_max) continue;
468  hit_group.push_back(hit);
469  }
470  // Check if any hit found in this plane
471  size_t v_nhits = hit_group.size() - u_nhits;
472  if (!(v_nhits) && !_debug_mode) return false;
473 
474  // Loop over hits in W-plane
475  if (use_wplane) {
476  for (auto const hit : _whits_v.at(windex)) {
477  if (hit->PeakTime() < trange_min) continue;
478  if (hit->PeakTime() > trange_max) continue;
479  hit_group.push_back(hit);
480  }
481  }
482  // Check if any hit found in this plane
483  size_t w_nhits = hit_group.size() - u_nhits - v_nhits;
484  if (!(w_nhits) && use_wplane && !_debug_mode) return false;
485 
486  // Run SpacePoint finder algo
487  if (u_nhits && v_nhits && (!use_wplane || (w_nhits && use_wplane))) {
489  _sps_algo->makeSpacePoints(clock_data, det_prop, hit_group, sps_v);
490  }
491 
492  size_t nsps = sps_v.size();
493  _u_nhits_v.push_back(u_nhits);
494  _v_nhits_v.push_back(v_nhits);
495  if (use_wplane) _w_nhits_v.push_back(w_nhits);
496  _nsps.push_back(nsps);
497 
498  if (nsps < _num_sps_cut) return false;
499  return true;
500  }
501 
502  std::vector<std::vector<unsigned int>>
504  {
505  std::vector<std::vector<unsigned int>> result;
506  result.push_back(_matched_uclusters_v);
507  result.push_back(_matched_vclusters_v);
508  result.push_back(_matched_wclusters_v);
509  return result;
510  }
511 
512  void
514  detinfo::DetectorPropertiesData const& detProp)
515  {
516  std::ostringstream msg;
517  msg << Form("Received (U,V,W) = (%zu,%zu,%zu) clusters...",
518  _uhits_v.size(),
519  _vhits_v.size(),
520  _whits_v.size())
521  << std::endl;
522  _tot_u = _ucluster_v.size();
523  _tot_v = _vcluster_v.size();
524  _tot_w = _wcluster_v.size();
525 
526  if (!(_tot_u + _tot_v + _tot_w)) {
527 
528  mf::LogError(__PRETTY_FUNCTION__)
529  << "No input cluster info found! Aborting the function call...";
530 
531  return;
532  }
533 
534  // Initialization
535  PrepareTTree();
537 
538  bool overlay_2d = false;
539  for (size_t uci_index = 0; uci_index < _ucluster_v.size(); ++uci_index) {
540 
541  for (size_t vci_index = 0; vci_index < _vcluster_v.size(); ++vci_index) {
542 
543  overlay_2d = true;
544 
545  // Apply cuts
546  // Rough z-position overlay cut
547  if (_match_methods[kRoughZ]) {
548 
549  if (Match_RoughZ(_ucluster_v.at(uci_index), _vcluster_v.at(vci_index), geo::kU, geo::kV))
550  _tot_pass_z++;
551  else if (!_debug_mode)
552  continue;
553  else
554  overlay_2d = false;
555  }
556 
557  // Sum charge cut
558  if (_match_methods[kSumCharge]) {
559 
560  if (Match_SumCharge(_ucluster_v.at(uci_index), _vcluster_v.at(vci_index)))
561  _tot_pass_qsum++;
562  else if (!_debug_mode)
563  continue;
564  else
565  overlay_2d = false;
566  }
567 
568  // Rough time overlap cut
569  if (_match_methods[kRoughT]) {
570 
571  if (Match_RoughTime(_ucluster_v.at(uci_index), _vcluster_v.at(vci_index)))
572  _tot_pass_t++;
573  else if (!_debug_mode)
574  continue;
575  else
576  overlay_2d = false;
577  }
578 
579  // SpacePoint cut
580  std::vector<recob::SpacePoint> sps_v;
582 
583  if (Match_SpacePoint(clockData, detProp, uci_index, vci_index, 0, sps_v))
584  _tot_pass_sps++;
585  else if (!_debug_mode)
586  continue;
587  else
588  overlay_2d = false;
589  }
590 
591  if (overlay_2d) {
592  _matched_uclusters_v.push_back((unsigned int)(_ucluster_v[uci_index].cluster_index));
593  _matched_vclusters_v.push_back((unsigned int)(_vcluster_v[vci_index].cluster_index));
594  if (_store_sps) _matched_sps_v.push_back(sps_v);
595  }
596  } // end of ... _vcluster_v loop
597  } // end of ... _ucluster_v loop
598 
599  // Report
600  msg << std::endl
601  << Form("Found %zu matched cluster pairs...", _matched_uclusters_v.size()) << std::endl;
602  for (size_t i = 0; i < _matched_uclusters_v.size(); ++i) {
603 
604  if (i == 0) msg << "Listing matched clusters (U,V)..." << std::endl;
605 
606  msg << Form("Pair %-2zu: (%-3d, %-3d)", i, _matched_uclusters_v[i], _matched_vclusters_v[i])
607  << std::endl;
608  }
609  msg << std::endl;
610  mf::LogWarning("ClusterMatchAlg") << msg.str();
611 
612  if (_match_tree) _match_tree->Fill();
613  if (_cluster_tree) _cluster_tree->Fill();
614 
615  // Clear input event data holders
617  /// Clear TTree variables
618  ClearTTreeInfo();
619  }
620 
621  void
623  detinfo::DetectorPropertiesData const& det_prop)
624  {
625  std::ostringstream msg;
626  msg << Form("Received (U,V,W) = (%zu,%zu,%zu) clusters...",
627  _uhits_v.size(),
628  _vhits_v.size(),
629  _whits_v.size())
630  << std::endl;
631  _tot_u = _ucluster_v.size();
632  _tot_v = _vcluster_v.size();
633  _tot_w = _wcluster_v.size();
634 
635  if (!(_tot_u + _tot_v + _tot_w)) {
636 
637  mf::LogError(__PRETTY_FUNCTION__)
638  << "No input cluster info found! Aborting the function call...";
639 
640  return;
641  }
642 
643  // Clear match information
644  PrepareTTree();
646 
647  bool overlay_2d = true;
648  bool overlay_3d = true;
649  // Loop over all possible u-v-w cluster combination
650  for (size_t uci_index = 0; uci_index < _ucluster_v.size(); ++uci_index) {
651 
652  for (size_t vci_index = 0; vci_index < _vcluster_v.size(); ++vci_index) {
653 
654  // Apply cuts that can be done with U&V planes here
655  overlay_2d = true;
656 
657  // Rough z-position overlay cut
658  if (_match_methods[kRoughZ]) {
659 
660  if (Match_RoughZ(_ucluster_v.at(uci_index), _vcluster_v.at(vci_index), geo::kU, geo::kV))
661  _tot_pass_z++;
662  else if (!_debug_mode)
663  continue;
664  else
665  overlay_2d = false;
666  }
667 
668  // Sum charge cut
669  if (_match_methods[kSumCharge]) {
670 
671  if (Match_SumCharge(_ucluster_v.at(uci_index), _vcluster_v.at(vci_index)))
672  _tot_pass_qsum++;
673  else if (!_debug_mode)
674  continue;
675  else
676  overlay_2d = false;
677  }
678 
679  for (size_t wci_index = 0; wci_index < _wcluster_v.size(); ++wci_index) {
680 
681  overlay_3d = overlay_2d;
682  // Apply cuts that requires 3 planes here
683 
684  // Rough time overlap cut
685  if (_match_methods[kRoughT]) {
686 
687  bool rough_time_match =
688  Match_RoughTime(_ucluster_v.at(uci_index), _vcluster_v.at(vci_index));
689  if (!_debug_mode && !rough_time_match) continue;
690 
691  rough_time_match =
692  (Match_RoughTime(_vcluster_v.at(vci_index), _wcluster_v.at(wci_index)) &&
693  rough_time_match);
694  if (!_debug_mode && !rough_time_match) continue;
695 
696  rough_time_match =
697  (Match_RoughTime(_wcluster_v.at(wci_index), _ucluster_v.at(uci_index)) &&
698  rough_time_match);
699 
700  overlay_3d = overlay_3d && rough_time_match;
701  if (rough_time_match)
702  _tot_pass_t++;
703  else if (!_debug_mode)
704  continue;
705  }
706 
707  // SpacePoint cut
708  std::vector<recob::SpacePoint> sps_v;
710 
711  if (Match_SpacePoint(clock_data, det_prop, uci_index, vci_index, wci_index, sps_v))
712  _tot_pass_sps++;
713  else if (!_debug_mode)
714  continue;
715  else
716  overlay_3d = false;
717  }
718 
719  if (overlay_3d) {
720  _matched_uclusters_v.push_back((unsigned int)(_ucluster_v[uci_index].cluster_index));
721  _matched_vclusters_v.push_back((unsigned int)(_vcluster_v[vci_index].cluster_index));
722  _matched_wclusters_v.push_back((unsigned int)(_wcluster_v[wci_index].cluster_index));
723  if (_store_sps) _matched_sps_v.push_back(sps_v);
724  }
725  } // end of ... _wcluster_v loop
726  } // end of ... _vcluster_v loop
727  } // end of ... _ucluster_v loop
728 
729  // Report
730  msg << std::endl
731  << Form("Found %zu matched cluster pairs...", _matched_uclusters_v.size()) << std::endl;
732  for (size_t i = 0; i < _matched_uclusters_v.size(); ++i) {
733 
734  if (i == 0) msg << "Listing matched clusters (U,V,W)..." << std::endl;
735 
736  msg << Form("Pair %-2zu: (%-3d, %-3d, %-3d)",
737  i,
741  << std::endl;
742  }
743  msg << std::endl;
744  mf::LogWarning("ClusterMatchAlg") << msg.str();
745 
746  if (_match_tree) _match_tree->Fill();
747  if (_cluster_tree) _cluster_tree->Fill();
748 
750  ClearTTreeInfo();
751  }
752 
753 } // namespace match
std::vector< uint16_t > _w_nhits_v
Use summed charge comparison ... see Match_SumCharge() description.
void ClearMatchInputInfo()
Method to clear input cluster information.
void reserve(size_type n)
Definition: PtrVector.h:337
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.
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.
static QCString result
double end_time_min
Minimum "end time" 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.
void msg(const char *fmt,...)
Definition: message.cpp:107
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.
Planes which measure V.
Definition: geo_types.h:130
double GetXTicksOffset(int p, int t, int c) const
double end_time_max
Maximum "end time" 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)
std::vector< double > _wu_tratio_v
std::vector< double > _qratio_v
void ClearEventInfo()
Method to clear event-wise information.
struct vector vector
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
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.
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
Cluster finding and building.
Particle class.
unsigned short _tot_pass_sps
art framework interface to geometry description
std::vector< double > _tpeak_max_v
std::vector< double > _charge_v
std::vector< uint16_t > _nhits_v
Planes which measure U.
Definition: geo_types.h:129
void PrepareTTree()
Internal method to create output TTree for quality checking of the algorithm.
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
const double e
std::vector< double > _tend_max_v
std::vector< double > _tend_min_v
std::vector< uint16_t > _nsps
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
double peak_time_min
Minimum "peak time" among all hits in this cluster.
void makeSpacePoints(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const art::PtrVector< recob::Hit > &hits, std::vector< recob::SpacePoint > &spts) const
T get(std::string const &key) const
Definition: ParameterSet.h:271
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.
double maxDT() const noexcept
Definition: SpacePointAlg.h:96
std::vector< uint16_t > _u_nhits_v
static int max(int a, int b)
void AppendClusterTreeVariables(const cluster_match_info &ci)
Internal method to fill cluster-info tree.
bool IntersectionPoint(geo::WireID const &wid1, geo::WireID const &wid2, double &y, double &z) const
Returns the intersection point of two wires.
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.
size_type size() const
Definition: PtrVector.h:302
Detector simulation of raw signals on wires.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
void ClearTTreeInfo()
Method to clear TTree variables.
geo::View_t View() const
Returns the view for this cluster.
Definition: Cluster.h:741
void clearHitMap() const
bool Match_SumCharge(const cluster_match_info &uc, const cluster_match_info &vc)
Declaration of signal hit object.
std::vector< uint16_t > _view_v
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
unsigned short _tot_pass_qsum
Contains all timing reference information for the detector.
ID_t ID() const
Identifier of this cluster.
Definition: Cluster.h:738
unsigned short wire_max
Maximum wire number in this cluster.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
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.
double start_time_max
Maximum "start time" among all hits in this cluster.
std::size_t getView(std::string const &moduleLabel, std::string const &productInstanceName, std::string const &processName, std::vector< ELEMENT const * > &result) const
Definition: DataViewImpl.h:500
Planes which measure W (third view for Bo, MicroBooNE, etc).
Definition: geo_types.h:131
bool _match_methods[kMATCH_METHOD_MAX]
Boolean list for enabled algorithms.
TCEvent evt
Definition: DataStructs.cxx:7
void MatchTwoPlanes(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop)
Two plane version of cluster matching method.
Event generator information.
Definition: MCTruth.h:32
void ReportConfig() const
Method to report the current configuration.
Algorithm for generating space points from hits.
LArSoft geometry interface.
Definition: ChannelGeo.h:16
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.
QTextStream & endl(QTextStream &s)
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.