]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDv1.cxx
HMPIDGeom.C renamed as Hgeom.C
[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
18#include "AliHMPIDParam.h" //CreateMaterials()
19#include "AliHMPIDHit.h" //Hits2SDigs(),StepManager()
20#include "AliHMPIDDigit.h" //CreateMaterials()
21#include "AliRawReader.h" //Raw2SDigits()
22#include <TParticle.h> //Hits2SDigits()
23#include <TRandom.h>
24#include <TVirtualMC.h> //StepManager() for gMC
25#include <TPDGCode.h> //StepHistory()
26#include <AliStack.h> //StepManager(),Hits2SDigits()
27#include <AliLoader.h> //Hits2SDigits()
28#include <AliRunLoader.h> //Hits2SDigits()
29#include <AliConst.h>
30#include <AliPDG.h>
31#include <AliMC.h> //StepManager()
32#include <AliRawDataHeader.h> //Digits2Raw()
33#include <AliDAQ.h> //Digits2Raw()
34#include <AliRun.h> //CreateMaterials()
35#include <AliMagF.h> //CreateMaterials()
36#include <TGeoManager.h> //CreateGeometry()
37#include <TMultiGraph.h> //Optics()
38#include <TGraph.h> //Optics()
39#include <TLegend.h> //Optics()
40#include <TCanvas.h> //Optics()
41#include <TF2.h> //CreateMaterials()
42#include <AliCDBManager.h> //CreateMaterials()
43#include <AliCDBEntry.h> //CreateMaterials()
44
45ClassImp(AliHMPIDv1)
46//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47void AliHMPIDv1::AddAlignableVolumes()const
48{
49// Associates the symbolic volume name with the corresponding volume path. Interface methode from AliModule ivoked from AliMC
50// Arguments: none
51// Returns: none
52 for(Int_t i=0;i<7;i++)
53 gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",i),Form("ALIC_1/HMPID_%i",i));
54}
55//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
56void AliHMPIDv1::CreateMaterials()
57{
58// Definition of available HMPID materials
59// Arguments: none
60// Returns: none
61 AliDebug(1,"Start v1 HMPID.");
62
63 Float_t emin=5.5,emax=8.5; //Photon energy range,[eV]
64
65 TF2 *pRaIF=new TF2("RidxRad","sqrt(1+0.554*(1239.84/x)^2/((1239.84/x)^2-5796)-0.0005*(y-20))" ,emin,emax,0,50); //DiMauro mail temp 0-50 degrees C
66 TF1 *pWiIF=new TF1("RidxWin","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
67 TF1 *pGaIF=new TF1("RidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)" ,emin,emax); //?????? from where
68
69 TF1 *pRaAF=new TF1("RabsRad","(x<7.8)*(gaus+gaus(3))+(x>=7.8)*0.0001" ,emin,emax); //fit from DiMauro data 28.10.03
70 pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
71 TF1 *pWiAF=new TF1("RabsWin","(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
72 TF1 *pGaAF=new TF1("RabsGap","(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
73
74 TF1 *pQeF =new TF1("Qe" ,"0+(x>6.07267)*0.344811*(1-exp(-1.29730*(x-6.07267)))" ,emin,emax); //fit from DiMauro data 28.10.03
75
76 const Int_t kNbins=30; //number of photon energy points
77 Float_t aEckov [kNbins];
78 Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
79 Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins];
80 Float_t aQeAll [kNbins], aQePc [kNbins];
81
82 for(Int_t i=0;i<kNbins;i++){
83 Float_t eV=emin+0.1*i; //Ckov energy in eV
84 aEckov [i] =1e-9*eV; //Ckov energy in GeV
85 aAbsRad[i]=pRaAF->Eval(eV); aIdxRad[i]=1.292;//pRaIF->Eval(eV,20); //Simulation for 20 degress C
86 aAbsWin[i]=pWiAF->Eval(eV); aIdxWin[i]=1.5787;//pWiIF->Eval(eV);
87 aAbsGap[i]=pGaAF->Eval(eV); aIdxGap[i]=1.0005;//pGaIF->Eval(eV); aQeAll[i] =1; //QE for all other materials except for PC must be 1.
88 aAbsMet[i] =0.0001; aIdxMet[i]=0; //metal ref idx must be 0 in order to reflect photon
89 aIdxPc [i]=1; aQePc [i]=pQeF->Eval(eV); //PC ref idx must be 1 in order to apply photon to QE conversion
90
91 }
92
93//data from PDG booklet 2002 density [gr/cm^3] rad len [cm] abs len [cm]
94 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
95 Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9} , wC6F14[2]={6 , 14} , dC6F14=1.68 ; Int_t nC6F14=-2;
96 Float_t aSiO2[2]={ 28.09 , 15.99} , zSiO2[2]={14 , 8} , wSiO2[2]={1 , 2} , dSiO2=2.64 ; Int_t nSiO2=-2;
97 Float_t aCH4[2]={ 12.01 , 1.01} , zCH4[2]={ 6 , 1} , wCH4[2]={1 , 4} , dCH4=7.17e-4 ; Int_t nCH4=-2;
98 Float_t aCsI[2]={132.90 ,126.90} , zCsI[2]={55 ,53} , wCsI[2]={1 , 1} , dCsI=0.1 ; Int_t nCsI=-2;
99 Float_t aRoha= 12.01 , zRoha= 6 , dRoha= 0.10 , radRoha= 18.80 , absRoha= 86.3/dRoha; //special material- quazi carbon
100 Float_t aCu= 63.55 , zCu= 29 , dCu= 8.96 , radCu= 1.43 , absCu= 134.9/dCu ;
101 Float_t aW=183.84 , zW= 74 , dW= 19.30 , radW= 0.35 , absW= 185.0/dW ;
102 Float_t aAl= 26.98 , zAl= 13 , dAl= 2.70 , radAl= 8.90 , absAl= 106.4/dAl ;
103
104 Int_t matId=0; //tmp material id number
105 Int_t unsens = 0, sens=1; //sensitive or unsensitive medium
106 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
107 Float_t maxfld = gAlice->Field()->Max(); //max field value
108 Float_t tmaxfd = -10.0; //max deflection angle due to magnetic field in one step
109 Float_t deemax = - 0.2; //max fractional energy loss in one step
110 Float_t stemax = - 0.1; //mas step allowed [cm]
111 Float_t epsil = 0.001; //abs tracking precision [cm]
112 Float_t stmin = - 0.001; //min step size [cm] in continius process transport, negative value: choose it automatically
113 AliMixture(++matId,"Air" ,aAir ,zAir ,dAir ,nAir ,wAir ); AliMedium(kAir ,"Air" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
114 AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
115 AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
116 AliMixture(++matId,"CH4" ,aCH4 ,zCH4 ,dCH4 ,nCH4 ,wCH4 ); AliMedium(kCH4 ,"CH4" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
117 AliMixture(++matId,"CsI" ,aCsI ,zCsI ,dCsI ,nCsI ,wCsI ); AliMedium(kCsI ,"CsI" ,matId, sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
118
119 AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha); AliMedium(kRoha,"Roha", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
120 AliMaterial(++matId,"Cu" ,aCu ,zCu ,dCu ,radCu ,absCu ); AliMedium(kCu ,"Cu" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
121 AliMaterial(++matId,"W" ,aW ,zW ,dW ,radW ,absW ); AliMedium(kW ,"W" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
122 AliMaterial(++matId,"Al" ,aAl ,zAl ,dAl ,radAl ,absAl ); AliMedium(kAl ,"Al" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
123
124 gMC->SetCerenkov((*fIdtmed)[kC6F14] , kNbins, aEckov, aAbsRad , aQeAll , aIdxRad );
125 gMC->SetCerenkov((*fIdtmed)[kSiO2] , kNbins, aEckov, aAbsWin , aQeAll , aIdxWin );
126 gMC->SetCerenkov((*fIdtmed)[kCH4] , kNbins, aEckov, aAbsGap , aQeAll , aIdxGap );
127 gMC->SetCerenkov((*fIdtmed)[kCu] , kNbins, aEckov, aAbsMet , aQeAll , aIdxMet );
128 gMC->SetCerenkov((*fIdtmed)[kW] , kNbins, aEckov, aAbsMet , aQeAll , aIdxMet ); //n=0 means reflect photons
129 gMC->SetCerenkov((*fIdtmed)[kCsI] , kNbins, aEckov, aAbsMet , aQePc , aIdxPc ); //n=1 means convert photons
130 gMC->SetCerenkov((*fIdtmed)[kAl] , kNbins, aEckov, aAbsMet , aQeAll , aIdxMet );
131
132 delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
133 AliDebug(1,"Stop v1 HMPID.");
134
135 TString ttl=GetTitle(); if(!ttl.Contains("ShowOptics")) return; //user didn't aks to plot optical curves
136
137 const Double_t kWidth=0.25,kHeight=0.2;
138 const Int_t kRadM=24 , kRadC=kRed;
139 const Int_t kWinM=26 , kWinC=kBlue;
140 const Int_t kGapM=25 , kGapC=kGreen;
141 const Int_t kPcM = 2 , kPcC =kMagenta;
142
143 Float_t aTraRad[kNbins],aTraWin[kNbins],aTraGap[kNbins],aTraTot[kNbins];
144 for(Int_t i=0;i<kNbins;i++){//calculate probability for photon to survive during transversing a volume of material with absorption length
145 aTraRad[i]=TMath::Exp(-AliHMPIDDigit::SizeRad()/ (aAbsRad[i]+0.0001)); //radiator
146 aTraWin[i]=TMath::Exp(-AliHMPIDDigit::SizeWin()/ (aAbsWin[i] +0.0001)); //window
147 aTraGap[i]=TMath::Exp(-AliHMPIDDigit::SizeGap()/ (aAbsGap[i] +0.0001)); //from window to PC
148 aTraTot[i]=aTraRad[i]*aTraWin[i]*aTraGap[i]*aQePc[i];
149 }
150
151 TGraph *pRaAG=new TGraph(kNbins,aEckov,aAbsRad);pRaAG->SetMarkerStyle(kRadM);pRaAG->SetMarkerColor(kRadC);
152 TGraph *pRaIG=new TGraph(kNbins,aEckov,aIdxRad);pRaIG->SetMarkerStyle(kRadM);pRaIG->SetMarkerColor(kRadC);
153 TGraph *pRaTG=new TGraph(kNbins,aEckov,aTraRad);pRaTG->SetMarkerStyle(kRadM);pRaTG->SetMarkerColor(kRadC);
154
155 TGraph *pWiAG=new TGraph(kNbins,aEckov,aAbsWin);pWiAG->SetMarkerStyle(kWinM);pWiAG->SetMarkerColor(kWinC);
156 TGraph *pWiIG=new TGraph(kNbins,aEckov,aIdxWin);pWiIG->SetMarkerStyle(kWinM);pWiIG->SetMarkerColor(kWinC);
157 TGraph *pWiTG=new TGraph(kNbins,aEckov,aTraWin);pWiTG->SetMarkerStyle(kWinM);pWiTG->SetMarkerColor(kWinC);
158
159 TGraph *pGaAG=new TGraph(kNbins,aEckov,aAbsGap);pGaAG->SetMarkerStyle(kGapM);pGaAG->SetMarkerColor(kGapC);
160 TGraph *pGaIG=new TGraph(kNbins,aEckov,aIdxGap);pGaIG->SetMarkerStyle(kGapM);pGaIG->SetMarkerColor(kGapC);
161 TGraph *pGaTG=new TGraph(kNbins,aEckov,aTraGap);pGaTG->SetMarkerStyle(kGapM);pGaTG->SetMarkerColor(kGapC);
162
163 TGraph *pQeG =new TGraph(kNbins,aEckov,aQePc); pQeG ->SetMarkerStyle(kPcM );pQeG->SetMarkerColor(kPcC);
164 TGraph *pToG =new TGraph(kNbins,aEckov,aTraTot);pToG ->SetMarkerStyle(30) ;pToG->SetMarkerColor(kYellow);
165
166 TMultiGraph *pIdxMG=new TMultiGraph("idx","Ref index;E_{#check{C}} [GeV]");
167 TMultiGraph *pAbsMG=new TMultiGraph("abs","Absorption [cm];E_{#check{C}} [GeV]");
168 TMultiGraph *pTraMG=new TMultiGraph("tra","Transmission;E_{#check{C}} [GeV]"); TLegend *pTraLe=new TLegend(0.2,0.4,0.2+kWidth,0.4+kHeight);
169 pAbsMG->Add(pRaAG); pIdxMG->Add(pRaIG); pTraMG->Add(pRaTG); pTraLe->AddEntry(pRaTG, "Rad", "p");
170 pAbsMG->Add(pWiAG); pIdxMG->Add(pWiIG); pTraMG->Add(pWiTG); pTraLe->AddEntry(pWiTG, "Win", "p");
171 pAbsMG->Add(pGaAG); pIdxMG->Add(pGaIG); pTraMG->Add(pGaTG); pTraLe->AddEntry(pGaTG, "Gap", "p");
172 pTraMG->Add(pToG); pTraLe->AddEntry(pToG, "Tot", "p");
173 pTraMG->Add(pQeG); pTraLe->AddEntry(pQeG, "QE" , "p");
174 TCanvas *pC=new TCanvas("c1","HMPID optics to check",1100,900); pC->Divide(2,2);
175 pC->cd(1); pIdxMG->Draw("AP");
176 pC->cd(2); gPad->SetLogy(); pAbsMG->Draw("AP");
177 pC->cd(3); pTraLe->Draw();
178 pC->cd(4); pTraMG->Draw("AP");
179}//void AliHMPID::CreateMaterials()
180//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
181void AliHMPIDv1::CreateGeometry()
182{
183//Creates detailed geometry simulation (currently GEANT volumes tree)
184 AliDebug(1,"Start main.");
185 if(!gMC->IsRootGeometrySupported()) return;
186
187 Double_t cm=1,mm=0.1*cm,mkm=0.001*mm,dx,dy,dz;//default is cm
188
189 TGeoVolume *pRich=gGeoManager->MakeBox("HMPID",gGeoManager->GetMedium("HMPID_CH4"),dx=(6*mm+1681*mm+6*mm)/2, //main HMPID volume
190 dy=(6*mm+1466*mm+6*mm)/2,
191 dz=(80*mm+40*mm)*2/2); //x,y taken from 2033P1 z from p84 TDR
192 const Double_t kAngHor=19.5; // horizontal angle between chambers 19.5 grad
193 const Double_t kAngVer=20; // vertical angle between chambers 20 grad
194 const Double_t kAngCom=30; // common HMPID rotation with respect to x axis 30 grad
195 const Double_t trans[3]={490,0,0}; //center of the chamber is on window-gap surface
196 for(Int_t iCh=0;iCh<7;iCh++){//place 7 chambers
197 TGeoHMatrix *pMatrix=new TGeoHMatrix;
198 pMatrix->RotateY(90); //rotate around y since initial position is in XY plane -> now in YZ plane
199 pMatrix->SetTranslation(trans); //now plane in YZ is shifted along x
200 switch(iCh){
201 case 0: pMatrix->RotateY(kAngHor); pMatrix->RotateZ(-kAngVer); break; //right and down
202 case 1: pMatrix->RotateZ(-kAngVer); break; //down
203 case 2: pMatrix->RotateY(kAngHor); break; //right
204 case 3: break; //no rotation
205 case 4: pMatrix->RotateY(-kAngHor); break; //left
206 case 5: pMatrix->RotateZ(kAngVer); break; //up
207 case 6: pMatrix->RotateY(-kAngHor); pMatrix->RotateZ(kAngVer); break; //left and up
208 }
209 pMatrix->RotateZ(kAngCom); //apply common rotation in XY plane
210 gGeoManager->GetVolume("ALIC")->AddNode(pRich,iCh,pMatrix);
211 }
212
213 Float_t par[3];
214 Int_t matrixIdReturn=0; //matrix id returned by AliMatrix
215//Pad Panel frame 6 sectors
216 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
217 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
218 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
219 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)
220
221 gMC->Gspos("Rppf",0,"HMPID", -335*mm, -433*mm, 8*cm+20*mm, 0,"ONLY");//F1 2040P1 z p.84 TDR
222 gMC->Gspos("Rppf",1,"HMPID", +335*mm, -433*mm, 8*cm+20*mm, 0,"ONLY");
223 gMC->Gspos("Rppf",2,"HMPID", -335*mm, 0*mm, 8*cm+20*mm, 0,"ONLY");
224 gMC->Gspos("Rppf",3,"HMPID", +335*mm, 0*mm, 8*cm+20*mm, 0,"ONLY");
225 gMC->Gspos("Rppf",4,"HMPID", -335*mm, +433*mm, 8*cm+20*mm, 0,"ONLY");
226 gMC->Gspos("Rppf",5,"HMPID", +335*mm, +433*mm, 8*cm+20*mm, 0,"ONLY");
227 gMC->Gspos("Rpc" ,1,"Rppf", 0*mm, 0*mm, -19.15*mm, 0,"ONLY");//PPF 2001P2
228 gMC->Gspos("RppfLarge",1,"Rppf", -224.5*mm, -151.875*mm, 0.85*mm, 0,"ONLY");
229 gMC->Gspos("RppfLarge",2,"Rppf", -224.5*mm, - 50.625*mm, 0.85*mm, 0,"ONLY");
230 gMC->Gspos("RppfLarge",3,"Rppf", -224.5*mm, + 50.625*mm, 0.85*mm, 0,"ONLY");
231 gMC->Gspos("RppfLarge",4,"Rppf", -224.5*mm, +151.875*mm, 0.85*mm, 0,"ONLY");
232 gMC->Gspos("RppfSmall",1,"Rppf", - 65.0*mm, -151.875*mm, 0.85*mm, 0,"ONLY");
233 gMC->Gspos("RppfSmall",2,"Rppf", - 65.0*mm, - 50.625*mm, 0.85*mm, 0,"ONLY");
234 gMC->Gspos("RppfSmall",3,"Rppf", - 65.0*mm, + 50.625*mm, 0.85*mm, 0,"ONLY");
235 gMC->Gspos("RppfSmall",4,"Rppf", - 65.0*mm, +151.875*mm, 0.85*mm, 0,"ONLY");
236 gMC->Gspos("RppfSmall",5,"Rppf", + 65.0*mm, -151.875*mm, 0.85*mm, 0,"ONLY");
237 gMC->Gspos("RppfSmall",6,"Rppf", + 65.0*mm, - 50.625*mm, 0.85*mm, 0,"ONLY");
238 gMC->Gspos("RppfSmall",7,"Rppf", + 65.0*mm, + 50.625*mm, 0.85*mm, 0,"ONLY");
239 gMC->Gspos("RppfSmall",8,"Rppf", + 65.0*mm, +151.875*mm, 0.85*mm, 0,"ONLY");
240 gMC->Gspos("RppfLarge",5,"Rppf", +224.5*mm, -151.875*mm, 0.85*mm, 0,"ONLY");
241 gMC->Gspos("RppfLarge",6,"Rppf", +224.5*mm, - 50.625*mm, 0.85*mm, 0,"ONLY");
242 gMC->Gspos("RppfLarge",7,"Rppf", +224.5*mm, + 50.625*mm, 0.85*mm, 0,"ONLY");
243 gMC->Gspos("RppfLarge",8,"Rppf", +224.5*mm, +151.875*mm, 0.85*mm, 0,"ONLY");
244//Gap - anod wires 6 copies to HMPID
245 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
246 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
247 AliMatrix(matrixIdReturn,180,0, 90,90, 90,0); //wires along x
248
249 gMC->Gspos("Rgap",0,"HMPID", -335*mm, -433*mm,8*cm-2.225*mm, 0,"ONLY"); //F1 2040P1 z WP 2099P1
250 gMC->Gspos("Rgap",1,"HMPID", +335*mm, -433*mm,8*cm-2.225*mm, 0,"ONLY");
251 gMC->Gspos("Rgap",2,"HMPID", -335*mm, 0*mm,8*cm-2.225*mm, 0,"ONLY");
252 gMC->Gspos("Rgap",3,"HMPID", +335*mm, 0*mm,8*cm-2.225*mm, 0,"ONLY");
253 gMC->Gspos("Rgap",4,"HMPID", -335*mm, +433*mm,8*cm-2.225*mm, 0,"ONLY");
254 gMC->Gspos("Rgap",5,"HMPID", +335*mm, +433*mm,8*cm-2.225*mm, 0,"ONLY");
255 for(int i=1;i<=96;i++)
256 gMC->Gspos("Rano",i,"Rgap", 0*mm, -411/2*mm+i*4*mm, 0.185*mm, matrixIdReturn,"ONLY"); //WP 2099P1
257//Defines radiators geometry
258 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
259 par[0]=1330*mm/2 ;par[1]= 413*mm/2 ;par[2]= 4*mm/2; gMC->Gsvolu("RradFront" ,"BOX ",(*fIdtmed)[kRoha] ,par,3); //front
260 par[0]=1330*mm/2 ;par[1]= 413*mm/2 ;par[2]= 5*mm/2; gMC->Gsvolu("RradWin" ,"BOX ",(*fIdtmed)[kSiO2] ,par,3); //window
261 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
262 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
263 par[0]= 0 ;par[1]= 10*mm/2 ;par[2]= 15*mm/2; gMC->Gsvolu("RradSpacer","TUBE",(*fIdtmed)[kSiO2] ,par,3); //spacer
264
265 gMC->Gspos("Rrad",1,"HMPID", 0*mm,-434*mm, -12*mm, 0,"ONLY"); //3 radiators to HMPID
266 gMC->Gspos("Rrad",2,"HMPID", 0*mm, 0*mm, -12*mm, 0,"ONLY");
267 gMC->Gspos("Rrad",3,"HMPID", 0*mm,+434*mm, -12*mm, 0,"ONLY");
268 gMC->Gspos("RradFront",1,"Rrad", 0*mm, 0*mm, -10.0*mm, 0,"ONLY"); //front cover
269 gMC->Gspos("RradWin" ,1,"Rrad", 0*mm, 0*mm, 9.5*mm, 0,"ONLY"); //quartz window (back cover)
270 gMC->Gspos("RradLong" ,1,"Rrad", 0*mm,-204*mm, -0.5*mm, 0,"ONLY"); //long side
271 gMC->Gspos("RradLong" ,2,"Rrad", 0*mm,+204*mm, -0.5*mm, 0,"ONLY"); //long side
272 gMC->Gspos("RradShort",1,"Rrad",-660*mm, 0*mm, -0.5*mm, 0,"ONLY"); //short side
273 gMC->Gspos("RradShort",2,"Rrad",+660*mm, 0*mm, -0.5*mm, 0,"ONLY"); //short side
274 for(int i=0;i<3;i++)
275 for(int j=0;j<10;j++)
276 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
277//Defines SandBox geometry
278 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
279 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
280 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
281
282 gMC->Gspos("Rsb",1,"HMPID", 0*mm, 0*mm, -73.75*mm, 0,"ONLY"); //p.84 TDR sandbox to rich
283 gMC->Gspos("RsbComb" ,1,"Rsb", 0*mm, 0*mm, 0*mm, 0,"ONLY"); //2072P1 honeycomv to sandbox
284 gMC->Gspos("RsbCover",1,"Rsb", 0*mm, 0*mm, +25*mm, 0,"ONLY"); //cover to sandbox
285 gMC->Gspos("RsbCover",2,"Rsb", 0*mm, 0*mm, -25*mm, 0,"ONLY"); //cover to sandbox
286 AliDebug(1,"Stop v1. HMPID option");
287}//CreateGeometry()
288//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
289void AliHMPIDv1::Init()
290{
291// This methode defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(gMC->CurrentVolID()==XXX) statements in StepManager()
292// Arguments: none
293// Returns: none
294 AliDebug(1,"Start v1 HMPID.");
295 fIdRad = gMC->VolId("Rrad");
296 fIdWin = gMC->VolId("RradWin");
297 fIdPc = gMC->VolId("Rpc");
298 fIdAmpGap = gMC->VolId("Rgap");
299 fIdProxGap = gMC->VolId("Rgap");
300 AliDebug(1,"Stop v1 HMPID.");
301}
302//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
303Bool_t AliHMPIDv1::IsLostByFresnel()
304{
305// Calculate probability for the photon to be lost by Fresnel reflection.
306 TLorentzVector p4;
307 Double_t mom[3],localMom[3];
308 gMC->TrackMomentum(p4); mom[0]=p4(1); mom[1]=p4(2); mom[2]=p4(3);
309 localMom[0]=0; localMom[1]=0; localMom[2]=0;
310 gMC->Gmtod(mom,localMom,2);
311 Double_t localTc = localMom[0]*localMom[0]+localMom[2]*localMom[2];
312 Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
313 Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
314 if(gMC->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
315 AliDebug(1,"Photon lost");
316 return kTRUE;
317 }else
318 return kFALSE;
319}//IsLostByFresnel()
320//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
321void AliHMPIDv1::GenFee(Int_t iCh,Float_t eloss)
322{
323// Generate FeedBack photons for the current particle. To be invoked from StepManager().
324// eloss=0 means photon so only pulse height distribution is to be analysed. This one is done in AliHMPIDParam::TotQdc()
325 TLorentzVector x4;
326 gMC->TrackPosition(x4);
327 Int_t iNphotons=gMC->GetRandom()->Poisson(0.02*200); eloss++; iCh++; //??????????????????????
328 AliDebug(1,Form("N photons=%i",iNphotons));
329 Int_t j;
330 Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
331//Generate photons
332 for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
333 Double_t ranf[2];
334 gMC->GetRandom()->RndmArray(2,ranf); //Sample direction
335 cthf=ranf[0]*2-1.0;
336 if(cthf<0) continue;
337 sthf = TMath::Sqrt((1 - cthf) * (1 + cthf));
338 phif = ranf[1] * 2 * TMath::Pi();
339
340 if(Double_t randomNumber=gMC->GetRandom()->Rndm()<=0.57)
341 enfp = 7.5e-9;
342 else if(randomNumber<=0.7)
343 enfp = 6.4e-9;
344 else
345 enfp = 7.9e-9;
346
347
348 dir[0] = sthf * TMath::Sin(phif); dir[1] = cthf; dir[2] = sthf * TMath::Cos(phif);
349 gMC->Gdtom(dir, mom, 2);
350 mom[0]*=enfp; mom[1]*=enfp; mom[2]*=enfp;
351 mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
352
353 // Polarisation
354 e1[0]= 0; e1[1]=-dir[2]; e1[2]= dir[1];
355 e2[0]=-dir[1]; e2[1]= dir[0]; e2[2]= 0;
356 e3[0]= dir[1]; e3[1]= 0; e3[2]=-dir[0];
357
358 vmod=0;
359 for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
360 if (!vmod) for(j=0;j<3;j++) {
361 uswop=e1[j];
362 e1[j]=e3[j];
363 e3[j]=uswop;
364 }
365 vmod=0;
366 for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
367 if (!vmod) for(j=0;j<3;j++) {
368 uswop=e2[j];
369 e2[j]=e3[j];
370 e3[j]=uswop;
371 }
372
373 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;
374 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;
375
376 phi = gMC->GetRandom()->Rndm()* 2 * TMath::Pi();
377 for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
378 gMC->Gdtom(pol, pol, 2);
379 Int_t outputNtracksStored;
380 gAlice->GetMCApp()->PushTrack(1, //transport
381 gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track
382 kFeedback, //PID
383 mom[0],mom[1],mom[2],mom[3], //track momentum
384 x4.X(),x4.Y(),x4.Z(),x4.T(), //track origin
385 pol[0],pol[1],pol[2], //polarization
386 kPFeedBackPhoton, //process ID
387 outputNtracksStored, //on return how many new photons stored on stack
388 1.0); //weight
389 }//feedbacks loop
390 AliDebug(1,"Stop.");
391}//GenerateFeedbacks()
392//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
393void AliHMPIDv1::Hits2SDigits()
394{
395// Interface methode ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
396// Arguments: none
397// Returns: none
398 AliDebug(1,"Start.");
399 for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){ //events loop
400 GetLoader()->GetRunLoader()->GetEvent(iEvt); //get next event
401
402 if(!GetLoader()->TreeH()) {GetLoader()->LoadHits(); }
403 if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
404
405 for(Int_t iPrimN=0;iPrimN<GetLoader()->TreeH()->GetEntries();iPrimN++){//prims loop
406 GetLoader()->TreeH()->GetEntry(iPrimN);
407 Hit2Sdi(Hits(),SdiLst());
408 }//prims loop
409 GetLoader()->TreeS()->Fill();
410 GetLoader()->WriteSDigits("OVERWRITE");
411 SdiReset();
412 }//events loop
413 GetLoader()->UnloadHits();
414 GetLoader()->UnloadSDigits();
415 AliDebug(1,"Stop.");
416}//Hits2SDigits()
417//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
418void AliHMPIDv1::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
419{
420// Converts list of hits to list of sdigits. For each hit in a loop the following steps are done:
421// - calcultion of the total charge induced by the hit
422// - determination of the pad contaning the hit and shifting hit y position to the nearest anod wire y
423// - defining a set of pads affected (up to 9 including the hitted pad)
424// - calculating charge induced to all those pads using integrated Mathieson distribution and creating sdigit
425// Arguments: pHitLst - list of hits provided not empty
426// pSDigLst - list of sdigits where to store the results
427// Returns: none
428 for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){ //hits loop
429 AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit); //get pointer to current hit
430 AliHMPIDDigit::Hit2Sdi(pHit,pSdiLst); //convert this hit to list of sdigits
431 }//hits loop loop
432}//Hits2SDigs() for TVector2
433//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
434void AliHMPIDv1::Digits2Raw()
435{
436// Creates raw data files in DDL format. Invoked by AliSimulation where loop over events is done
437// Arguments: none
438// Returns: none
439 AliDebug(1,"Start.");
440 GetLoader()->LoadDigits();
441 TTree * treeD = GetLoader()->TreeD();
442 if(!treeD) {
443 AliError("No digits tree!");
444 return;
445 }
446 treeD->GetEntry(0);
447
448 ofstream file[AliHMPIDDigit::kNddls]; //output streams
449 Int_t cnt[AliHMPIDDigit::kNddls]; //data words counters for DDLs
450 AliRawDataHeader header; //empty DDL header
451 UInt_t w32=0; //32 bits data word
452
453 for(Int_t i=0;i<AliHMPIDDigit::kNddls;i++){
454 file[i].open(AliDAQ::DdlFileName(GetName(),i)); //open all 14 DDL in parallel
455 file[i].write((char*)&header,sizeof(header)); //write dummy header as place holder, actual will be written later when total size of DDL is known
456 cnt[i]=0; //reset counters
457 }
458
459 for(Int_t iCh=0;iCh<7;iCh++)
460 for(Int_t iDig=0;iDig<DigLst(iCh)->GetEntriesFast();iDig++){//digits loop for a given chamber
461 AliHMPIDDigit *pDig=(AliHMPIDDigit*)DigLst(iCh)->At(iDig);
462 Int_t ddl=pDig->Raw(w32); //ddl is 0..13
463 file[ddl].write((char*)&w32,sizeof(w32)); cnt[ddl]++;//write formated digit to the propriate file (as decided in Dig2Raw) and increment corresponding counter
464 }//digits
465
466
467 for(Int_t i=0;i<AliHMPIDDigit::kNddls;i++){
468 header.fSize=sizeof(header)+cnt[i]*sizeof(w32); //now calculate total number of bytes for each DDL file
469 header.SetAttribute(0);
470 file[i].seekp(0); file[i].write((char*)&header,sizeof(header));//rewrite DDL header with fSize field properly set
471 file[i].close(); //close DDL file
472 }
473 GetLoader()->UnloadDigits();
474 AliDebug(1,"Stop.");
475}//Digits2Raw()
476//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
477Float_t AliHMPIDv1::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
478{
479// Correction for Fresnel ???????????
480// Arguments: ene - photon energy [GeV],
481// PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
482// Returns:
483 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,
484 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,
485 7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
486 Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
487 2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
488 2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
489 1.72,1.53};
490 Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
491 0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
492 0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
493 1.714,1.498};
494 Float_t xe=ene;
495 Int_t j=Int_t(xe*10)-49;
496 Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
497 Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
498
499 //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
500 //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
501
502 Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
503 Float_t tanin=sinin/pdoti;
504
505 Float_t c1=cn*cn-ck*ck-sinin*sinin;
506 Float_t c2=4*cn*cn*ck*ck;
507 Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
508 Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
509
510 Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
511 Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
512
513
514 //CORRECTION FACTOR FOR SURFACE ROUGHNESS
515 //B.J. STAGG APPLIED OPTICS, 30(1991),4113
516
517 Float_t sigraf=18.;
518 Float_t lamb=1240/ene;
519 Float_t fresn;
520
521 Float_t rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
522
523 if(pola)
524 {
525 Float_t pdotr=0.8; //DEGREE OF POLARIZATION : 1->P , -1->S
526 fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
527 }
528 else
529 fresn=0.5*(rp+rs);
530
531 fresn = fresn*rO;
532 return fresn;
533}//Fresnel()
534//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
535void AliHMPIDv1::Print(Option_t *option)const
536{
537// Debug printout
538 TObject::Print(option);
539}//void AliHMPID::Print(Option_t *option)const
540//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
541Bool_t AliHMPIDv1::Raw2SDigits(AliRawReader *pRR)
542{
543// Interface methode ivoked from AliSimulation to create a list of sdigits from raw digits. Events loop is done in AliSimulation
544// Arguments: pRR- raw reader
545// Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())
546 AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
547
548 if(!GetLoader()->TreeS()) {MakeTree("S"); MakeBranch("S");}
549
550 TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
551 pRR->Select("HMPID",0,13);//select all HMPID DDL files
552 UInt_t w32=0;
553 while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
554 UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
555 sdi.ReadRaw(ddl,w32);
556 new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
557 }//raw records loop
558 GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
559 SdiReset();
560 return kTRUE;
561}//Raw2SDigits
562//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
563void AliHMPIDv1::StepCount()
564{
565// Count number of ckovs created
566}
567//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
568void AliHMPIDv1::StepHistory()
569{
570// This methode is invoked from StepManager() in order to print out
571 static Int_t iStepN;
572 const char *sParticle;
573 switch(gMC->TrackPid()){
574 case kProton: sParticle="PROTON" ;break;
575 case kNeutron: sParticle="neutron" ;break;
576 case kGamma: sParticle="gamma" ;break;
577 case kCerenkov: sParticle="CKOV" ;break;
578 case kPi0: sParticle="Pi0" ;break;
579 case kPiPlus: sParticle="Pi+" ;break;
580 case kPiMinus: sParticle="Pi-" ;break;
581 case kElectron: sParticle="electron" ;break;
582 default: sParticle="not known" ;break;
583 }
584
585 TString flag="fanny combination";
586 if(gMC->IsTrackAlive())
587 if(gMC->IsTrackEntering()) flag="enters to";
588 else if(gMC->IsTrackExiting()) flag="exits from";
589 else if(gMC->IsTrackInside()) flag="inside";
590 else
591 if(gMC->IsTrackStop()) flag="stoped in";
592
593 Int_t vid=0,copy=0;
594 TString path=gMC->CurrentVolName(); path.Prepend("-");path.Prepend(gMC->CurrentVolOffName(1));//current volume and his mother are always there
595 vid=gMC->CurrentVolOffID(2,copy); if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
596 vid=gMC->CurrentVolOffID(3,copy); if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
597
598 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);
599
600 Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
601 iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
602 gMC->IsTrackAlive(), gMC->IsTrackDisappeared(),gMC->IsTrackEntering(), gMC->IsTrackExiting(),
603 gMC->IsTrackInside(),gMC->IsTrackOut(), gMC->IsTrackStop(), gMC->IsNewTrack());
604
605 Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
606 Int_t mid=gMC->CurrentMaterial(a,z,den,rad,abs);
607 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);
608 iStepN++;
609}//StepHistory()
610//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
611void AliHMPIDv1::StepManager()
612{
613// Full Step Manager.
614// Arguments: none
615// Returns: none
616// StepHistory(); return; //uncomment to print tracks history
617// StepCount(); return; //uncomment to count photons
618
619 Int_t copy; //volume copy aka node
620
621//Treat photons
622 if((gMC->TrackPid()==kCerenkov||gMC->TrackPid()==kFeedback)&&gMC->CurrentVolID(copy)==fIdPc){ //photon (Ckov or feedback) hit PC (fIdPc)
623 if(gMC->Edep()>0){ //photon survided QE test i.e. produces electron
624 if(IsLostByFresnel()){ gMC->StopTrack(); return;} //photon lost due to fersnel reflection on PC
625 gMC->CurrentVolOffID(2,copy); //current chamber since geomtry tree is HMPID-Rppf-Rpc
626 Int_t tid= gMC->GetStack()->GetCurrentTrackNumber(); //take TID
627 Int_t pid= gMC->TrackPid(); //take PID
628 Float_t etot= gMC->Etot(); //total hpoton energy, [GeV]
629 Double_t x[3]; gMC->TrackPosition(x[0],x[1],x[2]); //take MARS position at entrance to PC
630 Float_t xl,yl; AliHMPIDParam::Instance()->Mars2Lors(copy,x,xl,yl); //take LORS position
631 new((*fHits)[fNhits++])AliHMPIDHit(copy,etot,pid,tid,xl,yl,x); //HIT for photon, position at PC
632 GenFee(copy); //generate feedback photons
633 }//photon hit PC and DE >0
634 }//photon hit PC
635
636//Treat charged particles
637 static Double_t dEdX; //need to store mip parameters between different steps
638 static Double_t in[3];
639 if(gMC->TrackCharge() && gMC->CurrentVolID(copy)==fIdAmpGap){ //charged particle in amplification gap (fIdAmpGap)
640 if(gMC->IsTrackEntering()||gMC->IsNewTrack()) { //entering or newly created
641 dEdX=0; //reset dEdX collector
642 gMC->TrackPosition(in[0],in[1],in[2]); //take position at the entrance
643 }else if(gMC->IsTrackExiting()||gMC->IsTrackStop()||gMC->IsTrackDisappeared()){ //exiting or disappeared
644 dEdX +=gMC->Edep(); //take into account last step dEdX
645 gMC->CurrentVolOffID(1,copy); //take current chamber since geometry tree is HMPID-Rgap
646 Int_t tid= gMC->GetStack()->GetCurrentTrackNumber(); //take TID
647 Int_t pid= gMC->TrackPid(); //take PID
648 Double_t out[3]; gMC->TrackPosition(out[0],out[1],out[2]); //take MARS position at exit
649 out[0]=0.5*(out[0]+in[0]); out[1]=0.5*(out[1]+in[1]); out[1]=0.5*(out[1]+in[1]); //take hit position at the anod plane
650 Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(copy,out,xl,yl); //take LORS position
651 new((*fHits)[fNhits++])AliHMPIDHit(copy,dEdX,pid,tid,xl,yl,out); //HIT for MIP, position near anod plane
652 GenFee(copy,dEdX); //generate feedback photons
653 }else //just going inside
654 dEdX += gMC->Edep(); //collect this step dEdX
655 }//MIP in GAP
656}//StepManager()
657//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++