CMS 3D CMS Logo

SteppingAction.h
Go to the documentation of this file.
1 #ifndef SimG4Core_SteppingAction_H
2 #define SimG4Core_SteppingAction_H
3 
7 
8 #include "G4LogicalVolume.hh"
9 #include "G4Region.hh"
10 #include "G4UserSteppingAction.hh"
11 #include "G4VPhysicalVolume.hh"
12 #include "G4VTouchable.hh"
13 #include "G4Track.hh"
14 
15 #include <string>
16 #include <vector>
17 
18 class CMSSteppingVerbose;
19 
20 class SteppingAction : public G4UserSteppingAction {
21 public:
22  explicit SteppingAction(const CMSSteppingVerbose*, const edm::ParameterSet&, bool hasW);
23  ~SteppingAction() override = default;
24 
25  void UserSteppingAction(const G4Step* aStep) final;
26 
28 
29 private:
30  bool initPointer();
31 
32  inline bool isInsideDeadRegion(const G4Region* reg) const;
33  inline bool isOutOfTimeWindow(const G4Region* reg, const double& time) const;
34 
35  bool isLowEnergy(const G4LogicalVolume*, const G4Track*) const;
36  void PrintKilledTrack(const G4Track*, const TrackStatus&) const;
37 
38  const G4VPhysicalVolume* tracker{nullptr};
39  const G4VPhysicalVolume* calo{nullptr};
43  double maxTrackTime;
46 
47  unsigned int numberTimes;
48  unsigned int numberEkins;
49  unsigned int numberPart;
50  unsigned int ndeadRegions;
51  unsigned int nWarnings{0};
53 
54  bool initialized{false};
55  bool killBeamPipe{false};
56  bool hasWatcher;
57 
58  std::vector<double> maxTrackTimes, ekinMins;
59  std::vector<std::string> maxTimeNames, ekinNames, ekinParticles;
60  std::vector<std::string> deadRegionNames;
61  std::vector<const G4Region*> maxTimeRegions;
62  std::vector<const G4Region*> deadRegions;
63  std::vector<G4LogicalVolume*> ekinVolumes;
64  std::vector<int> ekinPDG;
65 };
66 
67 inline bool SteppingAction::isInsideDeadRegion(const G4Region* reg) const {
68  bool res = false;
69  for (auto& region : deadRegions) {
70  if (reg == region) {
71  res = true;
72  break;
73  }
74  }
75  return res;
76 }
77 
78 inline bool SteppingAction::isOutOfTimeWindow(const G4Region* reg, const double& time) const {
79  double tofM = maxTrackTime;
80  for (unsigned int i = 0; i < numberTimes; ++i) {
81  if (reg == maxTimeRegions[i]) {
82  tofM = maxTrackTimes[i];
83  break;
84  }
85  }
86  return (time > tofM);
87 }
88 
89 #endif
void UserSteppingAction(const G4Step *aStep) final
const G4VPhysicalVolume * tracker
std::vector< int > ekinPDG
SteppingAction(const CMSSteppingVerbose *, const edm::ParameterSet &, bool hasW)
std::vector< const G4Region * > deadRegions
double theCriticalDensity
SimActivityRegistry::G4StepSignal m_g4StepSignal
unsigned int numberPart
bool isOutOfTimeWindow(const G4Region *reg, const double &time) const
const CMSSteppingVerbose * steppingVerbose
unsigned int numberEkins
Definition: Electron.h:6
G4int maxNumberOfSteps
void PrintKilledTrack(const G4Track *, const TrackStatus &) const
std::vector< std::string > deadRegionNames
bool isLowEnergy(const G4LogicalVolume *, const G4Track *) const
unsigned int numberTimes
std::vector< double > maxTrackTimes
bool isInsideDeadRegion(const G4Region *reg) const
std::vector< G4LogicalVolume * > ekinVolumes
double theCriticalEnergyForVacuum
double maxZCentralCMS
std::vector< double > ekinMins
std::vector< std::string > maxTimeNames
std::vector< std::string > ekinNames
double maxTrackTimeForward
unsigned int ndeadRegions
~SteppingAction() override=default
std::vector< std::string > ekinParticles
TrackStatus
Definition: Common.h:9
std::vector< const G4Region * > maxTimeRegions
unsigned int nWarnings