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