Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BrachyPrimaryGeneratorActionI.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 // --------------------------------------------------------------
28 // GEANT 4 - Brachytherapy example
29 // --------------------------------------------------------------
30 //
31 // Code developed by:
32 // S.Guatelli
33 //
34 // ********************************************
35 // * *
36 // * BrachyPrimaryGeneratorActionI.cc *
37 // * *
38 // ********************************************
39 //
40 // $Id$
41 //
43 #include "globals.hh"
44 #include "Randomize.hh"
45 #include "G4SystemOfUnits.hh"
46 #include "G4ParticleTable.hh"
47 #include "G4Event.hh"
48 #include "G4ParticleGun.hh"
49 #include "G4UImanager.hh"
50 #include "G4RunManager.hh"
51 
53 {
54  G4int numberParticles = 1;
55 
56  particleGun = new G4ParticleGun(numberParticles);
57 
58  // Gamma energy spectrum ...
59  // Fill a vector with the energy probabilities
60  energySpectrum.push_back(0.783913);
61  energySpectrum.push_back(0.170416);
62  energySpectrum.push_back(0.045671);
63 }
64 
66 {
67  if(particleGun)
68  delete particleGun;
69 }
70 
72 {
73  // Define the primary particle type
75  G4String ParticleName = "gamma";
76  G4ParticleDefinition* particle = particleTable -> FindParticle(ParticleName);
77 
78  particleGun -> SetParticleDefinition(particle);
79 
80  // Random generation of gamma source point inside the Iodium core ...
81  G4double x,y,z;
82  G4double radiuMax = 0.30*mm;
83  G4double radiusMin = 0.085*mm;
84 
85  do{
86  x = (G4UniformRand()-0.5)*(radiuMax)/0.5;
87  y = (G4UniformRand()-0.5)*(radiuMax)/0.5;
88  }while(((x*x+y*y )> (radiuMax*radiuMax))||((x*x+y*y)<(radiusMin*radiusMin)));
89 
90  z = (G4UniformRand()-0.5)*1.75*mm/0.5 ;
91 
92  G4ThreeVector position(x,y,z);
93  particleGun -> SetParticlePosition(position);
94 
95  // Random generation of the impulse direction of primary particles ...
96  G4double a,b,c;
97  G4double n;
98  do{
99  a = (G4UniformRand()-0.5)/0.5;
100  b = (G4UniformRand()-0.5)/0.5;
101  c = (G4UniformRand()-0.5)/0.5;
102  n = a*a+b*b+c*c;
103  }while(n > 1 || n == 0.0);
104  n = std::sqrt(n);
105  a /= n;
106  b /= n;
107  c /= n;
108 
109  G4ThreeVector direction(a,b,c);
110  particleGun -> SetParticleMomentumDirection(direction);
111 
112  // Generate the primary particles with a defined energy spectrum
113  G4double random = G4UniformRand();
114  G4double sum = 0;
115  G4int i = 0;
116  while(sum<random){sum+=energySpectrum[i];
117  i++;}
118 
119  // energy spectrum
120  if(i==1){primaryParticleEnergy = 27.4*keV;}
121  else{
122  if(i==2){primaryParticleEnergy = 31.4*keV;}
123  else {primaryParticleEnergy = 35.5*keV;}}
124 
125  particleGun -> SetParticleEnergy(primaryParticleEnergy);
126 
127  // generate primary particle
128  particleGun->GeneratePrimaryVertex(anEvent);
129 }
130 
132 {
133  return primaryParticleEnergy;
134 }
135