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

//---------------------------------------------------------------//
//        Base class for handling the pid response               //
//        functions of all detectors                             //
//        and give access to the nsigmas                         //
//                                                               //
//   Origin: Jens Wiechula, Uni Tuebingen, jens.wiechula@cern.ch //
//---------------------------------------------------------------//

#include "AliVParticle.h"
#include "AliVTrack.h"

#include "AliITSPIDResponse.h"
#include "AliTPCPIDResponse.h"
#include "AliTRDPIDResponse.h"
#include "AliTOFPIDResponse.h"
#include "AliHMPIDPIDResponse.h"
#include "AliEMCALPIDResponse.h"
#include "AliPID.h"

#include "TNamed.h"

class TF1;
class TObjArray;
class TLinearFitter;

class AliVEvent;
class AliTRDPIDResponseObject;
class AliTRDdEdxParams;
class AliTOFPIDParams;
class AliHMPIDPIDParams;
class AliOADBContainer;

class AliPIDResponse : public TNamed {
public:
  AliPIDResponse(Bool_t isMC=kFALSE);
  virtual ~AliPIDResponse();

  enum EDetector {
    kITS=0,
    kTPC=1,
    kTRD=2,
    kTOF=3,
    kHMPID=4,
    kEMCAL=5,
    kPHOS=6,
    kNdetectors=7
  };
  
  enum EDetCode {
    kDetITS = 0x1,
    kDetTPC = 0x2,
    kDetTRD = 0x4,
    kDetTOF = 0x8,
    kDetHMPID = 0x10,
    kDetEMCAL = 0x20,
    kDetPHOS = 0x40
  };

  enum EBeamType {
    kPP = 0,
    kPPB,
    kPBPB
  };
 
  enum EStartTimeType_t {kFILL_T0,kTOF_T0, kT0_T0, kBest_T0};

  enum ITSPIDmethod { kITSTruncMean, kITSLikelihood };

  enum EDetPidStatus {
    kDetNoSignal=0,
    kDetPidOk=1,
    kDetMismatch=2,
    kDetNoParams=3
  };

  AliITSPIDResponse &GetITSResponse() {return fITSResponse;}
  AliTPCPIDResponse &GetTPCResponse() {return fTPCResponse;}
  AliTOFPIDResponse &GetTOFResponse() {return fTOFResponse;}
  AliTRDPIDResponse &GetTRDResponse() {return fTRDResponse;}
  AliEMCALPIDResponse &GetEMCALResponse() {return fEMCALResponse;}

  // -----------------------------------------
  // buffered getters
  //
  
  // Number of sigmas
  EDetPidStatus NumberOfSigmas(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type, Double_t &val) const;
  
  Float_t NumberOfSigmas(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type) const;
  
  virtual Float_t NumberOfSigmasITS  (const AliVParticle *track, AliPID::EParticleType type) const;
  virtual Float_t NumberOfSigmasTPC  (const AliVParticle *track, AliPID::EParticleType type) const;
  virtual Float_t NumberOfSigmasTPC  (const AliVParticle *track, AliPID::EParticleType type, AliTPCPIDResponse::ETPCdEdxSource dedxSource) const;
  virtual Float_t NumberOfSigmasTRD  (const AliVParticle *track, AliPID::EParticleType type) const;
  virtual Float_t NumberOfSigmasEMCAL(const AliVParticle *track, AliPID::EParticleType type, Double_t &eop, Double_t showershape[4]) const;
  virtual Float_t NumberOfSigmasTOF  (const AliVParticle *track, AliPID::EParticleType type) const;
  virtual Float_t NumberOfSigmasTOF  (const AliVParticle *track, AliPID::EParticleType type, Float_t /*timeZeroTOF*/) const { return NumberOfSigmasTOF(track,type); }
  virtual Float_t NumberOfSigmasHMPID(const AliVParticle *track, AliPID::EParticleType type) const;
  virtual Float_t NumberOfSigmasEMCAL(const AliVParticle *track, AliPID::EParticleType type) const;

  Bool_t IdentifiedAsElectronTRD(const AliVTrack *track, Double_t efficiencyLevel,Double_t centrality=-1,AliTRDPIDResponse::ETRDPIDMethod PIDmethod=AliTRDPIDResponse::kLQ1D) const;
  Bool_t IdentifiedAsElectronTRD(const AliVTrack *track, Int_t &ntracklets, Double_t efficiencyLevel,Double_t centrality=-1,AliTRDPIDResponse::ETRDPIDMethod PIDmethod=AliTRDPIDResponse::kLQ1D) const;
  

  // Signal delta
  EDetPidStatus GetSignalDelta(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type, Double_t &val, Bool_t ratio=kFALSE) const;
  Double_t GetSignalDelta(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type, Bool_t ratio=kFALSE) const;
  
  // Probabilities
  EDetPidStatus ComputePIDProbability  (EDetCode  detCode, const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  EDetPidStatus ComputePIDProbability  (EDetector detCode, const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  
  virtual EDetPidStatus ComputeITSProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  virtual EDetPidStatus ComputeTPCProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  virtual EDetPidStatus ComputeTOFProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  virtual EDetPidStatus ComputeTRDProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  virtual EDetPidStatus ComputeEMCALProbability(const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  virtual EDetPidStatus ComputePHOSProbability (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  virtual EDetPidStatus ComputeHMPIDProbability(const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;

  virtual EDetPidStatus ComputeTRDProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[],AliTRDPIDResponse::ETRDPIDMethod PIDmethod) const;
  
  // pid status
  EDetPidStatus CheckPIDStatus(EDetector detCode, const AliVTrack *track)  const;

  AliTOFPIDParams *GetTOFPIDParams() const {return fTOFPIDParams;}
  Float_t GetTOFMismatchProbability(const AliVTrack *track = NULL) const; // if empty argument return the value stored during TOF probability computation

  void SetITSPIDmethod(ITSPIDmethod pmeth) { fITSPIDmethod = pmeth; }
  
  void SetOADBPath(const char* path) {fOADBPath=path;}
  const char *GetOADBPath() const {return fOADBPath.Data();}

  void SetCustomTPCpidResponse(const char* tpcpid) { fCustomTPCpidResponse = tpcpid; }
  const char* GetCustomTPCpidResponse() const { return fCustomTPCpidResponse.Data(); }
  
  void SetCustomTPCetaMaps(const char* tpcEtaMaps) { fCustomTPCetaMaps = tpcEtaMaps; }
  const char* GetCustomTPCetaMaps() const { return fCustomTPCetaMaps.Data(); }
  
  void InitialiseEvent(AliVEvent *event, Int_t pass, Int_t run=-1);
  void SetCurrentFile(const char* file) { fCurrentFile=file; }
  
  void SetCurrentAliRootRev(Int_t alirootRev) { fCurrentAliRootRev = alirootRev; }
  Int_t GetCurrentAliRootRev() const { return fCurrentAliRootRev; }

  // cache PID in the track
  void SetCachePID(Bool_t cache)    { fCachePID=cache;  }
  Bool_t GetCachePID() const { return fCachePID; }
  void FillTrackDetectorPID(const AliVTrack *track, EDetector detector) const;
  void FillTrackDetectorPID();

  AliVEvent * GetCurrentEvent() const {return fCurrentEvent;}

  // User settings for the MC period and reco pass
  void SetMCperiod(const char *mcPeriod) {fMCperiodUser=mcPeriod;}
  void SetRecoPass(Int_t recoPass)       {fRecoPassUser=recoPass;}

  // event info
  Float_t GetCurrentCentrality() const {return fCurrCentrality;};
  void SetCurrentCentrality(Float_t centrality) {fCurrCentrality=centrality;fEMCALResponse.SetCentrality(fCurrCentrality);};
  // TPC setting
  void SetUseTPCEtaCorrection(Bool_t useEtaCorrection = kTRUE) { fUseTPCEtaCorrection = useEtaCorrection; };
  Bool_t UseTPCEtaCorrection() const { return fUseTPCEtaCorrection; };
  
  void SetUseTPCMultiplicityCorrection(Bool_t useMultiplicityCorrection = kTRUE) { fUseTPCMultiplicityCorrection = useMultiplicityCorrection; };
  Bool_t UseTPCMultiplicityCorrection() const { return fUseTPCMultiplicityCorrection; };
  
  // TOF setting
  void SetTOFtail(Float_t tail=0.9){if(tail > 0) fTOFtail=tail; else printf("TOF tail should be greater than 0 (nothing done)\n");};
  void SetTOFResponse(AliVEvent *vevent,EStartTimeType_t option);

  virtual Float_t GetTPCsignalTunedOnData(const AliVTrack *t) const {return t->GetTPCsignal();};
  virtual Float_t GetTOFsignalTunedOnData(const AliVTrack *t) const {return t->GetTOFsignal();};
  Bool_t IsTunedOnData() const {return fTuneMConData;};
  void SetTunedOnData(Bool_t flag=kTRUE,Int_t recoPass=0){fTuneMConData = flag; if(recoPass>0) fRecoPassUser = recoPass;};
  Int_t GetTunedOnDataMask() const {return fTuneMConDataMask;};
  void SetTunedOnDataMask(Int_t detMask) {fTuneMConDataMask = detMask;}

  // Utilities
  TString GetChecksum(const TObject* obj) const;
    
  AliPIDResponse(const AliPIDResponse &other);
  AliPIDResponse& operator=(const AliPIDResponse &other);

  EBeamType GetBeamType() const {return fBeamTypeNum;};

  void SetNoTOFmism(Bool_t value=kTRUE){fNoTOFmism=value;};

protected:
  AliITSPIDResponse   fITSResponse;    //PID response function of the ITS
  AliTPCPIDResponse   fTPCResponse;    //PID response function of the TPC
  AliTRDPIDResponse   fTRDResponse;    //PID response function of the TRD
  AliTOFPIDResponse   fTOFResponse;    //PID response function of the TOF
  AliHMPIDPIDResponse fHMPIDResponse;  //PID response function of the HMPID
  AliEMCALPIDResponse fEMCALResponse;  //PID response function of the EMCAL

  Float_t           fRange;          // nSigma max in likelihood
  ITSPIDmethod      fITSPIDmethod;   // 0 = trunc mean; 1 = likelihood

  //unbuffered PID calculation
  virtual Float_t GetNumberOfSigmasTOFold  (const AliVParticle */*track*/, AliPID::EParticleType /*type*/) const {return 0;}
  virtual Float_t GetSignalDeltaTOFold(const AliVParticle */*track*/, AliPID::EParticleType /*type*/, Bool_t /*ratio*/=kFALSE) const {return -9999.;}

  Int_t CalculateTRDResponse(const AliVTrack *track, Double_t p[],AliTRDPIDResponse::ETRDPIDMethod PIDmethod) const;
  EDetPidStatus GetComputeTRDProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[],AliTRDPIDResponse::ETRDPIDMethod PIDmethod=AliTRDPIDResponse::kLQ1D) const;
  EDetPidStatus GetTOFPIDStatus(const AliVTrack *track) const;

  Bool_t fTuneMConData;                // switch to force the MC to be similar to data
  Int_t fTuneMConDataMask;             // select for which detectors enable MC tuning on data


private:
  static Float_t fgTOFmismatchProb;    // TOF mismatch probability (Bayesian)

  Bool_t fIsMC;                        //  If we run on MC data
  Bool_t fCachePID;

  TString fOADBPath;                   // OADB path to use
  TString fCustomTPCpidResponse;       // Custom TPC Pid Response file for debugging purposes
  TString fCustomTPCetaMaps;           // Custom TPC eta map file for debugging purposes
  
  TString fBeamType;                   //! beam type (PP) or (PBPB)
  TString fLHCperiod;                  //! LHC period
  TString fMCperiodTPC;                //! corresponding MC period to use for the TPC splines
  TString fMCperiodUser;               //  MC prodution requested by the user
  TString fCurrentFile;                //! name of currently processed file
  Int_t   fCurrentAliRootRev;          //! Aliroot rev. used to reconstruct the data
  Int_t   fRecoPass;                   //! reconstruction pass
  Int_t   fRecoPassUser;               //  reconstruction pass explicitly set by the user
  Int_t   fRun;                        //! current run number
  Int_t   fOldRun;                     //! current run number
  Float_t fResT0A;                     //! T0A resolution in current run
  Float_t fResT0C;                     //! T0C resolution in current run
  Float_t fResT0AC;                    //! T0A.and.T0C resolution in current run
  
  TObjArray *fArrPidResponseMaster;     //!  TPC pid splines
  TF1       *fResolutionCorrection;     //! TPC resolution correction
  AliOADBContainer* fOADBvoltageMaps;   //! container with the voltage maps
  Bool_t fUseTPCEtaCorrection;          // Use TPC eta correction
  Bool_t fUseTPCMultiplicityCorrection; // Use TPC multiplicity correction

  AliTRDPIDResponseObject *fTRDPIDResponseObject; //! TRD PID Response Object
  AliTRDdEdxParams * fTRDdEdxParams; //! TRD dEdx Response for truncated mean signal

  Float_t fTOFtail;                    //! TOF tail effect used in TOF probability
  AliTOFPIDParams *fTOFPIDParams;      //! TOF PID Params - period depending (OADB loaded)
  
  AliHMPIDPIDParams *fHMPIDPIDParams;  //! HMPID PID Params (OADB loaded)

  TObjArray *fEMCALPIDParams;             //! EMCAL PID Params

  AliVEvent *fCurrentEvent;            //! event currently being processed

  Float_t fCurrCentrality;             //! current centrality

  EBeamType fBeamTypeNum;              //! beam type enum 
  
  Bool_t fNoTOFmism;                   //! flag to switch off the TOF mismatch in the TOF weights (to check with old aliroot version)

  void ExecNewRun();
  
  //
  //setup parametrisations
  //

  //ITS
  void SetITSParametrisation();
  
  //TPC
  void SetTPCEtaMaps(Double_t refineFactorMapX = 6.0, Double_t refineFactorMapY = 6.0, Double_t refineFactorSigmaMapX = 6.0,
                     Double_t refineFactorSigmaMapY = 6.0);
  void SetTPCPidResponseMaster();
  void SetTPCParametrisation();
  Double_t GetTPCMultiplicityBin(const AliVEvent * const event);
  
  // TPC helpers for the eta maps
  void AddPointToHyperplane(TH2D* h, TLinearFitter* linExtrapolation, Int_t binX, Int_t binY);
  TH2D* RefineHistoViaLinearInterpolation(TH2D* h, Double_t refineFactorX = 6.0, Double_t refineFactorY = 6.0);

  //TRD
  void SetTRDPidResponseMaster();
  void InitializeTRDResponse();
  void SetTRDSlices(UInt_t TRDslicesForPID[2],AliTRDPIDResponse::ETRDPIDMethod method) const;
  void SetTRDdEdxParams();

  //TOF
  void SetTOFPidResponseMaster();
  void InitializeTOFResponse();

  //HMPID
  void SetHMPIDPidResponseMaster();
  void InitializeHMPIDResponse();

  //EMCAL
  void SetEMCALPidResponseMaster();
  void InitializeEMCALResponse();

  //
  void SetRecoInfo();

  //-------------------------------------------------
  //unbuffered PID calculation
  //
  
  // Number of sigmas
  Float_t GetNumberOfSigmas(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type) const;
  Float_t GetNumberOfSigmasITS  (const AliVParticle *track, AliPID::EParticleType type) const;
  Float_t GetNumberOfSigmasTPC  (const AliVParticle *track, AliPID::EParticleType type) const;
  Float_t GetNumberOfSigmasTRD  (const AliVParticle *track, AliPID::EParticleType type) const;
  Float_t GetNumberOfSigmasTOF  (const AliVParticle *track, AliPID::EParticleType type) const;
  Float_t GetNumberOfSigmasHMPID(const AliVParticle *track, AliPID::EParticleType type) const;
  Float_t GetNumberOfSigmasEMCAL(const AliVParticle *track, AliPID::EParticleType type, Double_t &eop, Double_t showershape[4]) const;
  Float_t GetNumberOfSigmasEMCAL(const AliVParticle *track, AliPID::EParticleType type) const;

  Float_t GetBufferedNumberOfSigmas(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type) const;

  // Signal deltas
  EDetPidStatus GetSignalDeltaITS(const AliVParticle *track, AliPID::EParticleType type, Double_t &val, Bool_t ratio=kFALSE) const;
  EDetPidStatus GetSignalDeltaTPC(const AliVParticle *track, AliPID::EParticleType type, Double_t &val, Bool_t ratio=kFALSE) const;
  EDetPidStatus GetSignalDeltaTRD(const AliVParticle *track, AliPID::EParticleType type, Double_t &val, Bool_t ratio=kFALSE) const;
  EDetPidStatus GetSignalDeltaTOF(const AliVParticle *track, AliPID::EParticleType type, Double_t &val, Bool_t ratio=kFALSE) const;
  EDetPidStatus GetSignalDeltaHMPID(const AliVParticle *vtrack, AliPID::EParticleType type, Double_t &val, Bool_t ratio=kFALSE) const;
  
  // Probabilities
  EDetPidStatus GetComputePIDProbability  (EDetector detCode,  const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  EDetPidStatus GetComputeITSProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  EDetPidStatus GetComputeTPCProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  EDetPidStatus GetComputeTOFProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[],Bool_t kNoMism=kFALSE) const;
  EDetPidStatus GetComputeEMCALProbability(const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  EDetPidStatus GetComputePHOSProbability (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;
  EDetPidStatus GetComputeHMPIDProbability(const AliVTrack *track, Int_t nSpecies, Double_t p[]) const;

  // pid status
  EDetPidStatus GetPIDStatus(EDetector det, const AliVTrack *track) const;
  EDetPidStatus GetITSPIDStatus(const AliVTrack *track) const;
  EDetPidStatus GetTPCPIDStatus(const AliVTrack *track) const;
  EDetPidStatus GetTRDPIDStatus(const AliVTrack *track) const;
  EDetPidStatus GetHMPIDPIDStatus(const AliVTrack *track) const;
  EDetPidStatus GetPHOSPIDStatus(const AliVTrack *track) const;
  EDetPidStatus GetEMCALPIDStatus(const AliVTrack *track) const;

  ClassDef(AliPIDResponse, 14);  //PID response handling
};

#endif
 AliPIDResponse.h:1
 AliPIDResponse.h:2
 AliPIDResponse.h:3
 AliPIDResponse.h:4
 AliPIDResponse.h:5
 AliPIDResponse.h:6
 AliPIDResponse.h:7
 AliPIDResponse.h:8
 AliPIDResponse.h:9
 AliPIDResponse.h:10
 AliPIDResponse.h:11
 AliPIDResponse.h:12
 AliPIDResponse.h:13
 AliPIDResponse.h:14
 AliPIDResponse.h:15
 AliPIDResponse.h:16
 AliPIDResponse.h:17
 AliPIDResponse.h:18
 AliPIDResponse.h:19
 AliPIDResponse.h:20
 AliPIDResponse.h:21
 AliPIDResponse.h:22
 AliPIDResponse.h:23
 AliPIDResponse.h:24
 AliPIDResponse.h:25
 AliPIDResponse.h:26
 AliPIDResponse.h:27
 AliPIDResponse.h:28
 AliPIDResponse.h:29
 AliPIDResponse.h:30
 AliPIDResponse.h:31
 AliPIDResponse.h:32
 AliPIDResponse.h:33
 AliPIDResponse.h:34
 AliPIDResponse.h:35
 AliPIDResponse.h:36
 AliPIDResponse.h:37
 AliPIDResponse.h:38
 AliPIDResponse.h:39
 AliPIDResponse.h:40
 AliPIDResponse.h:41
 AliPIDResponse.h:42
 AliPIDResponse.h:43
 AliPIDResponse.h:44
 AliPIDResponse.h:45
 AliPIDResponse.h:46
 AliPIDResponse.h:47
 AliPIDResponse.h:48
 AliPIDResponse.h:49
 AliPIDResponse.h:50
 AliPIDResponse.h:51
 AliPIDResponse.h:52
 AliPIDResponse.h:53
 AliPIDResponse.h:54
 AliPIDResponse.h:55
 AliPIDResponse.h:56
 AliPIDResponse.h:57
 AliPIDResponse.h:58
 AliPIDResponse.h:59
 AliPIDResponse.h:60
 AliPIDResponse.h:61
 AliPIDResponse.h:62
 AliPIDResponse.h:63
 AliPIDResponse.h:64
 AliPIDResponse.h:65
 AliPIDResponse.h:66
 AliPIDResponse.h:67
 AliPIDResponse.h:68
 AliPIDResponse.h:69
 AliPIDResponse.h:70
 AliPIDResponse.h:71
 AliPIDResponse.h:72
 AliPIDResponse.h:73
 AliPIDResponse.h:74
 AliPIDResponse.h:75
 AliPIDResponse.h:76
 AliPIDResponse.h:77
 AliPIDResponse.h:78
 AliPIDResponse.h:79
 AliPIDResponse.h:80
 AliPIDResponse.h:81
 AliPIDResponse.h:82
 AliPIDResponse.h:83
 AliPIDResponse.h:84
 AliPIDResponse.h:85
 AliPIDResponse.h:86
 AliPIDResponse.h:87
 AliPIDResponse.h:88
 AliPIDResponse.h:89
 AliPIDResponse.h:90
 AliPIDResponse.h:91
 AliPIDResponse.h:92
 AliPIDResponse.h:93
 AliPIDResponse.h:94
 AliPIDResponse.h:95
 AliPIDResponse.h:96
 AliPIDResponse.h:97
 AliPIDResponse.h:98
 AliPIDResponse.h:99
 AliPIDResponse.h:100
 AliPIDResponse.h:101
 AliPIDResponse.h:102
 AliPIDResponse.h:103
 AliPIDResponse.h:104
 AliPIDResponse.h:105
 AliPIDResponse.h:106
 AliPIDResponse.h:107
 AliPIDResponse.h:108
 AliPIDResponse.h:109
 AliPIDResponse.h:110
 AliPIDResponse.h:111
 AliPIDResponse.h:112
 AliPIDResponse.h:113
 AliPIDResponse.h:114
 AliPIDResponse.h:115
 AliPIDResponse.h:116
 AliPIDResponse.h:117
 AliPIDResponse.h:118
 AliPIDResponse.h:119
 AliPIDResponse.h:120
 AliPIDResponse.h:121
 AliPIDResponse.h:122
 AliPIDResponse.h:123
 AliPIDResponse.h:124
 AliPIDResponse.h:125
 AliPIDResponse.h:126
 AliPIDResponse.h:127
 AliPIDResponse.h:128
 AliPIDResponse.h:129
 AliPIDResponse.h:130
 AliPIDResponse.h:131
 AliPIDResponse.h:132
 AliPIDResponse.h:133
 AliPIDResponse.h:134
 AliPIDResponse.h:135
 AliPIDResponse.h:136
 AliPIDResponse.h:137
 AliPIDResponse.h:138
 AliPIDResponse.h:139
 AliPIDResponse.h:140
 AliPIDResponse.h:141
 AliPIDResponse.h:142
 AliPIDResponse.h:143
 AliPIDResponse.h:144
 AliPIDResponse.h:145
 AliPIDResponse.h:146
 AliPIDResponse.h:147
 AliPIDResponse.h:148
 AliPIDResponse.h:149
 AliPIDResponse.h:150
 AliPIDResponse.h:151
 AliPIDResponse.h:152
 AliPIDResponse.h:153
 AliPIDResponse.h:154
 AliPIDResponse.h:155
 AliPIDResponse.h:156
 AliPIDResponse.h:157
 AliPIDResponse.h:158
 AliPIDResponse.h:159
 AliPIDResponse.h:160
 AliPIDResponse.h:161
 AliPIDResponse.h:162
 AliPIDResponse.h:163
 AliPIDResponse.h:164
 AliPIDResponse.h:165
 AliPIDResponse.h:166
 AliPIDResponse.h:167
 AliPIDResponse.h:168
 AliPIDResponse.h:169
 AliPIDResponse.h:170
 AliPIDResponse.h:171
 AliPIDResponse.h:172
 AliPIDResponse.h:173
 AliPIDResponse.h:174
 AliPIDResponse.h:175
 AliPIDResponse.h:176
 AliPIDResponse.h:177
 AliPIDResponse.h:178
 AliPIDResponse.h:179
 AliPIDResponse.h:180
 AliPIDResponse.h:181
 AliPIDResponse.h:182
 AliPIDResponse.h:183
 AliPIDResponse.h:184
 AliPIDResponse.h:185
 AliPIDResponse.h:186
 AliPIDResponse.h:187
 AliPIDResponse.h:188
 AliPIDResponse.h:189
 AliPIDResponse.h:190
 AliPIDResponse.h:191
 AliPIDResponse.h:192
 AliPIDResponse.h:193
 AliPIDResponse.h:194
 AliPIDResponse.h:195
 AliPIDResponse.h:196
 AliPIDResponse.h:197
 AliPIDResponse.h:198
 AliPIDResponse.h:199
 AliPIDResponse.h:200
 AliPIDResponse.h:201
 AliPIDResponse.h:202
 AliPIDResponse.h:203
 AliPIDResponse.h:204
 AliPIDResponse.h:205
 AliPIDResponse.h:206
 AliPIDResponse.h:207
 AliPIDResponse.h:208
 AliPIDResponse.h:209
 AliPIDResponse.h:210
 AliPIDResponse.h:211
 AliPIDResponse.h:212
 AliPIDResponse.h:213
 AliPIDResponse.h:214
 AliPIDResponse.h:215
 AliPIDResponse.h:216
 AliPIDResponse.h:217
 AliPIDResponse.h:218
 AliPIDResponse.h:219
 AliPIDResponse.h:220
 AliPIDResponse.h:221
 AliPIDResponse.h:222
 AliPIDResponse.h:223
 AliPIDResponse.h:224
 AliPIDResponse.h:225
 AliPIDResponse.h:226
 AliPIDResponse.h:227
 AliPIDResponse.h:228
 AliPIDResponse.h:229
 AliPIDResponse.h:230
 AliPIDResponse.h:231
 AliPIDResponse.h:232
 AliPIDResponse.h:233
 AliPIDResponse.h:234
 AliPIDResponse.h:235
 AliPIDResponse.h:236
 AliPIDResponse.h:237
 AliPIDResponse.h:238
 AliPIDResponse.h:239
 AliPIDResponse.h:240
 AliPIDResponse.h:241
 AliPIDResponse.h:242
 AliPIDResponse.h:243
 AliPIDResponse.h:244
 AliPIDResponse.h:245
 AliPIDResponse.h:246
 AliPIDResponse.h:247
 AliPIDResponse.h:248
 AliPIDResponse.h:249
 AliPIDResponse.h:250
 AliPIDResponse.h:251
 AliPIDResponse.h:252
 AliPIDResponse.h:253
 AliPIDResponse.h:254
 AliPIDResponse.h:255
 AliPIDResponse.h:256
 AliPIDResponse.h:257
 AliPIDResponse.h:258
 AliPIDResponse.h:259
 AliPIDResponse.h:260
 AliPIDResponse.h:261
 AliPIDResponse.h:262
 AliPIDResponse.h:263
 AliPIDResponse.h:264
 AliPIDResponse.h:265
 AliPIDResponse.h:266
 AliPIDResponse.h:267
 AliPIDResponse.h:268
 AliPIDResponse.h:269
 AliPIDResponse.h:270
 AliPIDResponse.h:271
 AliPIDResponse.h:272
 AliPIDResponse.h:273
 AliPIDResponse.h:274
 AliPIDResponse.h:275
 AliPIDResponse.h:276
 AliPIDResponse.h:277
 AliPIDResponse.h:278
 AliPIDResponse.h:279
 AliPIDResponse.h:280
 AliPIDResponse.h:281
 AliPIDResponse.h:282
 AliPIDResponse.h:283
 AliPIDResponse.h:284
 AliPIDResponse.h:285
 AliPIDResponse.h:286
 AliPIDResponse.h:287
 AliPIDResponse.h:288
 AliPIDResponse.h:289
 AliPIDResponse.h:290
 AliPIDResponse.h:291
 AliPIDResponse.h:292
 AliPIDResponse.h:293
 AliPIDResponse.h:294
 AliPIDResponse.h:295
 AliPIDResponse.h:296
 AliPIDResponse.h:297
 AliPIDResponse.h:298
 AliPIDResponse.h:299
 AliPIDResponse.h:300
 AliPIDResponse.h:301
 AliPIDResponse.h:302
 AliPIDResponse.h:303
 AliPIDResponse.h:304
 AliPIDResponse.h:305
 AliPIDResponse.h:306
 AliPIDResponse.h:307
 AliPIDResponse.h:308
 AliPIDResponse.h:309
 AliPIDResponse.h:310
 AliPIDResponse.h:311
 AliPIDResponse.h:312
 AliPIDResponse.h:313
 AliPIDResponse.h:314
 AliPIDResponse.h:315
 AliPIDResponse.h:316
 AliPIDResponse.h:317
 AliPIDResponse.h:318
 AliPIDResponse.h:319
 AliPIDResponse.h:320
 AliPIDResponse.h:321
 AliPIDResponse.h:322
 AliPIDResponse.h:323
 AliPIDResponse.h:324
 AliPIDResponse.h:325
 AliPIDResponse.h:326
 AliPIDResponse.h:327
 AliPIDResponse.h:328
 AliPIDResponse.h:329
 AliPIDResponse.h:330
 AliPIDResponse.h:331
 AliPIDResponse.h:332
 AliPIDResponse.h:333
 AliPIDResponse.h:334
 AliPIDResponse.h:335
 AliPIDResponse.h:336
 AliPIDResponse.h:337
 AliPIDResponse.h:338
 AliPIDResponse.h:339
 AliPIDResponse.h:340
 AliPIDResponse.h:341
 AliPIDResponse.h:342
 AliPIDResponse.h:343
 AliPIDResponse.h:344
 AliPIDResponse.h:345
 AliPIDResponse.h:346
 AliPIDResponse.h:347
 AliPIDResponse.h:348
 AliPIDResponse.h:349
 AliPIDResponse.h:350