Loading...
Searching...
No Matches
G4TAtomicHitsCollection.hh
Go to the documentation of this file.
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 parallel/ThreadsafeScorers/include/G4TAtomicHitsCollection.hh
27/// \brief Definition of the G4TAtomicHitsCollection class
28//
29//
30//
31//
32/// This is an implementation of G4THitsCollection<T> where the underlying
33/// type is G4atomic<T>, not just T. A static assert is provided to
34/// ensure that T is fundamental. This class should be used in lieu
35/// of G4THitsCollection<T> when memory is a concern. Atomics are
36/// thread-safe and *generally* faster that mutexes (as long as the
37/// STL implementation is lock-free) but the synchronization does
38/// not come without a cost. If performance is the primary concern,
39/// use G4THitsCollection<T> in thread-local instances.
40//
41//
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
45#ifndef G4TAtomicHitsCollection_h
46#define G4TAtomicHitsCollection_h 1
47
48#include "G4VHitsCollection.hh"
49#include "G4Allocator.hh"
50#include "globals.hh"
51#include "G4Threading.hh"
52#include "G4AutoLock.hh"
53#include "G4atomic.hh"
54
55#include <deque>
56#include <type_traits>
57
58// class description:
59//
60// This is a template class of hits collection and parametrized by
61// The concrete class of G4VHit. This is a uniform collection for
62// a particular concrete hit class objects.
63// An intermediate layer class G4HitsCollection appeared in this
64// header file is used just for G4Allocator, because G4Allocator
65// cannot be instansiated with a template class. Thus G4HitsCollection
66// class MUST NOT be directly used by the user.
67
68/*class G4HitsCollection : public G4VHitsCollection
69{
70 public:
71 G4HitsCollection();
72 G4HitsCollection(G4String detName,G4String colNam);
73 virtual ~G4HitsCollection();
74 G4bool operator==(const G4HitsCollection &right) const;
75
76 protected:
77 void* theCollection;
78};*/
79
80template <class T>
82{
83 protected:
84 static_assert(std::is_fundamental<T>::value,
85 "G4TAtomicHitsCollection must use fundamental type");
86
87 public:
88 typedef T base_type;
90 typedef typename std::deque<value_type*> container_type;
91
92 public:
94
95 public:
96 // with description
98
99 // constructor.
100 public:
101 virtual ~G4TAtomicHitsCollection();
102 G4bool operator==(const G4TAtomicHitsCollection<T>& right) const;
103
104 // inline void *operator new(size_t);
105 // inline void operator delete(void* anHC);
106
107 public: // with description
108 virtual void DrawAllHits();
109 virtual void PrintAllHits();
110 // These two methods invokes Draw() and Print() methods of all of
111 // hit objects stored in this collection, respectively.
112
113 public: // with description
114 inline value_type* operator[](size_t i) const { return (*theCollection)[i]; }
115 // Returns a pointer to a concrete hit object.
116 inline container_type* GetVector() const { return theCollection; }
117 // Returns a collection vector.
118 inline G4int insert(T* aHit)
119 {
120 G4AutoLock l(&fMutex);
121 theCollection->push_back(aHit);
122 return theCollection->size();
123 }
124 // Insert a hit object. Total number of hit objects stored in this
125 // collection is returned.
126 inline G4int entries() const
127 {
128 G4AutoLock l(&fMutex);
129 return theCollection->size();
130 }
131 // Returns the number of hit objects stored in this collection
132
133 public:
134 virtual G4VHit* GetHit(size_t i) const { return (*theCollection)[i]; }
135 virtual size_t GetSize() const
136 {
137 G4AutoLock l(&fMutex);
138 return theCollection->size();
139 }
140
141 protected:
143 G4Mutex fMutex;
144};
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147template <class T>
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152template <class T>
154 G4String colNam)
155 : G4VHitsCollection(detName, colNam)
156 , theCollection(new container_type)
157{}
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159template <class T>
161{
162 for(size_t i = 0; i < theCollection->size(); i++)
163 delete(*theCollection)[i];
164 theCollection->clear();
165 delete theCollection;
166}
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
168template <class T>
171{
172 return (collectionName == right.collectionName);
173}
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175template <class T>
177{
178 G4AutoLock l(&fMutex);
179 for(size_t i = 0; i < theCollection->size(); i++)
180 (*theCollection)[i]->Draw();
181}
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183template <class T>
185{
186 G4AutoLock l(&fMutex);
187 for(size_t i = 0; i < theCollection->size(); i++)
188 (*theCollection)[i]->Print();
189}
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191
192#endif
Definition of the G4atomic class.
This is an implementation of G4THitsCollection<T> where the underlying type is G4atomic<T>,...
value_type * operator[](size_t i) const
G4bool operator==(const G4TAtomicHitsCollection< T > &right) const
container_type * GetVector() const
virtual G4VHit * GetHit(size_t i) const
std::deque< value_type * > container_type

Applications | User Support | Publications | Collaboration