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