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