Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PassiveProtonBeamLine Class Reference

#include <PassiveProtonBeamLine.hh>

Inheritance diagram for PassiveProtonBeamLine:
Collaboration diagram for PassiveProtonBeamLine:

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 46 of file PassiveProtonBeamLine.hh.

Constructor & Destructor Documentation

PassiveProtonBeamLine::PassiveProtonBeamLine ( )

Definition at line 49 of file PassiveProtonBeamLine.cc.

49  :
50 modulator(0), physicalTreatmentRoom(0),hadrontherapyDetectorConstruction(0),
51 physiBeamLineSupport(0), physiBeamLineCover(0), physiBeamLineCover2(0),
52 firstScatteringFoil(0), physiFirstScatteringFoil(0), physiKaptonWindow(0),
53 solidStopper(0), physiStopper(0), secondScatteringFoil(0), physiSecondScatteringFoil(0),
54 physiFirstCollimator(0), solidRangeShifterBox(0), logicRangeShifterBox(0),
55 physiRangeShifterBox(0), physiSecondCollimator(0), physiFirstCollimatorModulatorBox(0),
56 physiHoleFirstCollimatorModulatorBox(0), physiSecondCollimatorModulatorBox(0),
57 physiHoleSecondCollimatorModulatorBox(0), physiMOPIMotherVolume(0),
58 physiFirstMonitorLayer1(0), physiFirstMonitorLayer2(0), physiFirstMonitorLayer3(0),
59 physiFirstMonitorLayer4(0), physiSecondMonitorLayer1(0), physiSecondMonitorLayer2(0),
60 physiSecondMonitorLayer3(0), physiSecondMonitorLayer4(0), physiNozzleSupport(0), //physiHoleNozzleSupport(0),
61 physiBrassTube(0), solidFinalCollimator(0), physiFinalCollimator(0)
62 {
63  // Messenger to change parameters of the passiveProtonBeamLine geometry
64  passiveMessenger = new PassiveProtonBeamLineMessenger(this);
65 
66  //***************************** PW ***************************************
67  static G4String ROGeometryName = "DetectorROGeometry";
68  RO = new HadrontherapyDetectorROGeometry(ROGeometryName);
69 
70  G4cout << "Going to register Parallel world...";
72  G4cout << "... done" << G4endl;
73  //***************************** PW ***************************************
74  /* // For the Faraday Cup activation
75  if (name)
76  {
77  doCalculation = true;
78  }*/
79 }
void RegisterParallelWorld(G4VUserParallelWorld *)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

PassiveProtonBeamLine::~PassiveProtonBeamLine ( )

Definition at line 81 of file PassiveProtonBeamLine.cc.

82 {
83  delete passiveMessenger;
84  delete hadrontherapyDetectorConstruction;
85  /* if (!pFaradayCup)
86  {
87  delete pFaradayCup;
88  }*/
89 }

Member Function Documentation

G4VPhysicalVolume * PassiveProtonBeamLine::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 92 of file PassiveProtonBeamLine.cc.

93 {
94  // Sets default geometry and materials
95  SetDefaultDimensions();
96 
97  // Construct the whole Passive Beam Line
98  ConstructPassiveProtonBeamLine();
99 
100  //***************************** PW ***************************************
101  if (!hadrontherapyDetectorConstruction)
102 
103  //***************************** PW ***************************************
104 
105  // HadrontherapyDetectorConstruction builds ONLY the phantom and the detector with its associated ROGeometry
106  hadrontherapyDetectorConstruction = new HadrontherapyDetectorConstruction(physicalTreatmentRoom);
107 
108 
109  //***************************** PW ***************************************
110 
111  hadrontherapyDetectorConstruction->InitializeDetectorROGeometry(RO,hadrontherapyDetectorConstruction->GetDetectorToWorldPosition());
112 
113  //***************************** PW ***************************************
114  /* if (doCalculation)
115  {
116  pFaradayCup = new FaradayCup(physicalTreatmentRoom);
117  G4cout << "Faraday Cup detector has been activated" << G4endl;
118  }*/
119  return physicalTreatmentRoom;
120 }
void InitializeDetectorROGeometry(HadrontherapyDetectorROGeometry *, G4ThreeVector detectorToWorldPosition)

Here is the call graph for this function:

static PassiveProtonBeamLine* PassiveProtonBeamLine::GetInstance ( )
static
void PassiveProtonBeamLine::HadrontherapyBeamCollimators ( )

Definition at line 811 of file PassiveProtonBeamLine.cc.

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

Here is the call graph for this function:

void PassiveProtonBeamLine::HadrontherapyBeamFinalCollimator ( )

Definition at line 1484 of file PassiveProtonBeamLine.cc.

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

Here is the call graph for this function:

void PassiveProtonBeamLine::HadrontherapyBeamLineSupport ( )

Definition at line 594 of file PassiveProtonBeamLine.cc.

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

Definition at line 1019 of file PassiveProtonBeamLine.cc.

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

Definition at line 1297 of file PassiveProtonBeamLine.cc.

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

Here is the call graph for this function:

void PassiveProtonBeamLine::HadrontherapyBeamScatteringFoils ( )

Definition at line 672 of file PassiveProtonBeamLine.cc.

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

Here is the call graph for this function:

void PassiveProtonBeamLine::HadrontherapyMOPIDetector ( )

Definition at line 1122 of file PassiveProtonBeamLine.cc.

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

Definition at line 784 of file PassiveProtonBeamLine.cc.

785 {
786  // ---------------------------- //
787  // THE RANGE SHIFTER //
788  // -----------------------------//
789  // It is a slab of PMMA acting as energy degreader of
790  // primary beam
791  solidRangeShifterBox = new G4Box("RangeShifterBox",
792  rangeShifterXSize,
793  rangeShifterYSize,
794  rangeShifterZSize);
795 
796  logicRangeShifterBox = new G4LogicalVolume(solidRangeShifterBox,
797  rangeShifterMaterial,
798  "RangeShifterBox");
799  physiRangeShifterBox = new G4PVPlacement(0,
800  G4ThreeVector(rangeShifterXPosition, 0., 0.),
801  "RangeShifterBox",
802  logicRangeShifterBox,
803  physicalTreatmentRoom,
804  false,
805  0);
806 
807 
808  logicRangeShifterBox -> SetVisAttributes(yellow);
809 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:64
void PassiveProtonBeamLine::SetFirstScatteringFoilXSize ( G4double  value)

Definition at line 1539 of file PassiveProtonBeamLine.cc.

1540 {
1541  firstScatteringFoil -> SetXHalfLength(value);
1543  G4cout <<"The X size of the first scattering foil is (mm):"<<
1544  ((firstScatteringFoil -> GetXHalfLength())*2.)/mm
1545  << G4endl;
1546 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void PassiveProtonBeamLine::SetInnerRadiusFinalCollimator ( G4double  value)

Definition at line 1569 of file PassiveProtonBeamLine.cc.

1570 {
1571  solidFinalCollimator -> SetInnerRadius(value);
1573  G4cout<<"Inner Radius of the final collimator is (mm):"
1574  << solidFinalCollimator -> GetInnerRadius()/mm
1575  << G4endl;
1576 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void PassiveProtonBeamLine::SetModulatorAngle ( G4double  angle)

Definition at line 1599 of file PassiveProtonBeamLine.cc.

1600 {
1601  modulator -> SetModulatorAngle(value);
1602  //G4RunManager::GetRunManager() -> GeometryHasBeenModified();
1603 }
void SetModulatorAngle(G4double angle)
const XML_Char int const XML_Char * value
Definition: expat.h:331
void PassiveProtonBeamLine::SetOuterRadiusStopper ( G4double  value)

Definition at line 1559 of file PassiveProtonBeamLine.cc.

1560 {
1561  solidStopper -> SetOuterRadius(value);
1563  G4cout << "OuterRadius od the Stopper is (mm):"
1564  << solidStopper -> GetOuterRadius()/mm
1565  << G4endl;
1566 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void PassiveProtonBeamLine::SetRangeShifterXPosition ( G4double  value)

Definition at line 1522 of file PassiveProtonBeamLine.cc.

1523 {
1524  physiRangeShifterBox -> SetTranslation(G4ThreeVector(value, 0., 0.));
1526  G4cout << "The Range Shifter is translated to"<< value/mm <<"mm along the X axis" <<G4endl;
1527 }
static constexpr double mm
Definition: G4SIunits.hh:115
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void PassiveProtonBeamLine::SetRangeShifterXSize ( G4double  halfSize)

Definition at line 1530 of file PassiveProtonBeamLine.cc.

1531 {
1532  solidRangeShifterBox -> SetXHalfLength(value) ;
1533  G4cout << "RangeShifter size X (mm): "<< ((solidRangeShifterBox -> GetXHalfLength())*2.)/mm
1534  << G4endl;
1536 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void PassiveProtonBeamLine::SetRSMaterial ( G4String  materialChoice)

Definition at line 1579 of file PassiveProtonBeamLine.cc.

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

Here is the call graph for this function:

void PassiveProtonBeamLine::SetSecondScatteringFoilXSize ( G4double  value)

Definition at line 1549 of file PassiveProtonBeamLine.cc.

1550 {
1551  secondScatteringFoil -> SetXHalfLength(value);
1553  G4cout <<"The X size of the second scattering foil is (mm):"<<
1554  ((secondScatteringFoil -> GetXHalfLength())*2.)/mm
1555  << G4endl;
1556 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:


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