ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, 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.                  *
 **************************************************************************/

// This class derives from AliEMCALClustrerizer

// --- Root ---
#include <TMath.h> 
#include <TMinuit.h>
#include <TTree.h> 
#include <TBenchmark.h>
#include <TBrowser.h>
#include <TROOT.h>
#include <TClonesArray.h>
#include <TH1I.h>

// --- AliRoot ---
#include "AliLog.h"
#include "AliEMCALRecPoint.h"
#include "AliEMCALDigit.h"
#include "AliEMCALGeometry.h"
#include "AliCaloCalibPedestal.h"
#include "AliEMCALCalibData.h"
#include "AliESDCaloCluster.h"
#include "AliEMCALUnfolding.h"

#include "AliEMCALClusterizerFixedWindow.h"

ClassImp(AliEMCALClusterizerFixedWindow)

//__________________________________________________________________________________________
AliEMCALClusterizerFixedWindow::AliEMCALClusterizerFixedWindow() :
  AliEMCALClusterizer(), 
  fNphi(4), 
  fNeta(4), 
  fShiftPhi(2), 
  fShiftEta(2),
  fTRUshift(0),
  fNEtaDigitsSupMod(0),
  fNPhiDigitsSupMod(0),
  fNTRUPhi(0),
  fNTRUEta(0),
  fNEtaDigits(0),
  fNPhiDigits(0),
  fMaxShiftPhi(0),
  fMaxShiftEta(0),
  fNDigitsCluster(0),
  fNClusEtaNoShift(0),
  fNClusPhiNoShift(0),
  fNClusters(0),
  fNTotalClus(0),
  fClustersArray(0),
  fInitialized(0)
{
  // Constructor
}

//__________________________________________________________________________________________
AliEMCALClusterizerFixedWindow::AliEMCALClusterizerFixedWindow(AliEMCALGeometry* geometry) :
  AliEMCALClusterizer(geometry), 
  fNphi(4), 
  fNeta(4), 
  fShiftPhi(2), 
  fShiftEta(2),
  fTRUshift(0),
  fNEtaDigitsSupMod(0),
  fNPhiDigitsSupMod(0),
  fNTRUPhi(0),
  fNTRUEta(0),
  fNEtaDigits(0),
  fNPhiDigits(0),
  fMaxShiftPhi(0),
  fMaxShiftEta(0),
  fNDigitsCluster(0),
  fNClusEtaNoShift(0),
  fNClusPhiNoShift(0),
  fNClusters(0),
  fNTotalClus(0),
  fClustersArray(0),
  fInitialized(0)
{
  // Constructor
}

//__________________________________________________________________________________________
AliEMCALClusterizerFixedWindow::AliEMCALClusterizerFixedWindow(AliEMCALGeometry* geometry, AliEMCALCalibData * calib, AliCaloCalibPedestal * caloped) :
  AliEMCALClusterizer(geometry, calib, caloped),
  fNphi(4), 
  fNeta(4), 
  fShiftPhi(2), 
  fShiftEta(2),
  fTRUshift(0),
  fNEtaDigitsSupMod(0),
  fNPhiDigitsSupMod(0),
  fNTRUPhi(0),
  fNTRUEta(0),
  fNEtaDigits(0),
  fNPhiDigits(0),
  fMaxShiftPhi(0),
  fMaxShiftEta(0),
  fNDigitsCluster(0),
  fNClusEtaNoShift(0),
  fNClusPhiNoShift(0),
  fNClusters(0),
  fNTotalClus(0),
  fClustersArray(0),
  fInitialized(0)
{
  // Constructor
}

//__________________________________________________________________________________________
AliEMCALClusterizerFixedWindow::~AliEMCALClusterizerFixedWindow()
{
  // Destructor

  if (fClustersArray) {
    for (Int_t i = 0; i < fNTotalClus; i++) {
      if (fClustersArray[i]) {
	delete[] fClustersArray[i];
	fClustersArray[i] = 0;
      }
    }
    delete[] fClustersArray;
    fClustersArray = 0;
  }

}

//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::SetNphi (Int_t n) 
{
  // Set fNphi; if clusterizer already initialized gives a warning and does nothing
  
  if (fInitialized != 0)
    AliWarning("Clusterizer already initialized. Unable to change the parameters.");
  else
    fNphi = n;
}

//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::SetNeta (Int_t n) 
{
  // Set fNeta; if clusterizer already initialized gives a warning and does nothing
  
  if (fInitialized != 0)
    AliWarning("Clusterizer already initialized. Unable to change the parameters.");
  else
    fNeta = n;
}

//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::SetShiftPhi (Int_t s) 
{
  // Set fShiftPhi; if clusterizer already initialized gives a warning and does nothing
  
  if (fInitialized != 0)
    AliWarning("Clusterizer already initialized. Unable to change the parameters.");
  else
    fShiftPhi = s;
}

//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::SetShiftEta (Int_t s) 
{
  // Set fShiftEta; if clusterizer already initialized gives a warning and does nothing
  
  if (fInitialized != 0)
    AliWarning("Clusterizer already initialized. Unable to change the parameters.");
  else
    fShiftEta = s;
}

//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::SetTRUshift(Bool_t b) 
{
  // Set fTRUshift; if clusterizer already initialized gives a warning and does nothing
  
  if (fInitialized != 0)
    AliWarning("Clusterizer already initialized. Unable to change the parameters.");
  else
    fTRUshift = b;
}


//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::Digits2Clusters(Option_t * option)
{
  // Steering method to perform clusterization for the current event 
 
  static Float_t cputime = 0;
  static Float_t realtime = 0;

  if (strstr(option,"tim"))
    gBenchmark->Start("EMCALClusterizer"); 
	
  if (strstr(option,"print"))
    Print(""); 
	
  //Get calibration parameters from file or digitizer default values.
  GetCalibrationParameters();
	
  //Get dead channel map from file or digitizer default values.
  GetCaloCalibPedestal();
	
  MakeClusters();  //only the real clusters
	
  if (fToUnfold) {
    fClusterUnfolding->SetInput(fNumberOfECAClusters,fRecPoints,fDigitsArr);
    fClusterUnfolding->MakeUnfolding();
  }
	
  //Evaluate position, dispersion and other RecPoint properties for EC section 
  for (Int_t index = 0; index < fRecPoints->GetEntries(); index++) { 
    AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index));
    if (rp) {
      rp->EvalAll(fECAW0,fDigitsArr,fJustClusters);
      AliDebug(5, Form("MAX INDEX %d ", rp->GetMaximalEnergyIndex()));
      //For each rec.point set the distance to the nearest bad crystal
      if (fCaloPed)
        rp->EvalDistanceToBadChannels(fCaloPed);
    }
  }
  
  fRecPoints->Sort();
	
  for (Int_t index = 0; index < fRecPoints->GetEntries(); index++) {
    AliEMCALRecPoint *rp = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index));
    if (rp) {
      rp->SetIndexInList(index);
    }
    else AliFatal("RecPoint NULL!!");
  }
	
  if (fTreeR)
    fTreeR->Fill();
	
  if (strstr(option,"deb") || strstr(option,"all"))  
    PrintRecPoints(option);
	
  AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
	
  if (strstr(option,"tim")) {
    gBenchmark->Stop("EMCALClusterizer");
    Printf("Exec took %f CPU time (%f real time) for clusterizing", 
           gBenchmark->GetCpuTime("EMCALClusterizer")-cputime,gBenchmark->GetRealTime("EMCALClusterizer")-realtime);
    cputime = gBenchmark->GetCpuTime("EMCALClusterizer");
    realtime = gBenchmark->GetRealTime("EMCALClusterizer");
  }    
}

//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::ExecOnce()
{
  // Initialize clusterizer.
  
  fInitialized = -1;

  if (!fGeom) {
    AliError("Did not get geometry!");
    return;
  }
	
  // Defining geometry and clusterization parameter
  fNEtaDigitsSupMod = fGeom->GetNEta() * fGeom->GetNETAdiv(); // always 48?;
  fNPhiDigitsSupMod = fGeom->GetNPhi() * fGeom->GetNPHIdiv(); // always 24?;
  
  fNTRUPhi = 1;
  fNTRUEta = 1;
  
  fNEtaDigits = fNEtaDigitsSupMod * fGeom->GetNumberOfSuperModules() / fGeom->GetNPhiSuperModule();
  fNPhiDigits = fNPhiDigitsSupMod * fGeom->GetNPhiSuperModule();    
  
  if (fTRUshift){
    fNTRUPhi = fGeom->GetNPhiSuperModule() * 3;
    fNTRUEta = fGeom->GetNumberOfSuperModules() / fGeom->GetNPhiSuperModule();
    fNEtaDigits /= fNTRUEta;
    fNPhiDigits /= fNTRUPhi;
  }

  // Check if clusterizer parameter are compatible with calorimeter geometry
  if (fNEtaDigits < fNeta){
    AliError(Form("Error: fNeta = %d is greater than nEtaDigits = %d.",fNeta,fNEtaDigits));
    return;
  }
  if (fNPhiDigits < fNphi){
    AliError(Form("Error: fNphi = %d is greater than nPhiDigits = %d.",fNphi,fNPhiDigits));
    return;
  }
  if (fNEtaDigits % fShiftEta != 0){
    AliError(Form("Error: fShiftEta = %d is such that clusters cannot slide the whole calorimeter (nEtaDigits = %d).",fShiftEta,fNEtaDigits));
    return;
  }
  if (fNPhiDigits % fShiftPhi != 0){
    AliError(Form("Error: fShiftPhi = %d is such that clusters cannot slide the whole calorimeter (nPhiDigits = %d).",fShiftPhi,fNPhiDigits));
    return;
  }
  if (fNeta % fShiftEta != 0){
    AliError(Form("Error: fShiftEta = %d is not divisor of fNeta = %d.",fShiftEta,fNeta));
    return;
  }
  if (fNphi % fShiftPhi != 0){
    AliError(Form("Error: fShiftPhi = %d is not divisor of fNphi = %d).",fShiftPhi,fNphi));
    return;
  }
  
  fMaxShiftPhi = fNphi / fShiftPhi;
  fMaxShiftEta = fNeta / fShiftEta;
  
  fNClusEtaNoShift = fNEtaDigits / fNeta;
  fNClusPhiNoShift = fNPhiDigits / fNphi;

  fNClusters = fNClusEtaNoShift * fNClusPhiNoShift * fNTRUEta * fNTRUPhi;
 
  fNTotalClus = fNClusters * fMaxShiftEta * fMaxShiftPhi;

  fNDigitsCluster = fNphi * fNeta;

  if (fClustersArray) {
    for (Int_t i = 0; i < fNTotalClus; i++) {
      if (fClustersArray[i]) {
	delete[] fClustersArray[i];
	fClustersArray[i] = 0;
      }
    }
    delete[] fClustersArray;
    fClustersArray = 0;
  }

  fClustersArray = new AliEMCALDigit**[fNTotalClus];
  for (Int_t i = 0; i < fNTotalClus; i++) {
    fClustersArray[i] = new AliEMCALDigit*[fNDigitsCluster];
    for (Int_t j = 0; j < fNDigitsCluster; j++) {
      fClustersArray[i][j] = 0;
    }
  }

  AliDebug(1,Form("****ExecOnce*****\n"
		  "fNphi = %d, fNeta = %d, fShiftPhi = %d, fShiftEta = %d, fTRUshift = %d\n"
		  "fNEtaDigitsSupMod = %d, fNPhiDigitsSupMod = %d, fNTRUPhi = %d, fNTRUEta = %d, fNEtaDigits = %d, fNPhiDigits = %d\n"
		  "fMaxShiftPhi = %d, fMaxShiftEta = %d, fNDigitsCluster = %d, fNClusEtaNoShift = %d, fNClusPhiNoShift = %d\n"
		  "fNClusters = %d, fNTotalClus = %d\n",
		  fNphi,fNeta,fShiftPhi,fShiftEta,fTRUshift,
		  fNEtaDigitsSupMod,fNPhiDigitsSupMod,fNTRUPhi,fNTRUEta,fNEtaDigits,fNPhiDigits,
		  fMaxShiftPhi,fMaxShiftEta,fNDigitsCluster,fNClusEtaNoShift,fNClusPhiNoShift,
		  fNClusters,fNTotalClus));

  fInitialized = 1;
}
  
//__________________________________________________________________________________________
void AliEMCALClusterizerFixedWindow::MakeClusters()
{
  // Make clusters.

  fNumberOfECAClusters = 0;
  fRecPoints->Delete();

  if (fInitialized == 0)
    ExecOnce();

  if (fInitialized == -1) {
    AliError(Form("%s: error initializing the clusterizer. No clusterization will be performed.",GetName()));
    return;
  }
  
  // Set up TObjArray with pointers to digits to work on calibrated digits 
  TObjArray *digitsC = new TObjArray();
  AliEMCALDigit *digit;
  Float_t dEnergyCalibrated = 0.0, ehs = 0.0, time = 0.0;
  TIter nextdigit(fDigitsArr);
  while ((digit = dynamic_cast<AliEMCALDigit*>(nextdigit()))) { // calibrate and clean up digits
    dEnergyCalibrated =  digit->GetAmplitude();
    time              =  digit->GetTime();
    Calibrate(dEnergyCalibrated, time, digit->GetId());
    digit->SetCalibAmp(dEnergyCalibrated);
    digit->SetTime(time);
    if (dEnergyCalibrated < fMinECut || time > fTimeMax || time < fTimeMin) {
      continue;
    }
    else if (!fGeom->CheckAbsCellId(digit->GetId())) {
      continue;
    }
    else {
      ehs += dEnergyCalibrated;
      digitsC->AddLast(digit);
    }
  } 
  
  AliDebug(1,Form("MakeClusters: Number of digits %d  -> (e %f), ehs %f\n",
                  fDigitsArr->GetEntries(),fMinECut,ehs));
   
  Int_t nSupMod=0, nModule=0, nIphi=0, nIeta=0;
  Int_t iphi=0, ieta=0;  // cell eta-phi indexes in SM

  for (Int_t ishiftPhi = 0; ishiftPhi < fMaxShiftPhi; ishiftPhi++){
    Int_t nClusPhi = (fNPhiDigits - fShiftPhi * ishiftPhi) / fNphi;
    
    for (Int_t ishiftEta = 0; ishiftEta < fMaxShiftEta; ishiftEta++) {
      
      Int_t nClusEta = (fNEtaDigits - fShiftEta * ishiftEta) / fNeta; 
      
      Int_t iTotalClus = fNClusters * (ishiftPhi * fMaxShiftEta + ishiftEta);
      
      TIter nextdigitC(digitsC);
      while ((digit = dynamic_cast<AliEMCALDigit*>(nextdigitC()))) { // scan over the list of digitsC
        
        fGeom->GetCellIndex (digit->GetId(), nSupMod, nModule, nIphi, nIeta);
        fGeom->GetCellPhiEtaIndexInSModule (nSupMod, nModule, nIphi, nIeta, iphi, ieta);
        
        Int_t iphi_eff = iphi - fShiftPhi * ishiftPhi + fNPhiDigitsSupMod * (nSupMod / 2); // N supermodules along phi
        
        Int_t iTRUphi = iphi_eff / fNPhiDigits;
        
        iphi_eff -= iTRUphi * fNPhiDigits;
        
        Int_t iClusPhi = iphi_eff / fNphi; 
        
        if (iphi_eff < 0 || iClusPhi >= nClusPhi) 
          continue;
        
        Int_t ieta_eff = ieta - fShiftEta * ishiftEta + fNEtaDigitsSupMod * (nSupMod % 2); // 2 supermodules along eta
        
        Int_t iTRUeta = ieta_eff / fNEtaDigits;
        
        ieta_eff -= iTRUeta * fNEtaDigits;
        
        Int_t iClusEta = ieta_eff / fNeta; 
        
        if (ieta_eff < 0 || iClusEta >= nClusEta) 
          continue;
        
        iphi_eff += iTRUphi * fNPhiDigits;
        iClusPhi = iphi_eff / fNphi; 
        
        ieta_eff += iTRUeta * fNEtaDigits;
        iClusEta = ieta_eff / fNeta; 
        
        Int_t iCluster = iClusPhi + iClusEta * fNClusPhiNoShift * fNTRUPhi; 
        Int_t iDigit = iphi_eff % fNphi + (ieta_eff % fNeta) * fNphi;

        if (iCluster >= fNClusters){
          AliError(Form("iCluster out of range! iCluster = %d, fNClusters = %d (should never happen...)", iCluster, fNClusters));
          return;
        }
        
        iCluster += iTotalClus;

        if (iCluster >= fNTotalClus){
          AliError(Form("iCluster out of range! iCluster = %d, fNTotalClus = %d (should never happen...)", iCluster, fNTotalClus));
          return;
        }

        if (iDigit >= fNDigitsCluster){
          AliError(Form("iDigit out of range! iDigit = %d, fNDigitsCluster = %d (should never happen...)", iDigit, fNDigitsCluster));
          return;
        }

        if (fClustersArray[iCluster][iDigit] != 0){
          AliError("Digit already added! (should never happen...)");
          return;
        }
        
        fClustersArray[iCluster][iDigit] = digit;
        
      } // loop on digit
      
    } // loop on eta shift
    
  } // loop on phi shift
  
  AliEMCALRecPoint *recPoint = 0;
  Bool_t recPointOk = kFALSE;
  for (Int_t iCluster = 0; iCluster < fNTotalClus; iCluster++) {

    if (!recPoint) {
      if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(fNumberOfECAClusters*2+1);

      (*fRecPoints)[fNumberOfECAClusters] = new AliEMCALRecPoint("");
      recPoint = static_cast<AliEMCALRecPoint*>(fRecPoints->At(fNumberOfECAClusters));
    }
		
    if (recPoint) {
      recPoint->SetClusterType(AliVCluster::kEMCALClusterv1);
      recPoint->SetUniqueID(iCluster);
      fNumberOfECAClusters++;

      for (Int_t iDigit = 0; iDigit < fNDigitsCluster; iDigit++) {
        if (fClustersArray[iCluster][iDigit] == 0) continue;
        digit = fClustersArray[iCluster][iDigit];
        recPoint->AddDigit(*digit, digit->GetCalibAmp(), kFALSE); //Time or TimeR?
        fClustersArray[iCluster][iDigit] = 0;
	recPointOk = kTRUE;
      }

      if (recPointOk) { // unset the pointer so that a new rec point will be allocated in the next iteration
	recPoint = 0;
	recPointOk = kFALSE;
      }
    }
    else {
      AliError("Error allocating rec points!");
      break;
    }
  }

  if (!recPointOk) {
    fRecPoints->RemoveLast();
    fNumberOfECAClusters--;
  }

  delete digitsC;
  AliDebug(1, Form("MakeClusters: Number of digits %d  -> (e %f)\n", fDigitsArr->GetEntries(),fMinECut));
  AliDebug(1, Form("total no of clusters %d from %d digits", fNumberOfECAClusters, fDigitsArr->GetEntriesFast())); 
}
 AliEMCALClusterizerFixedWindow.cxx:1
 AliEMCALClusterizerFixedWindow.cxx:2
 AliEMCALClusterizerFixedWindow.cxx:3
 AliEMCALClusterizerFixedWindow.cxx:4
 AliEMCALClusterizerFixedWindow.cxx:5
 AliEMCALClusterizerFixedWindow.cxx:6
 AliEMCALClusterizerFixedWindow.cxx:7
 AliEMCALClusterizerFixedWindow.cxx:8
 AliEMCALClusterizerFixedWindow.cxx:9
 AliEMCALClusterizerFixedWindow.cxx:10
 AliEMCALClusterizerFixedWindow.cxx:11
 AliEMCALClusterizerFixedWindow.cxx:12
 AliEMCALClusterizerFixedWindow.cxx:13
 AliEMCALClusterizerFixedWindow.cxx:14
 AliEMCALClusterizerFixedWindow.cxx:15
 AliEMCALClusterizerFixedWindow.cxx:16
 AliEMCALClusterizerFixedWindow.cxx:17
 AliEMCALClusterizerFixedWindow.cxx:18
 AliEMCALClusterizerFixedWindow.cxx:19
 AliEMCALClusterizerFixedWindow.cxx:20
 AliEMCALClusterizerFixedWindow.cxx:21
 AliEMCALClusterizerFixedWindow.cxx:22
 AliEMCALClusterizerFixedWindow.cxx:23
 AliEMCALClusterizerFixedWindow.cxx:24
 AliEMCALClusterizerFixedWindow.cxx:25
 AliEMCALClusterizerFixedWindow.cxx:26
 AliEMCALClusterizerFixedWindow.cxx:27
 AliEMCALClusterizerFixedWindow.cxx:28
 AliEMCALClusterizerFixedWindow.cxx:29
 AliEMCALClusterizerFixedWindow.cxx:30
 AliEMCALClusterizerFixedWindow.cxx:31
 AliEMCALClusterizerFixedWindow.cxx:32
 AliEMCALClusterizerFixedWindow.cxx:33
 AliEMCALClusterizerFixedWindow.cxx:34
 AliEMCALClusterizerFixedWindow.cxx:35
 AliEMCALClusterizerFixedWindow.cxx:36
 AliEMCALClusterizerFixedWindow.cxx:37
 AliEMCALClusterizerFixedWindow.cxx:38
 AliEMCALClusterizerFixedWindow.cxx:39
 AliEMCALClusterizerFixedWindow.cxx:40
 AliEMCALClusterizerFixedWindow.cxx:41
 AliEMCALClusterizerFixedWindow.cxx:42
 AliEMCALClusterizerFixedWindow.cxx:43
 AliEMCALClusterizerFixedWindow.cxx:44
 AliEMCALClusterizerFixedWindow.cxx:45
 AliEMCALClusterizerFixedWindow.cxx:46
 AliEMCALClusterizerFixedWindow.cxx:47
 AliEMCALClusterizerFixedWindow.cxx:48
 AliEMCALClusterizerFixedWindow.cxx:49
 AliEMCALClusterizerFixedWindow.cxx:50
 AliEMCALClusterizerFixedWindow.cxx:51
 AliEMCALClusterizerFixedWindow.cxx:52
 AliEMCALClusterizerFixedWindow.cxx:53
 AliEMCALClusterizerFixedWindow.cxx:54
 AliEMCALClusterizerFixedWindow.cxx:55
 AliEMCALClusterizerFixedWindow.cxx:56
 AliEMCALClusterizerFixedWindow.cxx:57
 AliEMCALClusterizerFixedWindow.cxx:58
 AliEMCALClusterizerFixedWindow.cxx:59
 AliEMCALClusterizerFixedWindow.cxx:60
 AliEMCALClusterizerFixedWindow.cxx:61
 AliEMCALClusterizerFixedWindow.cxx:62
 AliEMCALClusterizerFixedWindow.cxx:63
 AliEMCALClusterizerFixedWindow.cxx:64
 AliEMCALClusterizerFixedWindow.cxx:65
 AliEMCALClusterizerFixedWindow.cxx:66
 AliEMCALClusterizerFixedWindow.cxx:67
 AliEMCALClusterizerFixedWindow.cxx:68
 AliEMCALClusterizerFixedWindow.cxx:69
 AliEMCALClusterizerFixedWindow.cxx:70
 AliEMCALClusterizerFixedWindow.cxx:71
 AliEMCALClusterizerFixedWindow.cxx:72
 AliEMCALClusterizerFixedWindow.cxx:73
 AliEMCALClusterizerFixedWindow.cxx:74
 AliEMCALClusterizerFixedWindow.cxx:75
 AliEMCALClusterizerFixedWindow.cxx:76
 AliEMCALClusterizerFixedWindow.cxx:77
 AliEMCALClusterizerFixedWindow.cxx:78
 AliEMCALClusterizerFixedWindow.cxx:79
 AliEMCALClusterizerFixedWindow.cxx:80
 AliEMCALClusterizerFixedWindow.cxx:81
 AliEMCALClusterizerFixedWindow.cxx:82
 AliEMCALClusterizerFixedWindow.cxx:83
 AliEMCALClusterizerFixedWindow.cxx:84
 AliEMCALClusterizerFixedWindow.cxx:85
 AliEMCALClusterizerFixedWindow.cxx:86
 AliEMCALClusterizerFixedWindow.cxx:87
 AliEMCALClusterizerFixedWindow.cxx:88
 AliEMCALClusterizerFixedWindow.cxx:89
 AliEMCALClusterizerFixedWindow.cxx:90
 AliEMCALClusterizerFixedWindow.cxx:91
 AliEMCALClusterizerFixedWindow.cxx:92
 AliEMCALClusterizerFixedWindow.cxx:93
 AliEMCALClusterizerFixedWindow.cxx:94
 AliEMCALClusterizerFixedWindow.cxx:95
 AliEMCALClusterizerFixedWindow.cxx:96
 AliEMCALClusterizerFixedWindow.cxx:97
 AliEMCALClusterizerFixedWindow.cxx:98
 AliEMCALClusterizerFixedWindow.cxx:99
 AliEMCALClusterizerFixedWindow.cxx:100
 AliEMCALClusterizerFixedWindow.cxx:101
 AliEMCALClusterizerFixedWindow.cxx:102
 AliEMCALClusterizerFixedWindow.cxx:103
 AliEMCALClusterizerFixedWindow.cxx:104
 AliEMCALClusterizerFixedWindow.cxx:105
 AliEMCALClusterizerFixedWindow.cxx:106
 AliEMCALClusterizerFixedWindow.cxx:107
 AliEMCALClusterizerFixedWindow.cxx:108
 AliEMCALClusterizerFixedWindow.cxx:109
 AliEMCALClusterizerFixedWindow.cxx:110
 AliEMCALClusterizerFixedWindow.cxx:111
 AliEMCALClusterizerFixedWindow.cxx:112
 AliEMCALClusterizerFixedWindow.cxx:113
 AliEMCALClusterizerFixedWindow.cxx:114
 AliEMCALClusterizerFixedWindow.cxx:115
 AliEMCALClusterizerFixedWindow.cxx:116
 AliEMCALClusterizerFixedWindow.cxx:117
 AliEMCALClusterizerFixedWindow.cxx:118
 AliEMCALClusterizerFixedWindow.cxx:119
 AliEMCALClusterizerFixedWindow.cxx:120
 AliEMCALClusterizerFixedWindow.cxx:121
 AliEMCALClusterizerFixedWindow.cxx:122
 AliEMCALClusterizerFixedWindow.cxx:123
 AliEMCALClusterizerFixedWindow.cxx:124
 AliEMCALClusterizerFixedWindow.cxx:125
 AliEMCALClusterizerFixedWindow.cxx:126
 AliEMCALClusterizerFixedWindow.cxx:127
 AliEMCALClusterizerFixedWindow.cxx:128
 AliEMCALClusterizerFixedWindow.cxx:129
 AliEMCALClusterizerFixedWindow.cxx:130
 AliEMCALClusterizerFixedWindow.cxx:131
 AliEMCALClusterizerFixedWindow.cxx:132
 AliEMCALClusterizerFixedWindow.cxx:133
 AliEMCALClusterizerFixedWindow.cxx:134
 AliEMCALClusterizerFixedWindow.cxx:135
 AliEMCALClusterizerFixedWindow.cxx:136
 AliEMCALClusterizerFixedWindow.cxx:137
 AliEMCALClusterizerFixedWindow.cxx:138
 AliEMCALClusterizerFixedWindow.cxx:139
 AliEMCALClusterizerFixedWindow.cxx:140
 AliEMCALClusterizerFixedWindow.cxx:141
 AliEMCALClusterizerFixedWindow.cxx:142
 AliEMCALClusterizerFixedWindow.cxx:143
 AliEMCALClusterizerFixedWindow.cxx:144
 AliEMCALClusterizerFixedWindow.cxx:145
 AliEMCALClusterizerFixedWindow.cxx:146
 AliEMCALClusterizerFixedWindow.cxx:147
 AliEMCALClusterizerFixedWindow.cxx:148
 AliEMCALClusterizerFixedWindow.cxx:149
 AliEMCALClusterizerFixedWindow.cxx:150
 AliEMCALClusterizerFixedWindow.cxx:151
 AliEMCALClusterizerFixedWindow.cxx:152
 AliEMCALClusterizerFixedWindow.cxx:153
 AliEMCALClusterizerFixedWindow.cxx:154
 AliEMCALClusterizerFixedWindow.cxx:155
 AliEMCALClusterizerFixedWindow.cxx:156
 AliEMCALClusterizerFixedWindow.cxx:157
 AliEMCALClusterizerFixedWindow.cxx:158
 AliEMCALClusterizerFixedWindow.cxx:159
 AliEMCALClusterizerFixedWindow.cxx:160
 AliEMCALClusterizerFixedWindow.cxx:161
 AliEMCALClusterizerFixedWindow.cxx:162
 AliEMCALClusterizerFixedWindow.cxx:163
 AliEMCALClusterizerFixedWindow.cxx:164
 AliEMCALClusterizerFixedWindow.cxx:165
 AliEMCALClusterizerFixedWindow.cxx:166
 AliEMCALClusterizerFixedWindow.cxx:167
 AliEMCALClusterizerFixedWindow.cxx:168
 AliEMCALClusterizerFixedWindow.cxx:169
 AliEMCALClusterizerFixedWindow.cxx:170
 AliEMCALClusterizerFixedWindow.cxx:171
 AliEMCALClusterizerFixedWindow.cxx:172
 AliEMCALClusterizerFixedWindow.cxx:173
 AliEMCALClusterizerFixedWindow.cxx:174
 AliEMCALClusterizerFixedWindow.cxx:175
 AliEMCALClusterizerFixedWindow.cxx:176
 AliEMCALClusterizerFixedWindow.cxx:177
 AliEMCALClusterizerFixedWindow.cxx:178
 AliEMCALClusterizerFixedWindow.cxx:179
 AliEMCALClusterizerFixedWindow.cxx:180
 AliEMCALClusterizerFixedWindow.cxx:181
 AliEMCALClusterizerFixedWindow.cxx:182
 AliEMCALClusterizerFixedWindow.cxx:183
 AliEMCALClusterizerFixedWindow.cxx:184
 AliEMCALClusterizerFixedWindow.cxx:185
 AliEMCALClusterizerFixedWindow.cxx:186
 AliEMCALClusterizerFixedWindow.cxx:187
 AliEMCALClusterizerFixedWindow.cxx:188
 AliEMCALClusterizerFixedWindow.cxx:189
 AliEMCALClusterizerFixedWindow.cxx:190
 AliEMCALClusterizerFixedWindow.cxx:191
 AliEMCALClusterizerFixedWindow.cxx:192
 AliEMCALClusterizerFixedWindow.cxx:193
 AliEMCALClusterizerFixedWindow.cxx:194
 AliEMCALClusterizerFixedWindow.cxx:195
 AliEMCALClusterizerFixedWindow.cxx:196
 AliEMCALClusterizerFixedWindow.cxx:197
 AliEMCALClusterizerFixedWindow.cxx:198
 AliEMCALClusterizerFixedWindow.cxx:199
 AliEMCALClusterizerFixedWindow.cxx:200
 AliEMCALClusterizerFixedWindow.cxx:201
 AliEMCALClusterizerFixedWindow.cxx:202
 AliEMCALClusterizerFixedWindow.cxx:203
 AliEMCALClusterizerFixedWindow.cxx:204
 AliEMCALClusterizerFixedWindow.cxx:205
 AliEMCALClusterizerFixedWindow.cxx:206
 AliEMCALClusterizerFixedWindow.cxx:207
 AliEMCALClusterizerFixedWindow.cxx:208
 AliEMCALClusterizerFixedWindow.cxx:209
 AliEMCALClusterizerFixedWindow.cxx:210
 AliEMCALClusterizerFixedWindow.cxx:211
 AliEMCALClusterizerFixedWindow.cxx:212
 AliEMCALClusterizerFixedWindow.cxx:213
 AliEMCALClusterizerFixedWindow.cxx:214
 AliEMCALClusterizerFixedWindow.cxx:215
 AliEMCALClusterizerFixedWindow.cxx:216
 AliEMCALClusterizerFixedWindow.cxx:217
 AliEMCALClusterizerFixedWindow.cxx:218
 AliEMCALClusterizerFixedWindow.cxx:219
 AliEMCALClusterizerFixedWindow.cxx:220
 AliEMCALClusterizerFixedWindow.cxx:221
 AliEMCALClusterizerFixedWindow.cxx:222
 AliEMCALClusterizerFixedWindow.cxx:223
 AliEMCALClusterizerFixedWindow.cxx:224
 AliEMCALClusterizerFixedWindow.cxx:225
 AliEMCALClusterizerFixedWindow.cxx:226
 AliEMCALClusterizerFixedWindow.cxx:227
 AliEMCALClusterizerFixedWindow.cxx:228
 AliEMCALClusterizerFixedWindow.cxx:229
 AliEMCALClusterizerFixedWindow.cxx:230
 AliEMCALClusterizerFixedWindow.cxx:231
 AliEMCALClusterizerFixedWindow.cxx:232
 AliEMCALClusterizerFixedWindow.cxx:233
 AliEMCALClusterizerFixedWindow.cxx:234
 AliEMCALClusterizerFixedWindow.cxx:235
 AliEMCALClusterizerFixedWindow.cxx:236
 AliEMCALClusterizerFixedWindow.cxx:237
 AliEMCALClusterizerFixedWindow.cxx:238
 AliEMCALClusterizerFixedWindow.cxx:239
 AliEMCALClusterizerFixedWindow.cxx:240
 AliEMCALClusterizerFixedWindow.cxx:241
 AliEMCALClusterizerFixedWindow.cxx:242
 AliEMCALClusterizerFixedWindow.cxx:243
 AliEMCALClusterizerFixedWindow.cxx:244
 AliEMCALClusterizerFixedWindow.cxx:245
 AliEMCALClusterizerFixedWindow.cxx:246
 AliEMCALClusterizerFixedWindow.cxx:247
 AliEMCALClusterizerFixedWindow.cxx:248
 AliEMCALClusterizerFixedWindow.cxx:249
 AliEMCALClusterizerFixedWindow.cxx:250
 AliEMCALClusterizerFixedWindow.cxx:251
 AliEMCALClusterizerFixedWindow.cxx:252
 AliEMCALClusterizerFixedWindow.cxx:253
 AliEMCALClusterizerFixedWindow.cxx:254
 AliEMCALClusterizerFixedWindow.cxx:255
 AliEMCALClusterizerFixedWindow.cxx:256
 AliEMCALClusterizerFixedWindow.cxx:257
 AliEMCALClusterizerFixedWindow.cxx:258
 AliEMCALClusterizerFixedWindow.cxx:259
 AliEMCALClusterizerFixedWindow.cxx:260
 AliEMCALClusterizerFixedWindow.cxx:261
 AliEMCALClusterizerFixedWindow.cxx:262
 AliEMCALClusterizerFixedWindow.cxx:263
 AliEMCALClusterizerFixedWindow.cxx:264
 AliEMCALClusterizerFixedWindow.cxx:265
 AliEMCALClusterizerFixedWindow.cxx:266
 AliEMCALClusterizerFixedWindow.cxx:267
 AliEMCALClusterizerFixedWindow.cxx:268
 AliEMCALClusterizerFixedWindow.cxx:269
 AliEMCALClusterizerFixedWindow.cxx:270
 AliEMCALClusterizerFixedWindow.cxx:271
 AliEMCALClusterizerFixedWindow.cxx:272
 AliEMCALClusterizerFixedWindow.cxx:273
 AliEMCALClusterizerFixedWindow.cxx:274
 AliEMCALClusterizerFixedWindow.cxx:275
 AliEMCALClusterizerFixedWindow.cxx:276
 AliEMCALClusterizerFixedWindow.cxx:277
 AliEMCALClusterizerFixedWindow.cxx:278
 AliEMCALClusterizerFixedWindow.cxx:279
 AliEMCALClusterizerFixedWindow.cxx:280
 AliEMCALClusterizerFixedWindow.cxx:281
 AliEMCALClusterizerFixedWindow.cxx:282
 AliEMCALClusterizerFixedWindow.cxx:283
 AliEMCALClusterizerFixedWindow.cxx:284
 AliEMCALClusterizerFixedWindow.cxx:285
 AliEMCALClusterizerFixedWindow.cxx:286
 AliEMCALClusterizerFixedWindow.cxx:287
 AliEMCALClusterizerFixedWindow.cxx:288
 AliEMCALClusterizerFixedWindow.cxx:289
 AliEMCALClusterizerFixedWindow.cxx:290
 AliEMCALClusterizerFixedWindow.cxx:291
 AliEMCALClusterizerFixedWindow.cxx:292
 AliEMCALClusterizerFixedWindow.cxx:293
 AliEMCALClusterizerFixedWindow.cxx:294
 AliEMCALClusterizerFixedWindow.cxx:295
 AliEMCALClusterizerFixedWindow.cxx:296
 AliEMCALClusterizerFixedWindow.cxx:297
 AliEMCALClusterizerFixedWindow.cxx:298
 AliEMCALClusterizerFixedWindow.cxx:299
 AliEMCALClusterizerFixedWindow.cxx:300
 AliEMCALClusterizerFixedWindow.cxx:301
 AliEMCALClusterizerFixedWindow.cxx:302
 AliEMCALClusterizerFixedWindow.cxx:303
 AliEMCALClusterizerFixedWindow.cxx:304
 AliEMCALClusterizerFixedWindow.cxx:305
 AliEMCALClusterizerFixedWindow.cxx:306
 AliEMCALClusterizerFixedWindow.cxx:307
 AliEMCALClusterizerFixedWindow.cxx:308
 AliEMCALClusterizerFixedWindow.cxx:309
 AliEMCALClusterizerFixedWindow.cxx:310
 AliEMCALClusterizerFixedWindow.cxx:311
 AliEMCALClusterizerFixedWindow.cxx:312
 AliEMCALClusterizerFixedWindow.cxx:313
 AliEMCALClusterizerFixedWindow.cxx:314
 AliEMCALClusterizerFixedWindow.cxx:315
 AliEMCALClusterizerFixedWindow.cxx:316
 AliEMCALClusterizerFixedWindow.cxx:317
 AliEMCALClusterizerFixedWindow.cxx:318
 AliEMCALClusterizerFixedWindow.cxx:319
 AliEMCALClusterizerFixedWindow.cxx:320
 AliEMCALClusterizerFixedWindow.cxx:321
 AliEMCALClusterizerFixedWindow.cxx:322
 AliEMCALClusterizerFixedWindow.cxx:323
 AliEMCALClusterizerFixedWindow.cxx:324
 AliEMCALClusterizerFixedWindow.cxx:325
 AliEMCALClusterizerFixedWindow.cxx:326
 AliEMCALClusterizerFixedWindow.cxx:327
 AliEMCALClusterizerFixedWindow.cxx:328
 AliEMCALClusterizerFixedWindow.cxx:329
 AliEMCALClusterizerFixedWindow.cxx:330
 AliEMCALClusterizerFixedWindow.cxx:331
 AliEMCALClusterizerFixedWindow.cxx:332
 AliEMCALClusterizerFixedWindow.cxx:333
 AliEMCALClusterizerFixedWindow.cxx:334
 AliEMCALClusterizerFixedWindow.cxx:335
 AliEMCALClusterizerFixedWindow.cxx:336
 AliEMCALClusterizerFixedWindow.cxx:337
 AliEMCALClusterizerFixedWindow.cxx:338
 AliEMCALClusterizerFixedWindow.cxx:339
 AliEMCALClusterizerFixedWindow.cxx:340
 AliEMCALClusterizerFixedWindow.cxx:341
 AliEMCALClusterizerFixedWindow.cxx:342
 AliEMCALClusterizerFixedWindow.cxx:343
 AliEMCALClusterizerFixedWindow.cxx:344
 AliEMCALClusterizerFixedWindow.cxx:345
 AliEMCALClusterizerFixedWindow.cxx:346
 AliEMCALClusterizerFixedWindow.cxx:347
 AliEMCALClusterizerFixedWindow.cxx:348
 AliEMCALClusterizerFixedWindow.cxx:349
 AliEMCALClusterizerFixedWindow.cxx:350
 AliEMCALClusterizerFixedWindow.cxx:351
 AliEMCALClusterizerFixedWindow.cxx:352
 AliEMCALClusterizerFixedWindow.cxx:353
 AliEMCALClusterizerFixedWindow.cxx:354
 AliEMCALClusterizerFixedWindow.cxx:355
 AliEMCALClusterizerFixedWindow.cxx:356
 AliEMCALClusterizerFixedWindow.cxx:357
 AliEMCALClusterizerFixedWindow.cxx:358
 AliEMCALClusterizerFixedWindow.cxx:359
 AliEMCALClusterizerFixedWindow.cxx:360
 AliEMCALClusterizerFixedWindow.cxx:361
 AliEMCALClusterizerFixedWindow.cxx:362
 AliEMCALClusterizerFixedWindow.cxx:363
 AliEMCALClusterizerFixedWindow.cxx:364
 AliEMCALClusterizerFixedWindow.cxx:365
 AliEMCALClusterizerFixedWindow.cxx:366
 AliEMCALClusterizerFixedWindow.cxx:367
 AliEMCALClusterizerFixedWindow.cxx:368
 AliEMCALClusterizerFixedWindow.cxx:369
 AliEMCALClusterizerFixedWindow.cxx:370
 AliEMCALClusterizerFixedWindow.cxx:371
 AliEMCALClusterizerFixedWindow.cxx:372
 AliEMCALClusterizerFixedWindow.cxx:373
 AliEMCALClusterizerFixedWindow.cxx:374
 AliEMCALClusterizerFixedWindow.cxx:375
 AliEMCALClusterizerFixedWindow.cxx:376
 AliEMCALClusterizerFixedWindow.cxx:377
 AliEMCALClusterizerFixedWindow.cxx:378
 AliEMCALClusterizerFixedWindow.cxx:379
 AliEMCALClusterizerFixedWindow.cxx:380
 AliEMCALClusterizerFixedWindow.cxx:381
 AliEMCALClusterizerFixedWindow.cxx:382
 AliEMCALClusterizerFixedWindow.cxx:383
 AliEMCALClusterizerFixedWindow.cxx:384
 AliEMCALClusterizerFixedWindow.cxx:385
 AliEMCALClusterizerFixedWindow.cxx:386
 AliEMCALClusterizerFixedWindow.cxx:387
 AliEMCALClusterizerFixedWindow.cxx:388
 AliEMCALClusterizerFixedWindow.cxx:389
 AliEMCALClusterizerFixedWindow.cxx:390
 AliEMCALClusterizerFixedWindow.cxx:391
 AliEMCALClusterizerFixedWindow.cxx:392
 AliEMCALClusterizerFixedWindow.cxx:393
 AliEMCALClusterizerFixedWindow.cxx:394
 AliEMCALClusterizerFixedWindow.cxx:395
 AliEMCALClusterizerFixedWindow.cxx:396
 AliEMCALClusterizerFixedWindow.cxx:397
 AliEMCALClusterizerFixedWindow.cxx:398
 AliEMCALClusterizerFixedWindow.cxx:399
 AliEMCALClusterizerFixedWindow.cxx:400
 AliEMCALClusterizerFixedWindow.cxx:401
 AliEMCALClusterizerFixedWindow.cxx:402
 AliEMCALClusterizerFixedWindow.cxx:403
 AliEMCALClusterizerFixedWindow.cxx:404
 AliEMCALClusterizerFixedWindow.cxx:405
 AliEMCALClusterizerFixedWindow.cxx:406
 AliEMCALClusterizerFixedWindow.cxx:407
 AliEMCALClusterizerFixedWindow.cxx:408
 AliEMCALClusterizerFixedWindow.cxx:409
 AliEMCALClusterizerFixedWindow.cxx:410
 AliEMCALClusterizerFixedWindow.cxx:411
 AliEMCALClusterizerFixedWindow.cxx:412
 AliEMCALClusterizerFixedWindow.cxx:413
 AliEMCALClusterizerFixedWindow.cxx:414
 AliEMCALClusterizerFixedWindow.cxx:415
 AliEMCALClusterizerFixedWindow.cxx:416
 AliEMCALClusterizerFixedWindow.cxx:417
 AliEMCALClusterizerFixedWindow.cxx:418
 AliEMCALClusterizerFixedWindow.cxx:419
 AliEMCALClusterizerFixedWindow.cxx:420
 AliEMCALClusterizerFixedWindow.cxx:421
 AliEMCALClusterizerFixedWindow.cxx:422
 AliEMCALClusterizerFixedWindow.cxx:423
 AliEMCALClusterizerFixedWindow.cxx:424
 AliEMCALClusterizerFixedWindow.cxx:425
 AliEMCALClusterizerFixedWindow.cxx:426
 AliEMCALClusterizerFixedWindow.cxx:427
 AliEMCALClusterizerFixedWindow.cxx:428
 AliEMCALClusterizerFixedWindow.cxx:429
 AliEMCALClusterizerFixedWindow.cxx:430
 AliEMCALClusterizerFixedWindow.cxx:431
 AliEMCALClusterizerFixedWindow.cxx:432
 AliEMCALClusterizerFixedWindow.cxx:433
 AliEMCALClusterizerFixedWindow.cxx:434
 AliEMCALClusterizerFixedWindow.cxx:435
 AliEMCALClusterizerFixedWindow.cxx:436
 AliEMCALClusterizerFixedWindow.cxx:437
 AliEMCALClusterizerFixedWindow.cxx:438
 AliEMCALClusterizerFixedWindow.cxx:439
 AliEMCALClusterizerFixedWindow.cxx:440
 AliEMCALClusterizerFixedWindow.cxx:441
 AliEMCALClusterizerFixedWindow.cxx:442
 AliEMCALClusterizerFixedWindow.cxx:443
 AliEMCALClusterizerFixedWindow.cxx:444
 AliEMCALClusterizerFixedWindow.cxx:445
 AliEMCALClusterizerFixedWindow.cxx:446
 AliEMCALClusterizerFixedWindow.cxx:447
 AliEMCALClusterizerFixedWindow.cxx:448
 AliEMCALClusterizerFixedWindow.cxx:449
 AliEMCALClusterizerFixedWindow.cxx:450
 AliEMCALClusterizerFixedWindow.cxx:451
 AliEMCALClusterizerFixedWindow.cxx:452
 AliEMCALClusterizerFixedWindow.cxx:453
 AliEMCALClusterizerFixedWindow.cxx:454
 AliEMCALClusterizerFixedWindow.cxx:455
 AliEMCALClusterizerFixedWindow.cxx:456
 AliEMCALClusterizerFixedWindow.cxx:457
 AliEMCALClusterizerFixedWindow.cxx:458
 AliEMCALClusterizerFixedWindow.cxx:459
 AliEMCALClusterizerFixedWindow.cxx:460
 AliEMCALClusterizerFixedWindow.cxx:461
 AliEMCALClusterizerFixedWindow.cxx:462
 AliEMCALClusterizerFixedWindow.cxx:463
 AliEMCALClusterizerFixedWindow.cxx:464
 AliEMCALClusterizerFixedWindow.cxx:465
 AliEMCALClusterizerFixedWindow.cxx:466
 AliEMCALClusterizerFixedWindow.cxx:467
 AliEMCALClusterizerFixedWindow.cxx:468
 AliEMCALClusterizerFixedWindow.cxx:469
 AliEMCALClusterizerFixedWindow.cxx:470
 AliEMCALClusterizerFixedWindow.cxx:471
 AliEMCALClusterizerFixedWindow.cxx:472
 AliEMCALClusterizerFixedWindow.cxx:473
 AliEMCALClusterizerFixedWindow.cxx:474
 AliEMCALClusterizerFixedWindow.cxx:475
 AliEMCALClusterizerFixedWindow.cxx:476
 AliEMCALClusterizerFixedWindow.cxx:477
 AliEMCALClusterizerFixedWindow.cxx:478
 AliEMCALClusterizerFixedWindow.cxx:479
 AliEMCALClusterizerFixedWindow.cxx:480
 AliEMCALClusterizerFixedWindow.cxx:481
 AliEMCALClusterizerFixedWindow.cxx:482
 AliEMCALClusterizerFixedWindow.cxx:483
 AliEMCALClusterizerFixedWindow.cxx:484
 AliEMCALClusterizerFixedWindow.cxx:485
 AliEMCALClusterizerFixedWindow.cxx:486
 AliEMCALClusterizerFixedWindow.cxx:487
 AliEMCALClusterizerFixedWindow.cxx:488
 AliEMCALClusterizerFixedWindow.cxx:489
 AliEMCALClusterizerFixedWindow.cxx:490
 AliEMCALClusterizerFixedWindow.cxx:491
 AliEMCALClusterizerFixedWindow.cxx:492
 AliEMCALClusterizerFixedWindow.cxx:493
 AliEMCALClusterizerFixedWindow.cxx:494
 AliEMCALClusterizerFixedWindow.cxx:495
 AliEMCALClusterizerFixedWindow.cxx:496
 AliEMCALClusterizerFixedWindow.cxx:497
 AliEMCALClusterizerFixedWindow.cxx:498
 AliEMCALClusterizerFixedWindow.cxx:499
 AliEMCALClusterizerFixedWindow.cxx:500
 AliEMCALClusterizerFixedWindow.cxx:501
 AliEMCALClusterizerFixedWindow.cxx:502
 AliEMCALClusterizerFixedWindow.cxx:503
 AliEMCALClusterizerFixedWindow.cxx:504
 AliEMCALClusterizerFixedWindow.cxx:505
 AliEMCALClusterizerFixedWindow.cxx:506
 AliEMCALClusterizerFixedWindow.cxx:507
 AliEMCALClusterizerFixedWindow.cxx:508
 AliEMCALClusterizerFixedWindow.cxx:509
 AliEMCALClusterizerFixedWindow.cxx:510
 AliEMCALClusterizerFixedWindow.cxx:511
 AliEMCALClusterizerFixedWindow.cxx:512
 AliEMCALClusterizerFixedWindow.cxx:513
 AliEMCALClusterizerFixedWindow.cxx:514
 AliEMCALClusterizerFixedWindow.cxx:515
 AliEMCALClusterizerFixedWindow.cxx:516
 AliEMCALClusterizerFixedWindow.cxx:517
 AliEMCALClusterizerFixedWindow.cxx:518
 AliEMCALClusterizerFixedWindow.cxx:519
 AliEMCALClusterizerFixedWindow.cxx:520
 AliEMCALClusterizerFixedWindow.cxx:521
 AliEMCALClusterizerFixedWindow.cxx:522
 AliEMCALClusterizerFixedWindow.cxx:523
 AliEMCALClusterizerFixedWindow.cxx:524
 AliEMCALClusterizerFixedWindow.cxx:525
 AliEMCALClusterizerFixedWindow.cxx:526