ROOT logo
// $Id$
//
// Container with name, TClonesArray and cuts for jets
//
// Author: M. Verweij

#include <TClonesArray.h>

#include "AliEmcalJet.h"
#include "AliVEvent.h"
#include "AliLog.h"
#include "AliEMCALGeometry.h"
#include "AliParticleContainer.h"
#include "AliClusterContainer.h"
#include "AliLocalRhoParameter.h"
#include "AliPythiaInfo.h"

#include "AliJetContainer.h"

ClassImp(AliJetContainer)

//________________________________________________________________________
AliJetContainer::AliJetContainer():
  AliEmcalContainer("AliJetContainer"),
  fJetAcceptanceType(kUser),
  fJetRadius(0),
  fRhoName(),
  fLocalRhoName(),
  fRhoMassName(),
  fPythiaInfoName(),
  fFlavourSelection(0),
  fPtBiasJetTrack(0),
  fPtBiasJetClus(0),
  fJetPtCut(1),
  fJetAreaCut(-1),
  fAreaEmcCut(0),
  fJetMinEta(-0.9),
  fJetMaxEta(0.9),
  fJetMinPhi(-10),
  fJetMaxPhi(10),
  fMaxClusterPt(1000),
  fMaxTrackPt(100),
  fZLeadingEmcCut(10.),
  fZLeadingChCut(10.),
  fNEFMinCut(-10.),
  fNEFMaxCut(10.),
  fLeadingHadronType(0),
  fNLeadingJets(1),
  fJetBitMap(0),
  fJetTrigger(0),
  fTagStatus(-1),
  fParticleContainer(0),
  fClusterContainer(0),
  fRho(0),
  fLocalRho(0),
  fRhoMass(0),
  fPythiaInfo(0),
  fGeom(0),
  fRunNumber(0)
{
  // Default constructor.

  fClassName = "AliEmcalJet";
}

//________________________________________________________________________
AliJetContainer::AliJetContainer(const char *name):
  AliEmcalContainer(name),
  fJetAcceptanceType(kUser),
  fJetRadius(0),
  fRhoName(),
  fLocalRhoName(),
  fRhoMassName(),
  fPythiaInfoName(),
  fFlavourSelection(0),
  fPtBiasJetTrack(0),
  fPtBiasJetClus(0),
  fJetPtCut(1),
  fJetAreaCut(-1),
  fAreaEmcCut(0),
  fJetMinEta(-0.9),
  fJetMaxEta(0.9),
  fJetMinPhi(-10),
  fJetMaxPhi(10),
  fMaxClusterPt(1000),
  fMaxTrackPt(100),
  fZLeadingEmcCut(10.),
  fZLeadingChCut(10.),
  fNEFMinCut(-10.),
  fNEFMaxCut(10.),
  fLeadingHadronType(0),
  fNLeadingJets(1),
  fJetBitMap(0),
  fJetTrigger(0),
  fTagStatus(-1),
  fParticleContainer(0),
  fClusterContainer(0),
  fRho(0),
  fLocalRho(0),
  fRhoMass(0),
  fPythiaInfo(0),
  fGeom(0),
  fRunNumber(0)
{
  // Standard constructor.

  fClassName = "AliEmcalJet";
}

//________________________________________________________________________
void AliJetContainer::SetArray(AliVEvent *event) 
{
  // Set jet array

  AliEmcalContainer::SetArray(event);

  if(fJetAcceptanceType==kTPC) {
    AliDebug(2,Form("%s: set TPC acceptance cuts",GetName()));
    SetJetEtaPhiTPC();
  }
  else if(fJetAcceptanceType==kEMCAL) {
    AliDebug(2,Form("%s: set EMCAL acceptance cuts",GetName()));
    SetJetEtaPhiEMCAL();
 }
}


//________________________________________________________________________
void AliJetContainer::SetEMCALGeometry() {
  fGeom = AliEMCALGeometry::GetInstance();
  if (!fGeom) {
    AliError(Form("%s: Can not create geometry", GetName()));
    return;
  }
}

//________________________________________________________________________
void AliJetContainer::LoadRho(AliVEvent *event)
{
  // Load rho

  if (!fRhoName.IsNull() && !fRho) {
    fRho = dynamic_cast<AliRhoParameter*>(event->FindListObject(fRhoName));
    if (!fRho) {
      AliError(Form("%s: Could not retrieve rho %s!", GetName(), fRhoName.Data()));
      return;
    }
  }
}

//________________________________________________________________________
void AliJetContainer::LoadLocalRho(AliVEvent *event)
{
  // Load local rho

  if (!fLocalRhoName.IsNull() && !fLocalRho) {
    fLocalRho = dynamic_cast<AliLocalRhoParameter*>(event->FindListObject(fLocalRhoName));
    if (!fLocalRho) {
      AliError(Form("%s: Could not retrieve rho %s!", GetName(), fLocalRhoName.Data()));
      return;
    }
  }
}

//________________________________________________________________________
void AliJetContainer::LoadRhoMass(AliVEvent *event)
{
  // Load rho

  if (!fRhoMassName.IsNull() && !fRhoMass) {
    fRhoMass = dynamic_cast<AliRhoParameter*>(event->FindListObject(fRhoMassName));
    if (!fRhoMass) {
      AliError(Form("%s: Could not retrieve rho_mass %s!", GetName(), fRhoMassName.Data()));
      return;
    }
  }
}
//________________________________________________________________________
void AliJetContainer::LoadPythiaInfo(AliVEvent *event)
{
    // Load parton info
    
    if (!fPythiaInfoName.IsNull() && !fPythiaInfo) {
        fPythiaInfo = dynamic_cast<AliPythiaInfo*>(event->FindListObject(fPythiaInfoName));
        if (!fPythiaInfo) {
           AliError(Form("%s: Could not retrieve parton infos! %s!", GetName(), fPythiaInfoName.Data()));    
        return;
        }
    }
}



//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetLeadingJet(const char* opt)
{
  // Get the leading jet; if opt contains "rho" the sorting is according to pt-A*rho

  TString option(opt);
  option.ToLower();

  Int_t tempID = fCurrentID;

  AliEmcalJet *jetMax = GetNextAcceptJet(0);
  AliEmcalJet *jet = 0;

  if (option.Contains("rho")) {
    while ((jet = GetNextAcceptJet())) {
      if ( (jet->Pt()-jet->Area()*GetRhoVal()) > (jetMax->Pt()-jetMax->Area()*GetRhoVal()) )
	jetMax = jet;
    }
  }
  else {
    while ((jet = GetNextAcceptJet())) {
      if (jet->Pt() > jetMax->Pt()) jetMax = jet;
    }
  }

  fCurrentID = tempID;

  return jetMax;
}

//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetJet(Int_t i) const {

  //Get i^th jet in array

  if(i<0 || i>fClArray->GetEntriesFast()) return 0;
  AliEmcalJet *jet = static_cast<AliEmcalJet*>(fClArray->At(i));
  return jet;

}

//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetAcceptJet(Int_t i) {

  //Only return jet if is accepted

  AliEmcalJet *jet = GetJet(i);
  if(!AcceptJet(jet)) return 0;

  return jet;
}

//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetJetWithLabel(Int_t lab) const {

  //Get particle with label lab in array
  
  Int_t i = GetIndexFromLabel(lab);
  return GetJet(i);
}

//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetAcceptJetWithLabel(Int_t lab) {

  //Get particle with label lab in array
  
  Int_t i = GetIndexFromLabel(lab);
  return GetAcceptJet(i);
}

//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetNextAcceptJet(Int_t i) {

  //Get next accepted jet; if i >= 0 (re)start counter from i; return 0 if no accepted jet could be found

  if (i>=0) fCurrentID = i;

  const Int_t njets = GetNEntries();
  AliEmcalJet *jet = 0;
  while (fCurrentID < njets && !jet) { 
    jet = GetAcceptJet(fCurrentID);
    fCurrentID++;
  }

  return jet;
}

//________________________________________________________________________
AliEmcalJet* AliJetContainer::GetNextJet(Int_t i) {

  //Get next jet; if i >= 0 (re)start counter from i; return 0 if no jet could be found

  if (i>=0) fCurrentID = i;

  const Int_t njets = GetNEntries();
  AliEmcalJet *jet = 0;
  while (fCurrentID < njets && !jet) { 
    jet = GetJet(fCurrentID);
    fCurrentID++;
  }

  return jet;
}

//________________________________________________________________________
Double_t AliJetContainer::GetJetPtCorr(Int_t i) const {
  AliEmcalJet *jet = GetJet(i);

  return jet->Pt() - fRho->GetVal()*jet->Area();
}

//________________________________________________________________________
Double_t AliJetContainer::GetJetPtCorrLocal(Int_t i) const {
  AliEmcalJet *jet = GetJet(i);

  return jet->Pt() - fLocalRho->GetLocalVal(jet->Phi(), fJetRadius)*jet->Area();
}

//________________________________________________________________________
void AliJetContainer::GetMomentum(TLorentzVector &mom, Int_t i) const
{
  //Get momentum of the i^th jet in array

  AliEmcalJet *jet = GetJet(i);
  if(jet) jet->GetMom(mom);
}

//________________________________________________________________________
Bool_t AliJetContainer::AcceptBiasJet(const AliEmcalJet *jet)
{ 
  // Accept jet with a bias.

  if (fLeadingHadronType == 0) {
    if (jet->MaxTrackPt() < fPtBiasJetTrack) return kFALSE;
  }
  else if (fLeadingHadronType == 1) {
    if (jet->MaxClusterPt() < fPtBiasJetClus) return kFALSE;
  }
  else {
    if (jet->MaxTrackPt() < fPtBiasJetTrack && jet->MaxClusterPt() < fPtBiasJetClus) return kFALSE;
  }

  return kTRUE;
}

//________________________________________________________________________
Bool_t AliJetContainer::AcceptJet(const AliEmcalJet *jet)
{   
   // Return true if jet is accepted.

  fRejectionReason = 0;

  if (!jet) {
    AliDebug(11,"No jet found");
    fRejectionReason |= kNullObject;
    return kFALSE;
  }

  if (jet->Pt() <= fJetPtCut) {
    AliDebug(11,Form("Cut rejecting jet: JetPtCut %.1f",fJetPtCut));
    fRejectionReason |= kPtCut;
    return kFALSE;
  }

  Double_t jetPhi = jet->Phi();
  Double_t jetEta = jet->Eta();
   
  // if limits are given in (-pi, pi) range
  if (fJetMinPhi < 0) jetPhi -= TMath::Pi() * 2;
   
  if (jetEta < fJetMinEta || jetEta > fJetMaxEta || jetPhi < fJetMinPhi || jetPhi > fJetMaxPhi) {
    AliDebug(11,"Cut rejecting jet: Acceptance");
    fRejectionReason |= kAcceptanceCut;
    return kFALSE;
  }

  if (jet->TestBits(fJetBitMap) != (Int_t)fJetBitMap) {
    AliDebug(11,"Cut rejecting jet: Bit map");
    fRejectionReason |= kBitMapCut;
    return kFALSE;
  }

  if (jet->Area() <= fJetAreaCut)  {
    AliDebug(11,"Cut rejecting jet: Area");
    fRejectionReason |= kAreaCut;
    return kFALSE;
  }

  if (jet->AreaEmc() < fAreaEmcCut) {
    AliDebug(11,"Cut rejecting jet: AreaEmc");
    fRejectionReason |= kAreaEmcCut;
    return kFALSE;
  }
   
  if (fZLeadingChCut < 1 && GetZLeadingCharged(jet) > fZLeadingChCut) {
    AliDebug(11,"Cut rejecting jet: ZLeading");
    fRejectionReason |= kZLeadingChCut;
    return kFALSE;
  }
   
  if (fZLeadingEmcCut < 1 && GetZLeadingEmc(jet) > fZLeadingEmcCut) {
    AliDebug(11,"Cut rejecting jet: ZLeadEmc");
    fRejectionReason |= kZLeadingEmcCut;
    return kFALSE;
  }

  if (jet->NEF() < fNEFMinCut || jet->NEF() > fNEFMaxCut) {
    AliDebug(11,"Cut rejecting jet: NEF");
    fRejectionReason |= kNEFCut;
    return kFALSE;
  }
   
  if (!AcceptBiasJet(jet)) {
    AliDebug(11,"Cut rejecting jet: Bias");
    fRejectionReason |= kMinLeadPtCut;
    return kFALSE;
  }

  if (jet->MaxTrackPt() > fMaxTrackPt) {
    AliDebug(11,"Cut rejecting jet: MaxTrackPt");
    fRejectionReason |= kMaxTrackPtCut;
    return kFALSE;

  }

  if (jet->MaxClusterPt() > fMaxClusterPt) {
    AliDebug(11,"Cut rejecting jet: MaxClusPt");
    fRejectionReason |= kMaxClusterPtCut;
    return kFALSE;
  }
   
  if (fFlavourSelection != 0 && !jet->TestFlavourTag(fFlavourSelection)) {
    AliDebug(11,"Cut rejecting jet: Flavour");
    fRejectionReason |= kFlavourCut;
    return kFALSE;
  }
   
  if(fTagStatus>-1 && jet->GetTagStatus()!=fTagStatus) {
    AliDebug(11,"Cut rejecting jet: tag status");
    fRejectionReason |= kTagStatus;
    return kFALSE;
  }

  return kTRUE;
}

//________________________________________________________________________
Double_t AliJetContainer::GetLeadingHadronPt(const AliEmcalJet *jet) const
{
  if (fLeadingHadronType == 0)       // charged leading hadron
    return jet->MaxTrackPt();
  else if (fLeadingHadronType == 1)  // neutral leading hadron
    return jet->MaxClusterPt();
  else                               // charged or neutral
    return jet->MaxPartPt();
}

//________________________________________________________________________
void AliJetContainer::GetLeadingHadronMomentum(TLorentzVector &mom, const AliEmcalJet *jet) const
{
  Double_t maxClusterPt = 0;
  Double_t maxClusterEta = 0;
  Double_t maxClusterPhi = 0;
  
  Double_t maxTrackPt = 0;
  Double_t maxTrackEta = 0;
  Double_t maxTrackPhi = 0;
      
  if (fClusterContainer && fClusterContainer->GetArray() && (fLeadingHadronType == 1 || fLeadingHadronType == 2)) {
    AliVCluster *cluster = jet->GetLeadingCluster(fClusterContainer->GetArray());
    if (cluster) {
      TLorentzVector nPart;
      cluster->GetMomentum(nPart, const_cast<Double_t*>(fVertex));
      
      maxClusterEta = nPart.Eta();
      maxClusterPhi = nPart.Phi();
      maxClusterPt = nPart.Pt();
    }
  }
      
  if (fParticleContainer && fParticleContainer->GetArray() && (fLeadingHadronType == 0 || fLeadingHadronType == 2)) {
    AliVParticle *track = jet->GetLeadingTrack(fParticleContainer->GetArray());
    if (track) {
      maxTrackEta = track->Eta();
      maxTrackPhi = track->Phi();
      maxTrackPt = track->Pt();
    }
  }
      
  if (maxTrackPt > maxClusterPt) 
    mom.SetPtEtaPhiM(maxTrackPt,maxTrackEta,maxTrackPhi,0.139);
  else 
    mom.SetPtEtaPhiM(maxClusterPt,maxClusterEta,maxClusterPhi,0.139);
}

//________________________________________________________________________
Double_t AliJetContainer::GetZLeadingEmc(const AliEmcalJet *jet) const
{

  if (fClusterContainer && fClusterContainer->GetArray()) {
    TLorentzVector mom;
    
    AliVCluster *cluster = jet->GetLeadingCluster(fClusterContainer->GetArray());
    if (cluster) {
      cluster->GetMomentum(mom, const_cast<Double_t*>(fVertex));
      
      return GetZ(jet,mom);
    }
    else
      return -1;
  }
  else
    return -1;
}

//________________________________________________________________________
Double_t AliJetContainer::GetZLeadingCharged(const AliEmcalJet *jet) const
{

  if (fParticleContainer && fParticleContainer->GetArray() ) {
    TLorentzVector mom;
    
    AliVParticle *track = jet->GetLeadingTrack(fParticleContainer->GetArray());
    if (track) {
      mom.SetPtEtaPhiM(track->Pt(),track->Eta(),track->Phi(),0.139);
      
      return GetZ(jet,mom);
    }
    else
      return -1;
  }
  else
    return -1;
}

//________________________________________________________________________
Double_t AliJetContainer::GetZ(const AliEmcalJet *jet, TLorentzVector mom) const
{

  Double_t pJetSq = jet->Px()*jet->Px() + jet->Py()*jet->Py() + jet->Pz()*jet->Pz();

  if(pJetSq>1e-6)
    return (mom.Px()*jet->Px() + mom.Py()*jet->Py() + mom.Pz()*jet->Pz())/pJetSq;
  else {
    AliWarning(Form("%s: strange, pjet*pjet seems to be zero pJetSq: %f",GetName(), pJetSq));
    return -1;
  }

}

//________________________________________________________________________
void AliJetContainer::SetJetEtaPhiEMCAL()
{
  //Set default cuts for full jets

  if(!fGeom) SetEMCALGeometry();
  if(fGeom) {
    SetJetEtaLimits(fGeom->GetArm1EtaMin() + fJetRadius, fGeom->GetArm1EtaMax() - fJetRadius);

    if(fRunNumber>=177295 && fRunNumber<=197470) //small SM masked in 2012 and 2013
      SetJetPhiLimits(1.405+fJetRadius,3.135-fJetRadius);
    else
      SetJetPhiLimits(fGeom->GetArm1PhiMin() * TMath::DegToRad() + fJetRadius, fGeom->GetArm1PhiMax() * TMath::DegToRad() - fJetRadius);
  }
  else {
    AliWarning("Could not get instance of AliEMCALGeometry. Using manual settings for EMCAL year 2011!!");
    SetJetEtaLimits(-0.7+fJetRadius,0.7-fJetRadius);
    SetJetPhiLimits(1.405+fJetRadius,3.135-fJetRadius);
  }
}

//________________________________________________________________________
void AliJetContainer::SetJetEtaPhiTPC()
{
  //Set default cuts for charged jets

  SetJetEtaLimits(-0.9+fJetRadius, 0.9-fJetRadius);
  SetJetPhiLimits(-10, 10);
}

//________________________________________________________________________
void AliJetContainer::PrintCuts() 
{
  // Reset cuts to default values

  Printf("PtBiasJetTrack: %f",fPtBiasJetTrack);
  Printf("PtBiasJetClus: %f",fPtBiasJetClus);
  Printf("JetPtCut: %f", fJetPtCut);
  Printf("JetAreaCut: %f",fJetAreaCut);
  Printf("AreaEmcCut: %f",fAreaEmcCut);
  Printf("JetMinEta: %f", fJetMinEta);
  Printf("JetMaxEta: %f", fJetMaxEta);
  Printf("JetMinPhi: %f", fJetMinPhi);
  Printf("JetMaxPhi: %f", fJetMaxPhi);
  Printf("MaxClusterPt: %f",fMaxClusterPt);
  Printf("MaxTrackPt: %f",fMaxTrackPt);
  Printf("LeadingHadronType: %d",fLeadingHadronType);
  Printf("ZLeadingEmcCut: %f",fZLeadingEmcCut);
  Printf("ZLeadingChCut: %f",fZLeadingChCut);

}

//________________________________________________________________________
void AliJetContainer::ResetCuts() 
{
  // Reset cuts to default values

  fPtBiasJetTrack = 0;
  fPtBiasJetClus = 0;
  fJetPtCut = 1;
  fJetAreaCut = -1;
  fAreaEmcCut = 0;
  fJetMinEta = -0.9;
  fJetMaxEta = 0.9;
  fJetMinPhi = -10;
  fJetMaxPhi = 10;
  fMaxClusterPt = 1000;
  fMaxTrackPt = 100;
  fLeadingHadronType = 0;
  fZLeadingEmcCut = 10.;
  fZLeadingChCut = 10.;
}

//________________________________________________________________________
void AliJetContainer::SetClassName(const char *clname)
{
  // Set the class name

  TClass cls(clname);
  if (cls.InheritsFrom("AliEmcalJet")) fClassName = clname;
  else AliError(Form("Unable to set class name %s for a AliJetContainer, it must inherits from AliEmcalJet!",clname));
}

//________________________________________________________________________
Double_t AliJetContainer::GetFractionSharedPt(const AliEmcalJet *jet1) const
{
  //
  // Get fraction of shared pT between matched jets
  // Uses ClosestJet() jet pT as baseline: fraction = \Sum_{const,jet1} pT,const,i / pT,jet,closest
  // Only works if tracks array of both jets is the same
  //

  AliEmcalJet *jet2 = jet1->ClosestJet();
  if(!jet2) return -1;

  Double_t fraction = 0.;
  Double_t jetPt2 = jet2->Pt();
 
  if(jetPt2>0) {
    Double_t sumPt = 0.;
    AliVParticle *vpf = 0x0;
    Int_t iFound = 0;
    for(Int_t icc=0; icc<jet2->GetNumberOfTracks(); icc++) {
      Int_t idx = (Int_t)jet2->TrackAt(icc);
      iFound = 0;
      for(Int_t icf=0; icf<jet1->GetNumberOfTracks(); icf++) {
	if(idx == jet1->TrackAt(icf) && iFound==0 ) {
	  iFound=1;
	  vpf = static_cast<AliVParticle*>(jet1->TrackAt(icf, fParticleContainer->GetArray()));
	  if(vpf) sumPt += vpf->Pt();
	  continue;
	}
      }
    }
    fraction = sumPt/jetPt2;
  } else 
    fraction = -1;
  
  return fraction;
}

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