LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
BlurredClusteringAlg.h
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 
12 #ifndef BlurredClustering_h
13 #define BlurredClustering_h
14 
15 // Framework includes
19 
20 // LArSoft includes
23 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
24 namespace detinfo {
25  class DetectorProperties;
26 }
27 namespace fhicl {
28  class ParameterSet;
29 }
30 namespace lariov {
31  class ChannelStatusProvider;
32 }
33 namespace geo {
34  struct WireID;
35 }
36 
37 // ROOT
38 #include "TString.h"
39 class TCanvas;
40 class TH2F;
41 
42 // c++
43 #include <array>
44 #include <string>
45 #include <vector>
46 
47 namespace cluster {
48  class BlurredClusteringAlg;
49 }
50 
52 public:
55 
57  void CreateDebugPDF(int run, int subrun, int event);
58 
60  void ConvertBinsToClusters(std::vector<std::vector<double>> const& image,
61  std::vector<std::vector<int>> const& allClusterBins,
62  std::vector<art::PtrVector<recob::Hit>>& clusters) const;
63 
65  std::vector<std::vector<double>> ConvertRecobHitsToVector(
67  int readoutWindowSize);
68 
70  int FindClusters(std::vector<std::vector<double>> const& image,
71  std::vector<std::vector<int>>& allcluster) const;
72 
74  int GlobalWire(geo::WireID const& wireID) const;
75 
77  std::vector<std::vector<double>> GaussianBlur(
78  std::vector<std::vector<double>> const& image) const;
79 
81  unsigned int GetMinSize() const noexcept { return fMinSize; }
82 
84  TH2F* MakeHistogram(std::vector<std::vector<double>> const& image, TString name) const;
85 
88  void SaveImage(TH2F* image,
89  std::vector<art::PtrVector<recob::Hit>> const& allClusters,
90  int pad,
91  int tpc,
92  int plane);
93 
95  void SaveImage(TH2F* image, int pad, int tpc, int plane);
96 
99  void SaveImage(TH2F* image,
100  std::vector<std::vector<int>> const& allClusterBins,
101  int pad,
102  int tpc,
103  int plane);
104 
105 private:
107  art::PtrVector<recob::Hit> ConvertBinsToRecobHits(std::vector<std::vector<double>> const& image,
108  std::vector<int> const& bins) const;
109 
111  art::Ptr<recob::Hit> ConvertBinToRecobHit(std::vector<std::vector<double>> const& image,
112  int bin) const;
113 
115  int ConvertWireTickToBin(std::vector<std::vector<double>> const& image, int xbin, int ybin) const;
116 
118  double ConvertBinToCharge(std::vector<std::vector<double>> const& image, int bin) const;
119 
122  std::pair<int, int> DeadWireCount(int wire_bin, int width) const;
123 
125  std::array<int, 4> FindBlurringParameters() const;
126 
128  double GetTimeOfBin(std::vector<std::vector<double>> const& image, int bin) const;
129 
131  std::vector<std::vector<std::vector<double>>> MakeKernels() const;
132 
134  unsigned int NumNeighbours(int nx, std::vector<bool> const& used, int bin) const;
135 
137  bool PassesTimeCut(std::vector<double> const& times, double time) const;
138 
139  bool fDebug;
140  std::string fDetector;
141 
142  // Parameters used in the Blurred Clustering algorithm
143  int fBlurWire; // blur radius for Gauss kernel in the wire direction
144  int fBlurTick; // blur radius for Gauss kernel in the tick direction
145  double fSigmaWire; // sigma for Gaussian kernel in the wire direction
146  double fSigmaTick; // sigma for Gaussian kernel in the tick direction
147  int fMaxTickWidthBlur; // maximum distance to blur a hit based on its natural width in time
148  int fClusterWireDistance; // how far to cluster from seed in wire direction
149  int fClusterTickDistance; // how far to cluster from seed in tick direction
150  unsigned int fNeighboursThreshold; // min. number of neighbors to add to cluster
151  unsigned int fMinNeighbours; // minumum number of neighbors to keep in the cluster
152  unsigned int fMinSize; // minimum size for cluster
153  double fMinSeed; // minimum seed after blurring needed before clustering proceeds
154  double fTimeThreshold; // time threshold for clustering
155  double fChargeThreshold; // charge threshold for clustering
156 
157  // Blurring stuff
158  int fKernelWidth, fKernelHeight;
159  std::vector<std::vector<std::vector<double>>> fAllKernels;
160 
161  // Hit containers
162  std::vector<std::vector<art::Ptr<recob::Hit>>> fHitMap;
163  std::vector<bool> fDeadWires;
164 
165  int fLowerTick, fUpperTick;
166  int fLowerWire, fUpperWire;
167 
168  // For the debug pdf
169  TCanvas* fDebugCanvas{nullptr};
170  std::string fDebugPDFName{};
171 
172  // art service handles
174  lariov::ChannelStatusProvider const& fChanStatus{
176 };
177 
178 #endif
std::vector< std::vector< art::Ptr< recob::Hit > > > fHitMap
Declaration of signal hit object.
Cluster finding and building.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
void hits()
Definition: readHits.C:15
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
float DeadWireCount(const TCSlice &slc, const TrajPoint &tp1, const TrajPoint &tp2)
Definition: Utils.cxx:2094
parameter set interface
General LArSoft Utilities.
float bin[41]
Definition: plottest35.C:14
art::ServiceHandle< geo::Geometry const > fGeom
Filters for channels, events, etc.
std::vector< std::vector< std::vector< double > > > fAllKernels
Namespace collecting geometry-related classes utilities.
art framework interface to geometry description
unsigned int GetMinSize() const noexcept
Minimum size of cluster to save.
Event finding and building.