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