Geant4 Cross Reference |
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 // Hadrontherapy advanced example for Geant4 27 // See more at: https://twiki.cern.ch/twiki/bin/view/Geant4/AdvancedExamplesHadrontherapy 28 // Simulation of the "Zero degree" experimental beamline of INFN-LNS (Catania, Italy). 29 30 #include "G4Box.hh" 31 #include "G4Tubs.hh" 32 #include "G4VisAttributes.hh" 33 #include "G4Colour.hh" 34 #include "globals.hh" 35 #include "G4RunManager.hh" 36 #include "G4LogicalVolume.hh" 37 #include "G4PVPlacement.hh" 38 #include "G4RotationMatrix.hh" 39 #include "G4NistManager.hh" 40 #include "G4NistElementBuilder.hh" 41 #include "HadrontherapyDetectorConstruction.hh" 42 #include "HadrontherapyModulator.hh" 43 #include "PassiveCarbonBeamLine.hh" 44 #include "G4SystemOfUnits.hh" 45 #include "G4Trd.hh" 46 #include "PassiveCarbonBeamLineMessenger.hh" 47 #include "G4UnionSolid.hh" 48 49 using namespace std; 50 51 //G4bool PassiveCarbonBeamLine::doCalculation = false; 52 ///////////////////////////////////////////////////////////////////////////// 53 PassiveCarbonBeamLine::PassiveCarbonBeamLine(): 54 physicalTreatmentRoom(0),hadrontherapyDetectorConstruction(0), 55 physiBeamLineSupport(0), physiBeamLineCover(0), physiBeamLineCover2(0), 56 physiKaptonWindow(0),PhysiRippleFilter(0),PhysiRippleFilterBase(0),PhysiRippleFilterTrd(0), 57 physiFirstMonitorLayer1(0), physiFirstMonitorLayer2(0), 58 physiFirstMonitorLayer3(0), physiFirstMonitorLayer4(0), 59 physiNozzleSupport(0), physiHoleNozzleSupport(0) 60 { 61 62 // Messenger to change parameters of the passiveCarbonBeamLine geometry 63 PassiveCarbonMessenger = new PassiveCarbonBeamLineMessenger(this); 64 65 //***************************** PW *************************************** 66 67 static G4String ROGeometryName = "DetectorROGeometry"; 68 RO = new HadrontherapyDetectorROGeometry(ROGeometryName); 69 70 71 G4cout << "Going to register Parallel world..."; 72 RegisterParallelWorld(RO); 73 G4cout << "... done" << G4endl; 74 //***************************** PW *************************************** 75 } 76 77 ///////////////////////////////////////////////////////////////////////////// 78 PassiveCarbonBeamLine::~PassiveCarbonBeamLine() 79 { 80 delete hadrontherapyDetectorConstruction; 81 delete PassiveCarbonMessenger; 82 } 83 84 ///////////////////////////////////////////////////////////////////////////// 85 G4VPhysicalVolume* PassiveCarbonBeamLine::Construct() 86 { 87 // Sets default geometry and materials 88 SetDefaultDimensions(); 89 90 // Construct the whole CarbonPassive Beam Line 91 ConstructPassiveCarbonBeamLine(); 92 93 94 //***************************** PW *************************************** 95 if (!hadrontherapyDetectorConstruction) 96 97 //***************************** PW *************************************** 98 // HadrontherapyDetectorConstruction builds ONLY the phantom and the detector with its associated ROGeometry 99 hadrontherapyDetectorConstruction = new HadrontherapyDetectorConstruction(physicalTreatmentRoom); 100 101 //***************************** PW *************************************** 102 103 hadrontherapyDetectorConstruction->InitializeDetectorROGeometry(RO,hadrontherapyDetectorConstruction->GetDetectorToWorldPosition()); 104 105 //***************************** PW *************************************** 106 return physicalTreatmentRoom; 107 } 108 109 // In the following method the DEFAULTS used in the geometry of 110 // passive beam line are provided 111 // HERE THE USER CAN CHANGE THE GEOMETRY CHARACTERISTICS OF BEAM 112 // LINE ELEMENTS, ALTERNATIVELY HE/SHE CAN USE THE MACRO FILE (IF A 113 // MESSENGER IS PROVIDED) 114 // 115 // DEFAULT MATERIAL ARE ALSO PROVIDED 116 // and COLOURS ARE ALSO DEFINED 117 // ---------------------------------------------------------- 118 ///////////////////////////////////////////////////////////////////////////// 119 void PassiveCarbonBeamLine::SetDefaultDimensions() 120 { 121 // Set of coulors that can be used 122 white = new G4VisAttributes( G4Colour()); 123 white -> SetVisibility(true); 124 white -> SetForceSolid(true); 125 126 black = new G4VisAttributes( G4Colour(1., 1., 1.)); 127 black -> SetVisibility(true); 128 black -> SetForceSolid(true); 129 130 131 blue = new G4VisAttributes(G4Colour(0. ,0. ,1.)); 132 blue -> SetVisibility(true); 133 blue -> SetForceSolid(true); 134 135 gray = new G4VisAttributes( G4Colour(0.5, 0.5, 0.5 )); 136 gray-> SetVisibility(true); 137 gray-> SetForceSolid(true); 138 139 red = new G4VisAttributes(G4Colour(1. ,0. ,0.)); 140 red-> SetVisibility(true); 141 red-> SetForceSolid(true); 142 143 yellow = new G4VisAttributes(G4Colour(1., 1., 0. )); 144 yellow-> SetVisibility(true); 145 yellow-> SetForceSolid(true); 146 147 green = new G4VisAttributes( G4Colour(25/255. , 255/255. , 25/255. )); 148 green -> SetVisibility(true); 149 green -> SetForceSolid(true); 150 151 darkGreen = new G4VisAttributes( G4Colour(0/255. , 100/255. , 0/255. )); 152 darkGreen -> SetVisibility(true); 153 darkGreen -> SetForceSolid(true); 154 155 darkOrange3 = new G4VisAttributes( G4Colour(205/255. , 102/255. , 000/255. )); 156 darkOrange3 -> SetVisibility(true); 157 darkOrange3 -> SetForceSolid(true); 158 159 skyBlue = new G4VisAttributes( G4Colour(135/255. , 206/255. , 235/255. )); 160 skyBlue -> SetVisibility(true); 161 skyBlue -> SetForceSolid(true); 162 163 164 // FINAL COLLIMATOR: is the collimator giving the final transversal shape 165 // of the beam 166 167 G4double defaultinnerRadiusFinalCollimator = 12.5 *mm; 168 innerRadiusFinalCollimator = defaultinnerRadiusFinalCollimator; 169 170 // DEFAULT DEFINITION OF THE MATERIALS 171 // All elements and compound definition follows the NIST database 172 173 // ELEMENTS 174 G4bool isotopes = false; 175 aluminumNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Al", isotopes); 176 G4Element* zincNist = G4NistManager::Instance()->FindOrBuildElement("Zn"); 177 G4Element* copperNist = G4NistManager::Instance()->FindOrBuildElement("Cu"); 178 179 // MATERIAL (Including compounds) 180 copperNistMaterial = G4NistManager::Instance()->FindOrBuildMaterial("G4_Cu", isotopes); 181 airNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_AIR", isotopes); 182 kaptonNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_KAPTON", isotopes); 183 galacticNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Galactic", isotopes); 184 PMMANist = G4NistManager::Instance()->FindOrBuildMaterial("G4_PLEXIGLASS", isotopes); 185 tantalumNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Ta", isotopes); 186 187 G4double d; // Density 188 G4int nComponents;// Number of components 189 G4double fractionmass; // Fraction in mass of an element in a material 190 191 d = 8.40*g/cm3; 192 nComponents = 2; 193 brass = new G4Material("Brass", d, nComponents); 194 brass -> AddElement(zincNist, fractionmass = 30 *perCent); 195 brass -> AddElement(copperNist, fractionmass = 70 *perCent); 196 197 //***************************** PW *************************************** 198 199 // DetectorROGeometry Material 200 new G4Material("dummyMat", 1., 1.*g/mole, 1.*g/cm3); 201 202 //***************************** PW *************************************** 203 204 // MATERIAL ASSIGNMENT 205 // Support of the beam line 206 beamLineSupportMaterial = aluminumNist; 207 208 // Vacuum pipe 209 vacuumZoneMaterial = galacticNist; 210 firstScatteringFoilMaterial = tantalumNist; 211 212 // Material of kapton window 213 kaptonWindowMaterial = kaptonNist; 214 215 // Material of ripple filter 216 rippleFilterMaterial = PMMANist; 217 rippleFilterBoxMaterial = airNist; 218 219 // Materials of the monitor chamber 220 layer1MonitorChamberMaterial = kaptonNist; 221 layer2MonitorChamberMaterial = copperNistMaterial; 222 layer3MonitorChamberMaterial = airNist; 223 layer4MonitorChamberMaterial = copperNistMaterial; 224 225 // Material of the final nozzle 226 nozzleSupportMaterial = PMMANist; 227 holeNozzleSupportMaterial = airNist; 228 seconHoleNozzleSupportMaterial = airNist; 229 230 // Material of the final collimator 231 brassTubeMaterial = brassTube2Material = brassTube3Material = brass; 232 233 // Material of the final collimator 234 finalCollimatorMaterial = brass; 235 236 // Material of the PMMA collimator 237 PMMACollimatorMaterial = airNist; 238 239 240 G4Element* hydrogenNist = G4NistManager::Instance()->FindOrBuildElement("H"); 241 G4Element* oxygenNist = G4NistManager::Instance()->FindOrBuildElement("O"); 242 G4Element* carbonNist = G4NistManager::Instance()->FindOrBuildElement("C"); 243 244 G4Material* plastic = new G4Material("Plastic", 1.18*g/cm3, nComponents=3); 245 plastic -> AddElement(carbonNist, 21); 246 plastic -> AddElement(oxygenNist, 4); 247 plastic -> AddElement(hydrogenNist, 24); 248 PMMANist = plastic; 249 250 251 252 } 253 254 ///////////////////////////////////////////////////////////////////////////// 255 void PassiveCarbonBeamLine::ConstructPassiveCarbonBeamLine() 256 { 257 // ----------------------------- 258 // Treatment room - World volume 259 //------------------------------ 260 // Treatment room sizes 261 262 const G4double worldX = 400.0 *cm; 263 const G4double worldY = 400.0 *cm; 264 const G4double worldZ = 400.0 *cm; 265 G4bool isotopes = false; 266 267 airNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_AIR", isotopes); 268 treatmentRoom = new G4Box("TreatmentRoom",worldX,worldY,worldZ); 269 logicTreatmentRoom = new G4LogicalVolume(treatmentRoom, 270 airNist, 271 "logicTreatmentRoom", 272 0,0,0); 273 physicalTreatmentRoom = new G4PVPlacement(0, 274 G4ThreeVector(), 275 "physicalTreatmentRoom", 276 logicTreatmentRoom, 277 0,false,0); 278 279 280 // The treatment room is invisible in the Visualisation 281 logicTreatmentRoom -> SetVisAttributes (G4VisAttributes::GetInvisible()); 282 283 // Components of the Passive Carbon Beam Line 284 HadrontherapyBeamLineSupport(); 285 VacuumToAirInterface(); 286 HadrontherapyBeamMonitoring(); 287 HadrontherapyBeamNozzle(); 288 HadrontherapyBeamFinalCollimator(); 289 HadrontherapyPMMACollimator(); 290 // HadrontherapyRippleFilter(); 291 292 293 // StopperCostruction(); 294 295 296 HadrontherapyRidgeFilter(); 297 } 298 299 ///////////////////////////////////////////////////////////////////////////// 300 void PassiveCarbonBeamLine::HadrontherapyBeamLineSupport() 301 { 302 // ------------------// 303 // BEAM LINE SUPPORT // 304 //-------------------// 305 306 beamLineSupportXSize = 1.5*m; 307 beamLineSupportYSize = 20.*mm; 308 beamLineSupportZSize = 600.*mm; 309 310 beamLineSupportXPosition = -1745.09 *mm; 311 beamLineSupportYPosition = -230. *mm; 312 beamLineSupportZPosition = 0.*mm; 313 314 beamLineSupport = new G4Box("BeamLineSupport", 315 beamLineSupportXSize, 316 beamLineSupportYSize, 317 beamLineSupportZSize); 318 319 logicBeamLineSupport = new G4LogicalVolume(beamLineSupport, 320 beamLineSupportMaterial, 321 "BeamLineSupport"); 322 physiBeamLineSupport = new G4PVPlacement(0, G4ThreeVector(beamLineSupportXPosition, 323 beamLineSupportYPosition, 324 beamLineSupportZPosition), 325 "BeamLineSupport", 326 logicBeamLineSupport, 327 physicalTreatmentRoom, false, 0); 328 329 // Visualisation attributes of the beam line support 330 logicBeamLineSupport -> SetVisAttributes(gray); 331 332 //---------------------------------// 333 // Beam line cover 1 (left panel) // 334 //---------------------------------// 335 beamLineCoverXSize = 1.5*m; 336 beamLineCoverYSize = 750.*mm; 337 beamLineCoverZSize = 10.*mm; 338 339 beamLineCoverXPosition = -1745.09 *mm; 340 beamLineCoverYPosition = -1000.*mm; 341 beamLineCoverZPosition = 610.*mm; 342 343 beamLineCover = new G4Box("BeamLineCover", 344 beamLineCoverXSize, 345 beamLineCoverYSize, 346 beamLineCoverZSize); 347 348 logicBeamLineCover = new G4LogicalVolume(beamLineCover, 349 beamLineSupportMaterial, 350 "BeamLineCover"); 351 352 physiBeamLineCover = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition, 353 beamLineCoverYPosition, 354 beamLineCoverZPosition), 355 "BeamLineCover", 356 logicBeamLineCover, 357 physicalTreatmentRoom, 358 false, 359 0); 360 361 // ---------------------------------// 362 // Beam line cover 2 (rigth panel) // 363 // ---------------------------------// 364 // It has the same characteristic of beam line cover 1 but set in a different position 365 physiBeamLineCover2 = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition, 366 beamLineCoverYPosition, 367 - beamLineCoverZPosition), 368 "BeamLineCover2", 369 logicBeamLineCover, 370 physicalTreatmentRoom, 371 false, 372 0); 373 374 375 logicBeamLineCover -> SetVisAttributes(blue); 376 } 377 378 ///////////////////////////////////////////////////////////////////////////// 379 void PassiveCarbonBeamLine::VacuumToAirInterface() 380 { 381 // ------------// 382 // VACUUM PIPE // 383 //-------------// 384 // 385 // First track of the beam line is inside vacuum; 386 387 vacuumZoneXSize = 100 *mm; 388 vacuumZoneYSize = 52.5 *mm; 389 vacuumZoneZSize = 52.5 *mm; 390 vacuumPipeXPosition = -1708.0 *mm; 391 392 393 vacuumZone = new G4Box("VacuumZone", 394 vacuumZoneXSize/2, 395 vacuumZoneYSize/2, 396 vacuumZoneZSize/2); 397 398 logicVacuumZone = new G4LogicalVolume(vacuumZone, vacuumZoneMaterial,"VacuumZone"); 399 400 physiVacuumZone = new G4PVPlacement(0, G4ThreeVector(vacuumPipeXPosition, 0., 0.), 401 "VacuumZone", 402 logicVacuumZone, 403 physicalTreatmentRoom, 404 false, 405 0); 406 407 // --------------------------// 408 // THE FIRST SCATTERING FOIL // 409 // --------------------------// 410 // A thin foil performing a first scattering 411 // of the original beam 412 413 firstScatteringFoilXSize = 0.015 *mm; 414 firstScatteringFoilYSize = 52.5 *mm; 415 firstScatteringFoilZSize = 52.5 *mm; 416 firstScatteringFoilXPosition = 0.0 *mm; 417 418 firstScatteringFoil = new G4Box("FirstScatteringFoil", 419 firstScatteringFoilXSize/2, 420 firstScatteringFoilYSize/2, 421 firstScatteringFoilZSize/2); 422 423 logicFirstScatteringFoil = new G4LogicalVolume(firstScatteringFoil, 424 firstScatteringFoilMaterial, 425 "FirstScatteringFoil"); 426 427 physiFirstScatteringFoil = new G4PVPlacement(0, G4ThreeVector(firstScatteringFoilXPosition, 0.,0.), 428 "FirstScatteringFoil", logicFirstScatteringFoil, physiVacuumZone, 429 false, 0); 430 431 logicFirstScatteringFoil -> SetVisAttributes(skyBlue); 432 433 // -------------------// 434 // THE KAPTON WINDOWS // 435 //--------------------// 436 //It permits the passage of the beam from vacuum to air 437 438 // KAPTON WINDOW: it permits the passage of the beam from vacuum to air 439 kaptonWindowXSize = 0.050 *mm; 440 kaptonWindowYSize = 52.5 *mm; 441 kaptonWindowZSize = 52.5 *mm; 442 kaptonWindowXPosition = vacuumZoneXSize/2 - kaptonWindowXSize/2; 443 444 solidKaptonWindow = new G4Box("KaptonWindow", 445 kaptonWindowXSize/2, 446 kaptonWindowYSize/2, 447 kaptonWindowZSize/2); 448 449 logicKaptonWindow = new G4LogicalVolume(solidKaptonWindow, 450 kaptonWindowMaterial, 451 "KaptonWindow"); 452 453 physiKaptonWindow = new G4PVPlacement(0, G4ThreeVector(kaptonWindowXPosition, 0., 0.), 454 "KaptonWindow", logicKaptonWindow, 455 physiVacuumZone, false, 0); 456 457 logicKaptonWindow -> SetVisAttributes(darkOrange3); 458 } 459 ///////////////////////////////////////////////////////////////////////////// 460 void PassiveCarbonBeamLine::HadrontherapyRippleFilter() 461 { 462 463 G4double defaultRippleFilterXPosition = -1638.0*mm; 464 G4double ripple_position=(defaultRippleFilterXPosition); 465 G4double RF_x = 200.0 * mm; 466 G4double RF_y = 200.0 * mm; 467 G4double RF_z = 1.4 * mm; 468 G4double RFbase_z = 0.2 * mm; 469 G4double RFtrd_z = RF_z - RFbase_z; 470 G4double RFtrd_top = 1e-4 * mm; 471 G4double RFtrd_bottom = 1.5 * mm; 472 G4double distanceBetweenTrd = 0.1*mm; 473 474 475 476 477 G4double theta = -90. *deg; 478 // Matrix definition for a "theta" deg rotation with respect to Y axis 479 G4RotationMatrix rot; 480 rot.rotateY(theta); 481 482 483 SolidRippleFilter= new G4Box("RippleFilter", 484 RF_x/2 + 1*mm, 485 RF_y/2 + 1*mm, 486 RF_z/2 + 1*mm); 487 488 LogicRippleFilter = new G4LogicalVolume(SolidRippleFilter, 489 rippleFilterBoxMaterial, 490 "LogicRippleFilter", 491 0,0,0); 492 493 PhysiRippleFilter = new G4PVPlacement(G4Transform3D(rot,G4ThreeVector(ripple_position,0,0)), 494 "PhysiRippleFilter", 495 LogicRippleFilter, 496 physicalTreatmentRoom, 497 false, 498 1, 499 true); 500 501 PhysiRippleFilter = new G4PVPlacement(G4Transform3D(rot,G4ThreeVector(ripple_position + 10*cm,0,0)), 502 "PhysiRippleFilter", 503 LogicRippleFilter, 504 physicalTreatmentRoom, 505 false, 506 2, 507 true); 508 509 LogicRippleFilter -> SetVisAttributes(G4VisAttributes::GetInvisible()); 510 511 SolidRippleFilterBase = new G4Box("RippleFilterBase", 512 RF_x/2, 513 RF_y/2, 514 RFbase_z/2); 515 516 LogicRippleFilterBase = new G4LogicalVolume(SolidRippleFilterBase, 517 rippleFilterMaterial, 518 "LogicRippleFilterBase", 519 0,0,0); 520 521 LogicRippleFilterBase -> SetVisAttributes(green); 522 523 PhysiRippleFilterBase = new G4PVPlacement(0, 524 G4ThreeVector(0, 0, -RF_z/2 + RFbase_z/2), 525 "PhysiRippleFilter", 526 LogicRippleFilterBase, 527 PhysiRippleFilter, 528 false, 529 0, 530 false); 531 532 SolidRippleFilterTrd = new G4Trd("SolidRippleFilterTrd", 533 RF_x/2, 534 RF_x/2, 535 RFtrd_bottom/2, 536 RFtrd_top/2, 537 RFtrd_z/2); 538 539 LogicRippleFilterTrd = new G4LogicalVolume(SolidRippleFilterTrd, 540 rippleFilterMaterial, 541 "LogicRippleFilterTrd", 542 0,0,0); 543 544 LogicRippleFilterTrd -> SetVisAttributes(green); 545 546 G4int numberOfTrd = static_cast<int>(std::floor( RF_y / (RFtrd_bottom+distanceBetweenTrd) )); 547 548 G4int N = static_cast<int>( std::floor(numberOfTrd-1)/2 ); 549 550 G4int copyNumber = 0; 551 552 for( int i = -N; i <= N; i++ ) 553 { 554 PhysiRippleFilterTrd = new G4PVPlacement(0, 555 G4ThreeVector(0, 556 i*(RFtrd_bottom+distanceBetweenTrd), 557 -RF_z/2+RFbase_z+RFtrd_z/2), 558 "PhysiRippleFilterTrd", 559 LogicRippleFilterTrd, 560 PhysiRippleFilter, 561 false, 562 copyNumber, 563 false); 564 565 copyNumber++; 566 } 567 568 } 569 ///////////////////////////////////////////////////////////////////////////// 570 571 void PassiveCarbonBeamLine::StopperCostruction(){ 572 573 supportFoil= new G4Box("supportFoil",25/2*um,10*cm,10*cm); 574 LogicSupportFoil = new G4LogicalVolume(supportFoil,tantalumNist,"logicSupportFoil"); 575 PhysiSupportFoil = new G4PVPlacement(0,G4ThreeVector(-1398.*mm+25/2*um,0.*mm,0.*mm),"logicSupportFoil",LogicSupportFoil,physicalTreatmentRoom,false,0 ); 576 LogicSupportFoil -> SetVisAttributes(skyBlue); 577 578 579 580 stopper = new G4Tubs("Stopper",0*mm,3*mm,3.5*mm,0.*deg,360.*deg); 581 LogicStopper= new G4LogicalVolume(stopper,brass,"logicalStopper"); 582 583 G4double ti = -270. *deg; 584 G4RotationMatrix rt; 585 rt.rotateY(ti); 586 587 PhysicStopper= new G4PVPlacement(G4Transform3D(rt,G4ThreeVector(-1398.*mm-3.5*mm,0*mm,0*mm)),"logicalStopper",LogicStopper,physicalTreatmentRoom,false,0); 588 589 LogicStopper -> SetVisAttributes(red); 590 591 592 } 593 594 void PassiveCarbonBeamLine::HadrontherapyRidgeFilter(){ 595 596 597 G4double defaultRidgeXPosition= -1270.0*mm; 598 const G4double XBase = 0.5 *mm; 599 const G4double YBase =6.03*cm; 600 const G4double ZBase =6.03*cm; 601 602 603 604 SolidRidgeBase = new G4Box("BaseRidgeFilter", 605 XBase, 606 YBase/2, 607 ZBase/2); 608 609 LogicRidgeBase = new G4LogicalVolume(SolidRidgeBase, 610 PMMANist, 611 "BaseRidgeFilter"); 612 613 PhysiRidgeFilterBase = new G4PVPlacement(0, 614 G4ThreeVector( 615 -1270*mm, 616 0., 617 0.), 618 "BaseRidgeFilter", 619 LogicRidgeBase, 620 physicalTreatmentRoom, 621 false, 622 0); 623 624 LogicRidgeBase->SetVisAttributes(red); 625 626 627 SolidRidgeMother = new G4Box("MotherRidgeSOL",1.7*mm/2,1.7*mm/2,4.72*mm/2); 628 LogicRidgeMother = new G4LogicalVolume(SolidRidgeMother,airNist,"MotherRidge"); 629 630 G4Trd* trapp1=new G4Trd("Trapp1SOL",1.7*mm/2,1.68*mm/2,1.7*mm/2,1.68*mm/2,0.22*mm/2); 631 G4LogicalVolume* LogicTrapp1=new G4LogicalVolume(trapp1,PMMANist,"Trapp1"); 632 PhysiTrapp1=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22/2)*mm),LogicTrapp1,"Trapp1",LogicRidgeMother,false,0); 633 634 G4Trd* trapp2=new G4Trd("Trapp2SOL",1.68*mm/2,1.64*mm/2,1.68*mm/2,1.64*mm/2,0.08*mm/2); 635 G4LogicalVolume* LogicTrapp2=new G4LogicalVolume(trapp2,PMMANist,"Trapp2"); 636 PhysiTrapp2=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08/2)*mm),LogicTrapp2,"Trapp2",LogicRidgeMother,false,0); 637 638 G4Trd* trapp3=new G4Trd("Trapp3SOL",1.64*mm/2,1.58*mm/2,1.64*mm/2,1.58*mm/2,0.14*mm/2); 639 G4LogicalVolume* LogicTrapp3=new G4LogicalVolume(trapp3,PMMANist,"Trapp3"); 640 PhysiTrapp3=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14/2)*mm),LogicTrapp3,"Trapp3",LogicRidgeMother,false,0); 641 642 G4Trd* trapp4=new G4Trd("Trapp4SOL",1.58*mm/2,1.50*mm/2,1.58*mm/2,1.50*mm/2,0.22*mm/2); 643 G4LogicalVolume* LogicTrapp4=new G4LogicalVolume(trapp4,PMMANist,"Trapp4"); 644 PhysiTrapp4=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22/2)*mm),LogicTrapp4,"Trapp4",LogicRidgeMother,false,0); 645 646 G4Trd* trapp5=new G4Trd("Trapp5SOL",1.50*mm/2,1.40*mm/2,1.50*mm/2,1.40*mm/2,0.31*mm/2); 647 G4LogicalVolume* LogicTrapp5=new G4LogicalVolume(trapp5,PMMANist,"Trapp5"); 648 PhysiTrapp5=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31/2)*mm),LogicTrapp5,"Trapp5",LogicRidgeMother,false,0); 649 650 G4Trd* trapp6=new G4Trd("Trapp6SOL",1.40*mm/2,1.26*mm/2,1.40*mm/2,1.26*mm/2,0.48*mm/2); 651 G4LogicalVolume* LogicTrapp6=new G4LogicalVolume(trapp6,PMMANist,"Trapp6"); 652 PhysiTrapp6=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48/2)*mm),LogicTrapp6,"Trapp6",LogicRidgeMother,false,0); 653 654 G4Trd* trapp7=new G4Trd("Trapp7SOL",1.26*mm/2,0.94*mm/2,1.26*mm/2,0.94*mm/2,1.20*mm/2); 655 G4LogicalVolume* LogicTrapp7=new G4LogicalVolume(trapp7,PMMANist,"Trapp7"); 656 PhysiTrapp7=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20/2)*mm),LogicTrapp7,"Trapp7",LogicRidgeMother,false,0); 657 658 G4Trd* trapp8=new G4Trd("Trapp8SOL",0.94*mm/2,0.78*mm/2,0.94*mm/2,0.78*mm/2,0.57*mm/2); 659 G4LogicalVolume* LogicTrapp8=new G4LogicalVolume(trapp8,PMMANist,"Trapp8"); 660 PhysiTrapp8=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57/2)*mm),LogicTrapp8,"Trapp8",LogicRidgeMother,false,0); 661 662 G4Trd* trapp9=new G4Trd("Trapp9SOL",0.78*mm/2,0.66*mm/2,0.78*mm/2,0.66*mm/2,0.39*mm/2); 663 G4LogicalVolume* LogicTrapp9=new G4LogicalVolume(trapp9,PMMANist,"Trapp9"); 664 PhysiTrapp9=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39/2)*mm),LogicTrapp9,"Trapp9",LogicRidgeMother,false,0); 665 666 G4Trd* trapp10=new G4Trd("Trapp10SOL",0.66*mm/2,0.56*mm/2,0.66*mm/2,0.56*mm/2,0.29*mm/2); 667 G4LogicalVolume* LogicTrapp10=new G4LogicalVolume(trapp10,PMMANist,"Trapp10"); 668 PhysiTrapp10=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29/2)*mm),LogicTrapp10,"Trapp10",LogicRidgeMother,false,0); 669 670 G4Trd* trapp11=new G4Trd("Trapp11SOL",0.56*mm/2,0.46*mm/2,0.56*mm/2,0.46*mm/2,0.25*mm/2); 671 G4LogicalVolume* LogicTrapp11=new G4LogicalVolume(trapp11,PMMANist,"Trapp11"); 672 PhysiTrapp11=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29+0.25/2)*mm),LogicTrapp11,"Trapp11",LogicRidgeMother,false,0); 673 674 G4Trd* trapp12=new G4Trd("Trapp12SOL",0.46*mm/2,0.38*mm/2,0.46*mm/2,0.38*mm/2,0.17*mm/2); 675 G4LogicalVolume* LogicTrapp12=new G4LogicalVolume(trapp12,PMMANist,"Trapp12"); 676 PhysiTrapp12=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29+0.25+0.17/2)*mm),LogicTrapp12,"Trapp12",LogicRidgeMother,false,0); 677 678 G4Trd* trapp13=new G4Trd("Trapp13SOL",0.38*mm/2,0.30*mm/2,0.38*mm/2,0.30*mm/2,0.14*mm/2); 679 G4LogicalVolume* LogicTrapp13=new G4LogicalVolume(trapp13,PMMANist,"Trapp13"); 680 PhysiTrapp13=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29+0.25+0.17+0.14/2)*mm),LogicTrapp13,"Trapp13",LogicRidgeMother,false,0); 681 682 G4Trd* trapp14=new G4Trd("Trapp14SOL",0.30*mm/2,0.24*mm/2,0.30*mm/2,0.24*mm/2,0.09*mm/2); 683 G4LogicalVolume* LogicTrapp14=new G4LogicalVolume(trapp14,PMMANist,"Trapp13"); 684 PhysiTrapp14=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29+0.25+0.17+0.14+0.09/2)*mm),LogicTrapp14,"Trapp14",LogicRidgeMother,false,0); 685 686 G4Trd* trapp15=new G4Trd("Trapp14SOL",0.24*mm/2,0.18*mm/2,0.24*mm/2,0.18*mm/2,0.07*mm/2); 687 G4LogicalVolume* LogicTrapp15=new G4LogicalVolume(trapp15,PMMANist,"Trapp15"); 688 PhysiTrapp15=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29+0.25+0.17+0.14+0.09+0.07/2)*mm),LogicTrapp15,"Trapp15",LogicRidgeMother,false,0); 689 690 G4Trd* trapp16=new G4Trd("Trapp16SOL",0.18*mm/2,0.12*mm/2,0.18*mm/2,0.12*mm/2,0.10*mm/2); 691 G4LogicalVolume* LogicTrapp16=new G4LogicalVolume(trapp16,PMMANist,"Trapp16"); 692 PhysiTrapp16=new G4PVPlacement(0,G4ThreeVector(0.,0.,(-4.72/2+0.22+0.08+0.14+0.22+0.31+0.48+1.20+0.57+0.39+0.29+0.25+0.17+0.14+0.09+0.07+0.10/2)*mm),LogicTrapp16,"Trapp16",LogicRidgeMother,false,0); 693 694 LogicTrapp1->SetVisAttributes(green); 695 LogicTrapp2->SetVisAttributes(green); 696 LogicTrapp3->SetVisAttributes(green); 697 LogicTrapp4->SetVisAttributes(green); 698 LogicTrapp5->SetVisAttributes(green); 699 LogicTrapp6->SetVisAttributes(green); 700 LogicTrapp7->SetVisAttributes(green); 701 LogicTrapp8->SetVisAttributes(green); 702 LogicTrapp9->SetVisAttributes(green); 703 LogicTrapp10->SetVisAttributes(green); 704 LogicTrapp11->SetVisAttributes(green); 705 LogicTrapp12->SetVisAttributes(green); 706 LogicTrapp13->SetVisAttributes(green); 707 LogicTrapp14->SetVisAttributes(green); 708 LogicTrapp15->SetVisAttributes(green); 709 LogicTrapp16->SetVisAttributes(green); 710 G4VisAttributes* visAttr = new G4VisAttributes(); 711 visAttr->SetVisibility(false); 712 LogicRidgeMother->SetVisAttributes(visAttr); 713 714 715 716 717 718 719 G4int numberOfLayers = 30; 720 G4double minZ = 30.15*mm-0.3*mm-1.70/2*mm; 721 G4double minY = 30.15*mm-0.3*mm-1.70/2*mm; 722 G4double sum_space = 0.3*mm+1.70*mm; 723 724 725 std::vector<G4ThreeVector> singleTrapPositions; 726 727 for (int i = 0; i < numberOfLayers; i++) 728 { 729 for (int j = 0; j < numberOfLayers; j++) 730 { 731 singleTrapPositions.push_back({defaultRidgeXPosition-4.72/2*mm-0.5*mm, 732 minY - i*sum_space, 733 minZ - j*sum_space, 734 735 }); 736 } 737 } 738 739 G4double ti = - 90. *deg; 740 G4RotationMatrix rt; 741 rt.rotateY(ti); 742 743 744 G4int peaks = numberOfLayers*numberOfLayers; 745 for (int i = 0; i < peaks; i++) 746 { 747 748 std::ostringstream tName;tName << "singleTrap" << i; 749 new G4PVPlacement(G4Transform3D(rt, 750 singleTrapPositions[i]), 751 LogicRidgeMother, 752 "tName.str()", 753 logicTreatmentRoom, 754 0, 755 i); 756 757 } 758 759 } 760 761 762 763 764 765 766 ///////////////////////////////////////////////////////////////////////////// 767 void PassiveCarbonBeamLine::HadrontherapyPMMACollimator() 768 { 769 770 // ----------------------// 771 // PMMA COLLIMATOR // 772 // ----------------------// 773 PMMACollimatorSupportXSize = 25.0 *mm; 774 PMMACollimatorSupportYSize = 200. *mm; 775 PMMACollimatorSupportZSize = 200. *mm; 776 777 778 PMMACollimatorXPosition = -1257.0 *mm; 779 780 G4double phi = 90. *deg; 781 G4RotationMatrix rm; 782 rm.rotateY(phi); 783 784 solidPMMACollimatorSupport = new G4Box("PMMACollimatorSupport", 785 PMMACollimatorSupportXSize/2, 786 PMMACollimatorSupportYSize/2, 787 PMMACollimatorSupportZSize/2); 788 789 logicPMMACollimatorSupport = new G4LogicalVolume(solidPMMACollimatorSupport, 790 PMMACollimatorMaterial, 791 "PMMACollimatorSupport"); 792 793 physiPMMACollimatorSupport = new G4PVPlacement(0, G4ThreeVector(PMMACollimatorXPosition,0., 0.), 794 "PMMACollimatorSupport", 795 logicPMMACollimatorSupport, 796 physicalTreatmentRoom, 797 false, 798 0); 799 800 801 yellow = new G4VisAttributes(G4Colour(1., 1., 0. )); 802 yellow-> SetVisibility(true); 803 yellow-> SetForceWireframe(true); 804 805 logicPMMACollimatorSupport -> SetVisAttributes(yellow); 806 807 // ----------------------// 808 // PMMA COLLIMATOR // 809 //-----------------------// 810 innerRadiusPMMACollimator= 4.5 *cm; 811 outerRadiusPMMACollimator = 4.6 *cm; 812 hightPMMACollimator = PMMACollimatorSupportXSize; 813 startAnglePMMACollimator = 0.*deg; 814 spanningAnglePMMACollimator = 360.*deg; 815 816 817 solidPMMACollimator = new G4Tubs("PMMACollimator", 818 innerRadiusPMMACollimator, 819 outerRadiusPMMACollimator, 820 hightPMMACollimator/2, 821 startAnglePMMACollimator, 822 spanningAnglePMMACollimator); 823 824 logicPMMACollimator = new G4LogicalVolume(solidPMMACollimator, 825 PMMACollimatorMaterial, 826 "PMMACollimator", 827 0, 828 0, 829 0); 830 831 physiPMMACollimator = new G4PVPlacement(G4Transform3D(rm, 832 G4ThreeVector(0,0.,0.)), 833 "PMMACollimator", 834 logicPMMACollimator, 835 physiPMMACollimatorSupport, 836 false, 837 0); 838 839 logicPMMACollimator -> SetVisAttributes(yellow); 840 841 842 843 844 } 845 ///////////////////////////////////////////////////////////////////////////// 846 void PassiveCarbonBeamLine::HadrontherapyBeamMonitoring() 847 { 848 // ---------------------------- 849 // MONITOR CHAMBER 850 // ---------------------------- 851 // A monitor chamber is a free-air ionisation chamber 852 // able to measure do carbon fluence during the treatment. 853 // Here its responce is not simulated in terms of produced 854 // charge but only the energy losses are taked into account. 855 // Each chamber consist of 9 mm of air in a box 856 // that has two layers one of kapton and one 857 // of copper 858 859 monitor1XSize = 4.525022*mm; 860 monitor2XSize = 0.000011*mm; 861 monitor3XSize = 4.5*mm; 862 monitorYSize = 10.*cm; 863 monitorZSize = 10.*cm; 864 monitor1XPosition = -1239.974978 *mm; 865 monitor2XPosition = -4.500011*mm; 866 monitor4XPosition = 4.500011*mm; 867 868 solidFirstMonitorLayer1 = new G4Box("FirstMonitorLayer1", 869 monitor1XSize, 870 monitorYSize, 871 monitorZSize); 872 873 logicFirstMonitorLayer1 = new G4LogicalVolume(solidFirstMonitorLayer1, 874 layer1MonitorChamberMaterial, 875 "FirstMonitorLayer1"); 876 877 physiFirstMonitorLayer1 = new G4PVPlacement(0, 878 G4ThreeVector(monitor1XPosition,0.*cm,0.*cm), 879 "FirstMonitorLayer1", 880 logicFirstMonitorLayer1, 881 physicalTreatmentRoom, 882 false, 883 0); 884 885 solidFirstMonitorLayer2 = new G4Box("FirstMonitorLayer2", 886 monitor2XSize, 887 monitorYSize, 888 monitorZSize); 889 890 logicFirstMonitorLayer2 = new G4LogicalVolume(solidFirstMonitorLayer2, 891 layer2MonitorChamberMaterial, 892 "FirstMonitorLayer2"); 893 894 physiFirstMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector(monitor2XPosition,0.*cm,0.*cm), 895 "FirstMonitorLayer2", 896 logicFirstMonitorLayer2, 897 physiFirstMonitorLayer1, 898 false, 899 0); 900 901 solidFirstMonitorLayer3 = new G4Box("FirstMonitorLayer3", 902 monitor3XSize, 903 monitorYSize, 904 monitorZSize); 905 906 logicFirstMonitorLayer3 = new G4LogicalVolume(solidFirstMonitorLayer3, 907 layer3MonitorChamberMaterial, 908 "FirstMonitorLayer3"); 909 910 physiFirstMonitorLayer3 = new G4PVPlacement(0, 911 G4ThreeVector(0.*mm,0.*cm,0.*cm), 912 "MonitorLayer3", 913 logicFirstMonitorLayer3, 914 physiFirstMonitorLayer1, 915 false, 916 0); 917 918 solidFirstMonitorLayer4 = new G4Box("FirstMonitorLayer4", 919 monitor2XSize, 920 monitorYSize, 921 monitorZSize); 922 923 logicFirstMonitorLayer4 = new G4LogicalVolume(solidFirstMonitorLayer4, 924 layer4MonitorChamberMaterial, 925 "FirstMonitorLayer4"); 926 927 physiFirstMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm), 928 "FirstMonitorLayer4", 929 logicFirstMonitorLayer4, 930 physiFirstMonitorLayer1, false, 0); 931 932 logicFirstMonitorLayer3 -> SetVisAttributes(white); 933 934 } 935 936 ///////////////////////////////////////////////////////////////////////////// 937 ///////////////////////////////////////////////////////////////////////////// 938 void PassiveCarbonBeamLine::HadrontherapyBeamNozzle() 939 { 940 // ------------------------------// 941 // THE FINAL TUBE AND COLLIMATOR // 942 //-------------------------------// 943 // The last part of the transport beam line consists of 944 // a 59 mm thick PMMA slab (to stop all the diffused radiation), a 370 mm brass tube 945 // (to well collimate the proton beam) and a final collimator with 25 mm diameter 946 // aperture (that provide the final trasversal shape of the beam) 947 948 // -------------------// 949 // PMMA SUPPORT // 950 // -------------------// 951 952 nozzleSupportXSize = 50 *mm; 953 nozzleSupportYSize = 360. *mm; 954 nozzleSupportZSize = 360. *mm; 955 nozzleSupportXPosition = -423.0 *mm; 956 957 G4double phi = 90. *deg; 958 // Matrix definition for a 90 deg rotation. Also used for other volumes 959 G4RotationMatrix rm; 960 rm.rotateY(phi); 961 962 solidNozzleSupport = new G4Box("NozzleSupport", 963 nozzleSupportXSize/2, 964 nozzleSupportYSize/2, 965 nozzleSupportZSize/2); 966 967 logicNozzleSupport = new G4LogicalVolume(solidNozzleSupport, 968 nozzleSupportMaterial, 969 "NozzleSupport"); 970 971 physiNozzleSupport = new G4PVPlacement(0, G4ThreeVector(nozzleSupportXPosition,0., 0.), 972 "NozzleSupport", 973 logicNozzleSupport, 974 physicalTreatmentRoom, 975 false, 976 0); 977 978 yellow = new G4VisAttributes(G4Colour(1., 1., 0. )); 979 yellow-> SetVisibility(true); 980 yellow-> SetForceWireframe(true); 981 982 //------------------------------------// 983 // HOLE IN THE SUPPORT // 984 //------------------------------------// 985 innerRadiusHoleNozzleSupport = 0.*mm; 986 outerRadiusHoleNozzleSupport = 22.0*mm; 987 hightHoleNozzleSupport = nozzleSupportXSize; 988 startAngleHoleNozzleSupport = 0.*deg; 989 spanningAngleHoleNozzleSupport = 360.*deg; 990 991 solidHoleNozzleSupport = new G4Tubs("HoleNozzleSupport", 992 innerRadiusHoleNozzleSupport, 993 outerRadiusHoleNozzleSupport, 994 hightHoleNozzleSupport/2, 995 startAngleHoleNozzleSupport, 996 spanningAngleHoleNozzleSupport); 997 998 logicHoleNozzleSupport = new G4LogicalVolume(solidHoleNozzleSupport, 999 holeNozzleSupportMaterial, 1000 "HoleNozzleSupport", 1001 0, 1002 0, 1003 0); 1004 1005 physiHoleNozzleSupport = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()), 1006 "HoleNozzleSupport", 1007 logicHoleNozzleSupport, 1008 physiNozzleSupport, 1009 false, 0); 1010 1011 1012 // --------------------------------------// 1013 // BRASS TUBE 3 (beam line side) // 1014 // -------------------------------------// 1015 innerRadiusBrassTube3 = 13.5 *mm; 1016 outerRadiusBrassTube3 = 21.5 *mm; 1017 hightBrassTube3 = 20.0 *mm; 1018 startAngleBrassTube3 = 0.*deg; 1019 spanningAngleBrassTube3 = 360.*deg; 1020 brassTube3XPosition = -458.0 *mm; 1021 1022 solidBrassTube3 = new G4Tubs("BrassTube3", 1023 innerRadiusBrassTube3, 1024 outerRadiusBrassTube3, 1025 hightBrassTube3/2, 1026 startAngleBrassTube3, 1027 spanningAngleBrassTube3); 1028 1029 logicBrassTube3 = new G4LogicalVolume(solidBrassTube3, 1030 brassTube3Material, 1031 "BrassTube3", 1032 0, 0, 0); 1033 1034 physiBrassTube3 = new G4PVPlacement(G4Transform3D(rm, 1035 G4ThreeVector(brassTube3XPosition, 1036 0., 1037 0.)), 1038 "BrassTube3", 1039 logicBrassTube3, 1040 physicalTreatmentRoom, 1041 false, 1042 0); 1043 1044 logicBrassTube3 -> SetVisAttributes(darkOrange3); 1045 1046 // ----------------------------------------------// 1047 // BRASS TUBE 2 (inside the PMMA support) // 1048 // ----------------------------------------------// 1049 1050 innerRadiusBrassTube2 = 13.5 *mm; 1051 outerRadiusBrassTube2 = 21.5 *mm; 1052 hightBrassTube2 = nozzleSupportXSize; 1053 startAngleBrassTube2 = 0.*deg; 1054 spanningAngleBrassTube2 = 360.*deg; 1055 1056 1057 solidBrassTube2 = new G4Tubs("BrassTube2", 1058 innerRadiusBrassTube2, 1059 outerRadiusBrassTube2, 1060 hightBrassTube2/2, 1061 startAngleBrassTube2, 1062 spanningAngleBrassTube2); 1063 1064 logicBrassTube2 = new G4LogicalVolume(solidBrassTube2, 1065 brassTube2Material, 1066 "BrassTube2", 1067 0, 0, 0); 1068 physiBrassTube2 = new G4PVPlacement(0, 1069 G4ThreeVector(0,0.,0.), 1070 logicBrassTube2, 1071 "BrassTube2", 1072 logicHoleNozzleSupport, 1073 false, 1074 0); 1075 1076 logicBrassTube2 -> SetVisAttributes(darkOrange3); 1077 1078 // ---------------------------------// 1079 // BRASS TUBE 1 (phantom side) // 1080 // ---------------------------------// 1081 innerRadiusBrassTube= 18.*mm; 1082 outerRadiusBrassTube = 21.5 *mm; 1083 hightBrassTube = 208.0 *mm; 1084 startAngleBrassTube = 0.*deg; 1085 spanningAngleBrassTube = 360.*deg; 1086 brassTubeXPosition = -294 *mm; 1087 solidBrassTube = new G4Tubs("BrassTube", 1088 innerRadiusBrassTube, 1089 outerRadiusBrassTube, 1090 hightBrassTube/2, 1091 startAngleBrassTube, 1092 spanningAngleBrassTube); 1093 1094 logicBrassTube = new G4LogicalVolume(solidBrassTube, 1095 brassTubeMaterial, 1096 "BrassTube", 1097 0, 0, 0); 1098 1099 physiBrassTube = new G4PVPlacement(G4Transform3D(rm, 1100 G4ThreeVector(brassTubeXPosition, 1101 0., 1102 0.)), 1103 "BrassTube", 1104 logicBrassTube, 1105 physicalTreatmentRoom, 1106 false, 1107 0); 1108 1109 logicBrassTube -> SetVisAttributes(darkOrange3); 1110 } 1111 1112 ///////////////////////////////////////////////////////////////////////////// 1113 void PassiveCarbonBeamLine::HadrontherapyBeamFinalCollimator() 1114 { 1115 // -----------------------// 1116 // FINAL COLLIMATOR // 1117 //------------------------// 1118 outerRadiusFinalCollimator = 21.5*mm; 1119 hightFinalCollimator = 7.0 *mm; 1120 startAngleFinalCollimator = 0.*deg; 1121 spanningAngleFinalCollimator = 360.*deg; 1122 finalCollimatorXPosition = -186.5 *mm; 1123 1124 G4double phi = 90. *deg; 1125 1126 // Matrix definition for a 90 deg rotation. Also used for other volumes 1127 G4RotationMatrix rm; 1128 rm.rotateY(phi); 1129 1130 solidFinalCollimator = new G4Tubs("FinalCollimator", 1131 innerRadiusFinalCollimator, 1132 outerRadiusFinalCollimator, 1133 hightFinalCollimator/2, 1134 startAngleFinalCollimator, 1135 spanningAngleFinalCollimator); 1136 1137 logicFinalCollimator = new G4LogicalVolume(solidFinalCollimator, 1138 finalCollimatorMaterial, 1139 "FinalCollimator", 1140 0, 1141 0, 1142 0); 1143 1144 physiFinalCollimator = new G4PVPlacement(G4Transform3D(rm, 1145 G4ThreeVector(finalCollimatorXPosition,0.,0.)), 1146 "FinalCollimator", 1147 logicFinalCollimator, 1148 physicalTreatmentRoom, 1149 false, 1150 0); 1151 1152 logicFinalCollimator -> SetVisAttributes(yellow); 1153 } 1154 1155 1156 ///////////////////////////////////////////////////////////////////////////// 1157 /////////////////////////// MESSENGER /////////////////////////////////////// 1158 ///////////////////////////////////////////////////////////////////////////// 1159 1160 void PassiveCarbonBeamLine::SetRippleFilterXPosition(G4double value) 1161 { 1162 PhysiRippleFilter -> SetTranslation(G4ThreeVector(value, 0., 0.)); 1163 G4RunManager::GetRunManager() -> GeometryHasBeenModified(); 1164 G4cout << "The Ripple Filter is translated to"<< value/mm <<"mm along the X axis" <<G4endl; 1165 } 1166 1167 1168 ///////////////////////////////////////////////////////////////////////////// 1169 void PassiveCarbonBeamLine::SetInnerRadiusFinalCollimator(G4double value) 1170 { 1171 solidFinalCollimator -> SetInnerRadius(value); 1172 G4RunManager::GetRunManager() -> GeometryHasBeenModified(); 1173 G4cout<<"Inner Radius of the final collimator is (mm):" 1174 << solidFinalCollimator -> GetInnerRadius()/mm 1175 << G4endl; 1176 } 1177 1178 ///////////////////////////////////////////////////////////////////////////// 1179 void PassiveCarbonBeamLine::SetRippleFilterMaterial(G4String materialChoice) 1180 { 1181 if (G4Material* RFMaterial = G4NistManager::Instance()->FindOrBuildMaterial(materialChoice, false) ) 1182 { 1183 if (RFMaterial) 1184 { 1185 rippleFilterMaterial = RFMaterial; 1186 LogicRippleFilter -> SetMaterial(RFMaterial); 1187 LogicRippleFilterBase -> SetMaterial(RFMaterial); 1188 LogicRippleFilterTrd -> SetMaterial(RFMaterial); 1189 G4cout << "The material of the Ripple Filter has been changed to " << materialChoice << G4endl; 1190 } 1191 } 1192 else 1193 { 1194 G4cout << "WARNING: material \"" << materialChoice << "\" doesn't exist in NIST elements/materials" 1195 " table [located in $G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc]" << G4endl; 1196 G4cout << "Use command \"/parameter/nist\" to see full materials list!" << G4endl; 1197 } 1198 } 1199 1200 1201 1202