LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
hit::CCHitFinderAlg Class Reference

Hit finder algorithm designed to work with Cluster Crawler. More...

#include "CCHitFinderAlg.h"

Classes

struct  FitStats_t
 
class  HitChannelInfo_t
 exchange data about the originating wire More...
 

Public Member Functions

 CCHitFinderAlg (fhicl::ParameterSet const &pset)
 
virtual ~CCHitFinderAlg ()=default
 
virtual void reconfigure (fhicl::ParameterSet const &pset)
 
void RunCCHitFinder (std::vector< recob::Wire > const &Wires)
 
std::vector< recob::Hit > && YieldHits ()
 Returns (and loses) the collection of reconstructed hits. More...
 
template<typename Stream >
void PrintStats (Stream &out) const
 Print the fit statistics. More...
 

Public Attributes

std::vector< recob::Hitallhits
 

Private Member Functions

void FitNG (unsigned short nGaus, unsigned short npt, float *ticks, float *signl)
 
void MakeCrudeHit (unsigned short npt, float *ticks, float *signl)
 
void StoreHits (unsigned short TStart, unsigned short npt, HitChannelInfo_t info, float adcsum)
 
void StudyHits (unsigned short flag, unsigned short npt=0, float *ticks=0, float *signl=0, unsigned short tstart=0)
 

Static Private Member Functions

static bool FastGaussianFit (unsigned short npt, float const *ticks, float const *signl, std::array< double, 3 > &params, std::array< double, 3 > &paramerrors, float &chidof)
 Performs a "fast" fit. More...
 

Private Attributes

std::vector< float > fMinPeak
 
std::vector< float > fMinRMS
 
unsigned short fMaxBumps
 
unsigned short fMaxXtraHits
 
float fChiSplit
 Estimated noise error on the Signal. More...
 
std::vector< float > fChiNorms
 
std::vector< float > fTimeOffsets
 
std::vector< float > fChgNorms
 
raw::ChannelID_t theChannel
 
unsigned short theWireNum
 
unsigned short thePlane
 
float chinorm
 
bool fUseChannelFilter
 
art::ServiceHandle< geo::Geometrygeom
 
std::vector< double > par
 
std::vector< double > parerr
 
std::vector< double > parmin
 
std::vector< double > parmax
 
float chidof
 
int dof
 
std::vector< unsigned short > bumps
 
bool fStudyHits
 
std::vector< short > fUWireRange
 
std::vector< short > fUTickRange
 
std::vector< short > fVWireRange
 
std::vector< short > fVTickRange
 
std::vector< short > fWWireRange
 
std::vector< short > fWTickRange
 
std::vector< int > bumpCnt
 
std::vector< int > RATCnt
 
std::vector< float > bumpChi
 
std::vector< float > bumpRMS
 
std::vector< int > hitCnt
 
std::vector< float > hitRMS
 
std::vector< float > loWire
 
std::vector< float > loTime
 
std::vector< float > hiWire
 
std::vector< float > hiTime
 
bool SelRAT
 
bool fUseFastFit
 whether to attempt using a fast fit on single gauss. More...
 
std::unique_ptr< GausFitCacheFitCache
 a set of functions ready to be used More...
 
FitStats_t FinalFitStats
 counts of the good fits More...
 
FitStats_t TriedFitStats
 counts of the tried fits More...
 

Static Private Attributes

static constexpr float Sqrt2Pi = 2.5066
 
static constexpr float SqrtPi = 1.7725
 
static constexpr unsigned int MaxGaussians = 20
 

Detailed Description

Hit finder algorithm designed to work with Cluster Crawler.

This algorithm used to store hits in a proprietary CCHit data structure. It has now been changed to use recob::Hit class directly. It is possible to translate the former into the latter, with one exception, as follows:

// this is the original CCHit definition
struct CCHit {
  float Charge;            // recob::Hit::Integral()
  float ChargeErr;         // recob::Hit::SigmaIntegral()
  float Amplitude;         // recob::Hit::PeakAmplitude()
  float AmplitudeErr;      // recob::Hit::SigmaPeakAmplitude()
  float Time;              // recob::Hit::PeakTime()
  float TimeErr;           // recob::Hit::SigmaPeakTime()
  float RMS;               // recob::Hit::RMS()
  float RMSErr;            // dropped
  float ChiDOF;            // recob::Hit::GoodnessOfFit()
  int   DOF;               // recob::Hit::DegreesOfFreedom()
  float ADCSum;            // recob::Hit::SummedADC()
  unsigned short WireNum;  // recob::Hit::WireID().Wire
  unsigned short numHits;  // recob::Hit::Multiplicity()
  unsigned int LoHitID;    // see below
  float LoTime;            // recob::Hit::StartTick()
  float HiTime;            // recob::Hit::EndTick()
  short InClus;            // dropped; see below
  geo::WireID WirID;       // recob::Hit::WireID()
  recob::Wire const* Wire; // dropped; see below
};

The uncertainty on RMS has been dropped for good.

The LoHitID member used to mean the index of the first hit in the "hit train" (that is the set of hits extracted from the same region of interest). That is a concept that is not portable. If your hit list is still the original one as produced by this algorithm, or if at least the hits from the same train are stored sorted and contiguously, for a hit with index iHit, the equivalent value of LoHitID is iHit - hit.LocalIndex().

There is no pointer to the wire any more in recob::Hit. The wire can be obtained through associations, that are typically produced by the art module that runs CCHitFinderAlg (e.g. CCHitFinder). The channel ID is also directly available as recob::Hit::Channel().

Definition at line 82 of file CCHitFinderAlg.h.

Constructor & Destructor Documentation

hit::CCHitFinderAlg::CCHitFinderAlg ( fhicl::ParameterSet const &  pset)

Definition at line 51 of file CCHitFinderAlg.cxx.

References reconfigure().

51  :
52  FitCache(new
53  GausFitCache // run-time, on demand TFormula cache
54  // CompiledGausFitCache<MaxGaussians> // precompiled Gaussian set
55  // CompiledTruncatedGausFitCache<MaxGaussians, 4> // precompiled truncated Gaussian set
56  // CompiledTruncatedGausFitCache<MaxGaussians, 5> // precompiled truncated Gaussian set
57  ("GausFitCache_CCHitFinderAlg")
58  )
59  {
60  this->reconfigure(pset);
61  }
virtual void reconfigure(fhicl::ParameterSet const &pset)
std::unique_ptr< GausFitCache > FitCache
a set of functions ready to be used
virtual hit::CCHitFinderAlg::~CCHitFinderAlg ( )
virtualdefault

Member Function Documentation

bool hit::CCHitFinderAlg::FastGaussianFit ( unsigned short  npt,
float const *  ticks,
float const *  signl,
std::array< double, 3 > &  params,
std::array< double, 3 > &  paramerrors,
float &  chidof 
)
staticprivate

Performs a "fast" fit.

Parameters
nptnumber of points to be fitted
tickstick coordinates
signlsignal amplitude
paramsan array where the fit parameters will be stored
paramerrorsan array where the fit parameter errors will be stored
chidofa variable where to store chi^2 over degrees of freedom
Returns
whether the fit was successful or not

Note that the fit will bail out and rteurn false if any of the input signal amplitudes is zero or negative.

Also note that currently the chi^2 is not the one from comparing the Gaussian to the signal, but from comparing a fitted parabola to the logarithm of the signal.

Definition at line 282 of file CCHitFinderAlg.cxx.

References lar::util::GaussianFit< T >::add(), lar::util::GaussianFit< T >::ChiSquare(), lar::util::GaussianFit< T >::FillResults(), LOG_DEBUG, lar::util::GaussianFit< T >::NDF(), and par.

Referenced by FitNG().

287  {
288  // parameters: amplitude, mean, sigma
289 
290  lar::util::GaussianFit<double> fitter; // probably "double" is overdoing
291 
292  // apply a time shift so that the center of the time interval is 0
293  const float time_shift = (ticks[npt - 1] + ticks[0]) / 2.F;
294 
295  // fill the input data, no uncertainty
296  for (size_t i = 0; i < npt; ++i) {
297  if (signl[i] <= 0) {
298  LOG_DEBUG("CCHitFinderAlg")
299  << "Non-positive charge encountered. Backing up to ROOT fit.";
300  return false;
301  }
302  // we could freely add a Poisson uncertainty (as third parameter)
303  fitter.add(ticks[i] - time_shift, signl[i]);
304  } // for
305 
306  // we might have found that we don't want the fast fit after all...
307  if (!fitter.FillResults(params, paramerrors)) {
308  // something went wrong...
309  LOG_DEBUG("CCHitFinderAlg") << "Fast Gaussian fit failed.";
310  return false;
311  }
312 
313  // note that this is not the full chi^2, but it is the chi^2 of the
314  // parabolic fit underlying the Gaussian one
315  const double chi2 = fitter.ChiSquare();
316  chidof = chi2 / fitter.NDF();
317 
318  // remove the time shift
319  params[1] += time_shift; // mean
320 
321  // GP: inflate the uncertainties on the fit parameters according to chi2/NDF
322  // (not sure if this is in any way justified)
323  if (chidof > 1.)
324  for (double& par: paramerrors) par *= std::sqrt(chidof);
325 
326  return true;
327  } // FastGaussianFit()
bool add(Data_t x, Data_t y, Data_t sy=Data_t(1.0))
Clears all the input statistics.
Definition: SimpleFits.h:1790
virtual bool FillResults(FitParameters_t &params, FitMatrix_t &Xmat, Data_t &det, FitMatrix_t &Smat) const override
Fills the specified parameters.
Definition: SimpleFits.h:1888
"Fast" Gaussian fit
Definition: SimpleFits.h:1018
virtual int NDF() const override
Returns the degrees of freedom in the determination of the fit.
Definition: SimpleFits.h:1170
#define LOG_DEBUG(id)
virtual Data_t ChiSquare() const override
Returns the of the original fit.
Definition: SimpleFits.h:1161
std::vector< double > par
void hit::CCHitFinderAlg::FitNG ( unsigned short  nGaus,
unsigned short  npt,
float *  ticks,
float *  signl 
)
private

Definition at line 331 of file CCHitFinderAlg.cxx.

References hit::CCHitFinderAlg::FitStats_t::AddFast(), hit::CCHitFinderAlg::FitStats_t::AddMultiGaus(), bumps, chidof, chinorm, dof, FastGaussianFit(), FinalFitStats, fMinPeak, fMinRMS, fUseFastFit, par, parerr, thePlane, and TriedFitStats.

Referenced by RunCCHitFinder().

333  {
334  // Fit the signal to n Gaussians
335 
336  dof = npt - 3 * nGaus;
337 
338  chidof = 9999.;
339 
340  if(dof < 3) return;
341  if(bumps.size() == 0) return;
342 
343  // load the fit into a temp vector
344  std::vector<double> partmp;
345  std::vector<double> partmperr;
346 
347  //
348  // if it is possible, we try first with the quick single Gaussian fit
349  //
351 
352  bool bNeedROOTfit = (nGaus > 1) || !fUseFastFit;
353  if (!bNeedROOTfit) {
354  // so, we need only one puny Gaussian;
355  std::array<double, 3> params, paramerrors;
356 
358 
359  if (FastGaussianFit(npt, ticks, signl, params, paramerrors, chidof)) {
360  // success? copy the results in the proper structures
361  partmp.resize(3);
362  std::copy(params.begin(), params.end(), partmp.begin());
363  partmperr.resize(3);
364  std::copy(paramerrors.begin(), paramerrors.end(), partmperr.begin());
365  }
366  else bNeedROOTfit = true; // if we fail, let's schedule ROOT to back us up
367 
368  if (!bNeedROOTfit) FinalFitStats.AddFast();
369 
370  } // if we don't need ROOT to fit
371 
372  if (bNeedROOTfit) {
373  // we may land here either because the simple Gaussian fit did not work
374  // (either failed, or we chose not to trust it)
375  // or because the fit is multi-Gaussian
376 
377  // define the fit string to pass to TF1
378 
379  std::stringstream numConv;
380  std::string eqn = "gaus";
381  if(nGaus > 1) eqn = "gaus(0)";
382  for(unsigned short ii = 3; ii < nGaus*3; ii+=3){
383  eqn.append(" + gaus(");
384  numConv.str("");
385  numConv << ii;
386  eqn.append(numConv.str());
387  eqn.append(")");
388  }
389 
390  std::unique_ptr<TF1> Gn(new TF1("gn",eqn.c_str()));
391  /*
392  TF1* Gn = FitCache->Get(nGaus);
393  */
394  TGraph *fitn = new TGraph(npt, ticks, signl);
395  /*
396  if(prt) mf::LogVerbatim("CCHitFinder")
397  <<"FitNG nGaus "<<nGaus<<" nBumps "<<bumps.size();
398  */
399  // put in the bump parameters. Assume that nGaus >= bumps.size()
400  for(unsigned short ii = 0; ii < bumps.size(); ++ii) {
401  unsigned short index = ii * 3;
402  unsigned short bumptime = bumps[ii];
403  double amp = signl[bumptime];
404  Gn->SetParameter(index , amp);
405  Gn->SetParLimits(index, 0., 9999.);
406  Gn->SetParameter(index + 1, (double)bumptime);
407  Gn->SetParLimits(index + 1, 0, (double)npt);
408  Gn->SetParameter(index + 2, (double)fMinRMS[thePlane]);
409  Gn->SetParLimits(index + 2, 1., 3*(double)fMinRMS[thePlane]);
410  /*
411  if(prt) mf::LogVerbatim("CCHitFinder")<<"Bump params "<<ii<<" "<<(short)amp
412  <<" "<<(int)bumptime<<" "<<(int)fMinRMS[thePlane];
413  */
414  } // ii bumps
415 
416  // search for other bumps that may be hidden by the already found ones
417  for(unsigned short ii = bumps.size(); ii < nGaus; ++ii) {
418  // bump height must exceed fMinPeak
419  float big = fMinPeak[thePlane];
420  unsigned short imbig = 0;
421  for(unsigned short jj = 0; jj < npt; ++jj) {
422  float diff = signl[jj] - Gn->Eval((Double_t)jj, 0, 0, 0);
423  if(diff > big) {
424  big = diff;
425  imbig = jj;
426  }
427  } // jj
428  if(imbig > 0) {
429  /*
430  if(prt) mf::LogVerbatim("CCHitFinder")<<"Found bump "<<ii<<" "<<(short)big
431  <<" "<<imbig;
432  */
433  // set the parameters for the bump
434  unsigned short index = ii * 3;
435  Gn->SetParameter(index , (double)big);
436  Gn->SetParLimits(index, 0., 9999.);
437  Gn->SetParameter(index + 1, (double)imbig);
438  Gn->SetParLimits(index + 1, 0, (double)npt);
439  Gn->SetParameter(index + 2, (double)fMinRMS[thePlane]);
440  Gn->SetParLimits(index + 2, 1., 5*(double)fMinRMS[thePlane]);
441  } // imbig > 0
442  } // ii
443 
444  // W = set weights to 1, N = no drawing or storing, Q = quiet
445  // B = bounded parameters
446  fitn->Fit(&*Gn,"WNQB");
447 
448  for(unsigned short ipar = 0; ipar < 3 * nGaus; ++ipar) {
449  partmp.push_back(Gn->GetParameter(ipar));
450  partmperr.push_back(Gn->GetParError(ipar));
451  }
452  chidof = Gn->GetChisquare() / ( dof * chinorm);
453 
454  delete fitn;
455  // delete Gn;
456 
457  } // if ROOT fit
458 
459  // Sort by increasing time if necessary
460  if(nGaus > 1) {
461  std::vector< std::pair<unsigned short, unsigned short> > times;
462  // fill the sort vector
463  for(unsigned short ii = 0; ii < nGaus; ++ii) {
464  unsigned short index = ii * 3;
465  times.push_back(std::make_pair(partmp[index + 1],ii));
466  } // ii
467  std::sort(times.begin(), times.end());
468  // see if re-arranging is necessary
469  bool sortem = false;
470  for(unsigned short ii = 0; ii < nGaus; ++ii) {
471  if(times[ii].second != ii) {
472  sortem = true;
473  break;
474  }
475  } // ii
476  if(sortem) {
477  // temp temp vectors for putting things in the right time order
478  std::vector<double> partmpt;
479  std::vector<double> partmperrt;
480  for(unsigned short ii = 0; ii < nGaus; ++ii) {
481  unsigned short index = times[ii].second * 3;
482  partmpt.push_back(partmp[index]);
483  partmpt.push_back(partmp[index+1]);
484  partmpt.push_back(partmp[index+2]);
485  partmperrt.push_back(partmperr[index]);
486  partmperrt.push_back(partmperr[index+1]);
487  partmperrt.push_back(partmperr[index+2]);
488  } // ii
489  partmp = partmpt;
490  partmperr = partmperrt;
491  } // sortem
492  } // nGaus > 1
493 /*
494  if(prt) {
495  mf::LogVerbatim("CCHitFinder")<<"Fit "<<nGaus<<" chi "<<chidof
496  <<" npars "<<partmp.size();
497  mf::LogVerbatim("CCHitFinder")<<"pars errs ";
498  for(unsigned short ii = 0; ii < partmp.size(); ++ii) {
499  mf::LogVerbatim("CCHitFinder")<<ii<<" "<<partmp[ii]<<" "
500  <<partmperr[ii];
501  }
502  }
503 */
504  // ensure that the fit is reasonable
505  bool fitok = true;
506  for(unsigned short ii = 0; ii < nGaus; ++ii) {
507  unsigned short index = ii * 3;
508  // ensure that the fitted time is within the signal bounds
509  short fittime = partmp[index + 1];
510  if(fittime < 0 || fittime > npt - 1) {
511  fitok = false;
512  break;
513  }
514  // ensure that the signal peak is large enough
515  if(partmp[index] < fMinPeak[thePlane]) {
516  fitok = false;
517  break;
518  }
519  // ensure that the RMS is large enough but not too large
520  float rms = partmp[index + 2];
521  if(rms < 0.5 * fMinRMS[thePlane] || rms > 5 * fMinRMS[thePlane]) {
522  fitok = false;
523  break;
524  }
525  // ensure that the hits are not too similar in time (< 2 ticks)
526  for(unsigned short jj = 0; jj < nGaus; ++jj) {
527  if(jj == ii) continue;
528  unsigned short jndex = jj * 3;
529  float timediff = std::abs(partmp[jndex + 1] - partmp[index + 1]);
530  if(timediff < 2.) {
531  fitok = false;
532  break;
533  }
534  }
535  if(!fitok) break;
536  }
537 
538  if(fitok) {
539  par = partmp;
540  parerr = partmperr;
541  } else {
542  chidof = 9999.;
543  dof = -1;
544 // if(prt) mf::LogVerbatim("CCHitFinder")<<"Bad fit parameters";
545  }
546 
547  return;
548  } // FitNG
void AddMultiGaus(unsigned int nGaus)
std::vector< unsigned short > bumps
bool fUseFastFit
whether to attempt using a fast fit on single gauss.
FitStats_t TriedFitStats
counts of the tried fits
unsigned short thePlane
static bool FastGaussianFit(unsigned short npt, float const *ticks, float const *signl, std::array< double, 3 > &params, std::array< double, 3 > &paramerrors, float &chidof)
Performs a "fast" fit.
std::vector< float > fMinPeak
FitStats_t FinalFitStats
counts of the good fits
std::vector< double > parerr
std::vector< float > fMinRMS
std::vector< double > par
void hit::CCHitFinderAlg::MakeCrudeHit ( unsigned short  npt,
float *  ticks,
float *  signl 
)
private

Definition at line 551 of file CCHitFinderAlg.cxx.

References chidof, dof, pmtana::mean(), par, parerr, and Sqrt2Pi.

Referenced by RunCCHitFinder().

553  {
554  // make a single crude hit if fitting failed
555  float sumS = 0.;
556  float sumST = 0.;
557  for(unsigned short ii = 0; ii < npt; ++ii) {
558  sumS += signl[ii];
559  sumST += signl[ii] * ticks[ii];
560  }
561  float mean = sumST / sumS;
562  float rms = 0.;
563  for(unsigned short ii = 0; ii < npt; ++ii) {
564  float arg = ticks[ii] - mean;
565  rms += signl[ii] * arg * arg;
566  }
567  rms = std::sqrt(rms / sumS);
568  float amp = sumS / (Sqrt2Pi * rms);
569  par.clear();
570 /*
571  if(prt) mf::LogVerbatim("CCHitFinder")<<"Crude hit Amp "<<(int)amp<<" mean "
572  <<(int)mean<<" rms "<<rms;
573 */
574  par.push_back(amp);
575  par.push_back(mean);
576  par.push_back(rms);
577  // need to do the errors better
578  parerr.clear();
579  float amperr = npt;
580  float meanerr = std::sqrt(1/sumS);
581  float rmserr = 0.2 * rms;
582  parerr.push_back(amperr);
583  parerr.push_back(meanerr);
584  parerr.push_back(rmserr);
585 /*
586  if(prt) mf::LogVerbatim("CCHitFinder")<<" errors Amp "<<amperr<<" mean "
587  <<meanerr<<" rms "<<rmserr;
588 */
589  chidof = 9999.;
590  dof = -1;
591  } // MakeCrudeHit
double mean(const std::vector< short > &wf, size_t start, size_t nsample)
Definition: UtilFunc.cxx:15
static constexpr float Sqrt2Pi
std::vector< double > parerr
std::vector< double > par
template<typename Stream >
void hit::CCHitFinderAlg::PrintStats ( Stream &  out) const

Print the fit statistics.

Definition at line 235 of file CCHitFinderAlg.h.

References hit::CCHitFinderAlg::FitStats_t::FastFits, FinalFitStats, fUseFastFit, MaxGaussians, hit::CCHitFinderAlg::FitStats_t::MultiGausFits, and TriedFitStats.

Referenced by YieldHits().

235  {
236 
237  out << "CCHitFinderAlg fit statistics:";
238  if (fUseFastFit) {
239  out << "\n fast 1-Gaussian fits: " << FinalFitStats.FastFits
240  << " succeeded (" << TriedFitStats.FastFits << " tried)";
241  }
242  else
243  out << "\n fast 1-Gaussian fits: disabled";
244 
245  for (unsigned int nGaus = 1; nGaus < MaxGaussians; ++nGaus) {
246  if (TriedFitStats.MultiGausFits[nGaus-1] == 0) continue;
247  out << "\n " << nGaus << "-Gaussian fits: "
248  << FinalFitStats.MultiGausFits[nGaus-1]
249  << " accepted (" << TriedFitStats.MultiGausFits[nGaus-1] << " tried)";
250  } // for nGaus
251  if (TriedFitStats.MultiGausFits.back() > 0) {
252  out << "\n " << FinalFitStats.MultiGausFits.size()
253  << "-Gaussian fits or higher: " << FinalFitStats.MultiGausFits.back()
254  << " accepted (" << TriedFitStats.MultiGausFits.back() << " tried)";
255  }
256  out << std::endl;
257 
258 } // CCHitFinderAlg::FitStats_t::Print()
bool fUseFastFit
whether to attempt using a fast fit on single gauss.
FitStats_t TriedFitStats
counts of the tried fits
unsigned int FastFits
count of single-Gaussian fast fits
std::vector< unsigned int > MultiGausFits
multi-Gaussian stats
static constexpr unsigned int MaxGaussians
FitStats_t FinalFitStats
counts of the good fits
void hit::CCHitFinderAlg::reconfigure ( fhicl::ParameterSet const &  pset)
virtual

Definition at line 63 of file CCHitFinderAlg.cxx.

References art::errors::Configuration, fChiNorms, fChiSplit, FinalFitStats, fMaxBumps, fMaxXtraHits, fMinPeak, fMinRMS, fStudyHits, fUseChannelFilter, fUseFastFit, fUTickRange, fUWireRange, fVTickRange, fVWireRange, fWTickRange, fWWireRange, fhicl::ParameterSet::get(), fhicl::ParameterSet::has_key(), hit::CCHitFinderAlg::HitChannelInfo_t::HitChannelInfo_t(), LOG_WARNING, MaxGaussians, hit::CCHitFinderAlg::FitStats_t::Reset(), and TriedFitStats.

Referenced by CCHitFinderAlg(), and cluster::ClusterCrawler::reconfigure().

64  {
65  if(pset.has_key("MinSigInd")) throw art::Exception(art::errors::Configuration)
66  << "CCHitFinderAlg: Using no-longer-valid fcl input: MinSigInd, MinSigCol, etc";
67 
68  fMinPeak = pset.get<std::vector<float>>("MinPeak");
69  fMinRMS = pset.get<std::vector<float>>("MinRMS");
70  fMaxBumps = pset.get<unsigned short>("MaxBumps");
71  fMaxXtraHits = pset.get<unsigned short>("MaxXtraHits");
72  fChiSplit = pset.get<float>("ChiSplit");
73  fChiNorms = pset.get<std::vector< float > >("ChiNorms");
74  fUseFastFit = pset.get<bool>("UseFastFit", false);
75  fUseChannelFilter = pset.get<bool>("UseChannelFilter", true);
76  fStudyHits = pset.get<bool>("StudyHits", false);
77  // The following variables are only used in StudyHits mode
78  fUWireRange = pset.get< std::vector< short >>("UWireRange");
79  fUTickRange = pset.get< std::vector< short >>("UTickRange");
80  fVWireRange = pset.get< std::vector< short >>("VWireRange");
81  fVTickRange = pset.get< std::vector< short >>("VTickRange");
82  fWWireRange = pset.get< std::vector< short >>("WWireRange");
83  fWTickRange = pset.get< std::vector< short >>("WTickRange");
84 
85  if(fMinPeak.size() != fMinRMS.size()) {
86  mf::LogError("CCTF")<<"MinPeak size != MinRMS size";
87  return;
88  }
89 
90  if (fMaxBumps > MaxGaussians) {
91  // LOG_WARNING will point the user to this line of code.
92  // Any value of MaxGaussians can be used.
93  // That value is defined in the header file.
94  LOG_WARNING("CCHitFinderAlg")
95  << "CCHitFinder algorithm is currently hard-coded to support at most "
96  << MaxGaussians << " bumps per region of interest, but " << fMaxBumps
97  << " have been requested.\n"
98  << "We are forcing the parameter to " << MaxGaussians
99  << ". If this is not acceptable, increase CCHitFinderAlg::MaxGaussians"
100  << " value and recompile.";
101  fMaxBumps = MaxGaussians;
102  } // if too many gaussians
103 
106 
107  // sanity check for StudyHits mode
108  if(fStudyHits) {
109  if(fUWireRange.size() != 2 || fUTickRange.size() != 2 ||
110  fVWireRange.size() != 2 || fVTickRange.size() != 2 ||
111  fWWireRange.size() != 2 || fWTickRange.size() != 2) {
112  mf::LogError("CCHF")<<"Invalid vector size for StudyHits. Must be 2";
113  return;
114  }
115  } // fStudyHits
116 
117  }
unsigned short fMaxBumps
float fChiSplit
Estimated noise error on the Signal.
std::vector< short > fWTickRange
std::vector< short > fWWireRange
std::vector< float > fChiNorms
std::vector< short > fUWireRange
bool fUseFastFit
whether to attempt using a fast fit on single gauss.
std::vector< short > fVWireRange
std::vector< short > fVTickRange
FitStats_t TriedFitStats
counts of the tried fits
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
void Reset(unsigned int nGaus)
#define LOG_WARNING(category)
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
std::vector< float > fMinPeak
static constexpr unsigned int MaxGaussians
unsigned short fMaxXtraHits
std::vector< short > fUTickRange
FitStats_t FinalFitStats
counts of the good fits
std::vector< float > fMinRMS
void hit::CCHitFinderAlg::RunCCHitFinder ( std::vector< recob::Wire > const &  Wires)

Definition at line 128 of file CCHitFinderAlg.cxx.

References hit::CCHitFinderAlg::FitStats_t::AddMultiGaus(), allhits, bumps, recob::Wire::Channel(), geo::GeometryCore::ChannelToWire(), chidof, chinorm, dof, fChiNorms, fChiSplit, FinalFitStats, FitNG(), fMaxBumps, fMaxXtraHits, fMinPeak, fMinRMS, fStudyHits, MakeCrudeHit(), SelRAT, recob::Wire::Signal(), StoreHits(), StudyHits(), theChannel, thePlane, and theWireNum.

Referenced by cluster::ClusterCrawler::produce().

128  {
129 
130  allhits.clear();
131 
132  unsigned short maxticks = 1000;
133  float *ticks = new float[maxticks];
134  // define the ticks array used for fitting
135  for(unsigned short ii = 0; ii < maxticks; ++ii) {
136  ticks[ii] = ii;
137  }
138  float *signl = new float[maxticks];
139  float adcsum = 0;
140  // initialize the vectors for the hit study
141  if(fStudyHits) StudyHits(0);
142  bool first;
143 
144 // prt = false;
145  lariov::ChannelStatusProvider const& channelStatus
147 
148  for(size_t wireIter = 0; wireIter < Wires.size(); wireIter++){
149 
150  recob::Wire const& theWire = Wires[wireIter];
151  theChannel = theWire.Channel();
152  // ignore bad channels
153  if(channelStatus.IsBad(theChannel)) continue;
154 /*
155  geo::SigType_t SigType = geom->SignalType(theChannel);
156  minSig = 0.;
157  minRMS = 0.;
158  if(SigType == geo::kInduction){
159  minSig = fMinSigInd;
160  minRMS = fMinRMSInd;
161  }//<-- End if Induction Plane
162  else if(SigType == geo::kCollection){
163  minSig = fMinSigCol;
164  minRMS = fMinRMSCol;
165  }//<-- End if Collection Plane
166 */
167 
168  std::vector<geo::WireID> wids = geom->ChannelToWire(theChannel);
169  thePlane = wids[0].Plane;
170  if(thePlane > fMinPeak.size() - 1) {
171  mf::LogError("CCHF")<<"MinPeak vector too small for plane "<<thePlane;
172  return;
173  }
174  theWireNum = wids[0].Wire;
175  HitChannelInfo_t WireInfo(&theWire, wids[0], *geom);
176 
177  // minimum number of time samples
178  unsigned short minSamples = 2 * fMinRMS[thePlane];
179 
180  // factor used to normalize the chi/dof fits for each plane
182 
183  // edit this line to debug hit fitting on a particular plane/wire
184 // prt = (thePlane == 1 && theWireNum == 839);
185  std::vector<float> signal(theWire.Signal());
186 
187  unsigned short nabove = 0;
188  unsigned short tstart = 0;
189  unsigned short maxtime = signal.size() - 2;
190  // find the min time when the signal is below threshold
191  unsigned short mintime = 3;
192  for(unsigned short time = 3; time < maxtime; ++time) {
193  if(signal[time] < fMinPeak[thePlane]) {
194  mintime = time;
195  break;
196  }
197  }
198  for(unsigned short time = mintime; time < maxtime; ++time) {
199  if(signal[time] > fMinPeak[thePlane]) {
200  if(nabove == 0) tstart = time;
201  ++nabove;
202  } else {
203  // check for a wide enough signal above threshold
204  if(nabove > minSamples) {
205  // skip this wire if the RAT is too long
206  if(nabove > maxticks) mf::LogError("CCHitFinder")
207  <<"Long RAT "<<nabove<<" "<<maxticks
208  <<" No signal on wire "<<theWireNum<<" after time "<<time;
209  if(nabove > maxticks) break;
210  unsigned short npt = 0;
211  // look for bumps to inform the fit
212  bumps.clear();
213  adcsum = 0;
214  for(unsigned short ii = tstart; ii < time; ++ii) {
215  signl[npt] = signal[ii];
216  adcsum += signl[npt];
217  if(signal[ii ] > signal[ii - 1] &&
218  signal[ii - 1] > signal[ii - 2] &&
219  signal[ii ] > signal[ii + 1] &&
220  signal[ii + 1] > signal[ii + 2]) bumps.push_back(npt);
221 // if(prt) mf::LogVerbatim("CCHitFinder")<<"signl "<<ii<<" "<<signl[npt];
222  ++npt;
223  }
224  // decide if this RAT should be studied
225  if(fStudyHits) StudyHits(1, npt, ticks, signl, tstart);
226  // just make a crude hit if too many bumps
227  if(bumps.size() > fMaxBumps) {
228  MakeCrudeHit(npt, ticks, signl);
229  StoreHits(tstart, npt, WireInfo, adcsum);
230  nabove = 0;
231  continue;
232  }
233  // start looking for hits with the found bumps
234  unsigned short nHitsFit = bumps.size();
235  unsigned short nfit = 0;
236  chidof = 0.;
237  dof = -1;
238  bool HitStored = false;
239  unsigned short nMaxFit = bumps.size() + fMaxXtraHits;
240  // only used in StudyHits mode
241  first = true;
242  while(nHitsFit <= nMaxFit) {
243 
244  FitNG(nHitsFit, npt, ticks, signl);
245  if(fStudyHits && first && SelRAT) {
246  first = false;
247  StudyHits(2, npt, ticks, signl, tstart);
248  }
249  // good chisq so store it
250  if(chidof < fChiSplit) {
251  StoreHits(tstart, npt, WireInfo, adcsum);
252  HitStored = true;
253  break;
254  }
255  // the previous fit was better, so revert to it and
256  // store it
257  ++nHitsFit;
258  ++nfit;
259  } // nHitsFit < fMaxXtraHits
260  if( !HitStored && npt < maxticks) {
261  // failed all fitting. Make a crude hit
262  MakeCrudeHit(npt, ticks, signl);
263  StoreHits(tstart, npt, WireInfo, adcsum);
264  }
265  else if (nHitsFit > 0) FinalFitStats.AddMultiGaus(nHitsFit);
266  } // nabove > minSamples
267  nabove = 0;
268  } // signal < fMinPeak
269  } // time
270  } // wireIter
271 
272  // print out
273  if(fStudyHits) StudyHits(4);
274 
275  delete[] ticks;
276  delete[] signl;
277 
278  } //RunCCHitFinder
unsigned short fMaxBumps
void AddMultiGaus(unsigned int nGaus)
float fChiSplit
Estimated noise error on the Signal.
void FitNG(unsigned short nGaus, unsigned short npt, float *ticks, float *signl)
raw::ChannelID_t theChannel
void StoreHits(unsigned short TStart, unsigned short npt, HitChannelInfo_t info, float adcsum)
void StudyHits(unsigned short flag, unsigned short npt=0, float *ticks=0, float *signl=0, unsigned short tstart=0)
std::vector< unsigned short > bumps
std::vector< float > fChiNorms
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
std::vector< recob::Hit > allhits
art::ServiceHandle< geo::Geometry > geom
unsigned short thePlane
raw::ChannelID_t Channel() const
Returns the ID of the channel (or InvalidChannelID)
Definition: Wire.h:164
std::vector< float > Signal() const
Return a zero-padded full length vector filled with RoI signal.
Definition: Wire.cxx:47
unsigned short theWireNum
std::vector< float > fMinPeak
Class holding the deconvoluted signals from a channel.
Definition: Wire.h:80
unsigned short fMaxXtraHits
void MakeCrudeHit(unsigned short npt, float *ticks, float *signl)
FitStats_t FinalFitStats
counts of the good fits
std::vector< float > fMinRMS
void hit::CCHitFinderAlg::StoreHits ( unsigned short  TStart,
unsigned short  npt,
HitChannelInfo_t  info,
float  adcsum 
)
private

Definition at line 595 of file CCHitFinderAlg.cxx.

References allhits, recob::Wire::Channel(), chidof, dof, fStudyHits, hiTime, loTime, par, parerr, hit::CCHitFinderAlg::HitChannelInfo_t::sigType, Sqrt2Pi, SqrtPi, StudyHits(), recob::Wire::View(), hit::CCHitFinderAlg::HitChannelInfo_t::wire, and hit::CCHitFinderAlg::HitChannelInfo_t::wireID.

Referenced by RunCCHitFinder().

597  {
598  // store the hits in the struct
599  size_t nhits = par.size() / 3;
600 
601  if(allhits.max_size() - allhits.size() < nhits) {
602  mf::LogError("CCHitFinder")
603  << "Too many hits: existing " << allhits.size() << " plus new " << nhits
604  << " beyond the maximum " << allhits.max_size();
605  return;
606  }
607 
608  if(nhits == 0) return;
609 
610  // fill RMS for single hits
611  if(fStudyHits) StudyHits(3);
612 
613  const float loTime = TStart;
614  const float hiTime = TStart + npt;
615 
616  // Find sum of the areas of all Gaussians
617  float gsum = 0.;
618  for(size_t hit = 0; hit < nhits; ++hit) {
619  const unsigned short index = 3 * hit;
620  gsum += Sqrt2Pi * par[index] * par[index + 2];
621  }
622  for(size_t hit = 0; hit < nhits; ++hit) {
623  const size_t index = 3 * hit;
624  const float charge = Sqrt2Pi * par[index] * par[index + 2];
625  const float charge_err = SqrtPi
626  * (parerr[index] * par[index + 2] + par[index] * parerr[index + 2]);
627 
628  allhits.emplace_back(
629  info.wire->Channel(), // channel
630  loTime, // start_tick
631  hiTime, // end_tick
632  par[index + 1] + TStart, // peak_time
633  parerr[index + 1], // sigma_peak_time
634  par[index + 2], // rms
635  par[index], // peak_amplitude
636  parerr[index], // sigma_peak_amplitude
637  adcsum * charge / gsum, // summedADC
638  charge, // hit_integral
639  charge_err, // hit_sigma_integral
640  nhits, // multiplicity
641  hit, // local_index
642  chidof, // goodness_of_fit
643  dof, // dof
644  info.wire->View(), // view
645  info.sigType, // signal_type
646  info.wireID // wireID
647  );
648 /*
649  if(prt) {
650  mf::LogVerbatim("CCHitFinder")<<"W:T "<<allhits.back().WireID().Wire
651  <<":"<<(short)allhits.back().PeakTime()
652  <<" Chg "<<(short)allhits.back().Integral()
653  <<" RMS "<<allhits.back().RMS()
654  <<" lo ID "<<allhits.back().LocalIndex()
655  <<" numHits "<<allhits.back().Multiplicity()
656  <<" loTime "<<allhits.back().StartTick()<<" hiTime "<<allhits.back().EndTick()
657  <<" chidof "<<allhits.back().GoodnessOfFit()
658  << " DOF " << allhits.back().DegreesOfFreedom();
659  }
660 */
661  } // hit
662  } // StoreHits
static constexpr float SqrtPi
void StudyHits(unsigned short flag, unsigned short npt=0, float *ticks=0, float *signl=0, unsigned short tstart=0)
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
std::vector< recob::Hit > allhits
std::vector< float > hiTime
Detector simulation of raw signals on wires.
std::vector< float > loTime
static constexpr float Sqrt2Pi
std::vector< double > parerr
std::vector< double > par
void hit::CCHitFinderAlg::StudyHits ( unsigned short  flag,
unsigned short  npt = 0,
float *  ticks = 0,
float *  signl = 0,
unsigned short  tstart = 0 
)
private

Definition at line 666 of file CCHitFinderAlg.cxx.

References bumpChi, bumpCnt, bumpRMS, bumps, chidof, fChiNorms, fMinPeak, fUTickRange, fUWireRange, fVTickRange, fVWireRange, fWTickRange, fWWireRange, hitCnt, hiTime, hitRMS, hiWire, loTime, loWire, par, RATCnt, art::right(), SelRAT, thePlane, and theWireNum.

Referenced by RunCCHitFinder(), and StoreHits().

667  {
668  // study hits in user-selected ranges of wires and ticks in each plane. The user should identify
669  // a shallow-angle isolated track, e.g. using the event display, to determine the wire/tick ranges.
670  // One hit should be reconstructed on each wire when the hit finding fcl parameters are set correctly.
671  // The intent of this study is to determine the correct fcl parameters. The flag variable determines
672  // the operation performed.
673  // flag = 0: Initialize study vectors
674  // flag = 1: Set SelRat true if the Region Above Threshold resides within a wire/hit range
675  // flag = 2: Find the maximum signal and calculate the RMS. Also find the low and high ticks of signals
676  // in the wire range to allow a later calculation of the track angle. This isn't strictly
677  // necessary for the study and presumes that the user has selected compatible regions in each plane.
678  // flag = 3: Accumulate the RMS from the first Gaussian fit
679  // flag = 4: Calculate recommended fcl parameters and print the results to the screen
680 
681  // init
682  if(flag == 0) {
683  for(unsigned short ipl = 0; ipl < 3; ++ipl) {
684  // Average chisq of the first fit on a single bump in each plane
685  bumpChi.push_back(0.);
686  // Average RMS of the dump
687  bumpRMS.push_back(0.);
688  // The number of single bumps in each plane
689  bumpCnt.push_back(0.);
690  // number of RATs
691  RATCnt.push_back(0);
692  // The number of single hits found in each plane
693  hitCnt.push_back(0.);
694  // Average reconstructed hit RMS
695  hitRMS.push_back(0.);
696  // lo/hi wire/time
697  loWire.push_back(9999.);
698  loTime.push_back(0.);
699  hiWire.push_back(-1.);
700  hiTime.push_back(0.);
701  } // ii
702  return;
703  } // flag == 0
704 
705  if(flag == 1) {
706  SelRAT = false;
707  if(thePlane == 0) {
708  if(theWireNum > fUWireRange[0] && theWireNum < fUWireRange[1] &&
709  tstart > fUTickRange[0] && tstart < fUTickRange[1]) {
710  SelRAT = true;
711  RATCnt[thePlane] += 1;
712  }
713  return;
714  } // thePlane == 0
715  if(thePlane == 1) {
716  if(theWireNum > fVWireRange[0] && theWireNum < fVWireRange[1] &&
717  tstart > fVTickRange[0] && tstart < fVTickRange[1]) {
718  SelRAT = true;
719  RATCnt[thePlane] += 1;
720  }
721  return;
722  } // thePlane == 1
723  if(thePlane == 2) {
724  if(theWireNum > fWWireRange[0] && theWireNum < fWWireRange[1] &&
725  tstart > fWTickRange[0] && tstart < fWTickRange[1]) {
726  SelRAT = true;
727  RATCnt[thePlane] += 1;
728  }
729  return;
730  } // thePlane == 2
731  } // flag == 1
732 
733  if(flag == 2) {
734  if(!SelRAT) return;
735  // in this section we find the low/hi wire/time for a signal. This can be used to calculate
736  // the slope dT/dW to study hit width, fraction of crude hits, etc vs dT/dW
737  float big = 0.;
738  float imbig = 0.;
739  for(unsigned short ii = 0; ii < npt; ++ii) {
740  if(signl[ii] > big) {
741  big = signl[ii];
742  imbig = ii;
743  }
744  } // ii
745  // require a significant PH
746  if(big > fMinPeak[0]) {
747  // get the Lo info
748  if(theWireNum < loWire[thePlane]) {
750  loTime[thePlane] = tstart + imbig;
751  }
752  // get the Hi info
753  if(theWireNum > hiWire[thePlane]) {
755  hiTime[thePlane] = tstart + imbig;
756  }
757  } // big > fMinPeak[0]
758  if(bumps.size() == 1 && chidof < 9999.) {
759  bumpCnt[thePlane] += bumps.size();
760  bumpChi[thePlane] += chidof;
761  // calculate the average bin
762  float sumt = 0.;
763  float sum = 0.;
764  for(unsigned short ii = 0; ii < npt; ++ii) {
765  sum += signl[ii];
766  sumt += signl[ii] * ii;
767  } // ii
768  float aveb = sumt / sum;
769  // now calculate the RMS
770  sumt = 0.;
771  for(unsigned short ii = 0; ii < npt; ++ii) {
772  float dbin = (float)ii - aveb;
773  sumt += signl[ii] * dbin * dbin;
774  } // ii
775  bumpRMS[thePlane] += std::sqrt(sumt / sum);
776  } // bumps.size() == 1 && chidof < 9999.
777  return;
778  } // flag == 2
779 
780  // fill info for single hits
781  if(flag == 3) {
782  if(!SelRAT) return;
783  if(par.size() == 3) {
784  hitCnt[thePlane] += 1;
785  hitRMS[thePlane] += par[2];
786  }
787  return;
788  }
789 
790 
791  if(flag == 4) {
792  // The goal is to adjust the fcl inputs so that the number of single
793  // hits found is ~equal to the number of single bumps found for shallow
794  // angle tracks. The ChiNorm inputs should be adjusted so the average
795  // chisq/DOF is ~1 in each plane.
796  std::cout<<"Check lo and hi W/T for each plane"<<std::endl;
797  for(unsigned short ipl = 0; ipl < 3; ++ipl) {
798  std::cout<<ipl<<" lo "<<loWire[ipl]<<" "<<loTime[ipl]
799  <<" hi "<<hiWire[ipl]<<" "<<hiTime[ipl]<<std::endl;
800  }
801  std::cout<<" ipl nRAT bCnt bChi bRMS hCnt hRMS dT/dW New_ChiNorm"<<std::endl;
802  for(unsigned short ipl = 0; ipl < 3; ++ipl) {
803  if(bumpCnt[ipl] > 0) {
804  bumpChi[ipl] = bumpChi[ipl] / (float)bumpCnt[ipl];
805  bumpRMS[ipl] = bumpRMS[ipl] / (float)bumpCnt[ipl];
806  hitRMS[ipl] = hitRMS[ipl] / (float)hitCnt[ipl];
807  // calculate the slope
808  float dTdW = std::abs((hiTime[ipl] - loTime[ipl]) / (hiWire[ipl] - loWire[ipl]));
809  std::cout<<ipl<<std::right<<std::setw(5)<<RATCnt[ipl]
810  <<std::setw(5)<<bumpCnt[ipl]
811  <<std::setw(7)<<std::fixed<<std::setprecision(2)<<bumpChi[ipl]
812  <<std::setw(7)<<bumpRMS[ipl]
813  <<std::setw(7)<<hitCnt[ipl]
814  <<std::setw(7)<<std::setprecision(1)<<hitRMS[ipl]
815  <<std::setw(7)<<dTdW
816  <<std::setw(7)<<std::setprecision(2)
817  <<bumpChi[ipl]*fChiNorms[ipl]
818  <<std::endl;
819  } //
820  } // ipl
821  std::cout<<"nRAT is the number of Regions Above Threshold (RAT) used in the study.\n";
822  std::cout<<"bCnt is the number of single bumps that were successfully fitted \n";
823  std::cout<<"bChi is the average chisq/DOF of the first fit\n";
824  std::cout<<"bRMS is the average calculated RMS of the bumps\n";
825  std::cout<<"hCnt is the number of RATs that have a single hit\n";
826  std::cout<<"hRMS is the average RMS from the Gaussian fit -> use this value for fMinRMS[plane] in the fcl file\n";
827  std::cout<<"dTdW is the slope of the track\n";
828  std::cout<<"New_ChiNorm is the recommended values of ChiNorm that should be used in the fcl file\n";
829  bumpChi.clear();
830  bumpRMS.clear();
831  bumpCnt.clear();
832  RATCnt.clear();
833  hitRMS.clear();
834  hitCnt.clear();
835  loWire.clear();
836  loTime.clear();
837  hiWire.clear();
838  hiTime.clear();
839  }
840  } // StudyHits
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:112
std::vector< short > fWTickRange
std::vector< short > fWWireRange
std::vector< float > loWire
std::vector< float > hiWire
std::vector< unsigned short > bumps
std::vector< float > fChiNorms
std::vector< float > bumpRMS
std::vector< short > fUWireRange
std::vector< short > fVWireRange
std::vector< short > fVTickRange
unsigned short thePlane
std::vector< float > hiTime
std::vector< int > bumpCnt
std::vector< float > loTime
unsigned short theWireNum
std::vector< float > fMinPeak
std::vector< int > hitCnt
std::vector< short > fUTickRange
std::vector< float > hitRMS
std::vector< int > RATCnt
std::vector< float > bumpChi
std::vector< double > par
std::vector<recob::Hit>&& hit::CCHitFinderAlg::YieldHits ( )
inline

Returns (and loses) the collection of reconstructed hits.

Definition at line 96 of file CCHitFinderAlg.h.

References PrintStats().

Referenced by cluster::ClusterCrawler::produce().

96 { return std::move(allhits); }
std::vector< recob::Hit > allhits

Member Data Documentation

std::vector<recob::Hit> hit::CCHitFinderAlg::allhits

Definition at line 86 of file CCHitFinderAlg.h.

Referenced by RunCCHitFinder(), and StoreHits().

std::vector<float> hit::CCHitFinderAlg::bumpChi
private

Definition at line 169 of file CCHitFinderAlg.h.

Referenced by StudyHits().

std::vector<int> hit::CCHitFinderAlg::bumpCnt
private

Definition at line 167 of file CCHitFinderAlg.h.

Referenced by StudyHits().

std::vector<float> hit::CCHitFinderAlg::bumpRMS
private

Definition at line 170 of file CCHitFinderAlg.h.

Referenced by StudyHits().

std::vector<unsigned short> hit::CCHitFinderAlg::bumps
private

Definition at line 140 of file CCHitFinderAlg.h.

Referenced by FitNG(), RunCCHitFinder(), and StudyHits().

float hit::CCHitFinderAlg::chidof
private

Definition at line 138 of file CCHitFinderAlg.h.

Referenced by FitNG(), MakeCrudeHit(), RunCCHitFinder(), StoreHits(), and StudyHits().

float hit::CCHitFinderAlg::chinorm
private

Definition at line 119 of file CCHitFinderAlg.h.

Referenced by FitNG(), and RunCCHitFinder().

int hit::CCHitFinderAlg::dof
private

Definition at line 139 of file CCHitFinderAlg.h.

Referenced by FitNG(), MakeCrudeHit(), RunCCHitFinder(), and StoreHits().

std::vector<float> hit::CCHitFinderAlg::fChgNorms
private

Definition at line 113 of file CCHitFinderAlg.h.

std::vector<float> hit::CCHitFinderAlg::fChiNorms
private

Definition at line 111 of file CCHitFinderAlg.h.

Referenced by reconfigure(), RunCCHitFinder(), and StudyHits().

float hit::CCHitFinderAlg::fChiSplit
private

Estimated noise error on the Signal.

Definition at line 108 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and RunCCHitFinder().

FitStats_t hit::CCHitFinderAlg::FinalFitStats
private

counts of the good fits

Definition at line 197 of file CCHitFinderAlg.h.

Referenced by FitNG(), PrintStats(), reconfigure(), and RunCCHitFinder().

std::unique_ptr<GausFitCache> hit::CCHitFinderAlg::FitCache
private

a set of functions ready to be used

Definition at line 182 of file CCHitFinderAlg.h.

unsigned short hit::CCHitFinderAlg::fMaxBumps
private

Definition at line 106 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and RunCCHitFinder().

unsigned short hit::CCHitFinderAlg::fMaxXtraHits
private

Definition at line 107 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and RunCCHitFinder().

std::vector<float> hit::CCHitFinderAlg::fMinPeak
private

Definition at line 104 of file CCHitFinderAlg.h.

Referenced by FitNG(), reconfigure(), RunCCHitFinder(), and StudyHits().

std::vector<float> hit::CCHitFinderAlg::fMinRMS
private

Definition at line 105 of file CCHitFinderAlg.h.

Referenced by FitNG(), reconfigure(), and RunCCHitFinder().

bool hit::CCHitFinderAlg::fStudyHits
private

Definition at line 161 of file CCHitFinderAlg.h.

Referenced by reconfigure(), RunCCHitFinder(), and StoreHits().

std::vector<float> hit::CCHitFinderAlg::fTimeOffsets
private

Definition at line 112 of file CCHitFinderAlg.h.

bool hit::CCHitFinderAlg::fUseChannelFilter
private

Definition at line 124 of file CCHitFinderAlg.h.

Referenced by reconfigure().

bool hit::CCHitFinderAlg::fUseFastFit
private

whether to attempt using a fast fit on single gauss.

Definition at line 180 of file CCHitFinderAlg.h.

Referenced by FitNG(), PrintStats(), and reconfigure().

std::vector< short > hit::CCHitFinderAlg::fUTickRange
private

Definition at line 162 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and StudyHits().

std::vector< short > hit::CCHitFinderAlg::fUWireRange
private

Definition at line 162 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and StudyHits().

std::vector< short > hit::CCHitFinderAlg::fVTickRange
private

Definition at line 163 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and StudyHits().

std::vector< short > hit::CCHitFinderAlg::fVWireRange
private

Definition at line 163 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and StudyHits().

std::vector< short > hit::CCHitFinderAlg::fWTickRange
private

Definition at line 164 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and StudyHits().

std::vector< short > hit::CCHitFinderAlg::fWWireRange
private

Definition at line 164 of file CCHitFinderAlg.h.

Referenced by reconfigure(), and StudyHits().

art::ServiceHandle<geo::Geometry> hit::CCHitFinderAlg::geom
private

Definition at line 128 of file CCHitFinderAlg.h.

std::vector<int> hit::CCHitFinderAlg::hitCnt
private

Definition at line 171 of file CCHitFinderAlg.h.

Referenced by StudyHits().

std::vector<float> hit::CCHitFinderAlg::hiTime
private

Definition at line 177 of file CCHitFinderAlg.h.

Referenced by StoreHits(), and StudyHits().

std::vector<float> hit::CCHitFinderAlg::hitRMS
private

Definition at line 172 of file CCHitFinderAlg.h.

Referenced by StudyHits().

std::vector<float> hit::CCHitFinderAlg::hiWire
private

Definition at line 176 of file CCHitFinderAlg.h.

Referenced by StudyHits().

std::vector<float> hit::CCHitFinderAlg::loTime
private

Definition at line 175 of file CCHitFinderAlg.h.

Referenced by StoreHits(), and StudyHits().

std::vector<float> hit::CCHitFinderAlg::loWire
private

Definition at line 174 of file CCHitFinderAlg.h.

Referenced by StudyHits().

constexpr unsigned int hit::CCHitFinderAlg::MaxGaussians = 20
staticprivate

Definition at line 224 of file CCHitFinderAlg.h.

Referenced by PrintStats(), and reconfigure().

std::vector<double> hit::CCHitFinderAlg::par
private

Definition at line 134 of file CCHitFinderAlg.h.

Referenced by FastGaussianFit(), FitNG(), MakeCrudeHit(), StoreHits(), and StudyHits().

std::vector<double> hit::CCHitFinderAlg::parerr
private

Definition at line 135 of file CCHitFinderAlg.h.

Referenced by FitNG(), MakeCrudeHit(), and StoreHits().

std::vector<double> hit::CCHitFinderAlg::parmax
private

Definition at line 137 of file CCHitFinderAlg.h.

std::vector<double> hit::CCHitFinderAlg::parmin
private

Definition at line 136 of file CCHitFinderAlg.h.

std::vector<int> hit::CCHitFinderAlg::RATCnt
private

Definition at line 168 of file CCHitFinderAlg.h.

Referenced by StudyHits().

bool hit::CCHitFinderAlg::SelRAT
private

Definition at line 178 of file CCHitFinderAlg.h.

Referenced by RunCCHitFinder(), and StudyHits().

constexpr float hit::CCHitFinderAlg::Sqrt2Pi = 2.5066
staticprivate

Definition at line 121 of file CCHitFinderAlg.h.

Referenced by MakeCrudeHit(), and StoreHits().

constexpr float hit::CCHitFinderAlg::SqrtPi = 1.7725
staticprivate

Definition at line 122 of file CCHitFinderAlg.h.

Referenced by StoreHits().

raw::ChannelID_t hit::CCHitFinderAlg::theChannel
private

Definition at line 115 of file CCHitFinderAlg.h.

Referenced by RunCCHitFinder().

unsigned short hit::CCHitFinderAlg::thePlane
private

Definition at line 117 of file CCHitFinderAlg.h.

Referenced by FitNG(), RunCCHitFinder(), and StudyHits().

unsigned short hit::CCHitFinderAlg::theWireNum
private

Definition at line 116 of file CCHitFinderAlg.h.

Referenced by RunCCHitFinder(), and StudyHits().

FitStats_t hit::CCHitFinderAlg::TriedFitStats
private

counts of the tried fits

Definition at line 198 of file CCHitFinderAlg.h.

Referenced by FitNG(), PrintStats(), and reconfigure().


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