Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4TAtomicHitsCollection< T > Class Template Reference

#include <G4TAtomicHitsCollection.hh>

Inheritance diagram for G4TAtomicHitsCollection< T >:
Collaboration diagram for G4TAtomicHitsCollection< T >:

Public Types

typedef T base_type
 
typedef G4atomic< T > value_type
 
typedef std::deque< value_type * > container_type
 

Public Member Functions

 G4TAtomicHitsCollection ()
 
 G4TAtomicHitsCollection (G4String detName, G4String colNam)
 
virtual ~G4TAtomicHitsCollection ()
 
G4int operator== (const G4TAtomicHitsCollection< T > &right) const
 
virtual void DrawAllHits ()
 
virtual void PrintAllHits ()
 
value_typeoperator[] (size_t i) const
 
container_typeGetVector () const
 
G4int insert (T *aHit)
 
G4int entries () const
 
virtual G4VHitGetHit (size_t i) const
 
virtual size_t GetSize () const
 
- Public Member Functions inherited from G4VHitsCollection
 G4VHitsCollection ()
 
 G4VHitsCollection (G4String detName, G4String colNam)
 
virtual ~G4VHitsCollection ()
 
G4int operator== (const G4VHitsCollection &right) const
 
G4String GetName ()
 
G4String GetSDname ()
 
 G4VHitsCollection ()
 
 G4VHitsCollection (G4String detName, G4String colNam)
 
virtual ~G4VHitsCollection ()
 
G4int operator== (const G4VHitsCollection &right) const
 
G4String GetName ()
 
G4String GetSDname ()
 
 G4VHitsCollection ()
 
 G4VHitsCollection (G4String detName, G4String colNam)
 
virtual ~G4VHitsCollection ()
 
G4int operator== (const G4VHitsCollection &right) const
 
G4StringGetName ()
 
G4StringGetSDname ()
 
void SetColID (G4int i)
 
G4int GetColID () const
 

Protected Attributes

container_typetheCollection
 
G4Mutex fMutex
 
- Protected Attributes inherited from G4VHitsCollection
G4String collectionName
 
G4String SDname
 
G4int colID
 

Detailed Description

template<class T>
class G4TAtomicHitsCollection< T >

This is an implementation of G4THitsCollection<T> where the underlying type is G4atomic<T>, not just T. A static assert is provided to ensure that T is fundamental. This class should be used in lieu of G4THitsCollection<T> when memory is a concern. Atomics are thread-safe and generally faster that mutexes (as long as the STL implementation is lock-free) but the synchronization does not come without a cost. If performance is the primary concern, use G4THitsCollection<T> in thread-local instances.

Definition at line 86 of file G4TAtomicHitsCollection.hh.

Member Typedef Documentation

template<class T>
typedef T G4TAtomicHitsCollection< T >::base_type

Definition at line 90 of file G4TAtomicHitsCollection.hh.

template<class T>
typedef std::deque<value_type*> G4TAtomicHitsCollection< T >::container_type

Definition at line 95 of file G4TAtomicHitsCollection.hh.

template<class T>
typedef G4atomic<T> G4TAtomicHitsCollection< T >::value_type

Definition at line 94 of file G4TAtomicHitsCollection.hh.

Constructor & Destructor Documentation

template<class T >
G4TAtomicHitsCollection< T >::G4TAtomicHitsCollection ( )

Definition at line 163 of file G4TAtomicHitsCollection.hh.

165 { }
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
std::deque< value_type * > container_type
template<class T >
G4TAtomicHitsCollection< T >::G4TAtomicHitsCollection ( G4String  detName,
G4String  colNam 
)

Definition at line 168 of file G4TAtomicHitsCollection.hh.

170  : G4VHitsCollection(detName,colNam),
173 { }
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
std::deque< value_type * > container_type
template<class T >
G4TAtomicHitsCollection< T >::~G4TAtomicHitsCollection ( )
virtual

Definition at line 175 of file G4TAtomicHitsCollection.hh.

176 {
177  for(size_t i = 0; i < theCollection->size(); i++)
178  delete (*theCollection)[i];
179  theCollection->clear();
180  delete theCollection;
182 }
#define G4MUTEXDESTROY(mutex)
Definition: G4Threading.hh:178

Member Function Documentation

template<class T >
void G4TAtomicHitsCollection< T >::DrawAllHits ( )
virtual

Reimplemented from G4VHitsCollection.

Definition at line 192 of file G4TAtomicHitsCollection.hh.

193 {
194  G4AutoLock l(&fMutex);
195  for(size_t i = 0; i < theCollection->size(); i++)
196  (*theCollection)[i]->Draw();
197 }
template<class T>
G4int G4TAtomicHitsCollection< T >::entries ( ) const
inline

Definition at line 137 of file G4TAtomicHitsCollection.hh.

138  {
139  G4AutoLock l(&fMutex);
140  return theCollection->size();
141  }
template<class T>
virtual G4VHit* G4TAtomicHitsCollection< T >::GetHit ( size_t  i) const
inlinevirtual

Reimplemented from G4VHitsCollection.

Definition at line 145 of file G4TAtomicHitsCollection.hh.

146  {
147  return (*theCollection)[i];
148  }
template<class T>
virtual size_t G4TAtomicHitsCollection< T >::GetSize ( ) const
inlinevirtual

Reimplemented from G4VHitsCollection.

Definition at line 149 of file G4TAtomicHitsCollection.hh.

150  {
151  G4AutoLock l(&fMutex);
152  return theCollection->size();
153  }
template<class T>
container_type* G4TAtomicHitsCollection< T >::GetVector ( ) const
inline

Definition at line 124 of file G4TAtomicHitsCollection.hh.

125  {
126  return theCollection;
127  }
template<class T>
G4int G4TAtomicHitsCollection< T >::insert ( T *  aHit)
inline

Definition at line 129 of file G4TAtomicHitsCollection.hh.

130  {
131  G4AutoLock l(&fMutex);
132  theCollection->push_back(aHit);
133  return theCollection->size();
134  }
template<class T >
G4int G4TAtomicHitsCollection< T >::operator== ( const G4TAtomicHitsCollection< T > &  right) const

Definition at line 186 of file G4TAtomicHitsCollection.hh.

187 {
188  return (collectionName == right.collectionName);
189 }
template<class T>
value_type* G4TAtomicHitsCollection< T >::operator[] ( size_t  i) const
inline

Definition at line 119 of file G4TAtomicHitsCollection.hh.

120  {
121  return (*theCollection)[i];
122  }
template<class T >
void G4TAtomicHitsCollection< T >::PrintAllHits ( )
virtual

Reimplemented from G4VHitsCollection.

Definition at line 200 of file G4TAtomicHitsCollection.hh.

201 {
202  G4AutoLock l(&fMutex);
203  for(size_t i = 0; i < theCollection->size(); i++)
204  (*theCollection)[i]->Print();
205 }

Member Data Documentation

template<class T>
G4Mutex G4TAtomicHitsCollection< T >::fMutex
protected

Definition at line 157 of file G4TAtomicHitsCollection.hh.

template<class T>
container_type* G4TAtomicHitsCollection< T >::theCollection
protected

Definition at line 156 of file G4TAtomicHitsCollection.hh.


The documentation for this class was generated from the following file: