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 /// \file radiobiology/src/VoxelizedSensitiveDetector.cc 28 /// \brief Implementation of the RadioBio::VoxelizedSensitiveDetector class 29 30 #include "VoxelizedSensitiveDetector.hh" 31 32 #include "G4Box.hh" 33 #include "G4LogicalVolume.hh" 34 #include "G4PVReplica.hh" 35 #include "G4PhysicalVolumeStore.hh" 36 #include "G4RunManager.hh" 37 #include "G4SDManager.hh" 38 #include "G4SystemOfUnits.hh" 39 #include "G4TransportationManager.hh" 40 #include "G4VPhysicalVolume.hh" 41 42 #include "DetectorConstruction.hh" 43 #include "SD.hh" 44 #include "VoxelizedSensitiveDetectorMessenger.hh" 45 46 namespace RadioBio 47 { 48 49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 50 51 VoxelizedSensitiveDetector* VoxelizedSensitiveDetector::fInstance = nullptr; 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 VoxelizedSensitiveDetector* VoxelizedSensitiveDetector::CreateInstance(DetectorConstruction* det, 56 double xWidth, double yWidth, 57 double zWidth) 58 { 59 if (fInstance) { 60 delete fInstance; 61 G4Exception("VoxelizedSensitiveDetector::createInstance", "RecreatingVoxelization", 62 FatalException, "Creating another, new, instance of VoxelizedSensitiveDetector"); 63 } 64 fInstance = new VoxelizedSensitiveDetector(det, xWidth, yWidth, zWidth); 65 return fInstance; 66 } 67 68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 69 70 VoxelizedSensitiveDetector* VoxelizedSensitiveDetector::GetInstance() 71 { 72 return fInstance; 73 } 74 75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 76 77 VoxelizedSensitiveDetector::VoxelizedSensitiveDetector(DetectorConstruction* det, double xWidth, 78 double yWidth, double zWidth) 79 : fDetector(det), fVoxelWidthX(xWidth), fVoxelWidthY(yWidth), fVoxelWidthZ(zWidth) 80 { 81 fVoxelizedSensitiveDetectorMessenger = new VoxelizedSensitiveDetectorMessenger(this); 82 UpdateVoxelVolume(); 83 CalculateVoxelNumber(); 84 } 85 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 87 88 VoxelizedSensitiveDetector::~VoxelizedSensitiveDetector() 89 { 90 delete fVoxelizedSensitiveDetectorMessenger; 91 } 92 93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 94 95 void VoxelizedSensitiveDetector::UpdateVoxelVolume() 96 { 97 fVoxelVolume = fVoxelWidthX * fVoxelWidthY * fVoxelWidthZ; 98 fVoxelDensity = fDetector->GetMaterial()->GetDensity(); 99 fVoxelMass = fVoxelVolume * fVoxelDensity; 100 } 101 102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 103 104 void VoxelizedSensitiveDetector::SetVoxelWidth(G4ThreeVector voxWidth) 105 { 106 fVoxelWidthX = voxWidth.getX(); 107 fVoxelWidthY = voxWidth.getY(); 108 fVoxelWidthZ = voxWidth.getZ(); 109 CalculateVoxelNumber(); 110 } 111 112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 113 114 void VoxelizedSensitiveDetector::SetVoxelWidthX(G4double voxWidthX) 115 { 116 if (fVoxelWidthX == voxWidthX) return; 117 fVoxelWidthX = voxWidthX; 118 CalculateVoxelNumber(); 119 } 120 121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 122 123 void VoxelizedSensitiveDetector::SetVoxelWidthY(G4double voxWidthY) 124 { 125 if (fVoxelWidthY == voxWidthY) return; 126 fVoxelWidthY = voxWidthY; 127 CalculateVoxelNumber(); 128 } 129 130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 131 132 void VoxelizedSensitiveDetector::SetVoxelWidthZ(G4double voxWidthZ) 133 { 134 if (fVoxelWidthZ == voxWidthZ) return; 135 fVoxelWidthZ = voxWidthZ; 136 CalculateVoxelNumber(); 137 } 138 139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 140 141 // Calculte number of voxel approximating for an integer number of voxels 142 // Then recalculates voxels size according to approximations 143 void VoxelizedSensitiveDetector::CalculateVoxelNumber() 144 { 145 fVoxelNumberAlongX = G4int(fDetector->GetSizeX() / fVoxelWidthX); 146 fVoxelWidthX = fDetector->GetSizeX() / G4double(fVoxelNumberAlongX); 147 148 fVoxelNumberAlongY = G4int(fDetector->GetSizeY() / fVoxelWidthY); 149 fVoxelWidthY = fDetector->GetSizeY() / G4double(fVoxelNumberAlongY); 150 151 fVoxelNumberAlongZ = G4int(fDetector->GetSizeZ() / fVoxelWidthZ); 152 fVoxelWidthZ = fDetector->GetSizeZ() / G4double(fVoxelNumberAlongZ); 153 154 if (fVoxelNumberAlongY % 2 == 0) 155 G4Exception("VoxelizedSensitiveDetector::CalculateVoxelNumber", "VoxelNumberYEven", JustWarning, 156 "Trying to voxelize with an even number of voxels along the Y axis." 157 "Please select an odd number to prevent from warnings due to tracking"); 158 159 if (fVoxelNumberAlongZ % 2 == 0) 160 G4Exception("VoxelizedSensitiveDetector::CalculateVoxelNumber", "VoxelNumberZEven", JustWarning, 161 "Trying to voxelize with an even number of voxels along the Z axis." 162 "Please select an odd number to prevent from warnings due to tracking"); 163 164 fTotalVoxelNumber = fVoxelNumberAlongX * fVoxelNumberAlongY * fVoxelNumberAlongZ; 165 166 UpdateVoxelVolume(); 167 } 168 169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 170 171 void VoxelizedSensitiveDetector::ConstructXDivision() 172 { 173 if (fWorldLogical == nullptr) 174 G4Exception("VoxelizedSensitiveDetector::ConstructXDivision", "WorldNotInit", FatalException, 175 "Voxelizing without having a pointer to world logical volume!"); 176 177 if (!fDetector) 178 G4Exception("VoxelizedSensitiveDetector::ConstructXDivision", "DetConstInit", FatalException, 179 "Voxelizing without having a pointer to DetectorConstruction!"); 180 181 fVoxelizedDetectorXDivision = new G4Box("VoxelizedDetectorXDivision", fVoxelWidthX / 2, 182 fDetector->GetSizeY() / 2, fDetector->GetSizeZ() / 2); 183 184 fVoxelizedDetectorXDivisionLog = 185 new G4LogicalVolume(fVoxelizedDetectorXDivision, fWorldLogical->GetMaterial(), 186 "VoxelizedDetectorXDivisionLog", 0, 0, 0); 187 188 fVoxelizedDetectorXDivisionPhys = 189 new G4PVReplica("VoxelizedDetectorXDivisionPhys", fVoxelizedDetectorXDivisionLog, fWorldLogical, 190 kXAxis, fVoxelNumberAlongX, fVoxelWidthX); 191 } 192 193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 194 195 void VoxelizedSensitiveDetector::ConstructYDivision() 196 { 197 fVoxelizedDetectorYDivision = new G4Box("VoxelizedDetectorYDivision", fVoxelWidthX / 2, 198 fVoxelWidthY / 2, fDetector->GetSizeZ() / 2); 199 200 fVoxelizedDetectorYDivisionLog = 201 new G4LogicalVolume(fVoxelizedDetectorYDivision, fWorldLogical->GetMaterial(), 202 "VoxelizedDetectorYDivisionLog", 0, 0, 0); 203 204 fVoxelizedDetectorYDivisionPhys = 205 new G4PVReplica("VoxelizedDetectorYDivisionPhys", fVoxelizedDetectorYDivisionLog, 206 fVoxelizedDetectorXDivisionLog, kYAxis, fVoxelNumberAlongY, fVoxelWidthY); 207 } 208 209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 210 211 void VoxelizedSensitiveDetector::ConstructZDivision() 212 { 213 fVoxelizedDetectorZDivision = 214 new G4Box("VoxelizedDetectorZDivision", fVoxelWidthX / 2, fVoxelWidthY / 2, fVoxelWidthZ / 2); 215 216 fVoxelizedDetectorZDivisionLog = 217 new G4LogicalVolume(fVoxelizedDetectorZDivision, fWorldLogical->GetMaterial(), 218 "VoxelizedDetectorZDivisionLog", 0, 0, 0); 219 220 fVoxelizedDetectorZDivisionPhys = 221 new G4PVReplica("VoxelizedDetectorZDivisionPhys", fVoxelizedDetectorZDivisionLog, 222 fVoxelizedDetectorYDivisionPhys, kZAxis, fVoxelNumberAlongZ, fVoxelWidthZ); 223 224 fSensitiveLogicalVolume = fVoxelizedDetectorZDivisionLog; 225 } 226 227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 228 229 // First voxelize along X, then Y, then Z 230 G4bool VoxelizedSensitiveDetector::ConstructVoxelizedDetector() 231 { 232 // Creating X division 233 ConstructXDivision(); 234 235 // Creating Y division 236 ConstructYDivision(); 237 238 // Creating Z division 239 ConstructZDivision(); 240 241 // Set last, smallest volumes as sensitive 242 fSensitiveLogicalVolume = fVoxelizedDetectorZDivisionLog; 243 fIsBuilt = true; 244 245 return true; 246 } 247 248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 249 250 void VoxelizedSensitiveDetector::UpdateVoxelizedGeometry() 251 { 252 // Nothing happens if the voxelized geometry is not built. But parameters are properly set. 253 if (!fIsBuilt) { 254 return; 255 } 256 257 CalculateVoxelNumber(); 258 259 // Volume that will be deleted in order to update 260 G4VPhysicalVolume* myVol; 261 262 G4PhysicalVolumeStore* store = G4PhysicalVolumeStore::GetInstance(); 263 264 myVol = store->GetVolume("VoxelizedDetectorXDivisionPhys"); 265 store->DeRegister(myVol); 266 myVol = store->GetVolume("VoxelizedDetectorYDivisionPhys"); 267 store->DeRegister(myVol); 268 myVol = store->GetVolume("VoxelizedDetectorZDivisionPhys"); 269 store->DeRegister(myVol); 270 fVoxelizedDetectorXDivisionPhys = 271 new G4PVReplica("VoxelizedDetectorXDivisionPhys", fVoxelizedDetectorXDivisionLog, fWorldLogical, 272 kXAxis, fVoxelNumberAlongX, fVoxelWidthX); 273 274 fVoxelizedDetectorYDivisionPhys = 275 new G4PVReplica("VoxelizedDetectorYDivisionPhys", fVoxelizedDetectorYDivisionLog, 276 fVoxelizedDetectorXDivisionPhys, kYAxis, fVoxelNumberAlongY, fVoxelWidthY); 277 278 fVoxelizedDetectorZDivisionPhys = 279 new G4PVReplica("VoxelizedDetectorZDivisionPhys", fVoxelizedDetectorZDivisionLog, 280 fVoxelizedDetectorYDivisionPhys, kZAxis, fVoxelNumberAlongZ, fVoxelWidthZ); 281 282 G4RunManager::GetRunManager()->GeometryHasBeenModified(); 283 G4RunManager::GetRunManager()->PhysicsHasBeenModified(); 284 } 285 286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 287 288 void VoxelizedSensitiveDetector::ConstructSD() 289 { 290 G4String sensitiveDetectorName = "VoxelizedDetector"; 291 G4String HCname = "LETdata"; 292 293 SD* detectorSD = new SD(sensitiveDetectorName, HCname); 294 G4SDManager::GetSDMpointer()->AddNewDetector(detectorSD); 295 fSensitiveLogicalVolume->SetSensitiveDetector(detectorSD); 296 } 297 298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 299 300 void VoxelizedSensitiveDetector::Construct() 301 { 302 ConstructVoxelizedDetector(); 303 } 304 305 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 306 307 void VoxelizedSensitiveDetector::InitializeWorldPtr(G4VPhysicalVolume* pWorld) 308 { 309 if (pWorld == nullptr) 310 G4Exception("VoxelizedSensitiveDetector::InitializeWorldPtr", "WorldinitNull", FatalException, 311 "Initializing Voxelization Class with a Null Pointer to World!"); 312 fWorldLogical = pWorld->GetLogicalVolume(); 313 } 314 315 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 316 317 // Returns absolute voxel index given matrix indexes 318 G4int VoxelizedSensitiveDetector::GetThisVoxelNumber(G4int x, G4int y, G4int z) const 319 { 320 G4int nz = GetVoxelNumberAlongZ(); 321 G4int ny = GetVoxelNumberAlongY(); 322 323 return z + nz * (y + ny * (x)); 324 } 325 326 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 327 328 } // namespace RadioBio 329