4663d63d |
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 | $Log$ |
7cc664ae |
18 | Revision 1.1 2001/05/14 09:57:39 coppedis |
19 | A different geometry for the ZDCs |
20 | |
4663d63d |
21 | |
22 | */ |
23 | |
24 | /////////////////////////////////////////////////////////////////////////////// |
25 | // // |
26 | // Zero Degree Calorimeter // |
27 | // This class contains the basic functions for the ZDC // |
28 | // Functions specific to one particular geometry are // |
29 | // contained in the derived classes // |
30 | // // |
31 | /////////////////////////////////////////////////////////////////////////////// |
32 | |
33 | // --- Standard libraries |
34 | #include "stdio.h" |
35 | |
36 | // --- ROOT system |
37 | #include <TBRIK.h> |
38 | #include <TNode.h> |
39 | #include <TMath.h> |
40 | #include <TRandom.h> |
41 | #include <TSystem.h> |
42 | #include <TTree.h> |
43 | |
44 | |
45 | // --- AliRoot classes |
46 | #include "AliZDCv2.h" |
47 | #include "AliZDCHit.h" |
48 | #include "AliZDCDigit.h" |
49 | #include "AliRun.h" |
50 | #include "AliDetector.h" |
51 | #include "AliMagF.h" |
52 | #include "AliMC.h" |
53 | #include "AliCallf77.h" |
54 | #include "AliConst.h" |
55 | #include "AliPDG.h" |
56 | #include "TLorentzVector.h" |
57 | |
58 | |
59 | ClassImp(AliZDCv2) |
60 | |
61 | |
62 | /////////////////////////////////////////////////////////////////////////////// |
63 | // // |
64 | // Zero Degree Calorimeter version 2 // |
65 | // // |
66 | /////////////////////////////////////////////////////////////////////////////// |
67 | |
68 | //_____________________________________________________________________________ |
69 | AliZDCv2::AliZDCv2() : AliZDC() |
70 | { |
71 | // |
72 | // Default constructor for Zero Degree Calorimeter |
73 | // |
74 | |
75 | fMedSensF1 = 0; |
76 | fMedSensF2 = 0; |
77 | fMedSensZN = 0; |
78 | fMedSensZP = 0; |
79 | fMedSensZEM = 0; |
80 | fMedSensGR = 0; |
81 | // fMedSensPI = 0; |
82 | // fMedSensTDI = 0; |
83 | } |
84 | |
85 | //_____________________________________________________________________________ |
86 | AliZDCv2::AliZDCv2(const char *name, const char *title) |
87 | : AliZDC(name,title) |
88 | { |
89 | // |
90 | // Standard constructor for Zero Degree Calorimeter |
91 | // |
92 | // |
93 | // Check that DIPO, ABSO, DIPO and SHIL is there (otherwise tracking is wrong!!!) |
94 | |
95 | AliModule* PIPE=gAlice->GetModule("PIPE"); |
96 | AliModule* ABSO=gAlice->GetModule("ABSO"); |
97 | AliModule* DIPO=gAlice->GetModule("DIPO"); |
98 | AliModule* SHIL=gAlice->GetModule("SHIL"); |
99 | if((!PIPE) || (!ABSO) || (!DIPO) || (!SHIL)) { |
100 | Error("Constructor","ZDC needs PIPE, ABSO, DIPO and SHIL!!!\n"); |
101 | exit(1); |
102 | } |
103 | |
104 | fMedSensF1 = 0; |
105 | fMedSensF2 = 0; |
106 | fMedSensZN = 0; |
107 | fMedSensZP = 0; |
108 | fMedSensZEM = 0; |
109 | fMedSensGR = 0; |
110 | // fMedSensPI = 0; |
111 | // fMedSensTDI = 0; |
112 | |
113 | |
114 | // Parameters for light tables |
115 | fNalfan = 90; // Number of Alfa (neutrons) |
116 | fNalfap = 90; // Number of Alfa (protons) |
117 | fNben = 18; // Number of beta (neutrons) |
118 | fNbep = 28; // Number of beta (protons) |
119 | Int_t ip,jp,kp; |
120 | for(ip=0; ip<4; ip++){ |
121 | for(kp=0; kp<fNalfap; kp++){ |
122 | for(jp=0; jp<fNbep; jp++){ |
123 | fTablep[ip][kp][jp] = 0; |
124 | } |
125 | } |
126 | } |
127 | Int_t in,jn,kn; |
128 | for(in=0; in<4; in++){ |
129 | for(kn=0; kn<fNalfan; kn++){ |
130 | for(jn=0; jn<fNben; jn++){ |
131 | fTablen[in][kn][jn] = 0; |
132 | } |
133 | } |
134 | } |
135 | |
136 | // Parameters for hadronic calorimeters geometry |
137 | fDimZP[0] = 11.2; |
138 | fDimZP[1] = 6.; |
139 | fDimZP[2] = 75.; |
140 | fPosZN[0] = 0.; |
141 | fPosZN[1] = -1.2; |
142 | fPosZN[2] = 11650.; |
143 | fPosZP[0] = -24.; |
144 | fPosZP[1] = 0.; |
145 | fPosZP[2] = 11600.; |
146 | fFibZN[0] = 0.; |
147 | fFibZN[1] = 0.01825; |
148 | fFibZN[2] = 50.; |
149 | fFibZP[0] = 0.; |
150 | fFibZP[1] = 0.0275; |
151 | fFibZP[2] = 75.; |
152 | |
153 | // Parameters for EM calorimeter geometry |
154 | fPosZEM[0] = 8.5; |
155 | fPosZEM[1] = 0.; |
156 | fPosZEM[2] = -1000.; |
157 | |
158 | |
159 | fDigits = new TClonesArray("AliZDCDigit",1000); |
160 | } |
161 | |
162 | //_____________________________________________________________________________ |
163 | void AliZDCv2::CreateGeometry() |
164 | { |
165 | // |
166 | // Create the geometry for the Zero Degree Calorimeter version 1 |
167 | //* Initialize COMMON block ZDC_CGEOM |
168 | //* |
169 | |
170 | CreateBeamLine(); |
171 | CreateZDC(); |
172 | } |
173 | |
174 | //_____________________________________________________________________________ |
175 | void AliZDCv2::CreateBeamLine() |
176 | { |
177 | |
178 | Float_t zq, zd1, zd2; |
179 | Float_t conpar[9], tubpar[3], tubspar[5], boxpar[3]; |
180 | Int_t im1, im2; |
181 | |
182 | Int_t *idtmed = fIdtmed->GetArray(); |
183 | |
184 | // -- Mother of the ZDCs (Vacuum PCON) |
185 | |
186 | conpar[0] = 0.; |
187 | conpar[1] = 360.; |
188 | conpar[2] = 2.; |
189 | conpar[3] = -1100.; |
190 | conpar[4] = 0.; |
191 | conpar[5] = 155.; |
192 | conpar[6] = 13060.; |
193 | conpar[7] = 0.; |
194 | conpar[8] = 155.; |
195 | gMC->Gsvolu("ZDC ", "PCON", idtmed[11], conpar, 9); |
196 | gMC->Gspos("ZDC ", 1, "ALIC", 0., 0., 0., 0, "ONLY"); |
197 | |
198 | // -- FIRST SECTION OF THE BEAM PIPE (from compensator dipole to |
199 | // the beginning of D1) |
200 | |
201 | zd1 = 2000.; |
202 | |
203 | tubpar[0] = 6.3/2.; |
204 | tubpar[1] = 6.7/2.; |
205 | tubpar[2] = 3838.3/2.; |
206 | gMC->Gsvolu("QT01", "TUBE", idtmed[7], tubpar, 3); |
207 | gMC->Gspos("QT01", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
208 | |
209 | //-- SECOND SECTION OF THE BEAM PIPE (from the end of D1 to the |
210 | // beginning of D2) |
211 | |
212 | //-- FROM MAGNETIC BEGINNING OF D1 TO MAGNETIC END OF D1 + 13.5 cm |
213 | //-- Cylindrical pipe (r = 3.47) + conical flare |
214 | |
215 | // -> Beginning of D1 |
216 | zd1 += 2.*tubpar[2]; |
217 | |
218 | tubpar[0] = 3.47; |
219 | tubpar[1] = 3.47+0.2; |
220 | tubpar[2] = 958.5/2.; |
221 | gMC->Gsvolu("QT02", "TUBE", idtmed[7], tubpar, 3); |
222 | gMC->Gspos("QT02", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
223 | |
224 | zd1 += 2.*tubpar[2]; |
225 | |
226 | conpar[0] = 25./2.; |
227 | conpar[1] = 6.44/2.; |
228 | conpar[2] = 6.84/2.; |
229 | conpar[3] = 10./2.; |
230 | conpar[4] = 10.4/2.; |
231 | gMC->Gsvolu("QC01", "CONE", idtmed[7], conpar, 5); |
232 | gMC->Gspos("QC01", 1, "ZDC ", 0., 0., conpar[0] + zd1, 0, "ONLY"); |
233 | |
234 | zd1 += 2.*conpar[0]; |
235 | |
236 | tubpar[0] = 10./2.; |
237 | tubpar[1] = 10.4/2.; |
238 | tubpar[2] = 50./2.; |
239 | gMC->Gsvolu("QT03", "TUBE", idtmed[7], tubpar, 3); |
240 | gMC->Gspos("QT03", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
241 | |
242 | zd1 += tubpar[2]*2.; |
243 | |
244 | tubpar[0] = 10./2.; |
245 | tubpar[1] = 10.4/2.; |
246 | tubpar[2] = 10./2.; |
247 | gMC->Gsvolu("QT04", "TUBE", idtmed[7], tubpar, 3); |
248 | gMC->Gspos("QT04", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
249 | |
250 | zd1 += tubpar[2] * 2.; |
251 | |
252 | tubpar[0] = 10./2.; |
253 | tubpar[1] = 10.4/2.; |
254 | tubpar[2] = 3.16/2.; |
255 | gMC->Gsvolu("QT05", "TUBE", idtmed[7], tubpar, 3); |
256 | gMC->Gspos("QT05", 1, "ZDC ", 0., 0., tubpar[0] + zd1, 0, "ONLY"); |
257 | |
258 | zd1 += tubpar[2] * 2.; |
259 | |
260 | tubpar[0] = 10.0/2.; |
261 | tubpar[1] = 10.4/2; |
262 | tubpar[2] = 190./2.; |
263 | gMC->Gsvolu("QT06", "TUBE", idtmed[7], tubpar, 3); |
264 | gMC->Gspos("QT06", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
265 | |
266 | zd1 += tubpar[2] * 2.; |
267 | |
268 | conpar[0] = 30./2.; |
269 | conpar[1] = 10./2.; |
270 | conpar[2] = 10.4/2.; |
271 | conpar[3] = 20.6/2.; |
272 | conpar[4] = 21./2.; |
273 | gMC->Gsvolu("QC02", "CONE", idtmed[7], conpar, 5); |
274 | gMC->Gspos("QC02", 1, "ZDC ", 0., 0., conpar[0] + zd1, 0, "ONLY"); |
275 | |
276 | zd1 += conpar[0] * 2.; |
277 | |
278 | tubpar[0] = 20.6/2.; |
279 | tubpar[1] = 21./2.; |
280 | tubpar[2] = 450./2.; |
281 | gMC->Gsvolu("QT07", "TUBE", idtmed[7], tubpar, 3); |
282 | gMC->Gspos("QT07", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
283 | |
284 | zd1 += tubpar[2] * 2.; |
285 | |
286 | conpar[0] = 13.6/2.; |
287 | conpar[1] = 20.6/2.; |
288 | conpar[2] = 21./2.; |
289 | conpar[3] = 25.4/2.; |
290 | conpar[4] = 25.8/2.; |
291 | gMC->Gsvolu("QC03", "CONE", idtmed[7], conpar, 5); |
292 | gMC->Gspos("QC03", 1, "ZDC ", 0., 0., conpar[0] + zd1, 0, "ONLY"); |
293 | |
294 | zd1 += conpar[0] * 2.; |
295 | |
296 | tubpar[0] = 25.4/2.; |
297 | tubpar[1] = 25.8/2.; |
298 | tubpar[2] = 205.8/2.; |
299 | gMC->Gsvolu("QT08", "TUBE", idtmed[7], tubpar, 3); |
300 | gMC->Gspos("QT08", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
301 | |
302 | zd1 += tubpar[2] * 2.; |
303 | |
304 | tubpar[0] = 50./2.; |
305 | tubpar[1] = 50.4/2.; |
306 | // QT09 is 10 cm longer to accomodate TDI |
307 | tubpar[2] = 515.4/2.; |
308 | gMC->Gsvolu("QT09", "TUBE", idtmed[7], tubpar, 3); |
309 | gMC->Gspos("QT09", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
310 | |
311 | // --- Insert TDI (inside ZDC volume) |
312 | |
313 | boxpar[0] = 5.6; |
314 | boxpar[1] = 5.6; |
315 | boxpar[2] = 400./2.; |
316 | gMC->Gsvolu("QTD1", "BOX ", idtmed[7], boxpar, 3); |
317 | gMC->Gspos("QTD1", 1, "ZDC ", 0., 10.6, tubpar[2] + zd1 + 56.3, 0, "ONLY"); |
318 | gMC->Gspos("QTD1", 2, "ZDC ", 0., -10.6, tubpar[2] + zd1 + 56.3, 0, "ONLY"); |
319 | |
320 | boxpar[0] = 0.2/2.; |
321 | boxpar[1] = 5.6; |
322 | boxpar[2] = 400./2.; |
323 | gMC->Gsvolu("QTD2", "BOX ", idtmed[6], boxpar, 3); |
324 | gMC->Gspos("QTD2", 1, "ZDC ", 5.6+boxpar[0], 0., tubpar[2] + zd1 + 56.3, 0, "ONLY"); |
325 | |
326 | tubspar[0] = 6.2; |
327 | tubspar[1] = 6.4; |
328 | tubspar[2] = 400./2.; |
329 | tubspar[3] = 180.-62.5; |
330 | tubspar[4] = 180.+62.5; |
331 | gMC->Gsvolu("QTD3", "TUBS", idtmed[6], tubspar, 5); |
332 | gMC->Gspos("QTD3", 1, "ZDC ", -3., 0., tubpar[2] + zd1 + 56.3, 0, "ONLY"); |
333 | |
334 | zd1 += tubpar[2] * 2.; |
335 | |
336 | tubpar[0] = 50./2.; |
337 | tubpar[1] = 50.4/2.; |
338 | // QT10 is 10 cm shorter |
339 | tubpar[2] = 690./2.; |
340 | gMC->Gsvolu("QT10", "TUBE", idtmed[7], tubpar, 3); |
341 | gMC->Gspos("QT10", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
342 | |
343 | zd1 += tubpar[2] * 2.; |
344 | |
345 | tubpar[0] = 50./2.; |
346 | tubpar[1] = 50.4/2.; |
347 | tubpar[2] = 778.5/2.; |
348 | gMC->Gsvolu("QT11", "TUBE", idtmed[7], tubpar, 3); |
349 | gMC->Gspos("QT11", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
350 | |
351 | zd1 += tubpar[2] * 2.; |
352 | |
353 | conpar[0] = 14.18/2.; |
354 | conpar[1] = 50./2.; |
355 | conpar[2] = 50.4/2.; |
356 | conpar[3] = 55./2.; |
357 | conpar[4] = 55.4/2.; |
358 | gMC->Gsvolu("QC04", "CONE", idtmed[7], conpar, 5); |
359 | gMC->Gspos("QC04", 1, "ZDC ", 0., 0., conpar[0] + zd1, 0, "ONLY"); |
360 | |
361 | zd1 += conpar[0] * 2.; |
362 | |
363 | tubpar[0] = 55./2.; |
364 | tubpar[1] = 55.4/2.; |
365 | tubpar[2] = 730./2.; |
366 | gMC->Gsvolu("QT12", "TUBE", idtmed[7], tubpar, 3); |
367 | gMC->Gspos("QT12", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
368 | |
369 | zd1 += tubpar[2] * 2.; |
370 | |
371 | conpar[0] = 36.86/2.; |
372 | conpar[1] = 55./2.; |
373 | conpar[2] = 55.4/2.; |
374 | conpar[3] = 68./2.; |
375 | conpar[4] = 68.4/2.; |
376 | gMC->Gsvolu("QC05", "CONE", idtmed[7], conpar, 5); |
377 | gMC->Gspos("QC05", 1, "ZDC ", 0., 0., conpar[0] + zd1, 0, "ONLY"); |
378 | |
379 | zd1 += conpar[0] * 2.; |
380 | |
381 | tubpar[0] = 68./2.; |
382 | tubpar[1] = 68.4/2.; |
383 | tubpar[2] = 927.3/2.; |
384 | gMC->Gsvolu("QT13", "TUBE", idtmed[7], tubpar, 3); |
385 | gMC->Gspos("QT13", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
386 | |
387 | zd1 += tubpar[2] * 2.; |
388 | |
389 | tubpar[0] = 0./2.; |
390 | tubpar[1] = 68.4/2.; |
391 | tubpar[2] = 0.2/2.; |
392 | gMC->Gsvolu("QT14", "TUBE", idtmed[8], tubpar, 3); |
393 | gMC->Gspos("QT14", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
394 | |
395 | zd1 += tubpar[2] * 2.; |
396 | |
397 | tubpar[0] = 0./2.; |
398 | tubpar[1] = 6.4/2.; |
399 | tubpar[2] = 0.2/2.; |
400 | gMC->Gsvolu("QT15", "TUBE", idtmed[11], tubpar, 3); |
401 | |
402 | //-- Position QT15 inside QT14 |
403 | gMC->Gspos("QT15", 1, "QT14", -7.7, 0., 0., 0, "ONLY"); |
404 | |
405 | tubpar[0] = 0./2.; |
406 | tubpar[1] = 6.4/2.; |
407 | tubpar[2] = 0.2/2.; |
408 | gMC->Gsvolu("QT16", "TUBE", idtmed[11], tubpar, 3); |
409 | |
410 | //-- Position QT16 inside QT14 |
411 | gMC->Gspos("QT16", 1, "QT14", 7.7, 0., 0., 0, "ONLY"); |
412 | |
413 | |
414 | //-- BEAM PIPE BETWEEN END OF CONICAL PIPE AND BEGINNING OF D2 |
415 | |
416 | tubpar[0] = 6.4/2.; |
417 | tubpar[1] = 6.8/2.; |
418 | tubpar[2] = 680.8/2.; |
419 | gMC->Gsvolu("QT17", "TUBE", idtmed[7], tubpar, 3); |
420 | |
421 | tubpar[0] = 6.4/2.; |
422 | tubpar[1] = 6.8/2.; |
423 | tubpar[2] = 680.8/2.; |
424 | gMC->Gsvolu("QT18", "TUBE", idtmed[7], tubpar, 3); |
425 | |
426 | // -- ROTATE PIPES |
427 | |
428 | Float_t angle = 0.143*kDegrad; |
429 | |
430 | AliMatrix(im1, 90.-0.143, 0., 90., 90., 0.143, 180.); |
431 | gMC->Gspos("QT17", 1, "ZDC ", TMath::Sin(angle) * 680.8/ 2. - 9.4, |
432 | 0., tubpar[2] + zd1, im1, "ONLY"); |
433 | |
434 | AliMatrix(im2, 90.+0.143, 0., 90., 90., 0.143, 0.); |
435 | gMC->Gspos("QT18", 1, "ZDC ", 9.7 - TMath::Sin(angle) * 680.8 / 2., |
436 | 0., tubpar[2] + zd1, im2, "ONLY"); |
437 | |
438 | // -- BEAM PIPE ON THE OTHER SIDE OF I.P. TILL THE EM ZDC |
439 | |
440 | Float_t zb = -800.; // End of QBPM (from AliPIPEv0.cxx) |
441 | tubpar[0] = 8.0/2.; |
442 | tubpar[1] = 8.2/2.; |
7cc664ae |
443 | tubpar[2] = (1050+zb)/2.; // From the end of QBPM to z=1050. |
4663d63d |
444 | gMC->Gsvolu("QT19", "TUBE", idtmed[7], tubpar, 3); |
445 | gMC->Gspos("QT19", 1, "ZDC ", 0., 0., zb - tubpar[2], 0, "ONLY"); |
446 | |
447 | |
448 | // -- END OF BEAM PIPE VOLUME DEFINITION. |
449 | // ---------------------------------------------------------------- |
450 | |
451 | // -- MAGNET DEFINITION -> LHC OPTICS 6.2 (preliminary version) |
452 | |
453 | // ---------------------------------------------------------------- |
454 | // Replaced by the muon dipole |
455 | // ---------------------------------------------------------------- |
456 | // -- COMPENSATOR DIPOLE (MBXW) |
457 | // GAP (VACUUM WITH MAGNETIC FIELD) |
458 | |
459 | // tubpar[0] = 0.; |
460 | // tubpar[1] = 4.5; |
461 | // tubpar[2] = 340./2.; |
462 | // gMC->Gsvolu("MBXW", "TUBE", idtmed[11], tubpar, 3); |
463 | // gMC->Gspos("MBXW", 1, "ZDC ", 0., 0., tubpar[2] + 805., 0, "ONLY"); |
464 | |
465 | // -- YOKE (IRON WITHOUT MAGNETIC FIELD) |
466 | |
467 | // tubpar[0] = 4.5; |
468 | // tubpar[1] = 55.; |
469 | // tubpar[2] = 340./2.; |
470 | // gMC->Gsvolu("YMBX", "TUBE", idtmed[7], tubpar, 3); |
471 | // gMC->Gspos("YMBX", 1, "ZDC ", 0., 0., tubpar[2] + 805., 0, "ONLY"); |
472 | |
473 | // ---------------------------------------------------------------- |
474 | // Replaced by the second dipole |
475 | // ---------------------------------------------------------------- |
476 | // -- COMPENSATOR DIPOLE (MCBWA) |
477 | // GAP (VACUUM WITH MAGNETIC FIELD) |
478 | |
479 | // tubpar[0] = 0.; |
480 | // tubpar[1] = 4.5; |
481 | // tubpar[2] = 170./2.; |
482 | // gMC->Gsvolu("MCBW", "TUBE", idtmed[11], tubpar, 3); |
483 | // gMC->Gspos("MCBW", 1, "ZDC ", 0., 0., tubpar[2] + 1921.6, 0, "ONLY"); |
484 | |
485 | // -- YOKE (IRON WITHOUT MAGNETIC FIELD) |
486 | |
487 | // tubpar[0] = 4.5; |
488 | // tubpar[1] = 55.; |
489 | // tubpar[2] = 170./2.; |
490 | // gMC->Gsvolu("YMCB", "TUBE", idtmed[7], tubpar, 3); |
491 | // gMC->Gspos("YMCB", 1, "ZDC ", 0., 0., tubpar[2] + 1921.6, 0, "ONLY"); |
492 | |
493 | // -- INNER TRIPLET |
494 | |
495 | zq = 2296.5; |
496 | |
497 | // -- DEFINE MQXL AND MQX QUADRUPOLE ELEMENT |
498 | |
499 | // MQXL |
500 | // -- GAP (VACUUM WITH MAGNETIC FIELD) |
501 | |
502 | tubpar[0] = 0.; |
503 | tubpar[1] = 3.5; |
504 | tubpar[2] = 637./2.; |
505 | gMC->Gsvolu("MQXL", "TUBE", idtmed[11], tubpar, 3); |
506 | |
507 | // -- YOKE |
508 | |
509 | tubpar[0] = 3.5; |
510 | tubpar[1] = 22.; |
511 | tubpar[2] = 637./2.; |
512 | gMC->Gsvolu("YMQL", "TUBE", idtmed[7], tubpar, 3); |
513 | |
514 | gMC->Gspos("MQXL", 1, "ZDC ", 0., 0., tubpar[2] + zq, 0, "ONLY"); |
515 | gMC->Gspos("YMQL", 1, "ZDC ", 0., 0., tubpar[2] + zq, 0, "ONLY"); |
516 | |
517 | gMC->Gspos("MQXL", 2, "ZDC ", 0., 0., tubpar[2] + zq + 2430., 0, "ONLY"); |
518 | gMC->Gspos("YMQL", 2, "ZDC ", 0., 0., tubpar[2] + zq + 2430., 0, "ONLY"); |
519 | |
520 | // -- MQX |
521 | // -- GAP (VACUUM WITH MAGNETIC FIELD) |
522 | |
523 | tubpar[0] = 0.; |
524 | tubpar[1] = 3.5; |
525 | tubpar[2] = 550./2.; |
526 | gMC->Gsvolu("MQX ", "TUBE", idtmed[11], tubpar, 3); |
527 | |
528 | // -- YOKE |
529 | |
530 | tubpar[0] = 3.5; |
531 | tubpar[1] = 22.; |
532 | tubpar[2] = 550./2.; |
533 | gMC->Gsvolu("YMQ ", "TUBE", idtmed[7], tubpar, 3); |
534 | |
535 | gMC->Gspos("MQX ", 1, "ZDC ", 0., 0., tubpar[2] + zq + 883.5, 0, "ONLY"); |
536 | gMC->Gspos("YMQ ", 1, "ZDC ", 0., 0., tubpar[2] + zq + 883.5, 0, "ONLY"); |
537 | |
538 | gMC->Gspos("MQX ", 2, "ZDC ", 0., 0., tubpar[2] + zq + 1533.5, 0, "ONLY"); |
539 | gMC->Gspos("YMQ ", 2, "ZDC ", 0., 0., tubpar[2] + zq + 1533.5, 0, "ONLY"); |
540 | |
541 | // -- SEPARATOR DIPOLE D1 |
542 | |
543 | zd1 = 5838.3; |
544 | |
545 | // -- GAP (VACUUM WITH MAGNETIC FIELD) |
546 | |
547 | tubpar[0] = 0.; |
548 | tubpar[1] = 6.94/2.; |
549 | tubpar[2] = 945./2.; |
550 | gMC->Gsvolu("MD1 ", "TUBE", idtmed[11], tubpar, 3); |
551 | |
552 | // -- Insert horizontal Cu plates inside D1 |
553 | // -- (to simulate the vacuum chamber) |
554 | |
555 | boxpar[0] = TMath::Sqrt(tubpar[1]*tubpar[1]-(2.98+0.2)*(2.98+0.2)); |
556 | boxpar[1] = 0.2/2.; |
557 | boxpar[2] =945./2.; |
558 | gMC->Gsvolu("MD1V", "BOX ", idtmed[6], boxpar, 3); |
559 | gMC->Gspos("MD1V", 1, "MD1 ", 0., 2.98+boxpar[1], 0., 0, "ONLY"); |
560 | gMC->Gspos("MD1V", 2, "MD1 ", 0., -2.98-boxpar[1], 0., 0, "ONLY"); |
561 | |
562 | // -- YOKE |
563 | |
564 | tubpar[0] = 0.; |
565 | tubpar[1] = 110./2; |
566 | tubpar[2] = 945./2.; |
567 | gMC->Gsvolu("YD1 ", "TUBE", idtmed[7], tubpar, 3); |
568 | |
569 | gMC->Gspos("YD1 ", 1, "ZDC ", 0., 0., tubpar[2] + zd1, 0, "ONLY"); |
570 | gMC->Gspos("MD1 ", 1, "YD1 ", 0., 0., 0., 0, "ONLY"); |
571 | |
572 | // -- DIPOLE D2 |
573 | |
574 | zd2 = 12147.6; |
575 | |
576 | // -- GAP (VACUUM WITH MAGNETIC FIELD) |
577 | |
578 | tubpar[0] = 0.; |
579 | tubpar[1] = 7.5/2.; |
580 | tubpar[2] = 945./2.; |
581 | gMC->Gsvolu("MD2 ", "TUBE", idtmed[11], tubpar, 3); |
582 | |
583 | // -- YOKE |
584 | |
585 | tubpar[0] = 0.; |
586 | tubpar[1] = 55.; |
587 | tubpar[2] = 945./2.; |
588 | gMC->Gsvolu("YD2 ", "TUBE", idtmed[7], tubpar, 3); |
589 | |
590 | gMC->Gspos("YD2 ", 1, "ZDC ", 0., 0., tubpar[2] + zd2, 0, "ONLY"); |
591 | |
592 | gMC->Gspos("MD2 ", 1, "YD2 ", -9.4, 0., 0., 0, "ONLY"); |
593 | gMC->Gspos("MD2 ", 2, "YD2 ", 9.4, 0., 0., 0, "ONLY"); |
594 | |
595 | // -- END OF MAGNET DEFINITION |
596 | } |
597 | |
598 | //_____________________________________________________________________________ |
599 | void AliZDCv2::CreateZDC() |
600 | { |
601 | |
602 | Float_t DimPb[6], DimVoid[6]; |
603 | |
604 | Int_t *idtmed = fIdtmed->GetArray(); |
605 | |
606 | // Parameters for hadronic calorimeters geometry |
607 | // NB -> parameters used ONLY in CreateZDC() |
608 | Float_t fDimZN[3] = {3.52, 3.52, 50.}; // Dimensions of neutron detector |
609 | Float_t fGrvZN[3] = {0.03, 0.03, 50.}; // Grooves for neutron detector |
610 | Float_t fGrvZP[3] = {0.04, 0.04, 75.}; // Grooves for proton detector |
611 | Int_t fDivZN[3] = {11, 11, 0}; // Division for neutron detector |
612 | Int_t fDivZP[3] = {7, 15, 0}; // Division for proton detector |
613 | Int_t fTowZN[2] = {2, 2}; // Tower for neutron detector |
614 | Int_t fTowZP[2] = {4, 1}; // Tower for proton detector |
615 | |
616 | // Parameters for EM calorimeter geometry |
617 | // NB -> parameters used ONLY in CreateZDC() |
618 | Float_t fDimZEMPb = 0.15*(TMath::Sqrt(2.)); // z-dimension of the Pb slice |
619 | Float_t fDimZEMAir = 0.001; // scotch |
620 | Float_t fFibRadZEM = 0.0315; // External fiber radius (including cladding) |
621 | Int_t fDivZEM[3] = {92, 0, 20}; // Divisions for EM detector |
622 | Float_t fDimZEM0 = 2*fDivZEM[2]*(fDimZEMPb+fDimZEMAir+fFibRadZEM*(TMath::Sqrt(2.))); |
623 | Float_t fDimZEM[6] = {fDimZEM0, 3.5, 3.5, 45., 0., 0.}; // Dimensions of EM detector |
624 | Float_t fFibZEM2 = fDimZEM[2]/TMath::Sin(fDimZEM[3]*kDegrad)-fFibRadZEM; |
625 | Float_t fFibZEM[3] = {0., 0.0275, fFibZEM2}; // Fibers for EM calorimeter |
626 | |
627 | |
628 | //-- Create calorimeters geometry |
629 | |
630 | // ------------------------------------------------------------------------------- |
631 | //--> Neutron calorimeter (ZN) |
632 | |
633 | gMC->Gsvolu("ZNEU", "BOX ", idtmed[1], fDimZN, 3); // Passive material |
634 | gMC->Gsvolu("ZNF1", "TUBE", idtmed[3], fFibZN, 3); // Active material |
635 | gMC->Gsvolu("ZNF2", "TUBE", idtmed[4], fFibZN, 3); |
636 | gMC->Gsvolu("ZNF3", "TUBE", idtmed[4], fFibZN, 3); |
637 | gMC->Gsvolu("ZNF4", "TUBE", idtmed[3], fFibZN, 3); |
638 | gMC->Gsvolu("ZNG1", "BOX ", idtmed[12], fGrvZN, 3); // Empty grooves |
639 | gMC->Gsvolu("ZNG2", "BOX ", idtmed[12], fGrvZN, 3); |
640 | gMC->Gsvolu("ZNG3", "BOX ", idtmed[12], fGrvZN, 3); |
641 | gMC->Gsvolu("ZNG4", "BOX ", idtmed[12], fGrvZN, 3); |
642 | |
643 | // Divide ZNEU in towers (for hits purposes) |
644 | |
645 | gMC->Gsdvn("ZNTX", "ZNEU", fTowZN[0], 1); // x-tower |
646 | gMC->Gsdvn("ZN1 ", "ZNTX", fTowZN[1], 2); // y-tower |
647 | |
648 | //-- Divide ZN1 in minitowers |
649 | // fDivZN[0]= NUMBER OF FIBERS PER TOWER ALONG X-AXIS, |
650 | // fDivZN[1]= NUMBER OF FIBERS PER TOWER ALONG Y-AXIS |
651 | // (4 fibres per minitower) |
652 | |
653 | gMC->Gsdvn("ZNSL", "ZN1 ", fDivZN[1], 2); // Slices |
654 | gMC->Gsdvn("ZNST", "ZNSL", fDivZN[0], 1); // Sticks |
655 | |
656 | // --- Position the empty grooves in the sticks (4 grooves per stick) |
657 | Float_t dx = fDimZN[0] / fDivZN[0] / 4.; |
658 | Float_t dy = fDimZN[1] / fDivZN[1] / 4.; |
659 | |
660 | gMC->Gspos("ZNG1", 1, "ZNST", 0.-dx, 0.+dy, 0., 0, "ONLY"); |
661 | gMC->Gspos("ZNG2", 1, "ZNST", 0.+dx, 0.+dy, 0., 0, "ONLY"); |
662 | gMC->Gspos("ZNG3", 1, "ZNST", 0.-dx, 0.-dy, 0., 0, "ONLY"); |
663 | gMC->Gspos("ZNG4", 1, "ZNST", 0.+dx, 0.-dy, 0., 0, "ONLY"); |
664 | |
665 | // --- Position the fibers in the grooves |
666 | gMC->Gspos("ZNF1", 1, "ZNG1", 0., 0., 0., 0, "ONLY"); |
667 | gMC->Gspos("ZNF2", 1, "ZNG2", 0., 0., 0., 0, "ONLY"); |
668 | gMC->Gspos("ZNF3", 1, "ZNG3", 0., 0., 0., 0, "ONLY"); |
669 | gMC->Gspos("ZNF4", 1, "ZNG4", 0., 0., 0., 0, "ONLY"); |
670 | |
671 | // --- Position the neutron calorimeter in ZDC |
672 | gMC->Gspos("ZNEU", 1, "ZDC ", fPosZN[0], fPosZN[1], fPosZN[2] + fDimZN[2], 0, "ONLY"); |
673 | |
674 | |
675 | // ------------------------------------------------------------------------------- |
676 | //--> Proton calorimeter (ZP) |
677 | |
678 | gMC->Gsvolu("ZPRO", "BOX ", idtmed[2], fDimZP, 3); // Passive material |
679 | gMC->Gsvolu("ZPF1", "TUBE", idtmed[3], fFibZP, 3); // Active material |
680 | gMC->Gsvolu("ZPF2", "TUBE", idtmed[4], fFibZP, 3); |
681 | gMC->Gsvolu("ZPF3", "TUBE", idtmed[4], fFibZP, 3); |
682 | gMC->Gsvolu("ZPF4", "TUBE", idtmed[3], fFibZP, 3); |
683 | gMC->Gsvolu("ZPG1", "BOX ", idtmed[12], fGrvZP, 3); // Empty grooves |
684 | gMC->Gsvolu("ZPG2", "BOX ", idtmed[12], fGrvZP, 3); |
685 | gMC->Gsvolu("ZPG3", "BOX ", idtmed[12], fGrvZP, 3); |
686 | gMC->Gsvolu("ZPG4", "BOX ", idtmed[12], fGrvZP, 3); |
687 | |
688 | //-- Divide ZPRO in towers(for hits purposes) |
689 | |
690 | gMC->Gsdvn("ZPTX", "ZPRO", fTowZP[0], 1); // x-tower |
691 | gMC->Gsdvn("ZP1 ", "ZPTX", fTowZP[1], 2); // y-tower |
692 | |
693 | |
694 | //-- Divide ZP1 in minitowers |
695 | // fDivZP[0]= NUMBER OF FIBERS ALONG X-AXIS PER MINITOWER, |
696 | // fDivZP[1]= NUMBER OF FIBERS ALONG Y-AXIS PER MINITOWER |
697 | // (4 fiber per minitower) |
698 | |
699 | gMC->Gsdvn("ZPSL", "ZP1 ", fDivZP[1], 2); // Slices |
700 | gMC->Gsdvn("ZPST", "ZPSL", fDivZP[0], 1); // Sticks |
701 | |
702 | // --- Position the empty grooves in the sticks (4 grooves per stick) |
703 | dx = fDimZP[0] / fTowZP[0] / fDivZP[0] / 2.; |
704 | dy = fDimZP[1] / fTowZP[1] / fDivZP[1] / 2.; |
705 | |
706 | gMC->Gspos("ZPG1", 1, "ZPST", 0.-dx, 0.+dy, 0., 0, "ONLY"); |
707 | gMC->Gspos("ZPG2", 1, "ZPST", 0.+dx, 0.+dy, 0., 0, "ONLY"); |
708 | gMC->Gspos("ZPG3", 1, "ZPST", 0.-dx, 0.-dy, 0., 0, "ONLY"); |
709 | gMC->Gspos("ZPG4", 1, "ZPST", 0.+dx, 0.-dy, 0., 0, "ONLY"); |
710 | |
711 | // --- Position the fibers in the grooves |
712 | gMC->Gspos("ZPF1", 1, "ZPG1", 0., 0., 0., 0, "ONLY"); |
713 | gMC->Gspos("ZPF2", 1, "ZPG2", 0., 0., 0., 0, "ONLY"); |
714 | gMC->Gspos("ZPF3", 1, "ZPG3", 0., 0., 0., 0, "ONLY"); |
715 | gMC->Gspos("ZPF4", 1, "ZPG4", 0., 0., 0., 0, "ONLY"); |
716 | |
717 | |
718 | // --- Position the proton calorimeter in ZDC |
719 | gMC->Gspos("ZPRO", 1, "ZDC ", fPosZP[0], fPosZP[1], fPosZP[2] + fDimZP[2], 0, "ONLY"); |
720 | |
721 | |
722 | // ------------------------------------------------------------------------------- |
723 | // -> EM calorimeter (ZEM) |
724 | |
725 | gMC->Gsvolu("ZEM ", "PARA", idtmed[10], fDimZEM, 6); |
726 | |
727 | Int_t irot1, irot2; |
728 | |
729 | gMC->Matrix(irot1,180.,0.,90.,90.,90.,0.); // Rotation matrix 1 |
730 | gMC->Matrix(irot2,180.,0.,90.,fDimZEM[3]+90.,90.,fDimZEM[3]);// Rotation matrix 2 |
731 | // printf("irot1 = %d, irot2 = %d \n", irot1, irot2); |
732 | |
733 | gMC->Gsvolu("ZEMF", "TUBE", idtmed[3], fFibZEM, 3); // Active material |
734 | |
735 | gMC->Gsdvn("ZETR", "ZEM ", fDivZEM[2], 1); // Tranches |
736 | |
737 | DimPb[0] = fDimZEMPb; // Lead slices |
738 | DimPb[1] = fDimZEM[2]; |
739 | DimPb[2] = fDimZEM[1]; |
740 | DimPb[3] = 90.-fDimZEM[3]; |
741 | DimPb[4] = 0.; |
742 | DimPb[5] = 0.; |
743 | gMC->Gsvolu("ZEL0", "PARA", idtmed[5], DimPb, 6); |
744 | gMC->Gsvolu("ZEL1", "PARA", idtmed[5], DimPb, 6); |
745 | // gMC->Gsvolu("ZEL2", "PARA", idtmed[5], DimPb, 6); |
746 | |
747 | // --- Position the lead slices in the tranche |
748 | Float_t zTran = fDimZEM[0]/fDivZEM[2]; |
749 | Float_t zTrPb = -zTran+fDimZEMPb; |
750 | gMC->Gspos("ZEL0", 1, "ZETR", zTrPb, 0., 0., 0, "ONLY"); |
751 | gMC->Gspos("ZEL1", 1, "ZETR", fDimZEMPb, 0., 0., 0, "ONLY"); |
752 | |
753 | // --- Vacuum zone (to be filled with fibres) |
754 | DimVoid[0] = (zTran-2*fDimZEMPb)/2.; |
755 | DimVoid[1] = fDimZEM[2]; |
756 | DimVoid[2] = fDimZEM[1]; |
757 | DimVoid[3] = 90.-fDimZEM[3]; |
758 | DimVoid[4] = 0.; |
759 | DimVoid[5] = 0.; |
760 | gMC->Gsvolu("ZEV0", "PARA", idtmed[10], DimVoid,6); |
761 | gMC->Gsvolu("ZEV1", "PARA", idtmed[10], DimVoid,6); |
762 | |
763 | // --- Divide the vacuum slice into sticks along x axis |
764 | gMC->Gsdvn("ZES0", "ZEV0", fDivZEM[0], 3); |
765 | gMC->Gsdvn("ZES1", "ZEV1", fDivZEM[0], 3); |
766 | |
767 | // --- Positioning the fibers into the sticks |
768 | gMC->Gspos("ZEMF", 1,"ZES0", 0., 0., 0., irot2, "ONLY"); |
769 | gMC->Gspos("ZEMF", 1,"ZES1", 0., 0., 0., irot2, "ONLY"); |
770 | |
771 | // --- Positioning the vacuum slice into the tranche |
772 | Float_t DisplFib = fDimZEM[1]/fDivZEM[0]; |
773 | gMC->Gspos("ZEV0", 1,"ZETR", -DimVoid[0], 0., 0., 0, "ONLY"); |
774 | gMC->Gspos("ZEV1", 1,"ZETR", -DimVoid[0]+zTran, 0., DisplFib, 0, "ONLY"); |
775 | |
776 | // --- Positioning the ZEM into the ZDC - rotation for 90 degrees |
777 | gMC->Gspos("ZEM ", 1,"ZDC ", fPosZEM[0], fPosZEM[1], fPosZEM[2]+fDimZEM[0], irot1, "ONLY"); |
778 | |
779 | // --- Adding last slice at the end of the EM calorimeter |
780 | // Float_t zLastSlice = fPosZEM[2]+fDimZEMPb+fDimZEM[0]; |
781 | // gMC->Gspos("ZEL2", 1,"ZDC ", fPosZEM[0], fPosZEM[1], zLastSlice, irot1, "ONLY"); |
782 | |
783 | } |
784 | |
785 | //_____________________________________________________________________________ |
786 | void AliZDCv2::DrawModule() |
787 | { |
788 | // |
789 | // Draw a shaded view of the Zero Degree Calorimeter version 1 |
790 | // |
791 | |
792 | // Set everything unseen |
793 | gMC->Gsatt("*", "seen", -1); |
794 | // |
795 | // Set ALIC mother transparent |
796 | gMC->Gsatt("ALIC","SEEN",0); |
797 | // |
798 | // Set the volumes visible |
799 | gMC->Gsatt("ZDC ","SEEN",0); |
800 | gMC->Gsatt("QT01","SEEN",1); |
801 | gMC->Gsatt("QT02","SEEN",1); |
802 | gMC->Gsatt("QT03","SEEN",1); |
803 | gMC->Gsatt("QT04","SEEN",1); |
804 | gMC->Gsatt("QT05","SEEN",1); |
805 | gMC->Gsatt("QT06","SEEN",1); |
806 | gMC->Gsatt("QT07","SEEN",1); |
807 | gMC->Gsatt("QT08","SEEN",1); |
808 | gMC->Gsatt("QT09","SEEN",1); |
809 | gMC->Gsatt("QT10","SEEN",1); |
810 | gMC->Gsatt("QT11","SEEN",1); |
811 | gMC->Gsatt("QT12","SEEN",1); |
812 | gMC->Gsatt("QT13","SEEN",1); |
813 | gMC->Gsatt("QT14","SEEN",1); |
814 | gMC->Gsatt("QT15","SEEN",1); |
815 | gMC->Gsatt("QT16","SEEN",1); |
816 | gMC->Gsatt("QT17","SEEN",1); |
817 | gMC->Gsatt("QT18","SEEN",1); |
818 | gMC->Gsatt("QC01","SEEN",1); |
819 | gMC->Gsatt("QC02","SEEN",1); |
820 | gMC->Gsatt("QC03","SEEN",1); |
821 | gMC->Gsatt("QC04","SEEN",1); |
822 | gMC->Gsatt("QC05","SEEN",1); |
823 | gMC->Gsatt("QTD1","SEEN",1); |
824 | gMC->Gsatt("QTD2","SEEN",1); |
825 | gMC->Gsatt("QTD3","SEEN",1); |
826 | gMC->Gsatt("MQXL","SEEN",1); |
827 | gMC->Gsatt("YMQL","SEEN",1); |
828 | gMC->Gsatt("MQX ","SEEN",1); |
829 | gMC->Gsatt("YMQ ","SEEN",1); |
830 | gMC->Gsatt("ZQYX","SEEN",1); |
831 | gMC->Gsatt("MD1 ","SEEN",1); |
832 | gMC->Gsatt("MD1V","SEEN",1); |
833 | gMC->Gsatt("YD1 ","SEEN",1); |
834 | gMC->Gsatt("MD2 ","SEEN",1); |
835 | gMC->Gsatt("YD2 ","SEEN",1); |
836 | gMC->Gsatt("ZNEU","SEEN",0); |
837 | gMC->Gsatt("ZNF1","SEEN",0); |
838 | gMC->Gsatt("ZNF2","SEEN",0); |
839 | gMC->Gsatt("ZNF3","SEEN",0); |
840 | gMC->Gsatt("ZNF4","SEEN",0); |
841 | gMC->Gsatt("ZNG1","SEEN",0); |
842 | gMC->Gsatt("ZNG2","SEEN",0); |
843 | gMC->Gsatt("ZNG3","SEEN",0); |
844 | gMC->Gsatt("ZNG4","SEEN",0); |
845 | gMC->Gsatt("ZNTX","SEEN",0); |
846 | gMC->Gsatt("ZN1 ","COLO",4); |
847 | gMC->Gsatt("ZN1 ","SEEN",1); |
848 | gMC->Gsatt("ZNSL","SEEN",0); |
849 | gMC->Gsatt("ZNST","SEEN",0); |
850 | gMC->Gsatt("ZPRO","SEEN",0); |
851 | gMC->Gsatt("ZPF1","SEEN",0); |
852 | gMC->Gsatt("ZPF2","SEEN",0); |
853 | gMC->Gsatt("ZPF3","SEEN",0); |
854 | gMC->Gsatt("ZPF4","SEEN",0); |
855 | gMC->Gsatt("ZPG1","SEEN",0); |
856 | gMC->Gsatt("ZPG2","SEEN",0); |
857 | gMC->Gsatt("ZPG3","SEEN",0); |
858 | gMC->Gsatt("ZPG4","SEEN",0); |
859 | gMC->Gsatt("ZPTX","SEEN",0); |
860 | gMC->Gsatt("ZP1 ","COLO",6); |
861 | gMC->Gsatt("ZP1 ","SEEN",1); |
862 | gMC->Gsatt("ZPSL","SEEN",0); |
863 | gMC->Gsatt("ZPST","SEEN",0); |
864 | gMC->Gsatt("ZEM ","COLO",7); |
865 | gMC->Gsatt("ZEM ","SEEN",1); |
866 | gMC->Gsatt("ZEMF","SEEN",0); |
867 | gMC->Gsatt("ZETR","SEEN",0); |
868 | gMC->Gsatt("ZEL0","SEEN",0); |
869 | gMC->Gsatt("ZEL1","SEEN",0); |
870 | gMC->Gsatt("ZEL2","SEEN",0); |
871 | gMC->Gsatt("ZEV0","SEEN",0); |
872 | gMC->Gsatt("ZEV1","SEEN",0); |
873 | gMC->Gsatt("ZES0","SEEN",0); |
874 | gMC->Gsatt("ZES1","SEEN",0); |
875 | |
876 | // |
877 | gMC->Gdopt("hide", "on"); |
878 | gMC->Gdopt("shad", "on"); |
879 | gMC->Gsatt("*", "fill", 7); |
880 | gMC->SetClipBox("."); |
881 | gMC->SetClipBox("*", 0, 100, -100, 100, 12000, 16000); |
882 | gMC->DefaultRange(); |
883 | gMC->Gdraw("alic", 40, 30, 0, 488, 220, .07, .07); |
884 | gMC->Gdhead(1111, "Zero Degree Calorimeter Version 1"); |
885 | gMC->Gdman(18, 4, "MAN"); |
886 | } |
887 | |
888 | //_____________________________________________________________________________ |
889 | void AliZDCv2::CreateMaterials() |
890 | { |
891 | // |
892 | // Create Materials for the Zero Degree Calorimeter |
893 | // |
894 | |
895 | Int_t *idtmed = fIdtmed->GetArray(); |
896 | |
897 | Float_t dens, ubuf[1], wmat[2], a[2], z[2], deemax = -1; |
898 | Int_t i; |
899 | |
900 | // --- Store in UBUF r0 for nuclear radius calculation R=r0*A**1/3 |
901 | |
902 | // --- Tantalum -> ZN passive material |
903 | ubuf[0] = 1.1; |
904 | AliMaterial(1, "TANT", 180.95, 73., 16.65, .4, 11.9, ubuf, 1); |
905 | |
906 | // --- Tungsten |
907 | // ubuf[0] = 1.11; |
908 | // AliMaterial(1, "TUNG", 183.85, 74., 19.3, .35, 10.3, ubuf, 1); |
909 | |
910 | // --- Brass (CuZn) -> ZP passive material |
911 | dens = 8.48; |
912 | a[0] = 63.546; |
913 | a[1] = 65.39; |
914 | z[0] = 29.; |
915 | z[1] = 30.; |
916 | wmat[0] = .63; |
917 | wmat[1] = .37; |
918 | AliMixture(2, "BRASS ", a, z, dens, 2, wmat); |
919 | |
920 | // --- SiO2 |
921 | dens = 2.64; |
922 | a[0] = 28.086; |
923 | a[1] = 15.9994; |
924 | z[0] = 14.; |
925 | z[1] = 8.; |
926 | wmat[0] = 1.; |
927 | wmat[1] = 2.; |
928 | AliMixture(3, "SIO2 ", a, z, dens, -2, wmat); |
929 | |
930 | // --- Lead |
931 | ubuf[0] = 1.12; |
932 | AliMaterial(5, "LEAD", 207.19, 82., 11.35, .56, 18.5, ubuf, 1); |
933 | |
934 | // --- Copper |
935 | ubuf[0] = 1.10; |
936 | AliMaterial(6, "COPP", 63.54, 29., 8.96, 1.4, 0., ubuf, 1); |
937 | |
938 | // --- Iron (energy loss taken into account) |
939 | ubuf[0] = 1.1; |
940 | AliMaterial(7, "IRON", 55.85, 26., 7.87, 1.76, 0., ubuf, 1); |
941 | |
942 | // --- Iron (no energy loss) |
943 | ubuf[0] = 1.1; |
944 | AliMaterial(8, "IRON", 55.85, 26., 7.87, 1.76, 0., ubuf, 1); |
945 | |
946 | // --- Vacuum (no magnetic field) |
947 | AliMaterial(10, "VOID", 1e-16, 1e-16, 1e-16, 1e16, 1e16, ubuf,0); |
948 | |
949 | // --- Vacuum (with magnetic field) |
950 | AliMaterial(11, "VOIM", 1e-16, 1e-16, 1e-16, 1e16, 1e16, ubuf,0); |
951 | |
952 | // --- Air (no magnetic field) |
953 | AliMaterial(12, "Air $", 14.61, 7.3, .001205, 30420., 67500., ubuf, 0); |
954 | |
955 | // --- Definition of tracking media: |
956 | |
957 | // --- Tantalum = 1 ; |
958 | // --- Brass = 2 ; |
959 | // --- Fibers (SiO2) = 3 ; |
960 | // --- Fibers (SiO2) = 4 ; |
961 | // --- Lead = 5 ; |
962 | // --- Copper = 6 ; |
963 | // --- Iron (with energy loss) = 7 ; |
964 | // --- Iron (without energy loss) = 8 ; |
965 | // --- Vacuum (no field) = 10 |
966 | // --- Vacuum (with field) = 11 |
967 | // --- Air (no field) = 12 |
968 | |
969 | |
970 | // --- Tracking media parameters |
971 | Float_t epsil = .01, stmin=0.01, stemax = 1.; |
7cc664ae |
972 | // Int_t isxfld = gAlice->Field()->Integ(); |
4663d63d |
973 | Float_t fieldm = 0., tmaxfd = 0.; |
974 | Int_t ifield = 0, isvolActive = 1, isvol = 0, inofld = 0; |
975 | |
976 | AliMedium(1, "ZTANT", 1, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
977 | // AliMedium(1, "ZW", 1, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
978 | AliMedium(2, "ZBRASS",2, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
979 | AliMedium(3, "ZSIO2", 3, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
980 | AliMedium(4, "ZQUAR", 3, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
981 | AliMedium(5, "ZLEAD", 5, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
982 | // AliMedium(6, "ZCOPP", 6, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
983 | // AliMedium(7, "ZIRON", 7, isvolActive, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
984 | AliMedium(6, "ZCOPP", 6, isvol, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
985 | AliMedium(7, "ZIRON", 7, isvol, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
986 | AliMedium(8, "ZIRONN",8, isvol, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
987 | AliMedium(10,"ZVOID",10, isvol, inofld, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
988 | AliMedium(12,"ZAIR", 12, 0, inofld, fieldm, tmaxfd, stemax,deemax, epsil, stmin); |
989 | |
990 | ifield =2; |
991 | fieldm = 45.; |
7cc664ae |
992 | AliMedium(11, "ZVOIM", 11, isvol, ifield, fieldm, tmaxfd, stemax, deemax, epsil, stmin); |
4663d63d |
993 | |
994 | // Thresholds for showering in the ZDCs |
995 | i = 1; //tantalum |
996 | gMC->Gstpar(idtmed[i], "CUTGAM", .001); |
997 | gMC->Gstpar(idtmed[i], "CUTELE", .001); |
998 | gMC->Gstpar(idtmed[i], "CUTNEU", .01); |
999 | gMC->Gstpar(idtmed[i], "CUTHAD", .01); |
1000 | i = 2; //brass |
1001 | gMC->Gstpar(idtmed[i], "CUTGAM", .001); |
1002 | gMC->Gstpar(idtmed[i], "CUTELE", .001); |
1003 | gMC->Gstpar(idtmed[i], "CUTNEU", .01); |
1004 | gMC->Gstpar(idtmed[i], "CUTHAD", .01); |
1005 | i = 5; //lead |
1006 | gMC->Gstpar(idtmed[i], "CUTGAM", .001); |
1007 | gMC->Gstpar(idtmed[i], "CUTELE", .001); |
1008 | gMC->Gstpar(idtmed[i], "CUTNEU", .01); |
1009 | gMC->Gstpar(idtmed[i], "CUTHAD", .01); |
1010 | |
1011 | // Avoid too detailed showering in TDI |
1012 | i = 6; //copper |
1013 | gMC->Gstpar(idtmed[i], "CUTGAM", .1); |
1014 | gMC->Gstpar(idtmed[i], "CUTELE", .1); |
1015 | gMC->Gstpar(idtmed[i], "CUTNEU", 1.); |
1016 | gMC->Gstpar(idtmed[i], "CUTHAD", 1.); |
1017 | |
1018 | // Avoid too detailed showering along the beam line |
1019 | i = 7; //iron with energy loss (ZIRON) |
1020 | gMC->Gstpar(idtmed[i], "CUTGAM", .1); |
1021 | gMC->Gstpar(idtmed[i], "CUTELE", .1); |
1022 | gMC->Gstpar(idtmed[i], "CUTNEU", 1.); |
1023 | gMC->Gstpar(idtmed[i], "CUTHAD", 1.); |
1024 | |
1025 | // Avoid too detailed showering along the beam line |
1026 | i = 8; //iron with energy loss (ZIRONN) |
1027 | gMC->Gstpar(idtmed[i], "CUTGAM", .1); |
1028 | gMC->Gstpar(idtmed[i], "CUTELE", .1); |
1029 | gMC->Gstpar(idtmed[i], "CUTNEU", 1.); |
1030 | gMC->Gstpar(idtmed[i], "CUTHAD", 1.); |
1031 | |
1032 | // Avoid interaction in fibers (only energy loss allowed) |
1033 | i = 3; //fibers (ZSI02) |
1034 | gMC->Gstpar(idtmed[i], "DCAY", 0.); |
1035 | gMC->Gstpar(idtmed[i], "MULS", 0.); |
1036 | gMC->Gstpar(idtmed[i], "PFIS", 0.); |
1037 | gMC->Gstpar(idtmed[i], "MUNU", 0.); |
1038 | gMC->Gstpar(idtmed[i], "LOSS", 1.); |
1039 | gMC->Gstpar(idtmed[i], "PHOT", 0.); |
1040 | gMC->Gstpar(idtmed[i], "COMP", 0.); |
1041 | gMC->Gstpar(idtmed[i], "PAIR", 0.); |
1042 | gMC->Gstpar(idtmed[i], "BREM", 0.); |
1043 | gMC->Gstpar(idtmed[i], "DRAY", 0.); |
1044 | gMC->Gstpar(idtmed[i], "ANNI", 0.); |
1045 | gMC->Gstpar(idtmed[i], "HADR", 0.); |
1046 | i = 4; //fibers (ZQUAR) |
1047 | gMC->Gstpar(idtmed[i], "DCAY", 0.); |
1048 | gMC->Gstpar(idtmed[i], "MULS", 0.); |
1049 | gMC->Gstpar(idtmed[i], "PFIS", 0.); |
1050 | gMC->Gstpar(idtmed[i], "MUNU", 0.); |
1051 | gMC->Gstpar(idtmed[i], "LOSS", 1.); |
1052 | gMC->Gstpar(idtmed[i], "PHOT", 0.); |
1053 | gMC->Gstpar(idtmed[i], "COMP", 0.); |
1054 | gMC->Gstpar(idtmed[i], "PAIR", 0.); |
1055 | gMC->Gstpar(idtmed[i], "BREM", 0.); |
1056 | gMC->Gstpar(idtmed[i], "DRAY", 0.); |
1057 | gMC->Gstpar(idtmed[i], "ANNI", 0.); |
1058 | gMC->Gstpar(idtmed[i], "HADR", 0.); |
1059 | |
1060 | // Avoid interaction in void |
1061 | i = 11; //void with field |
1062 | gMC->Gstpar(idtmed[i], "DCAY", 0.); |
1063 | gMC->Gstpar(idtmed[i], "MULS", 0.); |
1064 | gMC->Gstpar(idtmed[i], "PFIS", 0.); |
1065 | gMC->Gstpar(idtmed[i], "MUNU", 0.); |
1066 | gMC->Gstpar(idtmed[i], "LOSS", 0.); |
1067 | gMC->Gstpar(idtmed[i], "PHOT", 0.); |
1068 | gMC->Gstpar(idtmed[i], "COMP", 0.); |
1069 | gMC->Gstpar(idtmed[i], "PAIR", 0.); |
1070 | gMC->Gstpar(idtmed[i], "BREM", 0.); |
1071 | gMC->Gstpar(idtmed[i], "DRAY", 0.); |
1072 | gMC->Gstpar(idtmed[i], "ANNI", 0.); |
1073 | gMC->Gstpar(idtmed[i], "HADR", 0.); |
1074 | |
1075 | // |
1076 | fMedSensZN = idtmed[1]; // Sensitive volume: ZN passive material |
1077 | fMedSensZP = idtmed[2]; // Sensitive volume: ZP passive material |
1078 | fMedSensF1 = idtmed[3]; // Sensitive volume: fibres type 1 |
1079 | fMedSensF2 = idtmed[4]; // Sensitive volume: fibres type 2 |
1080 | fMedSensZEM = idtmed[5]; // Sensitive volume: ZEM passive material |
1081 | // fMedSensTDI = idtmed[6]; // Sensitive volume: TDI Cu shield |
1082 | // fMedSensPI = idtmed[7]; // Sensitive volume: beam pipes |
1083 | fMedSensGR = idtmed[12]; // Sensitive volume: air into the grooves |
1084 | } |
1085 | |
1086 | //_____________________________________________________________________________ |
1087 | void AliZDCv2::Init() |
1088 | { |
1089 | InitTables(); |
1090 | } |
1091 | |
1092 | //_____________________________________________________________________________ |
1093 | void AliZDCv2::InitTables() |
1094 | { |
1095 | Int_t k, j; |
1096 | |
1097 | char *lightfName1,*lightfName2,*lightfName3,*lightfName4, |
1098 | *lightfName5,*lightfName6,*lightfName7,*lightfName8; |
1099 | FILE *fp1, *fp2, *fp3, *fp4, *fp5, *fp6, *fp7, *fp8; |
1100 | |
1101 | // --- Reading light tables for ZN |
1102 | lightfName1 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620362207s"); |
1103 | if((fp1 = fopen(lightfName1,"r")) == NULL){ |
1104 | printf("Cannot open file fp1 \n"); |
1105 | return; |
1106 | } |
1107 | lightfName2 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620362208s"); |
1108 | if((fp2 = fopen(lightfName2,"r")) == NULL){ |
1109 | printf("Cannot open file fp2 \n"); |
1110 | return; |
1111 | } |
1112 | lightfName3 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620362209s"); |
1113 | if((fp3 = fopen(lightfName3,"r")) == NULL){ |
1114 | printf("Cannot open file fp3 \n"); |
1115 | return; |
1116 | } |
1117 | lightfName4 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620362210s"); |
1118 | if((fp4 = fopen(lightfName4,"r")) == NULL){ |
1119 | printf("Cannot open file fp4 \n"); |
1120 | return; |
1121 | } |
1122 | |
1123 | for(k=0; k<fNalfan; k++){ |
1124 | for(j=0; j<fNben; j++){ |
1125 | fscanf(fp1,"%f",&fTablen[0][k][j]); |
1126 | fscanf(fp2,"%f",&fTablen[1][k][j]); |
1127 | fscanf(fp3,"%f",&fTablen[2][k][j]); |
1128 | fscanf(fp4,"%f",&fTablen[3][k][j]); |
1129 | } |
1130 | } |
1131 | fclose(fp1); |
1132 | fclose(fp2); |
1133 | fclose(fp3); |
1134 | fclose(fp4); |
1135 | |
1136 | // --- Reading light tables for ZP and ZEM |
1137 | lightfName5 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620552207s"); |
1138 | if((fp5 = fopen(lightfName5,"r")) == NULL){ |
1139 | printf("Cannot open file fp5 \n"); |
1140 | return; |
1141 | } |
1142 | lightfName6 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620552208s"); |
1143 | if((fp6 = fopen(lightfName6,"r")) == NULL){ |
1144 | printf("Cannot open file fp6 \n"); |
1145 | return; |
1146 | } |
1147 | lightfName7 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620552209s"); |
1148 | if((fp7 = fopen(lightfName7,"r")) == NULL){ |
1149 | printf("Cannot open file fp7 \n"); |
1150 | return; |
1151 | } |
1152 | lightfName8 = gSystem->ExpandPathName("$ALICE/$ALICE_LEVEL/ZDC/light22620552210s"); |
1153 | if((fp8 = fopen(lightfName8,"r")) == NULL){ |
1154 | printf("Cannot open file fp8 \n"); |
1155 | return; |
1156 | } |
1157 | |
1158 | for(k=0; k<fNalfap; k++){ |
1159 | for(j=0; j<fNbep; j++){ |
1160 | fscanf(fp5,"%f",&fTablep[0][k][j]); |
1161 | fscanf(fp6,"%f",&fTablep[1][k][j]); |
1162 | fscanf(fp7,"%f",&fTablep[2][k][j]); |
1163 | fscanf(fp8,"%f",&fTablep[3][k][j]); |
1164 | } |
1165 | } |
1166 | fclose(fp5); |
1167 | fclose(fp6); |
1168 | fclose(fp7); |
1169 | fclose(fp8); |
1170 | } |
1171 | |
1172 | //_____________________________________________________________________________ |
1173 | Int_t AliZDCv2::Digitize(Int_t Det, Int_t Quad, Int_t Light) |
1174 | { |
1175 | // Evaluation of the ADC channel corresponding to the light yield Light |
1176 | |
1177 | if(fDebug == 1){ |
1178 | printf("\n Digitize -> Det = %d, Quad = %d, Light = %d\n", Det, Quad, Light); |
1179 | } |
1180 | |
1181 | // Parameters for conversion of light yield in ADC channels |
1182 | Float_t fPMGain[3][5]; // PM gain |
1183 | Float_t fADCRes; // ADC conversion factor |
1184 | |
1185 | Int_t j,i; |
1186 | for(i=0; i<3; i++){ |
1187 | for(j=0; j<5; j++){ |
1188 | fPMGain[i][j] = 100000.; |
1189 | } |
1190 | } |
1191 | fADCRes = 0.00000064; // ADC Resolution: 250 fC/ADCch |
1192 | |
1193 | Int_t ADCch = Int_t(Light*fPMGain[Det-1][Quad]*fADCRes); |
1194 | |
1195 | return ADCch; |
1196 | } |
1197 | |
1198 | |
1199 | //_____________________________________________________________________________ |
1200 | void AliZDCv2::SDigits2Digits() |
1201 | { |
1202 | Hits2Digits(gAlice->GetNtrack()); |
1203 | } |
1204 | |
1205 | //_____________________________________________________________________________ |
1206 | void AliZDCv2::Hits2Digits(Int_t ntracks) |
1207 | { |
1208 | AliZDCDigit *newdigit; |
1209 | AliZDCHit *hit; |
1210 | |
1211 | Int_t PMCZN = 0, PMCZP = 0, PMQZN[4], PMQZP[4], PMZEM = 0; |
1212 | |
1213 | Int_t i; |
1214 | for(i=0; i<4; i++){ |
1215 | PMQZN[i] =0; |
1216 | PMQZP[i] =0; |
1217 | } |
1218 | |
1219 | Int_t itrack = 0; |
1220 | for(itrack=0; itrack<ntracks; itrack++){ |
1221 | gAlice->ResetHits(); |
1222 | gAlice->TreeH()->GetEvent(itrack); |
1223 | for(i=0; i<fHits->GetEntries(); i++){ |
1224 | hit = (AliZDCHit*)fHits->At(i); |
1225 | Int_t det = hit->GetVolume(0); |
1226 | Int_t quad = hit->GetVolume(1); |
1227 | Int_t lightQ = Int_t(hit->GetLightPMQ()); |
1228 | Int_t lightC = Int_t(hit->GetLightPMC()); |
1229 | if(fDebug == 1) |
1230 | printf(" \n itrack = %d, fNhits = %d, det = %d, quad = %d," |
1231 | "lightC = %d lightQ = %d\n", itrack, fNhits, det, quad, lightC, lightQ); |
1232 | |
1233 | if(det == 1){ //ZN |
1234 | PMCZN = PMCZN + lightC; |
1235 | PMQZN[quad-1] = PMQZN[quad-1] + lightQ; |
1236 | } |
1237 | |
1238 | if(det == 2){ //ZP |
1239 | PMCZP = PMCZP + lightC; |
1240 | PMQZP[quad-1] = PMQZP[quad-1] + lightQ; |
1241 | } |
1242 | |
1243 | if(det == 3){ //ZEM |
1244 | PMZEM = PMZEM + lightC; |
1245 | } |
1246 | } // Hits loop |
1247 | |
1248 | } // Tracks loop |
1249 | |
1250 | if(fDebug == 1){ |
1251 | printf("\n PMCZN = %d, PMQZN[0] = %d, PMQZN[1] = %d, PMQZN[2] = %d, PMQZN[3] = %d\n" |
1252 | , PMCZN, PMQZN[0], PMQZN[1], PMQZN[2], PMQZN[3]); |
1253 | printf("\n PMCZP = %d, PMQZP[0] = %d, PMQZP[1] = %d, PMQZP[2] = %d, PMQZP[3] = %d\n" |
1254 | , PMCZP, PMQZP[0], PMQZP[1], PMQZP[2], PMQZP[3]); |
1255 | printf("\n PMZEM = %d\n", PMZEM); |
1256 | } |
1257 | |
1258 | // ------------------------------------ Hits2Digits |
1259 | // Digits for ZN |
1260 | newdigit = new AliZDCDigit(1, 0, Digitize(1, 0, PMCZN)); |
1261 | new((*fDigits)[fNdigits]) AliZDCDigit(*newdigit); |
1262 | fNdigits++; |
1263 | delete newdigit; |
1264 | |
1265 | Int_t j; |
1266 | for(j=0; j<4; j++){ |
1267 | newdigit = new AliZDCDigit(1, j+1, Digitize(1, j+1, PMQZN[j])); |
1268 | new((*fDigits)[fNdigits]) AliZDCDigit(*newdigit); |
1269 | fNdigits++; |
1270 | delete newdigit; |
1271 | } |
1272 | |
1273 | // Digits for ZP |
1274 | newdigit = new AliZDCDigit(2, 0, Digitize(2, 0, PMCZP)); |
1275 | new((*fDigits)[fNdigits]) AliZDCDigit(*newdigit); |
1276 | fNdigits++; |
1277 | delete newdigit; |
1278 | |
1279 | Int_t k; |
1280 | for(k=0; k<4; k++){ |
1281 | newdigit = new AliZDCDigit(2, k+1, Digitize(2, k+1, PMQZP[k])); |
1282 | new((*fDigits)[fNdigits]) AliZDCDigit(*newdigit); |
1283 | fNdigits++; |
1284 | delete newdigit; |
1285 | } |
1286 | |
1287 | // Digits for ZEM |
1288 | newdigit = new AliZDCDigit(3, 0, Digitize(3, 0, PMZEM)); |
1289 | new((*fDigits)[fNdigits]) AliZDCDigit(*newdigit); |
1290 | fNdigits++; |
1291 | delete newdigit; |
1292 | |
1293 | |
1294 | gAlice->TreeD()->Fill(); |
1295 | gAlice->TreeD()->Write(0,TObject::kOverwrite); |
1296 | |
1297 | // if(fDebug == 1){ |
1298 | // printf("\n Event Digits -----------------------------------------------------\n"); |
1299 | // fDigits->Print(""); |
1300 | // } |
1301 | |
1302 | } |
1303 | //_____________________________________________________________________________ |
1304 | void AliZDCv2::MakeBranch(Option_t *opt, char *file) |
1305 | { |
1306 | // |
1307 | // Create a new branch in the current Root Tree |
1308 | // |
1309 | |
1310 | AliDetector::MakeBranch(opt); |
1311 | |
1312 | Char_t branchname[10]; |
1313 | sprintf(branchname,"%s",GetName()); |
1314 | const char *cD = strstr(opt,"D"); |
1315 | |
1316 | if (gAlice->TreeD() && cD) { |
1317 | |
1318 | // Creation of the digits from hits |
1319 | |
1320 | if(fDigits!=0) fDigits->Clear(); |
1321 | else fDigits = new TClonesArray ("AliZDCDigit",1000); |
1322 | char branchname[10]; |
1323 | sprintf(branchname,"%s",GetName()); |
1324 | gAlice->MakeBranchInTree(gAlice->TreeD(), |
1325 | branchname, &fDigits, fBufferSize, file) ; |
1326 | printf("* AliZDCv2::MakeBranch * Making Branch %s for digits\n\n",branchname); |
1327 | } |
1328 | |
1329 | } |
1330 | //_____________________________________________________________________________ |
1331 | void AliZDCv2::StepManager() |
1332 | { |
1333 | // |
1334 | // Routine called at every step in the Zero Degree Calorimeters |
1335 | // |
1336 | |
1337 | Int_t j, vol[2], ibeta=0, ialfa, ibe, nphe; |
1338 | Float_t x[3], xdet[3], destep, hits[10], m, ekin, um[3], ud[3], be, radius, out; |
1339 | TLorentzVector s, p; |
1340 | const char *knamed; |
1341 | |
1342 | for (j=0;j<10;j++) hits[j]=0; |
1343 | |
1344 | if((gMC->GetMedium() == fMedSensZN) || (gMC->GetMedium() == fMedSensZP) || |
1345 | (gMC->GetMedium() == fMedSensGR) || (gMC->GetMedium() == fMedSensF1) || |
1346 | (gMC->GetMedium() == fMedSensF2) || (gMC->GetMedium() == fMedSensZEM)){ |
1347 | // (gMC->GetMedium() == fMedSensPI) || (gMC->GetMedium() == fMedSensTDI)){ |
1348 | |
1349 | // If particle interacts with beam pipe -> return |
1350 | // if((gMC->GetMedium() == fMedSensPI) || (gMC->GetMedium() == fMedSensTDI)){ |
1351 | // If option NoShower is set -> StopTrack |
1352 | // if(fNoShower==1) { |
1353 | // if(gMC->GetMedium() == fMedSensPI) { |
1354 | // knamed = gMC->CurrentVolName(); |
1355 | // if((!strncmp(knamed,"MQ",2)) || (!strncmp(knamed,"YM",2))) fpLostIT += 1; |
1356 | // if((!strncmp(knamed,"MD1",3))|| (!strncmp(knamed,"YD1",2))) fpLostD1 += 1; |
1357 | // } |
1358 | // if(gMC->GetMedium() == fMedSensTDI) fpLostTDI += 1; |
1359 | // gMC->StopTrack(); |
1360 | // printf("\n # of p lost in Inner Triplet = %d\n",fpLostIT); |
1361 | // printf("\n # of p lost in D1 = %d\n",fpLostD1); |
1362 | // printf("\n # of p lost in TDI = %d\n",fpLostTDI); |
4663d63d |
1363 | // } |
7cc664ae |
1364 | // return; |
4663d63d |
1365 | // } |
1366 | |
1367 | //Particle coordinates |
1368 | gMC->TrackPosition(s); |
1369 | for(j=0; j<=2; j++){ |
1370 | x[j] = s[j]; |
1371 | } |
1372 | hits[0] = x[0]; |
1373 | hits[1] = x[1]; |
1374 | hits[2] = x[2]; |
1375 | |
1376 | // Determine in which ZDC the particle is |
1377 | knamed = gMC->CurrentVolName(); |
1378 | if(!strncmp(knamed,"ZN",2))vol[0]=1; |
1379 | if(!strncmp(knamed,"ZP",2))vol[0]=2; |
1380 | if(!strncmp(knamed,"ZE",2))vol[0]=3; |
1381 | |
1382 | // Determine in which quadrant the particle is |
1383 | |
1384 | //Quadrant in ZN |
1385 | if(vol[0]==1){ |
1386 | xdet[0] = x[0]-fPosZN[0]; |
1387 | xdet[1] = x[1]-fPosZN[1]; |
1388 | if((xdet[0]<=0.) && (xdet[1]>=0.)) vol[1]=1; |
1389 | if((xdet[0]>0.) && (xdet[1]>0.)) vol[1]=2; |
1390 | if((xdet[0]<0.) && (xdet[1]<0.)) vol[1]=3; |
1391 | if((xdet[0]>0.) && (xdet[1]<0.)) vol[1]=4; |
1392 | } |
1393 | |
1394 | //Quadrant in ZP |
1395 | if(vol[0]==2){ |
1396 | xdet[0] = x[0]-fPosZP[0]; |
1397 | xdet[1] = x[1]-fPosZP[1]; |
1398 | if(xdet[0]>fDimZP[0])xdet[0]=fDimZP[0]-0.01; |
1399 | if(xdet[0]<-fDimZP[0])xdet[0]=-fDimZP[0]+0.01; |
1400 | Float_t xqZP = xdet[0]/(fDimZP[0]/2); |
1401 | for(int i=1; i<=4; i++){ |
1402 | if(xqZP>=(i-3) && xqZP<(i-2)){ |
1403 | vol[1] = i; |
1404 | break; |
1405 | } |
1406 | } |
1407 | } |
1408 | |
1409 | //ZEM has only 1 quadrant |
1410 | if(vol[0] == 3){ |
1411 | vol[1] = 1; |
1412 | xdet[0] = x[0]-fPosZEM[0]; |
1413 | xdet[1] = x[1]-fPosZEM[1]; |
1414 | } |
1415 | |
1416 | // Store impact point and kinetic energy of the ENTERING particle |
1417 | |
1418 | // if(Curtrack==Prim){ |
1419 | if(gMC->IsTrackEntering()){ |
1420 | //Particle energy |
1421 | gMC->TrackMomentum(p); |
1422 | hits[3] = p[3]; |
1423 | // Impact point on ZDC |
1424 | hits[4] = xdet[0]; |
1425 | hits[5] = xdet[1]; |
1426 | hits[6] = 0; |
1427 | hits[7] = 0; |
1428 | hits[8] = 0; |
1429 | hits[9] = 0; |
1430 | |
1431 | // Int_t PcID = gMC->TrackPid(); |
1432 | // printf("Pc ID -> %d\n",PcID); |
1433 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1434 | |
1435 | if(fNoShower==1){ |
1436 | // fpDetected += 1; |
1437 | gMC->StopTrack(); |
1438 | // printf("\n # of detected p = %d\n",fpDetected); |
1439 | return; |
1440 | } |
1441 | } |
1442 | // } // Curtrack IF |
1443 | |
1444 | // Charged particles -> Energy loss |
1445 | if((destep=gMC->Edep())){ |
1446 | if(gMC->IsTrackStop()){ |
1447 | gMC->TrackMomentum(p); |
1448 | m = gMC->TrackMass(); |
1449 | ekin = p[3]-m; |
1450 | hits[9] = ekin; |
1451 | hits[7] = 0.; |
1452 | hits[8] = 0.; |
1453 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1454 | } |
1455 | else{ |
1456 | hits[9] = destep; |
1457 | hits[7] = 0.; |
1458 | hits[8] = 0.; |
1459 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1460 | } |
1461 | // printf(" Dep. E = %f \n",hits[9]); |
1462 | } |
1463 | }// NB -> Questa parentesi (chiude il primo IF) io la sposterei al fondo!??? |
1464 | |
1465 | |
1466 | // *** Light production in fibres |
1467 | if((gMC->GetMedium() == fMedSensF1) || (gMC->GetMedium() == fMedSensF2)){ |
1468 | |
1469 | //Select charged particles |
1470 | if((destep=gMC->Edep())){ |
1471 | |
1472 | // Particle velocity |
1473 | gMC->TrackMomentum(p); |
1474 | Float_t ptot=TMath::Sqrt(p[0]*p[0]+p[1]*p[1]+p[2]*p[2]); |
1475 | Float_t beta = ptot/p[3]; |
1476 | if(beta<0.67) return; |
1477 | if((beta>=0.67) && (beta<=0.75)) ibeta = 0; |
1478 | if((beta>0.75) && (beta<=0.85)) ibeta = 1; |
1479 | if((beta>0.85) && (beta<=0.95)) ibeta = 2; |
1480 | if(beta>0.95) ibeta = 3; |
1481 | |
1482 | // Angle between particle trajectory and fibre axis |
1483 | // 1 -> Momentum directions |
1484 | um[0] = p[0]/ptot; |
1485 | um[1] = p[1]/ptot; |
1486 | um[2] = p[2]/ptot; |
1487 | gMC->Gmtod(um,ud,2); |
1488 | // 2 -> Angle < limit angle |
1489 | Double_t alfar = TMath::ACos(ud[2]); |
1490 | Double_t alfa = alfar*kRaddeg; |
1491 | if(alfa>=110.) return; |
1492 | ialfa = Int_t(1.+alfa/2.); |
1493 | |
1494 | // Distance between particle trajectory and fibre axis |
1495 | gMC->TrackPosition(s); |
1496 | for(j=0; j<=2; j++){ |
1497 | x[j] = s[j]; |
1498 | } |
1499 | gMC->Gmtod(x,xdet,1); |
1500 | if(TMath::Abs(ud[0])>0.00001){ |
1501 | Float_t dcoeff = ud[1]/ud[0]; |
1502 | be = TMath::Abs((xdet[1]-dcoeff*xdet[0])/TMath::Sqrt(dcoeff*dcoeff+1.)); |
1503 | } |
1504 | else{ |
1505 | be = TMath::Abs(ud[0]); |
1506 | } |
1507 | |
1508 | if((vol[0]==1)) radius = fFibZN[1]; |
1509 | if((vol[0]==2)) radius = fFibZP[1]; |
1510 | ibe = Int_t(be*1000.+1); |
1511 | |
1512 | //Looking into the light tables |
1513 | Float_t charge = gMC->TrackCharge(); |
1514 | |
1515 | // (1) ZN |
1516 | if((vol[0]==1)) { |
1517 | if(ibe>fNben) ibe=fNben; |
1518 | out = charge*charge*fTablen[ibeta][ialfa][ibe]; |
1519 | nphe = gRandom->Poisson(out); |
1520 | // printf("ZN --- ibeta = %d, ialfa = %d, ibe = %d" |
1521 | // " -> out = %f, nphe = %d\n", ibeta, ialfa, ibe, out, nphe); |
1522 | if(gMC->GetMedium() == fMedSensF1){ |
1523 | hits[7] = nphe; //fLightPMQ |
1524 | hits[8] = 0; |
1525 | hits[9] = 0; |
1526 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1527 | } |
1528 | else{ |
1529 | hits[7] = 0; |
1530 | hits[8] = nphe; //fLightPMC |
1531 | hits[9] = 0; |
1532 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1533 | } |
1534 | } |
1535 | |
1536 | // (2) ZP |
1537 | if((vol[0]==2)) { |
1538 | if(ibe>fNbep) ibe=fNbep; |
1539 | out = charge*charge*fTablep[ibeta][ialfa][ibe]; |
1540 | nphe = gRandom->Poisson(out); |
1541 | // printf("ZP --- ibeta = %d, ialfa = %d, ibe = %d" |
1542 | // " -> out = %f, nphe = %d\n", ibeta, ialfa, ibe, out, nphe); |
1543 | if(gMC->GetMedium() == fMedSensF1){ |
1544 | hits[7] = nphe; //fLightPMQ |
1545 | hits[8] = 0; |
1546 | hits[9] = 0; |
1547 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1548 | } |
1549 | else{ |
1550 | hits[7] = 0; |
1551 | hits[8] = nphe; //fLightPMC |
1552 | hits[9] = 0; |
1553 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1554 | } |
1555 | } |
1556 | // (3) ZEM |
1557 | if((vol[0]==3)) { |
1558 | if(ibe>fNbep) ibe=fNbep; |
1559 | out = charge*charge*fTablep[ibeta][ialfa][ibe]; |
1560 | nphe = gRandom->Poisson(out); |
1561 | // printf("ZEM --- ibeta = %d, ialfa = %d, ibe = %d" |
1562 | // " -> out = %f, nphe = %d\n", ibeta, ialfa, ibe, out, nphe); |
1563 | hits[7] = 0; |
1564 | hits[8] = nphe; //fLightPMC |
1565 | hits[9] = 0; |
1566 | AddHit(gAlice->CurrentTrack(), vol, hits); |
1567 | } |
1568 | } |
1569 | } |
1570 | } |