Geant4.10
Main Page
Related Pages
Modules
Namespaces
Data Structures
Files
File List
Globals
All
Data Structures
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4.10.00.p01
examples
extended
electromagnetic
TestEm5
src
electromagnetic/TestEm5/src/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
//
26
/// \file electromagnetic/TestEm5/src/PhysListEmStandardWVI.cc
27
/// \brief Implementation of the PhysListEmStandardWVI class
28
//
29
// $Id: PhysListEmStandardWVI.cc 68585 2013-04-01 23:35:07Z adotti $
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 "
G4IonCoulombScatteringModel.hh
"
46
#include "
G4eIonisation.hh
"
47
#include "
G4eBremsstrahlung.hh
"
48
#include "
G4eplusAnnihilation.hh
"
49
50
#include "
G4MuMultipleScattering.hh
"
51
#include "
G4MuIonisation.hh
"
52
#include "
G4MuBremsstrahlung.hh
"
53
#include "
G4MuPairProduction.hh
"
54
55
#include "
G4hMultipleScattering.hh
"
56
#include "
G4hIonisation.hh
"
57
#include "
G4hBremsstrahlung.hh
"
58
#include "
G4hPairProduction.hh
"
59
60
#include "
G4ionIonisation.hh
"
61
#include "
G4IonParametrisedLossModel.hh
"
62
#include "
G4NuclearStopping.hh
"
63
64
#include "
G4EmProcessOptions.hh
"
65
#include "
G4MscStepLimitType.hh
"
66
67
#include "
G4PhysicalConstants.hh
"
68
#include "
G4SystemOfUnits.hh
"
69
70
71
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
73
PhysListEmStandardWVI::PhysListEmStandardWVI
(
const
G4String
&
name
)
74
:
G4VPhysicsConstructor
(name)
75
{}
76
77
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
79
PhysListEmStandardWVI::~PhysListEmStandardWVI
()
80
{}
81
82
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
84
void
PhysListEmStandardWVI::ConstructProcess
()
85
{
86
// Add standard EM Processes
87
//
88
89
aParticleIterator
->reset();
90
while
( (*
aParticleIterator
)() ){
91
G4ParticleDefinition
* particle =
aParticleIterator
->value();
92
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
93
G4String
particleName = particle->
GetParticleName
();
94
95
if
(particleName ==
"gamma"
) {
96
// gamma
97
pmanager->
AddDiscreteProcess
(
new
G4PhotoElectricEffect
);
98
pmanager->
AddDiscreteProcess
(
new
G4ComptonScattering
);
99
pmanager->
AddDiscreteProcess
(
new
G4GammaConversion
);
100
101
}
else
if
(particleName ==
"e-"
) {
102
//electron
103
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
104
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
105
pmanager->
AddProcess
(msc, -1, 1, 1);
106
pmanager->
AddProcess
(
new
G4eIonisation
, -1, 2, 2);
107
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1, 3, 3);
108
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 4);
109
110
}
else
if
(particleName ==
"e+"
) {
111
//positron
112
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
113
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
114
pmanager->
AddProcess
(msc, -1, 1, 1);
115
pmanager->
AddProcess
(
new
G4eIonisation
, -1, 2, 2);
116
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1, 3, 3);
117
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 4);
118
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 5);
119
120
}
else
if
(particleName ==
"mu+"
||
121
particleName ==
"mu-"
) {
122
//muon
123
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
124
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
125
pmanager->
AddProcess
(msc, -1, 1, 1);
126
pmanager->
AddProcess
(
new
G4MuIonisation
, -1, 2, 2);
127
pmanager->
AddProcess
(
new
G4MuBremsstrahlung
, -1, 3, 3);
128
pmanager->
AddProcess
(
new
G4MuPairProduction
, -1, 4, 4);
129
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 5);
130
131
}
else
if
( particleName ==
"proton"
||
132
particleName ==
"pi-"
||
133
particleName ==
"pi+"
) {
134
//proton
135
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
136
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
137
pmanager->
AddProcess
(msc, -1, 1, 1);
138
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
139
pmanager->
AddProcess
(
new
G4hBremsstrahlung
, -1, 3, 3);
140
pmanager->
AddProcess
(
new
G4hPairProduction
, -1, 4, 4);
141
142
}
else
if
( particleName ==
"alpha"
||
143
particleName ==
"He3"
) {
144
//alpha
145
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
146
pmanager->
AddProcess
(msc, -1, 1, 1);
147
pmanager->
AddProcess
(
new
G4ionIonisation
, -1, 2, 2);
148
pmanager->
AddProcess
(
new
G4NuclearStopping
, -1, 3,-1);
149
150
}
else
if
( particleName ==
"GenericIon"
) {
151
//Ions
152
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
153
pmanager->
AddProcess
(msc, -1, 1, 1);
154
G4ionIonisation
* ionIoni =
new
G4ionIonisation
();
155
ionIoni->
SetEmModel
(
new
G4IonParametrisedLossModel
());
156
pmanager->
AddProcess
(ionIoni, -1, 2, 2);
157
pmanager->
AddProcess
(
new
G4NuclearStopping
, -1, 3,-1);
158
159
}
else
if
((!particle->
IsShortLived
()) &&
160
(particle->
GetPDGCharge
() != 0.0) &&
161
(particle->
GetParticleName
() !=
"chargedgeantino"
)) {
162
//all others charged particles except geantino
163
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
164
pmanager->
AddProcess
(msc, -1, 1, 1);
165
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
166
}
167
}
168
169
// Em options
170
//
171
// Main options and setting parameters are shown here.
172
// Several of them have default values.
173
//
174
G4EmProcessOptions
emOptions;
175
176
//multiple coulomb scattering
177
//
178
emOptions.
SetPolarAngleLimit
(0.2);
179
}
180
181
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
G4hMultipleScattering.hh
G4eIonisation.hh
G4IonParametrisedLossModel.hh
G4MscStepLimitType.hh
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4EmProcessOptions
Definition:
G4EmProcessOptions.hh:63
G4ComptonScattering.hh
G4hBremsstrahlung.hh
G4hIonisation.hh
G4hPairProduction
Definition:
G4hPairProduction.hh:59
G4NuclearStopping.hh
G4MuPairProduction.hh
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4hPairProduction.hh
G4eMultipleScattering.hh
name
const XML_Char * name
Definition:
include/expat.h:151
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:
electromagnetic/TestEm12/src/PhysListEmStandardWVI.cc:73
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
G4IonCoulombScatteringModel.hh
G4NuclearStopping
Definition:
G4NuclearStopping.hh:66
G4ProcessManager.hh
G4ParticleDefinition.hh
G4PhotoElectricEffect.hh
G4ionIonisation
Definition:
G4ionIonisation.hh:78
G4MuBremsstrahlung
Definition:
G4MuBremsstrahlung.hh:78
G4PhysicalConstants.hh
G4MuIonisation
Definition:
G4MuIonisation.hh:85
PhysListEmStandardWVI::~PhysListEmStandardWVI
~PhysListEmStandardWVI()
Definition:
electromagnetic/TestEm12/src/PhysListEmStandardWVI.cc:79
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:144
G4VEnergyLossProcess::SetEmModel
void SetEmModel(G4VEmModel *, G4int index=1)
Definition:
G4VEnergyLossProcess.cc:362
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:
electromagnetic/TestEm12/src/PhysListEmStandardWVI.cc:84
G4SystemOfUnits.hh
G4ParticleDefinition::GetPDGCharge
G4double GetPDGCharge() const
Definition:
G4ParticleDefinition.hh:163
G4EmProcessOptions.hh
G4VPhysicsConstructor
Definition:
G4VPhysicsConstructor.hh:121
G4WentzelVIModel
Definition:
G4WentzelVIModel.hh:69
G4hIonisation
Definition:
G4hIonisation.hh:85
G4String
Definition:
examples/extended/parallel/TopC/ParN02/AnnotatedFiles/G4String.hh:45
G4EmProcessOptions::SetPolarAngleLimit
void SetPolarAngleLimit(G4double val)
Definition:
G4EmProcessOptions.cc:369
Generated on Wed Apr 30 2014 15:55:22 for Geant4.10 by
1.8.7