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