Geant4_10
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 
28 //
29 // MODULE: G4SingleParticleSource.hh
30 //
31 // Version: 1.0
32 // Date: 5/02/04
33 // Author: Fan Lei
34 // Organisation: QinetiQ ltd.
35 // Customer: ESA/ESTEC
36 //
38 //
39 // CHANGE HISTORY
40 // --------------
41 //
42 // Version 1.0, 05/02/2004, Fan Lei, Created.
43 // Based on the G4GeneralParticleSource class in Geant4 v6.0
44 //
46 //
47 #include <cmath>
48 
50 
51 #include "G4SystemOfUnits.hh"
52 #include "G4PrimaryParticle.hh"
53 #include "G4Event.hh"
54 #include "Randomize.hh"
55 #include "G4ParticleTable.hh"
56 #include "G4Geantino.hh"
57 #include "G4ParticleDefinition.hh"
58 #include "G4IonTable.hh"
59 #include "G4Ions.hh"
60 #include "G4TrackingManager.hh"
61 #include "G4Track.hh"
62 
64  // Initialise all variables
65  // Position distribution Variables
66 
67  NumberOfParticlesToBeGenerated = 1;
68  particle_definition = G4Geantino::GeantinoDefinition();
69  G4ThreeVector zero;
70  particle_momentum_direction = G4ParticleMomentum(1, 0, 0);
71  particle_energy = 1.0 * MeV;
72  particle_position = zero;
73  particle_time = 0.0;
74  particle_polarization = zero;
75  particle_charge = 0.0;
76  particle_weight = 1.0;
77 
78  biasRndm = new G4SPSRandomGenerator();
79  posGenerator = new G4SPSPosDistribution();
80  posGenerator->SetBiasRndm(biasRndm);
81  angGenerator = new G4SPSAngDistribution();
82  angGenerator->SetPosDistribution(posGenerator);
83  angGenerator->SetBiasRndm(biasRndm);
84  eneGenerator = new G4SPSEneDistribution();
85  eneGenerator->SetBiasRndm(biasRndm);
86 
87  // verbosity
88  verbosityLevel = 0;
89 
90 }
91 
93  delete biasRndm;
94  delete posGenerator;
95  delete angGenerator;
96  delete eneGenerator;
97 }
98 
100  verbosityLevel = vL;
101  posGenerator->SetVerbosity(vL);
102  angGenerator->SetVerbosity(vL);
103  eneGenerator->SetVerbosity(vL);
104  G4cout << "Verbosity Set to: " << verbosityLevel << G4endl;
105 }
106 
108  G4ParticleDefinition* aParticleDefinition) {
109  particle_definition = aParticleDefinition;
110  particle_charge = particle_definition->GetPDGCharge();
111 }
112 
114  if (particle_definition == NULL)
115  return;
116 
117  if (verbosityLevel > 1)
118  G4cout << " NumberOfParticlesToBeGenerated: "
119  << NumberOfParticlesToBeGenerated << G4endl;
120 
121  // Position stuff
122  particle_position = posGenerator->GenerateOne();
123 
124  // create a new vertex
125  G4PrimaryVertex* vertex = new G4PrimaryVertex(particle_position,
126  particle_time);
127 
128  for (G4int i = 0; i < NumberOfParticlesToBeGenerated; i++) {
129  // Angular stuff
130  particle_momentum_direction = angGenerator->GenerateOne();
131  // Energy stuff
132  particle_energy = eneGenerator->GenerateOne(particle_definition);
133 
134  if (verbosityLevel >= 2)
135  G4cout << "Creating primaries and assigning to vertex" << G4endl;
136  // create new primaries and set them to the vertex
137  G4double mass = particle_definition->GetPDGMass();
138  G4PrimaryParticle* particle =
139  new G4PrimaryParticle(particle_definition);
140  particle->SetKineticEnergy( particle_energy );
141  particle->SetMass( mass );
142  particle->SetMomentumDirection( particle_momentum_direction );
143  particle->SetCharge( particle_charge );
144  particle->SetPolarization(particle_polarization.x(),
145  particle_polarization.y(),
146  particle_polarization.z());
147  if (verbosityLevel > 1) {
148  G4cout << "Particle name: "
149  << particle_definition->GetParticleName() << G4endl;
150  G4cout << " Energy: " << particle_energy << G4endl;
151  G4cout << " Position: " << particle_position << G4endl;
152  G4cout << " Direction: " << particle_momentum_direction
153  << G4endl;
154  }
155  // Set bweight equal to the multiple of all non-zero weights
156  particle_weight = eneGenerator->GetWeight()*biasRndm->GetBiasWeight();
157  // pass it to primary particle
158  particle->SetWeight(particle_weight);
159 
160  vertex->SetPrimary(particle);
161 
162  }
163  // now pass the weight to the primary vertex. CANNOT be used here!
164  // vertex->SetWeight(particle_weight);
165  evt->AddPrimaryVertex(vertex);
166  if (verbosityLevel > 1)
167  G4cout << " Primary Vetex generated !" << G4endl;
168 }
169 
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
void SetBiasRndm(G4SPSRandomGenerator *a)
double x() const
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:143
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
double z() const
void SetKineticEnergy(G4double eKin)
G4GLOB_DLL std::ostream G4cout
G4ParticleMomentum GenerateOne()
void SetBiasRndm(G4SPSRandomGenerator *a)
void SetMass(G4double mas)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GenerateOne(G4ParticleDefinition *)
void SetBiasRndm(G4SPSRandomGenerator *a)
G4double GetPDGMass() const
void SetPosDistribution(G4SPSPosDistribution *a)
void SetMomentumDirection(const G4ThreeVector &p)
void SetWeight(G4double w)
void SetCharge(G4double chg)
double y() const
void SetPrimary(G4PrimaryParticle *pp)
#define G4endl
Definition: G4ios.hh:61
void GeneratePrimaryVertex(G4Event *evt)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4ThreeVector G4ParticleMomentum
void SetPolarization(const G4ThreeVector &pol)