Geant4-11
Public Member Functions | Private Attributes
G4TransportationLogger Class Reference

#include <G4TransportationLogger.hh>

Public Member Functions

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

Private Attributes

G4String fClassName
 
G4double fThldImportantEnergy
 
G4int fThldTrials
 
G4double fThldWarningEnergy
 
G4int fVerbose
 

Detailed Description

Definition at line 47 of file G4TransportationLogger.hh.

Constructor & Destructor Documentation

◆ G4TransportationLogger() [1/2]

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

◆ G4TransportationLogger() [2/2]

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

Definition at line 49 of file G4TransportationLogger.cc.

◆ ~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; }

References fThldImportantEnergy.

Referenced by ReportLooperThresholds().

◆ GetThresholdTrials()

G4double G4TransportationLogger::GetThresholdTrials ( ) const
inline

Definition at line 78 of file G4TransportationLogger.hh.

78{ return fThldTrials; }

References fThldTrials.

Referenced by ReportLooperThresholds(), and ReportLoopingTrack().

◆ GetThresholdWarningEnergy()

G4double G4TransportationLogger::GetThresholdWarningEnergy ( ) const
inline

Definition at line 76 of file G4TransportationLogger.hh.

76{ return fThldWarningEnergy; }

References fThldWarningEnergy.

Referenced by ReportLooperThresholds(), and ReportLoopingTrack().

◆ GetVerboseLevel()

G4int G4TransportationLogger::GetVerboseLevel ( ) const
inline

Definition at line 71 of file G4TransportationLogger.hh.

71{ return fVerbose; }

References 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:57
G4GLOB_DLL std::ostream G4cout
G4double GetThresholdImportantEnergy() const
G4double GetThresholdTrials() const
G4double GetThresholdWarningEnergy() const
static constexpr double MeV

References className, G4cout, G4endl, GetThresholdImportantEnergy(), GetThresholdTrials(), GetThresholdWarningEnergy(), and CLHEP::MeV.

Referenced by G4CoupledTransportation::ReportLooperThresholds(), and 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)
Definition: G4Exception.cc:35
static constexpr double centimeter
Definition: G4SIunits.hh:70
static constexpr double gram
Definition: G4SIunits.hh:163
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double MeV
Definition: G4SIunits.hh:200
int G4int
Definition: G4Types.hh:85
double mag() 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
string material
Definition: eplot.py:19

References centimeter, fClassName, G4endl, G4Exception(), G4Track::GetCurrentStepNumber(), G4StepPoint::GetKineticEnergy(), G4Track::GetKineticEnergy(), G4Track::GetMomentum(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), G4Step::GetPreStepPoint(), G4Step::GetStepLength(), GetThresholdTrials(), GetThresholdWarningEnergy(), G4Track::GetVolume(), gram, JustWarning, CLHEP::Hep3Vector::mag(), eplot::material, MeV, and mm.

Referenced by G4CoupledTransportation::AlongStepDoIt(), and G4Transportation::AlongStepDoIt().

◆ SetThresholdImportantEnergy()

void G4TransportationLogger::SetThresholdImportantEnergy ( G4double  val)
inline

Definition at line 81 of file G4TransportationLogger.hh.

References fThldImportantEnergy.

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)

References SetThresholdImportantEnergy(), SetThresholdTrials(), and SetThresholdWarningEnergy().

◆ SetThresholdTrials()

void G4TransportationLogger::SetThresholdTrials ( G4int  maxNoTrials)
inline

Definition at line 82 of file G4TransportationLogger.hh.

82{ fThldTrials = std::max( maxNoTrials, 1); }
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References fThldTrials, and G4INCL::Math::max().

Referenced by SetThresholds().

◆ SetThresholdWarningEnergy()

void G4TransportationLogger::SetThresholdWarningEnergy ( G4double  val)
inline

Definition at line 80 of file G4TransportationLogger.hh.

80{ fThldWarningEnergy= val; }

References fThldWarningEnergy.

Referenced by SetThresholds().

◆ SetVerboseLevel()

void G4TransportationLogger::SetVerboseLevel ( G4int  level)
inline

Definition at line 72 of file G4TransportationLogger.hh.

72{ fVerbose = level; }

References fVerbose.

Field Documentation

◆ fClassName

G4String G4TransportationLogger::fClassName
private

Definition at line 86 of file G4TransportationLogger.hh.

Referenced by ReportLoopingTrack().

◆ fThldImportantEnergy

G4double G4TransportationLogger::fThldImportantEnergy
private

◆ fThldTrials

G4int G4TransportationLogger::fThldTrials
private

Definition at line 94 of file G4TransportationLogger.hh.

Referenced by GetThresholdTrials(), and SetThresholdTrials().

◆ fThldWarningEnergy

G4double G4TransportationLogger::fThldWarningEnergy
private

◆ fVerbose

G4int G4TransportationLogger::fVerbose
private

Definition at line 87 of file G4TransportationLogger.hh.

Referenced by GetVerboseLevel(), and SetVerboseLevel().


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