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 parallel/ThreadsafeScorers/src/TSDetectorConstruction.cc 27 /// \brief Implementation of the TSDetectorConstruction class 28 // 29 // 30 // 31 // 32 /// Construction of a target material (default = boron) surrounded by a 33 /// casing material (default = water) and a vacuum world (default = 34 /// target and casing fill world). The target + casing is brick 35 /// geometry with fTargetSections defining the number of divisions 36 /// in each dimension. The end sections in each dimension 37 /// is set to the casing. So a fTargetSections = G4ThreeVector(3, 3, 3) 38 /// would be one section of boron and 8 sections of water. 39 /// The idea behind this geometry is just to create a simple geometry that 40 /// scatters and produces a lot neutrons with a minimal number of sections 41 /// (i.e. coarse meshing) such that the contention in operating on 42 /// the atomic hits maps is higher and round-off errors in the 43 /// thread-local hits maps are detectable (printed out in TSRunAction) 44 /// from the sheer number of floating point sum operations. 45 /// Two scorers are implemented: EnergyDeposit and Number of steps 46 /// The energy deposit is to (possibly) show the round-off error seen 47 /// with thread-local hits maps. The # of steps scorer is to verify 48 /// the thread-safe and thread-local hits maps provide the same results. 49 // 50 // 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 53 54 #include "TSDetectorConstruction.hh" 55 56 #include "G4Box.hh" 57 #include "G4Colour.hh" 58 #include "G4LogicalVolume.hh" 59 #include "G4Material.hh" 60 #include "G4MultiFunctionalDetector.hh" 61 #include "G4NistManager.hh" 62 #include "G4PSEnergyDeposit.hh" 63 #include "G4PSNofStep.hh" 64 #include "G4PVPlacement.hh" 65 #include "G4RunManager.hh" 66 #include "G4SDManager.hh" 67 #include "G4UnitsTable.hh" 68 #include "G4UserLimits.hh" 69 #include "G4VPhysicalVolume.hh" 70 #include "G4VisAttributes.hh" 71 72 using namespace CLHEP; 73 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 75 76 TSDetectorConstruction* TSDetectorConstruction::fgInstance = 0; 77 78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 79 80 TSDetectorConstruction* TSDetectorConstruction::Instance() 81 { 82 return fgInstance; 83 } 84 85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 86 87 TSDetectorConstruction::TSDetectorConstruction() 88 : fWorldPhys(0), 89 fWorldMaterialName("G4_Galactic"), 90 fTargetMaterialName("G4_B"), 91 fCasingMaterialName("G4_WATER"), 92 fWorldDim(G4ThreeVector(0.5 * m, 0.5 * m, 0.5 * m)), 93 fTargetDim(G4ThreeVector(0.5 * m, 0.5 * m, 0.5 * m)), 94 fTargetSections(G4ThreeVector(5, 5, 5)), 95 fMfdName("Target_MFD") 96 { 97 fgInstance = this; 98 } 99 100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 101 102 TSDetectorConstruction::~TSDetectorConstruction() 103 { 104 fgInstance = 0; 105 } 106 107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 108 109 G4VPhysicalVolume* TSDetectorConstruction::Construct() 110 { 111 return ConstructWorld(ConstructMaterials()); 112 } 113 114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 115 116 TSDetectorConstruction::MaterialCollection_t TSDetectorConstruction::ConstructMaterials() 117 { 118 MaterialCollection_t materials; 119 G4NistManager* nist = G4NistManager::Instance(); 120 121 materials["World"] = nist->FindOrBuildMaterial(fWorldMaterialName); 122 materials["Target"] = nist->FindOrBuildMaterial(fTargetMaterialName); 123 materials["Casing"] = nist->FindOrBuildMaterial(fCasingMaterialName); 124 125 return materials; 126 } 127 128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 129 130 G4VPhysicalVolume* TSDetectorConstruction::ConstructWorld(const MaterialCollection_t& materials) 131 { 132 G4UserLimits* steplimit = new G4UserLimits(0.1 * (fTargetDim.z() / fTargetSections.z())); 133 G4bool check_overlap = false; 134 135 G4Box* world_solid = 136 new G4Box("World", 0.5 * fWorldDim.x(), 0.5 * fWorldDim.y(), 0.5 * fWorldDim.z()); 137 G4LogicalVolume* world_log = 138 new G4LogicalVolume(world_solid, materials.find("World")->second, "World"); 139 fWorldPhys = 140 new G4PVPlacement(0, G4ThreeVector(0.), "World", world_log, 0, false, 0, check_overlap); 141 142 G4int nz = fTargetSections.z(); 143 G4int ny = fTargetSections.y(); 144 G4int nx = fTargetSections.x(); 145 146 // spacing between sections 147 G4double sx = fTargetDim.x() / fTargetSections.x(); 148 G4double sy = fTargetDim.y() / fTargetSections.y(); 149 G4double sz = fTargetDim.z() / fTargetSections.z(); 150 151 // G4cout << "World has dimensions : " 152 //<< G4BestUnit(fWorldDim, "Length") << G4endl; 153 154 //------------------------------------------------------------------------// 155 // Set Visual Attributes 156 //------------------------------------------------------------------------// 157 G4VisAttributes* red = new G4VisAttributes(G4Color(1., 0., 0., 1.0)); 158 G4VisAttributes* green = new G4VisAttributes(G4Color(0., 1., 0., 0.25)); 159 G4VisAttributes* blue = new G4VisAttributes(G4Color(0., 0., 1., 0.1)); 160 G4VisAttributes* white = new G4VisAttributes(G4Color(1., 1., 1., 1.)); 161 162 white->SetVisibility(true); 163 red->SetVisibility(true); 164 green->SetVisibility(true); 165 blue->SetVisibility(true); 166 167 white->SetForceWireframe(true); 168 red->SetForceSolid(true); 169 green->SetForceSolid(true); 170 blue->SetForceSolid(true); 171 172 world_log->SetVisAttributes(white); 173 174 for (G4int k = 0; k < nz; ++k) 175 for (G4int j = 0; j < ny; ++j) 176 for (G4int i = 0; i < nx; ++i) { 177 // displacement of section 178 G4double dx = 0.5 * sx + static_cast<G4double>(i) * sx - 0.5 * fWorldDim.x(); 179 G4double dy = 0.5 * sy + static_cast<G4double>(j) * sy - 0.5 * fWorldDim.y(); 180 G4double dz = 0.5 * sz + static_cast<G4double>(k) * sz - 0.5 * fWorldDim.z(); 181 G4ThreeVector td = G4ThreeVector(dx, dy, -dz); 182 // make unique name 183 std::stringstream ss_name; 184 ss_name << "Target_" << i << "_" << j << "_" << k; 185 186 G4Box* target_solid = new G4Box(ss_name.str(), 0.5 * sx, 0.5 * sy, 0.5 * sz); 187 188 G4Material* target_material = 0; 189 G4bool is_casing = true; 190 191 if (j == 0 || j + 1 == ny || i == 0 || i + 1 == nx || (nz > 1 && (k == 0 || k + 1 == nz))) 192 target_material = materials.find("Casing")->second; 193 else { 194 target_material = materials.find("Target")->second; 195 is_casing = false; 196 } 197 198 G4LogicalVolume* target_log = 199 new G4LogicalVolume(target_solid, target_material, ss_name.str()); 200 201 target_log->SetUserLimits(steplimit); 202 203 new G4PVPlacement(0, td, ss_name.str(), target_log, fWorldPhys, true, 204 k * nx * ny + j * nx + i, check_overlap); 205 206 fScoringVolumes.insert(target_log); 207 208 if (is_casing) 209 target_log->SetVisAttributes(blue); 210 else { 211 // making a checkerboard for kicks... 212 G4bool even_z = (k % 2 == 0) ? true : false; 213 G4bool even_y = (j % 2 == 0) ? true : false; 214 G4bool even_x = (i % 2 == 0) ? true : false; 215 216 G4VisAttributes* theColor = nullptr; 217 218 if ((even_z)) { 219 if ((even_y && even_x) || (!even_y && !even_x)) 220 theColor = red; 221 else 222 theColor = green; 223 } 224 else // ! even_z 225 { 226 if ((!even_y && even_x) || (even_y && !even_x)) 227 theColor = red; 228 else 229 theColor = green; 230 } 231 232 target_log->SetVisAttributes(theColor); 233 } 234 } 235 236 return fWorldPhys; 237 } 238 239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 240 241 void TSDetectorConstruction::ConstructSDandField() 242 { 243 //------------------------------------------------// 244 // MultiFunctionalDetector // 245 //------------------------------------------------// 246 // Define MultiFunctionalDetector with name. 247 G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(fMfdName); 248 G4SDManager::GetSDMpointer()->AddNewDetector(MFDet); 249 G4VPrimitiveScorer* edep = new G4PSEnergyDeposit("EnergyDeposit"); 250 MFDet->RegisterPrimitive(edep); 251 G4VPrimitiveScorer* nstep = new G4PSNofStep("NumberOfSteps"); 252 MFDet->RegisterPrimitive(nstep); 253 254 // add scoring volumes 255 for (auto ite : fScoringVolumes) { 256 SetSensitiveDetector(ite, MFDet); 257 } 258 } 259 260 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 261