]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRawUtils.cxx
Test beam raw data reading
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRawUtils.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line 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. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17 /* History of cvs commits:
18  *
19  * $Log$
20  * Revision 1.5  2007/11/01 01:20:33  mvl
21  * Further improvement of peak finding; more robust fit
22  *
23  * Revision 1.4  2007/10/31 17:15:24  mvl
24  * Fixed bug in raw data unpacking; Added pedestal to signal fit; Added logic to deal with high/low gain
25  *
26  * Revision 1.3  2007/09/27 08:36:46  mvl
27  * More robust setting of fit range in FitRawSignal (P. Hristov)
28  *
29  * Revision 1.2  2007/09/03 20:55:35  jklay
30  * EMCAL e-by-e reconstruction methods from Cvetan
31  *
32  * Revision 1.1  2007/03/17 19:56:38  mvl
33  * Moved signal shape routines from AliEMCAL to separate class AliEMCALRawUtils to streamline raw data reconstruction code.
34  * */
35
36 //*-- Author: Marco van Leeuwen (LBL)
37 #include "AliEMCALRawUtils.h"
38
39 #include "TF1.h"
40 #include "TGraph.h"
41 #include "TSystem.h"
42
43 #include "AliLog.h"
44 #include "AliRunLoader.h"
45 #include "AliCaloAltroMapping.h"
46 #include "AliAltroBuffer.h"
47 #include "AliRawReader.h"
48 #include "AliCaloRawStream.h"
49 #include "AliDAQ.h"
50
51 #include "AliEMCALLoader.h"
52 #include "AliEMCALGeometry.h"
53 #include "AliEMCALDigitizer.h"
54 #include "AliEMCALDigit.h"
55
56
57 ClassImp(AliEMCALRawUtils)
58
59 // Signal shape parameters
60 Int_t    AliEMCALRawUtils::fgOrder         = 2 ;      // Order of gamma function 
61 Double_t AliEMCALRawUtils::fgTimeBinWidth  = 100E-9 ; // each sample is 100 ns
62 Double_t AliEMCALRawUtils::fgTau         = 235E-9 ;   // 235 ns (from CERN testbeam; not very accurate)
63 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ;   // 15 time bins ~ 1.5 musec
64
65 // some digitization constants
66 Int_t    AliEMCALRawUtils::fgThreshold = 1;
67 Int_t    AliEMCALRawUtils::fgDDLPerSuperModule = 2;  // 2 ddls per SuperModule
68
69 AliEMCALRawUtils::AliEMCALRawUtils(): fHighLowGainFactor(0.) {
70   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
71 }
72 //____________________________________________________________________________
73 AliEMCALRawUtils::~AliEMCALRawUtils() {
74 }
75 //____________________________________________________________________________
76 void AliEMCALRawUtils::Digits2Raw()
77 {
78   // convert digits of the current event to raw data
79   
80   AliRunLoader *rl = AliRunLoader::GetRunLoader();
81   AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
82
83   // get the digits
84   loader->LoadDigits("EMCAL");
85   loader->GetEvent();
86   TClonesArray* digits = loader->Digits() ;
87   
88   if (!digits) {
89     Warning("Digits2Raw", "no digits found !");
90     return;
91   }
92     
93   // get the geometry
94   AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
95   if (!geom) {
96     AliError(Form("No geometry found !"));
97     return;
98   }
99   
100   static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
101   AliAltroBuffer* buffers[nDDL];
102   for (Int_t i=0; i < nDDL; i++)
103     buffers[i] = 0;
104
105   Int_t adcValuesLow[fgkTimeBins];
106   Int_t adcValuesHigh[fgkTimeBins];
107
108   //Load Mapping RCU files once
109   TString path = gSystem->Getenv("ALICE_ROOT");
110   path += "/EMCAL/mapping/RCU";
111   TString path0 = path+"0.data";//This file will change in future
112   TString path1 = path+"1.data";//This file will change in future
113   AliAltroMapping * mapping[2] ; // For the moment only 2
114   mapping[0] = new AliCaloAltroMapping(path0.Data());
115   mapping[1] = new AliCaloAltroMapping(path1.Data());
116
117   // loop over digits (assume ordered digits)
118   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
119     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
120     if (digit->GetAmp() < fgThreshold) 
121       continue;
122
123     //get cell indices
124     Int_t nSM = 0;
125     Int_t nIphi = 0;
126     Int_t nIeta = 0;
127     Int_t iphi = 0;
128     Int_t ieta = 0;
129     Int_t nModule = 0;
130     geom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
131     geom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
132     
133     //Check which is the RCU of the cell.
134     Int_t iRCU = -111;
135     //RCU0
136     if (0<=iphi&&iphi<8) iRCU=0; // first cable row
137     else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half; 
138     //second cable row
139     //RCU1
140     else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half; 
141     //second cable row
142     else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
143     
144     //Which DDL?
145     Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
146     if (iDDL >= nDDL)
147       Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
148
149     if (buffers[iDDL] == 0) {      
150       // open new file and write dummy header
151       TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
152       buffers[iDDL] = new AliAltroBuffer(fileName.Data(),mapping[iRCU]);
153       buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
154     }
155     
156     // out of time range signal (?)
157     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
158       AliInfo("Signal is out of time range.\n");
159       buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
160       buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() );  // time bin
161       buffers[iDDL]->FillBuffer(3);          // bunch length      
162       buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM);  // trailer
163       // calculate the time response function
164     } else {
165       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ; 
166       if (lowgain) 
167         buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
168       else 
169         buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
170     }
171   }
172   
173   // write headers and close files
174   for (Int_t i=0; i < nDDL; i++) {
175     if (buffers[i]) {
176       buffers[i]->Flush();
177       buffers[i]->WriteDataHeader(kFALSE, kFALSE);
178       delete buffers[i];
179     }
180   }
181   mapping[0]->Delete();
182   mapping[1]->Delete();
183   loader->UnloadDigits();
184 }
185
186 //____________________________________________________________________________
187 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr)
188 {
189   // convert raw data of the current event to digits
190   AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
191   if (!geom) {
192     AliError(Form("No geometry found !"));
193     return;
194   }
195
196   digitsArr->Clear(); 
197
198   if (!digitsArr) {
199     Error("Raw2Digits", "no digits found !");
200     return;
201   }
202   if (!reader) {
203     Error("Raw2Digits", "no raw reader found !");
204     return;
205   }
206
207   // Use AliAltroRawStream to read the ALTRO format.  No need to
208   // reinvent the wheel :-) 
209   AliCaloRawStream in(reader,"EMCAL");
210   // Select EMCAL DDL's;
211   reader->Select("EMCAL");
212   //in.SetOldRCUFormat(kTRUE); // Needed for testbeam data
213   
214   cout << "Stream set up" << endl;
215
216   // reading is from previously existing AliEMCALGetter.cxx
217   // ReadRaw method
218   TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
219   
220   Int_t id =  -1;
221   Float_t time = 0. ; 
222   Float_t amp = 0. ; 
223
224   TGraph * gSig = new TGraph(GetRawFormatTimeBins()) ; 
225
226   Int_t readOk = 1;
227   Int_t lowGain = 0;
228
229   while (readOk && in.GetModule() < 0) 
230     readOk = in.Next();  // Go to first digit
231
232   while (readOk) { 
233     id =  geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
234     lowGain = in.IsLowGain();
235     Int_t maxTime = in.GetTime();  // timebins come in reverse order
236     if (maxTime < 0 || maxTime >= GetRawFormatTimeBins()) {
237       AliWarning(Form("Invalid time bin %d",maxTime));
238       maxTime = GetRawFormatTimeBins();
239     }
240     gSig->Set(maxTime+1);
241     // There is some kind of zero-suppression in the raw data, 
242     // so set up the TGraph in advance
243     for (Int_t i=0; i < maxTime; i++) {
244       gSig->SetPoint(i, i * GetRawFormatTimeBinWidth(), 0);
245     }
246
247     Int_t iTime = 0;
248     do {
249       if (in.GetTime() >= gSig->GetN()) {
250           AliWarning("Too many time bins");
251           gSig->Set(in.GetTime());
252       }
253       gSig->SetPoint(in.GetTime(), 
254                    in.GetTime() * GetRawFormatTimeBinWidth(), 
255                    in.GetSignal()) ;
256       if (in.GetTime() > maxTime)
257         maxTime = in.GetTime();
258       iTime++;
259     } while ((readOk = in.Next()) && !in.IsNewHWAddress());
260     signalF->SetRange(0,(Float_t)maxTime*GetRawFormatTimeBinWidth());
261
262     FitRaw(gSig, signalF, amp, time) ; 
263     
264     if (amp > 0) {
265       AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
266       AddDigit(digitsArr, id, lowGain, (Int_t)amp, time);
267     }
268         
269     // Reset graph
270     for (Int_t index = 0; index < gSig->GetN(); index++) {
271       gSig->SetPoint(index, index * GetRawFormatTimeBinWidth(), 0) ;  
272     } 
273   }; // EMCAL entries loop
274   
275   delete signalF ; 
276   delete gSig;
277   
278   return ; 
279 }
280
281 //____________________________________________________________________________ 
282 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
283   //
284   // Add a new digit. 
285   // This routine checks whether a digit exists already for this tower 
286   // and then decides whether to use the high or low gain info
287   //
288   // Called by Raw2Digits
289   
290   AliEMCALDigit *digit = 0, *tmpdigit = 0;
291   
292   TIter nextdigit(digitsArr);
293   while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
294     if (tmpdigit->GetId() == id)
295       digit = tmpdigit;
296   }
297
298   if (!digit) { // no digit existed for this tower; create one
299     if (lowGain) 
300       amp = Int_t(fHighLowGainFactor * amp); 
301     Int_t idigit = digitsArr->GetEntries();
302     new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;   
303   }
304   else { // a digit already exists, check range 
305          // (use high gain if signal < 800, otherwise low gain)
306     if (lowGain) { // new digit is low gain
307       if (digit->GetAmp() > 800) {  // use if stored digit is out of range
308         digit->SetAmp(Int_t(fHighLowGainFactor * amp));
309         digit->SetTime(time);
310       }
311     }
312     else if (amp < 800) { // new digit is high gain; use if not out of range
313       digit->SetAmp(amp);
314       digit->SetTime(time);
315     }
316   }
317 }
318
319 //____________________________________________________________________________ 
320 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Float_t & amp, Float_t & time)
321 {
322   // Fits the raw signal time distribution; from AliEMCALGetter 
323
324   const Int_t kNoiseThreshold = 5;
325   const Int_t kNPedSamples = 10;
326   amp = time = 0. ; 
327   Double_t ped = 0;
328   Int_t nPed = 0;
329
330   for (Int_t index = 0; index < kNPedSamples; index++) {
331     Double_t ttime, signal;
332     gSig->GetPoint(index, ttime, signal) ; 
333     if (signal > 0) {
334       ped += signal;
335       nPed++;
336     }
337   }
338
339   if (nPed > 0)
340     ped /= nPed;
341   else {
342     AliWarning("Could determine pedestal");       
343     ped = 10; // put some small value as first guess
344   }
345
346   Int_t max_found = 0;
347   Int_t i_max = 0;
348   Float_t max = -1;
349   Float_t tmax = 0;
350   Float_t max_fit = gSig->GetN()*GetRawFormatTimeBinWidth();
351   Float_t min_after_sig = 9999;
352   Int_t imin_after_sig = gSig->GetN();
353   Float_t tmin_after_sig = gSig->GetN()*GetRawFormatTimeBinWidth();
354   Int_t n_ped_after_sig = 0;
355
356   for (Int_t i=kNPedSamples; i < gSig->GetN(); i++) {
357     Double_t ttime, signal;
358     gSig->GetPoint(i, ttime, signal) ; 
359     if (!max_found && signal > max) {
360       i_max = i;
361       tmax = ttime;
362       max = signal;
363     }
364     else if ( max > ped + kNoiseThreshold ) {
365       max_found = 1;
366       min_after_sig = signal;
367       imin_after_sig = i;
368       tmin_after_sig = ttime;
369     }
370     if (max_found) {
371       if ( signal < min_after_sig) {
372         min_after_sig = signal;
373         imin_after_sig = i;
374         tmin_after_sig = ttime;
375       }
376       if (i > tmin_after_sig + 5) {  // Two close peaks; end fit at minimum
377         max_fit = tmin_after_sig;
378         break;
379       }
380       if ( signal < ped + kNoiseThreshold)
381         n_ped_after_sig++;
382       if (n_ped_after_sig >= 5) {  // include 5 pedestal bins after peak
383         max_fit = ttime;
384         break;
385       }
386     }
387   }
388
389   if ( max - ped > kNoiseThreshold ) { // else its noise 
390     AliDebug(2,Form("Fitting max %d ped %d", max, ped));
391     signalF->SetParameter(0, ped) ; 
392     signalF->SetParameter(1, max - ped) ; 
393     signalF->SetParameter(2, tmax) ; 
394     signalF->SetParLimits(2, 0, max_fit) ; 
395     gSig->Fit(signalF, "QRON", "", 0., max_fit); //, "QRON") ; 
396     amp = signalF->GetParameter(1); 
397     time = signalF->GetParameter(2) - fgTimeTrigger;
398   }
399   return;
400 }
401 //__________________________________________________________________
402 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
403 {
404   // Shape of the electronics raw reponse:
405   // It is a semi-gaussian, 2nd order Gamma function of the general form
406   //
407   // t' = (t - t0 + tau) / tau
408   // F = A * t**N * exp( N * ( 1 - t) )   for t >= 0
409   // F = 0                                for t < 0 
410   //
411   // parameters:
412   // ped: par[0]
413   // A:   par[1]   // Amplitude = peak value
414   // t0:  par[2]
415   // tau: fgTau
416   // N:   fgOrder
417   //
418   Double_t signal ;
419   Double_t xx = ( x[0] - par[2] + fgTau ) / fgTau ; 
420
421   if (xx <= 0) 
422     signal = par[0] ;  
423   else {  
424     signal = par[0] + par[1] * TMath::Power(xx , fgOrder) * TMath::Exp(fgOrder * (1 - xx )) ; 
425
426   }
427   return signal ;  
428 }
429
430 //__________________________________________________________________
431 Bool_t AliEMCALRawUtils::RawSampledResponse(
432 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const 
433 {
434   // for a start time dtime and an amplitude damp given by digit, 
435   // calculates the raw sampled response AliEMCAL::RawResponseFunction
436
437   const Int_t kRawSignalOverflow = 0x3FF ; 
438   const Int_t pedVal = 32;
439   Bool_t lowGain = kFALSE ; 
440
441   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
442   signalF.SetParameter(0, pedVal) ; 
443   signalF.SetParameter(1, damp) ; 
444   signalF.SetParameter(2, dtime + fgTimeTrigger) ; 
445
446   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
447     Double_t time = iTime * GetRawFormatTimeBinWidth() ;
448     Double_t signal = signalF.Eval(time) ;     
449     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
450     if ( adcH[iTime] > kRawSignalOverflow ){  // larger than 10 bits 
451       adcH[iTime] = kRawSignalOverflow ;
452       lowGain = kTRUE ; 
453     }
454
455     signal /= fHighLowGainFactor;
456
457     adcL[iTime] =  static_cast<Int_t>(signal + 0.5) ;
458     if ( adcL[iTime] > kRawSignalOverflow)  // larger than 10 bits 
459       adcL[iTime] = kRawSignalOverflow ;
460   }
461   return lowGain ; 
462 }