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.                  *
**************************************************************************/
//
// HFE correction framework container
// Contains many single containers
// Extra fuctionality like appending added
//
// Author:
//   Markus Fasel <M.Fasel@gsi.de>
//
#include <iostream>
#include <TAxis.h>
#include <TClass.h>
#include <TCollection.h>
#include <THashList.h>
#include <THnSparse.h>
#include <TList.h>
#include <TObjArray.h>
#include <TObjString.h>
#include <TString.h>

#include "AliCFContainer.h"
#include "AliCFGridSparse.h"
#include "AliHFEcontainer.h"
#include "AliHFEtools.h"

ClassImp(AliHFEcontainer)
ClassImp(AliHFEcontainer::AliHFEvarInfo)

//__________________________________________________________________
AliHFEcontainer::AliHFEcontainer():
  TNamed("HFEcontainer", ""),
  fContainers(NULL),
  fCorrelationMatrices(NULL),
  fVariables(NULL),
  fNVars(0),
  fNEvents(0)
{
  //
  // Default constructor
  //
}

//__________________________________________________________________
AliHFEcontainer::AliHFEcontainer(const Char_t *name):
  TNamed(name, ""),
  fContainers(NULL),
  fCorrelationMatrices(NULL),
  fVariables(NULL),
  fNVars(0),
  fNEvents(0)
{
  //
  // Default constructor
  //
  fContainers = new THashList();
  fContainers->SetOwner();
}

//__________________________________________________________________
AliHFEcontainer::AliHFEcontainer(const Char_t *name, UInt_t nVar):
  TNamed(name, ""),
  fContainers(NULL),
  fCorrelationMatrices(NULL),
  fVariables(NULL),
  fNVars(0),
  fNEvents(0)
{
  //
  // Constructor
  // Setting Number of Variables too
  //
  fContainers = new THashList();
  fContainers->SetOwner();
  SetNumberOfVariables(nVar);
}

//__________________________________________________________________
AliHFEcontainer::AliHFEcontainer(const AliHFEcontainer &ref):
  TNamed(ref),
  fContainers(NULL),
  fCorrelationMatrices(NULL),
  fVariables(NULL),
  fNVars(ref.fNVars),
  fNEvents(ref.fNEvents)
{
  //
  // Copy constructor
  // creates a new object with new (empty) containers
  //
  if(fNVars){
    fVariables = new TObjArray(fNVars);
    AliHFEvarInfo *vtmp = NULL;
    for(UInt_t ivar = 0; ivar < fNVars; ivar++){
      vtmp = static_cast<AliHFEvarInfo *>(ref.fVariables->UncheckedAt(ivar));
      fVariables->AddAt(new AliHFEvarInfo(*vtmp), ivar);
    }
  }
  fContainers = new THashList;
  fContainers->SetOwner();
  AliCFContainer *ctmp = NULL;
  for(Int_t ien = 0; ien < ref.fContainers->GetEntries(); ien++){
    ctmp = static_cast<AliCFContainer *>(ref.fContainers->At(ien));
    CreateContainer(ctmp->GetName(), ctmp->GetTitle(), ctmp->GetNStep());
  }
  // Copy also correlation matrices
  if(ref.fCorrelationMatrices){
    THnSparseF *htmp = NULL;
    fCorrelationMatrices = new THashList;
    fCorrelationMatrices->SetOwner();
    for(Int_t ien = 0; ien < ref.fCorrelationMatrices->GetEntries(); ien++){
      htmp = static_cast<THnSparseF *>(ref.fCorrelationMatrices->At(ien));
      CreateCorrelationMatrix(htmp->GetName(), htmp->GetTitle());
    }
  }
}

//__________________________________________________________________
AliHFEcontainer &AliHFEcontainer::operator=(const AliHFEcontainer &ref){
  //
  // Assignment operator
  // Cleanup old object, create a new one with new containers inside
  //
  if(this == &ref) return *this;
  this->~AliHFEcontainer(); // cleanup old object before creating the new onwe
  TNamed::operator=(ref);
  fContainers = new THashList();
  fCorrelationMatrices = NULL;
  fNVars = ref.fNVars;
  if(fNVars){
    fVariables = new TObjArray(fNVars);
    AliHFEvarInfo *vtmp = NULL;
    for(UInt_t ivar = 0; ivar < fNVars; ivar++){
      vtmp = static_cast<AliHFEvarInfo *>(ref.fVariables->UncheckedAt(ivar));
      fVariables->AddAt(new AliHFEvarInfo(*vtmp), ivar);
    }
  } else {
    // No varible defined, do not try to copy anything
    fVariables = NULL;
    return *this;
  }

  // Reference contains content, try copying also the containers and the correlation matrices
  fContainers = new THashList();
  AliCFContainer *ctmp = NULL;
  for(Int_t ien = 0; ien < ref.fContainers->GetEntries(); ien++){
    ctmp = static_cast<AliCFContainer *>(ref.fContainers->At(ien));
    fContainers->Add(new AliCFContainer(*ctmp));
  }
  // Copy also correlation matrices
  if(ref.fCorrelationMatrices){
    THnSparseF *htmp = NULL;
    fCorrelationMatrices = new THashList;
    fCorrelationMatrices->SetOwner();
    for(Int_t ien = 0; ien < ref.fCorrelationMatrices->GetEntries(); ien++){
      htmp = static_cast<THnSparseF *>(ref.fCorrelationMatrices->At(ien));
      CreateCorrelationMatrix(htmp->GetName(), htmp->GetTitle());
    }
  }
  return *this;
}

//__________________________________________________________________
AliHFEcontainer::~AliHFEcontainer(){
  //
  // Destructor
  //
  delete fContainers;
  if(fCorrelationMatrices) delete fCorrelationMatrices;
  if(fVariables){
    fVariables->Delete();
    delete fVariables;
  }
}

//__________________________________________________________________
Long64_t AliHFEcontainer::Merge(TCollection *coll){
  //
  // Merge Container
  //
  if(!coll)
    return 0;
  if(coll->IsEmpty())
    return 1;

  TIter iter(coll);
  TObject *o = NULL;
  Long64_t count = 0;
  while((o = iter())){
    AliHFEcontainer *cont = dynamic_cast<AliHFEcontainer *>(o);
    if(!cont) continue;

    // Merge the two TObjArrays
    TList containers;
    containers.Add(cont->fContainers);
    fContainers->Merge(&containers);

    if(fCorrelationMatrices && cont->fCorrelationMatrices){
      containers.Clear();
      containers.Add(cont->fCorrelationMatrices);
      fCorrelationMatrices->Merge(&containers);
    }

    fNEvents += cont->GetNumberOfEvents();
    count++;
  }
  return count + 1;
}

//__________________________________________________________________
void AliHFEcontainer::SetNumberOfVariables(UInt_t nVar){
  //
  // Define the number of variables 
  // Initialize containers for the variable informations
  //
  if(fNVars) return;

  fNVars = nVar;
  fVariables = new TObjArray(nVar);
  for(UInt_t ivar = 0; ivar < nVar; ivar++)
    fVariables->AddAt(new AliHFEvarInfo, ivar);
}

//__________________________________________________________________
void AliHFEcontainer::CreateContainer(const Char_t *name, const Char_t *title, UInt_t nStep){
  //
  // Create a new Correction Framework Container and store it 
  //
  if(fContainers->FindObject(name)){
    AliError(Form("Container %s already exists. Cannot replace it!", name));
    return;
  }
  
  Int_t *nBins = new Int_t[fNVars];
  AliHFEvarInfo *var = NULL;
  for(UInt_t ivar = 0; ivar < fNVars; ivar++){ 
    var = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(ivar));
    nBins[ivar] = var ? var->GetNumberOfBins() : 0;
  }
  AliCFContainer *cont = new AliCFContainer(name, title, nStep, fNVars, nBins);
  for(UInt_t ivar = 0; ivar < fNVars; ivar++){
    var = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(ivar));
    if(var){
      cont->SetBinLimits(ivar, var->GetBinning());
      cont->SetVarTitle(ivar, var->GetVarName()->Data());
    }
  }
  delete[] nBins;
  fContainers->Add(cont);
  AliInfo(Form("Container %s created with %d cut steps", name, nStep));
}

//__________________________________________________________________
void AliHFEcontainer::CreateCorrelationMatrix(const Char_t *name, const Char_t *title){
  //
  // Create Correlation Matrix
  //
  if(!fCorrelationMatrices){
    fCorrelationMatrices = new THashList;
    fCorrelationMatrices->SetName("fCorrelationMatrices");
    fCorrelationMatrices->SetOwner();
  }

  Int_t *nBins = new Int_t[2*fNVars];
  AliHFEvarInfo *var = NULL;
  for(UInt_t ivar = 0; ivar < fNVars; ivar++){
    var = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(ivar));
    if(var){
      nBins[ivar] = var->GetNumberOfBins();
      nBins[ivar+fNVars] = var->GetNumberOfBins();
    }
  }

  THnSparseF * hTmp = new THnSparseF(name, title, 2*fNVars, nBins);
  for(UInt_t ivar = 0; ivar < fNVars; ivar++){
    var = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(ivar));
    if(var){
      hTmp->SetBinEdges(ivar,var->GetBinning());
      //hTmp->GetAxis(ivar)->Set(var->GetNumberOfBins(), var->GetBinning());
      hTmp->GetAxis(ivar)->SetTitle(var->GetVarName()->Data());
      //hTmp->GetAxis(ivar + fNVars)->Set(var->GetNumberOfBins(), var->GetBinning());
      hTmp->GetAxis(ivar + fNVars)->SetTitle(Form("%s_{MC}", var->GetVarName()->Data()));
      hTmp->SetBinEdges(ivar+fNVars,var->GetBinning());
    }
  }
  hTmp->Sumw2();
  fCorrelationMatrices->AddLast(hTmp);
}

//__________________________________________________________________
AliCFContainer *AliHFEcontainer::GetCFContainer(const Char_t *name) const{
  //
  // Find a given container 
  //
  return dynamic_cast<AliCFContainer *>(fContainers->FindObject(name));
}

//__________________________________________________________________
THnSparseF *AliHFEcontainer::GetCorrelationMatrix(const Char_t *name) const{
  //
  // Find Correlation Matrix
  //
  if(fCorrelationMatrices) return dynamic_cast<THnSparseF *>(fCorrelationMatrices->FindObject(name));
  else return 0x0;

}

//__________________________________________________________________
void AliHFEcontainer::FillCFContainer(const Char_t *name, UInt_t step, const Double_t * const content, Double_t weight) const {
  //
  // Fill container
  //
  AliCFContainer *cont = GetCFContainer(name);
  if(!cont) return;
  cont->Fill(content, step, weight);
}

//__________________________________________________________________
void AliHFEcontainer::FillCFContainerStepname(const Char_t *name, const Char_t *steptitle, const Double_t * const content, Double_t weight)const{
  //
  // Fill container
  //
  AliCFContainer *cont = GetCFContainer(name);
  if(!cont) return;
  // find the matching step title
  Int_t mystep = -1;
  for(Int_t istep = 0; istep < cont->GetNStep(); istep++){
    TString tstept = cont->GetStepTitle(istep);
    if(!tstept.CompareTo(steptitle)){
      mystep = istep;
      break;
    }
  }
  if(mystep < 0){
    // step not found
    AliDebug(1, Form("Step %s not found in container %s", steptitle, name));
    return;
  }
  AliDebug(1, Form("Filling step %s(%d) for container %s", steptitle, mystep, name));
  cont->Fill(content, mystep, weight);
}

//__________________________________________________________________
AliCFContainer *AliHFEcontainer::MakeMergedCFContainer(const Char_t *name, const Char_t *title, const Char_t* contnames) const {
  //
  // Merge CF Container out of several containers 
  // Container names are separated by :
  // returns a new object which has to be taken care of by the user
  //

  TObjArray *containers = TString(contnames).Tokenize(":");
  // we first need the size of the container to be merged
  Int_t nStepMerged = 0;
  AliCFContainer *ctemp = NULL;
  TObjString *cname = NULL;
  for(Int_t icont = 0; icont < containers->GetEntries(); icont++){
    cname = dynamic_cast<TObjString *>(containers->At(icont));
    ctemp = dynamic_cast<AliCFContainer *>(fContainers->FindObject(cname->String().Data()));
    if(!ctemp){
      AliWarning(Form("Container %s not found. It will be unprocessed", cname->String().Data()));
      continue;
    }
    nStepMerged += ctemp->GetNStep(); 
  }
  AliInfo("Please Ignore the messgae comming from AliCFContainer!");
  Int_t *dummyBinning = new Int_t[fNVars];
  for(UInt_t ibin = 0; ibin < fNVars; ibin++) dummyBinning[ibin] = 1;
  AliCFContainer *cmerged = new AliCFContainer(name, title, nStepMerged, fNVars, dummyBinning);
  delete[] dummyBinning;
  // Fill container with content
  AliInfo("Filling new container");
  Int_t cstep = 0;
  for(Int_t icont = 0; icont < containers->GetEntries(); icont++){
    cname = dynamic_cast<TObjString *>(containers->At(icont));
    ctemp = dynamic_cast<AliCFContainer *>(fContainers->FindObject(cname->String().Data()));
    if(!ctemp) continue;
    for(Int_t istep = 0; istep < ctemp->GetNStep(); istep++)
      cmerged->SetGrid(cstep++, new AliCFGridSparse(*ctemp->GetGrid(istep)));
  }
  delete containers;
  return cmerged;
}

//__________________________________________________________________
void AliHFEcontainer::SetStepTitle(const Char_t *contname, const Char_t *steptitle, UInt_t step){
  //
  // Set title for given analysis step in container with name contname
  //
  AliCFContainer *cont = GetCFContainer(contname);
  if(!cont) return;
  if(step >= static_cast<UInt_t>(cont->GetNStep())) return;
  cont->SetStepTitle(step, steptitle);
}

//__________________________________________________________________
void AliHFEcontainer::MakeLinearBinning(UInt_t var, UInt_t nBins, Double_t begin, Double_t end){
  //
  // Set Linear binning for the given container
  //
  AliHFEvarInfo *myvar = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(var));
  if(myvar) myvar->SetBinning(nBins, AliHFEtools::MakeLinearBinning(nBins, begin, end));
}

//__________________________________________________________________
void AliHFEcontainer::MakeLogarithmicBinning(UInt_t var, UInt_t nBins, Double_t begin, Double_t end){
  //
  // Set Logarithmic binning for the given container
  //
  AliHFEvarInfo *myvar = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(var));
  if(myvar) myvar->SetBinning(nBins, AliHFEtools::MakeLogarithmicBinning(nBins, begin, end));
}

//__________________________________________________________________
void AliHFEcontainer::MakeUserDefinedBinning(UInt_t var, UInt_t nBins, const Double_t *binning){
  //
  // Set User defined binning
  //
  AliHFEvarInfo *myvar = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(var));
  if(myvar) myvar->SetBinning(nBins, binning);
}

//__________________________________________________________________
void AliHFEcontainer::SetVariableName(UInt_t var, const Char_t *varname){
  //
  // Variable name
  // 
  AliHFEvarInfo *myvar = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(var));
  if(myvar) myvar->SetVarName(varname);
}

//__________________________________________________________________
Int_t AliHFEcontainer::GetNumberOfCFContainers() const{
  //
  // Get the number of entries
  //
  return fContainers->GetEntries();
}

//__________________________________________________________________
void AliHFEcontainer::Sumw2(const char *contname) const {
  // 
  // Set sums of weights for all steps of a given container
  //
  AliCFContainer *cont = GetCFContainer(contname);
  if(cont){
    for(Int_t istep = 0; istep < cont->GetNStep(); istep++)
      cont->GetGrid(istep)->SumW2();
  }
}

//__________________________________________________________________
void AliHFEcontainer::Print(const Option_t *)const{
  //
  // Print Container Status
  //
  std::cout << "Container status: " << std::endl;
  std::cout << "=====================================================\n";
  std::cout << "Number of variables: " << fNVars << std::endl;
  if(fNVars){
    UInt_t nVars = fVariables ? fVariables->GetEntriesFast() : 0;
    if(nVars != fNVars)
      std::cout << "Inconsistency in number of Variables [" << fNVars << "|" << nVars << "]" << std::endl;
    AliHFEvarInfo *var = NULL;
    if(fVariables){
      for(UInt_t ivar = 0; ivar < fNVars; ivar++){
        var = dynamic_cast<AliHFEvarInfo *>(fVariables->UncheckedAt(ivar));
        if(var)
          std::cout << "Variable " << ivar << ": Name: " << var->GetVarName()->Data() << ", Number of Bins: " << var->GetNumberOfBins() << std::endl;
      }
    }
  }
  std::cout << std::endl;

  // Print CF Containers:
  if(fContainers){
    std::cout << "Containers[" << fContainers->GetEntries() << "]: "<< std::endl;
    std::cout << "=====================================================\n";
    for(Int_t icont = 0; icont < fContainers->GetEntries(); icont++){
      AliCFContainer *c = dynamic_cast<AliCFContainer *>(fContainers->At(icont));
      if(c){
        std::cout << "Name: " << c->GetName() << ", Title: "  << c->GetTitle() << std::endl;
        for(Int_t istep = 0; istep < c->GetNStep(); istep++)
          std::cout << "Step " << istep << ": Title " << c->GetStepTitle(istep) << std::endl;
      }
      std::cout << "------------------------------------------------------\n";
    }
  }
  std::cout << "Number of Events: " << fNEvents << std::endl;
}

//------------------------------------ Content of class AliHFEvarInfo -----------------------------------
//__________________________________________________________________
AliHFEcontainer::AliHFEvarInfo::AliHFEvarInfo():
  TObject(),
  fVarName(NULL),
  fBinning(NULL)
{
  // Default constructor
  fBinning = new TArrayD;
  fVarName = new TString;
}

//__________________________________________________________________
AliHFEcontainer::AliHFEvarInfo::AliHFEvarInfo(const Char_t *name):
  TObject(),
  fVarName(NULL),
  fBinning(NULL)
{
  fBinning = new TArrayD;
  fVarName = new TString(name);
}

//__________________________________________________________________
AliHFEcontainer::AliHFEvarInfo::AliHFEvarInfo(const AliHFEvarInfo &ref):
  TObject(ref),
  fVarName(NULL),
  fBinning(NULL)
{
  //
  // Constructor
  //
  fVarName = new TString(*(ref.fVarName));
  fBinning = new TArrayD(*(ref.fBinning));
}

//__________________________________________________________________
AliHFEcontainer::AliHFEvarInfo &AliHFEcontainer::AliHFEvarInfo::operator=(const AliHFEvarInfo &ref){
  //
  // Assignment operator
  //
  TObject::operator=(ref);
  *fVarName = *(ref.fVarName);
  *fBinning = *(ref.fBinning);
  return *this;
}

//__________________________________________________________________
AliHFEcontainer::AliHFEvarInfo::~AliHFEvarInfo(){
  //
  // Destructor
  //
  delete fVarName;
  delete fBinning;
}

//__________________________________________________________________
void AliHFEcontainer::AliHFEvarInfo::SetVarName(const Char_t *name){
  //
  // Setter for var name
  //
  *fVarName = name;
}

//__________________________________________________________________
void AliHFEcontainer::AliHFEvarInfo::SetBinning(UInt_t nBins, const Double_t *content){
  // Setter for binning
  //
  fBinning->Set(nBins + 1, content);
}

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