source: trunk/source/processes/management/include/G4VContinuousDiscreteProcess.hh@ 1115

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

update to geant4.9.2

File size: 8.4 KB
RevLine 
[819]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//
[963]27// $Id: G4VContinuousDiscreteProcess.hh,v 1.9 2007/11/15 04:10:18 kurasige Exp $
[1007]28// GEANT4 tag $Name: geant4-09-02 $
[819]29//
30//
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34// Abstract class which defines the public behavior of
35//
36// Class Description
37// 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 G4VContinuousDiscreteProcess_h
49#define G4VContinuousDiscreteProcess_h 1
50
51#include "globals.hh"
52#include "G4ios.hh"
53
54#include "G4VProcess.hh"
55
56class G4VContinuousDiscreteProcess : public G4VProcess
57{
58 // Abstract class which defines the public behavior of
59 // discrete physics interactions.
60 public:
61
62 G4VContinuousDiscreteProcess(const G4String& ,
63 G4ProcessType aType = fNotDefined );
64 G4VContinuousDiscreteProcess(G4VContinuousDiscreteProcess &);
65
66 virtual ~G4VContinuousDiscreteProcess();
67
68 public :// with description
69 virtual G4double PostStepGetPhysicalInteractionLength(
70 const G4Track& track,
71 G4double previousStepSize,
72 G4ForceCondition* condition
73 );
74
75 virtual G4VParticleChange* PostStepDoIt(
76 const G4Track& ,
77 const G4Step&
78 );
79
80 virtual G4double AlongStepGetPhysicalInteractionLength(
81 const G4Track&,
82 G4double previousStepSize,
83 G4double currentMinimumStep,
84 G4double& currentSafety,
85 G4GPILSelection* selection
86 );
87
88 virtual G4VParticleChange* AlongStepDoIt(
89 const G4Track& ,
90 const G4Step&
91 );
92
93 // no operation in AtRestDoIt
94 virtual G4double AtRestGetPhysicalInteractionLength(
95 const G4Track& ,
96 G4ForceCondition*
97 ) { return -1.0; };
98
99 // no operation in AtRestDoIt
100 virtual G4VParticleChange* AtRestDoIt(
101 const G4Track& ,
102 const G4Step&
103 ) {return 0;};
104
105 protected:// with description
106 virtual G4double GetMeanFreePath(const G4Track& aTrack,
107 G4double previousStepSize,
108 G4ForceCondition* condition
109 )=0;
110 // Calculates from the macroscopic cross section a mean
111 // free path, the value is returned in units of distance.
112
113 protected:// with description
114 virtual G4double GetContinuousStepLimit(const G4Track& aTrack,
115 G4double previousStepSize,
116 G4double currentMinimumStep,
117 G4double& currentSafety
118 )=0;
119 private:
120 // this is the returnd value of G4GPILSelection in
121 // the arguments of AlongStepGPIL()
122 G4GPILSelection valueGPILSelection;
123
124 protected:// with description
125 // these two methods are set/get methods for valueGPILSelection
126 void SetGPILSelection(G4GPILSelection selection)
127 { valueGPILSelection = selection;};
128
129 G4GPILSelection GetGPILSelection() const{return valueGPILSelection;};
130
131 private:
132 // hide default constructor and assignment operator as private
133 G4VContinuousDiscreteProcess();
134 G4VContinuousDiscreteProcess & operator=(const G4VContinuousDiscreteProcess &right);
135
136};
137// -----------------------------------------
138// inlined function members implementation
139// -----------------------------------------
140#include "G4Step.hh"
141#include "G4Track.hh"
142#include "G4MaterialTable.hh"
143#include "G4VParticleChange.hh"
144
145inline G4double G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(
146 const G4Track& track,
147 G4double previousStepSize,
148 G4ForceCondition* condition
149 )
150{
151 if ( (previousStepSize <=0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
152 // beggining of tracking (or just after DoIt of this process)
153 ResetNumberOfInteractionLengthLeft();
154 } else if ( previousStepSize > 0.0) {
155 // subtract NumberOfInteractionLengthLeft
156 SubtractNumberOfInteractionLengthLeft(previousStepSize);
157 } else {
158 // zero step
159 // DO NOTHING
160 }
161
162 // condition is set to "Not Forced"
163 *condition = NotForced;
164
165 // get mean free path
166 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
167
168 G4double value;
169 if (currentInteractionLength <DBL_MAX) {
170 value = theNumberOfInteractionLengthLeft * currentInteractionLength;
171 } else {
172 value = DBL_MAX;
173 }
174#ifdef G4VERBOSE
175 if (verboseLevel>1){
176 G4cout << "G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength ";
177 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
178 track.GetDynamicParticle()->DumpInfo();
179 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
180 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
181 }
182#endif
183 return value;
184}
185
186
187inline G4VParticleChange* G4VContinuousDiscreteProcess::PostStepDoIt(
188 const G4Track& ,
189 const G4Step&
190 )
191{
192// clear NumberOfInteractionLengthLeft
193 ClearNumberOfInteractionLengthLeft();
194 return pParticleChange;
195}
196
197inline G4VParticleChange* G4VContinuousDiscreteProcess::AlongStepDoIt(
198 const G4Track& ,
199 const G4Step&
200 )
201{
202// clear NumberOfInteractionLengthLeft
203 ClearNumberOfInteractionLengthLeft();
204 return pParticleChange;
205}
206
207inline G4double G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(
208 const G4Track& track,
209 G4double previousStepSize,
210 G4double currentMinimumStep,
211 G4double& currentSafety,
212 G4GPILSelection* selection
213 )
214{
215 // GPILSelection is set to defaule value of CandidateForSelection
216 valueGPILSelection = CandidateForSelection;
217
218 // get Step limit proposed by the process
219 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
220
221 // set return value for G4GPILSelection
222 *selection = valueGPILSelection;
223
224#ifdef G4VERBOSE
225 if (verboseLevel>1){
226 G4cout << "G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
227 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
228 track.GetDynamicParticle()->DumpInfo();
229 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
230 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
231 }
232#endif
233 return steplength ;
234}
235
236#endif
237
238
239
240
241
242
Note: See TracBrowser for help on using the repository browser.