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