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