]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ACORDE/AliACORDEv0.cxx
Fixed warnings
[u/mrichter/AliRoot.git] / ACORDE / AliACORDEv0.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // ALICE Cosmic Ray Trigger                                                  //
21 //                                                                           //
22 //  This class contains the functions for version 0 of the ALICE Cosmic Ray  //
23 //  Trigger. This version will be used to simulation comic rays in alice with//
24 //  all the detectors. It include geometry and hits (posicion and momentum)  //
25 //                                                                           // 
26 //   Author: Enrique Gamez                                                   //
27 //                                                                           //
28 //                  Send comments to:                                        //
29 //      Arturo Fernandez <afernand@fcfm.buap.mx>                             //
30 //      Eleazar Cuautle  <ecuautle@nucleares.unam.mx>                        //
31 //                                                                           //
32 //      Last update: Nov. 17th. 2009                                         //
33 //      Mario Rodriguez Cahuantzi <mrodrigu@mail.cern.ch                     //
34 //      FCFM-BUAP, Puebla, Pue. Mexico                                       //
35 //                                                                           //
36 ///////////////////////////////////////////////////////////////////////////////
37
38
39 #include "AliACORDEv0.h"
40 #include <TClonesArray.h>
41 #include <TLorentzVector.h>
42 #include <TVirtualMC.h>
43 #include <TPDGCode.h>
44
45
46 #include "AliRun.h"
47 #include "AliConst.h"
48 #include "AliACORDEhit.h"
49 #include "AliACORDEConstants.h"
50 #include "AliMC.h"
51 #include "AliLog.h"
52
53 ClassImp(AliACORDEv0)
54  
55 //_____________________________________________________________________________
56 AliACORDEv0::AliACORDEv0()
57   : AliACORDE()
58 {
59   //
60   // Default constructor
61   fIshunt = 0;
62   fHits = 0;
63   //
64
65 //_____________________________________________________________________________
66 AliACORDEv0::AliACORDEv0(const char *name, const char *title)
67   : AliACORDE(name, title)
68 {
69   //
70   // Standard constructor
71   //
72   fIshunt = 1; // All hits are associated with primary particles 
73   fHits =  new TClonesArray("AliACORDEhit",400);
74   gAlice->GetMCApp()->AddHitList(fHits);
75 }
76 //_____________________________________________________________________________
77 AliACORDEv0::~AliACORDEv0()
78 {
79   //
80   // Default destructor
81   //
82 }
83
84 //_____________________________________________________________________________
85 void AliACORDEv0::CreateGeometry()
86 {
87   CreateAcorde();
88   if (GetCreateCavern()) CreateCavern();
89 }
90
91 void AliACORDEv0::CreateCavern()
92 {
93   Int_t* idtmed = fIdtmed->GetArray() - 1099 ;
94     // Create the mother volume, the one which will contain all the material
95   // above the hall.
96   Float_t pbox[3];
97   pbox[0] = AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad);
98   //pbox[0] = 12073;
99   pbox[1] = AliACORDEConstants::Instance()->Depth();
100   pbox[2] = pbox[0];
101   gMC->Gsvolu("ACORDE", "BOX", idtmed[1114], pbox, 3);
102   gMC->Gspos("ACORDE", 1, "ALIC", 0, 0, 0, 0, "ONLY");
103   CreateShafts();
104   CreateMolasse();
105 }
106
107 void AliACORDEv0::CreateShafts()
108
109 {
110
111   //
112   Int_t  idrotm[2499];    // The rotation matrix.
113   Int_t* idtmed = fIdtmed->GetArray() - 1099 ;
114
115   //
116   // Acces shafts
117   //
118   AliMatrix(idrotm[2001], 0, 0, 90, 0, 90, 90);
119
120
121   // Create a bing cilinder to hold the main structures in the shaft.
122   //   All the structures relative to the shaft will be put into
123   // this volume.
124   //   This shaft is composed by an open tube down in the hall, and
125   // a cilinder avobe the level of the ceiling.
126   Float_t ptube[3];
127   ptube[0] = 0;    // inner radius
128   ptube[1] = 1250; // outer radius
129   ptube[2] = 5150/2; // Half lenght in Z
130   gMC->Gsvolu("CSF1", "TUBE", idtmed[1114], ptube, 3);
131
132   Float_t ptubs[5];
133   // The open section of the PX24
134   ptubs[0] = 1150; // Inner radius
135   ptubs[1] = 1250; // Outer radius
136   ptubs[2] = 1300; // Half length
137   ptubs[3] = 180 + kRaddeg*TMath::ASin(1070/ptubs[0]); // starting angle
138   ptubs[4] = 180 -  kRaddeg*TMath::ASin(1070/ptubs[0]);
139   gMC->Gsvolu("CSF2", "TUBS", idtmed[1116], ptubs, 5);
140   gMC->Gspos("CSF2", 1, "CSF1", 0, 0, -ptube[2] + ptubs[2], 0, "MANY");
141
142   // The other part of the shaft.
143   ptube[0] = ptubs[0]; // Inner radius
144   ptube[1] = ptubs[1]; // Outer radius
145   ptube[2] = 5150/2 - ptubs[2]; // Half lenght
146   gMC->Gsvolu("CSF3", "TUBE", idtmed[1116], ptube, 3);
147   gMC->Gspos("CSF3", 1, "CSF1", 0, 0, 5150/2 - ptube[2], 0, "MANY");
148
149   Float_t pbox[3];
150   // Concrete walls along the shaft (next to the elevator.)
151   pbox[0] = 480/2;  // Half length in X
152   pbox[1] = 120/2;  // Half length in Y
153   pbox[2] = 5150/2; // Half length in Z
154   gMC->Gsvolu("CSW1", "BOX", idtmed[1116], pbox, 3);
155   gMC->Gspos("CSW1", 1, "CSF1", 820+pbox[0],  150+pbox[1], 0, 0, "MANY");
156   gMC->Gspos("CSW1", 2, "CSF1", 820+pbox[0], -300-pbox[1], 0, 0, "MANY");
157
158   //
159   pbox[0] = 120/2;  // Half length in X
160   pbox[1] = 750/2;  // Half length in Y
161   pbox[2] = 5150/2; // Half length in Z
162   gMC->Gsvolu("CSW2", "BOX", idtmed[1116], pbox, 3);
163   gMC->Gspos("CSW2", 1, "CSF1", 820-60, 150+pbox[1], 0, 0, "MANY");
164
165   //
166   pbox[0] = 120/2;  // Half length in X
167   pbox[1] = 600/2;  // Half lenght in Y
168   pbox[2] = 5150/2; // Half length in Z
169   gMC->Gsvolu("CSW3", "BOX", idtmed[1116], pbox, 3);
170   gMC->Gspos("CSW3", 1, "CSF1", 820-60, -300-pbox[1], 0, 0, "MANY");
171
172   // Material below the counting rooms.
173   pbox[0] = 400/2;
174   pbox[1] = 2300/2;
175   pbox[2] = 300/2;
176   gMC->Gsvolu("CSW4", "BOX", idtmed[1116], pbox, 3);
177   gMC->Gspos("CSW4",1,"CSF1",2300/2-pbox[0],0,3000-5150/2-pbox[2], 0, "MANY");
178
179   // Shielding plug.
180   pbox[0] = 1400/2;
181   pbox[1] = 2300/2;
182   pbox[2] = 170/2;
183   gMC->Gsvolu("CSW5", "BOX", idtmed[1116], pbox, 3);
184   gMC->Gspos("CSW5", 1, "CSF1", 0, 0, 3000-5150/2-130, 0, "MANY");
185
186   // The end of the support for the shielding plug.
187   pbox[0] = 170/2;
188   pbox[1] = 2300/2;
189   pbox[2] = 300/2;
190   gMC->Gsvolu("CSW6", "BOX", idtmed[1116], pbox, 3);
191   gMC->Gspos("CSW6",1,"CSF1",-1400/2-pbox[0],0,3000-5150/2-pbox[2],0,"MANY");
192
193   // ...
194   pbox[0] = 100/2;
195   pbox[1] = 2300/2;
196   pbox[2] = 450/2;
197   gMC->Gsvolu("CSW7", "BOX", idtmed[1116], pbox, 3);
198   gMC->Gspos("CSW7",1,"CSF1",-1400/2-170-pbox[0],0,3000-5150/2+pbox[2],0,"MANY");
199
200   // Material close to the pipe.
201   pbox[0] = 300/2;
202   pbox[1] = 2300/2;
203   pbox[2] = 170/2;
204   gMC->Gsvolu("CSW8", "BOX", idtmed[1116], pbox, 3);
205   gMC->Gspos("CSW8",1,"CSF1",-2300/2+pbox[0],0,2500-5150/2,0,"MANY");
206
207   // Now put the shaft into the mother volume.
208   gMC->Gspos("CSF1", 1, "ACORDE", 0, AliACORDEConstants::Instance()->Depth() - 5150/2, 2300, idrotm[2001], "MANY");
209
210   // PM25 Access Shaft
211   ptube[0] = 910/2;
212   ptube[1] = ptube[0] + 100;
213   ptube[2] = (5150 - 1166)/2;
214   gMC->Gsvolu("CSF4", "TUBE", idtmed[1116], ptube, 3);
215   gMC->Gspos("CSF4", 1, "ACORDE", 2100, AliACORDEConstants::Instance()->Depth()-ptube[2], 0, idrotm[2001], "MANY");
216
217   // PGC2 Access Shaft
218   ptube[0] = 1100/2;
219   ptube[1] = ptube[0] + 100;
220   ptube[2] = (5150 - 690)/2;
221   gMC->Gsvolu("CSF5", "TUBE", idtmed[1116], ptube, 3);
222   gMC->Gspos("CSF5", 1, "ACORDE", -375, AliACORDEConstants::Instance()->Depth()-ptube[2], -1900 - 2987.7, idrotm[2001], "MANY");
223
224 }
225
226
227 void AliACORDEv0::CreateMolasse()
228
229 {
230
231   //
232   Int_t  idrotm[2499];    // The rotation matrix.
233   Int_t* idtmed = fIdtmed->GetArray() - 1099 ;
234
235   Float_t px24radius = 2300/2;
236   Float_t px24X = 0;
237   //Float_t px24Y = ;
238   Float_t px24Z = 2300;
239
240   Float_t pm25radius = 910/2;
241   Float_t pm25X = 2100;
242   //Float_t pm25Y = ;
243   Float_t pm25Z = 0;
244
245   Float_t pgc2radius = 1100/2;
246   Float_t pgc2X = -375;
247   //Float_t pgc2Y = ;
248   Float_t pgc2Z = -(1900 + 2987.7);
249
250   Float_t concreteWidth = 100; // Standard width of the hall walls.
251
252
253   // Create a local mother volume.
254   Float_t pbox[3];
255   pbox[0] = AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad);
256   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
257   pbox[2] = pbox[0];
258   gMC->Gsvolu("CMO1", "BOX", idtmed[1114], pbox, 3);
259
260   // Now put the molasse exactly above the hall. OK
261   // Above the ceiling
262   Float_t ptubs[5];
263   ptubs[0] = 1170;
264   ptubs[1] = 2100 - pm25radius;
265   ptubs[2] = 1900/2 + px24radius;
266   ptubs[3] = 0;
267   ptubs[4] = 180;
268   gMC->Gsvolu("CMO2", "TUBS", idtmed[1123], ptubs, 5);
269   gMC->Gspos("CMO2", 1, "CMO1", 0, 500-AliACORDEConstants::Instance()->Depth()/2, ptubs[2]-1900, 0, "MANY");
270
271   // Molasse around the RB24/26 Wall. OK
272   ptubs[0] = 220 + 1600;
273   ptubs[1] = AliACORDEConstants::Instance()->Depth() - ptubs[0];
274   ptubs[2] = 2987.7/2 - 1100/4 - concreteWidth/2;
275   ptubs[3] = 0;
276   ptubs[4] = 180;
277   gMC->Gsvolu("CMO3", "TUBS", idtmed[1123], ptubs, 5);
278   gMC->Gspos("CMO3", 1, "CMO1", 70, 40-AliACORDEConstants::Instance()->Depth()/2, -1900 - ptubs[2], 0, "MANY");
279
280   // A big block above the RB24/26 wall. OK
281   pbox[0] = AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad);
282   pbox[1] = (AliACORDEConstants::Instance()->Depth() - 220 - 1600)/2;
283   pbox[2] = 2987.7/2 - 1100/4 - concreteWidth/2;
284   gMC->Gsvolu("CMO4", "BOX", idtmed[1123], pbox, 3);
285   gMC->Gspos("CMO4", 1, "CMO1", 0, AliACORDEConstants::Instance()->Depth()/2 - pbox[1], -1900 - pbox[2], 0, "MANY");
286   // Small blocks below the volume CMO4 on both sides of the wall RB24/26. OK
287   pbox[0] = (AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) -
288 ptubs[0])/2;
289   pbox[1] = AliACORDEConstants::Instance()->Depth()/2 - pbox[1];
290   gMC->Gsvolu("CM17", "BOX", idtmed[1123], pbox, 3);
291   gMC->Gspos("CM17", 1, "CMO1", AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) - pbox[0], -AliACORDEConstants::Instance()->Depth()/2 + pbox[1], -1900 - pbox[2], 0, "MANY");
292   gMC->Gspos("CM17", 2, "CMO1", -AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad)+ pbox[0], -AliACORDEConstants::Instance()->Depth()/2 + pbox[1], -1900 - pbox[2], 0, "MANY");
293
294   // And a big block of molasse above the hall up to the surface. OK
295   pbox[0] = pm25X - pm25radius;
296   pbox[1] = (AliACORDEConstants::Instance()->Depth()-500-1170)/2;
297   pbox[2] = (1900 + 1150)/2;
298   gMC->Gsvolu("CMO5", "BOX", idtmed[1123], pbox, 3);
299   gMC->Gspos("CMO5", 1, "CMO1", 0,AliACORDEConstants::Instance()->Depth()/2-pbox[1], pbox[2]-1900, 0, "MANY");
300   // Small blocks of molasse betwen the blocks CMO2, CMO5 and PM25. Ok
301   pbox[0] = (pm25X - pm25radius - 1170)/2;
302   pbox[1] = 1000;
303   gMC->Gsvolu("CM16", "BOX", idtmed[1123], pbox, 3);
304   gMC->Gspos("CM16", 1, "CMO1", 1170 + pbox[0], -AliACORDEConstants::Instance()->Depth()/2+pbox[1], pbox[2] - 1900, 0, "MANY");
305
306   // Molasse around the shafts.
307   AliMatrix(idrotm[2003], 0, 0, 90, 0, 90, 90);
308   // Around the PX24, the open section. OK
309   ptubs[0] = px24radius + concreteWidth;
310   ptubs[1] = ptubs[0] + 1000;
311   ptubs[2] = (2300 - (5150 - AliACORDEConstants::Instance()->Depth()))/2;
312   ptubs[3] = 180 + kRaddeg*TMath::ASin(1070/ptubs[0]);
313   ptubs[4] = 180 -  kRaddeg*TMath::ASin(1070/ptubs[0]);
314   gMC->Gsvolu("CMO6", "TUBS", idtmed[1123], ptubs, 5);
315   gMC->Gspos("CMO6", 1, "CMO1", px24X, ptubs[2] - AliACORDEConstants::Instance()->Depth()/2, px24Z, idrotm[2003], "MANY");
316   // Around the PX24, the closed section. OK
317   Float_t ptube[3];
318   ptube[0] = px24radius + concreteWidth;
319   ptube[1] = ptube[0] + 1000;
320   ptube[2] = (5150 - 2300)/2;
321   gMC->Gsvolu("CMO7", "TUBE", idtmed[1123], ptube, 3);
322   gMC->Gspos("CMO7", 1, "CMO1", px24X, AliACORDEConstants::Instance()->Depth()/2 - ptube[2], px24Z, idrotm[2003], "MANY");
323
324   // Around PM25. OK
325   ptube[0] = pm25radius + concreteWidth;
326   ptube[1] = ptube[0] + 400;
327   ptube[2] = AliACORDEConstants::Instance()->Depth()/2;
328   gMC->Gsvolu("CMO8", "TUBE", idtmed[1123], ptube, 3);
329   gMC->Gspos("CMO8", 1, "CMO1", pm25X, 0, pm25Z, idrotm[2003], "MANY");
330   // On both sides of the PM25 along the HALL.
331   pbox[0] = (2100 + pm25radius - 1170)/2;
332   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
333   pbox[2] = (3*px24radius - pm25radius)/2;
334   gMC->Gsvolu("CM18", "BOX", idtmed[1123], pbox, 3);
335   gMC->Gspos("CM18", 1, "CMO1", 2100, 0, pbox[2] + pm25radius, 0, "MANY");
336
337   pbox[2] = (1900 - pm25radius)/2;
338   gMC->Gsvolu("CM19", "BOX", idtmed[1123], pbox, 3);
339   gMC->Gspos("CM19", 1, "CMO1", 2100, 0, -pbox[2] - pm25radius, 0, "MANY");
340
341   // Around the PGC2. OK
342   ptube[0] = pgc2radius + concreteWidth;
343   ptube[1] = 2987.7 - 740;
344   ptube[2] = AliACORDEConstants::Instance()->Depth()/2;
345   gMC->Gsvolu("CMO9", "TUBE", idtmed[1123], ptube, 3);
346   gMC->Gspos("CMO9", 1, "CMO1", pgc2X, 0, pgc2Z, idrotm[2003], "MANY");
347
348   // On both sides of the PGC2.OK
349   pbox[0] = (AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) -
350 1100 - 375)/2;
351   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
352   pbox[2] = pgc2radius + concreteWidth;
353   gMC->Gsvolu("CM10", "BOX", idtmed[1123], pbox, 3);
354   gMC->Gspos("CM10", 1, "CMO1", AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) - pbox[0], 0, pgc2Z, 0, "MANY");
355   gMC->Gspos("CM10", 2, "CMO1", -AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) + pbox[0], 0, pgc2Z, 0, "MANY");
356
357   // big block of molasse behind the PX24. OK
358   pbox[0] = AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad);
359   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
360   pbox[2] = (pbox[0] - (2300 + 1150 + 100))/2;
361   gMC->Gsvolu("CM12", "BOX", idtmed[1123], pbox, 3);
362   gMC->Gspos("CM12", 1, "CMO1", px24X, 0, px24Z + px24radius + concreteWidth + pbox[2], 0, "MANY");
363
364   // big block of molasse in the opposite side of the PM25. OK
365   pbox[0] = (AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) -
366 1150)/2;
367   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
368   pbox[2] = (1900 + 2300 + 1150)/2;
369   gMC->Gsvolu("CM13", "BOX", idtmed[1123], pbox, 3);
370   gMC->Gspos("CM13", 1, "CMO1", -1150 - pbox[0], 0, pbox[2] - 1900, 0, "MANY");
371
372   // big block of molasse behind the PM25. OK
373   pbox[0] = (AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad) -
374 (2100 + 910/2 + 100))/2;
375   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
376   pbox[2] = (1900 + 2300 + 1150)/2;
377   gMC->Gsvolu("CM14", "BOX", idtmed[1123], pbox, 3);
378   gMC->Gspos("CM14", 1, "CMO1", pm25X + pm25radius + concreteWidth + pbox[0], 0, pbox[2] - 1900, 0, "MANY");
379
380   // big block of molasse behind the PGC2. OK
381   pbox[0] = AliACORDEConstants::Instance()->Depth()*TMath::Tan(67.5*kDegrad);
382   pbox[1] = AliACORDEConstants::Instance()->Depth()/2;
383   pbox[2] = (pbox[0] - (2987.7 + 1900 + 1100/2 + 100))/2;
384   gMC->Gsvolu("CM15", "BOX", idtmed[1123], pbox, 3);
385   gMC->Gspos("CM15", 1, "CMO1", 0, 0, -pbox[0] + pbox[2], 0, "MANY");
386
387   gMC->Gspos("CMO1",1,"ACORDE",0,AliACORDEConstants::Instance()->Depth()/2,0,0,"MANY");
388
389 }
390
391 void AliACORDEv0::CreateAcorde()
392 {
393   //
394   // Create geometry for the ACORDE array
395   // done in two main steps
396   //  1.- definition of the modules
397   //  2.- placement of the modules
398   //
399   Int_t  idrotm[2499];    // The rotation matrix.
400   Int_t* idtmed = fIdtmed->GetArray() - 1099;
401   AliACORDEConstants* constants = AliACORDEConstants::Instance();
402   Float_t box[3];
403   Float_t placed_at;
404   Float_t placed_at2;
405   Float_t small = 0.05; // to separate slightly some volumes
406                         // by half a mm so that they do not overlap
407
408
409   // 1.- Definition of a module
410   // *  ACORDE1 => volume filled with air, representing a module
411   //               it contains all other volumes defining the module
412   //               there are 60 copies of it
413   // *  ACORDE2 => volume defining one scintillator pad
414   //               there are 2 copies of it per module
415   // *  ACORDE3-6 => volumes representing the Al walls of box
416   //               surrounding the plastic
417   //               3: long wall, 2 copies (front, back)
418   //               4: end caps, 2 copies (left, right)
419   //               5: long stripe to model the profile 
420   //                  4 copies (upper front and back, lower)
421   //               6: short stripe to model the profile
422   //                  4 copies (upper left, right; lower)
423
424   // The full module volume.
425   // This volume will be ocupied by all the material of the module
426   // the scintillators, the aluminium frame, etc.
427   box[0] = constants->ModuleLength()/2;
428   box[1] = constants->ModuleHeight()/2;
429   box[2] = constants->ModuleWidth()/2;
430   gMC->Gsvolu("ACORDE1", "BOX", idtmed[1114], box, 3);
431
432   // The scintillators
433   box[0] = constants->PlasticLength()/2;
434   box[1] = constants->PlasticHeight()/2;
435   box[2] = constants->PlasticWidth()/2;
436   gMC->Gsvolu("ACORDE2", "BOX", idtmed[1112], box, 3);
437
438   // it is important to keep this order for easy assignment of 
439   // a volume to a physical module:
440   placed_at = box[1]+constants->ProfileThickness()
441     - constants->ModuleHeight()/2+small;
442   gMC->Gspos("ACORDE2", 1, "ACORDE1", 0, placed_at, 0, 0, "MANY");
443   placed_at = placed_at + 2.0*box[1]+small;
444   gMC->Gspos("ACORDE2", 2, "ACORDE1", 0, placed_at, 0, 0, "MANY");
445
446
447   // The metallic frame: long walls of box
448   // back,front,left,right, defined looking
449   // from the + z diraction into alice; i.e.
450   // back ==> z<0, front ==> z>0
451   // left ==> x<0, right ==> x>0
452   // up ==> increasing y, down ==> decreasing y
453   box[0] = constants->ModuleLength()/2;
454   box[1] = constants->ModuleHeight()/2;
455   box[2] = constants->ProfileThickness()/2.0; 
456   gMC->Gsvolu("ACORDE3", "BOX", idtmed[1108], box, 3);
457   // front wall
458   placed_at = constants->ModuleWidth()/2-constants->ProfileThickness()/2.0;
459   gMC->Gspos("ACORDE3", 1, "ACORDE1", 0, 0, placed_at, 0, "MANY");
460   // back wall
461   gMC->Gspos("ACORDE3", 2, "ACORDE1", 0, 0, -placed_at , 0, "MANY");
462
463   // The metallic frame: end caps
464   box[0] = constants->ProfileThickness()/2.0;
465   box[1] = constants->ModuleHeight()/2;
466   box[2] = constants->ModuleWidth()/2;
467   gMC->Gsvolu("ACORDE4", "BOX", idtmed[1108], box, 3);
468   // right cap
469   placed_at = constants->ModuleLength()/2-constants->ProfileThickness()/2.0;
470   gMC->Gspos("ACORDE4", 1, "ACORDE1", placed_at, 0, 0, 0, "MANY");
471   // left cap
472   gMC->Gspos("ACORDE4", 2, "ACORDE1", -placed_at, 0, 0, 0, "MANY");
473
474   // The metallic frame: the profile, long stripes
475   box[0] = constants->ModuleLength()/2.0;
476   box[1] = constants->ProfileThickness()/2;
477   box[2] = constants->ProfileWidth()/2;
478   gMC->Gsvolu("ACORDE5", "BOX", idtmed[1108], box, 3);
479   // upper front
480   placed_at = constants->ModuleHeight()/2-box[1];
481   placed_at2 = constants->ModuleWidth()/2-
482     constants->ProfileThickness()-box[2];
483   gMC->Gspos("ACORDE5", 1, "ACORDE1",0,placed_at,placed_at2, 0, "MANY");
484   // upper back
485   gMC->Gspos("ACORDE5", 2, "ACORDE1",0,placed_at,-placed_at2, 0, "MANY");
486   // lower front
487   gMC->Gspos("ACORDE5", 3, "ACORDE1",0,-placed_at,placed_at2, 0, "MANY");
488   // lower back
489   gMC->Gspos("ACORDE5", 4, "ACORDE1",0,-placed_at,-placed_at2, 0, "MANY");
490
491   // The metallic frame: the profile, long stripes
492   box[0] = constants->ProfileWidth()/2.0;
493   box[1] = constants->ProfileThickness()/2;
494   box[2] = constants->ModuleWidth()/2-constants->ProfileWidth();
495   gMC->Gsvolu("ACORDE6", "BOX", idtmed[1108], box, 3);
496   // upper right
497   placed_at = constants->ModuleHeight()/2-box[1];
498   placed_at2 = constants->ModuleLength()/2-
499     constants->ProfileThickness()-box[0];
500   gMC->Gspos("ACORDE6", 1, "ACORDE1",placed_at2,placed_at,0, 0, "MANY");
501   // upper left
502   gMC->Gspos("ACORDE6", 2, "ACORDE1",-placed_at2,placed_at,0, 0, "MANY");
503   // lower right
504   gMC->Gspos("ACORDE6", 3, "ACORDE1",placed_at2,-placed_at,0, 0, "MANY");
505   // lower left
506   gMC->Gspos("ACORDE6", 4, "ACORDE1",-placed_at2,-placed_at,0, 0, "MANY");
507
508   // End of MODULE definition
509
510   ////////////////////////////////////////////////////////////////////
511   ////////////////////////////////////////////////////////////////////
512
513   // 2.- placement of the module
514   // Now put all of them in the right position in 
515   // master volume ALIC
516
517   // rotation matrices (see Geant manual for conventions)
518   // for columns 4 and 5
519   AliMatrix(idrotm[231], 90, 45, 90, 135, 0, 0);
520   // for columns 0 and 1
521   AliMatrix(idrotm[232], 90, 315, 90, 45, 0, 0);
522
523   // place each one of the 6 columns in turn
524   // for the first and the last column the position
525   // of the two last modules depends on the value 
526   // of the fITSGeometry variable
527
528   // it is important to keep this order because
529   // the copy number defines the module!
530
531   // first column, except first and last  modules
532   for (Int_t copy = 2; copy < 10; copy++)
533     gMC->Gspos("ACORDE1",copy,"ALIC",
534                constants->OldModulePositionX(copy-1),
535                constants->OldModulePositionY(copy-1),
536                constants->OldModulePositionZ(copy-1),
537                idrotm[232], "MANY");
538   // second column
539   for (Int_t copy = 11; copy < 21; copy++)
540     gMC->Gspos("ACORDE1",copy,"ALIC",
541                constants->OldModulePositionX(copy-1),
542                constants->OldModulePositionY(copy-1),
543                constants->OldModulePositionZ(copy-1),
544                idrotm[232], "MANY");
545   // third and fourth columns
546   for (Int_t copy = 21; copy < 41; copy++)
547     gMC->Gspos("ACORDE1",copy,"ALIC",
548                constants->OldModulePositionX(copy-1),
549                constants->OldModulePositionY(copy-1),
550                constants->OldModulePositionZ(copy-1),
551                0, "MANY");
552   // fifth column
553   for (Int_t copy = 41; copy < 51; copy++)
554     gMC->Gspos("ACORDE1",copy,"ALIC",
555                constants->OldModulePositionX(copy-1),
556                constants->OldModulePositionY(copy-1),
557                constants->OldModulePositionZ(copy-1),
558                idrotm[231], "MANY");
559   // last column, except first and last  modules
560   for (Int_t copy = 52; copy < 60; copy++)
561     gMC->Gspos("ACORDE1",copy,"ALIC",
562                constants->OldModulePositionX(copy-1),
563                constants->OldModulePositionY(copy-1),
564                constants->OldModulePositionZ(copy-1),
565                idrotm[231], "MANY");
566   // the last four modules
567   if (Get4CentralModulesGeometry()) {
568     gMC->Gspos("ACORDE1",1,"ALIC",
569                constants->OldExtraModulePositionX(),
570                constants->OldExtraModulePositionY(),
571                constants->OldExtraModulePositionZ(0),
572                0, "MANY");  
573     gMC->Gspos("ACORDE1",10,"ALIC",
574                constants->OldExtraModulePositionX(),
575                constants->OldExtraModulePositionY(),
576                constants->OldExtraModulePositionZ(1),
577                0, "MANY");  
578     gMC->Gspos("ACORDE1",51,"ALIC",
579                constants->OldExtraModulePositionX(),
580                constants->OldExtraModulePositionY(),
581                constants->OldExtraModulePositionZ(2),
582                0, "MANY");  
583     gMC->Gspos("ACORDE1",60,"ALIC",
584                constants->OldExtraModulePositionX(),
585                constants->OldExtraModulePositionY(),
586                constants->OldExtraModulePositionZ(3),
587                0, "MANY");  
588   } else {
589     gMC->Gspos("ACORDE1",1,"ALIC",
590                constants->OldModulePositionX(0),
591                constants->OldModulePositionY(0),
592                constants->OldModulePositionZ(0),
593                idrotm[232], "MANY");
594     gMC->Gspos("ACORDE1",10,"ALIC",
595                constants->OldModulePositionX(9),
596                constants->OldModulePositionY(9),
597                constants->OldModulePositionZ(9),
598                idrotm[232], "MANY");
599     gMC->Gspos("ACORDE1",51,"ALIC",
600                constants->OldModulePositionX(50),
601                constants->OldModulePositionY(50),
602                constants->OldModulePositionZ(50),
603                idrotm[231], "MANY");
604     gMC->Gspos("ACORDE1",60,"ALIC",
605                constants->OldModulePositionX(59),
606                constants->OldModulePositionY(59),
607                constants->OldModulePositionZ(59),
608                idrotm[231], "MANY");
609   } // end if (fITSGeometry)
610
611 }
612 //_____________________________________________________________________________
613 void AliACORDEv0::DrawDetector() const
614 {
615
616   // not needed anymore
617
618 }
619
620 //____________________________________________________________________________
621
622 void AliACORDEv0::Init()
623 {
624   // Initialise L3 magnet after it has been built
625   Int_t i;
626   if(AliLog::GetGlobalDebugLevel()>0) {
627     printf("\n%s: ",ClassName());
628     for(i=0;i<35;i++) printf("*");
629     printf(" ACORDEv0_INIT ");
630     for(i=0;i<35;i++) printf("*");
631     printf("\n%s: ",ClassName());
632     // Here the ACORDEv initialisation code (if any!)
633     for(i=0;i<80;i++) printf("*");
634     printf("\n");
635   }
636  // AliACORDE::Init();  
637 }
638 //____________________________________________________________________________
639 void AliACORDEv0::StepManager()
640 {
641   //
642   // Called for every step in the Cosmic Ray Trigger
643   //
644
645
646   // volume: 
647   //  [0] = module number 1-60 (1==>(0-0), 60 (5-9)
648   //  [1] = Plastic number: 0 (down) to 1 (up)
649   static Int_t   vol[2]; 
650   //
651   // hit
652   // [0] = PID
653   // [1-3] = x, y, z 
654   // [4] = time 
655   // [5-7] = px, py, pz
656   // [8] = energy 
657   // [9] = energy loss
658   // [10] = length of track through plastic
659   static Float_t hits[11];
660
661   // local static variables
662   static Float_t eloss;
663   static Float_t step;
664   // scintillator volume
665   static Int_t idScint = gMC->VolId("ACORDE2");
666
667   // local variables
668   Int_t copy;
669   TLorentzVector pos;
670   TLorentzVector mom;
671
672   // only charged tracks
673   if ( !gMC->TrackCharge() || !gMC->IsTrackAlive() ) return;
674
675   // only in sensitive material
676   if (gMC->CurrentVolID(copy) == idScint) {
677     step  += gMC->TrackStep();
678     eloss += gMC->Edep();
679     // set all hit variables except eloss which is resetted
680     // set volume variables
681     if (gMC->IsTrackEntering()) {
682       eloss = 0.0;
683       step = 0.0;
684       gMC->TrackPosition(pos);
685       gMC->TrackMomentum(mom);
686       // hit
687       // [0] = PID
688       // [1-3] = x, y, z 
689       // [4] = time 
690       // [5-7] = px, py, pz
691       // [8] = energy 
692       // [9] = energy loss
693       hits[0]  = (Float_t ) gMC->TrackPid(); 
694       hits[1] = pos[0]; 
695       hits[2] = pos[1]; 
696       hits[3] = pos[2]; 
697       hits[4] = gMC->TrackTime();
698       hits[5] = mom[0]; 
699       hits[6] = mom[1]; 
700       hits[7] = mom[2]; 
701       hits[8] = gMC->Etot();
702       // volume: 
703       //  [0] = module number 1-60 (1==>(0-0), 60 (5-9)
704       //  [1] = Plastic number: 0 (down) to 1 (up)
705       Int_t copyPlastic; // plastic: down=1, up=2
706       Int_t copyModule; // module: 1-60
707       gMC->CurrentVolID(copyPlastic);
708       gMC->CurrentVolOffID(1, copyModule);
709       // module
710       vol[0] = copyModule;
711       // plastic: 0 = down, 1 = up
712       vol[1] = copyPlastic;
713     } // end if gMC->IsTrackEntering()
714
715     // set hit[9] = total energy loss and book hit
716     if( gMC->IsTrackExiting() || 
717         gMC->IsTrackStop() || 
718         gMC->IsTrackDisappeared()){
719       hits[9] = eloss;
720       hits[10] = step;
721       eloss = 0.0;
722       step = 0.0;
723       AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber(),vol, hits);
724     }
725   } // end if in scintillator
726
727 }
728
729 //_____________________________________________________________________________
730 void AliACORDEv0::AddHit(Int_t track, Int_t *vol, Float_t *hits)
731 {
732   //
733   // Add a ACORDE hit
734   //
735   TClonesArray &lhits = *fHits;
736   new(lhits[fNhits++]) AliACORDEhit(fIshunt,track,vol,hits);
737 }
738