Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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>
68{
69 public:
71
72 public: // with description
73 G4THitsCollection(G4String detName, G4String colNam);
74 // constructor.
75 public:
76 virtual ~G4THitsCollection();
77 G4bool operator==(const G4THitsCollection<T>& right) const;
78
79 inline void* operator new(size_t);
80 inline void operator delete(void* anHC);
81
82 public: // with description
83 virtual void DrawAllHits();
84 virtual void PrintAllHits();
85 // These two methods invokes Draw() and Print() methods of all of
86 // hit objects stored in this collection, respectively.
87
88 public: // with description
89 inline T* operator[](size_t i) const
90 {
91 return (*((std::vector<T*>*) theCollection))[i];
92 }
93 // Returns a pointer to a concrete hit object.
94 inline std::vector<T*>* GetVector() const
95 {
96 return (std::vector<T*>*) theCollection;
97 }
98 // Returns a collection vector.
99 inline size_t insert(T* aHit)
100 {
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 size_t entries() const
108 {
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 {
117 return (*((std::vector<T*>*) theCollection))[i];
118 }
119 virtual size_t GetSize() const
120 {
121 return ((std::vector<T*>*) theCollection)->size();
122 }
123};
124
125template <class T>
126inline void* G4THitsCollection<T>::operator new(size_t)
127{
128 if(anHCAllocator_G4MT_TLS_() == nullptr)
129 {
131 }
132 return (void*) anHCAllocator_G4MT_TLS_()->MallocSingle();
133}
134
135template <class T>
136inline void G4THitsCollection<T>::operator delete(void* anHC)
137{
138 anHCAllocator_G4MT_TLS_()->FreeSingle((G4HitsCollection*) anHC);
139}
140
141template <class T>
143{
144 std::vector<T*>* theHitsCollection = new std::vector<T*>;
145 theCollection = (void*) theHitsCollection;
146}
147
148template <class T>
150 : G4HitsCollection(detName, colNam)
151{
152 std::vector<T*>* theHitsCollection = new std::vector<T*>;
153 theCollection = (void*) theHitsCollection;
154}
155
156template <class T>
158{
159 std::vector<T*>* theHitsCollection = (std::vector<T*>*) theCollection;
160 // theHitsCollection->clearAndDestroy();
161 for(size_t i = 0; i < theHitsCollection->size(); ++i)
162 {
163 delete(*theHitsCollection)[i];
164 }
165 theHitsCollection->clear();
166 delete theHitsCollection;
167}
168
169template <class T>
171{
172 return (collectionName == right.collectionName);
173}
174
175template <class T>
177{
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 {
182 (*theHitsCollection)[i]->Draw();
183 }
184}
185
186template <class T>
188{
189 std::vector<T*>* theHitsCollection = (std::vector<T*>*) theCollection;
190 size_t n = theHitsCollection->size();
191 for(size_t i = 0; i < n; ++i)
192 {
193 (*theHitsCollection)[i]->Print();
194 }
195}
196
197#endif
G4DLLIMPORT G4Allocator< G4HitsCollection > *& anHCAllocator_G4MT_TLS_()
#define G4DLLIMPORT
Definition: G4Types.hh:69
#define G4DLLEXPORT
Definition: G4Types.hh:68
bool G4bool
Definition: G4Types.hh:86
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