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.                  *
 **************************************************************************/

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// AliHMPIDPid                                                          //
//                                                                      //
// HMPID class to perfom particle identification                        //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

#include "AliHMPIDPid.h"       //class header
#include "AliHMPIDParam.h"     //class header
#include "AliHMPIDRecon.h"     //class header
#include <AliESDtrack.h>       //FindPid()
#include <TRandom.h>           //Resolution()

//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
AliHMPIDPid::AliHMPIDPid():TNamed("HMPIDrec","HMPIDPid")
{
//..
//init of data members
//..
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDPid::FindPid(AliESDtrack *pTrk,Double_t nmean,Int_t nsp,Double_t *prob)
{
// Calculates probability to be a electron-muon-pion-kaon-proton with the "amplitude" method
// from the given Cerenkov angle and momentum assuming no initial particle composition
// (i.e. apriory probability to be the particle of the given sort is the same for all sorts)

  AliPID *pPid = new AliPID();
  
  Double_t thetaCerExp = -999.;                                                                           
  if(pTrk->GetHMPIDsignal()<=0) thetaCerExp = pTrk->GetHMPIDsignal();                                                                           
  else                          thetaCerExp = pTrk->GetHMPIDsignal() - (Int_t)pTrk->GetHMPIDsignal();     //  measured thetaCherenkov
  
  if(thetaCerExp<=0){                                                                                     //HMPID does not find anything reasonable for this track, assign 0.2 for all species
    for(Int_t iPart=0;iPart<nsp;iPart++) prob[iPart]=1.0/(Float_t)nsp;
    delete pPid ; pPid=0x0; return;
  } 
  
  Double_t p[3] = {0}, pmod = 0;
  if(pTrk->GetOuterHmpPxPyPz(p))  pmod = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]+p[2]*p[2]);  // Momentum of the charged particle
  
  else {                                         
    for(Int_t iPart=0;iPart<nsp;iPart++) prob[iPart]=1.0/(Float_t)nsp;
    delete pPid ; pPid=0x0; return;
  } 
  
  Double_t hTot=0;                               // Initialize the total height of the amplitude method
  Double_t *h = new Double_t [nsp];              // number of charged particles to be considered

  Bool_t desert = kTRUE;                                                                                                     //  Flag to evaluate if ThetaC is far ("desert") from the given Gaussians
  
  for(Int_t iPart=0;iPart<nsp;iPart++){                                                                                      //  for each particle
    
    h[iPart] = 0;                                                                                                            //  reset the height
    Double_t mass = pPid->ParticleMass(iPart);                                                                             //  with the given mass
    Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(nmean*pmod);                   //  evaluate the theor. Theta Cherenkov
    if(cosThetaTh>1) continue;                                                                                               //  no light emitted, zero height
    Double_t thetaCerTh = TMath::ACos(cosThetaTh);                                                                           //  theoretical Theta Cherenkov
    Double_t sigmaRing = Resolution(iPart,thetaCerTh,pTrk);
    
    if(sigmaRing==0) continue; 
   
    if(TMath::Abs(thetaCerExp-thetaCerTh)<4*sigmaRing) desert = kFALSE;                                                                //   
    h[iPart] =TMath::Gaus(thetaCerTh,thetaCerExp,sigmaRing,kTRUE);
    hTot    +=h[iPart]; //total height of all theoretical heights for normalization
    
  }//species loop

  for(Int_t iPart=0;iPart<nsp;iPart++) {//species loop to assign probabilities
    
    if(!desert) prob[iPart]=h[iPart]/hTot;
    else prob[iPart]=1.0/(Float_t)nsp;            //all theoretical values are far away from experemental one
    
  }
  
  delete [] h;
  delete pPid ; pPid=0x0;
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Double_t AliHMPIDPid::Resolution(Int_t iPart, Double_t thetaCerTh, AliESDtrack *pTrk)
{
  AliHMPIDParam *pParam = AliHMPIDParam::Instance();
      
  AliHMPIDRecon rec;
  Float_t xPc,yPc,thRa,phRa;
  Float_t x, y;
  Int_t q, nph;
  pTrk->GetHMPIDtrk(xPc,yPc,thRa,phRa);
  pTrk->GetHMPIDmip(x,y,q,nph);
  
  Double_t xRa = xPc - (pParam->RadThick()+pParam->WinThick()+pParam->GapThick())*TMath::Cos(phRa)*TMath::Tan(thRa); //just linear extrapolation back to RAD
  Double_t yRa = yPc - (pParam->RadThick()+pParam->WinThick()+pParam->GapThick())*TMath::Sin(phRa)*TMath::Tan(thRa); //just linear extrapolation back to RAD
  rec.SetTrack(xRa,yRa,thRa,phRa);
  
  Double_t occupancy = pTrk->GetHMPIDoccupancy();
  Double_t thetaMax = TMath::ACos(1./pParam->MeanIdxRad());
  Int_t    nPhotsTh = (Int_t)(12.*TMath::Sin(thetaCerTh)*TMath::Sin(thetaCerTh)/(TMath::Sin(thetaMax)*TMath::Sin(thetaMax))+0.01);

  Double_t sigmatot = 0;
  Int_t nTrks = 20;

  for(Int_t iTrk=0;iTrk<nTrks;iTrk++) {
    Double_t invSigma = 0;
    Int_t nPhotsAcc = 0;
    
    Int_t nPhots = 0; 
    if(nph<nPhotsTh+TMath::Sqrt(nPhotsTh) && nph>nPhotsTh-TMath::Sqrt(nPhotsTh)) nPhots = nph;
    else nPhots = gRandom->Poisson(nPhotsTh);
  
    for(Int_t j=0;j<nPhots;j++){
      Double_t phi = gRandom->Rndm()*TMath::TwoPi();
      TVector2 pos; pos=rec.TracePhot(thetaCerTh,phi);
      if(!pParam->IsInside(pos.X(),pos.Y())) continue;
      if(pParam->IsInDead(pos.X(),pos.Y()))  continue;
      Double_t sigma2 = pParam->Sigma2(thRa,phRa,thetaCerTh,phi);//photon candidate sigma^2
      if(sigma2!=0) {
        invSigma += 1./sigma2;
        nPhotsAcc++;
      }
    }      
    if(invSigma!=0) sigmatot += 1./TMath::Sqrt(invSigma);  
  }
  
  return (sigmatot/nTrks)*pParam->SigmaCorrFact(iPart,occupancy);
}
 AliHMPIDPid.cxx:1
 AliHMPIDPid.cxx:2
 AliHMPIDPid.cxx:3
 AliHMPIDPid.cxx:4
 AliHMPIDPid.cxx:5
 AliHMPIDPid.cxx:6
 AliHMPIDPid.cxx:7
 AliHMPIDPid.cxx:8
 AliHMPIDPid.cxx:9
 AliHMPIDPid.cxx:10
 AliHMPIDPid.cxx:11
 AliHMPIDPid.cxx:12
 AliHMPIDPid.cxx:13
 AliHMPIDPid.cxx:14
 AliHMPIDPid.cxx:15
 AliHMPIDPid.cxx:16
 AliHMPIDPid.cxx:17
 AliHMPIDPid.cxx:18
 AliHMPIDPid.cxx:19
 AliHMPIDPid.cxx:20
 AliHMPIDPid.cxx:21
 AliHMPIDPid.cxx:22
 AliHMPIDPid.cxx:23
 AliHMPIDPid.cxx:24
 AliHMPIDPid.cxx:25
 AliHMPIDPid.cxx:26
 AliHMPIDPid.cxx:27
 AliHMPIDPid.cxx:28
 AliHMPIDPid.cxx:29
 AliHMPIDPid.cxx:30
 AliHMPIDPid.cxx:31
 AliHMPIDPid.cxx:32
 AliHMPIDPid.cxx:33
 AliHMPIDPid.cxx:34
 AliHMPIDPid.cxx:35
 AliHMPIDPid.cxx:36
 AliHMPIDPid.cxx:37
 AliHMPIDPid.cxx:38
 AliHMPIDPid.cxx:39
 AliHMPIDPid.cxx:40
 AliHMPIDPid.cxx:41
 AliHMPIDPid.cxx:42
 AliHMPIDPid.cxx:43
 AliHMPIDPid.cxx:44
 AliHMPIDPid.cxx:45
 AliHMPIDPid.cxx:46
 AliHMPIDPid.cxx:47
 AliHMPIDPid.cxx:48
 AliHMPIDPid.cxx:49
 AliHMPIDPid.cxx:50
 AliHMPIDPid.cxx:51
 AliHMPIDPid.cxx:52
 AliHMPIDPid.cxx:53
 AliHMPIDPid.cxx:54
 AliHMPIDPid.cxx:55
 AliHMPIDPid.cxx:56
 AliHMPIDPid.cxx:57
 AliHMPIDPid.cxx:58
 AliHMPIDPid.cxx:59
 AliHMPIDPid.cxx:60
 AliHMPIDPid.cxx:61
 AliHMPIDPid.cxx:62
 AliHMPIDPid.cxx:63
 AliHMPIDPid.cxx:64
 AliHMPIDPid.cxx:65
 AliHMPIDPid.cxx:66
 AliHMPIDPid.cxx:67
 AliHMPIDPid.cxx:68
 AliHMPIDPid.cxx:69
 AliHMPIDPid.cxx:70
 AliHMPIDPid.cxx:71
 AliHMPIDPid.cxx:72
 AliHMPIDPid.cxx:73
 AliHMPIDPid.cxx:74
 AliHMPIDPid.cxx:75
 AliHMPIDPid.cxx:76
 AliHMPIDPid.cxx:77
 AliHMPIDPid.cxx:78
 AliHMPIDPid.cxx:79
 AliHMPIDPid.cxx:80
 AliHMPIDPid.cxx:81
 AliHMPIDPid.cxx:82
 AliHMPIDPid.cxx:83
 AliHMPIDPid.cxx:84
 AliHMPIDPid.cxx:85
 AliHMPIDPid.cxx:86
 AliHMPIDPid.cxx:87
 AliHMPIDPid.cxx:88
 AliHMPIDPid.cxx:89
 AliHMPIDPid.cxx:90
 AliHMPIDPid.cxx:91
 AliHMPIDPid.cxx:92
 AliHMPIDPid.cxx:93
 AliHMPIDPid.cxx:94
 AliHMPIDPid.cxx:95
 AliHMPIDPid.cxx:96
 AliHMPIDPid.cxx:97
 AliHMPIDPid.cxx:98
 AliHMPIDPid.cxx:99
 AliHMPIDPid.cxx:100
 AliHMPIDPid.cxx:101
 AliHMPIDPid.cxx:102
 AliHMPIDPid.cxx:103
 AliHMPIDPid.cxx:104
 AliHMPIDPid.cxx:105
 AliHMPIDPid.cxx:106
 AliHMPIDPid.cxx:107
 AliHMPIDPid.cxx:108
 AliHMPIDPid.cxx:109
 AliHMPIDPid.cxx:110
 AliHMPIDPid.cxx:111
 AliHMPIDPid.cxx:112
 AliHMPIDPid.cxx:113
 AliHMPIDPid.cxx:114
 AliHMPIDPid.cxx:115
 AliHMPIDPid.cxx:116
 AliHMPIDPid.cxx:117
 AliHMPIDPid.cxx:118
 AliHMPIDPid.cxx:119
 AliHMPIDPid.cxx:120
 AliHMPIDPid.cxx:121
 AliHMPIDPid.cxx:122
 AliHMPIDPid.cxx:123
 AliHMPIDPid.cxx:124
 AliHMPIDPid.cxx:125
 AliHMPIDPid.cxx:126
 AliHMPIDPid.cxx:127
 AliHMPIDPid.cxx:128
 AliHMPIDPid.cxx:129
 AliHMPIDPid.cxx:130
 AliHMPIDPid.cxx:131
 AliHMPIDPid.cxx:132
 AliHMPIDPid.cxx:133
 AliHMPIDPid.cxx:134
 AliHMPIDPid.cxx:135
 AliHMPIDPid.cxx:136
 AliHMPIDPid.cxx:137
 AliHMPIDPid.cxx:138
 AliHMPIDPid.cxx:139
 AliHMPIDPid.cxx:140
 AliHMPIDPid.cxx:141