BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSPTCOneTurnMap.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 BDSPTCONETURNMAP_H
20#define BDSPTCONETURNMAP_H
21
22#include "BDSParticleCoordsFullGlobal.hh"
23
24#include "globals.hh" // Geant4 typedefs
25#include "G4Track.hh"
26
27#include <vector>
28#include <set>
29
31
41{
42public:
44 {
45 G4double coefficient;
46 G4int nx;
47 G4int npx;
48 G4int ny;
49 G4int npy;
50 G4int ndeltaP;
51 };
52
53 BDSPTCOneTurnMap() = delete;
54 BDSPTCOneTurnMap(const BDSPTCOneTurnMap &other) = default;
55 BDSPTCOneTurnMap(BDSPTCOneTurnMap &&other) noexcept = default;
56 virtual ~BDSPTCOneTurnMap() {;}
57
59 BDSPTCOneTurnMap(const G4String& maptableFile,
60 const BDSParticleDefinition* designParticle);
61
63 G4bool ShouldApplyToPrimary(G4double momentum, G4int turnstaken);
64
67 G4bool offsetS0);
68
70 void UpdateCoordinates(G4ThreeVector localPosition,
71 G4ThreeVector localMomentum,
72 G4int turnstaken);
73
75 void GetThisTurn(G4double& x,
76 G4double& px,
77 G4double& y,
78 G4double& py,
79 G4double& pz,
80 G4int turnstaken);
81
82private:
83 G4double Evaluate(std::vector<PTCMapTerm>& terms,
84 G4double x,
85 G4double px,
86 G4double y,
87 G4double py,
88 G4double deltaP) const;
89
90 G4double initialPrimaryMomentum;
91 G4bool beamOffsetS0;
92 G4double referenceMomentum;
93 G4double mass;
94
95 std::set<G4int> turnsScattered;
96
97 G4int lastTurnNumber;
98 G4double xLastTurn;
99 G4double pxLastTurn;
100 G4double yLastTurn;
101 G4double pyLastTurn;
102 G4double deltaPLastTurn;
103
104 std::vector<PTCMapTerm> xTerms;
105 std::vector<PTCMapTerm> yTerms;
106 std::vector<PTCMapTerm> pxTerms;
107 std::vector<PTCMapTerm> pyTerms;
108 std::vector<PTCMapTerm> deltaPTerms;
109};
110
111#endif
Class to load and use PTC 1 turn map.
G4bool ShouldApplyToPrimary(G4double momentum, G4int turnstaken)
Decides whether or not this should be applied. Can add more.
virtual ~BDSPTCOneTurnMap()
Destructor.
BDSPTCOneTurnMap(BDSPTCOneTurnMap &&other) noexcept=default
Move constructor.
void GetThisTurn(G4double &x, G4double &px, G4double &y, G4double &py, G4double &pz, G4int turnstaken)
Return the coordinates for this turn.
BDSPTCOneTurnMap()=delete
Default constructor.
void UpdateCoordinates(G4ThreeVector localPosition, G4ThreeVector localMomentum, G4int turnstaken)
Update coordinates if the last turn is greater than the number of turns taken.
void SetInitialPrimaryCoordinates(const BDSParticleCoordsFullGlobal &coords, G4bool offsetS0)
Load initial coordinates where the beam started and convert to PTC coordinates.
BDSPTCOneTurnMap(const BDSPTCOneTurnMap &other)=default
Copy constructor.
A set of particle coordinates in both local and global.
Wrapper for particle definition.