]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRawUtils.cxx
77feaaac49d525d9e8a791fb0d88c567e6d36033
[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 class TSystem;
34   
35 class AliLog;
36 #include "AliRun.h"
37 #include "AliRunLoader.h"
38 class AliCaloAltroMapping;
39 #include "AliAltroBuffer.h"
40 #include "AliRawReader.h"
41 #include "AliCaloRawStreamV3.h"
42 #include "AliDAQ.h"
43   
44 #include "AliEMCALRecParam.h"
45 #include "AliEMCALLoader.h"
46 #include "AliEMCALGeometry.h"
47 class AliEMCALDigitizer;
48 #include "AliEMCALDigit.h"
49 #include "AliEMCAL.h"
50 #include "AliCaloCalibPedestal.h"  
51 #include "AliCaloFastAltroFitv0.h"
52
53 ClassImp(AliEMCALRawUtils)
54   
55 // Signal shape parameters
56 Int_t    AliEMCALRawUtils::fgTimeBins = 256; // number of sampling bins of the raw RO signal (we typically use 15-50; theoretical max is 1k+) 
57 Double_t AliEMCALRawUtils::fgTimeBinWidth  = 100E-9 ; // each sample is 100 ns
58 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ;   // 15 time bins ~ 1.5 musec
59
60 // some digitization constants
61 Int_t    AliEMCALRawUtils::fgThreshold = 1;
62 Int_t    AliEMCALRawUtils::fgDDLPerSuperModule = 2;  // 2 ddls per SuperModule
63 Int_t    AliEMCALRawUtils::fgPedestalValue = 32;     // pedestal value for digits2raw
64 Double_t AliEMCALRawUtils::fgFEENoise = 3.;          // 3 ADC channels of noise (sampled)
65
66 AliEMCALRawUtils::AliEMCALRawUtils()
67   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
68     fNPedSamples(0), fGeom(0), fOption(""),
69     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0)
70 {
71
72   //These are default parameters.  
73   //Can be re-set from without with setter functions
74   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
75   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
76   fOrder = 2;                         // order of gamma fn
77   fTau = 2.35;                        // in units of timebin, from CERN 2007 testbeam
78   fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
79   fNPedSamples = 4;    // less than this value => likely pedestal samples
80   fRemoveBadChannels = kTRUE; //Remove bad channels before fitting
81   fFittingAlgorithm  = kFastFit;//kStandard; // Use default minuit fitter
82         
83   //Get Mapping RCU files from the AliEMCALRecParam                                 
84   const TObjArray* maps = AliEMCALRecParam::GetMappings();
85   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
86
87   for(Int_t i = 0; i < 4; i++) {
88     fMapping[i] = (AliAltroMapping*)maps->At(i);
89   }
90
91   //To make sure we match with the geometry in a simulation file,
92   //let's try to get it first.  If not, take the default geometry
93   AliRunLoader *rl = AliRunLoader::Instance();
94   if(!rl) AliError("Cannot find RunLoader!");
95   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
96     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
97   } else {
98     AliInfo(Form("Using default geometry in raw reco"));
99     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
100   }
101
102   if(!fGeom) AliFatal(Form("Could not get geometry!"));
103
104 }
105
106 //____________________________________________________________________________
107 AliEMCALRawUtils::AliEMCALRawUtils(AliEMCALGeometry *pGeometry)
108   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
109     fNPedSamples(0), fGeom(pGeometry), fOption(""),
110     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0)
111 {
112   //
113   // Initialize with the given geometry - constructor required by HLT
114   // HLT does not use/support AliRunLoader(s) instances
115   // This is a minimum intervention solution
116   // Comment by MPloskon@lbl.gov
117   //
118
119   //These are default parameters. 
120   //Can be re-set from without with setter functions 
121   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
122   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits)
123   fOrder = 2;                         // order of gamma fn
124   fTau = 2.35;                        // in units of timebin, from CERN 2007 testbeam
125   fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
126   fNPedSamples = 4;    // less than this value => likely pedestal samples
127   fRemoveBadChannels = kTRUE; //Remove bad channels before fitting
128   fFittingAlgorithm  = kStandard; // Use default minuit fitter
129         
130   //Get Mapping RCU files from the AliEMCALRecParam
131   const TObjArray* maps = AliEMCALRecParam::GetMappings();
132   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
133
134   for(Int_t i = 0; i < 4; i++) {
135     fMapping[i] = (AliAltroMapping*)maps->At(i);
136   }
137
138   if(!fGeom) AliFatal(Form("Could not get geometry!"));
139
140 }
141
142 //____________________________________________________________________________
143 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
144   : TObject(),
145     fHighLowGainFactor(rawU.fHighLowGainFactor), 
146     fOrder(rawU.fOrder),
147     fTau(rawU.fTau),
148     fNoiseThreshold(rawU.fNoiseThreshold),
149     fNPedSamples(rawU.fNPedSamples),
150     fGeom(rawU.fGeom), 
151     fOption(rawU.fOption),
152     fRemoveBadChannels(rawU.fRemoveBadChannels),
153     fFittingAlgorithm(rawU.fFittingAlgorithm)
154 {
155   //copy ctor
156   fMapping[0] = rawU.fMapping[0];
157   fMapping[1] = rawU.fMapping[1];
158   fMapping[2] = rawU.fMapping[2];
159   fMapping[3] = rawU.fMapping[3];
160 }
161
162 //____________________________________________________________________________
163 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
164 {
165   //assignment operator
166
167   if(this != &rawU) {
168     fHighLowGainFactor = rawU.fHighLowGainFactor;
169     fOrder = rawU.fOrder;
170     fTau = rawU.fTau;
171     fNoiseThreshold = rawU.fNoiseThreshold;
172     fNPedSamples = rawU.fNPedSamples;
173     fGeom = rawU.fGeom;
174     fOption = rawU.fOption;
175     fRemoveBadChannels = rawU.fRemoveBadChannels;
176     fFittingAlgorithm  = rawU.fFittingAlgorithm;
177     fMapping[0] = rawU.fMapping[0];
178     fMapping[1] = rawU.fMapping[1];
179     fMapping[2] = rawU.fMapping[2];
180     fMapping[3] = rawU.fMapping[3];
181   }
182
183   return *this;
184
185 }
186
187 //____________________________________________________________________________
188 AliEMCALRawUtils::~AliEMCALRawUtils() {
189   //dtor
190
191 }
192
193 //____________________________________________________________________________
194 void AliEMCALRawUtils::Digits2Raw()
195 {
196   // convert digits of the current event to raw data
197   
198   AliRunLoader *rl = AliRunLoader::Instance();
199   AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
200
201   // get the digits
202   loader->LoadDigits("EMCAL");
203   loader->GetEvent();
204   TClonesArray* digits = loader->Digits() ;
205   
206   if (!digits) {
207     Warning("Digits2Raw", "no digits found !");
208     return;
209   }
210
211   static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
212   AliAltroBuffer* buffers[nDDL];
213   for (Int_t i=0; i < nDDL; i++)
214     buffers[i] = 0;
215
216   TArrayI adcValuesLow(fgTimeBins);
217   TArrayI adcValuesHigh(fgTimeBins);
218
219   // loop over digits (assume ordered digits)
220   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
221     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
222     if (digit->GetAmp() < fgThreshold) 
223       continue;
224
225     //get cell indices
226     Int_t nSM = 0;
227     Int_t nIphi = 0;
228     Int_t nIeta = 0;
229     Int_t iphi = 0;
230     Int_t ieta = 0;
231     Int_t nModule = 0;
232     fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
233     fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
234     
235     //Check which is the RCU, 0 or 1, of the cell.
236     Int_t iRCU = -111;
237     //RCU0
238     if (0<=iphi&&iphi<8) iRCU=0; // first cable row
239     else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half; 
240     //second cable row
241     //RCU1
242     else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half; 
243     //second cable row
244     else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
245
246     if (nSM%2==1) iRCU = 1 - iRCU; // swap for odd=C side, to allow us to cable both sides the same
247
248     if (iRCU<0) 
249       Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
250     
251     //Which DDL?
252     Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
253     if (iDDL >= nDDL)
254       Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
255
256     if (buffers[iDDL] == 0) {      
257       // open new file and write dummy header
258       TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
259       //Select mapping file RCU0A, RCU0C, RCU1A, RCU1C
260       Int_t iRCUside=iRCU+(nSM%2)*2;
261       //iRCU=0 and even (0) SM -> RCU0A.data   0
262       //iRCU=1 and even (0) SM -> RCU1A.data   1
263       //iRCU=0 and odd  (1) SM -> RCU0C.data   2
264       //iRCU=1 and odd  (1) SM -> RCU1C.data   3
265       //cout<<" nSM "<<nSM<<"; iRCU "<<iRCU<<"; iRCUside "<<iRCUside<<endl;
266       buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCUside]);
267       buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
268     }
269     
270     // out of time range signal (?)
271     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
272       AliInfo("Signal is out of time range.\n");
273       buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
274       buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() );  // time bin
275       buffers[iDDL]->FillBuffer(3);          // bunch length      
276       buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM);  // trailer
277       // calculate the time response function
278     } else {
279       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh.GetArray(), adcValuesLow.GetArray()) ; 
280       if (lowgain) 
281         buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow.GetArray(), fgThreshold);
282       else 
283         buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh.GetArray(), fgThreshold);
284     }
285   }
286   
287   // write headers and close files
288   for (Int_t i=0; i < nDDL; i++) {
289     if (buffers[i]) {
290       buffers[i]->Flush();
291       buffers[i]->WriteDataHeader(kFALSE, kFALSE);
292       delete buffers[i];
293     }
294   }
295
296   loader->UnloadDigits();
297 }
298
299 //____________________________________________________________________________
300 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr, AliCaloCalibPedestal* pedbadmap)
301 {
302   // convert raw data of the current event to digits                                                                                     
303
304   digitsArr->Clear(); 
305
306   if (!digitsArr) {
307     Error("Raw2Digits", "no digits found !");
308     return;
309   }
310   if (!reader) {
311     Error("Raw2Digits", "no raw reader found !");
312     return;
313   }
314
315   AliCaloRawStreamV3 in(reader,"EMCAL",fMapping);
316   // Select EMCAL DDL's;
317   reader->Select("EMCAL",0,43); // 43 = AliEMCALGeoParams::fgkLastAltroDDL
318
319   //Updated fitting routine from 2007 beam test takes into account
320   //possibility of two peaks in data and selects first one for fitting
321   //Also sets some of the starting parameters based on the shape of the
322   //given raw signal being fit
323
324   TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
325   signalF->SetParameters(10.,5.,fTau,fOrder,0.); //set all defaults once, just to be safe
326   signalF->SetParNames("amp","t0","tau","N","ped");
327   signalF->FixParameter(2,fTau); // tau in units of time bin
328   signalF->FixParameter(3,fOrder); // order
329   
330   Int_t id =  -1;
331   Float_t time = 0. ; 
332   Float_t amp = 0. ; 
333   Float_t ped = 0. ;
334   Float_t ampEstimate  = 0;
335   Float_t timeEstimate = 0;
336   Float_t pedEstimate = 0;
337   Int_t i = 0;
338   Int_t startBin = 0;
339
340   //Graph to hold data we will fit (should be converted to an array
341   //later to speed up processing
342   TGraph * gSig = new TGraph(GetRawFormatTimeBins()); 
343
344   Int_t lowGain = 0;
345   Int_t caloFlag = 0; // low, high gain, or TRU, or LED ref.
346
347   // start loop over input stream 
348   while (in.NextDDL()) {
349     while (in.NextChannel()) {
350
351       //Check if the signal  is high or low gain and then do the fit, 
352       //if it  is from TRU do not fit
353       caloFlag = in.GetCaloFlag();
354       if (caloFlag != 0 && caloFlag != 1) continue; 
355               
356       //Do not fit bad channels
357       if(fRemoveBadChannels && pedbadmap->IsBadChannel(in.GetModule(),in.GetColumn(),in.GetRow())) {
358         //printf("Tower from SM %d, column %d, row %d is BAD!!! Skip \n", in.GetModule(),in.GetColumn(),in.GetRow());
359         continue;
360       }  
361
362       // There can be zero-suppression in the raw data, 
363       // so set up the TGraph in advance
364       for (i=0; i < GetRawFormatTimeBins(); i++) {
365         gSig->SetPoint(i, i , -1); // init to out-of-range values
366       }
367
368       Int_t maxTimeBin = 0;
369       Int_t min = 0x3ff; // init to 10-bit max
370       Int_t max = 0; // init to 10-bit min
371       while (in.NextBunch()) {
372
373         const UShort_t *sig = in.GetSignals();
374         startBin = in.GetStartTimeBin();
375         if (maxTimeBin < startBin) {
376           maxTimeBin = startBin; // timebins come in reverse order
377         }       
378         if (maxTimeBin < 0 || maxTimeBin >= GetRawFormatTimeBins()) {
379           AliWarning(Form("Invalid time bin %d",maxTimeBin));
380           maxTimeBin = GetRawFormatTimeBins();
381         }
382         
383         for (i = 0; i < in.GetBunchLength(); i++) {
384           time = startBin--;
385           gSig->SetPoint((Int_t)time, time, (Double_t) sig[i]) ;
386           if (max < sig[i]) max = sig[i];
387           if (min > sig[i]) min = sig[i];
388           
389         }
390       } // loop over bunches
391
392       gSig->Set(maxTimeBin+1); // set actual max size of TGraph
393       
394       //Initialize the variables, do not keep previous values.
395       // not really necessary to reset all of them (only amp and time at the moment), but better safe than sorry
396       amp  = -1 ;
397       time = -1 ;
398       ped = -1;
399       ampEstimate  = -1 ;
400       timeEstimate = -1 ;
401       pedEstimate = -1;
402
403           if ( (max - min) > fNoiseThreshold) {
404                   switch(fFittingAlgorithm) 
405                   {
406                           case kStandard:
407                           {
408                                   //printf("Standard fitter \n");
409                                   FitRaw(gSig, signalF, maxTimeBin, amp, time, ped,
410                ampEstimate, timeEstimate, pedEstimate);
411                                   break;
412                           }       
413                           case kFastFit:
414                           {
415                                   //printf("FastFitter \n");
416                                   Double_t eSignal = 0;
417                                   Double_t dAmp = amp;
418                                   Double_t dTimeEstimate = timeEstimate;
419                                   Double_t eTimeEstimate = 0;
420                                   Double_t eAmp = 0;
421                                   Double_t chi2 = 0;
422
423                                   AliCaloFastAltroFitv0::FastFit(gSig->GetX(), gSig->GetY(), gSig->GetN(),
424                                                                                                  eSignal, fTau,
425                                                                                                  dAmp, eAmp, dTimeEstimate, eTimeEstimate, chi2);
426                                   amp=dAmp;
427                                   timeEstimate = dTimeEstimate;
428                                   //printf("FastFitter: Amp %f, time %f, eAmp %f, eTimeEstimate %f, chi2 %f\n",amp, timeEstimate,eAmp,eTimeEstimate,chi2);
429
430                                   break;
431                           }  
432                   }
433           }
434            
435       if ( amp>0 && amp<2000 && time>0 && time<(maxTimeBin*GetRawFormatTimeBinWidth()) ) {  //check both high and low end of amplitude result, and time
436         //2000 is somewhat arbitrary - not nice with magic numbers in the code..
437         id =  fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
438         lowGain = in.IsLowGain();
439
440         // check fit results: should be consistent with initial estimates
441         // more magic numbers, but very loose cuts, for now..
442         // We have checked that amp an time values are positive so division for assymmetry
443         // calculation should be OK/safe
444         Float_t ampAsymm = (amp - ampEstimate)/(amp + ampEstimate);
445         if ( (TMath::Abs(ampAsymm) > 0.1) ||
446              (TMath::Abs(time - timeEstimate) > 2*GetRawFormatTimeBinWidth()) ) {
447           AliDebug(2,Form("Fit results ped %f amp %f time %f not consistent with expectations ped %f max-ped %f time %d",
448                       ped, amp, time, pedEstimate, ampEstimate, timeEstimate));
449
450           // what should do we do then? skip this channel or assign the simple estimate? 
451           // for now just overwrite the fit results with the simple estimate
452           amp = ampEstimate;
453           time = timeEstimate; 
454         }
455
456         AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
457         // printf("Added tower: SM %d, row %d, column %d, amp %3.2f\n",in.GetModule(), in.GetRow(), in.GetColumn(),amp);
458         // round off amplitude value to nearest integer
459         AddDigit(digitsArr, id, lowGain, TMath::Nint(amp), time); 
460       }
461       
462       // Reset graph
463       for (Int_t index = 0; index < gSig->GetN(); index++) {
464         gSig->SetPoint(index, index, -1) ;  
465       } 
466       // Reset starting parameters for fit function
467       signalF->SetParameters(10.,5.,fTau,fOrder,0.); //reset all defaults just to be safe
468
469    } // end while over channel   
470   } //end while over DDL's, of input stream 
471   
472   delete signalF ; 
473   delete gSig;
474   
475   return ; 
476 }
477
478 //____________________________________________________________________________ 
479 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
480   //
481   // Add a new digit. 
482   // This routine checks whether a digit exists already for this tower 
483   // and then decides whether to use the high or low gain info
484   //
485   // Called by Raw2Digits
486   
487   AliEMCALDigit *digit = 0, *tmpdigit = 0;
488   TIter nextdigit(digitsArr);
489   while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
490     if (tmpdigit->GetId() == id)
491       digit = tmpdigit;
492   }
493
494   if (!digit) { // no digit existed for this tower; create one
495     if (lowGain && amp > fgkOverflowCut) 
496       amp = Int_t(fHighLowGainFactor * amp); 
497     Int_t idigit = digitsArr->GetEntries();
498     new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;   
499   }
500   else { // a digit already exists, check range 
501          // (use high gain if signal < cut value, otherwise low gain)
502     if (lowGain) { // new digit is low gain
503       if (digit->GetAmp() > fgkOverflowCut) {  // use if stored digit is out of range
504         digit->SetAmp(Int_t(fHighLowGainFactor * amp));
505         digit->SetTime(time);
506       }
507     }
508     else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
509       digit->SetAmp(amp);
510       digit->SetTime(time);
511     }
512   }
513 }
514
515 //____________________________________________________________________________ 
516 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, const Int_t lastTimeBin, Float_t & amp, Float_t & time, Float_t & ped, Float_t & ampEstimate, Float_t & timeEstimate, Float_t & pedEstimate, const Float_t cut) const 
517 {
518   // Fits the raw signal time distribution; from AliEMCALGetter 
519   // last argument: Float_t cut = 0.0; // indicating how much of amplitude w.r.t. max value fit should be above noise and pedestal 
520
521   // initialize return values
522   amp = 0; 
523   time = 0; 
524   ped = 0;
525   ampEstimate = 0;
526   timeEstimate = 0;
527   pedEstimate = 0;
528
529   // 0th step: remove plateau / overflow candidates
530   // before trying to estimate amplitude, search for maxima etc.
531   //
532   Int_t nOrig = gSig->GetN(); // number of samples before we remove any overflows
533   // Values for readback from input graph
534   Double_t ttime = 0;
535   Double_t signal = 0;
536
537   /*
538   // start: tmp dump of all values
539   for (Int_t i=0; i<gSig->GetN(); i++) {
540     gSig->GetPoint(i, ttime, signal) ; // get values
541     printf("orig: i %d, time %f, signal %f\n",i, ttime, signal);
542   }
543   // end: tmp dump of all values
544   */
545
546   // start from back of TGraph since RemovePoint will downshift indices
547   for (Int_t i=nOrig-1; i>=0; i--) {
548     gSig->GetPoint(i, ttime, signal) ; // get values
549     if (signal >= (pedEstimate + fgkOverflowCut) ) {
550       gSig->RemovePoint(i);
551     }
552   }
553
554   // 1st step: we try to estimate the pedestal value
555   Int_t nPed = 0;
556   for (Int_t index = 0; index < gSig->GetN(); index++) {
557     gSig->GetPoint(index, ttime, signal) ; 
558     // ttime < fNPedsamples used for pedestal estimate; 
559     // ttime >= fNPedSamples used for signal checks
560     if (signal >= 0 && ttime<fNPedSamples) { // valid value
561       pedEstimate += signal;
562       nPed++;
563     }
564   }
565
566   if (nPed > 0)
567     pedEstimate /= nPed;
568   else {
569     //AliWarning("Could not determine pedestal");         
570     AliDebug(1,"Could not determine pedestal");
571     pedEstimate = 0; // good estimate for ZeroSupp data (non ZS data should have no problem with pedestal estimate)
572   }
573
574   // 2nd step: we look through the rest of the time-bins/ADC values and
575   // see if we have something that looks like a signal.
576   // We look for a first local maxima, as well as for a global maxima 
577   Int_t locMaxFound = 0;
578   Int_t locMaxId = 0; // time-bin index at first local max
579   Float_t locMaxSig = -1; // actual local max value
580   Int_t globMaxId = 0; // time-bin index at global max
581   Float_t globMaxSig = -1; // actual global max value
582   // We will also look for any values that look like they are in overflow region
583   for (Int_t i=0; i<gSig->GetN(); i++) {
584     gSig->GetPoint(i, ttime, signal) ; // get values
585
586     // ttime < fNPedsamples used for pedestal estimate; 
587     // ttime >= fNPedSamples used for signal checks
588     if (ttime >= fNPedSamples) { 
589
590       // look for first local maximum signal=ADC value
591       if (!locMaxFound && signal > locMaxSig) {
592         locMaxId = i;
593         locMaxSig = signal;
594       }
595       else if ( locMaxSig > (pedEstimate + fNoiseThreshold) ) { 
596         // we enter this condition after signal<=max, but previous
597         // max value was large enough. I.e. at least a significant local 
598         // maxima has been found (just before)
599         locMaxFound = 1;
600       }
601
602       // also check for global maximum..
603       if (signal > globMaxSig) {
604         globMaxId = i;
605         globMaxSig = signal;
606       }
607     } // ttime check
608   } // end for-loop over samples after pedestal
609
610   // OK, we have looked through the signal spectra, let's see if we should try to make the fit
611   ampEstimate = locMaxSig - pedEstimate; // estimate using first local maxima 
612   if ( ampEstimate > fNoiseThreshold ) { // else it's just noise 
613
614     //Check that the local maximum we will use is not at the end or beginning of time sample range
615     Double_t timeMax = -1;
616     Int_t iMax = locMaxId;
617     gSig->GetPoint(locMaxId, timeMax, signal) ;
618     if (timeMax < 2 || timeMax > lastTimeBin-1) { // lastTimeBin is the lowest kept time-sample; current (Dec 2009) case
619       //    if (timeMax < 2 || timeMax > lastTimeBin-2) { // for when lastTimeBin is the lowest read-out time-sample, future (2010) case
620       AliDebug(1,Form("Skip fit, maximum of the sample close to the edges : timeMax %3.2f, ampEstimate %3.2f",timeMax, ampEstimate));
621       return;
622     }
623
624     // Check if the local and global maximum disagree
625     if (locMaxId != globMaxId) {
626       AliDebug(1,Form("Warning, local first maximum %d does not agree with global maximum %d\n", locMaxId, globMaxId));
627       return;
628     }
629     
630     // Get the maximum and find the lowest timebin (tailmin) where the ADC value is not 
631     // significantly different from the pedestal
632     // first lower times edge a.k.a. tailmin
633     Int_t tailMin = 0;
634     Double_t tmptime = 0;
635     for (Int_t i=iMax-1; i > 0; i--) {
636       gSig->GetPoint(i, tmptime, signal) ;
637       if((signal-pedEstimate) < fNoiseThreshold){
638         tailMin = i;
639         break;
640       }
641     }
642     // then same exercise for the higher times edge a.k.a. tailmax
643     Int_t tailMax = lastTimeBin;
644     for (Int_t i=iMax+1; i < gSig->GetN(); i++) {
645       gSig->GetPoint(i, tmptime, signal) ;
646       if ((signal-pedEstimate) <= (ampEstimate*cut + fNoiseThreshold)) { // stop fit at cut-fraction of amplitude above noise-threshold (cut>0 would mean avoid the pulse shape falling edge)
647         tailMax = i;
648         break;
649       }
650     }
651
652     // remove all points which are not in the distribution around maximum
653     // i.e. up to tailmin, and from tailmax
654     if ( tailMax != (gSig->GetN()-1) ){ // else nothing to remove
655       nOrig = gSig->GetN(); // can't use GetN call in for loop below since gSig size changes..
656       for(int j = tailMax; j < nOrig; j++) gSig->RemovePoint(tailMax);
657     }
658     for(int j = 0; j<=tailMin; j++) gSig->RemovePoint(0);
659
660     if(gSig->GetN() < 3) {
661       AliDebug(2,Form("Skip fit, number of entries in sample smaller than number of fitting parameters: in sample %d, fitting param 3", 
662                       gSig->GetN() ));
663       return;
664     }
665
666     timeEstimate = timeMax * GetRawFormatTimeBinWidth();
667
668     // determine what the valid fit range is
669     Double_t minFit = 9999;
670     Double_t maxFit = 0;
671     for (Int_t i=0; i < gSig->GetN(); i++) {
672       gSig->GetPoint(i, ttime, signal); 
673       if (minFit > ttime) minFit=ttime;
674       if (maxFit < ttime) maxFit=ttime;
675       //debug: printf("no tail: i %d, time %f, signal %f\n",i, ttime, signal); 
676     } 
677     signalF->SetRange(minFit, maxFit);
678
679     signalF->FixParameter(4, pedEstimate) ; 
680     signalF->SetParameter(1, timeMax);
681     signalF->SetParameter(0, ampEstimate);
682     
683     gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
684
685     // assign fit results
686     amp = signalF->GetParameter(0); 
687     time = signalF->GetParameter(1) * GetRawFormatTimeBinWidth(); // skip subtraction of fgTimeTrigger?
688     ped = signalF->GetParameter(4); 
689
690     //BEG YS alternative methods to calculate the amplitude
691     Double_t * ymx = gSig->GetX() ; 
692     Double_t * ymy = gSig->GetY() ; 
693     const Int_t kN = 3 ; 
694     Double_t ymMaxX[kN] = {0., 0., 0.} ; 
695     Double_t ymMaxY[kN] = {0., 0., 0.} ; 
696     Double_t ymax = 0. ; 
697       // find the maximum amplitude
698     Int_t ymiMax = 0 ;  
699     for (Int_t ymi = 0; ymi < gSig->GetN(); ymi++) {
700       if (ymy[ymi] > ymMaxY[0] ) {
701         ymMaxY[0] = ymy[ymi] ; //<========== This is the maximum amplitude
702         ymMaxX[0] = ymx[ymi] ;
703         ymiMax = ymi ; 
704       }
705     }
706       // find the maximum by fitting a parabola through the max and the two adjacent samples
707     if ( ymiMax < gSig->GetN()-1 && ymiMax > 0) {
708       ymMaxY[1] = ymy[ymiMax+1] ;
709       ymMaxY[2] = ymy[ymiMax-1] ; 
710       ymMaxX[1] = ymx[ymiMax+1] ;
711       ymMaxX[2] = ymx[ymiMax-1] ; 
712       if (ymMaxY[0]*ymMaxY[1]*ymMaxY[2] > 0) {
713           //fit a parabola through the 3 points y= a+bx+x*x*x
714         Double_t sy = 0 ; 
715         Double_t sx = 0 ; 
716         Double_t sx2 = 0 ; 
717         Double_t sx3 = 0 ; 
718         Double_t sx4 = 0 ; 
719         Double_t sxy = 0 ; 
720         Double_t sx2y = 0 ; 
721         for (Int_t i = 0; i < kN ; i++) {
722           sy += ymMaxY[i] ; 
723           sx += ymMaxX[i] ;             
724           sx2 += ymMaxX[i]*ymMaxX[i] ; 
725           sx3 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
726           sx4 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
727           sxy += ymMaxX[i]*ymMaxY[i] ; 
728           sx2y += ymMaxX[i]*ymMaxX[i]*ymMaxY[i] ; 
729         }
730         Double_t cN = (sx2y*kN-sy*sx2)*(sx3*sx-sx2*sx2)-(sx2y*sx-sxy*sx2)*(sx3*kN-sx*sx2); 
731         Double_t cD = (sx4*kN-sx2*sx2)*(sx3*sx-sx2*sx2)-(sx4*sx-sx3*sx2)*(sx3*kN-sx*sx2) ;
732         Double_t c  = cN / cD ; 
733         Double_t b  = ((sx2y*kN-sy*sx2)-c*(sx4*kN-sx2*sx2))/(sx3*kN-sx*sx2) ;
734         Double_t a  = (sy-b*sx-c*sx2)/kN  ;
735         Double_t xmax = -b/(2*c) ; 
736         ymax = a + b*xmax + c*xmax*xmax ;//<========== This is the maximum amplitude
737       }
738     }
739
740     Double_t diff = TMath::Abs(1-ymMaxY[0]/amp) ; 
741     if (diff > 0.1) 
742       amp = ymMaxY[0] ; 
743
744       //END YS
745
746   } // ampEstimate > fNoiseThreshold
747   return;
748 }
749 //__________________________________________________________________
750 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
751 {
752   // Matches version used in 2007 beam test
753   //
754   // Shape of the electronics raw reponse:
755   // It is a semi-gaussian, 2nd order Gamma function of the general form
756   //
757   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
758   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
759   // F = 0                                   for xx < 0 
760   //
761   // parameters:
762   // A:   par[0]   // Amplitude = peak value
763   // t0:  par[1]
764   // tau: par[2]
765   // N:   par[3]
766   // ped: par[4]
767   //
768   Double_t signal ;
769   Double_t tau =par[2];
770   Double_t n =par[3];
771   Double_t ped = par[4];
772   Double_t xx = ( x[0] - par[1] + tau ) / tau ;
773
774   if (xx <= 0) 
775     signal = ped ;  
776   else {  
777     signal = ped + par[0] * TMath::Power(xx , n) * TMath::Exp(n * (1 - xx )) ; 
778   }
779   return signal ;  
780 }
781
782 //__________________________________________________________________
783 Bool_t AliEMCALRawUtils::RawSampledResponse(
784 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const 
785 {
786   // for a start time dtime and an amplitude damp given by digit, 
787   // calculates the raw sampled response AliEMCAL::RawResponseFunction
788
789   Bool_t lowGain = kFALSE ; 
790
791   // A:   par[0]   // Amplitude = peak value
792   // t0:  par[1]                            
793   // tau: par[2]                            
794   // N:   par[3]                            
795   // ped: par[4]
796
797   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
798   signalF.SetParameter(0, damp) ; 
799   signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ; 
800   signalF.SetParameter(2, fTau) ; 
801   signalF.SetParameter(3, fOrder);
802   signalF.SetParameter(4, fgPedestalValue);
803
804   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
805     Double_t signal = signalF.Eval(iTime) ;     
806
807     // Next lines commeted for the moment but in principle it is not necessary to add
808     // extra noise since noise already added at the digits level.       
809
810     //According to Terry Awes, 13-Apr-2008
811     //add gaussian noise in quadrature to each sample
812     //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
813     //signal = sqrt(signal*signal + noise*noise);
814
815     // March 17,09 for fast fit simulations by Alexei Pavlinov.
816     // Get from PHOS analysis. In some sense it is open questions.
817     //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
818     //signal += noise; 
819
820     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
821     if ( adcH[iTime] > fgkRawSignalOverflow ){  // larger than 10 bits 
822       adcH[iTime] = fgkRawSignalOverflow ;
823       lowGain = kTRUE ; 
824     }
825
826     signal /= fHighLowGainFactor;
827
828     adcL[iTime] =  static_cast<Int_t>(signal + 0.5) ;
829     if ( adcL[iTime] > fgkRawSignalOverflow)  // larger than 10 bits 
830       adcL[iTime] = fgkRawSignalOverflow ;
831   }
832   return lowGain ; 
833 }