Geant4
10.02.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
PhysListEmStandardWVI.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: PhysListEmStandardWVI.cc 86064 2014-11-07 08:49:32Z gcosmo $
30
//
31
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34
#include "
PhysListEmStandardWVI.hh
"
35
#include "
G4ParticleDefinition.hh
"
36
#include "
G4ProcessManager.hh
"
37
38
#include "
G4ComptonScattering.hh
"
39
#include "
G4GammaConversion.hh
"
40
#include "
G4PhotoElectricEffect.hh
"
41
42
#include "
G4eMultipleScattering.hh
"
43
#include "
G4WentzelVIModel.hh
"
44
#include "
G4CoulombScattering.hh
"
45
#include "
G4eIonisation.hh
"
46
#include "
G4eBremsstrahlung.hh
"
47
#include "
G4eplusAnnihilation.hh
"
48
49
#include "
G4MuMultipleScattering.hh
"
50
#include "
G4MuIonisation.hh
"
51
#include "
G4MuBremsstrahlung.hh
"
52
#include "
G4MuPairProduction.hh
"
53
54
#include "
G4hMultipleScattering.hh
"
55
#include "
G4hIonisation.hh
"
56
#include "
G4hBremsstrahlung.hh
"
57
#include "
G4hPairProduction.hh
"
58
59
#include "
G4ionIonisation.hh
"
60
#include "
G4IonParametrisedLossModel.hh
"
61
#include "
G4NuclearStopping.hh
"
62
63
#include "
G4EmProcessOptions.hh
"
64
#include "
G4MscStepLimitType.hh
"
65
66
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68
PhysListEmStandardWVI::PhysListEmStandardWVI
(
const
G4String
&
name
)
69
:
G4VPhysicsConstructor
(name)
70
{}
71
72
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74
PhysListEmStandardWVI::~PhysListEmStandardWVI
()
75
{}
76
77
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
79
void
PhysListEmStandardWVI::ConstructProcess
()
80
{
81
// Add standard EM Processes
82
//
83
84
aParticleIterator
->reset();
85
while
( (*
aParticleIterator
)() ){
86
G4ParticleDefinition
* particle =
aParticleIterator
->value();
87
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
88
G4String
particleName = particle->
GetParticleName
();
89
90
if
(particleName ==
"gamma"
) {
91
// gamma
92
pmanager->
AddDiscreteProcess
(
new
G4PhotoElectricEffect
);
93
pmanager->
AddDiscreteProcess
(
new
G4ComptonScattering
);
94
pmanager->
AddDiscreteProcess
(
new
G4GammaConversion
);
95
96
}
else
if
(particleName ==
"e-"
) {
97
//electron
98
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
99
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
100
pmanager->
AddProcess
(msc, -1, 1, 1);
101
pmanager->
AddProcess
(
new
G4eIonisation
, -1, 2, 2);
102
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1, 3, 3);
103
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 4);
104
105
}
else
if
(particleName ==
"e+"
) {
106
//positron
107
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
108
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
109
pmanager->
AddProcess
(msc, -1, 1, 1);
110
pmanager->
AddProcess
(
new
G4eIonisation
, -1, 2, 2);
111
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1, 3, 3);
112
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 4);
113
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 5);
114
115
}
else
if
(particleName ==
"mu+"
||
116
particleName ==
"mu-"
) {
117
//muon
118
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
119
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
120
pmanager->
AddProcess
(msc, -1, 1, 1);
121
pmanager->
AddProcess
(
new
G4MuIonisation
, -1, 2, 2);
122
pmanager->
AddProcess
(
new
G4MuBremsstrahlung
, -1, 3, 3);
123
pmanager->
AddProcess
(
new
G4MuPairProduction
, -1, 4, 4);
124
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 5);
125
126
}
else
if
( particleName ==
"proton"
||
127
particleName ==
"pi-"
||
128
particleName ==
"pi+"
) {
129
//proton
130
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
131
pmanager->
AddProcess
(msc, -1, 1, 1);
132
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
133
pmanager->
AddProcess
(
new
G4hBremsstrahlung
, -1, 3, 3);
134
pmanager->
AddProcess
(
new
G4hPairProduction
, -1, 4, 4);
135
136
}
else
if
( particleName ==
"alpha"
||
137
particleName ==
"He3"
) {
138
//alpha
139
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
140
pmanager->
AddProcess
(msc, -1, 1, 1);
141
pmanager->
AddProcess
(
new
G4ionIonisation
, -1, 2, 2);
142
pmanager->
AddProcess
(
new
G4NuclearStopping
, -1, 3,-1);
143
144
}
else
if
( particleName ==
"GenericIon"
) {
145
//Ions
146
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
147
pmanager->
AddProcess
(msc, -1, 1, 1);
148
G4ionIonisation
* ionIoni =
new
G4ionIonisation
();
149
ionIoni->
SetEmModel
(
new
G4IonParametrisedLossModel
());
150
pmanager->
AddProcess
(ionIoni, -1, 2, 2);
151
pmanager->
AddProcess
(
new
G4NuclearStopping
, -1, 3,-1);
152
153
}
else
if
((!particle->
IsShortLived
()) &&
154
(particle->
GetPDGCharge
() != 0.0) &&
155
(particle->
GetParticleName
() !=
"chargedgeantino"
)) {
156
//all others charged particles except geantino
157
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
158
pmanager->
AddProcess
(msc, -1, 1, 1);
159
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
160
}
161
}
162
163
// Em options
164
//
165
// Main options and setting parameters are shown here.
166
// Several of them have default values.
167
//
168
G4EmProcessOptions
emOptions;
169
170
//multiple coulomb scattering
171
//
172
emOptions.
SetPolarAngleLimit
(0.2);
173
}
174
175
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176
G4hMultipleScattering.hh
G4eIonisation.hh
G4IonParametrisedLossModel.hh
G4MscStepLimitType.hh
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4EmProcessOptions
Definition:
G4EmProcessOptions.hh:64
G4ComptonScattering.hh
G4hBremsstrahlung.hh
G4hIonisation.hh
name
G4String name
Definition:
TRTMaterials.hh:40
G4hPairProduction
Definition:
G4hPairProduction.hh:59
G4NuclearStopping.hh
G4MuPairProduction.hh
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4hPairProduction.hh
G4eMultipleScattering.hh
G4ComptonScattering
Definition:
G4ComptonScattering.hh:71
G4MuIonisation.hh
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:111
G4MuBremsstrahlung.hh
G4ParticleDefinition::GetProcessManager
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition::GetParticleName
const G4String & GetParticleName() const
Definition:
G4ParticleDefinition.hh:159
G4CoulombScattering
Definition:
G4CoulombScattering.hh:55
G4ionIonisation.hh
G4CoulombScattering.hh
PhysListEmStandardWVI::PhysListEmStandardWVI
PhysListEmStandardWVI(const G4String &name="standardWVI")
Definition:
PhysListEmStandardWVI.cc:68
G4GammaConversion
Definition:
G4GammaConversion.hh:75
G4hBremsstrahlung
Definition:
G4hBremsstrahlung.hh:61
G4GammaConversion.hh
G4MuMultipleScattering
Definition:
G4MuMultipleScattering.hh:61
aParticleIterator
#define aParticleIterator
Definition:
G4VPhysicsConstructor.hh:119
G4ProcessManager::AddProcess
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
Definition:
G4ProcessManager.cc:410
G4IonParametrisedLossModel
Definition:
G4IonParametrisedLossModel.hh:93
G4NuclearStopping
Definition:
G4NuclearStopping.hh:66
G4ProcessManager.hh
G4ParticleDefinition.hh
G4PhotoElectricEffect.hh
G4ionIonisation
Definition:
G4ionIonisation.hh:78
G4MuBremsstrahlung
Definition:
G4MuBremsstrahlung.hh:78
G4MuIonisation
Definition:
G4MuIonisation.hh:85
PhysListEmStandardWVI::~PhysListEmStandardWVI
~PhysListEmStandardWVI()
Definition:
PhysListEmStandardWVI.cc:74
G4ParticleDefinition::IsShortLived
G4bool IsShortLived() const
Definition:
G4ParticleDefinition.hh:196
G4eBremsstrahlung
Definition:
G4eBremsstrahlung.hh:81
G4MuMultipleScattering.hh
G4WentzelVIModel.hh
G4eIonisation
Definition:
G4eIonisation.hh:80
G4eplusAnnihilation.hh
G4VMultipleScattering::AddEmModel
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
Definition:
G4VMultipleScattering.cc:146
G4VEnergyLossProcess::SetEmModel
void SetEmModel(G4VEmModel *, G4int index=1)
Definition:
G4VEnergyLossProcess.cc:411
G4MuPairProduction
Definition:
G4MuPairProduction.hh:74
G4eBremsstrahlung.hh
G4ProcessManager
Definition:
G4ProcessManager.hh:106
G4eplusAnnihilation
Definition:
G4eplusAnnihilation.hh:65
G4hMultipleScattering
Definition:
G4hMultipleScattering.hh:62
PhysListEmStandardWVI::ConstructProcess
virtual void ConstructProcess()
Definition:
PhysListEmStandardWVI.cc:79
G4ParticleDefinition::GetPDGCharge
G4double GetPDGCharge() const
Definition:
G4ParticleDefinition.hh:163
PhysListEmStandardWVI.hh
Definition of the PhysListEmStandardWVI class.
G4EmProcessOptions.hh
G4VPhysicsConstructor
Definition:
G4VPhysicsConstructor.hh:121
G4WentzelVIModel
Definition:
G4WentzelVIModel.hh:66
G4hIonisation
Definition:
G4hIonisation.hh:85
G4String
Definition:
G4String.hh:45
G4EmProcessOptions::SetPolarAngleLimit
void SetPolarAngleLimit(G4double val)
Definition:
G4EmProcessOptions.cc:236
geant4.10.02.p02
examples
extended
medical
electronScattering
src
PhysListEmStandardWVI.cc
Generated on Fri Jul 1 2016 14:10:52 for Geant4 by
1.8.8