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