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