Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
standard
src
G4ionIonisation.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
// $Id$
27
//
28
// -------------------------------------------------------------------
29
//
30
// GEANT4 Class file
31
//
32
//
33
// File name: G4ionIonisation
34
//
35
// Author: Vladimir Ivanchenko
36
//
37
// Creation date: 07.05.2002
38
//
39
// Modifications:
40
//
41
// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42
// 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
43
// 13-02-03 SubCutoff regime is assigned to a region (V.Ivanchenko)
44
// 18-04-03 Use IonFluctuations (V.Ivanchenko)
45
// 03-08-03 Add effective charge (V.Ivanchenko)
46
// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
47
// 27-05-04 Set integral to be a default regime (V.Ivanchenko)
48
// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
49
// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
50
// 10-01-06 SetStepLimits -> SetStepFunction (V.Ivantchenko)
51
// 10-05-06 Add a possibility to download user data (V.Ivantchenko)
52
// 13-05-06 Add data for light ion stopping in water (V.Ivantchenko)
53
// 14-01-07 use SetEmModel() and SetFluctModel() from G4VEnergyLossProcess (mma)
54
// 16-05-07 Add data for light ion stopping only for GenericIon (V.Ivantchenko)
55
// 07-11-07 Fill non-ionizing energy loss (V.Ivantchenko)
56
// 12-09-08 Removed InitialiseMassCharge and CorrectionsAlongStep (VI)
57
//
58
//
59
// -------------------------------------------------------------------
60
//
61
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
64
#include "
G4ionIonisation.hh
"
65
#include "
G4PhysicalConstants.hh
"
66
#include "
G4SystemOfUnits.hh
"
67
#include "
G4Electron.hh
"
68
#include "
G4Proton.hh
"
69
#include "
G4GenericIon.hh
"
70
#include "
G4BraggModel.hh
"
71
#include "
G4BraggIonModel.hh
"
72
#include "
G4BetheBlochModel.hh
"
73
#include "
G4UnitsTable.hh
"
74
#include "
G4LossTableManager.hh
"
75
#include "
G4WaterStopping.hh
"
76
#include "
G4EmCorrections.hh
"
77
#include "
G4IonFluctuations.hh
"
78
79
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
81
using namespace
std;
82
83
G4ionIonisation::G4ionIonisation
(
const
G4String
&
name
)
84
:
G4VEnergyLossProcess
(name),
85
theParticle(0),
86
isInitialised(false),
87
stopDataActive(true)
88
{
89
SetLinearLossLimit
(0.02);
90
SetStepFunction
(0.1, 0.01*
mm
);
91
SetProcessSubType
(
fIonisation
);
92
SetSecondaryParticle
(
G4Electron::Electron
());
93
corr =
G4LossTableManager::Instance
()->
EmCorrections
();
94
eth = 2*
MeV
;
95
}
96
97
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
99
G4ionIonisation::~G4ionIonisation
()
100
{}
101
102
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103
104
G4bool
G4ionIonisation::IsApplicable
(
const
G4ParticleDefinition
&
p
)
105
{
106
return
(p.
GetPDGCharge
() != 0.0 && !p.
IsShortLived
() &&
107
p.
GetParticleType
() ==
"nucleus"
);
108
}
109
110
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111
112
G4double
G4ionIonisation::MinPrimaryEnergy
(
const
G4ParticleDefinition
*
p
,
113
const
G4Material
*,
114
G4double
cut)
115
{
116
return
117
p->
GetPDGMass
()*(std::sqrt(1. + 0.5*cut/CLHEP::electron_mass_c2) - 1.0);
118
}
119
120
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
122
void
G4ionIonisation::InitialiseEnergyLossProcess
(
123
const
G4ParticleDefinition
*
part
,
124
const
G4ParticleDefinition
* bpart)
125
{
126
const
G4ParticleDefinition
*
ion
=
G4GenericIon::GenericIon
();
127
128
if
(!isInitialised) {
129
130
theParticle =
part
;
131
132
// define base particle
133
const
G4ParticleDefinition
* theBaseParticle = 0;
134
135
if
(part == ion) { theBaseParticle = 0; }
136
else
if
(bpart == 0) { theBaseParticle = ion; }
137
else
{ theBaseParticle = bpart; }
138
139
SetBaseParticle
(theBaseParticle);
140
141
if
(!
EmModel
(1)) {
SetEmModel
(
new
G4BraggIonModel
(), 1); }
142
EmModel
(1)->
SetLowEnergyLimit
(
MinKinEnergy
());
143
144
// model limit defined for protons
145
eth = (
EmModel
(1)->
HighEnergyLimit
())*part->
GetPDGMass
()/
proton_mass_c2
;
146
EmModel
(1)->
SetHighEnergyLimit
(eth);
147
148
if
(!
FluctModel
()) {
SetFluctModel
(
new
G4IonFluctuations
()); }
149
AddEmModel
(1,
EmModel
(1),
FluctModel
());
150
151
if
(!
EmModel
(2)) {
SetEmModel
(
new
G4BetheBlochModel
(),2); }
152
EmModel
(2)->
SetLowEnergyLimit
(eth);
153
EmModel
(2)->
SetHighEnergyLimit
(
MaxKinEnergy
());
154
AddEmModel
(2,
EmModel
(2),
FluctModel
());
155
156
// Add ion stoping tables for Generic Ion
157
if
(part == ion) {
158
G4WaterStopping
ws(corr);
159
corr->
SetIonisationModels
(
EmModel
(1),
EmModel
(2));
160
}
161
isInitialised =
true
;
162
}
163
// reinitialisation of corrections for the new run
164
if
(part == ion) { corr->
InitialiseForNewRun
(); }
165
}
166
167
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
169
void
G4ionIonisation::PrintInfo
()
170
{
171
if
(stopDataActive &&
G4GenericIon::GenericIon
() == theParticle) {
172
G4cout
<<
" Stopping Power data for "
173
<< corr->
GetNumberOfStoppingVectors
()
174
<<
" ion/material pairs "
175
<<
G4endl
;
176
}
177
}
178
179
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181
void
G4ionIonisation::AddStoppingData
(
G4int
Z
,
G4int
A,
182
const
G4String
& mname,
183
G4PhysicsVector
* dVector)
184
{
185
corr->
AddStoppingData
(Z, A, mname, dVector);
186
}
187
188
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Generated on Sat May 25 2013 14:33:37 for Geant4 by
1.8.4