Geant4 9.6.0
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 &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

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

Additional Inherited Members

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

Detailed Description

Definition at line 58 of file G4VDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4VDiscreteProcess() [1/2]

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

Definition at line 54 of file G4VDiscreteProcess.cc.

55 : G4VProcess(aName, aType)
56{
57 enableAtRestDoIt = false;
58 enableAlongStepDoIt = false;
59
60}
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351

Referenced by G4VDiscreteProcess().

◆ G4VDiscreteProcess() [2/2]

G4VDiscreteProcess::G4VDiscreteProcess ( G4VDiscreteProcess right)

Definition at line 66 of file G4VDiscreteProcess.cc.

67 : G4VProcess(right)
68{
69}

◆ ~G4VDiscreteProcess()

G4VDiscreteProcess::~G4VDiscreteProcess ( )
virtual

Definition at line 62 of file G4VDiscreteProcess.cc.

63{
64}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4NuclearStopping.

Definition at line 102 of file G4VDiscreteProcess.hh.

105 {return 0;};

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4NuclearStopping.

Definition at line 83 of file G4VDiscreteProcess.hh.

89 { return -1.0; };

◆ AtRestDoIt()

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

Implements G4VProcess.

Reimplemented in G4HadronStoppingProcess, G4eplusPolarizedAnnihilation, and G4eplusAnnihilation.

Definition at line 97 of file G4VDiscreteProcess.hh.

100 {return 0;};

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4eplusPolarizedAnnihilation, G4eplusAnnihilation, and G4HadronStoppingProcess.

Definition at line 91 of file G4VDiscreteProcess.hh.

94 { return -1.0; };

◆ GetMeanFreePath()

◆ PostStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4ErrorTrackLengthTarget, G4VErrorLimitProcess, G4VAdjointReverseReaction, G4VEmProcess, G4TransitionRadiation, G4NeutronKiller, G4UnknownDecay, G4AnnihiToMuPair, G4GammaConversionToMuons, G4ForwardXrayTR, G4VXTRenergyLoss, G4HadronicProcess, G4QAtomicElectronScattering, G4QCoherentChargeExchange, G4QDiffraction, G4QDiscProcessMixer, G4QElastic, G4QInelastic, G4QIonIonElastic, G4QLowEnergy, G4QNGamma, G4HadronElasticProcess, G4WHadronElasticProcess, G4OpAbsorption, G4OpBoundaryProcess, G4OpMieHG, G4OpRayleigh, G4OpWLS, G4SynchrotronRadiation, G4SynchrotronRadiationInMat, G4VTransitionRadiation, and G4QSynchRad.

Definition at line 112 of file G4VDiscreteProcess.cc.

116{
117// clear NumberOfInteractionLengthLeft
119
120 return pParticleChange;
121}
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Referenced by G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QDiscProcessMixer::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), and G4QSynchRad::PostStepDoIt().

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4ErrorMagFieldLimitProcess, G4ErrorStepLengthLimitProcess, G4ErrorTrackLengthTarget, G4UnknownDecay, G4eplusPolarizedAnnihilation, G4PolarizedCompton, G4VEmProcess, G4QDiscProcessMixer, G4HadronStoppingProcess, G4NeutronKiller, and G4VErrorLimitProcess.

Definition at line 71 of file G4VDiscreteProcess.cc.

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

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