Loading...
Searching...
No Matches
FFPrimaryGeneratorAction.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// =============== Begin Documentation Comments ===============
28//!
29//! \file FFPrimaryGeneratorAction.cc
30//! \author B. Wendt (brycen.linn.wendt@cern.ch)
31//! \date June 06, 2014
32//!
33//! \brief Implementation of the FFPrimaryGeneratorAction class
34//!
35//! \details Generates 4.5 MeV neutrons from the solid volume defined in
36//! FFDetectorConstruction as "NeturonSource" and shoots them off
37//! in a isotropically sampled direction
38//!
39// ================ End Documentation Comments ================
40//
41// Modified:
42//
43// 23-06-14 BWendt
44// Implemented "GetNeutronSourceCenter()" and added code to correctly sample
45// the neutron starting location
46//
47// -------------------------------------------------------------
48
49#include "globals.hh"
50
51#include "G4Event.hh"
52#include "G4LogicalVolume.hh"
53#include "G4LogicalVolumeStore.hh"
54#include "G4Neutron.hh"
55#include "G4ParticleGun.hh"
56#include "G4PhysicalVolumeStore.hh"
57#include "G4SystemOfUnits.hh"
58#include "G4Tubs.hh"
59#include "G4VPhysicalVolume.hh"
60#include "Randomize.hh"
61
63
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69#ifndef NDEBUG
70 fEventNumber(0),
71#endif // NDEBUG
72 fH2OPhysical(NULL),
73 fNeutronPhysical(NULL),
74 fNeutronSolid(NULL),
75 fParticleGun(new G4ParticleGun(1)),
76 fTankPhysical(NULL)
77{
78 fParticleGun->SetParticleDefinition(G4Neutron::Definition());
79 fParticleGun->SetParticleEnergy(4.5 * MeV);
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85{
86#ifndef NDEBUG
87 G4cout << "Shooting event " << ++fEventNumber << G4endl;
88#endif // NDEBUG
89 const G4ThreeVector sourceCenter = GetNeutronSourceCenter();
90
91 // Sample the neutron source location
92 const G4double radius = fNeutronSolid->GetOuterRadius();
93 const G4double z = fNeutronSolid->GetZHalfLength() * 2;
94 G4ThreeVector randomLocation;
95 randomLocation.setRThetaPhi(radius * std::sqrt(G4UniformRand()),
96 G4UniformRand() * 180 * deg,
97 0);
98 randomLocation.setZ(z * (G4UniformRand() - 0.5));
99 G4ThreeVector location(randomLocation.x() + sourceCenter.x(),
100 randomLocation.y() + sourceCenter.y(),
101 randomLocation.z() + sourceCenter.z());
102#ifndef NDEBUG
103 G4cout << "Emission Location: r: " << location << G4endl;
104#endif // NDEBUG
105
106 // Sample the neutron emission direction
107 G4ThreeVector direction;
108 direction.setRThetaPhi(1.0,
109 std::acos(G4UniformRand() * 2 - 1),
110 (G4UniformRand() * 2 - 1) * 180 * deg);
111#ifndef NDEBUG
112 G4cout << "Emission Direction: r: " << direction << G4endl;
113#endif // NDEBUG
114
115 // Load the event
116 fParticleGun->SetParticlePosition(location);
117 fParticleGun->SetParticleMomentumDirection(direction);
118 fParticleGun->GeneratePrimaryVertex(event);
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124{
125 // Get the dimensions of the neutron source
126 if(fNeutronSolid == NULL)
127 {
128 G4LogicalVolume* temp
129 = G4LogicalVolumeStore::GetInstance()->GetVolume("NeutronSource");
130 if(temp != NULL)
131 {
132 fNeutronSolid = dynamic_cast< G4Tubs* >(temp->GetSolid());
133 }
134
135 if(fNeutronSolid == NULL)
136 {
137 G4Exception("FFPrimaryGeneratorAction::"
138 "GeneratePrimaries(G4Event*)",
139 "Neutron source solid volume not found",
140 EventMustBeAborted,
141 "This run will be aborted");
142 }
143 }
144
145 // Get the position of the neutron source within the water
146 if(fNeutronPhysical == NULL)
147 {
148 fNeutronPhysical = G4PhysicalVolumeStore::GetInstance()
149 ->GetVolume("NeutronSource");
150 }
151
152 if(fNeutronPhysical == NULL)
153 {
154 G4Exception("FFPrimaryGeneratorAction::GetNeutronSourceCenter(void)",
155 "Neutron source physical volume not found",
156 EventMustBeAborted,
157 "This run will be aborted");
158 }
159
160 // Get the position of the water within the tank
161 if(fH2OPhysical == NULL)
162 {
163 fH2OPhysical = G4PhysicalVolumeStore::GetInstance()
164 ->GetVolume("Tank_H2O");
165
166 if(fH2OPhysical == NULL)
167 {
168 G4Exception("FFPrimaryGeneratorAction::"
169 "GetNeutronSourceCenter(void)",
170 "Tank H2O physical volume not found",
171 EventMustBeAborted,
172 "This run will be aborted");
173 }
174 }
175
176 // Get the position of the tank within the world
177 if(fTankPhysical == NULL)
178 {
179 fTankPhysical = G4PhysicalVolumeStore::GetInstance()
180 ->GetVolume("Tank_Wall");
181
182 if(fTankPhysical == NULL)
183 {
184 G4Exception("FFPrimaryGeneratorAction::"
185 "GetNeutronSourceCenter(void)",
186 "Tank physical volume not found",
187 EventMustBeAborted,
188 "This run will be aborted");
189 }
190 }
191
192 return fNeutronPhysical->GetTranslation()
193 + fH2OPhysical->GetTranslation()
194 + fTankPhysical->GetTranslation();
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203
204
Definition of the FFPrimaryGeneratorAction class.
G4ThreeVector GetNeutronSourceCenter(void)
virtual void GeneratePrimaries(G4Event *event)

Applications | User Support | Publications | Collaboration