source: trunk/source/processes/management/include/G4IVRestDiscreteProcess.hh@ 1111

Last change on this file since 1111 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 6.4 KB
Line 
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: G4IVRestDiscreteProcess.hh,v 1.8 2006/06/29 21:07:16 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30// $Id:
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34//
35// Class Description
36// Abstract class which defines the public behavior of
37// rest + discrete physics interactions.
38//
39// ------------------------------------------------------------
40// New Physics scheme 8 Mar. 1997 H.Kurahige
41// ------------------------------------------------------------
42
43
44#ifndef G4IVRestDiscreteProcess_h
45#define G4IVRestDiscreteProcess_h 1
46
47#include "globals.hh"
48#include "G4ios.hh"
49
50#include "G4VProcess.hh"
51
52
53class G4IVRestDiscreteProcess : public G4VProcess
54{
55 // Abstract class which defines the public behavior of
56 // rest + discrete physics interactions.
57
58 public: // with description
59
60 // constructors
61 G4IVRestDiscreteProcess(const G4String& ,
62 G4ProcessType aType = fNotDefined );
63 G4IVRestDiscreteProcess(G4IVRestDiscreteProcess &);
64
65 public:
66 virtual ~G4IVRestDiscreteProcess();
67
68 public: // with description
69 // GPIL and DoIt methods derived from G4VProcess
70 virtual G4double PostStepGetPhysicalInteractionLength(
71 const G4Track& track,
72 G4double previousStepSize,
73 G4ForceCondition* condition
74 );
75
76 virtual G4VParticleChange* PostStepDoIt(
77 const G4Track& ,
78 const G4Step&
79 );
80
81 virtual G4double AtRestGetPhysicalInteractionLength(
82 const G4Track& ,
83 G4ForceCondition*
84 );
85
86 G4VParticleChange* AtRestDoIt(
87 const G4Track& ,
88 const G4Step&
89 );
90
91 // no operation in AlongStepDoIt
92 virtual G4double AlongStepGetPhysicalInteractionLength(
93 const G4Track&,
94 G4double ,
95 G4double ,
96 G4double& ,
97 G4GPILSelection*
98 ){ return -1.0; }
99
100 // no operation in AlongStepDoIt
101 virtual G4VParticleChange* AlongStepDoIt(
102 const G4Track& ,
103 const G4Step&
104 ) {return 0;}
105
106 protected:// with description
107 virtual void SubtractNumberOfInteractionLengthLeft(
108 G4double previousStepSize) ;
109
110 virtual G4double GetMeanLifeTime(const G4Track& aTrack,G4ForceCondition* condition)=0;
111 // Calculates the mean life-time (i.e. for decays) of the
112 // particle at rest due to the occurence of the given process,
113 // or converts the probability of interaction (i.e. for
114 // annihilation) into the life-time of the particle for the
115 // occurence of the given process.
116
117 private:
118 // hide default constructor and assignment operator as private
119 G4IVRestDiscreteProcess();
120 G4IVRestDiscreteProcess & operator=(const G4IVRestDiscreteProcess &right);
121
122
123 protected:
124 G4PhysicsTable* theNlambdaTable ;
125 G4PhysicsTable* theInverseNlambdaTable ;
126
127 const G4double BIGSTEP;
128};
129
130// -----------------------------------------
131// inlined function members implementation
132// -----------------------------------------
133#include "Randomize.hh"
134#include "G4Step.hh"
135#include "G4Track.hh"
136#include "G4MaterialTable.hh"
137#include "G4VParticleChange.hh"
138
139inline
140 void G4IVRestDiscreteProcess::SubtractNumberOfInteractionLengthLeft(
141 G4double )
142 {
143 // dummy routine
144 ;
145 }
146
147inline G4VParticleChange* G4IVRestDiscreteProcess::PostStepDoIt(
148 const G4Track& ,
149 const G4Step&
150 )
151{
152// reset NumberOfInteractionLengthLeft
153 ClearNumberOfInteractionLengthLeft();
154
155 return pParticleChange;
156}
157
158inline G4double G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(
159 const G4Track& track,
160 G4ForceCondition* condition
161 )
162{
163 // beggining of tracking
164 ResetNumberOfInteractionLengthLeft();
165
166 // condition is set to "Not Forced"
167 *condition = NotForced;
168
169 // get mean life time
170 currentInteractionLength = GetMeanLifeTime(track, condition);
171
172#ifdef G4VERBOSE
173 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
174 G4cout << "G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
175 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
176 track.GetDynamicParticle()->DumpInfo();
177 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
178 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
179 }
180#endif
181
182 return theNumberOfInteractionLengthLeft * currentInteractionLength;
183}
184
185
186inline G4VParticleChange* G4IVRestDiscreteProcess::AtRestDoIt(
187 const G4Track&,
188 const G4Step&
189 )
190{
191// clear NumberOfInteractionLengthLeft
192 ClearNumberOfInteractionLengthLeft();
193
194 return pParticleChange;
195}
196
197
198
199#endif
200
201
202
203
204
205
206
207
208
209
210
211
212
Note: See TracBrowser for help on using the repository browser.