]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSRawDecoderv1.cxx
histo emum added; limits changed to fit cosmic
[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.SetOldRCUFormat(kTRUE);
27 //       dc.SubtractPedestals(kTRUE);
28 //       while ( dc.NextDigit() ) {
29 //         Int_t module = dc.GetModule();
30 //         Int_t column = dc.GetColumn();
31 //         Int_t row = dc.GetRow();
32 //         Double_t amplitude = dc.GetEnergy();
33 //         Double_t time = dc.GetTime();
34 //         Bool_t IsLowGain = dc.IsLowGain();
35 //            ..........
36 //       }
37 //     }
38
39 // Author: Dmitri Peressounko
40
41 // --- ROOT system ---
42 #include "TList.h"
43 #include "TMath.h"
44 #include "TMinuit.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 = (TCanvas*)gROOT->FindObjectAny("CSample") ;
152 //  if(!c)
153 //    c = new TCanvas("CSample","CSample") ;
154 // 
155 //  TH1D * h = (TH1D*)gROOT->FindObjectAny("hSample") ;
156 //  if(!h)
157 //    h=new TH1D("hSample","",200,0.,200.) ;
158 // 
159 //  TF1 * fff = (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   Int_t   nPed = 0;
171   Float_t baseLine = 1.0;
172   const Float_t nPreSamples = 10;
173   fQuality= 999. ;
174   const Float_t sampleMaxHG=102.332 ;  //maximal height of HG sample with given parameterization
175   const Float_t sampleMaxLG=277.196 ;  //maximal height of HG sample with given parameterization
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       Double_t pedestal =0. ;
197       if(fPedSubtract){ 
198         if (nPed > 0)
199           pedestal = (Double_t)(pedMean/nPed); 
200         else
201           return kFALSE;
202       }
203
204       //calculate time and energy
205       Int_t maxBin=0 ;
206       Int_t maxAmp=0 ;
207       Double_t aMean=0. ;                                                                                                                  
208       Double_t aRMS=0. ;                                                                                                                   
209       Int_t tStart = 0 ;                                                                                                                   
210       Int_t cnts=0 ;                                                                                                                       
211       for(Int_t i=iBin; i<fSamples->GetSize(); i++){
212         if(fSamples->At(i)>0){                                                                                                             
213           Double_t de=fSamples->At(i)-pedestal ;                                                                                           
214           aMean+=de ;                                                                                                                      
215           aRMS+=de*de ;                                                                                                                    
216           cnts++;                                                                                                                          
217           if(de>2 && tStart==0) 
218             tStart=i ;                                                                                                                     
219           if(maxAmp<fSamples->At(i)){
220             maxBin=i ;
221             maxAmp=fSamples->At(i) ;
222           }
223         }
224       }
225       if(maxBin==fSamples->GetSize()-1){//bad sample
226         fEnergy=0. ;
227         fTime=-999.;
228         fQuality= 999. ;
229         return kTRUE ;
230       }
231       fEnergy=Double_t(maxAmp)-pedestal ;
232       fOverflow =0 ;  //look for plato on the top of sample
233       if(fEnergy>500 &&  //this is not fluctuation of soft sample
234          maxBin<fSamples->GetSize()-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
235          fOverflow = kTRUE ;
236       }
237       
238       if(cnts>0){
239         aMean/=cnts; 
240         aRMS=aRMS/cnts-aMean*aMean;
241       }
242       
243 //Debug:=====Draw sample
244 //if(fEnergy>pedestal+10.){
245 //    c->cd() ;
246 //    h->Draw() ;
247 //    c->Update() ;
248 // printf("fEnergy=%f, cnts=%d, aMean=%f, aRMS=%f \n",fEnergy,cnts,aMean,aRMS) ;   
249 //}
250 //======================
251
252       //IF sample has reasonable mean and RMS, try to fit it with gamma2
253       if(fEnergy>2.&& cnts >20 && aMean>0. && aRMS>2.){ //more or less reasonable sample
254         
255         gMinuit->mncler();                     // Reset Minuit's list of paramters
256         gMinuit->SetPrintLevel(-1) ;           // No Printout
257         gMinuit->SetFCN(AliPHOSRawDecoderv1::UnfoldingChiSquare) ;  
258         // To set the address of the minimization function 
259         
260        fToFit->Clear() ;
261        if(fLowGainFlag){
262          fSampleParamsLow->AddAt(pedestal,4) ;
263          if(fOverflow)
264            fSampleParamsLow->AddAt(double(maxAmp),5) ;
265          else
266            fSampleParamsLow->AddAt(double(1023),5) ;
267          fSampleParamsLow->AddAt(double(iBin),6) ;
268          fToFit->AddFirst((TObject*)fSampleParamsLow) ; 
269         }
270         else{
271          fSampleParamsHigh->AddAt(pedestal,4) ;
272          if(fOverflow)
273            fSampleParamsHigh->AddAt(double(maxAmp),5) ;
274          else
275            fSampleParamsHigh->AddAt(double(1023),5);
276          fSampleParamsHigh->AddAt(double(iBin),6);
277          fToFit->AddFirst((TObject*)fSampleParamsHigh) ; 
278         }
279         fToFit->AddLast((TObject*)fSamples) ;
280         fToFit->AddLast((TObject*)fTimes) ;
281
282         gMinuit->SetObjectFit((TObject*)fToFit) ;         // To tranfer pointer to UnfoldingChiSquare
283         Int_t ierflg ;
284         gMinuit->mnparm(0, "t0",  1.*tStart, 0.1, 0, 0, ierflg) ;
285         if(ierflg != 0){
286 //        AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
287           fEnergy=0. ;
288           fTime=-999. ;
289           fQuality=999 ;
290           return kTRUE ; //will scan further
291         }
292         Double_t amp0; 
293         if(fLowGainFlag)
294           amp0=fEnergy/sampleMaxLG;
295         else
296           amp0=fEnergy/sampleMaxHG;
297
298         gMinuit->mnparm(1, "Energy", amp0 , 0.001*amp0, 0, 0, ierflg) ;
299         if(ierflg != 0){
300 //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
301           fEnergy=0. ;
302           fTime=-999. ;
303           fQuality=999 ;
304           return kTRUE ; //will scan further
305         }
306         
307         Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
308         //  depends on it. 
309         Double_t p1 = 1.0 ;
310         Double_t p2 = 0.0 ;
311         gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
312         gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
313         //      gMinuit->SetMaxIterations(100);
314         gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
315         
316         gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
317
318         
319         Double_t err,t0err ;
320         Double_t t0,efit ;
321         gMinuit->GetParameter(0,t0, t0err) ;    
322         gMinuit->GetParameter(1,efit, err) ;    
323
324         //Calculate total energy
325         if(fLowGainFlag)
326           efit*=sampleMaxLG;
327         else
328           efit*=sampleMaxHG;
329
330         if(efit<0. || efit > 10000.){                                                                          
331           fEnergy=0 ; //bad sample                                                    
332           fTime=-999.;                                                                
333           fQuality=999 ;                                                              
334           return kTRUE;
335         }                                                                             
336  
337         //evaluate fit quality
338         Double_t fmin,fedm,errdef ;
339         Int_t npari,nparx,istat;
340         gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
341         fQuality=fmin/(fSamples->GetSize()-iBin) ;
342         //compare quality with some parameterization
343         fQuality/=5.*TMath::Exp(0.0025*efit) ;
344
345 //Debug================
346 //        Double_t n,alpha,b,beta ;
347 //        if(fLowGainFlag){
348 //          n=fSampleParamsLow->At(0) ;
349 //          alpha=fSampleParamsLow->At(1) ;
350 //          b=fSampleParamsLow->At(2) ;
351 //          beta=fSampleParamsLow->At(3) ;
352 //        }
353 //        else{
354 //          n=fSampleParamsHigh->At(0) ;
355 //          alpha=fSampleParamsHigh->At(1) ;
356 //          b=fSampleParamsHigh->At(2) ;
357 //          beta=fSampleParamsHigh->At(3) ;
358 //        }
359 //
360 //    if( fQuality > 5.*TMath::Exp(0.0025*efit)){
361 //    if( fQuality > 1.){
362 //   printf("Col=%d, row=%d, qual=%f, E=%f \n",fColumn,fRow,fQuality,efit) ;
363 //    c->cd() ;
364 //    h->Draw() ;
365 //    if(fLowGainFlag){
366 //      fff->SetParameters(pedestal,efit/sampleMaxLG,t0,n,alpha,b,beta) ;
367 //    }
368 //    else{
369 //      fff->SetParameters(pedestal,efit/sampleMaxHG,t0,n,alpha,b,beta) ;
370 //    }
371 //    fff->Draw("same") ;
372 //    c->Update();
373 //    getchar() ;
374 //    }
375 //====================
376
377         fEnergy=efit ;
378         fTime=t0 ;
379         
380         fTime*=fPulseGenerator->GetRawFormatTimeTrigger() ; 
381
382         if (fEnergy < baseLine) fEnergy = 0;
383       }
384       else{ //bad sample
385         fEnergy=0. ;
386         fTime=-999. ;
387         fQuality=999.;
388       }
389       
390       return kTRUE;
391     }
392     
393     fLowGainFlag = in->IsLowGain();
394     fTime = fPulseGenerator->GetRawFormatTimeTrigger() * in->GetTime();
395     fModule = in->GetModule()+1;
396     fRow    = in->GetRow()   +1;
397     fColumn = in->GetColumn()+1;
398     
399     // Fill array with samples
400     iBin--;
401     if(fPedSubtract && (in->GetTime() < nPreSamples)) {
402       pedMean += in->GetSignal();
403       nPed++;
404     }
405     fSamples->AddAt(in->GetSignal(),iBin);
406     fTimes->AddAt(in->GetTime(),iBin);
407  
408 //Debug==============
409 //    h->SetBinContent(iBin,in->GetSignal()) ;
410     
411   } // in.Next()
412   
413   return kFALSE;
414 }
415 //_____________________________________________________________________________
416 void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
417 {
418   // Calculates the Chi square for the samples minimization
419   // Number of parameters, Gradient, Chi squared, parameters, what to do
420
421   TList * toFit= (TList*)gMinuit->GetObjectFit() ;
422   TArrayD * params=(TArrayD*)toFit->At(0) ; 
423   TArrayI * samples = (TArrayI*)toFit->At(1) ;
424   TArrayI * times = (TArrayI*)toFit->At(2) ;
425
426   fret = 0. ;     
427   if(iflag == 2)
428     for(Int_t iparam = 0 ; iparam < 2 ; iparam++)    
429       Grad[iparam] = 0 ; // Will evaluate gradient
430   
431   Double_t t0=x[0] ;
432   Double_t en=x[1] ;
433   Double_t n=params->At(0) ;
434   Double_t alpha=params->At(1) ;
435   Double_t b=params->At(2) ;
436   Double_t beta=params->At(3) ;
437   Int_t overflow=(Int_t)params->At(5) ;
438   Int_t iBin = (Int_t) params->At(6) ;
439   Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
440   
441   for(Int_t i = iBin ; i < nSamples ; i++) {
442     Int_t sample = samples->At(i) ;
443     if(sample==0 || sample==overflow) //zero or overflow - scip point
444       continue ;
445     Double_t dt=1.*times->At(i)-t0 ;
446     Double_t diff=float(sample)-Gamma2(dt,en,params) ;
447     if(iflag == 2){  // calculate gradient
448       if(dt>=0.){
449         Grad[0] +=  2.*en*diff*(TMath::Power(dt,n-1.)*(n-alpha*dt)*TMath::Exp(-alpha*dt)+
450                               b*dt*(2.-beta*dt)*TMath::Exp(-beta*dt))  ; //derivative over t0
451         Grad[1] += -2.*diff*(TMath::Power(dt,n)*TMath::Exp(-alpha*dt)+
452                      b*dt*dt*TMath::Exp(-dt*beta)) ;
453       }
454     }
455     fret += diff*diff ;
456   }
457   
458 }
459 //-----------------------------------------------------------------------------
460 Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,TArrayD * params){  //Function for fitting samples
461   //parameters:
462   //dt-time after start
463   //en-amplutude
464   //function parameters
465   
466   Double_t ped=params->At(4) ;
467   if(dt<0.)
468     return ped ; //pedestal
469   else
470     return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+params->At(2)*dt*dt*TMath::Exp(-dt*params->At(3))) ;
471 }
472