Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-31 09:22:19

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // Hadrontherapy advanced example for Geant4
0027 // See more at: https://twiki.cern.ch/twiki/bin/view/Geant4/AdvancedExamplesHadrontherapy
0028 
0029 #include "globals.hh"
0030 #include "G4SystemOfUnits.hh"
0031 #include "G4Box.hh"
0032 #include "G4Tubs.hh"
0033 #include "G4UnionSolid.hh"
0034 #include "G4Trd.hh"
0035 #include "G4AssemblyVolume.hh"
0036 #include "G4VisAttributes.hh"
0037 #include "G4Colour.hh"
0038 #include "G4RunManager.hh"
0039 #include "G4LogicalVolume.hh"
0040 #include "G4PVPlacement.hh"
0041 #include "G4PVReplica.hh"
0042 #include "G4RotationMatrix.hh"
0043 #include "G4NistManager.hh"
0044 #include "G4NistElementBuilder.hh"
0045 #include "HadrontherapyDetectorConstruction.hh"
0046 #include "HadrontherapyTIFPAPassiveProtonBeamLine.hh"
0047 #include "HadrontherapyTIFPAPassiveProtonBeamLineMessenger.hh"
0048 
0049 /////////////////////////////////////////////////////////////////////////////
0050 TrentoPassiveProtonBeamLine::TrentoPassiveProtonBeamLine():
0051  logicTreatmentRoom(0), physicalTreatmentRoom(0),hadrontherapyDetectorConstruction(0),
0052 physiBeamLineSupport(0), physiBeamLineCover(0), physiBeamLineCover2(0),
0053 physiMonitorLayer1(0), physiMonitorLayer2(0),
0054   ScatteringFoil(0), logicScatteringFoil(0), physiScatteringFoil(0), ridgeFilterPhys(0), preCollimator(0), physiPreCollimator(0), Collimator(0), physiCollimator(0), solidAirTube(0), solidAirPreTube(0), physiAirTube(0), physiAirPreTube(0)
0055 
0056 
0057 {
0058     // Messenger to change parameters of the passiveTrentoBeamLine geometry
0059     TrentoPassiveMessenger = new TrentoPassiveProtonBeamLineMessenger(this);
0060     
0061     //***************************** PW ***********************************
0062     static G4String ROGeometryName = "DetectorROGeometry";
0063     RO = new HadrontherapyDetectorROGeometry(ROGeometryName);
0064     
0065     G4cout << "Going to register Parallel world...";
0066     RegisterParallelWorld(RO);
0067     G4cout << "... done" << G4endl;
0068     //********************************************************************
0069     
0070 }
0071 /////////////////////////////////////////////////////////////////////////////
0072 TrentoPassiveProtonBeamLine::~TrentoPassiveProtonBeamLine()
0073 {
0074     delete TrentoPassiveMessenger;
0075     delete hadrontherapyDetectorConstruction;
0076 }
0077 
0078 using namespace std;
0079 
0080 /////////////////////////////////////////////////////////////////////////////
0081 G4VPhysicalVolume* TrentoPassiveProtonBeamLine::Construct()
0082 {
0083     // Sets default geometry and materials
0084     SetDefaultDimensions();
0085     
0086     // Construct the whole Passive Beam Line
0087     ConstructTrentoPassiveProtonBeamLine();
0088     
0089     //***************************** PW ***************************************
0090     if (!hadrontherapyDetectorConstruction)
0091         
0092         // HadrontherapyDetectorConstruction builds ONLY the phantom and the detector with its associated ROGeometry
0093         hadrontherapyDetectorConstruction = new HadrontherapyDetectorConstruction(physicalTreatmentRoom);
0094     
0095     
0096     //********************************************************************
0097     
0098     hadrontherapyDetectorConstruction->InitializeDetectorROGeometry(RO,hadrontherapyDetectorConstruction->GetDetectorToWorldPosition());
0099     
0100     return physicalTreatmentRoom;
0101 }
0102 
0103 // In the following method the DEFAULTS used in the geometry of
0104 // passive beam line are provided
0105 // HERE THE USER CAN CHANGE THE GEOMETRY CHARACTERISTICS OF BEAM
0106 // LINE ELEMENTS, ALTERNATIVELY HE/SHE CAN USE THE MACRO FILE (IF A
0107 // MESSENGER IS PROVIDED)
0108 //
0109 // DEFAULT MATERIAL ARE ALSO PROVIDED
0110 // and COLOURS ARE ALSO DEFINED
0111 // ----------------------------------------------------------
0112 /////////////////////////////////////////////////////////////////////////////
0113 void TrentoPassiveProtonBeamLine::SetDefaultDimensions()
0114 {
0115     // Set of coulors that can be used
0116     white = new G4VisAttributes( G4Colour());
0117     white -> SetVisibility(true);
0118     white -> SetForceSolid(true);
0119     
0120     blue = new G4VisAttributes(G4Colour(0. ,0. ,1.));
0121     blue -> SetVisibility(true);
0122     blue -> SetForceSolid(true);
0123     
0124     gray = new G4VisAttributes( G4Colour(0.5, 0.5, 0.5 ));
0125     gray-> SetVisibility(true);
0126     gray-> SetForceSolid(true);
0127     
0128     red = new G4VisAttributes(G4Colour(1. ,0. ,0.));
0129     red-> SetVisibility(true);
0130     red-> SetForceSolid(true);
0131     
0132     yellow = new G4VisAttributes(G4Colour(1., 1., 0. ));
0133     yellow-> SetVisibility(true);
0134     yellow-> SetForceSolid(true);
0135     
0136     green = new G4VisAttributes( G4Colour(25/255. , 255/255. ,  25/255. ));
0137     green -> SetVisibility(true);
0138     green -> SetForceSolid(true);
0139     
0140     darkGreen = new G4VisAttributes( G4Colour(0/255. , 100/255. ,  0/255. ));
0141     darkGreen -> SetVisibility(true);
0142     darkGreen -> SetForceSolid(true);
0143     
0144     darkOrange3 = new G4VisAttributes( G4Colour(205/255. , 102/255. ,  000/255. ));
0145     darkOrange3 -> SetVisibility(true);
0146     darkOrange3 -> SetForceSolid(false);
0147     
0148     skyBlue = new G4VisAttributes( G4Colour(135/255. , 206/255. ,  235/255. ));
0149     skyBlue -> SetVisibility(true);
0150     skyBlue -> SetForceSolid(true);
0151     
0152     
0153     
0154     // SCATTERING FOIL: it is a thin foil that provides the
0155     // final diffusion of the beam.
0156     G4double defaultScatteringFoilXSize = 0.75*cm;
0157     ScatteringFoilXSize = defaultScatteringFoilXSize;
0158     
0159     G4double defaultScatteringFoilYSize = 100*mm;
0160     ScatteringFoilYSize = defaultScatteringFoilYSize;
0161     
0162     G4double defaultScatteringFoilZSize = 100 *mm;
0163     ScatteringFoilZSize = defaultScatteringFoilZSize;
0164     
0165     G4double defaultScatteringFoilXPosition = -273.*cm;
0166     ScatteringFoilXPosition = defaultScatteringFoilXPosition;
0167     
0168     G4double defaultScatteringFoilYPosition =  0 *mm;
0169     ScatteringFoilYPosition = defaultScatteringFoilYPosition;
0170     
0171     G4double defaultScatteringFoilZPosition =  0 *mm;
0172     ScatteringFoilZPosition = defaultScatteringFoilZPosition;
0173 
0174     //PRE COLLIMATOR: it is a PMMA collimator
0175     G4double defaultPreCollimatorXHalfSide = 12.5 * cm;
0176     preCollimatorXHalfSide = defaultPreCollimatorXHalfSide;
0177 
0178     G4double defaultPreCollimatorXPosition = -50.*cm;
0179     preCollimatorXPosition = defaultPreCollimatorXPosition;
0180 
0181     //DEFAULT DEFINITION OF THE AIR TUBE
0182     G4double defaultYHalfSideAirTube= 5.*cm;
0183     YHalfSideAirTube = defaultYHalfSideAirTube;
0184 
0185     G4double defaultZHalfSideAirTube = 5.*cm;
0186     ZHalfSideAirTube = defaultZHalfSideAirTube;
0187 
0188     
0189     // DEFAULT DEFINITION OF THE MATERIALS
0190     // All elements and compound definition follows the NIST database
0191     
0192     // ELEMENTS
0193     G4bool isotopes = false;
0194     G4Material* aluminumNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Al", isotopes);
0195     G4Material* copperNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Cu", isotopes);
0196     G4Element* zincNist = G4NistManager::Instance()->FindOrBuildElement("Zn");
0197     G4Element* copper = G4NistManager::Instance()->FindOrBuildElement("Cu");
0198     G4Element* silicNist = G4NistManager::Instance()->FindOrBuildElement("Si");
0199     G4Element* hydrogenNist = G4NistManager::Instance()->FindOrBuildElement("H");
0200     G4Element* oxygenNist = G4NistManager::Instance()->FindOrBuildElement("O");
0201     G4Element* carbonNist = G4NistManager::Instance()->FindOrBuildElement("C");
0202     
0203     // COMPOUND
0204     G4Material* airNist =  G4NistManager::Instance()->FindOrBuildMaterial("G4_AIR", isotopes);
0205     G4Material* PMMANist = G4NistManager::Instance()->FindOrBuildMaterial("G4_PLEXIGLASS", isotopes);
0206     G4Material* MylarNist =  G4NistManager::Instance()->FindOrBuildMaterial("G4_MYLAR");
0207 
0208     G4int nComponents; // Number of components
0209     G4int nAtoms; //Number of atoms
0210     G4double fractionmass; // Fraction in mass of an element in a material
0211 
0212     //Quartz
0213     G4Material* SiO2 = new G4Material("quartz", 2.200*g/cm3, nComponents=2);
0214     SiO2->AddElement(silicNist, nAtoms=1);
0215     SiO2->AddElement(oxygenNist , nAtoms=2);
0216     
0217     //Epoxy (for FR4 )
0218     G4Material* Epoxy = new G4Material("Epoxy" , 1.2*g/cm3, nComponents=2);
0219     Epoxy->AddElement(hydrogenNist, nAtoms=2);
0220     Epoxy->AddElement(carbonNist, nAtoms=2);
0221     
0222     //FR4 (Glass + Epoxy)
0223     G4Material* FR4 = new G4Material("FR4"  , 1.86*g/cm3, nComponents=2);
0224     FR4->AddMaterial(SiO2, fractionmass=0.528);
0225     FR4->AddMaterial(Epoxy, fractionmass=0.472);
0226     
0227     //Brass
0228     G4Material* brass = new G4Material("Brass", 8.40*g/cm3, nComponents=2);
0229     brass -> AddElement(zincNist, fractionmass = 30 *perCent);
0230     brass -> AddElement(copper, fractionmass = 70 *perCent);
0231 
0232     //Plastic
0233     G4Material* plastic = new G4Material("Plastic", 1.40*g/cm3, nComponents=3);
0234     plastic -> AddElement(carbonNist, nAtoms = 21);
0235     plastic -> AddElement(oxygenNist, nAtoms = 4);
0236     plastic -> AddElement(hydrogenNist, nAtoms = 24);
0237      
0238     
0239     //***************************** PW ***************************************
0240     
0241     // DetectorROGeometry Material
0242     new G4Material("dummyMat", 1., 1.*g/mole, 1.*g/cm3);
0243     
0244     
0245     // MATERIAL ASSIGNMENT
0246     // Support of the beam line
0247     beamLineSupportMaterial = aluminumNist;
0248 
0249     // Matreial of the monitor chamber
0250     layerMonitorChamberMaterial = MylarNist;
0251     layerDefaultMaterial = airNist;
0252     internalStructureMaterial =  FR4;
0253     FoilMaterial = copperNist;
0254     airgapMaterial = airNist;
0255     
0256     // Material of the scattering foil
0257     ScatteringFoilMaterial = copperNist;
0258 
0259     //Material of the ridge filter
0260     singleTrapMaterial = plastic;
0261     
0262     // Materials of the collimators
0263     preCollimatorMaterial = PMMANist;
0264     CollimatorMaterial = brass;
0265 
0266     // Material of the final brass tube
0267     airTubeMaterial = airTube2Material = airTube3Material = airNist;
0268     
0269 }
0270 
0271 /////////////////////////////////////////////////////////////////////////////
0272 void TrentoPassiveProtonBeamLine::ConstructTrentoPassiveProtonBeamLine()
0273 {
0274     // -----------------------------
0275     // Treatment room - World volume
0276     //------------------------------
0277     // Treatment room sizes
0278     const G4double worldX = 400.0 *cm;
0279     const G4double worldY = 400.0 *cm;
0280     const G4double worldZ = 400.0 *cm;
0281     G4bool isotopes = false;
0282     
0283     G4Material* airNist =  G4NistManager::Instance()->FindOrBuildMaterial("G4_AIR", isotopes);
0284     
0285     G4Box* treatmentRoom = new G4Box("TreatmentRoom",
0286                                      worldX,
0287                                      worldY,
0288                                      worldZ);
0289     
0290     logicTreatmentRoom = new G4LogicalVolume(treatmentRoom,
0291                                              airNist,
0292                                              "logicTreatmentRoom",
0293                                              0,0,0);
0294     
0295     physicalTreatmentRoom = new G4PVPlacement(0,
0296                                               G4ThreeVector(),
0297                                               "physicalTreatmentRoom",
0298                                               logicTreatmentRoom,
0299                                               0,
0300                                               false,
0301                                               0);
0302     
0303     
0304     // The treatment room is invisible in the Visualisation
0305     logicTreatmentRoom -> SetVisAttributes (G4VisAttributes::GetInvisible());
0306     
0307     // Components of the Passive Proton Beam Line
0308     HadrontherapyBeamLineSupport();
0309     HadrontherapyBeamMonitoring();
0310     HadrontherapyBeamScatteringFoils();
0311     HadrontherapyRidgeFilter();
0312     HadrontherapyBeamCollimators();
0313     
0314 }
0315 
0316 /////////////////////////////////////////////////////////////////////////////
0317 void TrentoPassiveProtonBeamLine::HadrontherapyBeamLineSupport()
0318 {
0319     // ------------------//
0320     // BEAM LINE SUPPORT //
0321     //-------------------//
0322     const G4double beamLineSupportXSize = 1.2*m;
0323     const G4double beamLineSupportYSize = 20.*mm;
0324     const G4double beamLineSupportZSize = 600.*mm;
0325     
0326     const G4double beamLineSupportXPosition = -1745.09 *mm;
0327     const G4double beamLineSupportYPosition = -230. *mm;
0328     const G4double beamLineSupportZPosition = 0.*mm;
0329     
0330     G4Box* beamLineSupport = new G4Box("BeamLineSupport",
0331                                        beamLineSupportXSize,
0332                                        beamLineSupportYSize,
0333                                        beamLineSupportZSize);
0334     
0335     G4LogicalVolume* logicBeamLineSupport = new G4LogicalVolume(beamLineSupport,
0336                                                                 beamLineSupportMaterial,
0337                                                                 "BeamLineSupport");
0338     
0339     physiBeamLineSupport = new G4PVPlacement(0, G4ThreeVector(beamLineSupportXPosition,
0340                                                               beamLineSupportYPosition,
0341                                                               beamLineSupportZPosition),
0342                                              "BeamLineSupport",
0343                                              logicBeamLineSupport,
0344                                              physicalTreatmentRoom,
0345                                              false,
0346                                              0);
0347     
0348     // Visualisation attributes of the beam line support
0349     logicBeamLineSupport -> SetVisAttributes(gray);
0350     
0351     //---------------------------------//
0352     //  Beam line cover 1 (left panel) //
0353     //---------------------------------//
0354     const G4double beamLineCoverXSize = 1.2*m;
0355     const G4double beamLineCoverYSize = 750.*mm;
0356     const G4double beamLineCoverZSize = 10.*mm;
0357     
0358     const G4double beamLineCoverXPosition = -1745.09 *mm;
0359     const G4double beamLineCoverYPosition = -1000.*mm;
0360     const G4double beamLineCoverZPosition = 600.*mm;
0361     
0362     G4Box* beamLineCover = new G4Box("BeamLineCover",
0363                                      beamLineCoverXSize,
0364                                      beamLineCoverYSize,
0365                                      beamLineCoverZSize);
0366     
0367     G4LogicalVolume* logicBeamLineCover = new G4LogicalVolume(beamLineCover,
0368                                                               beamLineSupportMaterial,
0369                                                               "BeamLineCover");
0370     
0371     physiBeamLineCover = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition,
0372                                                             beamLineCoverYPosition,
0373                                                             beamLineCoverZPosition),
0374                                            "BeamLineCover",
0375                                            logicBeamLineCover,
0376                                            physicalTreatmentRoom,
0377                                            false,
0378                                            0);
0379     
0380     // ---------------------------------//
0381     //  Beam line cover 2 (rigth panel) //
0382     // ---------------------------------//
0383     // It has the same characteristic of beam line cover 1 but set in a different position
0384     physiBeamLineCover2 = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition,
0385                                                              beamLineCoverYPosition,
0386                                                              - beamLineCoverZPosition),
0387                                             "BeamLineCover2",
0388                                             logicBeamLineCover,
0389                                             physicalTreatmentRoom,
0390                                             false,
0391                                             0);
0392     
0393     logicBeamLineCover -> SetVisAttributes(blue);
0394 }
0395 
0396 /////////////////////////////////////////////////////////////////////////////
0397 void TrentoPassiveProtonBeamLine::HadrontherapyBeamMonitoring()
0398 {
0399   // ----------------------------
0400   //   MONITOR CHAMBER
0401   // ----------------------------
0402   // The  monitor chamber is a ionisation chamber
0403   // Each chamber consist in a sequence of 2 mm of air and
0404   // 8 microm copper layers which contain 800 microm of Fr4 in a box
0405   // that has the two layer of Mylar
0406 
0407     
0408 ////////////////////////////////////////////
0409 ///External Structure of the Monitor Chamber
0410 ////////////////////////////////////////////
0411     
0412   const G4double monitorXSize = 12.*um;
0413   const G4double monitorYSize = 12.7*cm;
0414   const G4double monitorZSize = 12.7*cm;
0415   
0416   const G4double shift = 17.*cm;
0417   
0418   const G4double monitorlayer1XPosition = -269.*cm +shift;
0419   const G4double monitorlayer2XPosition = -270.4104*cm +shift;
0420   
0421   
0422 // First Mylar layer of the monitor chamber
0423   G4Box* solidMonitorLayer1 = new G4Box("MonitorLayer1",
0424                     monitorXSize/2,
0425                     monitorYSize/2,
0426                     monitorZSize/2);
0427   
0428   G4LogicalVolume* logicMonitorLayer1 = new G4LogicalVolume(solidMonitorLayer1,
0429                                 layerMonitorChamberMaterial,
0430                                 "MonitorLayer1");
0431     
0432   physiMonitorLayer1 = new G4PVPlacement(0,
0433                                          G4ThreeVector(monitorlayer1XPosition, 0.*cm, 0.*cm),
0434                                          "MonitorLayer1",
0435                                          logicMonitorLayer1,
0436                                          physicalTreatmentRoom,
0437                                          false,
0438                                          0);
0439 
0440 // Second Mylar layer of the monitor chamber
0441   G4Box* solidMonitorLayer2 = new G4Box("MonitorLayer2",
0442                                         monitorXSize/2,
0443                                         monitorYSize/2,
0444                                         monitorZSize/2);
0445   
0446   G4LogicalVolume* logicMonitorLayer2 = new G4LogicalVolume(solidMonitorLayer2,
0447                                                             layerMonitorChamberMaterial,
0448                                                             "MonitorLayer2");
0449   
0450   physiMonitorLayer2 = new G4PVPlacement(0,
0451                      G4ThreeVector(monitorlayer2XPosition, 0.*cm, 0.*cm),
0452                      "MonitorLayer2",
0453                      logicMonitorLayer2,
0454                      physicalTreatmentRoom,
0455                      false,
0456                      0);    
0457   
0458   logicMonitorLayer1 -> SetVisAttributes(gray);
0459   logicMonitorLayer2 -> SetVisAttributes(gray);
0460   
0461     
0462 ///////////////////////////////////////////
0463 // Internal Structure of the Monitor Chamber
0464 ////////////////////////////////////////////
0465     
0466   const G4double layerThickness = 816*um;
0467   const G4double layerXposition = -270.2684*cm;
0468   const G4int nofLayers = 5;
0469   const G4double internalStructureThickness = 800.*um;
0470   const G4double airGapThickness = 2.*mm;
0471   const G4double foilThickness = 8.*um;
0472   const G4double FirstCoppperLayerXPosition = -404.*um;
0473   const G4double SecondCoppperLayerXPosition = 404.*um;
0474   const G4double InternalStructureXPosition = 0.;
0475   
0476 // Air Layer
0477   G4VSolid* SolidAirLayer= new G4Box("Layer",
0478                               layerThickness,
0479                               monitorYSize/2,
0480                               monitorZSize/2);
0481     
0482   new G4LogicalVolume(SolidAirLayer,
0483                       layerDefaultMaterial,
0484                       "Layer");
0485 
0486 
0487 
0488 // First Copper Layer
0489   G4VSolid* SolidFirstCoppperLayer = new G4Box("SolidFirstCoppperLayer",
0490                                   foilThickness/2,
0491                                   monitorYSize/2,
0492                                   monitorZSize/2);
0493   
0494   G4LogicalVolume* LogicFirstCoppperLayer = new G4LogicalVolume(
0495                                                      SolidFirstCoppperLayer,
0496                                                      FoilMaterial,
0497                                                      "SolidFirstCoppperLayer");
0498   
0499   
0500 // Fr4 Internal Layer
0501   G4VSolid* SolidInternalStructure = new G4Box("SolidInternalStructure",
0502                                           internalStructureThickness/2,
0503                                           monitorYSize/2,
0504                                           monitorZSize/2);
0505   
0506   G4LogicalVolume* LogicInternalStructure = new G4LogicalVolume(
0507                                                              SolidInternalStructure,
0508                                                              internalStructureMaterial,
0509                                                              "SolidInternalStructure");
0510 
0511     
0512     
0513 // Second Copper Layer
0514   G4VSolid* SolidSecondCoppperLayer = new G4Box("SolidSecondCoppperLayer",
0515                                    foilThickness/2,
0516                                    monitorYSize/2,
0517                                    monitorZSize/2); // its size
0518   
0519   G4LogicalVolume* LogicSecondCoppperLayer= new G4LogicalVolume(
0520                                                      SolidSecondCoppperLayer,
0521                                                      FoilMaterial,
0522                                                      "SolidSecondCoppperLayer");
0523 
0524   
0525   G4RotationMatrix Ra, Rm;
0526   G4ThreeVector Ta;
0527   G4Transform3D Tr;
0528   G4AssemblyVolume* assemblyMonitor = new G4AssemblyVolume();
0529 
0530  
0531     Ta.setX(FirstCoppperLayerXPosition); Ta.setY(0.); Ta.setZ(0.);
0532     Tr = G4Transform3D(Ra,Ta);
0533     assemblyMonitor->AddPlacedVolume(LogicFirstCoppperLayer, Tr);
0534     
0535     Ta.setX(InternalStructureXPosition); Ta.setY(0.); Ta.setZ(0.);
0536     Tr = G4Transform3D(Ra,Ta);
0537     assemblyMonitor->AddPlacedVolume(LogicInternalStructure, Tr );
0538     
0539     Ta.setX(SecondCoppperLayerXPosition); Ta.setY(0.); Ta.setZ(0.);
0540     Tr = G4Transform3D(Ra,Ta);
0541     
0542     assemblyMonitor->AddPlacedVolume(LogicSecondCoppperLayer, Tr );
0543    
0544   
0545     for( unsigned int i = 0; i<nofLayers; i++ )
0546    {
0547      G4ThreeVector Tm(shift+layerXposition + i*(airGapThickness), 0., 0.);
0548      Tr = G4Transform3D(Rm,Tm);
0549      assemblyMonitor -> MakeImprint(logicTreatmentRoom, Tr);
0550    }
0551     
0552 }
0553 
0554 
0555 void TrentoPassiveProtonBeamLine::HadrontherapyBeamScatteringFoils()
0556 {
0557   
0558   // ---------------------------//
0559   // SCATTERING FOIL            //
0560   // ---------------------------//
0561   // It is a metal foil and provides the
0562   // initial diffusion of the beam. 
0563   
0564   
0565   
0566   ScatteringFoil = new G4Box("ScatteringFoil",
0567                  ScatteringFoilXSize,
0568                  ScatteringFoilYSize,
0569                  ScatteringFoilZSize);
0570   
0571   
0572   logicScatteringFoil = new G4LogicalVolume(ScatteringFoil,
0573                         ScatteringFoilMaterial,
0574                         "ScatteringFoil");
0575   
0576   physiScatteringFoil = new G4PVPlacement(0, G4ThreeVector(ScatteringFoilXPosition,
0577                                ScatteringFoilYPosition,
0578                                ScatteringFoilZPosition),
0579                       "SeconScatteringFoil",
0580                                               logicScatteringFoil,
0581                       physicalTreatmentRoom,
0582                       false,
0583                       0);
0584 
0585   logicScatteringFoil -> SetVisAttributes(skyBlue);
0586 }
0587 
0588 void TrentoPassiveProtonBeamLine::HadrontherapyRidgeFilter()
0589 {
0590   // ---------------------------- //
0591   //         THE Ridge Filter    //
0592   // -----------------------------//
0593   // Energy degreader of
0594   // primary beam. Made of a series of pin-substructures.
0595   
0596     
0597     G4double ridgeXPosition = -200*cm;
0598     
0599     const G4double XBase = 0.1 *mm;
0600     const G4double YBase = 38.75*mm;
0601     
0602     G4VSolid* RidgeBase = new G4Box("Base_RidgeFilter",
0603                                     XBase,
0604                                     YBase,
0605                                     YBase);
0606     
0607     G4LogicalVolume* RidgeBaseLog = new G4LogicalVolume(RidgeBase,
0608                                                         singleTrapMaterial,
0609                                                         "Base_RidgeFilter");
0610     
0611     new G4PVPlacement(0,
0612                       G4ThreeVector(
0613                                     -199.7*cm,
0614                                     0.,
0615                                     0.),
0616                       "Base_RidgeFilter",
0617                       RidgeBaseLog,
0618                       physicalTreatmentRoom,
0619                       false,
0620                       0);
0621     
0622     RidgeBaseLog->SetVisAttributes(yellow);
0623     
0624     
0625 
0626   G4double x0 = 1.215*mm;
0627   G4double x1 = 1*mm;
0628   G4double x2 = 0.95*mm;
0629   G4double x3 = 0.87*mm;
0630   G4double x4 = 0.76*mm;
0631   G4double x5 = 0.71*mm;
0632   G4double x6 = 0.65*mm;
0633   G4double x7 = 0.56*mm;
0634   G4double x8 = 0.529*mm;
0635   G4double x9 = 0.258*mm;
0636   G4double x10 = 0.225*mm;
0637   G4double x11 = 0.144*mm;
0638   G4double x12 = 0.055*mm;
0639     
0640   G4double heigth = 0.13*mm;
0641   G4double heigth0 = 0.11*mm;
0642   G4double heigth1 = 2.3*mm;
0643   G4double heigth2 = 0.65*mm;
0644   G4double heigth3 = 1.9*mm;
0645   G4double heigth4 = 0.615*mm;
0646   G4double heigth5 = 2.23*mm;
0647   G4double heigth6 = 0.3*mm;
0648   G4double heigth7 = 4.1*mm;
0649   G4double heigth8 = 0.1*mm;
0650   G4double heigth9 = 0.105*mm;
0651   G4double heigth10 = 0.065*mm;
0652   
0653   G4VSolid* Trap00 = new G4Trd("singleTrap00", x0, x1, x0, x1, heigth);
0654   G4VSolid* Trap0 = new G4Trd("singleTrap0", x1, x2, x1, x2, heigth0);
0655   G4VSolid* Trap1 = new G4Trd("singleTrap1", x2, x3, x2, x3, heigth1);
0656   G4VSolid* Trap2 = new G4Trd("singleTrap2", x3, x4, x3, x4, heigth2);
0657   G4VSolid* Trap3 = new G4Trd("singleTrap3", x4, x5, x4, x5, heigth3);
0658   G4VSolid* Trap4 = new G4Trd("singleTrap4", x5, x6, x5, x6, heigth4);
0659   G4VSolid* Trap5 = new G4Trd("singleTrap5", x6, x7, x6, x7, heigth5);
0660   G4VSolid* Trap6 = new G4Trd("singleTrap6", x7, x8, x7, x8, heigth6);
0661   G4VSolid* Trap7 = new G4Trd("singleTrap7", x8, x9, x8, x9, heigth7);
0662   G4VSolid* Trap8 = new G4Trd("singleTrap8", x9, x10, x9, x10, heigth8);
0663   G4VSolid* Trap9 = new G4Trd("singleTrap9", x10, x11, x10, x11, heigth9);
0664   G4VSolid* Trap10 = new G4Trd("singleTrap10", x11, x12, x11, x12, heigth10);
0665   
0666  
0667   G4ThreeVector tr0(0., 0., 0.24);
0668   G4ThreeVector tr1(0., 0., 2.54);
0669   G4ThreeVector tr2(0., 0., 2.55);
0670   G4ThreeVector tr3(0., 0., 2.845);
0671   G4ThreeVector tr4(0., 0., 4.4);
0672   
0673   G4RotationMatrix* yRot = new G4RotationMatrix; 
0674   yRot->rotateY(0); 
0675   G4UnionSolid* unionTrap00 = new G4UnionSolid("Trap01", Trap00, Trap0, yRot, tr0);
0676   G4UnionSolid* unionTrap01 = new G4UnionSolid("Trap01", unionTrap00, Trap1, yRot, tr1);
0677   G4UnionSolid* unionTrap23 = new G4UnionSolid("Trap23", Trap2, Trap3, yRot, tr2);
0678   G4UnionSolid* unionTrap45 = new G4UnionSolid("Trap45", Trap4, Trap5, yRot, tr3);
0679   G4UnionSolid* unionTrap67 = new G4UnionSolid("Trap67", Trap6, Trap7, yRot, tr4);
0680   
0681   G4ThreeVector tr03(0., 0., 5.09);
0682   G4UnionSolid* unionTrap03 = new G4UnionSolid("unionTrap03", unionTrap01, unionTrap23, yRot, tr03);
0683   
0684   G4ThreeVector tr05(0., 0., 7.935);
0685   G4UnionSolid* unionTrap05 = new G4UnionSolid("unionTrap05", unionTrap03, unionTrap45, yRot, tr05);
0686     
0687   G4ThreeVector tr_blocco1(0., 0., 12.335);
0688   G4UnionSolid* unionTrap_blocco1 = new G4UnionSolid("unionTrap_blocco1", unionTrap05, unionTrap67, yRot, tr_blocco1);
0689   
0690   G4ThreeVector tr_blocco2( 0., 0., 12.435);
0691   G4UnionSolid* unionTrap_blocco2 = new G4UnionSolid("unionTrap_blocco2", unionTrap_blocco1, Trap8, yRot, tr_blocco2);
0692   
0693   G4ThreeVector tr_blocco3(0., 0., 12.54);
0694   G4UnionSolid* unionTrap_blocco3 = new G4UnionSolid("unionTrap_blocco3", unionTrap_blocco2, Trap9, yRot, tr_blocco3);
0695   
0696   G4ThreeVector tr_tot( 0., 0., 12.605);
0697   G4UnionSolid* unionTrap = new G4UnionSolid("unionTrap", unionTrap_blocco3, Trap10, yRot, tr_tot);
0698 
0699   
0700   G4LogicalVolume* singleTrapLog = new G4LogicalVolume(unionTrap, singleTrapMaterial, "singleTrap");
0701   
0702 
0703   singleTrapLog->SetVisAttributes(yellow);
0704   
0705   
0706     G4int numberOfLayers = 31;
0707     G4double minZ = -37.5*mm;
0708     G4double minY = -37.5*mm;
0709     G4double sum_space = 1.25*mm;
0710     
0711     vector<G4ThreeVector> singleTrapPositions;
0712     
0713     for (int i = 0; i < numberOfLayers; i++)
0714       {
0715         for (int j = 0; j < numberOfLayers; j++)
0716           {
0717         singleTrapPositions.push_back({-0.01*cm+ridgeXPosition,
0718               minY + 2*i*sum_space,
0719               minZ + 2*j*sum_space,});
0720           }
0721       }
0722     
0723     G4double ti = -  90. *deg;
0724     G4RotationMatrix rt;
0725     rt.rotateY(ti);
0726     
0727     G4int peaks = numberOfLayers*numberOfLayers;
0728     for (int i = 0; i < peaks; i++)
0729       {
0730     
0731         ostringstream tName;tName << "singleTrap" << i;
0732         new G4PVPlacement(G4Transform3D(rt,
0733                                         singleTrapPositions[i]),
0734                                         singleTrapLog,
0735                                         "tName.str()",
0736                                         logicTreatmentRoom,
0737                                         0,
0738                                         i);
0739     
0740       }
0741   
0742     
0743 }
0744 
0745 void TrentoPassiveProtonBeamLine::HadrontherapyBeamCollimators()
0746 {
0747   
0748   // ------------------------//
0749   //     PRE - COLLIMATOR    //
0750   // -----------------------//
0751   // It is a plastic collimator to limit neutron dose and reduce activation
0752   
0753   const G4double preCollimatorYHalfSide = 10.*cm;
0754   const G4double preCollimatorZHalfSide = 10.*cm;
0755   
0756   preCollimator = new G4Box("PreCollimator",
0757                  preCollimatorXHalfSide,
0758                 preCollimatorYHalfSide,
0759                 preCollimatorZHalfSide);
0760   
0761   
0762   G4LogicalVolume* logicPreCollimator = new G4LogicalVolume(preCollimator,
0763                                 preCollimatorMaterial,
0764                                 "PreCollimator");
0765   
0766   
0767   physiPreCollimator = new G4PVPlacement(0,
0768                                          G4ThreeVector(
0769                                                        -50.*cm,
0770                                                            0.,
0771                               0.),
0772                      "PreCollimator",
0773                      logicPreCollimator,
0774                                          physicalTreatmentRoom,
0775                                          false,
0776                                          0);
0777   
0778   
0779   
0780   logicPreCollimator -> SetVisAttributes(white);
0781   
0782   
0783   // -----------------//
0784   //    COLLIMATOR     //
0785   // -----------------//
0786   // It is a brass collimator
0787   
0788   
0789 G4double collimatorYHalfSide = 10.*cm;
0790 G4double collimatorZHalfSide = 10.*cm;
0791 G4double collimatorXHalfSide  = 3.25*cm;
0792   
0793 G4double CollimatorXPosition = -34.25*cm;
0794   
0795   Collimator = new G4Box("Collimator",
0796              collimatorXHalfSide,
0797              collimatorYHalfSide,
0798              collimatorZHalfSide);
0799   
0800   
0801   G4LogicalVolume* logicCollimator = new G4LogicalVolume(Collimator,
0802                              CollimatorMaterial,
0803                              "Collimator");
0804   
0805   
0806   
0807   physiCollimator = new G4PVPlacement(0, G4ThreeVector(CollimatorXPosition,
0808                                0.,
0809                                0.),
0810                       "Collimator",
0811                       logicCollimator,
0812                       physicalTreatmentRoom,
0813                       false,
0814                       0);
0815   
0816   
0817   
0818   logicCollimator -> SetVisAttributes(darkGreen);
0819   
0820 
0821     
0822     // ---------------------------------//
0823     //      AIR BOX Collimator          //
0824     // ---------------------------------//
0825     
0826 
0827     solidAirTube = new G4Box("AirTube",
0828                              collimatorXHalfSide,
0829                              YHalfSideAirTube,
0830                              ZHalfSideAirTube);
0831     
0832     G4LogicalVolume* logicAirTube = new G4LogicalVolume(solidAirTube,
0833                                                         airTubeMaterial,
0834                                                         "AirTube",
0835                                                         0, 0, 0);
0836     
0837     
0838     physiAirTube = new G4PVPlacement(0, G4ThreeVector(0,
0839                                                       0.,
0840                                                       0.),
0841                                      "AirTube",
0842                                      logicAirTube,
0843                                      physiCollimator,
0844                                      false,
0845                                      0);
0846     
0847     logicAirTube -> SetVisAttributes(darkOrange3);
0848 
0849     
0850     
0851     // // ---------------------------------//
0852     // //      AIR BOX PreCollimator       //
0853     // // ---------------------------------//
0854 
0855     
0856     solidAirPreTube = new G4Box("AirPreTube",
0857                                 preCollimatorXHalfSide,
0858                                 YHalfSideAirTube,
0859                                 ZHalfSideAirTube);
0860     
0861     G4LogicalVolume* logicAirPreTube = new G4LogicalVolume(solidAirPreTube,
0862                                                            airTubeMaterial,
0863                                                            "AirPreTube",
0864                                                            0, 0, 0);
0865     
0866     
0867     physiAirPreTube = new G4PVPlacement(0,
0868                                         G4ThreeVector(0,
0869                                                       0.,
0870                                                       0.),
0871                                         "AirPreTube",
0872                                         logicAirPreTube,
0873                                         physiPreCollimator,
0874                                         false,
0875                                         0);
0876     
0877     logicAirPreTube -> SetVisAttributes(darkOrange3);
0878 
0879 
0880 
0881 }
0882 
0883 
0884 /////////////////////////////////////////////////////////////////////////////
0885 /////////////////////////// MESSENGER ///////////////////////////////////////
0886 /////////////////////////////////////////////////////////////////////////////
0887 
0888 
0889 void TrentoPassiveProtonBeamLine::SetScatteringFoilXSize(G4double value)
0890 {
0891   ScatteringFoil -> SetXHalfLength(value);
0892   G4RunManager::GetRunManager() -> GeometryHasBeenModified();
0893     G4cout <<"The X size of the second scattering foil is (mm):"<<
0894       ((ScatteringFoil -> GetXHalfLength())*2.)/mm
0895        << G4endl;
0896 }
0897 
0898 
0899 void TrentoPassiveProtonBeamLine::SetPreCollimatorXSize(G4double value)
0900 {
0901   preCollimator -> SetXHalfLength(value);
0902   G4cout<<"The size of the pre collimator is (mm)"
0903     << ((preCollimator -> GetXHalfLength())*2)/mm << G4endl;
0904   G4RunManager::GetRunManager() -> GeometryHasBeenModified();
0905    
0906 }
0907 
0908 void TrentoPassiveProtonBeamLine::SetPreCollimatorXPosition(G4double value)
0909 {
0910   physiPreCollimator -> SetTranslation(G4ThreeVector(value, 0., 0.));
0911   G4RunManager::GetRunManager() -> GeometryHasBeenModified();
0912   G4cout  <<" The pre collimator is translated to" << value/mm <<"mm along the X axis" << G4endl;
0913    
0914 }
0915 
0916 void TrentoPassiveProtonBeamLine::SetAirTubeYSize(G4double value)
0917 {
0918   solidAirTube -> SetYHalfLength(value);
0919   G4cout<<"The y side of brass tube is (mm)"
0920     << ((preCollimator -> GetYHalfLength())*2) << G4endl;
0921   G4RunManager::GetRunManager() -> GeometryHasBeenModified();   
0922 }
0923 
0924 
0925 void TrentoPassiveProtonBeamLine::SetAirTubeZSize(G4double value)
0926 {
0927   solidAirTube -> SetZHalfLength(value);
0928   G4cout<<"The z side of the brass tube is (mm)"
0929     << ((Collimator -> GetZHalfLength())*2) << G4endl;
0930   G4RunManager::GetRunManager() -> GeometryHasBeenModified();
0931   
0932 }
0933 
0934 
0935 /////////////////////////////////////////////////////////////////////////////
0936 void TrentoPassiveProtonBeamLine::SetScattererMaterial(G4String materialChoice)
0937 {
0938   if (G4Material* pttoMaterial = G4NistManager::Instance()->FindOrBuildMaterial(materialChoice, false) )
0939     {
0940       if (pttoMaterial)
0941         {
0942       ScatteringFoilMaterial  = pttoMaterial;
0943       logicScatteringFoil -> SetMaterial(pttoMaterial);
0944       G4cout << "The material of the Scatterer has been changed to " << materialChoice << G4endl;
0945         }
0946     }
0947   else
0948     {
0949       G4cout << "WARNING: material \"" << materialChoice << "\" doesn't exist in NIST elements/materials"
0950         " table [located in $G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc]" << G4endl;
0951       G4cout << "Use command \"/parameter/nist\" to see full materials list!" << G4endl;
0952     }
0953 }