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

#include <G4StackChecker.hh>

+ Inheritance diagram for G4StackChecker:

Public Member Functions

 G4StackChecker ()=default
 
 ~G4StackChecker () override=default
 
G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *track) override
 
- Public Member Functions inherited from G4UserStackingAction
 G4UserStackingAction ()
 
virtual ~G4UserStackingAction ()=default
 
void SetStackManager (G4StackManager *value)
 
virtual void NewStage ()
 
virtual void PrepareNewEvent ()
 

Additional Inherited Members

- Protected Attributes inherited from G4UserStackingAction
G4StackManagerstackManager = nullptr
 

Detailed Description

Definition at line 43 of file G4StackChecker.hh.

Constructor & Destructor Documentation

◆ G4StackChecker()

G4StackChecker::G4StackChecker ( )
default

◆ ~G4StackChecker()

G4StackChecker::~G4StackChecker ( )
overridedefault

Member Function Documentation

◆ ClassifyNewTrack()

G4ClassificationOfNewTrack G4StackChecker::ClassifyNewTrack ( const G4Track * track)
overridevirtual

Reimplemented from G4UserStackingAction.

Definition at line 37 of file G4StackChecker.cc.

38{
40 G4double e = track->GetKineticEnergy();
41 if ( (!(e < 0.0) && !(e > 0.0) && !(e == 0.0))
42 || track->GetMomentumDirection() == nullDirection)
43 {
44 result = fKill;
45 G4String nam = track->GetDefinition()->GetParticleName();
46 G4cout << "### G4StackChecker: event# "
47 << (G4EventManager::GetEventManager())->GetConstCurrentEvent()->GetEventID()
48 << " unacceptable " << nam << " is killed in the stack" << G4endl;
49 G4cout << "### " << nam << " have been produced by the process "
50 << track->GetCreatorProcess()->GetProcessName()
51 << " trackID= " << track->GetTrackID()
52 << " parentID= " << track->GetParentID()
53 << G4endl;
54 G4cout << "### E= " << track->GetKineticEnergy()
55 << " position= " << track->GetPosition()
56 << " direction= " << track->GetMomentumDirection()
57 << " time= " << track->GetGlobalTime()
58 << G4endl;
59 }
60 return result;
61}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4EventManager * GetEventManager()
const G4String & GetParticleName() const
G4int GetTrackID() const
const G4VProcess * GetCreatorProcess() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetParentID() const
const G4String & GetProcessName() const

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