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


#include "AliHMPIDv1.h"     //class header
#include "AliHMPIDParam.h"  //StepManager()
#include "AliHMPIDHit.h"    //Hits2SDigs(),StepManager()
#include "AliHMPIDDigit.h"  //Digits2Raw(), Raw2SDigits()
#include "AliHMPIDRawStream.h"  //Digits2Raw(), Raw2SDigits()
#include "AliRawReader.h"  //Raw2SDigits()
#include <TVirtualMC.h>    //StepManager() for TVirtualMC::GetMC()
#include <TPDGCode.h>      //StepHistory() 
#include <AliStack.h>      //StepManager(),Hits2SDigits()
#include <AliLoader.h>        //Hits2SDigits()
#include <AliRunLoader.h>     //Hits2SDigits()
#include <AliMC.h>            //StepManager()      
#include <AliRun.h>           //CreateMaterials()    
#include <AliMagF.h>          //CreateMaterials()
//#include <TGeoManager.h>      //CreateGeometry()
#include <AliCDBEntry.h>      //CreateMaterials()
#include <AliCDBManager.h>    //CreateMaterials()
#include <TF1.h>              //DefineOpticalProperties()
#include <TF2.h>              //DefineOpticalProperties()
#include <TGeoGlobalMagField.h>
#include <TLorentzVector.h>   //IsLostByFresnel() 
#include <TTree.h>
 
ClassImp(AliHMPIDv1)    
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::AddAlignableVolumes()const
{
// Associates the symbolic volume name with the corresponding volume path. Interface method from AliModule invoked from AliMC
// Arguments: none
//   Returns: none   
  for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++)
    gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",i),Form("ALIC_1/HMPID_%i",i));
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::CreateMaterials()
{
// Definition of available HMPID materials  
// Arguments: none
//   Returns: none    
  AliDebug(1,"Start v1 HMPID.");
    
//data from PDG booklet 2002     density [gr/cm^3] rad len [cm] abs len [cm]    
  Float_t   aAir[4]={12,14,16,36}    ,   zAir[4]={6,7,8,18} ,   wAir[4]={0.000124,0.755267,0.231781,0.012827} , dAir=0.00120479; Int_t nAir=4;//mixture 0.9999999
  Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9}   , wC6F14[2]={6 , 14} , dC6F14=1.68    ; Int_t nC6F14=-2;
  Float_t  aSiO2[2]={ 28.09 , 15.99} ,  zSiO2[2]={14 , 8}   ,  wSiO2[2]={1 ,  2} ,  dSiO2=2.64    ; Int_t  nSiO2=-2; 
  Float_t   aCH4[2]={ 12.01 ,  1.01} ,   zCH4[2]={ 6 , 1}   ,   wCH4[2]={1 ,  4} ,   dCH4=7.17e-4 ; Int_t   nCH4=-2; 
  Float_t   aCsI[2]={132.90 ,126.90} ,   zCsI[2]={55 ,53}   ,   wCsI[2]={1 ,  1} ,   dCsI=0.1     ; Int_t   nCsI=-2; 
  Float_t     aRoha= 12.01 ,               zRoha=  6 ,                              dRoha=  0.10 , radRoha= 18.80 , absRoha=  86.3/dRoha; //special material- quazi carbon
  Float_t       aCu= 63.55 ,                 zCu= 29 ,                                dCu=  8.96 ,   radCu=  1.43 ,   absCu= 134.9/dCu  ;
  Float_t        aW=183.84 ,                  zW= 74 ,                                 dW= 19.30 ,    radW=  0.35 ,    absW= 185.0/dW   ;
  Float_t       aAl= 26.98 ,                 zAl= 13 ,                                dAl=  2.70 ,   radAl=  8.90 ,   absAl= 106.4/dAl  ;
  
    Int_t   matId=0;                           //tmp material id number
    Int_t   unsens =  0, sens=1;               //sensitive or unsensitive medium
    Int_t   itgfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ(); //type of field intergration 0 no field -1 user in guswim 1 Runge Kutta 2 helix 3 const field along z
    Float_t maxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();   //max field value
    Float_t tmaxfd = -10.0;                    //max deflection angle due to magnetic field in one step
    Float_t deemax = - 0.2;                    //max fractional energy loss in one step   
    Float_t stemax = - 0.1;                    //mas step allowed [cm]
    Float_t epsil  =   0.001;                  //abs tracking precision [cm]   
    Float_t stmin  = - 0.001;                  //min step size [cm] in continius process transport, negative value: choose it automatically
    AliMixture(++matId,"Air"  ,aAir  ,zAir  ,dAir  ,nAir  ,wAir  ); AliMedium(kAir  ,"Air"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
    AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);      
    AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);    
    AliMixture(++matId,"CH4"  ,aCH4  ,zCH4  ,dCH4  ,nCH4  ,wCH4  ); AliMedium(kCH4  ,"CH4"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);  
    AliMixture(++matId,"CsI"  ,aCsI  ,zCsI  ,dCsI  ,nCsI  ,wCsI  ); AliMedium(kCsI  ,"CsI"  ,matId,   sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
  
    AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha);  AliMedium(kRoha,"Roha", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
    AliMaterial(++matId,"Cu"  ,aCu  ,zCu  ,dCu  ,radCu  ,absCu  );  AliMedium(kCu  ,"Cu"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
    AliMaterial(++matId,"W"   ,aW   ,zW   ,dW   ,radW   ,absW   );  AliMedium(kW   ,"W"   , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
    AliMaterial(++matId,"Al"  ,aAl  ,zAl  ,dAl  ,radAl  ,absAl  );  AliMedium(kAl  ,"Al"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
    
//    DefineOpticalProperties(); // NOT TO BE CALLED BY USER CODE !!!
}//void AliHMPID::CreateMaterials()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::CreateGeometry()
{
//Creates detailed geometry simulation (currently GEANT volumes tree)         
  AliDebug(1,"Start main.");
  if(!TVirtualMC::GetMC()->IsRootGeometrySupported()) return;     
  
  Double_t cm=1,mm=0.1*cm,mkm=0.001*mm,dx,dy,dz;//default is cm
  
  TGeoVolume *pRich=gGeoManager->MakeBox("HMPID",gGeoManager->GetMedium("HMPID_CH4"),dx=(6*mm+1681*mm+6*mm)/2,  //main HMPID volume
                                                                                   dy=(6*mm+1466*mm+6*mm)/2,
                                                                                   dz=(80*mm+40*mm)*2/2);     //x,y taken from 2033P1  z from p84 TDR  
  for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//place 7 chambers
    TGeoHMatrix *pMatrix=new TGeoHMatrix;
    AliHMPIDParam::IdealPosition(iCh,pMatrix);
    gGeoManager->GetVolume("ALIC")->AddNode(pRich,iCh,pMatrix);
  }

  Float_t par[3];
  Int_t matrixIdReturn=0; //matrix id returned by AliMatrix
//Pad Panel frame  6 sectors
  par[0]=648*mm/2;par[1]=  411*mm/2;par[2]=40  *mm/2;TVirtualMC::GetMC()->Gsvolu("Rppf"     ,"BOX ",(*fIdtmed)[kAl]  ,par,3);//PPF 2001P2 inner size of the slab by 1mm more
  par[0]=181*mm/2;par[1]=89.25*mm/2;par[2]=38.3*mm/2;TVirtualMC::GetMC()->Gsvolu("RppfLarge","BOX ",(*fIdtmed)[kAir] ,par,3);//large whole
  par[0]=114*mm/2;par[1]=89.25*mm/2;par[2]=38.3*mm/2;TVirtualMC::GetMC()->Gsvolu("RppfSmall","BOX ",(*fIdtmed)[kAir] ,par,3);//small whole
  par[0]=644*mm/2;par[1]=  407*mm/2;par[2]= 1.7*mm/2;TVirtualMC::GetMC()->Gsvolu("Rpc"      ,"BOX ",(*fIdtmed)[kCsI] ,par,3);//by 0.2 mm more then actual size (PCB 2006P1)
  
  TVirtualMC::GetMC()->Gspos("Rppf",0,"HMPID",    -335*mm,      -433*mm,  8*cm+20*mm,  0,"ONLY");//F1 2040P1 z p.84 TDR
  TVirtualMC::GetMC()->Gspos("Rppf",1,"HMPID",    +335*mm,      -433*mm,  8*cm+20*mm,  0,"ONLY");
  TVirtualMC::GetMC()->Gspos("Rppf",2,"HMPID",    -335*mm,         0*mm,  8*cm+20*mm,  0,"ONLY");
  TVirtualMC::GetMC()->Gspos("Rppf",3,"HMPID",    +335*mm,         0*mm,  8*cm+20*mm,  0,"ONLY");
  TVirtualMC::GetMC()->Gspos("Rppf",4,"HMPID",    -335*mm,      +433*mm,  8*cm+20*mm,  0,"ONLY");
  TVirtualMC::GetMC()->Gspos("Rppf",5,"HMPID",    +335*mm,      +433*mm,  8*cm+20*mm,  0,"ONLY");  
    TVirtualMC::GetMC()->Gspos("Rpc"      ,1,"Rppf",       0*mm,         0*mm,   -19.15*mm,  0,"ONLY");//PPF 2001P2 
    TVirtualMC::GetMC()->Gspos("RppfLarge",1,"Rppf",  -224.5*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfLarge",2,"Rppf",  -224.5*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfLarge",3,"Rppf",  -224.5*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfLarge",4,"Rppf",  -224.5*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",1,"Rppf",  - 65.0*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",2,"Rppf",  - 65.0*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",3,"Rppf",  - 65.0*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",4,"Rppf",  - 65.0*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",5,"Rppf",  + 65.0*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",6,"Rppf",  + 65.0*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",7,"Rppf",  + 65.0*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfSmall",8,"Rppf",  + 65.0*mm,  +151.875*mm,     0.85*mm,  0,"ONLY"); 
    TVirtualMC::GetMC()->Gspos("RppfLarge",5,"Rppf",  +224.5*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfLarge",6,"Rppf",  +224.5*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfLarge",7,"Rppf",  +224.5*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
    TVirtualMC::GetMC()->Gspos("RppfLarge",8,"Rppf",  +224.5*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
//Gap - anod wires 6 copies to HMPID
  par[0]=648*mm/2;par[1]=  411*mm/2 ;par[2]=4.45*mm/2;TVirtualMC::GetMC()->Gsvolu("Rgap","BOX ",(*fIdtmed)[kCH4] ,par,3);//xy as PPF 2001P2 z WP 2099P1
  par[0]=  0*mm  ;par[1]=  20*mkm/2 ;par[2]= 648*mm/2;TVirtualMC::GetMC()->Gsvolu("Rano","TUBE",(*fIdtmed)[kW]   ,par,3);//WP 2099P1 z = gap x PPF 2001P2
  AliMatrix(matrixIdReturn,180,0, 90,90, 90,0); //wires along x
  
  TVirtualMC::GetMC()->Gspos("Rgap",0,"HMPID",    -335*mm,      -433*mm,8*cm-2.225*mm, 0,"ONLY"); //F1 2040P1 z WP 2099P1
  TVirtualMC::GetMC()->Gspos("Rgap",1,"HMPID",    +335*mm,      -433*mm,8*cm-2.225*mm, 0,"ONLY"); 
  TVirtualMC::GetMC()->Gspos("Rgap",2,"HMPID",    -335*mm,         0*mm,8*cm-2.225*mm, 0,"ONLY"); 
  TVirtualMC::GetMC()->Gspos("Rgap",3,"HMPID",    +335*mm,         0*mm,8*cm-2.225*mm, 0,"ONLY"); 
  TVirtualMC::GetMC()->Gspos("Rgap",4,"HMPID",    -335*mm,      +433*mm,8*cm-2.225*mm, 0,"ONLY"); 
  TVirtualMC::GetMC()->Gspos("Rgap",5,"HMPID",    +335*mm,      +433*mm,8*cm-2.225*mm, 0,"ONLY"); 
  for(int i=1;i<=96;i++)
    TVirtualMC::GetMC()->Gspos("Rano",i,"Rgap",     0*mm, -411/2*mm+i*4*mm, 0.185*mm, matrixIdReturn,"ONLY"); //WP 2099P1  
//Defines radiators geometry  
  par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=  24*mm/2;  TVirtualMC::GetMC()->Gsvolu("Rrad"      ,"BOX ",(*fIdtmed)[kC6F14]     ,par,3); // Rad 2011P1
  par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=   4*mm/2;  TVirtualMC::GetMC()->Gsvolu("RradFront" ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //front 
  par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=   5*mm/2;  TVirtualMC::GetMC()->Gsvolu("RradWin"   ,"BOX ",(*fIdtmed)[kSiO2]      ,par,3); //window
  par[0]=1330*mm/2 ;par[1]=   5*mm/2  ;par[2]=  15*mm/2;  TVirtualMC::GetMC()->Gsvolu("RradLong"  ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //long side  
  par[0]=  10*mm/2 ;par[1]= 403*mm/2  ;par[2]=  15*mm/2;  TVirtualMC::GetMC()->Gsvolu("RradShort" ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //short side 
  par[0]=   0      ;par[1]=  10*mm/2  ;par[2]=  15*mm/2;  TVirtualMC::GetMC()->Gsvolu("RradSpacer","TUBE",(*fIdtmed)[kSiO2]      ,par,3); //spacer        
    
  TVirtualMC::GetMC()->Gspos("Rrad",1,"HMPID",   0*mm,-434*mm,   -12*mm,  0,"ONLY"); //3 radiators to HMPID
  TVirtualMC::GetMC()->Gspos("Rrad",2,"HMPID",   0*mm,   0*mm,   -12*mm,  0,"ONLY"); 
  TVirtualMC::GetMC()->Gspos("Rrad",3,"HMPID",   0*mm,+434*mm,   -12*mm,  0,"ONLY"); 
    TVirtualMC::GetMC()->Gspos("RradFront",1,"Rrad",   0*mm,   0*mm, -10.0*mm,  0,"ONLY"); //front cover 
    TVirtualMC::GetMC()->Gspos("RradWin"  ,1,"Rrad",   0*mm,   0*mm,   9.5*mm,  0,"ONLY"); //quartz window (back cover)
    TVirtualMC::GetMC()->Gspos("RradLong" ,1,"Rrad",   0*mm,-204*mm,  -0.5*mm,  0,"ONLY"); //long side
    TVirtualMC::GetMC()->Gspos("RradLong" ,2,"Rrad",   0*mm,+204*mm,  -0.5*mm,  0,"ONLY"); //long side
    TVirtualMC::GetMC()->Gspos("RradShort",1,"Rrad",-660*mm,   0*mm,  -0.5*mm,  0,"ONLY"); //short side
    TVirtualMC::GetMC()->Gspos("RradShort",2,"Rrad",+660*mm,   0*mm,  -0.5*mm,  0,"ONLY"); //short side 
    for(int i=0;i<3;i++)
      for(int j=0;j<10;j++)
        TVirtualMC::GetMC()->Gspos("RradSpacer",10*i+j,"Rrad",-1330*mm/2+116*mm+j*122*mm,(i-1)*105*mm,-0.5*mm,0,"ONLY");//spacers
//Defines SandBox geometry
  par[0]=1419*mm/2 ;par[1]=1378*mm/2;par[2]=50.5*mm/2; TVirtualMC::GetMC()->Gsvolu("Rsb"     ,"BOX ",(*fIdtmed)[kAir]  ,par,3);  //2072P1   
  par[0]=1419*mm/2 ;par[1]=1378*mm/2;par[2]= 0.5*mm/2; TVirtualMC::GetMC()->Gsvolu("RsbCover","BOX ",(*fIdtmed)[kAl]   ,par,3);  //cover
  par[0]=1359*mm/2 ;par[1]=1318*mm/2;par[2]=49.5*mm/2; TVirtualMC::GetMC()->Gsvolu("RsbComb" ,"BOX ",(*fIdtmed)[kRoha] ,par,3);  //honeycomb structure 
  
  TVirtualMC::GetMC()->Gspos("Rsb",1,"HMPID",   0*mm, 0*mm, -73.75*mm, 0,"ONLY"); //p.84 TDR sandbox to rich
    TVirtualMC::GetMC()->Gspos("RsbComb" ,1,"Rsb", 0*mm, 0*mm,      0*mm, 0,"ONLY"); //2072P1 honeycomv to sandbox
    TVirtualMC::GetMC()->Gspos("RsbCover",1,"Rsb", 0*mm, 0*mm,    +25*mm, 0,"ONLY"); //cover to sandbox
    TVirtualMC::GetMC()->Gspos("RsbCover",2,"Rsb", 0*mm, 0*mm,    -25*mm, 0,"ONLY"); //cover to sandbox
  AliDebug(1,"Stop v1. HMPID option");  
}//CreateGeometry()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::Init()
{
// This methode defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(TVirtualMC::GetMC()->CurrentVolID()==XXX) statements in StepManager()
// Arguments: none
//   Returns: none      
  AliDebug(1,"Start v1 HMPID.");    
  fIdRad     = TVirtualMC::GetMC()->VolId("Rrad");
  fIdWin     = TVirtualMC::GetMC()->VolId("RradWin");
  fIdPc      = TVirtualMC::GetMC()->VolId("Rpc");
  fIdAmpGap  = TVirtualMC::GetMC()->VolId("Rgap");
  fIdProxGap = TVirtualMC::GetMC()->VolId("Rgap");

  AliDebug(1,"Stop v1 HMPID.");    
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::DefineOpticalProperties()
{
// Optical properties definition.
  const Int_t kNbins=30;       //number of photon energy points
  Float_t emin=5.5,emax=8.5;         //Photon energy range,[eV]
  Float_t aEckov [kNbins]; 
  Double_t dEckov [kNbins]; 
  Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
  Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins]; 
  Float_t                                                    aQeAll [kNbins], aQePc [kNbins];
  Double_t dReflMet[kNbins], dQePc[kNbins];

  TF2 *pRaIF=new TF2("HidxRad","sqrt(1+0.554*(1239.84/x)^2/((1239.84/x)^2-5769)-0.0005*(y-20))"                                       ,emin,emax,0,50); //DiMauro mail temp 0-50 degrees C
  TF1 *pWiIF=new TF1("HidxWin","sqrt(1+46.411/(10.666*10.666-x*x)+228.71/(18.125*18.125-x*x))"                                        ,emin,emax);      //SiO2 idx TDR p.35
  TF1 *pGaIF=new TF1("HidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)"                                                         ,emin,emax);      //?????? from where  

  TF1 *pRaAF=new TF1("HabsRad","(x<7.8)*(gaus+gaus(3))+(x>=7.8)*0.0001"                                                               ,emin,emax);  //fit from DiMauro data 28.10.03 
  pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
  TF1 *pWiAF=new TF1("HabsWin","(x<8.2)*(818.8638-301.0436*x+36.89642*x*x-1.507555*x*x*x)+(x>=8.2)*0.0001"                            ,emin,emax);  //fit from DiMauro data 28.10.03 
  TF1 *pGaAF=new TF1("HabsGap","(x<7.75)*6512.399+(x>=7.75)*3.90743e-2/(-1.655279e-1+6.307392e-2*x-8.011441e-3*x*x+3.392126e-4*x*x*x)",emin,emax);  //????? from where  
  
  TF1 *pQeF =new TF1("Hqe"    ,"0+(x>6.07267)*0.344811*(1-exp(-1.29730*(x-6.07267)))"                                                 ,emin,emax);  //fit from DiMauro data 28.10.03  
                            
  for(Int_t i=0;i<kNbins;i++){
    Float_t eV=emin+0.1*i;  //Ckov energy in eV
    aEckov [i] =1e-9*eV;    //Ckov energy in GeV
    dEckov [i] = aEckov[i];
    aAbsRad[i]=pRaAF->Eval(eV); aIdxRad[i]=1.292;//pRaIF->Eval(eV,20);      //Simulation for 20 degress C       
    aAbsWin[i]=pWiAF->Eval(eV); aIdxWin[i]=1.5787;//pWiIF->Eval(eV);
    aAbsGap[i]=pGaAF->Eval(eV); aIdxGap[i]=1.0005;//pGaIF->Eval(eV);   
    aQeAll[i] =1;                     //QE for all other materials except for PC must be 1.  
    aAbsMet[i] =0.0001;                aIdxMet[i]=0;                                             //metal ref idx must be 0 in order to reflect photon
                                       aIdxPc [i]=1;           aQePc [i]=pQeF->Eval(eV);         //PC ref idx must be 1 in order to apply photon to QE conversion 
    dQePc [i]=pQeF->Eval(eV);
    dReflMet[i] = 0.;     // no reflection on the surface of the pc (?)                                       
  }
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kC6F14]    , kNbins, aEckov, aAbsRad  , aQeAll , aIdxRad );    
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kSiO2]     , kNbins, aEckov, aAbsWin  , aQeAll , aIdxWin );    
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kCH4]      , kNbins, aEckov, aAbsGap  , aQeAll , aIdxGap );    
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kCu]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kW]        , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet ); //n=0 means reflect photons       
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kCsI]      , kNbins, aEckov, aAbsMet  , aQePc  , aIdxPc  ); //n=1 means convert photons    
  TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kAl]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    

  // Define a skin surface for the photocatode to enable 'detection' in G4
  TVirtualMC::GetMC()->DefineOpSurface("surfPc", kGlisur /*kUnified*/,kDielectric_metal,kPolished, 0.);
  TVirtualMC::GetMC()->SetMaterialProperty("surfPc", "EFFICIENCY", kNbins, dEckov, dQePc);
  TVirtualMC::GetMC()->SetMaterialProperty("surfPc", "REFLECTIVITY", kNbins, dEckov, dReflMet);
  TVirtualMC::GetMC()->SetSkinSurface("skinPc", "Rpc", "surfPc");

  delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Bool_t AliHMPIDv1::IsLostByFresnel()
{
// Calculate probability for the photon to be lost by Fresnel reflection.
  TLorentzVector p4;
  Double_t mom[3],localMom[3];
  TVirtualMC::GetMC()->TrackMomentum(p4);   mom[0]=p4(1);   mom[1]=p4(2);   mom[2]=p4(3);
  localMom[0]=0; localMom[1]=0; localMom[2]=0;
  TVirtualMC::GetMC()->Gmtod(mom,localMom,2);
  Double_t localTc    = localMom[0]*localMom[0]+localMom[2]*localMom[2];
  Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
  Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
  if(TVirtualMC::GetMC()->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
    AliDebug(1,"Photon lost");
    return kTRUE;
  }else
    return kFALSE;
}//IsLostByFresnel()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::GenFee(Float_t qtot)
{
// Generate FeedBack photons for the current particle. To be invoked from StepManager().
// eloss=0 means photon so only pulse height distribution is to be analysed.
  TLorentzVector x4;
  TVirtualMC::GetMC()->TrackPosition(x4); 
  Int_t iNphotons=TVirtualMC::GetMC()->GetRandom()->Poisson(0.02*qtot);  //# of feedback photons is proportional to the charge of hit
  AliDebug(1,Form("N photons=%i",iNphotons));
  Int_t j;
  Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
//Generate photons
  for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
    Double_t ranf[2];
    TVirtualMC::GetMC()->GetRandom()->RndmArray(2,ranf);    //Sample direction
    cthf=ranf[0]*2-1.0;
    if(cthf<0) continue;
    sthf = TMath::Sqrt((1. - cthf) * (1. + cthf));
    phif = ranf[1] * 2 * TMath::Pi();
    
    if(Double_t randomNumber=TVirtualMC::GetMC()->GetRandom()->Rndm()<=0.57)
      enfp = 7.5e-9;
    else if(randomNumber<=0.7)
      enfp = 6.4e-9;
    else
      enfp = 7.9e-9;
    

    dir[0] = sthf * TMath::Sin(phif);    dir[1] = cthf;    dir[2] = sthf * TMath::Cos(phif);
    TVirtualMC::GetMC()->Gdtom(dir, mom, 2);
    mom[0]*=enfp;    mom[1]*=enfp;    mom[2]*=enfp;
    mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
    
    // Polarisation
    e1[0]=      0;    e1[1]=-dir[2];    e1[2]= dir[1];
    e2[0]=-dir[1];    e2[1]= dir[0];    e2[2]=      0;
    e3[0]= dir[1];    e3[1]=      0;    e3[2]=-dir[0];
    
    vmod=0;
    for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
    if (!vmod) for(j=0;j<3;j++) {
      uswop=e1[j];
      e1[j]=e3[j];
      e3[j]=uswop;
    }
    vmod=0;
    for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
    if (!vmod) for(j=0;j<3;j++) {
      uswop=e2[j];
      e2[j]=e3[j];
      e3[j]=uswop;
    }
    
    vmod=0;  for(j=0;j<3;j++) vmod+=e1[j]*e1[j];  vmod=TMath::Sqrt(1/vmod);  for(j=0;j<3;j++) e1[j]*=vmod;    
    vmod=0;  for(j=0;j<3;j++) vmod+=e2[j]*e2[j];  vmod=TMath::Sqrt(1/vmod);  for(j=0;j<3;j++) e2[j]*=vmod;
    
    phi = TVirtualMC::GetMC()->GetRandom()->Rndm()* 2 * TMath::Pi();
    for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
    TVirtualMC::GetMC()->Gdtom(pol, pol, 2);
    Int_t outputNtracksStored;    
    gAlice->GetMCApp()->PushTrack(1,                             //transport
                     gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track 
                     50000051,                                   //PID
		     mom[0],mom[1],mom[2],mom[3],                //track momentum  
                     x4.X(),x4.Y(),x4.Z(),x4.T(),                //track origin 
                     pol[0],pol[1],pol[2],                       //polarization
		     kPFeedBackPhoton,                           //process ID   
                     outputNtracksStored,                        //on return how many new photons stored on stack
                     1.0);                                       //weight
  }//feedbacks loop
  AliDebug(1,"Stop.");
}//GenerateFeedbacks()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::Hits2SDigits()
{
// Interface method ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
// Arguments: none
//   Returns: none   
  AliDebug(1,"Start.");
  for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){                //events loop
    GetLoader()->GetRunLoader()->GetEvent(iEvt);                          //get next event
  
    if(!GetLoader()->TreeH()) {GetLoader()->LoadHits();                    }
    if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
          
    for(Int_t iEnt=0;iEnt<GetLoader()->TreeH()->GetEntries();iEnt++){//prims loop
      GetLoader()->TreeH()->GetEntry(iEnt);
      Hit2Sdi(Hits(),SdiLst());
    }//prims loop
    GetLoader()->TreeS()->Fill();
    GetLoader()->WriteSDigits("OVERWRITE");
    SdiReset();
  }//events loop  
  GetLoader()->UnloadHits();
  GetLoader()->UnloadSDigits();  
  AliDebug(1,"Stop.");
}//Hits2SDigits()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
{
// Converts list of hits to list of sdigits. 
// Arguments: pHitLst  - list of hits provided not empty
//            pSDigLst - list of sdigits where to store the results
//   Returns: none         
  for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){         //hits loop
    AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);         //get pointer to current hit   
    pHit->Hit2Sdi(pSdiLst);                                    //convert this hit to list of sdigits     
  }//hits loop loop
}//Hits2Sdi()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::Digits2Raw()
{
// Interface method invoked by AliSimulation to create raw data streams from digits. Events loop is done in AliSimulation
// Arguments: none
//   Returns: none    
  AliDebug(1,"Start.");
  GetLoader()->LoadDigits();
  TTree * treeD = GetLoader()->TreeD();
  if(!treeD) {
    AliError("No digits tree!");
    return;
  }
  treeD->GetEntry(0);
  
  //AliHMPIDDigit::WriteRaw(DigLst());
   AliHMPIDRawStream *pRS=0x0;
   pRS->WriteRaw(DigLst());
    
  GetLoader()->UnloadDigits();
  AliDebug(1,"Stop.");      
}//Digits2Raw()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Float_t AliHMPIDv1::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
{
// Correction for Fresnel   ???????????
// Arguments:   ene - photon energy [GeV],
//              PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
//   Returns:  
    Float_t en[36] = {5.0,5.1,5.2,5.3,5.4,5.5,5.6,5.7,5.8,5.9,6.0,6.1,6.2,
		      6.3,6.4,6.5,6.6,6.7,6.8,6.9,7.0,7.1,7.2,7.3,7.4,7.5,7.6,7.7,
		      7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
    Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
			2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
			2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
			1.72,1.53};
    Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
	 		0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
			0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
			1.714,1.498};
    Float_t xe=ene;
    Int_t  j=Int_t(xe*10)-49;
    Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
    Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);

    //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
    //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197

    Float_t sinin=TMath::Sqrt((1.-pdoti)*(1.+pdoti));
    Float_t tanin=sinin/pdoti;

    Float_t c1=cn*cn-ck*ck-sinin*sinin;
    Float_t c2=4*cn*cn*ck*ck;
    Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
    Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
    
    Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
    Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
    

    //CORRECTION FACTOR FOR SURFACE ROUGHNESS
    //B.J. STAGG  APPLIED OPTICS, 30(1991),4113

    Float_t sigraf=18.;
    Float_t lamb=1240/ene;
    Float_t fresn;
 
    Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));

    if(pola)
    {
	Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
	fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
    }
    else
	fresn=0.5*(rp+rs);
      
    fresn = fresn*rO;
    return fresn;
}//Fresnel()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::Print(Option_t *option)const
{
// Debug printout
  TObject::Print(option);
}//void AliHMPID::Print(Option_t *option)const
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Bool_t AliHMPIDv1::Raw2SDigits(AliRawReader *pRR)
{
// Interface methode ivoked from AliSimulation to create a list of sdigits from raw digits. Events loop is done in AliSimulation
// Arguments: pRR- raw reader 
//   Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())      
  //AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
  
  if(!GetLoader()->TreeS()) {MakeTree("S");  MakeBranch("S");}
    
  TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
  AliHMPIDRawStream stream(pRR);
  while(stream.Next())
  {
    for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
      AliHMPIDDigit sdi(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
      new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
    }
  }
  
  GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
  SdiReset();
  return kTRUE;
}//Raw2SDigits
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::StepCount()
{
// Count number of ckovs created  
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::StepHistory()
{
// This methode is invoked from StepManager() in order to print out 
  static Int_t iStepN;
  const char *sParticle;
  switch(TVirtualMC::GetMC()->TrackPid()){
    case kProton:      sParticle="PROTON"    ;break;
    case kNeutron:     sParticle="neutron"   ;break;
    case kGamma:       sParticle="gamma"     ;break;
    case 50000050:     sParticle="CKOV"      ;break;
    case kPi0:         sParticle="Pi0"       ;break;  
    case kPiPlus:      sParticle="Pi+"       ;break;  
    case kPiMinus:     sParticle="Pi-"       ;break;  
    case kElectron:    sParticle="electron"  ;break;  
    default:           sParticle="not known" ;break;
  }

  TString flag="fanny combination";
  if(TVirtualMC::GetMC()->IsTrackAlive()) {
    if(TVirtualMC::GetMC()->IsTrackEntering())      flag="enters to";
    else if(TVirtualMC::GetMC()->IsTrackExiting())  flag="exits from";
    else if(TVirtualMC::GetMC()->IsTrackInside())   flag="inside";
  } else {
    if(TVirtualMC::GetMC()->IsTrackStop())          flag="stopped in";        
  }
  
  Int_t vid=0,copy=0;
  TString path=TVirtualMC::GetMC()->CurrentVolName(); path.Prepend("-");path.Prepend(TVirtualMC::GetMC()->CurrentVolOffName(1));//current volume and his mother are always there
  vid=TVirtualMC::GetMC()->CurrentVolOffID(2,copy);  if(vid) {path.Prepend("-");path.Prepend(TVirtualMC::GetMC()->VolName(vid));}
  vid=TVirtualMC::GetMC()->CurrentVolOffID(3,copy);  if(vid) {path.Prepend("-");path.Prepend(TVirtualMC::GetMC()->VolName(vid));}
  
  Printf("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f",iStepN,sParticle,TVirtualMC::GetMC()->TrackPid(),flag.Data(),path.Data(),TVirtualMC::GetMC()->TrackMass(),TVirtualMC::GetMC()->TrackCharge(),TVirtualMC::GetMC()->Edep()*1e9);
  
  Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
                            iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
                            TVirtualMC::GetMC()->IsTrackAlive(), TVirtualMC::GetMC()->IsTrackDisappeared(),TVirtualMC::GetMC()->IsTrackEntering(), TVirtualMC::GetMC()->IsTrackExiting(),
                            TVirtualMC::GetMC()->IsTrackInside(),TVirtualMC::GetMC()->IsTrackOut(),        TVirtualMC::GetMC()->IsTrackStop(),     TVirtualMC::GetMC()->IsNewTrack());
  
  Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
  Int_t mid=TVirtualMC::GetMC()->CurrentMaterial(a,z,den,rad,abs);
  Printf("Step %i: id=%i a=%7.2f z=%7.2f den=%9.4f rad=%9.2f abs=%9.2f\n\n",iStepN,mid,a,z,den,rad,abs);
  iStepN++;
}//StepHistory()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDv1::StepManager()
{
// Full Step Manager.
// Arguments: none
//   Returns: none           
//  StepHistory(); return; //uncomment to print tracks history
//  StepCount(); return;     //uncomment to count photons
  
  Int_t   copy; //volume copy aka node
  
//Treat photons    
  if((TVirtualMC::GetMC()->TrackPid()==50000050||TVirtualMC::GetMC()->TrackPid()==50000051)&&TVirtualMC::GetMC()->CurrentVolID(copy)==fIdPc){  //photon (Ckov or feedback) hit PC (fIdPc)
    if(TVirtualMC::GetMC()->Edep()>0){                                                                           //photon survided QE test i.e. produces electron
      if(IsLostByFresnel()){ TVirtualMC::GetMC()->StopTrack(); return;}                                          //photon lost due to fersnel reflection on PC       
                       TVirtualMC::GetMC()->CurrentVolOffID(2,copy);                                             //current chamber since geomtry tree is HMPID-Rppf-Rpc
      Int_t   tid=     TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber();                                 //take TID
      Int_t   pid=     TVirtualMC::GetMC()->TrackPid();                                                          //take PID
      Float_t etot=    TVirtualMC::GetMC()->Etot();                                                              //total hpoton energy, [GeV] 
      Double_t x[3];   TVirtualMC::GetMC()->TrackPosition(x[0],x[1],x[2]);                                       //take MARS position at entrance to PC
      Float_t  hitTime=(Float_t)TVirtualMC::GetMC()->TrackTime();                                                //hit formation time
      Float_t xl,yl;   AliHMPIDParam::Instance()->Mars2Lors(copy,x,xl,yl);                       //take LORS position
      new((*fHits)[fNhits++])AliHMPIDHit(copy,etot,pid,tid,xl,yl,hitTime,x);                     //HIT for photon, position at P, etot will be set to Q
      GenFee(etot);                                                                              //generate feedback photons etot is modified in hit ctor to Q of hit
    }//photon hit PC and DE >0 
  }//photon hit PC
  
//Treat charged particles  
  static Float_t eloss;                                                                           //need to store mip parameters between different steps    
  static Double_t in[3];
  if(TVirtualMC::GetMC()->TrackCharge() && TVirtualMC::GetMC()->CurrentVolID(copy)==fIdAmpGap){                                   //charged particle in amplification gap (fIdAmpGap)
    if(TVirtualMC::GetMC()->IsTrackEntering()||TVirtualMC::GetMC()->IsNewTrack()) {                                               //entering or newly created
      eloss=0;                                                                                    //reset Eloss collector                         
      TVirtualMC::GetMC()->TrackPosition(in[0],in[1],in[2]);                                                      //take position at the entrance
    }else if(TVirtualMC::GetMC()->IsTrackExiting()||TVirtualMC::GetMC()->IsTrackStop()||TVirtualMC::GetMC()->IsTrackDisappeared()){               //exiting or disappeared
      eloss              +=TVirtualMC::GetMC()->Edep();                                                           //take into account last step Eloss
                          TVirtualMC::GetMC()->CurrentVolOffID(1,copy);                                           //take current chamber since geometry tree is HMPID-Rgap
      Int_t tid=          TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber();                               //take TID
      Int_t pid=          TVirtualMC::GetMC()->TrackPid();                                                        //take PID
      Double_t out[3];    TVirtualMC::GetMC()->TrackPosition(out[0],out[1],out[2]);                               //take MARS position at exit
      Float_t hitTime=    (Float_t)TVirtualMC::GetMC()->TrackTime();                                              //hit formation time
      out[0]=0.5*(out[0]+in[0]);                                                                  //>
      out[1]=0.5*(out[1]+in[1]);                                                                  //take hit position at the anod plane
      out[2]=0.5*(out[2]+in[2]);                                                                  //>
      Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(copy,out,xl,yl);                         //take LORS position
      new((*fHits)[fNhits++])AliHMPIDHit(copy,eloss,pid,tid,xl,yl,hitTime,out);                   //HIT for MIP, position near anod plane, eloss will be set to Q 
      GenFee(eloss);                                                                              //generate feedback photons 
    }else                                                                                         //just going inside
      eloss          += TVirtualMC::GetMC()->Edep();                                                              //collect this step eloss 
  }//MIP in GAP
}//StepManager()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 AliHMPIDv1.cxx:1
 AliHMPIDv1.cxx:2
 AliHMPIDv1.cxx:3
 AliHMPIDv1.cxx:4
 AliHMPIDv1.cxx:5
 AliHMPIDv1.cxx:6
 AliHMPIDv1.cxx:7
 AliHMPIDv1.cxx:8
 AliHMPIDv1.cxx:9
 AliHMPIDv1.cxx:10
 AliHMPIDv1.cxx:11
 AliHMPIDv1.cxx:12
 AliHMPIDv1.cxx:13
 AliHMPIDv1.cxx:14
 AliHMPIDv1.cxx:15
 AliHMPIDv1.cxx:16
 AliHMPIDv1.cxx:17
 AliHMPIDv1.cxx:18
 AliHMPIDv1.cxx:19
 AliHMPIDv1.cxx:20
 AliHMPIDv1.cxx:21
 AliHMPIDv1.cxx:22
 AliHMPIDv1.cxx:23
 AliHMPIDv1.cxx:24
 AliHMPIDv1.cxx:25
 AliHMPIDv1.cxx:26
 AliHMPIDv1.cxx:27
 AliHMPIDv1.cxx:28
 AliHMPIDv1.cxx:29
 AliHMPIDv1.cxx:30
 AliHMPIDv1.cxx:31
 AliHMPIDv1.cxx:32
 AliHMPIDv1.cxx:33
 AliHMPIDv1.cxx:34
 AliHMPIDv1.cxx:35
 AliHMPIDv1.cxx:36
 AliHMPIDv1.cxx:37
 AliHMPIDv1.cxx:38
 AliHMPIDv1.cxx:39
 AliHMPIDv1.cxx:40
 AliHMPIDv1.cxx:41
 AliHMPIDv1.cxx:42
 AliHMPIDv1.cxx:43
 AliHMPIDv1.cxx:44
 AliHMPIDv1.cxx:45
 AliHMPIDv1.cxx:46
 AliHMPIDv1.cxx:47
 AliHMPIDv1.cxx:48
 AliHMPIDv1.cxx:49
 AliHMPIDv1.cxx:50
 AliHMPIDv1.cxx:51
 AliHMPIDv1.cxx:52
 AliHMPIDv1.cxx:53
 AliHMPIDv1.cxx:54
 AliHMPIDv1.cxx:55
 AliHMPIDv1.cxx:56
 AliHMPIDv1.cxx:57
 AliHMPIDv1.cxx:58
 AliHMPIDv1.cxx:59
 AliHMPIDv1.cxx:60
 AliHMPIDv1.cxx:61
 AliHMPIDv1.cxx:62
 AliHMPIDv1.cxx:63
 AliHMPIDv1.cxx:64
 AliHMPIDv1.cxx:65
 AliHMPIDv1.cxx:66
 AliHMPIDv1.cxx:67
 AliHMPIDv1.cxx:68
 AliHMPIDv1.cxx:69
 AliHMPIDv1.cxx:70
 AliHMPIDv1.cxx:71
 AliHMPIDv1.cxx:72
 AliHMPIDv1.cxx:73
 AliHMPIDv1.cxx:74
 AliHMPIDv1.cxx:75
 AliHMPIDv1.cxx:76
 AliHMPIDv1.cxx:77
 AliHMPIDv1.cxx:78
 AliHMPIDv1.cxx:79
 AliHMPIDv1.cxx:80
 AliHMPIDv1.cxx:81
 AliHMPIDv1.cxx:82
 AliHMPIDv1.cxx:83
 AliHMPIDv1.cxx:84
 AliHMPIDv1.cxx:85
 AliHMPIDv1.cxx:86
 AliHMPIDv1.cxx:87
 AliHMPIDv1.cxx:88
 AliHMPIDv1.cxx:89
 AliHMPIDv1.cxx:90
 AliHMPIDv1.cxx:91
 AliHMPIDv1.cxx:92
 AliHMPIDv1.cxx:93
 AliHMPIDv1.cxx:94
 AliHMPIDv1.cxx:95
 AliHMPIDv1.cxx:96
 AliHMPIDv1.cxx:97
 AliHMPIDv1.cxx:98
 AliHMPIDv1.cxx:99
 AliHMPIDv1.cxx:100
 AliHMPIDv1.cxx:101
 AliHMPIDv1.cxx:102
 AliHMPIDv1.cxx:103
 AliHMPIDv1.cxx:104
 AliHMPIDv1.cxx:105
 AliHMPIDv1.cxx:106
 AliHMPIDv1.cxx:107
 AliHMPIDv1.cxx:108
 AliHMPIDv1.cxx:109
 AliHMPIDv1.cxx:110
 AliHMPIDv1.cxx:111
 AliHMPIDv1.cxx:112
 AliHMPIDv1.cxx:113
 AliHMPIDv1.cxx:114
 AliHMPIDv1.cxx:115
 AliHMPIDv1.cxx:116
 AliHMPIDv1.cxx:117
 AliHMPIDv1.cxx:118
 AliHMPIDv1.cxx:119
 AliHMPIDv1.cxx:120
 AliHMPIDv1.cxx:121
 AliHMPIDv1.cxx:122
 AliHMPIDv1.cxx:123
 AliHMPIDv1.cxx:124
 AliHMPIDv1.cxx:125
 AliHMPIDv1.cxx:126
 AliHMPIDv1.cxx:127
 AliHMPIDv1.cxx:128
 AliHMPIDv1.cxx:129
 AliHMPIDv1.cxx:130
 AliHMPIDv1.cxx:131
 AliHMPIDv1.cxx:132
 AliHMPIDv1.cxx:133
 AliHMPIDv1.cxx:134
 AliHMPIDv1.cxx:135
 AliHMPIDv1.cxx:136
 AliHMPIDv1.cxx:137
 AliHMPIDv1.cxx:138
 AliHMPIDv1.cxx:139
 AliHMPIDv1.cxx:140
 AliHMPIDv1.cxx:141
 AliHMPIDv1.cxx:142
 AliHMPIDv1.cxx:143
 AliHMPIDv1.cxx:144
 AliHMPIDv1.cxx:145
 AliHMPIDv1.cxx:146
 AliHMPIDv1.cxx:147
 AliHMPIDv1.cxx:148
 AliHMPIDv1.cxx:149
 AliHMPIDv1.cxx:150
 AliHMPIDv1.cxx:151
 AliHMPIDv1.cxx:152
 AliHMPIDv1.cxx:153
 AliHMPIDv1.cxx:154
 AliHMPIDv1.cxx:155
 AliHMPIDv1.cxx:156
 AliHMPIDv1.cxx:157
 AliHMPIDv1.cxx:158
 AliHMPIDv1.cxx:159
 AliHMPIDv1.cxx:160
 AliHMPIDv1.cxx:161
 AliHMPIDv1.cxx:162
 AliHMPIDv1.cxx:163
 AliHMPIDv1.cxx:164
 AliHMPIDv1.cxx:165
 AliHMPIDv1.cxx:166
 AliHMPIDv1.cxx:167
 AliHMPIDv1.cxx:168
 AliHMPIDv1.cxx:169
 AliHMPIDv1.cxx:170
 AliHMPIDv1.cxx:171
 AliHMPIDv1.cxx:172
 AliHMPIDv1.cxx:173
 AliHMPIDv1.cxx:174
 AliHMPIDv1.cxx:175
 AliHMPIDv1.cxx:176
 AliHMPIDv1.cxx:177
 AliHMPIDv1.cxx:178
 AliHMPIDv1.cxx:179
 AliHMPIDv1.cxx:180
 AliHMPIDv1.cxx:181
 AliHMPIDv1.cxx:182
 AliHMPIDv1.cxx:183
 AliHMPIDv1.cxx:184
 AliHMPIDv1.cxx:185
 AliHMPIDv1.cxx:186
 AliHMPIDv1.cxx:187
 AliHMPIDv1.cxx:188
 AliHMPIDv1.cxx:189
 AliHMPIDv1.cxx:190
 AliHMPIDv1.cxx:191
 AliHMPIDv1.cxx:192
 AliHMPIDv1.cxx:193
 AliHMPIDv1.cxx:194
 AliHMPIDv1.cxx:195
 AliHMPIDv1.cxx:196
 AliHMPIDv1.cxx:197
 AliHMPIDv1.cxx:198
 AliHMPIDv1.cxx:199
 AliHMPIDv1.cxx:200
 AliHMPIDv1.cxx:201
 AliHMPIDv1.cxx:202
 AliHMPIDv1.cxx:203
 AliHMPIDv1.cxx:204
 AliHMPIDv1.cxx:205
 AliHMPIDv1.cxx:206
 AliHMPIDv1.cxx:207
 AliHMPIDv1.cxx:208
 AliHMPIDv1.cxx:209
 AliHMPIDv1.cxx:210
 AliHMPIDv1.cxx:211
 AliHMPIDv1.cxx:212
 AliHMPIDv1.cxx:213
 AliHMPIDv1.cxx:214
 AliHMPIDv1.cxx:215
 AliHMPIDv1.cxx:216
 AliHMPIDv1.cxx:217
 AliHMPIDv1.cxx:218
 AliHMPIDv1.cxx:219
 AliHMPIDv1.cxx:220
 AliHMPIDv1.cxx:221
 AliHMPIDv1.cxx:222
 AliHMPIDv1.cxx:223
 AliHMPIDv1.cxx:224
 AliHMPIDv1.cxx:225
 AliHMPIDv1.cxx:226
 AliHMPIDv1.cxx:227
 AliHMPIDv1.cxx:228
 AliHMPIDv1.cxx:229
 AliHMPIDv1.cxx:230
 AliHMPIDv1.cxx:231
 AliHMPIDv1.cxx:232
 AliHMPIDv1.cxx:233
 AliHMPIDv1.cxx:234
 AliHMPIDv1.cxx:235
 AliHMPIDv1.cxx:236
 AliHMPIDv1.cxx:237
 AliHMPIDv1.cxx:238
 AliHMPIDv1.cxx:239
 AliHMPIDv1.cxx:240
 AliHMPIDv1.cxx:241
 AliHMPIDv1.cxx:242
 AliHMPIDv1.cxx:243
 AliHMPIDv1.cxx:244
 AliHMPIDv1.cxx:245
 AliHMPIDv1.cxx:246
 AliHMPIDv1.cxx:247
 AliHMPIDv1.cxx:248
 AliHMPIDv1.cxx:249
 AliHMPIDv1.cxx:250
 AliHMPIDv1.cxx:251
 AliHMPIDv1.cxx:252
 AliHMPIDv1.cxx:253
 AliHMPIDv1.cxx:254
 AliHMPIDv1.cxx:255
 AliHMPIDv1.cxx:256
 AliHMPIDv1.cxx:257
 AliHMPIDv1.cxx:258
 AliHMPIDv1.cxx:259
 AliHMPIDv1.cxx:260
 AliHMPIDv1.cxx:261
 AliHMPIDv1.cxx:262
 AliHMPIDv1.cxx:263
 AliHMPIDv1.cxx:264
 AliHMPIDv1.cxx:265
 AliHMPIDv1.cxx:266
 AliHMPIDv1.cxx:267
 AliHMPIDv1.cxx:268
 AliHMPIDv1.cxx:269
 AliHMPIDv1.cxx:270
 AliHMPIDv1.cxx:271
 AliHMPIDv1.cxx:272
 AliHMPIDv1.cxx:273
 AliHMPIDv1.cxx:274
 AliHMPIDv1.cxx:275
 AliHMPIDv1.cxx:276
 AliHMPIDv1.cxx:277
 AliHMPIDv1.cxx:278
 AliHMPIDv1.cxx:279
 AliHMPIDv1.cxx:280
 AliHMPIDv1.cxx:281
 AliHMPIDv1.cxx:282
 AliHMPIDv1.cxx:283
 AliHMPIDv1.cxx:284
 AliHMPIDv1.cxx:285
 AliHMPIDv1.cxx:286
 AliHMPIDv1.cxx:287
 AliHMPIDv1.cxx:288
 AliHMPIDv1.cxx:289
 AliHMPIDv1.cxx:290
 AliHMPIDv1.cxx:291
 AliHMPIDv1.cxx:292
 AliHMPIDv1.cxx:293
 AliHMPIDv1.cxx:294
 AliHMPIDv1.cxx:295
 AliHMPIDv1.cxx:296
 AliHMPIDv1.cxx:297
 AliHMPIDv1.cxx:298
 AliHMPIDv1.cxx:299
 AliHMPIDv1.cxx:300
 AliHMPIDv1.cxx:301
 AliHMPIDv1.cxx:302
 AliHMPIDv1.cxx:303
 AliHMPIDv1.cxx:304
 AliHMPIDv1.cxx:305
 AliHMPIDv1.cxx:306
 AliHMPIDv1.cxx:307
 AliHMPIDv1.cxx:308
 AliHMPIDv1.cxx:309
 AliHMPIDv1.cxx:310
 AliHMPIDv1.cxx:311
 AliHMPIDv1.cxx:312
 AliHMPIDv1.cxx:313
 AliHMPIDv1.cxx:314
 AliHMPIDv1.cxx:315
 AliHMPIDv1.cxx:316
 AliHMPIDv1.cxx:317
 AliHMPIDv1.cxx:318
 AliHMPIDv1.cxx:319
 AliHMPIDv1.cxx:320
 AliHMPIDv1.cxx:321
 AliHMPIDv1.cxx:322
 AliHMPIDv1.cxx:323
 AliHMPIDv1.cxx:324
 AliHMPIDv1.cxx:325
 AliHMPIDv1.cxx:326
 AliHMPIDv1.cxx:327
 AliHMPIDv1.cxx:328
 AliHMPIDv1.cxx:329
 AliHMPIDv1.cxx:330
 AliHMPIDv1.cxx:331
 AliHMPIDv1.cxx:332
 AliHMPIDv1.cxx:333
 AliHMPIDv1.cxx:334
 AliHMPIDv1.cxx:335
 AliHMPIDv1.cxx:336
 AliHMPIDv1.cxx:337
 AliHMPIDv1.cxx:338
 AliHMPIDv1.cxx:339
 AliHMPIDv1.cxx:340
 AliHMPIDv1.cxx:341
 AliHMPIDv1.cxx:342
 AliHMPIDv1.cxx:343
 AliHMPIDv1.cxx:344
 AliHMPIDv1.cxx:345
 AliHMPIDv1.cxx:346
 AliHMPIDv1.cxx:347
 AliHMPIDv1.cxx:348
 AliHMPIDv1.cxx:349
 AliHMPIDv1.cxx:350
 AliHMPIDv1.cxx:351
 AliHMPIDv1.cxx:352
 AliHMPIDv1.cxx:353
 AliHMPIDv1.cxx:354
 AliHMPIDv1.cxx:355
 AliHMPIDv1.cxx:356
 AliHMPIDv1.cxx:357
 AliHMPIDv1.cxx:358
 AliHMPIDv1.cxx:359
 AliHMPIDv1.cxx:360
 AliHMPIDv1.cxx:361
 AliHMPIDv1.cxx:362
 AliHMPIDv1.cxx:363
 AliHMPIDv1.cxx:364
 AliHMPIDv1.cxx:365
 AliHMPIDv1.cxx:366
 AliHMPIDv1.cxx:367
 AliHMPIDv1.cxx:368
 AliHMPIDv1.cxx:369
 AliHMPIDv1.cxx:370
 AliHMPIDv1.cxx:371
 AliHMPIDv1.cxx:372
 AliHMPIDv1.cxx:373
 AliHMPIDv1.cxx:374
 AliHMPIDv1.cxx:375
 AliHMPIDv1.cxx:376
 AliHMPIDv1.cxx:377
 AliHMPIDv1.cxx:378
 AliHMPIDv1.cxx:379
 AliHMPIDv1.cxx:380
 AliHMPIDv1.cxx:381
 AliHMPIDv1.cxx:382
 AliHMPIDv1.cxx:383
 AliHMPIDv1.cxx:384
 AliHMPIDv1.cxx:385
 AliHMPIDv1.cxx:386
 AliHMPIDv1.cxx:387
 AliHMPIDv1.cxx:388
 AliHMPIDv1.cxx:389
 AliHMPIDv1.cxx:390
 AliHMPIDv1.cxx:391
 AliHMPIDv1.cxx:392
 AliHMPIDv1.cxx:393
 AliHMPIDv1.cxx:394
 AliHMPIDv1.cxx:395
 AliHMPIDv1.cxx:396
 AliHMPIDv1.cxx:397
 AliHMPIDv1.cxx:398
 AliHMPIDv1.cxx:399
 AliHMPIDv1.cxx:400
 AliHMPIDv1.cxx:401
 AliHMPIDv1.cxx:402
 AliHMPIDv1.cxx:403
 AliHMPIDv1.cxx:404
 AliHMPIDv1.cxx:405
 AliHMPIDv1.cxx:406
 AliHMPIDv1.cxx:407
 AliHMPIDv1.cxx:408
 AliHMPIDv1.cxx:409
 AliHMPIDv1.cxx:410
 AliHMPIDv1.cxx:411
 AliHMPIDv1.cxx:412
 AliHMPIDv1.cxx:413
 AliHMPIDv1.cxx:414
 AliHMPIDv1.cxx:415
 AliHMPIDv1.cxx:416
 AliHMPIDv1.cxx:417
 AliHMPIDv1.cxx:418
 AliHMPIDv1.cxx:419
 AliHMPIDv1.cxx:420
 AliHMPIDv1.cxx:421
 AliHMPIDv1.cxx:422
 AliHMPIDv1.cxx:423
 AliHMPIDv1.cxx:424
 AliHMPIDv1.cxx:425
 AliHMPIDv1.cxx:426
 AliHMPIDv1.cxx:427
 AliHMPIDv1.cxx:428
 AliHMPIDv1.cxx:429
 AliHMPIDv1.cxx:430
 AliHMPIDv1.cxx:431
 AliHMPIDv1.cxx:432
 AliHMPIDv1.cxx:433
 AliHMPIDv1.cxx:434
 AliHMPIDv1.cxx:435
 AliHMPIDv1.cxx:436
 AliHMPIDv1.cxx:437
 AliHMPIDv1.cxx:438
 AliHMPIDv1.cxx:439
 AliHMPIDv1.cxx:440
 AliHMPIDv1.cxx:441
 AliHMPIDv1.cxx:442
 AliHMPIDv1.cxx:443
 AliHMPIDv1.cxx:444
 AliHMPIDv1.cxx:445
 AliHMPIDv1.cxx:446
 AliHMPIDv1.cxx:447
 AliHMPIDv1.cxx:448
 AliHMPIDv1.cxx:449
 AliHMPIDv1.cxx:450
 AliHMPIDv1.cxx:451
 AliHMPIDv1.cxx:452
 AliHMPIDv1.cxx:453
 AliHMPIDv1.cxx:454
 AliHMPIDv1.cxx:455
 AliHMPIDv1.cxx:456
 AliHMPIDv1.cxx:457
 AliHMPIDv1.cxx:458
 AliHMPIDv1.cxx:459
 AliHMPIDv1.cxx:460
 AliHMPIDv1.cxx:461
 AliHMPIDv1.cxx:462
 AliHMPIDv1.cxx:463
 AliHMPIDv1.cxx:464
 AliHMPIDv1.cxx:465
 AliHMPIDv1.cxx:466
 AliHMPIDv1.cxx:467
 AliHMPIDv1.cxx:468
 AliHMPIDv1.cxx:469
 AliHMPIDv1.cxx:470
 AliHMPIDv1.cxx:471
 AliHMPIDv1.cxx:472
 AliHMPIDv1.cxx:473
 AliHMPIDv1.cxx:474
 AliHMPIDv1.cxx:475
 AliHMPIDv1.cxx:476
 AliHMPIDv1.cxx:477
 AliHMPIDv1.cxx:478
 AliHMPIDv1.cxx:479
 AliHMPIDv1.cxx:480
 AliHMPIDv1.cxx:481
 AliHMPIDv1.cxx:482
 AliHMPIDv1.cxx:483
 AliHMPIDv1.cxx:484
 AliHMPIDv1.cxx:485
 AliHMPIDv1.cxx:486
 AliHMPIDv1.cxx:487
 AliHMPIDv1.cxx:488
 AliHMPIDv1.cxx:489
 AliHMPIDv1.cxx:490
 AliHMPIDv1.cxx:491
 AliHMPIDv1.cxx:492
 AliHMPIDv1.cxx:493
 AliHMPIDv1.cxx:494
 AliHMPIDv1.cxx:495
 AliHMPIDv1.cxx:496
 AliHMPIDv1.cxx:497
 AliHMPIDv1.cxx:498
 AliHMPIDv1.cxx:499
 AliHMPIDv1.cxx:500
 AliHMPIDv1.cxx:501
 AliHMPIDv1.cxx:502
 AliHMPIDv1.cxx:503
 AliHMPIDv1.cxx:504
 AliHMPIDv1.cxx:505
 AliHMPIDv1.cxx:506
 AliHMPIDv1.cxx:507
 AliHMPIDv1.cxx:508
 AliHMPIDv1.cxx:509
 AliHMPIDv1.cxx:510
 AliHMPIDv1.cxx:511
 AliHMPIDv1.cxx:512
 AliHMPIDv1.cxx:513
 AliHMPIDv1.cxx:514
 AliHMPIDv1.cxx:515
 AliHMPIDv1.cxx:516
 AliHMPIDv1.cxx:517
 AliHMPIDv1.cxx:518
 AliHMPIDv1.cxx:519
 AliHMPIDv1.cxx:520
 AliHMPIDv1.cxx:521
 AliHMPIDv1.cxx:522
 AliHMPIDv1.cxx:523
 AliHMPIDv1.cxx:524
 AliHMPIDv1.cxx:525
 AliHMPIDv1.cxx:526
 AliHMPIDv1.cxx:527
 AliHMPIDv1.cxx:528
 AliHMPIDv1.cxx:529
 AliHMPIDv1.cxx:530
 AliHMPIDv1.cxx:531
 AliHMPIDv1.cxx:532
 AliHMPIDv1.cxx:533
 AliHMPIDv1.cxx:534
 AliHMPIDv1.cxx:535
 AliHMPIDv1.cxx:536
 AliHMPIDv1.cxx:537
 AliHMPIDv1.cxx:538
 AliHMPIDv1.cxx:539
 AliHMPIDv1.cxx:540
 AliHMPIDv1.cxx:541
 AliHMPIDv1.cxx:542
 AliHMPIDv1.cxx:543
 AliHMPIDv1.cxx:544
 AliHMPIDv1.cxx:545
 AliHMPIDv1.cxx:546
 AliHMPIDv1.cxx:547
 AliHMPIDv1.cxx:548
 AliHMPIDv1.cxx:549
 AliHMPIDv1.cxx:550
 AliHMPIDv1.cxx:551
 AliHMPIDv1.cxx:552
 AliHMPIDv1.cxx:553
 AliHMPIDv1.cxx:554
 AliHMPIDv1.cxx:555
 AliHMPIDv1.cxx:556
 AliHMPIDv1.cxx:557
 AliHMPIDv1.cxx:558
 AliHMPIDv1.cxx:559
 AliHMPIDv1.cxx:560
 AliHMPIDv1.cxx:561
 AliHMPIDv1.cxx:562
 AliHMPIDv1.cxx:563
 AliHMPIDv1.cxx:564
 AliHMPIDv1.cxx:565
 AliHMPIDv1.cxx:566
 AliHMPIDv1.cxx:567
 AliHMPIDv1.cxx:568
 AliHMPIDv1.cxx:569
 AliHMPIDv1.cxx:570
 AliHMPIDv1.cxx:571
 AliHMPIDv1.cxx:572
 AliHMPIDv1.cxx:573
 AliHMPIDv1.cxx:574
 AliHMPIDv1.cxx:575
 AliHMPIDv1.cxx:576
 AliHMPIDv1.cxx:577
 AliHMPIDv1.cxx:578
 AliHMPIDv1.cxx:579
 AliHMPIDv1.cxx:580
 AliHMPIDv1.cxx:581
 AliHMPIDv1.cxx:582
 AliHMPIDv1.cxx:583
 AliHMPIDv1.cxx:584
 AliHMPIDv1.cxx:585
 AliHMPIDv1.cxx:586
 AliHMPIDv1.cxx:587
 AliHMPIDv1.cxx:588
 AliHMPIDv1.cxx:589