All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BlurredClusteringAlg.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////
2 // Implementation of the Blurred Clustering algorithm
3 //
4 // Converts a hit map into a 2D image of the hits before convoling
5 // with a Gaussian function to introduce a weighted blurring.
6 // Clustering proceeds on this blurred image to create more
7 // complete clusters.
8 //
9 // M Wallbank (m.wallbank@sheffield.ac.uk), May 2015
10 ////////////////////////////////////////////////////////////////////
11 
13 #include "cetlib/pow.h"
19 #include "messagefacility/MessageLogger/MessageLogger.h"
20 
21 #include "RtypesCore.h"
22 #include "TCanvas.h"
23 #include "TColor.h"
24 #include "TH2.h"
25 #include "TLatex.h"
26 #include "TMarker.h"
27 #include "TString.h"
28 #include "TStyle.h"
29 #include "TVector2.h"
30 #include "TVirtualPad.h"
31 
32 #include "range/v3/view.hpp"
33 
34 #include <cassert>
35 #include <cmath>
36 
38  : fDebug{pset.get<bool>("Debug", false)}
39  , fDetector{pset.get<std::string>("Detector", "dune35t")}
40  , fBlurWire{pset.get<int>("BlurWire")}
41  , fBlurTick{pset.get<int>("BlurTick")}
42  , fSigmaWire{pset.get<double>("SigmaWire")}
43  , fSigmaTick{pset.get<double>("SigmaTick")}
44  , fMaxTickWidthBlur{pset.get<int>("MaxTickWidthBlur")}
45  , fClusterWireDistance{pset.get<int>("ClusterWireDistance")}
46  , fClusterTickDistance{pset.get<int>("ClusterTickDistance")}
47  , fNeighboursThreshold{pset.get<unsigned int>("NeighboursThreshold")}
48  , fMinNeighbours{pset.get<unsigned int>("MinNeighbours")}
49  , fMinSize{pset.get<unsigned int>("MinSize")}
50  , fMinSeed{pset.get<double>("MinSeed")}
51  , fTimeThreshold{pset.get<double>("TimeThreshold")}
52  , fChargeThreshold{pset.get<double>("ChargeThreshold")}
53  , fKernelWidth{2 * fBlurWire + 1}
56 {}
57 
59 {
60  if (fDebugCanvas) {
61  std::string closeName = fDebugPDFName;
62  closeName.append("]");
63  fDebugCanvas->Print(closeName.c_str());
64  delete fDebugCanvas;
65  }
66 }
67 
68 void
69 cluster::BlurredClusteringAlg::CreateDebugPDF(int run, int subrun, int event)
70 {
71  if (!fDebugCanvas) {
72 
73  // Create the grayscale palette for the Z axis
74  Double_t Red[2] = {1.00, 0.00};
75  Double_t Green[2] = {1.00, 0.00};
76  Double_t Blue[2] = {1.00, 0.00};
77  Double_t Length[2] = {0.00, 1.00};
78  TColor::CreateGradientColorTable(2, Length, Red, Green, Blue, 1000);
79  gStyle->SetOptStat(110000);
80 
81  // Decide what to call this PDF
82  std::ostringstream oss;
83  oss << "BlurredImages_Run" << run << "_Subrun" << subrun;
84  fDebugPDFName = oss.str();
85  fDebugCanvas = new TCanvas(fDebugPDFName.c_str(), "Image canvas", 1000, 500);
86  fDebugPDFName.append(".pdf");
87 
88  std::string openName = fDebugPDFName;
89  openName.append("[");
90  fDebugCanvas->Print(openName.c_str());
91  fDebugCanvas->Divide(2, 2);
92  fDebugCanvas->SetGrid();
93  }
94 
95  // Clear the pads on the canvas
96  for (int i = 1; i <= 4; ++i) {
97  fDebugCanvas->GetPad(i)->Clear();
98  }
99 
100  std::ostringstream oss;
101  oss << "Event " << event;
102  fDebugCanvas->cd(1);
103  TLatex l;
104  l.SetTextSize(0.15);
105  l.DrawLatex(0.1, 0.1, oss.str().c_str());
106  fDebugCanvas->Print(fDebugPDFName.c_str());
107 }
108 
109 void
111  std::vector<std::vector<double>> const& image,
112  std::vector<std::vector<int>> const& allClusterBins,
113  std::vector<art::PtrVector<recob::Hit>>& clusters) const
114 {
115  // Loop through the clusters (each a vector of bins)
116  for (auto const& bins : allClusterBins) {
117  // Convert the clusters (vectors of bins) to hits in a vector of recob::Hits
118  art::PtrVector<recob::Hit> clusHits = ConvertBinsToRecobHits(image, bins);
119 
120  mf::LogInfo("BlurredClustering") << "Cluster made from " << bins.size() << " bins, of which "
121  << clusHits.size() << " were real hits";
122 
123  // Make sure the clusters are above the minimum cluster size
124  if (clusHits.size() < fMinSize) {
125  mf::LogVerbatim("BlurredClustering")
126  << "Cluster of size " << clusHits.size()
127  << " not saved since it is smaller than the minimum cluster size, set to " << fMinSize;
128  continue;
129  }
130 
131  clusters.push_back(clusHits);
132  }
133 }
134 
135 std::vector<std::vector<double>>
137  std::vector<art::Ptr<recob::Hit>> const& hits,
138  int const readoutWindowSize)
139 {
140  // Define the size of this particular plane -- dynamically to avoid huge histograms
141  int lowerTick = readoutWindowSize, upperTick{}, lowerWire = fGeom->MaxWires(), upperWire{};
142  using lar::to_element;
144  for (auto const& hit : hits | transform(to_element)) {
145  int histWire = GlobalWire(hit.WireID());
146  if (hit.PeakTime() < lowerTick) lowerTick = hit.PeakTime();
147  if (hit.PeakTime() > upperTick) upperTick = hit.PeakTime();
148  if (histWire < lowerWire) lowerWire = histWire;
149  if (histWire > upperWire) upperWire = histWire;
150  }
151  fLowerTick = lowerTick - 20;
152  fUpperTick = upperTick + 20;
153  fLowerWire = lowerWire - 20;
154  fUpperWire = upperWire + 20;
155 
156  // Use a map to keep a track of the real hits and their wire/ticks
157  fHitMap.clear();
158  fHitMap.resize(fUpperWire - fLowerWire,
159  std::vector<art::Ptr<recob::Hit>>(fUpperTick - fLowerTick));
160 
161  // Create a 2D vector
162  std::vector<std::vector<double>> image(fUpperWire - fLowerWire,
163  std::vector<double>(fUpperTick - fLowerTick));
164 
165  // Look through the hits
166  for (auto const& hit : hits) {
167  int const wire = GlobalWire(hit->WireID());
168  auto const tick = static_cast<int>(hit->PeakTime());
169  float const charge = hit->Integral();
170 
171  // Fill hit map and keep a note of all real hits for later
172  if (charge > image.at(wire - fLowerWire).at(tick - fLowerTick)) {
173  image.at(wire - fLowerWire).at(tick - fLowerTick) = charge;
174  fHitMap[wire - fLowerWire][tick - fLowerTick] = hit;
175  }
176  }
177 
178  // Keep a note of dead wires
179  fDeadWires = std::vector<bool>(fUpperWire - fLowerWire, false);
180  geo::PlaneID const planeID = hits.front()->WireID().planeID();
181 
182  for (int wire = fLowerWire; wire < fUpperWire; ++wire) {
183  raw::ChannelID_t const channel =
184  fGeom->PlaneWireToChannel(planeID.Plane, wire, planeID.TPC, planeID.Cryostat);
185  fDeadWires[wire - fLowerWire] = !fChanStatus.IsGood(channel);
186  }
187 
188  return image;
189 }
190 
191 int
192 cluster::BlurredClusteringAlg::FindClusters(std::vector<std::vector<double>> const& blurred,
193  std::vector<std::vector<int>>& allcluster) const
194 {
195  // Size of image in x and y
196  int const nbinsx = blurred.size();
197  int const nbinsy = blurred.at(0).size();
198  int const nbins = nbinsx * nbinsy;
199 
200  // Vectors to hold hit information
201  std::vector<bool> used(nbins);
202  std::vector<std::pair<double, int>> values;
203 
204  // Place the bin number and contents as a pair in the values vector
205  for (int xbin = 0; xbin < nbinsx; ++xbin) {
206  for (int ybin = 0; ybin < nbinsy; ++ybin) {
207  int const bin = ConvertWireTickToBin(blurred, xbin, ybin);
208  values.emplace_back(ConvertBinToCharge(blurred, bin), bin);
209  }
210  }
211 
212  // Sort the values into charge order
213  std::sort(values.rbegin(), values.rend());
214 
215  // Count the number of iterations of the cluster forming loop (== number of clusters)
216  int niter = 0;
217 
218  // Clustering loops
219  // First loop - considers highest charge hits in decreasing order, and puts them in a new cluster if they aren't already clustered (makes new cluster every iteration)
220  // Second loop - looks at the direct neighbours of this seed and clusters to this if above charge/time thresholds. Runs recursively over all hits in cluster (inc. new ones)
221  while (true) {
222 
223  // Start a new cluster each time loop is executed
224  std::vector<int> cluster;
225  std::vector<double> times;
226 
227  // Get the highest charge bin (go no further if below seed threshold)
228  if (double const blurred_binval = values[niter].first; blurred_binval < fMinSeed) break;
229 
230  // Iterate through the bins from highest charge down
231  int const bin = values[niter++].second;
232 
233  // Put this bin in used if not already there
234  if (used[bin]) continue;
235  used[bin] = true;
236 
237  // Start a new cluster
238  cluster.push_back(bin);
239 
240  // Get the time of this hit
241  if (double const time = GetTimeOfBin(blurred, bin); time > 0) times.push_back(time);
242 
243  // Now cluster neighbouring hits to this seed
244  while (true) {
245 
246  bool added_cluster{false};
247 
248  for (unsigned int clusBin = 0; clusBin < cluster.size(); ++clusBin) {
249 
250  // Get x and y values for bin (c++ returns a%b = a if a<b)
251  int const binx = cluster[clusBin] % nbinsx;
252  int const biny = ((cluster[clusBin] - binx) / nbinsx) % nbinsy;
253 
254  // Look for hits in the neighbouring x/y bins
255  for (int x = binx - fClusterWireDistance; x <= binx + fClusterWireDistance; x++) {
256  if (x >= nbinsx or x < 0) continue;
257  for (int y = biny - fClusterTickDistance; y <= biny + fClusterTickDistance; y++) {
258  if (y >= nbinsy or y < 0) continue;
259  if (x == binx and y == biny) continue;
260 
261  // Get this bin
262  auto const bin = ConvertWireTickToBin(blurred, x, y);
263  if (bin >= nbinsx * nbinsy or bin < 0) continue;
264  if (used[bin]) continue;
265 
266  // Get the blurred value and time for this bin
267  double const blurred_binval = ConvertBinToCharge(blurred, bin);
268  double const time =
269  GetTimeOfBin(blurred, bin); // NB for 'fake' hits, time is defaulted to -10000
270 
271  // Check real hits pass time cut (ignores fake hits)
272  if (time > 0 && times.size() > 0 && !PassesTimeCut(times, time)) continue;
273 
274  // Add to cluster if bin value is above threshold
275  if (blurred_binval > fChargeThreshold) {
276  used[bin] = true;
277  cluster.push_back(bin);
278  added_cluster = true;
279  if (time > 0) { times.push_back(time); }
280  } // End of adding blurred bin to cluster
281  }
282  } // End of looking at directly neighbouring bins
283 
284  } // End of looping over bins already in this cluster
285 
286  if (!added_cluster) break;
287 
288  } // End of adding hits to this cluster
289 
290  // Check this cluster is above minimum size
291  if (cluster.size() < fMinSize) {
292  for (auto const bin : cluster) {
293  assert(bin >= 0);
294  used[bin] = false;
295  }
296  continue;
297  }
298 
299  // Fill in holes in the cluster
300  for (unsigned int clusBin = 0; clusBin < cluster.size(); clusBin++) {
301 
302  // Looks at directly neighbouring bins (and not itself)
303  for (int x = -1; x <= 1; ++x) {
304  for (int y = -1; y <= 1; ++y) {
305  if (x == 0 && y == 0) continue;
306 
307  // Look at neighbouring bins to the clustered bin which are inside the cluster
308  int neighbouringBin = cluster[clusBin] + x + (y * nbinsx);
309  if (neighbouringBin < nbinsx || neighbouringBin % nbinsx == 0 ||
310  neighbouringBin % nbinsx == nbinsx - 1 || neighbouringBin >= nbinsx * (nbinsy - 1))
311  continue;
312 
313  double const time = GetTimeOfBin(blurred, neighbouringBin);
314 
315  // If not already clustered and passes neighbour/time thresholds, add to cluster
316  if (!used[neighbouringBin] &&
317  (NumNeighbours(nbinsx, used, neighbouringBin) > fNeighboursThreshold) &&
318  PassesTimeCut(times, time)) {
319  used[neighbouringBin] = true;
320  cluster.push_back(neighbouringBin);
321 
322  if (time > 0) { times.push_back(time); }
323  } // End of clustering neighbouring bin
324  }
325  } // End of looping over neighbouring bins
326 
327  } // End of looping over bins already in cluster
328 
329  mf::LogVerbatim("Blurred Clustering")
330  << "Size of cluster after filling in holes: " << cluster.size();
331 
332  // Remove peninsulas - hits which have too few neighbouring hits in the cluster (defined by fMinNeighbours)
333  while (true) {
334  bool removed_cluster{false};
335 
336  // Loop over all the bins in the cluster
337  for (int clusBin = cluster.size() - 1; clusBin >= 0; clusBin--) {
338  auto const bin = cluster[clusBin];
339 
340  // If bin is in cluster ignore
341  if (bin < nbinsx || bin % nbinsx == 0 || bin % nbinsx == nbinsx - 1 ||
342  bin >= nbinsx * (nbinsy - 1))
343  continue;
344 
345  // Remove hit if it has too few neighbouring hits
346  if (NumNeighbours(nbinsx, used, bin) < fMinNeighbours) {
347  used[bin] = false;
348  removed_cluster = true;
349  cluster.erase(cluster.begin() + clusBin);
350  }
351  }
352 
353  if (!removed_cluster) break;
354  }
355 
356  mf::LogVerbatim("Blurred Clustering")
357  << "Size of cluster after removing peninsulas: " << cluster.size();
358 
359  // Disregard cluster if not of minimum size
360  if (cluster.size() < fMinSize) {
361  for (auto const bin : cluster) {
362  assert(bin >= 0);
363  used[bin] = false;
364  }
365  continue;
366  }
367 
368  // Put this cluster in the vector of clusters
369  allcluster.push_back(cluster);
370 
371  } // End loop over this cluster
372 
373  // Return the number of clusters found in this hit map
374  return allcluster.size();
375 }
376 
377 int
379 {
380  double globalWire = -999;
381 
382  // Induction
383  if (fGeom->SignalType(wireID) == geo::kInduction) {
384  double wireCentre[3];
385  fGeom->WireIDToWireGeo(wireID).GetCenter(wireCentre);
386  if (wireID.TPC % 2 == 0)
387  globalWire =
388  fGeom->WireCoordinate(wireCentre[1], wireCentre[2], wireID.Plane, 0, wireID.Cryostat);
389  else
390  globalWire =
391  fGeom->WireCoordinate(wireCentre[1], wireCentre[2], wireID.Plane, 1, wireID.Cryostat);
392  }
393 
394  // Collection
395  else {
396  // FOR COLLECTION WIRES, HARD CODE THE GEOMETRY FOR GIVEN DETECTORS
397  // THIS _SHOULD_ BE TEMPORARY. GLOBAL WIRE SUPPORT IS BEING ADDED TO THE LARSOFT GEOMETRY AND SHOULD BE AVAILABLE SOON
398  if (fDetector == "dune35t") {
399  unsigned int nwires = fGeom->Nwires(wireID.Plane, 0, wireID.Cryostat);
400  if (wireID.TPC == 0 or wireID.TPC == 1)
401  globalWire = wireID.Wire;
402  else if (wireID.TPC == 2 or wireID.TPC == 3 or wireID.TPC == 4 or wireID.TPC == 5)
403  globalWire = nwires + wireID.Wire;
404  else if (wireID.TPC == 6 or wireID.TPC == 7)
405  globalWire = (2 * nwires) + wireID.Wire;
406  else
407  mf::LogError("BlurredClusterAlg")
408  << "Error when trying to find a global induction plane coordinate for TPC " << wireID.TPC
409  << " (geometry " << fDetector << ")";
410  }
411  else if (fDetector == "dune10kt") {
412  unsigned int nwires = fGeom->Nwires(wireID.Plane, 0, wireID.Cryostat);
413  // Detector geometry has four TPCs, two on top of each other, repeated along z...
414  int block = wireID.TPC / 4;
415  globalWire = (nwires * block) + wireID.Wire;
416  }
417  else {
418  double wireCentre[3];
419  fGeom->WireIDToWireGeo(wireID).GetCenter(wireCentre);
420  if (wireID.TPC % 2 == 0)
421  globalWire =
422  fGeom->WireCoordinate(wireCentre[1], wireCentre[2], wireID.Plane, 0, wireID.Cryostat);
423  else
424  globalWire =
425  fGeom->WireCoordinate(wireCentre[1], wireCentre[2], wireID.Plane, 1, wireID.Cryostat);
426  }
427  }
428 
429  return std::round(globalWire);
430 }
431 
432 std::vector<std::vector<double>>
433 cluster::BlurredClusteringAlg::GaussianBlur(std::vector<std::vector<double>> const& image) const
434 {
435  if (fSigmaWire == 0 and fSigmaTick == 0) return image;
436 
437  auto const [blur_wire, blur_tick, sigma_wire, sigma_tick] = FindBlurringParameters();
438 
439  // Convolve the Gaussian
440  int width = 2 * blur_wire + 1;
441  int height = 2 * blur_tick + 1;
442  int nbinsx = image.size();
443  int nbinsy = image.at(0).size();
444 
445  // Blurred histogram and normalisation for each bin
446  std::vector<std::vector<double>> copy(nbinsx, std::vector<double>(nbinsy, 0));
447 
448  // Loop through all the bins in the histogram to blur
449  for (int x = 0; x < nbinsx; ++x) {
450  for (int y = 0; y < nbinsy; ++y) {
451 
452  if (image[x][y] == 0) continue;
453 
454  // Scale the tick blurring based on the width of the hit
455  int tick_scale =
456  std::sqrt(cet::square(fHitMap[x][y]->RMS()) + cet::square(sigma_tick)) / (double)sigma_tick;
457  tick_scale = std::max(std::min(tick_scale, fMaxTickWidthBlur), 1);
458  auto const& correct_kernel = fAllKernels[sigma_wire][sigma_tick * tick_scale];
459 
460  // Find any dead wires in the potential blurring region
461  auto const [lower_bin_dead, upper_bin_dead] = DeadWireCount(x, width);
462 
463  // Note of how many dead wires we have passed whilst blurring in the wire direction
464  // If blurring below the seed hit, need to keep a note of how many dead wires to come
465  // If blurring above, need to keep a note of how many dead wires have passed
466  auto dead_wires_passed{lower_bin_dead};
467 
468  // Loop over the blurring region around this hit
469  for (int blurx = -(width / 2 + lower_bin_dead); blurx < (width + 1) / 2 + upper_bin_dead;
470  ++blurx) {
471  if (x + blurx < 0) continue;
472  for (int blury = -height / 2 * tick_scale;
473  blury < ((((height + 1) / 2) - 1) * tick_scale) + 1;
474  ++blury) {
475  if (blurx < 0 and fDeadWires[x + blurx]) dead_wires_passed -= 1;
476 
477  // Smear the charge of this hit
478  double const weight = correct_kernel[fKernelWidth * (fKernelHeight / 2 + blury) +
479  (fKernelWidth / 2 + (blurx - dead_wires_passed))];
480  if (x + blurx >= 0 and x + blurx < nbinsx and y + blury >= 0 and y + blury < nbinsy)
481  copy[x + blurx][y + blury] += weight * image[x][y];
482 
483  if (blurx > 0 and fDeadWires[x + blurx]) dead_wires_passed += 1;
484  }
485  } // blurring region
486  }
487  } // hits to blur
488 
489  // HAVE REMOVED NOMALISATION CODE
490  // WHEN USING DIFFERENT KERNELS, THERE'S NO EASY WAY OF DOING THIS...
491  // RECONSIDER...
492 
493  // Return the blurred histogram
494  return copy;
495 }
496 
497 TH2F*
499  TString const name) const
500 {
501  auto hist = new TH2F(name,
502  name,
504  fLowerWire - 0.5,
505  fUpperWire - 0.5,
507  fLowerTick - 0.5,
508  fUpperTick - 0.5);
509  hist->SetXTitle("Wire number");
510  hist->SetYTitle("Tick number");
511  hist->SetZTitle("Charge");
512 
513  for (unsigned int imageWireIt = 0; imageWireIt < image.size(); ++imageWireIt) {
514  int const wire = imageWireIt + fLowerWire;
515  for (unsigned int imageTickIt = 0; imageTickIt < image.at(imageWireIt).size(); ++imageTickIt) {
516  int const tick = imageTickIt + fLowerTick;
517  hist->Fill(wire, tick, image.at(imageWireIt).at(imageTickIt));
518  }
519  }
520 
521  return hist;
522 }
523 
524 void
526  std::vector<art::PtrVector<recob::Hit>> const& allClusters,
527  int const pad,
528  int const tpc,
529  int const plane)
530 {
531  // Make a vector of clusters
532  std::vector<std::vector<int>> allClusterBins;
533 
534  for (auto const& cluster : allClusters) {
535  if (cluster.empty()) continue;
536 
537  std::vector<int> clusterBins;
538 
539  for (auto const& hit : cluster) {
540  unsigned int const wire = GlobalWire(hit->WireID());
541  float const tick = hit->PeakTime();
542  int bin = image->GetBin((wire - fLowerWire) + 1, (tick - fLowerTick) + 1);
543  if (cluster.size() < fMinSize) bin *= -1;
544 
545  clusterBins.push_back(bin);
546  }
547 
548  allClusterBins.push_back(clusterBins);
549  }
550 
551  SaveImage(image, allClusterBins, pad, tpc, plane);
552 }
553 
554 void
555 cluster::BlurredClusteringAlg::SaveImage(TH2F* image, int const pad, int const tpc, int const plane)
556 {
557  std::vector<std::vector<int>> allClusterBins;
558  SaveImage(image, allClusterBins, pad, tpc, plane);
559 }
560 
561 void
563  std::vector<std::vector<int>> const& allClusterBins,
564  int const pad,
565  int const tpc,
566  int const plane)
567 {
568  fDebugCanvas->cd(pad);
569  std::string stage;
570 
571  switch (pad) {
572  case 1: stage = "Stage 1: Unblurred"; break;
573  case 2: stage = "Stage 2: Blurred"; break;
574  case 3: stage = "Stage 3: Blurred with clusters overlaid"; break;
575  case 4: stage = "Stage 4: Output clusters"; break;
576  default: stage = "Unknown stage"; break;
577  }
578 
579  std::stringstream title;
580  title << stage << " -- TPC " << tpc << ", Plane " << plane; // << " (Event " << fEvent << ")";
581 
582  image->SetName(title.str().c_str());
583  image->SetTitle(title.str().c_str());
584  image->DrawCopy("colz");
585 
586  // Draw the clustered hits on the histograms
587  int clusterNum = 2;
588  for (auto const& bins : allClusterBins) {
589  TMarker mark(0, 0, 20);
590  mark.SetMarkerColor(clusterNum);
591  mark.SetMarkerSize(0.1);
592 
593  for (auto bin : bins) {
594  // Hit from a cluster that we aren't going to save
595  if (bin < 0) {
596  bin *= -1;
597  mark.SetMarkerStyle(24);
598  }
599 
600  int wire, tick, z;
601  image->GetBinXYZ(bin, wire, tick, z);
602  mark.DrawMarker(wire + fLowerWire - 1, tick + fLowerTick - 1);
603  mark.SetMarkerStyle(20);
604  }
605  }
606 
607  if (pad == 4) {
608  fDebugCanvas->Print(fDebugPDFName.c_str());
609  fDebugCanvas->Clear("D");
610  }
611 }
612 
613 // Private member functions
614 
615 art::PtrVector<recob::Hit>
617  std::vector<int> const& bins) const
618 {
619  // Create the vector of hits to output
620  art::PtrVector<recob::Hit> hits;
621 
622  // Look through the hits in the cluster
623  for (auto const bin : bins) {
624  // Take each hit and convert it to a recob::Hit
625  art::Ptr<recob::Hit> const hit = ConvertBinToRecobHit(image, bin);
626 
627  // If this hit was a real hit put it in the hit selection
628  if (!hit.isNull()) hits.push_back(hit);
629  }
630 
631  // Return the vector of hits to make cluster
632  return hits;
633 }
634 
635 art::Ptr<recob::Hit>
637  int const bin) const
638 {
639  int const wire = bin % image.size();
640  int const tick = bin / image.size();
641  return fHitMap[wire][tick];
642 }
643 
644 int
646  int const xbin,
647  int const ybin) const
648 {
649  return ybin * image.size() + xbin;
650 }
651 
652 double
654  int const bin) const
655 {
656  int const x = bin % image.size();
657  int const y = bin / image.size();
658  return image.at(x).at(y);
659 }
660 
661 std::pair<int, int>
662 cluster::BlurredClusteringAlg::DeadWireCount(int const wire_bin, int const width) const
663 {
664  auto deadWires = std::make_pair(0, 0);
665 
666  int const lower_bin = width / 2;
667  int const upper_bin = (width + 1) / 2;
668 
669  auto const offset = wire_bin + fLowerWire;
670  for (int wire = std::max(offset - lower_bin, fLowerWire);
671  wire < std::min(offset + upper_bin, fUpperWire);
672  ++wire) {
673  if (!fDeadWires[wire - fLowerWire]) continue;
674 
675  if (wire < offset)
676  ++deadWires.first;
677  else if (wire > offset)
678  ++deadWires.second;
679  }
680 
681  return deadWires;
682 }
683 
684 std::array<int, 4>
686 {
687  // Calculate least squares slope
688  double nhits{}, sumx{}, sumy{}, sumx2{}, sumxy{};
689  for (unsigned int wireIt = 0; wireIt < fHitMap.size(); ++wireIt) {
690  for (unsigned int tickIt = 0; tickIt < fHitMap.at(wireIt).size(); ++tickIt) {
691  if (fHitMap[wireIt][tickIt].isNull()) continue;
692  ++nhits;
693  int const x = wireIt + fLowerWire;
694  int const y = tickIt + fLowerTick;
695  sumx += x;
696  sumy += y;
697  sumx2 += x * x;
698  sumxy += x * y;
699  }
700  }
701  double const gradient = (nhits * sumxy - sumx * sumy) / (nhits * sumx2 - sumx * sumx);
702 
703  // Get the rough unit vector for the trajectories, making sure to
704  // catch the vertical gradient.
705  auto const unit = std::isnan(gradient) ? TVector2{0, 1} : TVector2{1, gradient}.Unit();
706 
707  // Use this direction to scale the blurring radii and Gaussian sigma
708  int const blur_wire = std::max(std::abs(std::round(fBlurWire * unit.X())), 1.);
709  int const blur_tick = std::max(std::abs(std::round(fBlurTick * unit.Y())), 1.);
710 
711  int const sigma_wire = std::max(std::abs(std::round(fSigmaWire * unit.X())), 1.);
712  int const sigma_tick = std::max(std::abs(std::round(fSigmaTick * unit.Y())), 1.);
713  return {{blur_wire, blur_tick, sigma_wire, sigma_tick}};
714 }
715 
716 double
718  int const bin) const
719 {
720  auto const hit = ConvertBinToRecobHit(image, bin);
721  return hit.isNull() ? -10000. : hit->PeakTime();
722 }
723 
724 std::vector<std::vector<std::vector<double>>>
726 {
727  // Kernel size is the largest possible given the hit width rescaling
728  std::vector<std::vector<std::vector<double>>> allKernels(
729  fSigmaWire + 1,
730  std::vector<std::vector<double>>(fSigmaTick * fMaxTickWidthBlur + 1,
731  std::vector<double>(fKernelWidth * fKernelHeight)));
732 
733  // Ranges of kernels to make
734  // Complete range of sigmas possible after dynamic fixing and hit width convolution
735  for (int sigma_wire = 1; sigma_wire <= fSigmaWire; ++sigma_wire) {
736  for (int sigma_tick = 1; sigma_tick <= fSigmaTick * fMaxTickWidthBlur; ++sigma_tick) {
737 
738  // New kernel
739  std::vector<double> kernel(fKernelWidth * fKernelHeight, 0);
740 
741  // Smear out according to the blur radii in each direction
742  for (int i = -fBlurWire; i <= fBlurWire; i++) {
743  for (int j = -fBlurTick * fMaxTickWidthBlur; j <= fBlurTick * fMaxTickWidthBlur; j++) {
744 
745  // Fill kernel
746  double const sig2i = 2. * sigma_wire * sigma_wire;
747  double const sig2j = 2. * sigma_tick * sigma_tick;
748 
749  int const key = (fKernelWidth * (j + fBlurTick * fMaxTickWidthBlur)) + (i + fBlurWire);
750  double const value = 1. / std::sqrt(sig2i * M_PI) * std::exp(-i * i / sig2i) * 1. /
751  std::sqrt(sig2j * M_PI) * std::exp(-j * j / sig2j);
752  kernel.at(key) = value;
753  }
754  } // End loop over blurring region
755 
756  allKernels[sigma_wire][sigma_tick] = move(kernel);
757  }
758  }
759  return allKernels;
760 }
761 
762 unsigned int
764  std::vector<bool> const& used,
765  int const bin) const
766 {
767  unsigned int neighbours = 0;
768 
769  // Loop over all directly neighbouring hits (not itself)
770  for (int x = -1; x <= 1; x++) {
771  for (int y = -1; y <= 1; y++) {
772  if (x == 0 && y == 0) continue;
773 
774  // Determine bin
775  int neighbouringBin =
776  bin + x +
777  (y *
778  nbinsx); /// 2D hists can be considered a string of bins - the equation to convert between them is [bin = x + (nbinsx * y)]
779 
780  // If this bin is in the cluster, increase the neighbouring bin counter
781  if (used.at(neighbouringBin)) neighbours++;
782  }
783  }
784 
785  // Return the number of neighbours in the cluster of a particular hit
786  return neighbours;
787 }
788 
789 bool
790 cluster::BlurredClusteringAlg::PassesTimeCut(std::vector<double> const& times,
791  double const time) const
792 {
793  for (auto const t : times) {
794  if (std::abs(time - t) < fTimeThreshold) return true;
795  }
796  return false;
797 }
process_name opflash particleana ie ie ie z
float Length(const PFPStruct &pfp)
Definition: PFPUtils.cxx:3304
BlurredClusteringAlg(fhicl::ParameterSet const &pset)
BEGIN_PROLOG TPC Trig offset(g4 rise time) ProjectToHeight
Definition: CORSIKAGen.fcl:7
Utilities related to art service access.
TH2F * MakeHistogram(std::vector< std::vector< double >> const &image, TString name) const
Converts a 2D vector in a histogram for the debug pdf.
constexpr to_element_t to_element
Definition: ToElement.h:24
static constexpr Sample_t transform(Sample_t sample)
process_name opflash particleana ie x
std::vector< std::vector< art::Ptr< recob::Hit > > > fHitMap
std::vector< std::vector< std::vector< double > > > MakeKernels() const
Makes all the kernels which could be required given the tuned parameters.
process_name cluster
Definition: cheaterreco.fcl:51
double GetTimeOfBin(std::vector< std::vector< double >> const &image, int bin) const
Returns the hit time of a hit in a particular bin.
bool PassesTimeCut(std::vector< double > const &times, double time) const
Determine if a hit is within a time threshold of any other hits in a cluster.
The data type to uniquely identify a Plane.
Definition: geo_types.h:472
unsigned int NumNeighbours(int nx, std::vector< bool > const &used, int bin) const
Determines the number of clustered neighbours of a hit.
std::array< int, 4 > FindBlurringParameters() const
Dynamically find the blurring radii and Gaussian sigma in each dimension.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:580
process_name hit
Definition: cheaterreco.fcl:51
int FindClusters(std::vector< std::vector< double >> const &image, std::vector< std::vector< int >> &allcluster) const
Find clusters in the histogram.
constexpr details::BinObj< T > bin(T value)
Returns a wrapper to print the specified data in binary format.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
art::PtrVector< recob::Hit > ConvertBinsToRecobHits(std::vector< std::vector< double >> const &image, std::vector< int > const &bins) const
Converts a vector of bins into a hit selection - not all the hits in the bins vector are real hits...
T abs(T value)
void CreateDebugPDF(int run, int subrun, int event)
Create the PDF to save debug images.
std::pair< int, int > DeadWireCount(int wire_bin, int width) const
double ConvertBinToCharge(std::vector< std::vector< double >> const &image, int bin) const
Returns the charge stored in the global bin value.
process_name opflash particleana ie ie y
int ConvertWireTickToBin(std::vector< std::vector< double >> const &image, int xbin, int ybin) const
Converts an xbin and a ybin to a global bin number.
Signal from induction planes.
Definition: geo_types.h:145
art::Ptr< recob::Hit > ConvertBinToRecobHit(std::vector< std::vector< double >> const &image, int bin) const
Converts a bin into a recob::Hit (not all of these bins correspond to recob::Hits - some are fake hit...
int GlobalWire(geo::WireID const &wireID) const
Find the global wire position.
tick_as<> tick
Tick number, represented by std::ptrdiff_t.
Definition: electronics.h:75
return match has_match and(match.match_pdg==11 or match.match_pdg==-11)
virtual bool IsGood(raw::ChannelID_t channel) const
Returns whether the specified channel is physical and good.
std::vector< std::vector< double > > GaussianBlur(std::vector< std::vector< double >> const &image) const
Applies Gaussian blur to image.
void SaveImage(TH2F *image, std::vector< art::PtrVector< recob::Hit >> const &allClusters, int pad, int tpc, int plane)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
Definition of data types for geometry description.
art::ServiceHandle< geo::Geometry const > fGeom
Encapsulate the geometry of a wire.
void ConvertBinsToClusters(std::vector< std::vector< double >> const &image, std::vector< std::vector< int >> const &allClusterBins, std::vector< art::PtrVector< recob::Hit >> &clusters) const
Takes a vector of clusters (itself a vector of hits) and turns them into clusters using the initial h...
std::vector< std::vector< double > > ConvertRecobHitsToVector(std::vector< art::Ptr< recob::Hit >> const &hits, int readoutWindowSize)
Takes hit map and returns a 2D vector representing wire and tick, filled with the charge...
T copy(T const &v)
then echo fcl name
std::vector< std::vector< std::vector< double > > > fAllKernels
temporary value
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
lariov::ChannelStatusProvider const & fChanStatus