Geant4
10.00.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
PhysListEmStandard_WVI.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: PhysListEmStandard_WVI.cc 72961 2013-08-14 14:35:56Z gcosmo $
30
//
31
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34
#include "PhysListEmStandard_WVI.hh"
35
36
#include "
G4ParticleDefinition.hh
"
37
#include "
G4ProcessManager.hh
"
38
39
#include "
G4ComptonScattering.hh
"
40
#include "
G4GammaConversion.hh
"
41
#include "
G4PhotoElectricEffect.hh
"
42
43
#include "
G4eMultipleScattering.hh
"
44
#include "
G4MuMultipleScattering.hh
"
45
#include "
G4WentzelVIModel.hh
"
46
#include "
G4CoulombScattering.hh
"
47
48
#include "
G4eIonisation.hh
"
49
#include "MyMollerBhabhaModel.hh"
50
#include "
G4eBremsstrahlung.hh
"
51
#include "
G4eplusAnnihilation.hh
"
52
53
#include "
G4hIonisation.hh
"
54
#include "
G4hMultipleScattering.hh
"
55
56
#include "
G4EmProcessOptions.hh
"
57
#include "
G4MscStepLimitType.hh
"
58
59
#include "
G4SystemOfUnits.hh
"
60
61
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
63
PhysListEmStandard_WVI::PhysListEmStandard_WVI
(
const
G4String
&
name
)
64
:
G4VPhysicsConstructor
(name)
65
{}
66
67
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
69
PhysListEmStandard_WVI::~PhysListEmStandard_WVI
()
70
{}
71
72
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74
void
PhysListEmStandard_WVI::ConstructProcess
()
75
{
76
// Add standard EM Processes
77
//
78
79
aParticleIterator
->reset();
80
while
( (*
aParticleIterator
)() ){
81
G4ParticleDefinition
* particle =
aParticleIterator
->value();
82
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
83
G4String
particleName = particle->
GetParticleName
();
84
85
if
(particleName ==
"gamma"
) {
86
// gamma
87
88
pmanager->
AddDiscreteProcess
(
new
G4PhotoElectricEffect
);
89
pmanager->
AddDiscreteProcess
(
new
G4ComptonScattering
);
90
pmanager->
AddDiscreteProcess
(
new
G4GammaConversion
);
91
92
}
else
if
(particleName ==
"e-"
) {
93
//electron
94
95
G4MuMultipleScattering
* eMsc =
new
G4MuMultipleScattering
();
96
eMsc->
AddEmModel
(1,
new
G4WentzelVIModel
());
97
G4eIonisation
* eIoni =
new
G4eIonisation
();
98
eIoni->
SetEmModel
(
new
MyMollerBhabhaModel
);
99
100
pmanager->
AddProcess
(eMsc, -1, 1, 1);
101
pmanager->
AddProcess
(eIoni, -1, 2, 2);
103
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 3);
104
105
}
else
if
(particleName ==
"e+"
) {
106
//positron
107
108
G4MuMultipleScattering
* pMsc =
new
G4MuMultipleScattering
();
109
pMsc->
AddEmModel
(1,
new
G4WentzelVIModel
());
110
G4eIonisation
* pIoni =
new
G4eIonisation
();
111
pIoni->
SetEmModel
(
new
MyMollerBhabhaModel
);
112
113
pmanager->
AddProcess
(pMsc, -1, 1, 1);
114
pmanager->
AddProcess
(pIoni, -1, 2, 2);
116
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 3);
117
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 4);
118
119
}
else
if
( particleName ==
"proton"
) {
120
//proton
121
G4hMultipleScattering
* msc =
new
G4hMultipleScattering
();
122
msc->
AddEmModel
(1,
new
G4WentzelVIModel
());
123
pmanager->
AddProcess
(msc, -1, 1, 1);
124
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
125
pmanager->
AddProcess
(
new
G4CoulombScattering
, -1,-1, 3);
126
}
127
}
128
129
// Em options
130
//
131
// Main options and setting parameters are shown here.
132
// Several of them have default values.
133
//
134
G4EmProcessOptions
emOptions;
135
//physics tables
136
//
137
emOptions.
SetMinEnergy
(100*
eV
);
//default
138
emOptions.
SetMaxEnergy
(10*
GeV
);
//default
139
emOptions.
SetDEDXBinning
(8*20);
//default=8*7
140
emOptions.
SetLambdaBinning
(8*20);
//default=8*7
141
142
//multiple coulomb scattering
143
//
144
emOptions.
SetPolarAngleLimit
(0.2);
145
146
//energy loss
147
//
148
emOptions.
SetStepFunction
(0.2, 10*um);
//default=(0.2, 1*mm)
149
150
//build CSDA range
151
//
152
emOptions.
SetBuildCSDARange
(
true
);
//default=false
153
emOptions.
SetMaxEnergyForCSDARange
(10*
GeV
);
154
emOptions.
SetDEDXBinningForCSDARange
(8*7);
//default=8*7
155
156
//ionization
157
//
158
emOptions.
SetSubCutoff
(
false
);
//default
159
}
160
161
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162
G4hMultipleScattering.hh
G4eIonisation.hh
G4MscStepLimitType.hh
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4EmProcessOptions
Definition:
G4EmProcessOptions.hh:63
G4ComptonScattering.hh
G4hIonisation.hh
name
G4String name
Definition:
TRTMaterials.hh:40
G4EmProcessOptions::SetMinEnergy
void SetMinEnergy(G4double val)
Definition:
G4EmProcessOptions.cc:109
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4eMultipleScattering.hh
G4EmProcessOptions::SetStepFunction
void SetStepFunction(G4double v1, G4double v2)
Definition:
G4EmProcessOptions.cc:158
G4ComptonScattering
Definition:
G4ComptonScattering.hh:71
G4EmProcessOptions::SetDEDXBinningForCSDARange
void SetDEDXBinningForCSDARange(G4int val)
Definition:
G4EmProcessOptions.cc:144
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:111
G4ParticleDefinition::GetProcessManager
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition::GetParticleName
const G4String & GetParticleName() const
Definition:
G4ParticleDefinition.hh:159
G4EmProcessOptions::SetMaxEnergyForCSDARange
void SetMaxEnergyForCSDARange(G4double val)
Definition:
G4EmProcessOptions.cc:123
G4CoulombScattering
Definition:
G4CoulombScattering.hh:55
MyMollerBhabhaModel
Definition:
MyMollerBhabhaModel.hh:41
G4EmProcessOptions::SetDEDXBinning
void SetDEDXBinning(G4int val)
Definition:
G4EmProcessOptions.cc:137
G4CoulombScattering.hh
G4EmProcessOptions::SetLambdaBinning
void SetLambdaBinning(G4int val)
Definition:
G4EmProcessOptions.cc:151
G4GammaConversion
Definition:
G4GammaConversion.hh:75
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
GeV
static const double GeV
Definition:
G4SIunits.hh:196
G4ProcessManager.hh
G4ParticleDefinition.hh
G4PhotoElectricEffect.hh
G4EmProcessOptions::SetMaxEnergy
void SetMaxEnergy(G4double val)
Definition:
G4EmProcessOptions.cc:116
PhysListEmStandard_WVI::~PhysListEmStandard_WVI
~PhysListEmStandard_WVI()
Definition:
PhysListEmStandard_WVI.cc:72
eV
static const double eV
Definition:
G4SIunits.hh:194
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:145
G4VEnergyLossProcess::SetEmModel
void SetEmModel(G4VEmModel *, G4int index=1)
Definition:
G4VEnergyLossProcess.cc:389
PhysListEmStandard_WVI::PhysListEmStandard_WVI
PhysListEmStandard_WVI(const G4String &name, DetectorConstruction *det)
Definition:
PhysListEmStandard_WVI.cc:65
G4EmProcessOptions::SetBuildCSDARange
void SetBuildCSDARange(G4bool val)
Definition:
G4EmProcessOptions.cc:185
G4eBremsstrahlung.hh
G4ProcessManager
Definition:
G4ProcessManager.hh:106
G4eplusAnnihilation
Definition:
G4eplusAnnihilation.hh:65
G4hMultipleScattering
Definition:
G4hMultipleScattering.hh:62
G4SystemOfUnits.hh
PhysListEmStandard_WVI::ConstructProcess
virtual void ConstructProcess()
Definition:
PhysListEmStandard_WVI.cc:77
G4EmProcessOptions::SetSubCutoff
void SetSubCutoff(G4bool val, const G4Region *r=0)
Definition:
G4EmProcessOptions.cc:88
G4EmProcessOptions.hh
G4VPhysicsConstructor
Definition:
G4VPhysicsConstructor.hh:121
G4WentzelVIModel
Definition:
G4WentzelVIModel.hh:69
G4hIonisation
Definition:
G4hIonisation.hh:85
G4String
Definition:
G4String.hh:45
G4EmProcessOptions::SetPolarAngleLimit
void SetPolarAngleLimit(G4double val)
Definition:
G4EmProcessOptions.cc:369
geant4.10.00.p02
examples
extended
medical
fanoCavity2
src
PhysListEmStandard_WVI.cc
Generated on Thu Dec 31 2015 10:39:46 for Geant4 by
1.8.8