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