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