]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRawUtils.cxx
Remove QA reference from EMCAL directory, already available in QAref/EMCAL
[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
18 //_________________________________________________________________________
19 //  Utility Class for handling Raw data
20 //  Does all transitions from Digits to Raw and vice versa, 
21 //  for simu and reconstruction
22 //
23 //  Note: the current version is still simplified. Only 
24 //    one raw signal per digit is generated; either high-gain or low-gain
25 //    Need to add concurrent high and low-gain info in the future
26 //    No pedestal is added to the raw signal.
27 //*-- Author: Marco van Leeuwen (LBL)
28
29 #include "AliEMCALRawUtils.h"
30   
31 #include "TF1.h"
32 #include "TGraph.h"
33 #include <TRandom.h>
34 class TSystem;
35   
36 class AliLog;
37 #include "AliRun.h"
38 #include "AliRunLoader.h"
39 class AliCaloAltroMapping;
40 #include "AliAltroBuffer.h"
41 #include "AliRawReader.h"
42 #include "AliCaloRawStreamV3.h"
43 #include "AliDAQ.h"
44   
45 #include "AliEMCALRecParam.h"
46 #include "AliEMCALLoader.h"
47 #include "AliEMCALGeometry.h"
48 class AliEMCALDigitizer;
49 #include "AliEMCALDigit.h"
50 #include "AliEMCALRawDigit.h"
51 #include "AliEMCAL.h"
52 #include "AliCaloCalibPedestal.h"  
53 #include "AliCaloFastAltroFitv0.h"
54 #include "AliCaloNeuralFit.h"
55 #include "AliCaloBunchInfo.h"
56 #include "AliCaloFitResults.h"
57 #include "AliCaloRawAnalyzerFastFit.h"
58 #include "AliCaloRawAnalyzerNN.h"
59 #include "AliCaloRawAnalyzerLMS.h"
60 #include "AliCaloRawAnalyzerPeakFinder.h"
61 #include "AliCaloRawAnalyzerCrude.h"
62
63 ClassImp(AliEMCALRawUtils)
64   
65 // Signal shape parameters
66 Int_t    AliEMCALRawUtils::fgTimeBins = 256; // number of sampling bins of the raw RO signal (we typically use 15-50; theoretical max is 1k+) 
67 Double_t AliEMCALRawUtils::fgTimeBinWidth  = 100E-9 ; // each sample is 100 ns
68 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ;   // 15 time bins ~ 1.5 musec
69
70 // some digitization constants
71 Int_t    AliEMCALRawUtils::fgThreshold = 1;
72 Int_t    AliEMCALRawUtils::fgDDLPerSuperModule = 2;  // 2 ddls per SuperModule
73 Int_t    AliEMCALRawUtils::fgPedestalValue = 0;     // pedestal value for digits2raw, default generate ZS data
74 Double_t AliEMCALRawUtils::fgFEENoise = 3.;          // 3 ADC channels of noise (sampled)
75
76 AliEMCALRawUtils::AliEMCALRawUtils(fitAlgorithm fitAlgo)
77   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
78     fNPedSamples(0), fGeom(0), fOption(""),
79     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),fUseFALTRO(kFALSE),fRawAnalyzer(0)
80 {
81
82   //These are default parameters.  
83   //Can be re-set from without with setter functions
84   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
85   fHighLowGainFactor = 16. ;   // Adjusted for a low gain range of 82 GeV (10 bits) 
86   fOrder             = 2;      // Order of gamma fn
87   fTau               = 2.35;   // in units of timebin, from CERN 2007 testbeam
88   fNoiseThreshold    = 3;      // 3 ADC counts is approx. noise level
89   fNPedSamples       = 4;      // Less than this value => likely pedestal samples
90   fRemoveBadChannels = kFALSE; // Do not remove bad channels before fitting
91   fUseFALTRO         = kTRUE;  // Get the trigger FALTRO information and pass it to digits.
92   SetFittingAlgorithm(fitAlgo);
93
94   //Get Mapping RCU files from the AliEMCALRecParam                                 
95   const TObjArray* maps = AliEMCALRecParam::GetMappings();
96   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
97
98   for(Int_t i = 0; i < 4; i++) {
99     fMapping[i] = (AliAltroMapping*)maps->At(i);
100   }
101
102   //To make sure we match with the geometry in a simulation file,
103   //let's try to get it first.  If not, take the default geometry
104   AliRunLoader *rl = AliRunLoader::Instance();
105   if (rl && rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
106     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
107   } else {
108     AliInfo(Form("Using default geometry in raw reco"));
109     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
110   }
111
112   if(!fGeom) AliFatal(Form("Could not get geometry!"));
113
114 }
115
116 //____________________________________________________________________________
117 AliEMCALRawUtils::AliEMCALRawUtils(AliEMCALGeometry *pGeometry, fitAlgorithm fitAlgo)
118   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
119     fNPedSamples(0), fGeom(pGeometry), fOption(""),
120     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),fUseFALTRO(kFALSE),fRawAnalyzer()
121 {
122   //
123   // Initialize with the given geometry - constructor required by HLT
124   // HLT does not use/support AliRunLoader(s) instances
125   // This is a minimum intervention solution
126   // Comment by MPloskon@lbl.gov
127   //
128
129   //These are default parameters. 
130   //Can be re-set from without with setter functions 
131   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
132   fHighLowGainFactor = 16. ;   // adjusted for a low gain range of 82 GeV (10 bits)
133   fOrder             = 2;      // order of gamma fn
134   fTau               = 2.35;   // in units of timebin, from CERN 2007 testbeam
135   fNoiseThreshold    = 3;      // 3 ADC counts is approx. noise level
136   fNPedSamples       = 4;      // Less than this value => likely pedestal samples
137   fRemoveBadChannels = kFALSE; // Do not remove bad channels before fitting
138   fUseFALTRO         = kTRUE;  // Get the trigger FALTRO information and pass it to digits.
139   SetFittingAlgorithm(fitAlgo);
140
141   //Get Mapping RCU files from the AliEMCALRecParam
142   const TObjArray* maps = AliEMCALRecParam::GetMappings();
143   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
144
145   for(Int_t i = 0; i < 4; i++) {
146     fMapping[i] = (AliAltroMapping*)maps->At(i);
147   }
148
149   if(!fGeom) AliFatal(Form("Could not get geometry!"));
150
151 }
152
153 //____________________________________________________________________________
154 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
155   : TObject(),
156     fHighLowGainFactor(rawU.fHighLowGainFactor), 
157     fOrder(rawU.fOrder),
158     fTau(rawU.fTau),
159     fNoiseThreshold(rawU.fNoiseThreshold),
160     fNPedSamples(rawU.fNPedSamples),
161     fGeom(rawU.fGeom), 
162     fOption(rawU.fOption),
163     fRemoveBadChannels(rawU.fRemoveBadChannels),
164     fFittingAlgorithm(rawU.fFittingAlgorithm),
165         fUseFALTRO(rawU.fUseFALTRO),
166     fRawAnalyzer(rawU.fRawAnalyzer)
167 {
168   //copy ctor
169   fMapping[0] = rawU.fMapping[0];
170   fMapping[1] = rawU.fMapping[1];
171   fMapping[2] = rawU.fMapping[2];
172   fMapping[3] = rawU.fMapping[3];
173 }
174
175 //____________________________________________________________________________
176 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
177 {
178   //assignment operator
179
180   if(this != &rawU) {
181     fHighLowGainFactor = rawU.fHighLowGainFactor;
182     fOrder             = rawU.fOrder;
183     fTau               = rawU.fTau;
184     fNoiseThreshold    = rawU.fNoiseThreshold;
185     fNPedSamples       = rawU.fNPedSamples;
186     fGeom              = rawU.fGeom;
187     fOption            = rawU.fOption;
188     fRemoveBadChannels = rawU.fRemoveBadChannels;
189     fFittingAlgorithm  = rawU.fFittingAlgorithm;
190     fUseFALTRO         = rawU.fUseFALTRO;
191     fRawAnalyzer       = rawU.fRawAnalyzer;
192     fMapping[0]        = rawU.fMapping[0];
193     fMapping[1]        = rawU.fMapping[1];
194     fMapping[2]        = rawU.fMapping[2];
195     fMapping[3]        = rawU.fMapping[3];
196   }
197
198   return *this;
199
200 }
201
202 //____________________________________________________________________________
203 AliEMCALRawUtils::~AliEMCALRawUtils() {
204   //dtor
205
206 }
207
208 //____________________________________________________________________________
209 void AliEMCALRawUtils::Digits2Raw()
210 {
211   // convert digits of the current event to raw data
212   
213   AliRunLoader *rl = AliRunLoader::Instance();
214   AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
215
216   // get the digits
217   loader->LoadDigits("EMCAL");
218   loader->GetEvent();
219   TClonesArray* digits = loader->Digits() ;
220   
221   if (!digits) {
222     Warning("Digits2Raw", "no digits found !");
223     return;
224   }
225
226   static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
227   AliAltroBuffer* buffers[nDDL];
228   for (Int_t i=0; i < nDDL; i++)
229     buffers[i] = 0;
230
231   TArrayI adcValuesLow(fgTimeBins);
232   TArrayI adcValuesHigh(fgTimeBins);
233
234   // loop over digits (assume ordered digits)
235   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
236     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
237     if (digit->GetAmp() < fgThreshold) 
238       continue;
239
240     //get cell indices
241     Int_t nSM = 0;
242     Int_t nIphi = 0;
243     Int_t nIeta = 0;
244     Int_t iphi = 0;
245     Int_t ieta = 0;
246     Int_t nModule = 0;
247     fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
248     fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
249     
250     //Check which is the RCU, 0 or 1, of the cell.
251     Int_t iRCU = -111;
252     //RCU0
253     if (0<=iphi&&iphi<8) iRCU=0; // first cable row
254     else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half; 
255     //second cable row
256     //RCU1
257     else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half; 
258     //second cable row
259     else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
260
261     if (nSM%2==1) iRCU = 1 - iRCU; // swap for odd=C side, to allow us to cable both sides the same
262
263     if (iRCU<0) 
264       Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
265     
266     //Which DDL?
267     Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
268     if (iDDL >= nDDL)
269       Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
270
271     if (buffers[iDDL] == 0) {      
272       // open new file and write dummy header
273       TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
274       //Select mapping file RCU0A, RCU0C, RCU1A, RCU1C
275       Int_t iRCUside=iRCU+(nSM%2)*2;
276       //iRCU=0 and even (0) SM -> RCU0A.data   0
277       //iRCU=1 and even (0) SM -> RCU1A.data   1
278       //iRCU=0 and odd  (1) SM -> RCU0C.data   2
279       //iRCU=1 and odd  (1) SM -> RCU1C.data   3
280       //cout<<" nSM "<<nSM<<"; iRCU "<<iRCU<<"; iRCUside "<<iRCUside<<endl;
281       buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCUside]);
282       buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
283     }
284     
285     // out of time range signal (?)
286     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
287       AliInfo("Signal is out of time range.\n");
288       buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
289       buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() );  // time bin
290       buffers[iDDL]->FillBuffer(3);          // bunch length      
291       buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM);  // trailer
292       // calculate the time response function
293     } else {
294       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh.GetArray(), adcValuesLow.GetArray()) ; 
295       if (lowgain) 
296         buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow.GetArray(), fgThreshold);
297       else 
298         buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh.GetArray(), fgThreshold);
299     }
300   }
301   
302   // write headers and close files
303   for (Int_t i=0; i < nDDL; i++) {
304     if (buffers[i]) {
305       buffers[i]->Flush();
306       buffers[i]->WriteDataHeader(kFALSE, kFALSE);
307       delete buffers[i];
308     }
309   }
310
311   loader->UnloadDigits();
312 }
313
314 //____________________________________________________________________________
315 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr, const AliCaloCalibPedestal* pedbadmap, TClonesArray *digitsTRG)
316 {
317   // convert raw data of the current event to digits                                                                                     
318
319   digitsArr->Clear(); 
320
321   if (!digitsArr) {
322     Error("Raw2Digits", "no digits found !");
323     return;
324   }
325   if (!reader) {
326     Error("Raw2Digits", "no raw reader found !");
327     return;
328   }
329
330   AliCaloRawStreamV3 in(reader,"EMCAL",fMapping);
331   // Select EMCAL DDL's;
332   reader->Select("EMCAL",0,43); // 43 = AliEMCALGeoParams::fgkLastAltroDDL
333
334   // fRawAnalyzer setup
335   fRawAnalyzer->SetAmpCut(fNoiseThreshold);
336   fRawAnalyzer->SetFitArrayCut(fNoiseThreshold);
337   fRawAnalyzer->SetIsZeroSuppressed(true); // TMP - should use stream->IsZeroSuppressed(), or altro cfg registers later
338
339   // channel info parameters
340   Int_t lowGain = 0;
341   Int_t caloFlag = 0; // low, high gain, or TRU, or LED ref.
342
343   // start loop over input stream 
344   while (in.NextDDL()) {
345           
346 //    if ( in.GetDDLNumber() != 0 && in.GetDDLNumber() != 2 ) continue;
347
348     while (in.NextChannel()) {
349
350 /*
351           Int_t    hhwAdd    = in.GetHWAddress();
352           UShort_t iiBranch  = ( hhwAdd >> 11 ) & 0x1; // 0/1
353           UShort_t iiFEC     = ( hhwAdd >>  7 ) & 0xF;
354           UShort_t iiChip    = ( hhwAdd >>  4 ) & 0x7;
355           UShort_t iiChannel =   hhwAdd         & 0xF;
356                  
357           if ( !( iiBranch == 0 && iiFEC == 1 && iiChip == 3 && ( iiChannel >= 8 && iiChannel <= 15 ) ) && !( iiBranch == 1 && iiFEC == 0 && in.GetColumn() == 0 ) ) continue;
358 */
359                 
360       //Check if the signal  is high or low gain and then do the fit, 
361       //if it  is from TRU or LEDMon do not fit
362       caloFlag = in.GetCaloFlag();
363 //              if (caloFlag != 0 && caloFlag != 1) continue; 
364           if (caloFlag > 2) continue; // Work with ALTRO and FALTRO 
365                 
366       //Do not fit bad channels of ALTRO
367       if(caloFlag < 2 && fRemoveBadChannels && pedbadmap->IsBadChannel(in.GetModule(),in.GetColumn(),in.GetRow())) {
368         //printf("Tower from SM %d, column %d, row %d is BAD!!! Skip \n", in.GetModule(),in.GetColumn(),in.GetRow());
369         continue;
370       }  
371
372       vector<AliCaloBunchInfo> bunchlist; 
373       while (in.NextBunch()) {
374         bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
375       } // loop over bunches
376
377    
378     if ( caloFlag < 2 ){ // ALTRO
379                 
380           Float_t time = 0; 
381           Float_t amp = 0; 
382                 
383       if ( fFittingAlgorithm == kFastFit || fFittingAlgorithm == kNeuralNet || fFittingAlgorithm == kLMS || fFittingAlgorithm == kPeakFinder || fFittingAlgorithm == kCrude) {
384         // all functionality to determine amp and time etc is encapsulated inside the Evaluate call for these methods 
385         AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
386
387         amp = fitResults.GetAmp();
388         time = fitResults.GetTof();
389       }
390       else { // for the other methods we for now use the functionality of 
391         // AliCaloRawAnalyzer as well, to select samples and prepare for fits, 
392         // if it looks like there is something to fit
393
394         // parameters init.
395         Float_t ampEstimate  = 0;
396         short maxADC = 0;
397         short timeEstimate = 0;
398         Float_t pedEstimate = 0;
399         Int_t first = 0;
400         Int_t last = 0;
401         Int_t bunchIndex = 0;
402         //
403         // The PreFitEvaluateSamples + later call to FitRaw will hopefully 
404         // be replaced by a single Evaluate call or so soon, like for the other
405         // methods, but this should be good enough for evaluation of 
406         // the methods for now (Jan. 2010)
407         //
408         int nsamples = fRawAnalyzer->PreFitEvaluateSamples( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2(), bunchIndex, ampEstimate, maxADC, timeEstimate, pedEstimate, first, last); 
409         
410         if (ampEstimate > fNoiseThreshold) { // something worth looking at
411           
412           time = timeEstimate; 
413           amp = ampEstimate; 
414           
415           if ( nsamples > 1 ) { // possibly something to fit
416             FitRaw(first, last, amp, time);
417           }
418           
419           if ( amp>0 && time>0 ) { // brief sanity check of fit results
420             
421             // check fit results: should be consistent with initial estimates
422             // more magic numbers, but very loose cuts, for now..
423             // We have checked that amp and ampEstimate values are positive so division for assymmetry
424             // calculation should be OK/safe
425             Float_t ampAsymm = (amp - ampEstimate)/(amp + ampEstimate);
426             if ( (TMath::Abs(ampAsymm) > 0.1) ) {
427               AliDebug(2,Form("Fit results amp %f time %f not consistent with expectations ped %f max-ped %f time %d",
428                               amp, time, pedEstimate, ampEstimate, timeEstimate));
429               
430               // what should do we do then? skip this channel or assign the simple estimate? 
431               // for now just overwrite the fit results with the simple estimate
432               amp = ampEstimate;
433               time = timeEstimate; 
434             } // asymm check
435           } // amp & time check
436         } // ampEstimate check
437       } // method selection
438     
439       if (amp > fNoiseThreshold) { // something to be stored
440         Int_t id =  fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
441         lowGain = in.IsLowGain();
442
443         // go from time-bin units to physical time fgtimetrigger
444         time = time * GetRawFormatTimeBinWidth(); // skip subtraction of fgTimeTrigger?
445
446         AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
447         // printf("Added tower: SM %d, row %d, column %d, amp %3.2f\n",in.GetModule(), in.GetRow(), in.GetColumn(),amp);
448         // round off amplitude value to nearest integer
449         AddDigit(digitsArr, id, lowGain, TMath::Nint(amp), time); 
450       }
451       
452         }//ALTRO
453         else if(fUseFALTRO)
454         {// Fake ALTRO
455                 //              if (maxTimeBin && gSig->GetN() > maxTimeBin + 10) gSig->Set(maxTimeBin + 10); // set actual max size of TGraph
456                 Int_t    hwAdd    = in.GetHWAddress();
457                 UShort_t iRCU     = in.GetDDLNumber() % 2; // 0/1
458                 UShort_t iBranch  = ( hwAdd >> 11 ) & 0x1; // 0/1
459                 
460                 // Now find TRU number
461                 Int_t itru = 3 * in.GetModule() + ( (iRCU << 1) | iBranch ) - 1;
462                 
463                 AliDebug(1,Form("Found TRG digit in TRU: %2d ADC: %2d",itru,in.GetColumn()));
464                 
465                 Int_t idtrg;
466                 
467                 Bool_t isOK = fGeom->GetAbsFastORIndexFromTRU(itru, in.GetColumn(), idtrg);
468                 
469                 Int_t timeSamples[256]; for (Int_t j=0;j<256;j++) timeSamples[j] = 0;
470                 Int_t nSamples = 0;
471                 
472                 for (std::vector<AliCaloBunchInfo>::iterator itVectorData = bunchlist.begin(); itVectorData != bunchlist.end(); itVectorData++)
473                 {
474                         AliCaloBunchInfo bunch = *(itVectorData);
475                         
476                         const UShort_t* sig = bunch.GetData();
477                         Int_t startBin = bunch.GetStartBin();
478                         
479                         for (Int_t iS = 0; iS < bunch.GetLength(); iS++) 
480                         {
481                                 Int_t time = startBin--;
482                                 Int_t amp  = sig[iS];
483                                 
484                                 if ( amp ) timeSamples[nSamples++] = ( ( time << 12 ) & 0xFF000 ) | ( amp & 0xFFF );
485                         }
486                 }
487                 
488                 if (nSamples && isOK) AddDigit(digitsTRG, idtrg, timeSamples, nSamples);
489         }//Fake ALTRO
490    } // end while over channel   
491   } //end while over DDL's, of input stream 
492
493   return ; 
494 }
495
496 //____________________________________________________________________________ 
497 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t timeSamples[], Int_t nSamples) 
498 {
499         new((*digitsArr)[digitsArr->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);     
500         
501         //      Int_t idx = digitsArr->GetEntriesFast()-1;
502         //      AliEMCALRawDigit* d = (AliEMCALRawDigit*)digitsArr->At(idx);
503 }
504
505 //____________________________________________________________________________ 
506 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
507   //
508   // Add a new digit. 
509   // This routine checks whether a digit exists already for this tower 
510   // and then decides whether to use the high or low gain info
511   //
512   // Called by Raw2Digits
513   
514   AliEMCALDigit *digit = 0, *tmpdigit = 0;
515   TIter nextdigit(digitsArr);
516   while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
517     if (tmpdigit->GetId() == id)
518       digit = tmpdigit;
519   }
520
521   if (!digit) { // no digit existed for this tower; create one
522     if (lowGain && amp > fgkOverflowCut) 
523       amp = Int_t(fHighLowGainFactor * amp); 
524     Int_t idigit = digitsArr->GetEntries();
525     new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;   
526   }
527   else { // a digit already exists, check range 
528          // (use high gain if signal < cut value, otherwise low gain)
529     if (lowGain) { // new digit is low gain
530       if (digit->GetAmp() > fgkOverflowCut) {  // use if stored digit is out of range
531         digit->SetAmp(Int_t(fHighLowGainFactor * amp));
532         digit->SetTime(time);
533       }
534     }
535     else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
536       digit->SetAmp(amp);
537       digit->SetTime(time);
538     }
539   }
540 }
541
542 //____________________________________________________________________________ 
543 void AliEMCALRawUtils::FitRaw(const Int_t firstTimeBin, const Int_t lastTimeBin, Float_t & amp, Float_t & time) const 
544 { // Fits the raw signal time distribution
545   
546   //--------------------------------------------------
547   //Do the fit, different fitting algorithms available
548   //--------------------------------------------------
549   int nsamples = lastTimeBin - firstTimeBin + 1;
550
551   switch(fFittingAlgorithm) {
552   case kStandard:
553     {
554       if (nsamples < 3) { return; } // nothing much to fit
555       //printf("Standard fitter \n");
556
557       // Create Graph to hold data we will fit 
558       TGraph *gSig =  new TGraph( nsamples); 
559       for (int i=0; i<nsamples; i++) {
560         Int_t timebin = firstTimeBin + i;    
561         gSig->SetPoint(timebin, timebin, fRawAnalyzer->GetReversed(timebin)); 
562       }
563
564       TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
565       signalF->SetParameters(10.,5.,fTau,fOrder,0.); //set all defaults once, just to be safe
566       signalF->SetParNames("amp","t0","tau","N","ped");
567       signalF->FixParameter(2,fTau); // tau in units of time bin
568       signalF->FixParameter(3,fOrder); // order
569       signalF->FixParameter(4, 0); // pedestal should be subtracted when we get here 
570       signalF->SetParameter(1, time);
571       signalF->SetParameter(0, amp);
572                                 
573       gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
574                                 
575       // assign fit results
576       amp = signalF->GetParameter(0); 
577       time = signalF->GetParameter(1);
578
579       delete signalF;
580
581       // cross-check with ParabolaFit to see if the results make sense
582       FitParabola(gSig, amp); // amp is possibly updated
583
584       //printf("Std   : Amp %f, time %g\n",amp, time);
585       delete gSig; // delete TGraph
586                                 
587       break;
588     }//kStandard Fitter
589     //----------------------------
590   case kLogFit:
591     {
592       if (nsamples < 3) { return; } // nothing much to fit
593       //printf("LogFit \n");
594
595       // Create Graph to hold data we will fit 
596       TGraph *gSigLog =  new TGraph( nsamples); 
597       for (int i=0; i<nsamples; i++) {
598         Int_t timebin = firstTimeBin + i;    
599         gSigLog->SetPoint(timebin, timebin, TMath::Log(fRawAnalyzer->GetReversed(timebin) ) ); 
600       }
601
602       TF1 * signalFLog = new TF1("signalLog", RawResponseFunctionLog, 0, GetRawFormatTimeBins(), 5);
603       signalFLog->SetParameters(2.3, 5.,fTau,fOrder,0.); //set all defaults once, just to be safe
604       signalFLog->SetParNames("amplog","t0","tau","N","ped");
605       signalFLog->FixParameter(2,fTau); // tau in units of time bin
606       signalFLog->FixParameter(3,fOrder); // order
607       signalFLog->FixParameter(4, 0); // pedestal should be subtracted when we get here 
608       signalFLog->SetParameter(1, time);
609       if (amp>=1) {
610         signalFLog->SetParameter(0, TMath::Log(amp));
611       }
612         
613       gSigLog->Fit(signalFLog, "QROW"); // Note option 'W': equal errors on all points
614                                 
615       // assign fit results
616       Double_t amplog = signalFLog->GetParameter(0); //Not Amp, but Log of Amp
617       amp = TMath::Exp(amplog);
618       time = signalFLog->GetParameter(1);
619
620       delete signalFLog;
621       //printf("LogFit: Amp %f, time %g\n",amp, time);
622       delete gSigLog; 
623       break;
624     } //kLogFit 
625     //----------------------------      
626     
627     //----------------------------
628   }//switch fitting algorithms
629
630   return;
631 }
632
633 //__________________________________________________________________
634 void AliEMCALRawUtils::FitParabola(const TGraph *gSig, Float_t & amp) const 
635 {
636   //BEG YS alternative methods to calculate the amplitude
637   Double_t * ymx = gSig->GetX() ; 
638   Double_t * ymy = gSig->GetY() ; 
639   const Int_t kN = 3 ; 
640   Double_t ymMaxX[kN] = {0., 0., 0.} ; 
641   Double_t ymMaxY[kN] = {0., 0., 0.} ; 
642   Double_t ymax = 0. ; 
643   // find the maximum amplitude
644   Int_t ymiMax = 0 ;  
645   for (Int_t ymi = 0; ymi < gSig->GetN(); ymi++) {
646     if (ymy[ymi] > ymMaxY[0] ) {
647       ymMaxY[0] = ymy[ymi] ; //<========== This is the maximum amplitude
648       ymMaxX[0] = ymx[ymi] ;
649       ymiMax = ymi ; 
650     }
651   }
652   // find the maximum by fitting a parabola through the max and the two adjacent samples
653   if ( ymiMax < gSig->GetN()-1 && ymiMax > 0) {
654     ymMaxY[1] = ymy[ymiMax+1] ;
655     ymMaxY[2] = ymy[ymiMax-1] ; 
656     ymMaxX[1] = ymx[ymiMax+1] ;
657     ymMaxX[2] = ymx[ymiMax-1] ; 
658     if (ymMaxY[0]*ymMaxY[1]*ymMaxY[2] > 0) {
659       //fit a parabola through the 3 points y= a+bx+x*x*x
660       Double_t sy = 0 ; 
661       Double_t sx = 0 ; 
662       Double_t sx2 = 0 ; 
663       Double_t sx3 = 0 ; 
664       Double_t sx4 = 0 ; 
665       Double_t sxy = 0 ; 
666       Double_t sx2y = 0 ; 
667       for (Int_t i = 0; i < kN ; i++) {
668         sy += ymMaxY[i] ; 
669         sx += ymMaxX[i] ;               
670         sx2 += ymMaxX[i]*ymMaxX[i] ; 
671         sx3 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
672         sx4 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
673         sxy += ymMaxX[i]*ymMaxY[i] ; 
674         sx2y += ymMaxX[i]*ymMaxX[i]*ymMaxY[i] ; 
675       }
676       Double_t cN = (sx2y*kN-sy*sx2)*(sx3*sx-sx2*sx2)-(sx2y*sx-sxy*sx2)*(sx3*kN-sx*sx2); 
677       Double_t cD = (sx4*kN-sx2*sx2)*(sx3*sx-sx2*sx2)-(sx4*sx-sx3*sx2)*(sx3*kN-sx*sx2) ;
678       Double_t c  = cN / cD ; 
679       Double_t b  = ((sx2y*kN-sy*sx2)-c*(sx4*kN-sx2*sx2))/(sx3*kN-sx*sx2) ;
680       Double_t a  = (sy-b*sx-c*sx2)/kN  ;
681       Double_t xmax = -b/(2*c) ; 
682       ymax = a + b*xmax + c*xmax*xmax ;//<========== This is the maximum amplitude
683     }
684   }
685   
686   Double_t diff = TMath::Abs(1-ymMaxY[0]/amp) ; 
687   if (diff > 0.1) 
688     amp = ymMaxY[0] ; 
689   //printf("Yves   : Amp %f, time %g\n",amp, time);
690   //END YS
691   return;
692 }
693
694 //__________________________________________________________________
695 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
696 {
697   // Matches version used in 2007 beam test
698   //
699   // Shape of the electronics raw reponse:
700   // It is a semi-gaussian, 2nd order Gamma function of the general form
701   //
702   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
703   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
704   // F = 0                                   for xx < 0 
705   //
706   // parameters:
707   // A:   par[0]   // Amplitude = peak value
708   // t0:  par[1]
709   // tau: par[2]
710   // N:   par[3]
711   // ped: par[4]
712   //
713   Double_t signal ;
714   Double_t tau =par[2];
715   Double_t n =par[3];
716   Double_t ped = par[4];
717   Double_t xx = ( x[0] - par[1] + tau ) / tau ;
718
719   if (xx <= 0) 
720     signal = ped ;  
721   else {  
722     signal = ped + par[0] * TMath::Power(xx , n) * TMath::Exp(n * (1 - xx )) ; 
723   }
724   return signal ;  
725 }
726
727 //__________________________________________________________________
728 Double_t AliEMCALRawUtils::RawResponseFunctionLog(Double_t *x, Double_t *par)
729 {
730   // Matches version used in 2007 beam test
731   //
732   // Shape of the electronics raw reponse:
733   // It is a semi-gaussian, 2nd order Gamma function of the general form
734   //
735   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
736   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
737   // F = 0                                   for xx < 0 
738   //
739   // parameters:
740   // Log[A]:   par[0]   // Amplitude = peak value
741   // t0:  par[1]
742   // tau: par[2]
743   // N:   par[3]
744   // ped: par[4]
745   //
746   Double_t signal ;
747   Double_t tau =par[2];
748   Double_t n =par[3];
749   //Double_t ped = par[4]; // not used
750   Double_t xx = ( x[0] - par[1] + tau ) / tau ;
751
752   if (xx < 0) 
753     signal = par[0] - n*TMath::Log(TMath::Abs(xx)) + n * (1 - xx ) ;  
754   else {  
755     signal = par[0] + n*TMath::Log(xx) + n * (1 - xx ) ; 
756   }
757   return signal ;  
758 }
759
760 //__________________________________________________________________
761 Bool_t AliEMCALRawUtils::RawSampledResponse(const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL, const Int_t keyErr) const 
762 {
763   // for a start time dtime and an amplitude damp given by digit, 
764   // calculates the raw sampled response AliEMCAL::RawResponseFunction
765
766   Bool_t lowGain = kFALSE ; 
767
768   // A:   par[0]   // Amplitude = peak value
769   // t0:  par[1]                            
770   // tau: par[2]                            
771   // N:   par[3]                            
772   // ped: par[4]
773
774   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
775   signalF.SetParameter(0, damp) ; 
776   signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ; 
777   signalF.SetParameter(2, fTau) ; 
778   signalF.SetParameter(3, fOrder);
779   signalF.SetParameter(4, fgPedestalValue);
780         
781   Double_t signal=0.0, noise=0.0;
782   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
783     signal = signalF.Eval(iTime) ;  
784         
785     // Next lines commeted for the moment but in principle it is not necessary to add
786     // extra noise since noise already added at the digits level.       
787
788     //According to Terry Awes, 13-Apr-2008
789     //add gaussian noise in quadrature to each sample
790     //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
791     //signal = sqrt(signal*signal + noise*noise);
792
793     // March 17,09 for fast fit simulations by Alexei Pavlinov.
794     // Get from PHOS analysis. In some sense it is open questions.
795         if(keyErr>0) {
796                 noise = gRandom->Gaus(0.,fgFEENoise);
797                 signal += noise; 
798         }
799           
800     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
801     if ( adcH[iTime] > fgkRawSignalOverflow ){  // larger than 10 bits 
802       adcH[iTime] = fgkRawSignalOverflow ;
803       lowGain = kTRUE ; 
804     }
805
806     signal /= fHighLowGainFactor;
807
808     adcL[iTime] =  static_cast<Int_t>(signal + 0.5) ;
809     if ( adcL[iTime] > fgkRawSignalOverflow)  // larger than 10 bits 
810       adcL[iTime] = fgkRawSignalOverflow ;
811   }
812   return lowGain ; 
813 }
814
815 //__________________________________________________________________
816 void AliEMCALRawUtils::SetFittingAlgorithm(Int_t fitAlgo)              
817 {
818         //Set fitting algorithm and initialize it if this same algorithm was not set before.
819         //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
820
821         if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
822                 //Do nothing, this same algorithm already set before.
823                 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
824                 return;
825         }
826         //Initialize the requested algorithm
827         if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
828                 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
829                 
830                 fFittingAlgorithm = fitAlgo; 
831                 if (fRawAnalyzer) delete fRawAnalyzer;  // delete prev. analyzer if existed.
832                 
833                 if (fitAlgo == kFastFit) {
834                         fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
835                 }
836                 else if (fitAlgo == kNeuralNet) {
837                         fRawAnalyzer = new AliCaloRawAnalyzerNN();
838                 }
839                 else if (fitAlgo == kLMS) {
840                         fRawAnalyzer = new AliCaloRawAnalyzerLMS();
841                 }
842                 else if (fitAlgo == kPeakFinder) {
843                         fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
844                 }
845                 else if (fitAlgo == kCrude) {
846                         fRawAnalyzer = new AliCaloRawAnalyzerCrude();
847                 }
848                 else {
849                         fRawAnalyzer = new AliCaloRawAnalyzer();
850                 }
851         }
852         
853 }
854
855