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