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 SAXSSensitiveDetectorHit.hh 27 /// \brief Definition of the SAXSSensitiveDetectorHit class 28 // 29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 30 31 #ifndef SAXSSensitiveDetectorHit_h 32 #define SAXSSensitiveDetectorHit_h 1 33 34 #include "G4Allocator.hh" 35 #include "G4LogicalVolume.hh" 36 #include "G4RotationMatrix.hh" 37 #include "G4THitsCollection.hh" 38 #include "G4ThreeVector.hh" 39 #include "G4Transform3D.hh" 40 #include "G4VHit.hh" 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 44 /// Sensitive Detector hit 45 /// 46 /// It records: 47 /// - the spatial coordinate of the hit 48 /// - the time at which the hit occurred 49 /// - the type, ID, weight, momentum and energy of the impinging particle 50 /// - the ID of the impinging particle parent 51 52 class SAXSSensitiveDetectorHit : public G4VHit 53 { 54 public: 55 SAXSSensitiveDetectorHit(); 56 SAXSSensitiveDetectorHit(const SAXSSensitiveDetectorHit& right); 57 virtual ~SAXSSensitiveDetectorHit() { ; } 58 59 const SAXSSensitiveDetectorHit& operator=(const SAXSSensitiveDetectorHit& right); 60 61 int operator==(const SAXSSensitiveDetectorHit& right) const; 62 63 inline void* operator new(size_t); 64 inline void operator delete(void* aHit); 65 66 virtual void Draw(); 67 virtual void Print(); 68 69 private: 70 G4int fTrackID; 71 G4int fTrackIDP; 72 G4double fTime; 73 G4ThreeVector fPos; 74 G4ThreeVector fMom; 75 G4double fEnergy; 76 G4int fType; 77 G4double fWeight; 78 79 public: 80 inline void SetTrackID(G4int z) { fTrackID = z; } 81 inline G4int GetTrackID() const { return fTrackID; } 82 83 inline void SetTrackIDP(G4int z) { fTrackIDP = z; } 84 inline G4int GetTrackIDP() const { return fTrackIDP; } 85 86 inline void SetTime(G4double t) { fTime = t; } 87 inline G4double GetTime() const { return fTime; } 88 89 inline void SetPos(G4ThreeVector xyz) { fPos = xyz; } 90 inline G4ThreeVector GetPos() const { return fPos; } 91 92 inline void SetMom(G4ThreeVector xyz) { fMom = xyz; } 93 inline G4ThreeVector GetMom() const { return fMom; } 94 95 inline void SetEnergy(G4double energy) { fEnergy = energy; } 96 inline G4double GetEnergy() const { return fEnergy; } 97 98 inline void SetType(G4int type) { fType = type; } 99 inline G4int GetType() const { return fType; } 100 101 inline void SetWeight(G4double weight) { fWeight = weight; } 102 inline G4double GetWeight() const { return fWeight; } 103 }; 104 105 using SensitiveDetectorHitsCollection = G4THitsCollection<SAXSSensitiveDetectorHit>; 106 107 extern G4ThreadLocal G4Allocator<SAXSSensitiveDetectorHit>* hitAllocator; 108 109 inline void* SAXSSensitiveDetectorHit::operator new(size_t) 110 { 111 if (!hitAllocator) { 112 hitAllocator = new G4Allocator<SAXSSensitiveDetectorHit>; 113 } 114 return hitAllocator->MallocSingle(); 115 } 116 117 inline void SAXSSensitiveDetectorHit::operator delete(void* aHit) 118 { 119 if (!hitAllocator) { 120 hitAllocator = new G4Allocator<SAXSSensitiveDetectorHit>; 121 } 122 hitAllocator->FreeSingle((SAXSSensitiveDetectorHit*)aHit); 123 } 124 125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 126 127 #endif 128