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
advanced
human_phantom
src
G4HumanPhantomPhysicsList.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
// Code developed by:
28
// S. Guatelli, G. Guerrieri and M. G. Pia
29
//
30
//
31
//
32
// **********************************
33
// * *
34
// * G4HumanPhantomPhysicsList.cc *
35
// * *
36
// **********************************
37
//
38
// $Id$
39
//
40
#include "
G4HumanPhantomPhysicsList.hh
"
41
42
#include "
G4SystemOfUnits.hh
"
43
#include "
G4ParticleDefinition.hh
"
44
#include "
G4ProductionCutsTable.hh
"
45
#include "
G4ProcessManager.hh
"
46
#include "
G4ParticleTypes.hh
"
47
#include "
G4UnitsTable.hh
"
48
#include "
G4ios.hh
"
49
// gamma
50
#include "
G4PhotoElectricEffect.hh
"
51
#include "
G4LivermorePhotoElectricModel.hh
"
52
53
#include "
G4ComptonScattering.hh
"
54
#include "
G4LivermoreComptonModel.hh
"
55
56
#include "
G4GammaConversion.hh
"
57
#include "
G4LivermoreGammaConversionModel.hh
"
58
59
#include "
G4RayleighScattering.hh
"
60
#include "
G4LivermoreRayleighModel.hh
"
61
62
// e-
63
#include "
G4eMultipleScattering.hh
"
64
#include "
G4eIonisation.hh
"
65
#include "
G4LivermoreIonisationModel.hh
"
66
#include "
G4eBremsstrahlung.hh
"
67
#include "
G4LivermoreBremsstrahlungModel.hh
"
68
69
// e+
70
#include "
G4eIonisation.hh
"
71
#include "
G4eBremsstrahlung.hh
"
72
#include "
G4eplusAnnihilation.hh
"
73
74
G4HumanPhantomPhysicsList::G4HumanPhantomPhysicsList
():
G4VUserPhysicsList
()
75
{
76
SetVerboseLevel
(1);
77
}
78
79
G4HumanPhantomPhysicsList::~G4HumanPhantomPhysicsList
()
80
{
81
}
82
83
void
G4HumanPhantomPhysicsList::ConstructParticle
()
84
{
85
// In this method, static member functions should be called
86
// for all particles which you want to use.
87
// This ensures that objects of these particle types will be
88
// created in the program.
89
90
ConstructBosons
();
91
ConstructLeptons
();
92
93
}
94
95
void
G4HumanPhantomPhysicsList::ConstructBosons
()
96
{
97
// photons
98
G4Gamma::GammaDefinition
();
99
}
100
101
void
G4HumanPhantomPhysicsList::ConstructLeptons
()
102
{
103
// leptons
104
G4Electron::ElectronDefinition
();
105
G4Positron::PositronDefinition
();
106
}
107
108
void
G4HumanPhantomPhysicsList::ConstructProcess
()
109
{
110
AddTransportation
();
111
ConstructEM();
112
}
113
114
void
G4HumanPhantomPhysicsList::ConstructEM()
115
{
116
theParticleIterator
->
reset
();
117
118
while
( (*
theParticleIterator
)() ){
119
120
G4ParticleDefinition
* particle =
theParticleIterator
->
value
();
121
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
122
G4String
particleName = particle->
GetParticleName
();
123
124
// Processes
125
126
if
(particleName ==
"gamma"
) {
127
128
// Photon
129
G4RayleighScattering
* theRayleigh =
new
G4RayleighScattering
();
130
theRayleigh->
SetModel
(
new
G4LivermoreRayleighModel
());
//not strictly necessary
131
pmanager->
AddDiscreteProcess
(theRayleigh);
132
133
G4PhotoElectricEffect
* thePhotoElectricEffect =
new
G4PhotoElectricEffect
();
134
thePhotoElectricEffect->
SetModel
(
new
G4LivermorePhotoElectricModel
());
135
pmanager->
AddDiscreteProcess
(thePhotoElectricEffect);
136
137
G4ComptonScattering
* theComptonScattering =
new
G4ComptonScattering
();
138
theComptonScattering->
SetModel
(
new
G4LivermoreComptonModel
());
139
pmanager->
AddDiscreteProcess
(theComptonScattering);
140
141
G4GammaConversion
* theGammaConversion =
new
G4GammaConversion
();
142
theGammaConversion->
SetModel
(
new
G4LivermoreGammaConversionModel
());
143
pmanager->
AddDiscreteProcess
(theGammaConversion);
144
145
}
else
if
(particleName ==
"e-"
) {
146
// Electron
147
148
G4eMultipleScattering
* msc =
new
G4eMultipleScattering
();
149
msc->
SetStepLimitType
(
fUseDistanceToBoundary
);
150
pmanager->
AddProcess
(msc,-1, 1, 1);
151
152
// Ionisation
153
G4eIonisation
* eIonisation =
new
G4eIonisation
();
154
eIonisation->
SetEmModel
(
new
G4LivermoreIonisationModel
());
155
eIonisation->
SetStepFunction
(0.2, 100*um);
//improved precision in tracking
156
pmanager->
AddProcess
(eIonisation,-1, 2, 2);
157
158
// Bremsstrahlung
159
G4eBremsstrahlung
* eBremsstrahlung =
new
G4eBremsstrahlung
();
160
eBremsstrahlung->
SetEmModel
(
new
G4LivermoreBremsstrahlungModel
());
161
pmanager->
AddProcess
(eBremsstrahlung, -1,-3, 3);
162
163
}
else
if
(particleName ==
"e+"
) {
164
// Positron
165
G4eMultipleScattering
* msc =
new
G4eMultipleScattering
();
166
msc->
SetStepLimitType
(
fUseDistanceToBoundary
);
167
pmanager->
AddProcess
(msc,-1, 1, 1);
168
169
// Ionisation
170
G4eIonisation
* eIonisation =
new
G4eIonisation
();
171
eIonisation->
SetStepFunction
(0.2, 100*um);
//
172
pmanager->
AddProcess
(eIonisation, -1, 2, 2);
173
174
//Bremsstrahlung (use default, no low-energy available)
175
pmanager->
AddProcess
(
new
G4eBremsstrahlung
(), -1,-1, 3);
176
177
//Annihilation
178
pmanager->
AddProcess
(
new
G4eplusAnnihilation
(),0,-1, 4);
179
180
}
181
}
182
}
183
184
void
G4HumanPhantomPhysicsList::SetCuts
()
185
{
186
// The production threshold is fixed to 0.1 mm for all the particles
187
// Secondary particles with a range bigger than 0.1 mm
188
// are generated; otherwise their energy is considered deposited locally
189
190
defaultCutValue
= 0.1 *
mm
;
191
192
const
G4double
cutForGamma =
defaultCutValue
;
193
const
G4double
cutForElectron =
defaultCutValue
;
194
const
G4double
cutForPositron =
defaultCutValue
;
195
196
SetCutValue
(cutForGamma,
"gamma"
);
197
SetCutValue
(cutForElectron,
"e-"
);
198
SetCutValue
(cutForPositron,
"e+"
);
199
200
// Set the secondary production cut lower than 990. eV
201
// Very important for high precision of lowenergy processes at low energies
202
203
G4double
lowLimit = 250. *
eV
;
204
G4double
highLimit = 100. *
GeV
;
205
G4ProductionCutsTable::GetProductionCutsTable
()->
SetEnergyRange
(lowLimit, highLimit);
206
207
if
(
verboseLevel
>0)
DumpCutValuesTable
();
208
}
Generated on Sat May 25 2013 14:32:11 for Geant4 by
1.8.4