BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSSamplerHit.hh
1 /*
2 Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3 University of London 2001 - 2018.
4 
5 This file is part of BDSIM.
6 
7 BDSIM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published
9 by the Free Software Foundation version 3 of the License.
10 
11 BDSIM is distributed in the hope that it will be useful, but
12 WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15 
16 You should have received a copy of the GNU General Public License
17 along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 #ifndef BDSSAMPLERHIT_H
20 #define BDSSAMPLERHIT_H
21 
22 #include "G4VHit.hh"
23 #include "G4THitsCollection.hh"
24 #include "G4Allocator.hh"
25 
26 #include "BDSParticle.hh"
27 
28 #include <ostream>
29 
35 class BDSSamplerHit: public G4VHit
36 {
37 public:
38  BDSSamplerHit(G4String aName,
39  G4int samplerID,
40  BDSParticle init,
41  BDSParticle prod,
42  BDSParticle last_scat,
43  BDSParticle local,
44  BDSParticle global,
45  G4double s,
46  G4double weight,
47  G4int PDGtype,
48  G4int nEvent,
49  G4int ParentID,
50  G4int TrackID,
51  G4int TurnsTaken,
52  G4String process,
53  G4int beamlineIndex);
54 
55  virtual ~BDSSamplerHit();
56 
57  inline void* operator new(size_t) ;
58  inline void operator delete(void *aHit);
59 
60 private:
62  G4String itsName;
63 
66  G4int itsSamplerID;
67 
70 
73 
76 
80 
83 
85  G4double itsS;
86 
87  G4double itsWeight;
88  G4int itsPDGtype;
89  G4int itsEventNo;
90  G4int itsParentID;
91  G4int itsTrackID;
92  G4int itsTurnsTaken;
93 
95  G4String itsProcess;
96 
99 
100 public:
101  inline G4int GetSamplerID() const
102  {return itsSamplerID;}
103  inline G4double GetInitTotalEnergy() const
104  {return itsInit.GetTotalEnergy();}
105  inline BDSParticle GetInit() const
106  {return itsInit;}
107  inline G4double GetInitX() const
108  {return itsInit.GetX();}
109  inline G4double GetInitXPrime() const
110  {return itsInit.GetXp();}
111  inline G4double GetInitY() const
112  {return itsInit.GetY();}
113  inline G4double GetInitYPrime() const
114  {return itsInit.GetYp();}
115  inline G4double GetInitZ() const
116  {return itsInit.GetZ();}
117  inline G4double GetInitZPrime() const
118  {return itsInit.GetZp();}
119  inline G4double GetInitT() const
120  {return itsInit.GetT();}
121  inline BDSParticle GetProd() const
122  {return itsProd;}
123  inline G4double GetProdTotalEnergy() const
124  {return itsProd.GetTotalEnergy();}
125  inline G4double GetProdX() const
126  {return itsProd.GetX();}
127  inline G4double GetProdXPrime() const
128  {return itsProd.GetXp();}
129  inline G4double GetProdY() const
130  {return itsProd.GetY();}
131  inline G4double GetProdYPrime() const
132  {return itsProd.GetYp();}
133  inline G4double GetProdZ() const
134  {return itsProd.GetZ();}
135  inline G4double GetProdZPrime() const
136  {return itsProd.GetZp();}
137  inline G4double GetProdT() const
138  {return itsProd.GetT();}
139  inline BDSParticle GetLastScat() const
140  {return itsLastScat;}
141  inline G4double GetLastScatTotalEnergy() const
142  {return itsLastScat.GetTotalEnergy();}
143  inline G4double GetLastScatX() const
144  {return itsLastScat.GetX();}
145  inline G4double GetLastScatXPrime() const
146  {return itsLastScat.GetXp();}
147  inline G4double GetLastScatY() const
148  {return itsLastScat.GetY();}
149  inline G4double GetLastScatYPrime() const
150  {return itsLastScat.GetYp();}
151  inline G4double GetLastScatZ() const
152  {return itsLastScat.GetZ();}
153  inline G4double GetLastScatZPrime() const
154  {return itsLastScat.GetZp();}
155  inline G4double GetLastScatT() const
156  {return itsLastScat.GetT();}
157  inline BDSParticle GetLocal() const
158  {return itsLocal;}
159  inline G4double GetTotalEnergy() const
160  {return itsLocal.GetTotalEnergy();}
161  inline G4double GetX() const
162  {return itsLocal.GetX();}
163  inline G4double GetXPrime() const
164  {return itsLocal.GetXp();}
165  inline G4double GetY() const
166  {return itsLocal.GetY();}
167  inline G4double GetYPrime() const
168  {return itsLocal.GetYp();}
169  inline G4double GetZ() const
170  {return itsLocal.GetZ();}
171  inline G4double GetZPrime() const
172  {return itsLocal.GetZp();}
173  inline BDSParticle GetGlobal() const
174  {return itsGlobal;}
175  inline G4double GetGlobalX() const
176  {return itsGlobal.GetX();}
177  inline G4double GetGlobalXPrime() const
178  {return itsGlobal.GetXp();}
179  inline G4double GetGlobalY() const
180  {return itsGlobal.GetY();}
181  inline G4double GetGlobalYPrime() const
182  {return itsGlobal.GetYp();}
183  inline G4double GetGlobalZ() const
184  {return itsGlobal.GetZ();}
185  inline G4double GetGlobalZPrime() const
186  {return itsGlobal.GetZp();}
187  inline G4double GetT() const
188  {return itsLocal.GetT();}
189  inline G4double GetS() const
190  {return itsS;}
191  inline G4double GetWeight() const
192  {return itsWeight;}
193  inline G4String GetName() const
194  {return itsName;}
195  inline G4int GetEventNo() const
196  {return itsEventNo;}
197  inline G4int GetPDGtype() const
198  {return itsPDGtype;}
199  inline G4int GetParentID() const
200  {return itsParentID;}
201  inline G4int GetTrackID() const
202  {return itsTrackID;}
203  inline G4int GetTurnsTaken() const
204  {return itsTurnsTaken;}
205  inline G4String GetProcess() const
206  {return itsProcess;}
207  inline G4int GetBeamlineIndex() const
208  {return itsBeamlineIndex;}
209 };
210 
212 extern G4Allocator<BDSSamplerHit> BDSSamplerHitAllocator;
213 
214 inline void* BDSSamplerHit::operator new(size_t)
215 {
216  void* aHit;
217  aHit=(void*) BDSSamplerHitAllocator.MallocSingle();
218  return aHit;
219 }
220 
221 inline void BDSSamplerHit::operator delete(void *aHit)
222 {
223  BDSSamplerHitAllocator.FreeSingle((BDSSamplerHit*) aHit);
224 }
225 
226 #endif
G4int itsBeamlineIndex
simulation model component index
G4String itsName
Name of the sampler this hit was generated by.
G4String itsProcess
creating process
BDSParticle itsGlobal
actual position and momentum direction in GLOBAL coordinates
The information recorded from a particle impacting a sampler.
BDSParticle itsLastScat
point where the particle last scattered
G4double itsS
total current track length
BDSParticle itsProd
point where the particle was produced
a particle definition
Definition: BDSParticle.hh:36
BDSParticle itsLocal
BDSParticle itsInit
initial particle track in GLOBAL coordinates