Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 /// \file radiobiology/include/Hit.hh 27 /// \brief Definition of the RadioBio::Hit class 28 29 #ifndef RadiobiologyHit_h 30 #define RadiobiologyHit_h 1 31 32 #include "G4Allocator.hh" 33 #include "G4THitsCollection.hh" 34 #include "G4ThreeVector.hh" 35 #include "G4TouchableHandle.hh" 36 #include "G4VHit.hh" 37 #include "G4VPhysicalVolume.hh" 38 39 #include "tls.hh" 40 41 class G4ParticleDefinition; 42 43 namespace RadioBio 44 { 45 46 /// Detector hit class 47 /// 48 /// It defines data members to store the trackID, particle type, 49 /// mean kinetic energy, energy deposit, initial energy, track length, 50 /// electronic energy, physical volume and voxel index 51 52 class Hit : public G4VHit 53 { 54 public: 55 Hit(); 56 Hit(const Hit&); 57 ~Hit() override = default; 58 59 // Operators 60 const Hit& operator=(const Hit&); 61 G4int operator==(const Hit&) const; 62 63 inline void* operator new(size_t); 64 inline void operator delete(void*); 65 66 // Methods from base class 67 void Draw() override; 68 void Print() override; 69 70 // Set methods 71 void SetTrackID(G4int track) { fTrackID = track; } 72 void SetPartType(const G4ParticleDefinition* part) { fParticleType = part; } 73 void SetEkinMean(double EkinMean) { fEkinMean = EkinMean; } 74 void SetDeltaE(double DeltaE) { fDeltaE = DeltaE; } 75 void SetEinit(G4double e) { fEinit = e; } 76 void SetTrackLength(G4double x) { fTrackLength = x; } 77 void SetElectronEnergy(G4double elEnergy) { fElectronEnergy = elEnergy; } 78 void SetPhysicalVolume(G4VPhysicalVolume* PV) { fPhysicalVolume = PV; } 79 void SetVoxelIndexes(const G4TouchableHandle& TH); 80 81 // Get methods 82 G4int GetTrackID() const { return fTrackID; } 83 const G4ParticleDefinition* GetPartType() const { return fParticleType; } 84 G4double GetEkinMean() { return fEkinMean; } 85 G4double GetDeltaE() { return fDeltaE; } 86 G4double GetEinit() const { return fEinit; } 87 G4double GetTrackLength() const { return fTrackLength; } 88 G4double GetElectronEnergy() const { return fElectronEnergy; } 89 G4VPhysicalVolume* GetPhysicalVolume() const { return fPhysicalVolume; } 90 G4int GetXindex() { return fxIndex; } 91 G4int GetYindex() { return fyIndex; } 92 G4int GetZindex() { return fzIndex; } 93 94 private: 95 // Variables 96 G4int fTrackID = -1; 97 const G4ParticleDefinition* fParticleType = nullptr; 98 G4double fEkinMean = 0.; 99 G4double fDeltaE = 0.; 100 G4double fEinit = 0.; 101 G4double fTrackLength = 0.; 102 G4double fElectronEnergy = 0.; 103 G4VPhysicalVolume* fPhysicalVolume = nullptr; 104 G4int fxIndex = -1; 105 G4int fyIndex = -1; 106 G4int fzIndex = -1; 107 }; 108 109 // Definition of a HitColletion 110 typedef G4THitsCollection<Hit> RadioBioHitsCollection; 111 112 // Definition of the Allocator 113 extern G4ThreadLocal G4Allocator<Hit>* RadioBioHitAllocator; 114 115 inline void* Hit::operator new(size_t) 116 { 117 if (!RadioBioHitAllocator) RadioBioHitAllocator = new G4Allocator<Hit>; 118 return (void*)RadioBioHitAllocator->MallocSingle(); 119 } 120 121 inline void Hit::operator delete(void* hit) 122 { 123 RadioBioHitAllocator->FreeSingle((Hit*)hit); 124 } 125 126 } // namespace RadioBio 127 128 #endif // RadioBioHit_h 129