ROOT logo
#ifndef AliUEHistograms_H
#define AliUEHistograms_H

/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice                               */

/* $Id: AliUEHistograms.h 20164 2007-08-14 15:31:50Z morsch $ */

// encapsulates several AliUEHist objects for a full UE analysis plus additional control histograms

#include "TNamed.h"
#include "AliUEHist.h"
#include "TMath.h"
#include "THn.h" // in cxx file causes .../THn.h:257: error: conflicting declaration ‘typedef class THnT<float> THnF’

class AliVParticle;

class TList;
class TSeqCollection;
class TObjArray;
class TH1F;
class TH2F;
class TH3F;

class AliUEHistograms : public TNamed
{
 public:
  AliUEHistograms(const char* name = "AliUEHistograms", const char* histograms = "", const char* binning = 0);
  virtual ~AliUEHistograms();
  
  void Fill(Int_t eventType, Float_t zVtx, AliUEHist::CFStep step, AliVParticle* leading, TList* toward, TList* away, TList* min, TList* max);
  void FillCorrelations(Double_t centrality, Float_t zVtx, AliUEHist::CFStep step, TObjArray* particles, TObjArray* mixed = 0, Float_t weight = 1, Bool_t firstTime = kTRUE, Bool_t twoTrackEfficiencyCut = kFALSE, Float_t bSign = 0, Float_t twoTrackEfficiencyCutValue = 0.02, Bool_t applyEfficiency = kFALSE);
  void Fill(AliVParticle* leadingMC, AliVParticle* leadingReco);
  void FillEvent(Int_t eventType, Int_t step);
  void FillEvent(Double_t centrality, Int_t step);
  void FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, TObjArray* recoPrimPID, TObjArray* recoAllPID, TObjArray* fake, Int_t particleType, Double_t centrality = 0, Double_t zVtx = 0);
  void FillFakePt(TObjArray* fake, Double_t centrality);
 
  void CopyReconstructedData(AliUEHistograms* from);
  void DeepCopy(AliUEHistograms* from);
  
  AliUEHist* GetUEHist(Int_t id);
  
  AliUEHist* GetNumberDensitypT() { return fNumberDensitypT; }
  AliUEHist* GetSumpT() { return fSumpT; }
  AliUEHist* GetNumberDensityPhi() { return fNumberDensityPhi; }
  
  void SetNumberDensitypT(AliUEHist* obj) { fNumberDensitypT = obj; }
  void SetSumpT(AliUEHist* obj) { fSumpT = obj; }
  void SetNumberDensityPhi(AliUEHist* obj) { fNumberDensityPhi = obj; }
  
  void SetRunNumber(Long64_t runNumber) { fRunNumber = runNumber; }
  
  void SetEfficiencyCorrectionTriggers(THnF* hist)   { fEfficiencyCorrectionTriggers = hist;   }
  void SetEfficiencyCorrectionAssociated(THnF* hist) { fEfficiencyCorrectionAssociated = hist; }
  
  TH2F* GetCorrelationpT()  { return fCorrelationpT; }
  TH2F* GetCorrelationEta() { return fCorrelationEta; }
  TH2F* GetCorrelationPhi() { return fCorrelationPhi; }
  TH2F* GetCorrelationR()   { return fCorrelationR; }
  TH2F* GetCorrelationLeading2Phi() { return fCorrelationLeading2Phi; }
  TH2F* GetCorrelationMultiplicity() { return fCorrelationMultiplicity; }
  TH3F* GetYield() { return fYields; }
  TH2F* GetInvYield() { return fInvYield2; }
  
  TH2F* GetEventCount()     { return fEventCount; }
  TH3F* GetEventCountDifferential() { return fEventCountDifferential; }
  TH1F* GetVertexContributors() { return fVertexContributors; }
  TH1F* GetCentralityDistribution() { return fCentralityDistribution; }
  TH2F* GetCentralityCorrelation() { return fCentralityCorrelation; }
  Long64_t GetRunNumber() { return fRunNumber; }
  Int_t GetMergeCount() { return fMergeCount; }
  TH3F* GetTwoTrackDistance(Int_t i) { return fTwoTrackDistancePt[i]; }
  Bool_t GetWeightPerEvent() { return fWeightPerEvent; }
  
  void Correct(AliUEHistograms* corrections);
  
  void SetEtaRange(Float_t etaMin, Float_t etaMax);
  void SetPtRange(Float_t ptMin, Float_t ptMax);
  void SetPartSpecies(Int_t species);
  void SetZVtxRange(Float_t min, Float_t max);
  void SetContaminationEnhancement(TH1F* hist);
  void SetCombineMinMax(Bool_t flag);
  void SetTrackEtaCut(Float_t value);
  void SetWeightPerEvent(Bool_t flag);
  void SetSelectCharge(Int_t selectCharge) { fSelectCharge = selectCharge; }
  void SetSelectTriggerCharge(Int_t selectCharge) { fTriggerSelectCharge = selectCharge; }
  void SetSelectAssociatedCharge(Int_t selectCharge) { fAssociatedSelectCharge = selectCharge; }
  void SetTriggerRestrictEta(Float_t eta) { fTriggerRestrictEta = eta; }
  void SetEtaOrdering(Bool_t flag) { fEtaOrdering = flag; }
  void SetPairCuts(Bool_t conversions, Bool_t resonances) { fCutConversions = conversions; fCutResonances = resonances; }
  void SetRejectResonanceDaughters(Int_t value) { fRejectResonanceDaughters = value; }
  void SetOnlyOneEtaSide(Int_t flag)    { fOnlyOneEtaSide = flag; }
  void SetPtOrder(Bool_t flag) { fPtOrder = flag; }
  void SetTwoTrackCutMinRadius(Float_t min) { fTwoTrackCutMinRadius = min; }
  
  void ExtendTrackingEfficiency(Bool_t verbose = kFALSE);
  void Reset();

  AliUEHistograms(const AliUEHistograms &c);
  AliUEHistograms& operator=(const AliUEHistograms& c);
  virtual void Copy(TObject& c) const;

  virtual Long64_t Merge(TCollection* list);
  void Scale(Double_t factor);
  
protected:
  void FillRegion(AliUEHist::Region region, Float_t zVtx, AliUEHist::CFStep step, AliVParticle* leading, TList* list, Int_t multiplicity);
  Int_t CountParticles(TList* list, Float_t ptMin);
  void DeleteContainers();
  inline Float_t GetInvMassSquared(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2);
  inline Float_t GetInvMassSquaredCheap(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2);
  inline Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1, Float_t phi2, Float_t pt2, Float_t charge2, Float_t radius, Float_t bSign);
  
  static const Int_t fgkUEHists; // number of histograms

  AliUEHist* fNumberDensitypT;   // d^2N/dphideta vs pT,lead
  AliUEHist* fSumpT;             // d^2 sum(pT)/dphideta vs pT,lead
  AliUEHist* fNumberDensityPhi;  // d^2N/dphideta vs delta phi,lead (in pT,lead bins)
  
  TH2F* fCorrelationpT;         // pT,lead: true vs reco
  TH2F* fCorrelationEta;        // #eta,lead; true vs reco
  TH2F* fCorrelationPhi;        // #phi,lead; true vs reco
  TH2F* fCorrelationR;          // R = sqrt(delta eta^2 + delta phi^2) (true vs reco) vs pT,lead,MC
  TH2F* fCorrelationLeading2Phi;// delta phi (true vs reco) vs pT,lead,MC
  TH2F* fCorrelationMultiplicity; // number of mc particls vs reco particles (for pT > 0.5 GeV/c)
  TH3F* fYields;                // centrality vs pT vs eta
  TH2F* fInvYield2; 		// invariant yield as cross check of tracking
  
  TH2F* fEventCount;            // event count as function of step, (for pp: event type (plus additional step -1 for all events without vertex range even in MC)) (for PbPb: centrality)
  TH3F* fEventCountDifferential;// event count as function of leading pT, step, event type
  
  TH1F* fVertexContributors;    // number of contributors to the vertex
  TH1F* fCentralityDistribution; // distribution of the variable used for centrality selection
  TH2F* fCentralityCorrelation;  // centrality vs multiplicity
  
  TH3F* fITSClusterMap;          // its cluster map vs centrality vs pT
  
  TH3F* fTwoTrackDistancePt[2];    // control histograms for two-track efficiency study: dphi*_min vs deta (0 = before cut, 1 = after cut)
  TH2F* fControlConvResoncances; // control histograms for cuts on conversions and resonances
  
  THnF* fEfficiencyCorrectionTriggers;   // if non-0 this efficiency correction is applied on the fly to the filling for trigger particles. The factor is multiplicative, i.e. should contain 1/efficiency
  THnF* fEfficiencyCorrectionAssociated;   // if non-0 this efficiency correction is applied on the fly to the filling for associated particles. The factor is multiplicative, i.e. should contain 1/efficiency
  
  Int_t fSelectCharge;           // (un)like sign selection when building correlations: 0: no selection; 1: unlike sign; 2: like sign
  Int_t fTriggerSelectCharge;    // select charge of trigger particle
  Int_t fAssociatedSelectCharge; // select charge of associated particle
  Float_t fTriggerRestrictEta;   // restrict eta range for trigger particle (default: -1 [off])
  Bool_t fEtaOrdering;           // activate eta ordering to prevent shape distortions. see FillCorrelation for the details
  Bool_t fCutConversions;        // cut on conversions (inv mass)
  Bool_t fCutResonances;         // cut on resonances (inv mass)
  Int_t fRejectResonanceDaughters; // reject all daughters of all resonance candidates (1: test method (cut at m_inv=0.9); 2: k0; 3: lambda)
  Int_t fOnlyOneEtaSide;       // decides that only trigger particle from one eta side are considered (0 = all; -1 = negative, 1 = positive)
  Bool_t fWeightPerEvent;	// weight with the number of trigger particles per event
  Bool_t fPtOrder;		// apply pT,a < pT,t condition
  Float_t fTwoTrackCutMinRadius; // min radius for TTR cut
  
  Long64_t fRunNumber;           // run number that has been processed
  
  Int_t fMergeCount;		// counts how many objects have been merged together
  
  ClassDef(AliUEHistograms, 29)  // underlying event histogram container
};

Float_t AliUEHistograms::GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1, Float_t phi2, Float_t pt2, Float_t charge2, Float_t radius, Float_t bSign)
{ 
  //
  // calculates dphistar
  //
  
  Float_t dphistar = phi1 - phi2 - charge1 * bSign * TMath::ASin(0.075 * radius / pt1) + charge2 * bSign * TMath::ASin(0.075 * radius / pt2);
  
  static const Double_t kPi = TMath::Pi();
  
  // circularity
//   if (dphistar > 2 * kPi)
//     dphistar -= 2 * kPi;
//   if (dphistar < -2 * kPi)
//     dphistar += 2 * kPi;
  
  if (dphistar > kPi)
    dphistar = kPi * 2 - dphistar;
  if (dphistar < -kPi)
    dphistar = -kPi * 2 - dphistar;
  if (dphistar > kPi) // might look funny but is needed
    dphistar = kPi * 2 - dphistar;
  
  return dphistar;
}

Float_t AliUEHistograms::GetInvMassSquared(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2)
{
  // calculate inv mass squared
  // same can be achieved, but with more computing time with
  /*TLorentzVector photon, p1, p2;
  p1.SetPtEtaPhiM(triggerParticle->Pt(), triggerEta, triggerParticle->Phi(), 0.510e-3);
  p2.SetPtEtaPhiM(particle->Pt(), eta[j], particle->Phi(), 0.510e-3);
  photon = p1+p2;
  photon.M()*/
  
  Float_t tantheta1 = 1e10;
  
  if (eta1 < -1e-10 || eta1 > 1e-10)
  {
    Float_t expTmp = TMath::Exp(-eta1);
    tantheta1 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
  }
  
  Float_t tantheta2 = 1e10;
  if (eta2 < -1e-10 || eta2 > 1e-10)
  {
    Float_t expTmp = TMath::Exp(-eta2);
    tantheta2 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
  }
  
  Float_t e1squ = m0_1 * m0_1 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
  Float_t e2squ = m0_2 * m0_2 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
  
  Float_t mass2 = m0_1 * m0_1 + m0_2 * m0_2 + 2 * ( TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( TMath::Cos(phi1 - phi2) + 1.0 / tantheta1 / tantheta2 ) ) );
  
//   Printf(Form("%f %f %f %f %f %f %f %f %f", pt1, eta1, phi1, pt2, eta2, phi2, m0_1, m0_2, mass2));
  
  return mass2;
}

Float_t AliUEHistograms::GetInvMassSquaredCheap(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2)
{
  // calculate inv mass squared approximately
  
  Float_t tantheta1 = 1e10;
  
  if (eta1 < -1e-10 || eta1 > 1e-10)
  {
    Float_t expTmp = 1.0-eta1+eta1*eta1/2-eta1*eta1*eta1/6+eta1*eta1*eta1*eta1/24;
    tantheta1 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
  }
  
  Float_t tantheta2 = 1e10;
  if (eta2 < -1e-10 || eta2 > 1e-10)
  {
    Float_t expTmp = 1.0-eta2+eta2*eta2/2-eta2*eta2*eta2/6+eta2*eta2*eta2*eta2/24;
    tantheta2 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
  }
  
  Float_t e1squ = m0_1 * m0_1 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
  Float_t e2squ = m0_2 * m0_2 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
  
  // fold onto 0...pi
  Float_t deltaPhi = TMath::Abs(phi1 - phi2);
  while (deltaPhi > TMath::TwoPi())
    deltaPhi -= TMath::TwoPi();
  if (deltaPhi > TMath::Pi())
    deltaPhi = TMath::TwoPi() - deltaPhi;
  
  Float_t cosDeltaPhi = 0;
  if (deltaPhi < TMath::Pi()/3)
    cosDeltaPhi = 1.0 - deltaPhi*deltaPhi/2 + deltaPhi*deltaPhi*deltaPhi*deltaPhi/24;
  else if (deltaPhi < 2*TMath::Pi()/3)
    cosDeltaPhi = -(deltaPhi - TMath::Pi()/2) + 1.0/6 * TMath::Power((deltaPhi - TMath::Pi()/2), 3);
  else
    cosDeltaPhi = -1.0 + 1.0/2.0*(deltaPhi - TMath::Pi())*(deltaPhi - TMath::Pi()) - 1.0/24.0 * TMath::Power(deltaPhi - TMath::Pi(), 4);
  
  Float_t mass2 = m0_1 * m0_1 + m0_2 * m0_2 + 2 * ( TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( cosDeltaPhi + 1.0 / tantheta1 / tantheta2 ) ) );
  
//   Printf(Form("%f %f %f %f %f %f %f %f %f", pt1, eta1, phi1, pt2, eta2, phi2, m0_1, m0_2, mass2));
  
  return mass2;
}

#endif
 AliUEHistograms.h:1
 AliUEHistograms.h:2
 AliUEHistograms.h:3
 AliUEHistograms.h:4
 AliUEHistograms.h:5
 AliUEHistograms.h:6
 AliUEHistograms.h:7
 AliUEHistograms.h:8
 AliUEHistograms.h:9
 AliUEHistograms.h:10
 AliUEHistograms.h:11
 AliUEHistograms.h:12
 AliUEHistograms.h:13
 AliUEHistograms.h:14
 AliUEHistograms.h:15
 AliUEHistograms.h:16
 AliUEHistograms.h:17
 AliUEHistograms.h:18
 AliUEHistograms.h:19
 AliUEHistograms.h:20
 AliUEHistograms.h:21
 AliUEHistograms.h:22
 AliUEHistograms.h:23
 AliUEHistograms.h:24
 AliUEHistograms.h:25
 AliUEHistograms.h:26
 AliUEHistograms.h:27
 AliUEHistograms.h:28
 AliUEHistograms.h:29
 AliUEHistograms.h:30
 AliUEHistograms.h:31
 AliUEHistograms.h:32
 AliUEHistograms.h:33
 AliUEHistograms.h:34
 AliUEHistograms.h:35
 AliUEHistograms.h:36
 AliUEHistograms.h:37
 AliUEHistograms.h:38
 AliUEHistograms.h:39
 AliUEHistograms.h:40
 AliUEHistograms.h:41
 AliUEHistograms.h:42
 AliUEHistograms.h:43
 AliUEHistograms.h:44
 AliUEHistograms.h:45
 AliUEHistograms.h:46
 AliUEHistograms.h:47
 AliUEHistograms.h:48
 AliUEHistograms.h:49
 AliUEHistograms.h:50
 AliUEHistograms.h:51
 AliUEHistograms.h:52
 AliUEHistograms.h:53
 AliUEHistograms.h:54
 AliUEHistograms.h:55
 AliUEHistograms.h:56
 AliUEHistograms.h:57
 AliUEHistograms.h:58
 AliUEHistograms.h:59
 AliUEHistograms.h:60
 AliUEHistograms.h:61
 AliUEHistograms.h:62
 AliUEHistograms.h:63
 AliUEHistograms.h:64
 AliUEHistograms.h:65
 AliUEHistograms.h:66
 AliUEHistograms.h:67
 AliUEHistograms.h:68
 AliUEHistograms.h:69
 AliUEHistograms.h:70
 AliUEHistograms.h:71
 AliUEHistograms.h:72
 AliUEHistograms.h:73
 AliUEHistograms.h:74
 AliUEHistograms.h:75
 AliUEHistograms.h:76
 AliUEHistograms.h:77
 AliUEHistograms.h:78
 AliUEHistograms.h:79
 AliUEHistograms.h:80
 AliUEHistograms.h:81
 AliUEHistograms.h:82
 AliUEHistograms.h:83
 AliUEHistograms.h:84
 AliUEHistograms.h:85
 AliUEHistograms.h:86
 AliUEHistograms.h:87
 AliUEHistograms.h:88
 AliUEHistograms.h:89
 AliUEHistograms.h:90
 AliUEHistograms.h:91
 AliUEHistograms.h:92
 AliUEHistograms.h:93
 AliUEHistograms.h:94
 AliUEHistograms.h:95
 AliUEHistograms.h:96
 AliUEHistograms.h:97
 AliUEHistograms.h:98
 AliUEHistograms.h:99
 AliUEHistograms.h:100
 AliUEHistograms.h:101
 AliUEHistograms.h:102
 AliUEHistograms.h:103
 AliUEHistograms.h:104
 AliUEHistograms.h:105
 AliUEHistograms.h:106
 AliUEHistograms.h:107
 AliUEHistograms.h:108
 AliUEHistograms.h:109
 AliUEHistograms.h:110
 AliUEHistograms.h:111
 AliUEHistograms.h:112
 AliUEHistograms.h:113
 AliUEHistograms.h:114
 AliUEHistograms.h:115
 AliUEHistograms.h:116
 AliUEHistograms.h:117
 AliUEHistograms.h:118
 AliUEHistograms.h:119
 AliUEHistograms.h:120
 AliUEHistograms.h:121
 AliUEHistograms.h:122
 AliUEHistograms.h:123
 AliUEHistograms.h:124
 AliUEHistograms.h:125
 AliUEHistograms.h:126
 AliUEHistograms.h:127
 AliUEHistograms.h:128
 AliUEHistograms.h:129
 AliUEHistograms.h:130
 AliUEHistograms.h:131
 AliUEHistograms.h:132
 AliUEHistograms.h:133
 AliUEHistograms.h:134
 AliUEHistograms.h:135
 AliUEHistograms.h:136
 AliUEHistograms.h:137
 AliUEHistograms.h:138
 AliUEHistograms.h:139
 AliUEHistograms.h:140
 AliUEHistograms.h:141
 AliUEHistograms.h:142
 AliUEHistograms.h:143
 AliUEHistograms.h:144
 AliUEHistograms.h:145
 AliUEHistograms.h:146
 AliUEHistograms.h:147
 AliUEHistograms.h:148
 AliUEHistograms.h:149
 AliUEHistograms.h:150
 AliUEHistograms.h:151
 AliUEHistograms.h:152
 AliUEHistograms.h:153
 AliUEHistograms.h:154
 AliUEHistograms.h:155
 AliUEHistograms.h:156
 AliUEHistograms.h:157
 AliUEHistograms.h:158
 AliUEHistograms.h:159
 AliUEHistograms.h:160
 AliUEHistograms.h:161
 AliUEHistograms.h:162
 AliUEHistograms.h:163
 AliUEHistograms.h:164
 AliUEHistograms.h:165
 AliUEHistograms.h:166
 AliUEHistograms.h:167
 AliUEHistograms.h:168
 AliUEHistograms.h:169
 AliUEHistograms.h:170
 AliUEHistograms.h:171
 AliUEHistograms.h:172
 AliUEHistograms.h:173
 AliUEHistograms.h:174
 AliUEHistograms.h:175
 AliUEHistograms.h:176
 AliUEHistograms.h:177
 AliUEHistograms.h:178
 AliUEHistograms.h:179
 AliUEHistograms.h:180
 AliUEHistograms.h:181
 AliUEHistograms.h:182
 AliUEHistograms.h:183
 AliUEHistograms.h:184
 AliUEHistograms.h:185
 AliUEHistograms.h:186
 AliUEHistograms.h:187
 AliUEHistograms.h:188
 AliUEHistograms.h:189
 AliUEHistograms.h:190
 AliUEHistograms.h:191
 AliUEHistograms.h:192
 AliUEHistograms.h:193
 AliUEHistograms.h:194
 AliUEHistograms.h:195
 AliUEHistograms.h:196
 AliUEHistograms.h:197
 AliUEHistograms.h:198
 AliUEHistograms.h:199
 AliUEHistograms.h:200
 AliUEHistograms.h:201
 AliUEHistograms.h:202
 AliUEHistograms.h:203
 AliUEHistograms.h:204
 AliUEHistograms.h:205
 AliUEHistograms.h:206
 AliUEHistograms.h:207
 AliUEHistograms.h:208
 AliUEHistograms.h:209
 AliUEHistograms.h:210
 AliUEHistograms.h:211
 AliUEHistograms.h:212
 AliUEHistograms.h:213
 AliUEHistograms.h:214
 AliUEHistograms.h:215
 AliUEHistograms.h:216
 AliUEHistograms.h:217
 AliUEHistograms.h:218
 AliUEHistograms.h:219
 AliUEHistograms.h:220
 AliUEHistograms.h:221
 AliUEHistograms.h:222
 AliUEHistograms.h:223
 AliUEHistograms.h:224
 AliUEHistograms.h:225
 AliUEHistograms.h:226
 AliUEHistograms.h:227
 AliUEHistograms.h:228
 AliUEHistograms.h:229
 AliUEHistograms.h:230
 AliUEHistograms.h:231
 AliUEHistograms.h:232
 AliUEHistograms.h:233
 AliUEHistograms.h:234
 AliUEHistograms.h:235
 AliUEHistograms.h:236
 AliUEHistograms.h:237
 AliUEHistograms.h:238
 AliUEHistograms.h:239
 AliUEHistograms.h:240
 AliUEHistograms.h:241
 AliUEHistograms.h:242
 AliUEHistograms.h:243
 AliUEHistograms.h:244
 AliUEHistograms.h:245
 AliUEHistograms.h:246
 AliUEHistograms.h:247
 AliUEHistograms.h:248
 AliUEHistograms.h:249
 AliUEHistograms.h:250
 AliUEHistograms.h:251
 AliUEHistograms.h:252
 AliUEHistograms.h:253
 AliUEHistograms.h:254
 AliUEHistograms.h:255
 AliUEHistograms.h:256
 AliUEHistograms.h:257
 AliUEHistograms.h:258
 AliUEHistograms.h:259
 AliUEHistograms.h:260
 AliUEHistograms.h:261
 AliUEHistograms.h:262
 AliUEHistograms.h:263
 AliUEHistograms.h:264
 AliUEHistograms.h:265
 AliUEHistograms.h:266
 AliUEHistograms.h:267
 AliUEHistograms.h:268
 AliUEHistograms.h:269
 AliUEHistograms.h:270