#include <G4eBremsstrahlungSpectrum.hh>
|
| 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 |
|
| G4VEnergySpectrum () |
|
virtual | ~G4VEnergySpectrum () |
|
virtual G4double | Probability (G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0 |
|
virtual G4double | AverageEnergy (G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0 |
|
virtual G4double | SampleEnergy (G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0 |
|
virtual G4double | MaxEnergyOfSecondaries (G4double kineticEnergy, G4int Z=0, const G4ParticleDefinition *pd=0) const =0 |
|
virtual G4double | Excitation (G4int Z, G4double kineticEnergy) const =0 |
|
virtual void | PrintData () const =0 |
|
◆ G4eBremsstrahlungSpectrum()
◆ ~G4eBremsstrahlungSpectrum()
G4eBremsstrahlungSpectrum::~G4eBremsstrahlungSpectrum |
( |
| ) |
|
◆ AverageEnergy()
Implements G4VEnergySpectrum.
Definition at line 119 of file G4eBremsstrahlungSpectrum.cc.
125{
127 G4double t0 = std::max(tmin, lowestE);
128 if(t0 >= tm) return 0.0;
129
130 t0 /= e;
131 tm /= e;
132
134
136
137
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
147 if(zmin < t0) {
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
162 }
163 p.clear();
164
165 if(y > 0.0) x /= y;
166 else x = 0.0;
167
168
169 return x;
170}
G4DLLIMPORT std::ostream G4cout
G4double Parameter(G4int parameterIndex, G4int Z, G4double energy) const
G4double ParameterC(G4int index) const
◆ Excitation()
◆ MaxEnergyOfSecondaries()
◆ PrintData()
void G4eBremsstrahlungSpectrum::PrintData |
( |
| ) |
const |
|
virtual |
◆ Probability()
Implements G4VEnergySpectrum.
Definition at line 72 of file G4eBremsstrahlungSpectrum.cc.
78{
80 G4double t0 = std::max(tmin, lowestE);
81 if(t0 >= tm) return 0.0;
82
83 t0 /= e;
84 tm /= e;
85
88
89
90 for (size_t i=0; i<=length; i++) {
91 p.push_back(theBRparam->
Parameter(i, Z, e));
92 }
93
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
108 }
109 p.clear();
110
111 if(y > 0.0) x /= y;
112 else x = 0.0;
113
114
115 return x;
116}
◆ SampleEnergy()
Implements G4VEnergySpectrum.
Definition at line 173 of file G4eBremsstrahlungSpectrum.cc.
179{
181 G4double t0 = std::max(tmin, lowestE);
182 if(t0 >= tm) return 0.0;
183
184 t0 /= e;
185 tm /= e;
186
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
197
198 do {
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
208 }
209
211 } while (q > fun);
212
213 tgam *= e;
214
215 p.clear();
216
217 return tgam;
218}
The documentation for this class was generated from the following files: