Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
standard
src
G4eBremsstrahlung.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$
27
//
28
// -------------------------------------------------------------------
29
//
30
// GEANT4 Class file
31
//
32
//
33
// File name: G4eBremsstrahlung
34
//
35
// Author: Michel Maire
36
//
37
// Creation date: 26.06.1996
38
//
39
// Modifications:
40
//
41
// 26-09-96 extension of the total crosssection above 100 GeV, M.Maire
42
// 1-10-96 new type G4OrderedTable; ComputePartialSumSigma(), M.Maire
43
// 16-10-96 DoIt() call to the non static GetEnergyCuts(), L.Urban
44
// 13-12-96 Sign corrected in grejmax and greject
45
// error definition of screenvar, L.Urban
46
// 20-03-97 new energy loss+ionisation+brems scheme, L.Urban
47
// 07-04-98 remove 'tracking cut' of the diffracted particle, MMa
48
// 13-08-98 new methods SetBining() PrintInfo()
49
// 03-03-99 Bug fixed in LPM effect, L.Urban
50
// 10-02-00 modifications , new e.m. structure, L.Urban
51
// 07-08-00 new cross section/en.loss parametrisation, LPM flag , L.Urban
52
// 21-09-00 corrections in the LPM implementation, L.Urban
53
// 28-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
54
// 09-08-01 new methods Store/Retrieve PhysicsTable (mma)
55
// 17-09-01 migration of Materials to pure STL (mma)
56
// 21-09-01 completion of RetrievePhysicsTable() (mma)
57
// 29-10-01 all static functions no more inlined (mma)
58
// 08-11-01 particleMass becomes a local variable
59
// 30-04-02 V.Ivanchenko update to new design
60
// 23-12-02 Change interface in order to move to cut per region (VI)
61
// 26-12-02 Secondary production moved to derived classes (VI)
62
// 23-05-03 Define default integral + BohrFluctuations (V.Ivanchenko)
63
// 08-08-03 STD substitute standard (V.Ivanchenko)
64
// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
65
// 04-11-04 add gamma threshold (V.Ivanchenko)
66
// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
67
// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
68
// 22-05-06 Use gammaThreshold from manager (V.Ivantchenko)
69
// 15-01-07 use SetEmModel() from G4VEnergyLossProcess (mma)
70
// use RelEmModel above 1GeV (AS & VI)
71
// 13-11-08 reenable LPM switch (A.Schaelicke)
72
// -------------------------------------------------------------------
73
//
74
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77
#include "
G4eBremsstrahlung.hh
"
78
#include "
G4SystemOfUnits.hh
"
79
#include "
G4Gamma.hh
"
80
#include "
G4eBremsstrahlungModel.hh
"
81
#include "
G4SeltzerBergerModel.hh
"
82
#include "
G4eBremsstrahlungRelModel.hh
"
83
#include "
G4UnitsTable.hh
"
84
#include "
G4LossTableManager.hh
"
85
86
#include "
G4ProductionCutsTable.hh
"
87
#include "
G4MaterialCutsCouple.hh
"
88
89
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90
91
using namespace
std;
92
93
G4eBremsstrahlung::G4eBremsstrahlung
(
const
G4String
&
name
):
94
G4VEnergyLossProcess
(name),
95
isInitialised(false)
96
{
97
SetProcessSubType
(
fBremsstrahlung
);
98
SetSecondaryParticle
(
G4Gamma::Gamma
());
99
SetIonisation
(
false
);
100
}
101
102
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103
104
G4eBremsstrahlung::~G4eBremsstrahlung
()
105
{}
106
107
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108
109
G4bool
G4eBremsstrahlung::IsApplicable
(
const
G4ParticleDefinition
&
p
)
110
{
111
return
(&p ==
G4Electron::Electron
() || &p ==
G4Positron::Positron
());
112
}
113
114
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115
116
void
117
G4eBremsstrahlung::InitialiseEnergyLossProcess
(
const
G4ParticleDefinition
*,
118
const
G4ParticleDefinition
*)
119
{
120
if
(!
isInitialised
) {
121
122
// if (!EmModel(1)) { SetEmModel(new G4eBremsstrahlungModel(), 1); }
123
if
(!
EmModel
(1)) {
SetEmModel
(
new
G4SeltzerBergerModel
(), 1); }
124
if
(!
EmModel
(2)) {
SetEmModel
(
new
G4eBremsstrahlungRelModel
(), 2); }
125
126
G4double
energyLimit = 1*
GeV
;
127
128
EmModel
(1)->
SetLowEnergyLimit
(
MinKinEnergy
());
129
EmModel
(1)->
SetHighEnergyLimit
(energyLimit);
130
EmModel
(2)->
SetLowEnergyLimit
(energyLimit);
131
EmModel
(2)->
SetHighEnergyLimit
(
MaxKinEnergy
());
132
133
G4VEmFluctuationModel
* fm = 0;
134
AddEmModel
(1,
EmModel
(1), fm);
135
AddEmModel
(2,
EmModel
(2), fm);
136
isInitialised
=
true
;
137
}
138
G4LossTableManager
* man =
G4LossTableManager::Instance
();
139
G4double
eth = man->
BremsstrahlungTh
();
140
EmModel
(1)->
SetSecondaryThreshold
(eth);
141
EmModel
(2)->
SetSecondaryThreshold
(eth);
142
143
// Only high energy model LMP flag is ON/OFF
144
EmModel
(1)->
SetLPMFlag
(
false
);
145
EmModel
(2)->
SetLPMFlag
(man->
LPMFlag
());
146
}
147
148
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
150
void
G4eBremsstrahlung::PrintInfo
()
151
{
152
if
(
EmModel
(1)) {
153
G4LossTableManager
* man =
G4LossTableManager::Instance
();
154
G4double
eth = man->
BremsstrahlungTh
();
155
G4cout
<<
" LPM flag: "
<< man->
LPMFlag
() <<
" for E > "
156
<<
EmModel
(1)->
HighEnergyLimit
()/
GeV
<<
" GeV"
;
157
if
(eth <
DBL_MAX
) {
158
G4cout
<<
", HighEnergyThreshold(GeV)= "
<< eth/
GeV
;
159
}
160
G4cout
<<
G4endl
;
161
}
162
}
163
164
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Generated on Sat May 25 2013 14:33:37 for Geant4 by
1.8.4