BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSHitApertureImpact.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 BDSHITAPERTUREIMPACT_H
20#define BDSHITAPERTUREIMPACT_H
21
22#include "globals.hh"
23#include "G4VHit.hh"
24#include "G4ThreeVector.hh"
25#include "G4THitsCollection.hh"
26#include "G4Allocator.hh"
27
36class BDSHitApertureImpact: public G4VHit
37{
38public:
40 BDSHitApertureImpact(G4double totalEnergyIn,
41 G4double preStepKineticEnergyIn,
42 G4double SIn,
43 G4double xIn,
44 G4double yIn,
45 G4double pxIn,
46 G4double pyIn,
47 G4double weightIn,
48 G4double globalTimeIn,
49 G4int partIDIn,
50 G4int trackIDIn,
51 G4int parentIDIn,
52 G4int turnsTakenIn,
53 G4int beamlineIndexIn,
54 G4int nElectronsIn);
55
57 virtual ~BDSHitApertureImpact();
58
59 inline void* operator new(size_t);
60 inline void operator delete(void *aHit);
61
62 G4double totalEnergy;
63 G4double preStepKineticEnergy;
64 G4double S;
65 G4double x;
66 G4double y;
67 G4double xp;
68 G4double yp;
69 G4double weight;
70 G4double globalTime;
71 G4int partID;
72 G4int trackID;
73 G4int parentID;
74 G4int turnsTaken;
75 G4int beamlineIndex;
76 G4int nElectrons;
77};
78
80extern G4Allocator<BDSHitApertureImpact> BDSAllocatorApertureImpacts;
81
82inline void* BDSHitApertureImpact::operator new(size_t)
83{
84 void* aHit;
85 aHit=(void*) BDSAllocatorApertureImpacts.MallocSingle();
86 return aHit;
87}
88
89inline void BDSHitApertureImpact::operator delete(void *aHit)
90{
91 BDSAllocatorApertureImpacts.FreeSingle((BDSHitApertureImpact*) aHit);
92}
93
94#endif
Snapshot of information for particle passing through a collimator.
virtual ~BDSHitApertureImpact()
Note this should not be inline when we use a G4Allocator.