]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/icepack/IceMakeHits.cxx
06-jun-2007 NvE Explicit tests on null pointers for returned TObjArray* from e.g...
[u/mrichter/AliRoot.git] / RALICE / icepack / IceMakeHits.cxx
1 /*******************************************************************************
2  * Copyright(c) 2003, IceCube Experiment at the South Pole. All rights reserved.
3  *
4  * Author: The IceCube RALICE-based Offline Project.
5  * Contributors are mentioned in the code where appropriate.
6  *
7  * Permission to use, copy, modify and distribute this software and its
8  * documentation strictly for non-commercial purposes is hereby granted
9  * without fee, provided that the above copyright notice appears in all
10  * copies and that both the copyright notice and this permission notice
11  * appear in the supporting documentation.
12  * The authors make no claims about the suitability of this software for
13  * any purpose. It is provided "as is" without express or implied warranty.
14  *******************************************************************************/
15
16 // $Id$
17
18 ///////////////////////////////////////////////////////////////////////////
19 // Class IceMakeHits
20 // TTask derived class to perform hit extraction from waveforms.
21 //
22 // The code in this processor is based on the algorithms as developed by
23 // Nick van Eijndhoven and Garmt de Vries-Uiterweerd (Utrecht University, The Netherlands).
24 //
25 // Procedure applied for Amanda TWR data :
26 // ---------------------------------------
27 //
28 // 1) The waveform is fed to a TSpectrum object, and the peak locations 
29 //    are determined with the TSpectrum::Search() function.
30 //
31 // 2) The waveform is divided into regions corresponding to the peaks found by 
32 //    TSpectrum. The region boundary between two peaks is at the location of 
33 //    the minimum between the two peaks. 
34 //
35 // 3) For each region the "effective baseline" (used in the
36 //    evaluation of the leading edge value) is determined as :
37 //    effective baseline = fBasefracXXX * value at lower region boundary.
38 //    This takes into account the effect from the previous pulse.
39 //    For the first pulse, the effective baseline is equal to the overall 
40 //    baseline.
41 //
42 // 4) For each region, the point of steepest rise between the lower region
43 //    boundary and the peak location is determined. The tangent at this point
44 //    is extrapolated to the effective baseline. The point of intersection yields the
45 //    leading edge.
46 //
47 // 5) For each region the range of charge integration is determined as :
48 //    - Start of integration at the lower region boundary or at the leading edge,
49 //      whichever comes last;
50 //    - End of integration at the upper region boundary or at the point where the
51 //      signal drops below the overall baseline, whichever comes first.
52 //
53 // 6) For each region the integrated charge is determined as :
54 //    Sum over bins in integration range of (value in bin - overall baseline).
55 //
56 // 7) For each pulse the quality is evaluated by requiring that :
57 //    peak location - lower region boundary > lower region boundary - leading edge.
58 //    For a too shallow steepest rise, the leading edge value is unreliable, in 
59 //    which case the pulse is merged with the previous pulse. 
60 // 
61 // 8) Each pulse is checked for saturation and discarded if necessary.
62 //
63 // 9) TSpectrum needs a minimum number of bins for its Search function, otherwise
64 //    the clipping window is too large, which causes an error. If a waveform does
65 //    not contain enough bins, the following alternative approach is used : 
66 //    - A loop over all bins is performed.
67 //    - As soon as the signal exceeds a given threshold, a pulse is started.
68 //    - The pulse ends when the signal drops below the threshold again.
69 //    - While looping, the charge is integrated for each pulse.
70 //
71 // The defaults of the various parameters can be changed by the corresponding
72 // Set memberfunctions.
73 //
74 // Information about the actual parameter settings can be found in the event
75 // structure itself via the device named "IceMakeHits".
76 //
77 //--- Author: Nick van Eijndhoven and Garmt de Vries-Uiterweerd 15-jan-2007 Utrecht University
78 //- Modified: GdV $Date$ Utrecht University
79 ///////////////////////////////////////////////////////////////////////////
80  
81 #include "IceMakeHits.h"
82 #include "Riostream.h"
83
84 ClassImp(IceMakeHits) // Class implementation to enable ROOT I/O
85
86 IceMakeHits::IceMakeHits(const char* name,const char* title) : TTask(name,title)
87 {
88 // Default constructor.
89  fEvt=0;
90  fDaq=0;
91  fBasefracA=0.5;
92  fSigmaA=1.5;
93  fMaxPeaksA=10;
94  fMinPulseHeightA=50;
95  fThresholdA=0.2;
96 }
97 ///////////////////////////////////////////////////////////////////////////
98 IceMakeHits::~IceMakeHits()
99 {
100 // Default destructor.
101 }
102 ///////////////////////////////////////////////////////////////////////////
103 void IceMakeHits::SetBasefracA(Float_t val)
104 {
105 // Set baseline fractional update for Amanda TWR extraction.
106 // The default as set in the constructor of this class is 0.5.
107  fBasefracA=val;
108 }
109 ///////////////////////////////////////////////////////////////////////////
110 void IceMakeHits::SetSigmaA(Float_t val)
111 {
112 // Set clipping window width for Amanda TWR extraction.
113 // The default as set in the constructor of this class is 1.5.
114  fSigmaA=val;
115 }
116 ///////////////////////////////////////////////////////////////////////////
117 void IceMakeHits::SetMaxPeaksA(Int_t val)
118 {
119 // Set maximum number of peaks in a waveform for Amanda TWR extraction.
120 // The default as set in the constructor of this class is 10.
121  fMaxPeaksA=val;
122 }
123 ///////////////////////////////////////////////////////////////////////////
124 void IceMakeHits::SetMinPulseHeightA(Float_t val)
125 {
126 // Set minimum required pulse height for Amanda TWR extraction.
127 // This is used only for narrow pulses that cannot be handled with TSpectrum.
128 // The default as set in the constructor of this class is 50.
129  fMinPulseHeightA=val;
130 }
131 ///////////////////////////////////////////////////////////////////////////
132 void IceMakeHits::SetThresholdA(Float_t val)
133 {
134 // Set threshold for use in analysis of narrow pulses for Amanda TWR extraction.
135 // A peak is assumed to start when the signal rises above threshold*maxval,
136 // where maxval is the maximum value found in the waveform.
137 // The default as set in the constructor of this class is 0.2.
138  fThresholdA=val;
139 }
140 ///////////////////////////////////////////////////////////////////////////
141 void IceMakeHits::Exec(Option_t* opt)
142 {
143 // Implementation of the hit cleaning procedures.
144
145  TString name=opt;
146  AliJob* parent=(AliJob*)(gROOT->GetListOfTasks()->FindObject(name.Data()));
147
148  if (!parent) return;
149
150  fEvt=(IceEvent*)parent->GetObject("IceEvent");
151  if (!fEvt) return;
152
153  fDaq=(AliDevice*)fEvt->GetDevice("Daq");
154  if (!fDaq) return;
155
156  // Storage of the used parameters in the IceMakeHits device
157  AliSignal params;
158  params.SetNameTitle("IceMakeHits","IceMakeHits processor parameters");
159  params.SetSlotName("BasefracA",1);
160  params.SetSlotName("SigmaA",2);
161  params.SetSlotName("MaxPeaksA",3);
162  params.SetSignal(fBasefracA,1);
163  params.SetSignal(fSigmaA,2);
164  params.SetSignal(fMaxPeaksA,3);
165
166  fEvt->AddDevice(params);
167
168  Amanda();
169  InIce();
170  IceTop();
171
172 }
173 ///////////////////////////////////////////////////////////////////////////
174 void IceMakeHits::Amanda()
175 {
176 // Hit extraction from the Amanda TWR data.
177
178  // Arrays for storing info
179  Float_t* baseline=new Float_t[fMaxPeaksA];
180  Int_t* lowend=new Int_t[fMaxPeaksA];
181  Int_t* upend=new Int_t[fMaxPeaksA];
182  Int_t* startcharge=new Int_t[fMaxPeaksA];
183  Int_t* stopcharge=new Int_t[fMaxPeaksA];
184  Int_t* status=new Int_t[fMaxPeaksA]; // 0=OK, 1=rejected, 2=saturation
185  Float_t* leadingedge=new Float_t[fMaxPeaksA];
186  Float_t* charge=new Float_t[fMaxPeaksA];
187  Float_t* tot=new Float_t[fMaxPeaksA];
188
189  // Some objects and variables we will need
190  TH1F foo, diff;
191  TSpectrum spec(fMaxPeaksA);
192  Int_t nrIterations=(Int_t)(7*fSigmaA+0.5); // Number of iterations used in TSpectrum::SearchHighRes()
193  Int_t npeaks=0, ibin=0, lookforsteepestuntilbin=0, steep=0;
194  Float_t maxval=0, rise=0, rc=0, yyy=0;
195  Bool_t pulsegoingon=false;
196  Int_t* index=new Int_t[fMaxPeaksA];
197
198  // Update the DAQ device data for the nature of the hit info
199  fDaq->AddNamedSlot("TWR");
200  fDaq->SetSignal(1,"TWR");
201  Int_t idx=fDaq->GetSlotIndex("Muon");
202  if (idx) fDaq->SetSignal(0,idx);
203
204  // All Amanda OMs with a signal
205  TObjArray* aoms=fEvt->GetDevices("IceAOM");
206  if (!aoms) return;
207
208  // OM, waveform and hit
209  IceAOM* omx=0;
210  TH1F* wf=0;
211  AliSignal hit;
212  hit.SetSlotName("ADC",1);
213  hit.SetSlotName("LE",2);
214  hit.SetSlotName("TOT",3);
215
216  // Loop over all fired OMs and extract the hit info
217  for (Int_t iom=0; iom<aoms->GetEntries(); iom++)
218  {
219   omx=(IceAOM*)aoms->At(iom);
220   if (!omx) continue;
221   // Remove all existing hits of this OM 
222   omx->RemoveHits();
223   // Should we skip OMs that we know from the dbase to have problems ?
224 ////  if (omx->GetDeadValue("ADC") || omx->GetDeadValue("LE") || omx->GetDeadValue("TOT")) continue;
225
226   // Investigate all waveforms for this OM
227   for (Int_t iwf=1; iwf<=omx->GetNwaveforms(); iwf++)
228   {
229    wf=omx->GetWaveform(iwf);
230    if (!wf) continue; 
231    maxval=wf->GetMaximum();
232    // Check if clipping window is not too large
233    if(wf->GetNbinsX() > 2*nrIterations+1)
234    {
235     // Find peaks with TSpectrum
236     npeaks=spec.Search(wf,fSigmaA,"goff");
237     // Discard waveform if no or too many peaks found
238     if(npeaks<1 || npeaks>fMaxPeaksA) continue;
239
240     // Get differential of WF
241     diff=*wf;
242     for(ibin=2;ibin<diff.GetNbinsX();ibin++)
243     {
244      diff.SetBinContent(ibin,wf->GetBinContent(ibin)-wf->GetBinContent(ibin-1));
245     }
246     diff.SetBinContent(1,0);
247     // Set baseline and lower end for first peak,
248     baseline[0]=0;
249     lowend[0]=1;
250
251     // Sort peaks in time
252     TMath::Sort(npeaks,spec.GetPositionX(),index,false);
253     // For each of the peaks,
254     for(Int_t ipeak=0; ipeak<npeaks; ipeak++)
255     {
256      // Find baseline and region around peak
257      foo=*wf;
258      // (Second and later peaks: lower edge = upper edge previous peak,
259      // baseline is average of previous baseline and minimum value between two
260      // peaks)
261      if(ipeak>0)
262      {
263       lowend[ipeak]=upend[ipeak-1]+1;
264       baseline[ipeak]=fBasefracA*foo.GetBinContent(lowend[ipeak]);
265      }
266      // (Upper edge range is minimum between this and next peak)
267      if(ipeak<npeaks-1)
268      {
269       foo.SetAxisRange(spec.GetPositionX()[index[ipeak]],spec.GetPositionX()[index[ipeak+1]]);
270       upend[ipeak]=foo.GetMinimumBin()-1;
271      }
272      // (Last peak: upper edge is end of histo)
273      else
274      {
275       upend[ipeak]=wf->GetNbinsX();
276      }
277      // Find steepest rise
278      lookforsteepestuntilbin=wf->FindBin(spec.GetPositionX()[index[ipeak]]);
279      foo=diff;
280      // Look for steepest rise between lower edge and peak position
281      foo.SetAxisRange(wf->GetBinCenter(lowend[ipeak]),wf->GetBinCenter(lookforsteepestuntilbin));
282      steep=foo.GetMaximumBin();
283      rise=foo.GetBinContent(steep);
284
285      // Extrapolate tangent to find leading edge
286      yyy=wf->GetBinContent(steep)-baseline[ipeak];
287      rc=rise/foo.GetBinWidth(steep);
288      if(rc>0) leadingedge[ipeak]=wf->GetBinCenter(steep)-yyy/rc; else leadingedge[ipeak]=0;
289
290      // Determine peak status
291      status[ipeak]=0;
292      // Check for saturation
293      if(rc<0.1 && wf->GetBinContent(wf->FindBin(spec.GetPositionX()[index[ipeak]])) == maxval)
294      {
295       status[ipeak]=2;
296      }
297      // Check quality: LE should not be too far below lower edge
298      // Otherwise, ignore this peak and set baseline back to what it was
299      else if(wf->GetBinLowEdge(lowend[ipeak]) - leadingedge[ipeak] > spec.GetPositionX()[index[ipeak]] - wf->GetBinLowEdge(lowend[ipeak]))
300      {
301       status[ipeak]=1;
302       if(ipeak>0) baseline[ipeak]=baseline[ipeak-1];
303      }
304  
305      // Start charge integration at LE, or at lower edge of range
306      startcharge[ipeak]=wf->FindBin(leadingedge[ipeak]);
307      if(lowend[ipeak]>startcharge[ipeak]) startcharge[ipeak]=lowend[ipeak];
308  
309      // Integrate charge until pulse drop below baseline, or else until edge of range
310      stopcharge[ipeak]=upend[ipeak];
311      for(ibin=wf->FindBin(spec.GetPositionX()[index[ipeak]]); ibin<=upend[ipeak]; ibin++)
312      {
313       if(wf->GetBinContent(ibin)<0)
314       {
315        stopcharge[ipeak]=ibin-1;
316        break;
317       }
318      }
319
320      // Determine time over threshold
321      tot[ipeak]=wf->GetBinLowEdge(stopcharge[ipeak]+1)-wf->GetBinLowEdge(startcharge[ipeak]);
322  
323      // Determine charge
324      charge[ipeak]=0;
325      for(ibin=startcharge[ipeak]; ibin<=stopcharge[ipeak]; ibin++)
326      {
327       charge[ipeak]+=wf->GetBinContent(ibin);
328      }
329
330     } // end loop over peaks
331  
332     // Check all peaks, from latest to earliest
333     for(int ipeak=npeaks-1; ipeak>=0; ipeak--)
334     {
335
336      // If this peak was rejected, add charge and TOT to previous peak (if there is one)
337      if(status[ipeak]==1 && ipeak>0)
338      {
339       charge[ipeak-1]+=charge[ipeak];
340       charge[ipeak]=0;
341       tot[ipeak-1]+=tot[ipeak];
342       tot[ipeak]=0;
343      }
344
345      // If this peak is OK, add hit info
346      if(status[ipeak]==0)
347      {
348       hit.Reset();
349       hit.SetSignal(charge[ipeak],"ADC");
350       hit.SetSignal(leadingedge[ipeak],"LE");
351       hit.SetSignal(tot[ipeak],"TOT");
352       omx->AddHit(hit);
353      }
354
355     } // end loop over peaks
356    // If number of bins too small, use different method
357    }
358    else
359    {
360     // If maximum value high enough to suspect presence of peak,
361     if(maxval>fMinPulseHeightA)
362     {
363      // Loop over bins
364      pulsegoingon=false;
365      npeaks=0;
366      for(ibin=1; ibin<=wf->GetNbinsX(); ibin++)
367      {
368       // If bin content above threshold, start pulse
369       if(wf->GetBinContent(ibin)>fThresholdA*maxval){
370        if(!pulsegoingon)
371        {
372         // Pulse starts here
373         pulsegoingon=true;
374         leadingedge[npeaks]=wf->GetBinLowEdge(ibin);
375         charge[npeaks]=wf->GetBinContent(ibin);
376        }
377        else
378        {
379         // Pulse continues       
380         charge[npeaks]+=wf->GetBinContent(ibin);
381        }
382       }
383       else
384       {
385        if(pulsegoingon)
386        {
387         // Pulse ends here
388         tot[npeaks]=wf->GetBinLowEdge(ibin)-leadingedge[npeaks];
389
390         // Store pulse information
391         hit.Reset();
392         hit.SetSignal(charge[npeaks],"ADC");
393         hit.SetSignal(leadingedge[npeaks],"LE");
394         hit.SetSignal(tot[npeaks],"TOT");
395         omx->AddHit(hit);
396
397         // Get ready for next pulse
398         pulsegoingon=false;
399         npeaks++;
400        }
401       }
402
403      } // End of loop over bins
404     } 
405
406    } // End of alternative method for narrow pulses
407
408   } // End of WF loop
409  } // End of OM loop
410
411  // Clean up
412  delete[] baseline;
413  delete[] lowend;
414  delete[] upend;
415  delete[] startcharge;
416  delete[] stopcharge;
417  delete[] status;
418  delete[] leadingedge;
419  delete[] charge;
420  delete[] tot;
421  delete[] index;
422 }
423 ///////////////////////////////////////////////////////////////////////////
424 void IceMakeHits::InIce()
425 {
426 // Hit extraction from IceCube InIce ATWD data.
427 }
428 ///////////////////////////////////////////////////////////////////////////
429 void IceMakeHits::IceTop()
430 {
431 // Hit extraction from IceTop ATWD data.
432 }
433 ///////////////////////////////////////////////////////////////////////////