]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSRawDecoderv1.cxx
Updated histogram limits (PHOS energy)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSRawDecoderv1.cxx
1 /**************************************************************************
2  * Copyright(c) 2007, 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 // This class decodes the stream of ALTRO samples to extract
19 // the PHOS "digits" of current event. Uses fitting procedure
20 // to separate reasonable samples
21 // 
22 // Typical use case:
23 //     AliRawReader* rf = new AliRawReaderDate("2006run2211.raw");
24 //     AliPHOSRawDecoder dc(rf);
25 //     while (rf->NextEvent()) {
26 //       dc.SubtractPedestals(kTRUE);
27 //       while ( dc.NextDigit() ) {
28 //         Int_t module = dc.GetModule();
29 //         Int_t column = dc.GetColumn();
30 //         Int_t row = dc.GetRow();
31 //         Double_t amplitude = dc.GetEnergy();
32 //         Double_t time = dc.GetTime();
33 //         Bool_t IsLowGain = dc.IsLowGain();
34 //            ..........
35 //       }
36 //     }
37
38 // Author: Dmitri Peressounko
39
40 // --- ROOT system ---
41 #include "TList.h"
42 #include "TMath.h"
43 #include "TMinuit.h"
44 #include "AliPHOSCalibData.h"
45
46 #include "TCanvas.h"
47 #include "TH1.h"
48 #include "TH2.h"
49 #include "TF1.h"
50 #include "TROOT.h"
51
52 // --- AliRoot header files ---
53 //#include "AliLog.h"
54 #include "AliPHOSRawDecoderv1.h"
55 #include "AliPHOSPulseGenerator.h"
56
57
58 ClassImp(AliPHOSRawDecoderv1)
59
60 //-----------------------------------------------------------------------------
61   AliPHOSRawDecoderv1::AliPHOSRawDecoderv1():AliPHOSRawDecoder(),
62 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
63 {
64   //Default constructor.
65 }
66
67 //-----------------------------------------------------------------------------
68 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(AliRawReader* rawReader,  AliAltroMapping **mapping):
69   AliPHOSRawDecoder(rawReader,mapping),
70   fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
71 {
72   //Construct a decoder object.
73   //Is is user responsibility to provide next raw event 
74   //using AliRawReader::NextEvent().
75
76   if(!gMinuit) 
77     gMinuit = new TMinuit(100);
78   fSampleParamsHigh =new TArrayD(7) ;
79   fSampleParamsHigh->AddAt(2.174,0) ;
80   fSampleParamsHigh->AddAt(0.106,1) ;
81   fSampleParamsHigh->AddAt(0.173,2) ;
82   fSampleParamsHigh->AddAt(0.06106,3) ;
83   //last two parameters are pedestal and overflow
84   fSampleParamsLow=new TArrayD(7) ;
85   fSampleParamsLow->AddAt(2.456,0) ;
86   fSampleParamsLow->AddAt(0.137,1) ;
87   fSampleParamsLow->AddAt(2.276,2) ;
88   fSampleParamsLow->AddAt(0.08246,3) ;
89   fToFit = new TList() ;
90 }
91
92 //-----------------------------------------------------------------------------
93 AliPHOSRawDecoderv1::~AliPHOSRawDecoderv1()
94 {
95   //Destructor.
96   if(fSampleParamsLow){
97     delete fSampleParamsLow ; 
98     fSampleParamsLow=0 ;
99   }
100   if(fSampleParamsHigh){
101     delete fSampleParamsHigh ;
102     fSampleParamsHigh=0;
103   }
104   if(fToFit){
105     delete fToFit ;
106     fToFit=0 ;
107   }
108 }
109
110 //-----------------------------------------------------------------------------
111 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(const AliPHOSRawDecoderv1 &phosDecoder ):
112   AliPHOSRawDecoder(phosDecoder), 
113   fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
114 {
115   //Copy constructor.
116   fToFit = new TList() ;
117   fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
118   fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
119 }
120
121 //-----------------------------------------------------------------------------
122 AliPHOSRawDecoderv1& AliPHOSRawDecoderv1::operator = (const AliPHOSRawDecoderv1 &phosDecoder)
123 {
124   //Assignment operator.
125
126   //  if(this != &phosDecoder) {
127   //  }
128   fToFit = new TList() ;
129   if(fSampleParamsLow){
130     fSampleParamsLow = phosDecoder.fSampleParamsLow ;
131     fSampleParamsHigh= phosDecoder.fSampleParamsHigh ;
132   }
133   else{
134     fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ; 
135     fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
136   }
137   return *this;
138 }
139
140 //-----------------------------------------------------------------------------
141 Bool_t AliPHOSRawDecoderv1::NextDigit()
142 {
143   //Extract an energy deposited in the crystal,
144   //crystal' position (module,column,row),
145   //time and gain (high or low).
146   //First collects sample, then evaluates it and if it has
147   //reasonable shape, fits it with Gamma2 function and extracts 
148   //energy and time.
149
150 //Debug=====================
151 //  TCanvas * c = 0; //(TCanvas*)gROOT->FindObjectAny("CSample") ;
152 //  if(!c)
153 //    c = new TCanvas("CSample","CSample") ;
154 // 
155 //  TH1D * h = 0 ; //(TH1D*)gROOT->FindObjectAny("hSample") ;
156 //  if(!h)
157 //    h=new TH1D("hSample","",200,0.,200.) ;
158 // 
159 //  TF1 * fff = 0 ; //(TF1*)gROOT->FindObjectAny("fff") ;
160 //  if(!fff)
161 //    fff = new TF1("fff","[0]+[1]*((abs(x-[2]))^[3]*exp(-(x-[2])*[4])+[5]*(x-[2])*(x-[2])*exp(-(x-[2])*[6]))",0.,1000.) ;
162 //End debug===========
163   
164   AliCaloRawStream* in = fCaloStream;
165   
166   Int_t    iBin     = fSamples->GetSize() ;
167   Int_t    tLength  = 0;
168   fEnergy = -111;
169   Float_t pedMean = 0;
170   Float_t pedRMS = 0;
171   Int_t   nPed = 0;
172   Float_t baseLine = 1.0;
173   const Float_t nPreSamples = 10;
174   fQuality= 999. ;
175   const Float_t sampleMaxHG=102.332 ;  //maximal height of HG sample with given parameterization
176   const Float_t sampleMaxLG=277.196 ;  //maximal height of HG sample with given parameterization
177   const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
178
179   while ( in->Next() ) { 
180
181     if(!tLength) {
182       tLength = in->GetTimeLength();
183       if(tLength!=fSamples->GetSize()) {
184         delete fSamples ;
185         delete fTimes ;
186         fSamples = new TArrayI(tLength);
187         fTimes = new TArrayI(tLength);
188         iBin= fSamples->GetSize() ;
189       }
190       else{
191         fSamples->Reset() ;
192       }
193     }
194     
195     // Fit the full sample
196     if((in->IsNewHWAddress() && iBin != fSamples->GetSize()) //new HW address
197        ||(iBin<=0)) {  //or new signal in same address
198
199       //First remember new sample
200       fNewLowGainFlag = in->IsLowGain();                                                                                                        
201       fNewModule = in->GetModule()+1;                                                                                                           
202       fNewRow    = in->GetRow()   +1;                                                                                                           
203       fNewColumn = in->GetColumn()+1;                                                                                                           
204       fNewAmp = in->GetSignal() ;
205       fNewTime=in->GetTime() ;  
206   
207       //now handle already collected 
208       Double_t pedestal =0. ;
209       fPedestalRMS=0. ;
210       if(fPedSubtract){ 
211         if (nPed > 0){
212           pedestal = (Double_t)(pedMean/nPed); 
213           fPedestalRMS=pedRMS/nPed-pedestal*pedestal ;
214           if(fPedestalRMS>0.) fPedestalRMS=TMath::Sqrt(fPedestalRMS) ;
215         }
216         else
217           return kFALSE;
218       }
219       else{
220         //take pedestals from DB
221         pedestal = fAmpOffset ;
222         if(fCalibData){
223            Float_t truePed = fCalibData->GetADCpedestalEmc(fModule, fColumn, fRow) ;
224            Int_t   altroSettings = fCalibData->GetAltroOffsetEmc(fModule, fColumn, fRow) ;
225            pedestal += truePed - altroSettings ;
226          }
227          else{
228 //           printf("AliPHOSRawDecoderv1::NextDigit(): Can not read data from OCDB \n") ;
229          }
230       }
231
232       //calculate time and energy
233       Int_t maxBin=0 ;
234       Int_t maxAmp=0 ;
235       Double_t aMean=0. ;                                                                                                                  
236       Double_t aRMS=0. ;                                                                                                                   
237       Double_t wts=0 ;                                                                                                                       
238       Int_t tStart = 0 ;                                                                                                                   
239       for(Int_t i=iBin; i<fSamples->GetSize(); i++){
240         if(fSamples->At(i)>pedestal){                                                                                                             
241           Double_t de=fSamples->At(i)-pedestal ;                                                                                           
242           if(de>1.){
243             aMean+=de*i ;                                                                                                                      
244             aRMS+=de*i*i ;                                                                                                                    
245             wts+=de; 
246           }                                                                                                                         
247           if(de>2 && tStart==0) 
248             tStart=i ;                                                                                                                     
249           if(maxAmp<fSamples->At(i)){
250             maxBin=i ;
251             maxAmp=fSamples->At(i) ;
252           }
253         }
254       }
255       if(maxBin==fSamples->GetSize()-1){//bad "rising" sample
256         fEnergy=0. ;
257         fTime=-999.;
258         fQuality= 999. ;
259         return kTRUE ;
260       }
261       fEnergy=Double_t(maxAmp)-pedestal ;
262       fOverflow =0 ;  //look for plato on the top of sample
263       if(fEnergy>500 &&  //this is not fluctuation of soft sample
264          maxBin<fSamples->GetSize()-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
265          fOverflow = kTRUE ;
266       }
267       
268       if(wts>0){
269         aMean/=wts; 
270         aRMS=aRMS/wts-aMean*aMean;
271       }
272
273       //do not take too small energies
274       if(fEnergy < baseLine) 
275          fEnergy = 0;
276
277       //do not test quality of too soft samples
278       if(fEnergy<maxEtoFit){
279         fTime=fTimes->At(tStart);
280         if(aRMS<2.) //sigle peak
281           fQuality=999. ;
282         else
283           fQuality= 0. ;                                                                                                                   
284         return kTRUE ;                                                                                                                     
285       } 
286
287       
288 //Debug:=====Draw sample
289 //if(fEnergy>pedestal+10.){
290 //if(fLowGainFlag && fEnergy>2){
291 //  if(!c)
292 //    if(!fLowGainFlag && fRow==32 && fColumn==18){
293 //    TCanvas *c = new TCanvas("CSample","CSample") ;
294 //    c->cd() ;
295 //    h->Draw() ;
296 //    c->Update() ;
297 // printf("fEnergy=%f, aRMS=%f \n",fEnergy,aRMS) ;   
298 //getchar() ;
299 //}
300 //======================
301
302       //IF sample has reasonable mean and RMS, try to fit it with gamma2
303         
304         gMinuit->mncler();                     // Reset Minuit's list of paramters
305         gMinuit->SetPrintLevel(-1) ;           // No Printout
306         gMinuit->SetFCN(AliPHOSRawDecoderv1::UnfoldingChiSquare) ;  
307         // To set the address of the minimization function 
308         
309        fToFit->Clear("nodelete") ;
310        Double_t b,bmin,bmax ;
311        if(fLowGainFlag){
312          fSampleParamsLow->AddAt(pedestal,4) ;
313          if(fOverflow)
314            fSampleParamsLow->AddAt(double(maxAmp),5) ;
315          else
316            fSampleParamsLow->AddAt(double(1023),5) ;
317          fSampleParamsLow->AddAt(double(iBin),6) ;
318          fToFit->AddFirst((TObject*)fSampleParamsLow) ; 
319          b=fSampleParamsLow->At(2) ;
320          bmin=0.5 ;
321          bmax=10. ;
322        }
323        else{
324          fSampleParamsHigh->AddAt(pedestal,4) ;
325          if(fOverflow)
326            fSampleParamsHigh->AddAt(double(maxAmp),5) ;
327          else
328            fSampleParamsHigh->AddAt(double(1023),5);
329          fSampleParamsHigh->AddAt(double(iBin),6);
330          fToFit->AddFirst((TObject*)fSampleParamsHigh) ; 
331          b=fSampleParamsHigh->At(2) ;
332          bmin=0.05 ;
333          bmax=0.4 ;
334         }
335         fToFit->AddLast((TObject*)fSamples) ;
336         fToFit->AddLast((TObject*)fTimes) ;
337
338         gMinuit->SetObjectFit((TObject*)fToFit) ;         // To tranfer pointer to UnfoldingChiSquare
339         Int_t ierflg ;
340         gMinuit->mnparm(0, "t0",  1.*tStart, 0.01, -500., 500., ierflg) ;
341         if(ierflg != 0){
342 //        AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
343           fEnergy=0. ;
344           fTime=-999. ;
345           fQuality=999 ;
346           return kTRUE ; //will scan further
347         }
348         Double_t amp0; 
349         if(fLowGainFlag)
350           amp0=fEnergy/sampleMaxLG;
351         else
352           amp0=fEnergy/sampleMaxHG;
353
354         gMinuit->mnparm(1, "Energy", amp0 , 0.01*amp0, 0, 0, ierflg) ;
355         if(ierflg != 0){
356 //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
357           fEnergy=0. ;
358           fTime=-999. ;
359           fQuality=999 ;
360           return kTRUE ; //will scan further
361         }
362
363         gMinuit->mnparm(2, "p2", b, 0.01*b, bmin, bmax, ierflg) ;
364         if(ierflg != 0){                                         
365 //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;  
366           fEnergy=0. ;           
367           fTime=-999. ;         
368           fQuality=999 ;       
369           return kTRUE ; //will scan further  
370         }             
371  
372         
373         Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
374         //  depends on it. 
375         Double_t p1 = 1.0 ;
376         Double_t p2 = 0.0 ;
377         gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
378         gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
379         //      gMinuit->SetMaxIterations(100);
380         gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
381         
382         gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
383         
384         Double_t err,t0err ;
385         Double_t t0,efit ;
386         gMinuit->GetParameter(0,t0, t0err) ;    
387         gMinuit->GetParameter(1,efit, err) ;    
388
389         Double_t bfit, berr ;
390         gMinuit->GetParameter(2,bfit,berr) ;
391
392         //Calculate total energy
393         //this isparameterization of depetendence of pulse height on parameter b
394         if(fLowGainFlag)
395           efit*=99.54910 + 78.65038*bfit ;
396         else
397           efit*=80.33109+128.6433*bfit ;
398
399         if(efit<0. || efit > 10000.){                                                                          
400 //set energy to previously found max
401 //          fEnergy=0 ; //bad sample                                                    
402           fTime=-999.;                                                                
403           fQuality=999 ;                                                              
404           return kTRUE;
405         }                                                                             
406  
407         //evaluate fit quality
408         Double_t fmin,fedm,errdef ;
409         Int_t npari,nparx,istat;
410         gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
411         fQuality=fmin/(fSamples->GetSize()-iBin) ;
412         //compare quality with some parameterization
413         if(fLowGainFlag){
414           fQuality/=2.+0.002*fEnergy ;
415         }
416         else{
417           fQuality/=0.75+0.0025*fEnergy ;
418         }
419
420 //Debug================
421 //        Double_t n,alpha,beta ;
422 //        Double_t en ;
423 //       if(fLowGainFlag){
424 //          n=fSampleParamsLow->At(0) ;
425 //          alpha=fSampleParamsLow->At(1) ;
426 //          beta=fSampleParamsLow->At(3) ;
427 //          en=efit/(99.54910 + 78.65038*bfit) ;
428 //        }
429 //        else{
430 //          n=fSampleParamsHigh->At(0) ;
431 //          alpha=fSampleParamsHigh->At(1) ;
432 //          beta=fSampleParamsHigh->At(3) ;
433 //          en=efit/(80.33109+128.6433*bfit) ;
434 //        }
435 //
436 ////    if( fQuality > 1 && fEnergy > 20. && !fOverflow){
437 ////    if(!fLowGainFlag && fRow==32 && fColumn==18){
438 //{
439 //    printf("Col=%d, row=%d, qual=%f, E=%f, t0=%f, b=%f\n",fColumn,fRow,fQuality,efit,t0,bfit) ;
440 //    printf("    Energy = %f \n",fEnergy) ;
441 //    TCanvas * c = new TCanvas("samp") ;
442 //    c->cd() ;
443 //    h->Draw() ;
444 //    if(fLowGainFlag){
445 //      fff->SetParameters(pedestal,en,t0,n,alpha,bfit,beta) ;
446 //    }
447 //    else{
448 //     fff->SetParameters(pedestal,en,t0,n,alpha,bfit,beta) ;
449 //    }
450 //////    for(Int_t i=1;i<=h->GetNbinsX(); i++){
451 ////       Double_t x=h->GetBinCenter(i) ;
452 ////       h->SetBinContent(i,h->GetBinContent(i)-fff->Eval(x)) ;
453 ////    }
454 ////    h->SetMinimum(-15.) ;
455 ////    h->SetMaximum(15.) ;
456 //    h->Draw() ;
457 //    fff->Draw("same") ;
458 //    c->Update();
459 //    getchar() ;
460 //    }
461 //====================
462
463       fEnergy=efit ;
464       fTime=t0-4.024*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 70 samples
465 //      fTime=t0+2.8*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 50 samples
466 //      fQuality = bfit ;
467       
468       return kTRUE;
469     }
470     
471     fLowGainFlag = in->IsLowGain();
472     fModule = in->GetModule()+1;
473     fRow    = in->GetRow()   +1;
474     fColumn = in->GetColumn()+1;
475
476     //add previouly taken if coincides
477     if(fLowGainFlag==fNewLowGainFlag && fModule==fNewModule &&
478        fRow==fNewRow && fColumn==fNewColumn){
479        iBin--;                                                                                                                                
480        if(fPedSubtract && fNewTime < nPreSamples) {                                                                                    
481          pedMean += in->GetSignal();                                                                                                          
482          pedRMS += in->GetSignal()*in->GetSignal() ;
483          nPed++;                                                                                                                              
484        }                                                                                                                                      
485        fSamples->AddAt(fNewAmp,iBin);                                                                                                 
486        fTimes->AddAt(fNewTime,iBin);                                                                                                     
487     
488        //Mark that we already take it
489        fNewModule=-1 ;
490     }
491     
492     // Fill array with samples
493     iBin--;
494     if(fPedSubtract && (in->GetTime() < nPreSamples)) {
495       pedMean += in->GetSignal();
496       pedRMS += in->GetSignal()*in->GetSignal() ;
497       nPed++;
498     }
499     fSamples->AddAt(in->GetSignal()-10,iBin);
500     fTimes->AddAt(in->GetTime(),iBin);
501  
502 //Debug==============
503 //    h->SetBinContent(in->GetTime(),in->GetSignal()) ;
504 //EndDebug==============
505     
506   } // in.Next()
507   
508   return kFALSE;
509 }
510 //_____________________________________________________________________________
511 void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
512 {
513   // Number of parameters, Gradient, Chi squared, parameters, what to do
514
515   TList * toFit= (TList*)gMinuit->GetObjectFit() ;
516   TArrayD * params=(TArrayD*)toFit->At(0) ; 
517   TArrayI * samples = (TArrayI*)toFit->At(1) ;
518   TArrayI * times = (TArrayI*)toFit->At(2) ;
519
520   fret = 0. ;     
521   if(iflag == 2)
522     for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
523       Grad[iparam] = 0 ; // Will evaluate gradient
524   
525   Double_t t0=x[0] ;
526   Double_t en=x[1] ;
527   Double_t b=x[2] ;
528   Double_t n=params->At(0) ;
529   Double_t alpha=params->At(1) ;
530   Double_t beta=params->At(3) ;
531   Double_t ped=params->At(4) ;
532
533   Double_t overflow=params->At(5) ;
534   Int_t iBin = (Int_t) params->At(6) ;
535   Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
536   // iBin - first non-zero sample 
537   Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
538   Double_t ddt=times->At(iBin)-t0-tStep ;
539   Double_t exp1=TMath::Exp(-alpha*ddt) ;
540   Double_t exp2=TMath::Exp(-beta*ddt) ;
541   Double_t dexp1=TMath::Exp(-alpha*tStep) ;
542   Double_t dexp2=TMath::Exp(-beta*tStep) ;
543   for(Int_t i = iBin; i<nSamples ; i++) {
544     Double_t dt=double(times->At(i))-t0 ;
545     Double_t fsample = double(samples->At(i)) ;
546     if(fsample>=overflow)
547       continue ;
548     Double_t diff ;
549     exp1*=dexp1 ;
550     exp2*=dexp2 ;
551     if(dt<=0.){
552       diff=fsample - ped ; 
553       fret += diff*diff ;
554       continue ;
555     }
556     Double_t dtn=TMath::Power(dt,n) ;
557     Double_t dtnE=dtn*exp1 ;
558     Double_t dt2E=dt*dt*exp2 ;
559     Double_t fit=ped+en*(dtnE + b*dt2E) ;
560 //    if(fit>=overflow){
561 //      diff=fsample-overflow ;
562 //      fret += diff*diff ;
563 //      //zero gradient here
564 //    }
565 //    else{
566       diff = fsample - fit ;
567       fret += diff*diff ;
568       if(iflag == 2){  // calculate gradient
569         Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta))  ; //derivative over t0
570         Grad[1] -= diff*(dtnE+b*dt2E) ;
571         Grad[2] -= en*diff*dt2E ;
572       }
573 //    }
574   }
575   if(iflag == 2)
576     for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
577       Grad[iparam] *= 2. ; 
578 }
579 //-----------------------------------------------------------------------------
580 Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){  //Function for fitting samples
581   //parameters:
582   //dt-time after start
583   //en-amplutude
584   //function parameters
585   
586   Double_t ped=params->At(4) ;
587   if(dt<0.)
588     return ped ; //pedestal
589   else
590     return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;
591 }
592