All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RecoWireICARUS_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // $Id: RecoWireICARUS.cxx,v 1.36 2010/09/15 bpage Exp $
3 //
4 // RecoWireICARUS class
5 //
6 // brebel@fnal.gov
7 //
8 ////////////////////////////////////////////////////////////////////////
9 #ifndef RecoWireICARUS_H
10 #define RecoWireICARUS_H
11 
12 // ROOT includes
13 #include <TFile.h>
14 #include <TH2D.h>
15 #include <TH1D.h>
16 #include <TF1.h>
17 #include <TComplex.h>
18 
19 // Framework includes
20 #include "fhiclcpp/ParameterSet.h"
21 #include "messagefacility/MessageLogger/MessageLogger.h"
22 #include "cetlib_except/exception.h"
23 #include "cetlib/search_path.h"
24 #include "canvas/Persistency/Common/Ptr.h"
25 #include "canvas/Persistency/Common/Assns.h"
26 #include "art/Framework/Principal/Event.h"
27 #include "art/Framework/Principal/Handle.h"
28 #include "art/Framework/Services/Registry/ServiceHandle.h"
29 #include "art/Framework/Core/EDProducer.h" // include the proper bit of the framework
30 #include "art/Framework/Core/ModuleMacros.h"
31 
32 // LArSoft includes
33 #include "larcoreobj/SimpleTypesAndConstants/RawTypes.h" // raw::ChannelID_t
37 #include "lardataobj/RawData/raw.h"
41 
42 #include "icarus_signal_processing/Filters/ICARUSFFT.h"
43 
44 ///creation of calibrated signals on wires
45 namespace recowire {
46 
47  class RecoWireICARUS : public art::EDProducer {
48 
49  public:
50 
51  // create calibrated signals on wires. this class runs
52  // an fft to remove the electronics shaping.
53  explicit RecoWireICARUS(fhicl::ParameterSet const& pset);
54  virtual ~RecoWireICARUS();
55 
56  void produce(art::Event& evt);
57  void beginJob();
58  void endJob();
59  void reconfigure(fhicl::ParameterSet const& p);
60 
61  private:
62 
63  std::string fResponseFile; ///< response file containing transformed
64  ///< shape histograms and decay constants
65  int fExpEndBins; ///< number of end bins to consider for tail fit
66  int fPostsample; ///< number of postsample bins
67  std::string fDigitModuleLabel; ///< module that made digits
68 
69  std::vector<std::vector<std::complex<double>>> fKernelR; ///< holds transformed induction
70  ///< response function
71  std::vector<std::vector<std::complex<double>>> fKernelS; ///< holds transformed induction
72  ///< response function
73  std::vector<double> fDecayConstsR; ///< vector holding RC decay
74  ///< constants
75  std::vector<double> fDecayConstsS; ///< vector holding RC decay
76  ///< constants
77  std::vector<int> fKernMapR; ///< map telling which channels
78  ///< have which response functions
79  std::vector<int> fKernMapS; ///< map telling which channels
80  ///< have which response functions
81  protected:
82 
83  }; // class RecoWireICARUS
84 }
85 
86 namespace recowire{
87 
88  //-------------------------------------------------
89  RecoWireICARUS::RecoWireICARUS(fhicl::ParameterSet const& pset) : EDProducer{pset}
90  {
91  this->reconfigure(pset);
92 
93  produces< std::vector<recob::Wire> >();
94  produces<art::Assns<raw::RawDigit, recob::Wire>>();
95 
96  }
97 
98  //-------------------------------------------------
100  {
101  }
102 
103  //////////////////////////////////////////////////////
104  void RecoWireICARUS::reconfigure(fhicl::ParameterSet const& p)
105  {
106  fDigitModuleLabel = p.get< std::string >("DigitModuleLabel", "daq");
107  cet::search_path sp("FW_SEARCH_PATH");
108  sp.find_file(p.get<std::string>("ResponseFile"), fResponseFile);
109  fExpEndBins = p.get< int > ("ExponentialEndBins");
110  fPostsample = p.get< int > ("PostsampleBins");
111  }
112 
113  //-------------------------------------------------
115  {
116 
117  MF_LOG_DEBUG("RecoWireICARUS") << "RecoWireICARUS_plugin: Opening Electronics Response File: "
118  << fResponseFile.c_str();
119 
120  TFile f(fResponseFile.c_str());
121  if( f.IsZombie() )
122  mf::LogWarning("RecoWireICARUS") << "Cannot open response file "
123  << fResponseFile.c_str();
124 
125  TH2D *respRe = dynamic_cast<TH2D*>(f.Get("real/RespRe") );
126  TH2D *respIm = dynamic_cast<TH2D*>(f.Get("real/RespIm") );
127  TH1D *decayHist = dynamic_cast<TH1D*>(f.Get("real/decayHist"));
128  unsigned int wires = decayHist->GetNbinsX();
129  unsigned int bins = respRe->GetYaxis()->GetNbins();
130  unsigned int bin = 0;
131  unsigned int wire = 0;
132  fDecayConstsR.resize(wires);
133  fKernMapR.resize(wires);
134  fKernelR.resize(respRe->GetXaxis()->GetNbins());
135  const TArrayD *edges = respRe->GetXaxis()->GetXbins();
136  for(int i = 0; i < respRe->GetXaxis()->GetNbins(); ++i) {
137  fKernelR[i].resize(bins);
138  for(bin = 0; bin < bins; ++bin) {
139 
140  const TComplex a(respRe->GetBinContent(i+1,bin+1),
141  respIm->GetBinContent(i+1,bin+1));
142  fKernelR[i][bin]=a;
143  }
144  for(; wire < (*edges)[i+1]; ++wire) {
145  fKernMapR[wire]=i;
146  fDecayConstsR[wire]=decayHist->GetBinContent(wire+1);
147  }
148  }
149  respRe = dynamic_cast<TH2D*>(f.Get("sim/RespRe") );
150  respIm = dynamic_cast<TH2D*>(f.Get("sim/RespIm") );
151  decayHist = dynamic_cast<TH1D*>(f.Get("sim/decayHist"));
152  wires = decayHist->GetNbinsX();
153  bins = respRe->GetYaxis()->GetNbins();
154  fDecayConstsS.resize(wires);
155  fKernMapS.resize(wires);
156  fKernelS.resize(respRe->GetXaxis()->GetNbins());
157  const TArrayD *edges1 = respRe->GetXaxis()->GetXbins();
158  wire =0;
159  for(int i = 0; i < respRe->GetXaxis()->GetNbins(); ++i) {
160  fKernelS[i].resize(bins);
161  for(bin = 0; bin < bins; ++bin) {
162  const TComplex b(respRe->GetBinContent(i+1,bin+1),
163  respIm->GetBinContent(i+1,bin+1));
164  fKernelS[i][bin]=b;
165  }
166  for(; wire < (*edges1)[i+1]; ++wire) {
167  fKernMapS[wire]=i;
168  fDecayConstsS[wire]=decayHist->GetBinContent(wire+1);
169  }
170  }
171 
172  f.Close();
173  }
174 
175  //////////////////////////////////////////////////////
177  {
178  }
179 
180  //////////////////////////////////////////////////////
181  void RecoWireICARUS::produce(art::Event& evt)
182  {
183 
184 
185  // get the geometry
186  art::ServiceHandle<geo::Geometry> geom;
187 
188  std::vector<double> decayConsts;
189  std::vector<int> kernMap;
190  std::vector<std::vector<std::complex<double>>> kernel;
191  //Put correct response functions and decay constants in place
192  if(evt.isRealData()) {
193  decayConsts=fDecayConstsR;
194  kernMap=fKernMapR;
195  kernel=fKernelR;
196  }
197  else {
198  decayConsts=fDecayConstsS;
199  kernMap=fKernMapS;
200  kernel=fKernelS;
201  }
202 
203  // make a collection of Wires
204  std::unique_ptr<std::vector<recob::Wire> > wirecol(new std::vector<recob::Wire>);
205  // ... and an association set
206  std::unique_ptr<art::Assns<raw::RawDigit,recob::Wire> > WireDigitAssn
207  (new art::Assns<raw::RawDigit,recob::Wire>);
208 
209  // Read in the digit List object(s).
210  art::Handle< std::vector<raw::RawDigit> > digitVecHandle;
211  evt.getByLabel(fDigitModuleLabel, digitVecHandle);
212 
213  if (!digitVecHandle->size()) return;
214  mf::LogInfo("RecoWireICARUS") << "RecoWireICARUS:: digitVecHandle size is " << digitVecHandle->size();
215 
216  // Use the handle to get a particular (0th) element of collection.
217  art::Ptr<raw::RawDigit> digitVec0(digitVecHandle, 0);
218 
219  unsigned int dataSize = digitVec0->Samples(); //size of raw data vectors
220 
221  int transformSize = dataSize;
222  raw::ChannelID_t channel(raw::InvalidChannelID); // channel number
223  unsigned int bin(0); // time bin loop variable
224 
225  icarus_signal_processing::ICARUSFFT<double> fft(transformSize);
226 
227  double decayConst = 0.; // exponential decay constant of electronics shaping
228  double fitAmplitude = 0.; //This is the seed value for the amplitude in the exponential tail fit
229  std::vector<float> holder; // holds signal data
230  std::vector<float> shortADC; // holds signal data
231  std::vector<float> longADC; // holds signal data
232  std::vector<short> rawadc(transformSize); // vector holding uncompressed adc values
233  std::vector<std::complex<double>> freqHolder(transformSize+1); // temporary frequency data
234  wirecol->reserve(digitVecHandle->size());
235  // loop over all wires
236  for(unsigned int rdIter = 0; rdIter < digitVecHandle->size(); ++rdIter){ // ++ move
237  holder.clear();
238 
239  art::Ptr<raw::RawDigit> digitVec(digitVecHandle, rdIter);
240  channel = digitVec->Channel();
241  //std::cout<<channel<<std::endl;
242  holder.resize(transformSize);
243  shortADC.resize(transformSize);
244  longADC.resize(transformSize);
245 
246  // uncompress the data
247  raw::Uncompress(digitVec->ADCs(), rawadc, digitVec->Compression());
248 
249  for(bin = 0; bin < dataSize; ++bin)
250  holder[bin]=(rawadc[bin]-digitVec->GetPedestal());
251  // fExpEndBins only nonzero for detectors needing exponential tail fitting
252  geo::SigType_t sigtype = geom->SignalType(channel);
253  size_t k;
254  if(sigtype == geo::kInduction)
255  k = 0;
256  else if(sigtype == geo::kCollection)
257  k = 1;
258  else
259  throw cet::exception("RecoWireICARUS") << "Bad signal type = " << sigtype << "\n";
260  if (k >= kernel.size())
261  throw cet::exception("RecoWireICARUS") << "kernel size < " << k << "!\n";
262 
263  if(k==0)
264  {
265  for(bin = 1; bin < dataSize; ++bin)
266  {
267  shortADC[bin]=0;
268  }
269  for(bin = 0; bin < dataSize; ++bin)
270  {
271  if (holder[bin]>9) {
272  for (int ijk=0; ijk<30; ijk++) {
273  if ((int(bin)-ijk)>=0 && (bin+ijk)<4096) {
274  shortADC[bin-ijk]=1;
275  shortADC[bin+ijk]=1;
276  }
277  }
278  }
279  }
280  for(bin = 0; bin < dataSize; ++bin)
281  {
282  if(shortADC[bin]==0)holder[bin]=0;
283  }
284 
285  for(bin = 1; bin < dataSize; ++bin)
286  {
287  holder[bin]+= holder[bin-1];
288  }
289  for(bin = 0; bin < dataSize; ++bin)
290  {
291  if(shortADC[bin]==0)holder[bin]=0;
292  }
293 
294  /*
295  for(bin = 1; bin < dataSize; ++bin)
296  {
297  holder[bin]+= holder[bin-1];
298  longADC[bin]=0;shortADC[bin]=0;
299  }
300 
301  for(bin = 128; bin < dataSize; ++bin)
302  {
303  longADC[bin]=0;
304  for(int ijk=0;ijk<128;ijk++)
305  {
306  longADC[bin]+=(holder[bin-ijk])/128.;
307  }
308  }
309  for(bin = 8; bin < dataSize; ++bin)
310  {
311  shortADC[bin]=0;
312  for(int ijk=0;ijk<8;ijk++)
313  {
314  shortADC[bin]+=(holder[bin-ijk])/8.;
315  }
316  }
317 
318  //for(bin = 0; bin < dataSize; ++bin)
319  //holder[bin]= shortADC[bin]-longADC[bin];
320  for(bin = 0; bin < dataSize; ++bin)holder[bin]=holder[bin]*0.1;
321 
322  //fFFT->Convolute(holder,kernel[k]);
323  */
324  }
325 
326  if(k==2)
327  {
328  if(fExpEndBins && std::abs(decayConsts[channel]) > 0.0){
329 
330  TH1D expTailData("expTailData","Tail data for fit",
331  fExpEndBins,dataSize-fExpEndBins,dataSize);
332  TF1 expFit("expFit","[0]*exp([1]*x)");
333 
334  for(bin = 0; bin < (unsigned int)fExpEndBins; ++bin)
335  expTailData.Fill(dataSize-fExpEndBins+bin,holder[dataSize-fExpEndBins+bin]);
336  decayConst = decayConsts[channel];
337  fitAmplitude = holder[dataSize-fExpEndBins]/exp(decayConst*(dataSize-fExpEndBins));
338  expFit.FixParameter(1,decayConst);
339  expFit.SetParameter(0,fitAmplitude);
340  expTailData.Fit(&expFit,"QWN","",dataSize-fExpEndBins,dataSize);
341  expFit.SetRange(dataSize,transformSize);
342  for(bin = 0; bin < dataSize; ++bin)
343  holder[dataSize+bin]= expFit.Eval(bin+dataSize);
344  }
345  // This is actually deconvolution, by way of convolution with the inverted
346  // kernel. This code assumes the response function has already been
347  // been transformed and inverted. This way a complex multiplication, rather
348  // than a complex division is performed saving 2 multiplications and
349  // 2 divsions
350 
351  // the example below is for MicroBooNE, experiments should
352  // adapt as appropriate
353 
354  // Figure out which kernel to use (0=induction, 1=collection).
356  std::copy(holder.begin(),holder.end(),temp.begin());
357  fft.convolute(temp,kernel[k],0);
358  std::copy(temp.begin(),temp.end(),holder.begin());
359 
360  //This restores the DC component to signal removed by the deconvolution.
361  if(fPostsample) {
362  double average=0.0;
363  for(bin=0; bin < (unsigned int)fPostsample; ++bin)
364  average+=holder[holder.size()-1-bin]/(double)fPostsample;
365  for(bin = 0; bin < holder.size(); ++bin) holder[bin]-=average;
366  }
367  }
368  wirecol->push_back(recob::WireCreator(holder,*digitVec).move());
369  // add an association between the last object in wirecol
370  // (that we just inserted) and digitVec
371  if (!util::CreateAssn(*this, evt, *wirecol, digitVec, *WireDigitAssn)) {
372  throw art::Exception(art::errors::ProductRegistrationFailure)
373  << "Can't associate wire #" << (wirecol->size() - 1)
374  << " with raw digit #" << digitVec.key();
375  } // if failed to add association
376  } // for raw digits
377 
378  if(wirecol->size() == 0)
379  mf::LogWarning("RecoWireICARUS") << "No wires made for this event.";
380 
381  evt.put(std::move(wirecol));
382  evt.put(std::move(WireDigitAssn));
383 
384  return;
385  }
386 
387 } // end namespace recowire
388 
389 
390 namespace recowire{
391 
392  DEFINE_ART_MODULE(RecoWireICARUS)
393 
394 } // end namespace recowire
395 
396 
397 #endif // RecoWireICARUS_H
398 
std::vector< std::vector< std::complex< double > > > fKernelS
void reconfigure(fhicl::ParameterSet const &p)
Helper functions to create a wire.
pdgs p
Definition: selectors.fcl:22
std::vector< double > fDecayConstsR
Definition of basic raw digits.
std::vector< double > fDecayConstsS
for pfile in ack l reconfigure(.*) override"` do echo "checking $
int fExpEndBins
number of end bins to consider for tail fit
Class managing the creation of a new recob::Wire object.
Definition: WireCreator.h:53
RecoWireICARUS(fhicl::ParameterSet const &pset)
constexpr details::BinObj< T > bin(T value)
Returns a wrapper to print the specified data in binary format.
process_name gaushit a
T abs(T value)
std::string fDigitModuleLabel
module that made digits
constexpr ChannelID_t InvalidChannelID
ID of an invalid channel.
Definition: RawTypes.h:32
Signal from induction planes.
Definition: geo_types.h:145
enum geo::_plane_sigtype SigType_t
int fPostsample
number of postsample bins
Collect all the RawData header files together.
std::vector< SigProcPrecision > TimeVec
void produce(art::Event &evt)
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
Encapsulate the construction of a single detector plane.
std::vector< std::vector< std::complex< double > > > fKernelR
T copy(T const &v)
Declaration of basic channel signal object.
BEGIN_PROLOG recowire
TCEvent evt
Definition: DataStructs.cxx:8
pdgs k
Definition: selectors.fcl:22
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:776
art framework interface to geometry description
Signal from collection planes.
Definition: geo_types.h:146