Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RemSimPrimaryGeneratorAction.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 // $Id: RemSimPrimaryGeneratorAction.cc,v 1.17 2006-06-29 16:24:09 gunter Exp $// Author: Susanna Guatelli, guatelli@ge.infn.it
28 
29 #include <fstream>
30 #include <cmath>
31 
34 
35 #include "G4SystemOfUnits.hh"
36 #include "G4Event.hh"
37 #include "G4ParticleGun.hh"
38 #include "G4ParticleTable.hh"
39 #include "G4ParticleDefinition.hh"
40 #include "Randomize.hh"
41 #include "G4DataVector.hh"
42 
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 
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 
70 {
71  delete data;
72  delete energies;
73  delete [] cumulate;
74  delete [] energy;
75  delete messenger;
76  delete particleGun;
77 }
78 
80 {
81  G4double primaryParticleEnergy = particleGun -> GetParticleEnergy();
82  return primaryParticleEnergy;
83 }
84 
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 
131 particleGun -> GeneratePrimaryVertex(anEvent);
132 }
133 
135 {
136  // This method allows to read the .txt files containing the
137  // fluxes of the galactic cosmic rays derived form CREME96
138  readFile = true;
139 
140  std::ifstream file(fileName, std::ios::in);
141  std::filebuf* lsdp = file.rdbuf();
142 
143  if (! (lsdp->is_open()) )
144  {
145  G4String excep = "RemSimPrimaryGenerator - data file: "
146  + fileName + " not found";
147  G4Exception("RemSimPrimaryGeneratorAction::ReadProbability()",
148  "RadioP001",FatalException,excep);
149  }
150 
151  G4double a = 0;
152  G4int index = 0;
153  G4int k = 1;
154 
155  do
156  {
157  file >> a;
158  G4int nColumns = 2;
159  // The file is organized into two columns:
160  // column contains the probability distribution
161  // The file terminates with the pattern: -1
162  // -2
163  if (a == -1 || a == -2)
164  {
165 
166  }
167  else
168  {
169  if (k%nColumns != 0)
170  {
171  G4double ene = a *MeV;
172  energies -> push_back(ene);
173  // G4cout << "energy: " << ene << G4endl;
174  k++;
175  }
176  else if (k%nColumns == 0)
177  {
178  G4double data_value = a;
179  data -> push_back(data_value);
180  //G4cout<< "probability: "<< data_value << G4endl;
181  k=1;
182  index++;
183  }
184  }
185  } while (a != -2); // end of file
186 
187  file.close();
188  size = index -1 ;
189  cumulate = new G4double[size+1];
190  energy = new G4double[size+1];
191  cumulate[0] = 0;
192  energy [0] = (*energies)[0];
193 
194  for (G4int ptn = 0; ptn <size; ptn++ ) {
195 
196  G4double yy = ((*data)[ptn+1] + (*data)[ptn])/2.;
197  G4double weight = yy *((*energies)[ptn +1] - (*energies)[ptn]);
198  cumulate[ptn +1] = cumulate[ptn] + weight;
199  energy[ptn+1] = (*energies)[ptn+1];
200  }
201 
202  // Normalise the cumulate to 1 (that corresponds to the last value)
203  for ( G4int ptn= 0; ptn < size + 1; ptn++ )
204 {
205  cumulate[ptn] /= cumulate[size];
206 }
207 }
208 
209