Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 // Author: Haegin Han 28 // Reference: ICRP Publication 145. Ann. ICRP 49(3), 2020. 29 // Geant4 Contributors: J. Allison and S. Guatelli 30 // 31 32 #include "TETDetectorConstruction.hh" 33 34 #include "G4VisAttributes.hh" 35 36 TETDetectorConstruction::TETDetectorConstruction(TETModelImport* _tetData) 37 :fWorldPhysical(nullptr), fContainer_logic(nullptr), fTetData(_tetData), fTetLogic(nullptr) 38 { 39 // initialisation of the variables for phantom information 40 fPhantomSize = fTetData -> GetPhantomSize(); 41 fPhantomBoxMin = fTetData -> GetPhantomBoxMin(); 42 fPhantomBoxMax = fTetData -> GetPhantomBoxMax(); 43 fNOfTetrahedrons = fTetData -> GetNumTetrahedron(); 44 } 45 46 TETDetectorConstruction::~TETDetectorConstruction() 47 { 48 delete fTetData; 49 } 50 51 G4VPhysicalVolume* TETDetectorConstruction::Construct() 52 { 53 SetupWorldGeometry(); 54 ConstructPhantom(); 55 PrintPhantomInformation(); 56 return fWorldPhysical; 57 } 58 59 void TETDetectorConstruction::SetupWorldGeometry() 60 { 61 // Define the world box (size: 10*10*10 m3) 62 // 63 G4double worldXYZ = 10. * m; 64 G4Material* vacuum = G4NistManager::Instance()->FindOrBuildMaterial("G4_Galactic"); 65 66 G4VSolid* worldSolid = new G4Box("worldSolid", worldXYZ/2, worldXYZ/2, worldXYZ/2); 67 68 auto* worldLogical = new G4LogicalVolume(worldSolid,vacuum,"worldLogical"); 69 70 fWorldPhysical = new G4PVPlacement(nullptr,G4ThreeVector(), worldLogical,"worldPhysical", nullptr, false,0,false); 71 72 // Define the phantom container (10-cm margins from the bounding box of phantom) 73 // 74 auto* containerSolid = new G4Box("phantomBox", fPhantomSize.x()/2 + 10.*cm, 75 fPhantomSize.y()/2 + 10.*cm, 76 fPhantomSize.z()/2 + 10.*cm); 77 78 fContainer_logic = new G4LogicalVolume(containerSolid, vacuum, "phantomLogical"); 79 80 new G4PVPlacement(nullptr, G4ThreeVector(), fContainer_logic, "PhantomPhysical", 81 worldLogical, false, 0); 82 83 fContainer_logic->SetOptimisation(TRUE); 84 fContainer_logic->SetSmartless( 0.5 ); // for optimization (default=2) 85 } 86 87 void TETDetectorConstruction::ConstructPhantom() 88 { 89 // Define the tetrahedral mesh phantom as a parameterised geometry 90 // 91 // solid and logical volume to be used for parameterised geometry 92 G4VSolid* tetraSolid = new G4Tet("TetSolid", G4ThreeVector(), 93 G4ThreeVector(1.*cm,0,0), 94 G4ThreeVector(0,1.*cm,0), 95 G4ThreeVector(0,0,1.*cm)); 96 97 G4Material* vacuum = G4NistManager::Instance()->FindOrBuildMaterial("G4_Galactic"); 98 fTetLogic = new G4LogicalVolume(tetraSolid, vacuum, "TetLogic"); 99 100 // physical volume (phantom) constructed as parameterised geometry 101 new G4PVParameterised("wholePhantom",fTetLogic,fContainer_logic, 102 kUndefined, fTetData->GetNumTetrahedron(), 103 new TETParameterisation(fTetData)); 104 } 105 106 void TETDetectorConstruction::ConstructSDandField() 107 { 108 // Define detector (Phantom SD) and scorer (eDep) 109 // 110 G4SDManager* pSDman = G4SDManager::GetSDMpointer(); 111 G4String phantomSDname = "PhantomSD"; 112 113 // MultiFunctional detector 114 auto* MFDet = new G4MultiFunctionalDetector(phantomSDname); 115 pSDman->AddNewDetector( MFDet ); 116 117 // scorer for energy depositon in each organ 118 MFDet->RegisterPrimitive(new TETPSEnergyDeposit("eDep", fTetData)); 119 120 // attach the detector to logical volume for parameterised geometry (phantom geometry) 121 SetSensitiveDetector(fTetLogic, MFDet); 122 } 123 124 void TETDetectorConstruction::PrintPhantomInformation() 125 { 126 // print brief information on the imported phantom 127 G4cout<< G4endl; 128 G4cout.precision(3); 129 G4cout<<" Phantom name "<<fTetData->GetPhantomName() << " TET phantom"<<G4endl; 130 G4cout<<" Phantom size "<<fPhantomSize.x()<<" * "<<fPhantomSize.y()<<" * "<<fPhantomSize.z()<<" mm3"<<G4endl; 131 G4cout<<" Phantom box position (min) "<<fPhantomBoxMin.x()<<" mm, "<<fPhantomBoxMin.y()<<" mm, "<<fPhantomBoxMin.z()<<" mm"<<G4endl; 132 G4cout<<" Phantom box position (max) "<<fPhantomBoxMax.x()<<" mm, "<<fPhantomBoxMax.y()<<" mm, "<<fPhantomBoxMax.z()<<" mm"<<G4endl; 133 G4cout<<" Number of tetrahedrons "<<fNOfTetrahedrons<<G4endl<<G4endl; 134 } 135