Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleHPFissionSpectrum.hh
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 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
28 //
29 #ifndef G4ParticleHPFissionSpectrum_h
30 #define G4ParticleHPFissionSpectrum_h 1
31 
32 #include <fstream>
34 
35 #include "globals.hh"
36 #include "G4ios.hh"
37 #include "Randomize.hh"
38 #include "G4Exp.hh"
39 #include "G4ParticleHPVector.hh"
40 #include "G4VParticleHPEDis.hh"
41 
42 // we will need a List of these .... one per term.
43 
45 {
46  public:
48  {
49  expm1 = G4Exp(-1.);
50  }
52  {
53  }
54 
55  inline void Init(std::istream & aDataFile)
56  {
57  theFractionalProb.Init(aDataFile, CLHEP::eV);
58  theThetaDist.Init(aDataFile, CLHEP::eV);
59  }
60 
62  {
63  return theFractionalProb.GetY(anEnergy);
64  }
65 
66  inline G4double Sample(G4double anEnergy)
67  {
68  G4double theta = theThetaDist.GetY(anEnergy);
69  // here we need to sample Maxwells distribution, if
70  // need be.
71  G4double result, cut;
72  G4double range =50*CLHEP::MeV;
73  G4double max = Maxwell((theta*CLHEP::eV)/2., theta);
75  G4int icounter=0;
76  G4int icounter_max=1024;
77  do
78  {
79  icounter++;
80  if ( icounter > icounter_max ) {
81  G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
82  break;
83  }
84  result = range*G4UniformRand();
85  value = Maxwell(result, theta);
86  cut = G4UniformRand();
87  }
88  while(cut > value/max); // Loop checking, 11.05.2015, T. Koi
89  return result;
90  }
91 
92  private:
93 
94  // this is the function to sample from.
95  inline G4double Maxwell(G4double anEnergy, G4double theta)
96  {
97  G4double result = std::sqrt(anEnergy/CLHEP::eV)*G4Exp(-anEnergy/CLHEP::eV/theta);
98  return result;
99  }
100 
101  private:
102 
103  G4double expm1;
104 
105  G4ParticleHPVector theFractionalProb;
106 
107  G4ParticleHPVector theThetaDist;
108 
109 };
110 
111 #endif
G4double G4ParticleHPJENDLHEData::G4double result
int G4int
Definition: G4Types.hh:78
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static constexpr double MeV
G4double GetY(G4double x)
static constexpr double eV
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetFractionalProbability(G4double anEnergy)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define G4endl
Definition: G4ios.hh:61
void Init(std::istream &aDataFile)
G4double Sample(G4double anEnergy)
double G4double
Definition: G4Types.hh:76