1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 /* History of cvs commits:
20 * Revision 1.9 2007/12/06 02:19:51 jklay
21 * incorporated fitting procedure from testbeam analysis into AliRoot
23 * Revision 1.8 2007/12/05 02:30:51 jklay
24 * 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
26 * Revision 1.7 2007/11/14 15:51:46 gustavo
27 * Take out few unnecessary prints
29 * Revision 1.6 2007/11/01 01:23:51 mvl
30 * Removed call to SetOldRCUFormat, which is only needed for testbeam data
32 * Revision 1.5 2007/11/01 01:20:33 mvl
33 * Further improvement of peak finding; more robust fit
35 * Revision 1.4 2007/10/31 17:15:24 mvl
36 * Fixed bug in raw data unpacking; Added pedestal to signal fit; Added logic to deal with high/low gain
38 * Revision 1.3 2007/09/27 08:36:46 mvl
39 * More robust setting of fit range in FitRawSignal (P. Hristov)
41 * Revision 1.2 2007/09/03 20:55:35 jklay
42 * EMCAL e-by-e reconstruction methods from Cvetan
44 * Revision 1.1 2007/03/17 19:56:38 mvl
45 * Moved signal shape routines from AliEMCAL to separate class AliEMCALRawUtils to streamline raw data reconstruction code.
48 //*-- Author: Marco van Leeuwen (LBL)
49 #include "AliEMCALRawUtils.h"
56 #include "AliRunLoader.h"
57 #include "AliCaloAltroMapping.h"
58 #include "AliAltroBuffer.h"
59 #include "AliRawReader.h"
60 #include "AliCaloRawStream.h"
63 #include "AliEMCALRecParam.h"
64 #include "AliEMCALLoader.h"
65 #include "AliEMCALGeometry.h"
66 #include "AliEMCALDigitizer.h"
67 #include "AliEMCALDigit.h"
70 ClassImp(AliEMCALRawUtils)
72 // Signal shape parameters
73 Int_t AliEMCALRawUtils::fgOrder = 2 ; // Order of gamma function
74 Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ; // each sample is 100 ns
75 Double_t AliEMCALRawUtils::fgTau = 235E-9 ; // 235 ns (from CERN testbeam; not very accurate)
76 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
78 // some digitization constants
79 Int_t AliEMCALRawUtils::fgThreshold = 1;
80 Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
82 AliEMCALRawUtils::AliEMCALRawUtils()
83 : fHighLowGainFactor(0.), fOption("")
85 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
87 //____________________________________________________________________________
88 AliEMCALRawUtils::~AliEMCALRawUtils() {
90 //____________________________________________________________________________
91 void AliEMCALRawUtils::Digits2Raw(AliAltroMapping **mapping)
93 // convert digits of the current event to raw data
95 AliRunLoader *rl = AliRunLoader::GetRunLoader();
96 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
99 loader->LoadDigits("EMCAL");
101 TClonesArray* digits = loader->Digits() ;
104 Warning("Digits2Raw", "no digits found !");
109 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
111 AliError(Form("No geometry found !"));
115 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
116 AliAltroBuffer* buffers[nDDL];
117 for (Int_t i=0; i < nDDL; i++)
120 Int_t adcValuesLow[fgkTimeBins];
121 Int_t adcValuesHigh[fgkTimeBins];
123 // loop over digits (assume ordered digits)
124 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
125 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
126 if (digit->GetAmp() < fgThreshold)
136 geom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
137 geom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
139 //Check which is the RCU of the cell.
142 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
143 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
146 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
148 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
150 Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
153 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
155 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
157 if (buffers[iDDL] == 0) {
158 // open new file and write dummy header
159 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
160 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),mapping[iRCU]);
161 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
164 // out of time range signal (?)
165 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
166 AliInfo("Signal is out of time range.\n");
167 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
168 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
169 buffers[iDDL]->FillBuffer(3); // bunch length
170 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
171 // calculate the time response function
173 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ;
175 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
177 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
181 // write headers and close files
182 for (Int_t i=0; i < nDDL; i++) {
185 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
189 //PH mapping[0]->Delete();
190 //PH mapping[1]->Delete();
191 loader->UnloadDigits();
194 //____________________________________________________________________________
195 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr,
196 AliAltroMapping **mapping)
198 // convert raw data of the current event to digits
199 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
201 AliError(Form("No geometry found !"));
208 Error("Raw2Digits", "no digits found !");
212 Error("Raw2Digits", "no raw reader found !");
216 AliCaloRawStream in(reader,"EMCAL",mapping);
217 // Select EMCAL DDL's;
218 reader->Select("EMCAL");
220 TString option = GetOption();
221 if (option.Contains("OldRCUFormat"))
222 in.SetOldRCUFormat(kTRUE); // Needed for testbeam data
224 in.SetOldRCUFormat(kFALSE);
227 //Updated fitting routine from 2007 beam test takes into account
228 //possibility of two peaks in data and selects first one for fitting
229 //Also sets some of the starting parameters based on the shape of the
230 //given raw signal being fit
232 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
233 signalF->SetParNames("amp","t0","tau","N","ped");
234 signalF->SetParameter(2,2.35); // tau in units of time bin
235 signalF->SetParLimits(2,2,-1);
236 signalF->SetParameter(3,2); // order
237 signalF->SetParLimits(3,2,-1);
243 //Graph to hold data we will fit (should be converted to an array
244 //later to speed up processing
245 TGraph * gSig = new TGraph(GetRawFormatTimeBins());
250 while (readOk && in.GetModule() < 0)
251 readOk = in.Next(); // Go to first digit
257 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
258 lowGain = in.IsLowGain();
259 Int_t maxTime = in.GetTime(); // timebins come in reverse order
260 if (maxTime < 0 || maxTime >= GetRawFormatTimeBins()) {
261 AliWarning(Form("Invalid time bin %d",maxTime));
262 maxTime = GetRawFormatTimeBins();
264 gSig->Set(maxTime+1);
265 // There is some kind of zero-suppression in the raw data,
266 // so set up the TGraph in advance
267 for (Int_t i=0; i < maxTime; i++) {
268 gSig->SetPoint(i, i , 0);
273 if (in.GetTime() >= gSig->GetN()) {
274 AliWarning("Too many time bins");
275 gSig->Set(in.GetTime());
277 col = in.GetColumn();
280 gSig->SetPoint(in.GetTime(), in.GetTime(), in.GetSignal()) ;
281 if (in.GetTime() > maxTime)
282 maxTime = in.GetTime();
284 } while ((readOk = in.Next()) && !in.IsNewHWAddress());
286 FitRaw(gSig, signalF, amp, time) ;
289 AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
290 //cout << "col " << col-40 << " row " << row-8 << " lowGain " << lowGain << " amp " << amp << endl;
291 AddDigit(digitsArr, id, lowGain, (Int_t)amp, time);
295 for (Int_t index = 0; index < gSig->GetN(); index++) {
296 gSig->SetPoint(index, index, 0) ;
298 }; // EMCAL entries loop
306 //____________________________________________________________________________
307 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
310 // This routine checks whether a digit exists already for this tower
311 // and then decides whether to use the high or low gain info
313 // Called by Raw2Digits
315 AliEMCALDigit *digit = 0, *tmpdigit = 0;
317 TIter nextdigit(digitsArr);
318 while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
319 if (tmpdigit->GetId() == id)
323 if (!digit) { // no digit existed for this tower; create one
325 amp = Int_t(fHighLowGainFactor * amp);
326 Int_t idigit = digitsArr->GetEntries();
327 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;
329 else { // a digit already exists, check range
330 // (use high gain if signal < 800, otherwise low gain)
331 if (lowGain) { // new digit is low gain
332 if (digit->GetAmp() > 800) { // use if stored digit is out of range
333 digit->SetAmp(Int_t(fHighLowGainFactor * amp));
334 digit->SetTime(time);
337 else if (amp < 800) { // new digit is high gain; use if not out of range
339 digit->SetTime(time);
344 //____________________________________________________________________________
345 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Float_t & amp, Float_t & time)
347 // Fits the raw signal time distribution; from AliEMCALGetter
349 const Int_t kNoiseThreshold = 5;
350 const Int_t kNPedSamples = 5;
355 for (Int_t index = 0; index < kNPedSamples; index++) {
356 Double_t ttime, signal;
357 gSig->GetPoint(index, ttime, signal) ;
367 AliWarning("Could not determine pedestal");
368 ped = 10; // put some small value as first guess
374 Float_t max_fit = gSig->GetN();
375 Float_t min_after_sig = 9999;
376 Int_t tmin_after_sig = gSig->GetN();
377 Int_t n_ped_after_sig = 0;
379 for (Int_t i=kNPedSamples; i < gSig->GetN(); i++) {
380 Double_t ttime, signal;
381 gSig->GetPoint(i, ttime, signal) ;
382 if (!max_found && signal > max) {
386 else if ( max > ped + kNoiseThreshold ) {
388 min_after_sig = signal;
392 if ( signal < min_after_sig) {
393 min_after_sig = signal;
396 if (i > tmin_after_sig + 5) { // Two close peaks; end fit at minimum
397 max_fit = tmin_after_sig;
400 if ( signal < ped + kNoiseThreshold)
402 if (n_ped_after_sig >= 5) { // include 5 pedestal bins after peak
409 if ( max - ped > kNoiseThreshold ) { // else its noise
410 AliDebug(2,Form("Fitting max %d ped %d", max, ped));
411 signalF->SetRange(0,max_fit);
414 signalF->SetParLimits(2,1,3);
416 signalF->SetParameter(4, ped) ;
417 signalF->SetParameter(1, i_max);
418 signalF->SetParameter(0, max);
420 gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
421 amp = signalF->GetParameter(0);
422 time = signalF->GetParameter(1)*GetRawFormatTimeBinWidth() - fgTimeTrigger;
426 //__________________________________________________________________
427 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
429 // Matches version used in 2007 beam test
431 // Shape of the electronics raw reponse:
432 // It is a semi-gaussian, 2nd order Gamma function of the general form
434 // t' = (t - t0 + tau) / tau
435 // F = A * t**N * exp( N * ( 1 - t) ) for t >= 0
439 // A: par[0] // Amplitude = peak value
446 Double_t tau =par[2];
448 Double_t ped = par[4];
449 Double_t xx = ( x[0] - par[1] + tau ) / tau ;
454 signal = ped + par[0] * TMath::Power(xx , N) * TMath::Exp(N * (1 - xx )) ;
459 //__________________________________________________________________
460 Bool_t AliEMCALRawUtils::RawSampledResponse(
461 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
463 // for a start time dtime and an amplitude damp given by digit,
464 // calculates the raw sampled response AliEMCAL::RawResponseFunction
466 const Int_t kRawSignalOverflow = 0x3FF ;
467 const Int_t pedVal = 32;
468 Bool_t lowGain = kFALSE ;
470 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
471 signalF.SetParameter(0, pedVal) ;
472 signalF.SetParameter(1, damp) ;
473 signalF.SetParameter(2, dtime + fgTimeTrigger) ;
475 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
476 Double_t time = iTime * GetRawFormatTimeBinWidth() ;
477 Double_t signal = signalF.Eval(time) ;
478 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
479 if ( adcH[iTime] > kRawSignalOverflow ){ // larger than 10 bits
480 adcH[iTime] = kRawSignalOverflow ;
484 signal /= fHighLowGainFactor;
486 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
487 if ( adcL[iTime] > kRawSignalOverflow) // larger than 10 bits
488 adcL[iTime] = kRawSignalOverflow ;