LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
CandHitMorphological_tool.cc
Go to the documentation of this file.
1 
8 
13 #include "cetlib_except/exception.h"
16 
17 #include "TH1F.h"
18 #include "TProfile.h"
19 
20 #include <cmath>
21 #include <fstream>
22 
23 namespace reco_tool
24 {
25 
27 {
28 public:
29  explicit CandHitMorphological(const fhicl::ParameterSet& pset);
30 
32 
33  void configure(const fhicl::ParameterSet& pset) override;
34 
35  void findHitCandidates(const Waveform&,
36  size_t,
37  size_t,
38  size_t,
39  HitCandidateVec&) const override;
40 
41  void MergeHitCandidates(const Waveform&,
42  const HitCandidateVec&,
43  MergeHitCandidateVec&) const override;
44 
45 private:
46  // Internal functions
47  //< Top level hit finding using erosion/dilation vectors
51  size_t,
52  float,
53  HitCandidateVec&) const;
54 
55  //< Fine grain hit finding within candidate peak regions using derivative method
58  size_t,
59  int,
60  float,
61  HitCandidateVec&) const;
62 
63  // Finding the nearest maximum/minimum from current point
66 
67  // handle finding the "start" and "stop" of a candidate hit
70 
71  // some fhicl control variables
72  size_t fPlane; //< Identifies the plane this tool is meant to operate on
73  float fDilationThreshold; //< Dilation threshold
74  float fDilationFraction; //< Fraction of max dilation to set for min dilation
75  float fErosionFraction; //< Fraction of max dilation value to set min erosion
76  int fMinDeltaTicks; //< minimum ticks from max to min to consider
77  float fMinDeltaPeaks; //< minimum maximum to minimum peak distance
78  float fMinHitHeight; //< Drop candidate hits with height less than this
79  size_t fNumInterveningTicks; //< Number ticks between candidate hits to merge
80  int fStructuringElement; //< Window size for morphologcial filter
81  bool fOutputHistograms; //< If true will generate a very large file of hists!
82 
84 
85  // Global histograms
86  TH1F* fDStopStartHist; //< Basically keeps track of the length of hit regions
87  TH1F* fDMaxTickMinTickHist; //< This will be a measure of the width of candidate hits
88  TH1F* fDMaxDerivMinDerivHist; //< This is the difference peak to peak of derivative for cand hit
89 
90  mutable size_t fLastChannel; //< Kludge to keep track of last channel when histogramming in effect
91  mutable size_t fChannelCnt; //< Counts the number of times a channel is used (assumed in order)
92 
93  //< All of the real work is done in the waveform tool
94  std::unique_ptr<reco_tool::IWaveformTool> fWaveformTool;
95 
96  const geo::GeometryCore* fGeometry = lar::providerFrom<geo::Geometry>();
97 };
98 
99 //----------------------------------------------------------------------
100 // Constructor.
102 {
103  configure(pset);
104 }
105 
107 {
108 }
109 
111 {
112  // Recover our parameters
113  fPlane = pset.get< size_t >("Plane", 0);
114  fDilationThreshold = pset.get< float >("DilationThreshold", 4.);
115  fDilationFraction = pset.get< float >("DilationFraction", 0.75);
116  fErosionFraction = pset.get< float >("ErosionFraction", 0.2);
117  fMinDeltaTicks = pset.get< int >("MinDeltaTicks", 0);
118  fMinDeltaPeaks = pset.get< float >("MinDeltaPeaks", 0.025);
119  fMinHitHeight = pset.get< float >("MinHitHeight", 1.0);
120  fNumInterveningTicks = pset.get< size_t >("NumInterveningTicks", 6);
121  fStructuringElement = pset.get< int >("StructuringElement", 20);
122  fOutputHistograms = pset.get< bool >("OutputHistograms", false);
123 
124  // Recover the baseline tool
125  fWaveformTool = art::make_tool<reco_tool::IWaveformTool> (pset.get<fhicl::ParameterSet>("WaveformAlgs"));
126 
127  // Set the last channel to some nonsensical value
129  fChannelCnt = 0;
130 
131  // If asked, define the global histograms
132  if (fOutputHistograms)
133  {
134  // Access ART's TFileService, which will handle creating and writing
135  // histograms and n-tuples for us.
137 
138  fHistDirectory = tfs.get();
139 
140  // Make a directory for these histograms
141  art::TFileDirectory dir = fHistDirectory->mkdir(Form("HitPlane_%1zu",fPlane));
142 
143  fDStopStartHist = dir.make<TH1F>(Form("DStopStart_%1zu", fPlane), ";Delta Stop/Start;", 200, 0., 200.);
144  fDMaxTickMinTickHist = dir.make<TH1F>(Form("DMaxTMinT_%1zu", fPlane), ";Delta Max/Min Tick;", 200, 0., 200.);
145  fDMaxDerivMinDerivHist = dir.make<TH1F>(Form("DMaxDMinD_%1zu", fPlane), ";Delta Max/Min Deriv;", 200, 0., 200.);
146  }
147 
148  return;
149 }
150 
152  size_t roiStartTick,
153  size_t channel,
154  size_t eventCount,
155  HitCandidateVec& hitCandidateVec) const
156 {
157  // In this case we want to find hit candidates based on the derivative of of the input waveform
158  // We get this from our waveform algs too...
159  Waveform rawDerivativeVec;
160  Waveform derivativeVec;
161 
162  fWaveformTool->firstDerivative(waveform, rawDerivativeVec);
163  fWaveformTool->triangleSmooth(rawDerivativeVec, derivativeVec);
164 
165  // Now we get the erosion/dilation vectors
166  Waveform erosionVec;
167  Waveform dilationVec;
168  Waveform averageVec;
169  Waveform differenceVec;
170 
171  reco_tool::HistogramMap histogramMap;
172 
173  // Compute the morphological filter vectors
174  fWaveformTool->getErosionDilationAverageDifference(waveform, fStructuringElement, histogramMap, erosionVec, dilationVec, averageVec, differenceVec);
175 
176  // Now find the hits
177  findHitCandidates(derivativeVec.begin(), derivativeVec.end(),
178  erosionVec.begin(), erosionVec.end(),
179  dilationVec.begin(), dilationVec.end(),
180  roiStartTick,
182  hitCandidateVec);
183 
184  // Reset the hit height from the input waveform
185  for(auto& hitCandidate : hitCandidateVec)
186  {
187  size_t centerIdx = hitCandidate.hitCenter;
188 
189  hitCandidate.hitHeight = waveform.at(centerIdx);
190  }
191 
192  // Keep track of histograms if requested
193  if (fOutputHistograms)
194  {
195  // Recover the details...
196  std::vector<geo::WireID> wids = fGeometry->ChannelToWire(channel);
197  size_t plane = wids[0].Plane;
198  size_t cryo = wids[0].Cryostat;
199  size_t tpc = wids[0].TPC;
200  size_t wire = wids[0].Wire;
201 
202  if (channel != fLastChannel) fChannelCnt = 0;
203 
204  // Make a directory for these histograms
205  art::TFileDirectory dir = fHistDirectory->mkdir(Form("HitPlane_%1zu/ev%04zu/c%1zut%1zuwire_%05zu",plane,eventCount,cryo,tpc,wire));
206 
207  size_t waveformSize = waveform.size();
208  int waveStart = roiStartTick;
209  int waveStop = waveStart + waveformSize;
210 
211  TProfile* waveHist = dir.make<TProfile>(Form("HWfm_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Waveform", waveformSize, waveStart, waveStop, -500., 500.);
212  TProfile* derivHist = dir.make<TProfile>(Form("HDer_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Derivative", waveformSize, waveStart, waveStop, -500., 500.);
213  TProfile* erosionHist = dir.make<TProfile>(Form("HEro_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Erosion", waveformSize, waveStart, waveStop, -500., 500.);
214  TProfile* dilationHist = dir.make<TProfile>(Form("HDil_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Dilation", waveformSize, waveStart, waveStop, -500., 500.);
215  TProfile* candHitHist = dir.make<TProfile>(Form("HCan_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Cand Hits", waveformSize, waveStart, waveStop, -500., 500.);
216  TProfile* maxDerivHist = dir.make<TProfile>(Form("HMax_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Maxima", waveformSize, waveStart, waveStop, -500., 500.);
217  TProfile* strtStopHist = dir.make<TProfile>(Form("HSSS_%03zu_ctw%01zu-%01zu-%01zu-%05zu",fChannelCnt,cryo,tpc,plane,wire), "Start/Stop", waveformSize, waveStart, waveStop, -500., 500.);
218 
219  // Fill wave/derivative
220  for(size_t idx = 0; idx < waveform.size(); idx++)
221  {
222  waveHist->Fill(roiStartTick + idx, waveform.at(idx));
223  derivHist->Fill(roiStartTick + idx, derivativeVec.at(idx));
224  erosionHist->Fill(roiStartTick + idx, erosionVec.at(idx));
225  dilationHist->Fill(roiStartTick + idx, dilationVec.at(idx));
226  }
227 
228  // Fill hits
229  for(const auto& hitCandidate : hitCandidateVec)
230  {
231  candHitHist->Fill(hitCandidate.hitCenter, hitCandidate.hitHeight);
232  maxDerivHist->Fill(hitCandidate.maxTick, hitCandidate.maxDerivative);
233  maxDerivHist->Fill(hitCandidate.minTick, hitCandidate.minDerivative);
234  strtStopHist->Fill(hitCandidate.startTick, waveform.at(hitCandidate.startTick));
235  strtStopHist->Fill(hitCandidate.stopTick, waveform.at(hitCandidate.stopTick));
236 
237  fDStopStartHist->Fill(hitCandidate.stopTick - hitCandidate.startTick, 1.);
238  fDMaxTickMinTickHist->Fill(hitCandidate.minTick - hitCandidate.maxTick, 1.);
239  fDMaxDerivMinDerivHist->Fill(hitCandidate.maxDerivative - hitCandidate.minDerivative, 1.);
240  }
241 
242  fLastChannel = channel;
243  fChannelCnt++;
244  }
245 
246  return;
247 }
248 
250  Waveform::const_iterator erosionStartItr, Waveform::const_iterator erosionStopItr,
251  Waveform::const_iterator dilationStartItr, Waveform::const_iterator dilationStopItr,
252  size_t roiStartTick,
253  float dilationThreshold,
254  HitCandidateVec& hitCandidateVec) const
255 {
256  // This function aims to use the erosion/dilation vectors to find candidate hit regions
257  // Once armed with a region then the "standard" differential approach is used to return the candidate peaks
258 
259  // This function is recursive, we start by finding the largest element of the dilation vector
260  Waveform::const_iterator maxItr = std::max_element(dilationStartItr,dilationStopItr);
261  float maxVal = *maxItr;
262 
263  // Check that the peak is of reasonable height...
264  if (maxVal < dilationThreshold) return;
265 
266  int maxBin = std::distance(dilationStartItr,maxItr);
267 
268  // The candidate hit region we want lies between the two nearest minima to the maximum we just found
269  // subject to the condition that the erosion vector has return to less than zero
270  Waveform::const_iterator firstDilItr = maxItr;
271  Waveform::const_iterator erosionItr = erosionStartItr + maxBin;
272 
273  float firstDilationValue = *firstDilItr;
274  float dilationCutValue = fDilationFraction * maxVal;
275  float erosionCutValue = fErosionFraction * maxVal;
276 
277  // Search for starting point
278  while(firstDilItr != dilationStartItr)
279  {
280  // Look for the turnover point
281  if (*erosionItr-- < erosionCutValue && *firstDilItr < dilationCutValue && *firstDilItr >= firstDilationValue) break;
282 
283  firstDilationValue = *firstDilItr--;
284  }
285 
286  // Set the start bin
287  int hitRegionStart = std::distance(dilationStartItr,firstDilItr);
288 
289  // Now go the other way
290  Waveform::const_iterator lastDilItr = maxItr;
291 
292  // Reset the local variables
293  float lastDilationValue = *lastDilItr;
294 
295  erosionItr = erosionStartItr + maxBin;
296 
297  // Search for starting point
298  while(++lastDilItr != dilationStopItr)
299  {
300  if (*++erosionItr <= erosionCutValue && *lastDilItr < dilationCutValue && *lastDilItr >= lastDilationValue) break;
301 
302  lastDilationValue = *lastDilItr;
303  }
304 
305  // Set the stop bin
306  int hitRegionStop = std::distance(dilationStartItr,lastDilItr);
307 
308  // Recursive call to find any hits in front of where we are now
309  if (hitRegionStart > fMinDeltaTicks)
310  findHitCandidates(derivStartItr, derivStartItr + hitRegionStart,
311  erosionStartItr, erosionStartItr + hitRegionStart,
312  dilationStartItr, dilationStartItr + hitRegionStart,
313  roiStartTick,
314  2. * fDilationThreshold,
315  hitCandidateVec);
316 
317  // Call the differential hit finding to get the actual hits within the region
318  findHitCandidates(derivStartItr + hitRegionStart,
319  derivStartItr + hitRegionStop,
320  roiStartTick + hitRegionStart,
323  hitCandidateVec);
324 
325  // Now call ourselves again to find any hits trailing the region we just identified
326  if (std::distance(lastDilItr,dilationStopItr) > fMinDeltaTicks)
327  findHitCandidates(derivStartItr + hitRegionStop, derivStopItr,
328  erosionStartItr + hitRegionStop, erosionStopItr,
329  dilationStartItr + hitRegionStop, dilationStopItr,
330  roiStartTick + hitRegionStop,
331  2. * fDilationThreshold,
332  hitCandidateVec);
333 
334  return;
335 }
336 
338  Waveform::const_iterator stopItr,
339  size_t roiStartTick,
340  int dTicksThreshold,
341  float dPeakThreshold,
342  HitCandidateVec& hitCandidateVec) const
343 {
344  // Search for candidate hits...
345  // The idea will be to find the largest deviation in the input derivative waveform as the starting point. Depending
346  // on if a maximum or minimum, we search forward or backward to find the minimum or maximum that our extremum
347  // corresponds to.
348  std::pair<Waveform::const_iterator, Waveform::const_iterator> minMaxPair = std::minmax_element(startItr, stopItr);
349 
350  Waveform::const_iterator maxItr = minMaxPair.second;
351  Waveform::const_iterator minItr = minMaxPair.first;
352 
353  // Use the larger of the two as the starting point and recover the nearest max or min
354  if (std::fabs(*maxItr) > std::fabs(*minItr)) minItr = findNearestMin(maxItr, stopItr);
355  else maxItr = findNearestMax(minItr,startItr);
356 
357  int deltaTicks = std::distance(maxItr,minItr);
358  float range = *maxItr - *minItr;
359 
360  // At some point small rolling oscillations on the waveform need to be ignored...
361  if (deltaTicks >= dTicksThreshold && range > dPeakThreshold)
362  {
363  // Need to back up to find zero crossing, this will be the starting point of our
364  // candidate hit but also the endpoint of the pre sub-waveform we'll search next
365  Waveform::const_iterator newEndItr = findStartTick(maxItr, startItr);
366 
367  int startTick = std::distance(startItr,newEndItr);
368 
369  // Now need to go forward to again get close to zero, this will then be the end point
370  // of our candidate hit and the starting point for the post sub-waveform to search
371  Waveform::const_iterator newStartItr = findStopTick(minItr, stopItr);
372 
373  int stopTick = std::distance(startItr,newStartItr);
374 
375  // Find hits in the section of the waveform leading up to this candidate hit
376  if (startTick > 2)
377  {
378  // Special handling for merged hits
379  if (*(newEndItr-1) > 0.) {dTicksThreshold = 2; dPeakThreshold = 0.; }
380  else {dTicksThreshold = fMinDeltaTicks; dPeakThreshold = fMinDeltaPeaks;}
381 
382  findHitCandidates(startItr,newEndItr+1,roiStartTick,dTicksThreshold,dPeakThreshold,hitCandidateVec);
383  }
384 
385  // Create a new hit candidate and store away
386  HitCandidate_t hitCandidate;
387 
388  Waveform::const_iterator peakItr = std::min_element(maxItr,minItr,[](const auto& left, const auto& right){return std::fabs(left) < std::fabs(right);});
389 
390  // Check balance
391  if (2 * std::distance(peakItr,minItr) < std::distance(maxItr,peakItr)) peakItr--;
392  else if (2 * std::distance(maxItr,peakItr) < std::distance(peakItr,minItr)) peakItr++;
393 
394  // Special handling of the start tick for multiple hits
395  size_t hitCandidateStartTick = roiStartTick + startTick;
396 
397  if (!hitCandidateVec.empty())
398  {
399  int deltaTicks = hitCandidateStartTick - hitCandidateVec.back().stopTick;
400 
401  if (deltaTicks > 0)
402  {
403  hitCandidateStartTick -= deltaTicks / 2;
404  hitCandidateVec.back().stopTick += deltaTicks / 2;
405  }
406  }
407 
408  hitCandidate.startTick = hitCandidateStartTick;
409  hitCandidate.stopTick = roiStartTick + stopTick;
410  hitCandidate.maxTick = roiStartTick + std::distance(startItr,maxItr);
411  hitCandidate.minTick = roiStartTick + std::distance(startItr,minItr);
412  hitCandidate.maxDerivative = maxItr != stopItr ? *maxItr : 0.;
413  hitCandidate.minDerivative = minItr != stopItr ? *minItr : 0.;
414  hitCandidate.hitCenter = roiStartTick + std::distance(startItr,peakItr) + 0.5;
415  hitCandidate.hitSigma = 0.5 * float(hitCandidate.minTick - hitCandidate.maxTick);
416  hitCandidate.hitHeight = hitCandidate.hitSigma * (hitCandidate.maxDerivative - hitCandidate.minDerivative) / 1.2130;
417 
418  hitCandidateVec.push_back(hitCandidate);
419 
420  // Finally, search the section of the waveform following this candidate for more hits
421  if (std::distance(newStartItr,stopItr) > 2)
422  {
423  // Special handling for merged hits
424  if (*(newStartItr+1) < 0.) {dTicksThreshold = 2; dPeakThreshold = 0.; }
425  else {dTicksThreshold = fMinDeltaTicks; dPeakThreshold = fMinDeltaPeaks;}
426 
427  findHitCandidates(newStartItr,stopItr,roiStartTick + stopTick,dTicksThreshold,dPeakThreshold,hitCandidateVec);
428  }
429  }
430 
431  return;
432 }
433 
435  const HitCandidateVec& hitCandidateVec,
436  MergeHitCandidateVec& mergedHitsVec) const
437 {
438  // If nothing on the input end then nothing to do
439  if (hitCandidateVec.empty()) return;
440 
441  // The idea is to group hits that "touch" so they can be part of common fit, those that
442  // don't "touch" are fit independently. So here we build the output vector to achieve that
443  // Get a container for the hits...
444  HitCandidateVec groupedHitVec;
445 
446  // Initialize the end of the last hit which we'll set to the first input hit's stop
447  size_t lastStopTick = hitCandidateVec.front().stopTick;
448 
449  // Step through the input hit candidates and group them by proximity
450  for(const auto& hitCandidate : hitCandidateVec)
451  {
452  // Small pulse height hits should not be considered?
453  if (hitCandidate.hitHeight > fMinHitHeight)
454  {
455  // Check condition that we have a new grouping
456  if (hitCandidate.startTick > lastStopTick + fNumInterveningTicks && !groupedHitVec.empty())
457  {
458  mergedHitsVec.emplace_back(groupedHitVec);
459 
460  groupedHitVec.clear();
461  }
462 
463  // Add the current hit to the current group
464  groupedHitVec.emplace_back(hitCandidate);
465 
466  lastStopTick = hitCandidate.stopTick;
467  }
468  }
469 
470  // Check end condition
471  if (!groupedHitVec.empty()) mergedHitsVec.emplace_back(groupedHitVec);
472 
473  return;
474 }
475 
477 {
478  // reset the min iterator and search forward to find the nearest minimum
479  Waveform::const_iterator lastItr = maxItr;
480 
481  // The strategy is simple...
482  // We are at a maximum so we search forward until we find the lowest negative point
483  while((lastItr + 1) != stopItr)
484  {
485  if (*(lastItr + 1) > *lastItr) break;
486 
487  lastItr++;
488  }
489 
490  // The minimum will be the last iterator value...
491  return lastItr;
492 }
493 
495 {
496  // Set the internal loop variable...
497  Waveform::const_iterator lastItr = minItr;
498 
499  // One extra condition to watch for here, make sure we can actually back up!
500  if (std::distance(startItr,minItr) > 0)
501  {
502  // Similar to searching for a maximum, we loop backward over ticks looking for the waveform to start decreasing
503  while((lastItr - 1) != startItr)
504  {
505  if (*(lastItr - 1) < *lastItr) break;
506 
507  lastItr--;
508  }
509  }
510 
511  return lastItr;
512 }
513 
515 {
516  Waveform::const_iterator lastItr = maxItr;
517 
518  // If we can't back up then there is nothing to do
519  if (std::distance(startItr,lastItr) > 0)
520  {
521  // In theory, we are starting at a maximum and want to find the "start" of the candidate peak
522  // Ideally we would look to search backward to the point where the (derivative) waveform crosses zero again.
523  // However, the complication is that we need to watch for the case where two peaks are merged together and
524  // we might run through another peak before crossing zero...
525  // So... loop until we hit the startItr...
526  Waveform::const_iterator loopItr = lastItr - 1;
527 
528  while(loopItr != startItr)
529  {
530  // Ideal world case, we cross zero... but we might encounter a minimum... or an inflection point
531  if (*loopItr < 0. || !(*loopItr < *lastItr)) break;
532 
533  lastItr = loopItr--;
534  }
535  }
536  else lastItr = startItr;
537 
538  return lastItr;
539 }
540 
542 {
543  Waveform::const_iterator lastItr = minItr;
544 
545  // If we can't go forward then there is really nothing to do
546  if (std::distance(minItr,stopItr) > 1)
547  {
548  // Pretty much the same strategy as for finding the start tick...
549  Waveform::const_iterator loopItr = lastItr + 1;
550 
551  while(loopItr != stopItr)
552  {
553  // Ideal case that we have crossed zero coming from a minimum... but watch for a maximum as well
554  if (*loopItr > 0. || !(*loopItr > *lastItr)) break;
555 
556  lastItr = loopItr++;
557  }
558  }
559 
560  return lastItr;
561 }
562 
564 }
Waveform::const_iterator findNearestMin(Waveform::const_iterator, Waveform::const_iterator) const
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:112
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:45
T * get() const
Definition: ServiceHandle.h:71
This is the interface class for tools/algorithms that perform various operations on waveforms...
Waveform::const_iterator findStopTick(Waveform::const_iterator, Waveform::const_iterator) const
std::vector< HitCandidate_t > HitCandidateVec
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
void MergeHitCandidates(const Waveform &, const HitCandidateVec &, MergeHitCandidateVec &) const override
struct HitCandidate{size_t startTick HitCandidate_t
TFileDirectory mkdir(std::string const &dir, std::string const &descr="")
CandHitMorphological(const fhicl::ParameterSet &pset)
Int_t max
Definition: plot.C:27
Waveform::const_iterator findStartTick(Waveform::const_iterator, Waveform::const_iterator) const
intermediate_table::const_iterator const_iterator
std::unique_ptr< reco_tool::IWaveformTool > fWaveformTool
T get(std::string const &key) const
Definition: ParameterSet.h:231
Description of geometry of one entire detector.
void findHitCandidates(const Waveform &, size_t, size_t, size_t, HitCandidateVec &) const override
void configure(const fhicl::ParameterSet &pset) override
std::map< int, TProfile * > HistogramMap
Definition: IWaveformTool.h:37
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
T * make(ARGS...args) const
TDirectory * dir
Definition: macro.C:5
This provides an interface for tools which are tasked with finding candidate hits on input waveforms...
Waveform::const_iterator findNearestMax(Waveform::const_iterator, Waveform::const_iterator) const
std::vector< HitCandidateVec > MergeHitCandidateVec
art framework interface to geometry description