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