source: trunk/examples/advanced/hadrontherapy/src/HadrontherapyBeamLine.cc@ 1318

Last change on this file since 1318 was 807, checked in by garnier, 17 years ago

update

  • Property svn:executable set to *
File size: 33.2 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: HadrontherapyBeamLine.cc; May 2005
27// ----------------------------------------------------------------------------
28// GEANT 4 - Hadrontherapy example
29// ----------------------------------------------------------------------------
30// Code developed by:
31//
32// G.A.P. Cirrone(a)*, F. Di Rosa(a), S. Guatelli(b), G. Russo(a)
33//
34// (a) Laboratori Nazionali del Sud
35// of the INFN, Catania, Italy
36// (b) INFN Section of Genova, Genova, Italy
37//
38// * cirrone@lns.infn.it
39// ----------------------------------------------------------------------------
40
41#include "G4Box.hh"
42#include "G4Tubs.hh"
43#include "G4VisAttributes.hh"
44#include "G4Colour.hh"
45#include "HadrontherapyMaterial.hh"
46#include "globals.hh"
47#include "G4LogicalVolume.hh"
48#include "G4PVPlacement.hh"
49#include "G4RotationMatrix.hh"
50#include "HadrontherapyBeamLine.hh"
51#include "G4Material.hh"
52#include "G4SubtractionSolid.hh"
53
54HadrontherapyBeamLine::HadrontherapyBeamLine(G4VPhysicalVolume* motherVolume):
55 physiBeamLineSupport(0),
56 firstScatteringFoil(0), physiFirstScatteringFoil(0), physiKaptonWindow(0),
57 solidStopper(0), physiStopper(0),
58 secondScatteringFoil(0), physiSecondScatteringFoil(0),
59 physiFirstCollimator(0), solidRangeShifterBox(0), logicRangeShifterBox(0),
60 physiRangeShifterBox(0), physiSecondCollimator(0),
61 physiFirstCollimatorModulatorBox(0),
62 physiHoleFirstCollimatorModulatorBox(0),
63 physiSecondCollimatorModulatorBox(0),
64 physiHoleSecondCollimatorModulatorBox(0),
65 physiFirstMonitorLayer1(0), physiFirstMonitorLayer2(0),
66 physiFirstMonitorLayer3(0), physiFirstMonitorLayer4(0),
67 physiSecondMonitorLayer1(0), physiSecondMonitorLayer2(0),
68 physiSecondMonitorLayer3(0), physiSecondMonitorLayer4(0),
69 physiThirdMonitorLayer1(0), physiThirdMonitorLayer2(0),
70 physiThirdMonitorLayer3(0), physiThirdMonitorLayer4(0),
71 physiNozzleSupport(0), physiHoleNozzle(0),
72 solidFinalCollimator(0),
73 physiFinalCollimator(0)
74{
75 mother = motherVolume;
76 material = new HadrontherapyMaterial();
77
78 G4double defaultFirstScatteringFoilXSize = 0.0075 *mm;
79 firstScatteringFoilXSize = defaultFirstScatteringFoilXSize;
80
81 G4double defaultOuterRadiusStopper = 2 *mm;
82 outerRadiusStopper = defaultOuterRadiusStopper;
83
84 G4double defaultSecondScatteringFoilXSize = 0.0125 *mm;
85 secondScatteringFoilXSize = defaultSecondScatteringFoilXSize;
86
87 G4double defaultRangeShifterXSize = 5. *mm;
88 rangeShifterXSize = defaultRangeShifterXSize;
89
90 G4double defaultRangeShifterXPosition = -2530.5 *mm;
91 rangeShifterXPosition = defaultRangeShifterXPosition;
92
93 G4double defaultinnerRadiusFinalCollimator = 12.5 *mm;
94 innerRadiusFinalCollimator = defaultinnerRadiusFinalCollimator;
95
96 G4Material* air = material -> GetMat("Air") ;
97 RSMat = air;
98}
99
100HadrontherapyBeamLine::~HadrontherapyBeamLine()
101{
102 delete material;
103}
104
105void HadrontherapyBeamLine::HadrontherapyBeamLineSupport()
106{
107 // ------------------//
108 // Beam line support //
109 //-------------------//
110
111 const G4double beamLineSupportXSize = 1.5*m;
112 const G4double beamLineSupportYSize = 20.*mm;
113 const G4double beamLineSupportZSize = 600.*mm;
114
115 const G4double beamLineSupportXPosition = -1948.59 *mm;
116 const G4double beamLineSupportYPosition = -230. *mm;
117 const G4double beamLineSupportZPosition = 0.*mm;
118
119 G4Material* Al = material -> GetMat("MatAluminum");
120
121 G4Box* beamLineSupport = new G4Box("BeamLineSupport",
122 beamLineSupportXSize,
123 beamLineSupportYSize,
124 beamLineSupportZSize);
125
126 G4LogicalVolume* logicBeamLineSupport = new G4LogicalVolume(beamLineSupport,
127 Al,
128 "BeamLineSupport");
129 physiBeamLineSupport = new G4PVPlacement(0, G4ThreeVector(beamLineSupportXPosition,
130 beamLineSupportYPosition,
131 beamLineSupportZPosition),
132 "BeamLineSupport",
133 logicBeamLineSupport,
134 mother, false, 0);
135
136 // Visualisation attributes of the beam line support
137 G4VisAttributes * gray = new G4VisAttributes( G4Colour(0.5, 0.5, 0.5 ));
138 gray-> SetVisibility(true);
139 gray-> SetForceSolid(true);
140 logicBeamLineSupport -> SetVisAttributes(gray);
141
142
143}
144
145void HadrontherapyBeamLine::HadrontherapyBeamScatteringFoils()
146{
147 // ------------//
148 // Vacuum area //
149 //-------------//
150
151 const G4double vacuumZoneXSize = 60.5325 *mm;
152 const G4double vacuumZoneYSize = 52.5 *mm;
153 const G4double vacuumZoneZSize = 52.5 *mm;
154
155 const G4double vacuumZoneXPosition = -3188.05750 *mm;
156
157 G4Material* vacuum = material -> GetMat("Galactic");
158
159 G4Box* vacuumZone = new G4Box("VacuumZone", vacuumZoneXSize, vacuumZoneYSize, vacuumZoneZSize);
160
161 G4LogicalVolume* logicVacuumZone = new G4LogicalVolume(vacuumZone, vacuum, "VacuumZone");
162
163 G4VPhysicalVolume* physiVacuumZone = new G4PVPlacement(0, G4ThreeVector(vacuumZoneXPosition, 0., 0.),
164 "VacuumZone",
165 logicVacuumZone,
166 mother,
167 false,
168 0);
169 // ----------------------//
170 // First scattering foil //
171 // ----------------------//
172
173 const G4double firstScatteringFoilYSize = 52.5 *mm;
174 const G4double firstScatteringFoilZSize = 52.5 *mm;
175
176 const G4double firstScatteringFoilXPosition = -59.525 *mm;
177
178 G4Material* Ta = material -> GetMat("MatTantalum");
179
180 firstScatteringFoil = new G4Box("FirstScatteringFoil",
181 firstScatteringFoilXSize,
182 firstScatteringFoilYSize,
183 firstScatteringFoilZSize);
184
185 G4LogicalVolume* logicFirstScatteringFoil = new G4LogicalVolume(firstScatteringFoil,
186 Ta,
187 "FirstScatteringFoil");
188
189 physiFirstScatteringFoil = new G4PVPlacement(0, G4ThreeVector(firstScatteringFoilXPosition,
190 0.,
191 0.),
192 "FirstScatteringFoil",
193 logicFirstScatteringFoil,
194 physiVacuumZone,
195 false,
196 0);
197
198 // --------------//
199 // Kapton Window //
200 //---------------//
201
202 G4Material* kapton = material -> GetMat("Kapton") ;
203
204 const G4double kaptonWindowXSize = 0.025*mm;
205 const G4double kaptonWindowYSize = 5.25*cm;
206 const G4double kaptonWindowZSize = 5.25*cm;
207 const G4double kaptonWindowXPosition = 60.5075*mm;
208
209 G4Box* solidKaptonWindow = new G4Box("KaptonWindow",
210 kaptonWindowXSize,
211 kaptonWindowYSize,
212 kaptonWindowZSize);
213
214 G4LogicalVolume* logicKaptonWindow = new G4LogicalVolume(solidKaptonWindow,
215 kapton,
216 "KaptonWindow");
217
218 physiKaptonWindow = new G4PVPlacement(0, G4ThreeVector(kaptonWindowXPosition, 0., 0.),
219 "KaptonWindow",
220 logicKaptonWindow,
221 physiVacuumZone,
222 false,
223 0);
224
225 G4VisAttributes * white = new G4VisAttributes( G4Colour());
226 white -> SetVisibility(true);
227 white -> SetForceSolid(true);
228
229 logicKaptonWindow -> SetVisAttributes(white);
230
231 // --------//
232 // Stopper //
233 //---------//
234
235 G4double phi = 90. *deg;
236 // Matrix definition for a 90 deg rotation with respect to Y axis
237 G4RotationMatrix rm;
238 rm.rotateY(phi);
239
240 const G4double innerRadiusStopper = 0.*cm;
241 const G4double hightStopper = 3.5*mm;
242 const G4double startAngleStopper = 0.*deg;
243 const G4double spanningAngleStopper = 360.*deg;
244
245 const G4double stopperXPosition = -2956.04 *mm;
246 const G4double stopperYPosition = 0.*m;
247 const G4double stopperZPosition = 0.*m;
248
249 solidStopper = new G4Tubs("Stopper",
250 innerRadiusStopper,
251 outerRadiusStopper,
252 hightStopper,
253 startAngleStopper,
254 spanningAngleStopper);
255
256 G4Material* brass = material -> GetMat("Brass") ;
257
258 G4LogicalVolume* logicStopper = new G4LogicalVolume(solidStopper, brass, "Stopper", 0, 0, 0);
259
260 physiStopper = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(stopperXPosition,
261 stopperYPosition,
262 stopperZPosition)),
263 "Stopper",
264 logicStopper,
265 mother,
266 false,
267 0);
268
269 G4VisAttributes * red = new G4VisAttributes(G4Colour(1. ,0. ,0.));
270 red-> SetVisibility(true);
271 red-> SetForceSolid(true);
272 logicStopper -> SetVisAttributes(red);
273
274 // -----------------------//
275 // Second scattering foil //
276 // -----------------------//
277
278 const G4double secondScatteringFoilYSize = 52.5 *mm;
279 const G4double secondScatteringFoilZSize = 52.5 *mm;
280
281 const G4double secondScatteringFoilXPosition = -2952.52 *mm;
282 const G4double secondScatteringFoilYPosition = 0 *mm;
283 const G4double secondScatteringFoilZPosition = 0 *mm;
284
285 secondScatteringFoil = new G4Box("SecondScatteringFoil",
286 secondScatteringFoilXSize,
287 secondScatteringFoilYSize,
288 secondScatteringFoilZSize);
289
290 G4LogicalVolume* logicSecondScatteringFoil = new G4LogicalVolume(secondScatteringFoil,
291 Ta,
292 "SecondScatteringFoil");
293
294 physiSecondScatteringFoil = new G4PVPlacement(0, G4ThreeVector(secondScatteringFoilXPosition,
295 secondScatteringFoilYPosition,
296 secondScatteringFoilZPosition),
297 "SeconScatteringFoil",
298 logicSecondScatteringFoil,
299 mother,
300 false,
301 0);
302
303 logicSecondScatteringFoil -> SetVisAttributes(white);
304}
305
306void HadrontherapyBeamLine::HadrontherapyBeamCollimators()
307{
308 // -----------------//
309 // First collimator //
310 // -----------------//
311
312 const G4double firstCollimatorXSize = 20.*mm;
313 const G4double firstCollimatorYSize = 100.*mm;
314 const G4double firstCollimatorZSize = 100.*mm;
315
316 const G4double firstCollimatorXPosition = -2932.5*mm;
317 const G4double firstCollimatorYPosition = 0.*mm;
318 const G4double firstCollimatorZPosition = 0.*mm;
319
320 G4Material* PMMA = material -> GetMat("PMMA");
321
322 G4Box* solidFirstCollimator = new G4Box("FirstCollimator",
323 firstCollimatorXSize,
324 firstCollimatorYSize,
325 firstCollimatorZSize);
326
327
328 G4LogicalVolume* logicFirstCollimator = new G4LogicalVolume(solidFirstCollimator,
329 PMMA,
330 "FirstCollimator");
331
332 physiFirstCollimator = new G4PVPlacement(0, G4ThreeVector(firstCollimatorXPosition,
333 firstCollimatorYPosition,
334 firstCollimatorZPosition),
335 "FirstCollimator",
336 logicFirstCollimator,
337 mother,
338 false,
339 0);
340
341 // ----------------------------//
342 // Hole of the first collimator//
343 //-----------------------------//
344
345 G4double innerRadiusHoleFirstCollimator = 0.*mm;
346 G4double outerRadiusHoleFirstCollimator = 15.*mm;
347 G4double hightHoleFirstCollimator = 20.*mm;
348 G4double startAngleHoleFirstCollimator = 0.*deg;
349 G4double spanningAngleHoleFirstCollimator = 360.*deg;
350
351
352 G4Tubs* solidHoleFirstCollimator = new G4Tubs("HoleFirstCollimator",
353 innerRadiusHoleFirstCollimator,
354 outerRadiusHoleFirstCollimator,
355 hightHoleFirstCollimator,
356 startAngleHoleFirstCollimator,
357 spanningAngleHoleFirstCollimator);
358
359 G4Material* Air = material -> GetMat("Air") ;
360
361 G4LogicalVolume* logicHoleFirstCollimator = new G4LogicalVolume(solidHoleFirstCollimator,
362 Air,
363 "HoleFirstCollimator",
364 0, 0, 0);
365 G4double phi = 90. *deg;
366 // Matrix definition for a 90 deg rotation. Also used for other volumes
367 G4RotationMatrix rm;
368 rm.rotateY(phi);
369
370 physiHoleFirstCollimator = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
371 "HoleFirstCollimator",
372 logicHoleFirstCollimator,
373 physiFirstCollimator,
374 false,
375 0);
376 // --------------//
377 // Range shifter //
378 // --------------//
379
380 const G4double rangeShifterYSize = 176. *mm;
381 const G4double rangeShifterZSize = 176. *mm;
382
383 solidRangeShifterBox = new G4Box("RangeShifterBox",
384 rangeShifterXSize,
385 rangeShifterYSize,
386 rangeShifterZSize);
387
388 logicRangeShifterBox = new G4LogicalVolume(solidRangeShifterBox,
389 RSMat,
390 "RangeShifterBox");
391
392 physiRangeShifterBox = new G4PVPlacement(0,
393 G4ThreeVector(rangeShifterXPosition, 0., 0.),
394 "RangeShifterBox",
395 logicRangeShifterBox,
396 mother,
397 false,
398 0);
399
400 G4VisAttributes * yellow = new G4VisAttributes(G4Colour(1., 1., 0. ));
401 yellow-> SetVisibility(true);
402 yellow-> SetForceSolid(true);
403 logicRangeShifterBox -> SetVisAttributes(yellow);
404
405 // ------------------//
406 // Second collimator //
407 //-------------------.//
408
409 const G4double secondCollimatorXPosition = -2028.5*mm;
410 const G4double secondCollimatorYPosition = 0*mm;
411 const G4double secondCollimatorZPosition = 0*mm;
412
413 const G4double secondCollimatorXSize = 20.*mm;
414 const G4double secondCollimatorYSize = 100.*mm;
415 const G4double secondCollimatorZSize = 100.*mm;
416
417 G4Box* solidSecondCollimator = new G4Box("SecondCollimator",
418 secondCollimatorXSize,
419 secondCollimatorYSize,
420 secondCollimatorZSize);
421
422
423 G4LogicalVolume* logicSecondCollimator = new G4LogicalVolume(solidSecondCollimator,
424 PMMA,
425 "SecondCollimator");
426
427 physiSecondCollimator = new G4PVPlacement(0, G4ThreeVector(secondCollimatorXPosition,
428 secondCollimatorYPosition,
429 secondCollimatorZPosition),
430 "SecondCollimator",
431 logicSecondCollimator,
432 mother,
433 false,
434 0);
435
436
437
438 // ------------------------------//
439 // Hole of the second collimator //
440 // ------------------------------//
441
442 G4double innerRadiusHoleSecondCollimator = 0.*mm;
443 G4double outerRadiusHoleSecondCollimator = 15.*mm;
444 G4double hightHoleSecondCollimator = 20.*mm;
445 G4double startAngleHoleSecondCollimator = 0.*deg;
446 G4double spanningAngleHoleSecondCollimator = 360.*deg;
447
448
449 G4Tubs* solidHoleSecondCollimator = new G4Tubs("HoleSecondCollimator",
450 innerRadiusHoleSecondCollimator,
451 outerRadiusHoleSecondCollimator,
452 hightHoleSecondCollimator,
453 startAngleHoleSecondCollimator,
454 spanningAngleHoleSecondCollimator);
455
456
457 G4LogicalVolume* logicHoleSecondCollimator = new G4LogicalVolume(solidHoleSecondCollimator,
458 Air,
459 "HoleSecondCollimator",
460 0, 0, 0);
461 G4double phi2 = 90. *deg;
462 // Matrix definition for a 90 deg rotation. Also used for other volumes
463 G4RotationMatrix rm2;
464 rm2.rotateY(phi2);
465
466 physiHoleSecondCollimator = new G4PVPlacement(G4Transform3D(rm2, G4ThreeVector()),
467 "HoleSecondCollimator",
468 logicHoleSecondCollimator,
469 physiSecondCollimator,
470 false,
471 0);
472
473
474
475 // ---------------------------------//
476 // First Collimator modulator box //
477 // ---------------------------------//
478
479 const G4double firstCollimatorModulatorXSize = 10.*mm;
480 const G4double firstCollimatorModulatorYSize = 200.*mm;
481 const G4double firstCollimatorModulatorZSize = 200.*mm;
482
483 const G4double firstCollimatorModulatorXPosition = -2660.5*mm;
484 const G4double firstCollimatorModulatorYPosition = 0.*mm;
485 const G4double firstCollimatorModulatorZPosition = 0.*mm;
486
487 G4Material* Al = material -> GetMat("MatAluminum");
488
489 G4Box* solidFirstCollimatorModulatorBox = new G4Box("FirstCollimatorModulatorBox",
490 firstCollimatorModulatorXSize,
491 firstCollimatorModulatorYSize,
492 firstCollimatorModulatorZSize);
493
494 G4LogicalVolume* logicFirstCollimatorModulatorBox = new G4LogicalVolume(solidFirstCollimatorModulatorBox,
495 Al, "FirstCollimatorModulatorBox");
496
497 physiFirstCollimatorModulatorBox = new G4PVPlacement(0, G4ThreeVector(firstCollimatorModulatorXPosition,
498 firstCollimatorModulatorYPosition,
499 firstCollimatorModulatorZPosition),
500 "FirstCollimatorModulatorBox",
501 logicFirstCollimatorModulatorBox,
502 mother, false, 0);
503
504 // ----------------------------------//
505 // Hole of the first modulator box //
506 // ----------------------------------//
507 const G4double innerRadiusHoleFirstCollimatorModulatorBox = 0.*mm;
508 const G4double outerRadiusHoleFirstCollimatorModulatorBox = 31.*mm;
509 const G4double hightHoleFirstCollimatorModulatorBox = 10.*mm;
510 const G4double startAngleHoleFirstCollimatorModulatorBox = 0.*deg;
511 const G4double spanningAngleHoleFirstCollimatorModulatorBox = 360.*deg;
512
513 G4Tubs* solidHoleFirstCollimatorModulatorBox = new G4Tubs("HoleFirstCollimatorModulatorBox",
514 innerRadiusHoleFirstCollimatorModulatorBox,
515 outerRadiusHoleFirstCollimatorModulatorBox,
516 hightHoleFirstCollimatorModulatorBox ,
517 startAngleHoleFirstCollimatorModulatorBox,
518 spanningAngleHoleFirstCollimatorModulatorBox);
519
520 G4LogicalVolume* logicHoleFirstCollimatorModulatorBox = new G4LogicalVolume(solidHoleFirstCollimatorModulatorBox,
521 Air, "HoleFirstCollimatorModulatorBox", 0, 0, 0);
522
523 physiHoleFirstCollimatorModulatorBox = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
524 "HoleFirstCollimatorModulatorBox",
525 logicHoleFirstCollimatorModulatorBox,
526 physiFirstCollimatorModulatorBox, false, 0);
527 // -------------------------------------------//
528 // Second collimator modulator box //
529 // -------------------------------------------//
530
531
532 const G4double secondCollimatorModulatorXSize = 10.*mm;
533 const G4double secondCollimatorModulatorYSize = 200.*mm;
534 const G4double secondCollimatorModulatorZSize = 200.*mm;
535
536 const G4double secondCollimatorModulatorXPosition = -2090.5 *mm;
537 const G4double secondCollimatorModulatorYPosition = 0.*mm;
538 const G4double secondCollimatorModulatorZPosition = 0.*mm;
539
540
541 G4Box* solidSecondCollimatorModulatorBox = new G4Box("SecondCollimatorModulatorBox",
542 secondCollimatorModulatorXSize,
543 secondCollimatorModulatorYSize,
544 secondCollimatorModulatorZSize);
545
546 G4LogicalVolume* logicSecondCollimatorModulatorBox = new G4LogicalVolume(solidSecondCollimatorModulatorBox,
547 Al, "SecondCollimatorModulatorBox");
548
549 physiSecondCollimatorModulatorBox = new G4PVPlacement(0, G4ThreeVector(secondCollimatorModulatorXPosition,
550 secondCollimatorModulatorYPosition,
551 secondCollimatorModulatorZPosition),
552 "SecondCollimatorModulatorBox",
553 logicSecondCollimatorModulatorBox,
554 mother, false, 0);
555
556
557
558 // -------------------------------//
559 // Hole of the second collimator modulator box //
560 // -------------------------------//
561
562 const G4double innerRadiusHoleSecondCollimatorModulatorBox = 0.*mm;
563 const G4double outerRadiusHoleSecondCollimatorModulatorBox = 31.*mm;
564 const G4double hightHoleSecondCollimatorModulatorBox = 10.*mm;
565 const G4double startAngleHoleSecondCollimatorModulatorBox = 0.*deg;
566 const G4double spanningAngleHoleSecondCollimatorModulatorBox = 360.*deg;
567
568 G4Tubs* solidHoleSecondCollimatorModulatorBox = new G4Tubs("HoleSecondCollimatorModulatorBox",
569 innerRadiusHoleSecondCollimatorModulatorBox,
570 outerRadiusHoleSecondCollimatorModulatorBox,
571 hightHoleSecondCollimatorModulatorBox ,
572 startAngleHoleSecondCollimatorModulatorBox,
573 spanningAngleHoleSecondCollimatorModulatorBox);
574
575 G4LogicalVolume* logicHoleSecondCollimatorModulatorBox = new G4LogicalVolume(solidHoleSecondCollimatorModulatorBox,
576 Air, "HoleSecondCollimatorModulatorBox", 0, 0, 0);
577
578 physiHoleSecondCollimatorModulatorBox = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
579 "HoleSecondCollimatorModulatorBox",
580 logicHoleSecondCollimatorModulatorBox,
581 physiSecondCollimatorModulatorBox, false, 0);
582
583 G4VisAttributes * blue = new G4VisAttributes( G4Colour(0. ,0. ,1.));
584 blue -> SetVisibility(true);
585 blue -> SetForceSolid(true);
586
587 logicFirstCollimator -> SetVisAttributes(yellow);
588 logicFirstCollimatorModulatorBox -> SetVisAttributes(blue);
589 logicSecondCollimatorModulatorBox -> SetVisAttributes(blue);
590
591}
592
593void HadrontherapyBeamLine::HadrontherapyBeamMonitoring()
594{
595 // ----------------------------
596 // Firt monitor chamber
597 // ----------------------------
598
599 // Each chamber consist of 9 mm of air in a box
600 // that has two layers one of kapton and one
601 // of copper
602
603 G4Material* Kapton = material -> GetMat("Kapton") ;
604 G4Material* Cu = material -> GetMat("MatCopper");
605 G4Material* Air = material -> GetMat("Air") ;
606
607 const G4double monitor1XSize = 4.525022*mm;
608 const G4double monitor2XSize = 0.000011*mm;
609 const G4double monitor3XSize = 4.5*mm;
610 const G4double monitorYSize = 10.*cm;
611 const G4double monitorZSize = 10.*cm;
612
613 const G4double monitor1XPosition = -1765.97498 *mm;
614 const G4double monitor2XPosition = -4.500011*mm;
615 const G4double monitor4XPosition = 4.500011*mm;
616
617 G4Box* solidFirstMonitorLayer1 = new G4Box("FirstMonitorLayer1", monitor1XSize, monitorYSize, monitorZSize);
618
619 G4LogicalVolume* logicFirstMonitorLayer1 = new G4LogicalVolume(solidFirstMonitorLayer1, Kapton, "FirstMonitorLayer1");
620
621 physiFirstMonitorLayer1 = new G4PVPlacement(0,G4ThreeVector(monitor1XPosition,0.*cm,0.*cm),
622 "FirstMonitorLayer1", logicFirstMonitorLayer1, mother, false, 0);
623
624
625
626 G4Box* solidFirstMonitorLayer2 = new G4Box("FirstMonitorLayer2", monitor2XSize, monitorYSize, monitorZSize);
627
628 G4LogicalVolume* logicFirstMonitorLayer2 = new G4LogicalVolume(solidFirstMonitorLayer2, Cu, "FirstMonitorLayer2");
629
630 physiFirstMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector(monitor2XPosition,0.*cm,0.*cm),
631 "FirstMonitorLayer2", logicFirstMonitorLayer2, physiFirstMonitorLayer1,
632 false, 0);
633
634
635
636 G4Box* solidFirstMonitorLayer3 = new G4Box("FirstMonitorLayer3", monitor3XSize, monitorYSize, monitorZSize);
637
638 G4LogicalVolume* logicFirstMonitorLayer3 = new G4LogicalVolume(solidFirstMonitorLayer3, Air, "FirstMonitorLayer3");
639
640 physiFirstMonitorLayer3 = new G4PVPlacement(0, G4ThreeVector(0.*mm,0.*cm,0.*cm), "MonitorLayer3",
641 logicFirstMonitorLayer3, physiFirstMonitorLayer1, false, 0);
642
643
644
645 G4Box* solidFirstMonitorLayer4 = new G4Box("FirstMonitorLayer4", monitor2XSize, monitorYSize, monitorZSize);
646
647 G4LogicalVolume* logicFirstMonitorLayer4 = new G4LogicalVolume(solidFirstMonitorLayer4, Cu, "FirstMonitorLayer4");
648
649 physiFirstMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm),
650 "FirstMonitorLayer4", logicFirstMonitorLayer4, physiFirstMonitorLayer1, false, 0);
651
652
653
654 // ------------------------//
655 // Second monitor chamber //
656 // ------------------------//
657
658
659 G4Box* solidSecondMonitorLayer1 = new G4Box("SecondMonitorLayer1", monitor1XSize, monitorYSize, monitorZSize);
660
661 G4LogicalVolume* logicSecondMonitorLayer1 = new G4LogicalVolume(solidSecondMonitorLayer1, Kapton, "SecondMonitorLayer1");
662
663 physiSecondMonitorLayer1 = new G4PVPlacement(0, G4ThreeVector(-1634.92493 *mm,0.*cm,0.*cm),
664 "SecondMonitorLayer1", logicSecondMonitorLayer1,mother, false, 0);
665
666
667
668 G4Box* solidSecondMonitorLayer2 = new G4Box("SecondMonitorLayer2", monitor2XSize, monitorYSize, monitorZSize);
669
670 G4LogicalVolume* logicSecondMonitorLayer2 = new G4LogicalVolume(solidSecondMonitorLayer2, Cu, "SecondMonitorLayer2");
671
672 physiSecondMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector( monitor2XPosition,0.*cm,0.*cm), "SecondMonitorLayer2",
673 logicSecondMonitorLayer2, physiSecondMonitorLayer1, false, 0);
674
675
676
677 G4Box* solidSecondMonitorLayer3 = new G4Box("SecondMonitorLayer3", monitor3XSize, monitorYSize, monitorZSize);
678
679 G4LogicalVolume* logicSecondMonitorLayer3 = new G4LogicalVolume(solidSecondMonitorLayer3, Air, "SecondMonitorLayer3");
680
681 physiSecondMonitorLayer3 = new G4PVPlacement(0, G4ThreeVector(0.*mm,0.*cm,0.*cm), "MonitorLayer3",
682 logicSecondMonitorLayer3, physiSecondMonitorLayer1, false, 0);
683
684
685
686 G4Box* solidSecondMonitorLayer4 = new G4Box("SecondMonitorLayer4", monitor2XSize, monitorYSize, monitorZSize);
687
688 G4LogicalVolume* logicSecondMonitorLayer4 = new G4LogicalVolume(solidSecondMonitorLayer4, Cu, "SecondMonitorLayer4");
689
690 physiSecondMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm), "SecondMonitorLayer4",
691 logicSecondMonitorLayer4, physiSecondMonitorLayer1, false, 0);
692
693 // -----------------------//
694 // Third monitor chamber //
695 // -----------------------//
696
697
698 G4Box* solidThirdMonitorLayer1 = new G4Box("ThirdMonitorLayer1", monitor1XSize, monitorYSize, monitorZSize);
699
700 G4LogicalVolume* logicThirdMonitorLayer1 = new G4LogicalVolume(solidThirdMonitorLayer1, Kapton, "ThirdMonitorLayer1");
701
702 physiThirdMonitorLayer1 = new G4PVPlacement(0, G4ThreeVector(-1505.87489 *mm,0.*cm,0.*cm),
703 "ThirdMonitorLayer1", logicThirdMonitorLayer1, mother, false, 0);
704
705
706
707 G4Box* solidThirdMonitorLayer2 = new G4Box("ThirdMonitorLayer2", monitor2XSize, monitorYSize, monitorZSize);
708
709 G4LogicalVolume* logicThirdMonitorLayer2 = new G4LogicalVolume(solidThirdMonitorLayer2, Cu, "ThirdMonitorLayer2");
710
711 physiThirdMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector(monitor2XPosition, 0.*cm,0.*cm),
712 "ThirdMonitorLayer2",
713 logicThirdMonitorLayer2,
714 physiThirdMonitorLayer1,
715 false, 0);
716
717
718
719 G4Box* solidThirdMonitorLayer3 = new G4Box("ThirdMonitorLayer3", monitor3XSize, monitorYSize, monitorZSize);
720
721 G4LogicalVolume* logicThirdMonitorLayer3 = new G4LogicalVolume(solidThirdMonitorLayer3, Air, "ThirdMonitorLayer3");
722
723 physiThirdMonitorLayer3 = new G4PVPlacement(0, G4ThreeVector(0.*mm,0.*cm,0.*cm), "MonitorLayer3",
724 logicThirdMonitorLayer3, physiThirdMonitorLayer1, false, 0);
725
726
727
728
729
730
731 G4Box* solidThirdMonitorLayer4 = new G4Box("ThirdMonitorLayer4", monitor2XSize, monitorYSize, monitorZSize);
732
733 G4LogicalVolume* logicThirdMonitorLayer4 = new G4LogicalVolume(solidThirdMonitorLayer4, Cu, "ThirdMonitorLayer4");
734
735 physiThirdMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm), "ThirdMonitorLayer4",
736 logicThirdMonitorLayer4, physiThirdMonitorLayer1, false, 0);
737}
738
739
740
741void HadrontherapyBeamLine::HadrontherapyBeamNozzle()
742{
743 // ---------------//
744 // Nozzle support //
745 //----------------//
746
747 const G4double nozzleSupportXSize = 29.50 *mm;
748 const G4double nozzleSupportYSize = 180. *mm;
749 const G4double nozzleSupportZSize = 180. *mm;
750
751 const G4double nozzleSupportXPosition = -601.00 *mm;
752
753 G4Material* PMMA = material -> GetMat("PMMA");
754 G4Material* Brass = material -> GetMat("Brass") ;
755 // G4Material* Air = material -> GetMat("Air") ;
756
757 G4double phi = 90. *deg;
758
759 // Matrix definition for a 90 deg rotation. Also used for other volumes
760 G4RotationMatrix rm;
761 rm.rotateY(phi);
762
763 // G4Subtraction
764 G4Box* solidNozzleBox = new G4Box("NozzleBox", nozzleSupportXSize, nozzleSupportYSize, nozzleSupportZSize);
765
766
767 const G4double innerRadiusHoleNozzle = 0.*mm;
768 const G4double outerRadiusHoleNozzle = 21.5 *mm;
769 const G4double hightHoleNozzle = 29.5 *mm;
770 const G4double startAngleHoleNozzle = 0.*deg;
771 const G4double spanningAngleHoleNozzle = 360.*deg;
772
773 G4Tubs* solidHoleNozzle = new G4Tubs("HoleNozzle",
774 innerRadiusHoleNozzle,
775 outerRadiusHoleNozzle,
776 hightHoleNozzle,
777 startAngleHoleNozzle,
778 spanningAngleHoleNozzle);
779
780 G4SubtractionSolid* solidNozzleSupport = new G4SubtractionSolid("NozzleSupport",solidNozzleBox, solidHoleNozzle,
781 &rm, G4ThreeVector(0 ,0,0 ));
782
783 G4LogicalVolume* logicNozzleSupport = new G4LogicalVolume(solidNozzleSupport, PMMA, "NozzleSupport");
784
785 physiNozzleSupport = new G4PVPlacement(0, G4ThreeVector(nozzleSupportXPosition,0., 0.),
786 "NozzleSupport", logicNozzleSupport, mother, false, 0);
787
788
789 G4VisAttributes * blue = new G4VisAttributes( G4Colour(0. ,0. ,1.));
790 blue -> SetVisibility(true);
791 blue -> SetForceSolid(false);
792
793 logicNozzleSupport -> SetVisAttributes(blue);
794
795 //logicHoleNozzle -> SetVisAttributes(blue);
796
797 // ---------------------------------//
798 // First hole of the noozle support //
799 // ---------------------------------//
800
801 const G4double innerRadiusHoleNozzleSupport = 18.*mm;
802 const G4double outerRadiusHoleNozzleSupport = 21.5 *mm;
803 const G4double hightHoleNozzleSupport = 185.*mm;
804 const G4double startAngleHoleNozzleSupport = 0.*deg;
805 const G4double spanningAngleHoleNozzleSupport = 360.*deg;
806
807 const G4double holeNozzleSupportXPosition = -475.5 *mm;
808
809 G4Tubs* solidHoleNozzleSupport = new G4Tubs("HoleNozzleSupport",
810 innerRadiusHoleNozzleSupport,
811 outerRadiusHoleNozzleSupport,
812 hightHoleNozzleSupport,
813 startAngleHoleNozzleSupport,
814 spanningAngleHoleNozzleSupport);
815
816 G4LogicalVolume* logicHoleNozzleSupport = new G4LogicalVolume(solidHoleNozzleSupport,
817 Brass,
818 "HoleNozzleSupport",
819 0, 0, 0);
820
821 physiHoleNozzleSupport = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(holeNozzleSupportXPosition,
822 0., 0.)),
823 "HoleNozzleSupport", logicHoleNozzleSupport, mother, false, 0);
824
825 G4VisAttributes * yellow = new G4VisAttributes( G4Colour(1., 1., 0. ));
826 yellow-> SetVisibility(true);
827 yellow-> SetForceSolid(true);
828 logicHoleNozzleSupport -> SetVisAttributes(yellow);
829
830}
831
832void HadrontherapyBeamLine::HadrontherapyBeamFinalCollimator()
833{
834
835 // -----------------------//
836 // Final collimator //
837 //------------------------//
838
839 const G4double outerRadiusFinalCollimator = 21.5*mm;
840 const G4double hightFinalCollimator = 3.5*mm;
841 const G4double startAngleFinalCollimator = 0.*deg;
842 const G4double spanningAngleFinalCollimator = 360.*deg;
843 const G4double finalCollimatorXPosition = -287.0 *mm;
844
845 G4double phi = 90. *deg;
846
847 // Matrix definition for a 90 deg rotation. Also used for other volumes
848 G4RotationMatrix rm;
849 rm.rotateY(phi);
850
851 G4Material* Brass = material -> GetMat("Brass");
852
853 solidFinalCollimator = new G4Tubs("FinalCollimator", innerRadiusFinalCollimator,
854 outerRadiusFinalCollimator,
855 hightFinalCollimator,
856 startAngleFinalCollimator,
857 spanningAngleFinalCollimator);
858
859 G4LogicalVolume* logicFinalCollimator = new G4LogicalVolume(solidFinalCollimator,
860 Brass, "FinalCollimator", 0, 0, 0);
861
862 physiFinalCollimator = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(finalCollimatorXPosition,0.,0.)),
863 "FinalCollimator", logicFinalCollimator, mother, false, 0);
864
865 G4VisAttributes * yellow = new G4VisAttributes( G4Colour(1., 1., 0. ));
866 yellow-> SetVisibility(true);
867 yellow-> SetForceSolid(true);
868 logicFinalCollimator -> SetVisAttributes(yellow);
869
870
871}
872
873void HadrontherapyBeamLine::SetRangeShifterXPosition(G4double value)
874{
875 physiRangeShifterBox -> SetTranslation(G4ThreeVector(value, 0., 0.));
876 G4cout << "The Range Shifter is translated to"<< value/mm <<"mm along the X axis" <<G4endl;
877}
878
879void HadrontherapyBeamLine::SetRangeShifterXSize(G4double value)
880{
881 solidRangeShifterBox -> SetXHalfLength(value) ;
882 G4cout << "RangeShifter size X (mm): "<< ((solidRangeShifterBox -> GetXHalfLength())*2.)/mm
883 << G4endl;
884}
885
886void HadrontherapyBeamLine::SetFirstScatteringFoilXSize(G4double value)
887{
888 firstScatteringFoil -> SetXHalfLength(value);
889 G4cout <<"The X size of the first scattering foil is (mm):"<<
890 ((firstScatteringFoil -> GetXHalfLength())*2.)/mm
891 << G4endl;
892}
893
894void HadrontherapyBeamLine::SetSecondScatteringFoilXSize(G4double value)
895{
896 secondScatteringFoil -> SetXHalfLength(value);
897 G4cout <<"The X size of the second scattering foil is (mm):"<<
898 ((secondScatteringFoil -> GetXHalfLength())*2.)/mm
899 << G4endl;
900}
901
902void HadrontherapyBeamLine::SetOuterRadiusStopper(G4double value)
903{
904 solidStopper -> SetOuterRadius(value);
905 G4cout << "OuterRadius od the Stopper is (mm):"
906 << solidStopper -> GetOuterRadius()/mm
907 << G4endl;
908}
909
910void HadrontherapyBeamLine::SetInnerRadiusFinalCollimator(G4double value)
911{
912 solidFinalCollimator -> SetInnerRadius(value);
913 G4cout<<"Inner Radius of the final collimator is (mm):"
914 << solidFinalCollimator -> GetInnerRadius()/mm
915 << G4endl;
916}
917
918void HadrontherapyBeamLine::SetRSMaterial(G4String materialChoice)
919{
920 G4Material* pttoMaterial = G4Material::GetMaterial(materialChoice);
921
922 if (pttoMaterial)
923 {
924 RSMat = pttoMaterial;
925 logicRangeShifterBox -> SetMaterial(pttoMaterial);
926 }
927}
Note: See TracBrowser for help on using the repository browser.