DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
Geant4StepHandler.cpp
Go to the documentation of this file.
1 //==========================================================================
2 // AIDA Detector description implementation
3 //--------------------------------------------------------------------------
4 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
5 // All rights reserved.
6 //
7 // For the licensing terms see $DD4hepINSTALL/LICENSE.
8 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
9 //
10 // Author : M.Frank
11 //
12 //==========================================================================
13 
14 // Framework include files
15 #include <DDG4/Geant4StepHandler.h>
17 #include <DD4hep/DD4hepUnits.h>
18 #include <CLHEP/Units/SystemOfUnits.h>
19 
20 // Geant4 include files
21 #include <G4Version.hh>
22 
23 namespace units = dd4hep;
24 using namespace dd4hep::sim;
25 
27 const char* Geant4StepHandler::stepStatus(G4StepStatus status) {
28  switch (status) {
29  // Step reached the world boundary
30  case fWorldBoundary:
31  return "WorldBoundary";
32  // Step defined by a geometry boundary
33  case fGeomBoundary:
34  return "GeomBoundary";
35  // Step defined by a PreStepDoItVector
36  case fAtRestDoItProc:
37  return "AtRestDoItProc";
38  // Step defined by a AlongStepDoItVector
39  case fAlongStepDoItProc:
40  return "AlongStepDoItProc";
41  // Step defined by a PostStepDoItVector
42  case fPostStepDoItProc:
43  return "PostStepDoItProc";
44  // Step defined by the user Step limit in the logical volume
45  case fUserDefinedLimit:
46  return "UserDefinedLimit";
47  // Step defined by an exclusively forced PostStepDoIt process
48  case fExclusivelyForcedProc:
49  return "ExclusivelyForcedProc";
50  // Step not defined yet
51  case fUndefined:
52  default:
53  return "Undefined";
54  };
55 }
56 
58 const char* Geant4StepHandler::preStepStatus() const {
59  return stepStatus(pre ? pre->GetStepStatus() : fUndefined);
60 }
61 
63 const char* Geant4StepHandler::postStepStatus() const {
64  return stepStatus(post ? post->GetStepStatus() : fUndefined);
65 }
66 
69 #if G4VERSION_NUMBER >= 1001
70  static G4EmSaturation s_emSaturation(1);
71 #else
72  static G4EmSaturation s_emSaturation();
73  s_emSaturation.SetVerbose(1);
74 #endif
75 
76 #if G4VERSION_NUMBER >= 1030
77  static bool s_initialised = false;
78  if(not s_initialised) {
79  s_emSaturation.InitialiseG4Saturation();
80  s_initialised = true;
81  }
82 #endif
83 
84  double energyDeposition = step->GetTotalEnergyDeposit();
85  double length = step->GetStepLength();
86  double niel = step->GetNonIonizingEnergyDeposit();
87  const G4Track* trk = step->GetTrack();
88  const G4ParticleDefinition* particle = trk->GetDefinition();
89  const G4MaterialCutsCouple* couple = trk->GetMaterialCutsCouple();
90  double engyVis = s_emSaturation.VisibleEnergyDeposition(particle,
91  couple,
92  length,
93  energyDeposition,
94  niel);
95  return engyVis;
96 }
dd4hep::sim::Geant4StepHandler::pre
G4StepPoint * pre
Definition: Geant4StepHandler.h:49
Geant4StepHandler.h
dd4hep::sim::Geant4StepHandler::postStepStatus
const char * postStepStatus() const
Returns the post-step status in form of a string.
Definition: Geant4StepHandler.cpp:63
dd4hep::sim::Geant4StepHandler::post
G4StepPoint * post
Definition: Geant4StepHandler.h:50
dd4hep::sim::Geant4StepHandler::birkAttenuation
double birkAttenuation() const
Apply BirksLaw.
Definition: Geant4StepHandler.cpp:68
dd4hep::sim
Namespace for the Geant4 based simulation part of the AIDA detector description toolkit.
Definition: Geant4Output2EDM4hep.cpp:49
dd4hep::sim::Geant4StepHandler::step
const G4Step * step
Definition: Geant4StepHandler.h:48
dd4hep::sim::Geant4StepHandler::stepStatus
static const char * stepStatus(G4StepStatus status)
Returns the step status (argument) in form of a string.
Definition: Geant4StepHandler.cpp:27
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
dd4hep::sim::Geant4StepHandler::preStepStatus
const char * preStepStatus() const
Returns the pre-step status in form of a string.
Definition: Geant4StepHandler.cpp:58
Segmentation.h
DD4hepUnits.h