PeakFitterGaussian_tool.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file PeakFitterGaussian.cc
3 /// \author T. Usher
4 ////////////////////////////////////////////////////////////////////////
5 
7 #include "larreco/RecoAlg/GausFitCache.h" // hit::GausFitCache
8 
11 #include "art_root_io/TFileService.h"
13 
14 #include <cassert>
15 #include <fstream>
16 
17 #include "TF1.h"
18 #include "TH1F.h"
19 
20 namespace reco_tool
21 {
22 
23 /// Customized function cache for Gaussians with a baseline.
24 ///
25 /// The baseline parameter is always the last one.
27 
28 public:
29  /// Constructor (see base class constructor).
30  BaselinedGausFitCache(std::string const& new_name="BaselinedGausFitCache")
31  : hit::GausFitCache(new_name)
32  {}
33 
34 protected:
35 
36  /// Creates and returns the function with specified number of Gaussians.
37  ///
38  /// The formula is `gaus(0) + gaus(3) + ... + gaus(3*(nFunc-1)) + [nFunc*3]`.
39  virtual TF1* CreateFunction(size_t nFunc) const
40  {
41  // add the Gaussians first
42  std::string formula;
43  std::size_t iGaus = 0;
44  while (iGaus < nFunc) formula += "gaus(" + std::to_string(3 * (iGaus++)) + ") + ";
45  formula += "[" + std::to_string(3 * nFunc) + "]";
46 
47  std::string const func_name = FunctionName(nFunc);
48  auto* pF = new TF1(func_name.c_str(), formula.c_str());
49  pF->SetParName(iGaus * 3, "baseline");
50  return pF;
51  } // CreateFunction()
52 
53 }; // BaselinedGausFitCache
54 
55 
57 {
58 public:
59  explicit PeakFitterGaussian(const fhicl::ParameterSet& pset);
60 
61  void findPeakParameters(const std::vector<float>&,
64  double&,
65  int&) const override;
66 
67 private:
68  // Member variables from the fhicl file
69  const double fMinWidth; ///< minimum initial width for gaussian fit
70  const double fMaxWidthMult; ///< multiplier for max width for gaussian fit
71  const double fPeakRange; ///< set range limits for peak center
72  const double fAmpRange; ///< set range limit for peak amplitude
73  const bool fFloatBaseline; ///< Allow baseline to "float" away from zero
74  const bool fOutputHistograms; ///< If true will generate summary style histograms
75 
77  TH1F* fROISizeHist;
86 
87  mutable BaselinedGausFitCache fFitCache; ///< Preallocated ROOT functions for the fits.
88 
89  mutable TH1F fHistogram;
90 
91  const geo::GeometryCore* fGeometry = lar::providerFrom<geo::Geometry>();
92 };
93 
94 //----------------------------------------------------------------------
95 // Constructor.
97  fMinWidth(pset.get<double>("MinWidth", 0.5)),
98  fMaxWidthMult (pset.get<double>("MaxWidthMult", 3.)),
99  fPeakRange(pset.get<double>("PeakRangeFact", 2.)),
100  fAmpRange(pset.get<double>("PeakAmpRange", 2.)),
101  fFloatBaseline(pset.get< bool >("FloatBaseline", false)),
102  fOutputHistograms(pset.get< bool >("OutputHistograms", false))
103 {
104  fHistogram = TH1F("PeakFitterHitSignal","",500,0.,500.);
105 
106  fHistogram.Sumw2();
107 
108  std::string function = "Gaus(0)";
109 
110  // If asked, define the global histograms
111  if (fOutputHistograms)
112  {
113  // Access ART's TFileService, which will handle creating and writing
114  // histograms and n-tuples for us.
116 
117  // Make a directory for these histograms
118  art::TFileDirectory dir = tfs->mkdir("PeakFit");
119 
120  fNumCandHitsHist = dir.make<TH1F>("NumCandHits", "# Candidate Hits", 100, 0., 100.);
121  fROISizeHist = dir.make<TH1F>("ROISize", "ROI Size", 400, 0., 400.);
122  fCandPeakPositionHist = dir.make<TH1F>("CPeakPosition", "Peak Position", 200, 0., 400.);
123  fCandPeakWidHist = dir.make<TH1F>("CPeadWidth", "Peak Width", 100, 0., 25.);
124  fCandPeakAmpitudeHist = dir.make<TH1F>("CPeakAmplitude", "Peak Amplitude", 100, 0., 200.);
125  fCandBaselineHist = dir.make<TH1F>("CBaseline", "Baseline", 200, -25., 25.);
126  fFitPeakPositionHist = dir.make<TH1F>("FPeakPosition", "Peak Position", 200, 0., 400.);
127  fFitPeakWidHist = dir.make<TH1F>("FPeadWidth", "Peak Width", 100, 0., 25.);
128  fFitPeakAmpitudeHist = dir.make<TH1F>("FPeakAmplitude", "Peak Amplitude", 100, 0., 200.);
129  fFitBaselineHist = dir.make<TH1F>("FBaseline", "Baseline", 200, -25., 25.);
130  }
131 
132  return;
133 }
134 
135 // --------------------------------------------------------------------------------------------
136 void PeakFitterGaussian::findPeakParameters(const std::vector<float>& roiSignalVec,
137  const ICandidateHitFinder::HitCandidateVec& hitCandidateVec,
138  PeakParamsVec& peakParamsVec,
139  double& chi2PerNDF,
140  int& NDF) const
141 {
142  // The following is a translation of the original FitGaussians function in the original
143  // GausHitFinder module originally authored by Jonathan Asaadi
144  //
145  // *** NOTE: this algorithm assumes the reference time for input hit candidates is to
146  // the first tick of the input waveform (ie 0)
147  //
148  if (hitCandidateVec.empty()) return;
149 
150  // in case of a fit failure, set the chi-square to infinity
151  chi2PerNDF = std::numeric_limits<double>::infinity();
152 
153  int startTime = hitCandidateVec.front().startTick;
154  int endTime = hitCandidateVec.back().stopTick;
155  int roiSize = endTime - startTime;
156 
157  // Check to see if we need a bigger histogram for fitting
158  if (roiSize > fHistogram.GetNbinsX())
159  {
160  std::string histName = "PeakFitterHitSignal_" + std::to_string(roiSize);
161  fHistogram = TH1F(histName.c_str(),"",roiSize,0.,roiSize);
162  fHistogram.Sumw2();
163  }
164 
165  for(int idx = 0; idx < roiSize; idx++) fHistogram.SetBinContent(idx+1,roiSignalVec[startTime+idx]);
166 
167  // Build the string to describe the fit formula
168 #if 0
169  std::string equation = "gaus(0)";
170 
171  for(size_t idx = 1; idx < hitCandidateVec.size(); idx++) equation += "+gaus(" + std::to_string(3*idx) + ")";
172 
173  // Set the baseline if so desired
174  float baseline(0.);
175 
176  if (fFloatBaseline)
177  {
178  baseline = roiSignalVec[startTime];
179 
180  equation += "+" + std::to_string(baseline);
181  }
182 
183  // Now define the complete function to fit
184  TF1 Gaus("Gaus",equation.c_str(),0,roiSize,TF1::EAddToList::kNo);
185 #else
186  unsigned int const nGaus = hitCandidateVec.size();
187  assert(fFitCache.Get(nGaus));
188  TF1& Gaus = *(fFitCache.Get(nGaus));
189 
190  // Set the baseline if so desired
191  float baseline(0.);
192 
193  if (fFloatBaseline)
194  {
195  baseline = roiSignalVec[startTime];
196  Gaus.SetParameter(nGaus * 3, baseline);
197  Gaus.SetParLimits( nGaus * 3, baseline - 12., baseline + 12.);
198  }
199  else Gaus.FixParameter(nGaus * 3, baseline); // last parameter is the baseline
200 
201 #endif // 0
202 
203  if (fOutputHistograms)
204  {
205  fNumCandHitsHist->Fill(hitCandidateVec.size(), 1.);
206  fROISizeHist->Fill(roiSize, 1.);
207  fCandBaselineHist->Fill(baseline, 1.);
208  }
209 
210  // ### Setting the parameters for the Gaussian Fit ###
211  int parIdx{0};
212  for(auto const& candidateHit : hitCandidateVec)
213  {
214  double const peakMean = candidateHit.hitCenter - float(startTime);
215  double const peakWidth = candidateHit.hitSigma;
216  double const amplitude = candidateHit.hitHeight - baseline;
217  double const meanLowLim = std::max(peakMean - fPeakRange * peakWidth, 0.);
218  double const meanHiLim = std::min(peakMean + fPeakRange * peakWidth, double(roiSize));
219 
220  if (fOutputHistograms)
221  {
222  fCandPeakPositionHist->Fill(peakMean, 1.);
223  fCandPeakWidHist->Fill(peakWidth, 1.);
224  fCandPeakAmpitudeHist->Fill(amplitude, 1.);
225  }
226 
227  Gaus.SetParameter( parIdx, amplitude);
228  Gaus.SetParameter(1+parIdx, peakMean);
229  Gaus.SetParameter(2+parIdx, peakWidth);
230  Gaus.SetParLimits( parIdx, 0.1 * amplitude, fAmpRange * amplitude);
231  Gaus.SetParLimits(1+parIdx, meanLowLim, meanHiLim);
232  Gaus.SetParLimits(2+parIdx, std::max(fMinWidth, 0.1 * peakWidth), fMaxWidthMult * peakWidth);
233 
234  parIdx += 3;
235  }
236 
237  int fitResult{-1};
238 
239  try
240  { fitResult = fHistogram.Fit(&Gaus,"QNWB","", 0., roiSize);}
241  catch(...)
242  {mf::LogWarning("GausHitFinder") << "Fitter failed finding a hit";}
243 
244  // If the fit result is not zero there was an error
245  if (!fitResult)
246  {
247  // ##################################################
248  // ### Getting the fitted parameters from the fit ###
249  // ##################################################
250  chi2PerNDF = (Gaus.GetChisquare() / Gaus.GetNDF());
251  NDF = Gaus.GetNDF();
252 
253  int parIdx { 0 };
254  for(size_t idx = 0; idx < hitCandidateVec.size(); idx++)
255  {
256  PeakFitParams_t peakParams;
257 
258  peakParams.peakAmplitude = Gaus.GetParameter(parIdx);
259  peakParams.peakAmplitudeError = Gaus.GetParError( parIdx);
260  peakParams.peakCenter = Gaus.GetParameter(parIdx + 1) + float(startTime);
261  peakParams.peakCenterError = Gaus.GetParError( parIdx + 1);
262  peakParams.peakSigma = Gaus.GetParameter(parIdx + 2);
263  peakParams.peakSigmaError = Gaus.GetParError( parIdx + 2);
264 
265  if (fOutputHistograms)
266  {
267  fFitPeakPositionHist->Fill(peakParams.peakCenter, 1.);
268  fFitPeakWidHist->Fill(peakParams.peakSigma, 1.);
269  fFitPeakAmpitudeHist->Fill(peakParams.peakAmplitude, 1.);
270  }
271 
272  peakParamsVec.emplace_back(peakParams);
273 
274  parIdx += 3;
275  }
276 
277  if (fOutputHistograms) fFitBaselineHist->Fill(Gaus.GetParameter(3*nGaus), 1.);
278  }
279 #if 0
280  Gaus.Delete();
281 #endif // 0
282  return;
283 }
284 
286 }
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:42
std::string string
Definition: nybbler.cc:12
const double fPeakRange
set range limits for peak center
const double fMaxWidthMult
multiplier for max width for gaussian fit
const bool fOutputHistograms
If true will generate summary style histograms.
string dir
PeakFitterGaussian(const fhicl::ParameterSet &pset)
art framework interface to geometry description
A set of TF1 linear sum of base functions (Gaussians)
Definition: GausFitCache.h:46
const double fAmpRange
set range limit for peak amplitude
BaselinedGausFitCache fFitCache
Preallocated ROOT functions for the fits.
BaselinedGausFitCache(std::string const &new_name="BaselinedGausFitCache")
Constructor (see base class constructor).
virtual TF1 * Get(size_t nFunc)
Returns a function sum of nFunc base functions.
static int max(int a, int b)
Description of geometry of one entire detector.
GausFitCache(std::string new_name="GausFitCache")
Constructor; optionally set the name of the repository.
Definition: GausFitCache.h:49
Detector simulation of raw signals on wires.
virtual TF1 * CreateFunction(size_t nFunc) const
std::vector< PeakFitParams_t > PeakParamsVec
Definition: IPeakFitter.h:37
const bool fFloatBaseline
Allow baseline to "float" away from zero.
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
This provides an interface for tools which are tasked with fitting peaks on input waveforms...
virtual std::string FunctionName(size_t nFunc) const
Returns a name for the function with nFunc base functions.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
Provide caches for TF1 functions to be used with ROOT fitters.
auto const & get(AssnsNode< L, R, D > const &r)
Definition: AssnsNode.h:115
int bool
Definition: qglobal.h:345
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
void findPeakParameters(const std::vector< float > &, const ICandidateHitFinder::HitCandidateVec &, PeakParamsVec &, double &, int &) const override
const double fMinWidth
minimum initial width for gaussian fit
std::vector< HitCandidate > HitCandidateVec