1 /**************************************************************************
2 * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
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();
39 // Author: Dmitri Peressounko
41 // --- ROOT system ---
52 // --- AliRoot header files ---
54 #include "AliPHOSRawDecoderv1.h"
55 #include "AliPHOSPulseGenerator.h"
58 ClassImp(AliPHOSRawDecoderv1)
60 //-----------------------------------------------------------------------------
61 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1():AliPHOSRawDecoder(),
62 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
64 //Default constructor.
67 //-----------------------------------------------------------------------------
68 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(AliRawReader* rawReader, AliAltroMapping **mapping):
69 AliPHOSRawDecoder(rawReader,mapping),
70 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
72 //Construct a decoder object.
73 //Is is user responsibility to provide next raw event
74 //using AliRawReader::NextEvent().
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() ;
92 //-----------------------------------------------------------------------------
93 AliPHOSRawDecoderv1::~AliPHOSRawDecoderv1()
97 delete fSampleParamsLow ;
100 if(fSampleParamsHigh){
101 delete fSampleParamsHigh ;
110 //-----------------------------------------------------------------------------
111 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(const AliPHOSRawDecoderv1 &phosDecoder ):
112 AliPHOSRawDecoder(phosDecoder),
113 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
116 fToFit = new TList() ;
117 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
118 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
121 //-----------------------------------------------------------------------------
122 AliPHOSRawDecoderv1& AliPHOSRawDecoderv1::operator = (const AliPHOSRawDecoderv1 &phosDecoder)
124 //Assignment operator.
126 // if(this != &phosDecoder) {
128 fToFit = new TList() ;
129 if(fSampleParamsLow){
130 fSampleParamsLow = phosDecoder.fSampleParamsLow ;
131 fSampleParamsHigh= phosDecoder.fSampleParamsHigh ;
134 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
135 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
140 //-----------------------------------------------------------------------------
141 Bool_t AliPHOSRawDecoderv1::NextDigit()
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
150 //Debug=====================
151 // TCanvas * c = (TCanvas*)gROOT->FindObjectAny("CSample") ;
153 // c = new TCanvas("CSample","CSample") ;
155 // TH1D * h = (TH1D*)gROOT->FindObjectAny("hSample") ;
157 // h=new TH1D("hSample","",200,0.,200.) ;
159 // TF1 * fff = (TF1*)gROOT->FindObjectAny("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===========
164 AliCaloRawStream* in = fCaloStream;
166 Int_t iBin = fSamples->GetSize() ;
171 Float_t baseLine = 1.0;
172 const Float_t nPreSamples = 10;
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
177 while ( in->Next() ) {
180 tLength = in->GetTimeLength();
181 if(tLength!=fSamples->GetSize()) {
184 fSamples = new TArrayI(tLength);
185 fTimes = new TArrayI(tLength);
186 iBin= fSamples->GetSize() ;
193 // Fit the full sample
194 if(in->IsNewHWAddress() && iBin != fSamples->GetSize()) {
196 Double_t pedestal =0. ;
199 pedestal = (Double_t)(pedMean/nPed);
204 //calculate time and energy
211 for(Int_t i=iBin; i<fSamples->GetSize(); i++){
212 if(fSamples->At(i)>0){
213 Double_t de=fSamples->At(i)-pedestal ;
217 if(de>2 && tStart==0)
219 if(maxAmp<fSamples->At(i)){
221 maxAmp=fSamples->At(i) ;
225 if(maxBin==fSamples->GetSize()-1){//bad sample
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
240 aRMS=aRMS/cnts-aMean*aMean;
243 //Debug:=====Draw sample
244 //if(fEnergy>pedestal+10.){
248 // printf("fEnergy=%f, cnts=%d, aMean=%f, aRMS=%f \n",fEnergy,cnts,aMean,aRMS) ;
250 //======================
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
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
262 fSampleParamsLow->AddAt(pedestal,4) ;
264 fSampleParamsLow->AddAt(double(maxAmp),5) ;
266 fSampleParamsLow->AddAt(double(1023),5) ;
267 fSampleParamsLow->AddAt(double(iBin),6) ;
268 fToFit->AddFirst((TObject*)fSampleParamsLow) ;
271 fSampleParamsHigh->AddAt(pedestal,4) ;
273 fSampleParamsHigh->AddAt(double(maxAmp),5) ;
275 fSampleParamsHigh->AddAt(double(1023),5);
276 fSampleParamsHigh->AddAt(double(iBin),6);
277 fToFit->AddFirst((TObject*)fSampleParamsHigh) ;
279 fToFit->AddLast((TObject*)fSamples) ;
280 fToFit->AddLast((TObject*)fTimes) ;
282 gMinuit->SetObjectFit((TObject*)fToFit) ; // To tranfer pointer to UnfoldingChiSquare
284 gMinuit->mnparm(0, "t0", 1.*tStart, 0.1, 0, 0, ierflg) ;
286 // AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
290 return kTRUE ; //will scan further
294 amp0=fEnergy/sampleMaxLG;
296 amp0=fEnergy/sampleMaxHG;
298 gMinuit->mnparm(1, "Energy", amp0 , 0.001*amp0, 0, 0, ierflg) ;
300 // AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
304 return kTRUE ; //will scan further
307 Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
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
316 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
321 gMinuit->GetParameter(0,t0, t0err) ;
322 gMinuit->GetParameter(1,efit, err) ;
324 //Calculate total energy
330 if(efit<0. || efit > 10000.){
331 fEnergy=0 ; //bad sample
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) ;
345 //Debug================
346 // Double_t n,alpha,b,beta ;
348 // n=fSampleParamsLow->At(0) ;
349 // alpha=fSampleParamsLow->At(1) ;
350 // b=fSampleParamsLow->At(2) ;
351 // beta=fSampleParamsLow->At(3) ;
354 // n=fSampleParamsHigh->At(0) ;
355 // alpha=fSampleParamsHigh->At(1) ;
356 // b=fSampleParamsHigh->At(2) ;
357 // beta=fSampleParamsHigh->At(3) ;
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) ;
366 // fff->SetParameters(pedestal,efit/sampleMaxLG,t0,n,alpha,b,beta) ;
369 // fff->SetParameters(pedestal,efit/sampleMaxHG,t0,n,alpha,b,beta) ;
371 // fff->Draw("same") ;
375 //====================
380 fTime*=fPulseGenerator->GetRawFormatTimeTrigger() ;
382 if (fEnergy < baseLine) fEnergy = 0;
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;
399 // Fill array with samples
401 if(fPedSubtract && (in->GetTime() < nPreSamples)) {
402 pedMean += in->GetSignal();
405 fSamples->AddAt(in->GetSignal(),iBin);
406 fTimes->AddAt(in->GetTime(),iBin);
408 //Debug==============
409 // h->SetBinContent(iBin,in->GetSignal()) ;
415 //_____________________________________________________________________________
416 void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
418 // Calculates the Chi square for the samples minimization
419 // Number of parameters, Gradient, Chi squared, parameters, what to do
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) ;
428 for(Int_t iparam = 0 ; iparam < 2 ; iparam++)
429 Grad[iparam] = 0 ; // Will evaluate gradient
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)
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
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
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)) ;
459 //-----------------------------------------------------------------------------
460 Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,TArrayD * params){ //Function for fitting samples
462 //dt-time after start
464 //function parameters
466 Double_t ped=params->At(4) ;
468 return ped ; //pedestal
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))) ;