]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRawUtils.cxx
re-structured access to AliEMCALRawUtils in reconstruction
[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.), fGeom(0), 
87     fOption("")
88 {
89   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
90
91   //Get Mapping RCU files from the AliEMCALRecParam                                 
92   const TObjArray* maps = AliEMCALRecParam::GetMappings();
93   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
94
95   for(Int_t i = 0; i < 2; i++) {
96     fMapping[i] = (AliAltroMapping*)maps->At(i);
97   }
98
99
100   fGeom = AliEMCALGeometry::GetInstance();
101   if(!fGeom) {
102     fGeom = AliEMCALGeometry::GetInstance("","");
103     if(!fGeom) AliFatal(Form("Could not get geometry!!"));
104   }
105
106 }
107
108 //____________________________________________________________________________
109 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
110   : TObject(),
111     fHighLowGainFactor(rawU.fHighLowGainFactor), 
112     fGeom(rawU.fGeom), 
113     fOption(rawU.fOption)
114 {
115   //copy ctor
116   fMapping[0] = rawU.fMapping[0];
117   fMapping[1] = rawU.fMapping[1];
118 }
119
120 //____________________________________________________________________________
121 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
122 {
123   //assignment operator
124
125   if(this != &rawU) {
126     fHighLowGainFactor = rawU.fHighLowGainFactor;
127     fGeom = rawU.fGeom;
128     fOption = rawU.fOption;
129     fMapping[0] = rawU.fMapping[0];
130     fMapping[1] = rawU.fMapping[1];
131   }
132
133   return *this;
134
135 }
136
137 //____________________________________________________________________________
138 AliEMCALRawUtils::~AliEMCALRawUtils() {
139
140 }
141
142 //____________________________________________________________________________
143 void AliEMCALRawUtils::Digits2Raw()
144 {
145   // convert digits of the current event to raw data
146   
147   AliRunLoader *rl = AliRunLoader::GetRunLoader();
148   AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
149
150   // get the digits
151   loader->LoadDigits("EMCAL");
152   loader->GetEvent();
153   TClonesArray* digits = loader->Digits() ;
154   
155   if (!digits) {
156     Warning("Digits2Raw", "no digits found !");
157     return;
158   }
159
160   static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
161   AliAltroBuffer* buffers[nDDL];
162   for (Int_t i=0; i < nDDL; i++)
163     buffers[i] = 0;
164
165   Int_t adcValuesLow[fgkTimeBins];
166   Int_t adcValuesHigh[fgkTimeBins];
167
168   // loop over digits (assume ordered digits)
169   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
170     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
171     if (digit->GetAmp() < fgThreshold) 
172       continue;
173
174     //get cell indices
175     Int_t nSM = 0;
176     Int_t nIphi = 0;
177     Int_t nIeta = 0;
178     Int_t iphi = 0;
179     Int_t ieta = 0;
180     Int_t nModule = 0;
181     fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
182     fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
183     
184     //Check which is the RCU of the cell.
185     Int_t iRCU = -111;
186     //RCU0
187     if (0<=iphi&&iphi<8) iRCU=0; // first cable row
188     else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half; 
189     //second cable row
190     //RCU1
191     else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half; 
192     //second cable row
193     else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
194     if (iRCU<0) 
195       Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
196     
197     //Which DDL?
198     Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
199     if (iDDL >= nDDL)
200       Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
201
202     if (buffers[iDDL] == 0) {      
203       // open new file and write dummy header
204       TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
205       buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCU]);
206       buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
207     }
208     
209     // out of time range signal (?)
210     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
211       AliInfo("Signal is out of time range.\n");
212       buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
213       buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() );  // time bin
214       buffers[iDDL]->FillBuffer(3);          // bunch length      
215       buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM);  // trailer
216       // calculate the time response function
217     } else {
218       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ; 
219       if (lowgain) 
220         buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
221       else 
222         buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
223     }
224   }
225   
226   // write headers and close files
227   for (Int_t i=0; i < nDDL; i++) {
228     if (buffers[i]) {
229       buffers[i]->Flush();
230       buffers[i]->WriteDataHeader(kFALSE, kFALSE);
231       delete buffers[i];
232     }
233   }
234
235   loader->UnloadDigits();
236 }
237
238 //____________________________________________________________________________
239 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr)
240 {
241   // convert raw data of the current event to digits                                                                                     
242
243   digitsArr->Clear(); 
244
245   if (!digitsArr) {
246     Error("Raw2Digits", "no digits found !");
247     return;
248   }
249   if (!reader) {
250     Error("Raw2Digits", "no raw reader found !");
251     return;
252   }
253
254   AliCaloRawStream in(reader,"EMCAL",fMapping);
255   // Select EMCAL DDL's;
256   reader->Select("EMCAL");
257
258   TString option = GetOption();
259   if (option.Contains("OldRCUFormat"))
260     in.SetOldRCUFormat(kTRUE); // Needed for testbeam data
261   else
262     in.SetOldRCUFormat(kFALSE);
263
264
265   //Updated fitting routine from 2007 beam test takes into account
266   //possibility of two peaks in data and selects first one for fitting
267   //Also sets some of the starting parameters based on the shape of the
268   //given raw signal being fit
269
270   TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
271   signalF->SetParameters(10.,0.,2.35,2.,5.); //set all defaults once, just to be safe
272   signalF->SetParNames("amp","t0","tau","N","ped");
273   signalF->SetParameter(2,2.35); // tau in units of time bin
274   signalF->SetParLimits(2,2,-1);
275   signalF->SetParameter(3,2); // order
276   signalF->SetParLimits(3,2,-1);
277   
278   Int_t id =  -1;
279   Float_t time = 0. ; 
280   Float_t amp = 0. ; 
281
282   //Graph to hold data we will fit (should be converted to an array
283   //later to speed up processing
284   TGraph * gSig = new TGraph(GetRawFormatTimeBins()); 
285
286   Int_t readOk = 1;
287   Int_t lowGain = 0;
288
289   while (readOk && in.GetModule() < 0) 
290     readOk = in.Next();  // Go to first digit
291
292   Int_t col = 0;
293   Int_t row = 0;
294
295   while (readOk) { 
296     id =  fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
297     lowGain = in.IsLowGain();
298     Int_t maxTime = in.GetTime();  // timebins come in reverse order
299     if (maxTime < 0 || maxTime >= GetRawFormatTimeBins()) {
300       AliWarning(Form("Invalid time bin %d",maxTime));
301       maxTime = GetRawFormatTimeBins();
302     }
303     gSig->Set(maxTime+1);
304     // There is some kind of zero-suppression in the raw data, 
305     // so set up the TGraph in advance
306     for (Int_t i=0; i < maxTime; i++) {
307       gSig->SetPoint(i, i , 0);
308     }
309
310     Int_t iTime = 0;
311     do {
312       if (in.GetTime() >= gSig->GetN()) {
313           AliWarning("Too many time bins");
314           gSig->Set(in.GetTime());
315       }
316       col = in.GetColumn();
317       row = in.GetRow();
318       
319       gSig->SetPoint(in.GetTime(), in.GetTime(), in.GetSignal()) ;
320       if (in.GetTime() > maxTime)
321         maxTime = in.GetTime();
322       iTime++;
323     } while ((readOk = in.Next()) && !in.IsNewHWAddress());
324
325     FitRaw(gSig, signalF, amp, time) ; 
326     
327     if (amp > 0) {
328       AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
329       //cout << "col " << col-40 << " row " << row-8 << " lowGain " << lowGain << " amp " << amp << endl;
330       AddDigit(digitsArr, id, lowGain, (Int_t)amp, time);
331     }
332         
333     // Reset graph
334     for (Int_t index = 0; index < gSig->GetN(); index++) {
335       gSig->SetPoint(index, index, 0) ;  
336     } 
337     // Reset starting parameters for fit function
338     signalF->SetParameters(10.,0.,2.35,2.,5.); //reset all defaults just to be safe
339
340   }; // EMCAL entries loop
341   
342   delete signalF ; 
343   delete gSig;
344   
345   return ; 
346 }
347
348 //____________________________________________________________________________ 
349 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
350   //
351   // Add a new digit. 
352   // This routine checks whether a digit exists already for this tower 
353   // and then decides whether to use the high or low gain info
354   //
355   // Called by Raw2Digits
356   
357   AliEMCALDigit *digit = 0, *tmpdigit = 0;
358   
359   TIter nextdigit(digitsArr);
360   while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
361     if (tmpdigit->GetId() == id)
362       digit = tmpdigit;
363   }
364
365   if (!digit) { // no digit existed for this tower; create one
366     if (lowGain) 
367       amp = Int_t(fHighLowGainFactor * amp); 
368     Int_t idigit = digitsArr->GetEntries();
369     new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;   
370   }
371   else { // a digit already exists, check range 
372          // (use high gain if signal < 800, otherwise low gain)
373     if (lowGain) { // new digit is low gain
374       if (digit->GetAmp() > 800) {  // use if stored digit is out of range
375         digit->SetAmp(Int_t(fHighLowGainFactor * amp));
376         digit->SetTime(time);
377       }
378     }
379     else if (amp < 800) { // new digit is high gain; use if not out of range
380       digit->SetAmp(amp);
381       digit->SetTime(time);
382     }
383   }
384 }
385
386 //____________________________________________________________________________ 
387 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Float_t & amp, Float_t & time)
388 {
389   // Fits the raw signal time distribution; from AliEMCALGetter 
390
391   const Int_t kNoiseThreshold = 5;
392   const Int_t kNPedSamples = 5;
393   amp = time = 0. ; 
394   Double_t ped = 0;
395   Int_t nPed = 0;
396
397   for (Int_t index = 0; index < kNPedSamples; index++) {
398     Double_t ttime, signal;
399     gSig->GetPoint(index, ttime, signal) ; 
400     if (signal > 0) {
401       ped += signal;
402       nPed++;
403     }
404   }
405
406   if (nPed > 0)
407     ped /= nPed;
408   else {
409     AliWarning("Could not determine pedestal");   
410     ped = 10; // put some small value as first guess
411   }
412
413   Int_t max_found = 0;
414   Int_t i_max = 0;
415   Float_t max = -1;
416   Float_t max_fit = gSig->GetN();
417   Float_t min_after_sig = 9999;
418   Int_t tmin_after_sig = gSig->GetN();
419   Int_t n_ped_after_sig = 0;
420
421   for (Int_t i=kNPedSamples; i < gSig->GetN(); i++) {
422     Double_t ttime, signal;
423     gSig->GetPoint(i, ttime, signal) ; 
424     if (!max_found && signal > max) {
425       i_max = i;
426       max = signal;
427     }
428     else if ( max > ped + kNoiseThreshold ) {
429       max_found = 1;
430       min_after_sig = signal;
431       tmin_after_sig = i;
432     }
433     if (max_found) {
434       if ( signal < min_after_sig) {
435         min_after_sig = signal;
436         tmin_after_sig = i;
437       }
438       if (i > tmin_after_sig + 5) {  // Two close peaks; end fit at minimum
439         max_fit = tmin_after_sig;
440         break;
441       }
442       if ( signal < ped + kNoiseThreshold)
443         n_ped_after_sig++;
444       if (n_ped_after_sig >= 5) {  // include 5 pedestal bins after peak
445         max_fit = i;
446         break;
447       }
448     }
449   }
450
451   if ( max - ped > kNoiseThreshold ) { // else its noise 
452     AliDebug(2,Form("Fitting max %d ped %d", max, ped));
453     signalF->SetRange(0,max_fit);
454
455     if(max-ped > 50) 
456       signalF->SetParLimits(2,1,3);
457
458     signalF->SetParameter(4, ped) ; 
459     signalF->SetParameter(1, i_max);
460     signalF->SetParameter(0, max);
461     
462     gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
463     amp = signalF->GetParameter(0); 
464     time = signalF->GetParameter(1)*GetRawFormatTimeBinWidth() - fgTimeTrigger;
465   }
466   return;
467 }
468 //__________________________________________________________________
469 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
470 {
471   // Matches version used in 2007 beam test
472   //
473   // Shape of the electronics raw reponse:
474   // It is a semi-gaussian, 2nd order Gamma function of the general form
475   //
476   // t' = (t - t0 + tau) / tau
477   // F = A * t**N * exp( N * ( 1 - t) )   for t >= 0
478   // F = 0                                for t < 0 
479   //
480   // parameters:
481   // A:   par[0]   // Amplitude = peak value
482   // t0:  par[1]
483   // tau: par[2]
484   // N:   par[3]
485   // ped: par[4]
486   //
487   Double_t signal ;
488   Double_t tau =par[2];
489   Double_t N =par[3];
490   Double_t ped = par[4];
491   Double_t xx = ( x[0] - par[1] + tau ) / tau ;
492
493   if (xx <= 0) 
494     signal = ped ;  
495   else {  
496     signal = ped + par[0] * TMath::Power(xx , N) * TMath::Exp(N * (1 - xx )) ; 
497   }
498   return signal ;  
499 }
500
501 //__________________________________________________________________
502 Bool_t AliEMCALRawUtils::RawSampledResponse(
503 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const 
504 {
505   // for a start time dtime and an amplitude damp given by digit, 
506   // calculates the raw sampled response AliEMCAL::RawResponseFunction
507
508   const Int_t kRawSignalOverflow = 0x3FF ; 
509   const Int_t pedVal = 32;
510   Bool_t lowGain = kFALSE ; 
511
512   // A:   par[0]   // Amplitude = peak value
513   // t0:  par[1]                            
514   // tau: par[2]                            
515   // N:   par[3]                            
516   // ped: par[4]
517
518   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 5);
519   signalF.SetParameter(0, damp) ; 
520   signalF.SetParameter(1, dtime + fgTimeTrigger) ; 
521   signalF.SetParameter(2, fgTau) ; 
522   signalF.SetParameter(3, fgOrder);
523   signalF.SetParameter(4, pedVal);
524
525   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
526     Double_t time = iTime * GetRawFormatTimeBinWidth() ;
527     Double_t signal = signalF.Eval(time) ;     
528     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
529     if ( adcH[iTime] > kRawSignalOverflow ){  // larger than 10 bits 
530       adcH[iTime] = kRawSignalOverflow ;
531       lowGain = kTRUE ; 
532     }
533
534     signal /= fHighLowGainFactor;
535
536     adcL[iTime] =  static_cast<Int_t>(signal + 0.5) ;
537     if ( adcL[iTime] > kRawSignalOverflow)  // larger than 10 bits 
538       adcL[iTime] = kRawSignalOverflow ;
539   }
540   return lowGain ; 
541 }