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

#include <PassiveCarbonBeamLine.hh>

Inheritance diagram for PassiveCarbonBeamLine:
G4VUserDetectorConstruction

Public Member Functions

 PassiveCarbonBeamLine ()
 
 ~PassiveCarbonBeamLine ()
 
G4VPhysicalVolumeConstruct ()
 
void HadrontherapyBeamLineSupport ()
 
void ScatteringSystem ()
 
void VacuumToAirInterface ()
 
void HadrontherapyBeamMonitoring ()
 
void HadrontherapyBeamNozzle ()
 
void HadrontherapyBeamFinalCollimator ()
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void ConstructSDandField ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Data Fields

G4Materialkapton
 
G4VisAttributesredWire
 
G4VPhysicalVolumemother
 
G4double firstScatteringFoilXPosition
 
G4double firstScatteringFoilYPosition
 
G4double firstScatteringFoilZPosition
 

Additional Inherited Members

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

Detailed Description

Definition at line 44 of file PassiveCarbonBeamLine.hh.

Constructor & Destructor Documentation

PassiveCarbonBeamLine::PassiveCarbonBeamLine ( )

Definition at line 66 of file PassiveCarbonBeamLine.cc.

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

66  :
67 physicalTreatmentRoom(0), hadrontherapyDetectorConstruction(0),
68 physiBeamLineSupport(0), physiBeamLineCover(0), physiBeamLineCover2(0),
69 physiKaptonWindow(0),
70 physiFirstMonitorLayer1(0), physiFirstMonitorLayer2(0),
71 physiFirstMonitorLayer3(0), physiFirstMonitorLayer4(0),
72 physiNozzleSupport(0), physiHoleNozzleSupport(0),
73 physiSecondHoleNozzleSupport(0),
74 solidFinalCollimator(0),
75 physiFinalCollimator(0)
76 {
77  // Messenger to change parameters of the passiveProtonBeamLine geometry
78  //passiveMessenger = new PassiveProtonBeamLineMessenger(this);
79 
80 //***************************** PW ***************************************
81 
82  static G4String ROGeometryName = "DetectorROGeometry";
83  RO = new HadrontherapyDetectorROGeometry(ROGeometryName);
84 
85 
86 
87  G4cout << "Going to register Parallel world...";
89  G4cout << "... done" << G4endl;
90 //***************************** PW ***************************************
91 
92 
93 
94 }
void RegisterParallelWorld(G4VUserParallelWorld *)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
PassiveCarbonBeamLine::~PassiveCarbonBeamLine ( )

Definition at line 97 of file PassiveCarbonBeamLine.cc.

98 {
99  //delete passiveMessenger;
100  delete hadrontherapyDetectorConstruction;
101 }

Member Function Documentation

G4VPhysicalVolume * PassiveCarbonBeamLine::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 104 of file PassiveCarbonBeamLine.cc.

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

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

Definition at line 813 of file PassiveCarbonBeamLine.cc.

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

814 {
815  // -----------------------//
816  // FINAL COLLIMATOR //
817  //------------------------//
818  const G4double outerRadiusFinalCollimator = 21.5*mm;
819  const G4double hightFinalCollimator = 3.5*mm;
820  const G4double startAngleFinalCollimator = 0.*deg;
821  const G4double spanningAngleFinalCollimator = 360.*deg;
822  //XXX
823  const G4double finalCollimatorXPosition = -323.50 *mm;
824 
825  G4double phi = 90. *deg;
826 
827  // Matrix definition for a 90 deg rotation. Also used for other volumes
828  G4RotationMatrix rm;
829  rm.rotateY(phi);
830 
831  solidFinalCollimator = new G4Tubs("FinalCollimator",
832  innerRadiusFinalCollimator,
833  outerRadiusFinalCollimator,
834  hightFinalCollimator,
835  startAngleFinalCollimator,
836  spanningAngleFinalCollimator);
837 
838  G4LogicalVolume* logicFinalCollimator = new G4LogicalVolume(solidFinalCollimator,
839  finalCollimatorMaterial,
840  "FinalCollimator",
841  0,
842  0,
843  0);
844 
845  physiFinalCollimator = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(finalCollimatorXPosition,0.,0.)),
846  "FinalCollimator", logicFinalCollimator, physicalTreatmentRoom, false, 0);
847 
848  logicFinalCollimator -> SetVisAttributes(yellow);
849 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveCarbonBeamLine::HadrontherapyBeamLineSupport ( )

Definition at line 390 of file PassiveCarbonBeamLine.cc.

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

391 {
392  // ------------------//
393  // BEAM LINE SUPPORT //
394  //-------------------//
395  const G4double beamLineSupportXSize = 1.5*m;
396  const G4double beamLineSupportYSize = 20.*mm;
397  const G4double beamLineSupportZSize = 600.*mm;
398 
399  const G4double beamLineSupportXPosition = -1745.09 *mm;
400  const G4double beamLineSupportYPosition = -230. *mm;
401  const G4double beamLineSupportZPosition = 0.*mm;
402 
403  G4Box* beamLineSupport = new G4Box("BeamLineSupport",
404  beamLineSupportXSize,
405  beamLineSupportYSize,
406  beamLineSupportZSize);
407 
408  G4LogicalVolume* logicBeamLineSupport = new G4LogicalVolume(beamLineSupport,
409  beamLineSupportMaterial,
410  "BeamLineSupport");
411  physiBeamLineSupport = new G4PVPlacement(0, G4ThreeVector(beamLineSupportXPosition,
412  beamLineSupportYPosition,
413  beamLineSupportZPosition),
414  "BeamLineSupport",
415  logicBeamLineSupport,
416  physicalTreatmentRoom, false, 0);
417 
418  // Visualisation attributes of the beam line support
419 
420  logicBeamLineSupport -> SetVisAttributes(gray);
421 
422  //---------------------------------//
423  // Beam line cover 1 (left panel) //
424  //---------------------------------//
425  const G4double beamLineCoverXSize = 1.5*m;
426  const G4double beamLineCoverYSize = 750.*mm;
427  const G4double beamLineCoverZSize = 10.*mm;
428 
429  const G4double beamLineCoverXPosition = -1745.09 *mm;
430  const G4double beamLineCoverYPosition = -980.*mm;
431  const G4double beamLineCoverZPosition = 600.*mm;
432 
433  G4Box* beamLineCover = new G4Box("BeamLineCover",
434  beamLineCoverXSize,
435  beamLineCoverYSize,
436  beamLineCoverZSize);
437 
438  G4LogicalVolume* logicBeamLineCover = new G4LogicalVolume(beamLineCover,
439  beamLineSupportMaterial,
440  "BeamLineCover");
441 
442  physiBeamLineCover = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition,
443  beamLineCoverYPosition,
444  beamLineCoverZPosition),
445  "BeamLineCover",
446  logicBeamLineCover,
447  physicalTreatmentRoom,
448  false,
449  0);
450 
451  // ---------------------------------//
452  // Beam line cover 2 (rigth panel) //
453  // ---------------------------------//
454  // It has the same characteristic of beam line cover 1 but set in a different position
455  physiBeamLineCover2 = new G4PVPlacement(0, G4ThreeVector(beamLineCoverXPosition,
456  beamLineCoverYPosition,
457  - beamLineCoverZPosition),
458  "BeamLineCover2",
459  logicBeamLineCover,
460  physicalTreatmentRoom,
461  false,
462  0);
463 
464 
465  logicBeamLineCover -> SetVisAttributes(blue);
466 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
double G4double
Definition: G4Types.hh:76
void PassiveCarbonBeamLine::HadrontherapyBeamMonitoring ( )

Definition at line 613 of file PassiveCarbonBeamLine.cc.

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

614 {
615  // ----------------------------
616  // THE FIRST MONITOR CHAMBER
617  // ----------------------------
618  // A monitor chamber is a free-air ionisation chamber
619  // able to measure do carbon fluence during the treatment.
620  // Here its responce is not simulated in terms of produced
621  // charge but only the energy losses are taked into account.
622  // Each chamber consist of 9 mm of air in a box
623  // that has two layers one of kapton and one
624  // of copper
625  const G4double monitor1XSize = 4.525022*mm;
626  const G4double monitor2XSize = 0.000011*mm;
627  const G4double monitor3XSize = 4.5*mm;
628  const G4double monitorYSize = 10.*cm;
629  const G4double monitorZSize = 10.*cm;
630  // XXX (Camera monitor size = 9.050088 mm)
631  const G4double monitor1XPosition = -1450.474956 *mm;
632  const G4double monitor2XPosition = -4.500011*mm;
633  const G4double monitor4XPosition = 4.500011*mm;
634 
635  G4Box* solidFirstMonitorLayer1 = new G4Box("FirstMonitorLayer1",
636  monitor1XSize,
637  monitorYSize,
638  monitorZSize);
639 
640  G4LogicalVolume* logicFirstMonitorLayer1 = new G4LogicalVolume(solidFirstMonitorLayer1,
641  layer1MonitorChamberMaterial,
642  "FirstMonitorLayer1");
643 
644  physiFirstMonitorLayer1 = new G4PVPlacement(0,
645  G4ThreeVector(monitor1XPosition,0.*cm,0.*cm),
646  "FirstMonitorLayer1",
647  logicFirstMonitorLayer1,
648  physicalTreatmentRoom,
649  false,
650  0);
651 
652  G4Box* solidFirstMonitorLayer2 = new G4Box("FirstMonitorLayer2",
653  monitor2XSize,
654  monitorYSize,
655  monitorZSize);
656 
657  G4LogicalVolume* logicFirstMonitorLayer2 = new G4LogicalVolume(solidFirstMonitorLayer2,
658  layer2MonitorChamberMaterial,
659  "FirstMonitorLayer2");
660 
661  physiFirstMonitorLayer2 = new G4PVPlacement(0, G4ThreeVector(monitor2XPosition,0.*cm,0.*cm),
662  "FirstMonitorLayer2",
663  logicFirstMonitorLayer2,
664  physiFirstMonitorLayer1,
665  false,
666  0);
667 
668  G4Box* solidFirstMonitorLayer3 = new G4Box("FirstMonitorLayer3",
669  monitor3XSize,
670  monitorYSize,
671  monitorZSize);
672 
673  G4LogicalVolume* logicFirstMonitorLayer3 = new G4LogicalVolume(solidFirstMonitorLayer3,
674  layer3MonitorChamberMaterial,
675  "FirstMonitorLayer3");
676 
677  physiFirstMonitorLayer3 = new G4PVPlacement(0,
678  G4ThreeVector(0.*mm,0.*cm,0.*cm),
679  "MonitorLayer3",
680  logicFirstMonitorLayer3,
681  physiFirstMonitorLayer1,
682  false,
683  0);
684 
685  G4Box* solidFirstMonitorLayer4 = new G4Box("FirstMonitorLayer4",
686  monitor2XSize,
687  monitorYSize,
688  monitorZSize);
689 
690  G4LogicalVolume* logicFirstMonitorLayer4 = new G4LogicalVolume(solidFirstMonitorLayer4,
691  layer4MonitorChamberMaterial,
692  "FirstMonitorLayer4");
693 
694  physiFirstMonitorLayer4 = new G4PVPlacement(0, G4ThreeVector(monitor4XPosition,0.*cm,0.*cm),
695  "FirstMonitorLayer4",
696  logicFirstMonitorLayer4,
697  physiFirstMonitorLayer1, false, 0);
698 
699  logicFirstMonitorLayer3 -> SetVisAttributes(white);
700 
701 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
double G4double
Definition: G4Types.hh:76
void PassiveCarbonBeamLine::HadrontherapyBeamNozzle ( )

Definition at line 705 of file PassiveCarbonBeamLine.cc.

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

706 {
707  // ------------------------------//
708  // THE FINAL TUBE AND COLLIMATOR //
709  //-------------------------------//
710  // The last part of the transport beam line consists of
711  // a 59 mm thick PMMA slab (to stop all the diffused radiation), a 285 mm brass tube
712  // (to well collimate the carbon beam) and a final collimator with 25 mm diameter
713  // aperture (that provide the final trasversal shape of the beam)
714 
715  // -------------------//
716  // PMMA SUPPORT //
717  // -------------------//
718 
719  const G4double nozzleSupportXSize = 29.5 *mm;
720  const G4double nozzleSupportYSize = 180. *mm;
721  const G4double nozzleSupportZSize = 180. *mm;
722  //XXX Placed at
723  const G4double nozzleSupportXPosition = -558. *mm;
724 
725  G4double phi = 90. *deg;
726  // Matrix definition for a 90 deg rotation. Also used for other volumes
727  G4RotationMatrix rm;
728  rm.rotateY(phi);
729 
730  G4Box* solidNozzleSupport = new G4Box("NozzleSupport",
731  nozzleSupportXSize,
732  nozzleSupportYSize,
733  nozzleSupportZSize);
734 
735  G4LogicalVolume* logicNozzleSupport = new G4LogicalVolume(solidNozzleSupport,
736  nozzleSupportMaterial,
737  "NozzleSupport");
738 
739  physiNozzleSupport = new G4PVPlacement(0, G4ThreeVector(nozzleSupportXPosition,0., 0.),
740  "NozzleSupport",
741  logicNozzleSupport,
742  physicalTreatmentRoom,
743  false,
744  0);
745 
746  logicNozzleSupport -> SetVisAttributes(yellow);
747 
748  // -------------------//
749  // BRASS TUBE //
750  // -------------------//
751  const G4double innerRadiusHoleNozzleSupport = 18.*mm;
752  const G4double outerRadiusHoleNozzleSupport = 21.5 *mm;
753  //XXX h/2 = 142.5 mm
754  const G4double hightHoleNozzleSupport = 142.5*mm;
755  const G4double startAngleHoleNozzleSupport = 0.*deg;
756  const G4double spanningAngleHoleNozzleSupport = 360.*deg;
757  //XXX -(320+142.5)mm
758  const G4double holeNozzleSupportXPosition = -462.50 *mm;
759 
760  G4Tubs* solidHoleNozzleSupport = new G4Tubs("HoleNozzleSupport",
761  innerRadiusHoleNozzleSupport,
762  outerRadiusHoleNozzleSupport,
763  hightHoleNozzleSupport,
764  startAngleHoleNozzleSupport,
765  spanningAngleHoleNozzleSupport);
766 
767  G4LogicalVolume* logicHoleNozzleSupport = new G4LogicalVolume(solidHoleNozzleSupport,
768  holeNozzleSupportMaterial,
769  "HoleNozzleSupport",
770  0, 0, 0);
771 
772  physiHoleNozzleSupport = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(holeNozzleSupportXPosition, 0., 0.)),
773  "HoleNozzleSupport",
774  logicHoleNozzleSupport,
775  physicalTreatmentRoom, false, 0);
776 
777  logicHoleNozzleSupport -> SetVisAttributes(darkOrange3);
778 
779  //--------------------------------------------------------------//
780  // HOLE OF THE BRASS TUBE (otherwise we'll have PMMA) //
781  //--------------------------------------------------------------//
782  const G4double innerRadiusSecondHoleNozzleSupport = 0.*mm;
783  const G4double outerRadiusSecondHoleNozzleSupport = 18.*mm;
784  const G4double hightSecondHoleNozzleSupport = 29.5 *mm;
785  const G4double startAngleSecondHoleNozzleSupport = 0.*deg;
786  const G4double spanningAngleSecondHoleNozzleSupport = 360.*deg;
787 
788  G4Tubs* solidSecondHoleNozzleSupport = new G4Tubs("SecondHoleNozzleSupport",
789  innerRadiusSecondHoleNozzleSupport,
790  outerRadiusSecondHoleNozzleSupport,
791  hightSecondHoleNozzleSupport,
792  startAngleSecondHoleNozzleSupport,
793  spanningAngleSecondHoleNozzleSupport);
794 
795  G4LogicalVolume* logicSecondHoleNozzleSupport = new G4LogicalVolume(solidSecondHoleNozzleSupport,
796  seconHoleNozzleSupportMaterial,
797  "SecondHoleNozzleSupport",
798  0,
799  0,
800  0);
801 
802  physiSecondHoleNozzleSupport = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector()),
803  "SecondHoleNozzleSupport",
804  logicSecondHoleNozzleSupport,
805  physiNozzleSupport,
806  false, 0);
807 
808 
809  logicHoleNozzleSupport -> SetVisAttributes(darkOrange3);
810 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveCarbonBeamLine::ScatteringSystem ( )

Definition at line 546 of file PassiveCarbonBeamLine.cc.

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

547 {
548  // ------------//
549  // THE STOPPER //
550  //-------------//
551  // Is a small cylinder able to stop the central component
552  // of the beam (having a gaussian shape). It is connected to the SECON SCATTERING FOIL
553  // and represent the second element of the scattering system
554 
555  G4double phi = 90. *deg;
556  // Matrix definition for a 90 deg rotation with respect to Y axis
557  G4RotationMatrix rm;
558  rm.rotateY(phi);
559 
560  solidStopper = new G4Tubs("Stopper",
561  innerRadiusStopper,
562  outerRadiusStopper,
563  heightStopper/2,
564  startAngleStopper,
565  spanningAngleStopper);
566 
567  logicStopper = new G4LogicalVolume(solidStopper,
568  stopperMaterial,
569  "Stopper",
570  0, 0, 0);
571 
572  physiStopper = new G4PVPlacement(G4Transform3D(rm, G4ThreeVector(stopperXPosition,
573  stopperYPosition,
574  stopperZPosition)),
575  "Stopper",
576  logicStopper,
577  physicalTreatmentRoom,
578  false,
579  0);
580 
581  logicStopper -> SetVisAttributes(red);
582 
583  // ---------------------------//
584  // THE SECOND SCATTERING FOIL //
585  // ---------------------------//
586  // It is another thin foil and provides the
587  // final diffusion of the beam. It represents the third element of the scattering
588  // system;
589 
590  secondScatteringFoil = new G4Box("SecondScatteringFoil",
591  secondScatteringFoilXSize/2,
592  secondScatteringFoilYSize/2,
593  secondScatteringFoilZSize/2);
594 
595  G4LogicalVolume* logicSecondScatteringFoil = new G4LogicalVolume(secondScatteringFoil,
596  secondScatteringFoilMaterial,
597  "SecondScatteringFoil");
598 
599  physiSecondScatteringFoil = new G4PVPlacement(0, G4ThreeVector(secondScatteringFoilXPosition,
600  secondScatteringFoilYPosition,
601  secondScatteringFoilZPosition),
602  "SeconScatteringFoil",
603  logicSecondScatteringFoil,
604  physicalTreatmentRoom,
605  false,
606  0);
607 
608  logicSecondScatteringFoil -> SetVisAttributes(skyBlue);
609 
610 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:76
void PassiveCarbonBeamLine::VacuumToAirInterface ( )

Definition at line 469 of file PassiveCarbonBeamLine.cc.

References firstScatteringFoilXPosition.

470 {
471  // ------------//
472  // VACUUM PIPE //
473  //-------------//
474  //
475  // First track of the beam line is inside vacuum;
476  // The PIPE contains the FIRST SCATTERING FOIL and the KAPTON WINDOW
477  G4Box* vacuumZone = new G4Box("VacuumZone",
478  vacuumZoneXSize,
479  vacuumZoneYSize,
480  vacuumZoneZSize);
481 
482  G4LogicalVolume* logicVacuumZone = new G4LogicalVolume(vacuumZone,
483  vacuumZoneMaterial,
484  "VacuumZone");
485 
486  G4VPhysicalVolume* physiVacuumZone = new G4PVPlacement(0,
487  G4ThreeVector(vacuumZoneXPosition, 0., 0.),
488  "VacuumZone",
489  logicVacuumZone,
490  physicalTreatmentRoom,
491  false,
492  0);
493 
494 
495 
496 
497 
498  // --------------------------//
499  // THE FIRST SCATTERING FOIL //
500  // --------------------------//
501  // A thin foil performing a first scattering
502  // of the original beam
503 
504  firstScatteringFoil = new G4Box("FirstScatteringFoil",
505  firstScatteringFoilXSize/2,
506  firstScatteringFoilYSize/2,
507  firstScatteringFoilZSize/2);
508 
509  G4LogicalVolume* logicFirstScatteringFoil = new G4LogicalVolume(firstScatteringFoil,
510  firstScatteringFoilMaterial,
511  "FirstScatteringFoil");
512 
513  physiFirstScatteringFoil = new G4PVPlacement(0,
515  "FirstScatteringFoil",
516  logicFirstScatteringFoil,
517  physiVacuumZone,
518  false, 0);
519 
520  logicFirstScatteringFoil -> SetVisAttributes(skyBlue);
521 
522 
523 
524  // -------------------//
525  // THE KAPTON WINDOWS //
526  //--------------------//
527  //It permits the passage of the beam from vacuum to air
528 
529  G4Box* solidKaptonWindow = new G4Box("KaptonWindow",
530  kaptonWindowXSize,
531  kaptonWindowYSize,
532  kaptonWindowZSize);
533 
534  G4LogicalVolume* logicKaptonWindow = new G4LogicalVolume(solidKaptonWindow,
535  kaptonWindowMaterial,
536  "KaptonWindow");
537 
538  physiKaptonWindow = new G4PVPlacement(0, G4ThreeVector(kaptonWindowXPosition, 0., 0.),
539  "KaptonWindow", logicKaptonWindow,
540  physiVacuumZone, false, 0);
541 
542  logicKaptonWindow -> SetVisAttributes(darkOrange3);
543 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63

Field Documentation

G4double PassiveCarbonBeamLine::firstScatteringFoilXPosition

Definition at line 80 of file PassiveCarbonBeamLine.hh.

Referenced by VacuumToAirInterface().

G4double PassiveCarbonBeamLine::firstScatteringFoilYPosition

Definition at line 81 of file PassiveCarbonBeamLine.hh.

G4double PassiveCarbonBeamLine::firstScatteringFoilZPosition

Definition at line 82 of file PassiveCarbonBeamLine.hh.

G4Material* PassiveCarbonBeamLine::kapton

Definition at line 77 of file PassiveCarbonBeamLine.hh.

G4VPhysicalVolume* PassiveCarbonBeamLine::mother

Definition at line 79 of file PassiveCarbonBeamLine.hh.

G4VisAttributes* PassiveCarbonBeamLine::redWire

Definition at line 78 of file PassiveCarbonBeamLine.hh.


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