Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
source/digits_hits/hits/include/G4THitsCollection.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 //
27 // $Id: G4THitsCollection.hh 67992 2013-03-13 10:59:57Z gcosmo $
28 //
29 
30 #ifndef G4THitsCollection_h
31 #define G4THitsCollection_h 1
32 
33 #include "G4VHitsCollection.hh"
34 #include "G4Allocator.hh"
35 #include "globals.hh"
36 //#include "g4rw/tpordvec.h"
37 #include <vector>
38 
39 // class description:
40 //
41 // This is a template class of hits collection and parametrized by
42 // The concrete class of G4VHit. This is a uniform collection for
43 // a particular concrete hit class objects.
44 // An intermediate layer class G4HitsCollection appeared in this
45 // header file is used just for G4Allocator, because G4Allocator
46 // cannot be instansiated with a template class. Thus G4HitsCollection
47 // class MUST NOT be directly used by the user.
48 
50 {
51  public:
53  G4HitsCollection(G4String detName,G4String colNam);
54  virtual ~G4HitsCollection();
55  G4int operator==(const G4HitsCollection &right) const;
56 
57  protected:
58  void* theCollection;
59 };
60 
61 #if defined G4DIGI_ALLOC_EXPORT
63 #else
65 #endif
66 
67 template <class T> class G4THitsCollection : public G4HitsCollection
68 {
69  public:
71  public: // with description
72  G4THitsCollection(G4String detName,G4String colNam);
73  // constructor.
74  public:
75  virtual ~G4THitsCollection();
77 
78  inline void *operator new(size_t);
79  inline void operator delete(void* anHC);
80  public: // with description
81  virtual void DrawAllHits();
82  virtual void PrintAllHits();
83  // These two methods invokes Draw() and Print() methods of all of
84  // hit objects stored in this collection, respectively.
85 
86  public: // with description
87  inline T* operator[](size_t i) const
88  {
90  return (*((std::vector<T*>*)theCollection))[i];
91  }
92  // Returns a pointer to a concrete hit object.
93  inline std::vector<T*>* GetVector() const
94  {
96  return (std::vector<T*>*)theCollection; }
97  // Returns a collection vector.
98  inline G4int insert(T* aHit)
99  {
101  std::vector<T*>*theHitsCollection = (std::vector<T*>*)theCollection;
102  theHitsCollection->push_back(aHit);
103  return theHitsCollection->size();
104  }
105  // Insert a hit object. Total number of hit objects stored in this
106  // collection is returned.
107  inline G4int entries() const
108  {
110  std::vector<T*>*theHitsCollection = (std::vector<T*>*)theCollection;
111  return theHitsCollection->size();
112  }
113  // Returns the number of hit objects stored in this collection
114 
115  public:
116  virtual G4VHit* GetHit(size_t i) const
117  {
119  return (*((std::vector<T*>*)theCollection))[i];
120  }
121  virtual size_t GetSize() const
122  {
124  return ((std::vector<T*>*)theCollection)->size(); }
125 
126 };
127 
128 template <class T> inline void* G4THitsCollection<T>::operator new(size_t)
129 {
132  void* anHC;
133  anHC = (void*)anHCAllocator.MallocSingle();
134  return anHC;
135 }
136 
137 template <class T> inline void G4THitsCollection<T>::operator delete(void* anHC)
138 {
141  anHCAllocator.FreeSingle((G4HitsCollection*)anHC);
142 }
143 
144 template <class T> G4THitsCollection<T>::G4THitsCollection()
145 {
147  std::vector<T*> * theHitsCollection = new std::vector<T*>;
148  theCollection = (void*)theHitsCollection;
149 }
150 
151 template <class T> G4THitsCollection<T>::G4THitsCollection(G4String detName,G4String colNam)
152 : G4HitsCollection(detName,colNam)
153 {
155  std::vector<T*> * theHitsCollection = new std::vector<T*>;
156  theCollection = (void*)theHitsCollection;
157 }
158 
159 template <class T> G4THitsCollection<T>::~G4THitsCollection()
160 {
162  std::vector<T*> * theHitsCollection = (std::vector<T*>*)theCollection;
163  //theHitsCollection->clearAndDestroy();
164  for(size_t i=0;i<theHitsCollection->size();i++)
165  { delete (*theHitsCollection)[i]; }
166  theHitsCollection->clear();
167  delete theHitsCollection;
168 }
169 
170 template <class T> G4int G4THitsCollection<T>::operator==(const G4THitsCollection<T> &right) const
171 {
173  return (collectionName==right.collectionName);
174 }
175 
176 template <class T> void G4THitsCollection<T>::DrawAllHits()
177 {
179  std::vector<T*> * theHitsCollection = (std::vector<T*>*)theCollection;
180  size_t n = theHitsCollection->size();
181  for(size_t i=0;i<n;i++)
182  { (*theHitsCollection)[i]->Draw(); }
183 }
184 
185 template <class T> void G4THitsCollection<T>::PrintAllHits()
186 {
188  std::vector<T*> * theHitsCollection
189  = (std::vector<T*>*)theCollection;
190  size_t n = theHitsCollection->size();
191  for(size_t i=0;i<n;i++)
192  { (*theHitsCollection)[i]->Print(); }
193 }
194 
195 #endif
196 
Type * MallocSingle()
Definition: G4Allocator.hh:191
#define G4DLLEXPORT
Definition: G4Types.hh:62
G4int operator==(const G4THitsCollection< T > &right) const
#define G4DLLIMPORT
Definition: G4Types.hh:63
void FreeSingle(Type *anElement)
Definition: G4Allocator.hh:201
virtual ~G4HitsCollection()
Definition: G4VHit.hh:48
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
const G4int n
G4DLLIMPORT G4ThreadLocal G4Allocator< G4HitsCollection > * anHCAllocator_G4MT_TLS_
G4int operator==(const G4HitsCollection &right) const
G4DLLIMPORT G4Allocator< G4HitsCollection > anHCAllocator