Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TransportationLogger.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28//
29// class G4TransportationLogger Implementation
30//
31// Author: J. Apostolakis, June 2018
32//
33// --------------------------------------------------------------------
34
35#include <iomanip>
36
37#include "G4SystemOfUnits.hh"
39#include "G4Track.hh"
40#include "G4Step.hh"
41
43 G4int verbosity)
44 : fClassName(className), fVerbose(verbosity),
45 fThldWarningEnergy(0.), fThldImportantEnergy(0.), fThldTrials(0)
46{
47}
48
50 G4int verbosity)
51 : fClassName(className), fVerbose(verbosity),
52 fThldWarningEnergy(0.), fThldImportantEnergy(0.), fThldTrials(0)
53{
54}
55
57{
58}
59
60// ********************************************************************
61// SetThresholds
62// ********************************************************************
64SetThresholds( G4double newEnWarn, G4double importantEnergy,
65 G4int newMaxTrials )
66{
67 SetThresholdWarningEnergy( newEnWarn );
68 SetThresholdImportantEnergy( importantEnergy );
69 SetThresholdTrials(newMaxTrials );
70}
71
72/////////////////////////////////////////////////////////////////////////////
73//
74void
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}
87
88// ********************************************************************
89// ReportLoopingTrack
90// ********************************************************************
91//
93 const G4Step & stepData,
94 G4int numTrials,
95 G4long noCalls,
96 const char* methodName
97 ) const
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 << " ( Number of *calls* of Transport/AlongStepDoIt = " << noCalls << " )"
134 << G4endl;
135
136 const G4int numPrints= 5;
137 if( numAdviceExcessSteps++ < numPrints )
138 {
139 msg << " =============== Recommendations / advice ====================" << G4endl;
140 msg << " Recommendations to address this issue (Transport-001-ExcessSteps)" << G4endl;
141 msg << " This warning is controlled by the SetThresholdWarningEnergy "
142 << " method of G4Transportation. " << G4endl
143 << " Current value of 'warning' threshold= "
144 << GetThresholdWarningEnergy() / MeV << " MeV " << G4endl;
145 msg << " - If 'unimportant' particles (with energy low enough not to matter in your "
146 << " application, then increase its value. " << G4endl;
147 msg << " - If particles of high-enough energy to be important are being "
148 << " killed, you can " << G4endl
149 << " a) Increase the trial steps using the method SetThresholdTrials(). "
150 << " Particles above the 'important' threshold " << G4endl
151 << " will be given this many 'chances'."
152 << " The default value was 10, and the current value is " << GetThresholdTrials()
153 << G4endl
154 << " b) Increase the energy which you consider 'important' (above this they are"
155 << " killed only after extra trials), using the method SetThresholdImportantEnergy() " << G4endl
156 << " Note: this can incur a potentially high cost in extra simulation time "
157 << " if more tracks require very large number of integration steps . " << G4endl
158 << " c) investigate alternative integration methods " << G4endl
159 << " e.g. Helical methods for uniform or almost uniform fields"
160 << " or else higher order RK methods such as DormandPrince78 "
161 << G4endl;
162 msg << " This information is provided " << numPrints << " times. Current count: "
163 << numAdviceExcessSteps << " / " << numPrints << G4endl;
164 msg << " =============================================================" << G4endl;
165 }
166 const G4String fullMethodName= fClassName + "::" + methodName;
167 G4Exception( fullMethodName, "Transport-001-ExcessSteps", JustWarning, msg);
168}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag() const
G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:175
const G4String & GetParticleName() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
G4ThreeVector GetMomentum() const
G4double GetKineticEnergy() const
G4TransportationLogger(const G4String &className, G4int verbosity)
void SetThresholdWarningEnergy(G4double val)
void SetThresholdImportantEnergy(G4double val)
G4double GetThresholdImportantEnergy() const
void SetThresholdTrials(G4int maxNoTrials)
G4double GetThresholdTrials() const
G4double GetThresholdWarningEnergy() const
void SetThresholds(G4double newEnWarn, G4double importantEnergy, G4int newMaxTrials)
void ReportLoopingTrack(const G4Track &track, const G4Step &stepInfo, G4int numTrials, G4long noCalls, const char *methodName) const
void ReportLooperThresholds(const char *className)
G4LogicalVolume * GetLogicalVolume() const