Loading...
Searching...
No Matches
Hit.hh
Go to the documentation of this file.
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
42
43namespace 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
52class 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; }
79 void SetVoxelIndexes(const G4TouchableHandle& TH);
80
81 // Get methods
82 G4int GetTrackID() const { return fTrackID; }
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; }
90 G4int GetXindex() { return fxIndex; }
91 G4int GetYindex() { return fyIndex; }
92 G4int GetZindex() { return fzIndex; }
93
94 private:
95 // Variables
96 G4int fTrackID = -1;
98 G4double fEkinMean = 0.;
99 G4double fDeltaE = 0.;
100 G4double fEinit = 0.;
101 G4double fTrackLength = 0.;
102 G4double fElectronEnergy = 0.;
104 G4int fxIndex = -1;
105 G4int fyIndex = -1;
106 G4int fzIndex = -1;
107};
108
109// Definition of a HitColletion
111
112// Definition of the Allocator
113extern G4ThreadLocal G4Allocator<Hit>* RadioBioHitAllocator;
114
115inline void* Hit::operator new(size_t)
116{
118 return (void*)RadioBioHitAllocator->MallocSingle();
119}
120
121inline void Hit::operator delete(void* hit)
122{
123 RadioBioHitAllocator->FreeSingle((Hit*)hit);
124}
125
126} // namespace RadioBio
127
128#endif // RadioBioHit_h
Detector hit class.
G4int fzIndex
Definition Hit.hh:106
void SetTrackID(G4int track)
Definition Hit.hh:71
G4VPhysicalVolume * fPhysicalVolume
Definition Hit.hh:103
G4int operator==(const Hit &) const
Definition Hit.cc:88
void SetEinit(G4double e)
Definition Hit.hh:75
G4double GetElectronEnergy() const
Definition Hit.hh:88
void SetElectronEnergy(G4double elEnergy)
Definition Hit.hh:77
G4int fyIndex
Definition Hit.hh:105
void Print() override
Definition Hit.cc:102
G4int GetXindex()
Definition Hit.hh:90
void SetPartType(const G4ParticleDefinition *part)
Definition Hit.hh:72
const G4ParticleDefinition * GetPartType() const
Definition Hit.hh:83
G4double fEinit
Definition Hit.hh:100
G4double fElectronEnergy
Definition Hit.hh:102
const G4ParticleDefinition * fParticleType
Definition Hit.hh:97
void SetDeltaE(double DeltaE)
Definition Hit.hh:74
G4VPhysicalVolume * GetPhysicalVolume() const
Definition Hit.hh:89
G4double GetEkinMean()
Definition Hit.hh:84
G4int GetTrackID() const
Definition Hit.hh:82
G4double GetTrackLength() const
Definition Hit.hh:87
G4double GetDeltaE()
Definition Hit.hh:85
G4double fDeltaE
Definition Hit.hh:99
G4double fEkinMean
Definition Hit.hh:98
G4int fxIndex
Definition Hit.hh:104
const Hit & operator=(const Hit &)
Definition Hit.cc:70
G4int GetZindex()
Definition Hit.hh:92
G4int GetYindex()
Definition Hit.hh:91
~Hit() override=default
void SetEkinMean(double EkinMean)
Definition Hit.hh:73
void SetVoxelIndexes(const G4TouchableHandle &TH)
Definition Hit.cc:109
G4double GetEinit() const
Definition Hit.hh:86
G4int fTrackID
Definition Hit.hh:96
G4double fTrackLength
Definition Hit.hh:101
void SetTrackLength(G4double x)
Definition Hit.hh:76
void SetPhysicalVolume(G4VPhysicalVolume *PV)
Definition Hit.hh:78
void Draw() override
Definition Hit.cc:95
G4ThreadLocal G4Allocator< Hit > * RadioBioHitAllocator
Definition Hit.cc:48
G4THitsCollection< Hit > RadioBioHitsCollection
Definition Hit.hh:110

Applications | User Support | Publications | Collaboration