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