Loading...
Searching...
No Matches
TSDetectorConstruction.cc
Go to the documentation of this file.
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
55
56#include "G4RunManager.hh"
57#include "G4Box.hh"
58#include "G4LogicalVolume.hh"
59#include "G4VPhysicalVolume.hh"
60#include "G4Material.hh"
61#include "G4NistManager.hh"
62#include "G4PVPlacement.hh"
63#include "G4VisAttributes.hh"
64#include "G4Colour.hh"
65#include "G4UnitsTable.hh"
66#include "G4UserLimits.hh"
67
68#include "G4SDManager.hh"
69#include "G4MultiFunctionalDetector.hh"
70#include "G4PSEnergyDeposit.hh"
71#include "G4PSNofStep.hh"
72
73using namespace CLHEP;
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
89 : fWorldPhys(0)
90 , fWorldMaterialName("G4_Galactic")
91 , fTargetMaterialName("G4_B")
92 , fCasingMaterialName("G4_WATER")
93 , fWorldDim(G4ThreeVector(0.5 * m, 0.5 * m, 0.5 * m))
94 , fTargetDim(G4ThreeVector(0.5 * m, 0.5 * m, 0.5 * m))
95 , fTargetSections(G4ThreeVector(5, 5, 5))
96 , fMfdName("Target_MFD")
97{
98 fgInstance = this;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113
116{
117 MaterialCollection_t materials;
118 G4NistManager* nist = G4NistManager::Instance();
119
120 materials["World"] = nist->FindOrBuildMaterial(fWorldMaterialName);
121 materials["Target"] = nist->FindOrBuildMaterial(fTargetMaterialName);
122 materials["Casing"] = nist->FindOrBuildMaterial(fCasingMaterialName);
123
124 return materials;
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128
130 const MaterialCollection_t& materials)
131{
132 G4UserLimits* steplimit =
133 new G4UserLimits(0.1 * (fTargetDim.z() / fTargetSections.z()));
134 G4bool check_overlap = false;
135
136 G4Box* world_solid = new G4Box("World", 0.5 * fWorldDim.x(),
137 0.5 * fWorldDim.y(), 0.5 * fWorldDim.z());
138 G4LogicalVolume* world_log =
139 new G4LogicalVolume(world_solid, materials.find("World")->second, "World");
140 fWorldPhys = new G4PVPlacement(0, G4ThreeVector(0.), "World", world_log, 0,
141 false, 0, check_overlap);
142
143 G4int nz = fTargetSections.z();
144 G4int ny = fTargetSections.y();
145 G4int nx = fTargetSections.x();
146
147 // spacing between sections
148 G4double sx = fTargetDim.x() / fTargetSections.x();
149 G4double sy = fTargetDim.y() / fTargetSections.y();
150 G4double sz = fTargetDim.z() / fTargetSections.z();
151
152 // G4cout << "World has dimensions : "
153 //<< G4BestUnit(fWorldDim, "Length") << G4endl;
154
155 //------------------------------------------------------------------------//
156 // Set Visual Attributes
157 //------------------------------------------------------------------------//
158 G4VisAttributes* red = new G4VisAttributes(G4Color(1., 0., 0., 1.0));
159 G4VisAttributes* green = new G4VisAttributes(G4Color(0., 1., 0., 0.25));
160 G4VisAttributes* blue = new G4VisAttributes(G4Color(0., 0., 1., 0.1));
161 G4VisAttributes* white = new G4VisAttributes(G4Color(1., 1., 1., 1.));
162
163 white->SetVisibility(true);
164 red->SetVisibility(true);
165 green->SetVisibility(true);
166 blue->SetVisibility(true);
167
168 white->SetForceWireframe(true);
169 red->SetForceSolid(true);
170 green->SetForceSolid(true);
171 blue->SetForceSolid(true);
172
173 world_log->SetVisAttributes(white);
174
175 for(G4int k = 0; k < nz; ++k)
176 for(G4int j = 0; j < ny; ++j)
177 for(G4int i = 0; i < nx; ++i)
178 {
179 // displacement of section
180 G4double dx =
181 0.5 * sx + static_cast<G4double>(i) * sx - 0.5 * fWorldDim.x();
182 G4double dy =
183 0.5 * sy + static_cast<G4double>(j) * sy - 0.5 * fWorldDim.y();
184 G4double dz =
185 0.5 * sz + static_cast<G4double>(k) * sz - 0.5 * fWorldDim.z();
186 G4ThreeVector td = G4ThreeVector(dx, dy, -dz);
187 // make unique name
188 std::stringstream ss_name;
189 ss_name << "Target_" << i << "_" << j << "_" << k;
190
191 G4Box* target_solid =
192 new G4Box(ss_name.str(), 0.5 * sx, 0.5 * sy, 0.5 * sz);
193
194 G4Material* target_material = 0;
195 G4bool is_casing = true;
196
197 if(j == 0 || j + 1 == ny || i == 0 || i + 1 == nx ||
198 (nz > 1 && (k == 0 || k + 1 == nz)))
199 target_material = materials.find("Casing")->second;
200 else
201 {
202 target_material = materials.find("Target")->second;
203 is_casing = false;
204 }
205
206 G4LogicalVolume* target_log =
207 new G4LogicalVolume(target_solid, target_material, ss_name.str());
208
209 target_log->SetUserLimits(steplimit);
210
211 new G4PVPlacement(0, td, ss_name.str(), target_log, fWorldPhys, true,
212 k * nx * ny + j * nx + i, check_overlap);
213
214 fScoringVolumes.insert(target_log);
215
216 if(is_casing)
217 target_log->SetVisAttributes(blue);
218 else
219 {
220 // making a checkerboard for kicks...
221 G4bool even_z = (k % 2 == 0) ? true : false;
222 G4bool even_y = (j % 2 == 0) ? true : false;
223 G4bool even_x = (i % 2 == 0) ? true : false;
224
225 G4VisAttributes* theColor = nullptr;
226
227 if((even_z))
228 {
229 if((even_y && even_x) || (!even_y && !even_x))
230 theColor = red;
231 else
232 theColor = green;
233 }
234 else // ! even_z
235 {
236 if((!even_y && even_x) || (even_y && !even_x))
237 theColor = red;
238 else
239 theColor = green;
240 }
241
242 target_log->SetVisAttributes(theColor);
243 }
244 }
245
246 return fWorldPhys;
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250
252{
253 //------------------------------------------------//
254 // MultiFunctionalDetector //
255 //------------------------------------------------//
256 // Define MultiFunctionalDetector with name.
258 G4SDManager::GetSDMpointer()->AddNewDetector(MFDet);
259 G4VPrimitiveScorer* edep = new G4PSEnergyDeposit("EnergyDeposit");
260 MFDet->RegisterPrimitive(edep);
261 G4VPrimitiveScorer* nstep = new G4PSNofStep("NumberOfSteps");
262 MFDet->RegisterPrimitive(nstep);
263
264 // add scoring volumes
265 for(auto ite : fScoringVolumes)
266 {
267 SetSensitiveDetector(ite, MFDet);
268 }
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the TSDetectorConstruction class.
std::map< G4String, G4Material * > MaterialCollection_t
static TSDetectorConstruction * Instance()
static TSDetectorConstruction * fgInstance
virtual G4VPhysicalVolume * ConstructWorld(const MaterialCollection_t &)
G4VPhysicalVolume * Construct()
virtual MaterialCollection_t ConstructMaterials()

Applications | User Support | Publications | Collaboration