Geant4 10.7.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//
28
29#ifndef G4THitsCollection_h
30#define G4THitsCollection_h 1
31
32#include "G4VHitsCollection.hh"
33#include "G4Allocator.hh"
34#include "globals.hh"
35//#include "g4rw/tpordvec.h"
36#include <vector>
37
38// class description:
39//
40// This is a template class of hits collection and parametrized by
41// The concrete class of G4VHit. This is a uniform collection for
42// a particular concrete hit class objects.
43// An intermediate layer class G4HitsCollection appeared in this
44// header file is used just for G4Allocator, because G4Allocator
45// cannot be instansiated with a template class. Thus G4HitsCollection
46// class MUST NOT be directly used by the user.
47
49{
50 public:
52 G4HitsCollection(G4String detName,G4String colNam);
53 virtual ~G4HitsCollection();
54 G4bool operator==(const G4HitsCollection &right) const;
55
56 protected:
58};
59
60#if defined G4DIGI_ALLOC_EXPORT
62#else
64#endif
65
66template <class T> class G4THitsCollection : public G4HitsCollection
67{
68 public:
70 public: // with description
71 G4THitsCollection(G4String detName,G4String colNam);
72 // constructor.
73 public:
74 virtual ~G4THitsCollection();
75 G4bool operator==(const G4THitsCollection<T> &right) const;
76
77 inline void *operator new(size_t);
78 inline void operator delete(void* anHC);
79 public: // with description
80 virtual void DrawAllHits();
81 virtual void PrintAllHits();
82 // These two methods invokes Draw() and Print() methods of all of
83 // hit objects stored in this collection, respectively.
84
85 public: // with description
86 inline T* operator[](size_t i) const
87 {
89 return (*((std::vector<T*>*)theCollection))[i];
90 }
91 // Returns a pointer to a concrete hit object.
92 inline std::vector<T*>* GetVector() const
93 {
95 return (std::vector<T*>*)theCollection; }
96 // Returns a collection vector.
97 inline size_t insert(T* aHit)
98 {
100 std::vector<T*>*theHitsCollection = (std::vector<T*>*)theCollection;
101 theHitsCollection->push_back(aHit);
102 return theHitsCollection->size();
103 }
104 // Insert a hit object. Total number of hit objects stored in this
105 // collection is returned.
106 inline size_t entries() const
107 {
109 std::vector<T*>*theHitsCollection = (std::vector<T*>*)theCollection;
110 return theHitsCollection->size();
111 }
112 // Returns the number of hit objects stored in this collection
113
114 public:
115 virtual G4VHit* GetHit(size_t i) const
116 {
118 return (*((std::vector<T*>*)theCollection))[i];
119 }
120 virtual size_t GetSize() const
121 {
123 return ((std::vector<T*>*)theCollection)->size();
124 }
125};
126
127template <class T> inline void* G4THitsCollection<T>::operator new(size_t)
128{
131 void* anHC;
132 anHC = (void*)anHCAllocator.MallocSingle();
133 return anHC;
134}
135
136template <class T> inline void G4THitsCollection<T>::operator delete(void* anHC)
137{
140 anHCAllocator.FreeSingle((G4HitsCollection*)anHC);
141}
142
144{
146 std::vector<T*> * theHitsCollection = new std::vector<T*>;
147 theCollection = (void*)theHitsCollection;
148}
149
151: G4HitsCollection(detName,colNam)
152{
154 std::vector<T*> * theHitsCollection = new std::vector<T*>;
155 theCollection = (void*)theHitsCollection;
156}
157
159{
161 std::vector<T*> * theHitsCollection = (std::vector<T*>*)theCollection;
162 //theHitsCollection->clearAndDestroy();
163 for(size_t i=0;i<theHitsCollection->size();++i)
164 { delete (*theHitsCollection)[i]; }
165 theHitsCollection->clear();
166 delete theHitsCollection;
167}
168
170{
172 return (collectionName==right.collectionName);
173}
174
175template <class T> void G4THitsCollection<T>::DrawAllHits()
176{
178 std::vector<T*> * theHitsCollection = (std::vector<T*>*)theCollection;
179 size_t n = theHitsCollection->size();
180 for(size_t i=0;i<n;++i)
181 { (*theHitsCollection)[i]->Draw(); }
182}
183
184template <class T> void G4THitsCollection<T>::PrintAllHits()
185{
187 std::vector<T*> * theHitsCollection
188 = (std::vector<T*>*)theCollection;
189 size_t n = theHitsCollection->size();
190 for(size_t i=0;i<n;++i)
191 { (*theHitsCollection)[i]->Print(); }
192}
193
194#endif
195
G4DLLIMPORT G4Allocator< G4HitsCollection > *& anHCAllocator_G4MT_TLS_()
#define G4DLLIMPORT
Definition: G4Types.hh:69
#define G4DLLEXPORT
Definition: G4Types.hh:68
bool G4bool
Definition: G4Types.hh:86
void FreeSingle(Type *anElement)
Definition: G4Allocator.hh:206
Type * MallocSingle()
Definition: G4Allocator.hh:196
G4bool operator==(const G4HitsCollection &right) const
virtual ~G4HitsCollection()
size_t entries() const
virtual void DrawAllHits()
virtual size_t GetSize() const
std::vector< T * > * GetVector() const
G4bool operator==(const G4THitsCollection< T > &right) const
virtual G4VHit * GetHit(size_t i) const
size_t insert(T *aHit)
T * operator[](size_t i) const
virtual void PrintAllHits()
Definition: G4VHit.hh:48