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