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
pii
include
G4hRDEnergyLoss.hh
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
//
27
// $Id$
28
//
29
// $Id:
30
// ------------------------------------------------------------
31
// GEANT 4 class header file
32
//
33
// History: first implementation, based on object model of
34
// 2nd December 1995, G.Cosmo
35
// ---------- G4hEnergyLoss physics process -----------
36
// by Laszlo Urban, 30 May 1997
37
//
38
// ************************************************************
39
// It is the first implementation of the NEW UNIFIED ENERGY LOSS PROCESS.
40
// It calculates the continuous energy loss for charged hadrons.
41
// Processes giving contribution to the continuous loss :
42
// ionisation (= cont.ion.loss + delta ray production)
43
// can be added more easily ..........
44
// This class creates static proton/antiproton dE/dx and range tables ,
45
// which tables can be used by other processes.
46
// The energy loss for other charged hadrons is calculated from the p/pbar
47
// tables with scaled kinetic energy.
48
//
49
// 7/10/98 L.Urban some bugs fixed + some cleanup
50
// 22/10/98 L.Urban cleanup
51
// 02/02/99 L.Urban several bugs fixed
52
// 31/03/00 V.Ivanchenko rename to lowenergy as G4hLowEnergyLoss.hh
53
// 09/08/00 V.Ivanchenko remove GetContinuousStepLimit and IsApplicable
54
// 23/11/01 V.Ivanchenko Move static member-functions from header to source
55
// 22/01/03 V.Ivanchenko Cuts per region
56
// 18/04/03 V.Ivanchenko Make dRoverRange protected
57
//
58
// 31 Jul 2008 MGP Short term supply of energy loss of hadrons through clone of
59
// former G4hLowEnergyLoss (with some initial cleaning)
60
// To be replaced by reworked class to deal with condensed/discrete
61
// issues properly
62
//
63
// --------------------------------------------------------------
64
65
// Class description:
66
// Short term supply of energy loss of hadrons through clone of former G4hLowEnergyLoss
67
// (with some initial cleaning)
68
// To be replaced by reworked class to deal with condensed/discrete issues properly
69
70
// --------------------------------------------------------------
71
72
73
#ifndef G4HRDENERGYLOSS_HH
74
#define G4HRDENERGYLOSS_HH 1
75
76
#include "
G4ios.hh
"
77
#include "
globals.hh
"
78
#include "
Randomize.hh
"
79
#include "
G4VContinuousDiscreteProcess.hh
"
80
#include "
G4Material.hh
"
81
#include "
G4Element.hh
"
82
#include "
G4Proton.hh
"
83
#include "
G4AntiProton.hh
"
84
#include "
G4Electron.hh
"
85
#include "
G4VParticleChange.hh
"
86
#include "
G4Track.hh
"
87
#include "
G4Step.hh
"
88
#include "
G4PhysicsLogVector.hh
"
89
#include "
G4PhysicsLinearVector.hh
"
90
91
class
G4EnergyLossMessenger
;
92
93
class
G4hRDEnergyLoss
:
public
G4VContinuousDiscreteProcess
94
95
{
96
public
:
97
98
G4hRDEnergyLoss
(
const
G4String
& );
99
100
~G4hRDEnergyLoss
();
101
102
virtual
G4double
GetMeanFreePath
(
103
const
G4Track
& track,
104
G4double
previousStepSize,
105
enum
G4ForceCondition
*
condition
106
) = 0 ;
107
108
virtual
G4VParticleChange
*
PostStepDoIt
(
const
G4Track
& track,
109
const
G4Step
&
Step
) = 0 ;
110
111
// ---- MGP ---- All this static stuff is expected to disappear in a future
112
// development cycle
113
114
// get the number of processes contributing to the cont.energy loss
115
static
G4int
GetNumberOfProcesses
();
116
117
// set the number of processes contributing to the cont.energy loss
118
static
void
SetNumberOfProcesses
(
G4int
number);
119
120
// Increment the number of processes contributing to the cont.energy loss
121
static
void
PlusNumberOfProcesses
();
122
123
// decrement the number of processes contributing to the cont.energy loss
124
static
void
MinusNumberOfProcesses
();
125
126
static
void
SetdRoverRange
(
G4double
value
);
127
static
void
SetRndmStep
(
G4bool
value);
128
static
void
SetEnlossFluc
(
G4bool
value);
129
static
void
SetStepFunction
(
G4double
c1
,
G4double
c2);
130
131
protected
:
132
133
G4bool
CutsWhereModified
();
134
135
// G4Material *lastMaterial ;
136
const
G4double
MaxExcitationNumber
;
137
const
G4double
probLimFluct
;
138
const
long
nmaxDirectFluct
,
nmaxCont1
,
nmaxCont2
;
139
140
static
void
BuildDEDXTable
(
const
G4ParticleDefinition
& aParticleType);
141
142
protected
:
143
144
static
G4PhysicsTable
*
theDEDXpTable
;
145
static
G4PhysicsTable
*
theDEDXpbarTable
;
146
static
G4PhysicsTable
*
theRangepTable
;
147
static
G4PhysicsTable
*
theRangepbarTable
;
148
149
//inverse of the range tables
150
static
G4PhysicsTable
*
theInverseRangepTable
;
151
static
G4PhysicsTable
*
theInverseRangepbarTable
;
152
153
//lab and proper time tables
154
static
G4PhysicsTable
*
theLabTimepTable
;
155
static
G4PhysicsTable
*
theLabTimepbarTable
;
156
157
static
G4PhysicsTable
*
theProperTimepTable
;
158
static
G4PhysicsTable
*
theProperTimepbarTable
;
159
160
// processes inherited from G4hRDEnergyLoss
161
// register themselves in the static array Recorder
162
static
G4PhysicsTable
**
RecorderOfpProcess
;
163
static
G4PhysicsTable
**
RecorderOfpbarProcess
;
164
static
G4int
CounterOfpProcess
;
165
static
G4int
CounterOfpbarProcess
;
166
167
// particle mass
168
static
G4double
ParticleMass
;
169
170
// cut in range
171
static
G4double
ptableElectronCutInRange
;
172
static
G4double
pbartableElectronCutInRange
;
173
174
static
G4double
Charge
;
175
176
static
G4double
LowestKineticEnergy
;
177
static
G4double
HighestKineticEnergy
;
178
static
G4int
TotBin
;
// number of bins in table,
179
// calculated in BuildPhysicsTable
180
181
static
G4double
RTable
,
LOGRTable
;
// LOGRTable=std::log(HighestKineticEnergy
182
// /LowestKineticEnergy)/TotBin
183
// RTable = std::exp(LOGRTable)
184
185
G4PhysicsTable
*
theLossTable
;
186
187
G4double
linLossLimit
;
188
189
G4double
MinKineticEnergy
;
190
191
static
G4double
dRoverRange
;
// maximum allowed deltarange/range
192
// in one step
193
static
G4double
finalRange
;
// last step before stop
194
static
G4double
c1lim
,
c2lim
,
c3lim
;
// coeffs for computing steplimit
195
196
static
G4bool
rndmStepFlag
;
197
static
G4bool
EnlossFlucFlag
;
198
199
200
private
:
201
202
// hide assignment operator
203
204
G4hRDEnergyLoss
(
G4hRDEnergyLoss
&);
205
G4hRDEnergyLoss
& operator=(
const
G4hRDEnergyLoss
&
right
);
206
207
// variables for the integration routines
208
static
G4double
Mass,taulow,tauhigh,ltaulow,ltauhigh;
209
210
// ====================================================================
211
// static part of the class
212
213
static
void
BuildRangeTable(
const
G4ParticleDefinition
& aParticleType);
214
215
static
void
BuildInverseRangeTable(
const
G4ParticleDefinition
& aParticleType);
216
217
static
void
BuildTimeTables(
const
G4ParticleDefinition
& aParticleType);
218
219
static
void
BuildLabTimeVector(
G4int
materialIndex,
220
G4PhysicsLogVector
* rangeVector);
221
222
static
void
BuildProperTimeVector(
G4int
materialIndex,
223
G4PhysicsLogVector
* rangeVector);
224
225
static
void
InvertRangeVector(
G4int
materialIndex,
226
G4PhysicsLogVector
* rangeVector);
227
228
static
void
BuildRangeVector(
G4int
materialIndex,
229
G4PhysicsLogVector
* rangeVector);
230
231
static
G4double
LabTimeIntLog(
G4PhysicsVector
* physicsVector,
G4int
nbin);
232
233
static
G4double
ProperTimeIntLog(
G4PhysicsVector
* physicsVector,
G4int
nbin);
234
235
static
G4double
RangeIntLin(
G4PhysicsVector
* physicsVector,
G4int
nbin);
236
237
static
G4double
RangeIntLog(
G4PhysicsVector
* physicsVector,
G4int
nbin);
238
239
static
void
BuildRangeCoeffATable(
const
G4ParticleDefinition
& aParticleType);
240
static
void
BuildRangeCoeffBTable(
const
G4ParticleDefinition
& aParticleType);
241
static
void
BuildRangeCoeffCTable(
const
G4ParticleDefinition
& aParticleType);
242
243
// ====================================================================
244
245
static
G4PhysicsTable
* theDEDXTable;
246
247
static
G4PhysicsTable
* theRangeTable;
248
static
G4PhysicsTable
* theInverseRangeTable;
249
250
static
G4PhysicsTable
* theLabTimeTable;
251
static
G4PhysicsTable
* theProperTimeTable;
252
253
static
G4PhysicsTable
** RecorderOfProcess;
254
static
G4int
CounterOfProcess;
255
256
static
G4PhysicsTable
* thepRangeCoeffATable;
257
static
G4PhysicsTable
* thepRangeCoeffBTable;
258
static
G4PhysicsTable
* thepRangeCoeffCTable;
259
static
G4PhysicsTable
* thepbarRangeCoeffATable;
260
static
G4PhysicsTable
* thepbarRangeCoeffBTable;
261
static
G4PhysicsTable
* thepbarRangeCoeffCTable;
262
263
static
G4PhysicsTable
* theRangeCoeffATable;
264
static
G4PhysicsTable
* theRangeCoeffBTable;
265
static
G4PhysicsTable
* theRangeCoeffCTable;
266
static
G4int
NumberOfProcesses ;
267
268
};
269
270
#endif
271
272
273
Generated on Sat May 25 2013 14:33:35 for Geant4 by
1.8.4