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