source: trunk/examples/advanced/radioprotection/src/RemSimPrimaryGeneratorAction.cc@ 1305

Last change on this file since 1305 was 807, checked in by garnier, 17 years ago

update

File size: 6.7 KB
Line 
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// $Id: RemSimPrimaryGeneratorAction.cc,v 1.17 2006/06/29 16:24:09 gunter Exp $// Author: Susanna Guatelli, guatelli@ge.infn.it
28
29#include "RemSimPrimaryGeneratorAction.hh"
30#ifdef G4ANALYSIS_USE
31#include "RemSimAnalysisManager.hh"
32#endif
33#include "G4Event.hh"
34#include "G4ParticleGun.hh"
35#include "G4ParticleTable.hh"
36#include "G4ParticleDefinition.hh"
37#include "Randomize.hh"
38#include <fstream>
39#include "G4DataVector.hh"
40#include "RemSimPrimaryGeneratorMessenger.hh"
41#include <math.h>
42
43RemSimPrimaryGeneratorAction::RemSimPrimaryGeneratorAction()
44{
45 readFile = false;
46
47 particleGun = new G4ParticleGun(1);
48 messenger = new RemSimPrimaryGeneratorMessenger(this);
49
50 energies = new G4DataVector;
51 data = new G4DataVector;
52
53 size = 0;
54
55 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
56 G4String ParticleName = "proton";
57 G4ParticleDefinition* particle = particleTable->FindParticle(ParticleName);
58
59 particleGun->SetParticleDefinition(particle);
60
61 // Default values of the primary generator
62 G4ThreeVector position(0,0,-25. *m);
63 particleGun -> SetParticlePosition(position);
64 G4ThreeVector direction(0., 0., 1.);
65 particleGun -> SetParticleMomentumDirection(direction);
66 particleGun -> SetParticleEnergy(250. * MeV);
67}
68
69RemSimPrimaryGeneratorAction::~RemSimPrimaryGeneratorAction()
70{
71 delete data;
72 delete energies;
73 delete [] cumulate;
74 delete [] energy;
75 delete messenger;
76 delete particleGun;
77}
78
79G4double RemSimPrimaryGeneratorAction::GetInitialEnergy()
80{
81 G4double primaryParticleEnergy = particleGun -> GetParticleEnergy();
82 return primaryParticleEnergy;
83}
84
85void RemSimPrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
86{
87 G4int baryon = particleGun -> GetParticleDefinition() -> GetBaryonNumber();
88
89 // Generate the energy spectrum of primary particles according to
90 // the flux reported in the .txt file (firt column = Energy (MeV/nucl),
91 // second column = differnetila flux). The fluxes are derived from CREME96
92 if (readFile == true)
93 {
94 // Uniform number between 0. and 1.
95 G4double random = G4UniformRand();
96
97 G4int nbelow = 0; // largest k such that I[k] is known to be <= rand
98 G4int nabove = size; // largest k such that I[k] is known to be > rand
99 G4int middle;
100
101 while (nabove > nbelow+1) {
102 middle = (nabove + nbelow+1)>>1;
103 if (random >= cumulate[middle]) {
104 nbelow = middle;
105 } else {
106 nabove = middle;
107 }
108 }
109
110 // Linear interpolation
111
112 G4double energy_gun = 0. * MeV;
113
114 G4double binMeasure = cumulate[nabove] - cumulate[nbelow];
115
116 if ( binMeasure == 0 ) {
117
118 energy_gun = (energy[nbelow] + random*(energy[nabove] - energy[nbelow]))* baryon;
119 G4cout<<"BinMeasure is zero !!!!!"<<G4endl;
120 }
121
122 else
123 {
124 G4double binFraction = (random - cumulate[nbelow]) / binMeasure;
125 energy_gun = energy[nbelow] + binFraction *(energy[nabove] - energy[nbelow]);
126 }
127 particleGun -> SetParticleEnergy(energy_gun * baryon);
128 }
129
130#ifdef G4ANALYSIS_USE
131 G4double energy_plot = (particleGun -> GetParticleEnergy())/baryon;
132 // plot of MeV/nucl
133 RemSimAnalysisManager* analysis = RemSimAnalysisManager::getInstance();
134 // Plot the energy spectrum of primary particles
135 analysis -> primaryParticleEnergyDistribution(energy_plot/MeV);
136#endif
137
138particleGun -> GeneratePrimaryVertex(anEvent);
139}
140
141void RemSimPrimaryGeneratorAction::ReadProbability(G4String fileName)
142{
143 // This method allows to read the .txt files containing the
144 // fluxes of the galactic cosmic rays derived form CREME96
145 readFile = true;
146
147 std::ifstream file(fileName, std::ios::in);
148 std::filebuf* lsdp = file.rdbuf();
149
150 if (! (lsdp->is_open()) )
151 {
152 G4String excep = "RemSimPrimaryGenerator - data file: "
153 + fileName + " not found";
154 G4Exception(excep);
155 }
156
157 G4double a = 0;
158 G4int index = 0;
159 G4int k = 1;
160
161 do
162 {
163 file >> a;
164 G4int nColumns = 2;
165 // The file is organized into two columns:
166 // column contains the probability distribution
167 // The file terminates with the pattern: -1
168 // -2
169 if (a == -1 || a == -2)
170 {
171
172 }
173 else
174 {
175 if (k%nColumns != 0)
176 {
177 G4double energy = a *MeV;
178 energies -> push_back(energy);
179 // G4cout << "energy: " << energy << G4endl;
180 k++;
181 }
182 else if (k%nColumns == 0)
183 {
184 G4double data_value = a;
185 data -> push_back(data_value);
186 //G4cout<< "probability: "<< data_value << G4endl;
187 k=1;
188 index++;
189 }
190 }
191 } while (a != -2); // end of file
192
193 file.close();
194 size = index -1 ;
195 cumulate = new G4double[size+1];
196 energy = new G4double[size+1];
197 cumulate[0] = 0;
198 energy [0] = (*energies)[0];
199
200 for (G4int ptn = 0; ptn <size; ptn++ ) {
201
202 G4double yy = ((*data)[ptn+1] + (*data)[ptn])/2.;
203 G4double weight = yy *((*energies)[ptn +1] - (*energies)[ptn]);
204 cumulate[ptn +1] = cumulate[ptn] + weight;
205 energy[ptn+1] = (*energies)[ptn+1];
206 }
207
208 // Normalise the cumulate to 1 (that corresponds to the last value)
209 for ( G4int ptn= 0; ptn < size + 1; ptn++ )
210{
211 cumulate[ptn] /= cumulate[size];
212}
213}
214
215
Note: See TracBrowser for help on using the repository browser.