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 /// \file runAndEvent/RE02/src/RE02DetectorConstruction.cc 27 /// \brief Implementation of the RE02DetectorConstruction class 28 // 29 // 30 // 31 32 #include "RE02DetectorConstruction.hh" 33 34 #include "RE02NestedPhantomParameterisation.hh" 35 36 #include "G4Box.hh" 37 #include "G4Colour.hh" 38 #include "G4LogicalVolume.hh" 39 #include "G4Material.hh" 40 #include "G4NistManager.hh" 41 #include "G4PSCellFlux3D.hh" 42 #include "G4PSEnergyDeposit3D.hh" 43 #include "G4PSFlatSurfaceCurrent3D.hh" 44 #include "G4PSFlatSurfaceFlux3D.hh" 45 #include "G4PSNofStep3D.hh" 46 #include "G4PSPassageCellFlux3D.hh" 47 #include "G4PVParameterised.hh" 48 #include "G4PVPlacement.hh" 49 #include "G4SDChargedFilter.hh" 50 #include "G4SDManager.hh" 51 #include "G4SDParticleFilter.hh" 52 #include "G4SDParticleWithEnergyFilter.hh" 53 #include "G4SystemOfUnits.hh" 54 #include "G4VisAttributes.hh" 55 #include "G4ios.hh" 56 57 //======================================================================= 58 // RE02DetectorConstruction 59 // 60 // (Description) 61 // 62 // Detector construction for example RE02. 63 // 64 // [Geometry] 65 // The world volume is defined as 200 cm x 200 cm x 200 cm box with Air. 66 // Water phantom is defined as 200 mm x 200 mm x 400 mm box with Water. 67 // The water phantom is divided into 100 segments in x,y plane using 68 // replication, 69 // and then divided into 200 segments perpendicular to z axis using nested 70 // parameterised volume. 71 // These values are defined at constructor, 72 // e.g. the size of water phantom (fPhantomSize), and number of segmentation 73 // of water phantom (fNx, fNy, fNz). 74 // 75 // By default, lead plates are inserted into the position of even order 76 // segments. 77 // NIST database is used for materials. 78 // 79 // 80 // [Scorer] 81 // Assignment of G4MultiFunctionalDetector and G4PrimitiveScorer 82 // is demonstrated in this example. 83 // ------------------------------------------------- 84 // The collection names of defined Primitives are 85 // 0 PhantomSD/totalEDep 86 // 1 PhantomSD/protonEDep 87 // 2 PhantomSD/protonNStep 88 // 3 PhantomSD/chargedPassCellFlux 89 // 4 PhantomSD/chargedCellFlux 90 // 5 PhantomSD/chargedSurfFlux 91 // 6 PhantomSD/gammaSurfCurr000 92 // 7 PhantomSD/gammaSurfCurr001 93 // 9 PhantomSD/gammaSurdCurr002 94 // 10 PhantomSD/gammaSurdCurr003 95 // ------------------------------------------------- 96 // Please see README for detail description. 97 // 98 //======================================================================= 99 100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 101 RE02DetectorConstruction::RE02DetectorConstruction() : G4VUserDetectorConstruction() 102 { 103 // Default size of water phantom,and segmentation. 104 fPhantomSize.setX(200. * mm); 105 fPhantomSize.setY(200. * mm); 106 fPhantomSize.setZ(400. * mm); 107 fNx = fNy = fNz = 100; 108 fInsertLead = TRUE; 109 } 110 111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 112 RE02DetectorConstruction::~RE02DetectorConstruction() 113 { 114 ; 115 } 116 117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 118 G4VPhysicalVolume* RE02DetectorConstruction::Construct() 119 { 120 //===================== 121 // Material Definitions 122 //===================== 123 // 124 //-------- NIST Materials ---------------------------------------------------- 125 // Material Information imported from NIST database. 126 // 127 G4NistManager* NISTman = G4NistManager::Instance(); 128 G4Material* air = NISTman->FindOrBuildMaterial("G4_AIR"); 129 G4Material* water = NISTman->FindOrBuildMaterial("G4_WATER"); 130 G4Material* lead = NISTman->FindOrBuildMaterial("G4_Pb"); 131 132 // 133 // Print all the materials defined. 134 G4cout << G4endl << "The materials defined are : " << G4endl << G4endl; 135 G4cout << *(G4Material::GetMaterialTable()) << G4endl; 136 137 //============================================================================ 138 // Definitions of Solids, Logical Volumes, Physical Volumes 139 //============================================================================ 140 141 //------------- 142 // World Volume 143 //------------- 144 145 G4ThreeVector worldSize = G4ThreeVector(200 * cm, 200 * cm, 200 * cm); 146 147 G4Box* solidWorld = 148 new G4Box("world", worldSize.x() / 2., worldSize.y() / 2., worldSize.z() / 2.); 149 G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWorld, air, "World", 0, 0, 0); 150 151 // 152 // Must place the World Physical volume unrotated at (0,0,0). 153 G4VPhysicalVolume* physiWorld = new G4PVPlacement(0, // no rotation 154 G4ThreeVector(), // at (0,0,0) 155 logicWorld, // its logical volume 156 "World", // its name 157 0, // its mother volume 158 false, // no boolean operations 159 0); // copy number 160 161 //--------------- 162 // Water Phantom 163 //--------------- 164 165 //................................ 166 // Mother Volume of Water Phantom 167 //................................ 168 169 //-- Default size of water phantom is defined at constructor. 170 G4ThreeVector phantomSize = fPhantomSize; 171 172 G4Box* solidPhantom = 173 new G4Box("phantom", phantomSize.x() / 2., phantomSize.y() / 2., phantomSize.z() / 2.); 174 G4LogicalVolume* logicPhantom = new G4LogicalVolume(solidPhantom, water, "Phantom", 0, 0, 0); 175 176 G4RotationMatrix* rot = new G4RotationMatrix(); 177 // rot->rotateY(30.*deg); 178 G4ThreeVector positionPhantom; 179 // G4VPhysicalVolume * physiPhantom = 180 new G4PVPlacement(rot, // no rotation 181 positionPhantom, // at (x,y,z) 182 logicPhantom, // its logical volume 183 "Phantom", // its name 184 logicWorld, // its mother volume 185 false, // no boolean operations 186 0); // copy number 187 188 //.............................................. 189 // Phantom segmentation using Parameterisation 190 //.............................................. 191 // 192 G4cout << "<-- RE02DetectorConstruction::Construct-------" << G4endl; 193 G4cout << " Water Phantom Size " << fPhantomSize / mm << G4endl; 194 G4cout << " Segmentation (" << fNx << "," << fNy << "," << fNz << ")" << G4endl; 195 G4cout << " Lead plate at even copy # (0-False,1-True): " << IsLeadSegment() << G4endl; 196 G4cout << "<---------------------------------------------" << G4endl; 197 // Number of segmentation. 198 // - Default number of segmentation is defined at constructor. 199 G4int nxCells = fNx; 200 G4int nyCells = fNy; 201 G4int nzCells = fNz; 202 203 G4ThreeVector sensSize; 204 sensSize.setX(phantomSize.x() / (G4double)nxCells); 205 sensSize.setY(phantomSize.y() / (G4double)nyCells); 206 sensSize.setZ(phantomSize.z() / (G4double)nzCells); 207 // i.e Voxel size will be 2.0 x 2.0 x 2.0 mm3 cube by default. 208 // 209 210 // Replication of Water Phantom Volume. 211 // Y Slice 212 G4String yRepName("RepY"); 213 G4VSolid* solYRep = 214 new G4Box(yRepName, phantomSize.x() / 2., sensSize.y() / 2., phantomSize.z() / 2.); 215 G4LogicalVolume* logYRep = new G4LogicalVolume(solYRep, water, yRepName); 216 // G4PVReplica* yReplica = 217 new G4PVReplica(yRepName, logYRep, logicPhantom, kYAxis, fNy, sensSize.y()); 218 // X Slice 219 G4String xRepName("RepX"); 220 G4VSolid* solXRep = 221 new G4Box(xRepName, sensSize.x() / 2., sensSize.y() / 2., phantomSize.z() / 2.); 222 G4LogicalVolume* logXRep = new G4LogicalVolume(solXRep, water, xRepName); 223 // G4PVReplica* xReplica = 224 new G4PVReplica(xRepName, logXRep, logYRep, kXAxis, fNx, sensSize.x()); 225 226 // 227 //.................................. 228 // Voxel solid and logical volumes 229 //.................................. 230 // Z Slice 231 G4String zVoxName("phantomSens"); 232 G4VSolid* solVoxel = new G4Box(zVoxName, sensSize.x() / 2., sensSize.y() / 2., sensSize.z() / 2.); 233 fLVPhantomSens = new G4LogicalVolume(solVoxel, water, zVoxName); 234 // 235 // 236 std::vector<G4Material*> phantomMat(2, water); 237 if (IsLeadSegment()) phantomMat[1] = lead; 238 // 239 // Parameterisation for transformation of voxels. 240 // (voxel size is fixed in this example. 241 // e.g. nested parameterisation handles material and transfomation of voxels.) 242 RE02NestedPhantomParameterisation* paramPhantom = 243 new RE02NestedPhantomParameterisation(sensSize / 2., nzCells, phantomMat); 244 // G4VPhysicalVolume * physiPhantomSens = 245 new G4PVParameterised("PhantomSens", // their name 246 fLVPhantomSens, // their logical volume 247 logXRep, // Mother logical volume 248 kUndefined, // Are placed along this axis 249 nzCells, // Number of cells 250 paramPhantom); // Parameterisation. 251 // Optimization flag is avaiable for, 252 // kUndefined, kXAxis, kYAxis, kZAxis. 253 // 254 255 //=============================== 256 // Visualization attributes 257 //=============================== 258 259 G4VisAttributes* boxVisAtt = new G4VisAttributes(G4Colour(1.0, 1.0, 1.0)); 260 logicWorld->SetVisAttributes(boxVisAtt); 261 // logicWorld->SetVisAttributes(G4VisAttributes::GetInvisible()); 262 263 // Mother volume of WaterPhantom 264 G4VisAttributes* phantomVisAtt = new G4VisAttributes(G4Colour(1.0, 1.0, 0.0)); 265 logicPhantom->SetVisAttributes(phantomVisAtt); 266 267 // Replica 268 G4VisAttributes* yRepVisAtt = new G4VisAttributes(G4Colour(0.0, 1.0, 0.0)); 269 logYRep->SetVisAttributes(yRepVisAtt); 270 G4VisAttributes* xRepVisAtt = new G4VisAttributes(G4Colour(0.0, 1.0, 0.0)); 271 logXRep->SetVisAttributes(xRepVisAtt); 272 273 // Skip the visualization for those voxels. 274 fLVPhantomSens->SetVisAttributes(G4VisAttributes::GetInvisible()); 275 276 return physiWorld; 277 } 278 279 void RE02DetectorConstruction::ConstructSDandField() 280 { 281 //================================================ 282 // Sensitive detectors : MultiFunctionalDetector 283 //================================================ 284 // 285 // Sensitive Detector Manager. 286 G4SDManager* pSDman = G4SDManager::GetSDMpointer(); 287 // 288 // Sensitive Detector Name 289 G4String phantomSDname = "PhantomSD"; 290 291 //------------------------ 292 // MultiFunctionalDetector 293 //------------------------ 294 // 295 // Define MultiFunctionalDetector with name. 296 G4MultiFunctionalDetector* mFDet = new G4MultiFunctionalDetector(phantomSDname); 297 pSDman->AddNewDetector(mFDet); // Register SD to SDManager. 298 fLVPhantomSens->SetSensitiveDetector(mFDet); // Assign SD to the logical volume. 299 300 //--------------------------------------- 301 // SDFilter : Sensitive Detector Filters 302 //--------------------------------------- 303 // 304 // Particle Filter for Primitive Scorer with filter name(fltName) 305 // and particle name(particleName), 306 // or particle names are given by add("particle name"); method. 307 // 308 G4String fltName, particleName; 309 // 310 //-- proton filter 311 G4SDParticleFilter* protonFilter = 312 new G4SDParticleFilter(fltName = "protonFilter", particleName = "proton"); 313 // 314 //-- electron filter 315 G4SDParticleFilter* electronFilter = new G4SDParticleFilter(fltName = "electronFilter"); 316 electronFilter->add(particleName = "e+"); // accept electrons. 317 electronFilter->add(particleName = "e-"); // accept positorons. 318 // 319 //-- charged particle filter 320 G4SDChargedFilter* chargedFilter = new G4SDChargedFilter(fltName = "chargedFilter"); 321 322 //------------------------ 323 // PS : Primitive Scorers 324 //------------------------ 325 // Primitive Scorers are used with SDFilters according to your purpose. 326 // 327 // 328 //-- Primitive Scorer for Energy Deposit. 329 // Total, by protons, by electrons. 330 G4String psName; 331 G4PSEnergyDeposit3D* scorer0 = new G4PSEnergyDeposit3D(psName = "totalEDep", fNx, fNy, fNz); 332 G4PSEnergyDeposit3D* scorer1 = new G4PSEnergyDeposit3D(psName = "protonEDep", fNx, fNy, fNz); 333 scorer1->SetFilter(protonFilter); 334 335 // 336 //-- Number of Steps for protons 337 G4PSNofStep3D* scorer2 = new G4PSNofStep3D(psName = "protonNStep", fNx, fNy, fNz); 338 scorer2->SetFilter(protonFilter); 339 340 // 341 //-- CellFlux for charged particles 342 G4PSPassageCellFlux3D* scorer3 = 343 new G4PSPassageCellFlux3D(psName = "chargedPassCellFlux", fNx, fNy, fNz); 344 G4PSCellFlux3D* scorer4 = new G4PSCellFlux3D(psName = "chargedCellFlux", fNx, fNy, fNz); 345 G4PSFlatSurfaceFlux3D* scorer5 = 346 new G4PSFlatSurfaceFlux3D(psName = "chargedSurfFlux", fFlux_InOut, fNx, fNy, fNz); 347 scorer3->SetFilter(chargedFilter); 348 scorer4->SetFilter(chargedFilter); 349 scorer5->SetFilter(chargedFilter); 350 351 // 352 //------------------------------------------------------------ 353 // Register primitive scorers to MultiFunctionalDetector 354 //------------------------------------------------------------ 355 mFDet->RegisterPrimitive(scorer0); 356 mFDet->RegisterPrimitive(scorer1); 357 mFDet->RegisterPrimitive(scorer2); 358 mFDet->RegisterPrimitive(scorer3); 359 mFDet->RegisterPrimitive(scorer4); 360 mFDet->RegisterPrimitive(scorer5); 361 362 //======================== 363 // More additional Primitive Scoreres 364 //======================== 365 // 366 //--- Surface Current for gamma with energy bin. 367 // This example creates four primitive scorers. 368 // 4 bins with energy --- Primitive Scorer Name 369 // 1. to 10 KeV, gammaSurfCurr000 370 // 10 keV to 100 KeV, gammaSurfCurr001 371 // 100 keV to 1 MeV, gammaSurfCurr002 372 // 1 MeV to 10 MeV. gammaSurfCurr003 373 // 374 for (G4int i = 0; i < 4; i++) { 375 std::ostringstream name; 376 name << "gammaSurfCurr" << std::setfill('0') << std::setw(3) << i; 377 G4String psgName = name.str(); 378 G4double kmin = std::pow(10., (G4double)i) * keV; 379 G4double kmax = std::pow(10., (G4double)(i + 1)) * keV; 380 //-- Particle with kinetic energy filter. 381 G4SDParticleWithEnergyFilter* pkinEFilter = 382 new G4SDParticleWithEnergyFilter(fltName = "gammaE filter", kmin, kmax); 383 pkinEFilter->add("gamma"); // Accept only gamma. 384 pkinEFilter->show(); // Show accepting condition to stdout. 385 //-- Surface Current Scorer which scores number of tracks in unit area. 386 G4PSFlatSurfaceCurrent3D* scorer = 387 new G4PSFlatSurfaceCurrent3D(psgName, fCurrent_InOut, fNx, fNy, fNz); 388 scorer->SetFilter(pkinEFilter); // Assign filter. 389 mFDet->RegisterPrimitive(scorer); // Register it to MultiFunctionalDetector. 390 } 391 } 392