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

#include <G4TransportationLogger.hh>

Public Member Functions

 G4TransportationLogger (const G4String &className, G4int verbosity)
 
 G4TransportationLogger (const char *className, G4int verbosity)
 
 ~G4TransportationLogger ()
 
void ReportLoopingTrack (const G4Track &track, const G4Step &stepInfo, G4int numTrials, G4long noCalls, const char *methodName) const
 
void ReportLooperThresholds (const char *className)
 
void SetThresholds (G4double newEnWarn, G4double importantEnergy, G4int newMaxTrials)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int level)
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4double GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double val)
 
void SetThresholdImportantEnergy (G4double val)
 
void SetThresholdTrials (G4int maxNoTrials)
 

Detailed Description

Definition at line 47 of file G4TransportationLogger.hh.

Constructor & Destructor Documentation

◆ G4TransportationLogger() [1/2]

G4TransportationLogger::G4TransportationLogger ( const G4String & className,
G4int verbosity )

Definition at line 42 of file G4TransportationLogger.cc.

44 : fClassName(className), fVerbose(verbosity),
45 fThldWarningEnergy(0.), fThldImportantEnergy(0.), fThldTrials(0)
46{
47}

◆ G4TransportationLogger() [2/2]

G4TransportationLogger::G4TransportationLogger ( const char * className,
G4int verbosity )

Definition at line 49 of file G4TransportationLogger.cc.

51 : fClassName(className), fVerbose(verbosity),
52 fThldWarningEnergy(0.), fThldImportantEnergy(0.), fThldTrials(0)
53{
54}

◆ ~G4TransportationLogger()

G4TransportationLogger::~G4TransportationLogger ( )

Definition at line 56 of file G4TransportationLogger.cc.

57{
58}

Member Function Documentation

◆ GetThresholdImportantEnergy()

G4double G4TransportationLogger::GetThresholdImportantEnergy ( ) const
inline

Definition at line 77 of file G4TransportationLogger.hh.

77{ return fThldImportantEnergy; }

Referenced by ReportLooperThresholds().

◆ GetThresholdTrials()

G4double G4TransportationLogger::GetThresholdTrials ( ) const
inline

Definition at line 78 of file G4TransportationLogger.hh.

78{ return fThldTrials; }

Referenced by ReportLooperThresholds(), and ReportLoopingTrack().

◆ GetThresholdWarningEnergy()

G4double G4TransportationLogger::GetThresholdWarningEnergy ( ) const
inline

Definition at line 76 of file G4TransportationLogger.hh.

76{ return fThldWarningEnergy; }

Referenced by ReportLooperThresholds(), and ReportLoopingTrack().

◆ GetVerboseLevel()

G4int G4TransportationLogger::GetVerboseLevel ( ) const
inline

Definition at line 71 of file G4TransportationLogger.hh.

71{ return fVerbose; }

◆ ReportLooperThresholds()

void G4TransportationLogger::ReportLooperThresholds ( const char * className)

Definition at line 75 of file G4TransportationLogger.cc.

76{
77 G4cout << className << ": Current values for thresholds related to "
78 << " the killing of looping tracks: " << G4endl
79 << " Warning Energy = " << GetThresholdWarningEnergy() / CLHEP::MeV << " MeV "
80 << " ( below this tracks are killed without warning ) " << G4endl
81 << " Important Energy = " << GetThresholdImportantEnergy() / CLHEP::MeV
82 << " ( above this tracks are given multiple chances ) " << G4endl
83 << " Extra Trials = " << GetThresholdTrials()
84 << " 'important' tracks, i.e. those above 'important' energy "
85 << G4endl;
86}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetThresholdImportantEnergy() const
G4double GetThresholdWarningEnergy() const

Referenced by G4Transportation::ReportLooperThresholds().

◆ ReportLoopingTrack()

void G4TransportationLogger::ReportLoopingTrack ( const G4Track & track,
const G4Step & stepInfo,
G4int numTrials,
G4long noCalls,
const char * methodName ) const

Definition at line 92 of file G4TransportationLogger.cc.

98{
99 static std::atomic<unsigned int> numAdviceExcessSteps(0);
100 constexpr double gramPerCm3 = gram / ( centimeter * centimeter * centimeter ) ;
101 std::ostringstream msg;
102 auto preStepPt= stepData.GetPreStepPoint();
103 auto preStepEn= preStepPt ? preStepPt->GetKineticEnergy() / MeV : -1.0 ;
104 msg << " Transportation is killing track that is looping or stuck. " << G4endl
105 << " Track is "
107 << " and has " << track.GetKineticEnergy() / MeV
108 << " MeV energy ( pre-Step = " << preStepEn << " ) " << G4endl;
109 msg << " momentum = " << track.GetMomentum() << " mag= " << track.GetMomentum().mag()
110 << G4endl
111 << " position = " << track.GetPosition();
112 auto physVolume= track.GetVolume();
113 auto material= physVolume->GetLogicalVolume()->GetMaterial();
114 msg << " is in volume '" << physVolume->GetName() << "', ";
115 if( material )
116 {
117 msg << " its material is '" << material->GetName() << "'";
118 msg << " with density = " << material->GetDensity() / gramPerCm3
119 << " g/cm^3 ";
120 }
121 else
122 {
123 msg << " unable to obtain material information (including density.) ";
124 }
125 msg << G4endl;
126 msg << " Total number of Steps by this track: " << track.GetCurrentStepNumber()
127 << G4endl
128 << " Length of this step = " << stepData.GetStepLength() / mm << " mm "
129 << G4endl
130 << " Number of propagation trials = " << numTrials
131 << " ( vs maximum = " << GetThresholdTrials() << " for 'important' particles ) "
132 << G4endl;
133
134 if (noCalls)
135 msg << " ( Number of *calls* of Transport/AlongStepDoIt = " << noCalls << " )" << G4endl;
136
137 const G4int numPrints= 5;
138 if( numAdviceExcessSteps++ < numPrints )
139 {
140 msg << " =============== Recommendations / advice ====================" << G4endl;
141 msg << " Recommendations to address this issue (Transport-001-ExcessSteps)" << G4endl;
142 msg << " This warning is controlled by the SetThresholdWarningEnergy "
143 << " method of G4Transportation. " << G4endl
144 << " Current value of 'warning' threshold= "
145 << GetThresholdWarningEnergy() / MeV << " MeV " << G4endl;
146 msg << " - If 'unimportant' particles (with energy low enough not to matter in your "
147 << " application, then increase its value. " << G4endl;
148 msg << " - If particles of high-enough energy to be important are being "
149 << " killed, you can " << G4endl
150 << " a) Increase the trial steps using the method SetThresholdTrials(). "
151 << " Particles above the 'important' threshold " << G4endl
152 << " will be given this many 'chances'."
153 << " The default value was 10, and the current value is " << GetThresholdTrials()
154 << G4endl
155 << " b) Increase the energy which you consider 'important' (above this they are"
156 << " killed only after extra trials), using the method SetThresholdImportantEnergy() " << G4endl
157 << " Note: this can incur a potentially high cost in extra simulation time "
158 << " if more tracks require very large number of integration steps . " << G4endl
159 << " c) investigate alternative integration methods " << G4endl
160 << " e.g. Helical methods for uniform or almost uniform fields"
161 << " or else higher order RK methods such as DormandPrince78 "
162 << G4endl;
163 msg << " This information is provided " << numPrints << " times. Current count: "
164 << numAdviceExcessSteps << " / " << numPrints << G4endl;
165 msg << " =============================================================" << G4endl;
166 }
167 const G4String fullMethodName= fClassName + "::" + methodName;
168 G4Exception( fullMethodName, "Transport-001-ExcessSteps", JustWarning, msg);
169}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
int G4int
Definition G4Types.hh:85
double mag() const
G4Material * GetMaterial() const
const G4String & GetName() const
const G4String & GetParticleName() const
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
G4ThreeVector GetMomentum() const
G4double GetKineticEnergy() const
G4LogicalVolume * GetLogicalVolume() const

Referenced by G4Transportation::AlongStepDoIt().

◆ SetThresholdImportantEnergy()

void G4TransportationLogger::SetThresholdImportantEnergy ( G4double val)
inline

Definition at line 81 of file G4TransportationLogger.hh.

81{ fThldImportantEnergy= val; }

Referenced by SetThresholds().

◆ SetThresholds()

void G4TransportationLogger::SetThresholds ( G4double newEnWarn,
G4double importantEnergy,
G4int newMaxTrials )

Definition at line 63 of file G4TransportationLogger.cc.

66{
67 SetThresholdWarningEnergy( newEnWarn );
68 SetThresholdImportantEnergy( importantEnergy );
69 SetThresholdTrials(newMaxTrials );
70}
void SetThresholdWarningEnergy(G4double val)
void SetThresholdImportantEnergy(G4double val)
void SetThresholdTrials(G4int maxNoTrials)

◆ SetThresholdTrials()

void G4TransportationLogger::SetThresholdTrials ( G4int maxNoTrials)
inline

Definition at line 82 of file G4TransportationLogger.hh.

82{ fThldTrials = std::max( maxNoTrials, 1); }

Referenced by SetThresholds().

◆ SetThresholdWarningEnergy()

void G4TransportationLogger::SetThresholdWarningEnergy ( G4double val)
inline

Definition at line 80 of file G4TransportationLogger.hh.

80{ fThldWarningEnergy= val; }

Referenced by SetThresholds().

◆ SetVerboseLevel()

void G4TransportationLogger::SetVerboseLevel ( G4int level)
inline

Definition at line 72 of file G4TransportationLogger.hh.

72{ fVerbose = level; }

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