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