Geant4
10.00.p01
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
PhysListEmLivermore.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
//
30
// $Id: PhysListEmLivermore.cc 68585 2013-04-01 23:35:07Z adotti $
31
//
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
34
35
#include "PhysListEmLivermore.hh"
36
#include "
G4ParticleDefinition.hh
"
37
#include "
G4ProcessManager.hh
"
38
39
// gamma
40
41
#include "
G4PhotoElectricEffect.hh
"
42
#include "
G4LivermorePhotoElectricModel.hh
"
43
44
#include "
G4ComptonScattering.hh
"
45
#include "
G4LivermoreComptonModel.hh
"
46
47
#include "
G4GammaConversion.hh
"
48
#include "
G4LivermoreGammaConversionModel.hh
"
49
50
#include "
G4RayleighScattering.hh
"
51
#include "
G4LivermoreRayleighModel.hh
"
52
53
// e-
54
55
#include "
G4eIonisation.hh
"
56
#include "
G4LivermoreIonisationModel.hh
"
57
#include "
G4UniversalFluctuation.hh
"
58
59
#include "
G4eBremsstrahlung.hh
"
60
#include "
G4LivermoreBremsstrahlungModel.hh
"
61
62
// e+
63
64
#include "
G4eplusAnnihilation.hh
"
65
66
// mu
67
68
#include "
G4MuIonisation.hh
"
69
#include "
G4MuBremsstrahlung.hh
"
70
#include "
G4MuPairProduction.hh
"
71
72
// hadrons, ions
73
74
#include "
G4hIonisation.hh
"
75
#include "
G4ionIonisation.hh
"
76
77
#include "
G4EmProcessOptions.hh
"
78
#include "
G4LossTableManager.hh
"
79
#include "
G4UAtomicDeexcitation.hh
"
80
81
#include "
G4SystemOfUnits.hh
"
82
83
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
85
PhysListEmLivermore::PhysListEmLivermore
(
const
G4String
&
name
)
86
:
G4VPhysicsConstructor
(name)
87
{ }
88
89
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
91
PhysListEmLivermore::~PhysListEmLivermore
()
92
{ }
93
94
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
96
void
PhysListEmLivermore::ConstructProcess
()
97
{
98
// Add standard EM Processes
99
100
aParticleIterator
->reset();
101
while
( (*
aParticleIterator
)() ){
102
G4ParticleDefinition
* particle =
aParticleIterator
->value();
103
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
104
G4String
particleName = particle->
GetParticleName
();
105
106
//Applicability range for Livermore models
107
//for higher energies, the Standard models are used
108
G4double
highEnergyLimit = 1*
GeV
;
109
110
if
(particleName ==
"gamma"
) {
111
// gamma
112
113
G4PhotoElectricEffect
* phot =
new
G4PhotoElectricEffect
();
114
G4LivermorePhotoElectricModel
*
115
photModel =
new
G4LivermorePhotoElectricModel
();
116
photModel->
SetHighEnergyLimit
(highEnergyLimit);
117
phot->
AddEmModel
(0, photModel);
118
pmanager->
AddDiscreteProcess
(phot);
119
120
G4ComptonScattering
* compt =
new
G4ComptonScattering
();
121
G4LivermoreComptonModel
*
122
comptModel =
new
G4LivermoreComptonModel
();
123
comptModel->
SetHighEnergyLimit
(highEnergyLimit);
124
compt->
AddEmModel
(0, comptModel);
125
pmanager->
AddDiscreteProcess
(compt);
126
127
G4GammaConversion
* conv =
new
G4GammaConversion
();
128
G4LivermoreGammaConversionModel
*
129
convModel =
new
G4LivermoreGammaConversionModel
();
130
convModel->
SetHighEnergyLimit
(highEnergyLimit);
131
conv->
AddEmModel
(0, convModel);
132
pmanager->
AddDiscreteProcess
(conv);
133
134
G4RayleighScattering
* rayl =
new
G4RayleighScattering
();
135
G4LivermoreRayleighModel
*
136
raylModel =
new
G4LivermoreRayleighModel
();
137
raylModel->
SetHighEnergyLimit
(highEnergyLimit);
138
rayl->
AddEmModel
(0, raylModel);
139
pmanager->
AddDiscreteProcess
(rayl);
140
141
}
else
if
(particleName ==
"e-"
) {
142
//electron
143
144
G4eIonisation
* eIoni =
new
G4eIonisation
();
145
G4LivermoreIonisationModel
*
146
eIoniModel =
new
G4LivermoreIonisationModel
();
147
eIoniModel->
SetHighEnergyLimit
(highEnergyLimit);
148
eIoni->
AddEmModel
(0, eIoniModel,
new
G4UniversalFluctuation
() );
149
pmanager->
AddProcess
(eIoni, -1,-1, 1);
150
151
G4eBremsstrahlung
* eBrem =
new
G4eBremsstrahlung
();
152
G4LivermoreBremsstrahlungModel
*
153
eBremModel =
new
G4LivermoreBremsstrahlungModel
();
154
eBremModel->
SetHighEnergyLimit
(highEnergyLimit);
155
eBrem->
AddEmModel
(0, eBremModel);
156
pmanager->
AddProcess
(eBrem, -1,-1, 2);
157
158
}
else
if
(particleName ==
"e+"
) {
159
//positron
160
pmanager->
AddProcess
(
new
G4eIonisation
, -1,-1, 1);
161
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1,-1, 2);
162
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 3);
163
164
}
else
if
( particleName ==
"mu+"
||
165
particleName ==
"mu-"
) {
166
//muon
167
pmanager->
AddProcess
(
new
G4MuIonisation
, -1,-1, 1);
168
pmanager->
AddProcess
(
new
G4MuBremsstrahlung
, -1,-1, 2);
169
pmanager->
AddProcess
(
new
G4MuPairProduction
, -1,-1, 3);
170
171
}
else
if
( particleName ==
"alpha"
|| particleName ==
"GenericIon"
) {
172
pmanager->
AddProcess
(
new
G4ionIonisation
, -1,-1, 1);
173
174
}
else
if
((!particle->
IsShortLived
()) &&
175
(particle->
GetPDGCharge
() != 0.0) &&
176
(particle->
GetParticleName
() !=
"chargedgeantino"
)) {
177
//all others charged particles except geantino
178
pmanager->
AddProcess
(
new
G4hIonisation
, -1,-1, 1);
179
}
180
}
181
182
// Deexcitation
183
//
184
G4VAtomDeexcitation
* de =
new
G4UAtomicDeexcitation
();
185
de->
SetFluo
(
true
);
186
de->
SetAuger
(
false
);
187
de->
SetPIXE
(
false
);
188
G4LossTableManager::Instance
()->
SetAtomDeexcitation
(de);
189
}
190
191
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192
G4LivermoreBremsstrahlungModel
Definition:
G4LivermoreBremsstrahlungModel.hh:61
G4eIonisation.hh
G4LossTableManager::Instance
static G4LossTableManager * Instance()
Definition:
G4LossTableManager.cc:112
G4UniversalFluctuation.hh
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4VAtomDeexcitation::SetAuger
void SetAuger(G4bool)
Definition:
G4VAtomDeexcitation.hh:208
G4ComptonScattering.hh
G4hIonisation.hh
name
G4String name
Definition:
TRTMaterials.hh:40
G4LivermoreIonisationModel.hh
PhysListEmLivermore::ConstructProcess
virtual void ConstructProcess()
Definition:
PhysListEmLivermore.cc:92
G4MuPairProduction.hh
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4UAtomicDeexcitation
Definition:
G4UAtomicDeexcitation.hh:61
G4ComptonScattering
Definition:
G4ComptonScattering.hh:71
G4LossTableManager.hh
G4MuIonisation.hh
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:111
G4MuBremsstrahlung.hh
G4LivermoreComptonModel.hh
G4ParticleDefinition::GetProcessManager
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition::GetParticleName
const G4String & GetParticleName() const
Definition:
G4ParticleDefinition.hh:159
G4VEmModel::SetHighEnergyLimit
void SetHighEnergyLimit(G4double)
Definition:
G4VEmModel.hh:683
G4LivermoreRayleighModel
Definition:
G4LivermoreRayleighModel.hh:39
G4ionIonisation.hh
G4LivermoreIonisationModel
Definition:
G4LivermoreIonisationModel.hh:58
G4GammaConversion
Definition:
G4GammaConversion.hh:75
G4GammaConversion.hh
aParticleIterator
#define aParticleIterator
Definition:
G4VPhysicsConstructor.hh:119
G4VAtomDeexcitation
Definition:
G4VAtomDeexcitation.hh:63
G4ProcessManager::AddProcess
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
Definition:
G4ProcessManager.cc:410
G4LivermoreBremsstrahlungModel.hh
GeV
static const double GeV
Definition:
G4SIunits.hh:196
G4ProcessManager.hh
G4ParticleDefinition.hh
G4LivermoreComptonModel
Definition:
G4LivermoreComptonModel.hh:48
G4VEnergyLossProcess::AddEmModel
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
Definition:
G4VEnergyLossProcess.cc:344
G4PhotoElectricEffect.hh
G4ionIonisation
Definition:
G4ionIonisation.hh:78
G4MuBremsstrahlung
Definition:
G4MuBremsstrahlung.hh:78
G4LivermoreGammaConversionModel.hh
G4UAtomicDeexcitation.hh
G4MuIonisation
Definition:
G4MuIonisation.hh:85
G4LivermoreRayleighModel.hh
G4ParticleDefinition::IsShortLived
G4bool IsShortLived() const
Definition:
G4ParticleDefinition.hh:196
G4RayleighScattering.hh
G4eBremsstrahlung
Definition:
G4eBremsstrahlung.hh:81
G4eIonisation
Definition:
G4eIonisation.hh:80
G4VEmProcess::AddEmModel
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
Definition:
G4VEmProcess.cc:188
G4eplusAnnihilation.hh
G4LivermorePhotoElectricModel
Definition:
G4LivermorePhotoElectricModel.hh:50
G4MuPairProduction
Definition:
G4MuPairProduction.hh:74
G4UniversalFluctuation
Definition:
G4UniversalFluctuation.hh:62
G4eBremsstrahlung.hh
G4ProcessManager
Definition:
G4ProcessManager.hh:106
G4eplusAnnihilation
Definition:
G4eplusAnnihilation.hh:65
G4RayleighScattering
Definition:
G4RayleighScattering.hh:51
G4VAtomDeexcitation::SetPIXE
void SetPIXE(G4bool)
Definition:
G4VAtomDeexcitation.hh:219
G4LivermoreGammaConversionModel
Definition:
G4LivermoreGammaConversionModel.hh:41
G4double
double G4double
Definition:
G4Types.hh:76
G4SystemOfUnits.hh
G4ParticleDefinition::GetPDGCharge
G4double GetPDGCharge() const
Definition:
G4ParticleDefinition.hh:163
G4EmProcessOptions.hh
PhysListEmLivermore::~PhysListEmLivermore
~PhysListEmLivermore()
Definition:
PhysListEmLivermore.cc:87
G4LossTableManager::SetAtomDeexcitation
void SetAtomDeexcitation(G4VAtomDeexcitation *)
Definition:
G4LossTableManager.cc:1314
G4VPhysicsConstructor
Definition:
G4VPhysicsConstructor.hh:121
G4hIonisation
Definition:
G4hIonisation.hh:85
G4VAtomDeexcitation::SetFluo
void SetFluo(G4bool)
Definition:
G4VAtomDeexcitation.hh:198
G4LivermorePhotoElectricModel.hh
PhysListEmLivermore::PhysListEmLivermore
PhysListEmLivermore(const G4String &name="livermore")
Definition:
PhysListEmLivermore.cc:81
G4String
Definition:
G4String.hh:45
geant4.10.00.p01
examples
extended
electromagnetic
TestEm14
src
PhysListEmLivermore.cc
Generated on Thu Dec 31 2015 17:39:53 for Geant4 by
1.8.8