BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSHitThinThing.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 BDSHITTHINTHING_H
20#define BDSHITTHINTHING_H
21
22#include "G4Allocator.hh"
23#include "G4VHit.hh"
24#include "G4THitsCollection.hh"
25#include "G4Types.hh"
26
27#include <vector>
28
32
33class BDSHitThinThing; // forward declaration to allow typedef required for static function
35extern G4Allocator<BDSHitThinThing> BDSAllocatorThinThing;
36
45class BDSHitThinThing: public G4VHit
46{
47public:
48 BDSHitThinThing(G4int pdgIDIn,
49 G4int trackIDIn,
50 G4int parentIDIn,
51 G4int turnsTakenIn,
52 BDSTrajectoryPoint* hitIn);
53
55 virtual ~BDSHitThinThing();
56
57 inline void* operator new(size_t);
58 inline void operator delete(void* aHit);
59
62
63 G4int pdgID;
64 G4int trackID;
65 G4int parentID;
66 G4int turnsTaken;
68
70 static std::vector<const BDSTrajectoryPoint*> TrajectoryPointsFromHC(BDSHitsCollectionThinThing* hits);
71
72 static std::vector<const BDSTrajectoryPointHit*>
73 ResolvePossibleEarlierThinHits(const std::vector<const BDSTrajectoryPrimary*>& primaryTrajectoryHits,
74 const BDSHitsCollectionThinThing* thinThingHits);
75
76private:
77 BDSHitThinThing() = delete;
78};
79
80inline void* BDSHitThinThing::operator new(size_t)
81{
82 void* aHit;
83 aHit=(void*) BDSAllocatorThinThing.MallocSingle();
84 return aHit;
85}
86
87inline void BDSHitThinThing::operator delete(void* aHit)
88{
89 BDSAllocatorThinThing.FreeSingle((BDSHitThinThing*) aHit);
90}
91#endif
A hit if a particle lost energy in a thin object.
BDSHitThinThing()=delete
No default constructor.
BDSTrajectoryPointHit * GetNewTrajectoryPointHit() const
Allocate and return a new hit object.
virtual ~BDSHitThinThing()
Note this should not be inline when we use a G4Allocator.
static std::vector< const BDSTrajectoryPoint * > TrajectoryPointsFromHC(BDSHitsCollectionThinThing *hits)
Utility function to get a vector of trajectory points from the hits collection.
A summary trajectory object of a loss point.
A Point in a trajectory with extra information.
Trajectory information for only the primary.