LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
RawHitFinder_module.cc
Go to the documentation of this file.
1 #ifndef RAWHITFINDER_H
2 #define RAWHITFINDER_H
3 
5 //HIT FINDER THAT RUNS ON RAW SIGNALS INSTEAD OF DECONVOLUTED ONES.
6 //WRITTEN INITIALLY FOR DUNE 35T ONLINE FILTER.
7 //abooth@fnal.gov, mstancar@fnal.gov.
9 
10 // C/C++ standard libraries
11 #include <string>
12 #include <vector>
13 #include <fstream>
14 #include <set>
15 
16 //Framework
17 #include "fhiclcpp/ParameterSet.h"
26 
27 //LArSoft
32 #include "lardataobj/RawData/raw.h"
35 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
36 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
37 
38 //LArSoft From FFT
39 
45 
46 //ROOT from CalData
47 
48 #include "TComplex.h"
49 #include "TFile.h"
50 #include "TH2D.h"
51 #include "TF1.h"
52 
53 //ROOT From Gauss
54 
55 #include "TH1D.h"
56 #include "TDecompSVD.h"
57 #include "TMath.h"
58 
59 namespace hit {
60 
61  class RawHitFinder : public art::EDProducer {
62 
63  public:
64 
65  explicit RawHitFinder(fhicl::ParameterSet const& pset);
66  virtual ~RawHitFinder();
67 
68  void produce(art::Event& evt);
69  void beginJob();
70  void endJob();
71  void reconfigure(fhicl::ParameterSet const& p);
72 
73 
74  private:
75  unsigned int fDataSize; //SIZE OF RAW DATA ON ONE WIRE.
76  art::InputTag fDigitModuleLabel; //MODULE THAT MADE DIGITS.
77  std::string fSpillName; //NOMINAL SPILL IS AN EMPTY STRING.
78 
79  //FFT COPIED VARIABLES.
80  std::string fCalDataModuleLabel;
81  std::string fHitLabelName;
82  double fMinSigInd; //INDUCTION SIGNAL HEIGHT THRESHOLD.
83  double fMinSigCol; //COLLECTION SIGNAL HEIGHT THRESHOLD.
84  double fIndWidth; //INITIAL WIDTH FOR COLLECTION FIT.
85  double fColWidth; //INITIAL WIDTH FOR COLLECTION FIT.
86  double fIndMinWidth; //MINIMUM INDUCTION HIT WIDTH.
87  double fColMinWidth; //MINIMUM COLLECTION HIT WIDTH.
88  int fMaxMultiHit; //MAXIMUM HITS FOR MULTIFIT.
89  int fAreaMethod; //TYPE OF AREA CALCULATION.
90  std::vector<double> fAreaNorms; //FACTORS FOR CONVERTING AREA TO SAME UNIT AS PEAK HEIGHT.
91  bool fUncompressWithPed; //OPTION TO UNCOMPRESS WITH PEDESTAL.
92  bool fSkipInd; //OPTION TO SKIP INDUCTION PLANES
93  double fIncludeMoreTail; //FRACTION OF THE HIT WIDTH INCLUDED BELOW THRESHOLD IN
94  // THE CALCULATION OF CHARGE. COLLECTION PLANE ONLY
95  int fColMinWindow; // Minimum length of integration window for charge in ticks
96  int fIndCutoff; //MAX WIDTH FOR EARLY SIDE OF INDUCTION HIT IN TICKS
97 
98  protected:
99 
100  }; // class RawHitFinder
101 
102 
103  //-------------------------------------------------
105  {
106  this->reconfigure(pset);
107 
108  //LET HITCOLLECTIONCREATOR DECLARE THAT WE ARE GOING TO PRODUCE
109  //HITS AND ASSOCIATIONS TO RAW DIGITS BUT NOT ASSOCIATIONS TO WIRES
110  //(WITH NO PARTICULAR PRODUCT LABEL).
112  /*instance_name*/"",
113  /*doWireAssns*/false,
114  /*doRawDigitAssns*/true);
115  }
116 
117  //-------------------------------------------------
119  {
120  }
121 
123  {
124  fDigitModuleLabel = p.get< art::InputTag >("DigitModuleLabel", "daq");
125  fCalDataModuleLabel = p.get< std::string >("CalDataModuleLabel");
126  fMinSigInd = p.get< double >("MinSigInd");
127  fMinSigCol = p.get< double >("MinSigCol");
128  fIncludeMoreTail = p.get< double >("IncludeMoreTail", 0.);
129  fIndWidth = p.get< int >("IndWidth",20);
130  fColWidth = p.get< double >("ColWidth");
131  fIndMinWidth = p.get< double >("IndMinWidth");
132  fColMinWidth = p.get< double >("ColMinWidth",0.);
133  fMaxMultiHit = p.get< int >("MaxMultiHit");
134  fAreaMethod = p.get< int >("AreaMethod");
135  fAreaNorms = p.get< std::vector< double > >("AreaNorms");
136  fUncompressWithPed = p.get< bool >("UncompressWithPed", true);
137  fSkipInd = p.get< bool >("SkipInd", false);
138  fColMinWindow = p.get< int >("ColMinWindow",0);
139  fIndCutoff = p.get< int >("IndCutoff",20);
140  mf::LogInfo("RawHitFinder_module") << "fDigitModuleLabel: " << fDigitModuleLabel << std::endl;
141  }
142 
143  //-------------------------------------------------
145  {
146  }
147 
149  {
150  }
151 
152  //-------------------------------------------------
154  {
155  //GET THE GEOMETRY.
157 
158  //READ IN THE DIGIT LIST OBJECTS(S).
160 
161  bool retVal = evt.getByLabel(fDigitModuleLabel, digitVecHandle);
162  if(retVal==true)
163  mf::LogInfo("RawHitFinder_module") << "I got fDigitModuleLabel: " << fDigitModuleLabel << std::endl;
164  else
165  mf::LogWarning("RawHitFinder_module") << "Could not get fDigitModuleLabel: " << fDigitModuleLabel << std::endl;
166 
167  std::vector<float> holder; //HOLDS SIGNAL DATA.
168  std::vector<short> rawadc; //UNCOMPRESSED ADC VALUES.
169 
170  // ###############################################
171  // ### Making a ptr vector to put on the event ###
172  // ###############################################
173  // THIS CONTAINS THE HIT COLLECTION AND ITS ASSOCIATIONS TO WIRES AND RAW DIGITS.
174  recob::HitCollectionCreator hcol(*this, evt, false /* doWireAssns */, true /* doRawDigitAssns */);
175 
176  std::vector<float> startTimes; //STORES TIME OF WINDOW START.
177  std::vector<float> maxTimes; //STORES TIME OF LOCAL MAXIMUM.
178  std::vector<float> endTimes; //STORES TIME OF WINDOW END.
179  std::vector<float> peakHeight; //STORES ADC COUNT AT THE MAXIMUM.
180  std::vector<float> hitrms; //STORES CHARGE WEIGHTED RMS OF TIME ACROSS THE HIT.
181  std::vector<double> charge; //STORES THE TOTAL CHARGE ASSOCIATED WITH THE HIT.
182 
183  uint32_t channel = 0; //CHANNEL NUMBER.
184  double threshold = 0; //MINIMUM SIGNAL SIZE FOR ID'ING A HIT.
185  double totSig = 0;
186  double myrms = 0;
187  double mynorm = 0;
189  std::stringstream numConv;
190 
191  hcol.reserve(digitVecHandle->size());
192  for(size_t rdIter = 0; rdIter < digitVecHandle->size(); ++rdIter){
193  holder.clear();
194 
195  //GET THE REFERENCE TO THE CURRENT raw::RawDigit.
196  art::Ptr<raw::RawDigit> digitVec(digitVecHandle, rdIter);
197  channel = digitVec->Channel();
198  fDataSize = digitVec->Samples();
199 
200  rawadc.resize(fDataSize);
201  holder.resize(fDataSize);
202 
203  //UNCOMPRESS THE DATA.
204  if (fUncompressWithPed){
205  int pedestal = (int)digitVec->GetPedestal();
206  raw::Uncompress(digitVec->ADCs(), rawadc, pedestal, digitVec->Compression());
207  }
208  else{
209  raw::Uncompress(digitVec->ADCs(), rawadc, digitVec->Compression());
210  }
211 
212  //GET THE LIST OF BAD CHANNELS.
213  lariov::ChannelStatusProvider const& channelStatus
215 
216  lariov::ChannelStatusProvider::ChannelSet_t const BadChannels
217  = channelStatus.BadChannels();
218 
219  for(unsigned int bin = 0; bin < fDataSize; ++bin){
220  holder[bin]=(rawadc[bin]-digitVec->GetPedestal());
221  }
222 
223  sigType = geom->SignalType(channel);
224 
225  peakHeight.clear();
226  endTimes.clear();
227  startTimes.clear();
228  maxTimes.clear();
229  charge.clear();
230  hitrms.clear();
231 
232  bool channelSwitch = false;
233 
234  for(auto it = BadChannels.begin(); it != BadChannels.end(); it++)
235  {
236  if(channel==*it)
237  {
238  channelSwitch = true;
239  break;
240  }
241  }
242 
243  if(channelSwitch==false)
244  {
245  // ###############################################
246  // ### Induction Planes ###
247  // ###############################################
248 
249  //THE INDUCTION PLANE METHOD HAS NOT YET BEEN MODIFIED AND TESTED FOR REAL DATA.
250  // Or for detectors without a grid plane
251  //
252  if(sigType == geo::kInduction && !fSkipInd){
253  threshold = fMinSigInd;
254  // std::cout<< "Threshold is " << threshold << std::endl;
255  // fitWidth = fIndWidth;
256  // minWidth = fIndMinWidth;
257  // continue;
258  float negthr=-1.0*threshold;
259  unsigned int bin =1;
260  float minadc=0;
261 
262  // find the dips
263  while (bin<(fDataSize-1)) { // loop over ticks
264  float thisadc = holder[bin]; float nextadc = holder[bin+1];
265  if (thisadc<negthr && nextadc < negthr) { // new region, require two ticks above threshold
266  // std::cout << "new region" << bin << " " << thisadc << std::endl;
267  // step back to find zero crossing
268  unsigned int place = bin;
269  while (thisadc<=0 && bin>0) {
270  // std::cout << bin << " " << thisadc << std::endl;
271  bin--;
272  thisadc=holder[bin];
273  }
274  float hittime = bin+thisadc/(thisadc-holder[bin+1]);
275  maxTimes.push_back(hittime);
276 
277  // step back more to find the hit start time
278  uint32_t stop;
279  if (fIndCutoff<(int)bin) {stop=bin-fIndCutoff;} else {stop=0;}
280  while (thisadc<threshold && bin>stop) {
281  // std::cout << bin << " " << thisadc << std::endl;
282  bin--;
283  thisadc=holder[bin];
284  }
285  if (bin>=2) bin-=2;
286  while (thisadc>threshold && bin>stop) {
287  // std::cout << bin << " " << thisadc << std::endl;
288  bin--;
289  thisadc=holder[bin];
290  }
291  startTimes.push_back(bin+1);
292  // now step forward from hit time to find end time, area of dip
293  bin=place;
294  thisadc=holder[bin];
295  minadc=thisadc;
296  bin++;
297  totSig = fabs(thisadc);
298  while (thisadc<negthr && bin<fDataSize) {
299  totSig += fabs(thisadc);
300  thisadc=holder[bin];
301  if (thisadc<minadc) minadc=thisadc;
302  bin++;
303  }
304  endTimes.push_back(bin-1);
305  peakHeight.push_back(-1.0*minadc);
306  charge.push_back(totSig);
307  hitrms.push_back(5.0);
308  // std::cout << "TOTAL SIGNAL INDUCTION " << totSig << " 5.0" << std::endl;
309  // std::cout << "filled end times " << bin-1 << "peak height vector size " << peakHeight.size() << std::endl;
310 
311  // don't look for a new hit until it returns to baseline
312  while (thisadc<0 && bin<fDataSize) {
313  // std::cout << bin << " " << thisadc << std::endl;
314  bin++;
315  if (bin == fDataSize) break;
316  thisadc=holder[bin];
317  }
318  } // end region
319  bin++;
320  }// loop over ticks
321  }
322 
323  // ###############################################
324  // ### Collection Plane ###
325  // ###############################################
326 
327  else if(sigType == geo::kCollection)
328  {
329  threshold = fMinSigCol;
330 
331  float madc = threshold;
332  int ibin = 0;
333  int start = 0;
334  int end = 0;
335  unsigned int bin = 0;
336 
337  while (bin<fDataSize)
338  {
339  float thisadc = holder[bin];
340  madc = threshold;
341  ibin = 0;
342 
343  if (thisadc>madc)
344  {
345  start = bin;
346 
347  if(thisadc>threshold && bin<fDataSize)
348  {
349  while (thisadc>threshold && bin<fDataSize)
350  {
351  if (thisadc>madc)
352  {
353  ibin=bin;
354  madc=thisadc;
355  }
356  bin++;
357  if (bin == fDataSize) break;
358  thisadc=holder[bin];
359  }
360  }
361  else
362  {
363  bin++;
364  }
365 
366  end = bin-1;
367 
368  if(start!=end)
369  {
370  maxTimes.push_back(ibin);
371  peakHeight.push_back(madc);
372  startTimes.push_back(start);
373  endTimes.push_back(end);
374 
375  totSig = 0;
376  myrms = 0;
377  mynorm = 0;
378 
379  int moreTail = std::ceil(fIncludeMoreTail*(end-start));
380  if (moreTail<fColMinWindow) moreTail=fColMinWindow;
381 
382  for(int i = start-moreTail; i <= end+moreTail; i++)
383  {
384  if(i<(int)(holder.size()) && i>=0)
385  {
386  float temp = ibin-i;
387  myrms += temp*temp*holder[i];
388 
389  totSig += holder[i];
390  }
391  }
392 
393  charge.push_back(totSig);
394  mynorm = totSig;
395  myrms/=mynorm;
396  hitrms.push_back(sqrt(myrms));
397 
398  //PRE CHANGES MADE 04/14/16. A BOOTH, DUNE 35T.
399  /*
400  int moreTail = std::ceil(fIncludeMoreTail*(end-start));
401 
402  for(int i = start-moreTail; i <= end+moreTail; i++)
403  {
404  totSig += holder[i];
405  float temp2 = holder[i]*holder[i];
406  mynorm += temp2;
407  float temp = ibin-i;
408  myrms += temp*temp*temp2;
409  }
410 
411  charge.push_back(totSig);
412  myrms/=mynorm;
413  if((end-start+2*moreTail+1)!=0)
414  {
415  myrms/=(float)(end-start+2*moreTail+1);
416  hitrms.push_back(sqrt(myrms));
417  }
418  else
419  {
420  hitrms.push_back(sqrt(myrms));
421  }*/
422  }
423  }
424  start = 0;
425  end = 0;
426  bin++;
427  }
428  }
429  }
430 
431  int numHits(0); //NUMBER OF CONSECUTIVE HITS BEING FITTED.
432  int hitIndex(0); //INDEX OF CURRENT HIT IN SEQUENCE.
433  double amplitude(0), position(0); //FIT PARAMETERS.
434  double start(0), end(0);
435  double amplitudeErr(0), positionErr(0); //FIT ERRORS.
436  double goodnessOfFit(0), chargeErr(0); //CHI2/NDF and error on charge.
437  double hrms(0);
438 
439  numHits = maxTimes.size();
440  for (int i = 0; i < numHits; ++i)
441  {
442  amplitude = peakHeight[i];
443  position = maxTimes[i];
444  start = startTimes[i];
445  end = endTimes[i];
446  hrms = hitrms[i];
447  amplitudeErr = -1;
448  positionErr = 1.0;
449  goodnessOfFit = -1;
450  chargeErr = -1;
451  totSig = charge[i];
452 
453 
454  std::vector<geo::WireID> wids = geom->ChannelToWire(channel);
455  geo::WireID wid = wids[0];
456 
457  if (start>=end)
458  {
459  mf::LogWarning("RawHitFinder_module") << "Hit start " << start << " is >= hit end " << end;
460  continue;
461  }
462 
464  *digitVec, //RAW DIGIT REFERENCE.
465  wid, //WIRE ID.
466  start, //START TICK.
467  end, //END TICK.
468  hrms, //RMS.
469  position, //PEAK_TIME.
470  positionErr, //SIGMA_PEAK_TIME.
471  amplitude, //PEAK_AMPLITUDE.
472  amplitudeErr, //SIGMA_PEAK_AMPLITUDE.
473  totSig, //HIT_INTEGRAL.
474  chargeErr, //HIT_SIGMA_INTEGRAL.
475  std::accumulate(holder.begin() + (int) start, holder.begin() + (int) end, 0.), //SUMMED CHARGE.
476  1, //MULTIPLICITY.
477  -1, //LOCAL_INDEX.
478  goodnessOfFit, //WIRE ID.
479  int(end-start+1) //DEGREES OF FREEDOM.
480  );
481  hcol.emplace_back(hit.move(), digitVec);
482 
483  ++hitIndex;
484  }
485  }
486 
487  hcol.put_into(evt);
488  }
489 
491 
492 } // end namespace hit
493 
494 
495 #endif //RAWHITFINDER_H
float GetPedestal() const
Definition: RawDigit.h:213
const ADCvector_t & ADCs() const
Reference to the compressed ADC count vector.
Definition: RawDigit.h:209
Encapsulate the construction of a single cyostat.
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Declaration of signal hit object.
void reconfigure(fhicl::ParameterSet const &p)
ChannelID_t Channel() const
DAQ channel this raw data was read from.
Definition: RawDigit.h:211
std::vector< double > fAreaNorms
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.
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
void produce(art::Event &evt)
unsigned short Samples() const
Number of samples in the uncompressed ADC data.
Definition: RawDigit.h:212
Class managing the creation of a new recob::Hit object.
Definition: HitCreator.h:84
Helper functions to create a hit.
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
void reserve(size_t new_size)
Prepares the collection to host at least new_size hits.
Definition: HitCreator.h:651
enum geo::_plane_sigtype SigType_t
Enumerate the possible plane projections.
Collect all the RawData header files together.
T get(std::string const &key) const
Definition: ParameterSet.h:231
RawHitFinder(fhicl::ParameterSet const &pset)
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
float bin[41]
Definition: plottest35.C:14
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.
raw::Compress_t Compression() const
Compression algorithm used to store the ADC counts.
Definition: RawDigit.h:215
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
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
std::string fCalDataModuleLabel
art::InputTag fDigitModuleLabel
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
TCEvent evt
Definition: DataStructs.cxx:5
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:756
art framework interface to geometry description
Encapsulate the construction of a single detector plane.
Signal from collection planes.
Definition: geo_types.h:93