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