BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBunch.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 BDSBUNCH_H
20#define BDSBUNCH_H
21
22#include "BDSBunchType.hh"
23#include "BDSParticleCoords.hh"
24#include "BDSParticleCoordsFull.hh"
25#include "BDSParticleCoordsFullGlobal.hh"
26#include "BDSParticleDefinition.hh"
27
28#include "globals.hh"
29#include "G4Transform3D.hh"
30
31class BDSBeamline;
32
33namespace GMAD
34{
35 class Beam;
36}
37
47{
48public:
49 BDSBunch();
50 explicit BDSBunch(const G4String& nameIn);
51 virtual ~BDSBunch();
52
56
62 virtual void SetOptions(const BDSParticleDefinition* beamParticle,
63 const GMAD::Beam& beam,
64 const BDSBunchType& distrType,
65 G4Transform3D beamlineTransformIn = G4Transform3D::Identity,
66 const G4double beamlineS = 0);
67
70 virtual void CheckParameters();
71
73 virtual void Initialise();
74
78
80 virtual G4bool ExpectChangingParticleType() const {return false;}
81
87 virtual BDSParticleCoordsFullGlobal GetNextParticleValid(G4int maxTries = 100);
88
92 virtual void BeginOfRunAction(G4int numberOfEvents,
93 G4bool batchMode);
94
96 inline virtual const BDSParticleDefinition* ParticleDefinition() const {return particleDefinition;}
97
101 virtual void SetGeneratePrimariesOnly(G4bool generatePrimariesOnlyIn);
102
107
110
116 virtual void RecreateAdvanceToEvent(G4int eventOffset);
117
119 inline G4bool BeamParticleIsAnIon() const {return particleDefinition->IsAnIon();}
120
125 virtual void UpdateIonDefinition();
126
130
132 static G4double CalculateZp(G4double xp, G4double yp, G4double Zp0);
133
137 static void SetEmittances(const BDSParticleDefinition* beamParticle,
138 const GMAD::Beam& beam,
139 G4double& emittGeometricX,
140 G4double& emittGeometricY,
141 G4double& emittNormalisedX,
142 G4double& emittNormalisedY);
143
145 inline G4String Name() const {return name;}
146
148 inline G4int CurrentBunchIndex() const {return currentBunchIndex;}
149
151 void CalculateBunchIndex(G4int eventIndex);
152
153protected:
158
160 void ApplyTilt(BDSParticleCoordsFull& localIn) const;
161
164
167
168 G4String name;
169
171 G4double X0;
172 G4double Y0;
173 G4double Z0;
174 G4double S0;
175 G4double T0;
176 G4double Xp0;
177 G4double Yp0;
178 G4double Zp0;
179 G4double E0;
180 G4double P0;
181 G4double tilt;
182 G4double sigmaT;
183 G4double sigmaP;
184 G4double sigmaE;
185 G4double sigmaEk;
187
192 G4double bunchPeriod;
194
197
200
204
205 G4bool finiteTilt;
210
214
215private:
217 G4Transform3D beamlineTransform;
218
219 G4double beamlineS;
220 G4double mass2;
221
223 mutable const BDSBeamline* beamline;
224};
225
226#endif
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
The base class for bunch distribution generators.
Definition: BDSBunch.hh:47
BDSParticleCoordsFullGlobal ApplyCurvilinearTransform(const BDSParticleCoordsFull &localIn) const
Calculate the global coordinates from curvilinear coordinates of a beam line.
Definition: BDSBunch.cc:322
G4double mass2
Cache of mass squared as required to convert from p to E.
Definition: BDSBunch.hh:220
G4double sigmaE
Centre of distributions.
Definition: BDSBunch.hh:184
G4bool finiteTilt
Definition: BDSBunch.hh:205
virtual void BeginOfRunAction(G4int numberOfEvents, G4bool batchMode)
Definition: BDSBunch.cc:285
const BDSBeamline * beamline
A reference to the fully constructed beamline that's lazily instantiated.
Definition: BDSBunch.hh:223
G4bool ParticleDefinitionHasBeenUpdated() const
Definition: BDSBunch.hh:129
virtual G4bool ExpectChangingParticleType() const
A hint of whether we expect to require and extended particle set (ie pions, kaons,...
Definition: BDSBunch.hh:80
bool useBunchTiming
Bunch offset in time parameters.
Definition: BDSBunch.hh:189
G4double Yp0
Centre of distributions.
Definition: BDSBunch.hh:177
G4double T0
Centre of distributions.
Definition: BDSBunch.hh:175
G4bool finiteSigmaE
Flags to ignore random number generator in case of no finite E or T.
Definition: BDSBunch.hh:207
G4double tilt
Centre of distributions.
Definition: BDSBunch.hh:181
BDSParticleCoordsFullGlobal ApplyTransform(const BDSParticleCoordsFull &localIn) const
Definition: BDSBunch.cc:292
virtual void RecreateAdvanceToEvent(G4int eventOffset)
Definition: BDSBunch.cc:280
void ApplyTilt(BDSParticleCoordsFull &localIn) const
Apply a rotation about unitZ for the local coordinates according to member variable tilt.
Definition: BDSBunch.cc:304
G4double S0
Centre of distributions.
Definition: BDSBunch.hh:174
G4int currentBunchIndex
Bunch offset in time parameters.
Definition: BDSBunch.hh:190
G4double P0
central momentum
Definition: BDSBunch.hh:180
G4bool particleDefinitionHasBeenUpdated
Definition: BDSBunch.hh:203
G4bool generatePrimariesOnly
Definition: BDSBunch.hh:213
virtual void SetGeneratePrimariesOnly(G4bool generatePrimariesOnlyIn)
Definition: BDSBunch.cc:289
G4double sigmaP
Centre of distributions.
Definition: BDSBunch.hh:183
G4bool useCurvilinear
Whether to ignore z and use s and transform for curvilinear coordinates.
Definition: BDSBunch.hh:196
G4Transform3D beamlineTransform
Transform that beam line starts with that will also be applied to coordinates.
Definition: BDSBunch.hh:217
void ApplyBunchTiming(BDSParticleCoordsFullGlobal &localIn) const
Add on the offset in T for the current bunch number (i*bunchPeriod).
Definition: BDSBunch.cc:316
G4double sigmaT
Centre of distributions.
Definition: BDSBunch.hh:182
G4double Z0
Centre of distributions.
Definition: BDSBunch.hh:173
virtual const BDSParticleDefinition * ParticleDefinition() const
Access the beam particle definition.
Definition: BDSBunch.hh:96
void CalculateBunchIndex(G4int eventIndex)
Calculate which bunch index we should be at given an event index.
Definition: BDSBunch.cc:215
virtual BDSParticleCoordsFull GetNextParticleLocal()
Definition: BDSBunch.cc:272
G4double bunchPeriod
Bunch offset in time parameters.
Definition: BDSBunch.hh:192
G4int CurrentBunchIndex() const
Get the current bunch index for writing to output.
Definition: BDSBunch.hh:148
G4double X0
Centre of distributions.
Definition: BDSBunch.hh:171
G4double Zp0
Centre of distributions.
Definition: BDSBunch.hh:178
G4bool BeamParticleIsAnIon() const
Access whether the beam particle is an ion or not.
Definition: BDSBunch.hh:119
virtual void Initialise()
Any initialisation - to be used after SetOptions, then CheckParameters.
Definition: BDSBunch.cc:231
G4double Xp0
Centre of distributions.
Definition: BDSBunch.hh:176
G4double E0
Centre of distributions.
Definition: BDSBunch.hh:179
G4double sigmaEk
Centre of distributions.
Definition: BDSBunch.hh:185
static void SetEmittances(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, G4double &emittGeometricX, G4double &emittGeometricY, G4double &emittNormalisedX, G4double &emittNormalisedY)
Definition: BDSBunch.cc:175
G4double Y0
Centre of distributions.
Definition: BDSBunch.hh:172
virtual void SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
Definition: BDSBunch.cc:82
G4String Name() const
Distribution name.
Definition: BDSBunch.hh:145
virtual void UpdateIonDefinition()
Definition: BDSBunch.cc:382
G4double beamlineS
Beamline initial S position.
Definition: BDSBunch.hh:219
BDSParticleCoordsFullGlobal GetNextParticle()
Definition: BDSBunch.cc:262
G4bool finiteSigmaT
Flags to ignore random number generator in case of no finite E or T.
Definition: BDSBunch.hh:208
static G4double CalculateZp(G4double xp, G4double yp, G4double Zp0)
Calculate zp safely based on other components.
Definition: BDSBunch.cc:367
BDSParticleDefinition * particleDefinition
Particle definition for bunch - this class owns it.
Definition: BDSBunch.hh:199
virtual BDSParticleCoordsFullGlobal GetNextParticleValid(G4int maxTries=100)
Definition: BDSBunch.cc:234
G4String name
Name of distribution.
Definition: BDSBunch.hh:168
G4int eventsPerBunch
Bunch offset in time parameters.
Definition: BDSBunch.hh:191
G4bool UseCurvilinearTransform() const
Access whether there's a finite S offset and therefore we're using a CL transform.
Definition: BDSBunch.hh:109
virtual void CheckParameters()
Definition: BDSBunch.cc:223
A set of particle coordinates in both local and global.
A set of particle coordinates including energy and weight.
Wrapper for particle definition.
G4bool IsAnIon() const
Accessor.
Loader to use any HepMC3 compatible file.
Loader to read a specific sampler from a BDSIM ROOT output file.
Improve type-safety of native enum data type in C++.
Beam loader.
Definition: Beam.hh:37
Beam class.
Definition: beam.h:44
Parser namespace for GMAD language. Combination of Geant4 and MAD.