Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eBremsstrahlungSpectrum Class Reference

#include <G4eBremsstrahlungSpectrum.hh>

Inheritance diagram for G4eBremsstrahlungSpectrum:
Collaboration diagram for G4eBremsstrahlungSpectrum:

Public Member Functions

 G4eBremsstrahlungSpectrum (const G4DataVector &bins, const G4String &name)
 
 ~G4eBremsstrahlungSpectrum ()
 
G4double Probability (G4int Z, G4double tMin, G4double tMax, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const
 
G4double AverageEnergy (G4int Z, G4double tMin, G4double tMax, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const
 
G4double SampleEnergy (G4int Z, G4double tMin, G4double tMax, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const
 
G4double MaxEnergyOfSecondaries (G4double kineticEnergy, G4int Z=0, const G4ParticleDefinition *pd=0) const
 
G4double Excitation (G4int Z, G4double kineticEnergy) const
 
void PrintData () const
 
- Public Member Functions inherited from G4VEnergySpectrum
 G4VEnergySpectrum ()
 
virtual ~G4VEnergySpectrum ()
 

Detailed Description

Definition at line 65 of file G4eBremsstrahlungSpectrum.hh.

Constructor & Destructor Documentation

G4eBremsstrahlungSpectrum::G4eBremsstrahlungSpectrum ( const G4DataVector bins,
const G4String name 
)

Definition at line 55 of file G4eBremsstrahlungSpectrum.cc.

57  lowestE(0.1*eV),
58  xp(bins)
59 {
60  length = xp.size();
61  theBRparam = new G4BremsstrahlungParameters(name,length+1);
62 
63  verbose = 0;
64 }
static constexpr double eV
Definition: G4SIunits.hh:215
G4eBremsstrahlungSpectrum::~G4eBremsstrahlungSpectrum ( )

Definition at line 67 of file G4eBremsstrahlungSpectrum.cc.

68 {
69  delete theBRparam;
70 }

Member Function Documentation

G4double G4eBremsstrahlungSpectrum::AverageEnergy ( G4int  Z,
G4double  tMin,
G4double  tMax,
G4double  kineticEnergy,
G4int  shell = 0,
const G4ParticleDefinition pd = 0 
) const
virtual

Implements G4VEnergySpectrum.

Definition at line 120 of file G4eBremsstrahlungSpectrum.cc.

126 {
127  G4double tm = std::min(tmax, e);
128  G4double t0 = std::max(tmin, lowestE);
129  if(t0 >= tm) return 0.0;
130 
131  t0 /= e;
132  tm /= e;
133 
134  G4double z0 = lowestE/e;
135 
136  G4DataVector p;
137 
138  // Access parameters
139  for (size_t i=0; i<=length; i++) {
140  p.push_back(theBRparam->Parameter(i, Z, e));
141  }
142 
143  G4double x = AverageValue(t0, tm, p);
144  G4double y = IntSpectrum(z0, 1.0, p);
145 
146  // Add integrant over lowest energies
147  G4double zmin = tmin/e;
148  if(zmin < t0) {
149  G4double c = std::sqrt(theBRparam->ParameterC(Z));
150  x += p[0]*(t0 - zmin - c*(std::atan(t0/c) - std::atan(zmin/c)));
151  }
152  x *= e;
153 
154  if(1 < verbose) {
155  G4cout << "tcut(MeV)= " << tmin/MeV
156  << "; tMax(MeV)= " << tmax/MeV
157  << "; e(MeV)= " << e/MeV
158  << "; t0= " << t0
159  << "; tm= " << tm
160  << "; y= " << y
161  << "; x= " << x
162  << G4endl;
163  }
164  p.clear();
165 
166  if(y > 0.0) x /= y;
167  else x = 0.0;
168  // if(x < 0.0) x = 0.0;
169 
170  return x;
171 }
const char * p
Definition: xmltok.h:285
G4GLOB_DLL std::ostream G4cout
G4double Parameter(G4int parameterIndex, G4int Z, G4double energy) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
G4double ParameterC(G4int index) const

Here is the call graph for this function:

G4double G4eBremsstrahlungSpectrum::Excitation ( G4int  Z,
G4double  kineticEnergy 
) const
virtual

Implements G4VEnergySpectrum.

Definition at line 302 of file G4eBremsstrahlungSpectrum.cc.

303 {
304  return 0.0;
305 }
G4double G4eBremsstrahlungSpectrum::MaxEnergyOfSecondaries ( G4double  kineticEnergy,
G4int  Z = 0,
const G4ParticleDefinition pd = 0 
) const
virtual

Implements G4VEnergySpectrum.

Definition at line 307 of file G4eBremsstrahlungSpectrum.cc.

310 {
311  return kineticEnergy;
312 }
void G4eBremsstrahlungSpectrum::PrintData ( void  ) const
virtual

Implements G4VEnergySpectrum.

Definition at line 299 of file G4eBremsstrahlungSpectrum.cc.

300 { theBRparam->PrintData(); }

Here is the call graph for this function:

G4double G4eBremsstrahlungSpectrum::Probability ( G4int  Z,
G4double  tMin,
G4double  tMax,
G4double  kineticEnergy,
G4int  shell = 0,
const G4ParticleDefinition pd = 0 
) const
virtual

Implements G4VEnergySpectrum.

Definition at line 73 of file G4eBremsstrahlungSpectrum.cc.

79 {
80  G4double tm = std::min(tmax, e);
81  G4double t0 = std::max(tmin, lowestE);
82  if(t0 >= tm) return 0.0;
83 
84  t0 /= e;
85  tm /= e;
86 
87  G4double z0 = lowestE/e;
89 
90  // Access parameters
91  for (size_t i=0; i<=length; i++) {
92  p.push_back(theBRparam->Parameter(i, Z, e));
93  }
94 
95  G4double x = IntSpectrum(t0, tm, p);
96  G4double y = IntSpectrum(z0, 1.0, p);
97 
98 
99  if(1 < verbose) {
100  G4cout << "tcut(MeV)= " << tmin/MeV
101  << "; tMax(MeV)= " << tmax/MeV
102  << "; t0= " << t0
103  << "; tm= " << tm
104  << "; xp[0]= " << xp[0]
105  << "; z= " << z0
106  << "; val= " << x
107  << "; nor= " << y
108  << G4endl;
109  }
110  p.clear();
111 
112  if(y > 0.0) x /= y;
113  else x = 0.0;
114  // if(x < 0.0) x = 0.0;
115 
116  return x;
117 }
const char * p
Definition: xmltok.h:285
G4GLOB_DLL std::ostream G4cout
G4double Parameter(G4int parameterIndex, G4int Z, G4double energy) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4eBremsstrahlungSpectrum::SampleEnergy ( G4int  Z,
G4double  tMin,
G4double  tMax,
G4double  kineticEnergy,
G4int  shell = 0,
const G4ParticleDefinition pd = 0 
) const
virtual

Implements G4VEnergySpectrum.

Definition at line 174 of file G4eBremsstrahlungSpectrum.cc.

180 {
181  G4double tm = std::min(tmax, e);
182  G4double t0 = std::max(tmin, lowestE);
183  if(t0 >= tm) return 0.0;
184 
185  t0 /= e;
186  tm /= e;
187 
188  G4DataVector p;
189 
190  for (size_t i=0; i<=length; i++) {
191  p.push_back(theBRparam->Parameter(i, Z, e));
192  }
193  G4double amaj = std::max(p[length], 1. - (p[1] - p[0])*xp[0]/(xp[1] - xp[0]) );
194 
195  G4double amax = std::log(tm);
196  G4double amin = std::log(t0);
197  G4double tgam, q, fun;
198 
199  do {
200  G4double x = amin + G4UniformRand()*(amax - amin);
201  tgam = G4Exp(x);
202  fun = Function(tgam, p);
203 
204  if(fun > amaj) {
205  G4cout << "WARNING in G4eBremsstrahlungSpectrum::SampleEnergy:"
206  << " Majoranta " << amaj
207  << " < " << fun
208  << G4endl;
209  }
210 
211  q = amaj * G4UniformRand();
212  } while (q > fun);
213 
214  tgam *= e;
215 
216  p.clear();
217 
218  return tgam;
219 }
const char * p
Definition: xmltok.h:285
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double Parameter(G4int parameterIndex, G4int Z, G4double energy) const
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


The documentation for this class was generated from the following files: