source: trunk/source/processes/management/include/G4VRestContinuousProcess.hh@ 1159

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

update to geant4.9.2

File size: 8.0 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: G4VRestContinuousProcess.hh,v 1.6 2006/06/29 21:07:54 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30//
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34//
35// Class Description
36// Abstract class which defines the public behavior of
37// discrete physics interactions.
38//
39// ------------------------------------------------------------
40// New Physics scheme 8 Mar. 1997 H.Kurahige
41// ------------------------------------------------------------
42// modified 26 Mar. 1997 H.Kurahige
43// modified 16 Apr. 1997 L.Urban
44// modified 18 Sep. 1997 H.Kurashige
45// modified AlongStepGPIL etc. 17 Dec. 1997 H.Kurashige
46// fix bugs in GetGPILSelection() 24 Jan. 1998 H.Kurashige
47// modified for new ParticleChange 12 Mar. 1998 H.Kurashige
48
49#ifndef G4VRestContinuousProcess_h
50#define G4VRestContinuousProcess_h 1
51
52#include "globals.hh"
53#include "G4ios.hh"
54
55#include "G4VProcess.hh"
56
57class G4VRestContinuousProcess : public G4VProcess
58{
59 // Abstract class which defines the public behavior of
60 // discrete physics interactions.
61 public:
62
63 G4VRestContinuousProcess(const G4String& ,
64 G4ProcessType aType = fNotDefined );
65 G4VRestContinuousProcess(G4VRestContinuousProcess &);
66
67 virtual ~G4VRestContinuousProcess();
68
69 public: // with description
70
71 virtual G4double AtRestGetPhysicalInteractionLength(
72 const G4Track& ,
73 G4ForceCondition*
74 );
75
76 virtual G4VParticleChange* AtRestDoIt(
77 const G4Track& ,
78 const G4Step&
79 );
80
81 virtual G4double AlongStepGetPhysicalInteractionLength(
82 const G4Track& track,
83 G4double previousStepSize,
84 G4double currentMinimumStep,
85 G4double& currentSafety,
86 G4GPILSelection* selection
87 );
88
89 virtual G4VParticleChange* AlongStepDoIt(
90 const G4Track& ,
91 const G4Step&
92 );
93
94 // no operation in PostStepDoIt
95 virtual G4double PostStepGetPhysicalInteractionLength(
96 const G4Track& ,
97 G4double ,
98 G4ForceCondition*
99 ){ return -1.0; };
100
101 // no operation in PostStepDoIt
102 virtual G4VParticleChange* PostStepDoIt(
103 const G4Track& ,
104 const G4Step&
105 ) {return 0;};
106
107 protected: // with description
108 virtual G4double GetContinuousStepLimit(const G4Track& aTrack,
109 G4double previousStepSize,
110 G4double currentMinimumStep,
111 G4double& currentSafety
112 )=0;
113 // This pure virtual function is used to calculate step limit
114 // for AlongStep in the derived processes
115
116 private:
117 // this is the returnd value of G4GPILSelection in
118 // the arguments of AlongStepGPIL()
119 G4GPILSelection valueGPILSelection;
120
121 protected:// with description
122 // these two methods are set/get methods for valueGPILSelection
123 void SetGPILSelection(G4GPILSelection selection)
124 { valueGPILSelection = selection;};
125
126 G4GPILSelection GetGPILSelection() const{return valueGPILSelection;};
127
128 protected: // with description
129
130 virtual G4double GetMeanLifeTime(const G4Track& aTrack,G4ForceCondition* condition)=0;
131 // Calculates the mean life-time (i.e. for decays) of the
132 // particle at rest due to the occurence of the given process,
133 // or converts the probability of interaction (i.e. for
134 // annihilation) into the life-time of the particle for the
135 // occurence of the given process.
136
137 private:
138 // hide default constructor and assignment operator as private
139 G4VRestContinuousProcess();
140 G4VRestContinuousProcess & operator=(const G4VRestContinuousProcess &right);
141
142};
143
144// -----------------------------------------
145// inlined function members implementation
146// -----------------------------------------
147#include "G4Step.hh"
148#include "G4Track.hh"
149#include "G4MaterialTable.hh"
150#include "G4VParticleChange.hh"
151inline G4double G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(
152 const G4Track& track,
153 G4double previousStepSize,
154 G4double currentMinimumStep,
155 G4double& currentSafety,
156 G4GPILSelection* selection
157 )
158{
159 // GPILSelection is set to defaule value of CandidateForSelection
160 valueGPILSelection = CandidateForSelection;
161
162 // get Step limit proposed by the process
163 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
164
165 // set return value for G4GPILSelection
166 *selection = valueGPILSelection;
167#ifdef G4VERBOSE
168 if (verboseLevel>1){
169 G4cout << "G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength ";
170 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
171 track.GetDynamicParticle()->DumpInfo();
172 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
173 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
174 }
175#endif
176 return steplength ;
177}
178
179inline G4double G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(
180 const G4Track& track,
181 G4ForceCondition* condition
182 )
183{
184 // beggining of tracking
185 ResetNumberOfInteractionLengthLeft();
186
187 // condition is set to "Not Forced"
188 *condition = NotForced;
189
190 // get mean life time
191 currentInteractionLength = GetMeanLifeTime(track, condition);
192
193#ifdef G4VERBOSE
194 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
195 G4cout << "G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength ";
196 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
197 track.GetDynamicParticle()->DumpInfo();
198 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
199 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
200 }
201#endif
202
203 return theNumberOfInteractionLengthLeft * currentInteractionLength;
204}
205
206
207inline G4VParticleChange* G4VRestContinuousProcess::AtRestDoIt(
208 const G4Track&,
209 const G4Step&
210 )
211{
212// clear NumberOfInteractionLengthLeft
213 ClearNumberOfInteractionLengthLeft();
214
215 return pParticleChange;
216}
217
218inline G4VParticleChange* G4VRestContinuousProcess::AlongStepDoIt(
219 const G4Track& ,
220 const G4Step&
221 )
222{
223 return pParticleChange;
224}
225
226#endif
227
228
229
230
231
232
Note: See TracBrowser for help on using the repository browser.