Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4hBremsstrahlungModel.cc
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// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4hBremsstrahlungModel
33//
34// Author: Vladimir Ivanchenko on base of G4MuBremsstrahlungModel
35//
36// Creation date: 28.02.2008
37//
38// Modifications:
39//
40
41//
42// Class Description:
43//
44//
45// -------------------------------------------------------------------
46//
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49
52#include "G4SystemOfUnits.hh"
53#include "G4Log.hh"
54#include "G4NistManager.hh"
55
56using namespace std;
57
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
71 G4double tkin,
72 G4double Z,
73 G4double gammaEnergy)
74// differential cross section
75{
76 G4double dxsection = 0.;
77
78 if(gammaEnergy > tkin) return dxsection;
79 // G4cout << "G4hBremsstrahlungModel m= " << mass
80 // << " " << particle->GetParticleName() << G4endl;
81 G4double E = tkin + mass ;
82 G4double v = gammaEnergy/E ;
83 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy);
84 G4double rab0=delta*sqrte ;
85
86 G4int iz = std::max(G4lrint(Z), 1);
87
88 G4double z13 = 1.0/nist->GetZ13(iz);
89 G4double dn = mass*nist->GetA27(iz)/(70.*MeV);
90
91 G4double b = (1 == iz) ? bh : btf;
92
93 // nucleus contribution logarithm
94 G4double rab1=b*z13;
95 G4double fn=G4Log(rab1/(dn*(electron_mass_c2+rab0*rab1))*
96 (mass+delta*(dn*sqrte-2.))) ;
97 fn = std::max(fn, 0.0);
98
99 G4double x = 1.0 - v;
100 if(particle->GetPDGSpin() != 0) { x += 0.75*v*v; }
101
102 dxsection = coeff*x*Z*Z*fn/gammaEnergy;
103
104 return dxsection;
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4ParticleDefinition * particle
G4double GetA27(G4int Z) const
G4double GetZ13(G4double Z) const
G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy) override
G4hBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="hBrem")
int G4lrint(double ad)
Definition templates.hh:134