BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSHitCollimator.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 BDSHITCOLLIMATOR_H
20#define BDSHITCOLLIMATOR_H
21
22#include "globals.hh"
23#include "G4VHit.hh"
24#include "G4ThreeVector.hh"
25#include "G4THitsCollection.hh"
26#include "G4Allocator.hh"
27
28class BDSBeamline;
30
39class BDSHitCollimator: public G4VHit
40{
41public:
42 BDSHitCollimator() = delete;
43 BDSHitCollimator(const BDSBeamline* beamlineIn,
44 G4int collimatorIndexIn,
45 const G4ThreeVector& preStepPositionIn,
46 const G4ThreeVector& preStepMomentumIn,
47 G4double totalEnergyIn,
48 BDSHitEnergyDeposition* energyDepositionHitIn);
49
51 virtual ~BDSHitCollimator();
52
53 inline void* operator new(size_t);
54 inline void operator delete(void* aHit);
55
58 G4ThreeVector preStepPosition;
59 G4ThreeVector preStepMomentum;
60 G4double totalEnergy;
61
65};
66
68extern G4Allocator<BDSHitCollimator> BDSAllocatorCollimator;
69
70inline void* BDSHitCollimator::operator new(size_t)
71{
72 void* aHit;
73 aHit=(void*) BDSAllocatorCollimator.MallocSingle();
74 return aHit;
75}
76
77inline void BDSHitCollimator::operator delete(void* aHit)
78{
79 BDSAllocatorCollimator.FreeSingle((BDSHitCollimator*) aHit);
80}
81
82#endif
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
Snapshot of information for particle passing through a collimator.
G4double totalEnergy
Total energy of particle.
virtual ~BDSHitCollimator()
Note this should not be inline when we use a G4Allocator.
G4ThreeVector preStepMomentum
Local pre step point momentum.
G4int collimatorIndex
Index of collimator the hit is in.
const BDSBeamline * beamline
Which beam line the collimator is in.
BDSHitEnergyDeposition * energyDepositionHit
G4ThreeVector preStepPosition
Local pre step point (z from centre of object).
Information recorded for a single piece of energy deposition.