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
fanoCavity
src
src/PhysListEmStandard_option0.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/fanoCavity/src/PhysListEmStandard_option0.cc
27
/// \brief Implementation of the PhysListEmStandard_option0 class
28
//
29
// $Id: PhysListEmStandard_option0.cc 68525 2013-04-01 21:16:50Z adotti $
30
//
31
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34
#include "PhysListEmStandard_option0.hh"
35
#include "DetectorConstruction.hh"
36
37
#include "
G4ParticleDefinition.hh
"
38
#include "
G4ProcessManager.hh
"
39
40
#include "
G4ComptonScattering.hh
"
41
#include "
MyKleinNishinaCompton.hh
"
42
#include "
G4GammaConversion.hh
"
43
#include "
G4PhotoElectricEffect.hh
"
44
45
#include "
G4eMultipleScattering.hh
"
46
47
#include "
G4eIonisation.hh
"
48
#include "MyMollerBhabhaModel.hh"
49
#include "
G4eBremsstrahlung.hh
"
50
#include "
G4eplusAnnihilation.hh
"
51
52
#include "
G4hIonisation.hh
"
53
#include "
G4hMultipleScattering.hh
"
54
55
#include "
G4EmProcessOptions.hh
"
56
#include "
G4MscStepLimitType.hh
"
57
58
#include "
G4SystemOfUnits.hh
"
59
60
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62
PhysListEmStandard_option0::PhysListEmStandard_option0
(
const
G4String
&
name
,
63
DetectorConstruction
* det)
64
:
G4VPhysicsConstructor
(name), fDetector(det)
65
{}
66
67
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
69
PhysListEmStandard_option0::~PhysListEmStandard_option0
()
70
{}
71
72
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74
void
PhysListEmStandard_option0::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
G4ComptonScattering
* compton =
new
G4ComptonScattering
();
89
MyKleinNishinaCompton
* comptonModel =
new
MyKleinNishinaCompton
(fDetector);
90
comptonModel->
SetCSFactor
(1000.);
91
compton->
SetEmModel
(comptonModel );
92
93
pmanager->
AddDiscreteProcess
(
new
G4PhotoElectricEffect
);
94
pmanager->
AddDiscreteProcess
(compton);
95
pmanager->
AddDiscreteProcess
(
new
G4GammaConversion
);
96
97
}
else
if
(particleName ==
"e-"
) {
98
//electron
99
100
G4eIonisation
* eIoni =
new
G4eIonisation
();
101
eIoni->
SetEmModel
(
new
MyMollerBhabhaModel
);
102
103
pmanager->
AddProcess
(
new
G4eMultipleScattering
, -1, 1, 1);
104
pmanager->
AddProcess
(eIoni, -1, 2, 2);
105
/// pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);
106
107
}
else
if
(particleName ==
"e+"
) {
108
//positron
109
110
G4eIonisation
* pIoni =
new
G4eIonisation
();
111
pIoni->
SetEmModel
(
new
MyMollerBhabhaModel
);
112
113
pmanager->
AddProcess
(
new
G4eMultipleScattering
, -1, 1, 1);
114
pmanager->
AddProcess
(pIoni, -1, 2, 2);
115
/// pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);
116
pmanager->
AddProcess
(
new
G4eplusAnnihilation
, 0,-1, 3);
117
118
}
else
if
( particleName ==
"proton"
) {
119
//proton
120
pmanager->
AddProcess
(
new
G4hMultipleScattering
, -1, 1, 1);
121
pmanager->
AddProcess
(
new
G4hIonisation
, -1, 2, 2);
122
}
123
}
124
125
// Em options
126
//
127
// Main options and setting parameters are shown here.
128
// Several of them have default values.
129
//
130
G4EmProcessOptions
emOptions;
131
132
//build CSDA range
133
//
134
emOptions.
SetBuildCSDARange
(
true
);
//default=false
135
emOptions.
SetMaxEnergyForCSDARange
(10*
GeV
);
136
emOptions.
SetDEDXBinningForCSDARange
(8*7);
//default=8*7
137
}
138
139
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140
G4hMultipleScattering.hh
PhysListEmStandard_option0::~PhysListEmStandard_option0
~PhysListEmStandard_option0()
Definition:
src/PhysListEmStandard_option0.cc:69
G4eIonisation.hh
PhysListEmStandard_option0::ConstructProcess
virtual void ConstructProcess()
Definition:
src/PhysListEmStandard_option0.cc:74
python.hepunit.GeV
GeV
Definition:
hepunit.py:120
G4MscStepLimitType.hh
MyKleinNishinaCompton::SetCSFactor
void SetCSFactor(G4double factor)
Definition:
MyKleinNishinaCompton.hh:67
G4PhotoElectricEffect
Definition:
G4PhotoElectricEffect.hh:79
G4EmProcessOptions
Definition:
G4EmProcessOptions.hh:63
G4ComptonScattering.hh
G4hIonisation.hh
G4ProcessManager::AddDiscreteProcess
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4eMultipleScattering.hh
name
const XML_Char * name
Definition:
include/expat.h:151
G4ComptonScattering
Definition:
G4ComptonScattering.hh:71
G4EmProcessOptions::SetDEDXBinningForCSDARange
void SetDEDXBinningForCSDARange(G4int val)
Definition:
G4EmProcessOptions.cc:144
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:111
MyKleinNishinaCompton.hh
Definition of the MyKleinNishinaCompton class.
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
G4VEmProcess::SetEmModel
void SetEmModel(G4VEmModel *, G4int index=1)
Definition:
G4VEmProcess.cc:198
PhysListEmStandard_option0::PhysListEmStandard_option0
PhysListEmStandard_option0(const G4String &name, DetectorConstruction *det)
Definition:
src/PhysListEmStandard_option0.cc:62
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
G4ProcessManager.hh
G4ParticleDefinition.hh
G4PhotoElectricEffect.hh
MyKleinNishinaCompton
Definition:
MyKleinNishinaCompton.hh:43
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
DetectorConstruction
Definition:
environments/g4py/examples/demos/TestEm0/g4lib/DetectorConstruction.hh:46
G4eBremsstrahlung.hh
G4ProcessManager
Definition:
G4ProcessManager.hh:106
G4eplusAnnihilation
Definition:
G4eplusAnnihilation.hh:65
G4eMultipleScattering
Definition:
G4eMultipleScattering.hh:60
G4hMultipleScattering
Definition:
G4hMultipleScattering.hh:62
G4SystemOfUnits.hh
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