Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SteppingVerbose Class Reference

#include <G4SteppingVerbose.hh>

Inheritance diagram for G4SteppingVerbose:
Collaboration diagram for G4SteppingVerbose:

Public Member Functions

 G4SteppingVerbose ()
 
 ~G4SteppingVerbose ()
 
void NewStep ()
 
void AtRestDoItInvoked ()
 
void AlongStepDoItAllDone ()
 
void PostStepDoItAllDone ()
 
void AlongStepDoItOneByOne ()
 
void PostStepDoItOneByOne ()
 
void StepInfo ()
 
void TrackingStarted ()
 
void DPSLStarted ()
 
void DPSLUserLimit ()
 
void DPSLPostStep ()
 
void DPSLAlongStep ()
 
void VerboseTrack ()
 
void VerboseParticleChange ()
 
void ShowStep () const
 
- Public Member Functions inherited from G4VSteppingVerbose
virtual ~G4VSteppingVerbose ()
 
void CopyState ()
 
void SetManager (G4SteppingManager *const)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VSteppingVerbose
static void SetInstance (G4VSteppingVerbose *Instance)
 
static G4VSteppingVerboseGetInstance ()
 
static G4int GetSilent ()
 
static void SetSilent (G4int fSilent)
 
static G4int GetSilentStepInfo ()
 
static void SetSilentStepInfo (G4int fSilent)
 
- Protected Types inherited from G4VSteppingVerbose
typedef std::vector< G4intG4SelectedAtRestDoItVector
 
typedef std::vector< G4intG4SelectedAlongStepDoItVector
 
typedef std::vector< G4intG4SelectedPostStepDoItVector
 
- Protected Member Functions inherited from G4VSteppingVerbose
 G4VSteppingVerbose ()
 
- Protected Attributes inherited from G4VSteppingVerbose
G4SteppingManagerfManager
 
G4UserSteppingActionfUserSteppingAction
 
G4double PhysicalStep
 
G4double GeometricalStep
 
G4double CorrectedStep
 
G4bool PreStepPointIsGeom
 
G4bool FirstStep
 
G4StepStatus fStepStatus
 
G4double TempInitVelocity
 
G4double TempVelocity
 
G4double Mass
 
G4double sumEnergyChange
 
G4VParticleChangefParticleChange
 
G4TrackfTrack
 
G4TrackVectorfSecondary
 
G4StepfStep
 
G4StepPointfPreStepPoint
 
G4StepPointfPostStepPoint
 
G4VPhysicalVolumefCurrentVolume
 
G4VSensitiveDetectorfSensitive
 
G4VProcessfCurrentProcess
 
G4ProcessVectorfAtRestDoItVector
 
G4ProcessVectorfAlongStepDoItVector
 
G4ProcessVectorfPostStepDoItVector
 
G4ProcessVectorfAtRestGetPhysIntVector
 
G4ProcessVectorfAlongStepGetPhysIntVector
 
G4ProcessVectorfPostStepGetPhysIntVector
 
size_t MAXofAtRestLoops
 
size_t MAXofAlongStepLoops
 
size_t MAXofPostStepLoops
 
G4double currentMinimumStep
 
G4double numberOfInteractionLengthLeft
 
size_t fAtRestDoItProcTriggered
 
size_t fAlongStepDoItProcTriggered
 
size_t fPostStepDoItProcTriggered
 
G4int fN2ndariesAtRestDoIt
 
G4int fN2ndariesAlongStepDoIt
 
G4int fN2ndariesPostStepDoIt
 
G4NavigatorfNavigator
 
G4int verboseLevel
 
G4SelectedAtRestDoItVectorfSelectedAtRestDoItVector
 
G4SelectedAlongStepDoItVectorfSelectedAlongStepDoItVector
 
G4SelectedPostStepDoItVectorfSelectedPostStepDoItVector
 
G4double fPreviousStepSize
 
G4TouchableHandle fTouchableHandle
 
G4SteppingControl StepControlFlag
 
G4double physIntLength
 
G4ForceCondition fCondition
 
G4GPILSelection fGPILSelection
 
- Static Protected Attributes inherited from G4VSteppingVerbose
static G4ThreadLocal
G4VSteppingVerbose
fInstance = 0
 
static G4ThreadLocal G4int Silent = 0
 
static G4ThreadLocal G4int SilentStepInfo = 0
 

Detailed Description

Definition at line 51 of file G4SteppingVerbose.hh.

Constructor & Destructor Documentation

G4SteppingVerbose::G4SteppingVerbose ( )

Definition at line 57 of file G4SteppingVerbose.cc.

59 {
60 #ifdef G4_TRACKING_DEBUG
61  G4cout << "G4SteppingVerbose has instantiated" << G4endl;
62 #endif
63 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4SteppingVerbose::~G4SteppingVerbose ( )

Definition at line 66 of file G4SteppingVerbose.cc.

68 {
69 }

Member Function Documentation

void G4SteppingVerbose::AlongStepDoItAllDone ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 133 of file G4SteppingVerbose.cc.

135 {
136  if(Silent==1){ return; }
137 
138  G4VProcess* ptProcManager;
139 
140  CopyState();
141 
142  if(verboseLevel >= 3){
143  G4cout << G4endl;
144  G4cout << " >>AlongStepDoIt (after all invocations):" << G4endl;
145  G4cout << " ++List of invoked processes " << G4endl;
146 
147  for(size_t ci=0; ci<MAXofAlongStepLoops; ci++){
148  ptProcManager = (*fAlongStepDoItVector)(ci);
149  G4cout << " " << ci+1 << ") ";
150  if(ptProcManager != 0){
151  G4cout << ptProcManager->GetProcessName() << G4endl;
152  }
153  }
154 
155  ShowStep();
156  G4cout << G4endl;
157  G4cout << " ++List of secondaries generated "
158  << "(x,y,z,kE,t,PID):"
159  << " No. of secodaries = "
160  << (*fSecondary).size() << G4endl;
161 
162  if((*fSecondary).size()>0){
163  for(size_t lp1=0; lp1<(*fSecondary).size(); lp1++){
164  G4cout << " "
165  << std::setw( 9)
166  << G4BestUnit((*fSecondary)[lp1]->GetPosition().x(),"Length") << " "
167  << std::setw( 9)
168  << G4BestUnit((*fSecondary)[lp1]->GetPosition().y(),"Length") << " "
169  << std::setw( 9)
170  << G4BestUnit((*fSecondary)[lp1]->GetPosition().z(),"Length") << " "
171  << std::setw( 9)
172  << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy(),"Energy") << " "
173  << std::setw( 9)
174  << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime(),"Time") << " "
175  << std::setw(18)
176  << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
177  }
178  }
179  }
180 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4ThreadLocal G4int Silent
void ShowStep() const
G4GLOB_DLL std::ostream G4cout
G4TrackVector * fSecondary
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::AlongStepDoItOneByOne ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 467 of file G4SteppingVerbose.cc.

469 {
470  if(Silent==1){ return; }
471 
472  CopyState();
473 
474  if(verboseLevel >= 4){
475  G4cout << G4endl;
476  G4cout << " >>AlongStepDoIt (process by process): "
477  << " Process Name = "
479 
480  ShowStep();
481  G4cout << " "
482  << "!Note! Safety of PostStep is only valid "
483  << "after all DoIt invocations."
484  << G4endl;
485 
487  G4cout << G4endl;
488 
489  G4cout << " ++List of secondaries generated "
490  << "(x,y,z,kE,t,PID):"
491  << " No. of secodaries = "
493 
495  for(size_t lp1=(*fSecondary).size()-fN2ndariesAlongStepDoIt; lp1<(*fSecondary).size(); lp1++){
496  G4cout << " "
497  << std::setw( 9)
498  << G4BestUnit((*fSecondary)[lp1]->GetPosition().x() , "Length")<< " "
499  << std::setw( 9)
500  << G4BestUnit((*fSecondary)[lp1]->GetPosition().y() , "Length")<< " "
501  << std::setw( 9)
502  << G4BestUnit((*fSecondary)[lp1]->GetPosition().z() , "Length")<< " "
503  << std::setw( 9)
504  << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy() , "Energy")<< " "
505  << std::setw( 9)
506  << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime() , "Time")<< " "
507  << std::setw(18)
508  << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
509  }
510  }
511  }
512 }
G4VProcess * fCurrentProcess
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4ThreadLocal G4int Silent
void ShowStep() const
G4GLOB_DLL std::ostream G4cout
G4TrackVector * fSecondary
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::AtRestDoItInvoked ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 78 of file G4SteppingVerbose.cc.

80  {
81  if(Silent==1){ return; }
82 
83  G4VProcess* ptProcManager;
84  CopyState();
85 
86  if(verboseLevel >= 3 ){
87  G4int npt=0;
88  G4cout << " **List of AtRestDoIt invoked:" << G4endl;
89  for(size_t np=0; np < MAXofAtRestLoops; np++){
90  size_t npGPIL = MAXofAtRestLoops-np-1;
91  if( (*fSelectedAtRestDoItVector)[npGPIL] == 2 ){
92  npt++;
93  ptProcManager = (*fAtRestDoItVector)[np];
94  G4cout << " # " << npt << " : "
95  << ptProcManager->GetProcessName()
96  << " (Forced)" << G4endl;
97  } else if ( (*fSelectedAtRestDoItVector)[npGPIL] == 1 ){
98  npt++;
99  ptProcManager = (*fAtRestDoItVector)[np];
100  G4cout << " # " << npt << " : " << ptProcManager->GetProcessName() << G4endl;
101  }
102  }
103 
104  G4cout << " Generated secondries # : " << fN2ndariesAtRestDoIt << G4endl;
105 
106  if( fN2ndariesAtRestDoIt > 0 ){
107  G4cout << " -- List of secondaries generated : " << "(x,y,z,kE,t,PID) --" << G4endl;
108  for( size_t lp1=(*fSecondary).size()-fN2ndariesAtRestDoIt;
109  lp1<(*fSecondary).size(); lp1++) {
110  G4cout << " "
111  << std::setw( 9)
112  << G4BestUnit((*fSecondary)[lp1]->GetPosition().x(),"Length") << " "
113  << std::setw( 9)
114  << G4BestUnit((*fSecondary)[lp1]->GetPosition().y(),"Length") << " "
115  << std::setw( 9)
116  << G4BestUnit((*fSecondary)[lp1]->GetPosition().z(),"Length") << " "
117  << std::setw( 9)
118  << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy(),"Energy") << " "
119  << std::setw( 9)
120  << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime(),"Time") << " "
121  << std::setw(18)
122  << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
123  }
124  }
125  }
126 
127  if( verboseLevel >= 4 ){
128  ShowStep();
129  G4cout << G4endl;
130  }
131 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
static G4ThreadLocal G4int Silent
void ShowStep() const
G4GLOB_DLL std::ostream G4cout
G4TrackVector * fSecondary
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::DPSLAlongStep ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 388 of file G4SteppingVerbose.cc.

390 {
391  if(Silent==1){ return; }
392  CopyState();
393 
394  if( verboseLevel > 5 ){
395  G4cout << " ++ProposedStep(AlongStep) = "
396  << std::setw( 9) << G4BestUnit(physIntLength , "Length")
397  << " : ProcName = "
399  << " (";
401  G4cout << "CandidateForSelection)" << G4endl;
402  }
404  G4cout << "NotCandidateForSelection)" << G4endl;
405  }
406  else{
407  G4cout << "?!?)" << G4endl;
408  }
409  }
410 }
G4VProcess * fCurrentProcess
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4ThreadLocal G4int Silent
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
G4GPILSelection fGPILSelection

Here is the call graph for this function:

void G4SteppingVerbose::DPSLPostStep ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 361 of file G4SteppingVerbose.cc.

363 {
364  if(Silent==1){ return; }
365  CopyState();
366 
367  if( verboseLevel > 5 ){
368  G4cout << " ++ProposedStep(PostStep ) = " << std::setw( 9) << physIntLength
369  << " : ProcName = " << fCurrentProcess->GetProcessName() << " (";
371  G4cout << "ExclusivelyForced)" << G4endl;
372  }
373  else if(fCondition==StronglyForced){
374  G4cout << "StronglyForced)" << G4endl;
375  }
376  else if(fCondition==Conditionally){
377  G4cout << "Conditionally)" << G4endl;
378  }
379  else if(fCondition==Forced){
380  G4cout << "Forced)" << G4endl;
381  }
382  else{
383  G4cout << "No ForceCondition)" << G4endl;
384  }
385  }
386 }
G4VProcess * fCurrentProcess
static G4ThreadLocal G4int Silent
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4ForceCondition fCondition
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::DPSLStarted ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 336 of file G4SteppingVerbose.cc.

338 {
339  if(Silent==1){ return; }
340  CopyState();
341 
342  if( verboseLevel > 5 ){
343  G4cout << G4endl << " >>DefinePhysicalStepLength (List of proposed StepLengths): " << G4endl;
344  }
345 }
static G4ThreadLocal G4int Silent
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::DPSLUserLimit ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 347 of file G4SteppingVerbose.cc.

349 {
350  if(Silent==1){ return; }
351  CopyState();
352 
353  if( verboseLevel > 5 ){
354  G4cout << G4endl << G4endl;
355  G4cout << "=== Defined Physical Step Length (DPSL)" << G4endl;
356  G4cout << " ++ProposedStep(UserLimit) = " << std::setw( 9) << physIntLength
357  << " : ProcName = User defined maximum allowed Step" << G4endl;
358  }
359 }
static G4ThreadLocal G4int Silent
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::NewStep ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 72 of file G4SteppingVerbose.cc.

74 {
75 }
void G4SteppingVerbose::PostStepDoItAllDone ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 182 of file G4SteppingVerbose.cc.

184 {
185  if(Silent==1){ return; }
186 
187  G4VProcess* ptProcManager;
188 
189  CopyState();
190 
191  if( (fStepStatus == fPostStepDoItProc) |
192  (fCondition == Forced) |
195  (fCondition == StronglyForced) ){
196 
197  if(verboseLevel >= 3){
198  G4int npt=0;
199  G4cout << G4endl;
200  G4cout << " **PostStepDoIt (after all invocations):" << G4endl;
201  G4cout << " ++List of invoked processes " << G4endl;
202 
203  for(size_t np=0; np < MAXofPostStepLoops; np++){
204  size_t npGPIL = MAXofPostStepLoops-np-1;
205  if( (*fSelectedPostStepDoItVector)[npGPIL] == 2){
206  npt++;
207  ptProcManager = (*fPostStepDoItVector)[np];
208  G4cout << " " << npt << ") "
209  << ptProcManager->GetProcessName()
210  << " (Forced)" << G4endl;
211  } else if ( (*fSelectedPostStepDoItVector)[npGPIL] == 1){
212  npt++;
213  ptProcManager = (*fPostStepDoItVector)[np];
214  G4cout << " " << npt << ") " << ptProcManager->GetProcessName() << G4endl;
215  }
216  }
217 
218  ShowStep();
219  G4cout << G4endl;
220  G4cout << " ++List of secondaries generated "
221  << "(x,y,z,kE,t,PID):"
222  << " No. of secodaries = "
223  << (*fSecondary).size() << G4endl;
224  G4cout << " [Note]Secondaries from AlongStepDoIt included." << G4endl;
225 
226  if((*fSecondary).size()>0){
227  for(size_t lp1=0; lp1<(*fSecondary).size(); lp1++){
228  G4cout << " "
229  << std::setw( 9)
230  << G4BestUnit((*fSecondary)[lp1]->GetPosition().x() , "Length") << " "
231  << std::setw( 9)
232  << G4BestUnit((*fSecondary)[lp1]->GetPosition().y() , "Length") << " "
233  << std::setw( 9)
234  << G4BestUnit((*fSecondary)[lp1]->GetPosition().z() , "Length") << " "
235  << std::setw( 9)
236  << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy() , "Energy") << " "
237  << std::setw( 9)
238  << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime() , "Time") << " "
239  << std::setw(18)
240  << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
241  }
242  }
243  }
244  }
245 }
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
static G4ThreadLocal G4int Silent
void ShowStep() const
G4GLOB_DLL std::ostream G4cout
G4TrackVector * fSecondary
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4ForceCondition fCondition
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::PostStepDoItOneByOne ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 514 of file G4SteppingVerbose.cc.

516 {
517  if(Silent==1){ return; }
518 
519  CopyState();
520 
521  if(verboseLevel >= 4){
522  G4cout << G4endl;
523  G4cout << " >>PostStepDoIt (process by process): "
524  << " Process Name = "
526 
527  ShowStep();
528  G4cout << G4endl;
530  G4cout << G4endl;
531 
532  G4cout << " ++List of secondaries generated "
533  << "(x,y,z,kE,t,PID):"
534  << " No. of secodaries = "
536 
538  for(size_t lp1=(*fSecondary).size()-fN2ndariesPostStepDoIt; lp1<(*fSecondary).size(); lp1++){
539  G4cout << " "
540  << std::setw( 9)
541  << G4BestUnit((*fSecondary)[lp1]->GetPosition().x() , "Length")<< " "
542  << std::setw( 9)
543  << G4BestUnit((*fSecondary)[lp1]->GetPosition().y(), "Length") << " "
544  << std::setw( 9)
545  << G4BestUnit((*fSecondary)[lp1]->GetPosition().z(), "Length") << " "
546  << std::setw( 9)
547  << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy(), "Energy") << " "
548  << std::setw( 9)
549  << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime(), "Time") << " "
550  << std::setw(18)
551  << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
552  }
553  }
554  }
555 }
G4VProcess * fCurrentProcess
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4ThreadLocal G4int Silent
void ShowStep() const
G4GLOB_DLL std::ostream G4cout
G4TrackVector * fSecondary
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::ShowStep ( ) const

Definition at line 735 of file G4SteppingVerbose.cc.

737 {
738  if(Silent==1){ return; }
739  G4String volName;
740  G4int oldprc;
741 
742 // Show header
743  G4cout << G4endl;
744  G4cout << " ++G4Step Information " << G4endl;
745  oldprc = G4cout.precision(16);
746 
747 // Show G4Step specific information
748  G4cout << " Address of G4Track : " << fStep->GetTrack() << G4endl;
749  G4cout << " Step Length (mm) : " << fStep->GetTrack()->GetStepLength() << G4endl;
750  G4cout << " Energy Deposit (MeV) : " << fStep->GetTotalEnergyDeposit() << G4endl;
751 
752 // Show G4StepPoint specific information
753  G4cout << " -------------------------------------------------------"
754  << "----------------" << G4endl;
755  G4cout << " StepPoint Information " << std::setw(20) << "PreStep"
756  << std::setw(20) << "PostStep" << G4endl;
757  G4cout << " -------------------------------------------------------"
758  << "----------------" << G4endl;
759  G4cout << " Position - x (mm) : "
760  << std::setw(20) << fStep->GetPreStepPoint()->GetPosition().x()
761  << std::setw(20) << fStep->GetPostStepPoint()->GetPosition().x() << G4endl;
762  G4cout << " Position - y (mm) : "
763  << std::setw(20) << fStep->GetPreStepPoint()->GetPosition().y()
764  << std::setw(20) << fStep->GetPostStepPoint()->GetPosition().y() << G4endl;
765  G4cout << " Position - z (mm) : "
766  << std::setw(20) << fStep->GetPreStepPoint()->GetPosition().z()
767  << std::setw(20) << fStep->GetPostStepPoint()->GetPosition().z() << G4endl;
768  G4cout << " Global Time (ns) : "
769  << std::setw(20) << fStep->GetPreStepPoint()->GetGlobalTime()
770  << std::setw(20) << fStep->GetPostStepPoint()->GetGlobalTime() << G4endl;
771  G4cout << " Local Time (ns) : "
772  << std::setw(20) << fStep->GetPreStepPoint()->GetLocalTime()
773  << std::setw(20) << fStep->GetPostStepPoint()->GetLocalTime() << G4endl;
774  G4cout << " Proper Time (ns) : "
775  << std::setw(20) << fStep->GetPreStepPoint()->GetProperTime()
776  << std::setw(20) << fStep->GetPostStepPoint()->GetProperTime() << G4endl;
777  G4cout << " Momentum Direct - x : "
778  << std::setw(20) << fStep->GetPreStepPoint()->GetMomentumDirection().x()
779  << std::setw(20) << fStep->GetPostStepPoint()->GetMomentumDirection().x() << G4endl;
780  G4cout << " Momentum Direct - y : "
781  << std::setw(20) << fStep->GetPreStepPoint()->GetMomentumDirection().y()
782  << std::setw(20) << fStep->GetPostStepPoint()->GetMomentumDirection().y() << G4endl;
783  G4cout << " Momentum Direct - z : "
784  << std::setw(20) << fStep->GetPreStepPoint()->GetMomentumDirection().z()
785  << std::setw(20) << fStep->GetPostStepPoint()->GetMomentumDirection().z() << G4endl;
786  G4cout << " Momentum - x (MeV/c): "
787  << std::setw(20) << fStep->GetPreStepPoint()->GetMomentum().x()
788  << std::setw(20) << fStep->GetPostStepPoint()->GetMomentum().x() << G4endl;
789  G4cout << " Momentum - y (MeV/c): "
790  << std::setw(20) << fStep->GetPreStepPoint()->GetMomentum().y()
791  << std::setw(20) << fStep->GetPostStepPoint()->GetMomentum().y() << G4endl;
792  G4cout << " Momentum - z (MeV/c): "
793  << std::setw(20) << fStep->GetPreStepPoint()->GetMomentum().z()
794  << std::setw(20) << fStep->GetPostStepPoint()->GetMomentum().z() << G4endl;
795  G4cout << " Total Energy (MeV) : "
796  << std::setw(20) << fStep->GetPreStepPoint()->GetTotalEnergy()
797  << std::setw(20) << fStep->GetPostStepPoint()->GetTotalEnergy() << G4endl;
798  G4cout << " Kinetic Energy (MeV): "
799  << std::setw(20) << fStep->GetPreStepPoint()->GetKineticEnergy()
800  << std::setw(20) << fStep->GetPostStepPoint()->GetKineticEnergy() << G4endl;
801  G4cout << " Velocity (mm/ns) : "
802  << std::setw(20) << fStep->GetPreStepPoint()->GetVelocity()
803  << std::setw(20) << fStep->GetPostStepPoint()->GetVelocity() << G4endl;
804  G4cout << " Volume Name : "
805  << std::setw(20) << fStep->GetPreStepPoint()->GetPhysicalVolume()->GetName();
807  {
808  volName = fStep->GetPostStepPoint()->GetPhysicalVolume()->GetName();
809  }
810  else
811  {
812  volName = "OutOfWorld";
813  }
814  G4cout << std::setw(20) << volName << G4endl;
815  G4cout << " Safety (mm) : "
816  << std::setw(20) << fStep->GetPreStepPoint()->GetSafety()
817  << std::setw(20) << fStep->GetPostStepPoint()->GetSafety() << G4endl;
818  G4cout << " Polarization - x : "
819  << std::setw(20) << fStep->GetPreStepPoint()->GetPolarization().x()
820  << std::setw(20) << fStep->GetPostStepPoint()->GetPolarization().x() << G4endl;
821  G4cout << " Polarization - y : "
822  << std::setw(20) << fStep->GetPreStepPoint()->GetPolarization().y()
823  << std::setw(20) << fStep->GetPostStepPoint()->GetPolarization().y() << G4endl;
824  G4cout << " Polarization - Z : "
825  << std::setw(20) << fStep->GetPreStepPoint()->GetPolarization().z()
826  << std::setw(20) << fStep->GetPostStepPoint()->GetPolarization().z() << G4endl;
827  G4cout << " Weight : "
828  << std::setw(20) << fStep->GetPreStepPoint()->GetWeight()
829  << std::setw(20) << fStep->GetPostStepPoint()->GetWeight() << G4endl;
830  G4cout << " Step Status : " ;
831  G4StepStatus tStepStatus = fStep->GetPreStepPoint()->GetStepStatus();
832  if( tStepStatus == fGeomBoundary ){
833  G4cout << std::setw(20) << "Geom Limit";
834  } else if ( tStepStatus == fAlongStepDoItProc ){
835  G4cout << std::setw(20) << "AlongStep Proc.";
836  } else if ( tStepStatus == fPostStepDoItProc ){
837  G4cout << std::setw(20) << "PostStep Proc";
838  } else if ( tStepStatus == fAtRestDoItProc ){
839  G4cout << std::setw(20) << "AtRest Proc";
840  } else if ( tStepStatus == fUndefined ){
841  G4cout << std::setw(20) << "Undefined";
842  }
843 
844  tStepStatus = fStep->GetPostStepPoint()->GetStepStatus();
845  if( tStepStatus == fGeomBoundary ){
846  G4cout << std::setw(20) << "Geom Limit";
847  } else if ( tStepStatus == fAlongStepDoItProc ){
848  G4cout << std::setw(20) << "AlongStep Proc.";
849  } else if ( tStepStatus == fPostStepDoItProc ){
850  G4cout << std::setw(20) << "PostStep Proc";
851  } else if ( tStepStatus == fAtRestDoItProc ){
852  G4cout << std::setw(20) << "AtRest Proc";
853  } else if ( tStepStatus == fUndefined ){
854  G4cout << std::setw(20) << "Undefined";
855  }
856 
857  G4cout << G4endl;
858  G4cout << " Process defined Step: " ;
860  G4cout << std::setw(20) << "Undefined";
861  } else {
862  G4cout << std::setw(20) << fStep->GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName();
863  }
865  G4cout << std::setw(20) << "Undefined";
866  } else {
867  G4cout << std::setw(20) << fStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
868  }
869  G4cout.precision(oldprc);
870 
871  G4cout << G4endl;
872  G4cout << " -------------------------------------------------------"
873  << "----------------" << G4endl;
874 }
G4double GetTotalEnergy() const
G4double GetWeight() const
double x() const
G4StepStatus GetStepStatus() const
G4ThreeVector GetMomentum() const
G4double GetVelocity() const
int G4int
Definition: G4Types.hh:78
double z() const
G4StepStatus
Definition: G4StepStatus.hh:49
G4double GetLocalTime() const
static G4ThreadLocal G4int Silent
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
const G4ThreeVector & GetPosition() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTotalEnergyDeposit() const
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
double y() const
G4double GetSafety() const
G4double GetProperTime() const
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
const G4ThreeVector & GetPolarization() const
G4double GetStepLength() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4SteppingVerbose::StepInfo ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 248 of file G4SteppingVerbose.cc.

250 {
251  if(Silent==1){ return; }
252  if(SilentStepInfo==1){ return; }
253 
254  CopyState();
255  G4cout.precision(16);
256  G4int prec = G4cout.precision(3);
257 
258  if( verboseLevel >= 1 ){
259  if( verboseLevel >= 4 ) VerboseTrack();
260  if( verboseLevel >= 3 ){
261  G4cout << G4endl;
262 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
263  G4cout << std::setw( 5) << "#Step#" << " "
264  << std::setw( 8) << "X" << " " << std::setw( 8) << "Y" << " "
265  << std::setw( 8) << "Z" << " "
266  << std::setw( 9) << "KineE" << " " << std::setw( 8) << "dE" << " "
267  << std::setw(12) << "StepLeng" << " " << std::setw(12) << "TrackLeng" << " "
268  << std::setw(12) << "NextVolume" << " " << std::setw( 8) << "ProcName" << G4endl;
269 #else
270  G4cout << std::setw( 5) << "#Step#" << " "
271  << std::setw( 8) << "X(mm)" << " " << std::setw( 8) << "Y(mm)" << " "
272  << std::setw( 8) << "Z(mm)" << " "
273  << std::setw( 9) << "KinE(MeV)" << " " << std::setw( 8) << "dE(MeV)" << " "
274  << std::setw( 8) << "StepLeng" << " " << std::setw( 9) << "TrackLeng" << " "
275  << std::setw(11) << "NextVolume" << " " << std::setw( 8) << "ProcName" << G4endl;
276 #endif
277  }
278  G4cout << std::setw( 5) << fTrack->GetCurrentStepNumber() << " "
279  << std::setw( 8) << G4BestUnit(fTrack->GetPosition().x() , "Length") << " "
280  << std::setw( 8) << G4BestUnit(fTrack->GetPosition().y() , "Length") << " "
281  << std::setw( 8) << G4BestUnit(fTrack->GetPosition().z() , "Length") << " "
282  << std::setw( 9) << G4BestUnit(fTrack->GetKineticEnergy() , "Energy") << " "
283  << std::setw( 8) << G4BestUnit(fStep->GetTotalEnergyDeposit(), "Energy") << " "
284  << std::setw( 8) << G4BestUnit(fStep->GetStepLength() , "Length") << " "
285  << std::setw( 9) << G4BestUnit(fTrack->GetTrackLength() , "Length") << " ";
286 
287  // Put cut comment here
288  if( fTrack->GetNextVolume() != 0 ) {
289  G4cout << std::setw(11) << fTrack->GetNextVolume()->GetName() << " ";
290  } else {
291  G4cout << std::setw(11) << "OutOfWorld" << " ";
292  }
295  } else {
296  G4cout << "User Limit";
297  }
298  G4cout << G4endl;
299  if( verboseLevel == 2 )
300  {
302  if(tN2ndariesTot>0){
303  G4cout << " :----- List of 2ndaries - "
304  << "#SpawnInStep=" << std::setw(3) << tN2ndariesTot
305  << "(Rest=" << std::setw(2) << fN2ndariesAtRestDoIt
306  << ",Along=" << std::setw(2) << fN2ndariesAlongStepDoIt
307  << ",Post=" << std::setw(2) << fN2ndariesPostStepDoIt
308  << "), "
309  << "#SpawnTotal=" << std::setw(3) << (*fSecondary).size()
310  << " ---------------"
311  << G4endl;
312 
313  for(size_t lp1=(*fSecondary).size()-tN2ndariesTot; lp1<(*fSecondary).size(); lp1++){
314  G4cout << " : "
315  << std::setw( 9)
316  << G4BestUnit((*fSecondary)[lp1]->GetPosition().x() , "Length")<< " "
317  << std::setw( 9)
318  << G4BestUnit((*fSecondary)[lp1]->GetPosition().y() , "Length")<< " "
319  << std::setw( 9)
320  << G4BestUnit((*fSecondary)[lp1]->GetPosition().z() , "Length") << " "
321  << std::setw( 9)
322  << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy() , "Energy")<< " "
323  << std::setw(18)
324  << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
325  }
326  G4cout << " :-----------------------------" << "----------------------------------"
327  << "-- EndOf2ndaries Info ---------------" << G4endl;
328  }
329  }
330  }
331  G4cout.precision(prec);
332 }
G4double GetStepLength() const
double x() const
const G4ThreeVector & GetPosition() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
double z() const
G4VPhysicalVolume * GetNextVolume() const
static G4ThreadLocal G4int Silent
static G4ThreadLocal G4int SilentStepInfo
static const double prec
Definition: RanecuEngine.cc:58
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
G4TrackVector * fSecondary
const G4String & GetName() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTotalEnergyDeposit() const
G4double GetTrackLength() const
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
double y() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::TrackingStarted ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 414 of file G4SteppingVerbose.cc.

416 {
417  if(Silent==1){ return; }
418 
419  CopyState();
420 
421  G4int prec = G4cout.precision(3);
422  if( verboseLevel > 0 ){
423 
424 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
425  G4cout << std::setw( 5) << "Step#" << " "
426  << std::setw( 8) << "X" << " "
427  << std::setw( 8) << "Y" << " "
428  << std::setw( 8) << "Z" << " "
429  << std::setw( 9) << "KineE" << " "
430  << std::setw( 8) << "dE" << " "
431  << std::setw(12) << "StepLeng" << " "
432  << std::setw(12) << "TrackLeng" << " "
433  << std::setw(12) << "NextVolume" << " "
434  << std::setw( 8) << "ProcName" << G4endl;
435 #else
436  G4cout << std::setw( 5) << "Step#" << " "
437  << std::setw( 8) << "X(mm)" << " "
438  << std::setw( 8) << "Y(mm)" << " "
439  << std::setw( 8) << "Z(mm)" << " "
440  << std::setw( 9) << "KinE(MeV)" << " "
441  << std::setw( 8) << "dE(MeV)" << " "
442  << std::setw( 8) << "StepLeng" << " "
443  << std::setw( 9) << "TrackLeng" << " "
444  << std::setw(11) << "NextVolume" << " "
445  << std::setw( 8) << "ProcName" << G4endl;
446 #endif
447 
448  G4cout << std::setw( 5) << fTrack->GetCurrentStepNumber() << " "
449  << std::setw( 8) << G4BestUnit(fTrack->GetPosition().x(),"Length")<< " "
450  << std::setw( 8) << G4BestUnit(fTrack->GetPosition().y(),"Length") << " "
451  << std::setw( 8) << G4BestUnit(fTrack->GetPosition().z(),"Length")<< " "
452  << std::setw( 9) << G4BestUnit(fTrack->GetKineticEnergy(),"Energy")<< " "
453  << std::setw( 8) << G4BestUnit(fStep->GetTotalEnergyDeposit(),"Energy") << " "
454  << std::setw( 8) << G4BestUnit(fStep->GetStepLength(),"Length")<< " "
455  << std::setw( 9) << G4BestUnit(fTrack->GetTrackLength(),"Length") << " ";
456 
457  if(fTrack->GetNextVolume()){
458  G4cout << std::setw(11) << fTrack->GetNextVolume()->GetName() << " ";
459  } else {
460  G4cout << std::setw(11) << "OutOfWorld" << " ";
461  }
462  G4cout << "initStep" << G4endl;
463  }
464  G4cout.precision(prec);
465 }
G4double GetStepLength() const
double x() const
const G4ThreeVector & GetPosition() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
double z() const
G4VPhysicalVolume * GetNextVolume() const
static G4ThreadLocal G4int Silent
static const double prec
Definition: RanecuEngine.cc:58
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
const G4String & GetName() const
G4double GetTotalEnergyDeposit() const
G4double GetTrackLength() const
double y() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4SteppingVerbose::VerboseParticleChange ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 725 of file G4SteppingVerbose.cc.

727 {
728  if(Silent==1){ return; }
729 // Show header
730  G4cout << G4endl;
731  G4cout << " ++G4ParticleChange Information " << G4endl;
733 }
G4VParticleChange * fParticleChange
virtual void DumpInfo() const
static G4ThreadLocal G4int Silent
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4SteppingVerbose::VerboseTrack ( )
virtual

Implements G4VSteppingVerbose.

Definition at line 559 of file G4SteppingVerbose.cc.

561 {
562  if(Silent==1){ return; }
563 
564  CopyState();
565 // Show header
566  G4cout << G4endl;
567  G4cout << " ++G4Track Information " << G4endl;
568  G4int prec = G4cout.precision(3);
569 
570 
571  G4cout << " -----------------------------------------------"
572  << G4endl;
573  G4cout << " G4Track Information " << std::setw(20) << G4endl;
574  G4cout << " -----------------------------------------------"
575  << G4endl;
576 
577  G4cout << " Step number : "
578  << std::setw(20) << fTrack->GetCurrentStepNumber()
579  << G4endl;
580 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
581  G4cout << " Position - x : "
582  << std::setw(20) << G4BestUnit(fTrack->GetPosition().x(), "Length")
583  << G4endl;
584  G4cout << " Position - y : "
585  << std::setw(20) << G4BestUnit(fTrack->GetPosition().y(), "Length")
586  << G4endl;
587  G4cout << " Position - z : "
588  << std::setw(20) << G4BestUnit(fTrack->GetPosition().z(), "Length")
589  << G4endl;
590  G4cout << " Global Time : "
591  << std::setw(20) << G4BestUnit(fTrack->GetGlobalTime(), "Time")
592  << G4endl;
593  G4cout << " Local Time : "
594  << std::setw(20) << G4BestUnit(fTrack->GetLocalTime(), "Time")
595  << G4endl;
596 #else
597  G4cout << " Position - x (mm) : "
598  << std::setw(20) << fTrack->GetPosition().x() /mm
599  << G4endl;
600  G4cout << " Position - y (mm) : "
601  << std::setw(20) << fTrack->GetPosition().y() /mm
602  << G4endl;
603  G4cout << " Position - z (mm) : "
604  << std::setw(20) << fTrack->GetPosition().z() /mm
605  << G4endl;
606  G4cout << " Global Time (ns) : "
607  << std::setw(20) << fTrack->GetGlobalTime() /ns
608  << G4endl;
609  G4cout << " Local Time (ns) : "
610  << std::setw(20) << fTrack->GetLocalTime() /ns
611  << G4endl;
612 #endif
613  G4cout << " Momentum Direct - x : "
614  << std::setw(20) << fTrack->GetMomentumDirection().x()
615  << G4endl;
616  G4cout << " Momentum Direct - y : "
617  << std::setw(20) << fTrack->GetMomentumDirection().y()
618  << G4endl;
619  G4cout << " Momentum Direct - z : "
620  << std::setw(20) << fTrack->GetMomentumDirection().z()
621  << G4endl;
622 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
623  G4cout << " Kinetic Energy : "
624 #else
625  G4cout << " Kinetic Energy (MeV): "
626 #endif
627  << std::setw(20) << G4BestUnit(fTrack->GetKineticEnergy(), "Energy")
628  << G4endl;
629  G4cout << " Polarization - x : "
630  << std::setw(20) << fTrack->GetPolarization().x()
631  << G4endl;
632  G4cout << " Polarization - y : "
633  << std::setw(20) << fTrack->GetPolarization().y()
634  << G4endl;
635  G4cout << " Polarization - z : "
636  << std::setw(20) << fTrack->GetPolarization().z()
637  << G4endl;
638  G4cout << " Track Length : "
639  << std::setw(20) << G4BestUnit(fTrack->GetTrackLength(), "Length")
640  << G4endl;
641  G4cout << " Track ID # : "
642  << std::setw(20) << fTrack->GetTrackID()
643  << G4endl;
644  G4cout << " Parent Track ID # : "
645  << std::setw(20) << fTrack->GetParentID()
646  << G4endl;
647  G4cout << " Next Volume : "
648  << std::setw(20);
649  if( fTrack->GetNextVolume() != 0 ) {
650  G4cout << fTrack->GetNextVolume()->GetName() << " ";
651  } else {
652  G4cout << "OutOfWorld" << " ";
653  }
654  G4cout << G4endl;
655  G4cout << " Track Status : "
656  << std::setw(20);
657  if( fTrack->GetTrackStatus() == fAlive ){
658  G4cout << " Alive";
659  } else if( fTrack->GetTrackStatus() == fStopButAlive ){
660  G4cout << " StopButAlive";
661  } else if( fTrack->GetTrackStatus() == fStopAndKill ){
662  G4cout << " StopAndKill";
664  G4cout << " KillTrackAndSecondaries";
665  } else if( fTrack->GetTrackStatus() == fSuspend ){
666  G4cout << " Suspend";
667  } else if( fTrack->GetTrackStatus() == fPostponeToNextEvent ){
668  G4cout << " PostponeToNextEvent";
669  }
670  G4cout << G4endl;
671 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
672  G4cout << " Vertex - x : "
673  << std::setw(20) << G4BestUnit(fTrack->GetVertexPosition().x(),"Length")
674  << G4endl;
675  G4cout << " Vertex - y : "
676  << std::setw(20) << G4BestUnit(fTrack->GetVertexPosition().y(),"Length")
677  << G4endl;
678  G4cout << " Vertex - z : "
679  << std::setw(20) << G4BestUnit(fTrack->GetVertexPosition().z(),"Length")
680  << G4endl;
681 #else
682  G4cout << " Vertex - x (mm) : "
683  << std::setw(20) << fTrack->GetVertexPosition().x()/mm
684  << G4endl;
685  G4cout << " Vertex - y (mm) : "
686  << std::setw(20) << fTrack->GetVertexPosition().y()/mm
687  << G4endl;
688  G4cout << " Vertex - z (mm) : "
689  << std::setw(20) << fTrack->GetVertexPosition().z()/mm
690  << G4endl;
691 #endif
692  G4cout << " Vertex - Px (MomDir): "
693  << std::setw(20) << fTrack->GetVertexMomentumDirection().x()
694  << G4endl;
695  G4cout << " Vertex - Py (MomDir): "
696  << std::setw(20) << fTrack->GetVertexMomentumDirection().y()
697  << G4endl;
698  G4cout << " Vertex - Pz (MomDir): "
699  << std::setw(20) << fTrack->GetVertexMomentumDirection().z()
700  << G4endl;
701 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
702  G4cout << " Vertex - KineE : "
703 #else
704  G4cout << " Vertex - KineE (MeV): "
705 #endif
706  << std::setw(20) << G4BestUnit(fTrack->GetVertexKineticEnergy(),"Energy")
707  << G4endl;
708 
709  G4cout << " Creator Process : "
710  << std::setw(20);
711  if( fTrack->GetCreatorProcess() == 0){
712  G4cout << " Event Generator" << G4endl;
713  } else {
715  }
716 
717  G4cout << " -----------------------------------------------"
718  << G4endl;
719 
720  G4cout.precision(prec);
721 }
G4int GetParentID() const
const G4ThreeVector & GetPolarization() const
G4double GetLocalTime() const
static constexpr double mm
Definition: G4SIunits.hh:115
double x() const
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
double z() const
G4VPhysicalVolume * GetNextVolume() const
static G4ThreadLocal G4int Silent
const G4VProcess * GetCreatorProcess() const
static const double prec
Definition: RanecuEngine.cc:58
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
const G4String & GetName() const
G4double GetVertexKineticEnergy() const
G4int GetTrackID() const
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4ThreeVector & GetVertexPosition() const
G4double GetTrackLength() const
const G4ThreeVector & GetMomentumDirection() const
double y() const
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetVertexMomentumDirection() const
#define ns
Definition: xmlparse.cc:614

Here is the call graph for this function:

Here is the caller graph for this function:


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