Geant4
v4-10.4-release
메인 페이지
관련된 페이지
모듈
네임스페이스
클래스
파일들
파일 목록
파일 멤버
모두
클래스
네임스페이스들
파일들
함수
변수
타입정의
열거형 타입
열거형 멤버
Friends
매크로
그룹들
페이지들
examples
extended
medical
electronScattering
src
PhysListEmStandardWVI.cc
이 파일의 문서화 페이지로 가기
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 100262 2016-10-17 08:08:04Z gcosmo $
30
//
31
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34
#include "
PhysListEmStandardWVI.hh
"
35
36
#include "
G4BuilderType.hh
"
37
#include "
G4ParticleDefinition.hh
"
38
#include "
G4ProcessManager.hh
"
39
40
#include "
G4ComptonScattering.hh
"
41
#include "
G4GammaConversion.hh
"
42
#include "
G4PhotoElectricEffect.hh
"
43
44
#include "
G4eMultipleScattering.hh
"
45
#include "
G4WentzelVIModel.hh
"
46
#include "
G4CoulombScattering.hh
"
47
#include "
G4eIonisation.hh
"
48
#include "
G4eBremsstrahlung.hh
"
49
#include "
G4eplusAnnihilation.hh
"
50
51
#include "
G4MuMultipleScattering.hh
"
52
#include "
G4MuIonisation.hh
"
53
#include "
G4MuBremsstrahlung.hh
"
54
#include "
G4MuPairProduction.hh
"
55
56
#include "
G4hMultipleScattering.hh
"
57
#include "
G4hIonisation.hh
"
58
#include "
G4hBremsstrahlung.hh
"
59
#include "
G4hPairProduction.hh
"
60
61
#include "
G4ionIonisation.hh
"
62
#include "
G4IonParametrisedLossModel.hh
"
63
#include "
G4NuclearStopping.hh
"
64
65
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
67
PhysListEmStandardWVI::PhysListEmStandardWVI
(
const
G4String
&
name
)
68
:
G4VPhysicsConstructor
(name)
69
{
70
SetPhysicsType
(
bElectromagnetic
);
71
72
G4EmParameters
* param =
G4EmParameters::Instance
();
73
param->
SetDefaults
();
74
param->
SetMscThetaLimit
(0.2);
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
auto
particleIterator
=
GetParticleIterator
();
90
particleIterator
->reset();
91
while
( (*
particleIterator
)() ){
92
G4ParticleDefinition
* particle =
particleIterator
->value();
93
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
94
G4String
particleName = particle->
GetParticleName
();
95
96
if
(particleName ==
"gamma"
) {
97
// gamma
98
pmanager->
AddDiscreteProcess
(
new
G4PhotoElectricEffect
);
99
pmanager->
AddDiscreteProcess
(
new
G4ComptonScattering
);
100
pmanager->
AddDiscreteProcess
(
new
G4GammaConversion
);
101
102
}
else
if
(particleName ==
"e-"
) {
103
//electron
104
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
105
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
106
pmanager->
AddProcess
(msc, -1, 1, 1);
107
pmanager->
AddProcess
(
new
G4eIonisation
, -1, 2, 2);
108
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1, 3, 3);
109
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 4);
110
111
}
else
if
(particleName ==
"e+"
) {
112
//positron
113
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
114
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
115
pmanager->
AddProcess
(msc, -1, 1, 1);
116
pmanager->
AddProcess
(
new
G4eIonisation
, -1, 2, 2);
117
pmanager->
AddProcess
(
new
G4eBremsstrahlung
, -1, 3, 3);
118
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 4);
119
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 5);
120
121
}
else
if
(particleName ==
"mu+"
||
122
particleName ==
"mu-"
) {
123
//muon
124
G4MuMultipleScattering
* msc =
new
G4MuMultipleScattering
();
125
msc->
AddEmModel
(0,
new
G4WentzelVIModel
());
126
pmanager->
AddProcess
(msc, -1, 1, 1);
127
pmanager->
AddProcess
(
new
G4MuIonisation
, -1, 2, 2);
128
pmanager->
AddProcess
(
new
G4MuBremsstrahlung
, -1, 3, 3);
129
pmanager->
AddProcess
(
new
G4MuPairProduction
, -1, 4, 4);
130
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 5);
131
132
}
else
if
( particleName ==
"proton"
||
133
particleName ==
"pi-"
||
134
particleName ==
"pi+"
) {
135
//proton
136
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
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
170
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171
G4hMultipleScattering
Definition:
G4hMultipleScattering.hh:62
name
const XML_Char * name
Definition:
expat.h:151
G4ionIonisation
Definition:
G4ionIonisation.hh:78
G4MuBremsstrahlung.hh
G4MuIonisation
Definition:
G4MuIonisation.hh:85
G4PhotoElectricEffect.hh
G4VEnergyLossProcess::SetEmModel
void SetEmModel(G4VEmModel *, G4int index=0)
Definition:
G4VEnergyLossProcess.cc:409
G4hMultipleScattering.hh
G4eplusAnnihilation.hh
G4CoulombScattering.hh
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4IonParametrisedLossModel.hh
G4ParticleDefinition::GetParticleName
const G4String & GetParticleName() const
Definition:
G4ParticleDefinition.hh:121
G4String
Definition:
examples/extended/parallel/TopC/ParN02/AnnotatedFiles/G4String.hh:45
G4WentzelVIModel
Definition:
G4WentzelVIModel.hh:66
G4ParticleDefinition::GetPDGCharge
G4double GetPDGCharge() const
Definition:
G4ParticleDefinition.hh:125
bElectromagnetic
Definition:
G4BuilderType.hh:48
G4hPairProduction.hh
G4MuMultipleScattering
Definition:
G4MuMultipleScattering.hh:61
G4ComptonScattering
Definition:
G4ComptonScattering.hh:71
G4GammaConversion.hh
G4hIonisation
Definition:
G4hIonisation.hh:85
PhysListEmStandardWVI::PhysListEmStandardWVI
PhysListEmStandardWVI(const G4String &name="standardWVI")
Definition:
PhysListEmStandardWVI.cc:67
PhysListEmStandardWVI.hh
Definition of the PhysListEmStandardWVI class.
G4ProcessManager::AddProcess
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
Definition:
G4ProcessManager.cc:410
PhysListEmStandardWVI::ConstructProcess
virtual void ConstructProcess()
Definition:
PhysListEmStandardWVI.cc:84
G4EmParameters::SetMscThetaLimit
void SetMscThetaLimit(G4double val)
Definition:
G4EmParameters.cc:672
G4ProcessManager.hh
G4EmParameters
Definition:
G4EmParameters.hh:73
G4eplusAnnihilation
Definition:
G4eplusAnnihilation.hh:65
G4eMultipleScattering.hh
G4eBremsstrahlung.hh
G4eIonisation
Definition:
G4eIonisation.hh:80
G4MuPairProduction
Definition:
G4MuPairProduction.hh:74
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:73
PhysListEmStandardWVI::~PhysListEmStandardWVI
~PhysListEmStandardWVI()
Definition:
PhysListEmStandardWVI.cc:79
G4CoulombScattering
Definition:
G4CoulombScattering.hh:55
G4VPhysicsConstructor
Definition:
G4VPhysicsConstructor.hh:126
G4hPairProduction
Definition:
G4hPairProduction.hh:59
G4VPhysicsConstructor::SetPhysicsType
void SetPhysicsType(G4int)
Definition:
G4VPhysicsConstructor.hh:211
G4WentzelVIModel.hh
G4ionIonisation.hh
G4MuBremsstrahlung
Definition:
G4MuBremsstrahlung.hh:78
G4hBremsstrahlung
Definition:
G4hBremsstrahlung.hh:61
G4ParticleDefinition::GetProcessManager
G4ProcessManager * GetProcessManager() const
Definition:
G4ParticleDefinition.cc:258
G4NuclearStopping.hh
G4hBremsstrahlung.hh
G4GammaConversion
Definition:
G4GammaConversion.hh:75
particleIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition:
G4BigBanger.cc:65
G4eBremsstrahlung
Definition:
G4eBremsstrahlung.hh:81
G4VPhysicsConstructor::GetParticleIterator
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
Definition:
G4VPhysicsConstructor.cc:82
G4ComptonScattering.hh
G4MuPairProduction.hh
G4ParticleDefinition.hh
G4VMultipleScattering::AddEmModel
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
Definition:
G4VMultipleScattering.cc:144
G4NuclearStopping
Definition:
G4NuclearStopping.hh:64
G4MuIonisation.hh
G4ParticleDefinition::IsShortLived
G4bool IsShortLived() const
Definition:
G4ParticleDefinition.hh:158
G4hIonisation.hh
G4eIonisation.hh
G4IonParametrisedLossModel
Definition:
G4IonParametrisedLossModel.hh:93
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4ProcessManager
Definition:
G4ProcessManager.hh:106
G4EmParameters::Instance
static G4EmParameters * Instance()
Definition:
G4EmParameters.cc:70
G4EmParameters::SetDefaults
void SetDefaults()
Definition:
G4EmParameters.cc:107
G4BuilderType.hh
G4MuMultipleScattering.hh
다음에 의해 생성됨 :
1.8.5