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.2 2007/09/03 20:55:35 jklay
21 * EMCAL e-by-e reconstruction methods from Cvetan
23 * Revision 1.1 2007/03/17 19:56:38 mvl
24 * Moved signal shape routines from AliEMCAL to separate class AliEMCALRawUtils to streamline raw data reconstruction code.
27 //*-- Author: Marco van Leeuwen (LBL)
28 #include "AliEMCALRawUtils.h"
35 #include "AliRunLoader.h"
36 #include "AliCaloAltroMapping.h"
37 #include "AliAltroBuffer.h"
38 #include "AliRawReader.h"
39 #include "AliCaloRawStream.h"
42 #include "AliEMCALLoader.h"
43 #include "AliEMCALGeometry.h"
44 #include "AliEMCALDigitizer.h"
45 #include "AliEMCALDigit.h"
48 ClassImp(AliEMCALRawUtils)
50 // Signal shape parameters
51 Int_t AliEMCALRawUtils::fgOrder = 2 ; // Order of gamma function
52 Double_t AliEMCALRawUtils::fgTimeMax = 2.56E-5 ; // each sample is over 100 ns fTimeMax/fTimeBins
53 Double_t AliEMCALRawUtils::fgTau = 165E-9 ; // 165 ns (from testbeam; not very accurate)
54 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
56 // some digitization constants
57 Int_t AliEMCALRawUtils::fgThreshold = 1;
58 Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
60 AliEMCALRawUtils::AliEMCALRawUtils(): fHighLowGainFactor(0.) {
61 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
63 //____________________________________________________________________________
64 AliEMCALRawUtils::~AliEMCALRawUtils() {
66 //____________________________________________________________________________
67 void AliEMCALRawUtils::Digits2Raw()
69 // convert digits of the current event to raw data
71 AliRunLoader *rl = AliRunLoader::GetRunLoader();
72 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
75 loader->LoadDigits("EMCAL");
77 TClonesArray* digits = loader->Digits() ;
80 Warning("Digits2Raw", "no digits found !");
85 AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
87 AliError(Form("No geometry found !"));
91 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
92 AliAltroBuffer* buffers[nDDL];
93 for (Int_t i=0; i < nDDL; i++)
96 Int_t adcValuesLow[fgkTimeBins];
97 Int_t adcValuesHigh[fgkTimeBins];
99 //Load Mapping RCU files once
100 TString path = gSystem->Getenv("ALICE_ROOT");
101 path += "/EMCAL/mapping/RCU";
102 TString path0 = path+"0.data";//This file will change in future
103 TString path1 = path+"1.data";//This file will change in future
104 AliAltroMapping * mapping[2] ; // For the moment only 2
105 mapping[0] = new AliCaloAltroMapping(path0.Data());
106 mapping[1] = new AliCaloAltroMapping(path1.Data());
108 // loop over digits (assume ordered digits)
109 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
110 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
111 if (digit->GetAmp() < fgThreshold)
121 geom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
122 geom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
124 //Check which is the RCU of the cell.
127 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
128 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
131 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
133 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
136 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
138 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
140 if (buffers[iDDL] == 0) {
141 // open new file and write dummy header
142 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
143 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),mapping[iRCU]);
144 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
147 // out of time range signal (?)
148 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
149 AliInfo("Signal is out of time range.\n");
150 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
151 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
152 buffers[iDDL]->FillBuffer(3); // bunch length
153 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
154 // calculate the time response function
156 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ;
158 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
160 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
164 // write headers and close files
165 for (Int_t i=0; i < nDDL; i++) {
168 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
172 mapping[0]->Delete();
173 mapping[1]->Delete();
174 loader->UnloadDigits();
177 //____________________________________________________________________________
178 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr)
180 // convert raw data of the current event to digits
181 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
183 AliError(Form("No geometry found !"));
190 Error("Raw2Digits", "no digits found !");
194 Error("Raw2Digits", "no raw reader found !");
198 // Use AliAltroRawStream to read the ALTRO format. No need to
199 // reinvent the wheel :-)
200 AliCaloRawStream in(reader,"EMCAL");
201 // Select EMCAL DDL's;
202 reader->Select("EMCAL");
204 // reading is from previously existing AliEMCALGetter.cxx
206 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
207 signalF->SetParNames("Charge", "Gain", "Amplitude", "TimeZero");
214 TGraph * gSig = new TGraph(GetRawFormatTimeBins()) ;
216 Int_t eofReached = 0;
219 in.Next(); // Go to first digit
221 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
222 lowGain = in.IsLowGain();
223 gSig->SetPoint(in.GetTime(),
224 in.GetTime()* GetRawFormatTimeMax() / GetRawFormatTimeBins(),
232 gSig->SetPoint(in.GetTime(),
233 in.GetTime()* GetRawFormatTimeMax() / GetRawFormatTimeBins(),
237 } while (!eofReached && !in.IsNewRow() && !in.IsNewColumn() && !in.IsNewModule());
239 FitRaw(gSig, signalF, amp, time) ;
241 amp *= fHighLowGainFactor;
244 AliDebug(2,Form("id %d amp %g", id, amp));
245 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, (Int_t)amp, time, idigit) ;
251 for (index = 0; index < GetRawFormatTimeBins(); index++) {
252 gSig->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ;
254 } while (!eofReached); // EMCAL entries loop
262 //____________________________________________________________________________
263 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Double_t & amp, Double_t & time)
265 // Fits the raw signal time distribution; from AliEMCALGetter
267 const Int_t kNoiseThreshold = 0 ;
268 Double_t timezero1 = 0., timezero2 = 0., timemax = 0. ;
269 Double_t signal = 0., signalmax = 0. ;
272 timezero1 = signalmax = timemax = 0. ;
273 timezero1 = GetRawFormatTimeMax();
276 // Find the time,the value and the index of the maximum
279 for (Int_t index = 0; index < GetRawFormatTimeBins(); index++) {
280 gSig->GetPoint(index, time, signal) ;
281 if (signal >= signalmax) {
288 if ( timemax < GetRawFormatTimeMax() * 0.4 ) { // else its noise
290 // Find the start time of the signal;
291 for (Int_t index = indexmax; index >= 0; index--) {
292 gSig->GetPoint(index, time, signal) ;
294 if (signal>kNoiseThreshold) {
299 // Find the end time of the signal;
300 for (Int_t index = indexmax; index < GetRawFormatTimeBins(); index++) {
301 gSig->GetPoint(index, time, signal) ;
302 if (signal>kNoiseThreshold) {
308 signalF->SetParameter(0, signalmax) ;
309 // signalF->SetParameter(1, timemax) ;
310 signalF->SetParameter(1, 0.5*(timezero1+timezero2)) ;
311 gSig->Fit(signalF, "QRON", "", timezero1, timezero2); //, "QRON") ;
312 amp = signalF->GetParameter(0) ;
313 time = signalF->GetParameter(1) - fgTimeTrigger;
317 //__________________________________________________________________
318 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
320 // Shape of the electronics raw reponse:
321 // It is a semi-gaussian, 2nd order Gamma function of the general form
323 // t' = (t - t0 + tau) / tau
324 // F = A * t**N * exp( N * ( 1 - t) ) for t >= 0
328 // A: par[0] // Amplitude = peak value
334 Double_t xx = ( x[0] - par[1] + fgTau ) / fgTau ;
339 signal = par[0] * TMath::Power(xx , fgOrder) * TMath::Exp(fgOrder * (1 - xx )) ;
344 //__________________________________________________________________
345 Bool_t AliEMCALRawUtils::RawSampledResponse(
346 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
348 // for a start time dtime and an amplitude damp given by digit,
349 // calculates the raw sampled response AliEMCAL::RawResponseFunction
351 const Int_t kRawSignalOverflow = 0x3FF ;
352 Bool_t lowGain = kFALSE ;
354 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
355 signalF.SetParameter(0, damp) ;
356 signalF.SetParameter(1, dtime + fgTimeTrigger) ;
358 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
359 Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
360 Double_t signal = signalF.Eval(time) ;
361 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
362 if ( adcH[iTime] > kRawSignalOverflow ){ // larger than 10 bits
363 adcH[iTime] = kRawSignalOverflow ;
367 signal /= fHighLowGainFactor;
369 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
370 if ( adcL[iTime] > kRawSignalOverflow) // larger than 10 bits
371 adcL[iTime] = kRawSignalOverflow ;