Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SingleParticleSource.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// G4SingleParticleSource source implementation
27//
28// Author: Fan Lei, QinetiQ ltd. - 05/02/2004
29// Customer: ESA/ESTEC
30// Revision: Andrea Dotti, SLAC
31// --------------------------------------------------------------------
32
33#include <cmath>
34
36
37#include "G4SystemOfUnits.hh"
38#include "G4PrimaryParticle.hh"
39#include "G4Event.hh"
40#include "Randomize.hh"
41#include "G4ParticleTable.hh"
42#include "G4Geantino.hh"
44#include "G4IonTable.hh"
45#include "G4Ions.hh"
46#include "G4TrackingManager.hh"
47#include "G4Track.hh"
48#include "G4AutoLock.hh"
49
50G4SingleParticleSource::part_prop_t::part_prop_t()
51{
52 momentum_direction = G4ParticleMomentum(1,0,0);
53 energy = 1.*MeV;
55}
56
58{
59 // Initialise all variables
60 // Position distribution Variables
61
62 NumberOfParticlesToBeGenerated = 1;
63 definition = G4Geantino::GeantinoDefinition();
64
65 charge = 0.0;
66 time = 0;
67 polarization = G4ThreeVector();
68
69 biasRndm = new G4SPSRandomGenerator();
70 posGenerator = new G4SPSPosDistribution();
71 posGenerator->SetBiasRndm(biasRndm);
72 angGenerator = new G4SPSAngDistribution();
73 angGenerator->SetPosDistribution(posGenerator);
74 angGenerator->SetBiasRndm(biasRndm);
75 eneGenerator = new G4SPSEneDistribution();
76 eneGenerator->SetBiasRndm(biasRndm);
77
78 verbosityLevel = 0;
79
81}
82
84{
85 delete biasRndm;
86 delete posGenerator;
87 delete angGenerator;
88 delete eneGenerator;
89
91}
92
94{
95 G4AutoLock l(&mutex);
96 verbosityLevel = vL;
97 posGenerator->SetVerbosity(vL);
98 angGenerator->SetVerbosity(vL);
99 eneGenerator->SetVerbosity(vL);
100}
101
104{
105 definition = aParticleDefinition;
106 charge = aParticleDefinition->GetPDGCharge();
107}
108
110{
111 if (definition == nullptr)
112 {
113 // TODO: Should this rise an exception???
114 return;
115 }
116
117 if (verbosityLevel > 1)
118 {
119 G4cout << " NumberOfParticlesToBeGenerated: "
120 << NumberOfParticlesToBeGenerated << G4endl;
121 }
122
123 part_prop_t& pp = ParticleProperties.Get();
124
125 // Position stuff
126 pp.position = posGenerator->GenerateOne();
127
128 // Create a new vertex
129 auto* vertex = new G4PrimaryVertex(pp.position,time);
130
131 for (G4int i=0; i<NumberOfParticlesToBeGenerated; ++i)
132 {
133 // Angular stuff
134 pp.momentum_direction = angGenerator->GenerateOne();
135
136 // Energy stuff
137 pp.energy = eneGenerator->GenerateOne(definition);
138
139 if (verbosityLevel >= 2)
140 {
141 G4cout << "Creating primaries and assigning to vertex" << G4endl;
142 }
143
144 // Create new primaries and set them to the vertex
145 //
146 G4double mass = definition->GetPDGMass();
147 auto* particle = new G4PrimaryParticle(definition);
148 particle->SetKineticEnergy(pp.energy );
149 particle->SetMass( mass );
150 particle->SetMomentumDirection( pp.momentum_direction );
151 particle->SetCharge( charge );
152 particle->SetPolarization(polarization.x(),
153 polarization.y(),
154 polarization.z());
155 if (verbosityLevel > 1)
156 {
157 G4cout << "Particle name: " << definition->GetParticleName() << G4endl;
158 G4cout << " Energy: " << pp.energy << G4endl;
159 G4cout << " Position: " << pp.position << G4endl;
160 G4cout << " Direction: " << pp.momentum_direction << G4endl;
161 }
162
163 // Set bweight equal to the multiple of all non-zero weights
164 //
165 G4double weight = eneGenerator->GetWeight()*biasRndm->GetBiasWeight();
166
167 if(eneGenerator->IfApplyEnergyWeight())
168 {
169 weight *= eneGenerator->GetArbEneWeight(pp.energy);
170 }
171
172 // Pass it to primary particle
173 //
174 particle->SetWeight(weight);
175 vertex->SetPrimary(particle);
176 }
177
178 // Now pass the weight to the primary vertex. CANNOT be used here!
179 // vertex->SetWeight(particle_weight);
180 evt->AddPrimaryVertex(vertex);
181
182 if (verbosityLevel > 1)
183 {
184 G4cout << " Primary Vetex generated !" << G4endl;
185 }
186}
G4ThreeVector G4ParticleMomentum
#define G4MUTEXDESTROY(mutex)
#define G4MUTEXINIT(mutex)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
value_type & Get() const
Definition G4Cache.hh:315
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition G4Event.hh:126
static G4Geantino * GeantinoDefinition()
Definition G4Geantino.cc:76
const G4String & GetParticleName() const
void SetPosDistribution(G4SPSPosDistribution *a)
G4ParticleMomentum GenerateOne()
void SetBiasRndm(G4SPSRandomGenerator *a)
G4double GetArbEneWeight(G4double)
G4double GenerateOne(G4ParticleDefinition *)
void SetBiasRndm(G4SPSRandomGenerator *a)
void SetBiasRndm(G4SPSRandomGenerator *a)
void GeneratePrimaryVertex(G4Event *evt) override
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)