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
dna
management
include
G4VITRestProcess.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
// $Id: G4VITRestProcess.hh 64057 2012-10-30 15:04:49Z gcosmo $
27
//
29
//
30
// WARNING : This class is released as a prototype.
31
// It might strongly evolve or even disapear in the next releases.
32
//
33
// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
34
//
35
// History:
36
// -----------
37
// 10 Oct 2011 M.Karamitros created
38
//
39
// -------------------------------------------------------------------
40
41
42
#ifndef G4VITRestProcess_h
43
#define G4VITRestProcess_h 1
44
45
#include <
CLHEP/Units/SystemOfUnits.h
>
46
47
#include "
G4VITProcess.hh
"
48
53
class
G4VITRestProcess
:
public
G4VITProcess
54
{
55
// Abstract class which defines the public behavior of
56
// physics interactions at rest.
57
58
public
:
59
G4VITRestProcess
(
const
G4String
& ,
60
G4ProcessType
aType =
fNotDefined
);
61
G4VITRestProcess
(
const
G4VITRestProcess
& );
62
63
virtual
~G4VITRestProcess
();
64
65
public
:
// with description
66
virtual
G4double
AtRestGetPhysicalInteractionLength
(
67
const
G4Track
& track,
68
G4ForceCondition
*
condition
69
);
70
71
virtual
G4VParticleChange
*
AtRestDoIt
(
72
const
G4Track
& ,
73
const
G4Step
&
74
);
75
76
// no operation in PostStepDoIt and AlongStepDoIt
77
virtual
G4double
AlongStepGetPhysicalInteractionLength
(
78
const
G4Track
&,
79
G4double
,
80
G4double
,
81
G4double
& ,
82
G4GPILSelection
*
83
){
return
-1.0; }
84
85
virtual
G4double
PostStepGetPhysicalInteractionLength
(
86
const
G4Track
& ,
87
G4double
,
88
G4ForceCondition
*
89
) {
return
-1.0; }
90
91
// no operation in PostStepDoIt and AlongStepDoIt
92
virtual
G4VParticleChange
*
PostStepDoIt
(
93
const
G4Track
& ,
94
const
G4Step
&
95
) {
return
0;}
96
97
virtual
G4VParticleChange
*
AlongStepDoIt
(
98
const
G4Track
& ,
99
const
G4Step
&
100
) {
return
0;}
101
102
protected
:
// with description
103
104
virtual
G4double
GetMeanLifeTime
(
const
G4Track
& aTrack,
G4ForceCondition
* condition)=0;
105
// Calculates the mean life-time (i.e. for decays) of the
106
// particle at rest due to the occurence of the given process,
107
// or converts the probability of interaction (i.e. for
108
// annihilation) into the life-time of the particle for the
109
// occurence of the given process.
110
111
protected
:
112
// hide default constructor and assignment operator as private
113
G4VITRestProcess
();
114
G4VITRestProcess
&
operator=
(
const
G4VITRestProcess
&
right
);
115
};
116
117
// -----------------------------------------
118
// inlined function members implementation
119
// -----------------------------------------
120
inline
G4double
G4VITRestProcess::AtRestGetPhysicalInteractionLength
(
121
const
G4Track
& track,
122
G4ForceCondition
*
condition
123
)
124
{
125
// beggining of tracking
126
ResetNumberOfInteractionLengthLeft
();
127
128
// condition is set to "Not Forced"
129
*condition =
NotForced
;
130
131
// get mean life time
132
fpState
->
currentInteractionLength
=
GetMeanLifeTime
(track, condition);
133
134
#ifdef G4VERBOSE
135
if
((
fpState
->
currentInteractionLength
<0.0) || (
verboseLevel
>2)){
136
G4cout
<<
"G4VITRestProcess::AtRestGetPhysicalInteractionLength "
;
137
G4cout
<<
"[ "
<<
GetProcessName
() <<
"]"
<<
G4endl
;
138
track.
GetDynamicParticle
()->
DumpInfo
();
139
G4cout
<<
" in Material "
<< track.
GetMaterial
()->
GetName
() <<
G4endl
;
140
G4cout
<<
"MeanLifeTime = "
<<
fpState
->
currentInteractionLength
/CLHEP::ns <<
"[ns]"
<<
G4endl
;
141
}
142
#endif
143
144
return
(
fpState
->
theNumberOfInteractionLengthLeft
) * (
fpState
->
currentInteractionLength
);
145
}
146
147
148
inline
G4VParticleChange
*
G4VITRestProcess::AtRestDoIt
(
149
const
G4Track
&,
150
const
G4Step
&
151
)
152
{
153
// clear NumberOfInteractionLengthLeft
154
ClearNumberOfInteractionLengthLeft
();
155
156
return
pParticleChange
;
157
}
158
159
160
#endif
161
162
163
164
Generated on Sat May 25 2013 14:33:30 for Geant4 by
1.8.4