]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRawUtils.cxx
raw ana t0 calc. updates
[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 "AliEMCAL.h"
51 #include "AliCaloCalibPedestal.h"  
52 #include "AliCaloFastAltroFitv0.h"
53 #include "AliCaloNeuralFit.h"
54 #include "AliCaloBunchInfo.h"
55 #include "AliCaloFitResults.h"
56 #include "AliCaloRawAnalyzerFastFit.h"
57 #include "AliCaloRawAnalyzerNN.h"
58 #include "AliCaloRawAnalyzerLMS.h"
59 #include "AliCaloRawAnalyzerPeakFinder.h"
60 #include "AliCaloRawAnalyzerCrude.h"
61
62 ClassImp(AliEMCALRawUtils)
63   
64 // Signal shape parameters
65 Int_t    AliEMCALRawUtils::fgTimeBins = 256; // number of sampling bins of the raw RO signal (we typically use 15-50; theoretical max is 1k+) 
66 Double_t AliEMCALRawUtils::fgTimeBinWidth  = 100E-9 ; // each sample is 100 ns
67 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ;   // 15 time bins ~ 1.5 musec
68
69 // some digitization constants
70 Int_t    AliEMCALRawUtils::fgThreshold = 1;
71 Int_t    AliEMCALRawUtils::fgDDLPerSuperModule = 2;  // 2 ddls per SuperModule
72 Int_t    AliEMCALRawUtils::fgPedestalValue = 32;     // pedestal value for digits2raw
73 Double_t AliEMCALRawUtils::fgFEENoise = 3.;          // 3 ADC channels of noise (sampled)
74
75 AliEMCALRawUtils::AliEMCALRawUtils(fitAlgorithm fitAlgo)
76   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
77     fNPedSamples(0), fGeom(0), fOption(""),
78     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),fRawAnalyzer(0)
79 {
80
81   //These are default parameters.  
82   //Can be re-set from without with setter functions
83   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
84   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
85   fOrder = 2;                         // order of gamma fn
86   fTau = 2.35;                        // in units of timebin, from CERN 2007 testbeam
87   fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
88   fNPedSamples = 4;    // less than this value => likely pedestal samples
89   fRemoveBadChannels = kTRUE; //Remove bad channels before fitting
90   SetFittingAlgorithm(fitAlgo);
91
92   //Get Mapping RCU files from the AliEMCALRecParam                                 
93   const TObjArray* maps = AliEMCALRecParam::GetMappings();
94   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
95
96   for(Int_t i = 0; i < 4; i++) {
97     fMapping[i] = (AliAltroMapping*)maps->At(i);
98   }
99
100   //To make sure we match with the geometry in a simulation file,
101   //let's try to get it first.  If not, take the default geometry
102   AliRunLoader *rl = AliRunLoader::Instance();
103   if(!rl) AliError("Cannot find RunLoader!");
104   if (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)
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     while (in.NextChannel()) {
343
344       //Check if the signal  is high or low gain and then do the fit, 
345       //if it  is from TRU or LEDMon do not fit
346       caloFlag = in.GetCaloFlag();
347       if (caloFlag != 0 && caloFlag != 1) continue; 
348               
349       //Do not fit bad channels
350       if(fRemoveBadChannels && pedbadmap->IsBadChannel(in.GetModule(),in.GetColumn(),in.GetRow())) {
351         //printf("Tower from SM %d, column %d, row %d is BAD!!! Skip \n", in.GetModule(),in.GetColumn(),in.GetRow());
352         continue;
353       }  
354
355       vector<AliCaloBunchInfo> bunchlist; 
356       while (in.NextBunch()) {
357         bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
358       } // loop over bunches
359
360       Float_t time = 0; 
361       Float_t amp = 0; 
362
363       if ( fFittingAlgorithm == kFastFit || fFittingAlgorithm == kNeuralNet || fFittingAlgorithm == kLMS || fFittingAlgorithm == kPeakFinder || fFittingAlgorithm == kCrude) {
364         // all functionality to determine amp and time etc is encapsulated inside the Evaluate call for these methods 
365         AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
366
367         amp = fitResults.GetAmp();
368         time = fitResults.GetTof();
369       }
370       else { // for the other methods we for now use the functionality of 
371         // AliCaloRawAnalyzer as well, to select samples and prepare for fits, 
372         // if it looks like there is something to fit
373
374         // parameters init.
375         Float_t ampEstimate  = 0;
376         short maxADC = 0;
377         short timeEstimate = 0;
378         Float_t pedEstimate = 0;
379         Int_t first = 0;
380         Int_t last = 0;
381         Int_t bunchIndex = 0;
382         //
383         // The PreFitEvaluateSamples + later call to FitRaw will hopefully 
384         // be replaced by a single Evaluate call or so soon, like for the other
385         // methods, but this should be good enough for evaluation of 
386         // the methods for now (Jan. 2010)
387         //
388         int nsamples = fRawAnalyzer->PreFitEvaluateSamples( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2(), bunchIndex, ampEstimate, maxADC, timeEstimate, pedEstimate, first, last); 
389         
390         if (ampEstimate > fNoiseThreshold) { // something worth looking at
391           
392           time = timeEstimate; 
393           amp = ampEstimate; 
394           
395           if ( nsamples > 1 ) { // possibly something to fit
396             FitRaw(first, last, amp, time);
397           }
398           
399           if ( amp>0 && time>0 ) { // brief sanity check of fit results
400             
401             // check fit results: should be consistent with initial estimates
402             // more magic numbers, but very loose cuts, for now..
403             // We have checked that amp and ampEstimate values are positive so division for assymmetry
404             // calculation should be OK/safe
405             Float_t ampAsymm = (amp - ampEstimate)/(amp + ampEstimate);
406             if ( (TMath::Abs(ampAsymm) > 0.1) ) {
407               AliDebug(2,Form("Fit results amp %f time %f not consistent with expectations ped %f max-ped %f time %d",
408                               amp, time, pedEstimate, ampEstimate, timeEstimate));
409               
410               // what should do we do then? skip this channel or assign the simple estimate? 
411               // for now just overwrite the fit results with the simple estimate
412               amp = ampEstimate;
413               time = timeEstimate; 
414             } // asymm check
415           } // amp & time check
416         } // ampEstimate check
417       } // method selection
418     
419       if (amp > fNoiseThreshold) { // something to be stored
420         Int_t id =  fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
421         lowGain = in.IsLowGain();
422
423         // go from time-bin units to physical time fgtimetrigger
424         time = time * GetRawFormatTimeBinWidth(); // skip subtraction of fgTimeTrigger?
425
426         AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
427         // printf("Added tower: SM %d, row %d, column %d, amp %3.2f\n",in.GetModule(), in.GetRow(), in.GetColumn(),amp);
428         // round off amplitude value to nearest integer
429         AddDigit(digitsArr, id, lowGain, TMath::Nint(amp), time); 
430       }
431       
432    } // end while over channel   
433   } //end while over DDL's, of input stream 
434
435   return ; 
436 }
437
438 //____________________________________________________________________________ 
439 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
440   //
441   // Add a new digit. 
442   // This routine checks whether a digit exists already for this tower 
443   // and then decides whether to use the high or low gain info
444   //
445   // Called by Raw2Digits
446   
447   AliEMCALDigit *digit = 0, *tmpdigit = 0;
448   TIter nextdigit(digitsArr);
449   while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
450     if (tmpdigit->GetId() == id)
451       digit = tmpdigit;
452   }
453
454   if (!digit) { // no digit existed for this tower; create one
455     if (lowGain && amp > fgkOverflowCut) 
456       amp = Int_t(fHighLowGainFactor * amp); 
457     Int_t idigit = digitsArr->GetEntries();
458     new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;   
459   }
460   else { // a digit already exists, check range 
461          // (use high gain if signal < cut value, otherwise low gain)
462     if (lowGain) { // new digit is low gain
463       if (digit->GetAmp() > fgkOverflowCut) {  // use if stored digit is out of range
464         digit->SetAmp(Int_t(fHighLowGainFactor * amp));
465         digit->SetTime(time);
466       }
467     }
468     else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
469       digit->SetAmp(amp);
470       digit->SetTime(time);
471     }
472   }
473 }
474
475 //____________________________________________________________________________ 
476 void AliEMCALRawUtils::FitRaw(const Int_t firstTimeBin, const Int_t lastTimeBin, Float_t & amp, Float_t & time) const 
477 { // Fits the raw signal time distribution
478   
479   //--------------------------------------------------
480   //Do the fit, different fitting algorithms available
481   //--------------------------------------------------
482   int nsamples = lastTimeBin - firstTimeBin + 1;
483
484   switch(fFittingAlgorithm) {
485   case kStandard:
486     {
487       if (nsamples < 3) { return; } // nothing much to fit
488       //printf("Standard fitter \n");
489
490       // Create Graph to hold data we will fit 
491       TGraph *gSig =  new TGraph( nsamples); 
492       for (int i=0; i<nsamples; i++) {
493         Int_t timebin = firstTimeBin + i;    
494         gSig->SetPoint(timebin, timebin, fRawAnalyzer->GetReversed(timebin)); 
495       }
496
497       TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
498       signalF->SetParameters(10.,5.,fTau,fOrder,0.); //set all defaults once, just to be safe
499       signalF->SetParNames("amp","t0","tau","N","ped");
500       signalF->FixParameter(2,fTau); // tau in units of time bin
501       signalF->FixParameter(3,fOrder); // order
502       signalF->FixParameter(4, 0); // pedestal should be subtracted when we get here 
503       signalF->SetParameter(1, time);
504       signalF->SetParameter(0, amp);
505                                 
506       gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
507                                 
508       // assign fit results
509       amp = signalF->GetParameter(0); 
510       time = signalF->GetParameter(1);
511
512       delete signalF;
513
514       // cross-check with ParabolaFit to see if the results make sense
515       FitParabola(gSig, amp); // amp is possibly updated
516
517       //printf("Std   : Amp %f, time %g\n",amp, time);
518       delete gSig; // delete TGraph
519                                 
520       break;
521     }//kStandard Fitter
522     //----------------------------
523   case kLogFit:
524     {
525       if (nsamples < 3) { return; } // nothing much to fit
526       //printf("LogFit \n");
527
528       // Create Graph to hold data we will fit 
529       TGraph *gSigLog =  new TGraph( nsamples); 
530       for (int i=0; i<nsamples; i++) {
531         Int_t timebin = firstTimeBin + i;    
532         gSigLog->SetPoint(timebin, timebin, TMath::Log(fRawAnalyzer->GetReversed(timebin) ) ); 
533       }
534
535       TF1 * signalFLog = new TF1("signalLog", RawResponseFunctionLog, 0, GetRawFormatTimeBins(), 5);
536       signalFLog->SetParameters(2.3, 5.,fTau,fOrder,0.); //set all defaults once, just to be safe
537       signalFLog->SetParNames("amplog","t0","tau","N","ped");
538       signalFLog->FixParameter(2,fTau); // tau in units of time bin
539       signalFLog->FixParameter(3,fOrder); // order
540       signalFLog->FixParameter(4, 0); // pedestal should be subtracted when we get here 
541       signalFLog->SetParameter(1, time);
542       if (amp>=1) {
543         signalFLog->SetParameter(0, TMath::Log(amp));
544       }
545         
546       gSigLog->Fit(signalFLog, "QROW"); // Note option 'W': equal errors on all points
547                                 
548       // assign fit results
549       Double_t amplog = signalFLog->GetParameter(0); //Not Amp, but Log of Amp
550       amp = TMath::Exp(amplog);
551       time = signalFLog->GetParameter(1);
552
553       delete signalFLog;
554       //printf("LogFit: Amp %f, time %g\n",amp, time);
555       delete gSigLog; 
556       break;
557     } //kLogFit 
558     //----------------------------      
559     
560     //----------------------------
561   }//switch fitting algorithms
562
563   return;
564 }
565
566 //__________________________________________________________________
567 void AliEMCALRawUtils::FitParabola(const TGraph *gSig, Float_t & amp) const 
568 {
569   //BEG YS alternative methods to calculate the amplitude
570   Double_t * ymx = gSig->GetX() ; 
571   Double_t * ymy = gSig->GetY() ; 
572   const Int_t kN = 3 ; 
573   Double_t ymMaxX[kN] = {0., 0., 0.} ; 
574   Double_t ymMaxY[kN] = {0., 0., 0.} ; 
575   Double_t ymax = 0. ; 
576   // find the maximum amplitude
577   Int_t ymiMax = 0 ;  
578   for (Int_t ymi = 0; ymi < gSig->GetN(); ymi++) {
579     if (ymy[ymi] > ymMaxY[0] ) {
580       ymMaxY[0] = ymy[ymi] ; //<========== This is the maximum amplitude
581       ymMaxX[0] = ymx[ymi] ;
582       ymiMax = ymi ; 
583     }
584   }
585   // find the maximum by fitting a parabola through the max and the two adjacent samples
586   if ( ymiMax < gSig->GetN()-1 && ymiMax > 0) {
587     ymMaxY[1] = ymy[ymiMax+1] ;
588     ymMaxY[2] = ymy[ymiMax-1] ; 
589     ymMaxX[1] = ymx[ymiMax+1] ;
590     ymMaxX[2] = ymx[ymiMax-1] ; 
591     if (ymMaxY[0]*ymMaxY[1]*ymMaxY[2] > 0) {
592       //fit a parabola through the 3 points y= a+bx+x*x*x
593       Double_t sy = 0 ; 
594       Double_t sx = 0 ; 
595       Double_t sx2 = 0 ; 
596       Double_t sx3 = 0 ; 
597       Double_t sx4 = 0 ; 
598       Double_t sxy = 0 ; 
599       Double_t sx2y = 0 ; 
600       for (Int_t i = 0; i < kN ; i++) {
601         sy += ymMaxY[i] ; 
602         sx += ymMaxX[i] ;               
603         sx2 += ymMaxX[i]*ymMaxX[i] ; 
604         sx3 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
605         sx4 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
606         sxy += ymMaxX[i]*ymMaxY[i] ; 
607         sx2y += ymMaxX[i]*ymMaxX[i]*ymMaxY[i] ; 
608       }
609       Double_t cN = (sx2y*kN-sy*sx2)*(sx3*sx-sx2*sx2)-(sx2y*sx-sxy*sx2)*(sx3*kN-sx*sx2); 
610       Double_t cD = (sx4*kN-sx2*sx2)*(sx3*sx-sx2*sx2)-(sx4*sx-sx3*sx2)*(sx3*kN-sx*sx2) ;
611       Double_t c  = cN / cD ; 
612       Double_t b  = ((sx2y*kN-sy*sx2)-c*(sx4*kN-sx2*sx2))/(sx3*kN-sx*sx2) ;
613       Double_t a  = (sy-b*sx-c*sx2)/kN  ;
614       Double_t xmax = -b/(2*c) ; 
615       ymax = a + b*xmax + c*xmax*xmax ;//<========== This is the maximum amplitude
616     }
617   }
618   
619   Double_t diff = TMath::Abs(1-ymMaxY[0]/amp) ; 
620   if (diff > 0.1) 
621     amp = ymMaxY[0] ; 
622   //printf("Yves   : Amp %f, time %g\n",amp, time);
623   //END YS
624   return;
625 }
626
627 //__________________________________________________________________
628 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
629 {
630   // Matches version used in 2007 beam test
631   //
632   // Shape of the electronics raw reponse:
633   // It is a semi-gaussian, 2nd order Gamma function of the general form
634   //
635   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
636   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
637   // F = 0                                   for xx < 0 
638   //
639   // parameters:
640   // A:   par[0]   // Amplitude = peak value
641   // t0:  par[1]
642   // tau: par[2]
643   // N:   par[3]
644   // ped: par[4]
645   //
646   Double_t signal ;
647   Double_t tau =par[2];
648   Double_t n =par[3];
649   Double_t ped = par[4];
650   Double_t xx = ( x[0] - par[1] + tau ) / tau ;
651
652   if (xx <= 0) 
653     signal = ped ;  
654   else {  
655     signal = ped + par[0] * TMath::Power(xx , n) * TMath::Exp(n * (1 - xx )) ; 
656   }
657   return signal ;  
658 }
659
660 //__________________________________________________________________
661 Double_t AliEMCALRawUtils::RawResponseFunctionLog(Double_t *x, Double_t *par)
662 {
663   // Matches version used in 2007 beam test
664   //
665   // Shape of the electronics raw reponse:
666   // It is a semi-gaussian, 2nd order Gamma function of the general form
667   //
668   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
669   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
670   // F = 0                                   for xx < 0 
671   //
672   // parameters:
673   // Log[A]:   par[0]   // Amplitude = peak value
674   // t0:  par[1]
675   // tau: par[2]
676   // N:   par[3]
677   // ped: par[4]
678   //
679   Double_t signal ;
680   Double_t tau =par[2];
681   Double_t n =par[3];
682   //Double_t ped = par[4]; // not used
683   Double_t xx = ( x[0] - par[1] + tau ) / tau ;
684
685   if (xx < 0) 
686     signal = par[0] - n*TMath::Log(TMath::Abs(xx)) + n * (1 - xx ) ;  
687   else {  
688     signal = par[0] + n*TMath::Log(xx) + n * (1 - xx ) ; 
689   }
690   return signal ;  
691 }
692
693 //__________________________________________________________________
694 Bool_t AliEMCALRawUtils::RawSampledResponse(const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL, const Int_t keyErr) const 
695 {
696   // for a start time dtime and an amplitude damp given by digit, 
697   // calculates the raw sampled response AliEMCAL::RawResponseFunction
698
699   Bool_t lowGain = kFALSE ; 
700
701   // A:   par[0]   // Amplitude = peak value
702   // t0:  par[1]                            
703   // tau: par[2]                            
704   // N:   par[3]                            
705   // ped: par[4]
706
707   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
708   signalF.SetParameter(0, damp) ; 
709   signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ; 
710   signalF.SetParameter(2, fTau) ; 
711   signalF.SetParameter(3, fOrder);
712   signalF.SetParameter(4, fgPedestalValue);
713         
714   Double_t signal=0.0, noise=0.0;
715   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
716     signal = signalF.Eval(iTime) ;     
717
718     // Next lines commeted for the moment but in principle it is not necessary to add
719     // extra noise since noise already added at the digits level.       
720
721     //According to Terry Awes, 13-Apr-2008
722     //add gaussian noise in quadrature to each sample
723     //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
724     //signal = sqrt(signal*signal + noise*noise);
725
726     // March 17,09 for fast fit simulations by Alexei Pavlinov.
727     // Get from PHOS analysis. In some sense it is open questions.
728         if(keyErr>0) {
729                 noise = gRandom->Gaus(0.,fgFEENoise);
730                 signal += noise; 
731         }
732           
733     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
734     if ( adcH[iTime] > fgkRawSignalOverflow ){  // larger than 10 bits 
735       adcH[iTime] = fgkRawSignalOverflow ;
736       lowGain = kTRUE ; 
737     }
738
739     signal /= fHighLowGainFactor;
740
741     adcL[iTime] =  static_cast<Int_t>(signal + 0.5) ;
742     if ( adcL[iTime] > fgkRawSignalOverflow)  // larger than 10 bits 
743       adcL[iTime] = fgkRawSignalOverflow ;
744   }
745   return lowGain ; 
746 }
747
748 //__________________________________________________________________
749 void AliEMCALRawUtils::SetFittingAlgorithm(Int_t fitAlgo)              
750 {
751         //Set fitting algorithm and initialize it if this same algorithm was not set before.
752         
753         if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
754                 //Do nothing, this same algorithm already set before.
755                 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
756                 return;
757         }
758         //Initialize the requested algorithm
759         if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
760                 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
761                 
762                 fFittingAlgorithm = fitAlgo; 
763                 if (fRawAnalyzer) delete fRawAnalyzer;  // delete prev. analyzer if existed.
764                 
765                 if (fitAlgo == kFastFit) {
766                         fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
767                 }
768                 else if (fitAlgo == kNeuralNet) {
769                         fRawAnalyzer = new AliCaloRawAnalyzerNN();
770                 }
771                 else if (fitAlgo == kLMS) {
772                         fRawAnalyzer = new AliCaloRawAnalyzerLMS();
773                 }
774                 else if (fitAlgo == kPeakFinder) {
775                         fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
776                 }
777                 else if (fitAlgo == kCrude) {
778                         fRawAnalyzer = new AliCaloRawAnalyzerCrude();
779                 }
780                 else {
781                         fRawAnalyzer = new AliCaloRawAnalyzer();
782                 }
783         }
784         
785 }
786
787