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 // G4THitsCollection & G4HitsCollection 27 // 28 // Class description: 29 // 30 // This is a template class of hits collection and parametrised by 31 // The concrete class of G4VHit. This is a uniform collection for 32 // a particular concrete hit class objects. 33 // An intermediate layer class G4HitsCollection appeared in this 34 // header file is used just for G4Allocator, because G4Allocator 35 // cannot be instantiated with a template class. Thus G4HitsCollection 36 // class MUST NOT be directly used by the user. 37 // 38 // Author: Makoto Asai 39 // -------------------------------------------------------------------- 40 #ifndef G4THitsCollection_h 41 #define G4THitsCollection_h 1 42 43 #include "G4Allocator.hh" 44 #include "G4VHitsCollection.hh" 45 #include "globals.hh" 46 47 #include <vector> 48 49 class G4HitsCollection : public G4VHitsCollection 50 { 51 public: 52 53 using G4VHitsCollection::G4VHitsCollection; 54 ~G4HitsCollection() override = default; 55 G4bool operator==(const G4HitsCollection& right) const; 56 57 protected: 58 59 void* theCollection = nullptr; 60 }; 61 62 #if defined G4DIGI_ALLOC_EXPORT 63 extern G4DLLEXPORT G4Allocator<G4HitsCollection>*& anHCAllocator_G4MT_TLS_(); 64 #else 65 extern G4DLLIMPORT G4Allocator<G4HitsCollection>*& anHCAllocator_G4MT_TLS_(); 66 #endif 67 68 template <class T> 69 class G4THitsCollection : public G4HitsCollection 70 { 71 public: 72 73 G4THitsCollection(); 74 G4THitsCollection(const G4String& detName, const G4String& colNam); 75 ~G4THitsCollection() override; 76 77 G4bool operator==(const G4THitsCollection<T>& right) const; 78 79 inline void* operator new(std::size_t); 80 inline void operator delete(void* anHC); 81 82 // Invoke Draw() method on all hit objects in collection 83 void DrawAllHits() override; 84 85 // Invoke Print() method on all hit objects in collection 86 void PrintAllHits() override; 87 88 // Returns a pointer to the concrete hit object 89 // Returns pointer to the concrete hit object at given index 90 // Not bounds checked 91 inline T* operator[](std::size_t i) const { return (*((std::vector<T*>*)theCollection))[i]; } 92 93 // Return pointer to hit collection 94 inline std::vector<T*>* GetVector() const { return (std::vector<T*>*)theCollection; } 95 96 // Insert a hit object in the collection, taking onwenership 97 // Returns the total number of hit objects stored after insertion 98 inline std::size_t insert(T* aHit) 99 { 100 auto theHitsCollection = (std::vector<T*>*)theCollection; 101 theHitsCollection->push_back(aHit); 102 return theHitsCollection->size(); 103 } 104 105 // Returns the number of hit objects stored in this collection. 106 inline std::size_t entries() const 107 { 108 auto theHitsCollection = (std::vector<T*>*)theCollection; 109 return theHitsCollection->size(); 110 } 111 112 G4VHit* GetHit(std::size_t i) const override { return (*((std::vector<T*>*)theCollection))[i]; } 113 114 std::size_t GetSize() const override { return ((std::vector<T*>*)theCollection)->size(); } 115 }; 116 117 template <class T> 118 inline void* G4THitsCollection<T>::operator new(std::size_t) 119 { 120 if (anHCAllocator_G4MT_TLS_() == nullptr) { 121 anHCAllocator_G4MT_TLS_() = new G4Allocator<G4HitsCollection>; 122 } 123 return (void*)anHCAllocator_G4MT_TLS_()->MallocSingle(); 124 } 125 126 template <class T> 127 inline void G4THitsCollection<T>::operator delete(void* anHC) 128 { 129 anHCAllocator_G4MT_TLS_()->FreeSingle((G4HitsCollection*)anHC); 130 } 131 132 template <class T> 133 G4THitsCollection<T>::G4THitsCollection() 134 { 135 auto theHitsCollection = new std::vector<T*>; 136 theCollection = (void*)theHitsCollection; 137 } 138 139 template <class T> 140 G4THitsCollection<T>::G4THitsCollection(const G4String& detName, const G4String& colNam) 141 : G4HitsCollection(detName, colNam) 142 { 143 auto theHitsCollection = new std::vector<T*>; 144 theCollection = (void*)theHitsCollection; 145 } 146 147 template <class T> 148 G4THitsCollection<T>::~G4THitsCollection() 149 { 150 auto theHitsCollection = (std::vector<T*>*)theCollection; 151 for (const auto* hit : *theHitsCollection) { 152 delete hit; 153 } 154 theHitsCollection->clear(); 155 delete theHitsCollection; 156 } 157 158 template <class T> 159 G4bool G4THitsCollection<T>::operator==(const G4THitsCollection<T>& right) const 160 { 161 return (collectionName == right.collectionName); 162 } 163 164 template <class T> 165 void G4THitsCollection<T>::DrawAllHits() 166 { 167 auto theHitsCollection = (std::vector<T*>*)theCollection; 168 for (auto* hit : *theHitsCollection) { 169 hit->Draw(); 170 } 171 } 172 173 template <class T> 174 void G4THitsCollection<T>::PrintAllHits() 175 { 176 auto theHitsCollection = (std::vector<T*>*)theCollection; 177 for (auto* hit : *theHitsCollection) { 178 hit->Print(); 179 } 180 } 181 182 #endif 183