File indexing completed on 2025-01-31 09:22:12
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038 #ifndef GammaRayTelCalorimeterHit_h
0039 #define GammaRayTelCalorimeterHit_h 1
0040
0041 #include "G4Allocator.hh"
0042 #include "G4THitsCollection.hh"
0043 #include "G4ThreeVector.hh"
0044 #include "G4VHit.hh"
0045
0046 class GammaRayTelCalorimeterHit: public G4VHit {
0047 public:
0048 GammaRayTelCalorimeterHit();
0049
0050 ~GammaRayTelCalorimeterHit() override;
0051
0052 GammaRayTelCalorimeterHit(const GammaRayTelCalorimeterHit &right);
0053
0054 auto operator=(const GammaRayTelCalorimeterHit &right) -> const GammaRayTelCalorimeterHit&;
0055
0056 auto operator==(const GammaRayTelCalorimeterHit &right) const -> G4bool;
0057
0058 inline auto operator new(size_t) -> void*;
0059
0060 inline auto operator delete(void* hit) -> void;
0061
0062 void Draw() override;
0063
0064 void Print() override;
0065
0066 inline void AddEnergy(G4double value) {
0067 calDepositedEnergy += value;
0068 }
0069
0070 inline void SetCALBarNumber(const G4int &value) {
0071 calBarNumber = value;
0072 }
0073
0074 inline void SetCALPlaneNumber(const G4int &value) {
0075 calPlaneNumber = value;
0076 }
0077
0078 inline void SetCALType(const G4int &value) {
0079 isCALPlane = value;
0080 }
0081
0082 inline void SetPosition(const G4ThreeVector &vector) {
0083 position = vector;
0084 }
0085
0086 [[nodiscard]]
0087 inline auto GetCALDepositedEnergy() const -> G4double {
0088 return calDepositedEnergy;
0089 }
0090
0091 [[nodiscard]]
0092 inline auto GetCALBarNumber() const -> G4int {
0093 return calBarNumber;
0094 }
0095
0096 [[nodiscard]]
0097 inline auto GetCALPlaneNumber() const -> G4int {
0098 return calPlaneNumber;
0099 }
0100
0101 [[nodiscard]]
0102 inline auto GetCALType() const -> G4int {
0103 return isCALPlane;
0104 }
0105
0106 [[nodiscard]]
0107 inline auto GetPosition() const -> G4ThreeVector {
0108 return position;
0109 }
0110
0111 private:
0112 G4int calBarNumber{0};
0113
0114 G4int calPlaneNumber{0};
0115
0116 G4int isCALPlane{0};
0117
0118 G4double calDepositedEnergy{0.};
0119
0120 G4ThreeVector position{G4ThreeVector(0., 0., 0.)};
0121 };
0122
0123 using GammaRayTelCalorimeterHitsCollection = G4THitsCollection<GammaRayTelCalorimeterHit>;
0124 extern G4ThreadLocal G4Allocator<GammaRayTelCalorimeterHit> *calorimeterHitAllocator;
0125
0126 inline auto GammaRayTelCalorimeterHit::operator new(size_t) -> void* {
0127 if (calorimeterHitAllocator == nullptr) {
0128 calorimeterHitAllocator = new G4Allocator<GammaRayTelCalorimeterHit>;
0129 }
0130 return (void*) calorimeterHitAllocator->MallocSingle();
0131 }
0132
0133 inline auto GammaRayTelCalorimeterHit::operator delete(void* hit) -> void {
0134 calorimeterHitAllocator->FreeSingle((GammaRayTelCalorimeterHit*) hit);
0135 }
0136 #endif