]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSv0.cxx
Added protection. In case IROT=0 the address Q(LQ(JROTM-IROT)) should not
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv0.cxx
CommitLineData
fe4da5cc 1/////////////////////////////////////////////////////////
2// Manager and hits classes for set:PHOS version 0 //
3/////////////////////////////////////////////////////////
4
5// --- ROOT system ---
6#include "TH1.h"
7#include "TRandom.h"
8#include "TFile.h"
9#include "TTree.h"
10#include "TBRIK.h"
11#include "TNode.h"
12
13// --- galice header files ---
14#include "AliPHOSv0.h"
15#include "AliRun.h"
16#include "AliConst.h"
17#include "AliMC.h"
18
19ClassImp(AliPHOSv0)
20
21//______________________________________________________________________________
22
23
24AliPHOSv0::AliPHOSv0() : AliPHOS()
25{
26 fIdSens=0;
27}
28
29//______________________________________________________________________________
30
31AliPHOSv0::AliPHOSv0(const char *name, const char *title)
32 : AliPHOS(name, title)
33{
34 fIdSens=0;
35}
36
37//___________________________________________
38void AliPHOSv0::Init()
39{
40 AliMC* pMC = AliMC::GetMC();
41
42 fIdSens=pMC->VolId("PXTL");
43}
44
45//___________________________________________
46void AliPHOSv0::CreateGeometry()
47{
48// *** DEFINITION OF THE -0.25<y<0.25 TILTED GEOMETRY OF THE PHOS ***
49// ORIGIN : NICK VAN EIJNDHOVEN
50
51 AliMC* pMC = AliMC::GetMC();
52
53 Float_t pphi;
54 Float_t r, dptcb[3], dpair[3], dphos[3], dpucp[3], dpasp[3], dpcpv[3];
55 Float_t dpxtl[3];
56 Float_t yo;
57 Int_t idrotm[99];
58 Float_t xp1, yp1, xp2, yp2;
59
60 Int_t *idtmed = gAlice->Idtmed();
61
62// --- Dimensions of PbWO4 crystal ---
63 const Float_t XTL_X=2.2;
64 const Float_t XTL_Y=18.;
65 const Float_t XTL_Z=2.2;
66// --- Tyvek wrapper thickness
67 const Float_t PAP_THICK=0.01;
68// --- CPV thickness ---
69 const Float_t CPV_Y=0.5;
70// --- Polystyrene Foam Outer Cover dimensions ---
71 const Float_t FOC_X=214.6;
72 const Float_t FOC_Y=80.;
73 const Float_t FOC_Z=260.;
74// --- Inner AIR volume dimensions ---
75 const Float_t AIR_X=206.;
76 const Float_t AIR_Y=66.;
77 const Float_t AIR_Z=244.;
78// --- Tyvek Crystal Block dimensions ---
79 const Float_t TCB_X=198.;
80 const Float_t TCB_Y=25.;
81 const Float_t TCB_Z=234.;
82// --- Upper Cooling Plate thickness ---
83 const Float_t UCP_Y=0.06;
84// --- Al Support Plate thickness ---
85 const Float_t ASP_Y=10.;
86//--- Distance from IP to Foam Outer Cover top plate (needs to be 447.) ---
87 const Float_t FOC_R=467.;
88//--- Distance from IP to Crystal Block top Surface (needs to be 460.) ---
89 const Float_t CBS_R=480.;
90
91// --- Dimensions of volumes ---
92
93
94// --- Define PHOS box volume, fill with Polystyrene foam ---
95 dphos[0] = FOC_X/2.;
96 dphos[1] = FOC_Y/2.;
97 dphos[2] = FOC_Z/2.;
98 pMC->Gsvolu("PHOS", "BOX ", idtmed[703], dphos, 3);
99
100// --- Define air-filled box, place inside PHOS ---
101 dpair[0] = AIR_X/2.;
102 dpair[1] = AIR_Y/2.;
103 dpair[2] = AIR_Z/2.;
104 pMC->Gsvolu("PAIR", "BOX ", idtmed[798], dpair, 3);
105 pMC->Gspos("PAIR", 1, "PHOS", 0., 0., 0., 0, "ONLY");
106
107// --- Define Upper Cooling Panel ---
108// --- place it right behind upper foam plate ---
109 dpucp[0] = TCB_X/2.;
110 dpucp[1] = UCP_Y/2.;
111 dpucp[2] = TCB_Z/2.;
112 pMC->Gsvolu("PUCP", "BOX ", idtmed[701], dpucp, 3);
113 yo = (AIR_Y-UCP_Y)/2.;
114 pMC->Gspos("PUCP", 1, "PAIR", 0., yo, 0., 0, "ONLY");
115
116// --- Define Crystal Block, fill with Tyvek, position inside PAIR ---
117 dptcb[0] = TCB_X/2.;
118 dptcb[1] = TCB_Y/2.;
119 dptcb[2] = TCB_Z/2.;
120 pMC->Gsvolu("PTCB", "BOX ", idtmed[702], dptcb, 3);
121// --- Divide PTCB in X and Z directions --
122 pMC->Gsdvn("PSEC", "PTCB", 11, 1);
123 pMC->Gsdvn("PMOD", "PSEC", 13, 3);
124 pMC->Gsdvn("PSTR", "PMOD", 8, 1);
125 pMC->Gsdvn("PCEL", "PSTR", 8, 3);
126 yo = (FOC_Y-TCB_Y)/2. -(CBS_R-FOC_R);
127 pMC->Gspos("PTCB", 1, "PAIR", 0., yo, 0., 0, "ONLY");
128
129// --- Define PbWO4 crystal volume, place inside PCEL ---
130 dpxtl[0] = XTL_X/2.;
131 dpxtl[1] = XTL_Y/2.;
132 dpxtl[2] = XTL_Z/2.;
133 pMC->Gsvolu("PXTL", "BOX ", idtmed[699], dpxtl, 3);
134 yo = (TCB_Y-XTL_Y)/2. - PAP_THICK;
135 pMC->Gspos("PXTL", 1, "PCEL", 0., yo, 0., 0, "ONLY");
136
137// --- Define Al Support Plate, position it inside PAIR ---
138// --- right beneath PTCB ---
139 dpasp[0] = AIR_X/2.;
140 dpasp[1] = ASP_Y/2.;
141 dpasp[2] = AIR_Z/2.;
142 pMC->Gsvolu("PASP", "BOX ", idtmed[701], dpasp, 3);
143 yo = (FOC_Y-ASP_Y)/2. - (CBS_R-FOC_R+TCB_Y);
144 pMC->Gspos("PASP", 1, "PAIR", 0., yo, 0., 0, "ONLY");
145
146// --- Define CPV volume, DON'T PLACE IT YET ---
147 dpcpv[0] = TCB_X/2.;
148 dpcpv[1] = CPV_Y/2.;
149 dpcpv[2] = TCB_Z/2.;
150 pMC->Gsvolu("PCPV", "BOX ", idtmed[700], dpcpv, 3);
151// --- Divide in X and Z direction (same way as PTCB) ---
152 pMC->Gsdvn("PCSE", "PCPV", 11, 1);
153 pMC->Gsdvn("PCMO", "PCSE", 13, 3);
154 pMC->Gsdvn("PCST", "PCMO", 8, 1);
155 pMC->Gsdvn("PCCE", "PCST", 8, 3);
156
157// --- Position various PHOS units in ALICE setup ---
158// --- PHOS itself first ---
159 r = FOC_R+FOC_Y/2.;
160 pphi = TMath::ATan(FOC_X/(2.*FOC_R));
161 xp1 = -r * TMath::Sin(pphi * 3.);
162 yp1 = -r * TMath::Cos(pphi * 3.);
163 xp2 = -r * TMath::Sin(pphi);
164 yp2 = -r * TMath::Cos(pphi);
165 pphi *= kRaddeg;
166 AliMatrix(idrotm[0], 90.,-3*pphi, 90., 90-3*pphi, 0., 0.);
167 AliMatrix(idrotm[1], 90., -pphi, 90., 90-pphi, 0., 0.);
168 AliMatrix(idrotm[2], 90., pphi, 90., 90+pphi, 0., 0.);
169 AliMatrix(idrotm[3], 90., 3*pphi, 90., 90+3*pphi, 0., 0.);
170 pMC->Gspos("PHOS", 1, "ALIC", xp1, yp1, 0., idrotm[0], "ONLY");
171 pMC->Gspos("PHOS", 2, "ALIC", xp2, yp2, 0., idrotm[1], "ONLY");
172 pMC->Gspos("PHOS", 3, "ALIC",-xp2, yp2, 0., idrotm[2], "ONLY");
173 pMC->Gspos("PHOS", 4, "ALIC",-xp1, yp1, 0., idrotm[3], "ONLY");
174
175// --- Now position PCPV so that its plates are right on top of ---
176// --- corresponding PHOS supermodules (previously called cradles) ---
177 r = FOC_R-CPV_Y/2.;
178 pphi = TMath::ATan(FOC_X/(2.*FOC_R));
179 xp1 = -r * TMath::Sin(pphi * 3.);
180 yp1 = -r * TMath::Cos(pphi * 3.);
181 xp2 = -r * TMath::Sin(pphi);
182 yp2 = -r * TMath::Cos(pphi);
183 pMC->Gspos("PCPV", 1, "ALIC", xp1, yp1, 0., idrotm[0], "ONLY");
184 pMC->Gspos("PCPV", 2, "ALIC", xp2, yp2, 0., idrotm[1], "ONLY");
185 pMC->Gspos("PCPV", 3, "ALIC",-xp2, yp2, 0., idrotm[2], "ONLY");
186 pMC->Gspos("PCPV", 4, "ALIC",-xp1, yp1, 0., idrotm[3], "ONLY");
187
188// --- Set modules seen without tree for drawings ---
189 pMC->Gsatt("PMOD", "SEEN", -2);
190 pMC->Gsatt("PCMO", "SEEN", -2);
191}
192
193//___________________________________________
194void AliPHOSv0::CreateMaterials()
195{
196// *** DEFINITION OF AVAILABLE PHOS MATERIALS ***
197// ORIGIN : NICK VAN EIJNDHOVEN
198
199 AliMC* pMC = AliMC::GetMC();
200
201 Int_t ISXFLD = gAlice->Field()->Integ();
202 Float_t SXMGMX = gAlice->Field()->Max();
203
204// --- The PbWO4 crystals ---
205 Float_t ax[3] = { 207.19,183.85,16. };
206 Float_t zx[3] = { 82.,74.,8. };
207 Float_t wx[3] = { 1.,1.,4. };
208 Float_t dx = 8.28;
209// --- The polysterene scintillator (CH) ---
210 Float_t ap[2] = { 12.011,1.00794 };
211 Float_t zp[2] = { 6.,1. };
212 Float_t wp[2] = { 1.,1. };
213 Float_t dp = 1.032;
214// --- Tyvek (CnH2n)
215 Float_t at[2] = { 12.011,1.00794 };
216 Float_t zt[2] = { 6.,1. };
217 Float_t wt[2] = { 1.,2. };
218 Float_t dt = .331;
219// --- Polystyrene foam ---
220 Float_t af[2] = { 12.011,1.00794 };
221 Float_t zf[2] = { 6.,1. };
222 Float_t wf[2] = { 1.,1. };
223 Float_t df = .3;
224
225 Int_t *idtmed = gAlice->Idtmed();
226
227 AliMixture( 0, "PbWO4$", ax, zx, dx, -3, wx);
228 AliMixture( 1, "Polystyrene$", ap, zp, dp, -2, wp);
229 AliMaterial(2, "Al$", 26.98, 13., 2.7, 8.9, 999.);
230// --- Absorption length^ is ignored ---
231 AliMixture( 3, "Tyvek$", at, zt, dt, -2, wt);
232 AliMixture( 4, "Foam$", af, zf, df, -2, wf);
233 AliMaterial(9, "Air$", 14.61, 7.3, .001205, 30420., 67500);
234
235 AliMedium(700, "PHOS Xtal $", 0, 1, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
236 AliMedium(701, "CPV scint. $", 1, 1, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
237 AliMedium(702, "Al parts $", 2, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .001);
238 AliMedium(703, "Tyvek wrapper$", 3, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .001);
239 AliMedium(704, "Polyst. foam $", 4, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
240 AliMedium(799, "Air $", 9, 0, ISXFLD, SXMGMX, 10., 1., .1, .1, 10.);
241
242// --- Generate explicitly delta rays in aluminium parts ---
243 pMC->Gstpar(idtmed[701], "LOSS", 3.);
244 pMC->Gstpar(idtmed[701], "DRAY", 1.);
245}
246
247void AliPHOSv0::StepManager()
248{
249
250 AliMC* pMC = AliMC::GetMC();
251
252 TClonesArray &lhits = *fHits;
253 Int_t copy, i;
254 Int_t vol[5];
255 Float_t hits[4];
256 if(pMC->CurrentVol(0,copy) == fIdSens) {
257 //
258 //We are in the sensitive volume
259 for(i=0;i<4;i++) {
260 pMC->CurrentVolOff(i+1,0,copy);
261 vol[4-i]=copy;
262 }
263 pMC->CurrentVolOff(7,0,copy);
264 vol[0]=copy;
265 pMC->TrackPosition(hits);
266 hits[3]=pMC->Edep();
267 new(lhits[fNhits++]) AliPHOShit(fIshunt,gAlice->CurrentTrack(),vol,hits);
268 }
269}