source: trunk/source/processes/management/include/G4VRestContinuousDiscreteProcess.hh@ 819

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

import all except CVS

File size: 9.8 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: G4VRestContinuousDiscreteProcess.hh,v 1.8.4.1 2008/04/25 11:47:06 kurasige Exp $
28// GEANT4 tag $Name: geant4-09-01-patch-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// continuous and discrete physics interactions.
38//
39// ------------------------------------------------------------
40// New Physics scheme 8 Mar. 1997 H.Kurahige
41// ------------------------------------------------------------
42// fix bugs in GetGPILSelection() 24 Jan. 1998 H.Kurashige
43// modified for new ParticleChange 12 Mar. 1998 H.Kurashige
44// Fixed a bug in PostStepGetPhysicalInteractionLength
45// 15 Apr. 2002 H.Kurashige
46//
47
48#ifndef G4VRestContinuousDiscreteProcess_h
49#define G4VRestContinuousDiscreteProcess_h 1
50
51#include "globals.hh"
52#include "G4ios.hh"
53
54#include "G4VProcess.hh"
55
56class G4VRestContinuousDiscreteProcess : public G4VProcess
57{
58 // Abstract class which defines the public behavior of
59 // discrete physics interactions.
60 public:
61
62 G4VRestContinuousDiscreteProcess(const G4String& ,
63 G4ProcessType aType = fNotDefined );
64 G4VRestContinuousDiscreteProcess(G4VRestContinuousDiscreteProcess &);
65
66 virtual ~G4VRestContinuousDiscreteProcess();
67
68
69 public :// with description
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 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 virtual G4double AtRestGetPhysicalInteractionLength(
95 const G4Track& ,
96 G4ForceCondition*
97 );
98
99 virtual G4VParticleChange* AtRestDoIt(
100 const G4Track& ,
101 const G4Step&
102 );
103
104 protected: // with description
105 virtual G4double GetMeanLifeTime(const G4Track& aTrack,G4ForceCondition* condition)=0;
106 // Calculates the mean life-time (i.e. for decays) of the
107 // particle at rest due to the occurence of the given process,
108 // or converts the probability of interaction (i.e. for
109 // annihilation) into the life-time of the particle for the
110 // occurence of the given process.
111
112 virtual G4double GetContinuousStepLimit(const G4Track& aTrack,
113 G4double previousStepSize,
114 G4double currentMinimumStep,
115 G4double& currentSafety
116 )=0;
117 private:
118 // this is the returnd value of G4GPILSelection in
119 // the arguments of AlongStepGPIL()
120 G4GPILSelection valueGPILSelection;
121
122 protected: // with description
123 // these two methods are set/get methods for valueGPILSelection
124 void SetGPILSelection(G4GPILSelection selection)
125 { valueGPILSelection = selection;};
126
127 G4GPILSelection GetGPILSelection() const{return valueGPILSelection;};
128
129 protected: // with description
130 virtual G4double GetMeanFreePath(const G4Track& aTrack,
131 G4double previousStepSize,
132 G4ForceCondition* condition
133 )=0;
134 // Calculates from the macroscopic cross section a mean
135 // free path, the value is returned in units of distance.
136
137 private:
138 // hide default constructor and assignment operator as private
139 G4VRestContinuousDiscreteProcess();
140 G4VRestContinuousDiscreteProcess & operator=(const G4VRestContinuousDiscreteProcess &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"
151
152inline G4double G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(
153 const G4Track& track,
154 G4ForceCondition* condition
155 )
156{
157 // beggining of tracking
158 ResetNumberOfInteractionLengthLeft();
159
160 // condition is set to "Not Forced"
161 *condition = NotForced;
162
163 // get mean life time
164 currentInteractionLength = GetMeanLifeTime(track, condition);
165
166#ifdef G4VERBOSE
167 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
168 G4cout << "G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength ";
169 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
170 track.GetDynamicParticle()->DumpInfo();
171 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
172 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
173 }
174#endif
175
176 return theNumberOfInteractionLengthLeft * currentInteractionLength;
177}
178
179
180inline G4VParticleChange* G4VRestContinuousDiscreteProcess::AtRestDoIt(
181 const G4Track&,
182 const G4Step&
183 )
184{
185// clear NumberOfInteractionLengthLeft
186 ClearNumberOfInteractionLengthLeft();
187
188 return pParticleChange;
189}
190
191inline G4double G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(
192 const G4Track& track,
193 G4double previousStepSize,
194 G4double currentMinimumStep,
195 G4double& currentSafety,
196 G4GPILSelection* selection
197 )
198{
199 // GPILSelection is set to defaule value of CandidateForSelection
200 valueGPILSelection = CandidateForSelection;
201
202 // get Step limit proposed by the process
203 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
204
205 // set return value for G4GPILSelection
206 *selection = valueGPILSelection;
207
208#ifdef G4VERBOSE
209 if (verboseLevel>1){
210 G4cout << "G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
211 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
212 track.GetDynamicParticle()->DumpInfo();
213 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
214 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
215 }
216#endif
217 return steplength ;
218}
219
220inline G4VParticleChange* G4VRestContinuousDiscreteProcess::AlongStepDoIt(
221 const G4Track& ,
222 const G4Step&
223 )
224{
225 return pParticleChange;
226}
227
228inline G4double G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(
229 const G4Track& track,
230 G4double previousStepSize,
231 G4ForceCondition* condition
232 )
233{
234 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
235 // beggining of tracking (or just after DoIt of this process)
236 ResetNumberOfInteractionLengthLeft();
237 } else if ( previousStepSize > 0.0) {
238 // subtract NumberOfInteractionLengthLeft
239 SubtractNumberOfInteractionLengthLeft(previousStepSize);
240 } else {
241 // zero step
242 // DO NOTHING
243 }
244
245 // condition is set to "Not Forced"
246 *condition = NotForced;
247
248 // get mean free path
249 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
250
251
252 G4double value;
253 if (currentInteractionLength <DBL_MAX) {
254 value = theNumberOfInteractionLengthLeft * currentInteractionLength;
255 } else {
256 value = DBL_MAX;
257 }
258#ifdef G4VERBOSE
259 if (verboseLevel>1){
260 G4cout << "G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength ";
261 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
262 track.GetDynamicParticle()->DumpInfo();
263 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
264 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
265 }
266#endif
267 return value;
268}
269
270inline G4VParticleChange* G4VRestContinuousDiscreteProcess::PostStepDoIt(
271 const G4Track& ,
272 const G4Step&
273 )
274{
275// clear NumberOfInteractionLengthLeft
276 ClearNumberOfInteractionLengthLeft();
277
278 return pParticleChange;
279}
280
281
282#endif
283
Note: See TracBrowser for help on using the repository browser.