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.10 2007/12/06 13:58:11 hristov
21 * Additional pritection. Do not delete the mapping, it is owned by another class
23 * Revision 1.9 2007/12/06 02:19:51 jklay
24 * incorporated fitting procedure from testbeam analysis into AliRoot
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
29 * Revision 1.7 2007/11/14 15:51:46 gustavo
30 * Take out few unnecessary prints
32 * Revision 1.6 2007/11/01 01:23:51 mvl
33 * Removed call to SetOldRCUFormat, which is only needed for testbeam data
35 * Revision 1.5 2007/11/01 01:20:33 mvl
36 * Further improvement of peak finding; more robust fit
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
41 * Revision 1.3 2007/09/27 08:36:46 mvl
42 * More robust setting of fit range in FitRawSignal (P. Hristov)
44 * Revision 1.2 2007/09/03 20:55:35 jklay
45 * EMCAL e-by-e reconstruction methods from Cvetan
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.
51 //*-- Author: Marco van Leeuwen (LBL)
52 #include "AliEMCALRawUtils.h"
60 #include "AliRunLoader.h"
61 #include "AliCaloAltroMapping.h"
62 #include "AliAltroBuffer.h"
63 #include "AliRawReader.h"
64 #include "AliCaloRawStream.h"
67 #include "AliEMCALRecParam.h"
68 #include "AliEMCALLoader.h"
69 #include "AliEMCALGeometry.h"
70 #include "AliEMCALDigitizer.h"
71 #include "AliEMCALDigit.h"
74 ClassImp(AliEMCALRawUtils)
76 // Signal shape parameters
77 Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ; // each sample is 100 ns
78 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
80 // some digitization constants
81 Int_t AliEMCALRawUtils::fgThreshold = 1;
82 Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
83 Int_t AliEMCALRawUtils::fgPedestalValue = 32; // pedestal value for digits2raw
84 Double_t AliEMCALRawUtils::fgFEENoise = 3.; // 3 ADC channels of noise (sampled)
86 AliEMCALRawUtils::AliEMCALRawUtils()
87 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
88 fNPedSamples(0), fGeom(0), fOption("")
91 //These are default parameters.
92 //Can be re-set from without with setter functions
93 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
94 fOrder = 2; // order of gamma fn
95 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
99 //Get Mapping RCU files from the AliEMCALRecParam
100 const TObjArray* maps = AliEMCALRecParam::GetMappings();
101 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
103 for(Int_t i = 0; i < 4; i++) {
104 fMapping[i] = (AliAltroMapping*)maps->At(i);
107 //To make sure we match with the geometry in a simulation file,
108 //let's try to get it first. If not, take the default geometry
109 AliRunLoader *rl = AliRunLoader::Instance();
110 if(!rl) AliError("Cannot find RunLoader!");
111 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
112 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
114 AliInfo(Form("Using default geometry in raw reco"));
115 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
118 if(!fGeom) AliFatal(Form("Could not get geometry!"));
122 //____________________________________________________________________________
123 AliEMCALRawUtils::AliEMCALRawUtils(AliEMCALGeometry *pGeometry)
124 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
125 fNPedSamples(0), fGeom(pGeometry), fOption("")
128 // Initialize with the given geometry - constructor required by HLT
129 // HLT does not use/support AliRunLoader(s) instances
130 // This is a minimum intervention solution
131 // Comment by MPloskon@lbl.gov
134 //These are default parameters.
135 //Can be re-set from without with setter functions
136 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
137 fOrder = 2; // order of gamma fn
138 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
142 //Get Mapping RCU files from the AliEMCALRecParam
143 const TObjArray* maps = AliEMCALRecParam::GetMappings();
144 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
146 for(Int_t i = 0; i < 4; i++) {
147 fMapping[i] = (AliAltroMapping*)maps->At(i);
150 if(!fGeom) AliFatal(Form("Could not get geometry!"));
154 //____________________________________________________________________________
155 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
157 fHighLowGainFactor(rawU.fHighLowGainFactor),
160 fNoiseThreshold(rawU.fNoiseThreshold),
161 fNPedSamples(rawU.fNPedSamples),
163 fOption(rawU.fOption)
166 fMapping[0] = rawU.fMapping[0];
167 fMapping[1] = rawU.fMapping[1];
168 fMapping[2] = rawU.fMapping[2];
169 fMapping[3] = rawU.fMapping[3];
172 //____________________________________________________________________________
173 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
175 //assignment operator
178 fHighLowGainFactor = rawU.fHighLowGainFactor;
179 fOrder = rawU.fOrder;
181 fNoiseThreshold = rawU.fNoiseThreshold;
182 fNPedSamples = rawU.fNPedSamples;
184 fOption = rawU.fOption;
185 fMapping[0] = rawU.fMapping[0];
186 fMapping[1] = rawU.fMapping[1];
187 fMapping[2] = rawU.fMapping[2];
188 fMapping[3] = rawU.fMapping[3];
195 //____________________________________________________________________________
196 AliEMCALRawUtils::~AliEMCALRawUtils() {
200 //____________________________________________________________________________
201 void AliEMCALRawUtils::Digits2Raw()
203 // convert digits of the current event to raw data
205 AliRunLoader *rl = AliRunLoader::Instance();
206 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
209 loader->LoadDigits("EMCAL");
211 TClonesArray* digits = loader->Digits() ;
214 Warning("Digits2Raw", "no digits found !");
218 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
219 AliAltroBuffer* buffers[nDDL];
220 for (Int_t i=0; i < nDDL; i++)
223 Int_t adcValuesLow[fgkTimeBins];
224 Int_t adcValuesHigh[fgkTimeBins];
226 // loop over digits (assume ordered digits)
227 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
228 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
229 if (digit->GetAmp() < fgThreshold)
239 fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
240 fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
242 //Check which is the RCU, 0 or 1, of the cell.
245 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
246 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
249 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
251 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
253 if (nSM%2==1) iRCU = 1 - iRCU; // swap for odd=C side, to allow us to cable both sides the same
256 Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
259 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
261 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
263 if (buffers[iDDL] == 0) {
264 // open new file and write dummy header
265 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
266 //Select mapping file RCU0A, RCU0C, RCU1A, RCU1C
267 Int_t iRCUside=iRCU+(nSM%2)*2;
268 //iRCU=0 and even (0) SM -> RCU0A.data 0
269 //iRCU=1 and even (0) SM -> RCU1A.data 1
270 //iRCU=0 and odd (1) SM -> RCU0C.data 2
271 //iRCU=1 and odd (1) SM -> RCU1C.data 3
272 //cout<<" nSM "<<nSM<<"; iRCU "<<iRCU<<"; iRCUside "<<iRCUside<<endl;
273 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCUside]);
274 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
277 // out of time range signal (?)
278 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
279 AliInfo("Signal is out of time range.\n");
280 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
281 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
282 buffers[iDDL]->FillBuffer(3); // bunch length
283 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
284 // calculate the time response function
286 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ;
288 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
290 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
294 // write headers and close files
295 for (Int_t i=0; i < nDDL; i++) {
298 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
303 loader->UnloadDigits();
306 //____________________________________________________________________________
307 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr)
309 // convert raw data of the current event to digits
314 Error("Raw2Digits", "no digits found !");
318 Error("Raw2Digits", "no raw reader found !");
322 AliCaloRawStream in(reader,"EMCAL",fMapping);
323 // Select EMCAL DDL's;
324 reader->Select("EMCAL");
326 //Updated fitting routine from 2007 beam test takes into account
327 //possibility of two peaks in data and selects first one for fitting
328 //Also sets some of the starting parameters based on the shape of the
329 //given raw signal being fit
331 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
332 signalF->SetParameters(10.,0.,fTau,fOrder,5.); //set all defaults once, just to be safe
333 signalF->SetParNames("amp","t0","tau","N","ped");
334 signalF->SetParameter(2,fTau); // tau in units of time bin
335 signalF->SetParLimits(2,2,-1);
336 signalF->SetParameter(3,fOrder); // order
337 signalF->SetParLimits(3,2,-1);
343 //Graph to hold data we will fit (should be converted to an array
344 //later to speed up processing
345 TGraph * gSig = new TGraph(GetRawFormatTimeBins());
350 while (readOk && in.GetModule() < 0)
351 readOk = in.Next(); // Go to first digit
358 id = fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
359 lowGain = in.IsLowGain();
360 Int_t maxTime = in.GetTime(); // timebins come in reverse order
361 if (maxTime < 0 || maxTime >= GetRawFormatTimeBins()) {
362 AliWarning(Form("Invalid time bin %d",maxTime));
363 maxTime = GetRawFormatTimeBins();
365 gSig->Set(maxTime+1);
366 // There is some kind of zero-suppression in the raw data,
367 // so set up the TGraph in advance
368 for (Int_t i=0; i < maxTime; i++) {
369 gSig->SetPoint(i, i , 0);
374 if (in.GetTime() >= gSig->GetN()) {
375 AliWarning("Too many time bins");
376 gSig->Set(in.GetTime());
378 col = in.GetColumn();
381 gSig->SetPoint(in.GetTime(), in.GetTime(), in.GetSignal()) ;
383 if (in.GetTime() > maxTime)
384 maxTime = in.GetTime();
386 } while ((readOk = in.Next()) && !in.IsNewHWAddress());
388 FitRaw(gSig, signalF, amp, time) ;
390 if (amp > 0 && amp < 2000) { //check both high and low end of
391 //result, 2000 is somewhat arbitrary
392 AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
394 AddDigit(digitsArr, id, lowGain, (Int_t)amp, time);
398 for (Int_t index = 0; index < gSig->GetN(); index++) {
399 gSig->SetPoint(index, index, 0) ;
401 // Reset starting parameters for fit function
402 signalF->SetParameters(10.,0.,fTau,fOrder,5.); //reset all defaults just to be safe
404 }; // EMCAL entries loop
412 //____________________________________________________________________________
413 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
416 // This routine checks whether a digit exists already for this tower
417 // and then decides whether to use the high or low gain info
419 // Called by Raw2Digits
421 AliEMCALDigit *digit = 0, *tmpdigit = 0;
423 TIter nextdigit(digitsArr);
424 while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
425 if (tmpdigit->GetId() == id)
429 if (!digit) { // no digit existed for this tower; create one
431 amp = Int_t(fHighLowGainFactor * amp);
432 Int_t idigit = digitsArr->GetEntries();
433 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;
435 else { // a digit already exists, check range
436 // (use high gain if signal < cut value, otherwise low gain)
437 if (lowGain) { // new digit is low gain
438 if (digit->GetAmp() > fgkOverflowCut) { // use if stored digit is out of range
439 digit->SetAmp(Int_t(fHighLowGainFactor * amp));
440 digit->SetTime(time);
443 else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
445 digit->SetTime(time);
450 //____________________________________________________________________________
451 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Float_t & amp, Float_t & time)
453 // Fits the raw signal time distribution; from AliEMCALGetter
459 for (Int_t index = 0; index < fNPedSamples; index++) {
460 Double_t ttime, signal;
461 gSig->GetPoint(index, ttime, signal) ;
471 AliWarning("Could not determine pedestal");
472 ped = 10; // put some small value as first guess
478 Float_t max_fit = gSig->GetN();
479 Float_t min_after_sig = 9999;
480 Int_t tmin_after_sig = gSig->GetN();
481 Int_t n_ped_after_sig = 0;
482 Int_t plateau_width = 0;
483 Int_t plateau_start = 9999;
486 for (Int_t i=fNPedSamples; i < gSig->GetN(); i++) {
487 Double_t ttime, signal;
488 gSig->GetPoint(i, ttime, signal) ;
489 if (!max_found && signal > max) {
493 else if ( max > ped + fNoiseThreshold ) {
495 min_after_sig = signal;
499 if ( signal < min_after_sig) {
500 min_after_sig = signal;
503 if (i > tmin_after_sig + 5) { // Two close peaks; end fit at minimum
504 max_fit = tmin_after_sig;
507 if ( signal < Cut*max){ //stop fit at 30% amplitude(avoid the pulse shape falling edge)
511 if ( signal < ped + fNoiseThreshold)
513 if (n_ped_after_sig >= 5) { // include 5 pedestal bins after peak
518 //Add check on plateau
519 if (signal >= fgkRawSignalOverflow - fNoiseThreshold) {
520 if(plateau_width == 0) plateau_start = i;
525 if(plateau_width > 0) {
526 for(int j = 0; j < plateau_width; j++) {
527 //Note, have to remove the same point N times because after each
528 //remove, the positions of all subsequent points have shifted down
529 gSig->RemovePoint(plateau_start);
533 if ( max - ped > fNoiseThreshold ) { // else its noise
534 AliDebug(2,Form("Fitting max %d ped %d", max, ped));
535 signalF->SetRange(0,max_fit);
538 signalF->SetParLimits(2,1,3);
540 signalF->SetParameter(4, ped) ;
541 signalF->SetParameter(1, i_max);
542 signalF->SetParameter(0, max);
544 gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
545 amp = signalF->GetParameter(0);
546 time = signalF->GetParameter(1)*GetRawFormatTimeBinWidth() - fgTimeTrigger;
550 //__________________________________________________________________
551 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
553 // Matches version used in 2007 beam test
555 // Shape of the electronics raw reponse:
556 // It is a semi-gaussian, 2nd order Gamma function of the general form
558 // t' = (t - t0 + tau) / tau
559 // F = A * t**N * exp( N * ( 1 - t) ) for t >= 0
563 // A: par[0] // Amplitude = peak value
570 Double_t tau =par[2];
572 Double_t ped = par[4];
573 Double_t xx = ( x[0] - par[1] + tau ) / tau ;
578 signal = ped + par[0] * TMath::Power(xx , N) * TMath::Exp(N * (1 - xx )) ;
583 //__________________________________________________________________
584 Bool_t AliEMCALRawUtils::RawSampledResponse(
585 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
587 // for a start time dtime and an amplitude damp given by digit,
588 // calculates the raw sampled response AliEMCAL::RawResponseFunction
590 Bool_t lowGain = kFALSE ;
592 // A: par[0] // Amplitude = peak value
598 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
599 signalF.SetParameter(0, damp) ;
600 signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ;
601 signalF.SetParameter(2, fTau) ;
602 signalF.SetParameter(3, fOrder);
603 signalF.SetParameter(4, fgPedestalValue);
605 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
606 Double_t signal = signalF.Eval(iTime) ;
608 //According to Terry Awes, 13-Apr-2008
609 //add gaussian noise in quadrature to each sample
610 //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
611 //signal = sqrt(signal*signal + noise*noise);
613 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
614 if ( adcH[iTime] > fgkRawSignalOverflow ){ // larger than 10 bits
615 adcH[iTime] = fgkRawSignalOverflow ;
619 signal /= fHighLowGainFactor;
621 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
622 if ( adcL[iTime] > fgkRawSignalOverflow) // larger than 10 bits
623 adcL[iTime] = fgkRawSignalOverflow ;