BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSHitSampler.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef BDSHITSAMPLER_H
20#define BDSHITSAMPLER_H
21
22#include "BDSParticleCoordsFull.hh"
23
24#include "globals.hh"
25#include "G4VHit.hh"
26#include "G4THitsCollection.hh"
27#include "G4Allocator.hh"
28
35class BDSHitSampler: public G4VHit
36{
37public:
38 BDSHitSampler(G4int samplerIDIn,
39 const BDSParticleCoordsFull& coordsIn,
40 G4double momentumIn,
41 G4double massIn,
42 G4double chargeIn,
43 G4double rigidityIn,
44 G4int pdgIDIn,
45 G4int parentIDIn,
46 G4int trackIDIn,
47 G4int turnsTakenIn,
48 G4int beamlineIndexIn,
49 G4int nElectronsIn = 0);
50
52 virtual ~BDSHitSampler();
53
54 inline void* operator new(size_t);
55 inline void operator delete(void *aHit);
56
57 G4int samplerID;
59 G4double momentum;
60 G4double mass;
61 G4double charge;
62 G4double rigidity;
63 G4int pdgID;
64 G4int parentID;
65 G4int trackID;
66 G4int turnsTaken;
67 G4int beamlineIndex;
68 G4int nElectrons;
69
70private:
71 BDSHitSampler() = delete;
72};
73
75extern G4Allocator<BDSHitSampler> BDSAllocatorSampler;
76
77inline void* BDSHitSampler::operator new(size_t)
78{
79 void* aHit;
80 aHit=(void*) BDSAllocatorSampler.MallocSingle();
81 return aHit;
82}
83
84inline void BDSHitSampler::operator delete(void *aHit)
85{
86 BDSAllocatorSampler.FreeSingle((BDSHitSampler*) aHit);
87}
88
89#endif
The information recorded from a particle impacting a sampler.
virtual ~BDSHitSampler()
Note this should not be inline when we use a G4Allocator.
G4int nElectrons
Can only get this at inspection time so include here.
BDSHitSampler()=delete
No default constructor.
G4double charge
Double as g4 uses charge as a double.
A set of particle coordinates including energy and weight.