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
examples
extended
electromagnetic
TestEm7
src
PhysListEmStandardNR.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
//
28
//
29
// $Id: PhysListEmStandardNR.cc 66995 2013-01-29 14:46:45Z gcosmo $
30
//
31
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34
#include "
PhysListEmStandardNR.hh
"
35
36
#include "
G4SystemOfUnits.hh
"
37
#include "
G4ParticleDefinition.hh
"
38
#include "
G4LossTableManager.hh
"
39
#include "
G4EmProcessOptions.hh
"
40
41
#include "
G4ComptonScattering.hh
"
42
#include "
G4GammaConversion.hh
"
43
#include "
G4PhotoElectricEffect.hh
"
44
#include "
G4RayleighScattering.hh
"
45
#include "
G4PEEffectFluoModel.hh
"
46
#include "
G4KleinNishinaModel.hh
"
47
#include "
G4LowEPComptonModel.hh
"
48
#include "
G4PenelopeGammaConversionModel.hh
"
49
#include "
G4LivermorePhotoElectricModel.hh
"
50
51
#include "
G4eMultipleScattering.hh
"
52
#include "
G4MuMultipleScattering.hh
"
53
#include "
G4hMultipleScattering.hh
"
54
#include "
G4CoulombScattering.hh
"
55
#include "
G4eCoulombScatteringModel.hh
"
56
#include "
G4UrbanMscModel95.hh
"
57
58
#include "
G4eIonisation.hh
"
59
#include "
G4eBremsstrahlung.hh
"
60
#include "
G4Generator2BS.hh
"
61
#include "
G4SeltzerBergerModel.hh
"
62
#include "
G4PenelopeIonisationModel.hh
"
63
#include "
G4UniversalFluctuation.hh
"
64
65
#include "
G4eplusAnnihilation.hh
"
66
#include "
G4UAtomicDeexcitation.hh
"
67
68
#include "
G4MuIonisation.hh
"
69
#include "
G4MuBremsstrahlung.hh
"
70
#include "
G4MuPairProduction.hh
"
71
72
#include "
G4hIonisation.hh
"
73
#include "
G4ionIonisation.hh
"
74
#include "
G4IonParametrisedLossModel.hh
"
75
76
#include "
G4ScreenedNuclearRecoil.hh
"
77
78
#include "
G4PhysicsListHelper.hh
"
79
#include "
G4BuilderType.hh
"
80
81
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82
83
PhysListEmStandardNR::PhysListEmStandardNR
(
const
G4String
&
name
)
84
:
G4VPhysicsConstructor
(name)
85
{
86
G4LossTableManager::Instance
();
87
SetPhysicsType
(
bElectromagnetic
);
88
}
89
90
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91
92
PhysListEmStandardNR::~PhysListEmStandardNR
()
93
{}
94
95
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96
97
void
PhysListEmStandardNR::ConstructProcess
()
98
{
99
G4PhysicsListHelper
* ph =
G4PhysicsListHelper::GetPhysicsListHelper
();
100
101
// muon & hadron bremsstrahlung and pair production
102
G4MuBremsstrahlung
* mub =
new
G4MuBremsstrahlung
();
103
G4MuPairProduction
* mup =
new
G4MuPairProduction
();
104
105
G4ScreenedNuclearRecoil
* nucr =
new
G4ScreenedNuclearRecoil
();
106
G4double
energyLimit = 100.*
MeV
;
107
nucr->
SetMaxEnergyForScattering
(energyLimit);
108
G4eCoulombScatteringModel
* csm =
new
G4eCoulombScatteringModel
();
109
csm->
SetActivationLowEnergyLimit
(energyLimit);
110
111
theParticleIterator
->
reset
();
112
while
( (*
theParticleIterator
)() ){
113
G4ParticleDefinition
* particle =
theParticleIterator
->
value
();
114
G4String
particleName = particle->
GetParticleName
();
115
116
if
(particleName ==
"gamma"
) {
117
118
// Compton scattering
119
G4ComptonScattering
* cs =
new
G4ComptonScattering
;
120
cs->
SetEmModel
(
new
G4KleinNishinaModel
(),1);
121
ph->
RegisterProcess
(cs, particle);
122
123
// Photoelectric
124
G4PhotoElectricEffect
* pe =
new
G4PhotoElectricEffect
();
125
G4VEmModel
* theLivermorePEModel =
new
G4LivermorePhotoElectricModel
();
126
theLivermorePEModel->
SetHighEnergyLimit
(10*
GeV
);
127
pe->
SetEmModel
(theLivermorePEModel,1);
128
ph->
RegisterProcess
(pe, particle);
129
130
// Gamma conversion
131
G4GammaConversion
* gc =
new
G4GammaConversion
();
132
G4VEmModel
* thePenelopeGCModel =
new
G4PenelopeGammaConversionModel
();
133
thePenelopeGCModel->
SetHighEnergyLimit
(1*
GeV
);
134
gc->
SetEmModel
(thePenelopeGCModel,1);
135
ph->
RegisterProcess
(gc, particle);
136
137
// Rayleigh scattering
138
ph->
RegisterProcess
(
new
G4RayleighScattering
(), particle);
139
140
}
else
if
(particleName ==
"e-"
) {
141
142
// ionisation
143
G4eIonisation
* eIoni =
new
G4eIonisation
();
144
eIoni->
SetStepFunction
(0.2, 100*um);
145
//G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
146
//theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
147
//eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
148
149
// bremsstrahlung
150
G4eBremsstrahlung
* eBrem =
new
G4eBremsstrahlung
();
151
152
ph->
RegisterProcess
(
new
G4eMultipleScattering
(), particle);
153
ph->
RegisterProcess
(eIoni, particle);
154
ph->
RegisterProcess
(eBrem, particle);
155
156
}
else
if
(particleName ==
"e+"
) {
157
// ionisation
158
G4eIonisation
* eIoni =
new
G4eIonisation
();
159
eIoni->
SetStepFunction
(0.2, 100*um);
160
// G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
161
// theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
162
// eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
163
164
// bremsstrahlung
165
G4eBremsstrahlung
* eBrem =
new
G4eBremsstrahlung
();
166
167
ph->
RegisterProcess
(
new
G4eMultipleScattering
(), particle);
168
ph->
RegisterProcess
(eIoni, particle);
169
ph->
RegisterProcess
(eBrem, particle);
170
171
// annihilation at rest and in flight
172
ph->
RegisterProcess
(
new
G4eplusAnnihilation
(), particle);
173
174
}
else
if
(particleName ==
"mu+"
||
175
particleName ==
"mu-"
) {
176
177
G4MuIonisation
* muIoni =
new
G4MuIonisation
();
178
muIoni->
SetStepFunction
(0.2, 50*um);
179
180
ph->
RegisterProcess
(muIoni, particle);
181
ph->
RegisterProcess
(mub, particle);
182
ph->
RegisterProcess
(mup, particle);
183
ph->
RegisterProcess
(
new
G4CoulombScattering
(), particle);
184
185
}
else
if
(particleName ==
"alpha"
|| particleName ==
"He3"
) {
186
187
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
188
G4UrbanMscModel95
*
model
=
new
G4UrbanMscModel95
();
189
model->
SetActivationLowEnergyLimit
(energyLimit);
190
msc->
SetEmModel
(model, 1);
191
ph->
RegisterProcess
(msc, particle);
192
193
G4ionIonisation
* ionIoni =
new
G4ionIonisation
();
194
ionIoni->
SetStepFunction
(0.1, 10*um);
195
ph->
RegisterProcess
(ionIoni, particle);
196
197
ph->
RegisterProcess
(nucr, particle);
198
199
}
else
if
(particleName ==
"GenericIon"
) {
200
201
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
202
G4UrbanMscModel95
*
model
=
new
G4UrbanMscModel95
();
203
model->
SetActivationLowEnergyLimit
(energyLimit);
204
msc->
SetEmModel
(model, 1);
205
ph->
RegisterProcess
(msc, particle);
206
207
G4ionIonisation
* ionIoni =
new
G4ionIonisation
();
208
ionIoni->
SetEmModel
(
new
G4IonParametrisedLossModel
());
209
ionIoni->
SetStepFunction
(0.1, 1*um);
210
ph->
RegisterProcess
(ionIoni, particle);
211
212
ph->
RegisterProcess
(nucr, particle);
213
214
}
else
if
(particleName ==
"proton"
||
215
particleName ==
"deuteron"
||
216
particleName ==
"triton"
) {
217
218
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
219
G4UrbanMscModel95
*
model
=
new
G4UrbanMscModel95
();
220
model->
SetActivationLowEnergyLimit
(energyLimit);
221
msc->
SetEmModel
(model, 1);
222
ph->
RegisterProcess
(msc, particle);
223
224
G4hIonisation
* hIoni =
new
G4hIonisation
();
225
hIoni->
SetStepFunction
(0.05, 1*um);
226
ph->
RegisterProcess
(hIoni, particle);
227
228
ph->
RegisterProcess
(nucr, particle);
229
230
}
else
if
((!particle->
IsShortLived
()) &&
231
(particle->
GetPDGCharge
() != 0.0) &&
232
(particle->
GetParticleName
() !=
"chargedgeantino"
)) {
233
//all others charged particles except geantino
234
235
ph->
RegisterProcess
(
new
G4hMultipleScattering
(), particle);
236
ph->
RegisterProcess
(
new
G4hIonisation
(), particle);
237
}
238
}
239
240
// Em options
241
//
242
// Main options and setting parameters are shown here.
243
// Several of them have default values.
244
//
245
G4EmProcessOptions
emOptions;
246
247
//physics tables
248
//
249
emOptions.
SetMinEnergy
(10*
eV
);
250
emOptions.
SetMaxEnergy
(10*
TeV
);
251
emOptions.
SetDEDXBinning
(12*20);
252
emOptions.
SetLambdaBinning
(12*20);
253
254
// scattering
255
emOptions.
SetPolarAngleLimit
(0.0);
256
257
// Deexcitation
258
G4VAtomDeexcitation
* de =
new
G4UAtomicDeexcitation
();
259
G4LossTableManager::Instance
()->
SetAtomDeexcitation
(de);
260
de->
SetFluo
(
true
);
261
}
262
263
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264
Generated on Sat May 25 2013 14:32:21 for Geant4 by
1.8.4