Geant4 11.1.1
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 &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
 
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 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 &)
 
- Public Member Functions inherited from G4ErrorTarget
 G4ErrorTarget ()
 
virtual ~G4ErrorTarget ()=default
 
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 prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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
 
- Protected Attributes inherited from G4ErrorTarget
G4ErrorTargetType theType
 

Detailed Description

Definition at line 53 of file G4ErrorTrackLengthTarget.hh.

Constructor & Destructor Documentation

◆ G4ErrorTrackLengthTarget()

G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget ( const G4double  maxTrkLength)

Definition at line 45 of file G4ErrorTrackLengthTarget.cc.

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

◆ ~G4ErrorTrackLengthTarget()

virtual G4ErrorTrackLengthTarget::~G4ErrorTrackLengthTarget ( )
inlinevirtual

Definition at line 60 of file G4ErrorTrackLengthTarget.hh.

60{}

Member Function Documentation

◆ Dump()

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

Implements G4ErrorTarget.

Definition at line 126 of file G4ErrorTrackLengthTarget.cc.

127{
128 G4cout << msg << "G4ErrorTrackLengthTarget: max track length = "
129 << theMaximumTrackLength << G4endl;
130}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ GetDistanceFromPoint() [1/2]

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

Reimplemented from G4ErrorTarget.

Definition at line 70 of file G4ErrorTrackLengthTarget.hh.

71 {
72 return DBL_MAX;
73 }
#define DBL_MAX
Definition: templates.hh:62

◆ GetDistanceFromPoint() [2/2]

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

Reimplemented from G4ErrorTarget.

Definition at line 65 of file G4ErrorTrackLengthTarget.hh.

67 {
68 return DBL_MAX;
69 }

◆ GetMeanFreePath()

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

Definition at line 103 of file G4ErrorTrackLengthTarget.cc.

106{
107#ifdef G4VERBOSE
109 {
110 G4cout << " G4ErrorTrackLengthTarget::GetMeanFreePath "
111 << theMaximumTrackLength - track.GetTrackLength() << G4endl;
112 }
113#endif
114
115 return theMaximumTrackLength - track.GetTrackLength();
116}

Referenced by PostStepGetPhysicalInteractionLength().

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 118 of file G4ErrorTrackLengthTarget.cc.

120{
121 theParticleChange.Initialize(aTrack);
122 return &theParticleChange;
123}
virtual void Initialize(const G4Track &)

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 95 of file G4ErrorTrackLengthTarget.cc.

97{
99 return GetMeanFreePath(track, 0., condition);
100}
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: