Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4TransportationLogger.cc
이 파일의 문서화 페이지로 가기
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 // $Id: G4TransportationLogger.cc 110824 2018-06-15 13:25:52Z gcosmo $
28 //
29 //
30 // class G4TransportationLogger Implementation
31 //
32 // Author: J. Apostolakis, June 2018
33 //
34 // --------------------------------------------------------------------
35 
36 #include <iomanip>
37 
38 #include "G4SystemOfUnits.hh"
40 #include "G4Track.hh"
41 #include "G4Step.hh"
42 
44  G4int verbosity)
45  : fClassName(className), fVerbose(verbosity),
46  fThldWarningEnergy(0.), fThldImportantEnergy(0.), fThldTrials(0)
47 {
48 }
49 
51  G4int verbosity)
52  : fClassName(className), fVerbose(verbosity),
53  fThldWarningEnergy(0.), fThldImportantEnergy(0.), fThldTrials(0)
54 {
55 }
56 
58 {
59 }
60 
61 // ********************************************************************
62 // SetThresholds
63 // ********************************************************************
65 SetThresholds( G4double newEnWarn, G4double importantEnergy,
66  G4int newMaxTrials )
67 {
68  SetThresholdWarningEnergy( newEnWarn );
69  SetThresholdImportantEnergy( importantEnergy );
70  SetThresholdTrials(newMaxTrials );
71 }
72 
73 // ********************************************************************
74 // ReportLoopingTrack
75 // ********************************************************************
76 //
78  const G4Step & stepData,
79  G4int numTrials,
80  G4long noCalls,
81  const char* methodName
82  ) const
83 {
84  static std::atomic<unsigned int> numAdviceExcessSteps(0);
85  constexpr double gramPerCm3 = gram / ( centimeter * centimeter * centimeter ) ;
86  std::ostringstream msg;
87  auto preStepPt= stepData.GetPreStepPoint();
88  auto preStepEn= preStepPt ? preStepPt->GetKineticEnergy() / MeV : -1.0 ;
89  msg << " Transportation is killing track that is looping or stuck. " << G4endl
90  << " Track is "
92  << " and has " << track.GetKineticEnergy() / MeV
93  << " MeV energy ( pre-Step = " << preStepEn << " ) " << G4endl;
94  msg << " momentum = " << track.GetMomentum() << " mag= " << track.GetMomentum().mag()
95  << G4endl
96  << " position = " << track.GetPosition();
97  auto physVolume= track.GetVolume();
98  auto material= physVolume->GetLogicalVolume()->GetMaterial();
99  msg << " is in volume '" << physVolume->GetName() << "', ";
100  if( material )
101  {
102  msg << " its material is '" << material->GetName() << "'";
103  msg << " with density = " << material->GetDensity() / gramPerCm3
104  << " g/cm^3 ";
105  }
106  else
107  {
108  msg << " unable to obtain material information (including density.) ";
109  }
110  msg << G4endl;
111  msg << " Total number of Steps by this track: " << track.GetCurrentStepNumber()
112  << G4endl
113  << " Length of this step = " << stepData.GetStepLength() / mm << " mm "
114  << G4endl
115  << " Number of propagation trials = " << numTrials
116  << " ( vs maximum = " << GetThresholdTrials() << " for 'important' particles ) "
117  << G4endl
118  << " ( Number of *calls* of Transport/AlongStepDoIt = " << noCalls << " )"
119  << G4endl;
120 
121  const G4int numPrints= 5;
122  if( numAdviceExcessSteps++ < numPrints )
123  {
124  msg << " =============== Recommendations / advice ====================" << G4endl;
125  msg << " Recommendations to address this issue (Transport-001-ExcessSteps)" << G4endl;
126  msg << " This warning is controlled by the SetThresholdWarningEnergy "
127  << " method of G4Transportation. " << G4endl
128  << " Current value of 'warning' threshold= "
129  << GetThresholdWarningEnergy() / MeV << " MeV " << G4endl;
130  msg << " - If 'unimportant' particles (with energy low enough not to matter in your "
131  << " application, then increase its value. " << G4endl;
132  msg << " - If particles of high-enough energy to be important are being "
133  << " killed, you can " << G4endl
134  << " a) Increase the trial steps using the method SetThresholdTrials(). "
135  << " Particles above the 'important' threshold " << G4endl
136  << " will be given this many 'chances'."
137  << " The default value was 10, and the current value is " << GetThresholdTrials()
138  << G4endl
139  << " b) Increase the energy which you consider 'important' (above this they are"
140  << " killed only after extra trials), using the method SetThresholdImportantEnergy() " << G4endl
141  << " Note: this can incur a potentially high cost in extra simulation time "
142  << " if more tracks require very large number of integration steps . " << G4endl
143  << " c) investigate alternative integration methods " << G4endl
144  << " e.g. Helical methods for uniform or almost uniform fields"
145  << " or else higher order RK methods such as DormandPrince78 "
146  << G4endl;
147  msg << " This information is provided " << numPrints << " times. Current count: "
148  << numAdviceExcessSteps << " / " << numPrints << G4endl;
149  msg << " =============================================================" << G4endl;
150  }
151  const G4String fullMethodName= fClassName + "::" + methodName;
152  G4Exception( fullMethodName, "Transport-001-ExcessSteps", JustWarning, msg);
153 }
G4double GetKineticEnergy() const
void SetThresholdTrials(G4int maxNoTrials)
G4LogicalVolume * GetLogicalVolume() const
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetThresholds(G4double newEnWarn, G4double importantEnergy, G4int newMaxTrials)
G4StepPoint * GetPreStepPoint() const
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
G4TransportationLogger(const G4String &className, G4int verbosity)
const G4ParticleDefinition * GetParticleDefinition() const
G4int GetCurrentStepNumber() const
void SetThresholdImportantEnergy(G4double val)
G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:179
G4ThreeVector GetMomentum() const
double G4double
Definition: G4Types.hh:76
long G4long
Definition: G4Types.hh:80
G4double GetStepLength() const
const G4ThreeVector & GetPosition() const
G4double GetThresholdWarningEnergy() const
Definition: G4Step.hh:76
G4double GetKineticEnergy() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
double mag() const
int G4int
Definition: G4Types.hh:78
void SetThresholdWarningEnergy(G4double val)
static constexpr double centimeter
Definition: G4SIunits.hh:90
void ReportLoopingTrack(const G4Track &track, const G4Step &stepInfo, G4int numTrials, G4long noCalls, const char *methodName) const
G4double GetThresholdTrials() const
G4VPhysicalVolume * GetVolume() const
static constexpr double gram
Definition: G4SIunits.hh:178