![]() |
|
|||
File indexing completed on 2025-02-23 09:22:29
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 parallel/ThreadsafeScorers/include/G4TAtomicHitsCollection.hh 0027 /// \brief Definition of the G4TAtomicHitsCollection class 0028 // 0029 // 0030 // 0031 // 0032 /// This is an implementation of G4THitsCollection<T> where the underlying 0033 /// type is G4atomic<T>, not just T. A static assert is provided to 0034 /// ensure that T is fundamental. This class should be used in lieu 0035 /// of G4THitsCollection<T> when memory is a concern. Atomics are 0036 /// thread-safe and *generally* faster that mutexes (as long as the 0037 /// STL implementation is lock-free) but the synchronization does 0038 /// not come without a cost. If performance is the primary concern, 0039 /// use G4THitsCollection<T> in thread-local instances. 0040 // 0041 // 0042 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0043 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0044 0045 #ifndef G4TAtomicHitsCollection_h 0046 #define G4TAtomicHitsCollection_h 1 0047 0048 #include "G4Allocator.hh" 0049 #include "G4AutoLock.hh" 0050 #include "G4Threading.hh" 0051 #include "G4VHitsCollection.hh" 0052 #include "G4atomic.hh" 0053 #include "globals.hh" 0054 0055 #include <deque> 0056 #include <type_traits> 0057 0058 // class description: 0059 // 0060 // This is a template class of hits collection and parametrized by 0061 // The concrete class of G4VHit. This is a uniform collection for 0062 // a particular concrete hit class objects. 0063 // An intermediate layer class G4HitsCollection appeared in this 0064 // header file is used just for G4Allocator, because G4Allocator 0065 // cannot be instansiated with a template class. Thus G4HitsCollection 0066 // class MUST NOT be directly used by the user. 0067 0068 /*class G4HitsCollection : public G4VHitsCollection 0069 { 0070 public: 0071 G4HitsCollection(); 0072 G4HitsCollection(G4String detName,G4String colNam); 0073 virtual ~G4HitsCollection(); 0074 G4bool operator==(const G4HitsCollection &right) const; 0075 0076 protected: 0077 void* theCollection; 0078 };*/ 0079 0080 template<class T> 0081 class G4TAtomicHitsCollection : public G4VHitsCollection 0082 { 0083 protected: 0084 static_assert(std::is_fundamental<T>::value, 0085 "G4TAtomicHitsCollection must use fundamental type"); 0086 0087 public: 0088 typedef T base_type; 0089 typedef G4atomic<T> value_type; 0090 typedef typename std::deque<value_type*> container_type; 0091 0092 public: 0093 G4TAtomicHitsCollection(); 0094 0095 public: 0096 // with description 0097 G4TAtomicHitsCollection(G4String detName, G4String colNam); 0098 0099 // constructor. 0100 public: 0101 virtual ~G4TAtomicHitsCollection(); 0102 G4bool operator==(const G4TAtomicHitsCollection<T>& right) const; 0103 0104 // inline void *operator new(size_t); 0105 // inline void operator delete(void* anHC); 0106 0107 public: // with description 0108 virtual void DrawAllHits(); 0109 virtual void PrintAllHits(); 0110 // These two methods invokes Draw() and Print() methods of all of 0111 // hit objects stored in this collection, respectively. 0112 0113 public: // with description 0114 inline value_type* operator[](size_t i) const { return (*theCollection)[i]; } 0115 // Returns a pointer to a concrete hit object. 0116 inline container_type* GetVector() const { return theCollection; } 0117 // Returns a collection vector. 0118 inline G4int insert(T* aHit) 0119 { 0120 G4AutoLock l(&fMutex); 0121 theCollection->push_back(aHit); 0122 return theCollection->size(); 0123 } 0124 // Insert a hit object. Total number of hit objects stored in this 0125 // collection is returned. 0126 inline G4int entries() const 0127 { 0128 G4AutoLock l(&fMutex); 0129 return theCollection->size(); 0130 } 0131 // Returns the number of hit objects stored in this collection 0132 0133 public: 0134 virtual G4VHit* GetHit(size_t i) const { return (*theCollection)[i]; } 0135 virtual size_t GetSize() const 0136 { 0137 G4AutoLock l(&fMutex); 0138 return theCollection->size(); 0139 } 0140 0141 protected: 0142 container_type* theCollection; 0143 G4Mutex fMutex; 0144 }; 0145 0146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0147 template<class T> 0148 G4TAtomicHitsCollection<T>::G4TAtomicHitsCollection() : theCollection(new container_type) 0149 {} 0150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0151 template<class T> 0152 G4TAtomicHitsCollection<T>::G4TAtomicHitsCollection(G4String detName, G4String colNam) 0153 : G4VHitsCollection(detName, colNam), theCollection(new container_type) 0154 {} 0155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0156 template<class T> 0157 G4TAtomicHitsCollection<T>::~G4TAtomicHitsCollection() 0158 { 0159 for (size_t i = 0; i < theCollection->size(); i++) 0160 delete (*theCollection)[i]; 0161 theCollection->clear(); 0162 delete theCollection; 0163 } 0164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0165 template<class T> 0166 G4bool G4TAtomicHitsCollection<T>::operator==(const G4TAtomicHitsCollection<T>& right) const 0167 { 0168 return (collectionName == right.collectionName); 0169 } 0170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0171 template<class T> 0172 void G4TAtomicHitsCollection<T>::DrawAllHits() 0173 { 0174 G4AutoLock l(&fMutex); 0175 for (size_t i = 0; i < theCollection->size(); i++) 0176 (*theCollection)[i]->Draw(); 0177 } 0178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0179 template<class T> 0180 void G4TAtomicHitsCollection<T>::PrintAllHits() 0181 { 0182 G4AutoLock l(&fMutex); 0183 for (size_t i = 0; i < theCollection->size(); i++) 0184 (*theCollection)[i]->Print(); 0185 } 0186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0187 0188 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |