Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4eBremsstrahlungSpectrum Class Reference

#include <G4eBremsstrahlungSpectrum.hh>

Inheritance diagram for G4eBremsstrahlungSpectrum:
G4VEnergySpectrum

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 54 of file G4eBremsstrahlungSpectrum.cc.

56  lowestE(0.1*eV),
57  xp(bins)
58 {
59  length = xp.size();
60  theBRparam = new G4BremsstrahlungParameters(name,length+1);
61 
62  verbose = 0;
63 }
G4eBremsstrahlungSpectrum::~G4eBremsstrahlungSpectrum ( )

Definition at line 66 of file G4eBremsstrahlungSpectrum.cc.

67 {
68  delete theBRparam;
69 }

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 119 of file G4eBremsstrahlungSpectrum.cc.

References test::c, G4cout, G4endl, G4INCL::Math::max(), python.hepunit::MeV, G4INCL::Math::min(), G4BremsstrahlungParameters::Parameter(), G4BremsstrahlungParameters::ParameterC(), G4InuclParticleNames::tm, test::x, and G4InuclParticleNames::z0.

125 {
126  G4double tm = std::min(tmax, e);
127  G4double t0 = std::max(tmin, lowestE);
128  if(t0 >= tm) return 0.0;
129 
130  t0 /= e;
131  tm /= e;
132 
133  G4double z0 = lowestE/e;
134 
135  G4DataVector p;
136 
137  // Access parameters
138  for (size_t i=0; i<=length; i++) {
139  p.push_back(theBRparam->Parameter(i, Z, e));
140  }
141 
142  G4double x = AverageValue(t0, tm, p);
143  G4double y = IntSpectrum(z0, 1.0, p);
144 
145  // Add integrant over lowest energies
146  G4double zmin = tmin/e;
147  if(zmin < t0) {
148  G4double c = std::sqrt(theBRparam->ParameterC(Z));
149  x += p[0]*(t0 - zmin - c*(std::atan(t0/c) - std::atan(zmin/c)));
150  }
151  x *= e;
152 
153  if(1 < verbose) {
154  G4cout << "tcut(MeV)= " << tmin/MeV
155  << "; tMax(MeV)= " << tmax/MeV
156  << "; e(MeV)= " << e/MeV
157  << "; t0= " << t0
158  << "; tm= " << tm
159  << "; y= " << y
160  << "; x= " << x
161  << G4endl;
162  }
163  p.clear();
164 
165  if(y > 0.0) x /= y;
166  else x = 0.0;
167  // if(x < 0.0) x = 0.0;
168 
169  return x;
170 }
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
double G4double
Definition: G4Types.hh:76
G4double ParameterC(G4int index) const
G4double G4eBremsstrahlungSpectrum::Excitation ( G4int  Z,
G4double  kineticEnergy 
) const
virtual

Implements G4VEnergySpectrum.

Definition at line 301 of file G4eBremsstrahlungSpectrum.cc.

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

Implements G4VEnergySpectrum.

Definition at line 306 of file G4eBremsstrahlungSpectrum.cc.

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

Implements G4VEnergySpectrum.

Definition at line 298 of file G4eBremsstrahlungSpectrum.cc.

References G4BremsstrahlungParameters::PrintData().

299 { theBRparam->PrintData(); }
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 72 of file G4eBremsstrahlungSpectrum.cc.

References G4cout, G4endl, G4INCL::Math::max(), python.hepunit::MeV, G4INCL::Math::min(), G4BremsstrahlungParameters::Parameter(), G4InuclParticleNames::tm, test::x, and G4InuclParticleNames::z0.

78 {
79  G4double tm = std::min(tmax, e);
80  G4double t0 = std::max(tmin, lowestE);
81  if(t0 >= tm) return 0.0;
82 
83  t0 /= e;
84  tm /= e;
85 
86  G4double z0 = lowestE/e;
88 
89  // Access parameters
90  for (size_t i=0; i<=length; i++) {
91  p.push_back(theBRparam->Parameter(i, Z, e));
92  }
93 
94  G4double x = IntSpectrum(t0, tm, p);
95  G4double y = IntSpectrum(z0, 1.0, p);
96 
97 
98  if(1 < verbose) {
99  G4cout << "tcut(MeV)= " << tmin/MeV
100  << "; tMax(MeV)= " << tmax/MeV
101  << "; t0= " << t0
102  << "; tm= " << tm
103  << "; xp[0]= " << xp[0]
104  << "; z= " << z0
105  << "; val= " << x
106  << "; nor= " << y
107  << G4endl;
108  }
109  p.clear();
110 
111  if(y > 0.0) x /= y;
112  else x = 0.0;
113  // if(x < 0.0) x = 0.0;
114 
115  return x;
116 }
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
double G4double
Definition: G4Types.hh:76
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 173 of file G4eBremsstrahlungSpectrum.cc.

References G4cout, G4endl, G4UniformRand, G4INCL::Math::max(), G4INCL::Math::min(), G4BremsstrahlungParameters::Parameter(), G4InuclParticleNames::tm, and test::x.

179 {
180  G4double tm = std::min(tmax, e);
181  G4double t0 = std::max(tmin, lowestE);
182  if(t0 >= tm) return 0.0;
183 
184  t0 /= e;
185  tm /= e;
186 
187  G4DataVector p;
188 
189  for (size_t i=0; i<=length; i++) {
190  p.push_back(theBRparam->Parameter(i, Z, e));
191  }
192  G4double amaj = std::max(p[length], 1. - (p[1] - p[0])*xp[0]/(xp[1] - xp[0]) );
193 
194  G4double amax = std::log(tm);
195  G4double amin = std::log(t0);
196  G4double tgam, q, fun;
197 
198  do {
199  G4double x = amin + G4UniformRand()*(amax - amin);
200  tgam = std::exp(x);
201  fun = Function(tgam, p);
202 
203  if(fun > amaj) {
204  G4cout << "WARNING in G4eBremsstrahlungSpectrum::SampleEnergy:"
205  << " Majoranta " << amaj
206  << " < " << fun
207  << G4endl;
208  }
209 
210  q = amaj * G4UniformRand();
211  } while (q > fun);
212 
213  tgam *= e;
214 
215  p.clear();
216 
217  return tgam;
218 }
const char * p
Definition: xmltok.h:285
#define G4UniformRand()
Definition: Randomize.hh:87
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
double G4double
Definition: G4Types.hh:76

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