Geant4  10.02.p01
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 #include "G4AutoLock.hh"
63 
65  //definition = G4Geantino::GeantinoDefinition();
67  energy = 1.*MeV;
69 }
70 
72 // // Initialise all variables
73 // // Position distribution Variables
74 //
77 // G4ThreeVector zero;
78 // particle_momentum_direction = G4ParticleMomentum(1, 0, 0);
79 // particle_energy = 1.0 * MeV;
80 // particle_position = zero;
81 // particle_time = 0.0;
82 // particle_polarization = zero;
83  charge = 0.0;
84  time = 0;
86 
95 
96  // verbosity
97  verbosityLevel = 0;
98 
100 
101 }
102 
104  delete biasRndm;
105  delete posGenerator;
106  delete angGenerator;
107  delete eneGenerator;
109 }
110 
112  G4AutoLock l(&mutex);
113  verbosityLevel = vL;
117  //G4cout << "Verbosity Set to: " << verbosityLevel << G4endl;
118 }
119 
121  G4ParticleDefinition* aParticleDefinition) {
122  definition = aParticleDefinition;
123  charge = aParticleDefinition->GetPDGCharge();
124 }
125 
127 
128  //G4AutoLock l(&mutex);
129  //part_prop_t& pp = ParticleProperties.Get();
130  if (definition == NULL) {
131  //TODO: Should this rise an exception???
132  return ;
133  }
134  //return;
135 
136  if (verbosityLevel > 1)
137  G4cout << " NumberOfParticlesToBeGenerated: "
139 
141  // Position stuff
143 
144  // create a new vertex
145  G4PrimaryVertex* vertex = new G4PrimaryVertex(pp.position,time);
146 
147  for (G4int i = 0; i < NumberOfParticlesToBeGenerated; i++) {
148  // Angular stuff
150  // Energy stuff
152 
153  if (verbosityLevel >= 2)
154  G4cout << "Creating primaries and assigning to vertex" << G4endl;
155  // create new primaries and set them to the vertex
156  G4double mass = definition->GetPDGMass();
157  G4PrimaryParticle* particle =
159  particle->SetKineticEnergy(pp.energy );
160  particle->SetMass( mass );
162  particle->SetCharge( charge );
163  particle->SetPolarization(polarization.x(),
164  polarization.y(),
165  polarization.z());
166  if (verbosityLevel > 1) {
167  G4cout << "Particle name: "
169  G4cout << " Energy: " << pp.energy << G4endl;
170  G4cout << " Position: " << pp.position << G4endl;
171  G4cout << " Direction: " << pp.momentum_direction
172  << G4endl;
173  }
174  // Set bweight equal to the multiple of all non-zero weights
176  // pass it to primary particle
177  particle->SetWeight(weight);
178 
179  vertex->SetPrimary(particle);
180 
181  }
182  // now pass the weight to the primary vertex. CANNOT be used here!
183  // vertex->SetWeight(particle_weight);
184  evt->AddPrimaryVertex(vertex);
185  if (verbosityLevel > 1)
186  G4cout << " Primary Vetex generated !" << G4endl;
187 }
188 
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
G4Cache< part_prop_t > ParticleProperties
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * definition
void SetBiasRndm(G4SPSRandomGenerator *a)
CLHEP::Hep3Vector G4ThreeVector
Andrea Dotti Feb 2015 Important: This is a shared class between threads.
Andrea Dotti Feb 2015 Important: This is a shared class between threads.
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:154
#define G4MUTEXINIT(mutex)
Definition: G4Threading.hh:177
G4SPSEneDistribution * eneGenerator
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
Andrea Dotti Feb 2015 Important: This is a shared class between threads.
void SetKineticEnergy(G4double eKin)
G4SPSPosDistribution * posGenerator
G4GLOB_DLL std::ostream G4cout
G4ParticleMomentum GenerateOne()
void SetBiasRndm(G4SPSRandomGenerator *a)
G4SPSRandomGenerator * biasRndm
void SetMass(G4double mas)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GenerateOne(G4ParticleDefinition *)
void SetBiasRndm(G4SPSRandomGenerator *a)
Andrea Dotti Feb 2015 Important: This is a shared class between threads.
G4double GetPDGMass() const
void SetPosDistribution(G4SPSPosDistribution *a)
void SetMomentumDirection(const G4ThreeVector &p)
void SetWeight(G4double w)
void SetCharge(G4double chg)
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
#define G4MUTEXDESTROY(mutex)
Definition: G4Threading.hh:178
void SetPolarization(const G4ThreeVector &pol)
G4SPSAngDistribution * angGenerator