source: trunk/source/digits_hits/hits/include/G4THitsCollection.hh@ 1202

Last change on this file since 1202 was 1012, checked in by garnier, 17 years ago

update

File size: 6.2 KB
Line 
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,v 1.5 2006/06/29 18:06:34 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30
31#ifndef G4THitsCollection_h
32#define G4THitsCollection_h 1
33
34#include "G4VHitsCollection.hh"
35#include "G4Allocator.hh"
36#include "globals.hh"
37//#include "g4rw/tpordvec.h"
38#include <vector>
39
40// class description:
41//
42// This is a template class of hits collection and parametrized by
43// The concrete class of G4VHit. This is a uniform collection for
44// a particular concrete hit class objects.
45// An intermediate layer class G4HitsCollection appeared in this
46// header file is used just for G4Allocator, because G4Allocator
47// cannot be instansiated with a template class. Thus G4HitsCollection
48// class MUST NOT be directly used by the user.
49
50class G4HitsCollection : public G4VHitsCollection
51{
52 public:
53 G4HitsCollection();
54 G4HitsCollection(G4String detName,G4String colNam);
55 virtual ~G4HitsCollection();
56 G4int operator==(const G4HitsCollection &right) const;
57
58 protected:
59 void* theCollection;
60};
61
62#if defined G4DIGI_ALLOC_EXPORT
63 extern G4DLLEXPORT G4Allocator<G4HitsCollection> anHCAllocator;
64#else
65 extern G4DLLIMPORT G4Allocator<G4HitsCollection> anHCAllocator;
66#endif
67
68template <class T> class G4THitsCollection : public G4HitsCollection
69{
70 public:
71 G4THitsCollection();
72 public: // with description
73 G4THitsCollection(G4String detName,G4String colNam);
74 // constructor.
75 public:
76 virtual ~G4THitsCollection();
77 G4int operator==(const G4THitsCollection<T> &right) const;
78
79 inline void *operator new(size_t);
80 inline void operator delete(void* anHC);
81 public: // with description
82 virtual void DrawAllHits();
83 virtual void PrintAllHits();
84 // These two methods invokes Draw() and Print() methods of all of
85 // hit objects stored in this collection, respectively.
86
87 public: // with description
88 inline T* operator[](size_t i) const
89 { return (*((std::vector<T*>*)theCollection))[i]; }
90 // Returns a pointer to a concrete hit object.
91 inline std::vector<T*>* GetVector() const
92 { return (std::vector<T*>*)theCollection; }
93 // Returns a collection vector.
94 inline G4int insert(T* aHit)
95 {
96 std::vector<T*>*theHitsCollection
97 = (std::vector<T*>*)theCollection;
98 theHitsCollection->push_back(aHit);
99 return theHitsCollection->size();
100 }
101 // Insert a hit object. Total number of hit objects stored in this
102 // collection is returned.
103 inline G4int entries() const
104 {
105 std::vector<T*>*theHitsCollection
106 = (std::vector<T*>*)theCollection;
107 return theHitsCollection->size();
108 }
109 // Returns the number of hit objects stored in this collection
110
111 public:
112 virtual G4VHit* GetHit(size_t i) const
113 { return (*((std::vector<T*>*)theCollection))[i]; }
114 virtual size_t GetSize() const
115 { return ((std::vector<T*>*)theCollection)->size(); }
116
117};
118
119template <class T> inline void* G4THitsCollection<T>::operator new(size_t)
120{
121 void* anHC;
122 anHC = (void*)anHCAllocator.MallocSingle();
123 return anHC;
124}
125
126template <class T> inline void G4THitsCollection<T>::operator delete(void* anHC)
127{
128 anHCAllocator.FreeSingle((G4HitsCollection*)anHC);
129}
130
131template <class T> G4THitsCollection<T>::G4THitsCollection()
132{
133 std::vector<T*> * theHitsCollection
134 = new std::vector<T*>;
135 theCollection = (void*)theHitsCollection;
136}
137
138template <class T> G4THitsCollection<T>::G4THitsCollection(G4String detName,G4String colNam)
139: G4HitsCollection(detName,colNam)
140{
141 std::vector<T*> * theHitsCollection
142 = new std::vector<T*>;
143 theCollection = (void*)theHitsCollection;
144}
145
146template <class T> G4THitsCollection<T>::~G4THitsCollection()
147{
148 std::vector<T*> * theHitsCollection
149 = (std::vector<T*>*)theCollection;
150 //theHitsCollection->clearAndDestroy();
151 for(size_t i=0;i<theHitsCollection->size();i++)
152 { delete (*theHitsCollection)[i]; }
153 theHitsCollection->clear();
154 delete theHitsCollection;
155}
156
157template <class T> G4int G4THitsCollection<T>::operator==(const G4THitsCollection<T> &right) const
158{ return (collectionName==right.collectionName); }
159
160template <class T> void G4THitsCollection<T>::DrawAllHits()
161{
162 std::vector<T*> * theHitsCollection
163 = (std::vector<T*>*)theCollection;
164 size_t n = theHitsCollection->size();
165 for(size_t i=0;i<n;i++)
166 { (*theHitsCollection)[i]->Draw(); }
167}
168
169template <class T> void G4THitsCollection<T>::PrintAllHits()
170{
171 std::vector<T*> * theHitsCollection
172 = (std::vector<T*>*)theCollection;
173 size_t n = theHitsCollection->size();
174 for(size_t i=0;i<n;i++)
175 { (*theHitsCollection)[i]->Print(); }
176}
177
178#endif
179
Note: See TracBrowser for help on using the repository browser.