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.                  *
**************************************************************************/
//
// Class for TRD PID
// Class further contains TRD specific cuts and QA histograms 
//  
// Authors:
//   Jens Wiechula <Jens.Wiechula@cern.ch>
//   Markus Fasel <M.Fasel@gsi.de>  
// 
#include <TList.h>
#include <TString.h>

#include <AliVParticle.h>
#include <AliESDtrack.h>
#include <AliLog.h>
#include <AliPID.h>

#include "AliDielectronTRDpidCut.h"

ClassImp(AliDielectronTRDpidCut)

const Double_t AliDielectronTRDpidCut::fgkVerySmall = 1e-12;

//___________________________________________________________________
AliDielectronTRDpidCut::AliDielectronTRDpidCut() :
    AliAnalysisCuts()
  , fMinP(1.)
  , fElectronEfficiency(0.91)
  , fPIDMethod(kNN)
  , fRequirePIDbit(0)
{
  //
  // default  constructor
  // 
  memset(fThreshParams, 0, sizeof(Double_t) * kThreshParams);
  SetUseDefaultParameters();
}

//___________________________________________________________________
AliDielectronTRDpidCut::AliDielectronTRDpidCut(const char* name) :
    AliAnalysisCuts(name,name)
  , fMinP(1.)
  , fElectronEfficiency(0.91)
  , fPIDMethod(kNN)
  , fRequirePIDbit(0)
{
  //
  // default  constructor
  // 
  memset(fThreshParams, 0, sizeof(Double_t) * kThreshParams);
  SetUseDefaultParameters();
}

//___________________________________________________________________
AliDielectronTRDpidCut::AliDielectronTRDpidCut(const AliDielectronTRDpidCut &ref):
    AliAnalysisCuts(ref.GetName(),ref.GetTitle())
  , fMinP(ref.fMinP)
  , fElectronEfficiency(ref.fElectronEfficiency)
  , fPIDMethod(ref.fPIDMethod)
  , fRequirePIDbit(ref.fRequirePIDbit)
{
  //
  // Copy constructor
  //
  memset(fThreshParams, 0, sizeof(Double_t) * kThreshParams);
  ref.Copy(*this);
}

//___________________________________________________________________
AliDielectronTRDpidCut &AliDielectronTRDpidCut::operator=(const AliDielectronTRDpidCut &ref){
  //
  // Assignment operator
  //
  if(this != &ref){
    ref.Copy(*this);
  }
  return *this;
}

//___________________________________________________________________
void AliDielectronTRDpidCut::Copy(TObject &ref) const {
  //
  // Performs the copying of the object
  //
  AliDielectronTRDpidCut &target = dynamic_cast<AliDielectronTRDpidCut &>(ref);
  
  Bool_t defaultParameters = UseDefaultParameters();
  target.SetUseDefaultParameters(defaultParameters);
  target.fMinP = fMinP;
  target.fPIDMethod = fPIDMethod;
  target.fRequirePIDbit = fRequirePIDbit;
  target.fElectronEfficiency = fElectronEfficiency;
  memcpy(target.fThreshParams, fThreshParams, sizeof(Double_t) * kThreshParams);
  AliAnalysisCuts::Copy(ref);
}

//___________________________________________________________________
AliDielectronTRDpidCut::~AliDielectronTRDpidCut(){
  //
  // Destructor
  //
}

//______________________________________________________
Bool_t AliDielectronTRDpidCut::InitializePID(){
  //
  // InitializePID: Load TRD thresholds and create the electron efficiency axis
  // to navigate 
  //
  if(UseDefaultParameters()){
    if(fPIDMethod == kLQ)
      InitParameters1DLQ();
    else
      InitParameters();
  }
  return kTRUE;
}

//______________________________________________________
Bool_t AliDielectronTRDpidCut::IsSelected(TObject *track) {
  //
  // Does PID for TRD alone:
  // PID thresholds based on 90% Electron Efficiency level approximated by a linear 
  // step function
  //
  const AliESDtrack *part = dynamic_cast<const AliESDtrack *>(track);
  if (!part) return kFALSE;
  
  if (fRequirePIDbit==AliDielectronTRDpidCut::kRequire&&!(part->GetStatus()&AliESDtrack::kTRDpid)) return kFALSE;
  if (fRequirePIDbit==AliDielectronTRDpidCut::kIfAvailable&&!(part->GetStatus()&AliESDtrack::kTRDpid)) return kTRUE;
  
  Double_t p = GetP((AliVParticle*)track);
  if(p < fMinP){ 
    AliDebug(1, Form("Track momentum below %f", fMinP));
    if(fRequirePIDbit==AliDielectronTRDpidCut::kRequire) return kFALSE;
	if(fRequirePIDbit==AliDielectronTRDpidCut::kIfAvailable) return kTRUE;
  }

  Double_t electronLike = GetElectronLikelihood((AliVParticle*)track);
  Double_t threshold = GetTRDthresholds(fElectronEfficiency, p);
  AliDebug(1, Form("Threshold: %f\n", threshold));
  if(electronLike > threshold){
    return kTRUE;
  }
  return kFALSE;

}

//___________________________________________________________________
Double_t AliDielectronTRDpidCut::GetTRDthresholds(Double_t electronEff, Double_t p) const {
  //
  // Return momentum dependent and electron efficiency dependent TRD thresholds
  // 
  Double_t params[4];
  GetParameters(electronEff, params);
  Double_t threshold = 1. - params[0] - params[1] * p - params[2] * TMath::Exp(-params[3] * p);
  return TMath::Max(TMath::Min(threshold, 0.99), 0.2); // truncate the threshold upperwards to 0.999 and lowerwards to 0.2 and exclude unphysical values
}

//___________________________________________________________________
void AliDielectronTRDpidCut::SetThresholdParameters(Double_t electronEff, Double_t *params){
  //
  // Set threshold parameters for the given bin
  //
  if(electronEff >= 1. || electronEff < 0.7) return;
  Int_t effbin = static_cast<Int_t>((electronEff - 0.7)/0.05); 
  memcpy(&fThreshParams[effbin * 4], params, sizeof(Double_t) * 4); 
  SetUseDefaultParameters(kFALSE);
}

//___________________________________________________________________
void AliDielectronTRDpidCut::InitParameters(){
  //
  // Fill the Parameters into an array
  //

  AliDebug(2, "Loading threshold Parameter");
  // Parameters for 6 Layers
  fThreshParams[0] = -0.001839; // 0.7 electron eff
  fThreshParams[1] = 0.000276;
  fThreshParams[2] = 0.044902; 
  fThreshParams[3] = 1.726751;
  fThreshParams[4] = -0.002405; // 0.75 electron eff
  fThreshParams[5] = 0.000372;
  fThreshParams[6] = 0.061775;
  fThreshParams[7] = 1.739371;
  fThreshParams[8] = -0.003178; // 0.8 electron eff
  fThreshParams[9] = 0.000521;
  fThreshParams[10] = 0.087585;
  fThreshParams[11] = 1.749154;
  fThreshParams[12] = -0.004058; // 0.85 electron eff
  fThreshParams[13] = 0.000748;
  fThreshParams[14] = 0.129583;
  fThreshParams[15] = 1.782323;
  fThreshParams[16] = -0.004967; // 0.9 electron eff
  fThreshParams[17] = 0.001216;
  fThreshParams[18] = 0.210128;
  fThreshParams[19] = 1.807665;
  fThreshParams[20] = -0.000996; // 0.95 electron eff
  fThreshParams[21] = 0.002627;
  fThreshParams[22] = 0.409099;
  fThreshParams[23] = 1.787076;
}

//___________________________________________________________________
void AliDielectronTRDpidCut::InitParameters1DLQ(){
  //
  // Init Parameters for 1DLQ PID (M. Fasel, Sept. 6th, 2010)
  //

  // Parameters for 6 Layers
  AliDebug(2, Form("Loading threshold parameter for Method 1DLQ"));
  fThreshParams[0] = -0.02241; // 0.7 electron eff
  fThreshParams[1] = 0.05043;
  fThreshParams[2] = 0.7925; 
  fThreshParams[3] = 2.625;
  fThreshParams[4] = 0.07438; // 0.75 electron eff
  fThreshParams[5] = 0.05158;
  fThreshParams[6] = 2.864;
  fThreshParams[7] = 4.356;
  fThreshParams[8] = 0.1977; // 0.8 electron eff
  fThreshParams[9] = 0.05956;
  fThreshParams[10] = 2.853;
  fThreshParams[11] = 3.713;
  fThreshParams[12] = 0.5206; // 0.85 electron eff
  fThreshParams[13] = 0.03077;
  fThreshParams[14] = 2.966;
  fThreshParams[15] = 4.07;
  fThreshParams[16] = 0.8808; // 0.9 electron eff
  fThreshParams[17] = 0.002092;
  fThreshParams[18] = 1.17;
  fThreshParams[19] = 4.506;
  fThreshParams[20] = 1.; // 0.95 electron eff
  fThreshParams[21] = 0.;
  fThreshParams[22] = 0.;
  fThreshParams[23] = 0.;

}

//___________________________________________________________________
void AliDielectronTRDpidCut::RenormalizeElPi(const Double_t *likein, Double_t *likeout) const {
  //
  // Renormalize likelihoods for electrons and pions neglecting the 
  // likelihoods for protons, kaons and muons
  //
  memset(likeout, 0, sizeof(Double_t) * AliPID::kSPECIES);
  Double_t norm = likein[AliPID::kElectron] + likein[AliPID::kPion];
  if(norm == 0.) norm = 1.;   // Safety
  likeout[AliPID::kElectron] = likein[AliPID::kElectron] / norm;
  likeout[AliPID::kPion] = likein[AliPID::kPion] / norm;
}

//___________________________________________________________________
void AliDielectronTRDpidCut::GetParameters(Double_t electronEff, Double_t *parameters) const {
  //
  // return parameter set for the given efficiency bin
  //
  Int_t effbin = static_cast<Int_t>((electronEff - 0.7)/0.05);
  memcpy(parameters, fThreshParams + effbin * 4, sizeof(Double_t) * 4);
}

//___________________________________________________________________
Double_t AliDielectronTRDpidCut::GetElectronLikelihood(const AliVParticle *track) const {
  //
  // Get TRD likelihoods for ESD respectively AOD tracks
  //
  Double_t pidProbs[AliPID::kSPECIES]; memset(pidProbs, 0, sizeof(Double_t) * AliPID::kSPECIES);
  const AliESDtrack *esdtrack = dynamic_cast<const AliESDtrack *>(track);
  if(esdtrack) esdtrack->GetTRDpid(pidProbs);
  if(!IsRenormalizeElPi()) return pidProbs[AliPID::kElectron];
  Double_t probsNew[AliPID::kSPECIES];
  RenormalizeElPi(pidProbs, probsNew);
  return probsNew[AliPID::kElectron];
}

//___________________________________________________________________
Double_t AliDielectronTRDpidCut::GetP(const AliVParticle *track) const {
  //
  // Get the Momentum in the TRD
  //
  Double_t p = 0.;
  const AliESDtrack *esdtrack = dynamic_cast<const AliESDtrack *>(track);
  if(esdtrack) p = esdtrack->GetOuterParam() ? esdtrack->GetOuterParam()->P() : esdtrack->P();
  return p;
}

//___________________________________________________________________
Double_t AliDielectronTRDpidCut::GetChargeLayer(const AliVParticle *track, UInt_t layer) const {
  //
  // Get the Charge in a single TRD layer
  //
  if(layer >= 6) return 0.;
  Double_t charge = 0.;
  const AliESDtrack *esdtrack = dynamic_cast<const AliESDtrack *>(track);
  if(esdtrack)
    for(Int_t islice = 0; islice < esdtrack->GetNumberOfTRDslices(); islice++)
      charge += esdtrack->GetTRDslice(static_cast<UInt_t>(layer), islice);
  return charge;
}

//___________________________________________________________________
void AliDielectronTRDpidCut::GetTRDmomenta(const AliVParticle *track, Double_t *mom) const {
  //
  // Fill Array with momentum information at the TRD tracklet
  //
  const AliESDtrack *esdtrack = dynamic_cast<const AliESDtrack *>(track);
  if(esdtrack)
    for(Int_t itl = 0; itl < 6; itl++)
      mom[itl] = esdtrack->GetTRDmomentum(itl);
}

//___________________________________________________________________
Double_t AliDielectronTRDpidCut::GetTRDSignalV1(const AliESDtrack *track, Float_t truncation) const {
  //
  // Calculation of the TRD Signal via truncated mean
  // Method 1: Take all Slices available
  // cut out 0s
  // Order them in increasing order
  // Cut out the upper third
  // Calculate mean over the last 2/3 slices
  //
  const Int_t kNSlices = 48;
  const Int_t kSlicePerLayer = 7;
  // Weight the slice to equalize the MPV of the dQ/dl-distribution per slice to the one in the first slice
  // Pions are used as reference for the equalization
  const Double_t kWeightSlice[8] = {1., 2.122, 1.8, 1.635, 1.595, 1.614, 1.16, 7.0};
  AliDebug(3, Form("Number of Tracklets: %d\n", track->GetTRDntrackletsPID()));
  Double_t trdSlices[kNSlices], tmp[kNSlices];
  Int_t indices[48];
  Int_t icnt = 0;
  for(Int_t idet = 0; idet < 6; idet++)
    for(Int_t islice = 0; islice < kSlicePerLayer; islice++){
      AliDebug(2, Form("Chamber[%d], Slice[%d]: TRDSlice = %f", idet, islice, track->GetTRDslice(idet, islice)));
      if(TMath::Abs(track->GetTRDslice(idet, islice)) < fgkVerySmall) continue;;
      trdSlices[icnt++] = track->GetTRDslice(idet, islice) * kWeightSlice[islice];
    }
  AliDebug(1, Form("Number of Slices: %d\n", icnt));
  if(icnt < 6) return 0.;   // We need at least 6 Slices for the truncated mean
  TMath::Sort(icnt, trdSlices, indices, kFALSE);
  memcpy(tmp, trdSlices, sizeof(Double_t) * icnt);
  for(Int_t ien = 0; ien < icnt; ien++)
    trdSlices[ien] = tmp[indices[ien]];
  Double_t trdSignal = TMath::Mean(static_cast<Int_t>(static_cast<Float_t>(icnt) * truncation), trdSlices);
  Double_t mom = track->GetOuterParam() ? track->GetOuterParam()->P() : -1;
  AliDebug(3, Form("PID Meth. 1: p[%f], TRDSignal[%f]", mom, trdSignal));
  return trdSignal;
}

//___________________________________________________________________
Double_t AliDielectronTRDpidCut::GetTRDSignalV2(const AliESDtrack *track, Float_t truncation) const {
  //
  // Calculation of the TRD Signal via truncated mean
  // Method 2: Take only first 5 slices per chamber
  // Order them in increasing order
  // Cut out upper half 
  // Now do mean with the reamining 3 slices per chamber
  //
  const Double_t kWeightSlice[8] = {1., 2.122, 1.8, 1.635, 1.595, 1.614, 1.16, 7.0};
  const Int_t kLayers = 6;
  const Int_t kSlicesLow = 6;
  const Int_t kSlicesHigh = 1;
  Double_t trdSlicesLowTime[kLayers*kSlicesLow], trdSlicesRemaining[kLayers*(kSlicesHigh + kSlicesLow)];
  Int_t indices[kLayers*kSlicesLow];
  Int_t cntLowTime=0, cntRemaining = 0;
  for(Int_t idet = 0; idet < 6; idet++)
    for(Int_t islice = 0; islice < kSlicesLow+kSlicesHigh; islice++){
      if(TMath::Abs(track->GetTRDslice(idet, islice)) < fgkVerySmall) continue;;
      if(islice < kSlicesLow){
        AliDebug(3, Form("Part 1, Det[%d], Slice[%d], TRDSlice: %f", idet, islice, track->GetTRDslice(idet, islice)));
        trdSlicesLowTime[cntLowTime++] = track->GetTRDslice(idet, islice) * kWeightSlice[islice];
      } else{
        AliDebug(3, Form("Part 1, Det[%d], Slice[%d], TRDSlice: %f", idet, islice, track->GetTRDslice(idet, islice)));
        trdSlicesRemaining[cntRemaining++] = track->GetTRDslice(idet, islice) * kWeightSlice[islice];
      }
    }
  if(cntLowTime < 4 || cntRemaining < 2) return 0.; // Min. Number of Slices at high time is 2 (matches with 1 layer), for the truncated mean we need at least 4 Slices
  TMath::Sort(cntLowTime, trdSlicesLowTime, indices, kFALSE);
  // Fill the second array with the lower half of the first time bins
  for(Int_t ien = 0; ien < static_cast<Int_t>(static_cast<Float_t>(cntLowTime) * truncation); ien++)
    trdSlicesRemaining[cntRemaining++] = trdSlicesLowTime[indices[ien]];
  Double_t trdSignal = TMath::Mean(cntRemaining, trdSlicesRemaining);
  Double_t mom = track->GetOuterParam() ? track->GetOuterParam()->P() : -1;
  AliDebug(3, Form("PID Meth. 2: p[%f], TRDSignal[%f]", mom, trdSignal));
  return trdSignal;
}
 AliDielectronTRDpidCut.cxx:1
 AliDielectronTRDpidCut.cxx:2
 AliDielectronTRDpidCut.cxx:3
 AliDielectronTRDpidCut.cxx:4
 AliDielectronTRDpidCut.cxx:5
 AliDielectronTRDpidCut.cxx:6
 AliDielectronTRDpidCut.cxx:7
 AliDielectronTRDpidCut.cxx:8
 AliDielectronTRDpidCut.cxx:9
 AliDielectronTRDpidCut.cxx:10
 AliDielectronTRDpidCut.cxx:11
 AliDielectronTRDpidCut.cxx:12
 AliDielectronTRDpidCut.cxx:13
 AliDielectronTRDpidCut.cxx:14
 AliDielectronTRDpidCut.cxx:15
 AliDielectronTRDpidCut.cxx:16
 AliDielectronTRDpidCut.cxx:17
 AliDielectronTRDpidCut.cxx:18
 AliDielectronTRDpidCut.cxx:19
 AliDielectronTRDpidCut.cxx:20
 AliDielectronTRDpidCut.cxx:21
 AliDielectronTRDpidCut.cxx:22
 AliDielectronTRDpidCut.cxx:23
 AliDielectronTRDpidCut.cxx:24
 AliDielectronTRDpidCut.cxx:25
 AliDielectronTRDpidCut.cxx:26
 AliDielectronTRDpidCut.cxx:27
 AliDielectronTRDpidCut.cxx:28
 AliDielectronTRDpidCut.cxx:29
 AliDielectronTRDpidCut.cxx:30
 AliDielectronTRDpidCut.cxx:31
 AliDielectronTRDpidCut.cxx:32
 AliDielectronTRDpidCut.cxx:33
 AliDielectronTRDpidCut.cxx:34
 AliDielectronTRDpidCut.cxx:35
 AliDielectronTRDpidCut.cxx:36
 AliDielectronTRDpidCut.cxx:37
 AliDielectronTRDpidCut.cxx:38
 AliDielectronTRDpidCut.cxx:39
 AliDielectronTRDpidCut.cxx:40
 AliDielectronTRDpidCut.cxx:41
 AliDielectronTRDpidCut.cxx:42
 AliDielectronTRDpidCut.cxx:43
 AliDielectronTRDpidCut.cxx:44
 AliDielectronTRDpidCut.cxx:45
 AliDielectronTRDpidCut.cxx:46
 AliDielectronTRDpidCut.cxx:47
 AliDielectronTRDpidCut.cxx:48
 AliDielectronTRDpidCut.cxx:49
 AliDielectronTRDpidCut.cxx:50
 AliDielectronTRDpidCut.cxx:51
 AliDielectronTRDpidCut.cxx:52
 AliDielectronTRDpidCut.cxx:53
 AliDielectronTRDpidCut.cxx:54
 AliDielectronTRDpidCut.cxx:55
 AliDielectronTRDpidCut.cxx:56
 AliDielectronTRDpidCut.cxx:57
 AliDielectronTRDpidCut.cxx:58
 AliDielectronTRDpidCut.cxx:59
 AliDielectronTRDpidCut.cxx:60
 AliDielectronTRDpidCut.cxx:61
 AliDielectronTRDpidCut.cxx:62
 AliDielectronTRDpidCut.cxx:63
 AliDielectronTRDpidCut.cxx:64
 AliDielectronTRDpidCut.cxx:65
 AliDielectronTRDpidCut.cxx:66
 AliDielectronTRDpidCut.cxx:67
 AliDielectronTRDpidCut.cxx:68
 AliDielectronTRDpidCut.cxx:69
 AliDielectronTRDpidCut.cxx:70
 AliDielectronTRDpidCut.cxx:71
 AliDielectronTRDpidCut.cxx:72
 AliDielectronTRDpidCut.cxx:73
 AliDielectronTRDpidCut.cxx:74
 AliDielectronTRDpidCut.cxx:75
 AliDielectronTRDpidCut.cxx:76
 AliDielectronTRDpidCut.cxx:77
 AliDielectronTRDpidCut.cxx:78
 AliDielectronTRDpidCut.cxx:79
 AliDielectronTRDpidCut.cxx:80
 AliDielectronTRDpidCut.cxx:81
 AliDielectronTRDpidCut.cxx:82
 AliDielectronTRDpidCut.cxx:83
 AliDielectronTRDpidCut.cxx:84
 AliDielectronTRDpidCut.cxx:85
 AliDielectronTRDpidCut.cxx:86
 AliDielectronTRDpidCut.cxx:87
 AliDielectronTRDpidCut.cxx:88
 AliDielectronTRDpidCut.cxx:89
 AliDielectronTRDpidCut.cxx:90
 AliDielectronTRDpidCut.cxx:91
 AliDielectronTRDpidCut.cxx:92
 AliDielectronTRDpidCut.cxx:93
 AliDielectronTRDpidCut.cxx:94
 AliDielectronTRDpidCut.cxx:95
 AliDielectronTRDpidCut.cxx:96
 AliDielectronTRDpidCut.cxx:97
 AliDielectronTRDpidCut.cxx:98
 AliDielectronTRDpidCut.cxx:99
 AliDielectronTRDpidCut.cxx:100
 AliDielectronTRDpidCut.cxx:101
 AliDielectronTRDpidCut.cxx:102
 AliDielectronTRDpidCut.cxx:103
 AliDielectronTRDpidCut.cxx:104
 AliDielectronTRDpidCut.cxx:105
 AliDielectronTRDpidCut.cxx:106
 AliDielectronTRDpidCut.cxx:107
 AliDielectronTRDpidCut.cxx:108
 AliDielectronTRDpidCut.cxx:109
 AliDielectronTRDpidCut.cxx:110
 AliDielectronTRDpidCut.cxx:111
 AliDielectronTRDpidCut.cxx:112
 AliDielectronTRDpidCut.cxx:113
 AliDielectronTRDpidCut.cxx:114
 AliDielectronTRDpidCut.cxx:115
 AliDielectronTRDpidCut.cxx:116
 AliDielectronTRDpidCut.cxx:117
 AliDielectronTRDpidCut.cxx:118
 AliDielectronTRDpidCut.cxx:119
 AliDielectronTRDpidCut.cxx:120
 AliDielectronTRDpidCut.cxx:121
 AliDielectronTRDpidCut.cxx:122
 AliDielectronTRDpidCut.cxx:123
 AliDielectronTRDpidCut.cxx:124
 AliDielectronTRDpidCut.cxx:125
 AliDielectronTRDpidCut.cxx:126
 AliDielectronTRDpidCut.cxx:127
 AliDielectronTRDpidCut.cxx:128
 AliDielectronTRDpidCut.cxx:129
 AliDielectronTRDpidCut.cxx:130
 AliDielectronTRDpidCut.cxx:131
 AliDielectronTRDpidCut.cxx:132
 AliDielectronTRDpidCut.cxx:133
 AliDielectronTRDpidCut.cxx:134
 AliDielectronTRDpidCut.cxx:135
 AliDielectronTRDpidCut.cxx:136
 AliDielectronTRDpidCut.cxx:137
 AliDielectronTRDpidCut.cxx:138
 AliDielectronTRDpidCut.cxx:139
 AliDielectronTRDpidCut.cxx:140
 AliDielectronTRDpidCut.cxx:141
 AliDielectronTRDpidCut.cxx:142
 AliDielectronTRDpidCut.cxx:143
 AliDielectronTRDpidCut.cxx:144
 AliDielectronTRDpidCut.cxx:145
 AliDielectronTRDpidCut.cxx:146
 AliDielectronTRDpidCut.cxx:147
 AliDielectronTRDpidCut.cxx:148
 AliDielectronTRDpidCut.cxx:149
 AliDielectronTRDpidCut.cxx:150
 AliDielectronTRDpidCut.cxx:151
 AliDielectronTRDpidCut.cxx:152
 AliDielectronTRDpidCut.cxx:153
 AliDielectronTRDpidCut.cxx:154
 AliDielectronTRDpidCut.cxx:155
 AliDielectronTRDpidCut.cxx:156
 AliDielectronTRDpidCut.cxx:157
 AliDielectronTRDpidCut.cxx:158
 AliDielectronTRDpidCut.cxx:159
 AliDielectronTRDpidCut.cxx:160
 AliDielectronTRDpidCut.cxx:161
 AliDielectronTRDpidCut.cxx:162
 AliDielectronTRDpidCut.cxx:163
 AliDielectronTRDpidCut.cxx:164
 AliDielectronTRDpidCut.cxx:165
 AliDielectronTRDpidCut.cxx:166
 AliDielectronTRDpidCut.cxx:167
 AliDielectronTRDpidCut.cxx:168
 AliDielectronTRDpidCut.cxx:169
 AliDielectronTRDpidCut.cxx:170
 AliDielectronTRDpidCut.cxx:171
 AliDielectronTRDpidCut.cxx:172
 AliDielectronTRDpidCut.cxx:173
 AliDielectronTRDpidCut.cxx:174
 AliDielectronTRDpidCut.cxx:175
 AliDielectronTRDpidCut.cxx:176
 AliDielectronTRDpidCut.cxx:177
 AliDielectronTRDpidCut.cxx:178
 AliDielectronTRDpidCut.cxx:179
 AliDielectronTRDpidCut.cxx:180
 AliDielectronTRDpidCut.cxx:181
 AliDielectronTRDpidCut.cxx:182
 AliDielectronTRDpidCut.cxx:183
 AliDielectronTRDpidCut.cxx:184
 AliDielectronTRDpidCut.cxx:185
 AliDielectronTRDpidCut.cxx:186
 AliDielectronTRDpidCut.cxx:187
 AliDielectronTRDpidCut.cxx:188
 AliDielectronTRDpidCut.cxx:189
 AliDielectronTRDpidCut.cxx:190
 AliDielectronTRDpidCut.cxx:191
 AliDielectronTRDpidCut.cxx:192
 AliDielectronTRDpidCut.cxx:193
 AliDielectronTRDpidCut.cxx:194
 AliDielectronTRDpidCut.cxx:195
 AliDielectronTRDpidCut.cxx:196
 AliDielectronTRDpidCut.cxx:197
 AliDielectronTRDpidCut.cxx:198
 AliDielectronTRDpidCut.cxx:199
 AliDielectronTRDpidCut.cxx:200
 AliDielectronTRDpidCut.cxx:201
 AliDielectronTRDpidCut.cxx:202
 AliDielectronTRDpidCut.cxx:203
 AliDielectronTRDpidCut.cxx:204
 AliDielectronTRDpidCut.cxx:205
 AliDielectronTRDpidCut.cxx:206
 AliDielectronTRDpidCut.cxx:207
 AliDielectronTRDpidCut.cxx:208
 AliDielectronTRDpidCut.cxx:209
 AliDielectronTRDpidCut.cxx:210
 AliDielectronTRDpidCut.cxx:211
 AliDielectronTRDpidCut.cxx:212
 AliDielectronTRDpidCut.cxx:213
 AliDielectronTRDpidCut.cxx:214
 AliDielectronTRDpidCut.cxx:215
 AliDielectronTRDpidCut.cxx:216
 AliDielectronTRDpidCut.cxx:217
 AliDielectronTRDpidCut.cxx:218
 AliDielectronTRDpidCut.cxx:219
 AliDielectronTRDpidCut.cxx:220
 AliDielectronTRDpidCut.cxx:221
 AliDielectronTRDpidCut.cxx:222
 AliDielectronTRDpidCut.cxx:223
 AliDielectronTRDpidCut.cxx:224
 AliDielectronTRDpidCut.cxx:225
 AliDielectronTRDpidCut.cxx:226
 AliDielectronTRDpidCut.cxx:227
 AliDielectronTRDpidCut.cxx:228
 AliDielectronTRDpidCut.cxx:229
 AliDielectronTRDpidCut.cxx:230
 AliDielectronTRDpidCut.cxx:231
 AliDielectronTRDpidCut.cxx:232
 AliDielectronTRDpidCut.cxx:233
 AliDielectronTRDpidCut.cxx:234
 AliDielectronTRDpidCut.cxx:235
 AliDielectronTRDpidCut.cxx:236
 AliDielectronTRDpidCut.cxx:237
 AliDielectronTRDpidCut.cxx:238
 AliDielectronTRDpidCut.cxx:239
 AliDielectronTRDpidCut.cxx:240
 AliDielectronTRDpidCut.cxx:241
 AliDielectronTRDpidCut.cxx:242
 AliDielectronTRDpidCut.cxx:243
 AliDielectronTRDpidCut.cxx:244
 AliDielectronTRDpidCut.cxx:245
 AliDielectronTRDpidCut.cxx:246
 AliDielectronTRDpidCut.cxx:247
 AliDielectronTRDpidCut.cxx:248
 AliDielectronTRDpidCut.cxx:249
 AliDielectronTRDpidCut.cxx:250
 AliDielectronTRDpidCut.cxx:251
 AliDielectronTRDpidCut.cxx:252
 AliDielectronTRDpidCut.cxx:253
 AliDielectronTRDpidCut.cxx:254
 AliDielectronTRDpidCut.cxx:255
 AliDielectronTRDpidCut.cxx:256
 AliDielectronTRDpidCut.cxx:257
 AliDielectronTRDpidCut.cxx:258
 AliDielectronTRDpidCut.cxx:259
 AliDielectronTRDpidCut.cxx:260
 AliDielectronTRDpidCut.cxx:261
 AliDielectronTRDpidCut.cxx:262
 AliDielectronTRDpidCut.cxx:263
 AliDielectronTRDpidCut.cxx:264
 AliDielectronTRDpidCut.cxx:265
 AliDielectronTRDpidCut.cxx:266
 AliDielectronTRDpidCut.cxx:267
 AliDielectronTRDpidCut.cxx:268
 AliDielectronTRDpidCut.cxx:269
 AliDielectronTRDpidCut.cxx:270
 AliDielectronTRDpidCut.cxx:271
 AliDielectronTRDpidCut.cxx:272
 AliDielectronTRDpidCut.cxx:273
 AliDielectronTRDpidCut.cxx:274
 AliDielectronTRDpidCut.cxx:275
 AliDielectronTRDpidCut.cxx:276
 AliDielectronTRDpidCut.cxx:277
 AliDielectronTRDpidCut.cxx:278
 AliDielectronTRDpidCut.cxx:279
 AliDielectronTRDpidCut.cxx:280
 AliDielectronTRDpidCut.cxx:281
 AliDielectronTRDpidCut.cxx:282
 AliDielectronTRDpidCut.cxx:283
 AliDielectronTRDpidCut.cxx:284
 AliDielectronTRDpidCut.cxx:285
 AliDielectronTRDpidCut.cxx:286
 AliDielectronTRDpidCut.cxx:287
 AliDielectronTRDpidCut.cxx:288
 AliDielectronTRDpidCut.cxx:289
 AliDielectronTRDpidCut.cxx:290
 AliDielectronTRDpidCut.cxx:291
 AliDielectronTRDpidCut.cxx:292
 AliDielectronTRDpidCut.cxx:293
 AliDielectronTRDpidCut.cxx:294
 AliDielectronTRDpidCut.cxx:295
 AliDielectronTRDpidCut.cxx:296
 AliDielectronTRDpidCut.cxx:297
 AliDielectronTRDpidCut.cxx:298
 AliDielectronTRDpidCut.cxx:299
 AliDielectronTRDpidCut.cxx:300
 AliDielectronTRDpidCut.cxx:301
 AliDielectronTRDpidCut.cxx:302
 AliDielectronTRDpidCut.cxx:303
 AliDielectronTRDpidCut.cxx:304
 AliDielectronTRDpidCut.cxx:305
 AliDielectronTRDpidCut.cxx:306
 AliDielectronTRDpidCut.cxx:307
 AliDielectronTRDpidCut.cxx:308
 AliDielectronTRDpidCut.cxx:309
 AliDielectronTRDpidCut.cxx:310
 AliDielectronTRDpidCut.cxx:311
 AliDielectronTRDpidCut.cxx:312
 AliDielectronTRDpidCut.cxx:313
 AliDielectronTRDpidCut.cxx:314
 AliDielectronTRDpidCut.cxx:315
 AliDielectronTRDpidCut.cxx:316
 AliDielectronTRDpidCut.cxx:317
 AliDielectronTRDpidCut.cxx:318
 AliDielectronTRDpidCut.cxx:319
 AliDielectronTRDpidCut.cxx:320
 AliDielectronTRDpidCut.cxx:321
 AliDielectronTRDpidCut.cxx:322
 AliDielectronTRDpidCut.cxx:323
 AliDielectronTRDpidCut.cxx:324
 AliDielectronTRDpidCut.cxx:325
 AliDielectronTRDpidCut.cxx:326
 AliDielectronTRDpidCut.cxx:327
 AliDielectronTRDpidCut.cxx:328
 AliDielectronTRDpidCut.cxx:329
 AliDielectronTRDpidCut.cxx:330
 AliDielectronTRDpidCut.cxx:331
 AliDielectronTRDpidCut.cxx:332
 AliDielectronTRDpidCut.cxx:333
 AliDielectronTRDpidCut.cxx:334
 AliDielectronTRDpidCut.cxx:335
 AliDielectronTRDpidCut.cxx:336
 AliDielectronTRDpidCut.cxx:337
 AliDielectronTRDpidCut.cxx:338
 AliDielectronTRDpidCut.cxx:339
 AliDielectronTRDpidCut.cxx:340
 AliDielectronTRDpidCut.cxx:341
 AliDielectronTRDpidCut.cxx:342
 AliDielectronTRDpidCut.cxx:343
 AliDielectronTRDpidCut.cxx:344
 AliDielectronTRDpidCut.cxx:345
 AliDielectronTRDpidCut.cxx:346
 AliDielectronTRDpidCut.cxx:347
 AliDielectronTRDpidCut.cxx:348
 AliDielectronTRDpidCut.cxx:349
 AliDielectronTRDpidCut.cxx:350
 AliDielectronTRDpidCut.cxx:351
 AliDielectronTRDpidCut.cxx:352
 AliDielectronTRDpidCut.cxx:353
 AliDielectronTRDpidCut.cxx:354
 AliDielectronTRDpidCut.cxx:355
 AliDielectronTRDpidCut.cxx:356
 AliDielectronTRDpidCut.cxx:357
 AliDielectronTRDpidCut.cxx:358
 AliDielectronTRDpidCut.cxx:359
 AliDielectronTRDpidCut.cxx:360
 AliDielectronTRDpidCut.cxx:361
 AliDielectronTRDpidCut.cxx:362
 AliDielectronTRDpidCut.cxx:363
 AliDielectronTRDpidCut.cxx:364
 AliDielectronTRDpidCut.cxx:365
 AliDielectronTRDpidCut.cxx:366
 AliDielectronTRDpidCut.cxx:367
 AliDielectronTRDpidCut.cxx:368
 AliDielectronTRDpidCut.cxx:369
 AliDielectronTRDpidCut.cxx:370
 AliDielectronTRDpidCut.cxx:371
 AliDielectronTRDpidCut.cxx:372
 AliDielectronTRDpidCut.cxx:373
 AliDielectronTRDpidCut.cxx:374
 AliDielectronTRDpidCut.cxx:375
 AliDielectronTRDpidCut.cxx:376
 AliDielectronTRDpidCut.cxx:377
 AliDielectronTRDpidCut.cxx:378
 AliDielectronTRDpidCut.cxx:379
 AliDielectronTRDpidCut.cxx:380
 AliDielectronTRDpidCut.cxx:381
 AliDielectronTRDpidCut.cxx:382
 AliDielectronTRDpidCut.cxx:383
 AliDielectronTRDpidCut.cxx:384
 AliDielectronTRDpidCut.cxx:385
 AliDielectronTRDpidCut.cxx:386
 AliDielectronTRDpidCut.cxx:387
 AliDielectronTRDpidCut.cxx:388
 AliDielectronTRDpidCut.cxx:389
 AliDielectronTRDpidCut.cxx:390
 AliDielectronTRDpidCut.cxx:391
 AliDielectronTRDpidCut.cxx:392
 AliDielectronTRDpidCut.cxx:393
 AliDielectronTRDpidCut.cxx:394
 AliDielectronTRDpidCut.cxx:395
 AliDielectronTRDpidCut.cxx:396
 AliDielectronTRDpidCut.cxx:397
 AliDielectronTRDpidCut.cxx:398