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