CMS 3D CMS Logo

SteppingAction.h
Go to the documentation of this file.
1 #ifndef SimG4Core_SteppingAction_H
2 #define SimG4Core_SteppingAction_H
3 
6 
7 #include "G4LogicalVolume.hh"
8 #include "G4Region.hh"
9 #include "G4UserSteppingAction.hh"
10 #include "G4VPhysicalVolume.hh"
11 #include "G4VTouchable.hh"
12 #include "G4Track.hh"
13 
14 #include <string>
15 #include <vector>
16 
17 class EventAction;
18 class CMSSteppingVerbose;
19 
20 enum TrackStatus {
21  sAlive = 0,
24  sOutOfTime = 3,
25  sLowEnergy = 4,
28 };
29 
30 class SteppingAction: public G4UserSteppingAction {
31 
32 public:
33  explicit SteppingAction(EventAction * ea, const edm::ParameterSet & ps,
34  const CMSSteppingVerbose*, bool hasW);
35  ~SteppingAction() override;
36 
37  void UserSteppingAction(const G4Step * aStep) final;
38 
40 
41 private:
42 
43  bool initPointer();
44 
45  bool isInsideDeadRegion(const G4Region* reg) const;
46  bool isOutOfTimeWindow(G4Track* theTrack, const G4Region* reg) const;
47  bool isThisVolume(const G4VTouchable* touch, const G4VPhysicalVolume* pv) const;
48 
49  bool isLowEnergy(const G4Step * aStep) const;
50  void PrintKilledTrack(const G4Track*, const TrackStatus&) const;
51 
53  const G4VPhysicalVolume *tracker, *calo;
57  double maxTrackTime;
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  unsigned int numberTimes;
66  unsigned int numberEkins;
67  unsigned int numberPart;
68  unsigned int ndeadRegions;
69  unsigned int nWarnings;
70 
73  bool hasWatcher;
74 };
75 
76 inline bool SteppingAction::isInsideDeadRegion(const G4Region* reg) const
77 {
78  bool res = false;
79  for(unsigned int i=0; i<ndeadRegions; ++i) {
80  if(reg == deadRegions[i]) {
81  res = true;
82  break;
83  }
84  }
85  return res;
86 }
87 
88 inline bool
89 SteppingAction::isOutOfTimeWindow(G4Track* theTrack, const G4Region* reg) const
90 {
91  double tofM = maxTrackTime;
92  for (unsigned int i=0; i<numberTimes; ++i) {
93  if (reg == maxTimeRegions[i]) {
94  tofM = maxTrackTimes[i];
95  break;
96  }
97  }
98  return (theTrack->GetGlobalTime() > tofM) ? true : false;
99 }
100 
101 inline bool SteppingAction::isThisVolume(const G4VTouchable* touch,
102  const G4VPhysicalVolume* pv) const
103 {
104  int level = (touch->GetHistoryDepth())+1;
105  return (level >= 3) ? (touch->GetVolume(level - 3) == pv) : false;
106 }
107 
108 #endif
void UserSteppingAction(const G4Step *aStep) final
const G4VPhysicalVolume * tracker
std::vector< int > ekinPDG
std::vector< const G4Region * > deadRegions
double theCriticalDensity
SimActivityRegistry::G4StepSignal m_g4StepSignal
~SteppingAction() override
unsigned int numberPart
const CMSSteppingVerbose * steppingVerbose
bool isOutOfTimeWindow(G4Track *theTrack, const G4Region *reg) const
unsigned int numberEkins
Definition: Electron.h:6
TrackStatus
std::vector< std::string > deadRegionNames
const G4VPhysicalVolume * calo
def pv(vc)
Definition: MetAnalyzer.py:7
unsigned int numberTimes
std::vector< double > maxTrackTimes
std::vector< G4LogicalVolume * > ekinVolumes
double theCriticalEnergyForVacuum
EventAction * eventAction_
std::vector< double > ekinMins
std::vector< std::string > maxTimeNames
std::vector< std::string > ekinNames
bool isLowEnergy(const G4Step *aStep) const
bool isInsideDeadRegion(const G4Region *reg) const
SteppingAction(EventAction *ea, const edm::ParameterSet &ps, const CMSSteppingVerbose *, bool hasW)
void PrintKilledTrack(const G4Track *, const TrackStatus &) const
unsigned int ndeadRegions
std::vector< std::string > ekinParticles
bool isThisVolume(const G4VTouchable *touch, const G4VPhysicalVolume *pv) const
std::vector< const G4Region * > maxTimeRegions
unsigned int nWarnings