Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Generator2BN Class Reference

#include <G4Generator2BN.hh>

+ Inheritance diagram for G4Generator2BN:

Public Member Functions

 G4Generator2BN (const G4String &name="")
 
virtual ~G4Generator2BN ()
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) override
 
void PrintGeneratorInformation () const override
 
void SetInterpolationThetaIncrement (G4double increment)
 
G4double GetInterpolationThetaIncrement ()
 
void SetGammaCutValue (G4double cutValue)
 
G4double GetGammaCutValue ()
 
void ConstructMajorantSurface ()
 
G4Generator2BNoperator= (const G4Generator2BN &right)=delete
 
 G4Generator2BN (const G4Generator2BN &)=delete
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
 
virtual void PrintGeneratorInformation () const
 
const G4StringGetName () const
 
G4VEmAngularDistributionoperator= (const G4VEmAngularDistribution &right)=delete
 
 G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete
 

Protected Member Functions

G4double CalculateFkt (G4double k, G4double theta, G4double A, G4double c) const
 
G4double Calculatedsdkdt (G4double kout, G4double theta, G4double Eel) const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 61 of file G4Generator2BN.hh.

Constructor & Destructor Documentation

◆ G4Generator2BN() [1/2]

G4Generator2BN::G4Generator2BN ( const G4String name = "")
explicit

Definition at line 157 of file G4Generator2BN.cc.

158 : G4VEmAngularDistribution("AngularGen2BN")
159{
160 b = 1.2;
161 index_min = -300;
162 index_max = 319;
163
164 // Set parameters minimum limits Ekelectron = 250 eV and kphoton = 100 eV
165 kmin = 100*eV;
166 Ekmin = 250*eV;
167 kcut = 100*eV;
168
169 // Increment Theta value for surface interpolation
170 dtheta = 0.1*rad;
171
172 nwarn = 0;
173
174 // Construct Majorant Surface to 2BN cross-section
175 // (we dont need this. Pre-calculated values are used instead due to performance issues
176 // ConstructMajorantSurface();
177}

◆ ~G4Generator2BN()

virtual G4Generator2BN::~G4Generator2BN ( )
inlinevirtual

Definition at line 67 of file G4Generator2BN.hh.

67{;};

◆ G4Generator2BN() [2/2]

G4Generator2BN::G4Generator2BN ( const G4Generator2BN )
delete

Member Function Documentation

◆ Calculatedsdkdt()

G4double G4Generator2BN::Calculatedsdkdt ( G4double  kout,
G4double  theta,
G4double  Eel 
) const
protected

Definition at line 265 of file G4Generator2BN.cc.

266{
267 G4double dsdkdt_value = 0.;
268 G4double Z = 1;
269 // classic radius (in cm)
270 G4double r0 = 2.82E-13;
271 //squared classic radius (in barn)
272 G4double r02 = r0*r0*1.0E+24;
273
274 // Photon energy cannot be greater than electron kinetic energy
275 if(kout > (Eel-electron_mass_c2)){
276 dsdkdt_value = 0;
277 return dsdkdt_value;
278 }
279
280 G4double E0 = Eel/electron_mass_c2;
281 G4double k = kout/electron_mass_c2;
282 G4double E = E0 - k;
283
284 if(E <= 1*MeV ){ // Kinematic limit at 1 MeV !!
285 dsdkdt_value = 0;
286 return dsdkdt_value;
287 }
288
289 G4double p0 = std::sqrt(E0*E0-1);
290 G4double p = std::sqrt(E*E-1);
291 G4double LL = std::log((E*E0-1+p*p0)/(E*E0-1-p*p0));
292 G4double delta0 = E0 - p0*std::cos(theta);
293 G4double epsilon = std::log((E+p)/(E-p));
294 G4double Z2 = Z*Z;
295 G4double sintheta2 = std::sin(theta)*std::sin(theta);
296 G4double E02 = E0*E0;
297 G4double E2 = E*E;
298 G4double p02 = E0*E0-1;
299 G4double k2 = k*k;
300 G4double delta02 = delta0*delta0;
301 G4double delta04 = delta02* delta02;
302 G4double Q = std::sqrt(p02+k2-2*k*p0*std::cos(theta));
303 G4double Q2 = Q*Q;
304 G4double epsilonQ = std::log((Q+p)/(Q-p));
305
306
307 dsdkdt_value = Z2 * (r02/(8*pi*137)) * (1/k) * (p/p0) *
308 ( (8 * (sintheta2*(2*E02+1))/(p02*delta04)) -
309 ((2*(5*E02+2*E*E0+3))/(p02 * delta02)) -
310 ((2*(p02-k2))/((Q2*delta02))) +
311 ((4*E)/(p02*delta0)) +
312 (LL/(p*p0))*(
313 ((4*E0*sintheta2*(3*k-p02*E))/(p02*delta04)) +
314 ((4*E02*(E02+E2))/(p02*delta02)) +
315 ((2-2*(7*E02-3*E*E0+E2))/(p02*delta02)) +
316 (2*k*(E02+E*E0-1))/((p02*delta0))
317 ) -
318 ((4*epsilon)/(p*delta0)) +
319 ((epsilonQ)/(p*Q))*
320 (4/delta02-(6*k/delta0)-(2*k*(p02-k2))/(Q2*delta0))
321 );
322
323 dsdkdt_value = dsdkdt_value*std::sin(theta);
324 return dsdkdt_value;
325}
G4double epsilon(G4double density, G4double temperature)
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
const G4double pi

Referenced by ConstructMajorantSurface(), and SampleDirection().

◆ CalculateFkt()

G4double G4Generator2BN::CalculateFkt ( G4double  k,
G4double  theta,
G4double  A,
G4double  c 
) const
protected

Definition at line 256 of file G4Generator2BN.cc.

257{
258 G4double Fkt_value = 0;
259 Fkt_value = A*std::pow(k,-b)*theta/(1+c*theta*theta);
260 return Fkt_value;
261}
const G4double A[17]

Referenced by ConstructMajorantSurface().

◆ ConstructMajorantSurface()

void G4Generator2BN::ConstructMajorantSurface ( )

Definition at line 329 of file G4Generator2BN.cc.

330{
331 G4double Eel;
332 G4int vmax;
333 G4double Ek;
334 G4double k, theta, k0, theta0, thetamax;
335 G4double ds, df, dsmax, dk, dt;
336 G4double ratmin;
337 G4double ratio = 0;
338 G4double Vds, Vdf;
339 G4double A, c;
340
341 G4cout << "**** Constructing Majorant Surface for 2BN Distribution ****" << G4endl;
342
343 if(kcut > kmin) kmin = kcut;
344
345 G4int i = 0;
346 // Loop on energy index
347 for(G4int index = index_min; index < index_max; index++){
348
349 G4double fraction = index/100.;
350 Ek = std::pow(10.,fraction);
351 Eel = Ek + electron_mass_c2;
352
353 // find x-section maximum at k=kmin
354 dsmax = 0.;
355 thetamax = 0.;
356
357 for(theta = 0.; theta < pi; theta = theta + dtheta){
358
359 ds = Calculatedsdkdt(kmin, theta, Eel);
360
361 if(ds > dsmax){
362 dsmax = ds;
363 thetamax = theta;
364 }
365 }
366
367 // Compute surface paremeters at kmin
368 if(Ek < kmin || thetamax == 0){
369 c = 0;
370 A = 0;
371 }else{
372 c = 1/(thetamax*thetamax);
373 A = 2*std::sqrt(c)*dsmax/(std::pow(kmin,-b));
374 }
375
376 // look for correction factor to normalization at kmin
377 ratmin = 1.;
378
379 // Volume under surfaces
380 Vds = 0.;
381 Vdf = 0.;
382 k0 = 0.;
383 theta0 = 0.;
384
385 vmax = G4int(100.*std::log10(Ek/kmin));
386
387 for(G4int v = 0; v < vmax; v++){
388 G4double fractionLocal = (v/100.);
389 k = std::pow(10.,fractionLocal)*kmin;
390
391 for(theta = 0.; theta < pi; theta = theta + dtheta){
392 dk = k - k0;
393 dt = theta - theta0;
394 ds = Calculatedsdkdt(k,theta, Eel);
395 Vds = Vds + ds*dk*dt;
396 df = CalculateFkt(k,theta, A, c);
397 Vdf = Vdf + df*dk*dt;
398
399 if(ds != 0.){
400 if(df != 0.) ratio = df/ds;
401 }
402
403 if(ratio < ratmin && ratio != 0.){
404 ratmin = ratio;
405 }
406 }
407 }
408
409 // sampling efficiency
410 Atab[i] = A/ratmin * 1.04;
411 ctab[i] = c;
412
413// G4cout << Ek << " " << i << " " << index << " " << Atab[i] << " " << ctab[i] << " " << G4endl;
414 i++;
415 }
416}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double CalculateFkt(G4double k, G4double theta, G4double A, G4double c) const
G4double Calculatedsdkdt(G4double kout, G4double theta, G4double Eel) const

◆ GetGammaCutValue()

G4double G4Generator2BN::GetGammaCutValue ( )
inline

Definition at line 80 of file G4Generator2BN.hh.

80{return kcut;};

◆ GetInterpolationThetaIncrement()

G4double G4Generator2BN::GetInterpolationThetaIncrement ( )
inline

Definition at line 77 of file G4Generator2BN.hh.

77{return dtheta;};

◆ operator=()

G4Generator2BN & G4Generator2BN::operator= ( const G4Generator2BN right)
delete

◆ PrintGeneratorInformation()

void G4Generator2BN::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 420 of file G4Generator2BN.cc.

421{
422 G4cout << "\n" << G4endl;
423 G4cout << "Bremsstrahlung Angular Generator is 2BN Generator from 2BN Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
424 G4cout << "\n" << G4endl;
425}

◆ SampleDirection()

G4ThreeVector & G4Generator2BN::SampleDirection ( const G4DynamicParticle dp,
G4double  out_energy,
G4int  Z,
const G4Material mat = nullptr 
)
overridevirtual

Implements G4VEmAngularDistribution.

Definition at line 181 of file G4Generator2BN.cc.

185{
186 G4double Ek = dp->GetKineticEnergy();
187 G4double Eel = dp->GetTotalEnergy();
188 if(Eel > 2*MeV) {
189 return fGenerator2BS.SampleDirection(dp, out_energy, 0);
190 }
191
192 G4double k = Eel - out_energy;
193
194 G4double t;
195 G4double cte2;
196 G4double y, u;
197 G4double ds;
198 G4double dmax;
199
200 // find table index
201 G4int index = G4int(std::log10(Ek/MeV)*100) - index_min;
202 if(index > index_max) { index = index_max; }
203 else if(index < 0) { index = 0; }
204
205 G4double c = ctab[index];
206 G4double A = Atab[index];
207 if(index < index_max) { A = std::max(A, Atab[index+1]); }
208
209 do{
210 // generate k accordimg to std::pow(k,-b)
211
212 // normalization constant
213 // cte1 = (1-b)/(std::pow(kmax,(1-b))-std::pow(kmin2,(1-b)));
214 // y = G4UniformRand();
215 // k = std::pow(((1-b)/cte1*y+std::pow(kmin2,(1-b))),(1/(1-b)));
216
217 // generate theta accordimg to theta/(1+c*std::pow(theta,2)
218 // Normalization constant
219 cte2 = 2*c/std::log(1+c*pi2);
220
221 y = G4UniformRand();
222 t = std::sqrt((G4Exp(2*c*y/cte2)-1)/c);
223 u = G4UniformRand();
224
225 // point acceptance algorithm
226 //fk = std::pow(k,-b);
227 //ft = t/(1+c*t*t);
228 dmax = A*std::pow(k,-b)*t/(1+c*t*t);
229 ds = Calculatedsdkdt(k,t,Eel);
230
231 // violation point
232 if(ds > dmax && nwarn >= 20) {
233 ++nwarn;
234 G4cout << "### WARNING in G4Generator2BN: Ekin(MeV)= " << Ek/MeV
235 << " D(Ekin,k)/Dmax-1= " << (ds/dmax - 1)
236 << " results are not reliable!"
237 << G4endl;
238 if(20 == nwarn) {
239 G4cout << " WARNING in G4Generator2BN is closed" << G4endl;
240 }
241 }
242
243 } while(u*dmax > ds || t > CLHEP::pi);
244
245 G4double sint = std::sin(t);
246 G4double phi = twopi*G4UniformRand();
247
248 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi),std::cos(t));
250
251 return fLocalDirection;
252}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) override

◆ SetGammaCutValue()

void G4Generator2BN::SetGammaCutValue ( G4double  cutValue)
inline

Definition at line 79 of file G4Generator2BN.hh.

79{kcut = cutValue;};

◆ SetInterpolationThetaIncrement()

void G4Generator2BN::SetInterpolationThetaIncrement ( G4double  increment)
inline

Definition at line 76 of file G4Generator2BN.hh.

76{dtheta = increment;};

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