]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDv2.cxx
6034c4dc362d9e683e98721ef9b536680dc0aeef
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDv2.cxx
1
2 // **************************************************************************
3 // * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 // *                                                                        *
5 // * Author: The ALICE Off-line Project.                                    *
6 // * Contributors are mentioned in the code where appropriate.              *
7 // *                                                                        *
8 // * Permission to use, copy, modify and distribute this software and its   *
9 // * documentation strictly for non-commercial purposes is hereby granted   *
10 // * without fee, provided that the above copyright notice appears in all   *
11 // * copies and that both the copyright notice and this permission notice   *
12 // * appear in the supporting documentation. The authors make no claims     *
13 // * about the suitability of this software for any purpose. It is          *
14 // * provided "as is" without express or implied warranty.                  *
15 // **************************************************************************
16
17
18 #include "AliHMPIDv2.h"     //class header
19 #include "AliHMPIDParam.h"  //StepManager()
20 #include "AliHMPIDHit.h"    //Hits2SDigs(),StepManager()
21 #include "AliHMPIDDigit.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 <AliConst.h>
29 #include <AliPDG.h>
30 #include <AliMC.h>            //StepManager()      
31 #include <AliRun.h>           //CreateMaterials()    
32 #include <AliMagF.h>          //CreateMaterials()
33 #include <TGeoManager.h>      //CreateGeometry()
34 #include <TF1.h>              //DefineOpticalProperties()
35 #include <TF2.h>              //DefineOpticalProperties()
36 #include <TLorentzVector.h>   //IsLostByFresnel() 
37 #include <AliCDBManager.h>    //CreateMaterials()
38 #include <AliCDBEntry.h>      //CreateMaterials()
39  
40 ClassImp(AliHMPIDv2)    
41 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
42 void AliHMPIDv2::AddAlignableVolumes()const
43 {
44 // Associates the symbolic volume name with the corresponding volume path. Interface method from AliModule invoked from AliMC
45 // Arguments: none
46 //   Returns: none   
47   for(Int_t i=AliHMPIDDigit::kMinCh;i<=AliHMPIDDigit::kMaxCh;i++)
48     gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",i),Form("ALIC_1/Hmp_%i",i));           //clm ??? 
49 }
50 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
51 void AliHMPIDv2::CreateMaterials()
52 {
53 // Definition of available HMPID materials  
54 // Arguments: none
55 //   Returns: none    
56   AliDebug(1,"Start v2 HMPID.");
57     
58     //clm update material definition later on from Antonello
59     
60 //data from PDG booklet 2002     density [gr/cm^3] rad len [cm] abs len [cm]    
61   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
62   Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9}   , wC6F14[2]={6 , 14} , dC6F14=1.68    ; Int_t nC6F14=-2;
63   Float_t  aSiO2[2]={ 28.09 , 15.99} ,  zSiO2[2]={14 , 8}   ,  wSiO2[2]={1 ,  2} ,  dSiO2=2.64    ; Int_t  nSiO2=-2; 
64   Float_t   aCH4[2]={ 12.01 ,  1.01} ,   zCH4[2]={ 6 , 1}   ,   wCH4[2]={1 ,  4} ,   dCH4=7.17e-4 ; Int_t   nCH4=-2; 
65   Float_t   aCsI[2]={132.90 ,126.90} ,   zCsI[2]={55 ,53}   ,   wCsI[2]={1 ,  1} ,   dCsI=0.1     ; Int_t   nCsI=-2; 
66   
67   Float_t     aRoha= 12.01 ,               zRoha=  6 ,                              dRoha=  0.10 ,   radRoha= 18.80 , absRoha=  86.3/dRoha; //special material- quasi quartz
68   Float_t       aCu= 63.55 ,                 zCu= 29 ,                                dCu=  8.96 ,     radCu=  1.43 ,   absCu= 134.9/dCu  ;
69   Float_t        aW=183.84 ,                  zW= 74 ,                                 dW= 19.30 ,      radW=  0.35 ,    absW= 185.0/dW   ;
70   Float_t       aAl= 26.98 ,                 zAl= 13 ,                                dAl=  2.70 ,     radAl=  8.90 ,   absAl= 106.4/dAl  ;
71   Float_t       aAr= 39.94 ,                 zAr= 18 ,                                dAr=  1.396e-3,  radAr=  14.0 ,   absAr= 117.2/dAr  ;   
72            
73     Int_t   matId=0;                           //tmp material id number
74     Int_t   unsens =  0, sens=1;               //sensitive or unsensitive medium
75     Int_t   itgfld = gAlice->Field()->Integ(); //type of field intergration 0 no field -1 user in guswim 1 Runge Kutta 2 helix 3 const field along z
76     Float_t maxfld = gAlice->Field()->Max();   //max field value
77     Float_t tmaxfd = -10.0;                    //max deflection angle due to magnetic field in one step
78     Float_t deemax = - 0.2;                    //max fractional energy loss in one step   
79     Float_t stemax = - 0.1;                    //mas step allowed [cm]
80     Float_t epsil  =   0.001;                  //abs tracking precision [cm]   
81     Float_t stmin  = - 0.001;                  //min step size [cm] in continius process transport, negative value: choose it automatically
82     
83     AliMixture(++matId,"Air"  ,aAir  ,zAir  ,dAir  ,nAir  ,wAir  ); AliMedium(kAir  ,"Air"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
84     AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);      
85     AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);    
86     AliMixture(++matId,"CH4"  ,aCH4  ,zCH4  ,dCH4  ,nCH4  ,wCH4  ); AliMedium(kCH4  ,"CH4"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);  
87     AliMixture(++matId,"CsI"  ,aCsI  ,zCsI  ,dCsI  ,nCsI  ,wCsI  ); AliMedium(kCsI  ,"CsI"  ,matId,   sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
88   
89     AliMixture(++matId ,"Neo" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kNeo,"Neo" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin); //clm neoceram
90     AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha);  AliMedium(kRoha,"Roha", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin); //Roha->honeycomb
91
92     
93     AliMaterial(++matId,"Cu"  ,aCu  ,zCu  ,dCu  ,radCu  ,absCu  );  AliMedium(kCu  ,"Cu"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
94     AliMaterial(++matId,"W"   ,aW   ,zW   ,dW   ,radW   ,absW   );  AliMedium(kW   ,"W"   , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
95     AliMaterial(++matId,"Al"  ,aAl  ,zAl  ,dAl  ,radAl  ,absAl  );  AliMedium(kAl  ,"Al"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
96     AliMaterial(++matId,"Ar"  ,aAr  ,zAr  ,dAr  ,radAr  ,absAr  );  AliMedium(kAr  ,"Ar"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
97     
98     DefineOpticalProperties();
99 }//void AliHMPID::CreateMaterials()
100 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
101 void AliHMPIDv2::CreateGeometry()
102 {
103 //Creates detailed geometry simulation (currently GEANT volumes tree)         
104   AliDebug(1,"Start main.");
105   if(!gMC->IsRootGeometrySupported()) return;                
106  
107  Double_t cm=1,mm=0.1*cm,um=0.001*mm;//default is cm
108  
109   TGeoMedium *al   =gGeoManager->GetMedium("HMPID_Al");    
110   TGeoMedium *ch4  =gGeoManager->GetMedium("HMPID_CH4");    
111   TGeoMedium *roha =gGeoManager->GetMedium("HMPID_Roha");   
112   TGeoMedium *neoc =gGeoManager->GetMedium("HMPID_Neo");
113   TGeoMedium *c6f14=gGeoManager->GetMedium("HMPID_C6F14");  
114   TGeoMedium *sio2 =gGeoManager->GetMedium("HMPID_SiO2");   
115   TGeoMedium *cu   =gGeoManager->GetMedium("HMPID_Cu");     
116   TGeoMedium *w    =gGeoManager->GetMedium("HMPID_W");      
117   TGeoMedium *csi  =gGeoManager->GetMedium("HMPID_CsI");    
118   TGeoMedium *ar   =gGeoManager->GetMedium("HMPID_Ar");     
119
120   TGeoVolume *hmp=gGeoManager->MakeBox ("Hmp",ch4,1681*mm/2, 1466*mm/2,(2*80*mm+2*41*mm)/2);//2033P1  z from p84 TDR  
121
122   TString title=GetTitle();
123   if(title.Contains("TestBeam")  )
124       {
125         gGeoManager->GetVolume("ALIC")->AddNode(hmp,0,new TGeoTranslation(5.0*mm/2 ,  5.0*mm/2, 1000.0*mm));
126       }
127   else
128     {
129       for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++){//place 7 chambers
130       TGeoHMatrix *pMatrix=new TGeoHMatrix;
131       AliHMPIDParam::IdealPosition(iCh,pMatrix);
132       gGeoManager->GetVolume("ALIC")->AddNode(hmp,iCh,pMatrix);
133        }
134      }
135
136   TGeoRotation *rot=new TGeoRotation("HwireRot"); rot->RotateY(90); //rotate wires around Y to be along X (initially along Z)
137   TGeoVolume *sbo=gGeoManager->MakeBox ("Hsbo",ch4  , 1419*mm/2 , 1378.00*mm/2 ,   50.5*mm/2);//2072P1
138   TGeoVolume *cov=gGeoManager->MakeBox ("Hcov",al   , 1419*mm/2 , 1378.00*mm/2 ,    0.5*mm/2);  
139   TGeoVolume *hon=gGeoManager->MakeBox ("Hhon",roha , 1359*mm/2 , 1318.00*mm/2 ,   49.5*mm/2);  
140   TGeoVolume *rad=gGeoManager->MakeBox ("Hrad",c6f14, 1330*mm/2 ,  413.00*mm/2 ,   24.0*mm/2); //2011P1
141   TGeoVolume *neo=gGeoManager->MakeBox ("Hneo",neoc , 1330*mm/2 ,  413.00*mm/2 ,    4.0*mm/2); 
142   TGeoVolume *win=gGeoManager->MakeBox ("Hwin",sio2 , 1330*mm/2 ,  413.00*mm/2 ,    5.0*mm/2); 
143   TGeoVolume *si1=gGeoManager->MakeBox ("Hsi1",sio2 , 1330*mm/2 ,    5.00*mm/2 ,   15.0*mm/2);    
144   TGeoVolume *si2=gGeoManager->MakeBox ("Hsi2",neoc ,   10*mm/2 ,  403.00*mm/2 ,   15.0*mm/2);    
145   TGeoVolume *spa=gGeoManager->MakeTube("Hspa",sio2 ,    0*mm   ,    5.00*mm   ,   15.0*mm/2);         
146   TGeoVolume *fr4=gGeoManager->MakeBox ("Hfr4",ch4  , 1407*mm/2 , 1366.00*mm/2 ,   15.0*mm/2);//2043P1 
147   TGeoVolume *f4a=gGeoManager->MakeBox ("Hf4a",al   , 1407*mm/2 , 1366.00*mm/2 ,   10.0*mm/2); 
148   TGeoVolume *f4i=gGeoManager->MakeBox ("Hf4i",ch4  , 1323*mm/2 , 1296.00*mm/2 ,   10.0*mm/2); 
149   TGeoVolume *col=gGeoManager->MakeTube("Hcol",cu   ,    0*mm   ,  100.00*um   , 1323.0*mm/2);
150   TGeoVolume *sec=gGeoManager->MakeBox ("Hsec",ch4  ,  648*mm/2 ,  411.00*mm/2 ,   45.5*mm/2);//sec=gap+ppf
151   TGeoVolume *ppf=gGeoManager->MakeBox ("Hppf",al   ,  648*mm/2 ,  411.00*mm/2 ,   40.0*mm/2);//2001P2
152   TGeoVolume *lar=gGeoManager->MakeBox ("Hlar",ar   ,  181*mm/2 ,   89.25*mm/2 ,   38.3*mm/2);//2001P2
153   TGeoVolume *smo=gGeoManager->MakeBox ("Hsmo",ar   ,  114*mm/2 ,   89.25*mm/2 ,   38.3*mm/2);//2001P2
154   TGeoVolume *gap=gGeoManager->MakeBox ("Hgap",ch4  ,  640*mm/2 ,  403.20*mm/2 ,    5.5*mm/2);//gap=pad+ano+cat
155   TGeoVolume *cat=gGeoManager->MakeTube("Hcat",cu   ,    0*mm   ,   50.00*um   ,    8.0*mm/2); 
156   TGeoVolume *ano=gGeoManager->MakeTube("Hano",w    ,    0*mm   ,   20.00*um   ,    8.0*mm/2); 
157   TGeoVolume *pad=gGeoManager->MakeBox ("Hpad",csi  ,    8*mm/2 ,    8.40*mm/2 ,    1.0*mm/2);      
158 //
159 // ^ Y   z=         z=-12mm      z=98.25mm               ALIC->7xHmp (virtual)-->1xHsbo (virtual) --->2xHcov (real) 2072P1
160 // |  ____________________________________                                    |                   |-->1xHhon (real) 2072P1
161 // | |   ______     ____          ______  |                                   |
162 //   |  |      |   |    |   *    |      | |                                   |->3xHrad (virtual) --->1xHneo (real) 2011P1
163 //   |  |50.5mm|   |24mm|   *    |45.5mm| |                                   |                   |-->1xHwin (real) 2011P1
164 //   |  |      |   |    |   *    |      | |                                   |                   |-->2xHsi1 (real) 2011P1
165 //   |  |      |   |____|   *    |______| |                                   |                   |-->2xHsi2 (real) 2011P1
166 //   |  |      |    ____    *     ______  |                                   |                   |->30xHspa (real) 2011P1
167 //   |  |      |   |    |   *    |      | |                                   |
168 //   |  |      |   |    |   *    |      | |                                   |->1xHfr4 (vitual) --->1xHf4a (real)---->1xHf4i(virtual) 2043P1 
169 //   |  |  sb  |   | rad|   *    |      | |                                   |                  |-->322xHcol (real) 2043P1
170 //   |  |      |   |____|   *    |______| |                                   |
171 //   |  |      |    ____    *     ______  |                                   |->6xHsec (virtual) --> 1xHppf(real) ---->8xHlar (virtual) 2001P1
172 //   |  |      |   |    |   *    |      | |                                                                        |--->8xHsmo (virtual) 2001P1     
173 //   |  |      |   |    |   *    |      | |                                   |               
174 //   |  |      |   |    |   *    |      | |                                   |-> 1xHgap (virtual) --->48xHrow (virtual) -->80xHcel (virtual) -->4xHcat (real) from p84 TDR 
175 //   |  |______|   |____|   *    |______| |                                                                                                  |-->2xHano (real) from p84 TDR                                  
176 //   |____________________________________|                                                                                                  |-->1xHpad (real) from p84 TDR 
177 //                                                       --->Z 
178   hmp->AddNode(sbo ,1,new TGeoTranslation(   0*mm,   0*mm, -73.75*mm));                     //p.84 TDR
179      sbo->AddNode(hon ,1,new TGeoTranslation(  0*mm,0*mm,      0*mm)); //2072P1
180      sbo->AddNode(cov ,1,new TGeoTranslation(  0*mm,0*mm,    +25*mm)); 
181      sbo->AddNode(cov ,2,new TGeoTranslation(  0*mm,0*mm,    -25*mm)); 
182   hmp->AddNode(rad,2,new TGeoTranslation(   0*mm,+434*mm, -12.00*mm)); 
183   hmp->AddNode(rad,1,new TGeoTranslation(   0*mm,   0*mm, -12.00*mm)); 
184   hmp->AddNode(rad,0,new TGeoTranslation(   0*mm,-434*mm, -12.00*mm)); 
185     rad->AddNode(neo,1,new TGeoTranslation(   0*mm,   0*mm, -10.0*mm));
186     rad->AddNode(win,1,new TGeoTranslation(   0*mm,   0*mm,   9.5*mm));
187     rad->AddNode(si1,1,new TGeoTranslation(   0*mm,-204*mm,  -0.5*mm)); rad->AddNode(si1,2,new TGeoTranslation(   0*mm,+204*mm,  -0.5*mm));
188     rad->AddNode(si2,1,new TGeoTranslation(-660*mm,   0*mm,  -0.5*mm)); rad->AddNode(si2,2,new TGeoTranslation(+660*mm,   0*mm,  -0.5*mm));
189     for(Int_t i=0;i<3;i++) for(Int_t j=0;j<10;j++) rad->AddNode(spa,10*i+j,new TGeoTranslation(-1330*mm/2+116*mm+j*122*mm,(i-1)*105*mm,-0.5*mm));
190   hmp->AddNode(fr4,1,new TGeoTranslation(   0*mm,   0*mm,   9.00*mm));                     //p.84 TDR
191   for(int i=1;i<=322;i++)  fr4->AddNode(col,i,new TGeoCombiTrans( 0*mm, -1296/2*mm+i*4*mm,-5*mm,rot)); //F4 2043P1
192                            fr4->AddNode(f4a,1,new TGeoTranslation(   0*mm,0*mm, 2.5*mm));    
193                                         f4a->AddNode(f4i,1,new TGeoTranslation(   0*mm,0*mm,   0*mm));
194   hmp->AddNode(sec,4,new TGeoTranslation(-335*mm,+433*mm,  98.25*mm)); hmp->AddNode(sec,5,new TGeoTranslation(+335*mm,+433*mm,  98.25*mm));
195   hmp->AddNode(sec,2,new TGeoTranslation(-335*mm,   0*mm,  98.25*mm)); hmp->AddNode(sec,3,new TGeoTranslation(+335*mm,   0*mm,  98.25*mm));
196   hmp->AddNode(sec,0,new TGeoTranslation(-335*mm,-433*mm,  98.25*mm)); hmp->AddNode(sec,1,new TGeoTranslation(+335*mm,-433*mm,  98.25*mm));
197     sec->AddNode(gap,1,new TGeoTranslation(0,0,-20.00*mm));
198   TGeoVolume *row=          gap->Divide("Hrow",2,48,0,0);//along Y->48 rows
199   TGeoVolume *cel=          row->Divide("Hcel",1,80,0,0);//along X->80 cells
200       cel->AddNode(cat,1,new TGeoCombiTrans (0,  3.15*mm , -2.70*mm , rot)); //4 cathode wires
201       cel->AddNode(ano,1,new TGeoCombiTrans (0,  2.00*mm , -0.29*mm , rot)); //2 anod wires
202       cel->AddNode(cat,2,new TGeoCombiTrans (0,  1.05*mm , -2.70*mm , rot)); 
203       cel->AddNode(cat,3,new TGeoCombiTrans (0, -1.05*mm , -2.70*mm , rot)); 
204       cel->AddNode(ano,2,new TGeoCombiTrans (0, -2.00*mm , -0.29*mm , rot)); 
205       cel->AddNode(cat,4,new TGeoCombiTrans (0, -3.15*mm , -2.70*mm , rot));   
206       cel->AddNode(pad,1,new TGeoTranslation(0,  0.00*mm ,  2.25*mm));       //1 pad  
207     sec->AddNode(ppf,1,new TGeoTranslation(0,0,  2.75*mm));
208 // ^ Y  single cell                                                5.5mm CH4 = 1*mm CsI + 4.45*mm CsI x cath +0.05*mm safety margin         
209 // |      ______________________________           
210 // |     |                              |          ^                            ||     
211 //       |                              | 1.05mm                                ||     
212 // 2.2*mm| xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--              50um  x                || cat shift  x=0mm , y= 3.15mm , z=-2.70mm       
213 //       |                              |                                       ||     
214 //       |                              |                                       ||     
215 // __    |  ..........................  | 2.1mm                    20un .       ||  ano shift x=0mm , y= 2.00mm , z=-0.29mm   
216 //       |                              |                                       ||     
217 //       |                              |                                       ||     
218 //       | xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--                    x                ||  cat shift x=0mm , y= 1.05mm , z=-2.70mm   
219 //       |                              |                                       ||     
220 //       |                              |         8.4mm                         ||   
221 // 4*mm  |                              | 2.1mm                                 ||  pad shift x=0mm , y= 0.00mm , z=2.25*mm   
222 //       |                              |                                       ||  
223 //       |                              |                                       ||  
224 //       | xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--                    x                ||  cat shift x=0mm , y=-1.05mm , z=-2.70mm   
225 //       |                              |                                       ||  
226 //       |                              |                                       ||    
227 // __    |  ..........................  | 2.1mm                         . 2.04mm||  ano shift x=0mm , y=-2.00mm , z=-0.29mm   
228 //       |                              |                                       ||  
229 //       |                              |                                       ||  
230 //       | xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--                    x    4.45mm      ||  cat shift x=0mm , y=-3.15mm , z=-2.70mm   
231 // 2.2*mm|                              |                                       ||  
232 //       |                              | 1.05mm                                ||         
233 //       |______________________________|          v                            ||    
234 //       <             8 mm             >                          
235 //                                   ----->X                                 ----->Z
236   ppf->AddNode(lar,1,new TGeoTranslation(-224.5*mm,-151.875*mm,  0.85*mm));
237   ppf->AddNode(lar,2,new TGeoTranslation(-224.5*mm,- 50.625*mm,  0.85*mm));
238   ppf->AddNode(lar,3,new TGeoTranslation(-224.5*mm,+ 50.625*mm,  0.85*mm));
239   ppf->AddNode(lar,4,new TGeoTranslation(-224.5*mm,+151.875*mm,  0.85*mm));
240   ppf->AddNode(lar,5,new TGeoTranslation(+224.5*mm,-151.875*mm,  0.85*mm));
241   ppf->AddNode(lar,6,new TGeoTranslation(+224.5*mm,- 50.625*mm,  0.85*mm));
242   ppf->AddNode(lar,7,new TGeoTranslation(+224.5*mm,+ 50.625*mm,  0.85*mm));
243   ppf->AddNode(lar,8,new TGeoTranslation(+224.5*mm,+151.875*mm,  0.85*mm));
244   ppf->AddNode(smo,1,new TGeoTranslation(- 65.0*mm,-151.875*mm,  0.85*mm));
245   ppf->AddNode(smo,2,new TGeoTranslation(- 65.0*mm,- 50.625*mm,  0.85*mm));
246   ppf->AddNode(smo,3,new TGeoTranslation(- 65.0*mm,+ 50.625*mm,  0.85*mm));
247   ppf->AddNode(smo,4,new TGeoTranslation(- 65.0*mm,+151.875*mm,  0.85*mm));
248   ppf->AddNode(smo,5,new TGeoTranslation(+ 65.0*mm,-151.875*mm,  0.85*mm));
249   ppf->AddNode(smo,6,new TGeoTranslation(+ 65.0*mm,- 50.625*mm,  0.85*mm));
250   ppf->AddNode(smo,7,new TGeoTranslation(+ 65.0*mm,+ 50.625*mm,  0.85*mm));
251   ppf->AddNode(smo,8,new TGeoTranslation(+ 65.0*mm,+151.875*mm,  0.85*mm)); 
252
253
254
255   AliDebug(1,"Stop v2. HMPID option");  
256 }//CreateGeometry()
257 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
258 void AliHMPIDv2::Init()
259 {
260 // This methode defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(gMC->CurrentVolID()==XXX) statements in StepManager()
261 // Arguments: none
262 //   Returns: none      
263   AliDebug(1,"Start v2 HMPID.");    
264   fIdPad     = gMC->VolId("Hpad");
265   fIdCell    = gMC->VolId("Hcel");
266   AliDebug(1,"Stop v2 HMPID.");    
267 }
268 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
269 void AliHMPIDv2::DefineOpticalProperties()
270 {
271 // Optical properties definition.
272   const Int_t kNbins=30;       //number of photon energy points
273   Float_t emin=5.5,emax=8.5;         //Photon energy range,[eV]
274   Float_t aEckov [kNbins]; 
275   Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
276   Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins]; 
277   Float_t                                                    aQeAll [kNbins], aQePc [kNbins];
278
279   TF2 *pRaIF=new TF2("HidxRad","sqrt(1+0.554*(1239.84/x)^2/((1239.84/x)^2-5796)-0.0005*(y-20))"                                       ,emin,emax,0,50); //DiMauro mail temp 0-50 degrees C
280   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
281   TF1 *pGaIF=new TF1("HidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)"                                                         ,emin,emax);      //?????? from where  
282
283   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 
284   pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
285   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 
286   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  
287   
288   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  
289                             
290   for(Int_t i=0;i<kNbins;i++){
291     Float_t eV=emin+0.1*i;  //Ckov energy in eV
292     aEckov [i] =1e-9*eV;    //Ckov energy in GeV
293     aAbsRad[i]=pRaAF->Eval(eV); aIdxRad[i]=1.292;//pRaIF->Eval(eV,20);      //Simulation for 20 degress C       
294     aAbsWin[i]=pWiAF->Eval(eV); aIdxWin[i]=1.5787;//pWiIF->Eval(eV);
295     aAbsGap[i]=pGaAF->Eval(eV); aIdxGap[i]=1.0005;//pGaIF->Eval(eV);   
296     aQeAll[i] =1;                     //QE for all other materials except for PC must be 1.  
297     aAbsMet[i] =0.0001;                aIdxMet[i]=0;                                             //metal ref idx must be 0 in order to reflect photon
298                                        aIdxPc [i]=1;           aQePc [i]=pQeF->Eval(eV);         //PC ref idx must be 1 in order to apply photon to QE conversion 
299                                        
300   }
301   gMC->SetCerenkov((*fIdtmed)[kC6F14]    , kNbins, aEckov, aAbsRad  , aQeAll , aIdxRad );    
302   gMC->SetCerenkov((*fIdtmed)[kSiO2]     , kNbins, aEckov, aAbsWin  , aQeAll , aIdxWin );    
303   gMC->SetCerenkov((*fIdtmed)[kCH4]      , kNbins, aEckov, aAbsGap  , aQeAll , aIdxGap );    
304   gMC->SetCerenkov((*fIdtmed)[kCu]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
305   gMC->SetCerenkov((*fIdtmed)[kW]        , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet ); //n=0 means reflect photons       
306   gMC->SetCerenkov((*fIdtmed)[kCsI]      , kNbins, aEckov, aAbsMet  , aQePc  , aIdxPc  ); //n=1 means convert photons    
307   gMC->SetCerenkov((*fIdtmed)[kAl]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
308   delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
309 }
310 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
311 Bool_t AliHMPIDv2::IsLostByFresnel()
312 {
313 // Calculate probability for the photon to be lost by Fresnel reflection.
314   TLorentzVector p4;
315   Double_t mom[3],localMom[3];
316   gMC->TrackMomentum(p4);   mom[0]=p4(1);   mom[1]=p4(2);   mom[2]=p4(3);
317   localMom[0]=0; localMom[1]=0; localMom[2]=0;
318   gMC->Gmtod(mom,localMom,2);
319   Double_t localTc    = localMom[0]*localMom[0]+localMom[2]*localMom[2];
320   Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
321   Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
322   if(gMC->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
323     AliDebug(1,"Photon lost");
324     return kTRUE;
325   }else
326     return kFALSE;
327 }//IsLostByFresnel()
328 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
329 void AliHMPIDv2::GenFee(Float_t qtot)
330 {
331 // Generate FeedBack photons for the current particle. To be invoked from StepManager().
332 // eloss=0 means photon so only pulse height distribution is to be analysed.
333   TLorentzVector x4;
334   gMC->TrackPosition(x4); 
335   Int_t iNphotons=gMC->GetRandom()->Poisson(0.02*qtot);  //# of feedback photons is proportional to the charge of hit
336   AliDebug(1,Form("N photons=%i",iNphotons));
337   Int_t j;
338   Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
339 //Generate photons
340   for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
341     Double_t ranf[2];
342     gMC->GetRandom()->RndmArray(2,ranf);    //Sample direction
343     cthf=ranf[0]*2-1.0;
344     if(cthf<0) continue;
345     sthf = TMath::Sqrt((1 - cthf) * (1 + cthf));
346     phif = ranf[1] * 2 * TMath::Pi();
347     
348     if(Double_t randomNumber=gMC->GetRandom()->Rndm()<=0.57)
349       enfp = 7.5e-9;
350     else if(randomNumber<=0.7)
351       enfp = 6.4e-9;
352     else
353       enfp = 7.9e-9;
354     
355
356     dir[0] = sthf * TMath::Sin(phif);    dir[1] = cthf;    dir[2] = sthf * TMath::Cos(phif);
357     gMC->Gdtom(dir, mom, 2);
358     mom[0]*=enfp;    mom[1]*=enfp;    mom[2]*=enfp;
359     mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
360     
361     // Polarisation
362     e1[0]=      0;    e1[1]=-dir[2];    e1[2]= dir[1];
363     e2[0]=-dir[1];    e2[1]= dir[0];    e2[2]=      0;
364     e3[0]= dir[1];    e3[1]=      0;    e3[2]=-dir[0];
365     
366     vmod=0;
367     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
368     if (!vmod) for(j=0;j<3;j++) {
369       uswop=e1[j];
370       e1[j]=e3[j];
371       e3[j]=uswop;
372     }
373     vmod=0;
374     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
375     if (!vmod) for(j=0;j<3;j++) {
376       uswop=e2[j];
377       e2[j]=e3[j];
378       e3[j]=uswop;
379     }
380     
381     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;    
382     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;
383     
384     phi = gMC->GetRandom()->Rndm()* 2 * TMath::Pi();
385     for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
386     gMC->Gdtom(pol, pol, 2);
387     Int_t outputNtracksStored;    
388     gAlice->GetMCApp()->PushTrack(1,                             //transport
389                      gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track 
390                      50000051,                                   //PID
391                      mom[0],mom[1],mom[2],mom[3],                //track momentum  
392                      x4.X(),x4.Y(),x4.Z(),x4.T(),                //track origin 
393                      pol[0],pol[1],pol[2],                       //polarization
394                      kPFeedBackPhoton,                           //process ID   
395                      outputNtracksStored,                        //on return how many new photons stored on stack
396                      1.0);                                       //weight
397   }//feedbacks loop
398   AliDebug(1,"Stop.");
399 }//GenerateFeedbacks()
400 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
401 void AliHMPIDv2::Hits2SDigits()
402 {
403 // Interface method ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
404 // Arguments: none
405 //   Returns: none   
406   AliDebug(1,"Start.");
407   for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){                //events loop
408     GetLoader()->GetRunLoader()->GetEvent(iEvt);                          //get next event
409   
410     if(!GetLoader()->TreeH()) {GetLoader()->LoadHits();                    }
411     if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
412           
413     for(Int_t iEnt=0;iEnt<GetLoader()->TreeH()->GetEntries();iEnt++){//prims loop
414       GetLoader()->TreeH()->GetEntry(iEnt);
415       Hit2Sdi(Hits(),SdiLst());
416     }//prims loop
417     GetLoader()->TreeS()->Fill();
418     GetLoader()->WriteSDigits("OVERWRITE");
419     SdiReset();
420   }//events loop  
421   GetLoader()->UnloadHits();
422   GetLoader()->UnloadSDigits();  
423   AliDebug(1,"Stop.");
424 }//Hits2SDigits()
425 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
426 void AliHMPIDv2::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
427 {
428 // Converts list of hits to list of sdigits. 
429 // Arguments: pHitLst  - list of hits provided not empty
430 //            pSDigLst - list of sdigits where to store the results
431 //   Returns: none         
432   for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){         //hits loop
433     AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);         //get pointer to current hit   
434     pHit->Hit2Sdi(pSdiLst);                                    //convert this hit to list of sdigits     
435   }//hits loop loop
436 }//Hits2Sdi()
437 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
438 void AliHMPIDv2::Digits2Raw()
439 {
440 // Interface method invoked by AliSimulation to create raw data streams from digits. Events loop is done in AliSimulation
441 // Arguments: none
442 //   Returns: none    
443   AliDebug(1,"Start.");
444   GetLoader()->LoadDigits();
445   TTree * treeD = GetLoader()->TreeD();
446   if(!treeD) {
447     AliError("No digits tree!");
448     return;
449   }
450   treeD->GetEntry(0);
451   
452   AliHMPIDDigit::WriteRaw(DigLst());
453     
454   GetLoader()->UnloadDigits();
455   AliDebug(1,"Stop.");      
456 }//Digits2Raw()
457 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
458 Float_t AliHMPIDv2::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
459 {
460 // Correction for Fresnel   ???????????
461 // Arguments:   ene - photon energy [GeV],
462 //              PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
463 //   Returns:  
464     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,
465                       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,
466                       7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
467     Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
468                         2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
469                         2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
470                         1.72,1.53};
471     Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
472                         0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
473                         0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
474                         1.714,1.498};
475     Float_t xe=ene;
476     Int_t  j=Int_t(xe*10)-49;
477     Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
478     Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
479
480     //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
481     //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
482
483     Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
484     Float_t tanin=sinin/pdoti;
485
486     Float_t c1=cn*cn-ck*ck-sinin*sinin;
487     Float_t c2=4*cn*cn*ck*ck;
488     Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
489     Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
490     
491     Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
492     Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
493     
494
495     //CORRECTION FACTOR FOR SURFACE ROUGHNESS
496     //B.J. STAGG  APPLIED OPTICS, 30(1991),4113
497
498     Float_t sigraf=18.;
499     Float_t lamb=1240/ene;
500     Float_t fresn;
501  
502     Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
503
504     if(pola)
505     {
506         Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
507         fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
508     }
509     else
510         fresn=0.5*(rp+rs);
511       
512     fresn = fresn*rO;
513     return fresn;
514 }//Fresnel()
515 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
516 void AliHMPIDv2::Print(Option_t *option)const
517 {
518 // Debug printout
519   TObject::Print(option);
520 }//void AliHMPID::Print(Option_t *option)const
521 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
522 Bool_t AliHMPIDv2::Raw2SDigits(AliRawReader *pRR)
523 {
524 // Interface methode ivoked from AliSimulation to create a list of sdigits from raw digits. Events loop is done in AliSimulation
525 // Arguments: pRR- raw reader 
526 //   Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())      
527   AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
528   
529   if(!GetLoader()->TreeS()) {MakeTree("S");  MakeBranch("S");}
530     
531   TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
532   pRR->Select("HMPID",0,13);//select all HMPID DDL files
533   UInt_t w32=0;
534   while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
535     UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
536     sdi.Raw(ddl,w32);  
537     new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
538   }//raw records loop
539   GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
540   SdiReset();
541   return kTRUE;
542 }//Raw2SDigits
543 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
544 void AliHMPIDv2::StepCount()
545 {
546 // Count number of ckovs created  
547 }
548 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
549 void AliHMPIDv2::StepHistory()
550 {
551 // This methode is invoked from StepManager() in order to print out 
552   static Int_t iStepN;
553   const char *sParticle;
554   switch(gMC->TrackPid()){
555     case kProton:      sParticle="PROTON"    ;break;
556     case kNeutron:     sParticle="neutron"   ;break;
557     case kGamma:       sParticle="gamma"     ;break;
558     case 50000050:     sParticle="CKOV"      ;break;
559     case kPi0:         sParticle="Pi0"       ;break;  
560     case kPiPlus:      sParticle="Pi+"       ;break;  
561     case kPiMinus:     sParticle="Pi-"       ;break;  
562     case kElectron:    sParticle="electron"  ;break;  
563     default:           sParticle="not known" ;break;
564   }
565
566   TString flag="fanny combination";
567   if(gMC->IsTrackAlive())
568       if(gMC->IsTrackEntering())      flag="enters to";
569       else if(gMC->IsTrackExiting())  flag="exits from";
570       else if(gMC->IsTrackInside())   flag="inside";
571   else
572       if(gMC->IsTrackStop())          flag="stoped in";        
573   
574
575
576   
577   Int_t vid=0,copy=0;
578   TString path=gMC->CurrentVolName(); path.Prepend("-");path.Prepend(gMC->CurrentVolOffName(1));//current volume and his mother are always there
579   vid=gMC->CurrentVolOffID(2,copy);  if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
580   vid=gMC->CurrentVolOffID(3,copy);  if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
581  
582   
583   Printf("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f Etot=%.4f",iStepN,sParticle,gMC->TrackPid(),flag.Data(),path.Data(),gMC->TrackMass(),gMC->TrackCharge(),gMC->Edep()*1e9,gMC->Etot());
584   
585   Double_t gMcTrackPos[3]; gMC->TrackPosition(gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2]);
586   Double_t  gMcTrackPosLoc[3]; gMC->Gmtod(gMcTrackPos,gMcTrackPosLoc,1);
587   Printf("gMC Track Position (MARS) x: %5.3lf, y: %5.3lf, z: %5.3lf (r: %5.3lf) ---> (LOC) x: %5.3f, y: %5.3f, z: %5.3f",gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2],TMath::Sqrt(gMcTrackPos[0]*gMcTrackPos[0]+gMcTrackPos[1]*gMcTrackPos[1]+gMcTrackPos[2]*gMcTrackPos[2]),gMcTrackPosLoc[0],gMcTrackPosLoc[1],gMcTrackPosLoc[2]);
588   
589
590   
591   Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
592                             iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
593                             gMC->IsTrackAlive(), gMC->IsTrackDisappeared(),gMC->IsTrackEntering(), gMC->IsTrackExiting(),
594                             gMC->IsTrackInside(),gMC->IsTrackOut(),        gMC->IsTrackStop(),     gMC->IsNewTrack());
595   
596   Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
597   Int_t mid=gMC->CurrentMaterial(a,z,den,rad,abs);
598   Printf("Step %i: mid=%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);
599   
600   TArrayI proc;  gMC->StepProcesses(proc); 
601   Printf("Processes in this step:");
602   for ( int i = 0 ; i < proc.GetSize(); i++)
603   {
604     Printf("%s",TMCProcessName[proc.At(i)]);
605   }
606   Printf("End process list");
607   
608   iStepN++;
609 }//StepHistory()
610 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
611 void AliHMPIDv2::StepManager()
612 {
613 // Full Step Manager.
614 // Arguments: none
615 //   Returns: none           
616 //  StepHistory(); return; //uncomment to print tracks history
617  //  StepCount(); return;     //uncomment to count photons
618   
619   Int_t   copy; //volume copy aka node
620   
621 //Treat photons    
622   if((gMC->TrackPid()==50000050||gMC->TrackPid()==50000051)&&gMC->CurrentVolID(copy)==fIdPad){   //photon (Ckov or feedback) hit PC (fIdPad)
623     if(gMC->Edep()>0){                                                                           //photon survided QE test i.e. produces electron
624       if(IsLostByFresnel()){ gMC->StopTrack(); return;}                                          //photon lost due to fersnel reflection on PC       
625                        gMC->CurrentVolOffID(5,copy);                                             //current chamber since geomtry tree is Hmp-Hsec-Hgap-Hrow-Hcel-Hpad
626       Int_t   tid=     gMC->GetStack()->GetCurrentTrackNumber();                                 //take TID
627       Int_t   pid=     gMC->TrackPid();                                                          //take PID
628       Float_t etot=    gMC->Etot();                                                              //total hpoton energy, [GeV] 
629       Double_t x[3];   gMC->TrackPosition(x[0],x[1],x[2]);                                       //take MARS position at entrance to PC
630       Float_t xl,yl;   AliHMPIDParam::Instance()->Mars2Lors(copy,x,xl,yl);                       //take LORS position
631        if ( yl < 0  ) Printf("-------------------> SUPER PROBLEM PHOTON>>> Ch: %d, x[]: %f %f %f (MARS)-> xl: %f yl: %f",copy,x[0],x[1],x[2],xl,yl);
632       new((*fHits)[fNhits++])AliHMPIDHit(copy,etot,pid,tid,xl,yl,x);                             //HIT for photon, position at P, etot will be set to Q
633       GenFee(etot);                                                                              //generate feedback photons etot is modified in hit ctor to Q of hit
634     }//photon hit PC and DE >0 
635   }//photon hit PC
636   
637 //Treat charged particles  
638   static Float_t eloss;                                                                           //need to store mip parameters between different steps    
639   static Double_t in[3];
640   if(gMC->TrackCharge() && gMC->CurrentVolID(copy)==fIdCell){                                     //charged particle in amplification gap (fIdCell)
641     if(gMC->IsTrackEntering()||gMC->IsNewTrack()) {                                               //entering or newly created
642       eloss=0;                                                                                    //reset Eloss collector                         
643       gMC->TrackPosition(in[0],in[1],in[2]);                                                      //take position at the entrance
644     }else if(gMC->IsTrackExiting()||gMC->IsTrackStop()||gMC->IsTrackDisappeared()){               //exiting or disappeared
645       eloss              +=gMC->Edep();                                                           //take into account last step Eloss
646                           gMC->CurrentVolOffID(4,copy);                                           //take current chamber since geometry tree is Hmp-Hsec-Hgap-Hrow-Hcel
647       Int_t tid=          gMC->GetStack()->GetCurrentTrackNumber();                               //take TID
648       Int_t pid=          gMC->TrackPid();                                                        //take PID
649       Double_t out[3];    gMC->TrackPosition(out[0],out[1],out[2]);                               //take MARS position at exit
650       out[0]=0.5*(out[0]+in[0]);                                                                  //>
651       out[1]=0.5*(out[1]+in[1]);                                                                  //take hit position at the anod plane
652       out[2]=0.5*(out[2]+in[2]);                                                                  //>
653       Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(copy,out,xl,yl);                         //take LORS position
654        if ( yl < 0  ) Printf("-------------------> SUPER PROBLEM CHARGED>>> Ch: %d, x[]: %f %f %f (MARS)-> xl: %f yl: %f",copy,out[0],out[1],out[2],xl,yl);
655       new((*fHits)[fNhits++])AliHMPIDHit(copy,eloss,pid,tid,xl,yl,out);                           //HIT for MIP, position near anod plane, eloss will be set to Q 
656       GenFee(eloss);                                                                              //generate feedback photons 
657     }else                                                                                         //just going inside
658       eloss          += gMC->Edep();                                                              //collect this step eloss 
659   }//MIP in GAP
660 }//StepManager()
661 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++