]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CRT/AliCRTv0.cxx
cout substituted with Info
[u/mrichter/AliRoot.git] / CRT / AliCRTv0.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.9  2002/10/23 06:47:56  alibrary
19 Introducing Riostream.h
20
21 Revision 1.8  2002/10/14 14:55:34  hristov
22 Merging the VirtualMC branch to the main development branch (HEAD)
23
24 Revision 1.4.2.4  2002/10/10 14:40:31  hristov
25 Updating VirtualMC to v3-09-02
26
27 Revision 1.7  2002/10/07 11:13:25  gamez
28 Access shafts added
29
30 Revision 1.6  2002/07/26 06:21:12  gamez
31 CRT3 volume taken as sensitive volume
32
33 Revision 1.5  2002/07/25 12:52:34  morsch
34 AddHit call only if hit has been defined.
35
36 Revision 1.4  2002/07/12 12:57:29  gamez
37 Division of CRT1 corrected
38
39 Revision 1.3.2.1  2002/07/12 12:32:50  gamez
40 Division of CRT1 corrected
41
42 Revision 1.3  2002/07/10 15:57:04  gamez
43 CreateHall() removed, and new Molasse volumes
44
45 Revision 1.2  2002/07/09 08:45:35  hristov
46 Old style include files needed on HP (aCC)
47
48 Revision 1.1  2002/06/16 17:08:19  hristov
49 First version of CRT
50
51
52 */
53
54 ///////////////////////////////////////////////////////////////////////////////
55 //                                                                           //
56 // ALICE Cosmic Ray Trigger                                                  //
57 //                                                                           //
58 //  This class contains the functions for version 0 of the ALICE Cosmic Ray  //
59 //  Trigger. This version will be used to simulation comic rays in alice     //
60 //  with all the detectors.                                                  //
61 //
62 //   Authors:
63 //
64 //   Arturo Fernandez <afernand@fcfm.buap.mx>
65 //   Enrique Gamez    <egamez@fcfm.buap.mx>
66 //
67 //   Universidad Autonoma de Puebla
68 //
69 //
70 //Begin_Html
71 /*
72 <img src="picts/AliCRTv0Class.gif">
73 </pre>
74 <br clear=left>
75 <p>The responsible person for this module is
76 <a href="mailto:egamez@fcfm.buap.mx">Enrique Gamez</a>.
77 </font>
78 <pre>
79 */
80 //End_Html
81 //                                                                           //
82 ///////////////////////////////////////////////////////////////////////////////
83
84 #include <Riostream.h>
85
86 #include <TGeometry.h>
87 #include <TBRIK.h>
88 #include <TNode.h>
89 #include <TLorentzVector.h>
90
91 #include "AliRun.h"
92 #include "AliMC.h"
93 #include "AliMagF.h"
94 #include "AliConst.h"
95 #include "AliPDG.h"
96
97 #include "AliCRTv0.h"
98 #include "AliCRTConstants.h"
99
100 ClassImp(AliCRTv0)
101  
102 //_____________________________________________________________________________
103 AliCRTv0::AliCRTv0() : AliCRT()
104 {
105   //
106   // Default constructor for CRT v0
107   //
108 }
109  
110 //_____________________________________________________________________________
111 AliCRTv0::AliCRTv0(const char *name, const char *title)
112   : AliCRT(name,title)
113 {
114   //
115   // Standard constructor for CRT v0
116   //
117   //Begin_Html
118   /*
119     <img src="picts/AliCRTv0.gif">
120   */
121   //End_Html
122 }
123
124 //_____________________________________________________________________________
125 AliCRTv0::AliCRTv0(const AliCRTv0& crt)
126 {
127   //
128   // Copy ctor.
129   //
130   crt.Copy(*this);
131 }
132
133 //_____________________________________________________________________________
134 AliCRTv0& AliCRTv0::operator= (const AliCRTv0& crt)
135 {
136   //
137   // Asingment operator.
138   //
139   crt.Copy(*this);
140   return *this;
141 }
142
143 //_____________________________________________________________________________
144 void AliCRTv0::BuildGeometry()
145 {
146   //
147   // Create the ROOT TNode geometry for the CRT
148   //
149
150   TNode *node, *top;
151
152   const Int_t kColorCRT = kRed;
153
154   // Find the top node alice.
155   top = gAlice->GetGeometry()->GetNode("alice");
156
157   new TBRIK("S_CRT_A", "CRT box", "void", 
158             AliCRTConstants::fgActiveAreaLenght/2., 
159             AliCRTConstants::fgActiveAreaHeight/2., 
160             AliCRTConstants::fgActiveAreaWidth/2.);
161
162   
163   new TRotMatrix("Left", "Left", 90., 315., 90., 45., 0., 337.5);
164   new TRotMatrix("Right", "Right", 90., 45., 90., 315., 180., 202.5);
165   new TRotMatrix("Up", "Up", 90., 0., 90., 90., 0., 90.);
166   top->cd();
167
168   //
169   // Put 4 modules on the top of the magnet
170   Float_t box = AliCRTConstants::fgCageWidth/2.;
171   top->cd();
172   node = new TNode("upper1", "upper1", "S_CRT_A", 0., 790.,  3.*box, "Up");
173   node->SetLineColor(kColorCRT);
174   fNodes->Add(node);
175
176   top->cd();
177   node = new TNode("upper2", "upper2", "S_CRT_A", 0., 790.,    box, "Up");
178   node->SetLineColor(kColorCRT);
179   fNodes->Add(node);
180
181   top->cd();
182   node = new TNode("upper3", "upper3", "S_CRT_A", 0., 790., -1.*box, "Up");
183   node->SetLineColor(kColorCRT);
184   fNodes->Add(node);
185
186   top->cd();
187   node = new TNode("upper4", "upper4", "S_CRT_A", 0., 790., -3.*box, "Up");
188   node->SetLineColor(kColorCRT);
189   fNodes->Add(node);
190
191
192   // Modules on the left side.
193   Float_t xtragap = 10.;
194   Float_t initXside = (790.+xtragap)*TMath::Sin(2*22.5*kDegrad); //rigth side
195   Float_t initYside = (790.+xtragap)*TMath::Cos(2*22.5*kDegrad);
196   top->cd();
197   node = new TNode("upper5", "upper5", "S_CRT_A", initXside, initYside,  3.*box, "Left");
198   node->SetLineColor(kColorCRT);
199   fNodes->Add(node);
200
201   top->cd();
202   node = new TNode("upper6", "upper6", "S_CRT_A", initXside, initYside,    box, "Left");
203   node->SetLineColor(kColorCRT);
204   fNodes->Add(node);
205
206   top->cd();
207   node = new TNode("upper7", "upper7", "S_CRT_A", initXside, initYside, -1.*box, "Left");
208   node->SetLineColor(kColorCRT);
209   fNodes->Add(node);
210
211   top->cd();
212   node = new TNode("upper8", "upper8", "S_CRT_A", initXside, initYside, -3.*box, "Left");
213   node->SetLineColor(kColorCRT);
214   fNodes->Add(node);
215
216
217   // Modules on the right side.
218   top->cd();
219   node = new TNode("upper9", "upper9", "S_CRT_A", -initXside, initYside,  3.*box, "Right");
220   node->SetLineColor(kColorCRT);
221   fNodes->Add(node);
222
223   top->cd();
224   node = new TNode("upper10", "upper10", "S_CRT_A", -initXside, initYside,    box, "Right");
225   node->SetLineColor(kColorCRT);
226   fNodes->Add(node);
227
228   top->cd();
229   node = new TNode("upper11","upper11", "S_CRT_A", -initXside, initYside, -1.*box, "Right");
230   node->SetLineColor(kColorCRT);
231   fNodes->Add(node);
232
233   top->cd();
234   node = new TNode("upper12","upper12", "S_CRT_A", -initXside, initYside, -3.*box, "Right");
235   node->SetLineColor(kColorCRT);
236   fNodes->Add(node);
237
238
239 }
240
241 //_____________________________________________________________________________
242 void AliCRTv0::CreateGeometry()
243 {
244   //
245   // Create geometry for the CRT array
246   //
247   Int_t  idrotm[2499];    // The rotation matrix.
248
249   Int_t * idtmed = fIdtmed->GetArray() - 1099 ;
250
251   //
252   // Molasse
253   CreateMolasse();
254
255   //
256   // Scintillators
257
258   Float_t box[3];
259   box[0] = AliCRTConstants::fgCageLenght/2.; // Half Length of the box along the X axis, cm.
260   box[1] = AliCRTConstants::fgCageHeight/2.; // Half Length of the box along the Y axis, cm.
261   box[2] = AliCRTConstants::fgCageWidth/2.;  // Half Length of the box along the Z axis, cm.
262
263
264   // Define the Scintillators. as a big box.
265   Float_t scint[3];
266   scint[0] = AliCRTConstants::fgActiveAreaLenght/2.;       // Half Length in X
267   scint[1] = AliCRTConstants::fgActiveAreaHeight/2.;       // Half Length in Y
268   scint[2] = AliCRTConstants::fgActiveAreaWidth/2.;        // Half Length in Z
269   gMC->Gsvolu("CRT1", "BOX ", idtmed[1112], scint, 3);     // Scintillators
270
271   //
272   // Define the coordinates where the draw will begin.
273   //
274
275   //
276   // -- X axis.
277   // we'll start dawing from the center.
278   Float_t initX = 0.;
279
280   //
281   // -- Y axis
282   Float_t gapY   = 30.;        // 30 cms. above the barrel.
283   // For the height we staimate the from the center of the ceiling,
284   // if were a cilinder, must be about 280cm.
285   Float_t barrel = 790.; // Barrel radius.
286   Float_t height  = barrel + gapY - 30.;
287   Float_t initY = height;
288
289   //
290   // -- Z axis.
291   // we'll start dawing from the center.
292
293   //
294   // Put 4 modules on the top of the magnet
295   Int_t step = 4;
296   for ( Int_t i = 1 ; i <= 4 ; i++ ) {
297     gMC->Gspos("CRT1", i, "ALIC", initX, initY, (i-step)*box[2], 0, "ONLY");
298     step--;
299   }
300
301   // Modules on the barrel sides.
302   // Because the openenig angle for each face is 22.5, and if we want to
303   //    put the modules right in the middle
304   Float_t xtragap = 10.;
305   Float_t initXside = (height+xtragap)*TMath::Sin(2*22.5*kDegrad); //rigth side
306   Float_t initYside = (height+xtragap)*TMath::Cos(2*22.5*kDegrad);
307
308   // Put 4 modules on the left side of the magnet
309   // The rotation matrix parameters, for the left side.
310   AliMatrix(idrotm[232], 90., 315., 90., 45., 0., 337.5);
311   Int_t stepl = 4;
312   for ( Int_t i = 1 ; i <= 4 ; i++ ) {
313     gMC->Gspos("CRT1", i+4, "ALIC", initXside, initYside, (i-stepl)*box[2],
314                idrotm[232], "ONLY");
315     stepl--;
316   }
317
318   // Put 4 modules on the right side of the magnet
319   // The rotation matrix parameters for the right side.
320   AliMatrix(idrotm[231], 90., 45., 90., 315., 180., 202.5);
321   Int_t stepr = 4;
322   for ( Int_t i = 1 ; i <= 4 ; i++ ) {
323     gMC->Gspos("CRT1", i+8, "ALIC", -initXside, initYside, (i-stepr)*box[2],
324                idrotm[231], "ONLY");
325     stepr--;
326   }
327
328   // Divide the modules in 2 planes.
329   //gMC->Gsdvn("CRT2", "CRT1", 2, 2);
330   // Now divide each plane in 8 palettes
331   //gMC->Gsdvn("CRT3", "CRT2", 8, 3);
332
333 }
334
335 //_____________________________________________________________________________
336 void AliCRTv0::CreateMolasse()
337 {
338   Int_t  idrotm[2499];    // The rotation matrix.
339
340   Int_t * idtmed = fIdtmed->GetArray() - 1099 ;
341
342   //
343   // Molasse
344   //
345
346   // Exactly above the hall
347   Float_t tspar[5];
348   tspar[0] = 1170.;
349   tspar[1] = 1170. + 375.;
350   tspar[2] = (1900.+1150.)/2.+100.;
351   tspar[3] = 0.;
352   tspar[4] = 180.;
353   gMC->Gsvolu("CMO1", "TUBS", idtmed[1123], tspar, 5);
354   gMC->Gspos("CMO1", 1, "ALIC", 0., 500., 1900.-tspar[2]+400., 0, "MANY");
355
356   Float_t tbox[3];
357   tbox[0] = 1250.;
358   tbox[1] = (4420. - 1670.)/2.;
359   tbox[2] = (1900.+1150.)/2. + 200.;
360   gMC->Gsvolu("CM12", "BOX", idtmed[1123], tbox, 3);
361   gMC->Gspos("CM12", 1, "ALIC", 0., 4420. -tbox[1], 1900.-tbox[2]+400., 0, "MANY");
362
363   AliMatrix(idrotm[2003], 0., 0., 90., 0., 90., 90.);
364   // Along the PM25
365   Float_t tube[3];
366   tube[0] = 455. + 100.;
367   tube[1] = 555. + 375.;
368   tube[2] = (5150. - 1166.)/2.;
369   gMC->Gsvolu("CMO2", "TUBE", idtmed[1123], tube, 3);
370   gMC->Gspos("CMO2", 1, "ALIC", -2100., 4420.-tube[2], 0., idrotm[2003], "MANY");
371
372
373   // Along the PGC2
374   tube[0] = 650.;
375   tube[1] = 2987.7;
376   tube[2] = (5150. - 690.)/2.;
377   gMC->Gsvolu("CMO3", "TUBE", idtmed[1123], tube, 3);
378   gMC->Gspos("CMO3", 1, "ALIC", 375., 4420.-tube[2], 1900.+2987.7, idrotm[2003], "MANY");
379   // Behind the PGC2 up to the end of the M. volume.
380   tbox[0] = 12073.;
381   tbox[1] = 2575. + 95.;
382   tbox[2] = (12073. - 1900.-2987.7-650.)/2.;
383   gMC->Gsvolu("CMO7", "BOX", idtmed[1123], tbox, 3);
384   gMC->Gspos("CMO7", 1, "ALIC", 0., 4420.-tbox[1], 1900.+2987.7+650.+tbox[2], 0, "MANY");
385
386   // Along the PX24 , upper part.
387   tube[0] = 1250.;
388   tube[1] = 2300;
389   tube[2] = 2575. - 1300. + 95.;
390   gMC->Gsvolu("CMO4", "TUBE", idtmed[1123], tube, 3);
391   gMC->Gspos("CMO4", 1, "ALIC", 0., 404.+1300.+tube[2], -2300., idrotm[2003], "MANY");
392
393   // Along the PX24 , lower part
394   tspar[0] = 1250.;
395   tspar[1] = 2300;
396   tspar[2] = 1300.;
397   tspar[3] = kRaddeg*TMath::ASin(1070./1150.);
398   tspar[4] = 360. - tspar[3];
399   gMC->Gsvolu("CMO5", "TUBS", idtmed[1123], tspar, 5);
400   gMC->Gspos("CMO5", 1, "ALIC", 0., 404., -2300., idrotm[2003], "MANY");
401   // behind the PX24
402   tbox[0] = 12073.;
403   tbox[1] = 2575. + 95.;
404   tbox[2] = 8523./2.;
405   gMC->Gsvolu("CMO6", "BOX", idtmed[1123], tbox, 3);
406   gMC->Gspos("CMO6", 1, "ALIC", 0., 4420.-tbox[1], -3550.-tbox[2], 0, "MANY");
407
408
409   // On the right side of th hall
410   tbox[0] = (12073. - 1250.)/2.;
411   tbox[1] = 2575. + 95.;
412   tbox[2] = (8437.7+650.)/2.;
413   gMC->Gsvolu("CMO8", "BOX", idtmed[1123], tbox, 3);
414   gMC->Gspos("CMO8", 1, "ALIC", 1250.+tbox[0], 4420.-tbox[1], -3550.+tbox[2], 0, "MANY");
415
416   // on the left side of the hall, behind 
417   tbox[0] = (12073. - 2755.)/2.;
418   tbox[1] = 2575. + 95.;
419   tbox[2] = (8437.7+650.)/2.;
420   gMC->Gsvolu("CMO9", "BOX", idtmed[1123], tbox, 3);
421   gMC->Gspos("CMO9", 1, "ALIC", -2755.-tbox[0], 4420.-tbox[1], -3550.+tbox[2], 0, "MANY");
422
423
424   // Molasse betwen the PX24 & PM25 on the left side.
425   tbox[0] = (2755. - 1250.)/2.;
426   tbox[1] = 2575. + 95.;
427   tbox[2] = (3550. - 555.)/2.;
428   gMC->Gsvolu("CM10", "BOX", idtmed[1123], tbox, 3);
429   gMC->Gspos("CM10", 1, "ALIC", -1250.-tbox[0], 4420.-tbox[1], -tbox[2]-555., 0, "MANY");
430
431
432   // Molasse betwen the PGC2 & PM25 on the left side.
433   tbox[0] = (2755. - 1250.)/2.;
434   tbox[1] = 2575. + 95.;
435   tbox[2] = (1900.+2987.7 - 555. + 650.)/2.;
436   gMC->Gsvolu("CM11", "BOX", idtmed[1123], tbox, 3);
437   gMC->Gspos("CM11", 1, "ALIC", -1250.-tbox[0], 4420.-tbox[1], 555.+tbox[2], 0, "MANY");
438
439
440 }
441
442 //_____________________________________________________________________________
443 void AliCRTv0::CreateShafts()
444 {
445   //
446   //
447   //
448   Int_t  idrotm[2499];    // The rotation matrix.
449
450   Int_t * idtmed = fIdtmed->GetArray() - 1099 ;
451
452   // HAll ceiling
453   Float_t ptubs[5];
454   ptubs[0] = 1070.;
455   ptubs[1] = 1170.;
456   ptubs[2] = 1900.;
457   ptubs[3] = 0.;
458   ptubs[4] = 180.;
459   gMC->Gsvolu("CHC1", "TUBS", idtmed[1116], ptubs, 5);
460   gMC->Gspos("CHC1", 1, "ALIC", 0., 500., 0., 0, "ONLY");
461
462
463   //
464   // Acces shafts
465   //
466   AliMatrix(idrotm[2001], 0., 0., 90., 0., 90., 90.);
467   
468   // PX24
469   ptubs[0] = 1150.;
470   ptubs[1] = 1250.;
471   ptubs[2] = 1300.;
472   ptubs[3] = kRaddeg*TMath::ASin(1070./ptubs[0]);
473   ptubs[4] = 360 - ptubs[3];
474   gMC->Gsvolu("CSF1", "TUBS", idtmed[1116], ptubs, 5);
475   gMC->Gspos("CSF1", 1, "ALIC", 0., 404., -2300., idrotm[2001], "MANY");
476
477   Float_t ptube[3];
478   ptube[0] = ptubs[0];
479   ptube[1] = ptubs[1];
480   ptube[2] = 2575. - ptubs[2] + 95.;
481   gMC->Gsvolu("CSF2", "TUBE", idtmed[1116], ptube, 3);
482   gMC->Gspos("CSF2", 1, "ALIC", 0., 404.+ptubs[2]+ptube[2], -2300., idrotm[2001], "MANY");
483   
484   // Concrete walls along the shaft
485   Float_t pbox[3];
486   pbox[0] = 585./2.;
487   pbox[1] = 2575. + 95.;
488   pbox[2] = 20.;
489   gMC->Gsvolu("CSW1", "BOX", idtmed[1116], pbox, 3);
490   gMC->Gspos("CSW1", 1, "ALIC", -290-pbox[0], 404.-1300.+pbox[1], -3450.+210.*2, 0, "MANY");
491   
492   //
493   pbox[0] = 750./2.;
494   pbox[1] = 2575. + 95.;
495   pbox[2] = 20.;
496   gMC->Gsvolu("CSW3", "BOX", idtmed[1116], pbox, 3);
497   gMC->Gspos("CSW3", 1, "ALIC", 420.-290.+pbox[0], 404.-1300.+pbox[1], -3450.+210.*2, 0, "MANY");
498   
499   //
500   pbox[0] = 60.;
501   pbox[1] = 2575. + 95.;
502   pbox[2] = 210.;
503   gMC->Gsvolu("CSW2", "BOX", idtmed[1116], pbox, 3);
504   gMC->Gspos("CSW2", 1, "ALIC", -290-pbox[0], 404.-1300.+pbox[1], -3450.+pbox[2], 0, "MANY");
505   gMC->Gspos("CSW2", 2, "ALIC", 420.-290.+pbox[0], 404.-1300.+pbox[1], -3450.+pbox[2], 0, "MANY");
506   
507   
508   // 
509   pbox[0] = 1000.;
510   pbox[1] = 80.;
511   pbox[2] = 200.;
512   gMC->Gsvolu("CSP1", "BOX", idtmed[1116], pbox, 3);
513   gMC->Gspos("CSP1", 1, "ALIC", 0., 2600.-700., -1150-pbox[2], 0, "MANY");
514   
515   //
516   pbox[0] = 340.8;
517   pbox[1] = 300./2.;
518   pbox[2] = 460./2.;
519   gMC->Gsvolu("CSP2", "BOX", idtmed[1116], pbox, 3);
520   gMC->Gspos("CSP2", 1, "ALIC", 0., 2950.-700., -3450+pbox[2], 0, "MANY");
521   
522   //
523   pbox[0] = 600.;
524   pbox[1] = 150.;
525   pbox[2] = 75.;
526   gMC->Gsvolu("CSP3", "BOX", idtmed[1116], pbox, 3);
527   gMC->Gspos("CSP3", 1, "ALIC", 0., 2950.-700., -1150.-210.-pbox[2], 0, "MANY");
528   
529   //
530   pbox[0] = 600.;
531   pbox[1] = 250.;
532   pbox[2] = 38.;
533   gMC->Gsvolu("CSP4", "BOX", idtmed[1116], pbox, 3);
534   gMC->Gspos("CSP4", 1, "ALIC", 0., 2950.-700.+155.+pbox[1], -1150.-210.-pbox[2], 0, "MANY");
535   
536   
537   // Shielding plug
538   pbox[0] = 850.;
539   pbox[1] = 90.;
540   pbox[2] = 720.;
541   gMC->Gsvolu("CSP5", "BOX", idtmed[1116], pbox, 3);
542   gMC->Gspos("CSP5", 1, "ALIC", 0., 2950.-700., -3450.+460.+pbox[2], 0, "MANY");
543   
544   //
545   pbox[0] = 80.;
546   pbox[1] = 150.;
547   pbox[2] = 720.;
548   gMC->Gsvolu("CSP6", "BOX", idtmed[1116], pbox, 3);
549   gMC->Gspos("CSP6", 1, "ALIC", 1150.-600., 2950.-700., -3450.+460.+pbox[2], 0, "MANY");
550   gMC->Gspos("CSP6", 2, "ALIC", -1150.+600., 2950.-700., -3450.+460.+pbox[2], 0, "MANY");
551   
552   
553   //
554   pbox[0] = 130.;
555   pbox[1] = 60.;
556   pbox[2] = 750.;
557   gMC->Gsvolu("CSP7", "BOX", idtmed[1116], pbox, 3);
558   gMC->Gspos("CSP7", 1, "ALIC", 850.+pbox[0], 2950.-700.+100., -3450.+460.+pbox[2], 0, "MANY");
559   gMC->Gspos("CSP7", 2, "ALIC", -850.-pbox[0], 2950.-700.+100., -3450.+460.+pbox[2], 0, "MANY");
560   
561   
562   // PM25 Acces Shaft
563   ptube[0] = 910./2.;
564   ptube[1] = ptube[0] + 100.;
565   ptube[2] = (5150. - 1166.)/2.;
566   gMC->Gsvolu("CSF3", "TUBE", idtmed[1116], ptube, 3);
567   gMC->Gspos("CSF3", 1, "ALIC", -2100., AliCRTConstants::fgDepth-ptube[2], 0., idrotm[2001], "MANY");
568   
569   // PGC2 Access Shaft
570   ptube[0] = 1100./2.;
571   ptube[1] = ptube[0] + 100.;
572   ptube[2] = (5150. - 690.)/2.;
573   gMC->Gsvolu("CSF4", "TUBE", idtmed[1116], ptube, 3);
574   gMC->Gspos("CSF4", 1, "ALIC", 375., AliCRTConstants::fgDepth-ptube[2], 1900.+2987.7, idrotm[2001], "MANY");
575
576 }
577
578 //_____________________________________________________________________________
579
580 void AliCRTv0::CreateMaterials()
581 {
582   // Use the standard materials.
583   AliCRT::CreateMaterials();  
584 }
585
586
587 //_____________________________________________________________________________
588 void AliCRTv0::DrawDetector()
589 {
590   //
591   // Draw a shaded view of the L3 magnet
592   //
593    cout << "AliCRTv0::DrawModule() : Drawing the module" << endl;
594
595    gMC->Gsatt("*", "seen", -1);
596    gMC->Gsatt("alic", "seen", 0);
597
598    gMC->Gsatt("ALIC","seen",0);
599    gMC->Gsatt("L3MO","seen",1); // L3 Magnet
600    gMC->Gsatt("CRT1","seen",1); // Scintillators
601
602    // Draw the molasse volumes
603    gMC->Gsatt("CMO1","seen",0); // Exactly above the HALL
604    gMC->Gsatt("CMO2","seen",0); // Molasse, along the PM25
605    gMC->Gsatt("CMO3","seen",0); // molasse along the PGC2
606    gMC->Gsatt("CMO4","seen",0); // Molasse, behind the PX24 upper part
607    gMC->Gsatt("CMO5","seen",0); // molasse behind px24, lower part
608    gMC->Gsatt("CMO6","seen",0); // behind the PX24
609    gMC->Gsatt("CMO7","seen",0); // behind the PGC2
610    gMC->Gsatt("CMO8","seen",0); // on the right side.
611    gMC->Gsatt("CMO9","seen",0); // on the left side.
612    gMC->Gsatt("CM10","seen",0); // betwen PX24 & PM25.
613    gMC->Gsatt("CM11","seen",0); // betwen PGC2 & PM25.
614    gMC->Gsatt("CM12","seen",0); // box above the hall.
615
616    gMC->Gdopt("hide", "on");
617    gMC->Gdopt("edge","off");
618    gMC->Gdopt("shad", "on");
619    gMC->Gsatt("*", "fill", 7);
620    gMC->SetClipBox("ALIC", 0, 3000, -3000, 3000, -6000, 6000);
621    gMC->DefaultRange();
622    gMC->Gdraw("alic", 40, 30, 0, 10, 9.5, .009, .009);
623    gMC->Gdhead(1111, "View of CRT(ACORDE)");
624    gMC->Gdman(18, 4, "MAN");
625
626
627 }
628
629 //_____________________________________________________________________________
630 void AliCRTv0::Init()
631 {
632   //
633   // Initialise L3 magnet after it has been built
634   Int_t i;
635   //
636   if(fDebug) {
637     printf("\n%s: ",ClassName());
638     for(i=0;i<35;i++) printf("*");
639     printf(" CRTv0_INIT ");
640     for(i=0;i<35;i++) printf("*");
641     printf("\n%s: ",ClassName());
642     //
643     // Here the CRTv0 initialisation code (if any!)
644     for(i=0;i<80;i++) printf("*");
645     printf("\n");
646   }
647
648 }
649
650 //_____________________________________________________________________________
651 void AliCRTv0::StepManager()
652 {
653   //
654   // Called for every step in the Cosmic Ray Trigger
655   //
656   static Int_t   vol[5];
657   Int_t          copy;
658   Int_t          ipart;
659   TLorentzVector pos;
660   TLorentzVector mom;
661
662   static Float_t hits[14];
663   Int_t tracknumber = gAlice->CurrentTrack();
664
665   static Float_t eloss;
666   static Float_t tlength;
667   Float_t theta;
668   Float_t phi;
669
670   if ( !gMC->IsTrackAlive() ) return;
671
672   if (gMC->IsNewTrack()) {
673     // Reset the deposited energy
674     eloss = 0.;
675   }
676   
677   eloss += gMC->Edep(); // Store the energy loss along the trajectory.
678   tlength += gMC->TrackStep();
679
680   if (gMC->IsTrackEntering() && (strcmp(gMC->CurrentVolName(),"CM12") == 0) ) {
681
682   // Get current particle id (ipart), track position (pos) and momentum (mom)
683     gMC->TrackPosition(pos);
684     gMC->TrackMomentum(mom);
685     ipart = gMC->TrackPid();
686
687     Double_t tc   = mom[0]*mom[0]+mom[1]*mom[1];
688     Double_t pt   = TMath::Sqrt(tc);
689     theta   = Float_t(TMath::ATan2(pt,Double_t(mom[2])))*kRaddeg;
690     phi     = Float_t(TMath::ATan2(Double_t(mom[1]),Double_t(mom[0])))*kRaddeg;
691
692
693     vol[0]    = gMC->CurrentVolOffID(1, vol[1]);
694     vol[2]    = gMC->CurrentVolID(copy);
695     vol[3]    = copy;
696     
697     hits[0]  = 0.f; //                 (fnmou)
698     hits[1]  = (Float_t)ipart; //      (fId)
699
700     hits[2]  = pos[0]; // X coordinate (fX)
701     hits[3]  = pos[1]; // Y coordinate (fY)
702     hits[4]  = pos[2]; // Z coordinate (fZ)
703     hits[5]  = mom[0]; // Px           (fpxug)
704     hits[6]  = mom[1]; // Py           (fpyug)
705     hits[7]  = mom[2]; // Pz           (fpzug)
706     
707     hits[8]  = gMC->GetMedium();//layer(flay)
708     hits[9]  = theta;   // arrival angle
709     hits[10] = phi;     // 
710     hits[11] = eloss;   // Energy loss
711     hits[12] = tlength; // Trajectory lenght
712     hits[13] = (Float_t)tracknumber;
713
714     AddHit(gAlice->CurrentTrack(),vol, hits);
715
716   }
717
718 }
719