ROOT logo
/**************************************************************************
 * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved.      *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

/* $Id: $ */

// This class extracts the signal parameters (energy, time, quality)
// from ALTRO samples. Energy is in ADC counts, time is in time bin units.
// A fitting algorithm evaluates the energy and the time from Minuit minimization
// 
// Typical use case:
//     AliPHOSRawFitter *fitter=new AliPHOSRawFitter();
//     fitter->SetChannelGeo(module,cellX,cellZ,caloFlag);
//     fitter->SetCalibData(fgCalibData) ;
//     fitter->Eval(sig,sigStart,sigLength);
//     Double_t amplitude = fitter.GetEnergy();
//     Double_t time      = fitter.GetTime();
//     Bool_t   isLowGain = fitter.GetCaloFlag()==0;

// Author: Dmitri Peressounko (Oct.2008)
// Modified: Yuri Kharlov (Jul.2009)

// --- ROOT system ---
#include "TArrayI.h"
#include "TList.h"
#include "TMath.h"
#include "TMinuit.h"

// --- AliRoot header files ---
#include "AliLog.h"
#include "AliPHOSCalibData.h"
#include "AliPHOSRawFitterv1.h"
#include "AliPHOSPulseGenerator.h"

ClassImp(AliPHOSRawFitterv1)

//-----------------------------------------------------------------------------
AliPHOSRawFitterv1::AliPHOSRawFitterv1():
  AliPHOSRawFitterv0(),
  fSampleParamsLow(0x0),
  fSampleParamsHigh(0x0),
  fToFit(0x0)
{
  //Default constructor.
  if(!gMinuit) 
    gMinuit = new TMinuit(100);
  fSampleParamsHigh =new TArrayD(7) ;
  fSampleParamsHigh->AddAt(2.174,0) ;
  fSampleParamsHigh->AddAt(0.106,1) ;
  fSampleParamsHigh->AddAt(0.173,2) ;
  fSampleParamsHigh->AddAt(0.06106,3) ;
  //last two parameters are pedestal and overflow
  fSampleParamsLow=new TArrayD(7) ;
  fSampleParamsLow->AddAt(2.456,0) ;
  fSampleParamsLow->AddAt(0.137,1) ;
  fSampleParamsLow->AddAt(2.276,2) ;
  fSampleParamsLow->AddAt(0.08246,3) ;
  fToFit = new TList() ;
}

//-----------------------------------------------------------------------------
AliPHOSRawFitterv1::~AliPHOSRawFitterv1()
{
  //Destructor.
  //Destructor
  if(fSampleParamsLow){
    delete fSampleParamsLow ; 
    fSampleParamsLow=0 ;
  }
  if(fSampleParamsHigh){
    delete fSampleParamsHigh ;
    fSampleParamsHigh=0;
  }
  if(fToFit){
    delete fToFit ;
    fToFit=0 ;
  }
}

//-----------------------------------------------------------------------------
AliPHOSRawFitterv1::AliPHOSRawFitterv1(const AliPHOSRawFitterv1 &phosFitter ):
  AliPHOSRawFitterv0(phosFitter),
  fSampleParamsLow(0x0),
  fSampleParamsHigh(0x0),
  fToFit(0x0)
{
  //Copy constructor.
  fToFit = new TList() ;
  fSampleParamsLow =new TArrayD(*(phosFitter.fSampleParamsLow)) ;
  fSampleParamsHigh=new TArrayD(*(phosFitter.fSampleParamsHigh)) ;
}

//-----------------------------------------------------------------------------
AliPHOSRawFitterv1& AliPHOSRawFitterv1::operator = (const AliPHOSRawFitterv1 &phosFitter)
{
  //Assignment operator.
  if(this != &phosFitter) {
    fToFit = new TList() ;
    if(fSampleParamsLow){
      fSampleParamsLow = phosFitter.fSampleParamsLow ;
      fSampleParamsHigh= phosFitter.fSampleParamsHigh ;
    }
    else{
      fSampleParamsLow =new TArrayD(*(phosFitter.fSampleParamsLow)) ; 
      fSampleParamsHigh=new TArrayD(*(phosFitter.fSampleParamsHigh)) ;
    }
  }
  return *this;
}

//-----------------------------------------------------------------------------
Bool_t AliPHOSRawFitterv1::Eval(const UShort_t *signal, Int_t sigStart, Int_t sigLength)
{
  //Extract an energy deposited in the crystal,
  //crystal' position (module,column,row),
  //time and gain (high or low).
  //First collects sample, then evaluates it and if it has
  //reasonable shape, fits it with Gamma2 function and extracts 
  //energy and time.

  if (fCaloFlag == 2 || fNBunches > 1) {
    fQuality = 1000;
    return kTRUE;
  }

  Float_t pedMean = 0;
  Float_t pedRMS  = 0;
  Int_t   nPed    = 0;
  const Float_t kBaseLine   = 1.0;
  const Int_t   kPreSamples = 10;
  
  TArrayI *fSamples = new TArrayI(sigLength); // array of sample amplitudes
  TArrayI *fTimes   = new TArrayI(sigLength); // array of sample time stamps
  for (Int_t i=0; i<sigLength; i++) {
    if (i<kPreSamples) {
      nPed++;
      pedMean += signal[i];
      pedRMS  += signal[i]*signal[i] ;
    }
    fSamples->AddAt(signal[i],sigLength-i-1);
    fTimes  ->AddAt(i ,i);
  }

  fEnergy = -111;
  fQuality= 999. ;
  const Float_t sampleMaxHG=102.332 ;  //maximal height of HG sample with given parameterization
  const Float_t sampleMaxLG=277.196 ;  //maximal height of LG sample with given parameterization
  const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
  Double_t pedestal = 0;

  if (fPedSubtract) {
    if (nPed > 0) {
      fPedestalRMS=(pedRMS - pedMean*pedMean/nPed)/nPed ;
      if(fPedestalRMS > 0.) 
	fPedestalRMS = TMath::Sqrt(fPedestalRMS) ;
      fEnergy -= (Double_t)(pedMean/nPed); // pedestal subtraction
    }
    else
      return kFALSE;
  }
  else {
    //take pedestals from DB
    pedestal = (Double_t) fAmpOffset ;
    if (fCalibData) {
      Float_t truePed       = fCalibData->GetADCpedestalEmc(fModule, fCellZ, fCellX) ;
      Int_t   altroSettings = fCalibData->GetAltroOffsetEmc(fModule, fCellZ, fCellX) ;
      pedestal += truePed - altroSettings ;
    }
    else{
      AliWarning(Form("Can not read data from OCDB")) ;
    }
    fEnergy-=pedestal ;
  }

  if (fEnergy < kBaseLine) fEnergy = 0;
  //Evaluate time
  Int_t iStart = 0;
  while(iStart<sigLength && fSamples->At(iStart)-pedestal <kBaseLine) iStart++ ;
  fTime = sigStart-sigLength+iStart; 
  
  //calculate time and energy
  Int_t    maxBin=0 ;
  Int_t    maxAmp=0 ;
  Double_t aMean =0. ;
  Double_t aRMS  =0. ;
  Double_t wts   =0 ;
  Int_t    tStart=0 ;

  for (Int_t i=0; i<sigLength; i++){
    if(signal[i] > pedestal){
      Double_t de = signal[i] - pedestal ;
      if(de > 1.) {
	aMean += de*i ;
	aRMS  += de*i*i ;
	wts   += de; 
      }
      if(de > 2 && tStart==0) 
	tStart = i ;
      if(maxAmp < signal[i]){
	maxBin = i ;
	maxAmp = signal[i] ;
      }
    }
  }

  if (maxBin==sigLength-1){//bad "rising" sample
    fEnergy =    0. ;
    fTime   = -999. ;
    fQuality=  999. ;
    return kTRUE ;
  }

  fEnergy=Double_t(maxAmp)-pedestal ;
  fOverflow =0 ;  //look for plato on the top of sample
  if (fEnergy>500 &&  //this is not fluctuation of soft sample
     maxBin<sigLength-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
    fOverflow = kTRUE ;
  }
  
  if (wts > 0) {
    aMean /= wts; 
    aRMS   = aRMS/wts - aMean*aMean;
  }

  //do not take too small energies
  if (fEnergy < kBaseLine) 
    fEnergy = 0;
  
  //do not test quality of too soft samples
  if (fEnergy < maxEtoFit){
    fTime = tStart;
    if (aRMS < 2.) //sigle peak
      fQuality = 999. ;
    else
      fQuality =   0. ;
    return kTRUE ;
  }
      
  // if sample has reasonable mean and RMS, try to fit it with gamma2
  
  gMinuit->mncler();                     // Reset Minuit's list of paramters
  gMinuit->SetPrintLevel(-1) ;           // No Printout
  gMinuit->SetFCN(AliPHOSRawFitterv1::UnfoldingChiSquare) ;  

  // To set the address of the minimization function 
  
  fToFit->Clear("nodelete") ;
  Double_t b=0,bmin=0,bmax=0 ;
  if      (fCaloFlag == 0){ // Low gain
    fSampleParamsLow->AddAt(pedestal,4) ;
    if (fOverflow)
      fSampleParamsLow->AddAt(double(maxAmp),5) ;
    else
      fSampleParamsLow->AddAt(double(1023),5) ;
    fSampleParamsLow->AddAt(double(iStart),6) ;
    fToFit->AddFirst((TObject*)fSampleParamsLow) ; 
    b=fSampleParamsLow->At(2) ;
    bmin=0.5 ;
    bmax=10. ;
  }
  else if (fCaloFlag == 1){ // High gain
    fSampleParamsHigh->AddAt(pedestal,4) ;
    if (fOverflow)
      fSampleParamsHigh->AddAt(double(maxAmp),5) ;
    else
      fSampleParamsHigh->AddAt(double(1023),5);
    fSampleParamsHigh->AddAt(double(iStart),6);
    fToFit->AddFirst((TObject*)fSampleParamsHigh) ; 
    b=fSampleParamsHigh->At(2) ;
    bmin=0.05 ;
    bmax=0.4 ;
  }
  fToFit->AddLast((TObject*)fSamples) ;
  fToFit->AddLast((TObject*)fTimes) ;
  
  gMinuit->SetObjectFit((TObject*)fToFit) ;         // To tranfer pointer to UnfoldingChiSquare
  Int_t ierflg ;
  gMinuit->mnparm(0, "t0",  1.*tStart, 0.01, -500., 500., ierflg) ;
  if(ierflg != 0){
    //	  AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
   fEnergy =   0. ;
    fTime   =-999. ;
    fQuality= 999. ;
    return kTRUE ; //will scan further
  }
  Double_t amp0=0; 
  if      (fCaloFlag == 0) // Low gain
    amp0 = fEnergy/sampleMaxLG;
  else if (fCaloFlag == 1) // High gain
    amp0 = fEnergy/sampleMaxHG;
  
  gMinuit->mnparm(1, "Energy", amp0 , 0.01*amp0, 0, 0, ierflg) ;
  if(ierflg != 0){
    //	  AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
    fEnergy =   0. ;
    fTime   =-999. ;
    fQuality= 999. ;
    return kTRUE ; //will scan further
  }
  
  gMinuit->mnparm(2, "p2", b, 0.01*b, bmin, bmax, ierflg) ;
  if(ierflg != 0){                                         
    //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;  
    fEnergy =   0. ;
    fTime   =-999. ;
    fQuality= 999. ;
    return kTRUE ; //will scan further  
  }             
  
  Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
  //  depends on it. 
  Double_t p1 = 1.0 ;
  Double_t p2 = 0.0 ;
  gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
  gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
  //	gMinuit->SetMaxIterations(100);
  gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
  
  gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
  
  Double_t err=0.,t0err=0. ;
  Double_t t0=0.,efit=0. ;
  gMinuit->GetParameter(0,t0, t0err) ;
  gMinuit->GetParameter(1,efit, err) ;
  
  Double_t bfit=0., berr=0. ;
  gMinuit->GetParameter(2,bfit,berr) ;
  
  //Calculate total energy
  //this is parameterization of dependence of pulse height on parameter b
  if(fCaloFlag == 0) // Low gain
    efit *= 99.54910 + 78.65038*bfit ;
  else if(fCaloFlag == 1) // High gain
    efit *= 80.33109 + 128.6433*bfit ;
  
  if(efit < 0. || efit > 10000.){
    //set energy to previously found max
    fTime   =-999.;
    fQuality= 999 ;
    return kTRUE;
  }                                                                             
  
  //evaluate fit quality
  Double_t fmin=0.,fedm=0.,errdef=0. ;
  Int_t npari,nparx,istat;
  gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
  fQuality = fmin/sigLength ;
  //compare quality with some parameterization
  if      (fCaloFlag == 0) // Low gain
    fQuality /= 2.00 + 0.0020*fEnergy ;
  else if (fCaloFlag == 1) // High gain
    fQuality /= 0.75 + 0.0025*fEnergy ;
  
  fEnergy = efit ;
  fTime  += t0 - 4.024*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 70 samples
//  fTime  += sigStart;
  
  delete fSamples ;
  delete fTimes ;
  return kTRUE;
}
//-----------------------------------------------------------------------------
Double_t AliPHOSRawFitterv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){  //Function for fitting samples
  //parameters:
  //dt-time after start
  //en-amplutude
  //function parameters
  
  Double_t ped=params->At(4) ;
  if(dt<0.)
    return ped ; //pedestal
  else
    return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;
}
//_____________________________________________________________________________
void AliPHOSRawFitterv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
{
  // Number of parameters, Gradient, Chi squared, parameters, what to do

  TList * toFit= (TList*)gMinuit->GetObjectFit() ;
  TArrayD * params=(TArrayD*)toFit->At(0) ; 
  TArrayI * samples = (TArrayI*)toFit->At(1) ;
  TArrayI * times = (TArrayI*)toFit->At(2) ;

  fret = 0. ;     
  if(iflag == 2)
    for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
      Grad[iparam] = 0 ; // Will evaluate gradient
  
  Double_t t0=x[0] ;
  Double_t en=x[1] ;
  Double_t b=x[2] ;
  Double_t n=params->At(0) ;
  Double_t alpha=params->At(1) ;
  Double_t beta=params->At(3) ;
  //  Double_t ped=params->At(4) ;

  Double_t overflow=params->At(5) ;
  Int_t iBin = (Int_t) params->At(6) ;
  Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
  // iBin - first non-zero sample 
  Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
  Double_t ddt=times->At(iBin)-t0-tStep ;
  Double_t exp1=TMath::Exp(-alpha*ddt) ;
  Double_t exp2=TMath::Exp(-beta*ddt) ;
  Double_t dexp1=TMath::Exp(-alpha*tStep) ;
  Double_t dexp2=TMath::Exp(-beta*tStep) ;
  for(Int_t i = iBin; i<nSamples ; i++) {
    Double_t dt=double(times->At(i))-t0 ;
    Double_t fsample = double(samples->At(i)) ;
    Double_t diff=0. ;
    exp1*=dexp1 ;
    exp2*=dexp2 ;
//    if(fsample>=overflow)
//      continue ;    
    if(dt<=0.){
      diff=fsample ; 
      fret += diff*diff ;
      continue ;
    }
    
    Double_t dtn=TMath::Power(dt,n) ;
    Double_t dtnE=dtn*exp1 ;
    Double_t dt2E=dt*dt*exp2 ;
    Double_t fit=en*(dtnE + b*dt2E) ;
    if(fsample>=overflow && fit>=overflow)
      continue ;    

    diff = fsample - fit ;
    fret += diff*diff ;
    if(iflag == 2){  // calculate gradient
      Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta))  ; //derivative over t0
      Grad[1] -= diff*(dtnE+b*dt2E) ;
      Grad[2] -= en*diff*dt2E ;
    }
    
  }
  if(iflag == 2)
    for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
      Grad[iparam] *= 2. ; 
}
 AliPHOSRawFitterv1.cxx:1
 AliPHOSRawFitterv1.cxx:2
 AliPHOSRawFitterv1.cxx:3
 AliPHOSRawFitterv1.cxx:4
 AliPHOSRawFitterv1.cxx:5
 AliPHOSRawFitterv1.cxx:6
 AliPHOSRawFitterv1.cxx:7
 AliPHOSRawFitterv1.cxx:8
 AliPHOSRawFitterv1.cxx:9
 AliPHOSRawFitterv1.cxx:10
 AliPHOSRawFitterv1.cxx:11
 AliPHOSRawFitterv1.cxx:12
 AliPHOSRawFitterv1.cxx:13
 AliPHOSRawFitterv1.cxx:14
 AliPHOSRawFitterv1.cxx:15
 AliPHOSRawFitterv1.cxx:16
 AliPHOSRawFitterv1.cxx:17
 AliPHOSRawFitterv1.cxx:18
 AliPHOSRawFitterv1.cxx:19
 AliPHOSRawFitterv1.cxx:20
 AliPHOSRawFitterv1.cxx:21
 AliPHOSRawFitterv1.cxx:22
 AliPHOSRawFitterv1.cxx:23
 AliPHOSRawFitterv1.cxx:24
 AliPHOSRawFitterv1.cxx:25
 AliPHOSRawFitterv1.cxx:26
 AliPHOSRawFitterv1.cxx:27
 AliPHOSRawFitterv1.cxx:28
 AliPHOSRawFitterv1.cxx:29
 AliPHOSRawFitterv1.cxx:30
 AliPHOSRawFitterv1.cxx:31
 AliPHOSRawFitterv1.cxx:32
 AliPHOSRawFitterv1.cxx:33
 AliPHOSRawFitterv1.cxx:34
 AliPHOSRawFitterv1.cxx:35
 AliPHOSRawFitterv1.cxx:36
 AliPHOSRawFitterv1.cxx:37
 AliPHOSRawFitterv1.cxx:38
 AliPHOSRawFitterv1.cxx:39
 AliPHOSRawFitterv1.cxx:40
 AliPHOSRawFitterv1.cxx:41
 AliPHOSRawFitterv1.cxx:42
 AliPHOSRawFitterv1.cxx:43
 AliPHOSRawFitterv1.cxx:44
 AliPHOSRawFitterv1.cxx:45
 AliPHOSRawFitterv1.cxx:46
 AliPHOSRawFitterv1.cxx:47
 AliPHOSRawFitterv1.cxx:48
 AliPHOSRawFitterv1.cxx:49
 AliPHOSRawFitterv1.cxx:50
 AliPHOSRawFitterv1.cxx:51
 AliPHOSRawFitterv1.cxx:52
 AliPHOSRawFitterv1.cxx:53
 AliPHOSRawFitterv1.cxx:54
 AliPHOSRawFitterv1.cxx:55
 AliPHOSRawFitterv1.cxx:56
 AliPHOSRawFitterv1.cxx:57
 AliPHOSRawFitterv1.cxx:58
 AliPHOSRawFitterv1.cxx:59
 AliPHOSRawFitterv1.cxx:60
 AliPHOSRawFitterv1.cxx:61
 AliPHOSRawFitterv1.cxx:62
 AliPHOSRawFitterv1.cxx:63
 AliPHOSRawFitterv1.cxx:64
 AliPHOSRawFitterv1.cxx:65
 AliPHOSRawFitterv1.cxx:66
 AliPHOSRawFitterv1.cxx:67
 AliPHOSRawFitterv1.cxx:68
 AliPHOSRawFitterv1.cxx:69
 AliPHOSRawFitterv1.cxx:70
 AliPHOSRawFitterv1.cxx:71
 AliPHOSRawFitterv1.cxx:72
 AliPHOSRawFitterv1.cxx:73
 AliPHOSRawFitterv1.cxx:74
 AliPHOSRawFitterv1.cxx:75
 AliPHOSRawFitterv1.cxx:76
 AliPHOSRawFitterv1.cxx:77
 AliPHOSRawFitterv1.cxx:78
 AliPHOSRawFitterv1.cxx:79
 AliPHOSRawFitterv1.cxx:80
 AliPHOSRawFitterv1.cxx:81
 AliPHOSRawFitterv1.cxx:82
 AliPHOSRawFitterv1.cxx:83
 AliPHOSRawFitterv1.cxx:84
 AliPHOSRawFitterv1.cxx:85
 AliPHOSRawFitterv1.cxx:86
 AliPHOSRawFitterv1.cxx:87
 AliPHOSRawFitterv1.cxx:88
 AliPHOSRawFitterv1.cxx:89
 AliPHOSRawFitterv1.cxx:90
 AliPHOSRawFitterv1.cxx:91
 AliPHOSRawFitterv1.cxx:92
 AliPHOSRawFitterv1.cxx:93
 AliPHOSRawFitterv1.cxx:94
 AliPHOSRawFitterv1.cxx:95
 AliPHOSRawFitterv1.cxx:96
 AliPHOSRawFitterv1.cxx:97
 AliPHOSRawFitterv1.cxx:98
 AliPHOSRawFitterv1.cxx:99
 AliPHOSRawFitterv1.cxx:100
 AliPHOSRawFitterv1.cxx:101
 AliPHOSRawFitterv1.cxx:102
 AliPHOSRawFitterv1.cxx:103
 AliPHOSRawFitterv1.cxx:104
 AliPHOSRawFitterv1.cxx:105
 AliPHOSRawFitterv1.cxx:106
 AliPHOSRawFitterv1.cxx:107
 AliPHOSRawFitterv1.cxx:108
 AliPHOSRawFitterv1.cxx:109
 AliPHOSRawFitterv1.cxx:110
 AliPHOSRawFitterv1.cxx:111
 AliPHOSRawFitterv1.cxx:112
 AliPHOSRawFitterv1.cxx:113
 AliPHOSRawFitterv1.cxx:114
 AliPHOSRawFitterv1.cxx:115
 AliPHOSRawFitterv1.cxx:116
 AliPHOSRawFitterv1.cxx:117
 AliPHOSRawFitterv1.cxx:118
 AliPHOSRawFitterv1.cxx:119
 AliPHOSRawFitterv1.cxx:120
 AliPHOSRawFitterv1.cxx:121
 AliPHOSRawFitterv1.cxx:122
 AliPHOSRawFitterv1.cxx:123
 AliPHOSRawFitterv1.cxx:124
 AliPHOSRawFitterv1.cxx:125
 AliPHOSRawFitterv1.cxx:126
 AliPHOSRawFitterv1.cxx:127
 AliPHOSRawFitterv1.cxx:128
 AliPHOSRawFitterv1.cxx:129
 AliPHOSRawFitterv1.cxx:130
 AliPHOSRawFitterv1.cxx:131
 AliPHOSRawFitterv1.cxx:132
 AliPHOSRawFitterv1.cxx:133
 AliPHOSRawFitterv1.cxx:134
 AliPHOSRawFitterv1.cxx:135
 AliPHOSRawFitterv1.cxx:136
 AliPHOSRawFitterv1.cxx:137
 AliPHOSRawFitterv1.cxx:138
 AliPHOSRawFitterv1.cxx:139
 AliPHOSRawFitterv1.cxx:140
 AliPHOSRawFitterv1.cxx:141
 AliPHOSRawFitterv1.cxx:142
 AliPHOSRawFitterv1.cxx:143
 AliPHOSRawFitterv1.cxx:144
 AliPHOSRawFitterv1.cxx:145
 AliPHOSRawFitterv1.cxx:146
 AliPHOSRawFitterv1.cxx:147
 AliPHOSRawFitterv1.cxx:148
 AliPHOSRawFitterv1.cxx:149
 AliPHOSRawFitterv1.cxx:150
 AliPHOSRawFitterv1.cxx:151
 AliPHOSRawFitterv1.cxx:152
 AliPHOSRawFitterv1.cxx:153
 AliPHOSRawFitterv1.cxx:154
 AliPHOSRawFitterv1.cxx:155
 AliPHOSRawFitterv1.cxx:156
 AliPHOSRawFitterv1.cxx:157
 AliPHOSRawFitterv1.cxx:158
 AliPHOSRawFitterv1.cxx:159
 AliPHOSRawFitterv1.cxx:160
 AliPHOSRawFitterv1.cxx:161
 AliPHOSRawFitterv1.cxx:162
 AliPHOSRawFitterv1.cxx:163
 AliPHOSRawFitterv1.cxx:164
 AliPHOSRawFitterv1.cxx:165
 AliPHOSRawFitterv1.cxx:166
 AliPHOSRawFitterv1.cxx:167
 AliPHOSRawFitterv1.cxx:168
 AliPHOSRawFitterv1.cxx:169
 AliPHOSRawFitterv1.cxx:170
 AliPHOSRawFitterv1.cxx:171
 AliPHOSRawFitterv1.cxx:172
 AliPHOSRawFitterv1.cxx:173
 AliPHOSRawFitterv1.cxx:174
 AliPHOSRawFitterv1.cxx:175
 AliPHOSRawFitterv1.cxx:176
 AliPHOSRawFitterv1.cxx:177
 AliPHOSRawFitterv1.cxx:178
 AliPHOSRawFitterv1.cxx:179
 AliPHOSRawFitterv1.cxx:180
 AliPHOSRawFitterv1.cxx:181
 AliPHOSRawFitterv1.cxx:182
 AliPHOSRawFitterv1.cxx:183
 AliPHOSRawFitterv1.cxx:184
 AliPHOSRawFitterv1.cxx:185
 AliPHOSRawFitterv1.cxx:186
 AliPHOSRawFitterv1.cxx:187
 AliPHOSRawFitterv1.cxx:188
 AliPHOSRawFitterv1.cxx:189
 AliPHOSRawFitterv1.cxx:190
 AliPHOSRawFitterv1.cxx:191
 AliPHOSRawFitterv1.cxx:192
 AliPHOSRawFitterv1.cxx:193
 AliPHOSRawFitterv1.cxx:194
 AliPHOSRawFitterv1.cxx:195
 AliPHOSRawFitterv1.cxx:196
 AliPHOSRawFitterv1.cxx:197
 AliPHOSRawFitterv1.cxx:198
 AliPHOSRawFitterv1.cxx:199
 AliPHOSRawFitterv1.cxx:200
 AliPHOSRawFitterv1.cxx:201
 AliPHOSRawFitterv1.cxx:202
 AliPHOSRawFitterv1.cxx:203
 AliPHOSRawFitterv1.cxx:204
 AliPHOSRawFitterv1.cxx:205
 AliPHOSRawFitterv1.cxx:206
 AliPHOSRawFitterv1.cxx:207
 AliPHOSRawFitterv1.cxx:208
 AliPHOSRawFitterv1.cxx:209
 AliPHOSRawFitterv1.cxx:210
 AliPHOSRawFitterv1.cxx:211
 AliPHOSRawFitterv1.cxx:212
 AliPHOSRawFitterv1.cxx:213
 AliPHOSRawFitterv1.cxx:214
 AliPHOSRawFitterv1.cxx:215
 AliPHOSRawFitterv1.cxx:216
 AliPHOSRawFitterv1.cxx:217
 AliPHOSRawFitterv1.cxx:218
 AliPHOSRawFitterv1.cxx:219
 AliPHOSRawFitterv1.cxx:220
 AliPHOSRawFitterv1.cxx:221
 AliPHOSRawFitterv1.cxx:222
 AliPHOSRawFitterv1.cxx:223
 AliPHOSRawFitterv1.cxx:224
 AliPHOSRawFitterv1.cxx:225
 AliPHOSRawFitterv1.cxx:226
 AliPHOSRawFitterv1.cxx:227
 AliPHOSRawFitterv1.cxx:228
 AliPHOSRawFitterv1.cxx:229
 AliPHOSRawFitterv1.cxx:230
 AliPHOSRawFitterv1.cxx:231
 AliPHOSRawFitterv1.cxx:232
 AliPHOSRawFitterv1.cxx:233
 AliPHOSRawFitterv1.cxx:234
 AliPHOSRawFitterv1.cxx:235
 AliPHOSRawFitterv1.cxx:236
 AliPHOSRawFitterv1.cxx:237
 AliPHOSRawFitterv1.cxx:238
 AliPHOSRawFitterv1.cxx:239
 AliPHOSRawFitterv1.cxx:240
 AliPHOSRawFitterv1.cxx:241
 AliPHOSRawFitterv1.cxx:242
 AliPHOSRawFitterv1.cxx:243
 AliPHOSRawFitterv1.cxx:244
 AliPHOSRawFitterv1.cxx:245
 AliPHOSRawFitterv1.cxx:246
 AliPHOSRawFitterv1.cxx:247
 AliPHOSRawFitterv1.cxx:248
 AliPHOSRawFitterv1.cxx:249
 AliPHOSRawFitterv1.cxx:250
 AliPHOSRawFitterv1.cxx:251
 AliPHOSRawFitterv1.cxx:252
 AliPHOSRawFitterv1.cxx:253
 AliPHOSRawFitterv1.cxx:254
 AliPHOSRawFitterv1.cxx:255
 AliPHOSRawFitterv1.cxx:256
 AliPHOSRawFitterv1.cxx:257
 AliPHOSRawFitterv1.cxx:258
 AliPHOSRawFitterv1.cxx:259
 AliPHOSRawFitterv1.cxx:260
 AliPHOSRawFitterv1.cxx:261
 AliPHOSRawFitterv1.cxx:262
 AliPHOSRawFitterv1.cxx:263
 AliPHOSRawFitterv1.cxx:264
 AliPHOSRawFitterv1.cxx:265
 AliPHOSRawFitterv1.cxx:266
 AliPHOSRawFitterv1.cxx:267
 AliPHOSRawFitterv1.cxx:268
 AliPHOSRawFitterv1.cxx:269
 AliPHOSRawFitterv1.cxx:270
 AliPHOSRawFitterv1.cxx:271
 AliPHOSRawFitterv1.cxx:272
 AliPHOSRawFitterv1.cxx:273
 AliPHOSRawFitterv1.cxx:274
 AliPHOSRawFitterv1.cxx:275
 AliPHOSRawFitterv1.cxx:276
 AliPHOSRawFitterv1.cxx:277
 AliPHOSRawFitterv1.cxx:278
 AliPHOSRawFitterv1.cxx:279
 AliPHOSRawFitterv1.cxx:280
 AliPHOSRawFitterv1.cxx:281
 AliPHOSRawFitterv1.cxx:282
 AliPHOSRawFitterv1.cxx:283
 AliPHOSRawFitterv1.cxx:284
 AliPHOSRawFitterv1.cxx:285
 AliPHOSRawFitterv1.cxx:286
 AliPHOSRawFitterv1.cxx:287
 AliPHOSRawFitterv1.cxx:288
 AliPHOSRawFitterv1.cxx:289
 AliPHOSRawFitterv1.cxx:290
 AliPHOSRawFitterv1.cxx:291
 AliPHOSRawFitterv1.cxx:292
 AliPHOSRawFitterv1.cxx:293
 AliPHOSRawFitterv1.cxx:294
 AliPHOSRawFitterv1.cxx:295
 AliPHOSRawFitterv1.cxx:296
 AliPHOSRawFitterv1.cxx:297
 AliPHOSRawFitterv1.cxx:298
 AliPHOSRawFitterv1.cxx:299
 AliPHOSRawFitterv1.cxx:300
 AliPHOSRawFitterv1.cxx:301
 AliPHOSRawFitterv1.cxx:302
 AliPHOSRawFitterv1.cxx:303
 AliPHOSRawFitterv1.cxx:304
 AliPHOSRawFitterv1.cxx:305
 AliPHOSRawFitterv1.cxx:306
 AliPHOSRawFitterv1.cxx:307
 AliPHOSRawFitterv1.cxx:308
 AliPHOSRawFitterv1.cxx:309
 AliPHOSRawFitterv1.cxx:310
 AliPHOSRawFitterv1.cxx:311
 AliPHOSRawFitterv1.cxx:312
 AliPHOSRawFitterv1.cxx:313
 AliPHOSRawFitterv1.cxx:314
 AliPHOSRawFitterv1.cxx:315
 AliPHOSRawFitterv1.cxx:316
 AliPHOSRawFitterv1.cxx:317
 AliPHOSRawFitterv1.cxx:318
 AliPHOSRawFitterv1.cxx:319
 AliPHOSRawFitterv1.cxx:320
 AliPHOSRawFitterv1.cxx:321
 AliPHOSRawFitterv1.cxx:322
 AliPHOSRawFitterv1.cxx:323
 AliPHOSRawFitterv1.cxx:324
 AliPHOSRawFitterv1.cxx:325
 AliPHOSRawFitterv1.cxx:326
 AliPHOSRawFitterv1.cxx:327
 AliPHOSRawFitterv1.cxx:328
 AliPHOSRawFitterv1.cxx:329
 AliPHOSRawFitterv1.cxx:330
 AliPHOSRawFitterv1.cxx:331
 AliPHOSRawFitterv1.cxx:332
 AliPHOSRawFitterv1.cxx:333
 AliPHOSRawFitterv1.cxx:334
 AliPHOSRawFitterv1.cxx:335
 AliPHOSRawFitterv1.cxx:336
 AliPHOSRawFitterv1.cxx:337
 AliPHOSRawFitterv1.cxx:338
 AliPHOSRawFitterv1.cxx:339
 AliPHOSRawFitterv1.cxx:340
 AliPHOSRawFitterv1.cxx:341
 AliPHOSRawFitterv1.cxx:342
 AliPHOSRawFitterv1.cxx:343
 AliPHOSRawFitterv1.cxx:344
 AliPHOSRawFitterv1.cxx:345
 AliPHOSRawFitterv1.cxx:346
 AliPHOSRawFitterv1.cxx:347
 AliPHOSRawFitterv1.cxx:348
 AliPHOSRawFitterv1.cxx:349
 AliPHOSRawFitterv1.cxx:350
 AliPHOSRawFitterv1.cxx:351
 AliPHOSRawFitterv1.cxx:352
 AliPHOSRawFitterv1.cxx:353
 AliPHOSRawFitterv1.cxx:354
 AliPHOSRawFitterv1.cxx:355
 AliPHOSRawFitterv1.cxx:356
 AliPHOSRawFitterv1.cxx:357
 AliPHOSRawFitterv1.cxx:358
 AliPHOSRawFitterv1.cxx:359
 AliPHOSRawFitterv1.cxx:360
 AliPHOSRawFitterv1.cxx:361
 AliPHOSRawFitterv1.cxx:362
 AliPHOSRawFitterv1.cxx:363
 AliPHOSRawFitterv1.cxx:364
 AliPHOSRawFitterv1.cxx:365
 AliPHOSRawFitterv1.cxx:366
 AliPHOSRawFitterv1.cxx:367
 AliPHOSRawFitterv1.cxx:368
 AliPHOSRawFitterv1.cxx:369
 AliPHOSRawFitterv1.cxx:370
 AliPHOSRawFitterv1.cxx:371
 AliPHOSRawFitterv1.cxx:372
 AliPHOSRawFitterv1.cxx:373
 AliPHOSRawFitterv1.cxx:374
 AliPHOSRawFitterv1.cxx:375
 AliPHOSRawFitterv1.cxx:376
 AliPHOSRawFitterv1.cxx:377
 AliPHOSRawFitterv1.cxx:378
 AliPHOSRawFitterv1.cxx:379
 AliPHOSRawFitterv1.cxx:380
 AliPHOSRawFitterv1.cxx:381
 AliPHOSRawFitterv1.cxx:382
 AliPHOSRawFitterv1.cxx:383
 AliPHOSRawFitterv1.cxx:384
 AliPHOSRawFitterv1.cxx:385
 AliPHOSRawFitterv1.cxx:386
 AliPHOSRawFitterv1.cxx:387
 AliPHOSRawFitterv1.cxx:388
 AliPHOSRawFitterv1.cxx:389
 AliPHOSRawFitterv1.cxx:390
 AliPHOSRawFitterv1.cxx:391
 AliPHOSRawFitterv1.cxx:392
 AliPHOSRawFitterv1.cxx:393
 AliPHOSRawFitterv1.cxx:394
 AliPHOSRawFitterv1.cxx:395
 AliPHOSRawFitterv1.cxx:396
 AliPHOSRawFitterv1.cxx:397
 AliPHOSRawFitterv1.cxx:398
 AliPHOSRawFitterv1.cxx:399
 AliPHOSRawFitterv1.cxx:400
 AliPHOSRawFitterv1.cxx:401
 AliPHOSRawFitterv1.cxx:402
 AliPHOSRawFitterv1.cxx:403
 AliPHOSRawFitterv1.cxx:404
 AliPHOSRawFitterv1.cxx:405
 AliPHOSRawFitterv1.cxx:406
 AliPHOSRawFitterv1.cxx:407
 AliPHOSRawFitterv1.cxx:408
 AliPHOSRawFitterv1.cxx:409
 AliPHOSRawFitterv1.cxx:410
 AliPHOSRawFitterv1.cxx:411
 AliPHOSRawFitterv1.cxx:412
 AliPHOSRawFitterv1.cxx:413
 AliPHOSRawFitterv1.cxx:414
 AliPHOSRawFitterv1.cxx:415
 AliPHOSRawFitterv1.cxx:416
 AliPHOSRawFitterv1.cxx:417
 AliPHOSRawFitterv1.cxx:418
 AliPHOSRawFitterv1.cxx:419
 AliPHOSRawFitterv1.cxx:420
 AliPHOSRawFitterv1.cxx:421
 AliPHOSRawFitterv1.cxx:422
 AliPHOSRawFitterv1.cxx:423
 AliPHOSRawFitterv1.cxx:424
 AliPHOSRawFitterv1.cxx:425
 AliPHOSRawFitterv1.cxx:426
 AliPHOSRawFitterv1.cxx:427
 AliPHOSRawFitterv1.cxx:428
 AliPHOSRawFitterv1.cxx:429
 AliPHOSRawFitterv1.cxx:430
 AliPHOSRawFitterv1.cxx:431
 AliPHOSRawFitterv1.cxx:432
 AliPHOSRawFitterv1.cxx:433
 AliPHOSRawFitterv1.cxx:434
 AliPHOSRawFitterv1.cxx:435
 AliPHOSRawFitterv1.cxx:436
 AliPHOSRawFitterv1.cxx:437
 AliPHOSRawFitterv1.cxx:438
 AliPHOSRawFitterv1.cxx:439
 AliPHOSRawFitterv1.cxx:440
 AliPHOSRawFitterv1.cxx:441
 AliPHOSRawFitterv1.cxx:442
 AliPHOSRawFitterv1.cxx:443
 AliPHOSRawFitterv1.cxx:444
 AliPHOSRawFitterv1.cxx:445
 AliPHOSRawFitterv1.cxx:446
 AliPHOSRawFitterv1.cxx:447
 AliPHOSRawFitterv1.cxx:448
 AliPHOSRawFitterv1.cxx:449
 AliPHOSRawFitterv1.cxx:450
 AliPHOSRawFitterv1.cxx:451
 AliPHOSRawFitterv1.cxx:452