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