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.6 2007/11/01 01:23:51 mvl
21 * Removed call to SetOldRCUFormat, which is only needed for testbeam data
23 * Revision 1.5 2007/11/01 01:20:33 mvl
24 * Further improvement of peak finding; more robust fit
26 * Revision 1.4 2007/10/31 17:15:24 mvl
27 * Fixed bug in raw data unpacking; Added pedestal to signal fit; Added logic to deal with high/low gain
29 * Revision 1.3 2007/09/27 08:36:46 mvl
30 * More robust setting of fit range in FitRawSignal (P. Hristov)
32 * Revision 1.2 2007/09/03 20:55:35 jklay
33 * EMCAL e-by-e reconstruction methods from Cvetan
35 * Revision 1.1 2007/03/17 19:56:38 mvl
36 * Moved signal shape routines from AliEMCAL to separate class AliEMCALRawUtils to streamline raw data reconstruction code.
39 //*-- Author: Marco van Leeuwen (LBL)
40 #include "AliEMCALRawUtils.h"
47 #include "AliRunLoader.h"
48 #include "AliCaloAltroMapping.h"
49 #include "AliAltroBuffer.h"
50 #include "AliRawReader.h"
51 #include "AliCaloRawStream.h"
54 #include "AliEMCALLoader.h"
55 #include "AliEMCALGeometry.h"
56 #include "AliEMCALDigitizer.h"
57 #include "AliEMCALDigit.h"
60 ClassImp(AliEMCALRawUtils)
62 // Signal shape parameters
63 Int_t AliEMCALRawUtils::fgOrder = 2 ; // Order of gamma function
64 Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ; // each sample is 100 ns
65 Double_t AliEMCALRawUtils::fgTau = 235E-9 ; // 235 ns (from CERN testbeam; not very accurate)
66 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
68 // some digitization constants
69 Int_t AliEMCALRawUtils::fgThreshold = 1;
70 Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
72 AliEMCALRawUtils::AliEMCALRawUtils(): fHighLowGainFactor(0.) {
73 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
75 //____________________________________________________________________________
76 AliEMCALRawUtils::~AliEMCALRawUtils() {
78 //____________________________________________________________________________
79 void AliEMCALRawUtils::Digits2Raw()
81 // convert digits of the current event to raw data
83 AliRunLoader *rl = AliRunLoader::GetRunLoader();
84 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
87 loader->LoadDigits("EMCAL");
89 TClonesArray* digits = loader->Digits() ;
92 Warning("Digits2Raw", "no digits found !");
97 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
99 AliError(Form("No geometry found !"));
103 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
104 AliAltroBuffer* buffers[nDDL];
105 for (Int_t i=0; i < nDDL; i++)
108 Int_t adcValuesLow[fgkTimeBins];
109 Int_t adcValuesHigh[fgkTimeBins];
111 //Load Mapping RCU files once
112 TString path = gSystem->Getenv("ALICE_ROOT");
113 path += "/EMCAL/mapping/RCU";
114 TString path0 = path+"0.data";//This file will change in future
115 TString path1 = path+"1.data";//This file will change in future
116 AliAltroMapping * mapping[2] ; // For the moment only 2
117 mapping[0] = new AliCaloAltroMapping(path0.Data());
118 mapping[1] = new AliCaloAltroMapping(path1.Data());
120 // loop over digits (assume ordered digits)
121 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
122 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
123 if (digit->GetAmp() < fgThreshold)
133 geom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
134 geom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
136 //Check which is the RCU of the cell.
139 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
140 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
143 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
145 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
148 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
150 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
152 if (buffers[iDDL] == 0) {
153 // open new file and write dummy header
154 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
155 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),mapping[iRCU]);
156 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
159 // out of time range signal (?)
160 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
161 AliInfo("Signal is out of time range.\n");
162 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
163 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
164 buffers[iDDL]->FillBuffer(3); // bunch length
165 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
166 // calculate the time response function
168 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ;
170 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
172 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
176 // write headers and close files
177 for (Int_t i=0; i < nDDL; i++) {
180 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
184 mapping[0]->Delete();
185 mapping[1]->Delete();
186 loader->UnloadDigits();
189 //____________________________________________________________________________
190 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr)
192 // convert raw data of the current event to digits
193 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
195 AliError(Form("No geometry found !"));
202 Error("Raw2Digits", "no digits found !");
206 Error("Raw2Digits", "no raw reader found !");
210 // Use AliAltroRawStream to read the ALTRO format. No need to
211 // reinvent the wheel :-)
212 AliCaloRawStream in(reader,"EMCAL");
213 // Select EMCAL DDL's;
214 reader->Select("EMCAL");
215 //in.SetOldRCUFormat(kTRUE); // Needed for testbeam data
217 // reading is from previously existing AliEMCALGetter.cxx
219 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
225 TGraph * gSig = new TGraph(GetRawFormatTimeBins()) ;
230 while (readOk && in.GetModule() < 0)
231 readOk = in.Next(); // Go to first digit
234 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
235 lowGain = in.IsLowGain();
236 Int_t maxTime = in.GetTime(); // timebins come in reverse order
237 if (maxTime < 0 || maxTime >= GetRawFormatTimeBins()) {
238 AliWarning(Form("Invalid time bin %d",maxTime));
239 maxTime = GetRawFormatTimeBins();
241 gSig->Set(maxTime+1);
242 // There is some kind of zero-suppression in the raw data,
243 // so set up the TGraph in advance
244 for (Int_t i=0; i < maxTime; i++) {
245 gSig->SetPoint(i, i * GetRawFormatTimeBinWidth(), 0);
250 if (in.GetTime() >= gSig->GetN()) {
251 AliWarning("Too many time bins");
252 gSig->Set(in.GetTime());
254 gSig->SetPoint(in.GetTime(),
255 in.GetTime() * GetRawFormatTimeBinWidth(),
257 if (in.GetTime() > maxTime)
258 maxTime = in.GetTime();
260 } while ((readOk = in.Next()) && !in.IsNewHWAddress());
261 signalF->SetRange(0,(Float_t)maxTime*GetRawFormatTimeBinWidth());
263 FitRaw(gSig, signalF, amp, time) ;
266 AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
267 AddDigit(digitsArr, id, lowGain, (Int_t)amp, time);
271 for (Int_t index = 0; index < gSig->GetN(); index++) {
272 gSig->SetPoint(index, index * GetRawFormatTimeBinWidth(), 0) ;
274 }; // EMCAL entries loop
282 //____________________________________________________________________________
283 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
286 // This routine checks whether a digit exists already for this tower
287 // and then decides whether to use the high or low gain info
289 // Called by Raw2Digits
291 AliEMCALDigit *digit = 0, *tmpdigit = 0;
293 TIter nextdigit(digitsArr);
294 while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
295 if (tmpdigit->GetId() == id)
299 if (!digit) { // no digit existed for this tower; create one
301 amp = Int_t(fHighLowGainFactor * amp);
302 Int_t idigit = digitsArr->GetEntries();
303 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;
305 else { // a digit already exists, check range
306 // (use high gain if signal < 800, otherwise low gain)
307 if (lowGain) { // new digit is low gain
308 if (digit->GetAmp() > 800) { // use if stored digit is out of range
309 digit->SetAmp(Int_t(fHighLowGainFactor * amp));
310 digit->SetTime(time);
313 else if (amp < 800) { // new digit is high gain; use if not out of range
315 digit->SetTime(time);
320 //____________________________________________________________________________
321 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Float_t & amp, Float_t & time)
323 // Fits the raw signal time distribution; from AliEMCALGetter
325 const Int_t kNoiseThreshold = 5;
326 const Int_t kNPedSamples = 10;
331 for (Int_t index = 0; index < kNPedSamples; index++) {
332 Double_t ttime, signal;
333 gSig->GetPoint(index, ttime, signal) ;
343 AliWarning("Could determine pedestal");
344 ped = 10; // put some small value as first guess
351 Float_t max_fit = gSig->GetN()*GetRawFormatTimeBinWidth();
352 Float_t min_after_sig = 9999;
353 Int_t imin_after_sig = gSig->GetN();
354 Float_t tmin_after_sig = gSig->GetN()*GetRawFormatTimeBinWidth();
355 Int_t n_ped_after_sig = 0;
357 for (Int_t i=kNPedSamples; i < gSig->GetN(); i++) {
358 Double_t ttime, signal;
359 gSig->GetPoint(i, ttime, signal) ;
360 if (!max_found && signal > max) {
365 else if ( max > ped + kNoiseThreshold ) {
367 min_after_sig = signal;
369 tmin_after_sig = ttime;
372 if ( signal < min_after_sig) {
373 min_after_sig = signal;
375 tmin_after_sig = ttime;
377 if (i > tmin_after_sig + 5) { // Two close peaks; end fit at minimum
378 max_fit = tmin_after_sig;
381 if ( signal < ped + kNoiseThreshold)
383 if (n_ped_after_sig >= 5) { // include 5 pedestal bins after peak
390 if ( max - ped > kNoiseThreshold ) { // else its noise
391 AliDebug(2,Form("Fitting max %d ped %d", max, ped));
392 signalF->SetParameter(0, ped) ;
393 signalF->SetParameter(1, max - ped) ;
394 signalF->SetParameter(2, tmax) ;
395 signalF->SetParLimits(2, 0, max_fit) ;
396 gSig->Fit(signalF, "QRON", "", 0., max_fit); //, "QRON") ;
397 amp = signalF->GetParameter(1);
398 time = signalF->GetParameter(2) - fgTimeTrigger;
402 //__________________________________________________________________
403 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
405 // Shape of the electronics raw reponse:
406 // It is a semi-gaussian, 2nd order Gamma function of the general form
408 // t' = (t - t0 + tau) / tau
409 // F = A * t**N * exp( N * ( 1 - t) ) for t >= 0
414 // A: par[1] // Amplitude = peak value
420 Double_t xx = ( x[0] - par[2] + fgTau ) / fgTau ;
425 signal = par[0] + par[1] * TMath::Power(xx , fgOrder) * TMath::Exp(fgOrder * (1 - xx )) ;
431 //__________________________________________________________________
432 Bool_t AliEMCALRawUtils::RawSampledResponse(
433 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
435 // for a start time dtime and an amplitude damp given by digit,
436 // calculates the raw sampled response AliEMCAL::RawResponseFunction
438 const Int_t kRawSignalOverflow = 0x3FF ;
439 const Int_t pedVal = 32;
440 Bool_t lowGain = kFALSE ;
442 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
443 signalF.SetParameter(0, pedVal) ;
444 signalF.SetParameter(1, damp) ;
445 signalF.SetParameter(2, dtime + fgTimeTrigger) ;
447 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
448 Double_t time = iTime * GetRawFormatTimeBinWidth() ;
449 Double_t signal = signalF.Eval(time) ;
450 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
451 if ( adcH[iTime] > kRawSignalOverflow ){ // larger than 10 bits
452 adcH[iTime] = kRawSignalOverflow ;
456 signal /= fHighLowGainFactor;
458 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
459 if ( adcL[iTime] > kRawSignalOverflow) // larger than 10 bits
460 adcL[iTime] = kRawSignalOverflow ;