File indexing completed on 2025-04-04 08:03:50
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032 #include "TETDetectorConstruction.hh"
0033
0034 #include "G4VisAttributes.hh"
0035
0036 TETDetectorConstruction::TETDetectorConstruction(TETModelImport* _tetData)
0037 :fWorldPhysical(nullptr), fContainer_logic(nullptr), fTetData(_tetData), fTetLogic(nullptr)
0038 {
0039
0040 fPhantomSize = fTetData -> GetPhantomSize();
0041 fPhantomBoxMin = fTetData -> GetPhantomBoxMin();
0042 fPhantomBoxMax = fTetData -> GetPhantomBoxMax();
0043 fNOfTetrahedrons = fTetData -> GetNumTetrahedron();
0044 }
0045
0046 TETDetectorConstruction::~TETDetectorConstruction()
0047 {
0048 delete fTetData;
0049 }
0050
0051 G4VPhysicalVolume* TETDetectorConstruction::Construct()
0052 {
0053 SetupWorldGeometry();
0054 ConstructPhantom();
0055 PrintPhantomInformation();
0056 return fWorldPhysical;
0057 }
0058
0059 void TETDetectorConstruction::SetupWorldGeometry()
0060 {
0061
0062
0063 G4double worldXYZ = 10. * m;
0064 G4Material* vacuum = G4NistManager::Instance()->FindOrBuildMaterial("G4_Galactic");
0065
0066 G4VSolid* worldSolid = new G4Box("worldSolid", worldXYZ/2, worldXYZ/2, worldXYZ/2);
0067
0068 auto* worldLogical = new G4LogicalVolume(worldSolid,vacuum,"worldLogical");
0069
0070 fWorldPhysical = new G4PVPlacement(nullptr,G4ThreeVector(), worldLogical,"worldPhysical", nullptr, false,0,false);
0071
0072
0073
0074 auto* containerSolid = new G4Box("phantomBox", fPhantomSize.x()/2 + 10.*cm,
0075 fPhantomSize.y()/2 + 10.*cm,
0076 fPhantomSize.z()/2 + 10.*cm);
0077
0078 fContainer_logic = new G4LogicalVolume(containerSolid, vacuum, "phantomLogical");
0079
0080 new G4PVPlacement(nullptr, G4ThreeVector(), fContainer_logic, "PhantomPhysical",
0081 worldLogical, false, 0);
0082
0083 fContainer_logic->SetOptimisation(TRUE);
0084 fContainer_logic->SetSmartless( 0.5 );
0085 }
0086
0087 void TETDetectorConstruction::ConstructPhantom()
0088 {
0089
0090
0091
0092 G4VSolid* tetraSolid = new G4Tet("TetSolid", G4ThreeVector(),
0093 G4ThreeVector(1.*cm,0,0),
0094 G4ThreeVector(0,1.*cm,0),
0095 G4ThreeVector(0,0,1.*cm));
0096
0097 G4Material* vacuum = G4NistManager::Instance()->FindOrBuildMaterial("G4_Galactic");
0098 fTetLogic = new G4LogicalVolume(tetraSolid, vacuum, "TetLogic");
0099
0100
0101 new G4PVParameterised("wholePhantom",fTetLogic,fContainer_logic,
0102 kUndefined, fTetData->GetNumTetrahedron(),
0103 new TETParameterisation(fTetData));
0104 }
0105
0106 void TETDetectorConstruction::ConstructSDandField()
0107 {
0108
0109
0110 G4SDManager* pSDman = G4SDManager::GetSDMpointer();
0111 G4String phantomSDname = "PhantomSD";
0112
0113
0114 auto* MFDet = new G4MultiFunctionalDetector(phantomSDname);
0115 pSDman->AddNewDetector( MFDet );
0116
0117
0118 MFDet->RegisterPrimitive(new TETPSEnergyDeposit("eDep", fTetData));
0119
0120
0121 SetSensitiveDetector(fTetLogic, MFDet);
0122 }
0123
0124 void TETDetectorConstruction::PrintPhantomInformation()
0125 {
0126
0127 G4cout<< G4endl;
0128 G4cout.precision(3);
0129 G4cout<<" Phantom name "<<fTetData->GetPhantomName() << " TET phantom"<<G4endl;
0130 G4cout<<" Phantom size "<<fPhantomSize.x()<<" * "<<fPhantomSize.y()<<" * "<<fPhantomSize.z()<<" mm3"<<G4endl;
0131 G4cout<<" Phantom box position (min) "<<fPhantomBoxMin.x()<<" mm, "<<fPhantomBoxMin.y()<<" mm, "<<fPhantomBoxMin.z()<<" mm"<<G4endl;
0132 G4cout<<" Phantom box position (max) "<<fPhantomBoxMax.x()<<" mm, "<<fPhantomBoxMax.y()<<" mm, "<<fPhantomBoxMax.z()<<" mm"<<G4endl;
0133 G4cout<<" Number of tetrahedrons "<<fNOfTetrahedrons<<G4endl<<G4endl;
0134 }