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

#include <G4ErrorTrackLengthTarget.hh>

+ Inheritance diagram for G4ErrorTrackLengthTarget:

Public Member Functions

 G4ErrorTrackLengthTarget (const G4double maxTrkLength)
 
virtual ~G4ErrorTrackLengthTarget ()
 
virtual G4double GetDistanceFromPoint (const G4ThreeVector &, const G4ThreeVector &) const
 
virtual G4double GetDistanceFromPoint (const G4ThreeVector &) const
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetMeanFreePath (const class G4Track &track, G4double, G4ForceCondition *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual void Dump (const G4String &msg) const
 
- 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
 
- Public Member Functions inherited from G4ErrorTarget
 G4ErrorTarget ()
 
virtual ~G4ErrorTarget ()
 
virtual G4double GetDistanceFromPoint (const G4ThreeVector &, const G4ThreeVector &) const
 
virtual G4double GetDistanceFromPoint (const G4ThreeVector &) const
 
virtual G4bool TargetReached (const G4Step *)
 
virtual void Dump (const G4String &msg) const =0
 
G4ErrorTargetType GetType () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VDiscreteProcess
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
 
- Protected Attributes inherited from G4ErrorTarget
G4ErrorTargetType theType
 

Detailed Description

Definition at line 54 of file G4ErrorTrackLengthTarget.hh.

Constructor & Destructor Documentation

◆ G4ErrorTrackLengthTarget()

G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget ( const G4double  maxTrkLength)

Definition at line 46 of file G4ErrorTrackLengthTarget.cc.

48 : G4VDiscreteProcess ("G4ErrorTrackLengthTarget"),
49 theMaximumTrackLength( maxTrkLength )
50{
52
53 G4ParticleTable::G4PTblDicIterator* theParticleIterator
55
56 // loop over all particles in G4ParticleTable
57
58 theParticleIterator->reset();
59 while( (*theParticleIterator)() )
60 {
61 G4ParticleDefinition* particle = theParticleIterator->value();
62 G4ProcessManager* pmanager = particle->GetProcessManager();
63 if (!particle->IsShortLived())
64 {
65 // Add transportation process for all particles other than "shortlived"
66 if ( pmanager == 0)
67 {
68 // Error !! no process manager
69 G4String particleName = particle->GetParticleName();
70 G4Exception("G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget",
71 "No process manager", RunMustBeAborted, particleName );
72 }
73 else
74 {
75 G4ProcessVector* procvec = pmanager->GetProcessList();
76 size_t isiz = procvec->size();
77
78 for( size_t ii=0; ii < isiz; ii++ )
79 {
80 if( ((*procvec)[ii])->GetProcessName() == "G4ErrorTrackLengthTarget")
81 {
82 pmanager->RemoveProcess( (*procvec)[ii] );
83 }
84 }
85 pmanager ->AddDiscreteProcess(this,4);
86 isiz = procvec->size();
87 }
88 }
89 else
90 {
91 // shortlived particle case
92 }
93 }
94}
@ G4ErrorTarget_TrkL
@ RunMustBeAborted
G4ErrorTargetType theType
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4PTblDicIterator * GetIterator()
static G4ParticleTable * GetParticleTable()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ProcessVector * GetProcessList() const
G4VProcess * RemoveProcess(G4VProcess *aProcess)
G4int size() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4ErrorTrackLengthTarget()

virtual G4ErrorTrackLengthTarget::~G4ErrorTrackLengthTarget ( )
inlinevirtual

Definition at line 61 of file G4ErrorTrackLengthTarget.hh.

61{}

Member Function Documentation

◆ Dump()

void G4ErrorTrackLengthTarget::Dump ( const G4String msg) const
virtual

Implements G4ErrorTarget.

Definition at line 132 of file G4ErrorTrackLengthTarget.cc.

133{
134 G4cout << msg << "G4ErrorTrackLengthTarget: max track length = "
135 << theMaximumTrackLength << G4endl;
136}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ GetDistanceFromPoint() [1/2]

virtual G4double G4ErrorTrackLengthTarget::GetDistanceFromPoint ( const G4ThreeVector ) const
inlinevirtual

Reimplemented from G4ErrorTarget.

Definition at line 69 of file G4ErrorTrackLengthTarget.hh.

70 { return DBL_MAX; }
#define DBL_MAX
Definition: templates.hh:83

◆ GetDistanceFromPoint() [2/2]

virtual G4double G4ErrorTrackLengthTarget::GetDistanceFromPoint ( const G4ThreeVector ,
const G4ThreeVector  
) const
inlinevirtual

Reimplemented from G4ErrorTarget.

Definition at line 66 of file G4ErrorTrackLengthTarget.hh.

68 { return DBL_MAX; }

◆ GetMeanFreePath()

G4double G4ErrorTrackLengthTarget::GetMeanFreePath ( const class G4Track track,
G4double  ,
G4ForceCondition  
)
virtual

Definition at line 108 of file G4ErrorTrackLengthTarget.cc.

110{
111#ifdef G4VERBOSE
113 {
114 G4cout << " G4ErrorTrackLengthTarget::GetMeanFreePath "
115 << theMaximumTrackLength - track.GetTrackLength() << G4endl;
116 }
117#endif
118
119 return theMaximumTrackLength - track.GetTrackLength();
120}

Referenced by PostStepGetPhysicalInteractionLength().

◆ PostStepDoIt()

G4VParticleChange * G4ErrorTrackLengthTarget::PostStepDoIt ( const G4Track aTrack,
const G4Step  
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 123 of file G4ErrorTrackLengthTarget.cc.

125{
126 theParticleChange.Initialize(aTrack);
127 return &theParticleChange;
128}
virtual void Initialize(const G4Track &)

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 98 of file G4ErrorTrackLengthTarget.cc.

101{
103 return GetMeanFreePath( track, 0., condition );
104}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
virtual G4double GetMeanFreePath(const class G4Track &track, G4double, G4ForceCondition *)

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