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