Loading...
Searching...
No Matches
GB05DetectorConstruction.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//
27/// \file GB05DetectorConstruction.cc
28/// \brief Implementation of the GB05DetectorConstruction class
29
31#include "G4SystemOfUnits.hh"
32
33#include "G4Material.hh"
34#include "G4NistManager.hh"
35
36#include "G4Box.hh"
37#include "G4LogicalVolume.hh"
38#include "G4PVPlacement.hh"
39#include "G4LogicalVolumeStore.hh"
40
41#include "G4NistManager.hh"
42
44
45#include "GB05SD.hh"
46#include "G4SDManager.hh"
47
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
63{
64 // -- Collect a few materials from NIST database:
65 auto nistManager = G4NistManager::Instance();
66 G4Material* worldMaterial = nistManager->FindOrBuildMaterial("G4_Galactic");
67 G4Material* defaultMaterial = nistManager->FindOrBuildMaterial("G4_CONCRETE");
68
69
70 G4VSolid* solidWorld = new G4Box("World", 10*m, 10*m, 10*m );
71
72 G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWorld, // its solid
73 worldMaterial, // its material
74 "World"); // its name
75
76 G4PVPlacement* physiWorld = new G4PVPlacement(0, // no rotation
77 G4ThreeVector(), // at (0,0,0)
78 logicWorld, // its logical volume
79 "World", // its name
80 0, // its mother volume
81 false, // no boolean operation
82 0); // copy number
83
84 // -----------------------------------
85 // -- volume where biasing is applied:
86 // -----------------------------------
87 G4double halfXY = 5.0*m;
88 G4double halfZ = 1.0*m;
89 G4VSolid* solidShield = new G4Box("shield.solid", halfXY, halfXY, halfZ );
90
91 G4LogicalVolume* logicShield = new G4LogicalVolume( solidShield, // its solid
92 defaultMaterial, // its material
93 "shield.logical"); // its name
94
95 new G4PVPlacement(0, // no rotation
96 G4ThreeVector(0,0, halfZ), // volume entrance at (0,0,0)
97 logicShield, // its logical volume
98 "shield.phys", // its name
99 logicWorld, // its mother volume
100 false, // no boolean operation
101 0); // copy number
102
103 // --------------------------------------------------------------------
104 // -- dummy volume to print information on particle exiting the shield:
105 // --------------------------------------------------------------------
106 G4double halfz = 1*cm;
107 G4VSolid* solidMeasurement = new G4Box("meas.solid", halfXY, halfXY, halfz );
108
109 G4LogicalVolume* logicMeasurement = new G4LogicalVolume(solidMeasurement,// its solid
110 worldMaterial, // its material
111 "meas.logical"); // its name
112
113 new G4PVPlacement(0, // no rotation
114 G4ThreeVector(0,0, 2*halfZ + halfz), // volume entrance at (0,0,0)
115 logicMeasurement, // its logical volume
116 "meas.phys", // its name
117 logicWorld, // its mother volume
118 false, // no boolean operation
119 0); // copy number
120
121
122 return physiWorld;
123}
124
125
127{
128 // -- Fetch volume for biasing:
129 G4LogicalVolume* logicShield =
130 G4LogicalVolumeStore::GetInstance()->GetVolume("shield.logical");
131
132 // -------------------------------------------------------------
133 // -- operator creation, configuration and attachment to volume:
134 // -------------------------------------------------------------
135 GB05BOptrSplitAndKillByCrossSection* biasingOperator =
137 // -- Now, we declare to our biasing operator all the processes we
138 // -- their disapperance effect on neutrons to be counterbalanced
139 // -- by the splitting by cross-section :
140 biasingOperator->AddProcessToEquipoise("Decay");
141 biasingOperator->AddProcessToEquipoise("nCapture");
142 biasingOperator->AddProcessToEquipoise("neutronInelastic");
143
144 biasingOperator->AttachTo(logicShield);
145
146 G4cout << " Attaching biasing operator " << biasingOperator->GetName()
147 << " to logical volume " << biasingOperator->GetName()
148 << G4endl;
149
150 // ------------------------------------------------------------------------------------
151 // -- Attach a sensitive detector to print information on particles exiting the shield:
152 // ------------------------------------------------------------------------------------
153 // -- create and register sensitive detector code module:
154 G4SDManager* SDman = G4SDManager::GetSDMpointer();
155 G4VSensitiveDetector* sd = new GB05SD("Scorer");
156 SDman->AddNewDetector(sd);
157 // -- Fetch volume for sensitivity and attach sensitive module to it:
158 G4LogicalVolume* logicSD =
159 G4LogicalVolumeStore::GetInstance()->GetVolume("meas.logical");
160 logicSD->SetSensitiveDetector(sd);
161
162}
Definition of the GB05BOptrSplitAndKillByCrossSection class.
Definition of the GB05DetectorConstruction class.
Definition of the GB05SD class.
virtual G4VPhysicalVolume * Construct()

Applications | User Support | Publications | Collaboration