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