Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronKiller Class Reference

#include <G4NeutronKiller.hh>

+ Inheritance diagram for G4NeutronKiller:

Public Member Functions

 G4NeutronKiller (const G4String &processName="nKiller", G4ProcessType aType=fGeneral)
 
virtual ~G4NeutronKiller ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void SetTimeLimit (G4double)
 
void SetKinEnergyLimit (G4double)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *)
 
- Public Member Functions inherited from G4VDiscreteProcess
 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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 69 of file G4NeutronKiller.hh.

Constructor & Destructor Documentation

◆ G4NeutronKiller()

G4NeutronKiller::G4NeutronKiller ( const G4String processName = "nKiller",
G4ProcessType  aType = fGeneral 
)

Definition at line 49 of file G4NeutronKiller.cc.

50 : G4VDiscreteProcess(processName, aType)
51{
52 // set Process Sub Type
53 SetProcessSubType(static_cast<int>(NEUTRON_KILLER));
54
55 kinEnergyThreshold = 0.0;
56 timeThreshold = DBL_MAX;
57 pMess = new G4NeutronKillerMessenger(this);
58}
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
#define DBL_MAX
Definition: templates.hh:83

◆ ~G4NeutronKiller()

G4NeutronKiller::~G4NeutronKiller ( )
virtual

Definition at line 62 of file G4NeutronKiller.cc.

63{
64 delete pMess;
65}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4NeutronKiller::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition  
)
inlinevirtual

Implements G4VDiscreteProcess.

Definition at line 121 of file G4NeutronKiller.hh.

123{
124 return DBL_MAX;
125}

◆ IsApplicable()

G4bool G4NeutronKiller::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 69 of file G4NeutronKiller.cc.

70{
71 return (particle.GetParticleName() == "neutron");
72}
const G4String & GetParticleName() const

◆ PostStepDoIt()

G4VParticleChange * G4NeutronKiller::PostStepDoIt ( const G4Track aTrack,
const G4Step  
)
inlinevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 129 of file G4NeutronKiller.hh.

131{
134 return pParticleChange;
135}
@ fStopAndKill
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

◆ PostStepGetPhysicalInteractionLength()

G4double G4NeutronKiller::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
inlinevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 106 of file G4NeutronKiller.hh.

109{
110 // condition is set to "Not Forced"
112
113 G4double limit = DBL_MAX;
114 if(aTrack.GetGlobalTime() > timeThreshold ||
115 aTrack.GetKineticEnergy() < kinEnergyThreshold) limit = 0.0;
116 return limit;
117}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition: G4Types.hh:64

◆ SetKinEnergyLimit()

void G4NeutronKiller::SetKinEnergyLimit ( G4double  val)

Definition at line 86 of file G4NeutronKiller.cc.

87{
88 kinEnergyThreshold = val;
89 if(verboseLevel > 0)
90 G4cout << "### G4NeutronKiller: Tracking cut E(MeV) = "
91 << kinEnergyThreshold/MeV << G4endl;
92}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4int verboseLevel
Definition: G4VProcess.hh:368

Referenced by G4NeutronTrackingCut::ConstructProcess(), and G4NeutronKillerMessenger::SetNewValue().

◆ SetTimeLimit()

void G4NeutronKiller::SetTimeLimit ( G4double  val)

Definition at line 76 of file G4NeutronKiller.cc.

77{
78 timeThreshold = val;
79 if(verboseLevel > 0)
80 G4cout << "### G4NeutronKiller: timeLimit(ns) = "
81 << timeThreshold/ns << G4endl;
82}
#define ns
Definition: xmlparse.cc:597

Referenced by G4NeutronTrackingCut::ConstructProcess(), and G4NeutronKillerMessenger::SetNewValue().


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