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