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

Applications | User Support | Publications | Collaboration