Monte Carlo LYSO
Geant4 simulation for the LYSO calorimeter prototype
hit.hh
Go to the documentation of this file.
1 
5 #ifndef HIT_HH
6 #define HIT_HH
7 
8 #include "G4VHit.hh"
9 #include "G4THitsCollection.hh"
10 #include "G4Allocator.hh"
11 #include "G4ThreeVector.hh"
12 #include "tls.hh"
13 #include "G4UnitsTable.hh"
14 
19 class MyHit : public G4VHit
20 {
21 public:
22  MyHit() = default;
23  ~MyHit() override = default;
25  // Operators
26  MyHit& operator=(const MyHit&) = default;
27  G4bool operator==(const MyHit&) const;
29  // Operator overloads for memory management
30  inline void* operator new(size_t);
31  inline void operator delete(void*);
33  // Set methods
34  inline void SetDetectionTime(G4double t) { fDetectionTime = t; }
35  inline void SetDetectorPosition(G4ThreeVector xyz) { fDetectorPosition = xyz; }
36  inline void SetDetectorChannel(G4int ch) { fCh = ch; }
38  // Get methods
39  inline G4double GetDetectionTime() const { return fDetectionTime; }
40  inline G4ThreeVector GetDetectorPosition() const { return fDetectorPosition; }
41  inline G4int GetDetectorChannel() const { return fCh; }
43 private:
44  G4double fDetectionTime;
45  G4ThreeVector fDetectorPosition;
46  G4int fCh;
47 };
48 
54 using MyHitsCollection = G4THitsCollection<MyHit>;
55 
60 extern G4ThreadLocal G4Allocator<MyHit>* MyHitAllocator;
61 
62 
63 // Implementation of the new operator for MyHit by leveraging MyHitAllocator
64 inline void* MyHit::operator new(size_t)
65 {
66  // Create a new G4Allocator if it does not exist
67  if (!MyHitAllocator)
68  MyHitAllocator = new G4Allocator<MyHit>;
69 
70  // Allocate memory for a single MyHit object
71  return (void*)MyHitAllocator->MallocSingle();
72 }
73 
74 // Implementation of the delete operator for MyHit by leveraging MyHitAllocator
75 inline void MyHit::operator delete(void* hit)
76 {
77  // Free the memory for a single MyHit object
78  MyHitAllocator->FreeSingle(static_cast<MyHit*>(hit));
79 }
80 
81 #endif // HIT_HH
Concrete class of G4VHit, representing a hit in the MySensitiveDetector.
Definition: hit.hh:20
G4int fCh
Channel of the SiPM hit by the optical photon.
Definition: hit.hh:46
G4bool operator==(const MyHit &) const
Equality operator, compares two MyHit objects for equality.
Definition: hit.cc:11
MyHit & operator=(const MyHit &)=default
Assignment operator, default implementation.
void SetDetectorChannel(G4int ch)
Set the channel of the SiPM hit by the optical photon.
Definition: hit.hh:36
G4double GetDetectionTime() const
Get the detection time of the optical photon.
Definition: hit.hh:39
G4int GetDetectorChannel() const
Get the channel of the SiPM hit by the optical photon.
Definition: hit.hh:41
G4double fDetectionTime
Time of detection of the optical photon.
Definition: hit.hh:44
G4ThreeVector fDetectorPosition
Position (center) of the SiPM hit by the optical photon.
Definition: hit.hh:45
G4ThreeVector GetDetectorPosition() const
Get the position (center) of the SiPM hit by the optical photon.
Definition: hit.hh:40
void SetDetectionTime(G4double t)
Set the detection time of the optical photon.
Definition: hit.hh:34
MyHit()=default
Constructor of the class.
void SetDetectorPosition(G4ThreeVector xyz)
Set the position (center) of the SiPM hit by the optical photon.
Definition: hit.hh:35
~MyHit() override=default
Destructor of the class.
G4ThreadLocal G4Allocator< MyHit > * MyHitAllocator
G4Allocator for MyHit objects.
Definition: hit.cc:8
G4THitsCollection< MyHit > MyHitsCollection
Concrete hit collection class for MyHit.
Definition: hit.hh:54