ROOT logo
//
// This class calculates the exclusive charged particle density
// in each for the 5 FMD rings. 
//
#include "AliFMDCorrector.h"
#include <AliESDFMD.h>
#include <TAxis.h>
#include <TList.h>
#include <TMath.h>
#include "AliForwardCorrectionManager.h"
#include "AliFMDCorrSecondaryMap.h"
#include "AliFMDCorrVertexBias.h"
#include "AliFMDCorrMergingEfficiency.h"
#include "AliFMDCorrAcceptance.h"
#include "AliLog.h"
#include <TH2D.h>
#include <TROOT.h>
#include <THStack.h>
#include <iostream>
#include <iomanip>

ClassImp(AliFMDCorrector)
#if 0
; // For Emacs
#endif 

//____________________________________________________________________
AliFMDCorrector::AliFMDCorrector()
  : TNamed(), 
    fRingHistos(),
    fUseSecondaryMap(true),
    fUseVertexBias(false),
    fUseAcceptance(false),
    fUseMergingEfficiency(false),
    fDebug(0)
{
  // Constructor
  DGUARD(fDebug, 3, "Default CTOR of AliFMDCorrector");
}

//____________________________________________________________________
AliFMDCorrector::AliFMDCorrector(const char* title)
  : TNamed("fmdCorrector", title), 
    fRingHistos(), 
    fUseSecondaryMap(true),
    fUseVertexBias(false),
    fUseAcceptance(false),
    fUseMergingEfficiency(false),
    fDebug(0)
{
  // Constructor 
  // 
  // Parameters: 
  //   title   Title
  DGUARD(fDebug, 3, "Named CTOR of AliFMDCorrector: %s", title);
  fRingHistos.SetName(GetName());
  fRingHistos.Add(new RingHistos(1, 'I'));
  fRingHistos.Add(new RingHistos(2, 'I'));
  fRingHistos.Add(new RingHistos(2, 'O'));
  fRingHistos.Add(new RingHistos(3, 'I'));
  fRingHistos.Add(new RingHistos(3, 'O'));
}

//____________________________________________________________________
AliFMDCorrector::AliFMDCorrector(const AliFMDCorrector& o)
  : TNamed(o), 
    fRingHistos(), 
    fUseSecondaryMap(o.fUseSecondaryMap),
    fUseVertexBias(o.fUseVertexBias),
    fUseAcceptance(o.fUseAcceptance),
    fUseMergingEfficiency(o.fUseMergingEfficiency),
    fDebug(o.fDebug)
{
  // Copy constructor 
  // 
  // Parameters: 
  //  o  Object to copy from 
  DGUARD(fDebug, 3, "Copy CTOR of AliFMDCorrector");
  TIter    next(&o.fRingHistos);
  TObject* obj = 0;
  while ((obj = next())) fRingHistos.Add(obj);
}

//____________________________________________________________________
AliFMDCorrector::~AliFMDCorrector()
{
  // Destructor 
  // 
  // 
  DGUARD(fDebug, 3, "DTOR of AliFMDCorrector");
  fRingHistos.Delete();
}

//____________________________________________________________________
AliFMDCorrector&
AliFMDCorrector::operator=(const AliFMDCorrector& o)
{
  // Assignment operator 
  // 
  // Parameters:
  //   o   Object to assign from 
  DGUARD(fDebug, 3, "Assignment of AliFMDCorrector");
  if (&o == this) return *this; 
  TNamed::operator=(o);

  fDebug   = o.fDebug;
  fRingHistos.Delete();
  fUseSecondaryMap = o.fUseSecondaryMap;
  fUseVertexBias = o.fUseVertexBias;
  fUseAcceptance = o.fUseAcceptance;
  fUseMergingEfficiency = o.fUseMergingEfficiency;
  TIter    next(&o.fRingHistos);
  TObject* obj = 0;
  while ((obj = next())) fRingHistos.Add(obj);
  
  return *this;
}

//____________________________________________________________________
void
AliFMDCorrector::SetupForData(const TAxis&)
{
  //
  // Initialize this object
  //
  // Parameters:
  //    etaAxis Eta axis to use
  //
  DGUARD(fDebug, 1, "Initialization of AliFMDCorrector");
  if (!fUseSecondaryMap)
    AliWarning("Secondary maps not used - BE CAREFUL");
  if (!fUseVertexBias)
    AliWarning("Vertex bias not used");
  if (!fUseAcceptance)
    AliWarning("Acceptance from dead-channels not used");
}

//____________________________________________________________________
AliFMDCorrector::RingHistos*
AliFMDCorrector::GetRingHistos(UShort_t d, Char_t r) const
{
  // 
  // Get the ring histogram container 
  // Parameters:
  //    d Detector
  //    r Ring 
  // 
  // Return:
  //    Ring histogram container 
  //
  Int_t idx = -1;
  switch (d) { 
  case 1: idx = 0; break;
  case 2: idx = 1 + (r == 'I' || r == 'i' ? 0 : 1); break;
  case 3: idx = 3 + (r == 'I' || r == 'i' ? 0 : 1); break;
  }
  if (idx < 0 || idx >= fRingHistos.GetEntries()) return 0;
  
  return static_cast<RingHistos*>(fRingHistos.At(idx));
}
    
//____________________________________________________________________
void
AliFMDCorrector::DivideMap(TH2* num, const TH2* denom,
			   Bool_t alsoUnderOver) const
{
  // 
  // Implement TH1::Divide but 
  // - Assume compatible histograms 
  // - Unless third argument is true, do not divide over/under flow bins
  // 
  if (!num || !denom) return;

  Int_t first = (alsoUnderOver ? 0 : 1);
  Int_t lastX = num->GetNbinsX() + (alsoUnderOver ? 1 : 0);
  Int_t lastY = num->GetNbinsY() + (alsoUnderOver ? 1 : 0);
  
  for (Int_t ix = first; ix <= lastX; ix++) {
    for (Int_t iy = first; iy <= lastY; iy++) { 
      Int_t    bin = num->GetBin(ix,iy);
      Double_t c0  = num->GetBinContent(bin);
      Double_t c1  = denom->GetBinContent(bin);
      if (!c1) { 
	num->SetBinContent(bin,0);
	num->SetBinError(bin, 0);
	continue;
      }
      Double_t w   = c0 / c1;
      Double_t e0  = num->GetBinError(bin);
      Double_t e1  = denom->GetBinError(bin);
      Double_t c12 = c1*c1;
      Double_t e2  = (e0*e0*c1*c1 + e1*e1*c0*c0)/(c12*c12);
      
      num->SetBinContent(bin, w);
      num->SetBinError(bin, TMath::Sqrt(e2));
    }
  }
}
//____________________________________________________________________
Bool_t
AliFMDCorrector::Correct(AliForwardUtil::Histos& hists,
			 UShort_t                vtxbin)
{
  // 
  // Do the calculations 
  // Parameters:
  //    hists    Cache of histograms 
  //    vtxBin   Vertex bin 
  // 
  // Return:
  //    true on successs 
  //
  DGUARD(fDebug, 1, "Correct histograms of AliFMDCorrector");
  AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();

  UShort_t uvb = vtxbin;
  for (UShort_t d=1; d<=3; d++) { 
    UShort_t nr = (d == 1 ? 1 : 2);
    for (UShort_t q=0; q<nr; q++) { 
      Char_t      r  = (q == 0 ? 'I' : 'O');
      TH2D*       h  = hists.Get(d,r);
      RingHistos* rh = GetRingHistos(d,r);

      if (fUseSecondaryMap) {
        TH2D*  bg = fcm.GetSecondaryMap()->GetCorrection(d, r, uvb);
        if (!bg) {
          AliWarning(Form("No secondary correction for FMDM%d%c "
			  "in vertex bin %d", d, r, uvb));
          continue;
        }
        // Divide by primary/total ratio
	DivideMap(h, bg, false);
      }
      if (fUseVertexBias) {
        TH2D*  ef = fcm.GetVertexBias()->GetCorrection(r, uvb);
        if (!ef) {
          AliWarning(Form("No event %s vertex bias correction in vertex bin %d",
                          (r == 'I' || r == 'i' ? "inner" : "outer"), uvb));
          continue;
        }
        // Divide by the event selection efficiency
	DivideMap(h, ef, false);
      }
      if (fUseAcceptance) {
        TH2D*  ac = fcm.GetAcceptance()->GetCorrection(d, r, uvb);
        if (!ac) {
          AliWarning(Form("No acceptance correction for FMD%d%c in "
			  "vertex bin %d", d, r, uvb));
          continue;
        }
	// Fill overflow bin with ones 
	for (Int_t i = 1; i <= h->GetNbinsX(); i++) 
	  h->SetBinContent(i, h->GetNbinsY()+1, 1);

        // Divide by the acceptance correction - 
	DivideMap(h, ac, fcm.GetAcceptance()->HasOverflow());
      }

      if (fUseMergingEfficiency) {
	if (!fcm.GetMergingEfficiency()) { 
	  AliWarning("No merging efficiencies");
	  continue;
	}
	TH1D* sf = fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb);
	if (!fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb)) { 
	  AliWarning(Form("No merging efficiency for FMD%d%c at vertex bin %d",
			  d, r, uvb));
	  continue;
	}
	
      
	for (Int_t ieta = 1; ieta <= h->GetNbinsX(); ieta++) {
	  Float_t c  = sf->GetBinContent(ieta);
	  Float_t ec = sf->GetBinError(ieta);
	  	  
	  for (Int_t iphi = 1; iphi <= h->GetNbinsY(); iphi++) { 
	    if (c == 0) {
	      h->SetBinContent(ieta,iphi,0);
	      h->SetBinError(ieta,iphi,0);
	      continue;
	    }

	    Double_t m  = h->GetBinContent(ieta, iphi) / c;
	    Double_t em = h->GetBinError(ieta, iphi);
	  
	    Double_t e  = TMath::Sqrt(em * em + (m * ec) * (m * ec));
	    
	    h->SetBinContent(ieta,iphi,m);
	    h->SetBinError(ieta,iphi,e);
	  }
	}
      }
      
      rh->fDensity->Add(h);
    }
  }
  
  return kTRUE;
}

//____________________________________________________________________
void
AliFMDCorrector::Terminate(const TList* dir, TList* output, Int_t nEvents)
{
  // 
  // Scale the histograms to the total number of events 
  // Parameters:
  //    dir     Where the output is stored
  //    nEvents Number of events 
  //
  DGUARD(fDebug, 1, "Scale histograms of AliFMDCorrector");
  if (nEvents <= 0) return;
  TList* d = static_cast<TList*>(dir->FindObject(GetName()));
  if (!d) return;

  TList* out = new TList;
  out->SetName(d->GetName());
  out->SetOwner();

  TIter    next(&fRingHistos);
  RingHistos* o = 0;
  THStack* sums = new THStack("sums", "Sums of ring results");
  while ((o = static_cast<RingHistos*>(next()))) {
    o->Terminate(d, nEvents);
    if (!o->fDensity) { 
      Warning("Terminate", "No density from %s", o->GetName());
      continue;
    }
    TH1D* sum = o->fDensity->ProjectionX(o->GetName(), 1, 
					 o->fDensity->GetNbinsY(),"e");
    sum->Scale(1., "width");
    sum->SetTitle(o->GetName());
    sum->SetDirectory(0);
    sum->SetYTitle("#sum N_{ch,primary}");
    sums->Add(sum);
  }
  out->Add(sums);
  output->Add(out);
}
//____________________________________________________________________
void
AliFMDCorrector::CreateOutputObjects(TList* dir)
{
  // 
  // Output diagnostic histograms to directory 
  // 
  // Parameters:
  //    dir List to write in
  //  
  DGUARD(fDebug, 1, "Define output of AliFMDCorrector");
  TList* d = new TList;
  d->SetName(GetName());
  dir->Add(d);

  d->Add(AliForwardUtil::MakeParameter("secondary", fUseSecondaryMap));
  d->Add(AliForwardUtil::MakeParameter("vertexBias", fUseVertexBias));
  d->Add(AliForwardUtil::MakeParameter("acceptance", fUseAcceptance));
  d->Add(AliForwardUtil::MakeParameter("merging", fUseMergingEfficiency));
  

  TIter    next(&fRingHistos);
  RingHistos* o = 0;
  while ((o = static_cast<RingHistos*>(next()))) {
    o->CreateOutputObjects(d);
  }
}

#define PF(N,V,...)					\
  AliForwardUtil::PrintField(N,V, ## __VA_ARGS__)
#define PFB(N,FLAG)				\
  do {									\
    AliForwardUtil::PrintName(N);					\
    std::cout << std::boolalpha << (FLAG) << std::noboolalpha << std::endl; \
  } while(false)
#define PFV(N,VALUE)					\
  do {							\
    AliForwardUtil::PrintName(N);			\
    std::cout << (VALUE) << std::endl; } while(false)

//____________________________________________________________________
void
AliFMDCorrector::Print(Option_t* /* option */) const
{
  // 
  // Print information
  // Parameters:
  //    option Not used 
  //
  AliForwardUtil::PrintTask(*this);
  gROOT->IncreaseDirLevel();
  PFB("Use secondary maps",		fUseSecondaryMap );
  PFB("Use vertex bias",		fUseVertexBias );
  PFB("Use acceptance",			fUseAcceptance );
  PFB("Use merging efficiency",		fUseMergingEfficiency);
  gROOT->DecreaseDirLevel();  
}

//====================================================================
AliFMDCorrector::RingHistos::RingHistos()
  : AliForwardUtil::RingHistos(), 
    fDensity(0)
{
  // Constructor 
  //
  // 
}

//____________________________________________________________________
AliFMDCorrector::RingHistos::RingHistos(UShort_t d, Char_t r)
  : AliForwardUtil::RingHistos(d,r), 
    fDensity(0)
{
  // 
  // Constructor
  // Parameters:
  //    d detector
  //    r ring 
  //
  fDensity = new TH2D("primaryDensity", 
		      "#sum N_{ch,primary}/(#Delta#eta#Delta#phi)", 
		      200, -4, 6, (r == 'I' || r == 'i' ? 20 : 40), 
		      0, 2*TMath::Pi());
  fDensity->SetDirectory(0);
  fDensity->SetMarkerColor(Color());
  fDensity->Sumw2();
  fDensity->SetXTitle("#eta");
  fDensity->SetYTitle("#phi [radians]");
  fDensity->SetZTitle("Primary N_{ch} density");
}
//____________________________________________________________________
AliFMDCorrector::RingHistos::RingHistos(const RingHistos& o)
  : AliForwardUtil::RingHistos(o), 
    fDensity(o.fDensity)
{
  // 
  // Copy constructor 
  // Parameters:
  //    o Object to copy from 
  //
}

//____________________________________________________________________
AliFMDCorrector::RingHistos&
AliFMDCorrector::RingHistos::operator=(const RingHistos& o)
{
  // 
  // Assignment operator 
  // Parameters:
  //    o Object to assign from 
  // 
  // Return:
  //    Reference to this 
  //
  if (&o == this) return *this; 
  AliForwardUtil::RingHistos::operator=(o);
  
  if (fDensity) delete fDensity;
  
  fDensity = static_cast<TH2D*>(o.fDensity->Clone());
  
  return *this;
}
//____________________________________________________________________
AliFMDCorrector::RingHistos::~RingHistos()
{
  // 
  // Destructor 
  //
  // if (fDensity) delete fDensity;
}

//____________________________________________________________________
void
AliFMDCorrector::RingHistos::CreateOutputObjects(TList* dir)
{
  // 
  // Make output 
  // Parameters:
  //    dir Where to put it 
  //
  TList* d = DefineOutputList(dir);
  d->Add(fDensity);
}

//____________________________________________________________________
void
AliFMDCorrector::RingHistos::Terminate(TList* dir, Int_t nEvents)
{ 
  // 
  // Scale the histograms to the total number of events 
  // Parameters:
  //    dir     where the output is stored
  //    nEvents Number of events 
  //
  TList* l = GetOutputList(dir);
  if (!l) return; 

  TH2D* density = static_cast<TH2D*>(GetOutputHist(l,"primaryDensity"));
  if (density) density->Scale(1./nEvents);
  fDensity = density;
}

//____________________________________________________________________
//
// EOF
//
	  


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