GGS(GenericGEANT4Simulation)Software  2.6.0
 All Data Structures Namespaces Files Functions Variables Typedefs Macros
GGSIntHit.h
Go to the documentation of this file.
1 /*
2  * GGSIntHit.h
3  *
4  * Created on: 2010-10-17
5  * Authors: Franz Longo , Elena Vannuccini
6  */
7 
10 #ifndef GGSINTHIT_H
11 #define GGSINTHIT_H
12 
13 #include "G4VHit.hh"
14 #include "G4THitsCollection.hh"
15 #include "G4Allocator.hh"
16 #include "G4ThreeVector.hh"
17 #include "G4VPhysicalVolume.hh"
18 
20 #include "utils/GGSFactory.h"
21 
22 #include <numeric>
23 
31 class GGSIntHit: public G4VHit {
32 public:
33 
38  GGSIntHit(const std::vector<G4double> &timeBins);
39 
41  ~GGSIntHit();
42 
47  GGSIntHit(const GGSIntHit& right);
48 
54  const GGSIntHit& operator=(const GGSIntHit& right);
55 
61  int operator==(const GGSIntHit& right) const;
62 
70  inline void* operator new(size_t size);
71 
79  inline void operator delete(void* aHit);
80 
93  virtual void AddStep(const G4Step &step);
94 
99  void SetVolume(const G4VPhysicalVolume *volume) {
100  _volume = volume;
101  }
102 
116  void SetPosHitsStorage(bool flag);
117 
130  void SetPartHitsStorage(bool flag);
131 
136  void SetAbsTranslation(const G4ThreeVector &pos) {
137  _absTranslation = pos;
138  }
139 
170  void SetID(G4int id) {
171  _id = id;
172  }
173 
178  G4double GetTotalEnergyDep() const {
179  return accumulate(_eDep.begin(), _eDep.end(), 0.);
180  }
181 
186  const std::vector<G4double>& GetEnergyDep() const {
187  return _eDep;
188  }
189 
194  G4double GetTime() const {
195  return _time;
196  }
197 
208  const std::vector<G4double>& GetTimeBins() const {
209  return _timeBins;
210  }
211 
216  G4double GetPathLength() const {
217  return _pathLength;
218  }
219 
224  const G4VPhysicalVolume *GetVolume() {
225  return _volume;
226  }
227 
235  return _partHits;
236  }
237 
242  const G4ThreeVector& GetAbsTranslation() {
243  return _absTranslation;
244  }
245 
250  G4int GetID() {
251  return _id;
252  }
253 
254 protected:
255 
256  std::vector<G4double> _eDep; // Deposited energy
257  std::vector<G4double> _timeBins; // Time bins
258  G4double _time; // Time of first energy release
259  G4double _pathLength; // Path length of the Hit
260  const G4VPhysicalVolume *_volume;
261  G4ThreeVector _absTranslation;
262  G4int _id;
263 
264  bool _isReset;
265  GGSPartHitsCollection *_partHits;
266  bool _storePosHits;
267 
268 };
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
272 typedef G4THitsCollection<GGSIntHit> GGSIntHitsCollection;
273 
275 extern G4Allocator<GGSIntHit> *GGSIntHitAllocator;
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278 
279 inline void* GGSIntHit::operator new(size_t) {
280  void* aHit;
281  aHit = (void*) GGSIntHitAllocator->MallocSingle();
282  return aHit;
283 }
284 
285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
286 
287 inline void GGSIntHit::operator delete(void* aHit) {
288  GGSIntHitAllocator->FreeSingle((GGSIntHit*) aHit);
289 }
290 
291 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
292 
293 #include "utils/GGSFactory.h"
294 
297 
298 #define RegisterIntHit(className) \
299  GGSIntHit* className##Builder(const std::vector<G4double> &timeBins){ \
300  return new className(timeBins); \
301  }\
302  class className##Proxy{ \
303  public: \
304  className##Proxy(){ \
305  GGSIntHitFactory::GetInstance().RegisterBuilder(#className, className##Builder); \
306  } \
307  };\
308  className##Proxy proxyFor##className
309 #endif /* GGSINTHIT_H */
310 
G4double GetTotalEnergyDep() const
Energy release getter.
Definition: GGSIntHit.h:178
virtual void AddStep(const G4Step &step)
Adds a step to the particle hit.
Definition: GGSIntHit.cpp:104
const G4ThreeVector & GetAbsTranslation()
Getter of absolute position of sensitive element.
Definition: GGSIntHit.h:242
G4Allocator< GGSIntHit > * GGSIntHitAllocator
Alias for G4 template memory allocator for GGSIntHit.
Definition: GGSIntHit.cpp:7
~GGSIntHit()
Destructor.
Definition: GGSIntHit.cpp:20
GGSPartHitsCollection * GetPartHits()
Getter of container of particle hits.
Definition: GGSIntHit.h:234
const std::vector< G4double > & GetTimeBins() const
Time bins getter.
Definition: GGSIntHit.h:208
void SetAbsTranslation(const G4ThreeVector &pos)
Set the position in global coordinates of the sensitive element.
Definition: GGSIntHit.h:136
void SetPosHitsStorage(bool flag)
Turn on or off the storage of position hits.
Definition: GGSIntHit.cpp:125
const G4VPhysicalVolume * GetVolume()
Getter for volume.
Definition: GGSIntHit.h:224
G4double GetPathLength() const
Path length getter.
Definition: GGSIntHit.h:216
const GGSIntHit & operator=(const GGSIntHit &right)
Assignment operator.
Definition: GGSIntHit.cpp:47
Template factory class.
Definition: GGSFactory.h:30
G4int GetID()
Getter for volume ID.
Definition: GGSIntHit.h:250
void SetID(G4int id)
Setter for volume ID.
Definition: GGSIntHit.h:170
G4double GetTime() const
Time getter.
Definition: GGSIntHit.h:194
GGSIntHit(const std::vector< G4double > &timeBins)
Constructor.
Definition: GGSIntHit.cpp:12
void SetVolume(const G4VPhysicalVolume *volume)
Setter for volume.
Definition: GGSIntHit.h:99
G4THitsCollection< GGSIntHit > GGSIntHitsCollection
Alias for G4 template hits collection for GGSIntHit.
Definition: GGSIntHit.h:272
G4THitsCollection< GGSPartHit > GGSPartHitsCollection
Alias for G4 template hits collection for GGSPartHit.
Definition: GGSPartHit.h:259
Definition of GGS Integrated Hit.
Definition: GGSIntHit.h:31
int operator==(const GGSIntHit &right) const
Comparison operator.
Definition: GGSIntHit.cpp:70
const std::vector< G4double > & GetEnergyDep() const
Energy release getter.
Definition: GGSIntHit.h:186
void SetPartHitsStorage(bool flag)
Turn on or off the storage of particle hits.
Definition: GGSIntHit.cpp:144