LArSoft  v07_13_02
Liquid Argon Software toolkit - http://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
23 
24 // LArSoft includes
25 #include "larcore/CoreUtils/ServiceUtil.h" // lar::providerFrom<>()
27 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
28 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
35 
36 // ROOT
37 #include <TTree.h>
38 #include <TH2F.h>
39 #include <TH2.h>
40 #include <TCanvas.h>
41 #include <TCutG.h>
42 #include <TString.h>
43 #include <TMarker.h>
44 #include <TColor.h>
45 #include <TCanvas.h>
46 #include <TStyle.h>
47 #include <TVirtualPad.h>
48 #include <TLatex.h>
49 #include <TGraph.h>
50 #include <TF1.h>
51 #include <TLine.h>
52 #include <TPrincipal.h>
53 #include <TMath.h>
54 #include <TVector.h>
55 #include <TVectorD.h>
56 #include <TVector2.h>
57 
58 // c++
59 #include <array>
60 #include <string>
61 #include <vector>
62 #include <map>
63 #include <sstream>
64 
65 
66 namespace cluster {
67  class BlurredClusteringAlg;
68 }
69 
71 public:
72 
75 
77  void CreateDebugPDF(int run, int subrun, int event);
78 
80  void ConvertBinsToClusters(std::vector<std::vector<double>> const& image,
81  std::vector<std::vector<int>> const& allClusterBins,
82  std::vector<art::PtrVector<recob::Hit>>& clusters) const;
83 
85  std::vector<std::vector<double>> ConvertRecobHitsToVector(std::vector<art::Ptr<recob::Hit>> const& hits);
86 
88  int FindClusters(std::vector<std::vector<double>> const& image, std::vector<std::vector<int>>& allcluster);
89 
91  int GlobalWire(geo::WireID const& wireID);
92 
94  std::vector<std::vector<double>> GaussianBlur(std::vector<std::vector<double>> const& image);
95 
97  unsigned int GetMinSize() const noexcept { return fMinSize; }
98 
100  TH2F* MakeHistogram(std::vector<std::vector<double>> const& image, TString name);
101 
104  void SaveImage(TH2F* image, std::vector<art::PtrVector<recob::Hit>> const& allClusters, int pad, int tpc, int plane);
105 
107  void SaveImage(TH2F* image, int pad, int tpc, int plane);
108 
111  void SaveImage(TH2F* image, std::vector<std::vector<int>> const& allClusterBins, int pad, int tpc, int plane);
112 
113 private:
114 
116  art::PtrVector<recob::Hit> ConvertBinsToRecobHits(std::vector<std::vector<double>> const& image, std::vector<int> const& bins) const;
117 
119  art::Ptr<recob::Hit> ConvertBinToRecobHit(std::vector<std::vector<double>> const& image, int bin) const;
120 
122  int ConvertWireTickToBin(std::vector<std::vector<double>> const& image, int xbin, int ybin) const;
123 
125  double ConvertBinToCharge(std::vector<std::vector<double>> const& image, int bin) const;
126 
129  std::pair<int, int> DeadWireCount(int wire_bin, int width) const;
130 
132  std::array<int, 4> FindBlurringParameters() const;
133 
135  double GetTimeOfBin(std::vector<std::vector<double>> const& image, int bin) const;
136 
138  std::vector<std::vector<std::vector<double>>> MakeKernels() const;
139 
141  unsigned int NumNeighbours(int nx, std::vector<bool> const& used, int bin) const;
142 
144  bool PassesTimeCut(std::vector<double> const& times, double time) const;
145 
146  bool fDebug;
147  std::string fDetector;
148 
149  // Parameters used in the Blurred Clustering algorithm
150  int fBlurWire; // blur radius for Gauss kernel in the wire direction
151  int fBlurTick; // blur radius for Gauss kernel in the tick direction
152  double fSigmaWire; // sigma for Gaussian kernel in the wire direction
153  double fSigmaTick; // sigma for Gaussian kernel in the tick direction
154  int fMaxTickWidthBlur; // maximum distance to blur a hit based on its natural width in time
155  int fClusterWireDistance; // how far to cluster from seed in wire direction
156  int fClusterTickDistance; // how far to cluster from seed in tick direction
157  unsigned int fNeighboursThreshold; // min. number of neighbors to add to cluster
158  unsigned int fMinNeighbours; // minumum number of neighbors to keep in the cluster
159  unsigned int fMinSize; // minimum size for cluster
160  double fMinSeed; // minimum seed after blurring needed before clustering proceeds
161  double fTimeThreshold; // time threshold for clustering
162  double fChargeThreshold; // charge threshold for clustering
163 
164  // Blurring stuff
166  std::vector<std::vector<std::vector<double>>> fAllKernels;
167 
168  // Hit containers
169  std::vector<std::vector<art::Ptr<recob::Hit>>> fHitMap;
170  std::vector<bool> fDeadWires;
171 
174 
175  // For the debug pdf
176  TCanvas* fDebugCanvas{nullptr};
177  std::string fDebugPDFName{};
178 
179  // art service handles
182  lariov::ChannelStatusProvider const& fChanStatus{art::ServiceHandle<lariov::ChannelStatusService>()->GetProvider()};
183 
184 };
185 
186 #endif
BlurredClusteringAlg(fhicl::ParameterSet const &pset)
Utilities related to art service access.
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.
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.
Declaration of signal hit object.
unsigned int NumNeighbours(int nx, std::vector< bool > const &used, int bin) const
Determines the number of clustered neighbours of a hit.
std::vector< std::vector< double > > ConvertRecobHitsToVector(std::vector< art::Ptr< recob::Hit >> const &hits)
Takes hit map and returns a 2D vector representing wire and tick, filled with the charge...
std::array< int, 4 > FindBlurringParameters() const
Dynamically find the blurring radii and Gaussian sigma in each dimension.
Cluster finding and building.
art::ServiceHandle< geo::Geometry > fGeom
int FindClusters(std::vector< std::vector< double >> const &image, std::vector< std::vector< int >> &allcluster)
Find clusters in the histogram.
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...
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
Provides recob::Track data product.
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.
TH2F * MakeHistogram(std::vector< std::vector< double >> const &image, TString name)
Converts a 2D vector in a histogram for the debug pdf.
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...
void SaveImage(TH2F *image, std::vector< art::PtrVector< recob::Hit >> const &allClusters, int pad, int tpc, int plane)
float bin[41]
Definition: plottest35.C:14
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...
Encapsulate the construction of a single detector plane.
std::vector< std::vector< double > > GaussianBlur(std::vector< std::vector< double >> const &image)
Applies Gaussian blur to image.
std::vector< std::vector< std::vector< double > > > fAllKernels
lariov::ChannelStatusProvider const & fChanStatus
detinfo::DetectorProperties const * fDetProp
art framework interface to geometry description
unsigned int GetMinSize() const noexcept
Minimum size of cluster to save.
Event finding and building.
int GlobalWire(geo::WireID const &wireID)
Find the global wire position.