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 "G4Sphere.hh"
0034 #include "G4NistManager.hh"
0035 #include "G4NistElementBuilder.hh"
0036 #include "G4VisAttributes.hh"
0037 #include "G4Colour.hh"
0038 #include "G4RunManager.hh"
0039 #include "G4LogicalVolume.hh"
0040 #include "G4PVPlacement.hh"
0041 #include "G4RotationMatrix.hh"
0042 #include "HadrontherapyDetectorConstruction.hh"
0043 #include "LaserDrivenBeamLine.hh"
0044 #include "LaserDrivenBeamLineMessenger.hh"
0045 //
0046 #include "G4PhysicalConstants.hh"
0047 #include "G4ThreeVector.hh"
0048 #include "G4Material.hh"
0049 //
0050 #include "G4FieldManager.hh"
0051 #include "G4MagIntegratorStepper.hh"
0052 #include "G4Mag_UsualEqRhs.hh"
0053 #include "G4ExplicitEuler.hh"
0054 #include "G4ChordFinder.hh"
0055 //#include "G4TransportationManager.hh"
0056 #include "G4EqMagElectricField.hh"
0057 #include "G4UniformMagField.hh"
0058 #include "G4PropagatorInField.hh"
0059 #include "G4VisCommandsViewer.hh"
0060 #include "G4UImanager.hh"
0061 #include "G4ExplicitEuler.hh"
0062 #include "G4ImplicitEuler.hh"
0063 #include "G4SimpleRunge.hh"
0064 #include "G4SimpleHeum.hh"
0065 #include "G4ClassicalRK4.hh"
0066 #include "G4HelixExplicitEuler.hh"
0067 #include "G4HelixImplicitEuler.hh"
0068 #include "G4HelixSimpleRunge.hh"
0069 #include "G4CashKarpRKF45.hh"
0070 #include "G4RKG3_Stepper.hh"
0071 #include "G4MagIntegratorDriver.hh"
0072 #include "G4SubtractionSolid.hh"
0073 
0074 //
0075 #include "G4UniformElectricField.hh"
0076 #include "G4ElectricField.hh"
0077 #include "HadrontherapyElectricTabulatedField3D.hh"
0078 
0079 #include "HadrontherapyMagneticField3D.hh"
0080 //
0081 //G4bool LaserDrivenBeamLine::doCalculation = false;
0082 /////////////////////////////////////////////////////////////////////////////
0083 LaserDrivenBeamLine::LaserDrivenBeamLine():
0084 hadrontherapydetectorconstruction(0), physicTreatmentRoom(0),
0085 PFirstTriplet(0),PSecondTriplet(0),PThirdTriplet(0),PFourthTriplet(0), physicFirstQuad(0),physicSecondQuad(0),physicThirdQuad(0),physicFourthQuad(0),
0086 solidExternalChamber(0),logicExternalChamber(0),physicExternalChamber(0),
0087 solidInternalChamber(0),logicInternalChamber(0),physicInternalChamber(0),
0088 solidCollimator(0),logicCollimator(0),physicCollimator(0),
0089 solidCollimatorHole(0),logicCollimatorHole(0),physicCollimatorHole(0),
0090 solidFinalCollimator(0), logicFinalCollimator(0),physicFinalCollimator(0),
0091 solidFinalCollimatorHole(0),logicFinalCollimatorHole(0),physicFinalCollimatorHole(0),
0092 solidExternalMagnet_1(0),logicExternalMagnet_1(0),physicExternalMagnet_1(0), physicExternalMagnet_1Down(0),
0093 solidMagnet_1(0),logicMagnet_1(0),physicMagnet_1Right(0),physicMagnet_1Left(0), solidExternalMagnet_2(0),logicExternalMagnet_2(0),
0094 physicExternalMagnet_2(0),physicExternalMagnet_2Down(0),solidMagnet_2(0),logicMagnet_2(0),physicMagnet_2Right(0),physicMagnet_2Left(0), solidExternalMagnet_3(0),logicExternalMagnet_3(0),physicExternalMagnet_3(0),physicExternalMagnet_3Down(0),
0095 solidMagnet_3(0),logicMagnet_3(0),physicMagnet_3Right(0),physicMagnet_3Left(0),
0096 solidExternalMagnet_4(0),logicExternalMagnet_4(0),physicExternalMagnet_4(0),physicExternalMagnet_4Down(0),
0097 solidMagnet_4(0),logicMagnet_4(0),physicMagnet_4Right(0),physicMagnet_4Left(0),
0098 solidExternalSlit(0), logicExternalSlit(0), physicExternalSlit(0),
0099 solidInternalSlit(0),logicInternalSlit(0),physicInternalSlit(0),
0100 physicExitPipe(0),physicExitWindow(0),physicExithole(0),physicEntrancePipe(0),physicEntrancehole(0)
0101 {
0102     laserDrivenMessenger = new LaserDrivenBeamLineMessenger(this);
0103     
0104     //***************************** PW ***************************************
0105     
0106     static G4String ROGeometryName = "DetectorROGeometry";
0107     RO = new HadrontherapyDetectorROGeometry(ROGeometryName);
0108     
0109     G4cout << "Going to register Parallel world...";
0110     RegisterParallelWorld(RO);
0111     G4cout << "... done" << G4endl;
0112     //***************************** PW ***************************************
0113 }
0114 
0115 /////////////////////////////////////////////////////////////////////////////
0116 LaserDrivenBeamLine::~LaserDrivenBeamLine()
0117 {
0118     //delete laserDrivenMessenger;
0119     delete hadrontherapydetectorconstruction;
0120 }
0121 
0122 /////////////////////////////////////////////////////////////////////////////
0123 G4VPhysicalVolume* LaserDrivenBeamLine::Construct()
0124 {
0125     // Sets default geometry and materials
0126     SetDefaultDimensions();
0127     
0128     // Construct the energyselector (magnetic part and slit) and detector plane
0129     ConstructLaserDrivenBeamLine();
0130     
0131     //***************************** PW ***************************************
0132     if (!hadrontherapydetectorconstruction)
0133         
0134         //***************************** PW ***************************************
0135         
0136         // HadrontherapyDetectorConstruction builds ONLY the phantom and the detector with its associated ROGeometry
0137         hadrontherapydetectorconstruction = new HadrontherapyDetectorConstruction(physicTreatmentRoom);
0138     G4cout<<"HadrontherapyDetectorConstruction"<<G4endl;
0139     //***************************** PW ***************************************
0140     
0141     hadrontherapydetectorconstruction->InitializeDetectorROGeometry(RO,hadrontherapydetectorconstruction->GetDetectorToWorldPosition());
0142     
0143     //***************************** PW ***************************************
0144     return physicTreatmentRoom;
0145 }
0146 /////////////////////////////////////////////////////////////////////////////
0147 void LaserDrivenBeamLine::SetDefaultDimensions()
0148 {
0149     ///////////////////////////////////////////////////////////////////////
0150     // Definition of the colour sets
0151     white = new G4VisAttributes( G4Colour(1.,1.,1., 0.2));
0152     white -> SetVisibility(true);
0153     white -> SetForceSolid(true);
0154     white -> SetForceWireframe(true);
0155     
0156     blue = new G4VisAttributes(G4Colour(0. ,0. ,1.));
0157     blue -> SetVisibility(true);
0158     //blue -> SetForceSolid(true);
0159     
0160     gray = new G4VisAttributes( G4Colour(0.5, 0.5, 0.5, 0.5 ));
0161     gray-> SetVisibility(true);
0162     gray-> SetForceSolid(true);
0163     
0164     red = new G4VisAttributes(G4Colour(1. ,0. ,0., 0.2));
0165     red-> SetVisibility(true);
0166     red-> SetForceSolid(true);
0167     //red -> SetForceWireframe(true);
0168     
0169     yellow = new G4VisAttributes(G4Colour(1., 1., 0., 0.2));
0170     yellow-> SetVisibility(true);
0171     yellow-> SetForceSolid(true);
0172     
0173     green = new G4VisAttributes( G4Colour(25/255. , 255/255. ,  25/255., 0.4));
0174     green -> SetVisibility(true);
0175     green -> SetForceWireframe(true);
0176     green -> SetForceSolid(true);
0177     
0178     black = new G4VisAttributes( G4Colour(255/255. , 255/255.,  255/255.));
0179     black -> SetVisibility(true);
0180     black -> SetForceSolid(true);
0181     
0182     darkGreen = new G4VisAttributes( G4Colour(0/255. , 100/255. ,  0/255.));
0183     darkGreen -> SetVisibility(true);
0184     darkGreen -> SetForceSolid(true);
0185     
0186     darkOrange3 = new G4VisAttributes( G4Colour(205/255. , 102/255. ,  000/255., 0.7));
0187     darkOrange3 -> SetVisibility(true);
0188     darkOrange3 -> SetForceSolid(true);
0189     
0190     skyBlue = new G4VisAttributes( G4Colour(135/255. , 206/255. ,  235/255., 0.1));
0191     skyBlue -> SetVisibility(true);
0192     skyBlue -> SetForceSolid(true);
0193     
0194     // DEFAULT DIMENSIONS AND POSITIONS ARE PROVIDED HERE.
0195     /////////////////////// Exit Window ///////////////////////////////////////////////
0196     G4double defaultInnerRadiusExitWindow=0. *mm;
0197     InnerRadiusExitWindow=defaultInnerRadiusExitWindow;
0198     
0199     G4double defaultExternalRadiusExitWindow=55*mm;
0200     ExternalRadiusExitWindow=defaultExternalRadiusExitWindow;
0201     
0202     G4double defaultExitWindowThickness=25 *um;
0203     ExitWindowThickness=defaultExitWindowThickness;
0204     
0205     G4double defaultExitWindowXPosition=-ExitWindowThickness/2.;
0206     ExitWindowXPosition=defaultExitWindowXPosition;
0207     
0208     G4double defaultExitWindowYPosition=0.;
0209     ExitWindowYPosition=defaultExitWindowYPosition;
0210     
0211     G4double defaultExitWindowZPosition=0.0*mm;
0212     ExitWindowZPosition=defaultExitWindowZPosition;
0213     
0214     G4double defaultStartAngleExitWindow = 0.0 *deg;
0215     startAngleExitWindow = defaultStartAngleExitWindow;
0216     
0217     G4double defaultSpanningAngleExitWindow = 360.*deg;
0218     spanningAngleExitWindow = defaultSpanningAngleExitWindow;
0219     ////////////////////////////// Exit pipe ////////////////////////////////
0220     G4double defaultExitPipeheight=105. *mm;
0221     ExitPipeheight=defaultExitPipeheight;
0222     
0223     G4double defaultInnerRadiusExitPipe=50. *mm;
0224     InnerRadiusExitPipe=defaultInnerRadiusExitPipe;
0225     
0226     G4double defaultExternalRadiusExitPipe=55 *mm;
0227     ExternalRadiusExitPipe=defaultExternalRadiusExitPipe;
0228     
0229     G4double defaultExitPipeXPosition=-ExitPipeheight/2-ExitWindowThickness;
0230     ExitPipeXPosition=defaultExitPipeXPosition;
0231     
0232     G4double defaultExitPipeYPosition=0;
0233     ExitPipeYPosition=defaultExitPipeYPosition;
0234     
0235     G4double defaultExitPipeZPosition=0.0*mm;
0236     ExitPipeZPosition=defaultExitPipeZPosition;
0237     
0238     G4double defaultStartAngleExitPipe = 0.0 *deg;
0239     startAngleExitPipe = defaultStartAngleExitPipe;
0240     
0241     G4double defaultSpanningAngleExitPipe = 360.*deg;
0242     spanningAngleExitPipe = defaultSpanningAngleExitPipe;
0243     //////////////////////////////////////////////// Vacuum chamber //////////////////////////////
0244     G4double defaultExternalChamberXSize = 79.6*cm;
0245     externalChamberXSize = defaultExternalChamberXSize;
0246     
0247     G4double defaultExternalChamberYSize = 50. *cm;
0248     externalChamberYSize = defaultExternalChamberYSize;
0249     
0250     G4double defaultExternalChamberZSize = 50. *cm;
0251     externalChamberZSize = defaultExternalChamberZSize;
0252     
0253     G4double defaultExternalChamberXPosition = -(externalChamberXSize/2.+ExitPipeheight/2.)+ ExitPipeXPosition;
0254     externalChamberXPosition = defaultExternalChamberXPosition;
0255     
0256     G4double defaultExternalChamberYPosition = 0.0 *mm;
0257     externalChamberYPosition = defaultExternalChamberYPosition;
0258     
0259     G4double defaultExternalChamberZPosition = 0.0 *mm;
0260     externalChamberZPosition = defaultExternalChamberZPosition;
0261     
0262     // Defaults of the internal chamber dimensions
0263     // The position of its center is in the center
0264     // of the internal chamber while the dimension are
0265     // authomatically calculated respect to the external chamber ones
0266     G4double defaultVaccumChamberWallThickness=5 *mm;
0267     VaccumChamberWallThickness=defaultVaccumChamberWallThickness;
0268     
0269     G4double defaultInternalChamberXSize =externalChamberXSize - 2*VaccumChamberWallThickness;
0270     internalChamberXSize = defaultInternalChamberXSize;
0271     
0272     G4double defaultInternalChamberYSize =externalChamberYSize - 2*VaccumChamberWallThickness;
0273     internalChamberYSize = defaultInternalChamberYSize;
0274     
0275     G4double defaultInternalChamberZSize = externalChamberZSize - 2*VaccumChamberWallThickness;
0276     internalChamberZSize = defaultInternalChamberZSize;
0277     /////////////////////// Exit hole in vessel ///////////////////////////////////////////////
0278     G4double defaultInnerRadiusExithole=0.*mm;
0279     InnerRadiusExithole=defaultInnerRadiusExithole;
0280     
0281     G4double defaultExternalRadiusExithole=50.*mm;
0282     ExternalRadiusExithole=defaultExternalRadiusExithole;
0283     
0284     G4double defaultExitholeThickness=VaccumChamberWallThickness;
0285     ExitholeThickness=defaultExitholeThickness;
0286     
0287     G4double defaultExitholeXPosition=(externalChamberXSize/2.-ExitholeThickness/2.);
0288     ExitholeXPosition=defaultExitholeXPosition;
0289     
0290     G4double defaultExitholeYPosition=0.;
0291     ExitholeYPosition=defaultExitholeYPosition;
0292     
0293     G4double defaultExitholeZPosition=0.*mm;
0294     ExitholeZPosition=defaultExitholeZPosition;
0295     
0296     G4double defaultStartAngleExithole = 0.0 *deg;
0297     startAngleExithole= defaultStartAngleExithole;
0298     
0299     G4double defaultSpanningAngleExithole = 360.*deg;
0300     spanningAngleExithole = defaultSpanningAngleExithole;
0301     /////////////////////////////////Final collimator //////////////////////////////
0302     // The Final Collimator is located after  the 4th magnet
0303     G4double defaultExitholeToFinalCollimator=70 *mm;
0304     ExitholeToFinalCollimator=defaultExitholeToFinalCollimator;
0305     
0306     defaultInnerRadiusFinalCollimator = 0.0 *mm;
0307     innerRadiusFinalCollimator = defaultInnerRadiusFinalCollimator;
0308     
0309     defaultOuterRadiusFinalCollimator = 2.50 *mm;
0310     outerRadiusFinalCollimator = defaultOuterRadiusFinalCollimator;
0311     
0312     defaultFinalCollimatorThickness = 3.0 *mm;
0313     FinalCollimatorThickness = defaultFinalCollimatorThickness;
0314     
0315     defaultStartAngleFinalCollimator = 0.0 *deg;
0316     startAngleFinalCollimator = defaultStartAngleFinalCollimator;
0317     
0318     defaultSpanningAngleFinalCollimator = 360.*deg;
0319     spanningAngleFinalCollimator = defaultSpanningAngleFinalCollimator;
0320     
0321     defaultFinalCollimatorXPosition = internalChamberXSize/2.-ExitholeToFinalCollimator-FinalCollimatorThickness/2.;
0322     collimatorFinalBox_XPosition=defaultFinalCollimatorXPosition;
0323     FinalcollimatorXPosition = 0.0*mm; //HOLE IN THE FINAL COLLIMATOR
0324     
0325     defaultFinalCollimatorYPosition = 0.0*mm;
0326     collimatorFinalBox_YPosition=defaultFinalCollimatorYPosition;
0327     FinalcollimatorYPosition = defaultFinalCollimatorYPosition;
0328     
0329     defaultFinalCollimatorZPosition = 0.0*mm;
0330     collimatorFinalBox_ZPosition=0.0*mm;
0331     FinalcollimatorZPosition =defaultFinalCollimatorZPosition;
0332     
0333     defaultThicknessCollimator =3.0 *mm;
0334     collimatorFinalBoxXSize=defaultFinalCollimatorThickness;
0335     collimatorFinalBoxYSize=82.0*mm;
0336     collimatorFinalBoxZSize=210.0*mm;
0337     //////////////////ooooooooooOOOOOOOO000000000000OOOOOOOOOOOOooooooooooo/////////////////
0338     
0339     //Magnet characteristics
0340     G4double defaultExternalMagnet_XSize = 88.0*mm;
0341     G4double defaultExternalMagnet_YSizeTotal=87.*mm;
0342     G4double defaultInternalMagnet_YSize = 10. *mm;
0343     G4double defaultExternalMagnet_YSize =(defaultExternalMagnet_YSizeTotal-defaultInternalMagnet_YSize)/2.;
0344     G4double defaultExternalMagnet_ZSize = 104 *mm;
0345     
0346     G4double defaultExternalMagnet_YPosition =defaultInternalMagnet_YSize/2.+defaultExternalMagnet_YSize/2.;
0347     G4double defaultExternalMagnet_ZPosition = 0.0 *mm;
0348     
0349     G4double defaultMagnet_XSize=defaultExternalMagnet_XSize;
0350     G4double defaultMagnet_YSize=defaultExternalMagnet_YSizeTotal;
0351     G4double defaultMagnet_ZSize=19*mm;
0352     
0353     // Defaults of the external part of the magnet 4:
0354     G4double defaultFinalCollimatorToMagnet4=25.*mm;
0355     FinalCollimatorToMagnet4=defaultFinalCollimatorToMagnet4;
0356     
0357     externalMagnet_4XSize = defaultExternalMagnet_XSize;
0358     externalMagnet_4YSize = defaultExternalMagnet_YSize;
0359     externalMagnet_4ZSize = defaultExternalMagnet_ZSize;
0360     
0361     Magnet_4XSize=defaultMagnet_XSize;
0362     Magnet_4YSize=defaultMagnet_YSize;
0363     Magnet_4ZSize=defaultMagnet_ZSize;
0364     
0365     G4double defaultExternalMagnet_4XPosition = -(FinalCollimatorThickness/2.+FinalCollimatorToMagnet4+defaultExternalMagnet_XSize/2.)+ collimatorFinalBox_XPosition;
0366     externalMagnet_4XPosition = defaultExternalMagnet_4XPosition;
0367     
0368     externalMagnet_4YPosition = defaultExternalMagnet_YPosition;
0369     externalMagnet_4ZPosition = defaultExternalMagnet_ZPosition;
0370     
0371     Magnet_4XPosition=externalMagnet_4XPosition;
0372     Magnet_4YPosition=0.0*mm;
0373     Magnet_4ZPosition=(defaultExternalMagnet_ZSize+defaultMagnet_ZSize)/2.;
0374     //////////////////ooooooooooOOOOOOOO000000000000OOOOOOOOOOOOooooooooooo/////////////////
0375     // Defaults of the external part of the magnet 3:
0376     externalMagnet_3XSize = defaultExternalMagnet_XSize;
0377     externalMagnet_3YSize = defaultExternalMagnet_YSize;
0378     externalMagnet_3ZSize = defaultExternalMagnet_ZSize;
0379     
0380     Magnet_3XSize=defaultMagnet_XSize;
0381     Magnet_3YSize=defaultMagnet_YSize;
0382     Magnet_3ZSize=defaultMagnet_ZSize;
0383     
0384     G4double defaultMagnet4ToMagnet3=65.*mm; //85.*mm ANTONELLA
0385     Magnet4ToMagnet3=defaultMagnet4ToMagnet3;
0386     
0387     G4double defaultExternalMagnet_3XPosition =-(Magnet4ToMagnet3+defaultExternalMagnet_XSize/2.+defaultExternalMagnet_XSize/2.)+externalMagnet_4XPosition;
0388     externalMagnet_3XPosition = defaultExternalMagnet_3XPosition;
0389     
0390     externalMagnet_3YPosition =defaultExternalMagnet_YPosition;
0391     externalMagnet_3ZPosition = defaultExternalMagnet_ZPosition;
0392     
0393     
0394     Magnet_3XPosition=externalMagnet_3XPosition;
0395     Magnet_3YPosition=0.0*mm;
0396     Magnet_3ZPosition=(defaultExternalMagnet_ZSize+defaultMagnet_ZSize)/2.;
0397     //////////////////ooooooooooOOOOOOOO000000000000OOOOOOOOOOOOooooooooooo/////////////////
0398     // Defaults of the external part of the magnet 2:
0399     externalMagnet_2XSize = defaultExternalMagnet_XSize;
0400     externalMagnet_2YSize = defaultExternalMagnet_YSize;
0401     externalMagnet_2ZSize = defaultExternalMagnet_ZSize;
0402     
0403     Magnet_2XSize=defaultMagnet_XSize;
0404     Magnet_2YSize=defaultMagnet_YSize;
0405     Magnet_2ZSize=defaultMagnet_ZSize;
0406     
0407     G4double defaultMagnet3ToMagnet2=10 *mm;
0408     Magnet3ToMagnet2=defaultMagnet3ToMagnet2;
0409     
0410     G4double defaultExternalMagnet_2XPosition =-(Magnet3ToMagnet2+defaultExternalMagnet_XSize/2.+defaultExternalMagnet_XSize/2.)+externalMagnet_3XPosition;
0411     externalMagnet_2XPosition = defaultExternalMagnet_2XPosition;
0412     
0413     externalMagnet_2YPosition = defaultExternalMagnet_YPosition;
0414     externalMagnet_2ZPosition = defaultExternalMagnet_ZPosition;
0415     
0416     Magnet_2XPosition=externalMagnet_2XPosition;
0417     Magnet_2YPosition=0.0*mm;
0418     Magnet_2ZPosition=(defaultExternalMagnet_ZSize+defaultMagnet_ZSize)/2.;
0419     //////////////////ooooooooooOOOOOOOO000000000000OOOOOOOOOOOOooooooooooo/////////////////
0420     // Defaults of the external part of the magnet 1:
0421     externalMagnet_1XSize=defaultExternalMagnet_XSize;
0422     externalMagnet_1YSize = defaultExternalMagnet_YSize;
0423     externalMagnet_1ZSize = defaultExternalMagnet_ZSize;
0424     
0425     Magnet_1XSize=defaultMagnet_XSize;
0426     Magnet_1YSize=defaultMagnet_YSize;
0427     Magnet_1ZSize=defaultMagnet_ZSize;
0428     
0429     G4double defaultMagnet2ToMagnet1=85 *mm;
0430     Magnet2ToMagnet1=defaultMagnet2ToMagnet1;
0431     
0432     G4double defaultExternalMagnet_1XPosition = -(Magnet2ToMagnet1+defaultExternalMagnet_XSize/2.+defaultExternalMagnet_XSize/2.)+externalMagnet_2XPosition;
0433     externalMagnet_1XPosition = defaultExternalMagnet_1XPosition;
0434     
0435     externalMagnet_1YPosition = defaultExternalMagnet_YPosition;
0436     externalMagnet_1ZPosition = defaultExternalMagnet_ZPosition;
0437     
0438     Magnet_1XPosition=defaultExternalMagnet_1XPosition;
0439     Magnet_1YPosition=0.0*mm;
0440     Magnet_1ZPosition=(defaultExternalMagnet_ZSize+defaultMagnet_ZSize)/2.;
0441     
0442     // Defaults of the external part of the Slit
0443     G4double defaultExternalSlitXSize = 8.0 *mm;
0444     externalSlitXSize = defaultExternalSlitXSize;
0445     
0446     G4double defaultExternalSlitYSize = 82. *mm;
0447     externalSlitYSize = defaultExternalSlitYSize;
0448     
0449     G4double defaultExternalSlitZSize = 210. *mm;
0450     externalSlitZSize = defaultExternalSlitZSize;
0451     
0452     G4double defaultExternalSlitXPosition = -(Magnet3ToMagnet2/2.+defaultExternalMagnet_XSize/2.)+externalMagnet_3XPosition;
0453     externalSlitXPosition = defaultExternalSlitXPosition;
0454     
0455     G4double defaultExternalSlitYPosition = 0.0 *mm;
0456     externalSlitYPosition = defaultExternalSlitYPosition;
0457     
0458     G4double defaultExternalSlitZPosition = 0.0 *mm;
0459     externalSlitZPosition = defaultExternalSlitZPosition;
0460     
0461     // Defaults of the internal part of the Slit:
0462     internalSlitXSize = defaultExternalSlitXSize;
0463     
0464     G4double defaultInternalSlitYSize = 3 *mm;
0465     internalSlitYSize = defaultInternalSlitYSize;
0466     
0467     G4double defaultInternalSlitZSize = 3 *mm;
0468     internalSlitZSize = defaultInternalSlitZSize;
0469     
0470     G4double defaultInternalSlitXPosition = 0.0 *mm;
0471     internalSlitXPosition = defaultInternalSlitXPosition;
0472     
0473     G4double defaultInternalSlitYPosition = 0.0 *mm;
0474     internalSlitYPosition = defaultInternalSlitYPosition;
0475     
0476     G4double defaultInternalSlitZPosition = 40.0 *mm;
0477     internalSlitZPosition = defaultInternalSlitZPosition;
0478     
0479     // Defaults of the particle collimator (First collimator).
0480     // The Collimator should be located before the 1st magnet
0481     //
0482     defaultInnerRadiusCollimator = 0.0 *mm;
0483     innerRadiusCollimator = defaultInnerRadiusCollimator;
0484     
0485     defaultOuterRadiusCollimator = 2.5 *mm;
0486     outerRadiusCollimator = defaultOuterRadiusCollimator;
0487     
0488     thicknessCollimator = defaultThicknessCollimator;
0489     
0490     defaultStartAngleCollimator = 0.0 *deg;
0491     startAngleCollimator = defaultStartAngleCollimator;
0492     
0493     defaultSpanningAngleCollimator = 360.*deg;
0494     spanningAngleCollimator = defaultSpanningAngleCollimator;
0495     
0496     G4double defultMagnet1ToFirstCollimator=25.*mm;
0497     Magnet1ToFirstCollimator=defultMagnet1ToFirstCollimator;
0498     
0499     defaultCollimatorXPosition = -(thicknessCollimator/2.+Magnet1ToFirstCollimator+defaultExternalMagnet_XSize/2.)+externalMagnet_1XPosition;
0500     collimatorBox_XPosition=defaultCollimatorXPosition;
0501     collimatorXPosition = 0.0*mm;
0502     
0503     defaultCollimatorYPosition = 0.0*mm;
0504     collimatorBox_YPosition=defaultCollimatorYPosition;
0505     collimatorYPosition = 0.0*mm;
0506     
0507     defaultCollimatorZPosition = 0.0*mm;
0508     collimatorBox_ZPosition=defaultCollimatorZPosition;
0509     collimatorZPosition = 0.*mm;
0510     
0511     collimatorBoxYSize=82.0* mm;
0512     collimatorBoxZSize=210.0* mm;
0513     
0514     //////////////////// Entrance Hole //////////////////////////////////
0515     G4double defaultInnerRadiusEntrancehole=0. *mm;
0516     InnerRadiusEntrancehole=defaultInnerRadiusEntrancehole;
0517     
0518     G4double defaultExternalRadiusEntrancehole=50.*mm;
0519     ExternalRadiusEntrancehole=defaultExternalRadiusEntrancehole;
0520     
0521     G4double defaultEntranceholeThickness=VaccumChamberWallThickness;
0522     EntranceholeThickness=defaultEntranceholeThickness;
0523     
0524     G4double defaultEntranceholeXPosition=-(externalChamberXSize/2.-EntranceholeThickness/2.);
0525     EntranceholeXPosition=defaultEntranceholeXPosition;
0526     
0527     G4double defaultEntranceholeQuadXPosition=+(externalChamberXSize/2.-EntranceholeThickness/2.);
0528     EntranceholeQuadXPosition=defaultEntranceholeQuadXPosition;
0529     
0530     G4double defaultEntranceholeYPosition=0.;
0531     EntranceholeYPosition=defaultEntranceholeYPosition;
0532     
0533     G4double defaultEntranceholeZPosition=0.0*mm;
0534     EntranceholeZPosition=defaultEntranceholeZPosition;
0535     
0536     G4double defaultStartAngleEntrancehole= 0.0 *deg;
0537     startAngleEntrancehole= defaultStartAngleEntrancehole;
0538     
0539     G4double defaultSpanningAngleEntrancehole= 360.*deg;
0540     spanningAngleEntrancehole=defaultSpanningAngleEntrancehole;
0541     
0542     ///////////////// Entrance Pipe/////////////////////////////////////////////
0543     
0544     G4double defaultEntrancePipeheight=105. *mm;
0545     EntrancePipeheight=defaultEntrancePipeheight;
0546     
0547     G4double defaultInnerRadiusEntrancePipe=50. *mm;
0548     InnerRadiusEntrancePipe=defaultInnerRadiusEntrancePipe;
0549     
0550     G4double defaultExternalRadiusEntrancePipe=55 *mm;
0551     ExternalRadiusEntrancePipe=defaultExternalRadiusEntrancePipe;
0552     
0553     G4double defaultEntrancePipeXPosition=-EntrancePipeheight/2-externalChamberXSize/2+externalChamberXPosition;
0554     EntrancePipeXPosition=defaultEntrancePipeXPosition;
0555     
0556     G4double defaultEntrancePipeYPosition=0;
0557     EntrancePipeYPosition=defaultEntrancePipeYPosition;
0558     
0559     G4double defaultEntrancePipeZPosition=0.0*mm;
0560     EntrancePipeZPosition=defaultEntrancePipeZPosition;
0561     
0562     G4double defaultStartAngleEntrancePipe= 0.0 *deg;
0563     startAngleEntrancePipe= defaultStartAngleEntrancePipe;
0564     
0565     G4double defaultSpanningAngleEntrancePipe= 360.*deg;
0566     spanningAngleEntrancePipe=defaultSpanningAngleEntrancePipe;
0567     
0568     /////////////////////////////////////Quadrupole//////////////////////////////////
0569     G4double defaultQuadChamberWallPosX=-(externalChamberXSize/2.)-EntrancePipeheight/2.+EntrancePipeXPosition;
0570     QuadChamberWallPosX=defaultQuadChamberWallPosX;
0571     G4double defaultQuadChamberWallPosY=0.0*cm;
0572     QuadChamberWallPosY=defaultQuadChamberWallPosY;
0573     G4double defaultQuadChamberWallPosZ=0.0*cm;
0574     QuadChamberWallPosZ=defaultQuadChamberWallPosZ;
0575     
0576     G4double defaultInnerRadiusQuad=10.0*mm;
0577     InnerRadiusQuad=defaultInnerRadiusQuad;
0578     
0579     G4double defaultInnerRadiusTriplet=0.0*mm;
0580     InnerRadiusTriplet=defaultInnerRadiusTriplet;
0581     
0582     G4double defaultExternalRadiusQuad=30.0*mm;
0583     ExternalRadiusQuad=defaultExternalRadiusQuad;
0584     
0585     G4double defaultFirstQuadThickness=80.0*mm;
0586     FirstQuadThickness=defaultFirstQuadThickness;
0587     G4double defaultSecondQuadThickness=40.0*mm;
0588     SecondQuadThickness=defaultSecondQuadThickness;
0589     G4double defaultThirdQuadThickness=40.0*mm;
0590     ThirdQuadThickness=defaultThirdQuadThickness;
0591     G4double defaultFourthQuadThickness=80.0*mm;
0592     FourthQuadThickness=defaultFourthQuadThickness;
0593     
0594     G4double defaultStartAngleQuad = 0.0 *deg;
0595     startAngleQuad = defaultStartAngleQuad;
0596     
0597     G4double defaultSpanningAngleQuad = 360.*deg;
0598     spanningAngleQuad = defaultSpanningAngleQuad;
0599     
0600     G4double distancefromQuadChamber=100.0*mm;
0601     G4double defaultFourthQuadXPosition= internalChamberXSize/2.-distancefromQuadChamber-FourthQuadThickness/2.;
0602     FourthQuadXPosition=defaultFourthQuadXPosition;
0603     FourthQXPosition=0.0*mm;
0604     
0605     G4double distanceFQuadTQuad=100.0*mm;
0606     G4double defaultThirdQuadXPosition=-ThirdQuadThickness/2.-distanceFQuadTQuad-FourthQuadThickness/2.+FourthQuadXPosition;
0607     ThirdQuadXPosition=defaultThirdQuadXPosition;
0608     ThirdQXPosition=0.0*mm;
0609     
0610     G4double distanceTQuadSQuad=100.0*mm;
0611     G4double defaultSecondQuadXPosition=-SecondQuadThickness/2.-distanceTQuadSQuad-ThirdQuadThickness/2.+ThirdQuadXPosition;
0612     SecondQuadXPosition=defaultSecondQuadXPosition;
0613     SecondQXPosition=0.0*mm;
0614     
0615     G4double distanceSQuadFQuad=100.0*mm;
0616     G4double defaultFirstQuadXPosition=-FirstQuadThickness/2.-distanceSQuadFQuad-SecondQuadThickness/2.+SecondQuadXPosition;
0617     FirstQuadXPosition=defaultFirstQuadXPosition;
0618     FirstQXPosition=0.0*mm;
0619     
0620     G4double defaultQuadYPosition=0.0*mm;
0621     QuadYPosition=defaultQuadYPosition;
0622     QYPosition=defaultQuadYPosition;
0623     
0624     G4double defaultQuadTZPosition= 0.*mm;
0625     QuadZPosition=defaultQuadTZPosition;
0626     G4double defaultQuadZPosition=0.0*mm;
0627     QZPosition=defaultQuadZPosition;
0628     
0629     // DEFAULT DEFINITION OF THE MATERIALS
0630     // All elements and compound definition follows the NIST database
0631     
0632     //ELEMENTS
0633     G4bool isotopes = false;
0634     G4Element* zincNist = G4NistManager::Instance()->FindOrBuildElement("Zn");
0635     G4Element* copperNist = G4NistManager::Instance()->FindOrBuildElement("Cu");
0636     
0637     //COMPOUNDS
0638     G4Material* ironNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Fe", isotopes);
0639     G4Material* aluminiumNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_Al");
0640     G4Material* kaptonNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_KAPTON", isotopes);
0641     //G4Material* waterNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER", isotopes);
0642     G4Material* stainless_steelNist = G4NistManager::Instance()->FindOrBuildMaterial("G4_STAINLESS-STEEL", isotopes);
0643     
0644     // Elements and compunds not pre-defined in Geant4
0645     G4double d; // Density
0646     G4int nComponents;// Number of components
0647     G4double fractionmass; // Fraction in mass of an element in a material
0648     d = 8.40*g/cm3;
0649     nComponents = 2;
0650     G4Material* brass = new G4Material("Brass", d, nComponents);
0651     brass -> AddElement(zincNist, fractionmass = 30 *perCent);
0652     brass -> AddElement(copperNist, fractionmass = 70 *perCent);
0653     
0654     G4double atomicNumber = 1.;
0655     G4double massOfMole = 1.008*g/mole;
0656     d = 1.e-25*g/cm3;
0657     G4double temperature = 2.73*kelvin;
0658     G4double pressure = 3.e-18*pascal;
0659     G4Material* vacuum = new G4Material("interGalactic", atomicNumber,massOfMole, d, kStateGas,temperature, pressure);
0660     
0661     //***************************** PW ***************************************
0662     
0663     // DetectorROGeometry Material
0664     new G4Material("dummyMat", 1., 1.*g/mole, 1.*g/cm3);
0665     
0666     //***************************** PW ***************************************
0667     
0668     // MATERIAL ASSIGNMENT
0669     MotherMaterial=vacuum;
0670     QuadMaterial=ironNist;
0671     externalChamberMaterial = stainless_steelNist;
0672     internalChamberMaterial = vacuum;
0673     collimatorMaterial = aluminiumNist;
0674     collimatorHoleMaterial=vacuum;
0675     FinalcollimatorMaterial=aluminiumNist;
0676     FinalcollimatorHoleMaterial=vacuum;
0677     WindowMaterial=kaptonNist;
0678     PipeMaterial=stainless_steelNist;
0679     
0680     externalMagnet_1Material = ironNist;
0681     externalMagnet_2Material = ironNist;
0682     externalMagnet_3Material = ironNist;
0683     externalMagnet_4Material = ironNist;
0684     
0685     externalSlitMaterial = brass;
0686     internalSlitMaterial =vacuum;
0687     
0688     //FC Material
0689     
0690     KaptonEntranceWindowMaterial=kaptonNist;
0691     GuardRingMaterial=stainless_steelNist;
0692     FaradayCupBottomMaterial=aluminiumNist;
0693     CupMaterial=FaradayCupBottomMaterial;
0694     MassRingMaterial=GuardRingMaterial;
0695     
0696 }
0697 
0698 /////////////////////////////////////////////////////////////////////////////
0699 void LaserDrivenBeamLine::ConstructLaserDrivenBeamLine()
0700 {
0701     // -----------------------------
0702     // Treatment room - World volume
0703     //------------------------------
0704     
0705     const G4double worldX = 800.0 *cm;
0706     const G4double worldY = 400.0 *cm;
0707     const G4double worldZ = 400.0 *cm;
0708     
0709     solidTreatmentRoom = new G4Box("TreatmentRoom",
0710                                    worldX,
0711                                    worldY,
0712                                    worldZ);
0713     
0714     logicTreatmentRoom = new G4LogicalVolume(solidTreatmentRoom,
0715                                              MotherMaterial,
0716                                              "logicTreatmentRoom",
0717                                              0,
0718                                              0,
0719                                              0);
0720     
0721     physicTreatmentRoom = new G4PVPlacement(0,
0722                                             G4ThreeVector(),
0723                                             "physicalTreatmentRoom",
0724                                             logicTreatmentRoom,
0725                                             0,
0726                                             false,
0727                                             0);
0728     
0729     
0730     // The treatment room is invisible in the Visualisation
0731     logicTreatmentRoom -> SetVisAttributes (G4VisAttributes::GetInvisible());
0732     
0733     // The various components of the energyselector are constructed calling
0734     // the following methods
0735     
0736     // This method constructs the chamber where the energyselector is located
0737     EnergySelectorChamber();
0738     // This method construct the exit window
0739     ExitWindow();
0740     // This method construct the exit pipe
0741     ExitPipe();
0742     // This method construct the exit hole
0743     Exithole();
0744     
0745     // This method constructs a circular collimator of variable thickness and
0746     // aperture. It is placed befor the magnet to collimate particles caming from the
0747     // plasma;
0748     Collimator();
0749     
0750     // This method constructs the magnet 1 and its associated magnetic field
0751     Magnet_1();
0752     
0753     // This method constructs the magnet 2 and its associated magnetic field
0754     Magnet_2();
0755     
0756     // This method constructs the magnet 3 and its associated magnetic field
0757     Magnet_3();
0758     
0759     // This method constructs the magnet 4 and its associated magnetic field
0760     Magnet_4();
0761     
0762     // The selection slit is a square hole moveable inside a metallic plate
0763     Slit();
0764     
0765     FinalCollimator();
0766     
0767     // This method construct the quadrupoles
0768     Quadrupole();
0769     // This method construct the entrance hole
0770     Entrancehole();
0771     // This method construct the entrance pipe
0772     EntrancePipe();
0773     
0774     FaradayCup();
0775     
0776 }
0777 
0778 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0779 
0780 void LaserDrivenBeamLine::ConstructSDandField()
0781 {
0782     G4double minEps=1.0e-5;  //   Minimum & value for smallest steps
0783     G4double maxEps=1.0e-4;
0784     G4bool allLocal = true;
0785     // G4int nvar = 8;  For pure magnetic field, the number of integration variables is the default!
0786     
0787     //....oooOO0OOooo..........ENERGY SELECTOR SYSTEM FIELD..........oooOO0OOooo....
0788     if(logicInternalChamber){G4double xOffset =(internalChamberXSize/2.0)+externalSlitXPosition;
0789         PurgMagField = new HadrontherapyMagneticField3D("field/ESSMagneticField.TABLE", xOffset);
0790         pFieldMgr =new G4FieldManager();
0791         pFieldMgr -> SetDetectorField(PurgMagField);
0792         G4cout << "DeltaStep "<< pFieldMgr -> GetDeltaOneStep()/mm <<"mm" <<G4endl;
0793         pFieldMgr -> CreateChordFinder(PurgMagField);
0794         fEquation = new G4Mag_UsualEqRhs(PurgMagField);
0795         fstepper = new G4ClassicalRK4(fEquation);
0796         //////fstepper = new G4HelixImplicitEuler(fEquation);
0797         pIntgrDriver = new G4MagInt_Driver(1*mm,fstepper,fstepper-> GetNumberOfVariables());
0798         //the first parameter is the minimum step
0799         pChordFinder = new G4ChordFinder(pIntgrDriver);
0800         pFieldMgr->SetChordFinder(pChordFinder);
0801         pFieldMgr->SetMinimumEpsilonStep(minEps);
0802         pFieldMgr->SetMaximumEpsilonStep(maxEps);
0803         pFieldMgr->SetDeltaOneStep(0.5e-3*mm);//default value of DeltaChord is 0.25 mm
0804         logicInternalChamber -> SetFieldManager(pFieldMgr, allLocal);}
0805     //....oooOO0OOooo..........QUADS FIELDS..........oooOO0OOooo....
0806     //....oooOO0OOooo..........FOURTH QUAD FIELD..........oooOO0OOooo....
0807     if(LFourthTriplet){G4double xOffsetFQ =-(QuadChamberWallPosX+FourthQuadXPosition);
0808         PurgMagFieldQuadFourth = new HadrontherapyMagneticField3D("field/Quad80MagneticField.TABLE", xOffsetFQ);
0809         pFieldMgrQuadFourth =  new G4FieldManager();
0810         pFieldMgrQuadFourth -> SetDetectorField(PurgMagFieldQuadFourth);
0811         
0812         pFieldMgrQuadFourth -> CreateChordFinder(PurgMagFieldQuadFourth);
0813         fEquationQuadFourth = new G4Mag_UsualEqRhs(PurgMagFieldQuadFourth);
0814         fstepperQuadFourth = new G4ClassicalRK4(fEquationQuadFourth);
0815         pIntgrDriverQuadFourth = new G4MagInt_Driver(1*mm,fstepperQuadFourth,fstepperQuadFourth-> GetNumberOfVariables());
0816         //the first parameter is the minimum step
0817         pChordFinderQuadFourth = new G4ChordFinder(pIntgrDriverQuadFourth);
0818         pFieldMgrQuadFourth->SetChordFinder(pChordFinderQuadFourth);
0819         pFieldMgrQuadFourth->SetMinimumEpsilonStep(minEps);
0820         pFieldMgrQuadFourth->SetMaximumEpsilonStep(maxEps);
0821         pFieldMgrQuadFourth->SetDeltaOneStep(0.5e-3*mm);//default value of DeltaChord is 0.25 mm
0822         LFourthTriplet -> SetFieldManager(pFieldMgrQuadFourth, allLocal);}
0823     //....oooOO0OOooo..........THIRD QUAD FIELD..........oooOO0OOooo....
0824     if(LThirdTriplet){ G4double xOffsetTQ =-(QuadChamberWallPosX+ThirdQuadXPosition);
0825         PurgMagFieldQuadThird = new HadrontherapyMagneticField3D("field/Quad40MagneticField.TABLE", xOffsetTQ);
0826         pFieldMgrQuadThird =  new G4FieldManager();
0827         pFieldMgrQuadThird -> SetDetectorField(PurgMagFieldQuadThird);
0828         pFieldMgrQuadThird -> CreateChordFinder(PurgMagFieldQuadThird);
0829         fEquationQuadThird = new G4Mag_UsualEqRhs(PurgMagFieldQuadThird);
0830         fstepperQuadThird = new G4ClassicalRK4(fEquationQuadThird);
0831         pIntgrDriverQuadThird = new G4MagInt_Driver(1*mm,fstepperQuadThird,fstepperQuadThird-> GetNumberOfVariables());
0832         //the first parameter is the minimum step
0833         pChordFinderQuadThird = new G4ChordFinder(pIntgrDriverQuadThird);
0834         pFieldMgrQuadThird->SetChordFinder(pChordFinderQuadThird);
0835         pFieldMgrQuadThird->SetMinimumEpsilonStep(minEps);
0836         pFieldMgrQuadThird->SetMaximumEpsilonStep(maxEps);
0837         pFieldMgrQuadThird->SetDeltaOneStep(0.5e-3*mm);//default value of DeltaChord is 0.25 mm
0838         LThirdTriplet -> SetFieldManager(pFieldMgrQuadThird, allLocal);}
0839     //....oooOO0OOooo..........SECOND QUAD FIELD..........oooOO0OOooo....
0840     if(LSecondTriplet){G4double xOffsetSQ =-(QuadChamberWallPosX+SecondQuadXPosition);
0841         PurgMagFieldQuadSecond = new HadrontherapyMagneticField3D("field/Quad40MagneticField.TABLE", xOffsetSQ);
0842         pFieldMgrQuadSecond =  new G4FieldManager();
0843         pFieldMgrQuadSecond -> SetDetectorField(PurgMagFieldQuadSecond);
0844         pFieldMgrQuadSecond -> CreateChordFinder(PurgMagFieldQuadSecond);
0845         fEquationQuadSecond = new G4Mag_UsualEqRhs(PurgMagFieldQuadSecond);
0846         fstepperQuadSecond = new G4ClassicalRK4(fEquationQuadSecond);
0847         pIntgrDriverQuadSecond = new G4MagInt_Driver(1*mm,fstepperQuadSecond,fstepperQuadSecond-> GetNumberOfVariables());
0848         //the first parameter is the minimum step
0849         pChordFinderQuadSecond = new G4ChordFinder(pIntgrDriverQuadSecond);
0850         pFieldMgrQuadSecond->SetChordFinder(pChordFinderQuadSecond);
0851         pFieldMgrQuadSecond->SetMinimumEpsilonStep(minEps);
0852         pFieldMgrQuadSecond->SetMaximumEpsilonStep(maxEps);
0853         pFieldMgrQuadSecond->SetDeltaOneStep(0.5e-3*mm);//default value of DeltaChord is 0.25 mm
0854         LSecondTriplet -> SetFieldManager(pFieldMgrQuadSecond, allLocal);}
0855     //....oooOO0OOooo..........FIRST QUAD FIELD..........oooOO0OOooo....
0856     if(LFirstTriplet) {G4double xOffsetFirstQ =-(QuadChamberWallPosX+FirstQuadXPosition);
0857         PurgMagFieldQuadFirst = new HadrontherapyMagneticField3D("field/Quad80MagneticField.TABLE", xOffsetFirstQ);
0858         pFieldMgrQuadFirst =  new G4FieldManager();
0859         pFieldMgrQuadFirst -> SetDetectorField(PurgMagFieldQuadFirst);
0860         pFieldMgrQuadFirst -> CreateChordFinder(PurgMagFieldQuadFirst);
0861         fEquationQuadFirst = new G4Mag_UsualEqRhs(PurgMagFieldQuadFirst);
0862         fstepperQuadFirst = new G4ClassicalRK4(fEquationQuadFirst);
0863         pIntgrDriverQuadFirst = new G4MagInt_Driver(1*mm,fstepperQuadFirst,fstepperQuadFirst-> GetNumberOfVariables());
0864         //the first parameter is the minimum step
0865         pChordFinderQuadFirst = new G4ChordFinder(pIntgrDriverQuadFirst);
0866         pFieldMgrQuadFirst->SetChordFinder(pChordFinderQuadFirst);
0867         pFieldMgrQuadFirst->SetMinimumEpsilonStep(minEps);
0868         pFieldMgrQuadFirst->SetMaximumEpsilonStep(maxEps);
0869         pFieldMgrQuadFirst->SetDeltaOneStep(0.5e-3*mm);//default value of DeltaChord is 0.25 mm
0870         LFirstTriplet -> SetFieldManager(pFieldMgrQuadFirst, allLocal);}
0871     //....oooOO0OOooo..........FARADAY CUP FIELD..........oooOO0OOooo....
0872     if(logicVirtualMag) {G4double exOffset= -20*cm;
0873         G4double eyOffset= 0*cm;
0874         G4double ezOffset= 0*cm;
0875         G4FieldManager *pEFieldmanager = new G4FieldManager();
0876         G4ElectricField *ElectricField = new HadrontherapyElectricTabulatedField3D("field/ElectricFieldFC-600V.TABLE", exOffset, eyOffset, ezOffset);
0877         // UNIFORM FIELD
0878         // G4ElectroMagneticField* ElectricField = new G4UniformElectricField(G4ThreeVector(0.0, 10.0*volt/m, 0.0)); //G4UniformElectricField
0879         // The following is only for global field in the whole geometry
0880         //pEFieldmanager = G4TransportationManager::GetTransportationManager() -> GetFieldManager();
0881         
0882         const G4int nvarElectric=8;  // The Equation of motion for Electric (or combined Electric/Magnetic)
0883         // field requires 8 integration variables
0884         
0885         G4EqMagElectricField *fLocalEquation = new G4EqMagElectricField(ElectricField);
0886         G4MagIntegratorStepper* fLocalStepper = new G4ClassicalRK4(fLocalEquation, nvarElectric);
0887         G4MagInt_Driver  *pIntgrDriver_E = new G4MagInt_Driver(0.02*mm, fLocalStepper, fLocalStepper -> GetNumberOfVariables() );
0888         G4ChordFinder *fLocalChordFinder = new G4ChordFinder(pIntgrDriver_E);
0889         pEFieldmanager -> SetDetectorField(ElectricField);
0890         pEFieldmanager -> SetChordFinder(fLocalChordFinder);
0891         //G4double deltainter=0.0001*mm;
0892         //G4double missdist=0.1*mm;
0893         //pEFieldmanager->SetDeltaIntersection(deltainter);
0894         //fLocalChordFinder->SetDeltaChord(missdist);
0895         pEFieldmanager->SetMinimumEpsilonStep(minEps);
0896         pEFieldmanager->SetMaximumEpsilonStep(maxEps);
0897         pEFieldmanager->SetDeltaOneStep( 0.5e-3 * mm );
0898         //pEFieldmanager -> SetFieldChangesEnergy(true);
0899         logicVirtualMag -> SetFieldManager(pEFieldmanager, allLocal);}
0900     //....oooOO0OOooo....................oooOO0OOooo....
0901     G4cout<<" //....oooOO0OOooo.......... FIELDS HAVE BEEN IMPLEMENTED..........oooOO0OOooo...."<<G4endl;
0902     return;
0903 }
0904 
0905 /////////////////////////////////////////////////////////////////////////////
0906 void LaserDrivenBeamLine::FaradayCup()
0907 {
0908     /// FC sizes ///
0909     
0910     G4double InnerRadiusFC=25*mm;
0911     G4double OuterRadiusFC=45*mm;
0912     G4double MassRingThickness=5*mm;
0913     G4double GuardRingThickness=180*mm;
0914     G4double FaradayCupBottomThickness=120*mm;
0915     G4double CupThickness=10*cm;
0916     G4double KaptonEntranceWindowThickness=25*um;
0917     
0918     /// Virtual Volumes ///
0919     
0920     G4double VirtualWindowThickness=1.*um ;
0921     G4double VirtualMiddleThickness= 1.*um ;
0922     G4double VirtualBottomThickness= 1. *um ;
0923     G4double VirtualOverBottomThickness=1. *um ;
0924     G4double VirtualLateralLength=FaradayCupBottomThickness+CupThickness+VirtualBottomThickness;
0925     
0926     
0927     //// Position ////
0928     
0929     G4double virtualMagPosX=31*cm;
0930     G4double FC_XOffset=20*cm;
0931     G4double KaptonEntranceWindowPosX=-virtualMagPosX+KaptonEntranceWindowThickness/2+FC_XOffset;
0932     G4double MassRingPosX=KaptonEntranceWindowPosX+KaptonEntranceWindowThickness/2+MassRingThickness/2;
0933     G4double VirtualWindowPosX=MassRingPosX+MassRingThickness/2+VirtualWindowThickness/2;
0934     G4double GuardRingPosX=MassRingPosX+MassRingThickness/2+GuardRingThickness/2+2*mm;
0935     G4double VirtualMiddlePosX=GuardRingPosX+GuardRingThickness/2+VirtualMiddleThickness/2;
0936     G4double FaradayCupBottomPosX=GuardRingPosX+GuardRingThickness/2+FaradayCupBottomThickness/2+1*cm;
0937     G4double VirtualBottomPosX=FaradayCupBottomPosX+FaradayCupBottomThickness/2+VirtualBottomThickness/2;
0938     G4double CupPosX=VirtualBottomPosX+VirtualBottomThickness/2+CupThickness/2;
0939     G4double VirtualOverBottomPosX=CupPosX+CupThickness/2+VirtualOverBottomThickness/2;
0940     G4double VirtualLateralPosX=GuardRingPosX+GuardRingThickness/2+1*cm+(FaradayCupBottomThickness+CupThickness+VirtualBottomThickness)/2;
0941     G4double phi = 90. *deg;
0942     G4RotationMatrix rm;
0943     rm.rotateY(phi);
0944     
0945     virtualMag= new G4Box("virtualMag", 31.*cm, 6*cm, 6*cm );
0946     
0947     logicVirtualMag= new G4LogicalVolume( virtualMag,
0948                                          internalChamberMaterial,
0949                                          "LVirtualMag",
0950                                          0,0,0);
0951     physicVirtualMag = new G4PVPlacement(0,
0952                                          G4ThreeVector(virtualMagPosX, 0.*cm, 0*mm),
0953                                          "PVirtualMag",
0954                                          logicVirtualMag,
0955                                          physicTreatmentRoom,
0956                                          true, 0);
0957     
0958     
0959     logicVirtualMag -> SetVisAttributes(blue);
0960     
0961     //// BeveledCylinder ////
0962     
0963     G4RotationMatrix *Rot= new G4RotationMatrix;
0964     Rot->rotateX(14*deg);
0965     G4ThreeVector trans(0.,22.5*mm,-15*mm);
0966     Cylinder= new G4Tubs("cylinder",20*mm,22.5*mm,90*mm,0.,2*pi);
0967     Box= new G4Box("Box",22.5*mm,22.5*mm,90*mm);
0968     
0969     G4SubtractionSolid* BeveledCylinder=new G4SubtractionSolid("Cylinder-Box",
0970                                                                Cylinder,
0971                                                                Box,
0972                                                                Rot,
0973                                                                trans);
0974     
0975     logicBeveledCylinder= new G4LogicalVolume                   (BeveledCylinder,
0976                                                                  GuardRingMaterial,
0977                                                                  "LBeveledCylinder",
0978                                                                  0,0,0);
0979     
0980     physicBeveledCylinder =new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(GuardRingPosX,0,0)),
0981                                              "physicBeveledCylinder",
0982                                              logicBeveledCylinder,
0983                                              physicVirtualMag,
0984                                              true,0);
0985     
0986     logicBeveledCylinder->SetVisAttributes(green);
0987     
0988     
0989     ///// KaptonEntranceWindow /////
0990     
0991     KaptonEntranceWindow= new G4Tubs("KaptonEntranceWindow",
0992                                      0,
0993                                      OuterRadiusFC,
0994                                      KaptonEntranceWindowThickness/2,
0995                                      0*deg,360*deg);
0996     
0997     logicKaptonEntranceWindow=new G4LogicalVolume(          KaptonEntranceWindow,
0998                                                   // internalChamberMaterial, for track control
0999                                                   KaptonEntranceWindowMaterial,
1000                                                   "LKaptonEntranceWindow",
1001                                                   0,0,0);
1002     
1003     physicKaptonEntranceWindow=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(KaptonEntranceWindowPosX,0,0)),
1004                                                  "PhysicEntranceWindow",
1005                                                  logicKaptonEntranceWindow,
1006                                                  physicVirtualMag,true,0);
1007     logicKaptonEntranceWindow -> SetVisAttributes(gray);
1008     
1009     ////// MassRing /////
1010     
1011     MassRing=new G4Tubs ("MassRing",
1012                          InnerRadiusFC,
1013                          OuterRadiusFC,
1014                          MassRingThickness/2,
1015                          0*deg,360*deg);
1016     
1017     logicMassRing=new G4LogicalVolume(                      MassRing,
1018                                       MassRingMaterial,
1019                                       "logicMassRing",
1020                                       0,0,0);
1021     
1022     physicMassRing=new G4PVPlacement(                       G4Transform3D(rm,G4ThreeVector(MassRingPosX,0,0)),
1023                                      
1024                                      "PhysicMassRing",logicMassRing,
1025                                      
1026                                      physicVirtualMag,
1027                                      true,0);
1028     logicMassRing -> SetVisAttributes(green);
1029     
1030     
1031     
1032     
1033     ///// VirtualWindow /////
1034     
1035     
1036     VirtualWindow=new G4Tubs("VirtualWindow",
1037                              0,
1038                              OuterRadiusFC,
1039                              VirtualWindowThickness/2,
1040                              0*deg,360*deg);
1041     
1042     logicVirtualWindow=new G4LogicalVolume(                 VirtualWindow,
1043                                            internalChamberMaterial,
1044                                            "logicVirtualWindow",
1045                                            0,0,0);
1046     
1047     physicVirtualWindow=new G4PVPlacement(                   G4Transform3D(rm,G4ThreeVector(VirtualWindowPosX,0,0)),
1048                                           
1049                                           "PhysicVirtualWindow",
1050                                           logicVirtualWindow,
1051                                           physicVirtualMag,
1052                                           true,0);
1053     logicVirtualWindow->SetVisAttributes (G4VisAttributes::GetInvisible());
1054     
1055     ///// GuardRing /////
1056     
1057     GuardRing=new G4Tubs ("GuardRing",
1058                           InnerRadiusFC,
1059                           OuterRadiusFC,
1060                           GuardRingThickness/2,
1061                           0*deg,360*deg);
1062     
1063     logicGuardRing=new G4LogicalVolume(                      GuardRing,
1064                                        GuardRingMaterial,
1065                                        "logicGuardRing",
1066                                        0,0,0);
1067     
1068     physicGuardRing=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(GuardRingPosX,0,0)),
1069                                       
1070                                       "PhysicGuardRing", logicGuardRing,
1071                                       
1072                                       physicVirtualMag,
1073                                       true,0);
1074     logicGuardRing -> SetVisAttributes(red);
1075     
1076     
1077     /////VirtualMiddle /////
1078     
1079     
1080     VirtualMiddle=new G4Tubs ("VirtualMiddle",
1081                               0,
1082                               OuterRadiusFC,
1083                               VirtualMiddleThickness/2,
1084                               0*deg,360*deg);
1085     
1086     logicVirtualMiddle=new G4LogicalVolume(                      VirtualMiddle,
1087                                            internalChamberMaterial,
1088                                            "logicVirtualMiddle",
1089                                            0,0,0);
1090     
1091     physicVirtualMiddle=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(VirtualMiddlePosX,0,0)),
1092                                           
1093                                           "PhysicVirtualMiddle", logicVirtualMiddle,
1094                                           
1095                                           physicVirtualMag,
1096                                           true,0);
1097     
1098     logicVirtualMiddle->SetVisAttributes (G4VisAttributes::GetInvisible());
1099     
1100     ///// FaradayCupBottom /////
1101     
1102     FaradayCupBottom=new G4Tubs ("FaradayCupBottom",
1103                                  InnerRadiusFC,
1104                                  OuterRadiusFC,
1105                                  FaradayCupBottomThickness/2,
1106                                  0*deg,360*deg);
1107     
1108     logicFaradayCupBottom=new G4LogicalVolume(             FaradayCupBottom,
1109                                               FaradayCupBottomMaterial,
1110                                               "logicFaradayCupBottom",
1111                                               0,0,0);
1112     
1113     physicFaradayCupBottom=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(FaradayCupBottomPosX,0,0)),
1114                                              "PhysicFaradayCupBottom",logicFaradayCupBottom,
1115                                              physicVirtualMag,
1116                                              true,0);
1117     logicFaradayCupBottom -> SetVisAttributes(yellow);
1118     
1119     
1120     ///// Virtual Bottom //////
1121     
1122     VirtualBottom=new G4Tubs ("VirtualBottom",
1123                               0,
1124                               OuterRadiusFC,
1125                               VirtualBottomThickness/2,
1126                               0*deg,360*deg);
1127     
1128     logicVirtualBottom=new G4LogicalVolume(                  VirtualBottom,
1129                                            internalChamberMaterial,
1130                                            "logicVirtualBottom",
1131                                            0,0,0);
1132     
1133     physicVirtualBottom=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(VirtualBottomPosX,0,0)),
1134                                           "PhysicVirtualBottom",
1135                                           logicVirtualBottom,
1136                                           physicVirtualMag,
1137                                           true,0);
1138     
1139     logicVirtualBottom->SetVisAttributes (G4VisAttributes::GetInvisible());
1140     
1141     ///// Cup /////
1142     
1143     Cup=new G4Tubs ("Cup",
1144                     0,
1145                     OuterRadiusFC,
1146                     CupThickness/2,
1147                     0*deg,360*deg);
1148     
1149     logicCup=new G4LogicalVolume(                           Cup,
1150                                  CupMaterial,
1151                                  "logicCup",
1152                                  0,0,0);
1153     
1154     physicCup=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(CupPosX,0,0)),
1155                                 "PhysicCup", logicCup,
1156                                 
1157                                 physicVirtualMag,
1158                                 true,0);
1159     
1160     logicCup -> SetVisAttributes(darkGreen);
1161     
1162     
1163     ///// Virtual OverBottom /////
1164     
1165     VirtualOverBottom=new G4Tubs ("VirtualOverBottom",
1166                                   0,
1167                                   OuterRadiusFC,
1168                                   VirtualOverBottomThickness/2,
1169                                   0*deg,360*deg);
1170     
1171     logicVirtualOverBottom=new G4LogicalVolume(             VirtualOverBottom,
1172                                                internalChamberMaterial,
1173                                                "logicVirtualOverBottom",
1174                                                0,0,0);
1175     
1176     physicVirtualOverBottom=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(VirtualOverBottomPosX,0,0)),
1177                                               "PhysicVirtualOverBottom",logicVirtualOverBottom,
1178                                               
1179                                               physicVirtualMag,
1180                                               true,0);
1181     logicVirtualOverBottom->SetVisAttributes (G4VisAttributes::GetInvisible());
1182     
1183     
1184     ///// Virtual Lateral /////
1185     
1186     
1187     VirtualLateral=new G4Tubs ("VirtualLateral",
1188                                OuterRadiusFC,
1189                                OuterRadiusFC+1*um,// the VirtualLateralThickness is 1*um
1190                                VirtualLateralLength/2,
1191                                0*deg,360*deg);
1192     
1193     logicVirtualLateral=new G4LogicalVolume(                  VirtualLateral,
1194                                             internalChamberMaterial,
1195                                             "logicVirtualLateral",
1196                                             0,0,0);
1197     
1198     physicVirtualLateral=new G4PVPlacement(G4Transform3D(rm,G4ThreeVector(VirtualLateralPosX,0,0)),
1199                                            "VirtualLateral",logicVirtualLateral,
1200                                            
1201                                            physicVirtualMag,
1202                                            true,0);
1203     
1204     
1205     
1206     logicVirtualLateral->SetVisAttributes (G4VisAttributes::GetInvisible());
1207 }
1208 
1209 /////////////////////////////////////////////////////////////////////////////
1210 void LaserDrivenBeamLine::Quadrupole()
1211 {
1212     // To rotate the quadrupoles putting their axis (along X direction) parallel to the beam axis
1213     G4double phi = 90. *deg;
1214     G4RotationMatrix rm;
1215     rm.rotateY(phi);
1216     
1217     SQuadChamberWall = new G4Box("solidQuadChamberWall",externalChamberXSize/2., externalChamberYSize/2.,externalChamberZSize/2.);
1218     
1219     LQuadChamberWall = new G4LogicalVolume(SQuadChamberWall, externalChamberMaterial,"logicQuadChamberWall");
1220     
1221     PQuadChamberWall = new G4PVPlacement(0, G4ThreeVector(QuadChamberWallPosX, QuadChamberWallPosY, QuadChamberWallPosZ),
1222                                          "physQuadChamberWall", LQuadChamberWall,physicTreatmentRoom, false, 0);
1223     
1224     
1225     SQuadChamber = new G4Box("solidQuadChamber", internalChamberXSize/2., internalChamberYSize/2.,internalChamberZSize/2.);
1226     
1227     LQuadChamber = new G4LogicalVolume(SQuadChamber, internalChamberMaterial,"logicQuadChamber");
1228     
1229     PQuadChamber = new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, 0.0),
1230                                      "physQuadChamber", LQuadChamber,PQuadChamberWall, false, 0);
1231     
1232     LQuadChamberWall -> SetVisAttributes(red);
1233     LQuadChamber -> SetVisAttributes(white);
1234     ///////////----------------------------Fourth Quadrupole----------------------------/////////
1235     SFourthTriplet = new G4Tubs("SolidTQuad", InnerRadiusTriplet, ExternalRadiusQuad,((FourthQuadThickness/2.)+1*mm),
1236                                 startAngleQuad, spanningAngleQuad);
1237     
1238     LFourthTriplet = new G4LogicalVolume(SFourthTriplet, internalChamberMaterial,"LogicTQuad", 0, 0, 0);
1239     
1240     PFourthTriplet = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(FourthQuadXPosition, QuadYPosition, QuadZPosition)),
1241                                        "PhysFourthTQuad", LFourthTriplet, PQuadChamber, false, 0);
1242     
1243     solidFourthQuad = new G4Tubs("SolidQuad", InnerRadiusQuad, ExternalRadiusQuad, FourthQuadThickness/2.,
1244                                  startAngleQuad, spanningAngleQuad);
1245     
1246     logicFourthQuad = new G4LogicalVolume(solidFourthQuad, QuadMaterial, "LogicQuad", 0, 0, 0);
1247     
1248     physicFourthQuad = new G4PVPlacement(0, G4ThreeVector(FourthQXPosition, QYPosition, QZPosition),
1249                                          "PhysFourthQuad",logicFourthQuad, PFourthTriplet, false, 0);
1250     
1251     LFourthTriplet -> SetVisAttributes(yellow);
1252     logicFourthQuad -> SetVisAttributes(green);
1253     ///////////----------------------------Third Quadrupole----------------------------/////////
1254     SThirdTriplet = new G4Tubs("SolidTTQuad", InnerRadiusTriplet, ExternalRadiusQuad, (ThirdQuadThickness/2.)+1*mm,
1255                                startAngleQuad, spanningAngleQuad);
1256     
1257     LThirdTriplet = new G4LogicalVolume(SThirdTriplet, internalChamberMaterial,"LogicTTQuad", 0, 0, 0);
1258     
1259     PThirdTriplet = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(ThirdQuadXPosition, QuadYPosition, QuadZPosition)),
1260                                       "PhysThirdTQuad",LThirdTriplet,PQuadChamber, false, 0);
1261     
1262     solidThirdQuad = new G4Tubs("SolidTQuad", InnerRadiusQuad, ExternalRadiusQuad, ThirdQuadThickness/2.,
1263                                 startAngleQuad, spanningAngleQuad);
1264     
1265     logicThirdQuad = new G4LogicalVolume(solidThirdQuad, QuadMaterial, "LogicTQuad", 0, 0, 0);
1266     
1267     physicThirdQuad = new G4PVPlacement(0, G4ThreeVector(ThirdQXPosition, QYPosition, QZPosition),
1268                                         "PhysThirdQuad",logicThirdQuad, PThirdTriplet, false, 0);
1269     
1270     LThirdTriplet -> SetVisAttributes(yellow);
1271     logicThirdQuad -> SetVisAttributes(green);
1272     ///////////----------------------------Second Quadrupole----------------------------/////////
1273     SSecondTriplet = new G4Tubs("SolidTSQuad", InnerRadiusTriplet, ExternalRadiusQuad, (SecondQuadThickness/2.)+1*mm,
1274                                 startAngleQuad, spanningAngleQuad);
1275     
1276     LSecondTriplet = new G4LogicalVolume(SSecondTriplet, internalChamberMaterial,"LogicTSQuad", 0, 0, 0);
1277     
1278     PSecondTriplet = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(SecondQuadXPosition, QuadYPosition, QuadZPosition)),
1279                                        "PhysSecondTQuad", LSecondTriplet, PQuadChamber, false, 0);
1280     
1281     solidSecondQuad = new G4Tubs("SolidSQuad", InnerRadiusQuad, ExternalRadiusQuad, SecondQuadThickness/2.,
1282                                  startAngleQuad, spanningAngleQuad);
1283     
1284     logicSecondQuad = new G4LogicalVolume(solidSecondQuad, QuadMaterial, "LogicSQuad", 0, 0, 0);
1285     
1286     physicSecondQuad = new G4PVPlacement(0, G4ThreeVector(SecondQXPosition, QYPosition, QZPosition),
1287                                          "PhysSecondQuad", logicSecondQuad, PSecondTriplet, false, 0);
1288     
1289     LSecondTriplet -> SetVisAttributes(yellow);
1290     logicSecondQuad -> SetVisAttributes(green);
1291     ///////////----------------------------First Quadrupole----------------------------/////////
1292     SFirstTriplet = new G4Tubs("SolidTQuad", InnerRadiusTriplet, ExternalRadiusQuad, (FirstQuadThickness/2.)+1*mm,
1293                                startAngleQuad, spanningAngleQuad);
1294     
1295     LFirstTriplet = new G4LogicalVolume(SFirstTriplet, internalChamberMaterial,"LogicTQuad", 0, 0, 0);
1296     
1297     PFirstTriplet = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(FirstQuadXPosition, QuadYPosition, QuadZPosition)),
1298                                       "PhysFirstTQuad", LFirstTriplet, PQuadChamber, false, 0);
1299     
1300     solidFirstQuad = new G4Tubs("SolidQuad", InnerRadiusQuad, ExternalRadiusQuad, FirstQuadThickness/2.,
1301                                 startAngleQuad, spanningAngleQuad);
1302     
1303     logicFirstQuad = new G4LogicalVolume(solidFirstQuad, QuadMaterial, "LogicQuad", 0, 0, 0);
1304     
1305     physicFirstQuad = new G4PVPlacement(0, G4ThreeVector(FirstQXPosition, QYPosition, QZPosition),
1306                                         "PhysFirstQuad",logicFirstQuad, PFirstTriplet, false, 0);
1307     
1308     LFirstTriplet -> SetVisAttributes(yellow);
1309     logicFirstQuad -> SetVisAttributes(green);
1310 }
1311 
1312 /////////////////////////////////////////////////////////////////////////////
1313 void LaserDrivenBeamLine::EnergySelectorChamber()
1314 {
1315     // The whole energyselector is mounted inside a
1316     // a vacuum chamber  (called 'ExternalChamber')
1317     // inside which a vacuum box is inserted.
1318     
1319     solidExternalChamber = new G4Box("ExternalChamber",
1320                                      externalChamberXSize/2.0,
1321                                      externalChamberYSize/2.0,
1322                                      externalChamberZSize/2.0);
1323     
1324     logicExternalChamber = new G4LogicalVolume(solidExternalChamber,
1325                                                externalChamberMaterial,
1326                                                "ExternalChamber");
1327     
1328     physicExternalChamber = new G4PVPlacement(0,
1329                                               G4ThreeVector(externalChamberXPosition,
1330                                                             externalChamberYPosition,
1331                                                             externalChamberZPosition),
1332                                               "ExternalChamber",
1333                                               logicExternalChamber,
1334                                               physicTreatmentRoom,
1335                                               false,
1336                                               0);
1337     
1338     // Visualisation of the External part
1339     logicExternalChamber -> SetVisAttributes(red);
1340     
1341     // This is a vacuum box inside the steel box
1342     solidInternalChamber = new G4Box("SInternalChamber",
1343                                      internalChamberXSize/2.0,
1344                                      internalChamberYSize/2.0,
1345                                      internalChamberZSize/2.0);
1346     
1347     logicInternalChamber = new G4LogicalVolume(solidInternalChamber,
1348                                                internalChamberMaterial,
1349                                                "LInternalChamber");
1350     
1351     physicInternalChamber = new G4PVPlacement(0,
1352                                               G4ThreeVector(0,0,0),
1353                                               "InternalChamber",
1354                                               logicInternalChamber,
1355                                               physicExternalChamber,
1356                                               false,
1357                                               0);
1358     logicInternalChamber -> SetVisAttributes(white);
1359 }
1360 
1361 //////////////////////////////////////////////////// Entrance pipe ///////////////////////
1362 void LaserDrivenBeamLine::EntrancePipe()
1363 {
1364     // To rotate the EntrancePipe putting its axis (along X direction) parallel to the beam axis
1365     G4double phi = 90. *deg;
1366     G4RotationMatrix rm;
1367     rm.rotateY(phi);
1368     
1369     solidEntrancePipe = new G4Tubs("EntrancePipe",
1370                                    InnerRadiusEntrancePipe,
1371                                    ExternalRadiusEntrancePipe,
1372                                    EntrancePipeheight/2.,
1373                                    startAngleEntrancePipe,
1374                                    spanningAngleEntrancePipe);
1375     
1376     logicEntrancePipe = new G4LogicalVolume(solidEntrancePipe,
1377                                             PipeMaterial,
1378                                             "EntrancePipe",
1379                                             0,
1380                                             0,
1381                                             0);
1382     
1383     physicEntrancePipe = new G4PVPlacement(G4Transform3D(rm,
1384                                                          G4ThreeVector(EntrancePipeXPosition,
1385                                                                        EntrancePipeYPosition,
1386                                                                        EntrancePipeZPosition)),
1387                                            "EntrancePipe",
1388                                            logicEntrancePipe,
1389                                            physicTreatmentRoom,
1390                                            false,
1391                                            0);
1392     
1393     logicEntrancePipe -> SetVisAttributes(red);
1394     
1395 }
1396 
1397 //////////////////////////////////////////////////// Entrance hole ///////////////////////
1398 void LaserDrivenBeamLine::Entrancehole()
1399 {
1400     // To rotate the ExitPipe putting its axis (along X direction) parallel to the beam axis
1401     G4double phi = 90. *deg;
1402     G4RotationMatrix rm;
1403     rm.rotateY(phi);
1404     
1405     solidEntrancehole = new G4Tubs("Entrancehole",
1406                                    InnerRadiusEntrancehole,
1407                                    ExternalRadiusEntrancehole,
1408                                    EntranceholeThickness/2.,
1409                                    startAngleEntrancehole,
1410                                    spanningAngleEntrancehole);
1411     
1412     logicEntrancehole = new G4LogicalVolume(solidEntrancehole,
1413                                             internalChamberMaterial,
1414                                             "Entrancehole",
1415                                             0,
1416                                             0,
1417                                             0);
1418     //the hole in the energy selector chamber
1419     physicEntranceholeESSChamber = new G4PVPlacement(G4Transform3D(rm,
1420                                                                    G4ThreeVector(EntranceholeXPosition,
1421                                                                                  EntranceholeYPosition,
1422                                                                                  EntranceholeZPosition)),
1423                                                      "Entrancehole",
1424                                                      logicEntrancehole,
1425                                                      physicExternalChamber,
1426                                                      false,
1427                                                      0);
1428     //the hole in the quadrupoles chamber
1429     physicEntrancehole = new G4PVPlacement(G4Transform3D(rm,
1430                                                          G4ThreeVector(EntranceholeQuadXPosition,
1431                                                                        EntranceholeYPosition,
1432                                                                        EntranceholeZPosition)),
1433                                            "EntranceholeQuad",
1434                                            logicEntrancehole,
1435                                            PQuadChamberWall,
1436                                            false,
1437                                            0);
1438     
1439     logicEntrancehole -> SetVisAttributes(skyBlue);
1440     
1441     
1442 }
1443 /////////////////////////////////////////////////////////////////////////////
1444 void LaserDrivenBeamLine::Collimator()
1445 {
1446     // To rotate the collimator putting its axis (along X direction) parallel to the beam axis
1447     G4double phi = 90. *deg;
1448     G4RotationMatrix rm;
1449     rm.rotateY(phi);
1450     //8x82x210 mm are the collimator default dimensions
1451     solidCollimator = new G4Box("collimator",
1452                                 thicknessCollimator/2.0,
1453                                 collimatorBoxYSize/2.0,
1454                                 collimatorBoxZSize/2.0);
1455     
1456     logicCollimator = new G4LogicalVolume(solidCollimator,
1457                                           collimatorMaterial,
1458                                           "collimator");
1459     
1460     physicCollimator = new G4PVPlacement(0,
1461                                          G4ThreeVector(collimatorBox_XPosition,
1462                                                        collimatorBox_YPosition,
1463                                                        collimatorBox_ZPosition),
1464                                          "collimator",
1465                                          logicCollimator,
1466                                          physicInternalChamber,
1467                                          false,
1468                                          0);
1469     
1470     logicCollimator -> SetVisAttributes(darkOrange3);
1471     
1472     solidCollimatorHole = new G4Tubs("CollimatorHole",
1473                                      innerRadiusCollimator,
1474                                      outerRadiusCollimator,
1475                                      thicknessCollimator/2.,
1476                                      startAngleCollimator,
1477                                      spanningAngleCollimator);
1478     
1479     logicCollimatorHole = new G4LogicalVolume(solidCollimatorHole,
1480                                               collimatorHoleMaterial,
1481                                               "CollimatorHole",
1482                                               0,
1483                                               0,
1484                                               0);
1485     
1486     physicCollimatorHole = new G4PVPlacement(G4Transform3D(rm,
1487                                                            G4ThreeVector(collimatorXPosition,
1488                                                                          collimatorYPosition,
1489                                                                          collimatorZPosition)),
1490                                              "CollimatorHole",
1491                                              logicCollimatorHole,
1492                                              physicCollimator,
1493                                              false,
1494                                              0);
1495     
1496     logicCollimatorHole -> SetVisAttributes(skyBlue);
1497 }
1498 
1499 /////////////////////////////////////////////////////////////////////////////
1500 // Magnet number 1
1501 void LaserDrivenBeamLine::Magnet_1()
1502 {  // The positions of the external and internal partes are given as respect the external chamber.
1503     solidExternalMagnet_1 = new G4Box("SolidExternalMagnet_1",
1504                                       externalMagnet_1XSize/2.0,
1505                                       externalMagnet_1YSize/2.0,
1506                                       externalMagnet_1ZSize/2.0);
1507     
1508     logicExternalMagnet_1 = new G4LogicalVolume(solidExternalMagnet_1,
1509                                                 externalMagnet_1Material,
1510                                                 "LogicExternalMagnet_1");
1511     
1512     physicExternalMagnet_1 = new G4PVPlacement(0,
1513                                                G4ThreeVector(externalMagnet_1XPosition,
1514                                                              externalMagnet_2YPosition,
1515                                                              externalMagnet_2ZPosition),
1516                                                "PhysicExternalMagnet_1",
1517                                                logicExternalMagnet_1,
1518                                                physicInternalChamber,
1519                                                false,
1520                                                0);
1521     physicExternalMagnet_1Down = new G4PVPlacement(0,
1522                                                    G4ThreeVector(externalMagnet_1XPosition,
1523                                                                  -externalMagnet_2YPosition,
1524                                                                  externalMagnet_2ZPosition),
1525                                                    "PhysicExternalMagnet_1Down",
1526                                                    logicExternalMagnet_1,
1527                                                    physicInternalChamber,
1528                                                    false,
1529                                                    0);
1530     
1531     
1532     logicExternalMagnet_1 -> SetVisAttributes(gray);
1533     
1534     // The right and left part of the magnet
1535     solidMagnet_1 = new G4Box("SolidMagnet_1",
1536                               Magnet_1XSize/2.0,
1537                               Magnet_1YSize/2.0,
1538                               Magnet_1ZSize/2.0);
1539     
1540     logicMagnet_1 = new G4LogicalVolume(solidMagnet_1,
1541                                         externalMagnet_1Material,
1542                                         "LogicMagnet_1");
1543     
1544     physicMagnet_1Right = new G4PVPlacement(0,
1545                                             G4ThreeVector(Magnet_1XPosition,Magnet_1YPosition,
1546                                                           Magnet_1ZPosition),
1547                                             "PhysicMagnet_1Right",
1548                                             logicMagnet_1,
1549                                             physicInternalChamber,
1550                                             false,
1551                                             0);
1552     physicMagnet_1Left = new G4PVPlacement(0,
1553                                            G4ThreeVector(Magnet_1XPosition,Magnet_1YPosition,
1554                                                          -Magnet_1ZPosition),
1555                                            "PhysicMagnet_1Left",
1556                                            logicMagnet_1,
1557                                            physicInternalChamber,
1558                                            false,
1559                                            0);
1560     
1561     logicMagnet_1 -> SetVisAttributes(gray);
1562 }
1563 
1564 /////////////////////////////////////////////////////////////////////////////
1565 // Magnet number 2
1566 void LaserDrivenBeamLine::Magnet_2()
1567 { // The position of the external part are given as respect the external chamber.
1568     
1569     solidExternalMagnet_2 = new G4Box("SolidExternalMagnet_2",
1570                                       externalMagnet_2XSize/2.0,
1571                                       externalMagnet_2YSize/2.0,
1572                                       externalMagnet_2ZSize/2.0);
1573     
1574     logicExternalMagnet_2 = new G4LogicalVolume(solidExternalMagnet_2,
1575                                                 externalMagnet_2Material,
1576                                                 "LogicExternalMagnet_2");
1577     
1578     physicExternalMagnet_2 = new G4PVPlacement(0,
1579                                                G4ThreeVector(externalMagnet_2XPosition,
1580                                                              externalMagnet_2YPosition,
1581                                                              (externalMagnet_2ZPosition+32*mm)),
1582                                                "PhysicExternalMagnet_2",
1583                                                logicExternalMagnet_2,
1584                                                physicInternalChamber,
1585                                                false,
1586                                                0);
1587     
1588     physicExternalMagnet_2Down = new G4PVPlacement(0,
1589                                                    G4ThreeVector(externalMagnet_2XPosition,
1590                                                                  -externalMagnet_2YPosition,
1591                                                                  (externalMagnet_2ZPosition+32*mm)),
1592                                                    "PhysicExternalMagnet_2Down",
1593                                                    logicExternalMagnet_2,
1594                                                    physicInternalChamber,
1595                                                    false,
1596                                                    0);
1597     
1598     
1599     logicExternalMagnet_2 -> SetVisAttributes(gray);
1600     
1601     // The right and left part of the magnet
1602     solidMagnet_2 = new G4Box("SolidMagnet_2",
1603                               Magnet_2XSize/2.0,
1604                               Magnet_2YSize/2.0,
1605                               Magnet_2ZSize/2.0);
1606     
1607     logicMagnet_2 = new G4LogicalVolume(solidMagnet_2,
1608                                         externalMagnet_2Material,
1609                                         "LogicMagnet_2");
1610     
1611     physicMagnet_2Right = new G4PVPlacement(0,
1612                                             G4ThreeVector(Magnet_2XPosition,Magnet_2YPosition,
1613                                                           (Magnet_2ZPosition)+32*mm),
1614                                             "PhysicMagnet_2Right",
1615                                             logicMagnet_2,
1616                                             physicInternalChamber,
1617                                             false,
1618                                             0);
1619     physicMagnet_2Left = new G4PVPlacement(0,
1620                                            G4ThreeVector(Magnet_2XPosition,Magnet_2YPosition,
1621                                                          (-(Magnet_2ZPosition)+32*mm)),
1622                                            "PhysicMagnet_2Left",
1623                                            logicMagnet_2,
1624                                            physicInternalChamber,
1625                                            false,
1626                                            0);
1627     logicMagnet_2 -> SetVisAttributes(gray);
1628 }
1629 
1630 /////////////////////////////////////////////////////////////////////////////
1631 // Magnet number 3
1632 void LaserDrivenBeamLine::Magnet_3()
1633 { // The position of the external part are given as respect the external chamber.
1634     
1635     solidExternalMagnet_3 = new G4Box("SolidExternalMagnet_3",
1636                                       externalMagnet_3XSize/2.0,
1637                                       externalMagnet_3YSize/2.0,
1638                                       externalMagnet_3ZSize/2.0);
1639     
1640     logicExternalMagnet_3 = new G4LogicalVolume(solidExternalMagnet_3,
1641                                                 externalMagnet_3Material,
1642                                                 "LogicExternalMagnet_3");
1643     
1644     physicExternalMagnet_3 = new G4PVPlacement(0,
1645                                                G4ThreeVector((externalMagnet_3XPosition),
1646                                                              externalMagnet_3YPosition,
1647                                                              (externalMagnet_3ZPosition+32*mm)),
1648                                                "PhysicExternalMagnet_3",
1649                                                logicExternalMagnet_3,
1650                                                physicInternalChamber,
1651                                                false,
1652                                                0);
1653     
1654     physicExternalMagnet_3Down = new G4PVPlacement(0,
1655                                                    G4ThreeVector((externalMagnet_3XPosition),
1656                                                                  -externalMagnet_3YPosition,
1657                                                                  (externalMagnet_3ZPosition+32*mm)),
1658                                                    "PhysicExternalMagnet_3Down",
1659                                                    logicExternalMagnet_3,
1660                                                    physicInternalChamber,
1661                                                    false,
1662                                                    0);
1663     
1664     logicExternalMagnet_3 -> SetVisAttributes(gray);
1665     
1666     // The right and left part of the magnet
1667     solidMagnet_3 = new G4Box("SolidMagnet_3",
1668                               Magnet_3XSize/2.0,
1669                               Magnet_3YSize/2.0,
1670                               Magnet_3ZSize/2.0);
1671     
1672     logicMagnet_3 = new G4LogicalVolume(solidMagnet_3,
1673                                         externalMagnet_3Material,
1674                                         "LogicMagnet_3");
1675     
1676     physicMagnet_3Right = new G4PVPlacement(0,
1677                                             G4ThreeVector(Magnet_3XPosition,Magnet_3YPosition,
1678                                                           (Magnet_3ZPosition+32*mm)),
1679                                             "PhysicMagnet_3Right",
1680                                             logicMagnet_3,
1681                                             physicInternalChamber,
1682                                             false,
1683                                             0);
1684     physicMagnet_3Left = new G4PVPlacement(0,
1685                                            G4ThreeVector(Magnet_3XPosition,Magnet_3YPosition,
1686                                                          (-(Magnet_3ZPosition)+32*mm)),
1687                                            "PhysicMagnet_3Left",
1688                                            logicMagnet_3,
1689                                            physicInternalChamber,
1690                                            false,
1691                                            0);
1692     logicMagnet_3 -> SetVisAttributes(gray);
1693     
1694 }
1695 
1696 /////////////////////////////////////////////////////////////////////////////
1697 // Magnet number 4
1698 void LaserDrivenBeamLine::Magnet_4()
1699 { // The position of the external part are given as respect the external chamber.
1700     
1701     solidExternalMagnet_4 = new G4Box("SolidExternalMagnet_4",
1702                                       externalMagnet_4XSize/2.0,
1703                                       externalMagnet_4YSize/2.0,
1704                                       externalMagnet_4ZSize/2.0);
1705     
1706     logicExternalMagnet_4 = new G4LogicalVolume(solidExternalMagnet_4,
1707                                                 externalMagnet_4Material,
1708                                                 "LogicExternalMagnet_4");
1709     
1710     physicExternalMagnet_4 = new G4PVPlacement(0,
1711                                                G4ThreeVector(externalMagnet_4XPosition,
1712                                                              externalMagnet_4YPosition,
1713                                                              externalMagnet_4ZPosition),
1714                                                "PhysicExternalMagnet_4",
1715                                                logicExternalMagnet_4,
1716                                                physicInternalChamber,
1717                                                false,
1718                                                0);
1719     
1720     physicExternalMagnet_4Down = new G4PVPlacement(0,
1721                                                    G4ThreeVector(externalMagnet_4XPosition,
1722                                                                  -externalMagnet_4YPosition,
1723                                                                  externalMagnet_4ZPosition),
1724                                                    "PhysicExternalMagnet_4Down",
1725                                                    logicExternalMagnet_4,
1726                                                    physicInternalChamber,
1727                                                    false,
1728                                                    0);
1729     
1730     logicExternalMagnet_4 -> SetVisAttributes(gray);
1731     
1732     // The right and left part of the magnet
1733     solidMagnet_4 = new G4Box("SolidMagnet_4",
1734                               Magnet_4XSize/2.0,
1735                               Magnet_4YSize/2.0,
1736                               Magnet_4ZSize/2.0);
1737     
1738     logicMagnet_4 = new G4LogicalVolume(solidMagnet_4,
1739                                         externalMagnet_4Material,
1740                                         "LogicMagnet_4");
1741     
1742     physicMagnet_4Right = new G4PVPlacement(0,
1743                                             G4ThreeVector(Magnet_4XPosition,Magnet_4YPosition,
1744                                                           Magnet_4ZPosition),
1745                                             "PhysicMagnet_4Right",
1746                                             logicMagnet_4,
1747                                             physicInternalChamber,
1748                                             false,
1749                                             0);
1750     physicMagnet_4Left = new G4PVPlacement(0,
1751                                            G4ThreeVector(Magnet_4XPosition,Magnet_4YPosition,
1752                                                          -Magnet_4ZPosition),
1753                                            "PhysicMagnet_4Left",
1754                                            logicMagnet_4,
1755                                            physicInternalChamber,
1756                                            false,
1757                                            0);
1758     logicMagnet_4 -> SetVisAttributes(gray);
1759 }
1760 
1761 /////////////////////////////////////////////////////////////////////////////
1762 // Slit
1763 void LaserDrivenBeamLine::Slit()
1764 {
1765     solidExternalSlit = new G4Box("ExternalSlit",
1766                                   externalSlitXSize/2.0,
1767                                   externalSlitYSize/2.0,
1768                                   externalSlitZSize/2.0);
1769     
1770     logicExternalSlit = new G4LogicalVolume(solidExternalSlit,
1771                                             externalSlitMaterial,
1772                                             "ExternalSlit");
1773     
1774     physicExternalSlit = new G4PVPlacement(0,
1775                                            G4ThreeVector(externalSlitXPosition,
1776                                                          externalSlitYPosition,
1777                                                          externalSlitZPosition),
1778                                            "ExternalSlit",
1779                                            logicExternalSlit,
1780                                            physicInternalChamber,
1781                                            false,
1782                                            0);
1783     
1784     logicExternalSlit -> SetVisAttributes(green);
1785     // The hole
1786     solidInternalSlit = new G4Box("InternalSlit",
1787                                   internalSlitXSize/2.0,
1788                                   internalSlitYSize/2.0,
1789                                   internalSlitZSize/2.0);
1790     
1791     logicInternalSlit = new G4LogicalVolume(solidInternalSlit,
1792                                             internalSlitMaterial,
1793                                             "InternalSlit");
1794     
1795     physicInternalSlit = new G4PVPlacement(0,
1796                                            G4ThreeVector(internalSlitXPosition,
1797                                                          internalSlitYPosition,
1798                                                          internalSlitZPosition),
1799                                            "InternalSlit",
1800                                            logicInternalSlit,
1801                                            physicExternalSlit,
1802                                            false,
1803                                            0);
1804     
1805     logicInternalSlit -> SetVisAttributes(skyBlue);
1806     
1807 }
1808 ////////////////////////////////////// Final collimator ////////////////////////////////////////////
1809 void LaserDrivenBeamLine::FinalCollimator()
1810 {
1811     // To rotate the collimator putting its axis (along X direction) parallel to the beam axis
1812     G4double phi = 90. *deg;
1813     G4RotationMatrix rm;
1814     rm.rotateY(phi);
1815     
1816     solidFinalCollimator = new G4Box("collimatorFinal",
1817                                      collimatorFinalBoxXSize/2.0,
1818                                      collimatorFinalBoxYSize/2.0,
1819                                      collimatorFinalBoxZSize/2.0);
1820     
1821     logicFinalCollimator = new G4LogicalVolume(solidFinalCollimator,
1822                                                FinalcollimatorMaterial,
1823                                                "collimatorFinal");
1824     
1825     physicFinalCollimator = new G4PVPlacement(0,
1826                                               G4ThreeVector(collimatorFinalBox_XPosition,
1827                                                             collimatorFinalBox_YPosition,
1828                                                             collimatorFinalBox_ZPosition),
1829                                               "collimatorFinal",
1830                                               logicFinalCollimator,
1831                                               physicInternalChamber,
1832                                               false,
1833                                               0);
1834     logicFinalCollimator -> SetVisAttributes(darkOrange3);
1835     
1836     solidFinalCollimatorHole= new G4Tubs("FinalCollimatorHole",
1837                                          innerRadiusFinalCollimator,
1838                                          outerRadiusFinalCollimator,
1839                                          FinalCollimatorThickness/2.,
1840                                          startAngleFinalCollimator,
1841                                          spanningAngleFinalCollimator);
1842     
1843     logicFinalCollimatorHole = new G4LogicalVolume(solidFinalCollimatorHole,
1844                                                    FinalcollimatorHoleMaterial,
1845                                                    "FinalCollimatorHole",
1846                                                    0,
1847                                                    0,
1848                                                    0);
1849     
1850     physicFinalCollimatorHole = new G4PVPlacement(G4Transform3D(rm,
1851                                                                 G4ThreeVector(FinalcollimatorXPosition,
1852                                                                               FinalcollimatorYPosition,
1853                                                                               FinalcollimatorZPosition)),
1854                                                   "FinalCollimatorHole",
1855                                                   logicFinalCollimatorHole,
1856                                                   physicFinalCollimator,
1857                                                   false,
1858                                                   0);
1859     logicFinalCollimatorHole -> SetVisAttributes(skyBlue);
1860 }
1861 //////////////////////////// Exit Window ////////////////////////////////////////////
1862 void LaserDrivenBeamLine::ExitWindow()
1863 {
1864     // To rotate the ExitWindow putting its axis (along X direction) parallel to the beam axis
1865     G4double phi = 90. *deg;
1866     G4RotationMatrix rm;
1867     rm.rotateY(phi);
1868     
1869     solidExitWindow = new G4Tubs("ExitWindow",
1870                                  InnerRadiusExitWindow,
1871                                  ExternalRadiusExitWindow,
1872                                  ExitWindowThickness/2.,
1873                                  startAngleExitWindow,
1874                                  spanningAngleExitWindow);
1875     
1876     logicExitWindow = new G4LogicalVolume(solidExitWindow,
1877                                           WindowMaterial,
1878                                           "ExitWindow",
1879                                           0,
1880                                           0,
1881                                           0);
1882     
1883     physicExitWindow = new G4PVPlacement(G4Transform3D(rm,
1884                                                        G4ThreeVector(ExitWindowXPosition,
1885                                                                      ExitWindowYPosition,
1886                                                                      ExitWindowZPosition)),
1887                                          "ExitWindow",
1888                                          logicExitWindow,
1889                                          physicTreatmentRoom,
1890                                          false,
1891                                          0);
1892     
1893     logicExitWindow -> SetVisAttributes(skyBlue);
1894     
1895 }
1896 
1897 //////////////////////////////////////////////////// Exit pipe ///////////////////////
1898 void LaserDrivenBeamLine::ExitPipe()
1899 {
1900     // To rotate the ExitPipe putting its axis (along X direction) parallel to the beam axis
1901     G4double phi = 90. *deg;
1902     G4RotationMatrix rm;
1903     rm.rotateY(phi);
1904     
1905     solidExitPipe = new G4Tubs("ExitPipe",
1906                                InnerRadiusExitPipe,
1907                                ExternalRadiusExitPipe,
1908                                ExitPipeheight/2.,
1909                                startAngleExitPipe,
1910                                spanningAngleExitPipe);
1911     
1912     logicExitPipe = new G4LogicalVolume(solidExitPipe,
1913                                         PipeMaterial,
1914                                         "ExitPipe",
1915                                         0,
1916                                         0,
1917                                         0);
1918     
1919     physicExitPipe = new G4PVPlacement(G4Transform3D(rm,
1920                                                      G4ThreeVector(ExitPipeXPosition,
1921                                                                    ExitPipeYPosition,
1922                                                                    ExitPipeZPosition)),
1923                                        "ExitPipe",
1924                                        logicExitPipe,
1925                                        physicTreatmentRoom,
1926                                        false,
1927                                        0);
1928     
1929     logicExitPipe -> SetVisAttributes(red);
1930     
1931 }
1932 
1933 ///////////////////////////////////// Exit hole ///////////////////////
1934 void LaserDrivenBeamLine::Exithole()
1935 {
1936     // To rotate the ExitPipe putting its axis (along X direction) parallel to the beam axis
1937     G4double phi = 90. *deg;
1938     G4RotationMatrix rm;
1939     rm.rotateY(phi);
1940     
1941     solidExithole = new G4Tubs("Exithole",
1942                                InnerRadiusExithole,
1943                                ExternalRadiusExithole,
1944                                ExitholeThickness/2.,
1945                                startAngleExithole,
1946                                spanningAngleExithole);
1947     
1948     logicExithole = new G4LogicalVolume(solidExithole,
1949                                         internalChamberMaterial,
1950                                         "Exithole",
1951                                         0,
1952                                         0,
1953                                         0);
1954     
1955     physicExithole = new G4PVPlacement(G4Transform3D(rm,
1956                                                      G4ThreeVector(ExitholeXPosition,
1957                                                                    ExitholeYPosition,
1958                                                                    ExitholeZPosition)),
1959                                        "Exithole",
1960                                        logicExithole,
1961                                        physicExternalChamber,
1962                                        false,
1963                                        0);
1964     
1965     logicExithole -> SetVisAttributes(skyBlue);
1966     
1967 }
1968 /////////////////////////// MESSENGER ///////////////////////////////////////
1969 /////////////////////////////////////////////////////////////////////////////
1970 // Disable via external macro command the Energy Selector System
1971 void LaserDrivenBeamLine::RemoveESS()
1972 {
1973     if(physicMagnet_1Left) {delete physicMagnet_1Left; delete physicMagnet_1Right; delete logicMagnet_1; delete solidMagnet_1;}
1974     if(physicExternalMagnet_1Down){delete physicExternalMagnet_1Down; delete physicExternalMagnet_1; delete logicExternalMagnet_1; delete solidExternalMagnet_1;}
1975     if(physicMagnet_2Left){delete physicMagnet_2Left; delete physicMagnet_2Right; delete logicMagnet_2; delete solidMagnet_2;}
1976     if(physicExternalMagnet_2Down){ delete physicExternalMagnet_2Down; delete physicExternalMagnet_2; delete logicExternalMagnet_2; delete solidExternalMagnet_2; }
1977     if(physicMagnet_3Left){delete physicMagnet_3Left; delete physicMagnet_3Right; delete logicMagnet_3; delete solidMagnet_3; }
1978     if(physicExternalMagnet_3Down){delete physicExternalMagnet_3Down; delete physicExternalMagnet_3; delete logicExternalMagnet_3; delete solidExternalMagnet_3; }
1979     if(physicMagnet_4Left) {delete physicMagnet_4Left; delete physicMagnet_4Right; delete logicMagnet_4; delete solidMagnet_4; }
1980     if(physicExternalMagnet_4Down){delete physicExternalMagnet_4Down; delete physicExternalMagnet_4; delete logicExternalMagnet_4; delete solidExternalMagnet_4; }
1981     if(physicCollimatorHole){delete physicCollimatorHole; delete logicCollimatorHole; delete solidCollimatorHole; }
1982     if(physicCollimator) {delete physicCollimator; delete logicCollimator; delete solidCollimator; }
1983     if(physicFinalCollimatorHole) {delete physicFinalCollimatorHole; delete logicFinalCollimatorHole; delete solidFinalCollimatorHole; }
1984     if(physicFinalCollimator){delete physicFinalCollimator; delete logicFinalCollimator; delete solidFinalCollimator; }
1985     if(physicInternalSlit){ delete physicInternalSlit; delete logicInternalSlit; delete solidInternalSlit; }
1986     if(physicExternalSlit){delete physicExternalSlit; delete logicExternalSlit; delete solidExternalSlit; }
1987     if(physicExithole){delete physicExithole; delete logicExithole; delete solidExithole;}
1988     if(physicExitWindow){delete physicExitWindow; delete logicExitWindow; delete solidExitWindow;}
1989     if(physicExitPipe){delete physicExitPipe; delete logicExitPipe; delete solidExitPipe;}
1990     if(physicEntranceholeESSChamber){delete physicEntranceholeESSChamber;}
1991     if(physicInternalChamber){delete physicInternalChamber; delete logicInternalChamber; delete solidInternalChamber;}
1992     if(physicExternalChamber) {delete physicExternalChamber; delete logicExternalChamber; delete solidExternalChamber;}
1993     if(pFieldMgr) {delete pFieldMgr;}
1994     
1995     
1996     
1997     G4cout << "****************************************************" << G4endl;
1998     G4cout << "************ The ESS has been disabled *************"  << G4endl;
1999     G4cout << "****************************************************" << G4endl;
2000     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2001 
2002     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2003 
2004 }
2005 // Change via external macro command the diameter of the first collimator
2006 void LaserDrivenBeamLine::SetFirstCollimatorRadius(G4double valueR)
2007 {
2008     G4double radius = valueR;
2009     solidCollimatorHole -> SetOuterRadius(radius);
2010     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2011 
2012     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2013 
2014     G4cout << "The first collimator aperture has been modified to "<< valueR/mm <<"mm in diameter" << G4endl;
2015 }
2016 /////////////////////////////////////////////////////////////////////////////
2017 // Change via external macro command the thickness of the first collimator
2018 void LaserDrivenBeamLine::SetFirstCollimatorThickness(G4double valueC)
2019 {
2020     G4double thickness = valueC/2;
2021     solidCollimator -> SetXHalfLength(thickness);
2022     solidCollimatorHole -> SetZHalfLength(thickness);
2023     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2024 
2025     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2026 
2027     G4cout << "The first collimator thickness has been modified to "<< valueC/mm <<" mm in thickness" << G4endl;
2028 }
2029 
2030 // Change via external macro command the Z position of the first collimator hole
2031 void LaserDrivenBeamLine::SetFirstCollimatorPositionZ(G4double valueQ)
2032 {
2033     physicCollimatorHole -> SetTranslation(G4ThreeVector(0., 0., valueQ));
2034     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2035 
2036     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2037 
2038     G4cout << "The first collimator has been translated to "<< valueQ/mm <<"mm (along the z axis)" << G4endl;
2039 }
2040 
2041 // Change via external macro command the diameter of the second collimator
2042 void LaserDrivenBeamLine::SetSecondCollimatorRadius(G4double value)
2043 {
2044     G4double radius = value;
2045     solidFinalCollimatorHole -> SetOuterRadius(radius);
2046     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2047 
2048     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2049 
2050     G4cout << "The second collimator aperture has been modified to "<< value/mm <<"mm in diameter" << G4endl;
2051 }
2052 
2053 /////////////////////////////////////////////////////////////////////////////
2054 // Change via external macro command the thickness of the second collimator
2055 void LaserDrivenBeamLine::SetSecondCollimatorThickness(G4double value)
2056 {
2057     G4double thickness = value/2;
2058     solidFinalCollimator -> SetXHalfLength(thickness);
2059     solidFinalCollimatorHole -> SetZHalfLength(thickness);
2060     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2061 
2062     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2063 
2064     G4cout << "The second collimator thickness has been modified to "<< value/mm <<" mm in thickness" << G4endl;
2065 }
2066 
2067 // Change via external macro command the Z position of the second collimator hole
2068 void LaserDrivenBeamLine::SetSecondCollimatorPositionZ(G4double value)
2069 {
2070     physicFinalCollimatorHole -> SetTranslation(G4ThreeVector(0., 0., value));
2071     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2072 
2073     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2074 
2075     G4cout << "The second collimator has been translated to "<< value/mm <<"mm (along the z axis)" << G4endl;
2076 }
2077 // THE SLIT MESSENGERS
2078 /////////////////////////////////////////////////////////////////////////////
2079 // Change the thickness of the Slit
2080 void LaserDrivenBeamLine::SetThicknessSlit(G4double value)
2081 {
2082     if (value >(10.0*mm)) {
2083         G4cout <<"***************************************"<< G4endl;
2084         G4cout <<"******This is a warning messenger******"<< G4endl;
2085         G4cout <<"***************************************"<< G4endl;
2086         G4cout <<"The maximum value of the thickness of the slit is 10 mm, your value is >10 mm." << G4endl;
2087         G4cout <<"The default thickness value is used, it is: " << ((solidExternalSlit -> GetXHalfLength())*2.)/mm
2088         << G4endl;
2089         G4cout <<"***************************************"<< G4endl;
2090         
2091     }
2092     else  {
2093         G4double dimension = value/2;
2094         solidExternalSlit -> SetXHalfLength(dimension);
2095         solidInternalSlit -> SetXHalfLength(dimension);
2096         G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2097 
2098         G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2099 
2100         G4cout <<"The thickness of the slit is:" << ((solidExternalSlit -> GetXHalfLength())*2.)/mm
2101         << G4endl;
2102     }
2103 }
2104 /////////////////////////////////////////////////////////////////////////////
2105 // Change the hole size (in Y direction) of the Slit
2106 void LaserDrivenBeamLine::SetSlitHoleDimensionY(G4double value)
2107 {
2108     G4double hole = value/2;
2109     solidInternalSlit -> SetYHalfLength(hole);
2110     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2111 
2112     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2113 
2114     G4cout << "The hole of the Slit has been changed in the Y direction to "<< value/mm <<" mm" <<G4endl;
2115 }
2116 
2117 /////////////////////////////////////////////////////////////////////////////
2118 // Change the hole size (in Z direction) of the Slit
2119 void LaserDrivenBeamLine::SetSlitHoleDimensionZ(G4double value)
2120 {
2121     G4double hole = value/2;
2122     solidInternalSlit -> SetZHalfLength(hole);
2123     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2124 
2125     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2126 
2127     G4cout << "The hole of the Slit has been changed in the Z direction to "<< value/mm <<" mm" <<G4endl;
2128 }
2129 /////////////////////////////////////////////////////////////////////////////
2130 // Change the Z position of the hole of the Slit
2131 void LaserDrivenBeamLine::SetSlitHolePositionZ(G4double value)
2132 {
2133     physicInternalSlit -> SetTranslation(G4ThreeVector(0., 0., value));
2134     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2135 
2136     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2137 
2138     G4cout << "The hole of the slit has been translated to "<< value/mm <<" mm (along the Z axis)" <<G4endl;
2139 }
2140 
2141 // QUADRUPOLES
2142 
2143 // Disable via external macro command all quadrupoles
2144 void LaserDrivenBeamLine::RemoveQuads()
2145 {
2146     if(physicFirstQuad)
2147     {delete solidFirstQuad; delete logicFirstQuad; delete physicFirstQuad;delete SFirstTriplet; delete LFirstTriplet; delete PFirstTriplet;}
2148     if(physicSecondQuad)
2149     {delete solidSecondQuad; delete logicSecondQuad; delete physicSecondQuad;delete SSecondTriplet; delete LSecondTriplet;  delete PSecondTriplet;}
2150     if(physicThirdQuad)
2151     {delete solidThirdQuad; delete logicThirdQuad; delete physicThirdQuad;delete SThirdTriplet; delete LThirdTriplet;   delete PThirdTriplet;}
2152     if(physicFourthQuad)
2153     {delete solidFourthQuad; delete logicFourthQuad; delete physicFourthQuad;delete SFourthTriplet; delete LFourthTriplet;  delete PFourthTriplet;}
2154     if(pFieldMgrQuadFourth) {delete pFieldMgrQuadFourth;}
2155     if(pFieldMgrQuadThird) {delete pFieldMgrQuadThird;}
2156     if(pFieldMgrQuadSecond) {delete pFieldMgrQuadSecond;}
2157     if(pFieldMgrQuadFirst) {delete pFieldMgrQuadFirst;}
2158     
2159     
2160     G4cout << "******************************************************************" << G4endl;
2161     G4cout << "************ The Quadrupoles system has been disabled *************"  << G4endl;
2162     G4cout << "******************************************************************" << G4endl;
2163     G4RunManager::GetRunManager() -> GeometryHasBeenModified();
2164 
2165     G4UImanager::GetUIpointer() -> ApplyCommand("/vis/viewer/flush");
2166 
2167 }
2168