Geant4-11
Public Member Functions | Protected Types | Protected Attributes | Private Member Functions | Private Attributes
G4VITSteppingVerbose Class Referenceabstract

#include <G4VITSteppingVerbose.hh>

Inheritance diagram for G4VITSteppingVerbose:
G4UImessenger G4ITSteppingVerbose

Public Member Functions

virtual void AlongStepDoItAllDone ()=0
 
virtual void AlongStepDoItOneByOne ()=0
 
virtual void AtRestDoItInvoked ()=0
 
virtual void AtRestDoItOneByOne ()=0
 
void CopyState ()
 
virtual void DoItStarted ()=0
 
virtual void DPSLAlongStep ()=0
 
virtual void DPSLPostStep ()=0
 
virtual void DPSLStarted ()=0
 
virtual void DPSLUserLimit ()=0
 
 G4VITSteppingVerbose ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4int GetVerbose ()
 
virtual void NewStep ()=0
 
virtual void PostStepDoItAllDone ()=0
 
virtual void PostStepDoItOneByOne ()=0
 
virtual void PostStepVerbose (G4Track *)=0
 
virtual void PreStepVerbose (G4Track *)=0
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
void SetStepProcessor (const G4ITStepProcessor *stepProcessor)
 
void SetVerbose (int flag)
 
virtual void StepInfo ()=0
 
virtual void StepInfoForLeadingTrack ()=0
 
void TrackBanner (G4Track *track, const G4String &message)
 
virtual void TrackingEnded (G4Track *track)
 
virtual void TrackingStarted (G4Track *track)
 
virtual void VerboseParticleChange ()=0
 
virtual void VerboseTrack ()=0
 
virtual ~G4VITSteppingVerbose ()
 

Protected Types

typedef std::vector< G4intG4SelectedAlongStepDoItVector
 
typedef std::vector< G4intG4SelectedAtRestDoItVector
 
typedef std::vector< G4intG4SelectedPostStepDoItVector
 

Protected Attributes

G4ProcessVectorfAlongStepDoItVector
 
G4ProcessVectorfAlongStepGetPhysIntVector
 
size_t fAtRestDoItProcTriggered
 
G4ProcessVectorfAtRestDoItVector
 
G4ProcessVectorfAtRestGetPhysIntVector
 
G4ForceCondition fCondition
 
const G4VITProcessfCurrentProcess
 
const G4VPhysicalVolumefCurrentVolume
 
G4GPILSelection fGPILSelection
 
G4int fN2ndariesAlongStepDoIt
 
G4int fN2ndariesAtRestDoIt
 
G4int fN2ndariesPostStepDoIt
 
const G4VParticleChangefParticleChange
 
size_t fPostStepDoItProcTriggered
 
G4ProcessVectorfPostStepDoItVector
 
G4ProcessVectorfPostStepGetPhysIntVector
 
G4StepPointfPostStepPoint
 
const ProcessGeneralInfofpProcessGeneralInfo
 
G4StepPointfPreStepPoint
 
G4double fPreviousStepSize
 
G4ITStepProcessorStatefpState
 
const G4ITStepProcessorfpStepProcessor
 
G4UIcmdWithAnIntegerfpVerboseUI
 
const G4TrackVectorfSecondary
 
G4SelectedAtRestDoItVectorfSelectedAtRestDoItVector
 
G4SelectedPostStepDoItVectorfSelectedPostStepDoItVector
 
const G4StepfStep
 
G4StepStatus fStepStatus
 
G4TouchableHandle fTouchableHandle
 
const G4TrackfTrack
 
G4int fVerboseLevel
 
size_t MAXofAlongStepLoops
 
size_t MAXofAtRestLoops
 
size_t MAXofPostStepLoops
 
G4double PhysicalStep
 
G4double physIntLength
 

Private Member Functions

void AddUIcommand (G4UIcommand *newCommand)
 
G4String BtoS (G4bool b)
 
G4bool CommandsShouldBeInMaster () const
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
G4String DtoS (G4double a)
 
G4String ItoS (G4int i)
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool StoB (G4String s)
 
G4double StoD (G4String s)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 

Private Attributes

G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 67 of file G4VITSteppingVerbose.hh.

Member Typedef Documentation

◆ G4SelectedAlongStepDoItVector

Definition at line 183 of file G4VITSteppingVerbose.hh.

◆ G4SelectedAtRestDoItVector

Definition at line 182 of file G4VITSteppingVerbose.hh.

◆ G4SelectedPostStepDoItVector

Definition at line 184 of file G4VITSteppingVerbose.hh.

Constructor & Destructor Documentation

◆ G4VITSteppingVerbose()

G4VITSteppingVerbose::G4VITSteppingVerbose ( )

Definition at line 38 of file G4VITSteppingVerbose.cc.

39{
41 fpState = 0;
43
44 PhysicalStep = -1;
46
48 fTrack = 0;
49 fSecondary = 0;
50 fStep = 0;
51 fPreStepPoint = 0;
53
55 // fSensitive = fpStepProcessor->GetfSensitive();
57
61
65
69
72
76
77 // fNavigator = fpStepProcessor->GetfNavigator();
78
79 fVerboseLevel = 0;
80 fpVerboseUI = new G4UIcmdWithAnInteger("/chem/tracking/verbose", this);
81
84
86
88
89 // StepControlFlag = fpStepProcessor->GetStepControlFlag();
90
91 physIntLength = 0;
94
95}
@ InActivated
@ NotCandidateForSelection
@ fUndefined
Definition: G4StepStatus.hh:55
G4ITStepProcessorState * fpState
const G4ITStepProcessor * fpStepProcessor
G4ProcessVector * fAtRestGetPhysIntVector
G4ProcessVector * fPostStepDoItVector
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector
G4UIcmdWithAnInteger * fpVerboseUI
G4ProcessVector * fPostStepGetPhysIntVector
G4ProcessVector * fAtRestDoItVector
G4TouchableHandle fTouchableHandle
const G4VPhysicalVolume * fCurrentVolume
const ProcessGeneralInfo * fpProcessGeneralInfo
G4ProcessVector * fAlongStepGetPhysIntVector
const G4VITProcess * fCurrentProcess
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector
const G4VParticleChange * fParticleChange
const G4TrackVector * fSecondary
G4GPILSelection fGPILSelection
G4ProcessVector * fAlongStepDoItVector

References fAlongStepDoItVector, fAlongStepGetPhysIntVector, fAtRestDoItProcTriggered, fAtRestDoItVector, fAtRestGetPhysIntVector, fCondition, fCurrentProcess, fCurrentVolume, fGPILSelection, fN2ndariesAlongStepDoIt, fN2ndariesAtRestDoIt, fN2ndariesPostStepDoIt, fParticleChange, fPostStepDoItProcTriggered, fPostStepDoItVector, fPostStepGetPhysIntVector, fPostStepPoint, fpProcessGeneralInfo, fPreStepPoint, fPreviousStepSize, fpState, fpStepProcessor, fpVerboseUI, fSecondary, fSelectedAtRestDoItVector, fSelectedPostStepDoItVector, fStep, fStepStatus, fTouchableHandle, fTrack, fUndefined, fVerboseLevel, InActivated, MAXofAlongStepLoops, MAXofAtRestLoops, MAXofPostStepLoops, NotCandidateForSelection, PhysicalStep, and physIntLength.

◆ ~G4VITSteppingVerbose()

G4VITSteppingVerbose::~G4VITSteppingVerbose ( )
virtual

Definition at line 99 of file G4VITSteppingVerbose.cc.

100{
101 if(fpVerboseUI) delete fpVerboseUI;
102}

References fpVerboseUI.

Member Function Documentation

◆ AddUIcommand()

void G4UImessenger::AddUIcommand ( G4UIcommand newCommand)
protectedinherited

Definition at line 149 of file G4UImessenger.cc.

150{
151 G4cerr << "Warning : Old style definition of G4UIcommand <"
152 << newCommand->GetCommandPath() << ">." << G4endl;
153}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136

References G4cerr, G4endl, and G4UIcommand::GetCommandPath().

◆ AlongStepDoItAllDone()

virtual void G4VITSteppingVerbose::AlongStepDoItAllDone ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ AlongStepDoItOneByOne()

virtual void G4VITSteppingVerbose::AlongStepDoItOneByOne ( )
pure virtual

◆ AtRestDoItInvoked()

virtual void G4VITSteppingVerbose::AtRestDoItInvoked ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ AtRestDoItOneByOne()

virtual void G4VITSteppingVerbose::AtRestDoItOneByOne ( )
pure virtual

Implemented in G4ITSteppingVerbose.

◆ BtoS()

G4String G4UImessenger::BtoS ( G4bool  b)
protectedinherited

Definition at line 98 of file G4UImessenger.cc.

99{
100 G4String vl = "0";
101 if(b)
102 vl = "true";
103 return vl;
104}

◆ CommandsShouldBeInMaster()

G4bool G4UImessenger::CommandsShouldBeInMaster ( ) const
inlineinherited

Definition at line 77 of file G4UImessenger.hh.

78 {
80 }
G4bool commandsShouldBeInMaster

References G4UImessenger::commandsShouldBeInMaster.

Referenced by G4UIcommand::G4UIcommandCommonConstructorCode().

◆ CopyState()

void G4VITSteppingVerbose::CopyState ( )

Definition at line 106 of file G4VITSteppingVerbose.cc.

107{
108
110 else
111 {
113 }
114
116
119
126
128// fSensitive = fpStepProcessor->GetfSensitive();
130
134
139
143
146
150
151// fNavigator = fpStepProcessor->GetfNavigator();
152
155
157
159
160// StepControlFlag = fpStepProcessor->GetStepControlFlag();
161
165}
G4TouchableHandle fTouchableHandle
G4SelectedPostStepDoItVector fSelectedPostStepDoItVector
G4SelectedAtRestDoItVector fSelectedAtRestDoItVector
G4int GetN2ndariesAtRestDoIt() const
size_t GetAtRestDoItProcTriggered() const
size_t GetPostStepDoItProcTriggered() const
G4int GetN2ndariesAlongStepDoIt() const
const ProcessGeneralInfo * GetCurrentProcessInfo() const
G4GPILSelection GetGPILSelection() const
G4ForceCondition GetCondition() const
G4int GetN2ndariesPostStepDoIt() const
G4TrackVector * GetSecondaries() const
G4double GetPhysIntLength() const
const G4ITStepProcessorState * GetProcessorState() const
const G4VPhysicalVolume * GetCurrentVolume() const
const G4VParticleChange * GetParticleChange() const
const G4VITProcess * GetCurrentProcess() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4ProcessVector * fpPostStepDoItVector
G4ProcessVector * fpPostStepGetPhysIntVector
G4ProcessVector * fpAlongStepGetPhysIntVector
G4ProcessVector * fpAlongStepDoItVector
G4ProcessVector * fpAtRestGetPhysIntVector
G4ProcessVector * fpAtRestDoItVector

References fAlongStepDoItVector, fAlongStepGetPhysIntVector, fAtRestDoItProcTriggered, fAtRestDoItVector, fAtRestGetPhysIntVector, fCondition, fCurrentProcess, fCurrentVolume, fGPILSelection, fN2ndariesAlongStepDoIt, fN2ndariesAtRestDoIt, fN2ndariesPostStepDoIt, ProcessGeneralInfo::fpAlongStepDoItVector, ProcessGeneralInfo::fpAlongStepGetPhysIntVector, fParticleChange, ProcessGeneralInfo::fpAtRestDoItVector, ProcessGeneralInfo::fpAtRestGetPhysIntVector, fPostStepDoItProcTriggered, fPostStepDoItVector, fPostStepGetPhysIntVector, fPostStepPoint, ProcessGeneralInfo::fpPostStepDoItVector, ProcessGeneralInfo::fpPostStepGetPhysIntVector, fpProcessGeneralInfo, fPreStepPoint, G4ITStepProcessorState::fPreviousStepSize, fPreviousStepSize, fpState, fpStepProcessor, fSecondary, G4ITStepProcessorState::fSelectedAtRestDoItVector, fSelectedAtRestDoItVector, G4ITStepProcessorState::fSelectedPostStepDoItVector, fSelectedPostStepDoItVector, fStep, G4ITStepProcessorState::fStepStatus, fStepStatus, G4ITStepProcessorState::fTouchableHandle, fTouchableHandle, fTrack, G4ITStepProcessor::GetAtRestDoItProcTriggered(), G4ITStepProcessor::GetCondition(), G4ITStepProcessor::GetCurrentProcess(), G4ITStepProcessor::GetCurrentProcessInfo(), G4ITStepProcessor::GetCurrentVolume(), G4ITStepProcessor::GetGPILSelection(), G4ITStepProcessor::GetN2ndariesAlongStepDoIt(), G4ITStepProcessor::GetN2ndariesAtRestDoIt(), G4ITStepProcessor::GetN2ndariesPostStepDoIt(), G4ITStepProcessor::GetParticleChange(), G4ITStepProcessor::GetPhysIntLength(), G4ITStepProcessor::GetPostStepDoItProcTriggered(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4ITStepProcessor::GetProcessorState(), G4ITStepProcessor::GetSecondaries(), G4ITStepProcessor::GetStep(), G4ITStepProcessor::GetTrack(), ProcessGeneralInfo::MAXofAlongStepLoops, MAXofAlongStepLoops, ProcessGeneralInfo::MAXofAtRestLoops, MAXofAtRestLoops, ProcessGeneralInfo::MAXofPostStepLoops, MAXofPostStepLoops, PhysicalStep, and physIntLength.

Referenced by G4ITSteppingVerbose::AlongStepDoItAllDone(), G4ITSteppingVerbose::AlongStepDoItOneByOne(), G4ITSteppingVerbose::AtRestDoItInvoked(), G4ITSteppingVerbose::AtRestDoItOneByOne(), G4ITSteppingVerbose::DPSLAlongStep(), G4ITSteppingVerbose::DPSLPostStep(), G4ITSteppingVerbose::DPSLStarted(), G4ITSteppingVerbose::DPSLUserLimit(), G4ITSteppingVerbose::PostStepDoItAllDone(), G4ITSteppingVerbose::PostStepDoItOneByOne(), G4ITSteppingVerbose::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), and G4ITSteppingVerbose::VerboseTrack().

◆ CreateCommand()

template<typename T >
T * G4UImessenger::CreateCommand ( const G4String cname,
const G4String dsc 
)
protectedinherited

Definition at line 110 of file G4UImessenger.hh.

111{
112 G4String path;
113 if(cname[0] != '/')
114 {
115 path = baseDirName + cname;
116 if(path[0] != '/')
117 path = "/" + path;
118 }
119
120 T* command = new T(path.c_str(), this);
121 command->SetGuidance(dsc.c_str());
122
123 return command;
124}
G4String baseDirName

References G4UImessenger::baseDirName.

◆ CreateDirectory()

void G4UImessenger::CreateDirectory ( const G4String path,
const G4String dsc,
G4bool  commandsToBeBroadcasted = true 
)
protectedinherited

Definition at line 156 of file G4UImessenger.cc.

158{
160
161 G4String fullpath = path;
162 if(fullpath.back() != '/')
163 fullpath.append("/");
164
165 G4UIcommandTree* tree = ui->GetTree()->FindCommandTree(fullpath.c_str());
166 if(tree != nullptr)
167 {
168 baseDirName = tree->GetPathName();
169 }
170 else
171 {
172 baseDir = new G4UIdirectory(fullpath.c_str(), commandsToBeBroadcasted);
173 baseDirName = fullpath;
174 baseDir->SetGuidance(dsc.c_str());
175 }
176}
const G4String & GetPathName() const
G4UIcommandTree * FindCommandTree(const char *commandPath)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4UIdirectory * baseDir

References G4UImessenger::baseDir, G4UImessenger::baseDirName, G4UIcommandTree::FindCommandTree(), G4UIcommandTree::GetPathName(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), and G4UIcommand::SetGuidance().

Referenced by G4MoleculeShootMessenger::G4MoleculeShootMessenger(), and G4UImessenger::G4UImessenger().

◆ DoItStarted()

virtual void G4VITSteppingVerbose::DoItStarted ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoIt().

◆ DPSLAlongStep()

virtual void G4VITSteppingVerbose::DPSLAlongStep ( )
pure virtual

◆ DPSLPostStep()

virtual void G4VITSteppingVerbose::DPSLPostStep ( )
pure virtual

◆ DPSLStarted()

virtual void G4VITSteppingVerbose::DPSLStarted ( )
pure virtual

◆ DPSLUserLimit()

virtual void G4VITSteppingVerbose::DPSLUserLimit ( )
pure virtual

Implemented in G4ITSteppingVerbose.

◆ DtoS()

G4String G4UImessenger::DtoS ( G4double  a)
protectedinherited

Definition at line 90 of file G4UImessenger.cc.

91{
92 std::ostringstream os;
93 os << a;
94 return G4String(os.str());
95}

Referenced by G4ScoreQuantityMessenger::FilterCommands(), and G4UIcontrolMessenger::SetNewValue().

◆ GetCurrentValue()

G4String G4VITSteppingVerbose::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 179 of file G4VITSteppingVerbose.cc.

180{
181 return command->ConvertToString(fVerboseLevel);
182}
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:445

References G4UIcommand::ConvertToString(), and fVerboseLevel.

◆ GetVerbose()

G4int G4VITSteppingVerbose::GetVerbose ( )
inline

Definition at line 112 of file G4VITSteppingVerbose.hh.

113 {
114 return fVerboseLevel;
115 }

References fVerboseLevel.

Referenced by G4ITTrackingInteractivity::GetSteppingVerboseLevel().

◆ ItoS()

G4String G4UImessenger::ItoS ( G4int  i)
protectedinherited

Definition at line 82 of file G4UImessenger.cc.

83{
84 std::ostringstream os;
85 os << i;
86 return G4String(os.str());
87}

Referenced by G4GenericMessenger::DeclareMethod(), and G4ParticleGunMessenger::GetCurrentValue().

◆ NewStep()

virtual void G4VITSteppingVerbose::NewStep ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ operator!=()

G4bool G4UImessenger::operator!= ( const G4UImessenger messenger) const
inherited

Definition at line 76 of file G4UImessenger.cc.

77{
78 return this != &messenger;
79}

◆ operator==()

G4bool G4UImessenger::operator== ( const G4UImessenger messenger) const
inherited

Definition at line 70 of file G4UImessenger.cc.

71{
72 return this == &messenger;
73}

◆ PostStepDoItAllDone()

virtual void G4VITSteppingVerbose::PostStepDoItAllDone ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ PostStepDoItOneByOne()

virtual void G4VITSteppingVerbose::PostStepDoItOneByOne ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::InvokePSDIP().

◆ PostStepVerbose()

virtual void G4VITSteppingVerbose::PostStepVerbose ( G4Track )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ PreStepVerbose()

virtual void G4VITSteppingVerbose::PreStepVerbose ( G4Track )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ SetNewValue()

void G4VITSteppingVerbose::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 169 of file G4VITSteppingVerbose.cc.

170{
171 if(command == fpVerboseUI)
172 {
174 }
175}
static G4int GetNewIntValue(const char *paramString)

References fpVerboseUI, fVerboseLevel, and G4UIcmdWithAnInteger::GetNewIntValue().

◆ SetStepProcessor()

void G4VITSteppingVerbose::SetStepProcessor ( const G4ITStepProcessor stepProcessor)
inline

Definition at line 127 of file G4VITSteppingVerbose.hh.

128 {
129 this->fpStepProcessor = stepProcessor;
130 }

References fpStepProcessor.

Referenced by G4ITStepProcessor::Initialize().

◆ SetVerbose()

void G4VITSteppingVerbose::SetVerbose ( int  flag)
inline

Definition at line 107 of file G4VITSteppingVerbose.hh.

108 {
109 fVerboseLevel = flag;
110 }

References fVerboseLevel.

Referenced by G4ITTrackingInteractivity::SetSteppingVerboseLevel().

◆ StepInfo()

virtual void G4VITSteppingVerbose::StepInfo ( )
pure virtual

Implemented in G4ITSteppingVerbose.

◆ StepInfoForLeadingTrack()

virtual void G4VITSteppingVerbose::StepInfoForLeadingTrack ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Referenced by G4ITStepProcessor::DoStepping().

◆ StoB()

G4bool G4UImessenger::StoB ( G4String  s)
protectedinherited

Definition at line 137 of file G4UImessenger.cc.

138{
140 G4bool vl = false;
141 if(v == "Y" || v == "YES" || v == "1" || v == "T" || v == "TRUE")
142 {
143 vl = true;
144 }
145 return vl;
146}
bool G4bool
Definition: G4Types.hh:86
G4String to_upper_copy(G4String str)
Return uppercase copy of string.

References G4StrUtil::to_upper_copy().

Referenced by G4LocalThreadCoutMessenger::SetNewValue(), G4CascadeParamMessenger::SetNewValue(), G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ StoD()

G4double G4UImessenger::StoD ( G4String  s)
protectedinherited

◆ StoI()

G4int G4UImessenger::StoI ( G4String  s)
protectedinherited

◆ StoL()

G4long G4UImessenger::StoL ( G4String  s)
protectedinherited

Definition at line 117 of file G4UImessenger.cc.

118{
119 G4long vl;
120 const char* t = str;
121 std::istringstream is(t);
122 is >> vl;
123 return vl;
124}
long G4long
Definition: G4Types.hh:87

Referenced by G4RunMessenger::SetNewValue().

◆ TrackBanner()

void G4VITSteppingVerbose::TrackBanner ( G4Track track,
const G4String message 
)

Definition at line 219 of file G4VITSteppingVerbose.cc.

220{
221 G4cout << G4endl;
222 G4cout << "*******************************************************"
223 << "**************************************************"
224 << G4endl;
225 if(message != "")
226 {
227 G4cout << message;
228 }
229 G4cout << " * G4Track Information: "
230 << " Particle : " << track->GetDefinition()->GetParticleName()
231 << ","
232 << " Track ID : " << track->GetTrackID()
233 << ","
234 << " Parent ID : " << track->GetParentID()
235 << G4endl;
236 G4cout << "*******************************************************"
237 << "**************************************************"
238 << G4endl;
239 G4cout << G4endl;
240}
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleName() const
G4int GetTrackID() const
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const

References G4cout, G4endl, G4Track::GetDefinition(), G4Track::GetParentID(), G4ParticleDefinition::GetParticleName(), and G4Track::GetTrackID().

Referenced by TrackingEnded(), and TrackingStarted().

◆ TrackingEnded()

void G4VITSteppingVerbose::TrackingEnded ( G4Track track)
virtual

Reimplemented in G4ITSteppingVerbose.

Definition at line 203 of file G4VITSteppingVerbose.cc.

208{
209#ifdef G4VERBOSE
210 if(fVerboseLevel > 0)
211 {
212 TrackBanner(track, "G4ITTrackingManager::EndTracking : ");
213 }
214#endif
215}
void TrackBanner(G4Track *track, const G4String &message)

References fVerboseLevel, and TrackBanner().

Referenced by G4ITTrackingManager::EndTracking(), and G4ITTrackingManager::EndTrackingWOKill().

◆ TrackingStarted()

void G4VITSteppingVerbose::TrackingStarted ( G4Track track)
virtual

Reimplemented in G4ITSteppingVerbose.

Definition at line 186 of file G4VITSteppingVerbose.cc.

191{
192#ifdef G4VERBOSE
193 if(fVerboseLevel > 0)
194 {
195 TrackBanner(track, "G4ITTrackingManager::StartTracking : ");
196 }
197#endif
198
199}

References fVerboseLevel, and TrackBanner().

Referenced by G4ITTrackingManager::StartTracking().

◆ VerboseParticleChange()

virtual void G4VITSteppingVerbose::VerboseParticleChange ( )
pure virtual

Implemented in G4ITSteppingVerbose.

◆ VerboseTrack()

virtual void G4VITSteppingVerbose::VerboseTrack ( )
pure virtual

Implemented in G4ITSteppingVerbose.

Field Documentation

◆ baseDir

G4UIdirectory* G4UImessenger::baseDir = nullptr
protectedinherited

◆ baseDirName

G4String G4UImessenger::baseDirName = ""
protectedinherited

◆ commandsShouldBeInMaster

G4bool G4UImessenger::commandsShouldBeInMaster = false
protectedinherited

◆ fAlongStepDoItVector

G4ProcessVector* G4VITSteppingVerbose::fAlongStepDoItVector
protected

Definition at line 158 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fAlongStepGetPhysIntVector

G4ProcessVector* G4VITSteppingVerbose::fAlongStepGetPhysIntVector
protected

Definition at line 162 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fAtRestDoItProcTriggered

size_t G4VITSteppingVerbose::fAtRestDoItProcTriggered
protected

Definition at line 169 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fAtRestDoItVector

G4ProcessVector* G4VITSteppingVerbose::fAtRestDoItVector
protected

Definition at line 157 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fAtRestGetPhysIntVector

G4ProcessVector* G4VITSteppingVerbose::fAtRestGetPhysIntVector
protected

Definition at line 161 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fCondition

G4ForceCondition G4VITSteppingVerbose::fCondition
protected

◆ fCurrentProcess

const G4VITProcess* G4VITSteppingVerbose::fCurrentProcess
protected

◆ fCurrentVolume

const G4VPhysicalVolume* G4VITSteppingVerbose::fCurrentVolume
protected

Definition at line 151 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fGPILSelection

G4GPILSelection G4VITSteppingVerbose::fGPILSelection
protected

◆ fN2ndariesAlongStepDoIt

G4int G4VITSteppingVerbose::fN2ndariesAlongStepDoIt
protected

◆ fN2ndariesAtRestDoIt

G4int G4VITSteppingVerbose::fN2ndariesAtRestDoIt
protected

◆ fN2ndariesPostStepDoIt

G4int G4VITSteppingVerbose::fN2ndariesPostStepDoIt
protected

◆ fParticleChange

const G4VParticleChange* G4VITSteppingVerbose::fParticleChange
protected

◆ fPostStepDoItProcTriggered

size_t G4VITSteppingVerbose::fPostStepDoItProcTriggered
protected

Definition at line 170 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fPostStepDoItVector

G4ProcessVector* G4VITSteppingVerbose::fPostStepDoItVector
protected

Definition at line 159 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fPostStepGetPhysIntVector

G4ProcessVector* G4VITSteppingVerbose::fPostStepGetPhysIntVector
protected

Definition at line 163 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fPostStepPoint

G4StepPoint* G4VITSteppingVerbose::fPostStepPoint
protected

Definition at line 149 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fpProcessGeneralInfo

const ProcessGeneralInfo* G4VITSteppingVerbose::fpProcessGeneralInfo
protected

Definition at line 139 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fPreStepPoint

G4StepPoint* G4VITSteppingVerbose::fPreStepPoint
protected

Definition at line 148 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fPreviousStepSize

G4double G4VITSteppingVerbose::fPreviousStepSize
protected

Definition at line 188 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fpState

G4ITStepProcessorState* G4VITSteppingVerbose::fpState
protected

Definition at line 138 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fpStepProcessor

const G4ITStepProcessor* G4VITSteppingVerbose::fpStepProcessor
protected

Definition at line 135 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), G4VITSteppingVerbose(), and SetStepProcessor().

◆ fpVerboseUI

G4UIcmdWithAnInteger* G4VITSteppingVerbose::fpVerboseUI
protected

◆ fSecondary

const G4TrackVector* G4VITSteppingVerbose::fSecondary
protected

◆ fSelectedAtRestDoItVector

G4SelectedAtRestDoItVector* G4VITSteppingVerbose::fSelectedAtRestDoItVector
protected

◆ fSelectedPostStepDoItVector

G4SelectedPostStepDoItVector* G4VITSteppingVerbose::fSelectedPostStepDoItVector
protected

◆ fStep

const G4Step* G4VITSteppingVerbose::fStep
protected

◆ fStepStatus

G4StepStatus G4VITSteppingVerbose::fStepStatus
protected

◆ fTouchableHandle

G4TouchableHandle G4VITSteppingVerbose::fTouchableHandle
protected

Definition at line 190 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ fTrack

const G4Track* G4VITSteppingVerbose::fTrack
protected

◆ fVerboseLevel

G4int G4VITSteppingVerbose::fVerboseLevel
protected

◆ MAXofAlongStepLoops

size_t G4VITSteppingVerbose::MAXofAlongStepLoops
protected

◆ MAXofAtRestLoops

size_t G4VITSteppingVerbose::MAXofAtRestLoops
protected

◆ MAXofPostStepLoops

size_t G4VITSteppingVerbose::MAXofPostStepLoops
protected

◆ PhysicalStep

G4double G4VITSteppingVerbose::PhysicalStep
protected

Definition at line 141 of file G4VITSteppingVerbose.hh.

Referenced by CopyState(), and G4VITSteppingVerbose().

◆ physIntLength

G4double G4VITSteppingVerbose::physIntLength
protected

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