File indexing completed on 2025-01-18 09:55:24
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013 #ifndef DDG4_GEANT4STEPHANDLER_H
0014 #define DDG4_GEANT4STEPHANDLER_H
0015
0016
0017 #include <DDG4/Geant4HitHandler.h>
0018
0019
0020 #include <G4Step.hh>
0021 #include <G4StepPoint.hh>
0022 #include <G4VTouchable.hh>
0023 #include <G4VSensitiveDetector.hh>
0024 #include <G4EmSaturation.hh>
0025 #include <G4Version.hh>
0026
0027
0028 namespace dd4hep {
0029
0030
0031 namespace sim {
0032
0033
0034 class Geant4StepHandler;
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046 class Geant4StepHandler : public Geant4HitHandler {
0047 public:
0048 const G4Step* step;
0049 G4StepPoint* pre;
0050 G4StepPoint* post;
0051 bool applyBirksLaw;
0052
0053 Geant4StepHandler() = delete;
0054
0055 Geant4StepHandler(const G4Step* s)
0056 : Geant4HitHandler(s->GetTrack(), (s->GetPreStepPoint()->GetTouchableHandle())()),
0057 step(s), pre(s->GetPreStepPoint()), post(s->GetPostStepPoint())
0058 {
0059 applyBirksLaw = false;
0060 }
0061
0062 Geant4StepHandler(const Geant4StepHandler& copy) = delete;
0063
0064 Geant4StepHandler(Geant4StepHandler&& copy) = delete;
0065
0066 Geant4StepHandler& operator=(const Geant4StepHandler& copy) = delete;
0067
0068 Geant4StepHandler& operator=(Geant4StepHandler&& copy) = delete;
0069
0070
0071 static const char* stepStatus(G4StepStatus status);
0072
0073 const char* preStepStatus() const;
0074
0075 const char* postStepStatus() const;
0076
0077 double totalEnergy() const {
0078 if(applyBirksLaw == true)
0079 return this->birkAttenuation();
0080 else
0081 return step->GetTotalEnergyDeposit();
0082 }
0083
0084 Position prePos() const {
0085 const G4ThreeVector& p = pre->GetPosition();
0086 return Position(p.x(), p.y(), p.z());
0087 }
0088
0089 const G4ThreeVector& prePosG4() const {
0090 return pre->GetPosition();
0091 }
0092
0093 Position postPos() const {
0094 const G4ThreeVector& p = post->GetPosition();
0095 return Position(p.x(), p.y(), p.z());
0096 }
0097
0098 const G4ThreeVector& postPosG4() const {
0099 return post->GetPosition();
0100 }
0101
0102 G4ThreeVector avgPositionG4() const {
0103 const G4ThreeVector& p1 = pre->GetPosition();
0104 const G4ThreeVector& p2 = post->GetPosition();
0105 G4ThreeVector r = 0.5*(p2+p1);
0106 return r;
0107 }
0108
0109 Position avgPosition() const {
0110 const G4ThreeVector& p1 = pre->GetPosition();
0111 const G4ThreeVector& p2 = post->GetPosition();
0112 G4ThreeVector r = 0.5*(p2+p1);
0113 return Position(r.x(),r.y(),r.z());
0114 }
0115
0116 Position direction() const {
0117 const G4ThreeVector& p1 = pre->GetPosition();
0118 const G4ThreeVector& p2 = post->GetPosition();
0119 G4ThreeVector r = (p2-p1);
0120 double len = r.r();
0121 if ( len > 1e-15 ) {
0122 r /= len;
0123 return Position(r.x(),r.y(),r.z());
0124 }
0125 return Position();
0126 }
0127 Momentum preMom() const {
0128 const G4ThreeVector& p = pre->GetMomentum();
0129 return Momentum(p.x(), p.y(), p.z());
0130 }
0131 Momentum postMom() const {
0132 const G4ThreeVector& p = post->GetMomentum();
0133 return Momentum(p.x(), p.y(), p.z());
0134 }
0135 double deposit() const {
0136 return step->GetTotalEnergyDeposit();
0137 }
0138 double stepLength() const {
0139 return step->GetStepLength();
0140 }
0141 const G4VTouchable* preTouchable() const {
0142 return pre->GetTouchable();
0143 }
0144 const G4VTouchable* postTouchable() const {
0145 return post->GetTouchable();
0146 }
0147 const char* volName(const G4StepPoint* p, const char* undefined = "") const {
0148 G4VPhysicalVolume* v = volume(p);
0149 return v ? v->GetName().c_str() : undefined;
0150 }
0151 G4VPhysicalVolume* volume(const G4StepPoint* p) const {
0152 return p->GetTouchableHandle()->GetVolume();
0153 }
0154 G4VSolid* solid(const G4StepPoint* p) const {
0155 return p->GetTouchableHandle()->GetSolid();
0156 }
0157 G4VPhysicalVolume* physvol(const G4StepPoint* p) const {
0158 return p->GetPhysicalVolume();
0159 }
0160 G4LogicalVolume* logvol(const G4StepPoint* p) const {
0161 G4VPhysicalVolume* pv = physvol(p);
0162 return pv ? pv->GetLogicalVolume() : 0;
0163 }
0164 G4VSensitiveDetector* sd(const G4StepPoint* p) const {
0165 G4LogicalVolume* lv = logvol(p);
0166 return lv ? lv->GetSensitiveDetector() : 0;
0167 }
0168 std::string sdName(const G4StepPoint* p, const std::string& undefined = "") const {
0169 G4VSensitiveDetector* s = sd(p);
0170 return s ? s->GetName().c_str() : undefined;
0171 }
0172 bool isSensitive(const G4StepPoint* point) const {
0173 return point ? this->Geant4HitHandler::isSensitive(volume(point)) : false;
0174 }
0175 G4VPhysicalVolume* preVolume() const {
0176 return volume(pre);
0177 }
0178 G4VSensitiveDetector* preSD() const {
0179 return sd(pre);
0180 }
0181 G4VPhysicalVolume* postVolume() const {
0182 return volume(post);
0183 }
0184 G4VSensitiveDetector* postSD() const {
0185 return sd(post);
0186 }
0187 bool firstInVolume() const {
0188 return step->IsFirstStepInVolume();
0189 }
0190 bool lastInVolume() const {
0191 return step->IsLastStepInVolume();
0192 }
0193
0194 double birkAttenuation() const;
0195
0196 void doApplyBirksLaw(void) {
0197 applyBirksLaw = true;
0198 }
0199 };
0200
0201 }
0202 }
0203
0204 #endif