Geant4
10.03.p01
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
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
// $Id: G4hBremsstrahlungModel.cc 74020 2013-09-19 13:38:38Z gcosmo $
27
//
28
// -------------------------------------------------------------------
29
//
30
// GEANT4 Class file
31
//
32
//
33
// File name: G4hBremsstrahlungModel
34
//
35
// Author: Vladimir Ivanchenko on base of G4MuBremsstrahlungModel
36
//
37
// Creation date: 28.02.2008
38
//
39
// Modifications:
40
//
41
42
//
43
// Class Description:
44
//
45
//
46
// -------------------------------------------------------------------
47
//
48
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
51
#include "
G4hBremsstrahlungModel.hh
"
52
#include "
G4PhysicalConstants.hh
"
53
#include "
G4SystemOfUnits.hh
"
54
#include "
G4Log.hh
"
55
56
using namespace
std;
57
58
G4hBremsstrahlungModel::G4hBremsstrahlungModel
(
const
G4ParticleDefinition
*
p
,
59
const
G4String
& nam)
60
:
G4MuBremsstrahlungModel
(p, nam)
61
{}
62
63
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65
G4hBremsstrahlungModel::~G4hBremsstrahlungModel
()
66
{}
67
68
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70
G4double
G4hBremsstrahlungModel::ComputeDMicroscopicCrossSection
(
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 =
G4int
(Z);
87
if
(iz < 1) { iz = 1; }
88
89
G4double
z13 = 1.0/
nist
->
GetZ13
(iz);
90
G4double
dn = mass*
nist
->
GetA27
(iz)/(70.*
MeV
);
91
92
G4double
b =
btf
;
93
if
(1 == iz) b =
bh
;
94
95
// nucleus contribution logarithm
96
G4double
rab1=b*z13;
97
G4double
fn=
G4Log
(rab1/(dn*(
electron_mass_c2
+rab0*rab1))*
98
(mass+delta*(dn*sqrte-2.))) ;
99
if
(fn <0.) fn = 0. ;
100
101
G4double
x = 1.0 - v;
102
if
(
particle
->
GetPDGSpin
() != 0) { x += 0.75*v*v; }
103
104
dxsection =
coeff
*x*Z*Z*fn/gammaEnergy;
105
106
return
dxsection;
107
}
108
109
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4MuBremsstrahlungModel::btf
G4double btf
Definition:
G4MuBremsstrahlungModel.hh:147
p
const char * p
Definition:
xmltok.h:285
G4NistManager::GetA27
G4double GetA27(G4int Z) const
Definition:
G4NistManager.hh:586
G4hBremsstrahlungModel::ComputeDMicroscopicCrossSection
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy) override
Definition:
G4hBremsstrahlungModel.cc:70
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:72
G4int
int G4int
Definition:
G4Types.hh:78
G4MuBremsstrahlungModel::sqrte
G4double sqrte
Definition:
G4MuBremsstrahlungModel.hh:144
CLHEP::electron_mass_c2
static constexpr double electron_mass_c2
Definition:
PhysicalConstants.h:79
G4hBremsstrahlungModel.hh
G4NistManager::GetZ13
G4double GetZ13(G4double Z) const
Definition:
G4NistManager.hh:572
G4MuBremsstrahlungModel::particle
const G4ParticleDefinition * particle
Definition:
G4MuBremsstrahlungModel.hh:138
G4MuBremsstrahlungModel::coeff
G4double coeff
Definition:
G4MuBremsstrahlungModel.hh:143
G4hBremsstrahlungModel::G4hBremsstrahlungModel
G4hBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="hBrem")
Definition:
G4hBremsstrahlungModel.cc:58
G4MuBremsstrahlungModel::nist
G4NistManager * nist
Definition:
G4MuBremsstrahlungModel.hh:139
G4Log.hh
G4PhysicalConstants.hh
G4Log
G4double G4Log(G4double x)
Definition:
G4Log.hh:230
G4hBremsstrahlungModel::~G4hBremsstrahlungModel
virtual ~G4hBremsstrahlungModel()
Definition:
G4hBremsstrahlungModel.cc:65
G4ParticleDefinition::GetPDGSpin
G4double GetPDGSpin() const
Definition:
G4ParticleDefinition.hh:126
MeV
static constexpr double MeV
Definition:
G4SIunits.hh:214
G4MuBremsstrahlungModel
Definition:
G4MuBremsstrahlungModel.hh:72
G4double
double G4double
Definition:
G4Types.hh:76
G4SystemOfUnits.hh
G4MuBremsstrahlungModel::bh
G4double bh
Definition:
G4MuBremsstrahlungModel.hh:145
G4MuBremsstrahlungModel::mass
G4double mass
Definition:
G4MuBremsstrahlungModel.hh:140
Z
G4int Z
Definition:
G4ParticleHPJENDLHEData.cc:262
G4String
Definition:
G4String.hh:60
geant4.10.03.p01
source
processes
electromagnetic
highenergy
src
G4hBremsstrahlungModel.cc
Generated on Thu Mar 16 2017 22:38:02 for Geant4 by
1.8.5