LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
FFTHitFinder_module.cc
Go to the documentation of this file.
1 #ifndef FFTHITFINDER_H
2 #define FFTHITFINDER_H
3 //
5 // FFTHitFinder class
6 //
7 // pagebri3@msu.edu
8 //
9 // This algorithm is designed to find hits on wires after deconvolution
10 // with an average shape used as the input response.
12 
13 // C/C++ standard library
14 #include <string>
15 #include <vector>
16 #include <utility> // std::move
17 #include <algorithm> // std::accumulate
18 
19 // Framework includes
25 
26 // LArSoft Includes
28 #include "larcoreobj/SimpleTypesAndConstants/RawTypes.h" // raw::ChannelID_t
37 
38 // ROOT Includes
39 #include "TH1D.h"
40 #include "TDecompSVD.h"
41 #include "TMath.h"
42 #include "TF1.h"
43 
44 namespace hit{
45 
46  class FFTHitFinder : public art::EDProducer {
47 
48  public:
49 
50  explicit FFTHitFinder(fhicl::ParameterSet const& pset);
51  virtual ~FFTHitFinder();
52 
53  void produce(art::Event& evt);
54  void beginJob();
55  void endJob();
56  void reconfigure(fhicl::ParameterSet const& p);
57 
58  private:
59 
60  std::string fCalDataModuleLabel;
61  double fMinSigInd;
62  double fMinSigCol;
63  double fIndWidth;
64  double fColWidth;
65  double fIndMinWidth;
66  double fColMinWidth;
69  std::vector<double> fAreaNorms;
70  protected:
71 
72 
73  }; // class FFTHitFinder
74 
75  //-------------------------------------------------
77  {
78  this->reconfigure(pset);
79 
80  // let HitCollectionCreator declare that we are going to produce
81  // hits and associations with wires and raw digits
82  // (with no particular product label)
84  }
85 
86 
87  //-------------------------------------------------
89  {
90  }
91 
92  //-------------------------------------------------
94  {
95  fCalDataModuleLabel = p.get< std::string >("CalDataModuleLabel");
96  fMinSigInd = p.get< double >("MinSigInd");
97  fMinSigCol = p.get< double >("MinSigCol");
98  fIndWidth = p.get< double >("IndWidth");
99  fColWidth = p.get< double >("ColWidth");
100  fIndMinWidth = p.get< double >("IndMinWidth");
101  fColMinWidth = p.get< double >("ColMinWidth");
102  fMaxMultiHit = p.get< int >("MaxMultiHit");
103  fAreaMethod = p.get< int >("AreaMethod");
104  fAreaNorms = p.get< std::vector< double > >("AreaNorms");
105  }
106 
107  //-------------------------------------------------
109  {
110  }
111 
112  //-------------------------------------------------
114  {
115  }
116 
117  // This algorithm uses the fact that deconvolved signals are very smooth
118  // and looks for hits as areas between local minima that have signal above
119  // threshold.
120  //-------------------------------------------------
122  {
123 
124  // this object contains the hit collection
125  // and its associations to wires and raw digits:
126  recob::HitCollectionCreator hcol(*this, evt);
127 
128  // Read in the wire List object(s).
130  evt.getByLabel(fCalDataModuleLabel,wireVecHandle);
132 
133  // also get the raw digits associated with wires
134  art::FindOneP<raw::RawDigit> WireToRawDigits
135  (wireVecHandle, evt, fCalDataModuleLabel);
136 
137  std::vector<int> startTimes; // stores time of 1st local minimum
138  std::vector<int> maxTimes; // stores time of local maximum
139  std::vector<int> endTimes; // stores time of 2nd local minimum
140  int time = 0; // current time bin
141  int minTimeHolder = 0; // current start time
142  raw::ChannelID_t channel = raw::InvalidChannelID; // channel number
143  bool maxFound = false; // Flag for whether a peak > threshold has been found
144  double threshold = 0.; // minimum signal size for id'ing a hit
145  double fitWidth = 0.; // hit fit width initial value
146  double minWidth = 0.; // minimum hit width
147  std::string eqn = "gaus(0)"; // string for equation to fit
148  geo::SigType_t sigType = geo::kInduction;// type of plane we are looking at
149  std::stringstream numConv;
150 
151 
152  //loop over wires
153  for(unsigned int wireIter = 0; wireIter < wireVecHandle->size(); wireIter++) {
154 
155  art::Ptr<recob::Wire> wire(wireVecHandle, wireIter);
156  startTimes.clear();
157  maxTimes.clear();
158  endTimes.clear();
159  std::vector<float> signal(wire->Signal());
160  std::vector<float>::iterator timeIter; // iterator for time bins
161  time = 0;
162  minTimeHolder = 0;
163  maxFound = false;
164  channel = wire->Channel();
165  sigType = geom->SignalType(channel);
166 
167  //Set the appropriate signal widths and thresholds
168  if(sigType == geo::kInduction){
169  threshold = fMinSigInd;
170  fitWidth = fIndWidth;
171  minWidth = fIndMinWidth;
172  }
173  else if(sigType == geo::kCollection){
174  threshold = fMinSigCol;
175  fitWidth = fColWidth;
176  minWidth = fColMinWidth;
177  }
178  // loop over signal
179  for(timeIter = signal.begin(); timeIter+2 < signal.end(); timeIter++){
180  //test if timeIter+1 is a local minimum
181  if(*timeIter > *(timeIter+1) && *(timeIter+1) < *(timeIter+2)){
182  //only add points if already found a local max above threshold.
183  if(maxFound) {
184  endTimes.push_back(time+1);
185  maxFound = false;
186  //keep these in case new hit starts right away
187  minTimeHolder = time+2;
188  }
189  else minTimeHolder = time+1;
190  }
191  //if not a minimum, test if we are at a local maximum
192  //if so, and the max value is above threshold, add it and proceed.
193  else if(*timeIter < *(timeIter+1) &&
194  *(timeIter+1) > *(timeIter+2) &&
195  *(timeIter+1) > threshold){
196  maxFound = true;
197  maxTimes.push_back(time+1);
198  startTimes.push_back(minTimeHolder);
199  }
200  time++;
201  }//end loop over signal vec
202 
203 
204  //if no inflection found before end, but peak found add end point
205  while(maxTimes.size()>endTimes.size())
206  endTimes.push_back(signal.size()-1);
207  if(startTimes.size() == 0) continue;
208 
209  //All code below does the fitting, adding of hits
210  //to the hit vector and when all wires are complete
211  //saving them
212  double totSig(0); // stores the total hit signal
213  double startT(0); // stores the start time
214  double endT(0); // stores the end time
215  int numHits(0); // number of consecutive hits being fitted
216  int size(0); // size of data vector for fit
217  int hitIndex(0); // index of current hit in sequence
218  double amplitude(0), position(0), width(0); //fit parameters
219  double amplitudeErr(0), positionErr(0), widthErr(0); //fit errors
220  double goodnessOfFit(0), chargeErr(0); //Chi2/NDF and error on charge
221  double minPeakHeight(0); //lowest peak height in multi-hit fit
222 
223  //stores gaussian paramters first index is the hit number
224  //the second refers to height, position, and width respectively
225  std::vector<double> hitSig;
226 
227  //add found hits to hit vector
228  while(hitIndex < (signed)startTimes.size()) {
229 
230  startT = endT = 0;
231  numHits = 1;
232  minPeakHeight = signal[maxTimes[hitIndex]];
233 
234  //consider adding pulse to group of consecutive hits if:
235  //1 less than max consecutive hits
236  //2 we are not at the last point in the signal vector
237  //3 the height of the dip between the two is greater than threshold/2
238  //4 and there is no gap between them
239  while(numHits < fMaxMultiHit &&
240  numHits+hitIndex < (signed)endTimes.size() &&
241  signal[endTimes[hitIndex+numHits-1]] >threshold/2.0 &&
242  startTimes[hitIndex+numHits] - endTimes[hitIndex+numHits-1] < 2){
243 
244  if(signal[maxTimes[hitIndex+numHits]] < minPeakHeight)
245  minPeakHeight = signal[maxTimes[hitIndex+numHits]];
246 
247  ++numHits;
248  }
249 
250  //finds the first point > 1/2 the smallest peak
251  startT = startTimes[hitIndex];
252 
253  while(signal[(int)startT] < minPeakHeight/2.0) ++startT;
254 
255  //finds the first point from the end > 1/2 the smallest peak
256  endT = endTimes[hitIndex+numHits-1];
257 
258  while(signal[(int)endT] <minPeakHeight/2.0) --endT;
259  size = (int)(endT-startT);
260  TH1D hitSignal("hitSignal","",size,startT,endT);
261  for(int i = (int)startT; i < (int)endT; ++i)
262  hitSignal.Fill(i,signal[i]);
263 
264  //build the TFormula
265  eqn = "gaus(0)";
266 
267  for(int i = 3; i < numHits*3; i+=3) {
268  eqn.append("+gaus(");
269  numConv.str("");
270  numConv << i;
271  eqn.append(numConv.str());
272  eqn.append(")");
273  }
274 
275  TF1 gSum("gSum",eqn.c_str(),0,size);
276 
277  if(numHits > 1) {
278  TArrayD data(numHits*numHits);
279  TVectorD amps(numHits);
280  for(int i = 0; i < numHits; ++i) {
281  amps[i] = signal[maxTimes[hitIndex+i]];
282  for(int j = 0; j < numHits;j++)
283  data[i+numHits*j] = TMath::Gaus(maxTimes[hitIndex+j],
284  maxTimes[hitIndex+i],
285  fitWidth);
286  }//end loop over hits
287 
288  //This section uses a linear approximation in order to get an
289  //initial value of the individual hit amplitudes
290  try{
291  TMatrixD h(numHits,numHits);
292  h.Use(numHits,numHits,data.GetArray());
293  TDecompSVD a(h);
294  a.Solve(amps);
295  }
296  catch(...){
297  mf::LogInfo("FFTHitFinder")<<"TDcompSVD failed";
298  hitIndex += numHits;
299  continue;
300  }
301 
302  for(int i = 0; i < numHits; ++i) {
303  //if the approximation makes a peak vanish
304  //set initial height as average of threshold and
305  //raw peak height
306  if(amps[i] > 0 ) amplitude = amps[i];
307  else amplitude = 0.5*(threshold+signal[maxTimes[hitIndex+i]]);
308  gSum.SetParameter(3*i,amplitude);
309  gSum.SetParameter(1+3*i, maxTimes[hitIndex+i]);
310  gSum.SetParameter(2+3*i, fitWidth);
311  gSum.SetParLimits(3*i, 0.0, 3.0*amplitude);
312  gSum.SetParLimits(1+3*i, startT , endT);
313  gSum.SetParLimits(2+3*i, 0.0, 10.0*fitWidth);
314  }//end loop over hits
315  }//end if numHits > 1
316  else {
317  gSum.SetParameters(signal[maxTimes[hitIndex]],maxTimes[hitIndex],fitWidth);
318  gSum.SetParLimits(0,0.0,1.5*signal[maxTimes[hitIndex]]);
319  gSum.SetParLimits(1, startT , endT);
320  gSum.SetParLimits(2,0.0,10.0*fitWidth);
321  }
322 
324  hitSignal.Fit(&gSum,"QNRW","", startT, endT);
325  for(int hitNumber = 0; hitNumber < numHits; ++hitNumber) {
326  totSig = 0;
327  if(gSum.GetParameter(3*hitNumber) > threshold/2.0 &&
328  gSum.GetParameter(3*hitNumber+2) > minWidth) {
329  amplitude = gSum.GetParameter(3*hitNumber);
330  position = gSum.GetParameter(3*hitNumber+1);
331  width = gSum.GetParameter(3*hitNumber+2);
332  amplitudeErr = gSum.GetParError(3*hitNumber);
333  positionErr = gSum.GetParError(3*hitNumber+1);
334  widthErr = gSum.GetParError(3*hitNumber+2);
335  goodnessOfFit = gSum.GetChisquare()/(double)gSum.GetNDF();
336  int DoF = gSum.GetNDF();
337 
338  //estimate error from area of Gaussian
339  chargeErr = std::sqrt(TMath::Pi())*(amplitudeErr*width+widthErr*amplitude);
340 
341  hitSig.resize(size);
342 
343  for(int sigPos = 0; sigPos < size; ++sigPos){
344  hitSig[sigPos] = amplitude*TMath::Gaus(sigPos+startT,position, width);
345  totSig += hitSig[(int)sigPos];
346  }
347 
348  if(fAreaMethod)
349  totSig = std::sqrt(2*TMath::Pi())*amplitude*width/fAreaNorms[(size_t)sigType];
350 
351  // get the WireID for this hit
352  std::vector<geo::WireID> wids = geom->ChannelToWire(channel);
354  // for now, just take the first option returned from ChannelToWire
355  geo::WireID wid = wids[0];
356 
357  // make the hit
359  *wire, // wire
360  wid, // wireID
361  (int) startT, // start_tick
362  (int) endT, // end_tick
363  width, // rms
364  position, // peak_time
365  positionErr, // sigma_peak_time
366  amplitude, // peak_amplitude
367  amplitudeErr, // sigma_peak_amplitude
368  totSig, // hit_integral
369  chargeErr, // hit_sigma_integral
370  std::accumulate // summedADC
371  (signal.begin() + (int) startT, signal.begin() + (int) endT, 0.),
372  1, // multiplicity
373  -1, // local_index
375  goodnessOfFit, // goodness_of_fit
376  DoF // dof
377  );
378 
379  // get the object associated with the original hit
380  art::Ptr<raw::RawDigit> rawdigits = WireToRawDigits.at(wireIter);
381 
382  hcol.emplace_back(hit.move(), wire, rawdigits);
383  }//end if over threshold
384  }//end loop over hits
385  hitIndex += numHits;
386  } // end while on hitIndex<(signed)startTimes.size()
387 
388  } // while on Wires
389 
390  // put the hit collection and associations into the event
391  hcol.put_into(evt);
392 
393  } // End of produce()
394 
395 
396 
398 
399 } // end of hit namespace
400 
401 
402 #endif // FFTHITFINDER_H
double fIndMinWidth
Minimum induction hit width.
void reconfigure(fhicl::ParameterSet const &p)
FFTHitFinder(fhicl::ParameterSet const &pset)
Encapsulate the construction of a single cyostat.
intermediate_table::iterator iterator
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
double fMinSigInd
Induction signal height threshold.
Declaration of signal hit object.
double fIndWidth
Initial width for induction fit.
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
Definition of basic raw digits.
int fMaxMultiHit
maximum hits for multi fit
double fColWidth
Initial width for collection fit.
SigType_t SignalType(geo::PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
static void declare_products(ModuleType &producer, std::string instance_name="", bool doWireAssns=true, bool doRawDigitAssns=true)
Declares the hit products we are going to fill.
Definition: HitCreator.h:1117
std::vector< double > fAreaNorms
factors for converting area to same units as peak height
void produce(art::Event &evt)
Class managing the creation of a new recob::Hit object.
Definition: HitCreator.h:84
Helper functions to create a hit.
std::string fCalDataModuleLabel
constexpr ChannelID_t InvalidChannelID
ID of an invalid channel.
Definition: RawTypes.h:31
A class handling a collection of hits and its associations.
Definition: HitCreator.h:513
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
Signal from induction planes.
Definition: geo_types.h:92
raw::ChannelID_t Channel() const
Returns the ID of the channel (or InvalidChannelID)
Definition: Wire.h:164
enum geo::_plane_sigtype SigType_t
Enumerate the possible plane projections.
double fColMinWidth
Minimum collection hit width.
T get(std::string const &key) const
Definition: ParameterSet.h:231
void emplace_back(recob::Hit &&hit, art::Ptr< recob::Wire > const &wire=art::Ptr< recob::Wire >(), art::Ptr< raw::RawDigit > const &digits=art::Ptr< raw::RawDigit >())
Adds the specified hit to the data collection.
Definition: HitCreator.cxx:245
Definition of data types for geometry description.
void put_into(art::Event &)
Moves the data into an event.
Definition: HitCreator.h:663
Detector simulation of raw signals on wires.
std::vector< float > Signal() const
Return a zero-padded full length vector filled with RoI signal.
Definition: Wire.cxx:47
Encapsulate the construction of a single detector plane.
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
double fMinSigCol
Collection signal height threshold.
Declaration of basic channel signal object.
int fAreaMethod
Type of area calculation.
TCEvent evt
Definition: DataStructs.cxx:5
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:27
Definition: fwd.h:25
art framework interface to geometry description
Encapsulate the construction of a single detector plane.
Signal from collection planes.
Definition: geo_types.h:93