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
medical
fanoCavity2
src
/src/PhysListEmStandard_GS.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 medical/fanoCavity2/src/PhysListEmStandard_GS.cc
27
/// \brief Implementation of the PhysListEmStandard_GS class
28
//
29
// $Id: PhysListEmStandard_GS.cc 72961 2013-08-14 14:35:56Z gcosmo $
30
//
31
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34
#include "PhysListEmStandard_GS.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 "
G4GoudsmitSaundersonMscModel.hh
"
45
46
#include "
G4eIonisation.hh
"
47
#include "MyMollerBhabhaModel.hh"
48
#include "
G4eBremsstrahlung.hh
"
49
#include "
G4eplusAnnihilation.hh
"
50
51
#include "
G4hIonisation.hh
"
52
#include "
G4hMultipleScattering.hh
"
53
54
#include "
G4EmProcessOptions.hh
"
55
#include "
G4MscStepLimitType.hh
"
56
57
#include "
G4SystemOfUnits.hh
"
58
59
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61
PhysListEmStandard_GS::PhysListEmStandard_GS
(
const
G4String
&
name
)
62
:
G4VPhysicsConstructor
(name)
63
{}
64
65
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
67
PhysListEmStandard_GS::~PhysListEmStandard_GS
()
68
{}
69
70
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
72
void
PhysListEmStandard_GS::ConstructProcess
()
73
{
74
// Add standard EM Processes
75
//
76
77
aParticleIterator
->reset();
78
while
( (*
aParticleIterator
)() ){
79
G4ParticleDefinition
* particle =
aParticleIterator
->value();
80
G4ProcessManager
* pmanager = particle->
GetProcessManager
();
81
G4String
particleName = particle->
GetParticleName
();
82
83
if
(particleName ==
"gamma"
) {
84
// gamma
85
86
pmanager->
AddDiscreteProcess
(
new
G4PhotoElectricEffect
);
87
pmanager->
AddDiscreteProcess
(
new
G4ComptonScattering
);
88
pmanager->
AddDiscreteProcess
(
new
G4GammaConversion
);
89
90
}
else
if
(particleName ==
"e-"
) {
91
//electron
92
93
G4eMultipleScattering
* eMsc =
new
G4eMultipleScattering
();
94
eMsc->
AddEmModel
(1,
new
G4GoudsmitSaundersonMscModel
);
95
96
G4eIonisation
* eIoni =
new
G4eIonisation
();
97
eIoni->
SetEmModel
(
new
MyMollerBhabhaModel
);
98
99
pmanager->
AddProcess
(eMsc, -1, 1, 1);
100
pmanager->
AddProcess
(eIoni, -1, 2, 2);
101
/// pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);
102
103
}
else
if
(particleName ==
"e+"
) {
104
//positron
105
106
G4eMultipleScattering
* pMsc =
new
G4eMultipleScattering
();
107
pMsc->
AddEmModel
(1,
new
G4GoudsmitSaundersonMscModel
);
108
109
G4eIonisation
* pIoni =
new
G4eIonisation
();
110
pIoni->
SetEmModel
(
new
MyMollerBhabhaModel
);
111
112
pmanager->
AddProcess
(pMsc, -1, 1, 1);
113
pmanager->
AddProcess
(pIoni, -1, 2, 2);
114
/// pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);
115
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 3);
116
117
}
else
if
( particleName ==
"proton"
) {
118
//proton
119
pmanager->
AddProcess
(
new
G4hMultipleScattering
, -1, 1, 1);
120
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
121
}
122
}
123
124
// Em options
125
//
126
// Main options and setting parameters are shown here.
127
// Several of them have default values.
128
//
129
G4EmProcessOptions
emOptions;
130
131
//physics tables
132
//
133
emOptions.
SetMinEnergy
(100*
eV
);
//default
134
emOptions.
SetMaxEnergy
(10*
GeV
);
//default
135
emOptions.
SetDEDXBinning
(8*20);
//default=8*7
136
emOptions.
SetLambdaBinning
(8*20);
//default=8*7
137
138
//multiple coulomb scattering
139
//
140
emOptions.
SetMscStepLimitation
(
fUseDistanceToBoundary
);
//default=fUseSafety
141
142
//energy loss
143
//
144
emOptions.
SetStepFunction
(0.2, 10*um);
//default=(0.2, 1*mm)
145
146
//build CSDA range
147
//
148
emOptions.
SetBuildCSDARange
(
true
);
//default=false
149
emOptions.
SetMaxEnergyForCSDARange
(10*
GeV
);
150
emOptions.
SetDEDXBinningForCSDARange
(8*20);
//default=8*7
151
}
152
153
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
154
G4hMultipleScattering.hh
G4eIonisation.hh
python.hepunit.GeV
GeV
Definition:
hepunit.py:120
G4MscStepLimitType.hh
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4GoudsmitSaundersonMscModel
Definition:
G4GoudsmitSaundersonMscModel.hh:73
G4GoudsmitSaundersonMscModel.hh
G4EmProcessOptions
Definition:
G4EmProcessOptions.hh:63
G4ComptonScattering.hh
PhysListEmStandard_GS::ConstructProcess
virtual void ConstructProcess()
Definition:
src/PhysListEmStandard_GS.cc:75
G4hIonisation.hh
G4EmProcessOptions::SetMinEnergy
void SetMinEnergy(G4double val)
Definition:
G4EmProcessOptions.cc:109
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4eMultipleScattering.hh
name
const XML_Char * name
Definition:
include/expat.h:151
G4EmProcessOptions::SetStepFunction
void SetStepFunction(G4double v1, G4double v2)
Definition:
G4EmProcessOptions.cc:158
G4ComptonScattering
Definition:
G4ComptonScattering.hh:71
fUseDistanceToBoundary
Definition:
G4MscStepLimitType.hh:51
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
MyMollerBhabhaModel
Definition:
include/MyMollerBhabhaModel.hh:41
G4EmProcessOptions::SetDEDXBinning
void SetDEDXBinning(G4int val)
Definition:
G4EmProcessOptions.cc:137
G4EmProcessOptions::SetLambdaBinning
void SetLambdaBinning(G4int val)
Definition:
G4EmProcessOptions.cc:151
G4GammaConversion
Definition:
G4GammaConversion.hh:75
G4GammaConversion.hh
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
python.hepunit.eV
eV
Definition:
hepunit.py:118
G4ProcessManager.hh
G4ParticleDefinition.hh
G4PhotoElectricEffect.hh
G4EmProcessOptions::SetMaxEnergy
void SetMaxEnergy(G4double val)
Definition:
G4EmProcessOptions.cc:116
PhysListEmStandard_GS::PhysListEmStandard_GS
PhysListEmStandard_GS(const G4String &name, DetectorConstruction *det)
Definition:
src/PhysListEmStandard_GS.cc:63
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
G4EmProcessOptions::SetBuildCSDARange
void SetBuildCSDARange(G4bool val)
Definition:
G4EmProcessOptions.cc:185
G4eBremsstrahlung.hh
G4ProcessManager
Definition:
G4ProcessManager.hh:106
G4eplusAnnihilation
Definition:
G4eplusAnnihilation.hh:65
G4eMultipleScattering
Definition:
G4eMultipleScattering.hh:60
G4hMultipleScattering
Definition:
G4hMultipleScattering.hh:62
G4EmProcessOptions::SetMscStepLimitation
void SetMscStepLimitation(G4MscStepLimitType val)
Definition:
G4EmProcessOptions.cc:306
G4SystemOfUnits.hh
G4EmProcessOptions.hh
G4VPhysicsConstructor
Definition:
G4VPhysicsConstructor.hh:121
G4hIonisation
Definition:
G4hIonisation.hh:85
PhysListEmStandard_GS::~PhysListEmStandard_GS
~PhysListEmStandard_GS()
Definition:
src/PhysListEmStandard_GS.cc:70
G4String
Definition:
examples/extended/parallel/TopC/ParN02/AnnotatedFiles/G4String.hh:45
Generated on Wed Apr 30 2014 15:55:21 for Geant4.10 by
1.8.7