Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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:
59};
60
61#if defined G4DIGI_ALLOC_EXPORT
63#else
65#endif
66
67template <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();
76 G4int operator==(const G4THitsCollection<T> &right) const;
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 { return (*((std::vector<T*>*)theCollection))[i]; }
89 // Returns a pointer to a concrete hit object.
90 inline std::vector<T*>* GetVector() const
91 { return (std::vector<T*>*)theCollection; }
92 // Returns a collection vector.
93 inline G4int insert(T* aHit)
94 {
95 std::vector<T*>*theHitsCollection
96 = (std::vector<T*>*)theCollection;
97 theHitsCollection->push_back(aHit);
98 return theHitsCollection->size();
99 }
100 // Insert a hit object. Total number of hit objects stored in this
101 // collection is returned.
102 inline G4int entries() const
103 {
104 std::vector<T*>*theHitsCollection
105 = (std::vector<T*>*)theCollection;
106 return theHitsCollection->size();
107 }
108 // Returns the number of hit objects stored in this collection
109
110 public:
111 virtual G4VHit* GetHit(size_t i) const
112 { return (*((std::vector<T*>*)theCollection))[i]; }
113 virtual size_t GetSize() const
114 { return ((std::vector<T*>*)theCollection)->size(); }
115
116};
117
118template <class T> inline void* G4THitsCollection<T>::operator new(size_t)
119{
120 void* anHC;
121 anHC = (void*)anHCAllocator.MallocSingle();
122 return anHC;
123}
124
125template <class T> inline void G4THitsCollection<T>::operator delete(void* anHC)
126{
127 anHCAllocator.FreeSingle((G4HitsCollection*)anHC);
128}
129
131{
132 std::vector<T*> * theHitsCollection
133 = new std::vector<T*>;
134 theCollection = (void*)theHitsCollection;
135}
136
138: G4HitsCollection(detName,colNam)
139{
140 std::vector<T*> * theHitsCollection
141 = new std::vector<T*>;
142 theCollection = (void*)theHitsCollection;
143}
144
146{
147 std::vector<T*> * theHitsCollection
148 = (std::vector<T*>*)theCollection;
149 //theHitsCollection->clearAndDestroy();
150 for(size_t i=0;i<theHitsCollection->size();i++)
151 { delete (*theHitsCollection)[i]; }
152 theHitsCollection->clear();
153 delete theHitsCollection;
154}
155
156template <class T> G4int G4THitsCollection<T>::operator==(const G4THitsCollection<T> &right) const
157{ return (collectionName==right.collectionName); }
158
159template <class T> void G4THitsCollection<T>::DrawAllHits()
160{
161 std::vector<T*> * theHitsCollection
162 = (std::vector<T*>*)theCollection;
163 size_t n = theHitsCollection->size();
164 for(size_t i=0;i<n;i++)
165 { (*theHitsCollection)[i]->Draw(); }
166}
167
168template <class T> void G4THitsCollection<T>::PrintAllHits()
169{
170 std::vector<T*> * theHitsCollection
171 = (std::vector<T*>*)theCollection;
172 size_t n = theHitsCollection->size();
173 for(size_t i=0;i<n;i++)
174 { (*theHitsCollection)[i]->Print(); }
175}
176
177#endif
178
G4DLLIMPORT G4Allocator< G4HitsCollection > anHCAllocator
#define G4DLLIMPORT
Definition: G4Types.hh:56
#define G4DLLEXPORT
Definition: G4Types.hh:55
int G4int
Definition: G4Types.hh:66
G4int operator==(const G4HitsCollection &right) const
virtual ~G4HitsCollection()
G4int operator==(const G4THitsCollection< T > &right) const
virtual void DrawAllHits()
virtual size_t GetSize() const
std::vector< T * > * GetVector() const
virtual G4VHit * GetHit(size_t i) const
G4int insert(T *aHit)
T * operator[](size_t i) const
virtual void PrintAllHits()
Definition: G4VHit.hh:49