Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VDiscreteProcess Class Referenceabstract

#include <G4VDiscreteProcess.hh>

+ Inheritance diagram for G4VDiscreteProcess:

Public Member Functions

 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 45 of file G4VDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4VDiscreteProcess() [1/2]

G4VDiscreteProcess::G4VDiscreteProcess ( const G4String & aName,
G4ProcessType aType = fNotDefined )

Definition at line 50 of file G4VDiscreteProcess.cc.

52 : G4VProcess(aName, aType)
53{
54 enableAtRestDoIt = false;
55 enableAlongStepDoIt = false;
56}
G4bool enableAtRestDoIt
G4bool enableAlongStepDoIt

Referenced by G4VDiscreteProcess().

◆ G4VDiscreteProcess() [2/2]

G4VDiscreteProcess::G4VDiscreteProcess ( G4VDiscreteProcess & right)

Definition at line 64 of file G4VDiscreteProcess.cc.

65 : G4VProcess(right)
66{
67}

◆ ~G4VDiscreteProcess()

G4VDiscreteProcess::~G4VDiscreteProcess ( )
virtual

Definition at line 59 of file G4VDiscreteProcess.cc.

60{
61}

Member Function Documentation

◆ AlongStepDoIt()

virtual G4VParticleChange * G4VDiscreteProcess::AlongStepDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Reimplemented in G4NuclearStopping.

Definition at line 90 of file G4VDiscreteProcess.hh.

93 { return nullptr; }

◆ AlongStepGetPhysicalInteractionLength()

virtual G4double G4VDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4double ,
G4double & ,
G4GPILSelection *  )
inlinevirtual

Implements G4VProcess.

Reimplemented in G4NuclearStopping.

Definition at line 70 of file G4VDiscreteProcess.hh.

76 { return -1.0; }

◆ AtRestDoIt()

virtual G4VParticleChange * G4VDiscreteProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Reimplemented in G4eplusAnnihilation, and G4HadronStoppingProcess.

Definition at line 85 of file G4VDiscreteProcess.hh.

88 { return nullptr; }

◆ AtRestGetPhysicalInteractionLength()

virtual G4double G4VDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlinevirtual

Implements G4VProcess.

Reimplemented in G4eplusAnnihilation, and G4HadronStoppingProcess.

Definition at line 78 of file G4VDiscreteProcess.hh.

81 { return -1.0; }

◆ GetCrossSection()

virtual G4double G4VDiscreteProcess::GetCrossSection ( const G4double ,
const G4MaterialCutsCouple *  )
inlinevirtual

Reimplemented in G4VEmProcess.

Definition at line 97 of file G4VDiscreteProcess.hh.

100 { return 0.0; }

Referenced by G4EmUtility::FindCrossSectionMax().

◆ GetMeanFreePath()

◆ MinPrimaryEnergy()

virtual G4double G4VDiscreteProcess::MinPrimaryEnergy ( const G4ParticleDefinition * ,
const G4Material *  )
inlinevirtual

Reimplemented in G4CoulombScattering, and G4GammaConversion.

Definition at line 103 of file G4VDiscreteProcess.hh.

106 { return 0.0; }

Referenced by G4EmTableUtil::BuildLambdaTable().

◆ operator=()

G4VDiscreteProcess & G4VDiscreteProcess::operator= ( const G4VDiscreteProcess & )
delete

◆ PostStepDoIt()

G4VParticleChange * G4VDiscreteProcess::PostStepDoIt ( const G4Track & ,
const G4Step &  )
virtual

Implements G4VProcess.

Reimplemented in G4AnnihiToMuPair, G4Channeling, G4ElNeutrinoNucleusProcess, G4ErrorTrackLengthTarget, G4ForwardXrayTR, G4GammaConversionToMuons, G4GammaGeneralProcess, G4HadronElasticProcess, G4HadronicProcess, G4LowECapture, G4MicroElecCapture, G4MicroElecSurface, G4MuNeutrinoNucleusProcess, G4NeutrinoElectronProcess, G4NeutronGeneralProcess, G4NeutronKiller, G4NuVacOscProcess, G4OpAbsorption, G4OpBoundaryProcess, G4OpMieHG, G4OpRayleigh, G4OpWLS2, G4OpWLS, G4PhononDownconversion, G4PhononReflection, G4PhononScattering, G4SynchrotronRadiation, G4SynchrotronRadiationInMat, G4TauNeutrinoNucleusProcess, G4TransitionRadiation, G4UCNAbsorption, G4UCNBoundaryProcess, G4UCNLoss, G4UCNMultiScattering, G4UnknownDecay, G4VAdjointReverseReaction, G4VEmProcess, G4VErrorLimitProcess, G4VTransitionRadiation, G4VXTRenergyLoss, and G4XrayReflection.

Definition at line 120 of file G4VDiscreteProcess.cc.

122{
123 // clear NumberOfInteractionLengthLeft
125
126 return pParticleChange;
127}
void ClearNumberOfInteractionLengthLeft()
G4VParticleChange * pParticleChange

Referenced by G4AnnihiToMuPair::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), and G4VXTRenergyLoss::PostStepDoIt().

◆ PostStepGetPhysicalInteractionLength()

G4double G4VDiscreteProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Reimplemented in G4ElNeutrinoNucleusProcess, G4ErrorMagFieldLimitProcess, G4ErrorStepLengthLimitProcess, G4ErrorTrackLengthTarget, G4GammaGeneralProcess, G4HadronicProcess, G4HadronStoppingProcess, G4LowECapture, G4MuNeutrinoNucleusProcess, G4NeutrinoElectronProcess, G4NeutronGeneralProcess, G4NeutronKiller, G4PolarizedAnnihilation, G4PolarizedCompton, G4TauNeutrinoNucleusProcess, G4UnknownDecay, G4VEmProcess, and G4VErrorLimitProcess.

Definition at line 70 of file G4VDiscreteProcess.cc.

74{
75 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0))
76 {
77 // beginning of tracking (or just after DoIt() of this process)
79 }
80 else if ( previousStepSize > 0.0)
81 {
82 // subtract NumberOfInteractionLengthLeft
84 }
85 else
86 {
87 // zero step
88 // DO NOTHING
89 }
90
91 // condition is set to "Not Forced"
93
94 // get mean free path
95 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
96
97 G4double value;
99 {
101 }
102 else
103 {
104 value = DBL_MAX;
105 }
106#ifdef G4VERBOSE
107 if (verboseLevel>1)
108 {
109 G4cout << "G4VDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
110 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
111 track.GetDynamicParticle()->DumpInfo();
112 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
113 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
114 }
115#endif
116 return value;
117}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
G4double currentInteractionLength
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80
G4int verboseLevel
G4double theNumberOfInteractionLengthLeft
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62

Referenced by G4ElNeutrinoNucleusProcess::PostStepGetPhysicalInteractionLength(), G4MuNeutrinoNucleusProcess::PostStepGetPhysicalInteractionLength(), G4NeutrinoElectronProcess::PostStepGetPhysicalInteractionLength(), and G4TauNeutrinoNucleusProcess::PostStepGetPhysicalInteractionLength().


The documentation for this class was generated from the following files: