]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROv7.cxx
Revision of V0A by Lizardo Valencia and Andreas
[u/mrichter/AliRoot.git] / VZERO / AliVZEROv7.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 //  (V-zero) detector  version 7 as designed by the Lyon and         //
21 //   Mexico groups and Carlos Perez Lara from Pontificia Universidad //
22 //   Catolica del Peru                                               // 
23 //   All comments should be sent to Brigitte CHEYNIS:                //
24 //                     b.cheynis@ipnl.in2p3.fr                       // 
25 //   Geometry of April 2006 done with ROOT geometrical modeler       //
26 //   V0R (now V0C) sits between Z values  -89.5 and  -84.8 cm        //
27 //   V0L (now V0A) sits between Z values +338.5 and +342.5 cm        // 
28 //   New coordinate system has been implemented in october 2003      //
29 //   Revision of the V0A part by Lizardo Valencia  in July 2008      //
30 //                                                                   //
31 /////////////////////////////////////////////////////////////////////// 
32
33 // --- Standard libraries ---
34 #include <Riostream.h>
35
36 // --- ROOT libraries ---
37 #include <TClonesArray.h>
38 #include <TMath.h>
39 #include <TVirtualMC.h>
40 #include <TParticle.h>
41
42 #include <TGeoManager.h>
43 #include <TGeoMatrix.h>
44 #include <TGeoMaterial.h>
45 #include <TGeoMedium.h>
46 #include <TGeoVolume.h>
47 #include "TGeoTube.h"
48 #include "TGeoArb8.h"
49 #include "TGeoCompositeShape.h"
50
51 // --- AliRoot header files ---
52 #include "AliRun.h"
53 #include "AliMC.h"
54 #include "AliMagF.h"
55 #include "AliVZEROLoader.h"
56 #include "AliVZEROdigit.h"
57 #include "AliVZEROhit.h"
58 #include "AliVZEROv7.h"
59 #include "AliLog.h"
60  
61 ClassImp(AliVZEROv7)
62
63 //_____________________________________________________________________________
64 AliVZEROv7:: AliVZEROv7():AliVZERO(),
65    fCellId(0),
66    fTrackPosition(),
67    fTrackMomentum(), 
68    fV0CHeight1(2.5), 
69    fV0CHeight2(4.4), 
70    fV0CHeight3(7.4), 
71    fV0CHeight4(12.5),
72    fV0CRMin(4.6), 
73    fV0CRBox(38.0),
74    fV0CLidThickness(0.30),
75    fV0CCellThickness(2.00),
76    fV0CBoxThickness(4.70),
77    fV0COffsetFibers(1.0),
78    fV0CLightYield(93.75),
79    fV0CLightAttenuation(0.05),
80    fV0CnMeters(15.0),
81    fV0CFibToPhot(0.3),
82    fV0AR0(4.2), 
83    fV0AR1(7.6), 
84    fV0AR2(13.8), 
85    fV0AR3(22.7),
86    fV0AR4(41.3), 
87    fV0AR5(43.3), 
88    fV0AR6(72.6),
89    fV0AR7(92.0), // Distance from origin to outtermost intersection sector7 and sector8 
90    fV0ASciWd(2.5), 
91    fV0APlaWd(0.5), 
92    fV0APlaAl(0.06), 
93    fV0AOctWd(0.75), 
94    fV0AFraWd(0.2),
95    fV0AOctH1(1.0), 
96    fV0AOctH2(2.0), 
97    fV0ABasHt(2.0),
98    fV0AFibRd(0.1),
99    fV0APlaEx(4.4),
100    fV0APMBWd(24.6), 
101    fV0APMBHt(22.0), 
102    fV0APMBTh(7.1), 
103    fV0APMBWdW(0.3), 
104    fV0APMBHtW(1.0),
105    fV0APMBAng(30.0), 
106    fV0APMBThW(0.3), 
107    fV0APMTR1(2.44), 
108    fV0APMTR2(2.54), 
109    fV0APMTR3(2.54),
110    fV0APMTR4(2.70), 
111    fV0APMTH(10.0), 
112    fV0APMTB(1.0),
113    fV0AnMeters(fV0AR6*0.01),
114    fV0ALightYield(93.75),
115    fV0ALightAttenuation(0.05),
116    fV0AFibToPhot(0.3),
117    fVersion(7)
118 {
119 // Standard default constructor 
120 }
121
122 //_____________________________________________________________________________
123 AliVZEROv7::AliVZEROv7(const char *name, const char *title):AliVZERO(name,title),
124    fCellId(0),
125    fTrackPosition(),
126    fTrackMomentum(), 
127    fV0CHeight1(2.5), 
128    fV0CHeight2(4.4), 
129    fV0CHeight3(7.4), 
130    fV0CHeight4(12.5),
131    fV0CRMin(4.6), 
132    fV0CRBox(38.0),
133    fV0CLidThickness(0.30),
134    fV0CCellThickness(2.00),
135    fV0CBoxThickness(4.70),
136    fV0COffsetFibers(1.0),
137    fV0CLightYield(93.75),
138    fV0CLightAttenuation(0.05),
139    fV0CnMeters(15.0),
140    fV0CFibToPhot(0.3),
141    fV0AR0(4.2), 
142    fV0AR1(7.6), 
143    fV0AR2(13.8), 
144    fV0AR3(22.7),
145    fV0AR4(41.3), 
146    fV0AR5(43.3), 
147    fV0AR6(72.6),
148    fV0AR7(92.0), // Distance from origin to outtermost intersection of sector7 and sector8 
149    fV0ASciWd(2.5), 
150    fV0APlaWd(0.5), 
151    fV0APlaAl(0.06), 
152    fV0AOctWd(0.75), 
153    fV0AFraWd(0.2),
154    fV0AOctH1(1.0), 
155    fV0AOctH2(2.0), 
156    fV0ABasHt(2.0),
157    fV0AFibRd(0.1),
158    fV0APlaEx(4.4),
159    fV0APMBWd(24.6), 
160    fV0APMBHt(22.0), 
161    fV0APMBTh(7.1), 
162    fV0APMBWdW(0.3), 
163    fV0APMBHtW(1.0),
164    fV0APMBAng(30.0), 
165    fV0APMBThW(0.3), 
166    fV0APMTR1(2.44), 
167    fV0APMTR2(2.54), 
168    fV0APMTR3(2.54),
169    fV0APMTR4(2.70), 
170    fV0APMTH(10.0), 
171    fV0APMTB(1.0),
172    fV0AnMeters(fV0AR6*0.01),
173    fV0ALightYield(93.75),
174    fV0ALightAttenuation(0.05),
175    fV0AFibToPhot(0.3),
176    fVersion(7)
177
178
179 {
180 // Standard constructor for V-zero Detector  version 7
181
182   AliDebug(2,"Create VZERO object ");
183
184 //  fVersion            =     7;  // version number
185
186 //   // V0C Parameters related to geometry: All in cm
187 //   fV0CHeight1         =    2.5; // height of cell 1
188 //   fV0CHeight2         =    4.4; // height of cell 2
189 //   fV0CHeight3         =    7.4; // height of cell 3
190 //   fV0CHeight4         =   12.5; // height of cell 4
191 //   fV0CRMin            =    4.6; // inner radius of box
192 //   fV0CRBox            =   38.0; // outer radius of box
193 //   fV0CLidThickness    =   0.30; // thickness of Carbon lid
194 //   fV0CCellThickness   =   2.00; // thickness of elementary cell
195 //   fV0CBoxThickness    =   4.70; // thickness of V0C Box
196 //   fV0COffsetFibers    =    1.0; // offset to output fibers
197 //   // V0C Parameters related to light output
198 //   fV0CLightYield         =  93.75; // Light yield in BC408 (93.75 eV per photon)
199 //   fV0CLightAttenuation   =   0.05; // Light attenuation in fiber (0.05 per meter)
200 //   fV0CnMeters            =   15.0; // Number of meters of clear fibers to PM
201 //   fV0CFibToPhot          =    0.3; // Attenuation at fiber-photocathode interface
202 // 
203 //   // V0A Parameters related to geometry: All in cm
204 //   fV0AR0     =  4.2;  // Radius of hole
205 //   fV0AR1     =  7.6;  // Maximun radius of 1st cell
206 //   fV0AR2     = 13.8; // Maximun radius of 2nd cell
207 //   fV0AR3     = 22.7; // Maximun radius of 3rd cell
208 //   fV0AR4     = 41.3; // Maximun radius of 4th cell
209 //   fV0AR5     = 43.3; // Radius circunscrite to innermost octagon
210 //   fV0AR6     = 68.0; // Radius circunscrite to outtermost octagon
211 //   fV0ASciWd  =  2.5;  // Scintillator thickness 
212 //   fV0APlaWd  =  0.5;  // Plates thinckness
213 //   fV0APlaAl  = 0.06; // Plates AlMg3 thinckness
214 //   fV0AOctWd  = 0.75; // Innermost octagon thickness
215 //   fV0AOctH1  =  1.0;  // Height of innermost octagon
216 //   fV0AOctH2  =  2.0;  // Height of outtermost octagon
217 //   fV0AFibRd  =  0.1;  // Radius of Fiber
218 //   fV0AFraWd  =  0.2;  // Support Frame thickness
219 //   fV0APMBWd  = 24.6;  // Width of PM Box
220 //   fV0APMBHt  = 22.0;  // Height of PM Box
221 //   fV0APMBTh  =  7.1;  // Thickness of PM Box
222 //   fV0APMBWdW =  0.3;  // Thickness of PM Box Side1 Wall
223 //   fV0APMBHtW =  1.0;  // Thickness of PM Box Side2 Wall
224 //   fV0APMBThW =  0.3;  // Thickness of PM Box Top Wall
225 //   fV0APMBAng = 30.0;  // Angle between PM Box and Support
226 //   fV0APMTR1  = 2.44;  // PMT Glass
227 //   fV0APMTR2  = 2.54;  // PMT Glass
228 //   fV0APMTR3  = 2.54;  // PMT Cover
229 //   fV0APMTR4  = 2.70;  // PMT Cover
230 //   fV0APMTH   = 10.0;  // PMT Height
231 //   fV0APMTB   =  1.0;  // PMT Basis
232 //   fV0APlaEx  =  4.4;  // Plates Extension height
233 //   fV0ABasHt  =  2.0;  // Basis Height
234 //   // V0A Parameters related to light output
235 //   fV0ALightYield         =  93.75;      // Light yield in BC404
236 //   fV0ALightAttenuation   =   0.05;      // Light attenuation in WLS fiber, per meter
237 //   fV0AnMeters            = fV0AR6*0.01; // Tentative value, in meters
238 //   fV0AFibToPhot          =    0.3;      // Attenuation at fiber-photocathode interface
239 }
240 //_____________________________________________________________________________
241
242 void AliVZEROv7::BuildGeometry()
243
244 }
245             
246 //_____________________________________________________________________________
247 void AliVZEROv7::CreateGeometry()
248 {
249 // Constructs TGeo geometry 
250
251   AliDebug(2,"VZERO ConstructGeometry");
252   TGeoVolume *top = gGeoManager->GetVolume("ALIC");
253
254   ///////////////////////////////////////////////////////////////////////////
255   // Construct the geometry of V0C Detector. Brigitte CHEYNIS
256   
257     const int kColorVZERO  = kGreen;
258     TGeoMedium *medV0CAlu = gGeoManager->GetMedium("VZERO_V0CAlu");
259     TGeoMedium *medV0CCar = gGeoManager->GetMedium("VZERO_V0CCar");
260     TGeoMedium *medV0CSci = gGeoManager->GetMedium("VZERO_V0CSci");
261     TGeoVolume *v0RI = new TGeoVolumeAssembly("V0RI");
262     Float_t heightRight, r4Right;
263     Float_t zdet = 90.0 - 0.5 - fV0CBoxThickness/2.0;
264     heightRight  = fV0CHeight1 + fV0CHeight2 + fV0CHeight3 + fV0CHeight4;
265     r4Right      = fV0CRMin + heightRight + 3.0*0.2; // 3 spacings of 2mm between rings
266
267     // Creation of  carbon lids (3.0 mm thick) to keep V0C box shut :
268     Float_t   partube[3];
269     partube[0] =   fV0CRMin;
270     partube[1] =   fV0CRBox;
271     partube[2] =   fV0CLidThickness/2.0;
272     TGeoTube   *sV0CA = new TGeoTube("V0CA", partube[0], partube[1], partube[2]);
273     TGeoVolume *v0CA  = new TGeoVolume("V0CA",sV0CA,medV0CCar);
274     TGeoTranslation *tr2 = new TGeoTranslation(0.,0., fV0CBoxThickness/2.0-partube[2]);
275     TGeoTranslation *tr3 = new TGeoTranslation(0.,0.,-fV0CBoxThickness/2.0+partube[2]);
276     v0RI->AddNode(v0CA,1,tr2);
277     v0RI->AddNode(v0CA,2,tr3);
278     v0CA->SetLineColor(kYellow);
279
280     // Creation of aluminum rings 3.0 mm thick to maintain the v0RI pieces : 
281     partube[0] =   fV0CRMin - 0.3;
282     partube[1] =   fV0CRMin;
283     partube[2] =   fV0CBoxThickness/2.0;
284     TGeoTube   *sV0IR = new TGeoTube("V0IR", partube[0], partube[1], partube[2]);
285     TGeoVolume *v0IR  = new TGeoVolume("V0IR",sV0IR,medV0CAlu);
286     v0RI->AddNode(v0IR,1,0);
287     v0IR->SetLineColor(kYellow);
288     partube[0] =   fV0CRBox;
289     partube[1] =   fV0CRBox + 0.3; 
290     partube[2] =   fV0CBoxThickness/2.0;
291     TGeoTube   *sV0ER = new TGeoTube("V0ER", partube[0], partube[1], partube[2]);
292     TGeoVolume *v0ER  = new TGeoVolume("V0ER",sV0ER,medV0CAlu);
293     v0RI->AddNode(v0ER,1,0);
294     v0ER->SetLineColor(kYellow);
295
296     // Creation of assembly V0R0 of scintillator cells within one sector
297     TGeoVolume *v0R0 = new TGeoVolumeAssembly("V0R0");                                            
298
299     // Elementary cell of ring 1  - right part - :
300     // (cells of ring 1 will be shifted by 2.0 cm backwards to output fibers)
301     Float_t   r1Right =  fV0CRMin + fV0CHeight1;
302     Float_t   offset  = fV0CBoxThickness/2.0 - fV0CLidThickness - fV0CCellThickness/2.0;   
303     Float_t   partubs[5];   
304     partubs[0]     =  fV0CRMin;
305     partubs[1]     =  r1Right;
306     partubs[2]     =  fV0CCellThickness/2.0;
307     partubs[3]     =  90.0-22.5;
308     partubs[4]     = 135.0-22.5;
309     TGeoTubeSeg *sV0R1 = new TGeoTubeSeg("V0R1", partubs[0], partubs[1], partubs[2],
310                                          partubs[3], partubs[4]);
311     TGeoVolume  *v0R1  = new TGeoVolume("V0R1",sV0R1,medV0CSci);                                       
312     TGeoTranslation *tr4 = new TGeoTranslation(0.,0.,-offset);
313     v0R0->AddNode(v0R1,1,tr4);
314     v0R1->SetLineColor(kColorVZERO);
315
316     // Elementary cell of ring 2 - right part - :
317     // (cells of ring 2 will be shifted by 1.0 cm backwards to output fibers)
318     Float_t   r2Right  =  r1Right + fV0CHeight2;  
319     partubs[0]     =  r1Right;  //  must be equal to 7.1
320     partubs[1]     =  r2Right;  //  must be equal to 11.5
321     TGeoTubeSeg *sV0R2 = new TGeoTubeSeg("V0R2", partubs[0], partubs[1], partubs[2],
322                                          partubs[3], partubs[4]);
323     TGeoVolume  *v0R2  = new TGeoVolume("V0R2",sV0R2,medV0CSci);
324     TGeoTranslation *tr5 = new TGeoTranslation(0.0,0.2,-offset + fV0COffsetFibers);
325     v0R0->AddNode(v0R2,1,tr5);
326     v0R2->SetLineColor(kColorVZERO);
327
328     // Ring 3 - right part -  :
329     r2Right  =  r2Right + 0.2;
330     Float_t   r3Right  =  r2Right + fV0CHeight3;     
331     partubs[0]     =  r2Right;  //  must be equal to 11.7
332     partubs[1]     =  r3Right;  //  must be equal to 19.1
333     partubs[3]     =  90.0-22.5;
334     partubs[4]     = 112.5-22.5;
335     TGeoTubeSeg *sV0R3 = new TGeoTubeSeg("V0R3", partubs[0], partubs[1], partubs[2],
336                                          partubs[3], partubs[4]);
337     TGeoVolume  *v0R3  = new TGeoVolume("V0R3",sV0R3,medV0CSci);
338     TGeoTranslation *tr6 = new TGeoTranslation(0.,0.2,-offset + 2.0*fV0COffsetFibers);
339     v0R0->AddNode(v0R3,1,tr6);
340     v0R3->SetLineColor(kColorVZERO);
341     partubs[3]     = 112.5-22.5;
342     partubs[4]     = 135.0-22.5;
343     TGeoTubeSeg *sV0R4 = new TGeoTubeSeg("V0R4", partubs[0], partubs[1], partubs[2],
344                                          partubs[3], partubs[4]);
345     TGeoVolume  *v0R4  = new TGeoVolume("V0R4",sV0R4,medV0CSci);
346     v0R0->AddNode(v0R4,1,tr6);
347     v0R4->SetLineColor(kColorVZERO);
348   
349     // Ring 4 - right part -  : 
350     Float_t x = TMath::ATan(3.5/257.5) * ((180./TMath::Pi()));
351     r3Right = r3Right + 0.2 + 0.2;   // + 0.2 because no shift in translation here !!
352     partubs[0]     =  r3Right;  //  must be equal to 19.5
353     partubs[1]     =  r4Right;  //  must be equal to 32.0
354     partubs[3]     =  90.0-22.5+x;
355     partubs[4]     = 112.5-22.5-x;
356     TGeoTubeSeg *sV0R5 = new TGeoTubeSeg("V0R5", partubs[0], partubs[1], partubs[2],
357                                          partubs[3], partubs[4]);
358     TGeoVolume  *v0R5  = new TGeoVolume("V0R5",sV0R5,medV0CSci);
359     TGeoTranslation *tr7 = new TGeoTranslation(0.,0.0,-offset + 2.0*fV0COffsetFibers);                                        
360     v0R0->AddNode(v0R5,1,tr7);
361     v0R5->SetLineColor(kColorVZERO);
362     partubs[3]     = 112.5-22.5+x;
363     partubs[4]     = 135.0-22.5-x;
364     TGeoTubeSeg *sV0R6 = new TGeoTubeSeg("V0R6", partubs[0], partubs[1], partubs[2],
365                                          partubs[3], partubs[4]);
366     TGeoVolume  *v0R6  = new TGeoVolume("V0R6",sV0R6,medV0CSci);
367     v0R0->AddNode(v0R6,1,tr7);
368     v0R6->SetLineColor(kColorVZERO);
369     Float_t  phi;
370     Float_t  phiDeg= 180./4.;
371     Int_t    nsecR = 1;     // number of sectors in right part of V0
372     for (phi = 22.5; phi < 360.0; phi = phi + phiDeg) {
373       TGeoRotation  *rot1 = new TGeoRotation("rot1", 90.0, +phi, 90., 90.+phi, 0.0, 0.0 ); 
374       v0RI->AddNode(v0R0,nsecR,rot1);    
375       nsecR++;        
376     }
377
378   ///////////////////////////////////////////////////////////////////////////
379   // Construct the geometry of V0A Detector. Carlos PEREZ, PUCP
380   // Revision by Lizardo VALENCIA in July 2008
381
382     const int kV0AColorSci   = 5;
383     const int kV0AColorPlaIn = 3;
384     const int kV0AColorPlaOu = 41;
385     const int kV0AColorOct   = 7;
386     const int kV0AColorFra   = 6;
387     const int kV0AColorFib   = 11;
388     const int kV0AColorPMG   = 1;
389     const int kV0AColorPMA   = 2;
390     TGeoMedium *medV0ASci = gGeoManager->GetMedium("VZERO_V0ASci");
391     TGeoMedium *medV0APlaIn = gGeoManager->GetMedium("VZERO_V0APlaIn");
392     TGeoMedium *medV0APlaOu = gGeoManager->GetMedium("VZERO_V0APlaOu");
393     TGeoMedium *medV0ASup = gGeoManager->GetMedium("VZERO_V0ALuc");
394     TGeoMedium *medV0AFra = gGeoManager->GetMedium("VZERO_V0ALuc");
395     TGeoMedium *medV0AFib = gGeoManager->GetMedium("VZERO_V0AFib");
396     TGeoMedium *medV0APMGlass = gGeoManager->GetMedium("VZERO_V0APMG");
397     TGeoMedium *medV0APMAlum = gGeoManager->GetMedium("VZERO_V0APMA");
398     double pi = TMath::Pi();
399     double sin225   = TMath::Sin(pi/8.);
400     double cos225   = TMath::Cos(pi/8.);
401     double sin45    = TMath::Sin(pi/4.); // lucky: Sin45=Cos45
402     double cos45    = TMath::Cos(pi/4.); 
403     double v0APts[16];
404     //double sin6645   = TMath::Sin(1.16);
405     //double cos6645   = TMath::Cos(1.16);
406     double sin654   = TMath::Sin(1.14);
407     double cos654   = TMath::Cos(1.14);
408     double sin65   = TMath::Sin(1.13);//65
409     double cos65   = TMath::Cos(1.13);
410     double sin665   = TMath::Sin(1.16);
411     double cos665   = TMath::Cos(1.16);//66.5
412
413     ////////////////////////////
414     /// Definition sector 1
415     TGeoVolume *v0ASec = new TGeoVolumeAssembly("V0ASec");
416     //TGeoCompositeShape *sV0AHole = new TGeoCompositeShape("sV0AHole");
417     
418     /// For boolean sustraction
419     double preShape = 0.2;
420     for (int i=0;i<2;i++) {
421       v0APts[0+8*i] = fV0AR0-fV0AFraWd/2.-preShape;  v0APts[1+8*i] = -preShape;
422       v0APts[2+8*i] = fV0AR0-fV0AFraWd/2.-preShape;  v0APts[3+8*i] = fV0AFraWd/2.;
423       v0APts[4+8*i] = fV0AR4+fV0AFraWd/2.+preShape;  v0APts[5+8*i] = fV0AFraWd/2.;
424       v0APts[6+8*i] = fV0AR4+fV0AFraWd/2.+preShape;  v0APts[7+8*i] = -preShape;
425     }
426     new TGeoArb8("sV0ACha1",fV0ASciWd/1.5,v0APts);
427     for (int i=0;i<2;i++) {
428       v0APts[0+8*i] = fV0AR0*sin45-preShape;
429       v0APts[1+8*i] = (fV0AR0-fV0AFraWd)*sin45-preShape;
430       v0APts[2+8*i] = (fV0AR0-fV0AFraWd/2.)*sin45-preShape;
431       v0APts[3+8*i] = (fV0AR0-fV0AFraWd/2.)*sin45;
432       v0APts[4+8*i] = (fV0AR4+fV0AFraWd/2.)*sin45+preShape;
433       v0APts[5+8*i] = (fV0AR4+fV0AFraWd/2.)*sin45+2.*preShape;
434       v0APts[6+8*i] = (fV0AR4+fV0AFraWd)*sin45+preShape;
435       v0APts[7+8*i] = fV0AR4*sin45+preShape;
436     }
437     new TGeoArb8("sV0ACha2", fV0ASciWd/2.+2.*preShape, v0APts);
438     new TGeoCompositeShape("sV0ACha12","sV0ACha1+sV0ACha2");
439     new TGeoTube("sV0ANail1SciHole", 0.0, 0.4, 1.65);
440     TGeoTranslation *pos1 = new TGeoTranslation("pos1", 42.9, 0.51, 0.0);
441     pos1->RegisterYourself();
442     new TGeoTube("sV0ANail2SciHole", 0.0, 0.4, 1.65);
443     TGeoTranslation *pos2 = new TGeoTranslation("pos2", 30.8,30.04,0.0);
444     pos2->RegisterYourself();
445     new TGeoCompositeShape("sV0ANailsSciHoles","sV0ANail1SciHole:pos1+sV0ANail2SciHole:pos2");
446     new TGeoCompositeShape("sV0ACha","sV0ACha12+sV0ANailsSciHoles");
447     
448     
449     /// Frame
450     TGeoVolume *v0AFra = new TGeoVolumeAssembly("V0AFra");
451     for (int i=0;i<2;i++) {
452       v0APts[0+8*i] = fV0AR0-fV0AFraWd/2.;  v0APts[1+8*i] = 0.;
453       v0APts[2+8*i] = fV0AR0-fV0AFraWd/2.;  v0APts[3+8*i] = fV0AFraWd/2.;
454       v0APts[4+8*i] = fV0AR4+fV0AFraWd/2.;  v0APts[5+8*i] = fV0AFraWd/2.;
455       v0APts[6+8*i] = fV0AR4+fV0AFraWd/2.;  v0APts[7+8*i] = 0.;
456     }
457     TGeoArb8 *sV0AFraB1 = new TGeoArb8("sV0AFraB1",fV0ASciWd/2.,v0APts);
458     TGeoVolume *v0AFraB1 = new TGeoVolume("V0AFraB1",sV0AFraB1,medV0AFra);
459     for (int i=0;i<2;i++) {
460       v0APts[0+8*i] = fV0AR0*sin45;
461       v0APts[1+8*i] = (fV0AR0-fV0AFraWd)*sin45;
462       v0APts[2+8*i] = (fV0AR0-fV0AFraWd/2.)*sin45;
463       v0APts[3+8*i] = (fV0AR0-fV0AFraWd/2.)*sin45;
464       v0APts[4+8*i] = (fV0AR4+fV0AFraWd/2.)*sin45;
465       v0APts[5+8*i] = (fV0AR4+fV0AFraWd/2.)*sin45;
466       v0APts[6+8*i] = (fV0AR4+fV0AFraWd)*sin45;
467       v0APts[7+8*i] = fV0AR4*sin45;
468     }
469     TGeoArb8 *sV0AFraB2 = new TGeoArb8("sV0AFraB2", fV0ASciWd/2., v0APts);
470     TGeoVolume *v0AFraB2 = new TGeoVolume("V0AFraB2",sV0AFraB2,medV0AFra);
471     v0AFraB1->SetLineColor(kV0AColorFra); v0AFraB2->SetLineColor(kV0AColorFra);
472     v0AFra->AddNode(v0AFraB1,1);
473     v0AFra->AddNode(v0AFraB2,1);  // Prefer 2 GeoObjects insted of 3 GeoMovements
474     new TGeoTubeSeg( "sV0AFraR1b", fV0AR0-fV0AFraWd/2.,
475                      fV0AR0+fV0AFraWd/2., fV0ASciWd/2., 0, 45);
476     new TGeoTubeSeg( "sV0AFraR2b", fV0AR1-fV0AFraWd/2.,
477                      fV0AR1+fV0AFraWd/2., fV0ASciWd/2., 0, 45);
478     new TGeoTubeSeg( "sV0AFraR3b", fV0AR2-fV0AFraWd/2.,
479                      fV0AR2+fV0AFraWd/2., fV0ASciWd/2., 0, 45);
480     new TGeoTubeSeg( "sV0AFraR4b", fV0AR3-fV0AFraWd/2.,
481                      fV0AR3+fV0AFraWd/2., fV0ASciWd/2., 0, 45);
482     new TGeoTubeSeg( "sV0AFraR5b", fV0AR4-fV0AFraWd/2.,
483                      fV0AR4+fV0AFraWd/2., fV0ASciWd/2., 0, 45);
484     TGeoCompositeShape *sV0AFraR1 = new TGeoCompositeShape("sV0AFraR1","sV0AFraR1b-sV0ACha");
485     TGeoCompositeShape *sV0AFraR2 = new TGeoCompositeShape("sV0AFraR2","sV0AFraR2b-sV0ACha");
486     TGeoCompositeShape *sV0AFraR3 = new TGeoCompositeShape("sV0AFraR3","sV0AFraR3b-sV0ACha");
487     TGeoCompositeShape *sV0AFraR4 = new TGeoCompositeShape("sV0AFraR4","sV0AFraR4b-sV0ACha");
488     TGeoCompositeShape *sV0AFraR5 = new TGeoCompositeShape("sV0AFraR5","sV0AFraR5b-sV0ACha");
489     TGeoVolume *v0AFraR1 = new TGeoVolume("V0AFraR1",sV0AFraR1,medV0AFra);
490     TGeoVolume *v0AFraR2 = new TGeoVolume("V0AFraR2",sV0AFraR2,medV0AFra);
491     TGeoVolume *v0AFraR3 = new TGeoVolume("V0AFraR3",sV0AFraR3,medV0AFra);
492     TGeoVolume *v0AFraR4 = new TGeoVolume("V0AFraR4",sV0AFraR4,medV0AFra);
493     TGeoVolume *v0AFraR5 = new TGeoVolume("V0AFraR5",sV0AFraR5,medV0AFra);
494     v0AFraR1->SetLineColor(kV0AColorFra); v0AFraR2->SetLineColor(kV0AColorFra);
495     v0AFraR3->SetLineColor(kV0AColorFra); v0AFraR4->SetLineColor(kV0AColorFra);
496     v0AFraR5->SetLineColor(kV0AColorFra);
497     v0AFra->AddNode(v0AFraR1,1);
498     v0AFra->AddNode(v0AFraR2,1);
499     v0AFra->AddNode(v0AFraR3,1); 
500     v0AFra->AddNode(v0AFraR4,1);
501     v0AFra->AddNode(v0AFraR5,1);
502     v0ASec->AddNode(v0AFra,1);
503
504     /// Sensitive scintilator
505     TGeoVolume *v0ASci = new TGeoVolumeAssembly("V0ASci");
506     new TGeoTubeSeg( "sV0AR1b", fV0AR0+fV0AFraWd/2.,
507                      fV0AR1-fV0AFraWd/2., fV0ASciWd/2., 0, 45);
508     new TGeoTubeSeg( "sV0AR2b", fV0AR1+fV0AFraWd/2.,
509                      fV0AR2-fV0AFraWd/2., fV0ASciWd/2., 0, 45);
510     new TGeoTubeSeg( "sV0AR3b", fV0AR2+fV0AFraWd/2.,
511                      fV0AR3-fV0AFraWd/2., fV0ASciWd/2., 0, 45);
512     new TGeoTubeSeg( "sV0AR4b", fV0AR3+fV0AFraWd/2.,
513                      fV0AR4-fV0AFraWd/2., fV0ASciWd/2., 0, 45);
514     TGeoCompositeShape *sV0AR1 = new TGeoCompositeShape("sV0AR1","sV0AR1b-sV0ACha");
515     TGeoCompositeShape *sV0AR2 = new TGeoCompositeShape("sV0AR2","sV0AR2b-sV0ACha");
516     TGeoCompositeShape *sV0AR3 = new TGeoCompositeShape("sV0AR3","sV0AR3b-sV0ACha");
517     TGeoCompositeShape *sV0AR4 = new TGeoCompositeShape("sV0AR4","sV0AR4b-sV0ACha");
518     TGeoVolume *v0L1 = new TGeoVolume("V0L1",sV0AR1,medV0ASci);
519     TGeoVolume *v0L2 = new TGeoVolume("V0L2",sV0AR2,medV0ASci);
520     TGeoVolume *v0L3 = new TGeoVolume("V0L3",sV0AR3,medV0ASci);
521     TGeoVolume *v0L4 = new TGeoVolume("V0L4",sV0AR4,medV0ASci);
522     v0L1->SetLineColor(kV0AColorSci); v0L2->SetLineColor(kV0AColorSci);
523     v0L3->SetLineColor(kV0AColorSci); v0L4->SetLineColor(kV0AColorSci);
524     v0ASec->AddNode(v0L1,1);
525     v0ASec->AddNode(v0L2,1);
526     v0ASec->AddNode(v0L3,1);
527     v0ASec->AddNode(v0L4,1);
528     
529
530     /// Non-sensitive scintilator
531     new TGeoTubeSeg("sV0AR5S2", fV0AR4+fV0AFraWd/2., fV0AR4 + fV0AR0, fV0ASciWd/2.+2*preShape, 0, 45);
532     TGeoCompositeShape *sV0AR5 = new TGeoCompositeShape("V0AR5","sV0AR5S2 - sV0ACha");
533     TGeoVolume *v0AR5 = new TGeoVolume("V0AR5",sV0AR5,medV0ASci);
534     v0AR5->SetLineColor(kV0AColorSci);
535     v0ASci->AddNode(v0AR5,1);
536     v0ASec->AddNode(v0ASci,1); 
537
538     /// Segment of innermost octagon
539     TGeoVolume *v0ASup = new TGeoVolumeAssembly("V0ASup");
540     
541     /// Segment of outtermost octagon
542     for (int i=0;i<2;i++) {
543     v0APts[0+8*i] =  fV0AR6-fV0AOctH2;  v0APts[1+8*i] = 0.;
544       v0APts[2+8*i] = (fV0AR6-fV0AOctH2)*sin45;  v0APts[3+8*i] = (fV0AR6-fV0AOctH2)*sin45;
545       v0APts[4+8*i] = fV0AR6*sin45;             v0APts[5+8*i] = fV0AR6*sin45;
546       v0APts[6+8*i] = fV0AR6;                   v0APts[7+8*i] = 0.;
547     }
548     TGeoArb8 *sV0AOct2 = new TGeoArb8("sV0AOct2", (fV0ASciWd+2*fV0AOctWd)/2., v0APts);
549     TGeoVolume *v0AOct2 = new TGeoVolume("V0AOct2", sV0AOct2,medV0ASup);
550     v0AOct2->SetLineColor(kV0AColorOct);
551     v0ASup->AddNode(v0AOct2,1);
552     v0ASec->AddNode(v0ASup,1);
553
554
555     //Bunch of fibers
556     v0APts[ 0] = v0APts[ 2] = -14.0;
557     v0APts[ 1] = v0APts[ 7] = (fV0ASciWd+fV0AOctWd)/2.-0.01;
558     v0APts[ 3] = v0APts[ 5] = (fV0ASciWd+fV0AOctWd)/2.+0.01;
559     v0APts[ 4] = v0APts[ 6] = +14.0;
560     v0APts[ 8] = v0APts[10] = -10.0;
561     v0APts[ 9] = v0APts[15] = 0.;
562     v0APts[11] = v0APts[13] = 0.25;
563     v0APts[12] = v0APts[14] = +10.0;
564     TGeoArb8 *sV0AFib = new TGeoArb8("sV0AFib", 9.0, v0APts);
565     TGeoVolume *v0AFib1 = new TGeoVolume("V0AFib1",sV0AFib,medV0AFib);
566     TGeoVolume *v0AFib = new TGeoVolumeAssembly("V0AFib");
567     TGeoRotation *rot = new TGeoRotation("rot");
568     rot->RotateX(-90);
569     rot->RotateZ(-90.+22.5);
570     v0AFib->AddNode(v0AFib1,1,rot);
571     rot = new TGeoRotation("rot");
572     rot->RotateX(-90);
573     rot->RotateY(180);
574     rot->RotateZ(-90.+22.5);
575     v0AFib->SetLineColor(kV0AColorFib);
576     v0AFib->AddNode(v0AFib1,2,rot);
577     v0ASec->AddNode(v0AFib,1,new TGeoTranslation((fV0AR6-fV0AOctH2+fV0AR5)*cos225/2. - 1.0, (fV0AR6-fV0AOctH2+fV0AR5)*sin225/2., 0));
578     
579   
580     /// Plates
581     new TGeoTube("sV0ANail1PlaInHole", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
582     new TGeoTube("sV0ANail2PlaInHole", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
583     new TGeoCompositeShape("sV0ANailsPlaInHoles","sV0ANail1PlaInHole:pos1+sV0ANail2PlaInHole:pos2");
584     new TGeoTube("sV0ANail1PlaOuHole", 0.0, 0.4, (fV0APlaAl)/2.);
585     new TGeoTube("sV0ANail2PlaOuHole", 0.0, 0.4, (fV0APlaAl)/2.);
586     new TGeoCompositeShape("sV0ANailsPlaOuHoles","sV0ANail1PlaOuHole:pos1+sV0ANail2PlaOuHole:pos2");
587     for (int i=0;i<2;i++) {
588       v0APts[0+8*i] = fV0AR0;                   v0APts[1+8*i] = 0.;
589       v0APts[2+8*i] = fV0AR0*sin45;             v0APts[3+8*i] = fV0AR0*sin45;
590       v0APts[4+8*i] = fV0AR6 * sin45;   v0APts[5+8*i] = fV0AR6*sin45;
591       v0APts[6+8*i] = fV0AR6;           v0APts[7+8*i] = 0.;
592     }
593     new TGeoArb8("sV0APlaIn", (fV0APlaWd-2*fV0APlaAl)/2., v0APts);
594     TGeoCompositeShape *sV0APlaInNailsHoles = new TGeoCompositeShape("sV0APlaInNailsHoles","sV0APlaIn-sV0ANailsPlaInHoles");
595     TGeoVolume *v0APlaInNailsHoles = new TGeoVolume("V0APlaInNailsHoles", sV0APlaInNailsHoles, medV0APlaIn);
596     new TGeoArb8("sV0APlaOu", fV0APlaAl/2., v0APts);
597     TGeoCompositeShape *sV0APlaOuNailsHoles = new TGeoCompositeShape("sV0APlaOuNailsHoles","sV0APlaOu-sV0ANailsPlaOuHoles"); 
598     TGeoVolume *v0APlaOuNailsHoles = new TGeoVolume("V0APlaOuNailsHoles", sV0APlaOuNailsHoles, medV0APlaOu);
599     v0APlaInNailsHoles->SetLineColor(kV0AColorPlaIn); v0APlaOuNailsHoles->SetLineColor(kV0AColorPlaOu);
600     TGeoVolume *v0APla = new TGeoVolumeAssembly("V0APla");
601     v0APla->AddNode(v0APlaInNailsHoles,1);
602     v0APla->AddNode(v0APlaOuNailsHoles,1,new TGeoTranslation(0,0,(fV0APlaWd-fV0APlaAl)/2.));
603     v0APla->AddNode(v0APlaOuNailsHoles,2,new TGeoTranslation(0,0,-(fV0APlaWd-fV0APlaAl)/2.));
604     v0ASec->AddNode(v0APla,1,new TGeoTranslation(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
605     v0ASec->AddNode(v0APla,2,new TGeoTranslation(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
606
607     /// PMBox
608     TGeoVolume* v0APM = new TGeoVolumeAssembly("V0APM");
609     new TGeoBBox("sV0APMB1", fV0APMBWd/2., fV0APMBHt/2., fV0APMBTh/2.);
610     new TGeoBBox("sV0APMB2", fV0APMBWd/2.-fV0APMBWdW, fV0APMBHt/2.-fV0APMBHtW, fV0APMBTh/2.-fV0APMBThW);
611     TGeoCompositeShape *sV0APMB = new TGeoCompositeShape("sV0APMB","sV0APMB1-sV0APMB2");
612     TGeoVolume *v0APMB = new TGeoVolume("V0APMB",sV0APMB, medV0APMAlum);
613     v0APMB->SetLineColor(kV0AColorPMA);
614     v0APM->AddNode(v0APMB,1);
615
616     /// PMTubes
617     TGeoTube *sV0APMT1 = new TGeoTube("sV0APMT1", fV0APMTR1, fV0APMTR2, fV0APMTH/2.);
618     TGeoVolume *v0APMT1 = new TGeoVolume("V0APMT1", sV0APMT1, medV0APMGlass);
619     TGeoTube *sV0APMT2 = new TGeoTube("sV0APMT2", fV0APMTR3, fV0APMTR4, fV0APMTH/2.);
620     TGeoVolume *v0APMT2 = new TGeoVolume("V0APMT2", sV0APMT2, medV0APMAlum);
621     TGeoVolume *v0APMT = new TGeoVolumeAssembly("V0APMT");
622     TGeoTube *sV0APMTT = new TGeoTube("sV0APMTT", 0., fV0APMTR4, fV0APMTB/2.);
623     TGeoVolume *v0APMTT = new TGeoVolume("V0APMTT", sV0APMTT, medV0APMAlum);
624     v0APMT1->SetLineColor(kV0AColorPMG);
625     v0APMT2->SetLineColor(kV0AColorPMA);
626     v0APMTT->SetLineColor(kV0AColorPMA);
627     rot = new TGeoRotation("rot", 90, 0, 180, 0, 90, 90);
628     v0APMT->AddNode(v0APMT1,1,rot);
629     v0APMT->AddNode(v0APMT2,1,rot);
630     v0APMT->AddNode(v0APMTT,1,new TGeoCombiTrans(0,-(fV0APMTH+fV0APMTB)/2.,0,rot));
631     double autoShift = (fV0APMBWd-2*fV0APMBWdW)/4.;
632     v0APM->AddNode(v0APMT, 1, new TGeoTranslation(-1.5*autoShift, 0, 0));
633     v0APM->AddNode(v0APMT, 2, new TGeoTranslation(-0.5*autoShift, 0, 0));
634     v0APM->AddNode(v0APMT, 3, new TGeoTranslation(+0.5*autoShift, 0, 0));
635     v0APM->AddNode(v0APMT, 4, new TGeoTranslation(+1.5*autoShift, 0, 0));
636
637     /// PM
638     rot = new TGeoRotation("rot");
639     rot->RotateX(90-fV0APMBAng);
640     rot->RotateZ(-90.+22.5);
641     double cosAngPMB = TMath::Cos(fV0APMBAng*TMath::DegToRad());
642     double sinAngPMB = TMath::Sin(fV0APMBAng*TMath::DegToRad());
643     double shiftZ = fV0APMBHt/2. * cosAngPMB
644       -   ( fV0ASciWd + 2 * fV0AOctWd + 2 * fV0APlaWd )/2.   -   fV0APMBTh/2. * sinAngPMB;
645     double shiftR = fV0AR6  +  fV0AOctH2 + fV0APlaAl;
646     v0ASec->AddNode(v0APM,1, new TGeoCombiTrans( shiftR*cos225+1.07, shiftR*sin225, shiftZ, rot));
647     
648     // Aluminium nails 
649     TGeoTube *sV0ANail1 = new TGeoTube("sV0ANail1", 0.0, 0.4, 5.09/2.);
650     TGeoVolume *v0ANail1 = new TGeoVolume("V0ANail1", sV0ANail1, medV0APMAlum);
651     v0ANail1->SetLineColor(kV0AColorPMA);// this is the color for aluminium
652     v0ASec->AddNode(v0ANail1,1,new TGeoTranslation(42.9, 0.51, 0.0));
653     TGeoTube *sV0ANail2 = new TGeoTube("sV0ANail2", 0.0, 0.4, 5.09/2.);
654     TGeoVolume *v0ANail2 = new TGeoVolume("V0ANail2", sV0ANail2, medV0APMAlum);
655     v0ANail2->SetLineColor(kV0AColorPMA);
656     v0ASec->AddNode(v0ANail2,1,new TGeoTranslation(30.8,30.04,0.0));
657     
658     
659     /// Replicate sectors
660     TGeoVolume *v0LE = new TGeoVolumeAssembly("V0LE");
661     for(int i=0; i<4; i++) {
662       TGeoRotation *rot = new TGeoRotation("rot", 90., i*45., 90., 90.+i*45., 0., 0.);
663       v0LE->AddNode(v0ASec,i+1,rot);  /// modificacion +1 anhadido
664     }
665    
666     //Upper supports
667     for (int i=0;i<2;i++){
668     v0APts[0+8*i] = 0.2;            v0APts[1+8*i] = 45.5;  
669     v0APts[2+8*i] = 0.2;          v0APts[3+8*i] = 70.5;   //70.6
670     v0APts[4+8*i] = 4.0;            v0APts[5+8*i] = 68.9;
671     v0APts[6+8*i] = 4.0;            v0APts[7+8*i] = 45.5;  
672     }
673     TGeoArb8 *sV0ASuppur = new TGeoArb8("sV0ASuppur", 2.0, v0APts);    
674     TGeoVolume *v0ASuppur = new TGeoVolume("V0ASuppur", sV0ASuppur, medV0ASup);
675     v0ASuppur->SetLineColor(kV0AColorOct);
676     v0LE->AddNode(v0ASuppur,1);
677     for (int i=0;i<2;i++){
678     v0APts[0+8*i] = -0.2;           v0APts[1+8*i] = 45.3;
679     v0APts[2+8*i] = -0.2;            v0APts[3+8*i] = 70.6;
680     v0APts[4+8*i] = -4.0;           v0APts[5+8*i] = 68.9;
681     v0APts[6+8*i] = -4.0;           v0APts[7+8*i] = 45.3;
682     }
683     TGeoArb8 *sV0ASuppul = new TGeoArb8("sV0ASuppul", 2.0, v0APts);    
684     TGeoVolume *v0ASuppul = new TGeoVolume("V0ASuppul", sV0ASuppul, medV0ASup);
685     v0ASuppul->SetLineColor(kV0AColorOct);
686     v0LE->AddNode(v0ASuppul,1);
687   
688
689       //Definition of sector 5
690    
691        TGeoVolume *v0ASec5 = new TGeoVolumeAssembly("V0ASec5"); 
692
693  /// For boolean sustraction
694     double preShape5 = 0.2;
695     for (int i=0;i<2;i++) {
696       v0APts[0+8*i] = -fV0AR0+fV0AFraWd/2.-preShape5;  v0APts[1+8*i] = -preShape5;
697       v0APts[2+8*i] = -fV0AR0+fV0AFraWd/2.-preShape5;  v0APts[3+8*i] = fV0AFraWd/2.;
698       v0APts[4+8*i] = -fV0AR4-fV0AFraWd/2.+preShape5;  v0APts[5+8*i] = fV0AFraWd/2.;
699       v0APts[6+8*i] = -fV0AR4-fV0AFraWd/2.+preShape5;  v0APts[7+8*i] = -preShape5;
700     }
701     new TGeoArb8("sV0ACha15",fV0ASciWd/1.5,v0APts);
702     for (int i=0;i<2;i++) {
703       v0APts[0+8*i] = -fV0AR0*cos65;
704       v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin65;
705       v0APts[2+8*i] = -(fV0AR0-fV0AFraWd/2.)*cos65;
706       v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin65;
707       v0APts[4+8*i] = -(fV0AR4+fV0AFraWd/2.)*cos65;
708       v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin65;
709       v0APts[6+8*i] = -(fV0AR4+fV0AFraWd)*cos65;
710       v0APts[7+8*i] = -fV0AR4*sin65;
711     }
712     new TGeoArb8("sV0ACha25", fV0ASciWd/2.+2.*preShape5, v0APts);
713     new TGeoCompositeShape("sV0ACha125","sV0ACha15+sV0ACha25");
714     new TGeoTube("sV0ANail15Hole", 0.0, 0.4, 1.65);
715     TGeoTranslation *pos15 = new TGeoTranslation("pos15", -42.9, -0.51, 0.0);
716     pos15->RegisterYourself();
717     new TGeoTube("sV0ANail25Hole", 0.0, 0.4, 1.65);
718     TGeoTranslation *pos25 = new TGeoTranslation("pos25",-30.8,-30.04,0.0);
719     pos25->RegisterYourself();
720     new TGeoTube("sV0ANail35Hole", 0.0, 0.4, 1.65);
721     TGeoTranslation *pos35 = new TGeoTranslation("pos35", -30.05,-30.79,0.0);
722     pos35->RegisterYourself();
723     new TGeoCompositeShape("sV0ANailsHoles5","sV0ANail15Hole:pos15+sV0ANail25Hole:pos25+sV0ANail35Hole:pos35");
724     new TGeoCompositeShape("sV0ACha5","sV0ACha125+sV0ANailsHoles5");
725
726     /// Frame
727     TGeoVolume *v0AFra5 = new TGeoVolumeAssembly("V0AFra5");
728     for (int i=0;i<2;i++) {
729       v0APts[0+8*i] = -fV0AR0+fV0AFraWd/2.;  v0APts[1+8*i] = 0.;
730       v0APts[2+8*i] = -fV0AR0+fV0AFraWd/2.;  v0APts[3+8*i] = fV0AFraWd/2.;
731       v0APts[4+8*i] = -fV0AR4-fV0AFraWd/2.;  v0APts[5+8*i] = fV0AFraWd/2.;
732       v0APts[6+8*i] = -fV0AR4-fV0AFraWd/2.;  v0APts[7+8*i] = 0.;
733     }
734     TGeoArb8 *sV0AFraB15 = new TGeoArb8("sV0AFraB15",fV0ASciWd/2.,v0APts);
735     TGeoVolume *v0AFraB15 = new TGeoVolume("V0AFraB15",sV0AFraB15,medV0AFra);
736     for (int i=0;i<2;i++) {
737       v0APts[0+8*i] = -fV0AR0*cos65;
738       v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin65;
739       v0APts[2+8*i] = -(fV0AR0-fV0AFraWd/2.)*cos65;
740       v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin65;
741       v0APts[4+8*i] = -(fV0AR4+fV0AFraWd/2.)*cos65;
742       v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin65;
743       v0APts[6+8*i] = -(fV0AR4+fV0AFraWd)*cos65;
744       v0APts[7+8*i] = -fV0AR4*sin65;
745     }
746     TGeoArb8 *sV0AFraB25 = new TGeoArb8("sV0AFraB25", fV0ASciWd/2., v0APts);
747     TGeoVolume *v0AFraB25 = new TGeoVolume("V0AFraB25",sV0AFraB25,medV0AFra);
748     v0AFraB15->SetLineColor(kV0AColorFra); v0AFraB25->SetLineColor(kV0AColorFra);
749     v0AFra5->AddNode(v0AFraB15,1);
750     v0AFra5->AddNode(v0AFraB25,1);  // Prefer 2 GeoObjects insted of 3 GeoMovements
751     new TGeoTubeSeg( "sV0AFraR1b5", fV0AR0-fV0AFraWd/2.,
752                      fV0AR0+fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
753     new TGeoTubeSeg( "sV0AFraR2b5", fV0AR1-fV0AFraWd/2.,
754                      fV0AR1+fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
755     new TGeoTubeSeg( "sV0AFraR3b5", fV0AR2-fV0AFraWd/2.,
756                      fV0AR2+fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
757     new TGeoTubeSeg( "sV0AFraR4b5", fV0AR3-fV0AFraWd/2.,
758                      fV0AR3+fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
759     new TGeoTubeSeg( "sV0AFraR5b5", fV0AR4-fV0AFraWd/2.,
760                      fV0AR4+fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
761     TGeoCompositeShape *sV0AFraR15 = new TGeoCompositeShape("sV0AFraR15","sV0AFraR1b5-sV0ACha5");
762     TGeoCompositeShape *sV0AFraR25 = new TGeoCompositeShape("sV0AFraR25","sV0AFraR2b5-sV0ACha5");
763     TGeoCompositeShape *sV0AFraR35 = new TGeoCompositeShape("sV0AFraR35","sV0AFraR3b5-sV0ACha5");
764     TGeoCompositeShape *sV0AFraR45 = new TGeoCompositeShape("sV0AFraR45","sV0AFraR4b5-sV0ACha5");
765     TGeoCompositeShape *sV0AFraR55 = new TGeoCompositeShape("sV0AFraR55","sV0AFraR5b5-sV0ACha5");
766     TGeoVolume *v0AFraR15 = new TGeoVolume("V0AFraR15",sV0AFraR15,medV0AFra);
767     TGeoVolume *v0AFraR25 = new TGeoVolume("V0AFraR25",sV0AFraR25,medV0AFra);
768     TGeoVolume *v0AFraR35 = new TGeoVolume("V0AFraR35",sV0AFraR35,medV0AFra);
769     TGeoVolume *v0AFraR45 = new TGeoVolume("V0AFraR45",sV0AFraR45,medV0AFra);
770     TGeoVolume *v0AFraR55 = new TGeoVolume("V0AFraR55",sV0AFraR55,medV0AFra);
771     v0AFraR15->SetLineColor(kV0AColorFra); v0AFraR25->SetLineColor(kV0AColorFra);
772     v0AFraR35->SetLineColor(kV0AColorFra); v0AFraR45->SetLineColor(kV0AColorFra);
773     v0AFraR55->SetLineColor(kV0AColorFra);
774     v0AFra5->AddNode(v0AFraR15,1);
775     v0AFra5->AddNode(v0AFraR25,1);
776     v0AFra5->AddNode(v0AFraR35,1);
777     v0AFra5->AddNode(v0AFraR45,1);
778     v0AFra5->AddNode(v0AFraR55,1);
779     v0ASec5->AddNode(v0AFra5,1);
780
781     /// Sensitive scintilator
782     TGeoVolume *v0ASci5 = new TGeoVolumeAssembly("V0ASci5");
783     new TGeoTubeSeg( "sV0AR1b5", fV0AR0+fV0AFraWd/2.,
784                      fV0AR1-fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
785     new TGeoTubeSeg( "sV0AR2b5", fV0AR1+fV0AFraWd/2.,
786                      fV0AR2-fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
787     new TGeoTubeSeg( "sV0AR3b5", fV0AR2+fV0AFraWd/2.,
788                      fV0AR3-fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
789     new TGeoTubeSeg( "sV0AR4b5", fV0AR3+fV0AFraWd/2.,
790                      fV0AR4-fV0AFraWd/2., fV0ASciWd/2., 180.0, 245.4);
791     TGeoCompositeShape *sV0AR15 = new TGeoCompositeShape("sV0AR15","sV0AR1b5-sV0ACha5");
792     TGeoCompositeShape *sV0AR25 = new TGeoCompositeShape("sV0AR25","sV0AR2b5-sV0ACha5");
793     TGeoCompositeShape *sV0AR35 = new TGeoCompositeShape("sV0AR35","sV0AR3b5-sV0ACha5");
794     TGeoCompositeShape *sV0AR45 = new TGeoCompositeShape("sV0AR45","sV0AR4b5-sV0ACha5");
795     TGeoVolume *v0L15 = new TGeoVolume("V0L15",sV0AR15,medV0ASci);
796     TGeoVolume *v0L25 = new TGeoVolume("V0L25",sV0AR25,medV0ASci);
797     TGeoVolume *v0L35 = new TGeoVolume("V0L35",sV0AR35,medV0ASci);
798     TGeoVolume *v0L45 = new TGeoVolume("V0L45",sV0AR45,medV0ASci);
799     v0L15->SetLineColor(kV0AColorSci); v0L25->SetLineColor(kV0AColorSci);
800     v0L35->SetLineColor(kV0AColorSci); v0L45->SetLineColor(kV0AColorSci);
801     v0ASci5->AddNode(v0L15,1);
802     v0ASci5->AddNode(v0L25,1);
803     v0ASci5->AddNode(v0L35,1);
804     v0ASci5->AddNode(v0L45,1);
805
806     /// Non-sensitive scintilator
807     new TGeoTubeSeg("sV0AR5S25", fV0AR4+fV0AFraWd/2., fV0AR4 + fV0AR0, fV0ASciWd/2.+2*preShape, 180.0, 245.4);
808     TGeoCompositeShape *sV0AR55 = new TGeoCompositeShape("V0AR55","sV0AR5S25 - sV0ACha5");
809     TGeoVolume *v0AR55 = new TGeoVolume("V0AR55",sV0AR55,medV0ASci);
810     v0AR55->SetLineColor(kV0AColorSci);
811     v0ASci5->AddNode(v0AR55,1);
812     v0ASec5->AddNode(v0ASci5,1);
813
814     /// Segment of innermost octagon
815     TGeoVolume *v0ASup5 = new TGeoVolumeAssembly("V0ASup5");
816
817 /// Segment of outtermost octagon   
818      for (int i=0;i<2;i++) {
819     v0APts[0+8*i] = -fV0AR6+fV0AOctH2;           v0APts[1+8*i] = 0.;
820     v0APts[2+8*i] = -(fV0AR7-fV0AOctH2)*cos654;   v0APts[3+8*i] = -(fV0AR7-fV0AOctH2)*sin654;
821     v0APts[4+8*i] = -fV0AR7*cos654;              v0APts[5+8*i] = -fV0AR7*sin654;
822     v0APts[6+8*i] = -fV0AR6;                     v0APts[7+8*i] = 0.;
823     }
824     TGeoArb8 *sV0AOct25 = new TGeoArb8("sV0AOct25", (fV0ASciWd+2*fV0AOctWd)/2., v0APts);
825     TGeoVolume *v0AOct25 = new TGeoVolume("V0AOct25", sV0AOct25,medV0ASup);
826     v0AOct25->SetLineColor(kV0AColorOct);
827     v0ASup5->AddNode(v0AOct25,1);
828     v0ASec5->AddNode(v0ASup5,1);
829
830     //Bunch of fibers
831     v0APts[ 0] = v0APts[ 2] = -14.0;
832     v0APts[ 1] = v0APts[ 7] = (fV0ASciWd+fV0AOctWd)/2.-0.01;
833     v0APts[ 3] = v0APts[ 5] = (fV0ASciWd+fV0AOctWd)/2.+0.01;
834     v0APts[ 4] = v0APts[ 6] = +14.0;
835     v0APts[ 8] = v0APts[10] = -10.0;
836     v0APts[ 9] = v0APts[15] = 0.;
837     v0APts[11] = v0APts[13] = 0.25;
838     v0APts[12] = v0APts[14] = +10.0;
839     TGeoArb8 *sV0AFiba5 = new TGeoArb8("sV0AFiba5", 9.0, v0APts);
840     TGeoVolume *v0AFiba15 = new TGeoVolume("V0AFiba15",sV0AFiba5,medV0AFib);
841     TGeoVolume *v0AFiba5 = new TGeoVolumeAssembly("V0AFiba5");
842     rot = new TGeoRotation("rot");
843     rot->RotateX(-90);
844     rot->RotateZ(90+22.5);
845     v0AFiba5->AddNode(v0AFiba15,1,rot);
846     rot = new TGeoRotation("rot");
847     rot->RotateX(-90);
848     rot->RotateY(180);
849     rot->RotateZ(90+22.5);
850     v0AFiba5->SetLineColor(kV0AColorFib);
851     v0AFiba5->AddNode(v0AFiba15,2,rot);
852     v0ASec5->AddNode(v0AFiba5,1,new TGeoTranslation(-(fV0AR6-fV0AOctH2+fV0AR5)*cos225/2. + 1.5, -(fV0AR6-fV0AOctH2+fV0AR5)*sin225/2. - 1.0, 0));
853     for (int i=0;i<2;i++) {
854     v0APts[0+8*i] = -33.47;           v0APts[1+8*i] = -35.13;
855     v0APts[2+8*i] = -50.48;            v0APts[3+8*i] = -46.97;
856     v0APts[4+8*i] = -44.66;            v0APts[5+8*i] = -61.92;
857     v0APts[6+8*i] = -7.9;             v0APts[7+8*i] = -47.18;
858     }
859     TGeoArb8 *sV0AFibb5 = new TGeoArb8("sV0AFibb8", 1.25, v0APts);
860     TGeoVolume *v0AFibb15 = new TGeoVolume("V0AFibb15",sV0AFibb5,medV0AFib);
861     TGeoVolume *v0AFibb5 = new TGeoVolumeAssembly("V0AFibb5");
862     v0AFibb5->AddNode(v0AFibb15,1);
863     v0AFibb5->SetLineColor(kV0AColorFib);
864     v0ASec5->AddNode(v0AFibb5,1);    
865
866     /// Plates
867     new TGeoTube("sV0ANail1PlaInHole5", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
868     new TGeoTube("sV0ANail2PlaInHole5", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
869     new TGeoTube("sV0ANail3PlaInHole5", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
870     new TGeoCompositeShape("sV0ANailsPlaInHoles5","sV0ANail1PlaInHole5:pos15+sV0ANail2PlaInHole5:pos25+sV0ANail3PlaInHole5:pos35");
871     new TGeoTube("sV0ANail1PlaOuHole5", 0.0, 0.4, (fV0APlaAl)/2.);
872     new TGeoTube("sV0ANail2PlaOuHole5", 0.0, 0.4, (fV0APlaAl)/2.);
873     new TGeoTube("sV0ANail3PlaOuHole5", 0.0, 0.4, (fV0APlaAl)/2.);
874     new TGeoCompositeShape("sV0ANailsPlaOuHoles5","sV0ANail1PlaOuHole5:pos15+sV0ANail2PlaOuHole5:pos25+sV0ANail3PlaOuHole5:pos35");
875     for (int i=0;i<2;i++) {
876       v0APts[0+8*i] = -fV0AR0;                  v0APts[1+8*i] = 0.;
877       v0APts[2+8*i] = -fV0AR0*cos654;           v0APts[3+8*i] = -fV0AR0*sin654;
878       v0APts[4+8*i] = -fV0AR7*cos654;           v0APts[5+8*i] = -fV0AR7*sin654;
879       v0APts[6+8*i] = -fV0AR6;                  v0APts[7+8*i] = 0.;
880     }
881     new TGeoArb8("sV0APlaIn5", (fV0APlaWd-2*fV0APlaAl)/2., v0APts);
882     TGeoCompositeShape *sV0APlaInNailsHoles5 = new TGeoCompositeShape("sV0APlaInNailsHoles5","sV0APlaIn5-sV0ANailsPlaInHoles5");
883     TGeoVolume *v0APlaInNailsHoles5 = new TGeoVolume("V0APlaInNailsHoles5", sV0APlaInNailsHoles5, medV0APlaIn);
884     new TGeoArb8("sV0APlaOu5", fV0APlaAl/2., v0APts);
885     TGeoCompositeShape *sV0APlaOuNailsHoles5 = new TGeoCompositeShape("sV0APlaOuNailsHoles5","sV0APlaOu5-sV0ANailsPlaOuHoles5"); 
886     TGeoVolume *v0APlaOuNailsHoles5 = new TGeoVolume("V0APlaOuNailsHoles5", sV0APlaOuNailsHoles5, medV0APlaOu);
887     v0APlaInNailsHoles5->SetLineColor(kV0AColorPlaIn); v0APlaOuNailsHoles5->SetLineColor(kV0AColorPlaOu);
888     TGeoVolume *v0APla5 = new TGeoVolumeAssembly("V0APla5");
889     v0APla5->AddNode(v0APlaInNailsHoles5,1);
890     v0APla5->AddNode(v0APlaOuNailsHoles5,1,new TGeoTranslation(0,0,(fV0APlaWd-fV0APlaAl)/2.));
891     v0APla5->AddNode(v0APlaOuNailsHoles5,2,new TGeoTranslation(0,0,-(fV0APlaWd-fV0APlaAl)/2.));
892     v0ASec5->AddNode(v0APla5,1,new TGeoTranslation(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
893     v0ASec5->AddNode(v0APla5,2,new TGeoTranslation(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
894     
895     /// PMBox 1ero
896     TGeoVolume* v0APM5 = new TGeoVolumeAssembly("V0APM5");
897     new TGeoBBox("sV0APMB15", fV0APMBWd/2., fV0APMBHt/2., fV0APMBTh/2.);
898     new TGeoBBox("sV0APMB25", fV0APMBWd/2.-fV0APMBWdW, fV0APMBHt/2.-fV0APMBHtW, fV0APMBTh/2.-fV0APMBThW);
899     TGeoCompositeShape *sV0APMB5 = new TGeoCompositeShape("sV0APMB5","sV0APMB15-sV0APMB25");
900     TGeoVolume *v0APMB5 = new TGeoVolume("V0APMB5",sV0APMB5, medV0APMAlum);
901     v0APMB5->SetLineColor(kV0AColorPMA);
902     v0APM5->AddNode(v0APMB5,1);
903
904     /// PMTubes 1ero
905     TGeoTube *sV0APMT15 = new TGeoTube("sV0APMT15", fV0APMTR1, fV0APMTR2, fV0APMTH/2.);
906     TGeoVolume *v0APMT15 = new TGeoVolume("V0APMT15", sV0APMT15, medV0APMGlass);
907     TGeoTube *sV0APMT25 = new TGeoTube("sV0APMT25", fV0APMTR3, fV0APMTR4, fV0APMTH/2.);
908     TGeoVolume *v0APMT25 = new TGeoVolume("V0APMT25", sV0APMT25, medV0APMAlum);
909     TGeoVolume *v0APMT5 = new TGeoVolumeAssembly("V0APMT5");
910     TGeoTube *sV0APMTT5 = new TGeoTube("sV0APMTT5", 0., fV0APMTR4, fV0APMTB/2.);
911     TGeoVolume *v0APMTT5 = new TGeoVolume("V0APMT15", sV0APMTT5, medV0APMAlum);
912     v0APMT5->SetLineColor(kV0AColorPMG);
913     v0APMT25->SetLineColor(kV0AColorPMA);
914     v0APMTT5->SetLineColor(kV0AColorPMA);
915     rot = new TGeoRotation("rot", 90, 0, 180, 0, 90, 90);
916     v0APMT5->AddNode(v0APMT15,1,rot);
917     v0APMT5->AddNode(v0APMT25,1,rot);
918     v0APMT5->AddNode(v0APMTT5,1,new TGeoCombiTrans(0,(fV0APMTH+fV0APMTB)/2.,0,rot));
919     double autoShift5 = (fV0APMBWd-2*fV0APMBWdW)/4.;
920     v0APM5->AddNode(v0APMT5, 1, new TGeoTranslation(-1.5*autoShift5, 0, 0));
921     v0APM5->AddNode(v0APMT5, 2, new TGeoTranslation(-0.5*autoShift5, 0, 0));
922     v0APM5->AddNode(v0APMT5, 3, new TGeoTranslation(+0.5*autoShift5, 0, 0));
923     v0APM5->AddNode(v0APMT5, 4, new TGeoTranslation(+1.5*autoShift5, 0, 0));
924
925     /// PM 1ero
926     rot = new TGeoRotation("rot");
927     rot->RotateX(-90+30);
928     rot->RotateY(0); 
929     rot->RotateZ(-65-3);
930     double cosAngPMB5 = TMath::Cos(fV0APMBAng*TMath::DegToRad());
931     double sinAngPMB5 = TMath::Sin(fV0APMBAng*TMath::DegToRad());
932     double shiftZ5 = fV0APMBHt/2. * cosAngPMB5
933       -   ( fV0ASciWd + 2 * fV0AOctWd + 2 * fV0APlaWd )/2.   -   fV0APMBTh/2. * sinAngPMB5;
934     double shiftR5 = fV0AR6  +  fV0AOctH2 + fV0APlaAl;
935     v0ASec5->AddNode(v0APM5,1, new TGeoCombiTrans( -shiftR5*cos225-1.3, -shiftR5*sin225, shiftZ5, rot));
936
937
938     /// PMBox 2do 
939     TGeoVolume* v0APM52 = new TGeoVolumeAssembly("V0APM52");
940     new TGeoBBox("sV0APMB151", fV0APMBWd/2., fV0APMBHt/2., fV0APMBTh/2.);
941     new TGeoBBox("sV0APMB252", fV0APMBWd/2.-fV0APMBWdW, fV0APMBHt/2.-fV0APMBHtW, fV0APMBTh/2.-fV0APMBThW);
942     TGeoCompositeShape *sV0APMB52 = new TGeoCompositeShape("sV0APMB52","sV0APMB151-sV0APMB252");
943     TGeoVolume *v0APMB52 = new TGeoVolume("V0APMB52",sV0APMB52, medV0APMAlum);
944     v0APMB52->SetLineColor(kV0AColorPMA);
945     v0APM52->AddNode(v0APMB52,1);
946
947     /// PMTubes 2ndo
948     TGeoTube *sV0APMT152 = new TGeoTube("sV0APMT152", fV0APMTR1, fV0APMTR2, fV0APMTH/2.);
949     TGeoVolume *v0APMT152 = new TGeoVolume("V0APMT152", sV0APMT152, medV0APMGlass);
950     TGeoTube *sV0APMT252 = new TGeoTube("sV0APMT252", fV0APMTR3, fV0APMTR4, fV0APMTH/2.);
951     TGeoVolume *v0APMT252 = new TGeoVolume("V0APMT252", sV0APMT252, medV0APMAlum);
952     TGeoVolume *v0APMT52 = new TGeoVolumeAssembly("V0APMT52"); // pk si no se confunde con la 752 o con la 794
953     TGeoTube *sV0APMTT52 = new TGeoTube("sV0APMTT52", 0., fV0APMTR4, fV0APMTB/2.);
954     TGeoVolume *v0APMTT52 = new TGeoVolume("V0APMT52", sV0APMTT52, medV0APMAlum);
955     v0APMT52->SetLineColor(kV0AColorPMG);
956     v0APMT252->SetLineColor(kV0AColorPMA);
957     v0APMTT52->SetLineColor(kV0AColorPMA);
958     rot = new TGeoRotation("rot", 90, 0, 180, 0, 90, 90);
959     v0APMT52->AddNode(v0APMT152,1,rot);
960     v0APMT52->AddNode(v0APMT252,1,rot);
961     v0APMT52->AddNode(v0APMTT52,1,new TGeoCombiTrans(0,(fV0APMTH+fV0APMTB)/2.,0,rot));
962     double autoShift52 = (fV0APMBWd-2*fV0APMBWdW)/4.;
963     v0APM52->AddNode(v0APMT52, 1, new TGeoTranslation(-1.5*autoShift52, 0, 0));
964     v0APM52->AddNode(v0APMT52, 2, new TGeoTranslation(-0.5*autoShift52, 0, 0));
965     v0APM52->AddNode(v0APMT52, 3, new TGeoTranslation(+0.5*autoShift52, 0, 0));
966     v0APM52->AddNode(v0APMT52, 4, new TGeoTranslation(+1.5*autoShift52, 0, 0));
967
968     /// PM 2ndo
969     rot = new TGeoRotation("rot");
970     rot->RotateX(-90+30);
971     rot->RotateY(0);
972     rot->RotateZ(-65-3);
973     double cosAngPMB52 = TMath::Cos(fV0APMBAng*TMath::DegToRad());
974     double sinAngPMB52 = TMath::Sin(fV0APMBAng*TMath::DegToRad());
975     double shiftZ52 = fV0APMBHt/2. * cosAngPMB52
976       -   ( fV0ASciWd + 2 * fV0AOctWd + 2 * fV0APlaWd )/2.   -   fV0APMBTh/2. * sinAngPMB52;
977     double shiftR52 = fV0AR6  + fV0AR1 + fV0AOctWd + fV0APlaAl/3.;
978     v0ASec5->AddNode(v0APM52,1, new TGeoCombiTrans( -shiftR52*cos45-1.3, -shiftR52*sin45, shiftZ52, rot));    
979     
980
981     // Replicate sectors
982     //for(int i=0; i<1; i++) {
983     //TGeoRotation *rot = new TGeoRotation("rot", 90., i*45., 90., 90.+i*45., 0., 0.);
984     v0LE->AddNode(v0ASec5, 1);   //i + 1,rot);
985     //}
986
987
988      //Definition of  sector 6
989    
990     TGeoVolume *v0ASec6 = new TGeoVolumeAssembly("V0ASec6");
991
992     /// For boolean sustraction
993     double preShape6 = 0.2;
994     for (int i=0;i<2;i++) {
995     v0APts[0+8*i] = -preShape6;                      v0APts[1+8*i] = -fV0AR0+fV0AFraWd/2.-preShape6;
996     v0APts[2+8*i] = -fV0AFraWd/2.;                    v0APts[3+8*i] = -fV0AR0+fV0AFraWd/2.-preShape6;
997     v0APts[4+8*i] = -fV0AFraWd/2.;                    v0APts[5+8*i] = -fV0AR4-fV0AFraWd/2.+preShape6;
998     v0APts[6+8*i] = -preShape6;                      v0APts[7+8*i] = -fV0AR4-fV0AFraWd/2.+preShape6;
999     }
1000     new TGeoArb8("sV0ACha16",fV0ASciWd/1.5,v0APts);
1001     for (int i=0;i<2;i++) {
1002     v0APts[0+8*i] = -fV0AR0*cos665;
1003     v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin665;
1004     v0APts[2+8*i] = -(fV0AR0-fV0AFraWd/2.)*cos665;
1005     v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin665;
1006     v0APts[4+8*i] = -(fV0AR4+fV0AFraWd/2.)*cos665;
1007     v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin665;
1008     v0APts[6+8*i] = -(fV0AR4+fV0AFraWd)*cos665;
1009     v0APts[7+8*i] = -fV0AR4*sin665;
1010     }
1011     new TGeoArb8("sV0ACha26", fV0ASciWd/2.+2.*preShape, v0APts);
1012     new TGeoCompositeShape("sV0ACha126","sV0ACha16+sV0ACha26");
1013     new TGeoTube("sV0ANail16Hole", 0.0, 0.4, 1.65);
1014     TGeoTranslation *pos16 = new TGeoTranslation("pos16",-0.51,-42.9,0.0);
1015     pos16->RegisterYourself();
1016     new TGeoCompositeShape("sV0ACha6","sV0ACha126+sV0ANail16Hole:pos16");
1017
1018     /// Frame
1019     TGeoVolume *v0AFra6 = new TGeoVolumeAssembly("V0AFra6");
1020     for (int i=0;i<2;i++) {
1021     v0APts[0+8*i] = 0.;                              v0APts[1+8*i] = -fV0AR0+fV0AFraWd/2.;
1022     v0APts[2+8*i] = -fV0AFraWd/2.;                    v0APts[3+8*i] = -fV0AR0+fV0AFraWd/2.;
1023     v0APts[4+8*i] = -fV0AFraWd/2.;                    v0APts[5+8*i] = -fV0AR4-fV0AFraWd/2.;
1024     v0APts[6+8*i] = 0.;                              v0APts[7+8*i] = -fV0AR4-fV0AFraWd/2.;
1025     }
1026     TGeoArb8 *sV0AFraB16 = new TGeoArb8("sV0AFraB16",fV0ASciWd/2.,v0APts);
1027     TGeoVolume *v0AFraB16 = new TGeoVolume("V0AFraB16",sV0AFraB16,medV0AFra);
1028     for (int i=0;i<2;i++) {
1029     v0APts[0+8*i] = -fV0AR0*cos665;
1030     v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin665;
1031     v0APts[2+8*i] = -(fV0AR0-fV0AFraWd/2.)*cos665;
1032     v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin665;
1033     v0APts[4+8*i] = -(fV0AR4+fV0AFraWd/2.)*cos665;
1034     v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin665;
1035     v0APts[6+8*i] = -(fV0AR4+fV0AFraWd)*cos665;
1036     v0APts[7+8*i] = -fV0AR4*sin665;
1037     }
1038     TGeoArb8 *sV0AFraB26 = new TGeoArb8("sV0AFraB26", fV0ASciWd/2., v0APts);
1039     TGeoVolume *v0AFraB26 = new TGeoVolume("V0AFraB26",sV0AFraB26,medV0AFra);
1040     v0AFraB16->SetLineColor(kV0AColorFra); v0AFraB26->SetLineColor(kV0AColorFra);
1041     v0AFra6->AddNode(v0AFraB16,1);
1042     v0AFra6->AddNode(v0AFraB26,1);  // Prefer 2 GeoObjects insted of 3 GeoMovements
1043     new TGeoTubeSeg( "sV0AFraR1b6", fV0AR0-fV0AFraWd/2.,
1044              fV0AR0+fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1045     new TGeoTubeSeg( "sV0AFraR2b6", fV0AR1-fV0AFraWd/2.,
1046              fV0AR1+fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1047     new TGeoTubeSeg( "sV0AFraR3b6", fV0AR2-fV0AFraWd/2.,
1048              fV0AR2+fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1049     new TGeoTubeSeg( "sV0AFraR4b6", fV0AR3-fV0AFraWd/2.,
1050              fV0AR3+fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1051     new TGeoTubeSeg( "sV0AFraR5b6", fV0AR4-fV0AFraWd/2.,
1052              fV0AR4+fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1053     TGeoCompositeShape *sV0AFraR16 = new TGeoCompositeShape("sV0AFraR16","sV0AFraR1b6-sV0ACha6");
1054     TGeoCompositeShape *sV0AFraR26 = new TGeoCompositeShape("sV0AFraR26","sV0AFraR2b6-sV0ACha6");
1055     TGeoCompositeShape *sV0AFraR36 = new TGeoCompositeShape("sV0AFraR36","sV0AFraR3b6-sV0ACha6");
1056     TGeoCompositeShape *sV0AFraR46 = new TGeoCompositeShape("sV0AFraR46","sV0AFraR4b6-sV0ACha6");
1057     TGeoCompositeShape *sV0AFraR56 = new TGeoCompositeShape("sV0AFraR56","sV0AFraR5b6-sV0ACha6");
1058     TGeoVolume *v0AFraR16 = new TGeoVolume("V0AFraR16",sV0AFraR16,medV0AFra);
1059     TGeoVolume *v0AFraR26 = new TGeoVolume("V0AFraR26",sV0AFraR26,medV0AFra);
1060     TGeoVolume *v0AFraR36 = new TGeoVolume("V0AFraR36",sV0AFraR36,medV0AFra);
1061     TGeoVolume *v0AFraR46 = new TGeoVolume("V0AFraR46",sV0AFraR46,medV0AFra);
1062     TGeoVolume *v0AFraR56 = new TGeoVolume("V0AFraR56",sV0AFraR56,medV0AFra);
1063     v0AFraR16->SetLineColor(kV0AColorFra); v0AFraR26->SetLineColor(kV0AColorFra);
1064     v0AFraR36->SetLineColor(kV0AColorFra); v0AFraR46->SetLineColor(kV0AColorFra);
1065     v0AFraR56->SetLineColor(kV0AColorFra);
1066     v0AFra6->AddNode(v0AFraR16,1);
1067     v0AFra6->AddNode(v0AFraR26,1);
1068     v0AFra6->AddNode(v0AFraR36,1);
1069     v0AFra6->AddNode(v0AFraR46,1);
1070     v0AFra6->AddNode(v0AFraR56,1);
1071     v0ASec6->AddNode(v0AFra6,1);
1072
1073     /// Sensitive scintilator
1074     TGeoVolume *v0ASci6 = new TGeoVolumeAssembly("V0ASci6");
1075     new TGeoTubeSeg( "sV0AR1b6", fV0AR0+fV0AFraWd/2.,
1076                      fV0AR1-fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1077     new TGeoTubeSeg( "sV0AR2b6", fV0AR1+fV0AFraWd/2.,
1078              fV0AR2-fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1079     new TGeoTubeSeg( "sV0AR3b6", fV0AR2+fV0AFraWd/2.,
1080              fV0AR3-fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1081     new TGeoTubeSeg( "sV0AR4b6", fV0AR3+fV0AFraWd/2.,
1082              fV0AR4-fV0AFraWd/2., fV0ASciWd/2., 245.4, 270.0);
1083     TGeoCompositeShape *sV0AR16 = new TGeoCompositeShape("sV0AR16","sV0AR1b6-sV0ACha6");
1084     TGeoCompositeShape *sV0AR26 = new TGeoCompositeShape("sV0AR26","sV0AR2b6-sV0ACha6");
1085     TGeoCompositeShape *sV0AR36 = new TGeoCompositeShape("sV0AR36","sV0AR3b6-sV0ACha6");
1086     TGeoCompositeShape *sV0AR46 = new TGeoCompositeShape("sV0AR46","sV0AR4b6-sV0ACha6");
1087     TGeoVolume *v0L16 = new TGeoVolume("V0L16",sV0AR16,medV0ASci);
1088     TGeoVolume *v0L26 = new TGeoVolume("V0L26",sV0AR26,medV0ASci);
1089     TGeoVolume *v0L36 = new TGeoVolume("V0L36",sV0AR36,medV0ASci);
1090     TGeoVolume *v0L46 = new TGeoVolume("V0L46",sV0AR46,medV0ASci);
1091     v0L16->SetLineColor(kV0AColorSci); v0L26->SetLineColor(kV0AColorSci);
1092     v0L36->SetLineColor(kV0AColorSci); v0L46->SetLineColor(kV0AColorSci);
1093     v0ASci6->AddNode(v0L16,1);
1094     v0ASci6->AddNode(v0L26,1);
1095     v0ASci6->AddNode(v0L36,1);
1096     v0ASci6->AddNode(v0L46,1);
1097
1098     /// Non-sensitive scintilator
1099     new TGeoTubeSeg("sV0AR5S26", fV0AR4+fV0AFraWd/2., fV0AR4 + fV0AR0, fV0ASciWd/2.+2*preShape, 245.4, 270.0);
1100     TGeoCompositeShape *sV0AR56 = new TGeoCompositeShape("V0AR56","sV0AR5S26 - sV0ACha6");
1101     TGeoVolume *v0AR56 = new TGeoVolume("V0AR56",sV0AR56,medV0ASci);
1102     v0AR56->SetLineColor(kV0AColorSci);
1103     v0ASci6->AddNode(v0AR56,1);
1104     v0ASec6->AddNode(v0ASci6,1);
1105
1106     /// Segment of innermost octagon
1107     TGeoVolume *v0ASup6 = new TGeoVolumeAssembly("V0ASup6");
1108
1109     /// Segment of outtermost octagon   
1110     for (int i=0;i<2;i++) {
1111     v0APts[0+8*i] = 0.;                          v0APts[1+8*i] = -(fV0AR7-fV0AOctH2)*sin654;
1112     v0APts[2+8*i] = 0.;                          v0APts[3+8*i] = -fV0AR7*sin654;
1113     v0APts[4+8*i] = -fV0AR7*cos654;              v0APts[5+8*i] = -fV0AR7*sin654;
1114     v0APts[6+8*i] = -(fV0AR7-fV0AOctH2)*cos654;  v0APts[7+8*i] = -(fV0AR7-fV0AOctH2)*sin654;
1115     }
1116     TGeoArb8 *sV0AOct26 = new TGeoArb8("sV0AOct26", (fV0ASciWd+2*fV0AOctWd)/2., v0APts);
1117     TGeoVolume *v0AOct26 = new TGeoVolume("V0AOct26", sV0AOct26,medV0ASup);
1118     v0AOct26->SetLineColor(kV0AColorOct);
1119     v0ASup6->AddNode(v0AOct26,1);
1120     v0ASec6->AddNode(v0ASup6,1);
1121
1122
1123     /// Plates
1124     new TGeoTube("sV0ANail1PlaInHole6", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);    
1125     new TGeoTube("sV0ANail1PlaOuHole6", 0.0, 0.4, (fV0APlaAl)/2.);
1126     for (int i=0;i<2;i++) {
1127     v0APts[0+8*i] = 0.;                 v0APts[1+8*i] = -fV0AR0;
1128     v0APts[2+8*i] = 0.;                 v0APts[3+8*i] = -fV0AR7*sin654;
1129     v0APts[4+8*i] = -fV0AR7*cos654;     v0APts[5+8*i] = -fV0AR7*sin654;
1130     v0APts[6+8*i] = -fV0AR0*cos654;             v0APts[7+8*i] = -fV0AR0*sin654;
1131     }
1132     new TGeoArb8("sV0APlaIn6", (fV0APlaWd-2*fV0APlaAl)/2., v0APts);
1133     TGeoCompositeShape *sV0APlaInNailsHoles6 = new TGeoCompositeShape("sV0APlaInNailsHoles6","sV0APlaIn6-sV0ANail1PlaInHole6:pos16");
1134     TGeoVolume *v0APlaInNailsHoles6 = new TGeoVolume("V0APlaInNailsHoles6", sV0APlaInNailsHoles6, medV0APlaIn);
1135     new TGeoArb8("sV0APlaOu6", fV0APlaAl/2., v0APts);
1136     TGeoCompositeShape *sV0APlaOuNailsHoles6 = new TGeoCompositeShape("sV0APlaOuNailsHoles6","sV0APlaOu6-sV0ANail1PlaOuHole6:pos16"); 
1137     TGeoVolume *v0APlaOuNailsHoles6 = new TGeoVolume("V0APlaOuNailsHoles6", sV0APlaOuNailsHoles6, medV0APlaOu);
1138     v0APlaInNailsHoles6->SetLineColor(kV0AColorPlaIn); v0APlaOuNailsHoles6->SetLineColor(kV0AColorPlaOu);
1139     TGeoVolume *v0APla6 = new TGeoVolumeAssembly("V0APla6");
1140     v0APla6->AddNode(v0APlaInNailsHoles6,1);
1141     v0APla6->AddNode(v0APlaOuNailsHoles6,1,new TGeoTranslation(0,0,(fV0APlaWd-fV0APlaAl)/2.));
1142     v0APla6->AddNode(v0APlaOuNailsHoles6,2,new TGeoTranslation(0,0,-(fV0APlaWd-fV0APlaAl)/2.));
1143     v0ASec6->AddNode(v0APla6,1,new TGeoTranslation(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1144     v0ASec6->AddNode(v0APla6,2,new TGeoTranslation(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1145     
1146     /// Support
1147     TGeoBBox *sV0ASuppbl = new TGeoBBox("sV0ASuppbl", 2.0, 18.137, 2.0);       
1148     TGeoVolume *v0ASuppbl = new TGeoVolume("V0ASuppbl", sV0ASuppbl, medV0ASup);
1149     v0ASuppbl->SetLineColor(kV0AColorOct);
1150     v0ASec6->AddNode(v0ASuppbl,1,new TGeoTranslation(-2.0,-63.635,0.0));
1151
1152
1153     // Replicate sectors
1154     //for(int i=0; i<1; i++) {
1155     //TGeoRotation *rot = new TGeoRotation("rot", 90., i*45., 90., 90.+i*45., 0., 0.);
1156     v0LE->AddNode(v0ASec6, 1);   //i +1, rot);
1157     //}
1158
1159
1160     //Definition of sector 7
1161    
1162     TGeoVolume *v0ASec7 = new TGeoVolumeAssembly("V0ASec7");
1163
1164     /// For boolean sustraction
1165     double preShape7 = 0.2;
1166     for (int i=0;i<2;i++) {
1167     v0APts[0+8*i] = 0.0;                      v0APts[1+8*i] = -fV0AR0+fV0AFraWd/2.-preShape7;
1168     v0APts[2+8*i] = fV0AFraWd/2.;                    v0APts[3+8*i] = -fV0AR0+fV0AFraWd/2.-preShape7;
1169     v0APts[4+8*i] = fV0AFraWd/2.;                    v0APts[5+8*i] = -fV0AR4-fV0AFraWd/2.+preShape7;
1170     v0APts[6+8*i] = 0.0;                      v0APts[7+8*i] = -fV0AR4-fV0AFraWd/2.+preShape7;
1171     }
1172     new TGeoArb8("sV0ACha17",fV0ASciWd/1.5,v0APts);
1173     for (int i=0;i<2;i++) {
1174     v0APts[0+8*i] = fV0AR0*cos654-preShape7;
1175     v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin654-preShape7;
1176     v0APts[2+8*i] = (fV0AR0-fV0AFraWd/2.)*cos654-preShape7;
1177     v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin654;
1178     v0APts[4+8*i] = (fV0AR4+fV0AFraWd/2.)*cos654+preShape7;
1179     v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin654+2.*preShape7;
1180     v0APts[6+8*i] = (fV0AR4+fV0AFraWd)*cos654+preShape7;
1181     v0APts[7+8*i] = -fV0AR4*sin654+preShape7;
1182     }
1183     new TGeoArb8("sV0ACha27", fV0ASciWd/2.+2.*preShape, v0APts);
1184     new TGeoCompositeShape("sV0ACha127","sV0ACha17+sV0ACha27");
1185     new TGeoTube("sV0ANail17Hole", 0.0, 0.4, 1.65);
1186     TGeoTranslation *pos17 = new TGeoTranslation("pos17",0.51,-42.9,0.0);
1187     pos17->RegisterYourself();
1188     new TGeoCompositeShape("sV0ACha7","sV0ACha127+sV0ANail17Hole:pos17");
1189
1190     /// Frame
1191     TGeoVolume *v0AFra7 = new TGeoVolumeAssembly("V0AFra7");
1192     for (int i=0;i<2;i++) {
1193     v0APts[0+8*i] = 0.;                              v0APts[1+8*i] = -fV0AR0-fV0AFraWd/2.;
1194     v0APts[2+8*i] = fV0AFraWd/2.;                    v0APts[3+8*i] = -fV0AR0-fV0AFraWd/2.;
1195     v0APts[4+8*i] = fV0AFraWd/2.;                    v0APts[5+8*i] = -fV0AR4+fV0AFraWd/2.;
1196     v0APts[6+8*i] = 0.;                              v0APts[7+8*i] = -fV0AR4+fV0AFraWd/2.;
1197     }
1198     TGeoArb8 *sV0AFraB17 = new TGeoArb8("sV0AFraB17",fV0ASciWd/2.,v0APts);
1199     TGeoVolume *v0AFraB17 = new TGeoVolume("V0AFraB17",sV0AFraB17,medV0AFra);
1200     for (int i=0;i<2;i++) {
1201     v0APts[0+8*i] = fV0AR0*cos654;
1202     v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin654;
1203     v0APts[2+8*i] = (fV0AR0-fV0AFraWd/2.)*cos654;
1204     v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin654;
1205     v0APts[4+8*i] = (fV0AR4+fV0AFraWd/2.)*cos654;
1206     v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin654;
1207     v0APts[6+8*i] = (fV0AR4+fV0AFraWd)*cos654;
1208     v0APts[7+8*i] = -fV0AR4*sin654;
1209     }
1210     TGeoArb8 *sV0AFraB27 = new TGeoArb8("sV0AFraB27", fV0ASciWd/2., v0APts);
1211     TGeoVolume *v0AFraB27 = new TGeoVolume("V0AFraB27",sV0AFraB27,medV0AFra);
1212     v0AFraB17->SetLineColor(kV0AColorFra); v0AFraB27->SetLineColor(kV0AColorFra);
1213     v0AFra7->AddNode(v0AFraB17,1);
1214     v0AFra7->AddNode(v0AFraB27,1);  // Prefer 2 GeoObjects insted of 3 GeoMovements
1215     new TGeoTubeSeg( "sV0AFraR1b7", fV0AR0-fV0AFraWd/2.,
1216              fV0AR0+fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1217     new TGeoTubeSeg( "sV0AFraR2b7", fV0AR1-fV0AFraWd/2.,
1218              fV0AR1+fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1219     new TGeoTubeSeg( "sV0AFraR3b7", fV0AR2-fV0AFraWd/2.,
1220              fV0AR2+fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1221     new TGeoTubeSeg( "sV0AFraR4b7", fV0AR3-fV0AFraWd/2.,
1222              fV0AR3+fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1223     new TGeoTubeSeg( "sV0AFraR5b7", fV0AR4-fV0AFraWd/2.,
1224              fV0AR4+fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1225     TGeoCompositeShape *sV0AFraR17 = new TGeoCompositeShape("sV0AFraR17","sV0AFraR1b7-sV0ACha7");
1226     TGeoCompositeShape *sV0AFraR27 = new TGeoCompositeShape("sV0AFraR27","sV0AFraR2b7-sV0ACha7");
1227     TGeoCompositeShape *sV0AFraR37 = new TGeoCompositeShape("sV0AFraR37","sV0AFraR3b7-sV0ACha7");
1228     TGeoCompositeShape *sV0AFraR47 = new TGeoCompositeShape("sV0AFraR47","sV0AFraR4b7-sV0ACha7");
1229     TGeoCompositeShape *sV0AFraR57 = new TGeoCompositeShape("sV0AFraR57","sV0AFraR5b7-sV0ACha7");
1230     TGeoVolume *v0AFraR17 = new TGeoVolume("V0AFraR17",sV0AFraR17,medV0AFra);
1231     TGeoVolume *v0AFraR27 = new TGeoVolume("V0AFraR27",sV0AFraR27,medV0AFra);
1232     TGeoVolume *v0AFraR37 = new TGeoVolume("V0AFraR37",sV0AFraR37,medV0AFra);
1233     TGeoVolume *v0AFraR47 = new TGeoVolume("V0AFraR47",sV0AFraR47,medV0AFra);
1234     TGeoVolume *v0AFraR57 = new TGeoVolume("V0AFraR57",sV0AFraR57,medV0AFra);
1235     v0AFraR17->SetLineColor(kV0AColorFra); v0AFraR27->SetLineColor(kV0AColorFra);
1236     v0AFraR37->SetLineColor(kV0AColorFra); v0AFraR47->SetLineColor(kV0AColorFra);
1237     v0AFraR57->SetLineColor(kV0AColorFra);
1238     v0AFra7->AddNode(v0AFraR17,1);
1239     v0AFra7->AddNode(v0AFraR27,1);
1240     v0AFra7->AddNode(v0AFraR37,1);
1241     v0AFra7->AddNode(v0AFraR47,1);
1242     v0AFra7->AddNode(v0AFraR57,1);
1243     v0ASec7->AddNode(v0AFra7,1);
1244
1245     /// Sensitive scintilator
1246     TGeoVolume *v0ASci7 = new TGeoVolumeAssembly("V0ASci7");
1247     new TGeoTubeSeg( "sV0AR1b7", fV0AR0+fV0AFraWd/2.,
1248              fV0AR1-fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1249     new TGeoTubeSeg( "sV0AR2b7", fV0AR1+fV0AFraWd/2.,
1250              fV0AR2-fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1251     new TGeoTubeSeg( "sV0AR3b7", fV0AR2+fV0AFraWd/2.,
1252              fV0AR3-fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1253     new TGeoTubeSeg( "sV0AR4b7", fV0AR3+fV0AFraWd/2.,
1254              fV0AR4-fV0AFraWd/2., fV0ASciWd/2., 270.0, 294.52);
1255     TGeoCompositeShape *sV0AR17 = new TGeoCompositeShape("sV0AR17","sV0AR1b7-sV0ACha7");
1256     TGeoCompositeShape *sV0AR27 = new TGeoCompositeShape("sV0AR27","sV0AR2b7-sV0ACha7");
1257     TGeoCompositeShape *sV0AR37 = new TGeoCompositeShape("sV0AR37","sV0AR3b7-sV0ACha7");
1258     TGeoCompositeShape *sV0AR47 = new TGeoCompositeShape("sV0AR47","sV0AR4b7-sV0ACha7");
1259     TGeoVolume *v0L17 = new TGeoVolume("V0L17",sV0AR17,medV0ASci);
1260     TGeoVolume *v0L27 = new TGeoVolume("V0L27",sV0AR27,medV0ASci);
1261     TGeoVolume *v0L37 = new TGeoVolume("V0L37",sV0AR37,medV0ASci);
1262     TGeoVolume *v0L47 = new TGeoVolume("V0L47",sV0AR47,medV0ASci);
1263     v0L17->SetLineColor(kV0AColorSci); v0L27->SetLineColor(kV0AColorSci);
1264     v0L37->SetLineColor(kV0AColorSci); v0L47->SetLineColor(kV0AColorSci);
1265     v0ASci7->AddNode(v0L17,1);
1266     v0ASci7->AddNode(v0L27,1);
1267     v0ASci7->AddNode(v0L37,1);
1268     v0ASci7->AddNode(v0L47,1);
1269
1270     /// Non-sensitive scintilator
1271     new TGeoTubeSeg("sV0AR5S27", fV0AR4+fV0AFraWd/2., fV0AR4 + fV0AR0, fV0ASciWd/2.+2*preShape, 270.0, 294.52);
1272     TGeoCompositeShape *sV0AR57 = new TGeoCompositeShape("V0AR57","sV0AR5S27 - sV0ACha7");
1273     TGeoVolume *v0AR57 = new TGeoVolume("V0AR57",sV0AR57,medV0ASci);
1274     v0AR57->SetLineColor(kV0AColorSci);
1275     v0ASci7->AddNode(v0AR57,1);
1276     v0ASec7->AddNode(v0ASci7,1);
1277
1278     /// Segment of innermost octagon
1279     TGeoVolume *v0ASup7 = new TGeoVolumeAssembly("V0ASup7");
1280
1281     /// Segment of outtermost octagon   
1282     for (int i=0;i<2;i++) {
1283     v0APts[0+8*i] = 0.;                          v0APts[1+8*i] = -(fV0AR7-fV0AOctH2)*sin654;
1284     v0APts[2+8*i] = 0.;                          v0APts[3+8*i] = -fV0AR7*sin654;
1285     v0APts[4+8*i] = fV0AR7*cos654;               v0APts[5+8*i] = -fV0AR7*sin654;
1286     v0APts[6+8*i] = (fV0AR7-fV0AOctH2)*cos654;   v0APts[7+8*i] = -(fV0AR7-fV0AOctH2)*sin654;
1287     }
1288     TGeoArb8 *sV0AOct27 = new TGeoArb8("sV0AOct27", (fV0ASciWd+2*fV0AOctWd)/2., v0APts);
1289     TGeoVolume *v0AOct27 = new TGeoVolume("V0AOct27", sV0AOct27,medV0ASup);
1290     v0AOct27->SetLineColor(kV0AColorOct);
1291     v0ASup7->AddNode(v0AOct27,1);
1292     v0ASec7->AddNode(v0ASup7,1);
1293
1294
1295     /// Plates
1296     new TGeoTube("sV0ANail1PlaInHole7", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);    
1297     new TGeoTube("sV0ANail1PlaOuHole7", 0.0, 0.4, (fV0APlaAl)/2.);    
1298     for (int i=0;i<2;i++) {
1299     v0APts[0+8*i] = 0.;                 v0APts[1+8*i] = -fV0AR0;
1300     v0APts[2+8*i] = 0.;                 v0APts[3+8*i] = -fV0AR7*sin654;
1301     v0APts[4+8*i] = fV0AR7*cos654;      v0APts[5+8*i] = -fV0AR7*sin654;
1302     v0APts[6+8*i] = fV0AR0*cos654;              v0APts[7+8*i] = -fV0AR0*sin654;
1303     }
1304     new TGeoArb8("sV0APlaIn7", (fV0APlaWd-2*fV0APlaAl)/2., v0APts);
1305     TGeoCompositeShape *sV0APlaInNailsHoles7 = new TGeoCompositeShape("sV0APlaInNailsHoles7","sV0APlaIn7-sV0ANail1PlaInHole7:pos17");
1306     TGeoVolume *v0APlaInNailsHoles7 = new TGeoVolume("V0APlaInNailsHoles7", sV0APlaInNailsHoles7, medV0APlaIn);
1307     new TGeoArb8("sV0APlaOu7", fV0APlaAl/2., v0APts);
1308     TGeoCompositeShape *sV0APlaOuNailsHoles7 = new TGeoCompositeShape("sV0APlaOuNailsHoles7","sV0APlaOu7-sV0ANail1PlaOuHole7:pos17"); 
1309     TGeoVolume *v0APlaOuNailsHoles7 = new TGeoVolume("V0APlaOuNailsHoles7", sV0APlaOuNailsHoles7, medV0APlaOu);
1310     v0APlaInNailsHoles7->SetLineColor(kV0AColorPlaIn); v0APlaOuNailsHoles7->SetLineColor(kV0AColorPlaOu);
1311     TGeoVolume *v0APla7 = new TGeoVolumeAssembly("V0APla7");
1312     v0APla7->AddNode(v0APlaInNailsHoles7,1);
1313     v0APla7->AddNode(v0APlaOuNailsHoles7,1,new TGeoTranslation(0,0,(fV0APlaWd-fV0APlaAl)/2.));
1314     v0APla7->AddNode(v0APlaOuNailsHoles7,2,new TGeoTranslation(0,0,-(fV0APlaWd-fV0APlaAl)/2.));
1315     v0ASec7->AddNode(v0APla7,1,new TGeoTranslation(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1316     v0ASec7->AddNode(v0APla7,2,new TGeoTranslation(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1317
1318      /// Support
1319     TGeoBBox *sV0ASuppbr = new TGeoBBox("sV0ASuppbr", 2.0, 18.138, 2.0);      
1320     TGeoVolume *v0ASuppbr = new TGeoVolume("V0ASuppbr", sV0ASuppbr, medV0ASup);
1321     v0ASuppbr->SetLineColor(kV0AColorOct);
1322     v0ASec7->AddNode(v0ASuppbr,1,new TGeoTranslation(2.0,-63.639,0.0));
1323
1324     // Replicate sectors
1325     for(int i=0; i<1; i++) {
1326     TGeoRotation *rot = new TGeoRotation("rot", 90., i*45., 90., 90.+i*45., 0., 0.);
1327     v0LE->AddNode(v0ASec7,i + 1,rot);
1328     }
1329
1330
1331
1332     
1333     //Definition of sector 8
1334    
1335        TGeoVolume *v0ASec8 = new TGeoVolumeAssembly("V0ASec8");
1336
1337  /// For boolean sustraction
1338       for (int i=0;i<2;i++) {
1339       v0APts[0+8*i] = fV0AR0-fV0AFraWd/2.-preShape;  v0APts[1+8*i] = -preShape;
1340       v0APts[2+8*i] = fV0AR0-fV0AFraWd/2.-preShape;  v0APts[3+8*i] = fV0AFraWd/2.;
1341       v0APts[4+8*i] = fV0AR4+fV0AFraWd/2.+preShape;  v0APts[5+8*i] = fV0AFraWd/2.;
1342       v0APts[6+8*i] = fV0AR4+fV0AFraWd/2.+preShape;  v0APts[7+8*i] = -preShape;
1343     }
1344     new TGeoArb8("sV0ACha18",fV0ASciWd/1.5,v0APts);
1345     for (int i=0;i<2;i++) {
1346       v0APts[0+8*i] = fV0AR0*cos654-preShape;
1347       v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin654-preShape;
1348       v0APts[2+8*i] = (fV0AR0-fV0AFraWd/2.)*cos654-preShape;
1349       v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin654;
1350       v0APts[4+8*i] = (fV0AR4+fV0AFraWd/2.)*cos654+preShape;
1351       v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin654+2.*preShape;
1352       v0APts[6+8*i] = (fV0AR4+fV0AFraWd)*cos654+preShape;
1353       v0APts[7+8*i] = -fV0AR4*sin654+preShape;
1354     }
1355     new TGeoArb8("sV0ACha28", fV0ASciWd/2.+2.*preShape, v0APts);
1356     new TGeoCompositeShape("sV0ACha128","sV0ACha18+sV0ACha28");
1357     new TGeoTube("sV0ANail18Hole", 0.0, 0.4, 1.65);
1358     TGeoTranslation *pos18 = new TGeoTranslation("pos18",42.9,-.51,0.0);
1359     pos18->RegisterYourself();
1360     new TGeoTube("sV0ANail28Hole", 0.0, 0.4, 1.65);
1361     TGeoTranslation *pos28 = new TGeoTranslation("pos28", 30.8,-30.04,0.0);
1362     pos28->RegisterYourself();
1363     new TGeoTube("sV0ANail38Hole", 0.0, 0.4, 1.65);
1364     TGeoTranslation *pos38 = new TGeoTranslation("pos38",29.8,-31.04,0.0); 
1365     pos38->RegisterYourself();
1366     new TGeoCompositeShape("sV0ANailsHoles8","sV0ANail18Hole:pos18+sV0ANail28Hole:pos28+sV0ANail38Hole:pos38");
1367     new TGeoCompositeShape("sV0ACha8","sV0ACha128+sV0ANailsHoles8");
1368
1369     /// Frame
1370     TGeoVolume *v0AFra8 = new TGeoVolumeAssembly("V0AFra8");
1371     for (int i=0;i<2;i++) {
1372       v0APts[0+8*i] = fV0AR0-fV0AFraWd/2.;  v0APts[1+8*i] = preShape;
1373       v0APts[2+8*i] = fV0AR0-fV0AFraWd/2.;  v0APts[3+8*i] = fV0AFraWd/2.;
1374       v0APts[4+8*i] = fV0AR4+fV0AFraWd/2.;  v0APts[5+8*i] = fV0AFraWd/2.;
1375       v0APts[6+8*i] = fV0AR4+fV0AFraWd/2.;  v0APts[7+8*i] = preShape;
1376     }
1377     TGeoArb8 *sV0AFraB18 = new TGeoArb8("sV0AFraB18",fV0ASciWd/2.,v0APts);
1378     TGeoVolume *v0AFraB18 = new TGeoVolume("V0AFraB18",sV0AFraB18,medV0AFra);
1379     for (int i=0;i<2;i++) {
1380       v0APts[0+8*i] = fV0AR0*cos654;
1381       v0APts[1+8*i] = -(fV0AR0-fV0AFraWd)*sin654;
1382       v0APts[2+8*i] = (fV0AR0-fV0AFraWd/2.)*cos654;
1383       v0APts[3+8*i] = -(fV0AR0-fV0AFraWd/2.)*sin654;
1384       v0APts[4+8*i] = (fV0AR4+fV0AFraWd/2.)*cos654;
1385       v0APts[5+8*i] = -(fV0AR4+fV0AFraWd/2.)*sin654;
1386       v0APts[6+8*i] = (fV0AR4+fV0AFraWd)*cos654;
1387       v0APts[7+8*i] = -fV0AR4*sin654;
1388     }
1389     TGeoArb8 *sV0AFraB28 = new TGeoArb8("sV0AFraB28", fV0ASciWd/2., v0APts);
1390     TGeoVolume *v0AFraB28 = new TGeoVolume("V0AFraB28",sV0AFraB28,medV0AFra);
1391     v0AFraB18->SetLineColor(kV0AColorFra); v0AFraB28->SetLineColor(kV0AColorFra);
1392     v0AFra8->AddNode(v0AFraB18,1);
1393     v0AFra8->AddNode(v0AFraB28,1);  // Prefer 2 GeoObjects insted of 3 GeoMovements
1394     new TGeoTubeSeg( "sV0AFraR1b8", fV0AR0-fV0AFraWd/2.,
1395                      fV0AR0+fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1396     new TGeoTubeSeg( "sV0AFraR2b8", fV0AR1-fV0AFraWd/2.,
1397                      fV0AR1+fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1398     new TGeoTubeSeg( "sV0AFraR3b8", fV0AR2-fV0AFraWd/2.,
1399                      fV0AR2+fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1400     new TGeoTubeSeg( "sV0AFraR4b8", fV0AR3-fV0AFraWd/2.,
1401                      fV0AR3+fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1402     new TGeoTubeSeg( "sV0AFraR5b8", fV0AR4-fV0AFraWd/2.,
1403                      fV0AR4+fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1404     TGeoCompositeShape *sV0AFraR18 = new TGeoCompositeShape("sV0AFraR18","sV0AFraR1b8-sV0ACha8");
1405     TGeoCompositeShape *sV0AFraR28 = new TGeoCompositeShape("sV0AFraR28","sV0AFraR2b8-sV0ACha8");
1406     TGeoCompositeShape *sV0AFraR38 = new TGeoCompositeShape("sV0AFraR38","sV0AFraR3b8-sV0ACha8");
1407     TGeoCompositeShape *sV0AFraR48 = new TGeoCompositeShape("sV0AFraR48","sV0AFraR4b8-sV0ACha8");
1408     TGeoCompositeShape *sV0AFraR58 = new TGeoCompositeShape("sV0AFraR58","sV0AFraR5b8-sV0ACha8");
1409     TGeoVolume *v0AFraR18 = new TGeoVolume("V0AFraR18",sV0AFraR18,medV0AFra);
1410     TGeoVolume *v0AFraR28 = new TGeoVolume("V0AFraR28",sV0AFraR28,medV0AFra);
1411     TGeoVolume *v0AFraR38 = new TGeoVolume("V0AFraR38",sV0AFraR38,medV0AFra);
1412     TGeoVolume *v0AFraR48 = new TGeoVolume("V0AFraR48",sV0AFraR48,medV0AFra);
1413     TGeoVolume *v0AFraR58 = new TGeoVolume("V0AFraR58",sV0AFraR58,medV0AFra);
1414     v0AFraR18->SetLineColor(kV0AColorFra); v0AFraR28->SetLineColor(kV0AColorFra);
1415     v0AFraR38->SetLineColor(kV0AColorFra); v0AFraR48->SetLineColor(kV0AColorFra);
1416     v0AFraR58->SetLineColor(kV0AColorFra);
1417     v0AFra8->AddNode(v0AFraR18,1);
1418     v0AFra8->AddNode(v0AFraR28,1);
1419     v0AFra8->AddNode(v0AFraR38,1);
1420     v0AFra8->AddNode(v0AFraR48,1);
1421     v0AFra8->AddNode(v0AFraR58,1);
1422     v0ASec8->AddNode(v0AFra8,1);
1423
1424     /// Sensitive scintilator
1425     TGeoVolume *v0ASci8 = new TGeoVolumeAssembly("V0ASci8");
1426     new TGeoTubeSeg( "sV0AR1b8", fV0AR0+fV0AFraWd/2.,
1427                      fV0AR1-fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1428     new TGeoTubeSeg( "sV0AR2b8", fV0AR1+fV0AFraWd/2.,
1429                      fV0AR2-fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1430     new TGeoTubeSeg( "sV0AR3b8", fV0AR2+fV0AFraWd/2.,
1431                      fV0AR3-fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1432     new TGeoTubeSeg( "sV0AR4b8", fV0AR3+fV0AFraWd/2.,
1433                      fV0AR4-fV0AFraWd/2., fV0ASciWd/2., 294.52, 359.0);
1434     TGeoCompositeShape *sV0AR18 = new TGeoCompositeShape("sV0AR18","sV0AR1b8-sV0ACha8");
1435     TGeoCompositeShape *sV0AR28 = new TGeoCompositeShape("sV0AR28","sV0AR2b8-sV0ACha8");
1436     TGeoCompositeShape *sV0AR38 = new TGeoCompositeShape("sV0AR38","sV0AR3b8-sV0ACha8");
1437     TGeoCompositeShape *sV0AR48 = new TGeoCompositeShape("sV0AR48","sV0AR4b8-sV0ACha8");
1438     TGeoVolume *v0L18 = new TGeoVolume("V0L18",sV0AR18,medV0ASci);
1439     TGeoVolume *v0L28 = new TGeoVolume("V0L28",sV0AR28,medV0ASci);
1440     TGeoVolume *v0L38 = new TGeoVolume("V0L38",sV0AR38,medV0ASci);
1441     TGeoVolume *v0L48 = new TGeoVolume("V0L48",sV0AR48,medV0ASci);
1442     v0L18->SetLineColor(kV0AColorSci); v0L28->SetLineColor(kV0AColorSci);
1443     v0L38->SetLineColor(kV0AColorSci); v0L48->SetLineColor(kV0AColorSci);
1444     v0ASci8->AddNode(v0L18,1);
1445     v0ASci8->AddNode(v0L28,1);
1446     v0ASci8->AddNode(v0L38,1);
1447     v0ASci8->AddNode(v0L48,1);
1448
1449     /// Non-sensitive scintilator
1450     new TGeoTubeSeg("sV0AR5S28", fV0AR4+fV0AFraWd/2., fV0AR4 + fV0AR0, fV0ASciWd/2.+2*preShape, 294.52, 359.0);
1451     TGeoCompositeShape *sV0AR58 = new TGeoCompositeShape("V0AR58","sV0AR5S28 - sV0ACha8");
1452     TGeoVolume *v0AR58 = new TGeoVolume("V0AR58",sV0AR58,medV0ASci);
1453     v0AR58->SetLineColor(kV0AColorSci);
1454     v0ASci8->AddNode(v0AR58,1);
1455     v0ASec8->AddNode(v0ASci8,1);
1456
1457     /// Segment of innermost octagon
1458     TGeoVolume *v0ASup8 = new TGeoVolumeAssembly("V0ASup8");
1459
1460 /// Segment of outtermost octagon   
1461      for (int i=0;i<2;i++) {
1462     v0APts[0+8*i] = fV0AR6-fV0AOctH2;            v0APts[1+8*i] = 0.;
1463     v0APts[2+8*i] = (fV0AR7-fV0AOctH2)*cos654;   v0APts[3+8*i] = -(fV0AR7-fV0AOctH2)*sin654;
1464     v0APts[4+8*i] = fV0AR7*cos654;               v0APts[5+8*i] = -fV0AR7*sin654;
1465     v0APts[6+8*i] = fV0AR6;                      v0APts[7+8*i] = 0.;
1466     }
1467     TGeoArb8 *sV0AOct28 = new TGeoArb8("sV0AOct28", (fV0ASciWd+2*fV0AOctWd)/2., v0APts);
1468     TGeoVolume *v0AOct28 = new TGeoVolume("V0AOct28", sV0AOct28,medV0ASup);
1469     v0AOct28->SetLineColor(kV0AColorOct);
1470     v0ASup8->AddNode(v0AOct28,1);
1471     v0ASec8->AddNode(v0ASup8,1);
1472
1473     //Bunch of fibers
1474     v0APts[ 0] = v0APts[ 2] = -14.0;
1475     v0APts[ 1] = v0APts[ 7] = (fV0ASciWd+fV0AOctWd)/2.-0.01;
1476     v0APts[ 3] = v0APts[ 5] = (fV0ASciWd+fV0AOctWd)/2.+0.01;
1477     v0APts[ 4] = v0APts[ 6] = +14.0;
1478     v0APts[ 8] = v0APts[10] = -10.0;
1479     v0APts[ 9] = v0APts[15] = 0.;
1480     v0APts[11] = v0APts[13] = 0.25;
1481     v0APts[12] = v0APts[14] = +10.0;
1482     TGeoArb8 *sV0AFiba8 = new TGeoArb8("sV0AFiba8", 9.0, v0APts);
1483     TGeoVolume *v0AFiba18 = new TGeoVolume("V0AFiba18",sV0AFiba8,medV0AFib);
1484     TGeoVolume *v0AFiba8 = new TGeoVolumeAssembly("V0AFiba8");
1485     rot = new TGeoRotation("rot");
1486     rot->RotateX(-90);
1487     rot->RotateZ(-90-22.5);
1488     v0AFiba8->AddNode(v0AFiba18,1,rot);
1489     rot = new TGeoRotation("rot");
1490     rot->RotateX(-90);
1491     rot->RotateY(180);
1492     rot->RotateZ(-90-22.5);
1493     v0AFiba8->SetLineColor(kV0AColorFib);
1494     v0AFiba8->AddNode(v0AFiba18,2,rot);
1495     v0ASec8->AddNode(v0AFiba8,1,new TGeoTranslation((fV0AR6-fV0AOctH2+fV0AR5)*cos225/2. - 1.5, -(fV0AR6-fV0AOctH2+fV0AR5)*sin225/2. - 1.0, 0));
1496     for (int i=0;i<2;i++) {
1497     v0APts[0+8*i] = 33.47;            v0APts[1+8*i] = -35.13;
1498     v0APts[2+8*i] = 50.48;            v0APts[3+8*i] = -46.97;
1499     v0APts[4+8*i] = 44.66;            v0APts[5+8*i] = -61.92;
1500     v0APts[6+8*i] = 7.9;              v0APts[7+8*i] = -47.18;
1501     }
1502     TGeoArb8 *sV0AFibb8 = new TGeoArb8("sV0AFibb8", 1.25, v0APts);
1503     TGeoVolume *v0AFibb18 = new TGeoVolume("V0AFibb18",sV0AFibb8,medV0AFib);
1504     TGeoVolume *v0AFibb8 = new TGeoVolumeAssembly("V0AFibb8");
1505     //rot = new TGeoRotation("rot");
1506     //rot->RotateX(-90);
1507     //rot->RotateZ(-90-40);
1508     v0AFibb8->AddNode(v0AFibb18,1);//,rot);
1509     //rot = new TGeoRotation("rot");
1510     //rot->RotateX(-90);
1511     //rot->RotateY(180);
1512     //rot->RotateZ(-90-40);
1513     v0AFibb8->SetLineColor(kV0AColorFib);
1514     //v0AFibb8->AddNode(v0AFibb18,2,rot);
1515     v0ASec8->AddNode(v0AFibb8,1);//,new TGeoTranslation((fV0AR6-fV0AOctH2+fV0AR5)*cos225/2. - 17.0, -(fV0AR6-fV0AOctH2+fV0AR5)*sin225/2. - 30.0, 0));
1516     
1517
1518     /// Plates
1519     new TGeoTube("sV0ANail1PlaInHole8", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
1520     new TGeoTube("sV0ANail2PlaInHole8", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
1521     new TGeoTube("sV0ANail3PlaInHole8", 0.0, 0.4, (fV0APlaWd-2*fV0APlaAl)/2.);
1522     new TGeoCompositeShape("sV0ANailsPlaInHoles8","sV0ANail1PlaInHole8:pos18+sV0ANail2PlaInHole8:pos28+sV0ANail3PlaInHole8:pos38");
1523     new TGeoTube("sV0ANail1PlaOuHole8", 0.0, 0.4, (fV0APlaAl)/2.);
1524     new TGeoTube("sV0ANail2PlaOuHole8", 0.0, 0.4, (fV0APlaAl)/2.);
1525     new TGeoTube("sV0ANail3PlaOuHole8", 0.0, 0.4, (fV0APlaAl)/2.);
1526     new TGeoCompositeShape("sV0ANailsPlaOuHoles8","sV0ANail1PlaOuHole8:pos18+sV0ANail2PlaOuHole8:pos28+sV0ANail3PlaInHole8:pos38");
1527     for (int i=0;i<2;i++) {
1528       v0APts[0+8*i] = fV0AR0;                   v0APts[1+8*i] = 0.;
1529       v0APts[2+8*i] = fV0AR0*cos654;            v0APts[3+8*i] = -fV0AR0*sin654;
1530       v0APts[4+8*i] = fV0AR7*cos654;            v0APts[5+8*i] = -fV0AR7*sin654;
1531       v0APts[6+8*i] = fV0AR6;                   v0APts[7+8*i] = 0.;
1532     }
1533     new TGeoArb8("sV0APlaIn8", (fV0APlaWd-2*fV0APlaAl)/2., v0APts);
1534     TGeoCompositeShape *sV0APlaInNailsHoles8 = new TGeoCompositeShape("sV0APlaInNailsHoles8","sV0APlaIn8-sV0ANailsPlaInHoles8");
1535     TGeoVolume *v0APlaInNailsHoles8 = new TGeoVolume("V0APlaInNailsHoles8", sV0APlaInNailsHoles8, medV0APlaIn);
1536     new TGeoArb8("sV0APlaOu8", fV0APlaAl/2., v0APts);
1537     TGeoCompositeShape *sV0APlaOuNailsHoles8 = new TGeoCompositeShape("sV0APlaOuNailsHoles8","sV0APlaOu8-sV0ANailsPlaOuHoles8"); 
1538     TGeoVolume *v0APlaOuNailsHoles8 = new TGeoVolume("V0APlaOuNailsHoles8", sV0APlaOuNailsHoles8, medV0APlaOu);
1539     v0APlaInNailsHoles8->SetLineColor(kV0AColorPlaIn); v0APlaOuNailsHoles8->SetLineColor(kV0AColorPlaOu);
1540     TGeoVolume *v0APla8 = new TGeoVolumeAssembly("V0APla8");
1541     v0APla8->AddNode(v0APlaInNailsHoles8,1);
1542     v0APla8->AddNode(v0APlaOuNailsHoles8,1,new TGeoTranslation(0,0,(fV0APlaWd-fV0APlaAl)/2.));
1543     v0APla8->AddNode(v0APlaOuNailsHoles8,2,new TGeoTranslation(0,0,-(fV0APlaWd-fV0APlaAl)/2.));
1544     v0ASec8->AddNode(v0APla8,1,new TGeoTranslation(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1545     v0ASec8->AddNode(v0APla8,2,new TGeoTranslation(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1546
1547     /// PMBox 1ero
1548     TGeoVolume* v0APM8 = new TGeoVolumeAssembly("V0APM1");
1549     new TGeoBBox("sV0APMB18", fV0APMBWd/2., fV0APMBHt/2., fV0APMBTh/2.);
1550     new TGeoBBox("sV0APMB28", fV0APMBWd/2.-fV0APMBWdW, fV0APMBHt/2.-fV0APMBHtW, fV0APMBTh/2.-fV0APMBThW);
1551     TGeoCompositeShape *sV0APMB8 = new TGeoCompositeShape("sV0APMB8","sV0APMB18-sV0APMB28");
1552     TGeoVolume *v0APMB8 = new TGeoVolume("V0APMB8",sV0APMB8, medV0APMAlum);
1553     v0APMB8->SetLineColor(kV0AColorPMA);
1554     v0APM8->AddNode(v0APMB8,1);
1555
1556     /// PMTubes 1ero
1557     TGeoTube *sV0APMT18 = new TGeoTube("sV0APMT18", fV0APMTR1, fV0APMTR2, fV0APMTH/2.);
1558     TGeoVolume *v0APMT18 = new TGeoVolume("V0APMT18", sV0APMT18, medV0APMGlass);
1559     TGeoTube *sV0APMT28 = new TGeoTube("sV0APMT28", fV0APMTR3, fV0APMTR4, fV0APMTH/2.);
1560     TGeoVolume *v0APMT28 = new TGeoVolume("V0APMT28", sV0APMT28, medV0APMAlum);
1561     TGeoVolume *v0APMT8 = new TGeoVolumeAssembly("V0APMT8");
1562     TGeoTube *sV0APMTT8 = new TGeoTube("sV0APMTT8", 0., fV0APMTR4, fV0APMTB/2.);
1563     TGeoVolume *v0APMTT8 = new TGeoVolume("V0APMT18", sV0APMTT8, medV0APMAlum);
1564     v0APMT8->SetLineColor(kV0AColorPMG);
1565     v0APMT28->SetLineColor(kV0AColorPMA);
1566     v0APMTT8->SetLineColor(kV0AColorPMA);
1567     rot = new TGeoRotation("rot", 90, 0, 180, 0, 90, 90);
1568     v0APMT8->AddNode(v0APMT18,1,rot);
1569     v0APMT8->AddNode(v0APMT28,1,rot);
1570     v0APMT8->AddNode(v0APMTT8,1,new TGeoCombiTrans(0,(fV0APMTH+fV0APMTB)/2.,0,rot));
1571     double autoShift8 = (fV0APMBWd-2*fV0APMBWdW)/4.;
1572     v0APM8->AddNode(v0APMT8, 1, new TGeoTranslation(-1.5*autoShift8, 0, 0));
1573     v0APM8->AddNode(v0APMT8, 2, new TGeoTranslation(-0.5*autoShift8, 0, 0));
1574     v0APM8->AddNode(v0APMT8, 3, new TGeoTranslation(+0.5*autoShift8, 0, 0));
1575     v0APM8->AddNode(v0APMT8, 4, new TGeoTranslation(+1.5*autoShift8, 0, 0));
1576
1577     /// PM 1ero
1578     rot = new TGeoRotation("rot");
1579     rot->RotateX(-90+30);
1580     rot->RotateY(0);
1581     rot->RotateZ(65+3);
1582     double shiftZ8 = fV0APMBHt/2. * cosAngPMB
1583       -   ( fV0ASciWd + 2 * fV0AOctWd + 2 * fV0APlaWd )/2.   -   fV0APMBTh/2. * sinAngPMB;
1584     double shiftR8 = fV0AR6  +  fV0AOctH2 + fV0APlaAl;
1585     v0ASec8->AddNode(v0APM8,1, new TGeoCombiTrans( shiftR8*cos225, -shiftR8*sin225, shiftZ8, rot));
1586
1587     /// PMBox 2do 
1588     TGeoVolume* v0APM82 = new TGeoVolumeAssembly("V0APM82");
1589     new TGeoBBox("sV0APMB82", fV0APMBWd/2., fV0APMBHt/2., fV0APMBTh/2.);
1590     new TGeoBBox("sV0APMB82", fV0APMBWd/2.-fV0APMBWdW, fV0APMBHt/2.-fV0APMBHtW, fV0APMBTh/2.-fV0APMBThW);
1591     TGeoCompositeShape *sV0APMB82 = new TGeoCompositeShape("sV0APMB82","sV0APMB82-sV0APMB82");
1592     TGeoVolume *v0APMB82 = new TGeoVolume("V0APMB82",sV0APMB82, medV0APMAlum);
1593     v0APMB82->SetLineColor(kV0AColorPMA);
1594     v0APM82->AddNode(v0APMB82,1);
1595
1596     /// PMTubes 2ndo
1597     TGeoTube *sV0APMT182 = new TGeoTube("sV0APMT182", fV0APMTR1, fV0APMTR2, fV0APMTH/2.);
1598     TGeoVolume *v0APMT182 = new TGeoVolume("V0APMT182", sV0APMT182, medV0APMGlass);
1599     TGeoTube *sV0APMT282 = new TGeoTube("sV0APMT282", fV0APMTR3, fV0APMTR4, fV0APMTH/2.);
1600     TGeoVolume *v0APMT282 = new TGeoVolume("V0APMT282", sV0APMT282, medV0APMAlum);
1601     TGeoVolume *v0APMT82 = new TGeoVolumeAssembly("V0APMT82"); // pk si no choca con la 752 o con la 794
1602     TGeoTube *sV0APMTT82 = new TGeoTube("sV0APMTT82", 0., fV0APMTR4, fV0APMTB/2.);
1603     TGeoVolume *v0APMTT82 = new TGeoVolume("V0APMT82", sV0APMTT82, medV0APMAlum);
1604     v0APMT82->SetLineColor(kV0AColorPMG);
1605     v0APMT282->SetLineColor(kV0AColorPMA);
1606     v0APMTT82->SetLineColor(kV0AColorPMA);
1607     rot = new TGeoRotation("rot", 90, 0, 180, 0, 90, 90);
1608     v0APMT82->AddNode(v0APMT182,1,rot);
1609     v0APMT82->AddNode(v0APMT282,1,rot);
1610     v0APMT82->AddNode(v0APMTT82,1,new TGeoCombiTrans(0,(fV0APMTH+fV0APMTB)/2.,0,rot));
1611     v0APM82->AddNode(v0APMT82, 1, new TGeoTranslation(-1.5*autoShift, 0, 0));
1612     v0APM82->AddNode(v0APMT82, 2, new TGeoTranslation(-0.5*autoShift, 0, 0));
1613     v0APM82->AddNode(v0APMT82, 3, new TGeoTranslation(+0.5*autoShift, 0, 0));
1614     v0APM82->AddNode(v0APMT82, 4, new TGeoTranslation(+1.5*autoShift, 0, 0));
1615
1616     /// PM 2ndo
1617     rot = new TGeoRotation("rot");
1618     rot->RotateX(-90+30);
1619     rot->RotateY(0);
1620     rot->RotateZ(65+3);
1621     double shiftZ82 = fV0APMBHt/2. * cosAngPMB
1622       -   ( fV0ASciWd + 2 * fV0AOctWd + 2 * fV0APlaWd )/2.   -   fV0APMBTh/2. * sinAngPMB;
1623     double shiftR82 = fV0AR6  + fV0AR1 + fV0AOctWd + fV0APlaAl/3.;
1624     v0ASec8->AddNode(v0APM82,1, new TGeoCombiTrans( shiftR82*cos45, -shiftR82*sin45, shiftZ82, rot)); 
1625    
1626     // Replicate sectors
1627     for(int i=0; i<1; i++) {
1628     TGeoRotation *rot = new TGeoRotation("rot", 90., i*45., 90., 90.+i*45., 0., 0.);    
1629     v0LE->AddNode(v0ASec8,i + 1,rot);
1630     }
1631     
1632     
1633     ///Aluminium nails for sectors 5, 6, 7 and 8
1634     TGeoTube *sV0ANail51 = new TGeoTube("sV0ANail51", 0.0, 0.4, 5.09/2.);
1635     TGeoVolume *v0ANail51 = new TGeoVolume("V0ANail51", sV0ANail51, medV0APMAlum);
1636     v0ANail51->SetLineColor(kV0AColorPMA);// this is the color for aluminium
1637     v0LE->AddNode(v0ANail51,1,new TGeoTranslation(-42.9,-0.51,0.0));
1638     TGeoTube *sV0ANail52 = new TGeoTube("sV0ANail52", 0.0, 0.4, 5.09/2.);
1639     TGeoVolume *v0ANail52 = new TGeoVolume("V0ANail52", sV0ANail52, medV0APMAlum);
1640     v0ANail52->SetLineColor(kV0AColorPMA);
1641     v0LE->AddNode(v0ANail52,1,new TGeoTranslation(-30.8,-30.04,0.0)); 
1642     TGeoTube *sV0ANail61 = new TGeoTube("sV0ANail61", 0.0, 0.4, 5.09/2.);
1643     TGeoVolume *v0ANail61 = new TGeoVolume("V0ANail61", sV0ANail61, medV0APMAlum);
1644     v0ANail61->SetLineColor(kV0AColorPMA);// this is the color for aluminium
1645     v0LE->AddNode(v0ANail61,1,new TGeoTranslation(-0.51,-42.9,0.0));  
1646     TGeoTube *sV0ANail53 = new TGeoTube("sV0ANail53", 0.0, 0.4, 5.09/2.);
1647     TGeoVolume *v0ANail53 = new TGeoVolume("V0ANail53", sV0ANail53, medV0APMAlum);
1648     v0ANail53->SetLineColor(kV0AColorPMA);
1649     v0LE->AddNode(v0ANail53,1,new TGeoTranslation(-30.05,-30.79,0.0)); //     -29.8,-31.04,0.0
1650     TGeoTube *sV0ANail71 = new TGeoTube("sV0ANail71", 0.0, 0.4, 5.09/2.);
1651     TGeoVolume *v0ANail71 = new TGeoVolume("V0ANail71", sV0ANail71, medV0APMAlum);
1652     v0ANail71->SetLineColor(kV0AColorPMA);// this is the color for aluminium
1653     v0LE->AddNode(v0ANail71,1,new TGeoTranslation(0.51,-42.9,0.0));
1654     TGeoTube *sV0ANail83 = new TGeoTube("sV0ANail83", 0.0, 0.4, 5.09/2.);
1655     TGeoVolume *v0ANail83 = new TGeoVolume("V0ANail83", sV0ANail83, medV0APMAlum);
1656     v0ANail83->SetLineColor(kV0AColorPMA);
1657     v0LE->AddNode(v0ANail83,1,new TGeoTranslation(29.8,-31.04,0.0)); 
1658     TGeoTube *sV0ANail81 = new TGeoTube("sV0ANail81", 0.0, 0.4, 5.09/2.);
1659     TGeoVolume *v0ANail81 = new TGeoVolume("V0ANail81", sV0ANail81, medV0APMAlum);
1660     v0ANail81->SetLineColor(kV0AColorPMA);// this is the color for aluminium
1661     v0LE->AddNode(v0ANail81,1,new TGeoTranslation(42.9,-.51,0.0));
1662     TGeoTube *sV0ANail82 = new TGeoTube("sV0ANail82", 0.0, 0.4, 5.09/2.);
1663     TGeoVolume *v0ANail82 = new TGeoVolume("V0ANail82", sV0ANail82, medV0APMAlum); 
1664     v0ANail82->SetLineColor(kV0AColorPMA);
1665     v0LE->AddNode(v0ANail82,1,new TGeoTranslation(30.8,-30.04,0.0)); 
1666
1667
1668
1669    /* /// Basis Construction
1670     rot = new TGeoRotation("rot"); rot->RotateX(90-fV0APMBAng); rot->RotateZ(-22.5);
1671     TGeoCombiTrans *pos1 = new TGeoCombiTrans("pos1", shiftR*sin225, shiftR*cos225, shiftZ, rot);
1672     pos1->RegisterYourself();
1673     for (int i=0;i<2;i++) {
1674     v0APts[0+8*i] = fV0AR6/cos225*sin45;  v0APts[1+8*i] = fV0AR6/cos225*sin45;
1675     v0APts[2+8*i] = 0;              v0APts[3+8*i] = fV0AR6/cos225;
1676     v0APts[4+8*i] = 0;              v0APts[5+8*i] = fV0AR6/cos225+fV0APlaEx;
1677     v0APts[6+8*i] = fV0AR6/cos225-(fV0AR6/cos225+fV0APlaEx)/ctg225;
1678     v0APts[7+8*i] = fV0AR6/cos225+fV0APlaEx;
1679     }
1680     new TGeoArb8("sV0APlaExIn1", (fV0APlaWd-2*fV0APlaAl)/2., v0APts);
1681     new TGeoArb8("sV0APlaExOu1", fV0APlaAl/2., v0APts);
1682     TGeoCompositeShape *sV0APlaExIn = new TGeoCompositeShape("sV0APlaExIn","sV0APlaExIn1-sV0APMB1:pos1");
1683     TGeoVolume *v0APlaExIn = new TGeoVolume("V0APlaExIn", sV0APlaExIn, medV0APlaIn);
1684     TGeoCompositeShape *sV0APlaExOu = new TGeoCompositeShape("sV0APlaExOu","sV0APlaExOu1-sV0APMB1:pos1");
1685     TGeoVolume *v0APlaExOu = new TGeoVolume("V0APlaExOu", sV0APlaExOu, medV0APlaOu);
1686     v0APlaExIn->SetLineColor(kV0AColorPlaIn); v0APlaExOu->SetLineColor(kV0AColorPlaOu);
1687     TGeoVolume *v0APlaEx = new TGeoVolumeAssembly("V0APlaEx");
1688     v0APlaEx->AddNode(v0APlaExIn,1);
1689     v0APlaEx->AddNode(v0APlaExOu,1,new TGeoTranslation(0,0,(fV0APlaWd-fV0APlaAl)/2.));
1690     v0APlaEx->AddNode(v0APlaExOu,2,new TGeoTranslation(0,0,-(fV0APlaWd-fV0APlaAl)/2.));
1691     for (int i=0;i<2;i++) {
1692     v0APts[0+8*i] = fV0AR6/cos225-(fV0AR6/cos225+fV0APlaEx)/ctg225-fV0ABasHt*sin45;
1693     v0APts[1+8*i] = fV0AR6/cos225+fV0APlaEx-fV0ABasHt*sin45;
1694     v0APts[2+8*i] = 0;  v0APts[3+8*i] = fV0AR6/cos225+fV0APlaEx-fV0ABasHt;
1695     v0APts[4+8*i] = 0;  v0APts[5+8*i] = fV0AR6/cos225+fV0APlaEx;
1696     v0APts[6+8*i] = fV0AR6/cos225-(fV0AR6/cos225+fV0APlaEx)/ctg225;
1697     v0APts[7+8*i] = fV0AR6/cos225+fV0APlaEx;
1698     }
1699     new TGeoArb8("sV0ABas1", (fV0ASciWd+2*fV0AOctWd)/2., v0APts);
1700     TGeoCompositeShape *sV0ABas = new TGeoCompositeShape("sV0ABas","sV0ABas1-sV0APMB1:pos1");
1701     TGeoVolume *v0ABas = new TGeoVolume("V0ABas", sV0ABas, medV0ABas);
1702     v0ABas->SetLineColor(kV0AColorBas);
1703     TGeoVolume *v0ABasis = new TGeoVolumeAssembly("V0ABasis");
1704     rot = new TGeoRotation("rot",90.,180.,90.,90.,0.,0.);
1705     v0ABasis->AddNode(v0APlaEx,1, new TGeoTranslation(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1706     v0ABasis->AddNode(v0APlaEx,2, new TGeoTranslation(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.));
1707     v0ABasis->AddNode(v0APlaEx,3, new TGeoCombiTrans(0,0,(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.,rot));
1708     v0ABasis->AddNode(v0APlaEx,4, new TGeoCombiTrans(0,0,-(fV0ASciWd+2*fV0AOctWd+fV0APlaWd)/2.,rot));
1709     v0ABasis->AddNode(v0ABas,1);
1710     v0ABasis->AddNode(v0ABas,2,rot);
1711     rot = new TGeoRotation("rot");
1712     rot->RotateZ(180);
1713     v0LE->AddNode(v0ABasis,1,rot); */
1714
1715     // Adding detectors to top volume
1716     TGeoVolume *vZERO = new TGeoVolumeAssembly("VZERO");
1717     vZERO->AddNode(v0RI,1,new TGeoTranslation(0, 0, -zdet));
1718     // V0A position according to TB decision 13/12/2005 
1719     //rot=new TGeoRotation("rot");
1720     //rot->RotateX(90);
1721     //rot->RotateY(180);
1722     //rot->RotateZ(270);
1723     vZERO->AddNode(v0LE,1,new TGeoTranslation(0, 0, +327.5));//,rot));
1724     top->AddNode(vZERO,1);
1725 }
1726
1727 //_____________________________________________________________________________
1728 void AliVZEROv7::AddAlignableVolumes() const
1729 {
1730   //
1731   // Create entries for alignable volumes associating the symbolic volume
1732   // name with the corresponding volume path. Needs to be syncronized with
1733   // eventual changes in the geometry.
1734   // 
1735   TString vpC = "/ALIC_1/VZERO_1/V0RI_1";
1736   TString vpA = "/ALIC_1/VZERO_1/V0LE_1";
1737   TString snC = "VZERO/V0C";
1738   TString snA = "VZERO/V0A";
1739   
1740   if(!gGeoManager->SetAlignableEntry(snC.Data(),vpC.Data()))
1741     AliFatal(Form("Alignable entry %s not created. Volume path %s not valid", snC.Data(),vpC.Data()));
1742   if(!gGeoManager->SetAlignableEntry(snA.Data(),vpA.Data()))
1743     AliFatal(Form("Alignable entry %s not created. Volume path %s not valid", snA.Data(),vpA.Data()));
1744
1745
1746
1747 //_____________________________________________________________________________
1748 void AliVZEROv7::CreateMaterials()
1749 {
1750
1751 // Creates materials used for geometry 
1752
1753   AliDebug(2,"Create materials");
1754   // Parameters for simulation scope
1755   Int_t     fieldType       = gAlice->Field()->Integ();     // Field type 
1756   Double_t  maxField        = gAlice->Field()->Max();       // Field max.
1757   Double_t  maxBending      = 10;    // Max Angle
1758   Double_t  maxStepSize     = 0.01;  // Max step size 
1759   Double_t  maxEnergyLoss   = 1;     // Max Delta E
1760   Double_t  precision       = 0.003; // Precision
1761   Double_t  minStepSize     = 0.003; // Minimum step size 
1762
1763   Int_t    id;
1764   Double_t a, z, radLength, absLength;
1765   Float_t density, as[4], zs[4], ws[4];
1766
1767 // Parameters  for V0CPrePlates: Aluminium
1768    a = 26.98; 
1769    z = 13.00;
1770    density     = 2.7;
1771    radLength   = 8.9;
1772    absLength   = 37.2;
1773    id = 2;
1774    AliMaterial( id, "V0CAlu", a, z, density, radLength, absLength, 0, 0);
1775    AliMedium(id, "V0CAlu", id, 1, fieldType, maxField, maxBending, maxStepSize,
1776              maxEnergyLoss, precision, minStepSize);
1777                     
1778 // Parameters  for V0CPlates: Carbon 
1779    a = 12.01; 
1780    z =  6.00;
1781    density   = 2.265;
1782    radLength = 18.8;
1783    absLength = 49.9;
1784    id = 3;
1785    AliMaterial(id, "V0CCar",  a, z, density, radLength, absLength, 0, 0);
1786    AliMedium(id, "V0CCar", id, 1, fieldType, maxField, maxBending, maxStepSize,
1787              maxEnergyLoss, precision, minStepSize);
1788             
1789 // Parameters  for V0Cscintillator: BC408
1790    as[0] = 1.00794;     as[1] = 12.011;
1791    zs[0] = 1.;          zs[1] = 6.;
1792    ws[0] = 1.;          ws[1] = 1.;
1793    density      = 1.032;
1794    id           = 4;
1795    AliMixture(id, "V0CSci", as, zs, density, -2, ws);
1796    AliMedium(id,"V0CSci", id, 1, fieldType, maxField, maxBending, maxStepSize,
1797              maxEnergyLoss, precision, minStepSize);
1798
1799 // Parameters for V0Ascintilator: BC404
1800    as[0] = 1.00794;     as[1] = 12.011;
1801    zs[0] = 1.;          zs[1] = 6.;
1802    ws[0] = 5.21;        ws[1] = 4.74;
1803    density      = 1.032;
1804    id           = 5;
1805    AliMixture(id, "V0ASci", as, zs, density, -2, ws);
1806    AliMedium(id,  "V0ASci", id, 1, fieldType, maxField, maxBending, maxStepSize,
1807              maxEnergyLoss, precision, minStepSize);
1808
1809 // Parameters for V0ALuc: Lucita but for the simulation BC404
1810    as[0] = 1.00794;     as[1] = 12.011;
1811    zs[0] = 1.;          zs[1] = 6.;
1812    ws[0] = 5.21;        ws[1] = 4.74;
1813    density      = 1.032;
1814    id           = 6;
1815    AliMixture(id, "V0ALuc", as, zs, density, -2, ws);
1816    AliMedium(id, "V0ALuc", id, 1, fieldType, maxField, maxBending, maxStepSize,
1817              maxEnergyLoss, precision, minStepSize);
1818
1819 // Parameters for V0Aplate: EuroComposite - EC-PI 626 PS - AlMg3
1820    as[0] = 26.982;      as[1] = 24.305;
1821    zs[0] = 13.;         zs[1] = 12.;
1822    ws[0] = 1.;          ws[1] = 3.;
1823    density      = 3.034;
1824    id           = 7;
1825    AliMixture(id, "V0APlaOu", as, zs, density, -2, ws);
1826    AliMedium(id, "V0APlaOu", id, 1, fieldType, maxField, maxBending, maxStepSize,
1827              maxEnergyLoss, precision, minStepSize);
1828
1829 // Parameters for V0Aplate: EuroComposite - EC-PI 626 PS - EC-PI 6.4-42
1830    as[0] = 1.00794;     as[1] = 12.011;
1831    zs[0] = 1.;          zs[1] = 6.;
1832    ws[0] = 5.21;        ws[1] = 4.74;
1833    density      = 0.042;
1834    id           = 8;
1835    AliMixture(id, "V0APlaIn", as, zs, density, -2, ws);
1836    AliMedium(id, "V0APlaIn", id, 1, fieldType, maxField, maxBending, maxStepSize,
1837              maxEnergyLoss, precision, minStepSize);
1838
1839 // Parameters for V0Afiber: BC9929AMC Plastic Scintillating Fiber from Saint-Gobain
1840    as[0] = 1.00794;     as[1] = 12.011;
1841    zs[0] = 1.;          zs[1] = 6.;
1842    ws[0] = 4.82;        ws[1] = 4.85;
1843    density      = 1.05;
1844    id           = 9;
1845    AliMixture(id, "V0AFib", as, zs, density, -2, ws);
1846    AliMedium(id, "V0AFib", id, 1, fieldType, maxField, maxBending, maxStepSize,
1847              maxEnergyLoss, precision, minStepSize);
1848
1849 // Parameters for V0APMA: Aluminium
1850    a = 26.98; 
1851    z = 13.00;
1852    density     = 2.7;
1853    radLength   = 8.9;
1854    absLength   = 37.2;
1855    id = 10;
1856    AliMaterial(id, "V0APMA",  a, z, density, radLength, absLength, 0, 0);
1857    AliMedium(id, "V0APMA", id, 1, fieldType, maxField, maxBending, maxStepSize,
1858              maxEnergyLoss, precision, minStepSize);
1859
1860 // Parameters for V0APMG: Glass for the simulation Aluminium
1861    a = 26.98; 
1862    z = 13.00;
1863    density   = 2.7;
1864    radLength = 8.9;
1865    absLength = 37.2;
1866    id = 11;
1867    AliMaterial(id, "V0APMG",  a, z, density, radLength, absLength, 0, 0);
1868    AliMedium(id, "V0APMG", id, 1, fieldType, maxField, maxBending, maxStepSize,
1869              maxEnergyLoss, precision, minStepSize);
1870 }
1871
1872 //_____________________________________________________________________________
1873 void AliVZEROv7::DrawModule() const
1874 {
1875 //  Drawing is done in DrawVZERO.C
1876
1877    AliDebug(2,"DrawModule");
1878 }
1879
1880
1881 //_____________________________________________________________________________
1882 void AliVZEROv7::DrawGeometry() 
1883 {
1884 //  Drawing of V0 geometry done in DrawV0.C
1885
1886    AliDebug(2,"DrawGeometry");
1887 }
1888
1889 //_____________________________________________________________________________
1890 void AliVZEROv7::Init()
1891 {
1892 // Initialises version of the VZERO Detector given in Config
1893 // Just prints an information message
1894
1895 //   AliInfo(Form("VZERO version %d initialized \n",IsVersion()));
1896    
1897    AliDebug(1,"VZERO version 7 initialized");
1898    AliVZERO::Init();  
1899 }
1900
1901 //_____________________________________________________________________________
1902 void AliVZEROv7::StepManager()
1903 {
1904   // Step Manager, called at each step 
1905
1906   Int_t     copy;
1907   static    Int_t   vol[4];
1908   static    Float_t hits[21];
1909   static    Float_t eloss, tlength;
1910   static    Int_t   nPhotonsInStep = 0;
1911   static    Int_t   nPhotons = 0; 
1912   static    Int_t   numStep = 0;
1913   Int_t     ringNumber;
1914   Float_t   destep, step;
1915   numStep += 1; 
1916
1917   //   We keep only charged tracks :
1918   if ( !gMC->TrackCharge() || !gMC->IsTrackAlive() ) return;
1919
1920   vol[0]    = gMC->CurrentVolOffID(1, vol[1]);
1921   vol[2]    = gMC->CurrentVolID(copy);
1922   vol[3]    = copy;
1923   static Int_t idV0R1 = gMC->VolId("V0R1");
1924   static Int_t idV0L1 = gMC->VolId("V0L1");
1925   static Int_t idV0R2 = gMC->VolId("V0R2");
1926   static Int_t idV0L2 = gMC->VolId("V0L2");
1927   static Int_t idV0R3 = gMC->VolId("V0R3");
1928   static Int_t idV0L3 = gMC->VolId("V0L3");
1929   static Int_t idV0R4 = gMC->VolId("V0R4");
1930   static Int_t idV0L4 = gMC->VolId("V0L4");
1931   static Int_t idV0R5 = gMC->VolId("V0R5");
1932   static Int_t idV0R6 = gMC->VolId("V0R6");
1933   bool   hitOnV0C = true;
1934   double lightYield;
1935   double lightAttenuation;
1936   double nMeters;
1937   double fibToPhot;
1938   if      ( gMC->CurrentVolID(copy) == idV0R1 || gMC->CurrentVolID(copy) == idV0L1 )
1939     ringNumber = 1;
1940   else if ( gMC->CurrentVolID(copy) == idV0R2 || gMC->CurrentVolID(copy) == idV0L2 )
1941     ringNumber = 2;  
1942   else if ( gMC->CurrentVolID(copy) == idV0R3 || gMC->CurrentVolID(copy) == idV0R4
1943             || gMC->CurrentVolID(copy) == idV0L3 ) ringNumber = 3;
1944   else if ( gMC->CurrentVolID(copy) == idV0R5 || gMC->CurrentVolID(copy) == idV0R6
1945             || gMC->CurrentVolID(copy) == idV0L4 ) ringNumber = 4;             
1946   else ringNumber = 0;
1947   if  (ringNumber) {
1948     if (gMC->CurrentVolID(copy) == idV0L1 || gMC->CurrentVolID(copy) == idV0L2 ||
1949         gMC->CurrentVolID(copy) == idV0L3 || gMC->CurrentVolID(copy) == idV0L4)
1950       hitOnV0C = false;
1951     destep = gMC->Edep();
1952     step   = gMC->TrackStep();
1953     if (hitOnV0C) {
1954       lightYield = fV0CLightYield;
1955       lightAttenuation = fV0CLightAttenuation;
1956       nMeters = fV0CnMeters;
1957       fibToPhot = fV0CFibToPhot;
1958     } else {
1959       lightYield = fV0ALightYield;
1960       lightAttenuation = fV0ALightAttenuation;
1961       nMeters = fV0AnMeters;
1962       fibToPhot = fV0AFibToPhot;
1963     }
1964     nPhotonsInStep  = Int_t(destep / (lightYield *1e-9) );      
1965     nPhotonsInStep  = gRandom->Poisson(nPhotonsInStep);
1966     eloss    += destep;
1967     tlength  += step;    
1968     if ( gMC->IsTrackEntering() ) { 
1969       nPhotons  =  nPhotonsInStep;
1970       gMC->TrackPosition(fTrackPosition);
1971       gMC->TrackMomentum(fTrackMomentum);
1972       Float_t pt  = TMath::Sqrt( fTrackMomentum.Px() * fTrackMomentum.Px()
1973                                  + fTrackMomentum.Py() * fTrackMomentum.Py() );
1974       TParticle *par = gAlice->GetMCApp()->Particle(gAlice->GetMCApp()->GetCurrentTrackNumber());
1975       hits[0]  = fTrackPosition.X();
1976       hits[1]  = fTrackPosition.Y();
1977       hits[2]  = fTrackPosition.Z();             
1978       hits[3]  = Float_t (gMC->TrackPid()); 
1979       hits[4]  = gMC->TrackTime();
1980       hits[5]  = gMC->TrackCharge();
1981       hits[6]  = fTrackMomentum.Theta()*TMath::RadToDeg();
1982       hits[7]  = fTrackMomentum.Phi()*TMath::RadToDeg();
1983       hits[8]  = ringNumber;
1984       hits[9]  = pt;
1985       hits[10] = fTrackMomentum.P();
1986       hits[11] = fTrackMomentum.Px();
1987       hits[12] = fTrackMomentum.Py();
1988       hits[13] = fTrackMomentum.Pz();
1989       hits[14] = par->Vx();
1990       hits[15] = par->Vy();
1991       hits[16] = par->Vz();
1992       tlength  = 0.0;
1993       eloss    = 0.0;       
1994
1995       //////////////////////////
1996       ///// Display V0A geometry
1997       //      if (!hitOnV0C) {
1998       //        FILE *of;
1999       //        of = fopen("V0A.out", "a");
2000       //        // x, y, z, ringnumber, cellid
2001       //        fprintf( of, "%f %f %f %f %d \n",  hits[0], hits[1], hits[2], hits[8], GetCellId (vol, hits) );
2002       //        fclose(of);
2003       //      }
2004       //////////////////////////
2005     }
2006     nPhotons  = nPhotons + nPhotonsInStep;
2007     if( gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()){
2008       nPhotons = nPhotons - Int_t((Float_t(nPhotons) * lightAttenuation * nMeters));
2009       nPhotons = nPhotons - Int_t( Float_t(nPhotons) * fibToPhot);
2010       hits[17] = eloss;
2011       hits[18] = tlength;
2012       hits[19] = nPhotons;
2013       hits[20] = GetCellId (vol, hits);
2014       AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber(), vol, hits);
2015       tlength         = 0.0;
2016       eloss           = 0.0; 
2017       nPhotons        = 0;
2018       nPhotonsInStep  = 0;
2019       numStep         = 0;  
2020     }
2021   }
2022 }
2023
2024 //_____________________________________________________________________________
2025 void AliVZEROv7::AddHit(Int_t track, Int_t *vol, Float_t *hits)
2026 {
2027 //  Adds a VZERO hit
2028
2029   TClonesArray &lhits = *fHits;
2030   new(lhits[fNhits++]) AliVZEROhit(fIshunt,track,vol,hits);
2031 }
2032
2033 //_____________________________________________________________________________
2034 void AliVZEROv7::AddDigits(Int_t *tracks, Int_t* digits) 
2035 {
2036 //  Adds a VZERO digit
2037
2038    TClonesArray  &ldigits = *fDigits;
2039    new(ldigits[fNdigits++]) AliVZEROdigit(tracks, digits);
2040 }
2041
2042 //_____________________________________________________________________________
2043 void AliVZEROv7::MakeBranch(Option_t *option)
2044 {
2045 // Creates new branches in the current Root Tree
2046     
2047   char branchname[10];
2048   sprintf(branchname,"%s",GetName());
2049   AliDebug(2,Form("fBufferSize = %d",fBufferSize));
2050   const char *cH = strstr(option,"H");
2051   if (fHits   && TreeH() && cH) {
2052     TreeH()->Branch(branchname,&fHits, fBufferSize);
2053     AliDebug(2,Form("Making Branch %s for hits",branchname));
2054   }     
2055   const char *cD = strstr(option,"D");
2056   if (fDigits   && fLoader->TreeD() && cD) {
2057     fLoader->TreeD()->Branch(branchname,&fDigits, fBufferSize);
2058     AliDebug(2,Form("Making Branch %s for digits",branchname));
2059   }  
2060 }
2061
2062 //_____________________________________________________________________________
2063 Int_t AliVZEROv7::GetCellId(Int_t *vol, Float_t *hits) 
2064 {
2065   //   Returns Id of scintillator cell
2066   //   Right side from  0 to 47 
2067   //   Left  side from 48 to 79
2068   //   hits[8] = ring number (1 to 4)
2069   //   vol[1]  = copy number (1 to 8)
2070
2071   Int_t index      = vol[1];
2072   Int_t ringNumber = Int_t(hits[8]);
2073   fCellId          = 0;
2074
2075   Float_t phi = Float_t(TMath::ATan2(Double_t(hits[1]),Double_t(hits[0])) ); 
2076   Float_t kRaddeg = 180.0/TMath::Pi();
2077   phi = kRaddeg * phi;
2078
2079   if (index < 7) index = index + 8;
2080
2081   if (hits[2] < 0.0) {
2082     if(ringNumber < 3) {
2083       index = (index - 7) + ( ( ringNumber - 1 ) * 8);
2084     } else if (ringNumber >= 3) { 
2085       if ( gMC->CurrentVolID(vol[1]) == gMC->VolId("V0R3") || gMC->CurrentVolID(vol[1])
2086            == gMC->VolId("V0R5") )  index = (index*2-14)+((ringNumber-2)*16);
2087       if ( gMC->CurrentVolID(vol[1]) == gMC->VolId("V0R4") || gMC->CurrentVolID(vol[1])
2088            == gMC->VolId("V0R6") )  index = (index*2-13)+((ringNumber-2)*16);
2089     }
2090     fCellId   = index;           
2091   } else if (hits[2] > 0.0) {
2092     index = (index - 7 + 48) + ( ( ringNumber - 1 ) * 8);
2093     fCellId   = index;
2094   }
2095
2096   return fCellId;
2097 }