LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
reco_tool::PeakFitterGaussian Class Reference
Inheritance diagram for reco_tool::PeakFitterGaussian:
reco_tool::IPeakFitter

Public Member Functions

 PeakFitterGaussian (const fhicl::ParameterSet &pset)
 
 ~PeakFitterGaussian ()
 
void configure (const fhicl::ParameterSet &pset) override
 
void findPeakParameters (const std::vector< float > &, const ICandidateHitFinder::HitCandidateVec &, PeakParamsVec &, double &, int &) const override
 

Private Types

using PeakFitParams_t = struct PeakFitParams{float peakCenter
 
using PeakParamsVec = std::vector< PeakFitParams_t >
 

Private Attributes

double fMinWidth
 minimum initial width for gaussian fit More...
 
double fMaxWidthMult
 multiplier for max width for gaussian fit More...
 
double fPeakRange
 set range limits for peak center More...
 
double fAmpRange
 set range limit for peak amplitude More...
 
bool fFloatBaseline
 Allow baseline to "float" away from zero. More...
 
BaselinedGausFitCache fFitCache
 Preallocated ROOT functions for the fits. More...
 
TH1F fHistogram
 
const geo::GeometryCorefGeometry = lar::providerFrom<geo::Geometry>()
 
float peakCenterError
 
float peakSigma
 
float peakSigmaError
 
float peakAmplitude
 
float peakAmplitudeError
 

Detailed Description

Definition at line 59 of file PeakFitterGaussian_tool.cc.

Member Typedef Documentation

using reco_tool::IPeakFitter::PeakFitParams_t = struct PeakFitParams { float peakCenter
inherited

Definition at line 31 of file IPeakFitter.h.

Definition at line 39 of file IPeakFitter.h.

Constructor & Destructor Documentation

reco_tool::PeakFitterGaussian::PeakFitterGaussian ( const fhicl::ParameterSet pset)
explicit

Definition at line 91 of file PeakFitterGaussian_tool.cc.

92 {
93  configure(pset);
94 }
void configure(const fhicl::ParameterSet &pset) override
reco_tool::PeakFitterGaussian::~PeakFitterGaussian ( )

Definition at line 96 of file PeakFitterGaussian_tool.cc.

97 {
98 }

Member Function Documentation

void reco_tool::PeakFitterGaussian::configure ( const fhicl::ParameterSet pset)
overridevirtual

Implements reco_tool::IPeakFitter.

Definition at line 100 of file PeakFitterGaussian_tool.cc.

References fhicl::ParameterSet::get().

101 {
102  // Start by recovering the parameters
103  fMinWidth = pset.get<double>("MinWidth", 0.5);
104  fMaxWidthMult = pset.get<double>("MaxWidthMult", 3.);
105  fPeakRange = pset.get<double>("PeakRangeFact", 2.);
106  fAmpRange = pset.get<double>("PeakAmpRange", 2.);
107  fFloatBaseline = pset.get< bool >("FloatBaseline", false);
108 
109  fHistogram = TH1F("PeakFitterHitSignal","",500,0.,500.);
110 
111  fHistogram.Sumw2();
112 
113  std::string function = "Gaus(0)";
114 
115  return;
116 }
bool fFloatBaseline
Allow baseline to "float" away from zero.
double fMaxWidthMult
multiplier for max width for gaussian fit
double fAmpRange
set range limit for peak amplitude
double fPeakRange
set range limits for peak center
T get(std::string const &key) const
Definition: ParameterSet.h:231
double fMinWidth
minimum initial width for gaussian fit
void reco_tool::PeakFitterGaussian::findPeakParameters ( const std::vector< float > &  roiSignalVec,
const ICandidateHitFinder::HitCandidateVec hitCandidateVec,
PeakParamsVec peakParamsVec,
double &  chi2PerNDF,
int &  NDF 
) const
overridevirtual

Implements reco_tool::IPeakFitter.

Definition at line 119 of file PeakFitterGaussian_tool.cc.

References DEFINE_ART_CLASS_TOOL, max, min, and util::flags::to_string().

124 {
125  // The following is a translation of the original FitGaussians function in the original
126  // GausHitFinder module originally authored by Jonathan Asaadi
127  //
128  // *** NOTE: this algorithm assumes the reference time for input hit candidates is to
129  // the first tick of the input waveform (ie 0)
130  //
131  if (hitCandidateVec.empty()) return;
132 
133  // in case of a fit failure, set the chi-square to infinity
134  chi2PerNDF = std::numeric_limits<double>::infinity();
135 
136  int startTime = hitCandidateVec.front().startTick;
137  int endTime = hitCandidateVec.back().stopTick;
138  int roiSize = endTime - startTime;
139 
140  // Check to see if we need a bigger histogram for fitting
141  if (roiSize > fHistogram.GetNbinsX())
142  {
143  std::string histName = "PeakFitterHitSignal_" + std::to_string(roiSize);
144  fHistogram = TH1F(histName.c_str(),"",roiSize,0.,roiSize);
145  fHistogram.Sumw2();
146  }
147 
148  for(int idx = 0; idx < roiSize; idx++) fHistogram.SetBinContent(idx+1,roiSignalVec.at(startTime+idx));
149 
150  // Build the string to describe the fit formula
151 #if 0
152  std::string equation = "gaus(0)";
153 
154  for(size_t idx = 1; idx < hitCandidateVec.size(); idx++) equation += "+gaus(" + std::to_string(3*idx) + ")";
155 
156  // Set the baseline if so desired
157  float baseline(0.);
158 
159  if (fFloatBaseline)
160  {
161  baseline = roiSignalVec.at(startTime);
162 
163  equation += "+" + std::to_string(baseline);
164  }
165 
166  // Now define the complete function to fit
167  TF1 Gaus("Gaus",equation.c_str(),0,roiSize,TF1::EAddToList::kNo);
168 #else
169  unsigned int const nGaus = hitCandidateVec.size();
170  assert(fFitCache.Get(nGaus));
171  TF1& Gaus = *(fFitCache.Get(nGaus));
172 
173  // Set the baseline if so desired
174  float const baseline = fFloatBaseline? roiSignalVec.at(startTime): 0.0;
175 
176  Gaus.FixParameter(nGaus * 3, baseline); // last parameter is the baseline
177 
178 #endif // 0
179 
180  // ### Setting the parameters for the Gaussian Fit ###
181  int parIdx{0};
182  for(auto const& candidateHit : hitCandidateVec)
183  {
184  double const peakMean = candidateHit.hitCenter - float(startTime);
185  double const peakWidth = candidateHit.hitSigma;
186  double const amplitude = candidateHit.hitHeight - baseline;
187  double const meanLowLim = std::max(peakMean - fPeakRange * peakWidth, 0.);
188  double const meanHiLim = std::min(peakMean + fPeakRange * peakWidth, double(roiSize));
189 
190  Gaus.SetParameter( parIdx, amplitude);
191  Gaus.SetParameter(1+parIdx, peakMean);
192  Gaus.SetParameter(2+parIdx, peakWidth);
193  Gaus.SetParLimits( parIdx, 0.1 * amplitude, fAmpRange * amplitude);
194  Gaus.SetParLimits(1+parIdx, meanLowLim, meanHiLim);
195  Gaus.SetParLimits(2+parIdx, std::max(fMinWidth, 0.1 * peakWidth), fMaxWidthMult * peakWidth);
196 
197  parIdx += 3;
198  }
199 
200  int fitResult{-1};
201 
202  try
203  { fitResult = fHistogram.Fit(&Gaus,"QNWB","", 0., roiSize);}
204  catch(...)
205  {mf::LogWarning("GausHitFinder") << "Fitter failed finding a hit";}
206 
207  // If the fit result is not zero there was an error
208  if (!fitResult)
209  {
210  // ##################################################
211  // ### Getting the fitted parameters from the fit ###
212  // ##################################################
213  chi2PerNDF = (Gaus.GetChisquare() / Gaus.GetNDF());
214  NDF = Gaus.GetNDF();
215 
216  int parIdx { 0 };
217  for(size_t idx = 0; idx < hitCandidateVec.size(); idx++)
218  {
219  PeakFitParams_t peakParams;
220 
221  peakParams.peakAmplitude = Gaus.GetParameter(parIdx);
222  peakParams.peakAmplitudeError = Gaus.GetParError( parIdx);
223  peakParams.peakCenter = Gaus.GetParameter(parIdx + 1) + float(startTime);
224  peakParams.peakCenterError = Gaus.GetParError( parIdx + 1);
225  peakParams.peakSigma = Gaus.GetParameter(parIdx + 2);
226  peakParams.peakSigmaError = Gaus.GetParError( parIdx + 2);
227 
228  peakParamsVec.emplace_back(peakParams);
229 
230  parIdx += 3;
231  }
232  }
233 #if 0
234  Gaus.Delete();
235 #endif // 0
236  return;
237 }
bool fFloatBaseline
Allow baseline to "float" away from zero.
double fMaxWidthMult
multiplier for max width for gaussian fit
double fAmpRange
set range limit for peak amplitude
double fPeakRange
set range limits for peak center
struct PeakFitParams{float peakCenter PeakFitParams_t
Definition: IPeakFitter.h:31
BaselinedGausFitCache fFitCache
Preallocated ROOT functions for the fits.
Int_t max
Definition: plot.C:27
virtual TF1 * Get(size_t nFunc)
Returns a function sum of nFunc base functions.
std::string to_string(Flag_t< Storage > const flag)
Convert a flag into a stream (shows its index).
Definition: BitMask.h:187
Int_t min
Definition: plot.C:26
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
double fMinWidth
minimum initial width for gaussian fit

Member Data Documentation

double reco_tool::PeakFitterGaussian::fAmpRange
private

set range limit for peak amplitude

Definition at line 79 of file PeakFitterGaussian_tool.cc.

BaselinedGausFitCache reco_tool::PeakFitterGaussian::fFitCache
mutableprivate

Preallocated ROOT functions for the fits.

Definition at line 82 of file PeakFitterGaussian_tool.cc.

bool reco_tool::PeakFitterGaussian::fFloatBaseline
private

Allow baseline to "float" away from zero.

Definition at line 80 of file PeakFitterGaussian_tool.cc.

const geo::GeometryCore* reco_tool::PeakFitterGaussian::fGeometry = lar::providerFrom<geo::Geometry>()
private

Definition at line 86 of file PeakFitterGaussian_tool.cc.

TH1F reco_tool::PeakFitterGaussian::fHistogram
mutableprivate

Definition at line 84 of file PeakFitterGaussian_tool.cc.

double reco_tool::PeakFitterGaussian::fMaxWidthMult
private

multiplier for max width for gaussian fit

Definition at line 77 of file PeakFitterGaussian_tool.cc.

double reco_tool::PeakFitterGaussian::fMinWidth
private

minimum initial width for gaussian fit

Definition at line 76 of file PeakFitterGaussian_tool.cc.

double reco_tool::PeakFitterGaussian::fPeakRange
private

set range limits for peak center

Definition at line 78 of file PeakFitterGaussian_tool.cc.

float reco_tool::IPeakFitter::peakAmplitude
inherited

Definition at line 35 of file IPeakFitter.h.

float reco_tool::IPeakFitter::peakAmplitudeError
inherited

Definition at line 36 of file IPeakFitter.h.

float reco_tool::IPeakFitter::peakCenterError
inherited

Definition at line 32 of file IPeakFitter.h.

float reco_tool::IPeakFitter::peakSigma
inherited

Definition at line 33 of file IPeakFitter.h.

float reco_tool::IPeakFitter::peakSigmaError
inherited

Definition at line 34 of file IPeakFitter.h.


The documentation for this class was generated from the following file: