|
|
|||
File indexing completed on 2026-04-16 07:41:59
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 /// \file G4TAtomicHitsCollection.hh 0027 /// \brief Definition of the G4TAtomicHitsCollection class 0028 /// 0029 /// This is an implementation of G4THitsCollection<T> where the underlying 0030 /// type is G4atomic<T>, not just T. A static assert is provided to 0031 /// ensure that T is fundamental. This class should be used in lieu 0032 /// of G4THitsCollection<T> when memory is a concern. Atomics are 0033 /// thread-safe and *generally* faster that mutexes (as long as the 0034 /// STL implementation is lock-free) but the synchronization does 0035 /// not come without a cost. If performance is the primary concern, 0036 /// use G4THitsCollection<T> in thread-local instances. 0037 // 0038 // 0039 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0040 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0041 0042 #ifndef G4TAtomicHitsCollection_h 0043 #define G4TAtomicHitsCollection_h 1 0044 0045 #include "G4Allocator.hh" 0046 #include "G4AutoLock.hh" 0047 #include "G4Threading.hh" 0048 #include "G4VHitsCollection.hh" 0049 #include "G4atomic.hh" 0050 #include "globals.hh" 0051 0052 #include <deque> 0053 #include <type_traits> 0054 0055 /// class description: 0056 /// 0057 /// This is a template class of hits collection and parametrized by 0058 /// The concrete class of G4VHit. This is a uniform collection for 0059 /// a particular concrete hit class objects. 0060 /// An intermediate layer class G4HitsCollection appeared in this 0061 /// header file is used just for G4Allocator, because G4Allocator 0062 /// cannot be instansiated with a template class. Thus G4HitsCollection 0063 /// class MUST NOT be directly used by the user. 0064 0065 /*class G4HitsCollection : public G4VHitsCollection 0066 { 0067 public: 0068 G4HitsCollection(); 0069 G4HitsCollection(G4String detName,G4String colNam); 0070 virtual ~G4HitsCollection(); 0071 G4bool operator==(const G4HitsCollection &right) const; 0072 0073 protected: 0074 void* theCollection; 0075 };*/ 0076 0077 template<class T> 0078 class G4TAtomicHitsCollection : public G4VHitsCollection 0079 { 0080 protected: 0081 static_assert(std::is_fundamental<T>::value, 0082 "G4TAtomicHitsCollection must use fundamental type"); 0083 0084 public: 0085 typedef T base_type; 0086 typedef G4atomic<T> value_type; 0087 typedef typename std::deque<value_type*> container_type; 0088 0089 public: 0090 G4TAtomicHitsCollection(); 0091 0092 public: 0093 // with description 0094 G4TAtomicHitsCollection(G4String detName, G4String colNam); 0095 0096 // constructor. 0097 public: 0098 virtual ~G4TAtomicHitsCollection(); 0099 G4bool operator==(const G4TAtomicHitsCollection<T>& right) const; 0100 0101 // inline void *operator new(size_t); 0102 // inline void operator delete(void* anHC); 0103 0104 public: // with description 0105 virtual void DrawAllHits(); 0106 virtual void PrintAllHits(); 0107 // These two methods invokes Draw() and Print() methods of all of 0108 // hit objects stored in this collection, respectively. 0109 0110 public: // with description 0111 inline value_type* operator[](size_t i) const { return (*theCollection)[i]; } 0112 // Returns a pointer to a concrete hit object. 0113 inline container_type* GetVector() const { return theCollection; } 0114 // Returns a collection vector. 0115 inline G4int insert(T* aHit) 0116 { 0117 G4AutoLock l(&fMutex); 0118 theCollection->push_back(aHit); 0119 return theCollection->size(); 0120 } 0121 // Insert a hit object. Total number of hit objects stored in this 0122 // collection is returned. 0123 inline G4int entries() const 0124 { 0125 G4AutoLock l(&fMutex); 0126 return theCollection->size(); 0127 } 0128 // Returns the number of hit objects stored in this collection 0129 0130 public: 0131 virtual G4VHit* GetHit(size_t i) const { return (*theCollection)[i]; } 0132 virtual size_t GetSize() const 0133 { 0134 G4AutoLock l(&fMutex); 0135 return theCollection->size(); 0136 } 0137 0138 protected: 0139 container_type* theCollection; 0140 G4Mutex fMutex; 0141 }; 0142 0143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0144 template<class T> 0145 G4TAtomicHitsCollection<T>::G4TAtomicHitsCollection() : theCollection(new container_type) 0146 {} 0147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0148 template<class T> 0149 G4TAtomicHitsCollection<T>::G4TAtomicHitsCollection(G4String detName, G4String colNam) 0150 : G4VHitsCollection(detName, colNam), theCollection(new container_type) 0151 {} 0152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0153 template<class T> 0154 G4TAtomicHitsCollection<T>::~G4TAtomicHitsCollection() 0155 { 0156 for (size_t i = 0; i < theCollection->size(); i++) 0157 delete (*theCollection)[i]; 0158 theCollection->clear(); 0159 delete theCollection; 0160 } 0161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0162 template<class T> 0163 G4bool G4TAtomicHitsCollection<T>::operator==(const G4TAtomicHitsCollection<T>& right) const 0164 { 0165 return (collectionName == right.collectionName); 0166 } 0167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0168 template<class T> 0169 void G4TAtomicHitsCollection<T>::DrawAllHits() 0170 { 0171 G4AutoLock l(&fMutex); 0172 for (size_t i = 0; i < theCollection->size(); i++) 0173 (*theCollection)[i]->Draw(); 0174 } 0175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0176 template<class T> 0177 void G4TAtomicHitsCollection<T>::PrintAllHits() 0178 { 0179 G4AutoLock l(&fMutex); 0180 for (size_t i = 0; i < theCollection->size(); i++) 0181 (*theCollection)[i]->Print(); 0182 } 0183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0184 0185 #endif
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|