source: trunk/source/processes/management/include/G4VProcess.hh@ 1103

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

update to geant4.9.2

File size: 16.9 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: G4VProcess.hh,v 1.25 2007/11/15 04:09:58 kurasige Exp $
[1007]28// GEANT4 tag $Name: geant4-09-02 $
[819]29//
30//
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34// History: first implementation, based on object model of
35// 2nd December 1995, G.Cosmo
36//
37// Class Description
38// This class is the virtual class for physics process objects.
39// It defines public methods which describe the behavior of
40// a physics process.
41//
42// ------------------------------------------------------------
43// New Physics scheme 18 Dec. 1996 H.Kurahige
44// ------------------------------------------------------------
45// change DoIt/GetPIL arguments type 20 Mar. 1997 H.Kurashige
46// modified AlongStepGPIL 17 Dec. 1997 H.Kurashige
47// modified for new ParticleChange 12 Mar. 1998 H.Kurashige
48// Add process trype 27 Mar. 1998 H.Kurashige
49// Remove thePhysicsTable 2 Aug. 1998 H.Kurashige
50// Add PILfactor and GPIL 3 Nov. 2000 H.Kurashige
51// Add Store/RetrievePhysicsTable 8 Nov. 2000 H.Kurashige
52// Modify Store/RetrievePhysicsTable methods 9 Mar. 2001 H.Kurashige
53// Added PreparePhysicsTable 20 Aug. 2004 H.Kurashige
54// Added isXXXXDoItIsEnabled 2 Oct. 2007 H.Kurashige
[963]55// Added ProcessSubType 15 Nov. 2007 H.Kurashige
[819]56
57#ifndef G4VProcess_h
58#define G4VProcess_h 1
59
60#include "globals.hh"
61#include <cmath>
62#include "G4ios.hh"
63
64class G4ParticleDefinition;
65class G4DynamicParticle;
66class G4Track;
67class G4Step;
68
69#include "G4PhysicsTable.hh"
70#include "G4VParticleChange.hh"
71#include "G4ForceCondition.hh"
72#include "G4GPILSelection.hh"
73#include "G4ParticleChange.hh"
74#include "G4ProcessType.hh"
75
76class G4VProcess
77{
78 // A virtual class for physics process objects. It defines
79 // public methods which describe the behavior of a
80 // physics process.
81
82 private:
83 // hide default constructor and assignment operator as private
84 // do not hide default constructor for alpha version
85 // G4VProcess G4VProcess();
86 G4VProcess & operator=(const G4VProcess &right);
87
88 public: // with description
89 // constructor requires the process name and type
90 G4VProcess(const G4String& aName = "NoName",
91 G4ProcessType aType = fNotDefined );
92
93 // copy constructor copys the name but does not copy the
94 // physics table (0 pointer is assigned)
95 G4VProcess(const G4VProcess &right);
96
97 public:
98 // destructor
99 virtual ~G4VProcess();
100
101 // equal opperators
102 G4int operator==(const G4VProcess &right) const;
103 G4int operator!=(const G4VProcess &right) const;
104
105 public: // with description
106 ////////////////////////////
107 // DoIt /////////////////
108 ///////////////////////////
109 virtual G4VParticleChange* PostStepDoIt(
110 const G4Track& track,
111 const G4Step& stepData
112 ) = 0;
113
114 virtual G4VParticleChange* AlongStepDoIt(
115 const G4Track& track,
116 const G4Step& stepData
117 ) = 0;
118 virtual G4VParticleChange* AtRestDoIt(
119 const G4Track& track,
120 const G4Step& stepData
121 ) = 0;
122 // A virtual base class function that has to be overridden
123 // by any subclass. The DoIt method actually performs the
124 // physics process and determines either momentum change
125 // of the production of secondaries etc.
126 // arguments
127 // const G4Track& track:
128 // reference to the current G4Track information
129 // const G4Step& stepData:
130 // reference to the current G4Step information
131
132 //////////////////////////
133 // GPIL //////////////
134 /////////////////////////
135 virtual G4double AlongStepGetPhysicalInteractionLength(
136 const G4Track& track,
137 G4double previousStepSize,
138 G4double currentMinimumStep,
139 G4double& proposedSafety,
140 G4GPILSelection* selection) = 0;
141
142 virtual G4double AtRestGetPhysicalInteractionLength(
143 const G4Track& track,
144 G4ForceCondition* condition
145 ) = 0;
146
147 virtual G4double PostStepGetPhysicalInteractionLength(
148 const G4Track& track,
149 G4double previousStepSize,
150 G4ForceCondition* condition
151 ) = 0;
152
153 // Returns the Step-size (actual length) which is allowed
154 // by "this" process. (for AtRestGetPhysicalInteractionLength,
155 // return value is Step-time) The NumberOfInteractionLengthLeft is
156 // recalculated by using previousStepSize and the Step-size is
157 // calucalted accoding to the resultant NumberOfInteractionLengthLeft.
158 // using NumberOfInteractionLengthLeft, which is recalculated at
159 // arguments
160 // const G4Track& track:
161 // reference to the current G4Track information
162 // G4double* previousStepSize:
163 // the Step-size (actual length) of the previous Step
164 // of this track. Negative calue indicates that
165 // NumberOfInteractionLengthLeft must be reset.
166 // the current physical interaction legth of this process
167 // G4ForceCondition* condition:
168 // the flag indicates DoIt of this process is forced
169 // to be called
170 // Forced: Corresponding DoIt is forced
171 // NotForced: Corresponding DoIt is called
172 // if the Step size of this Step is determined
173 // by this process
174 // !! AlongStepDoIt is always called !!
175 // G4double& currentMinimumStep:
176 // this value is used for transformation of
177 // true path length to geometrical path length
178
179 G4double GetCurrentInteractionLength() const;
180 // Returns currentInteractionLength
181
182 ////////// PIL factor ////////
183 void SetPILfactor(G4double value);
184 G4double GetPILfactor() const;
185 // Set/Get factor for PhysicsInteractionLength
186 // which is passed to G4SteppingManager for both AtRest and PostStep
187
188 // These three GPIL methods are used by Stepping Manager.
189 // They invoke virtual GPIL methods listed above.
190 // As for AtRest and PostStep the returned value is multipled by thePILfactor
191 //
192 G4double AlongStepGPIL( const G4Track& track,
193 G4double previousStepSize,
194 G4double currentMinimumStep,
195 G4double& proposedSafety,
196 G4GPILSelection* selection );
197
198 G4double AtRestGPIL( const G4Track& track,
199 G4ForceCondition* condition );
200
201 G4double PostStepGPIL( const G4Track& track,
202 G4double previousStepSize,
203 G4ForceCondition* condition );
204
205 //////////////////////
206 virtual G4bool IsApplicable(const G4ParticleDefinition&){return true;};
207 // Returns true if this process object is applicable to
208 // the particle type
209 // Process will not be registered to a particle if IsApplicable is false
210
211 virtual void BuildPhysicsTable(const G4ParticleDefinition&){};
212 // Messaged by the Particle definition (via the Process manager)
213 // whenever cross section tables have to be rebuilt (i.e. if new
214 // materials have been defined).
215 // It is overloaded by individual processes when they need physics
216 // tables.
217
218 virtual void PreparePhysicsTable(const G4ParticleDefinition&){};
219 // Messaged by the Particle definition (via the Process manager)
220 // whenever cross section tables have to be prepare for rebuilt
221 // (i.e. if new materials have been defined).
222 // It is overloaded by individual processes when they need physics
223 // tables.
224
225 // Processes which Build physics tables independent of cuts
226 // (for example in their constructors)
227 // should preferably use private
228 // void BuildThePhysicsTable() and void PreparePhysicsTable().
229 // Not another BuildPhysicsTable, please.
230
231
232 virtual G4bool StorePhysicsTable(const G4ParticleDefinition* ,
233 const G4String&,
234 G4bool ascii = false)
235 {ascii=false; return true;}
236 // Store PhysicsTable in a file.
237 // (return false in case of failure at I/O )
238
239 virtual G4bool RetrievePhysicsTable( const G4ParticleDefinition* ,
240 const G4String&,
241 G4bool ascii = false)
242 {ascii=false; return false;}
243 // Retrieve Physics from a file.
244 // (return true if the Physics Table can be build by using file)
245 // (return false if the process has no functionality or in case of failure)
246 // File name should be defined by each process
247 // and the file should be placed under the directory specifed by the argument.
248 const G4String& GetPhysicsTableFileName(const G4ParticleDefinition* ,
249 const G4String& directory,
250 const G4String& tableName,
251 G4bool ascii =false);
252 // this method is utility for Store/RetreivePhysicsTable
253
254 ////////////////////////////
255 const G4String& GetProcessName() const;
256 // Returns the name of the process.
257
258 G4ProcessType GetProcessType() const;
259 // Returns the process type.
260
261 void SetProcessType(G4ProcessType );
262 // Set the process type.
263
[963]264 G4int GetProcessSubType() const;
265 // Returns the process sub type.
266
267 void SetProcessSubType(G4int );
268 // Set the process sub type.
269
[819]270 static const G4String& GetProcessTypeName(G4ProcessType );
271 // Returns the process type name
272
273 virtual void StartTracking(G4Track*);
274 virtual void EndTracking();
275 // inform Start/End of tracking for each track to the physics process
276
277 public:
278 virtual void SetProcessManager(const G4ProcessManager*);
279 // A process manager set its own pointer when the process is registered
280 // the process Manager
281 virtual const G4ProcessManager* GetProcessManager();
282 // Get the process manager which the process belongs to
283
284 protected:
285 const G4ProcessManager* aProcessManager;
286
287 protected:
288 G4VParticleChange* pParticleChange;
289 // The pointer to G4VParticleChange object
290 // which is modified and returned by address by the DoIt() method.
291 // This pointer should be set in each physics process
292 // after construction of derived class object.
293
294 G4ParticleChange aParticleChange;
295 // This object is kept for compatibility with old scheme
296 // This will be removed in future
297
298 G4double theNumberOfInteractionLengthLeft;
299 // The flight length left for the current tracking particle
300 // in unit of "Interaction length".
301
302 G4double currentInteractionLength;
303 // The InteractionLength in the current material
304
305 public: // with description
306 virtual void ResetNumberOfInteractionLengthLeft();
307 // reset (determine the value of)NumberOfInteractionLengthLeft
308
309 protected: // with description
310 virtual void SubtractNumberOfInteractionLengthLeft(
311 G4double previousStepSize
312 );
313 // subtract NumberOfInteractionLengthLeft by the value corresponding to
314 // previousStepSize
315
316 virtual void ClearNumberOfInteractionLengthLeft();
317 // clear NumberOfInteractionLengthLeft
318 // !!! This method should be at the end of PostStepDoIt()
319 // !!! and AtRestDoIt
320
321 public: // with description
322 // These methods indicate which DoIt is enabled
323 // These methods are used by G4ProcessManager to check
324 // that ordering parameters are set properly
325 G4bool isAtRestDoItIsEnabled() const;
326 G4bool isAlongStepDoItIsEnabled() const;
327 G4bool isPostStepDoItIsEnabled() const;
328
329 protected:
330 G4String theProcessName;
331 // The name of the process
332
333 G4String thePhysicsTableFileName;
334
335 G4ProcessType theProcessType;
336 // The type of the process
337
[963]338 G4int theProcessSubType;
339 // The sub type of the process
340
[819]341 G4double thePILfactor;
342 // factor for PhysicsInteractionLength
343 // which is passed to G4SteppingManager
344
345 G4bool enableAtRestDoIt;
346 G4bool enableAlongStepDoIt;
347 G4bool enablePostStepDoIt;
348
349 public: // with description
350 virtual void DumpInfo() const;
351 // dump out process information
352
353 public: // with description
354 void SetVerboseLevel(G4int value);
355 G4int GetVerboseLevel() const;
356 // set/get controle flag for output message
357 // 0: Silent
358 // 1: Warning message
359 // 2: More
360
361
362 protected:
363 G4int verboseLevel;
364 // controle flag for output message
365
366};
367
368// -----------------------------------------
369// inlined function members implementation
370// -----------------------------------------
371#include "Randomize.hh"
372
373inline
374 const G4String& G4VProcess::GetProcessName() const
375{
376 return theProcessName;
377}
378
379inline
380 G4ProcessType G4VProcess::GetProcessType() const
381{
382 return theProcessType;
383}
384
385inline
386 void G4VProcess::SetProcessType(G4ProcessType aType)
387{
388 theProcessType = aType;
389}
390
[963]391inline
392 G4int G4VProcess::GetProcessSubType() const
393{
394 return theProcessSubType;
395}
396
397inline
398 void G4VProcess::SetProcessSubType(G4int value)
399{
400 theProcessSubType = value;
401}
402
[819]403inline void G4VProcess::SetVerboseLevel(G4int value)
404{
405 verboseLevel = value;
406}
407
408inline G4int G4VProcess::GetVerboseLevel() const
409{
410 return verboseLevel;
411}
412
413inline void G4VProcess::ResetNumberOfInteractionLengthLeft()
414{
415 theNumberOfInteractionLengthLeft = -std::log( G4UniformRand() );
416}
417
418inline void G4VProcess::ClearNumberOfInteractionLengthLeft()
419{
420 theNumberOfInteractionLengthLeft = -1.0;
421}
422
423inline G4double G4VProcess::GetCurrentInteractionLength() const
424{
425 return currentInteractionLength;
426}
427
428inline void G4VProcess::SetPILfactor(G4double value)
429{
430 if (value>0.) {
431 thePILfactor = value;
432 }
433}
434
435inline G4double G4VProcess::GetPILfactor() const
436{
437 return thePILfactor;
438}
439
440inline G4double G4VProcess::AlongStepGPIL( const G4Track& track,
441 G4double previousStepSize,
442 G4double currentMinimumStep,
443 G4double& proposedSafety,
444 G4GPILSelection* selection )
445{
446 G4double value
447 =AlongStepGetPhysicalInteractionLength(track, previousStepSize, currentMinimumStep, proposedSafety, selection);
448 return value;
449}
450
451inline G4double G4VProcess::AtRestGPIL( const G4Track& track,
452 G4ForceCondition* condition )
453{
454 G4double value
455 =AtRestGetPhysicalInteractionLength(track, condition);
456 return thePILfactor*value;
457}
458
459inline G4double G4VProcess::PostStepGPIL( const G4Track& track,
460 G4double previousStepSize,
461 G4ForceCondition* condition )
462{
463 G4double value
464 =PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
465 return thePILfactor*value;
466}
467
468inline
469 void G4VProcess::SetProcessManager(const G4ProcessManager* procMan)
470{
471 aProcessManager = procMan;
472}
473
474inline
475 const G4ProcessManager* G4VProcess::GetProcessManager()
476{
477 return aProcessManager;
478}
479
480inline
481 G4bool G4VProcess::isAtRestDoItIsEnabled() const
482{
483 return enableAtRestDoIt;
484}
485
486inline
487 G4bool G4VProcess::isAlongStepDoItIsEnabled() const
488{
489 return enableAlongStepDoIt;
490}
491
492inline
493 G4bool G4VProcess::isPostStepDoItIsEnabled() const
494{
495 return enablePostStepDoIt;
496}
497
498#endif
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
Note: See TracBrowser for help on using the repository browser.