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