]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDv1.cxx
- fixing use of centrality class
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDv1.cxx
1 // **************************************************************************
2 // * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 // *                                                                        *
4 // * Author: The ALICE Off-line Project.                                    *
5 // * Contributors are mentioned in the code where appropriate.              *
6 // *                                                                        *
7 // * Permission to use, copy, modify and distribute this software and its   *
8 // * documentation strictly for non-commercial purposes is hereby granted   *
9 // * without fee, provided that the above copyright notice appears in all   *
10 // * copies and that both the copyright notice and this permission notice   *
11 // * appear in the supporting documentation. The authors make no claims     *
12 // * about the suitability of this software for any purpose. It is          *
13 // * provided "as is" without express or implied warranty.                  *
14 // **************************************************************************
15
16
17 #include "AliHMPIDv1.h"     //class header
18 #include "AliHMPIDParam.h"  //StepManager()
19 #include "AliHMPIDHit.h"    //Hits2SDigs(),StepManager()
20 #include "AliHMPIDDigit.h"  //Digits2Raw(), Raw2SDigits()
21 #include "AliHMPIDRawStream.h"  //Digits2Raw(), Raw2SDigits()
22 #include "AliRawReader.h"  //Raw2SDigits()
23 #include <TVirtualMC.h>    //StepManager() for gMC
24 #include <TPDGCode.h>      //StepHistory() 
25 #include <AliStack.h>      //StepManager(),Hits2SDigits()
26 #include <AliLoader.h>        //Hits2SDigits()
27 #include <AliRunLoader.h>     //Hits2SDigits()
28 #include <AliMC.h>            //StepManager()      
29 #include <AliRun.h>           //CreateMaterials()    
30 #include <AliMagF.h>          //CreateMaterials()
31 //#include <TGeoManager.h>      //CreateGeometry()
32 #include <AliCDBEntry.h>      //CreateMaterials()
33 #include <AliCDBManager.h>    //CreateMaterials()
34 #include <TF1.h>              //DefineOpticalProperties()
35 #include <TF2.h>              //DefineOpticalProperties()
36 #include <TGeoGlobalMagField.h>
37 #include <TLorentzVector.h>   //IsLostByFresnel() 
38  
39 ClassImp(AliHMPIDv1)    
40 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
41 void AliHMPIDv1::AddAlignableVolumes()const
42 {
43 // Associates the symbolic volume name with the corresponding volume path. Interface method from AliModule invoked from AliMC
44 // Arguments: none
45 //   Returns: none   
46   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++)
47     gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",i),Form("ALIC_1/HMPID_%i",i));
48 }
49 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
50 void AliHMPIDv1::CreateMaterials()
51 {
52 // Definition of available HMPID materials  
53 // Arguments: none
54 //   Returns: none    
55   AliDebug(1,"Start v1 HMPID.");
56     
57 //data from PDG booklet 2002     density [gr/cm^3] rad len [cm] abs len [cm]    
58   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
59   Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9}   , wC6F14[2]={6 , 14} , dC6F14=1.68    ; Int_t nC6F14=-2;
60   Float_t  aSiO2[2]={ 28.09 , 15.99} ,  zSiO2[2]={14 , 8}   ,  wSiO2[2]={1 ,  2} ,  dSiO2=2.64    ; Int_t  nSiO2=-2; 
61   Float_t   aCH4[2]={ 12.01 ,  1.01} ,   zCH4[2]={ 6 , 1}   ,   wCH4[2]={1 ,  4} ,   dCH4=7.17e-4 ; Int_t   nCH4=-2; 
62   Float_t   aCsI[2]={132.90 ,126.90} ,   zCsI[2]={55 ,53}   ,   wCsI[2]={1 ,  1} ,   dCsI=0.1     ; Int_t   nCsI=-2; 
63   Float_t     aRoha= 12.01 ,               zRoha=  6 ,                              dRoha=  0.10 , radRoha= 18.80 , absRoha=  86.3/dRoha; //special material- quazi carbon
64   Float_t       aCu= 63.55 ,                 zCu= 29 ,                                dCu=  8.96 ,   radCu=  1.43 ,   absCu= 134.9/dCu  ;
65   Float_t        aW=183.84 ,                  zW= 74 ,                                 dW= 19.30 ,    radW=  0.35 ,    absW= 185.0/dW   ;
66   Float_t       aAl= 26.98 ,                 zAl= 13 ,                                dAl=  2.70 ,   radAl=  8.90 ,   absAl= 106.4/dAl  ;
67   
68     Int_t   matId=0;                           //tmp material id number
69     Int_t   unsens =  0, sens=1;               //sensitive or unsensitive medium
70     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
71     Float_t maxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();   //max field value
72     Float_t tmaxfd = -10.0;                    //max deflection angle due to magnetic field in one step
73     Float_t deemax = - 0.2;                    //max fractional energy loss in one step   
74     Float_t stemax = - 0.1;                    //mas step allowed [cm]
75     Float_t epsil  =   0.001;                  //abs tracking precision [cm]   
76     Float_t stmin  = - 0.001;                  //min step size [cm] in continius process transport, negative value: choose it automatically
77     AliMixture(++matId,"Air"  ,aAir  ,zAir  ,dAir  ,nAir  ,wAir  ); AliMedium(kAir  ,"Air"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
78     AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);      
79     AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);    
80     AliMixture(++matId,"CH4"  ,aCH4  ,zCH4  ,dCH4  ,nCH4  ,wCH4  ); AliMedium(kCH4  ,"CH4"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);  
81     AliMixture(++matId,"CsI"  ,aCsI  ,zCsI  ,dCsI  ,nCsI  ,wCsI  ); AliMedium(kCsI  ,"CsI"  ,matId,   sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
82   
83     AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha);  AliMedium(kRoha,"Roha", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
84     AliMaterial(++matId,"Cu"  ,aCu  ,zCu  ,dCu  ,radCu  ,absCu  );  AliMedium(kCu  ,"Cu"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
85     AliMaterial(++matId,"W"   ,aW   ,zW   ,dW   ,radW   ,absW   );  AliMedium(kW   ,"W"   , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
86     AliMaterial(++matId,"Al"  ,aAl  ,zAl  ,dAl  ,radAl  ,absAl  );  AliMedium(kAl  ,"Al"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
87     
88 //    DefineOpticalProperties(); // NOT TO BE CALLED BY USER CODE !!!
89 }//void AliHMPID::CreateMaterials()
90 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
91 void AliHMPIDv1::CreateGeometry()
92 {
93 //Creates detailed geometry simulation (currently GEANT volumes tree)         
94   AliDebug(1,"Start main.");
95   if(!gMC->IsRootGeometrySupported()) return;     
96   
97   Double_t cm=1,mm=0.1*cm,mkm=0.001*mm,dx,dy,dz;//default is cm
98   
99   TGeoVolume *pRich=gGeoManager->MakeBox("HMPID",gGeoManager->GetMedium("HMPID_CH4"),dx=(6*mm+1681*mm+6*mm)/2,  //main HMPID volume
100                                                                                    dy=(6*mm+1466*mm+6*mm)/2,
101                                                                                    dz=(80*mm+40*mm)*2/2);     //x,y taken from 2033P1  z from p84 TDR  
102   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//place 7 chambers
103     TGeoHMatrix *pMatrix=new TGeoHMatrix;
104     AliHMPIDParam::IdealPosition(iCh,pMatrix);
105     gGeoManager->GetVolume("ALIC")->AddNode(pRich,iCh,pMatrix);
106   }
107
108   Float_t par[3];
109   Int_t matrixIdReturn=0; //matrix id returned by AliMatrix
110 //Pad Panel frame  6 sectors
111   par[0]=648*mm/2;par[1]=  411*mm/2;par[2]=40  *mm/2;gMC->Gsvolu("Rppf"     ,"BOX ",(*fIdtmed)[kAl]  ,par,3);//PPF 2001P2 inner size of the slab by 1mm more
112   par[0]=181*mm/2;par[1]=89.25*mm/2;par[2]=38.3*mm/2;gMC->Gsvolu("RppfLarge","BOX ",(*fIdtmed)[kAir] ,par,3);//large whole
113   par[0]=114*mm/2;par[1]=89.25*mm/2;par[2]=38.3*mm/2;gMC->Gsvolu("RppfSmall","BOX ",(*fIdtmed)[kAir] ,par,3);//small whole
114   par[0]=644*mm/2;par[1]=  407*mm/2;par[2]= 1.7*mm/2;gMC->Gsvolu("Rpc"      ,"BOX ",(*fIdtmed)[kCsI] ,par,3);//by 0.2 mm more then actual size (PCB 2006P1)
115   
116   gMC->Gspos("Rppf",0,"HMPID",    -335*mm,      -433*mm,  8*cm+20*mm,  0,"ONLY");//F1 2040P1 z p.84 TDR
117   gMC->Gspos("Rppf",1,"HMPID",    +335*mm,      -433*mm,  8*cm+20*mm,  0,"ONLY");
118   gMC->Gspos("Rppf",2,"HMPID",    -335*mm,         0*mm,  8*cm+20*mm,  0,"ONLY");
119   gMC->Gspos("Rppf",3,"HMPID",    +335*mm,         0*mm,  8*cm+20*mm,  0,"ONLY");
120   gMC->Gspos("Rppf",4,"HMPID",    -335*mm,      +433*mm,  8*cm+20*mm,  0,"ONLY");
121   gMC->Gspos("Rppf",5,"HMPID",    +335*mm,      +433*mm,  8*cm+20*mm,  0,"ONLY");  
122     gMC->Gspos("Rpc"      ,1,"Rppf",       0*mm,         0*mm,   -19.15*mm,  0,"ONLY");//PPF 2001P2 
123     gMC->Gspos("RppfLarge",1,"Rppf",  -224.5*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
124     gMC->Gspos("RppfLarge",2,"Rppf",  -224.5*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
125     gMC->Gspos("RppfLarge",3,"Rppf",  -224.5*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
126     gMC->Gspos("RppfLarge",4,"Rppf",  -224.5*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
127     gMC->Gspos("RppfSmall",1,"Rppf",  - 65.0*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
128     gMC->Gspos("RppfSmall",2,"Rppf",  - 65.0*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
129     gMC->Gspos("RppfSmall",3,"Rppf",  - 65.0*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
130     gMC->Gspos("RppfSmall",4,"Rppf",  - 65.0*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
131     gMC->Gspos("RppfSmall",5,"Rppf",  + 65.0*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
132     gMC->Gspos("RppfSmall",6,"Rppf",  + 65.0*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
133     gMC->Gspos("RppfSmall",7,"Rppf",  + 65.0*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
134     gMC->Gspos("RppfSmall",8,"Rppf",  + 65.0*mm,  +151.875*mm,     0.85*mm,  0,"ONLY"); 
135     gMC->Gspos("RppfLarge",5,"Rppf",  +224.5*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
136     gMC->Gspos("RppfLarge",6,"Rppf",  +224.5*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
137     gMC->Gspos("RppfLarge",7,"Rppf",  +224.5*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
138     gMC->Gspos("RppfLarge",8,"Rppf",  +224.5*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
139 //Gap - anod wires 6 copies to HMPID
140   par[0]=648*mm/2;par[1]=  411*mm/2 ;par[2]=4.45*mm/2;gMC->Gsvolu("Rgap","BOX ",(*fIdtmed)[kCH4] ,par,3);//xy as PPF 2001P2 z WP 2099P1
141   par[0]=  0*mm  ;par[1]=  20*mkm/2 ;par[2]= 648*mm/2;gMC->Gsvolu("Rano","TUBE",(*fIdtmed)[kW]   ,par,3);//WP 2099P1 z = gap x PPF 2001P2
142   AliMatrix(matrixIdReturn,180,0, 90,90, 90,0); //wires along x
143   
144   gMC->Gspos("Rgap",0,"HMPID",    -335*mm,      -433*mm,8*cm-2.225*mm, 0,"ONLY"); //F1 2040P1 z WP 2099P1
145   gMC->Gspos("Rgap",1,"HMPID",    +335*mm,      -433*mm,8*cm-2.225*mm, 0,"ONLY"); 
146   gMC->Gspos("Rgap",2,"HMPID",    -335*mm,         0*mm,8*cm-2.225*mm, 0,"ONLY"); 
147   gMC->Gspos("Rgap",3,"HMPID",    +335*mm,         0*mm,8*cm-2.225*mm, 0,"ONLY"); 
148   gMC->Gspos("Rgap",4,"HMPID",    -335*mm,      +433*mm,8*cm-2.225*mm, 0,"ONLY"); 
149   gMC->Gspos("Rgap",5,"HMPID",    +335*mm,      +433*mm,8*cm-2.225*mm, 0,"ONLY"); 
150   for(int i=1;i<=96;i++)
151     gMC->Gspos("Rano",i,"Rgap",     0*mm, -411/2*mm+i*4*mm, 0.185*mm, matrixIdReturn,"ONLY"); //WP 2099P1  
152 //Defines radiators geometry  
153   par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=  24*mm/2;  gMC->Gsvolu("Rrad"      ,"BOX ",(*fIdtmed)[kC6F14]     ,par,3); // Rad 2011P1
154   par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=   4*mm/2;  gMC->Gsvolu("RradFront" ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //front 
155   par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=   5*mm/2;  gMC->Gsvolu("RradWin"   ,"BOX ",(*fIdtmed)[kSiO2]      ,par,3); //window
156   par[0]=1330*mm/2 ;par[1]=   5*mm/2  ;par[2]=  15*mm/2;  gMC->Gsvolu("RradLong"  ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //long side  
157   par[0]=  10*mm/2 ;par[1]= 403*mm/2  ;par[2]=  15*mm/2;  gMC->Gsvolu("RradShort" ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //short side 
158   par[0]=   0      ;par[1]=  10*mm/2  ;par[2]=  15*mm/2;  gMC->Gsvolu("RradSpacer","TUBE",(*fIdtmed)[kSiO2]      ,par,3); //spacer        
159     
160   gMC->Gspos("Rrad",1,"HMPID",   0*mm,-434*mm,   -12*mm,  0,"ONLY"); //3 radiators to HMPID
161   gMC->Gspos("Rrad",2,"HMPID",   0*mm,   0*mm,   -12*mm,  0,"ONLY"); 
162   gMC->Gspos("Rrad",3,"HMPID",   0*mm,+434*mm,   -12*mm,  0,"ONLY"); 
163     gMC->Gspos("RradFront",1,"Rrad",   0*mm,   0*mm, -10.0*mm,  0,"ONLY"); //front cover 
164     gMC->Gspos("RradWin"  ,1,"Rrad",   0*mm,   0*mm,   9.5*mm,  0,"ONLY"); //quartz window (back cover)
165     gMC->Gspos("RradLong" ,1,"Rrad",   0*mm,-204*mm,  -0.5*mm,  0,"ONLY"); //long side
166     gMC->Gspos("RradLong" ,2,"Rrad",   0*mm,+204*mm,  -0.5*mm,  0,"ONLY"); //long side
167     gMC->Gspos("RradShort",1,"Rrad",-660*mm,   0*mm,  -0.5*mm,  0,"ONLY"); //short side
168     gMC->Gspos("RradShort",2,"Rrad",+660*mm,   0*mm,  -0.5*mm,  0,"ONLY"); //short side 
169     for(int i=0;i<3;i++)
170       for(int j=0;j<10;j++)
171         gMC->Gspos("RradSpacer",10*i+j,"Rrad",-1330*mm/2+116*mm+j*122*mm,(i-1)*105*mm,-0.5*mm,0,"ONLY");//spacers
172 //Defines SandBox geometry
173   par[0]=1419*mm/2 ;par[1]=1378*mm/2;par[2]=50.5*mm/2; gMC->Gsvolu("Rsb"     ,"BOX ",(*fIdtmed)[kAir]  ,par,3);  //2072P1   
174   par[0]=1419*mm/2 ;par[1]=1378*mm/2;par[2]= 0.5*mm/2; gMC->Gsvolu("RsbCover","BOX ",(*fIdtmed)[kAl]   ,par,3);  //cover
175   par[0]=1359*mm/2 ;par[1]=1318*mm/2;par[2]=49.5*mm/2; gMC->Gsvolu("RsbComb" ,"BOX ",(*fIdtmed)[kRoha] ,par,3);  //honeycomb structure 
176   
177   gMC->Gspos("Rsb",1,"HMPID",   0*mm, 0*mm, -73.75*mm, 0,"ONLY"); //p.84 TDR sandbox to rich
178     gMC->Gspos("RsbComb" ,1,"Rsb", 0*mm, 0*mm,      0*mm, 0,"ONLY"); //2072P1 honeycomv to sandbox
179     gMC->Gspos("RsbCover",1,"Rsb", 0*mm, 0*mm,    +25*mm, 0,"ONLY"); //cover to sandbox
180     gMC->Gspos("RsbCover",2,"Rsb", 0*mm, 0*mm,    -25*mm, 0,"ONLY"); //cover to sandbox
181   AliDebug(1,"Stop v1. HMPID option");  
182 }//CreateGeometry()
183 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
184 void AliHMPIDv1::Init()
185 {
186 // This methode defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(gMC->CurrentVolID()==XXX) statements in StepManager()
187 // Arguments: none
188 //   Returns: none      
189   AliDebug(1,"Start v1 HMPID.");    
190   fIdRad     = gMC->VolId("Rrad");
191   fIdWin     = gMC->VolId("RradWin");
192   fIdPc      = gMC->VolId("Rpc");
193   fIdAmpGap  = gMC->VolId("Rgap");
194   fIdProxGap = gMC->VolId("Rgap");
195
196   AliDebug(1,"Stop v1 HMPID.");    
197 }
198 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
199 void AliHMPIDv1::DefineOpticalProperties()
200 {
201 // Optical properties definition.
202   const Int_t kNbins=30;       //number of photon energy points
203   Float_t emin=5.5,emax=8.5;         //Photon energy range,[eV]
204   Float_t aEckov [kNbins]; 
205   Double_t dEckov [kNbins]; 
206   Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
207   Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins]; 
208   Float_t                                                    aQeAll [kNbins], aQePc [kNbins];
209   Double_t dReflMet[kNbins], dQePc[kNbins];
210
211   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
212   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
213   TF1 *pGaIF=new TF1("HidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)"                                                         ,emin,emax);      //?????? from where  
214
215   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 
216   pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
217   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 
218   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  
219   
220   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  
221                             
222   for(Int_t i=0;i<kNbins;i++){
223     Float_t eV=emin+0.1*i;  //Ckov energy in eV
224     aEckov [i] =1e-9*eV;    //Ckov energy in GeV
225     dEckov [i] = aEckov[i];
226     aAbsRad[i]=pRaAF->Eval(eV); aIdxRad[i]=1.292;//pRaIF->Eval(eV,20);      //Simulation for 20 degress C       
227     aAbsWin[i]=pWiAF->Eval(eV); aIdxWin[i]=1.5787;//pWiIF->Eval(eV);
228     aAbsGap[i]=pGaAF->Eval(eV); aIdxGap[i]=1.0005;//pGaIF->Eval(eV);   
229     aQeAll[i] =1;                     //QE for all other materials except for PC must be 1.  
230     aAbsMet[i] =0.0001;                aIdxMet[i]=0;                                             //metal ref idx must be 0 in order to reflect photon
231                                        aIdxPc [i]=1;           aQePc [i]=pQeF->Eval(eV);         //PC ref idx must be 1 in order to apply photon to QE conversion 
232     dQePc [i]=pQeF->Eval(eV);
233     dReflMet[i] = 0.;     // no reflection on the surface of the pc (?)                                       
234   }
235   gMC->SetCerenkov((*fIdtmed)[kC6F14]    , kNbins, aEckov, aAbsRad  , aQeAll , aIdxRad );    
236   gMC->SetCerenkov((*fIdtmed)[kSiO2]     , kNbins, aEckov, aAbsWin  , aQeAll , aIdxWin );    
237   gMC->SetCerenkov((*fIdtmed)[kCH4]      , kNbins, aEckov, aAbsGap  , aQeAll , aIdxGap );    
238   gMC->SetCerenkov((*fIdtmed)[kCu]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
239   gMC->SetCerenkov((*fIdtmed)[kW]        , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet ); //n=0 means reflect photons       
240   gMC->SetCerenkov((*fIdtmed)[kCsI]      , kNbins, aEckov, aAbsMet  , aQePc  , aIdxPc  ); //n=1 means convert photons    
241   gMC->SetCerenkov((*fIdtmed)[kAl]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
242
243   // Define a skin surface for the photocatode to enable 'detection' in G4
244   gMC->DefineOpSurface("surfPc", kGlisur /*kUnified*/,kDielectric_metal,kPolished, 0.);
245   gMC->SetMaterialProperty("surfPc", "EFFICIENCY", kNbins, dEckov, dQePc);
246   gMC->SetMaterialProperty("surfPc", "REFLECTIVITY", kNbins, dEckov, dReflMet);
247   gMC->SetSkinSurface("skinPc", "Rpc", "surfPc");
248
249   delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
250 }
251 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
252 Bool_t AliHMPIDv1::IsLostByFresnel()
253 {
254 // Calculate probability for the photon to be lost by Fresnel reflection.
255   TLorentzVector p4;
256   Double_t mom[3],localMom[3];
257   gMC->TrackMomentum(p4);   mom[0]=p4(1);   mom[1]=p4(2);   mom[2]=p4(3);
258   localMom[0]=0; localMom[1]=0; localMom[2]=0;
259   gMC->Gmtod(mom,localMom,2);
260   Double_t localTc    = localMom[0]*localMom[0]+localMom[2]*localMom[2];
261   Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
262   Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
263   if(gMC->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
264     AliDebug(1,"Photon lost");
265     return kTRUE;
266   }else
267     return kFALSE;
268 }//IsLostByFresnel()
269 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
270 void AliHMPIDv1::GenFee(Float_t qtot)
271 {
272 // Generate FeedBack photons for the current particle. To be invoked from StepManager().
273 // eloss=0 means photon so only pulse height distribution is to be analysed.
274   TLorentzVector x4;
275   gMC->TrackPosition(x4); 
276   Int_t iNphotons=gMC->GetRandom()->Poisson(0.02*qtot);  //# of feedback photons is proportional to the charge of hit
277   AliDebug(1,Form("N photons=%i",iNphotons));
278   Int_t j;
279   Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
280 //Generate photons
281   for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
282     Double_t ranf[2];
283     gMC->GetRandom()->RndmArray(2,ranf);    //Sample direction
284     cthf=ranf[0]*2-1.0;
285     if(cthf<0) continue;
286     sthf = TMath::Sqrt((1. - cthf) * (1. + cthf));
287     phif = ranf[1] * 2 * TMath::Pi();
288     
289     if(Double_t randomNumber=gMC->GetRandom()->Rndm()<=0.57)
290       enfp = 7.5e-9;
291     else if(randomNumber<=0.7)
292       enfp = 6.4e-9;
293     else
294       enfp = 7.9e-9;
295     
296
297     dir[0] = sthf * TMath::Sin(phif);    dir[1] = cthf;    dir[2] = sthf * TMath::Cos(phif);
298     gMC->Gdtom(dir, mom, 2);
299     mom[0]*=enfp;    mom[1]*=enfp;    mom[2]*=enfp;
300     mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
301     
302     // Polarisation
303     e1[0]=      0;    e1[1]=-dir[2];    e1[2]= dir[1];
304     e2[0]=-dir[1];    e2[1]= dir[0];    e2[2]=      0;
305     e3[0]= dir[1];    e3[1]=      0;    e3[2]=-dir[0];
306     
307     vmod=0;
308     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
309     if (!vmod) for(j=0;j<3;j++) {
310       uswop=e1[j];
311       e1[j]=e3[j];
312       e3[j]=uswop;
313     }
314     vmod=0;
315     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
316     if (!vmod) for(j=0;j<3;j++) {
317       uswop=e2[j];
318       e2[j]=e3[j];
319       e3[j]=uswop;
320     }
321     
322     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;    
323     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;
324     
325     phi = gMC->GetRandom()->Rndm()* 2 * TMath::Pi();
326     for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
327     gMC->Gdtom(pol, pol, 2);
328     Int_t outputNtracksStored;    
329     gAlice->GetMCApp()->PushTrack(1,                             //transport
330                      gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track 
331                      50000051,                                   //PID
332                      mom[0],mom[1],mom[2],mom[3],                //track momentum  
333                      x4.X(),x4.Y(),x4.Z(),x4.T(),                //track origin 
334                      pol[0],pol[1],pol[2],                       //polarization
335                      kPFeedBackPhoton,                           //process ID   
336                      outputNtracksStored,                        //on return how many new photons stored on stack
337                      1.0);                                       //weight
338   }//feedbacks loop
339   AliDebug(1,"Stop.");
340 }//GenerateFeedbacks()
341 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
342 void AliHMPIDv1::Hits2SDigits()
343 {
344 // Interface method ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
345 // Arguments: none
346 //   Returns: none   
347   AliDebug(1,"Start.");
348   for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){                //events loop
349     GetLoader()->GetRunLoader()->GetEvent(iEvt);                          //get next event
350   
351     if(!GetLoader()->TreeH()) {GetLoader()->LoadHits();                    }
352     if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
353           
354     for(Int_t iEnt=0;iEnt<GetLoader()->TreeH()->GetEntries();iEnt++){//prims loop
355       GetLoader()->TreeH()->GetEntry(iEnt);
356       Hit2Sdi(Hits(),SdiLst());
357     }//prims loop
358     GetLoader()->TreeS()->Fill();
359     GetLoader()->WriteSDigits("OVERWRITE");
360     SdiReset();
361   }//events loop  
362   GetLoader()->UnloadHits();
363   GetLoader()->UnloadSDigits();  
364   AliDebug(1,"Stop.");
365 }//Hits2SDigits()
366 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
367 void AliHMPIDv1::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
368 {
369 // Converts list of hits to list of sdigits. 
370 // Arguments: pHitLst  - list of hits provided not empty
371 //            pSDigLst - list of sdigits where to store the results
372 //   Returns: none         
373   for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){         //hits loop
374     AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);         //get pointer to current hit   
375     pHit->Hit2Sdi(pSdiLst);                                    //convert this hit to list of sdigits     
376   }//hits loop loop
377 }//Hits2Sdi()
378 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
379 void AliHMPIDv1::Digits2Raw()
380 {
381 // Interface method invoked by AliSimulation to create raw data streams from digits. Events loop is done in AliSimulation
382 // Arguments: none
383 //   Returns: none    
384   AliDebug(1,"Start.");
385   GetLoader()->LoadDigits();
386   TTree * treeD = GetLoader()->TreeD();
387   if(!treeD) {
388     AliError("No digits tree!");
389     return;
390   }
391   treeD->GetEntry(0);
392   
393   //AliHMPIDDigit::WriteRaw(DigLst());
394    AliHMPIDRawStream *pRS=0x0;
395    pRS->WriteRaw(DigLst());
396     
397   GetLoader()->UnloadDigits();
398   AliDebug(1,"Stop.");      
399 }//Digits2Raw()
400 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
401 Float_t AliHMPIDv1::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
402 {
403 // Correction for Fresnel   ???????????
404 // Arguments:   ene - photon energy [GeV],
405 //              PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
406 //   Returns:  
407     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,
408                       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,
409                       7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
410     Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
411                         2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
412                         2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
413                         1.72,1.53};
414     Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
415                         0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
416                         0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
417                         1.714,1.498};
418     Float_t xe=ene;
419     Int_t  j=Int_t(xe*10)-49;
420     Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
421     Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
422
423     //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
424     //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
425
426     Float_t sinin=TMath::Sqrt((1.-pdoti)*(1.+pdoti));
427     Float_t tanin=sinin/pdoti;
428
429     Float_t c1=cn*cn-ck*ck-sinin*sinin;
430     Float_t c2=4*cn*cn*ck*ck;
431     Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
432     Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
433     
434     Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
435     Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
436     
437
438     //CORRECTION FACTOR FOR SURFACE ROUGHNESS
439     //B.J. STAGG  APPLIED OPTICS, 30(1991),4113
440
441     Float_t sigraf=18.;
442     Float_t lamb=1240/ene;
443     Float_t fresn;
444  
445     Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
446
447     if(pola)
448     {
449         Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
450         fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
451     }
452     else
453         fresn=0.5*(rp+rs);
454       
455     fresn = fresn*rO;
456     return fresn;
457 }//Fresnel()
458 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
459 void AliHMPIDv1::Print(Option_t *option)const
460 {
461 // Debug printout
462   TObject::Print(option);
463 }//void AliHMPID::Print(Option_t *option)const
464 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
465 Bool_t AliHMPIDv1::Raw2SDigits(AliRawReader *pRR)
466 {
467 // Interface methode ivoked from AliSimulation to create a list of sdigits from raw digits. Events loop is done in AliSimulation
468 // Arguments: pRR- raw reader 
469 //   Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())      
470   //AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
471   
472   if(!GetLoader()->TreeS()) {MakeTree("S");  MakeBranch("S");}
473     
474   TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
475   AliHMPIDRawStream stream(pRR);
476   while(stream.Next())
477   {
478     for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
479       AliHMPIDDigit sdi(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
480       new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
481     }
482   }
483   
484   GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
485   SdiReset();
486   return kTRUE;
487 }//Raw2SDigits
488 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
489 void AliHMPIDv1::StepCount()
490 {
491 // Count number of ckovs created  
492 }
493 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
494 void AliHMPIDv1::StepHistory()
495 {
496 // This methode is invoked from StepManager() in order to print out 
497   static Int_t iStepN;
498   const char *sParticle;
499   switch(gMC->TrackPid()){
500     case kProton:      sParticle="PROTON"    ;break;
501     case kNeutron:     sParticle="neutron"   ;break;
502     case kGamma:       sParticle="gamma"     ;break;
503     case 50000050:     sParticle="CKOV"      ;break;
504     case kPi0:         sParticle="Pi0"       ;break;  
505     case kPiPlus:      sParticle="Pi+"       ;break;  
506     case kPiMinus:     sParticle="Pi-"       ;break;  
507     case kElectron:    sParticle="electron"  ;break;  
508     default:           sParticle="not known" ;break;
509   }
510
511   TString flag="fanny combination";
512   if(gMC->IsTrackAlive()) {
513     if(gMC->IsTrackEntering())      flag="enters to";
514     else if(gMC->IsTrackExiting())  flag="exits from";
515     else if(gMC->IsTrackInside())   flag="inside";
516   } else {
517     if(gMC->IsTrackStop())          flag="stopped in";        
518   }
519   
520   Int_t vid=0,copy=0;
521   TString path=gMC->CurrentVolName(); path.Prepend("-");path.Prepend(gMC->CurrentVolOffName(1));//current volume and his mother are always there
522   vid=gMC->CurrentVolOffID(2,copy);  if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
523   vid=gMC->CurrentVolOffID(3,copy);  if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
524   
525   Printf("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f",iStepN,sParticle,gMC->TrackPid(),flag.Data(),path.Data(),gMC->TrackMass(),gMC->TrackCharge(),gMC->Edep()*1e9);
526   
527   Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
528                             iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
529                             gMC->IsTrackAlive(), gMC->IsTrackDisappeared(),gMC->IsTrackEntering(), gMC->IsTrackExiting(),
530                             gMC->IsTrackInside(),gMC->IsTrackOut(),        gMC->IsTrackStop(),     gMC->IsNewTrack());
531   
532   Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
533   Int_t mid=gMC->CurrentMaterial(a,z,den,rad,abs);
534   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);
535   iStepN++;
536 }//StepHistory()
537 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
538 void AliHMPIDv1::StepManager()
539 {
540 // Full Step Manager.
541 // Arguments: none
542 //   Returns: none           
543 //  StepHistory(); return; //uncomment to print tracks history
544 //  StepCount(); return;     //uncomment to count photons
545   
546   Int_t   copy; //volume copy aka node
547   
548 //Treat photons    
549   if((gMC->TrackPid()==50000050||gMC->TrackPid()==50000051)&&gMC->CurrentVolID(copy)==fIdPc){  //photon (Ckov or feedback) hit PC (fIdPc)
550     if(gMC->Edep()>0){                                                                           //photon survided QE test i.e. produces electron
551       if(IsLostByFresnel()){ gMC->StopTrack(); return;}                                          //photon lost due to fersnel reflection on PC       
552                        gMC->CurrentVolOffID(2,copy);                                             //current chamber since geomtry tree is HMPID-Rppf-Rpc
553       Int_t   tid=     gMC->GetStack()->GetCurrentTrackNumber();                                 //take TID
554       Int_t   pid=     gMC->TrackPid();                                                          //take PID
555       Float_t etot=    gMC->Etot();                                                              //total hpoton energy, [GeV] 
556       Double_t x[3];   gMC->TrackPosition(x[0],x[1],x[2]);                                       //take MARS position at entrance to PC
557       Float_t  hitTime=(Float_t)gMC->TrackTime();                                                //hit formation time
558       Float_t xl,yl;   AliHMPIDParam::Instance()->Mars2Lors(copy,x,xl,yl);                       //take LORS position
559       new((*fHits)[fNhits++])AliHMPIDHit(copy,etot,pid,tid,xl,yl,hitTime,x);                     //HIT for photon, position at P, etot will be set to Q
560       GenFee(etot);                                                                              //generate feedback photons etot is modified in hit ctor to Q of hit
561     }//photon hit PC and DE >0 
562   }//photon hit PC
563   
564 //Treat charged particles  
565   static Float_t eloss;                                                                           //need to store mip parameters between different steps    
566   static Double_t in[3];
567   if(gMC->TrackCharge() && gMC->CurrentVolID(copy)==fIdAmpGap){                                   //charged particle in amplification gap (fIdAmpGap)
568     if(gMC->IsTrackEntering()||gMC->IsNewTrack()) {                                               //entering or newly created
569       eloss=0;                                                                                    //reset Eloss collector                         
570       gMC->TrackPosition(in[0],in[1],in[2]);                                                      //take position at the entrance
571     }else if(gMC->IsTrackExiting()||gMC->IsTrackStop()||gMC->IsTrackDisappeared()){               //exiting or disappeared
572       eloss              +=gMC->Edep();                                                           //take into account last step Eloss
573                           gMC->CurrentVolOffID(1,copy);                                           //take current chamber since geometry tree is HMPID-Rgap
574       Int_t tid=          gMC->GetStack()->GetCurrentTrackNumber();                               //take TID
575       Int_t pid=          gMC->TrackPid();                                                        //take PID
576       Double_t out[3];    gMC->TrackPosition(out[0],out[1],out[2]);                               //take MARS position at exit
577       Float_t hitTime=    (Float_t)gMC->TrackTime();                                              //hit formation time
578       out[0]=0.5*(out[0]+in[0]);                                                                  //>
579       out[1]=0.5*(out[1]+in[1]);                                                                  //take hit position at the anod plane
580       out[2]=0.5*(out[2]+in[2]);                                                                  //>
581       Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(copy,out,xl,yl);                         //take LORS position
582       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 
583       GenFee(eloss);                                                                              //generate feedback photons 
584     }else                                                                                         //just going inside
585       eloss          += gMC->Edep();                                                              //collect this step eloss 
586   }//MIP in GAP
587 }//StepManager()
588 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++