Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions
PassiveProtonBeamLine Class Reference

#include <PassiveProtonBeamLine.hh>

Inheritance diagram for PassiveProtonBeamLine:
G4VUserDetectorConstruction

Public Member Functions

 PassiveProtonBeamLine ()
 
 ~PassiveProtonBeamLine ()
 
G4VPhysicalVolumeConstruct ()
 
void HadrontherapyBeamLineSupport ()
 
void HadrontherapyBeamScatteringFoils ()
 
void HadrontherapyRangeShifter ()
 
void HadrontherapyBeamCollimators ()
 
void HadrontherapyBeamMonitoring ()
 
void HadrontherapyMOPIDetector ()
 
void HadrontherapyBeamNozzle ()
 
void HadrontherapyBeamFinalCollimator ()
 
void SetRangeShifterXPosition (G4double value)
 
void SetRangeShifterXSize (G4double halfSize)
 
void SetFirstScatteringFoilXSize (G4double)
 
void SetSecondScatteringFoilXSize (G4double)
 
void SetOuterRadiusStopper (G4double)
 
void SetInnerRadiusFinalCollimator (G4double)
 
void SetRSMaterial (G4String)
 
void SetModulatorAngle (G4double angle)
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void ConstructSDandField ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Static Public Member Functions

static PassiveProtonBeamLineGetInstance ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VUserDetectorConstruction
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 

Detailed Description

Definition at line 51 of file PassiveProtonBeamLine.hh.

Constructor & Destructor Documentation

PassiveProtonBeamLine::PassiveProtonBeamLine ( )

Definition at line 55 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, and G4VUserDetectorConstruction::RegisterParallelWorld().

55  :
56  modulator(0), physicalTreatmentRoom(0),hadrontherapyDetectorConstruction(0),
57  physiBeamLineSupport(0), physiBeamLineCover(0), physiBeamLineCover2(0),
58  firstScatteringFoil(0), physiFirstScatteringFoil(0), physiKaptonWindow(0),
59  solidStopper(0), physiStopper(0),
60  secondScatteringFoil(0), physiSecondScatteringFoil(0),
61  physiFirstCollimator(0), solidRangeShifterBox(0), logicRangeShifterBox(0),
62  physiRangeShifterBox(0), physiSecondCollimator(0),
63  physiFirstCollimatorModulatorBox(0),
64  physiHoleFirstCollimatorModulatorBox(0),
65  physiSecondCollimatorModulatorBox(0),
66  physiHoleSecondCollimatorModulatorBox(0),
67  physiMOPIMotherVolume(0),
68  physiFirstMonitorLayer1(0), physiFirstMonitorLayer2(0),
69  physiFirstMonitorLayer3(0), physiFirstMonitorLayer4(0),
70  physiSecondMonitorLayer1(0), physiSecondMonitorLayer2(0),
71  physiSecondMonitorLayer3(0), physiSecondMonitorLayer4(0),
72  physiNozzleSupport(0), //physiHoleNozzleSupport(0),
73  physiBrassTube(0),
74  solidFinalCollimator(0),
75  physiFinalCollimator(0)
76 {
77 
78  // Messenger to change parameters of the passiveProtonBeamLine geometry
79  passiveMessenger = new PassiveProtonBeamLineMessenger(this);
80 
81 //***************************** PW ***************************************
82 
83  static G4String ROGeometryName = "DetectorROGeometry";
84  RO = new HadrontherapyDetectorROGeometry(ROGeometryName);
85 
86 
87 
88  G4cout << "Going to register Parallel world...";
90  G4cout << "... done" << G4endl;
91 //***************************** PW ***************************************
92 
93 }
void RegisterParallelWorld(G4VUserParallelWorld *)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
PassiveProtonBeamLine::~PassiveProtonBeamLine ( )

Definition at line 95 of file PassiveProtonBeamLine.cc.

96 {
97  delete passiveMessenger;
98  delete hadrontherapyDetectorConstruction;
99 }

Member Function Documentation

G4VPhysicalVolume * PassiveProtonBeamLine::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 101 of file PassiveProtonBeamLine.cc.

References HadrontherapyDetectorConstruction::GetDetectorToWorldPosition(), and HadrontherapyDetectorConstruction::InitializeDetectorROGeometry().

102 {
103  // Sets default geometry and materials
104  SetDefaultDimensions();
105 
106  // Construct the whole Passive Beam Line
107  ConstructPassiveProtonBeamLine();
108 
109 
110 //***************************** PW ***************************************
111  if (!hadrontherapyDetectorConstruction)
112 
113 //***************************** PW ***************************************
114 
115 
116  // HadrontherapyDetectorConstruction builds ONLY the phantom and the detector with its associated ROGeometry
117  hadrontherapyDetectorConstruction = new HadrontherapyDetectorConstruction(physicalTreatmentRoom);
118 
119 
120  //***************************** PW ***************************************
121 
122  hadrontherapyDetectorConstruction->InitializeDetectorROGeometry(RO,hadrontherapyDetectorConstruction->GetDetectorToWorldPosition());
123 
124 //***************************** PW ***************************************
125 
126 
127 
128  return physicalTreatmentRoom;
129 }
void InitializeDetectorROGeometry(HadrontherapyDetectorROGeometry *, G4ThreeVector detectorToWorldPosition)
static PassiveProtonBeamLine* PassiveProtonBeamLine::GetInstance ( )
static
void PassiveProtonBeamLine::HadrontherapyBeamCollimators ( )

Definition at line 820 of file PassiveProtonBeamLine.cc.

References python.hepunit::deg, python.hepunit::mm, and CLHEP::HepRotation::rotateY().

821 {
822  // -----------------//
823  // FIRST COLLIMATOR //
824  // -----------------//
825  // It is a slab of PMMA with an hole in its center
826  const G4double firstCollimatorXSize = 20.*mm;
827  const G4double firstCollimatorYSize = 100.*mm;
828  const G4double firstCollimatorZSize = 100.*mm;
829 
830  const G4double firstCollimatorXPosition = -2673.00*mm;
831  const G4double firstCollimatorYPosition = 0.*mm;
832  const G4double firstCollimatorZPosition = 0.*mm;
833 
834  G4Box* solidFirstCollimator = new G4Box("FirstCollimator",
835  firstCollimatorXSize,
836  firstCollimatorYSize,
837  firstCollimatorZSize);
838 
839  G4LogicalVolume* logicFirstCollimator = new G4LogicalVolume(solidFirstCollimator,
840  firstCollimatorMaterial,
841  "FirstCollimator");
842 
843  physiFirstCollimator = new G4PVPlacement(0, G4ThreeVector(firstCollimatorXPosition,
844  firstCollimatorYPosition,
845  firstCollimatorZPosition),
846  "FirstCollimator",
847  logicFirstCollimator,
848  physicalTreatmentRoom,
849  false,
850  0);
851  // ----------------------------//
852  // Hole of the first collimator//
853  //-----------------------------//
854  G4double innerRadiusHoleFirstCollimator = 0.*mm;
855  G4double outerRadiusHoleFirstCollimator = 15.*mm;
856  G4double hightHoleFirstCollimator = 20.*mm;
857  G4double startAngleHoleFirstCollimator = 0.*deg;
858  G4double spanningAngleHoleFirstCollimator = 360.*deg;
859 
860  G4Tubs* solidHoleFirstCollimator = new G4Tubs("HoleFirstCollimator",
861  innerRadiusHoleFirstCollimator,
862  outerRadiusHoleFirstCollimator,
863  hightHoleFirstCollimator,
864  startAngleHoleFirstCollimator,
865  spanningAngleHoleFirstCollimator);
866 
867  G4LogicalVolume* logicHoleFirstCollimator = new G4LogicalVolume(solidHoleFirstCollimator,
868  holeFirstCollimatorMaterial,
869  "HoleFirstCollimator",
870  0, 0, 0);
871  G4double phi = 90. *deg;
872  // Matrix definition for a 90 deg rotation. Also used for other volumes
873  G4RotationMatrix rm;
874  rm.rotateY(phi);
875 
876  physiHoleFirstCollimator = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
877  "HoleFirstCollimator",
878  logicHoleFirstCollimator,
879  physiFirstCollimator,
880  false,
881  0);
882  // ------------------//
883  // SECOND COLLIMATOR //
884  //-------------------//
885  // It is a slab of PMMA with an hole in its center
886  const G4double secondCollimatorXPosition = -1900.00*mm;
887  const G4double secondCollimatorYPosition = 0*mm;
888  const G4double secondCollimatorZPosition = 0*mm;
889 
890  physiSecondCollimator = new G4PVPlacement(0, G4ThreeVector(secondCollimatorXPosition,
891  secondCollimatorYPosition,
892  secondCollimatorZPosition),
893  "SecondCollimator",
894  logicFirstCollimator,
895  physicalTreatmentRoom,
896  false,
897  0);
898 
899  // ------------------------------//
900  // Hole of the second collimator //
901  // ------------------------------//
902  physiHoleSecondCollimator = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
903  "HoleSecondCollimator",
904  logicHoleFirstCollimator,
905  physiSecondCollimator,
906  false,
907  0);
908 
909  // --------------------------------------//
910  // FIRST SIDE OF THE MODULATOR BOX //
911  // --------------------------------------//
912  // The modulator box is an aluminum box in which
913  // the range shifter and the energy modulator are located
914  // In this example only the entrance and exit
915  // faces of the box are simulated.
916  // Each face is an aluminum slab with an hole in its center
917 
918  const G4double firstCollimatorModulatorXSize = 10.*mm;
919  const G4double firstCollimatorModulatorYSize = 200.*mm;
920  const G4double firstCollimatorModulatorZSize = 200.*mm;
921 
922  const G4double firstCollimatorModulatorXPosition = -2523.00*mm;
923  const G4double firstCollimatorModulatorYPosition = 0.*mm;
924  const G4double firstCollimatorModulatorZPosition = 0.*mm;
925 
926  G4Box* solidFirstCollimatorModulatorBox = new G4Box("FirstCollimatorModulatorBox",
927  firstCollimatorModulatorXSize,
928  firstCollimatorModulatorYSize,
929  firstCollimatorModulatorZSize);
930 
931  G4LogicalVolume* logicFirstCollimatorModulatorBox = new G4LogicalVolume(solidFirstCollimatorModulatorBox,
932  modulatorBoxMaterial,
933  "FirstCollimatorModulatorBox");
934 
935  physiFirstCollimatorModulatorBox = new G4PVPlacement(0, G4ThreeVector(firstCollimatorModulatorXPosition,
936  firstCollimatorModulatorYPosition,
937  firstCollimatorModulatorZPosition),
938  "FirstCollimatorModulatorBox",
939  logicFirstCollimatorModulatorBox,
940  physicalTreatmentRoom, false, 0);
941 
942  // ----------------------------------------------------//
943  // Hole of the first collimator of the modulator box //
944  // ----------------------------------------------------//
945  const G4double innerRadiusHoleFirstCollimatorModulatorBox = 0.*mm;
946  const G4double outerRadiusHoleFirstCollimatorModulatorBox = 31.*mm;
947  const G4double hightHoleFirstCollimatorModulatorBox = 10.*mm;
948  const G4double startAngleHoleFirstCollimatorModulatorBox = 0.*deg;
949  const G4double spanningAngleHoleFirstCollimatorModulatorBox = 360.*deg;
950 
951  G4Tubs* solidHoleFirstCollimatorModulatorBox = new G4Tubs("HoleFirstCollimatorModulatorBox",
952  innerRadiusHoleFirstCollimatorModulatorBox,
953  outerRadiusHoleFirstCollimatorModulatorBox,
954  hightHoleFirstCollimatorModulatorBox ,
955  startAngleHoleFirstCollimatorModulatorBox,
956  spanningAngleHoleFirstCollimatorModulatorBox);
957 
958  G4LogicalVolume* logicHoleFirstCollimatorModulatorBox = new G4LogicalVolume(solidHoleFirstCollimatorModulatorBox,
959  holeModulatorBoxMaterial,
960  "HoleFirstCollimatorModulatorBox",
961  0, 0, 0);
962 
963  physiHoleFirstCollimatorModulatorBox = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
964  "HoleFirstCollimatorModulatorBox",
965  logicHoleFirstCollimatorModulatorBox,
966  physiFirstCollimatorModulatorBox, false, 0);
967 
968  // --------------------------------------------------//
969  // SECOND SIDE OF THE MODULATOR BOX //
970  // --------------------------------------------------//
971  const G4double secondCollimatorModulatorXSize = 10.*mm;
972  const G4double secondCollimatorModulatorYSize = 200.*mm;
973  const G4double secondCollimatorModulatorZSize = 200.*mm;
974 
975  const G4double secondCollimatorModulatorXPosition = -1953.00 *mm;
976 
977  const G4double secondCollimatorModulatorYPosition = 0.*mm;
978  const G4double secondCollimatorModulatorZPosition = 0.*mm;
979 
980  G4Box* solidSecondCollimatorModulatorBox = new G4Box("SecondCollimatorModulatorBox",
981  secondCollimatorModulatorXSize,
982  secondCollimatorModulatorYSize,
983  secondCollimatorModulatorZSize);
984 
985  G4LogicalVolume* logicSecondCollimatorModulatorBox = new G4LogicalVolume(solidSecondCollimatorModulatorBox,
986  modulatorBoxMaterial,
987  "SecondCollimatorModulatorBox");
988 
989  physiSecondCollimatorModulatorBox = new G4PVPlacement(0, G4ThreeVector(secondCollimatorModulatorXPosition,
990  secondCollimatorModulatorYPosition,
991  secondCollimatorModulatorZPosition),
992  "SecondCollimatorModulatorBox",
993  logicSecondCollimatorModulatorBox,
994  physicalTreatmentRoom, false, 0);
995 
996  // ----------------------------------------------//
997  // Hole of the second collimator modulator box //
998  // ----------------------------------------------//
999  const G4double innerRadiusHoleSecondCollimatorModulatorBox = 0.*mm;
1000  const G4double outerRadiusHoleSecondCollimatorModulatorBox = 31.*mm;
1001  const G4double hightHoleSecondCollimatorModulatorBox = 10.*mm;
1002  const G4double startAngleHoleSecondCollimatorModulatorBox = 0.*deg;
1003  const G4double spanningAngleHoleSecondCollimatorModulatorBox = 360.*deg;
1004 
1005  G4Tubs* solidHoleSecondCollimatorModulatorBox = new G4Tubs("HoleSecondCollimatorModulatorBox",
1006  innerRadiusHoleSecondCollimatorModulatorBox,
1007  outerRadiusHoleSecondCollimatorModulatorBox,
1008  hightHoleSecondCollimatorModulatorBox ,
1009  startAngleHoleSecondCollimatorModulatorBox,
1010  spanningAngleHoleSecondCollimatorModulatorBox);
1011 
1012  G4LogicalVolume* logicHoleSecondCollimatorModulatorBox = new G4LogicalVolume(solidHoleSecondCollimatorModulatorBox,
1013  holeModulatorBoxMaterial,
1014  "HoleSecondCollimatorModulatorBox",
1015  0, 0, 0);
1016 
1017  physiHoleSecondCollimatorModulatorBox = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
1018  "HoleSecondCollimatorModulatorBox",
1019  logicHoleSecondCollimatorModulatorBox,
1020  physiSecondCollimatorModulatorBox, false, 0);
1021 
1022  logicFirstCollimator -> SetVisAttributes(yellow);
1023  logicFirstCollimatorModulatorBox -> SetVisAttributes(blue);
1024  logicSecondCollimatorModulatorBox -> SetVisAttributes(blue);
1025 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveProtonBeamLine::HadrontherapyBeamFinalCollimator ( )

Definition at line 1490 of file PassiveProtonBeamLine.cc.

References python.hepunit::deg, python.hepunit::mm, and CLHEP::HepRotation::rotateY().

1491 {
1492  // -----------------------//
1493  // FINAL COLLIMATOR //
1494  //------------------------//
1495  const G4double outerRadiusFinalCollimator = 21.5*mm;
1496  const G4double hightFinalCollimator = 3.5*mm;
1497  const G4double startAngleFinalCollimator = 0.*deg;
1498  const G4double spanningAngleFinalCollimator = 360.*deg;
1499  const G4double finalCollimatorXPosition = -83.5 *mm;
1500 
1501  G4double phi = 90. *deg;
1502 
1503  // Matrix definition for a 90 deg rotation. Also used for other volumes
1504  G4RotationMatrix rm;
1505  rm.rotateY(phi);
1506 
1507  solidFinalCollimator = new G4Tubs("FinalCollimator",
1508  innerRadiusFinalCollimator,
1509  outerRadiusFinalCollimator,
1510  hightFinalCollimator,
1511  startAngleFinalCollimator,
1512  spanningAngleFinalCollimator);
1513 
1514  G4LogicalVolume* logicFinalCollimator = new G4LogicalVolume(solidFinalCollimator,
1515  finalCollimatorMaterial,
1516  "FinalCollimator",
1517  0,
1518  0,
1519  0);
1520 
1521  physiFinalCollimator = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(finalCollimatorXPosition,0.,0.)),
1522  "FinalCollimator", logicFinalCollimator, physicalTreatmentRoom, false, 0);
1523 
1524  logicFinalCollimator -> SetVisAttributes(yellow);
1525 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveProtonBeamLine::HadrontherapyBeamLineSupport ( )

Definition at line 603 of file PassiveProtonBeamLine.cc.

References python.hepunit::m, and python.hepunit::mm.

604 {
605  // ------------------//
606  // BEAM LINE SUPPORT //
607  //-------------------//
608  const G4double beamLineSupportXSize = 1.5*m;
609  const G4double beamLineSupportYSize = 20.*mm;
610  const G4double beamLineSupportZSize = 600.*mm;
611 
612  const G4double beamLineSupportXPosition = -1745.09 *mm;
613  const G4double beamLineSupportYPosition = -230. *mm;
614  const G4double beamLineSupportZPosition = 0.*mm;
615 
616  G4Box* beamLineSupport = new G4Box("BeamLineSupport",
617  beamLineSupportXSize,
618  beamLineSupportYSize,
619  beamLineSupportZSize);
620 
621  G4LogicalVolume* logicBeamLineSupport = new G4LogicalVolume(beamLineSupport,
622  beamLineSupportMaterial,
623  "BeamLineSupport");
624  physiBeamLineSupport = new G4PVPlacement(0, G4ThreeVector(beamLineSupportXPosition,
625  beamLineSupportYPosition,
626  beamLineSupportZPosition),
627  "BeamLineSupport",
628  logicBeamLineSupport,
629  physicalTreatmentRoom, false, 0);
630 
631  // Visualisation attributes of the beam line support
632 
633  logicBeamLineSupport -> SetVisAttributes(gray);
634 
635  //---------------------------------//
636  // Beam line cover 1 (left panel) //
637  //---------------------------------//
638  const G4double beamLineCoverXSize = 1.5*m;
639  const G4double beamLineCoverYSize = 750.*mm;
640  const G4double beamLineCoverZSize = 10.*mm;
641 
642  const G4double beamLineCoverXPosition = -1745.09 *mm;
643  const G4double beamLineCoverYPosition = -1000.*mm;
644  const G4double beamLineCoverZPosition = 600.*mm;
645 
646  G4Box* beamLineCover = new G4Box("BeamLineCover",
647  beamLineCoverXSize,
648  beamLineCoverYSize,
649  beamLineCoverZSize);
650 
651  G4LogicalVolume* logicBeamLineCover = new G4LogicalVolume(beamLineCover,
652  beamLineSupportMaterial,
653  "BeamLineCover");
654 
655  physiBeamLineCover = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition,
656  beamLineCoverYPosition,
657  beamLineCoverZPosition),
658  "BeamLineCover",
659  logicBeamLineCover,
660  physicalTreatmentRoom,
661  false,
662  0);
663 
664  // ---------------------------------//
665  // Beam line cover 2 (rigth panel) //
666  // ---------------------------------//
667  // It has the same characteristic of beam line cover 1 but set in a different position
668  physiBeamLineCover2 = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition,
669  beamLineCoverYPosition,
670  - beamLineCoverZPosition),
671  "BeamLineCover2",
672  logicBeamLineCover,
673  physicalTreatmentRoom,
674  false,
675  0);
676 
677  logicBeamLineCover -> SetVisAttributes(blue);
678 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
double G4double
Definition: G4Types.hh:76
void PassiveProtonBeamLine::HadrontherapyBeamMonitoring ( )

Definition at line 1028 of file PassiveProtonBeamLine.cc.

References python.hepunit::cm, and python.hepunit::mm.

1029 {
1030  // ----------------------------
1031  // THE FIRST MONITOR CHAMBER
1032  // ----------------------------
1033  // A monitor chamber is a free-air ionisation chamber
1034  // able to measure do proton fluence during the treatment.
1035  // Here its responce is not simulated in terms of produced
1036  // charge but only the energy losses are taked into account.
1037  // Each chamber consist of 9 mm of air in a box
1038  // that has two layers one of kapton and one
1039  // of copper
1040  const G4double monitor1XSize = 4.525022*mm;
1041  const G4double monitor2XSize = 0.000011*mm;
1042  const G4double monitor3XSize = 4.5*mm;
1043  const G4double monitorYSize = 10.*cm;
1044  const G4double monitorZSize = 10.*cm;
1045  const G4double monitor1XPosition = -1262.47498 *mm;
1046  const G4double monitor2XPosition = -4.500011*mm;
1047  const G4double monitor4XPosition = 4.500011*mm;
1048 
1049  G4Box* solidFirstMonitorLayer1 = new G4Box("FirstMonitorLayer1",
1050  monitor1XSize,
1051  monitorYSize,
1052  monitorZSize);
1053 
1054  G4LogicalVolume* logicFirstMonitorLayer1 = new G4LogicalVolume(solidFirstMonitorLayer1,
1055  layer1MonitorChamberMaterial,
1056  "FirstMonitorLayer1");
1057 
1058  physiFirstMonitorLayer1 = new G4PVPlacement(0,
1059  G4ThreeVector(monitor1XPosition,0.*cm,0.*cm),
1060  "FirstMonitorLayer1",
1061  logicFirstMonitorLayer1,
1062  physicalTreatmentRoom,
1063  false,
1064  0);
1065 
1066  G4Box* solidFirstMonitorLayer2 = new G4Box("FirstMonitorLayer2",
1067  monitor2XSize,
1068  monitorYSize,
1069  monitorZSize);
1070 
1071  G4LogicalVolume* logicFirstMonitorLayer2 = new G4LogicalVolume(solidFirstMonitorLayer2,
1072  layer2MonitorChamberMaterial,
1073  "FirstMonitorLayer2");
1074 
1075  physiFirstMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector(monitor2XPosition,0.*cm,0.*cm),
1076  "FirstMonitorLayer2",
1077  logicFirstMonitorLayer2,
1078  physiFirstMonitorLayer1,
1079  false,
1080  0);
1081 
1082  G4Box* solidFirstMonitorLayer3 = new G4Box("FirstMonitorLayer3",
1083  monitor3XSize,
1084  monitorYSize,
1085  monitorZSize);
1086 
1087  G4LogicalVolume* logicFirstMonitorLayer3 = new G4LogicalVolume(solidFirstMonitorLayer3,
1088  layer3MonitorChamberMaterial,
1089  "FirstMonitorLayer3");
1090 
1091  physiFirstMonitorLayer3 = new G4PVPlacement(0,
1092  G4ThreeVector(0.*mm,0.*cm,0.*cm),
1093  "MonitorLayer3",
1094  logicFirstMonitorLayer3,
1095  physiFirstMonitorLayer1,
1096  false,
1097  0);
1098 
1099  G4Box* solidFirstMonitorLayer4 = new G4Box("FirstMonitorLayer4",
1100  monitor2XSize,
1101  monitorYSize,
1102  monitorZSize);
1103 
1104  G4LogicalVolume* logicFirstMonitorLayer4 = new G4LogicalVolume(solidFirstMonitorLayer4,
1105  layer4MonitorChamberMaterial,
1106  "FirstMonitorLayer4");
1107 
1108  physiFirstMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm),
1109  "FirstMonitorLayer4",
1110  logicFirstMonitorLayer4,
1111  physiFirstMonitorLayer1, false, 0);
1112  // ----------------------------//
1113  // THE SECOND MONITOR CHAMBER //
1114  // ----------------------------//
1115  physiSecondMonitorLayer1 = new G4PVPlacement(0, G4ThreeVector(-1131.42493 *mm,0.*cm,0.*cm),
1116  "SecondMonitorLayer1", logicFirstMonitorLayer1,physicalTreatmentRoom, false, 0);
1117 
1118  physiSecondMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector( monitor2XPosition,0.*cm,0.*cm), "SecondMonitorLayer2",
1119  logicFirstMonitorLayer2, physiSecondMonitorLayer1, false, 0);
1120 
1121  physiSecondMonitorLayer3 = new G4PVPlacement(0, G4ThreeVector(0.*mm,0.*cm,0.*cm), "MonitorLayer3",
1122  logicFirstMonitorLayer3, physiSecondMonitorLayer1, false, 0);
1123 
1124  physiSecondMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm), "SecondMonitorLayer4",
1125  logicFirstMonitorLayer4, physiSecondMonitorLayer1, false, 0);
1126 
1127  logicFirstMonitorLayer3 -> SetVisAttributes(white);
1128 
1129 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
double G4double
Definition: G4Types.hh:76
void PassiveProtonBeamLine::HadrontherapyBeamNozzle ( )

Definition at line 1306 of file PassiveProtonBeamLine.cc.

References python.hepunit::deg, python.hepunit::mm, and CLHEP::HepRotation::rotateY().

1307 {
1308  // ------------------------------//
1309  // THE FINAL TUBE AND COLLIMATOR //
1310  //-------------------------------//
1311  // The last part of the transport beam line consists of
1312  // a 59 mm thick PMMA slab (to stop all the diffused radiation), a 370 mm brass tube
1313  // (to well collimate the proton beam) and a final collimator with 25 mm diameter
1314  // aperture (that provide the final trasversal shape of the beam)
1315 
1316  // -------------------//
1317  // PMMA SUPPORT //
1318  // -------------------//
1319  const G4double nozzleSupportXSize = 29.5 *mm;
1320  const G4double nozzleSupportYSize = 180. *mm;
1321  const G4double nozzleSupportZSize = 180. *mm;
1322 
1323  const G4double nozzleSupportXPosition = -397.50 *mm;
1324 
1325  G4double phi = 90. *deg;
1326  // Matrix definition for a 90 deg rotation. Also used for other volumes
1327  G4RotationMatrix rm;
1328  rm.rotateY(phi);
1329 
1330  G4Box* solidNozzleSupport = new G4Box("NozzleSupport",
1331  nozzleSupportXSize,
1332  nozzleSupportYSize,
1333  nozzleSupportZSize);
1334 
1335  G4LogicalVolume* logicNozzleSupport = new G4LogicalVolume(solidNozzleSupport,
1336  nozzleSupportMaterial,
1337  "NozzleSupport");
1338 
1339  physiNozzleSupport = new G4PVPlacement(0, G4ThreeVector(nozzleSupportXPosition,0., 0.),
1340  "NozzleSupport",
1341  logicNozzleSupport,
1342  physicalTreatmentRoom,
1343  false,
1344  0);
1345 
1346  logicNozzleSupport -> SetVisAttributes(yellow);
1347 
1348 
1349 
1350  //------------------------------------//
1351  // HOLE IN THE SUPPORT //
1352  //------------------------------------//
1353  const G4double innerRadiusHoleNozzleSupport = 0.*mm;
1354  const G4double outerRadiusHoleNozzleSupport = 21.5*mm;
1355  const G4double hightHoleNozzleSupport = 29.5 *mm;
1356  const G4double startAngleHoleNozzleSupport = 0.*deg;
1357  const G4double spanningAngleHoleNozzleSupport = 360.*deg;
1358 
1359  G4Tubs* solidHoleNozzleSupport = new G4Tubs("HoleNozzleSupport",
1360  innerRadiusHoleNozzleSupport,
1361  outerRadiusHoleNozzleSupport,
1362  hightHoleNozzleSupport,
1363  startAngleHoleNozzleSupport,
1364  spanningAngleHoleNozzleSupport);
1365 
1366  G4LogicalVolume* logicHoleNozzleSupport = new G4LogicalVolume(solidHoleNozzleSupport,
1367  holeNozzleSupportMaterial,
1368  "HoleNozzleSupport",
1369  0,
1370  0,
1371  0);
1372 
1373 
1374  physiHoleNozzleSupport = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
1375  "HoleNozzleSupport",
1376  logicHoleNozzleSupport,
1377  physiNozzleSupport,
1378  false, 0);
1379 
1380  logicHoleNozzleSupport -> SetVisAttributes(darkOrange3);
1381 
1382  // ---------------------------------//
1383  // BRASS TUBE 1 (phantom side) //
1384  // ---------------------------------//
1385  const G4double innerRadiusBrassTube= 18.*mm;
1386  const G4double outerRadiusBrassTube = 21.5 *mm;
1387  const G4double hightBrassTube = 140.5*mm;
1388  const G4double startAngleBrassTube = 0.*deg;
1389  const G4double spanningAngleBrassTube = 360.*deg;
1390 
1391  const G4double brassTubeXPosition = -227.5 *mm;
1392 
1393  G4Tubs* solidBrassTube = new G4Tubs("BrassTube",
1394  innerRadiusBrassTube,
1395  outerRadiusBrassTube,
1396  hightBrassTube,
1397  startAngleBrassTube,
1398  spanningAngleBrassTube);
1399 
1400  G4LogicalVolume* logicBrassTube = new G4LogicalVolume(solidBrassTube,
1401  brassTubeMaterial,
1402  "BrassTube",
1403  0, 0, 0);
1404 
1405  physiBrassTube = new G4PVPlacement(G4Transform3D(rm,
1406  G4ThreeVector(brassTubeXPosition,
1407  0.,
1408  0.)),
1409  "BrassTube",
1410  logicBrassTube,
1411  physicalTreatmentRoom,
1412  false,
1413  0);
1414 
1415  logicBrassTube -> SetVisAttributes(darkOrange3);
1416 
1417  // ----------------------------------------------//
1418  // BRASS TUBE 2 (inside the PMMA support) //
1419  // ----------------------------------------------//
1420  const G4double innerRadiusBrassTube2= 18.*mm;
1421  const G4double outerRadiusBrassTube2 = 21.5 *mm;
1422  const G4double hightBrassTube2 = 29.5*mm;
1423  const G4double startAngleBrassTube2 = 0.*deg;
1424  const G4double spanningAngleBrassTube2 = 360.*deg;
1425 
1426  // const G4double brassTube2XPosition = -227.5 *mm;
1427 
1428  G4Tubs* solidBrassTube2 = new G4Tubs("BrassTube2",
1429  innerRadiusBrassTube2,
1430  outerRadiusBrassTube2,
1431  hightBrassTube2,
1432  startAngleBrassTube2,
1433  spanningAngleBrassTube2);
1434 
1435  G4LogicalVolume* logicBrassTube2 = new G4LogicalVolume(solidBrassTube2,
1436  brassTube2Material,
1437  "BrassTube2",
1438  0, 0, 0);
1439 
1440  physiBrassTube2 = new G4PVPlacement(G4Transform3D(rm,
1441  G4ThreeVector(0,
1442  0.,
1443  0.)),
1444  "BrassTube2",
1445  logicBrassTube2,
1446  physiNozzleSupport,
1447  false,
1448  0);
1449 
1450  logicBrassTube2 -> SetVisAttributes(darkOrange3);
1451 
1452 
1453  // --------------------------------------//
1454  // BRASS TUBE 3 (beam line side) //
1455  // -------------------------------------//
1456  const G4double innerRadiusBrassTube3= 18.*mm;
1457  const G4double outerRadiusBrassTube3 = 21.5 *mm;
1458  const G4double hightBrassTube3 = 10.0 *mm;
1459  const G4double startAngleBrassTube3 = 0.*deg;
1460  const G4double spanningAngleBrassTube3 = 360.*deg;
1461 
1462  const G4double brassTube3XPosition = -437 *mm;
1463 
1464  G4Tubs* solidBrassTube3 = new G4Tubs("BrassTube3",
1465  innerRadiusBrassTube3,
1466  outerRadiusBrassTube3,
1467  hightBrassTube3,
1468  startAngleBrassTube3,
1469  spanningAngleBrassTube3);
1470 
1471  G4LogicalVolume* logicBrassTube3 = new G4LogicalVolume(solidBrassTube3,
1472  brassTube3Material,
1473  "BrassTube3",
1474  0, 0, 0);
1475 
1476  physiBrassTube3 = new G4PVPlacement(G4Transform3D(rm,
1477  G4ThreeVector(brassTube3XPosition,
1478  0.,
1479  0.)),
1480  "BrassTube3",
1481  logicBrassTube3,
1482  physicalTreatmentRoom,
1483  false,
1484  0);
1485 
1486  logicBrassTube3 -> SetVisAttributes(darkOrange3);
1487 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveProtonBeamLine::HadrontherapyBeamScatteringFoils ( )

Definition at line 681 of file PassiveProtonBeamLine.cc.

References python.hepunit::deg, and CLHEP::HepRotation::rotateY().

682 {
683  // ------------//
684  // VACUUM PIPE //
685  //-------------//
686  //
687  // First track of the beam line is inside vacuum;
688  // The PIPE contains the FIRST SCATTERING FOIL and the KAPTON WINDOW
689  G4Box* vacuumZone = new G4Box("VacuumZone", vacuumZoneXSize, vacuumZoneYSize, vacuumZoneZSize);
690  G4LogicalVolume* logicVacuumZone = new G4LogicalVolume(vacuumZone, vacuumZoneMaterial, "VacuumZone");
691  G4VPhysicalVolume* physiVacuumZone = new G4PVPlacement(0, G4ThreeVector(vacuumZoneXPosition, 0., 0.),
692  "VacuumZone", logicVacuumZone, physicalTreatmentRoom, false, 0);
693  // --------------------------//
694  // THE FIRST SCATTERING FOIL //
695  // --------------------------//
696  // A thin foil performing a first scattering
697  // of the original beam
698  firstScatteringFoil = new G4Box("FirstScatteringFoil",
699  firstScatteringFoilXSize,
700  firstScatteringFoilYSize,
701  firstScatteringFoilZSize);
702 
703  G4LogicalVolume* logicFirstScatteringFoil = new G4LogicalVolume(firstScatteringFoil,
704  firstScatteringFoilMaterial,
705  "FirstScatteringFoil");
706 
707  physiFirstScatteringFoil = new G4PVPlacement(0, G4ThreeVector(firstScatteringFoilXPosition, 0.,0.),
708  "FirstScatteringFoil", logicFirstScatteringFoil, physiVacuumZone,
709  false, 0);
710 
711  logicFirstScatteringFoil -> SetVisAttributes(skyBlue);
712  // -------------------//
713  // THE KAPTON WINDOWS //
714  //--------------------//
715  //It prmits the passage of the beam from vacuum to air
716  G4Box* solidKaptonWindow = new G4Box("KaptonWindow",
717  kaptonWindowXSize,
718  kaptonWindowYSize,
719  kaptonWindowZSize);
720 
721  G4LogicalVolume* logicKaptonWindow = new G4LogicalVolume(solidKaptonWindow,
722  kaptonWindowMaterial,
723  "KaptonWindow");
724 
725  physiKaptonWindow = new G4PVPlacement(0, G4ThreeVector(kaptonWindowXPosition, 0., 0.),
726  "KaptonWindow", logicKaptonWindow,
727  physiVacuumZone, false, 0);
728 
729  logicKaptonWindow -> SetVisAttributes(darkOrange3);
730 
731  // ------------//
732  // THE STOPPER //
733  //-------------//
734  // Is a small cylinder able to stop the central component
735  // of the beam (having a gaussian shape). It is connected to the SECON SCATTERING FOIL
736  // and represent the second element of the scattering system
737  G4double phi = 90. *deg;
738  // Matrix definition for a 90 deg rotation with respect to Y axis
739  G4RotationMatrix rm;
740  rm.rotateY(phi);
741 
742  solidStopper = new G4Tubs("Stopper",
743  innerRadiusStopper,
744  outerRadiusStopper,
745  heightStopper,
746  startAngleStopper,
747  spanningAngleStopper);
748 
749  logicStopper = new G4LogicalVolume(solidStopper,
750  stopperMaterial,
751  "Stopper",
752  0, 0, 0);
753 
754  physiStopper = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(stopperXPosition,
755  stopperYPosition,
756  stopperZPosition)),
757  "Stopper",
758  logicStopper,
759  physicalTreatmentRoom,
760  false,
761  0);
762 
763  logicStopper -> SetVisAttributes(red);
764 
765  // ---------------------------//
766  // THE SECOND SCATTERING FOIL //
767  // ---------------------------//
768  // It is another thin foil and provides the
769  // final diffusion of the beam. It represents the third element of the scattering
770  // system;
771 
772  secondScatteringFoil = new G4Box("SecondScatteringFoil",
773  secondScatteringFoilXSize,
774  secondScatteringFoilYSize,
775  secondScatteringFoilZSize);
776 
777  G4LogicalVolume* logicSecondScatteringFoil = new G4LogicalVolume(secondScatteringFoil,
778  secondScatteringFoilMaterial,
779  "SecondScatteringFoil");
780 
781  physiSecondScatteringFoil = new G4PVPlacement(0, G4ThreeVector(secondScatteringFoilXPosition,
782  secondScatteringFoilYPosition,
783  secondScatteringFoilZPosition),
784  "SeconScatteringFoil",
785  logicSecondScatteringFoil,
786  physicalTreatmentRoom,
787  false,
788  0);
789 
790  logicSecondScatteringFoil -> SetVisAttributes(skyBlue);
791 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveProtonBeamLine::HadrontherapyMOPIDetector ( )

Definition at line 1131 of file PassiveProtonBeamLine.cc.

1132 {
1133  // --------------------------------//
1134  // THE MOPI DETECTOR //
1135  // --------------------------------//
1136  // MOPI DETECTOR: two orthogonal microstrip gas detectors developed
1137  // by the INFN Section of Turin in collaboration with some
1138  // of the author of this example. It permits the
1139  // on-line check of the beam simmetry via the signal
1140  // integration of the collected charge for each strip.
1141  //
1142  // In this example it is simulated as:
1143  // 1. First anode: 35 mu of kapton + 15 mu of aluminum,
1144  // 2. First air gap: 6 mm of air,
1145  // 3. The cathode: 1 mu Al + 25 mu mylar + 1 mu Al
1146  // (in common with the two air gap),
1147  // 4. Second air gap: 6 mm of air,
1148  // 5 Second anode: 15 mu Al + 35 mu kapton
1149  // Color used in the graphical output
1150 
1151 
1152  // Mother volume
1153  solidMOPIMotherVolume = new G4Box("MOPIMotherVolume",
1154  MOPIMotherVolumeXSize/2,
1155  MOPIMotherVolumeYSize/2,
1156  MOPIMotherVolumeYSize/2);
1157 
1158  logicMOPIMotherVolume = new G4LogicalVolume(solidMOPIMotherVolume,
1159  MOPIMotherVolumeMaterial,
1160  "MOPIMotherVolume");
1161  physiMOPIMotherVolume = new G4PVPlacement(0,
1162  G4ThreeVector(MOPIMotherVolumeXPosition,
1163  MOPIMotherVolumeYPosition,
1164  MOPIMotherVolumeZPosition),
1165  "MOPIMotherVolume",
1166  logicMOPIMotherVolume,
1167  physicalTreatmentRoom,
1168  false,
1169  0);
1170 
1171  // First Kapton layer
1172  solidMOPIFirstKaptonLayer = new G4Box("MOPIFirstKaptonLayer",
1173  MOPIFirstKaptonLayerXSize/2,
1174  MOPIFirstKaptonLayerYSize/2 ,
1175  MOPIFirstKaptonLayerZSize/2);
1176 
1177  logicMOPIFirstKaptonLayer = new G4LogicalVolume(solidMOPIFirstKaptonLayer,
1178  MOPIFirstKaptonLayerMaterial,
1179  "MOPIFirstKaptonLayer");
1180 
1181  physiMOPIFirstKaptonLayer = new G4PVPlacement(0,
1182  G4ThreeVector(MOPIFirstKaptonLayerXPosition,
1183  MOPIFirstKaptonLayerYPosition ,
1184  MOPIFirstKaptonLayerZPosition),
1185  "MOPIFirstKaptonLayer",
1186  logicMOPIFirstKaptonLayer,
1187  physiMOPIMotherVolume,
1188  false,
1189  0);
1190 
1191  // First Aluminum layer
1192  solidMOPIFirstAluminumLayer = new G4Box("MOPIFirstAluminumLayer",
1193  MOPIFirstAluminumLayerXSize/2,
1194  MOPIFirstAluminumLayerYSize/2 ,
1195  MOPIFirstAluminumLayerZSize/2);
1196 
1197  logicMOPIFirstAluminumLayer = new G4LogicalVolume(solidMOPIFirstAluminumLayer,
1198  MOPIFirstAluminumLayerMaterial,
1199  "MOPIFirstAluminumLayer");
1200 
1201  physiMOPIFirstAluminumLayer = new G4PVPlacement(0,
1202  G4ThreeVector(MOPIFirstAluminumLayerXPosition,
1203  MOPIFirstAluminumLayerYPosition ,
1204  MOPIFirstAluminumLayerZPosition),
1205  "MOPIFirstAluminumLayer",
1206  logicMOPIFirstAluminumLayer, physiMOPIMotherVolume, false, 0);
1207 
1208  // First Air GAP
1209  solidMOPIFirstAirGap = new G4Box("MOPIFirstAirGap",
1210  MOPIFirstAirGapXSize/2,
1211  MOPIFirstAirGapYSize/2,
1212  MOPIFirstAirGapZSize/2);
1213 
1214  logicMOPIFirstAirGap = new G4LogicalVolume(solidMOPIFirstAirGap,
1215  MOPIFirstAirGapMaterial,
1216  "MOPIFirstAirgap");
1217 
1218  physiMOPIFirstAirGap = new G4PVPlacement(0,
1219  G4ThreeVector(MOPIFirstAirGapXPosition,
1220  MOPIFirstAirGapYPosition ,
1221  MOPIFirstAirGapZPosition),
1222  "MOPIFirstAirGap",
1223  logicMOPIFirstAirGap, physiMOPIMotherVolume, false, 0);
1224 
1225 
1226  // The Cathode
1227  solidMOPICathode = new G4Box("MOPICathode",
1228  MOPICathodeXSize/2,
1229  MOPICathodeYSize/2,
1230  MOPICathodeZSize/2);
1231 
1232  logicMOPICathode = new G4LogicalVolume(solidMOPICathode,
1233  MOPICathodeMaterial,
1234  "MOPICathode");
1235 
1236  physiMOPICathode = new G4PVPlacement(0,
1237  G4ThreeVector(MOPICathodeXPosition,
1238  MOPICathodeYPosition ,
1239  MOPICathodeZPosition),
1240  "MOPICathode",
1241  logicMOPICathode,
1242  physiMOPIMotherVolume, false, 0);
1243 
1244  // Second Air GAP
1245  solidMOPISecondAirGap = new G4Box("MOPISecondAirGap",
1246  MOPISecondAirGapXSize/2,
1247  MOPISecondAirGapYSize/2,
1248  MOPISecondAirGapZSize/2);
1249 
1250  logicMOPISecondAirGap = new G4LogicalVolume(solidMOPISecondAirGap,
1251  MOPISecondAirGapMaterial,
1252  "MOPISecondAirgap");
1253 
1254  physiMOPISecondAirGap = new G4PVPlacement(0,
1255  G4ThreeVector(MOPISecondAirGapXPosition,
1256  MOPISecondAirGapYPosition ,
1257  MOPISecondAirGapZPosition),
1258  "MOPISecondAirGap",
1259  logicMOPISecondAirGap, physiMOPIMotherVolume, false, 0);
1260 
1261  // Second Aluminum layer
1262  solidMOPISecondAluminumLayer = new G4Box("MOPISecondAluminumLayer",
1263  MOPISecondAluminumLayerXSize/2,
1264  MOPISecondAluminumLayerYSize/2 ,
1265  MOPISecondAluminumLayerZSize/2);
1266 
1267  logicMOPISecondAluminumLayer = new G4LogicalVolume(solidMOPISecondAluminumLayer,
1268  MOPISecondAluminumLayerMaterial,
1269  "MOPISecondAluminumLayer");
1270 
1271  physiMOPISecondAluminumLayer = new G4PVPlacement(0,
1272  G4ThreeVector(MOPISecondAluminumLayerXPosition,
1273  MOPISecondAluminumLayerYPosition ,
1274  MOPISecondAluminumLayerZPosition),
1275  "MOPISecondAluminumLayer",
1276  logicMOPISecondAluminumLayer,
1277  physiMOPIMotherVolume,
1278  false,
1279  0);
1280 
1281  // Second Kapton layer
1282  solidMOPISecondKaptonLayer = new G4Box("MOPISecondKaptonLayer",
1283  MOPISecondKaptonLayerXSize/2,
1284  MOPISecondKaptonLayerYSize/2 ,
1285  MOPISecondKaptonLayerZSize/2);
1286 
1287  logicMOPISecondKaptonLayer = new G4LogicalVolume(solidMOPISecondKaptonLayer,
1288  MOPIFirstKaptonLayerMaterial,
1289  "MOPISecondKaptonLayer");
1290 
1291  physiMOPISecondKaptonLayer = new G4PVPlacement(0,
1292  G4ThreeVector(MOPISecondKaptonLayerXPosition,
1293  MOPISecondKaptonLayerYPosition ,
1294  MOPISecondKaptonLayerZPosition),
1295  "MOPISecondKaptonLayer",
1296  logicMOPISecondKaptonLayer,
1297  physiMOPIMotherVolume,
1298  false,
1299  0);
1300 
1301  logicMOPIFirstAirGap -> SetVisAttributes(darkGreen);
1302  logicMOPISecondAirGap -> SetVisAttributes(darkGreen);
1303 
1304 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
void PassiveProtonBeamLine::HadrontherapyRangeShifter ( )

Definition at line 793 of file PassiveProtonBeamLine.cc.

794 {
795  // ---------------------------- //
796  // THE RANGE SHIFTER //
797  // -----------------------------//
798  // It is a slab of PMMA acting as energy degreader of
799  // primary beam
800  solidRangeShifterBox = new G4Box("RangeShifterBox",
801  rangeShifterXSize,
802  rangeShifterYSize,
803  rangeShifterZSize);
804 
805  logicRangeShifterBox = new G4LogicalVolume(solidRangeShifterBox,
806  rangeShifterMaterial,
807  "RangeShifterBox");
808  physiRangeShifterBox = new G4PVPlacement(0,
809  G4ThreeVector(rangeShifterXPosition, 0., 0.),
810  "RangeShifterBox",
811  logicRangeShifterBox,
812  physicalTreatmentRoom,
813  false,
814  0);
815 
816 
817  logicRangeShifterBox -> SetVisAttributes(yellow);
818 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
void PassiveProtonBeamLine::SetFirstScatteringFoilXSize ( G4double  value)

Definition at line 1545 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), and python.hepunit::mm.

1546 {
1547  firstScatteringFoil -> SetXHalfLength(value);
1548  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1549  G4cout <<"The X size of the first scattering foil is (mm):"<<
1550  ((firstScatteringFoil -> GetXHalfLength())*2.)/mm
1551  << G4endl;
1552 }
G4GLOB_DLL std::ostream G4cout
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61
void PassiveProtonBeamLine::SetInnerRadiusFinalCollimator ( G4double  value)

Definition at line 1575 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), and python.hepunit::mm.

1576 {
1577  solidFinalCollimator -> SetInnerRadius(value);
1578  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1579  G4cout<<"Inner Radius of the final collimator is (mm):"
1580  << solidFinalCollimator -> GetInnerRadius()/mm
1581  << G4endl;
1582 }
G4GLOB_DLL std::ostream G4cout
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61
void PassiveProtonBeamLine::SetModulatorAngle ( G4double  angle)

Definition at line 1605 of file PassiveProtonBeamLine.cc.

1606 {
1607  modulator -> SetModulatorAngle(value);
1608  //G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1609 }
void SetModulatorAngle(G4double angle)
const XML_Char int const XML_Char * value
void PassiveProtonBeamLine::SetOuterRadiusStopper ( G4double  value)

Definition at line 1565 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), and python.hepunit::mm.

1566 {
1567  solidStopper -> SetOuterRadius(value);
1568  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1569  G4cout << "OuterRadius od the Stopper is (mm):"
1570  << solidStopper -> GetOuterRadius()/mm
1571  << G4endl;
1572 }
G4GLOB_DLL std::ostream G4cout
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61
void PassiveProtonBeamLine::SetRangeShifterXPosition ( G4double  value)

Definition at line 1528 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), and python.hepunit::mm.

1529 {
1530  physiRangeShifterBox -> SetTranslation(G4ThreeVector(value, 0., 0.));
1531  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1532  G4cout << "The Range Shifter is translated to"<< value/mm <<"mm along the X axis" <<G4endl;
1533 }
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cout
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61
void PassiveProtonBeamLine::SetRangeShifterXSize ( G4double  halfSize)

Definition at line 1536 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), and python.hepunit::mm.

1537 {
1538  solidRangeShifterBox -> SetXHalfLength(value) ;
1539  G4cout << "RangeShifter size X (mm): "<< ((solidRangeShifterBox -> GetXHalfLength())*2.)/mm
1540  << G4endl;
1541  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1542 }
G4GLOB_DLL std::ostream G4cout
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61
void PassiveProtonBeamLine::SetRSMaterial ( G4String  materialChoice)

Definition at line 1585 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4NistManager::Instance(), and EmPlot::SetMaterial().

1586 {
1587  if (G4Material* pttoMaterial = G4NistManager::Instance()->FindOrBuildMaterial(materialChoice, false) )
1588  {
1589  if (pttoMaterial)
1590  {
1591  rangeShifterMaterial = pttoMaterial;
1592  logicRangeShifterBox -> SetMaterial(pttoMaterial);
1593  G4cout << "The material of the Range Shifter has been changed to " << materialChoice << G4endl;
1594  }
1595  }
1596  else
1597  {
1598  G4cout << "WARNING: material \"" << materialChoice << "\" doesn't exist in NIST elements/materials"
1599  " table [located in $G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc]" << G4endl;
1600  G4cout << "Use command \"/parameter/nist\" to see full materials list!" << G4endl;
1601  }
1602 }
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
def SetMaterial
Definition: EmPlot.py:25
#define G4endl
Definition: G4ios.hh:61
void PassiveProtonBeamLine::SetSecondScatteringFoilXSize ( G4double  value)

Definition at line 1555 of file PassiveProtonBeamLine.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), and python.hepunit::mm.

1556 {
1557  secondScatteringFoil -> SetXHalfLength(value);
1558  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1559  G4cout <<"The X size of the second scattering foil is (mm):"<<
1560  ((secondScatteringFoil -> GetXHalfLength())*2.)/mm
1561  << G4endl;
1562 }
G4GLOB_DLL std::ostream G4cout
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61

The documentation for this class was generated from the following files: