CMS 3D CMS Logo

ExceptionHandler.cc
Go to the documentation of this file.
2 
5 
6 #include "G4EventManager.hh"
7 #include "G4TrackingManager.hh"
8 #include "G4Track.hh"
9 #include "globals.hh"
10 #include <sstream>
11 
12 ExceptionHandler::ExceptionHandler(double th, bool tr) : m_eth(th), m_trace(tr) {}
13 
15 
16 bool ExceptionHandler::Notify(const char* exceptionOrigin,
17  const char* exceptionCode,
18  G4ExceptionSeverity severity,
19  const char* description) {
20  static const G4String es_banner = "\n-------- EEEE ------- G4Exception-START -------- EEEE -------\n";
21  static const G4String ee_banner = "\n-------- EEEE -------- G4Exception-END --------- EEEE -------\n";
22  static const G4String ws_banner = "\n-------- WWWW ------- G4Exception-START -------- WWWW -------\n";
23  static const G4String we_banner = "\n-------- WWWW -------- G4Exception-END --------- WWWW -------\n";
24 
25  G4Track* track = G4EventManager::GetEventManager()->GetTrackingManager()->GetTrack();
26  double ekin = m_eth;
27 
28  std::stringstream message;
29  message << "*** G4Exception : " << exceptionCode << "\n"
30  << " issued by : " << exceptionOrigin << "\n"
31  << description;
32 
33  // part of exception happens outside tracking loop
34  if (nullptr != track) {
35  ekin = track->GetKineticEnergy();
36  message << "\n CMS info: "
37  << "TrackID=" << track->GetTrackID() << " ParentID=" << track->GetParentID() << " "
38  << track->GetParticleDefinition()->GetParticleName() << "; Ekin(MeV)=" << ekin / CLHEP::MeV
39  << "; time(ns)=" << track->GetGlobalTime() / CLHEP::ns << "; status=" << track->GetTrackStatus()
40  << "\n position(mm): " << track->GetPosition() << "; direction: " << track->GetMomentumDirection();
41  const G4VPhysicalVolume* vol = track->GetVolume();
42  if (nullptr != vol) {
43  message << "\n PhysicalVolume: " << vol->GetName() << ";";
44  const G4LogicalVolume* lv = vol->GetLogicalVolume();
45  if (nullptr != lv) {
46  message << " material: " << lv->GetMaterial()->GetName();
47  }
48  }
49  message << "\n stepNumber=" << track->GetCurrentStepNumber()
50  << "; stepLength(mm)=" << track->GetStepLength() / CLHEP::mm << "; weight=" << track->GetWeight();
51  const G4VProcess* proc = track->GetCreatorProcess();
52  if (nullptr != proc) {
53  message << "; creatorProcess: " << proc->GetProcessName() << "; modelID=" << track->GetCreatorModelID();
54  }
55  }
56  message << " \n";
57 
58  G4ExceptionSeverity localSeverity = severity;
59  std::stringstream mescode;
60  mescode << exceptionCode << "\n";
61  G4String code;
62  mescode >> code;
63 
64  // track is killed
65  if (ekin < m_eth && (code == "GeomNav0003" || code == "GeomField0003")) {
66  localSeverity = JustWarning;
67  track->SetTrackStatus(fStopAndKill);
68  ++m_number;
69  }
70 
71  bool res = false;
72  switch (localSeverity) {
73  case FatalException:
74  case FatalErrorInArgument:
75  case RunMustBeAborted:
76  case EventMustBeAborted:
77  edm::LogWarning("SimG4CoreApplication") << es_banner << message.str() << ee_banner;
78  throw cms::Exception("Geant4 fatal exception");
79  res = m_trace;
80  break;
81 
82  case JustWarning:
83  if (m_number < 20)
84  edm::LogWarning("SimG4CoreApplication")
85  << ws_banner << message.str() << "*** This is just a warning message. ***" << we_banner;
86  break;
87  }
88  return res;
89 }
~ExceptionHandler() override
Definition: Electron.h:6
bool Notify(const char *exceptionOrigin, const char *exceptionCode, G4ExceptionSeverity severity, const char *description) override
ExceptionHandler(double th, bool tr)
Log< level::Warning, false > LogWarning