BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBunch.cc
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#include "BDSAcceleratorModel.hh"
20#include "BDSBeamline.hh"
21#include "BDSBunch.hh"
22#include "BDSDebug.hh"
23#include "BDSException.hh"
24#include "BDSGlobalConstants.hh"
25#include "BDSIonDefinition.hh"
26#include "BDSParticleCoords.hh"
27#include "BDSParticleCoordsFull.hh"
28#include "BDSParticleCoordsFullGlobal.hh"
29#include "BDSParticleDefinition.hh"
30#include "BDSPhysicsUtilities.hh"
31#include "BDSUtilities.hh"
32
33#include "parser/beam.h"
34
35#include "G4IonTable.hh"
36#include "G4ParticleTable.hh"
37#include "G4ThreeVector.hh"
38#include "G4Transform3D.hh"
39#include "G4TwoVector.hh"
40#include "G4Version.hh"
41
42#include "CLHEP/Geometry/Point3D.h"
43
44#include <cmath>
45#include <limits>
46#include <set>
47#include <string>
48
49
50BDSBunch::BDSBunch():
51 BDSBunch("reference")
52{;}
53
54BDSBunch::BDSBunch(const G4String& nameIn):
55 name(nameIn),
56 X0(0.0), Y0(0.0), Z0(0.0), S0(0.0), T0(0.0),
57 Xp0(0.0), Yp0(0.0), Zp0(0.0), E0(0.0), P0(0.0),
58 tilt(0.0),
59 sigmaT(0.0), sigmaP(0.0), sigmaE(0.0), sigmaEk(0.0),
60 useCurvilinear(false),
61 particleDefinition(nullptr),
62 particleDefinitionHasBeenUpdated(false),
63 finiteTilt(false),
64 finiteSigmaE(true),
65 finiteSigmaT(true),
66 generatePrimariesOnly(false),
67 beamlineTransform(G4Transform3D()),
68 beamlineS(0),
69 mass2(0.0),
70 beamline(nullptr)
71{;}
72
73BDSBunch::~BDSBunch()
74{
75 delete particleDefinition;
76}
77
79 const GMAD::Beam& beam,
80 const BDSBunchType& /*distrType*/,
81 G4Transform3D beamlineTransformIn,
82 G4double beamlineSIn)
83{
84 particleDefinition = new BDSParticleDefinition(*beamParticle); // copy it so this instance owns it
85
86 // back the starting point up by length safety to avoid starting on a boundary
87 G4ThreeVector unitZBeamline = G4ThreeVector(0,0,-1).transform(beamlineTransformIn.getRotation());
88 G4ThreeVector translation = BDSGlobalConstants::Instance()->LengthSafety() * unitZBeamline;
89 beamlineTransform = G4Transform3D(beamlineTransformIn.getRotation(), beamlineTransformIn.getTranslation()+translation);
90
91 beamlineS = beamlineSIn;
92
93 X0 = beam.X0 * CLHEP::m;
94 Y0 = beam.Y0 * CLHEP::m;
95 Z0 = beam.Z0 * CLHEP::m;
96 S0 = beam.S0 * CLHEP::m;
97 T0 = beam.T0 * CLHEP::s;
98 Xp0 = beam.Xp0 * CLHEP::rad;
99 Yp0 = beam.Yp0 * CLHEP::rad;
100 E0 = particleDefinition->TotalEnergy(); // already calculated and set earlier depending on available parameters
102 tilt = beam.tilt * CLHEP::rad;
103 sigmaT = beam.sigmaT;
104 sigmaP = beam.sigmaP;
105 sigmaE = beam.sigmaE;
106 sigmaEk = beam.sigmaEk;
107
111 G4bool finiteSigmaP = BDS::IsFinite(sigmaP);
112 G4bool finiteSigmaEk = BDS::IsFinite(sigmaEk);
113
114 std::set<std::string> keysDesign = {"sigmaE", "sigmaEk", "sigmaP"};
115 G4int nSetDesign = BDS::NBeamParametersSet(beam, keysDesign);
116 BDS::ConflictingParametersSet(beam, keysDesign, nSetDesign, false, "GeV");// warn only if too many set
117 if (finiteSigmaE)
118 {
119 sigmaP = (1./std::pow(beamParticle->Beta(),2)) * sigmaE; // dE/E = (beta^2) dP/P
120 sigmaEk = (beamParticle->TotalEnergy() / beamParticle->KineticEnergy()) * sigmaE;
121 }
122 else if (finiteSigmaP)
123 {
124 sigmaE = std::pow(beamParticle->Beta(),2) * sigmaP;
125 sigmaEk = (beamParticle->TotalEnergy() / beamParticle->KineticEnergy()) * sigmaE;
126 }
127 else if (finiteSigmaEk)
128 {
129 sigmaE = sigmaEk * (beamParticle->KineticEnergy() / beamParticle->TotalEnergy());
130 sigmaP = (1./std::pow(beamParticle->Beta(),2)) * sigmaE; // dE/E = (beta^2) dP/P
131 }
132 // else they'll all be 0 - no need for a calculation
133
134 finiteSigmaE = finiteSigmaE || finiteSigmaP || finiteSigmaEk; // finiteSigmaE used to know whether any variation in other classes
135 if (finiteSigmaE)
136 {
137 G4cout << "Beam> sigmaP: " << sigmaP << G4endl;
138 G4cout << "Beam> sigmaE: " << sigmaE << G4endl;
139 G4cout << "Beam> sigmaEk: " << sigmaEk << G4endl;
140 }
141
142 Zp0 = CalculateZp(Xp0,Yp0,beam.Zp0);
143
144 if (S0 > beamlineS)
145 {
146#ifdef BDSDEBUG
147 G4cout << __METHOD_NAME__ << "using curvilinear transform" << G4endl;
148#endif
149 if (BDS::IsFinite(Z0))
150 {throw BDSException(__METHOD_NAME__, "both Z0 and S0 are defined - please define only one!");}
151 useCurvilinear = true;
152 }
153}
154
156 const GMAD::Beam& beam,
157 G4double& emittGeometricX,
158 G4double& emittGeometricY,
159 G4double& emittNormalisedX,
160 G4double& emittNormalisedY)
161{
162 std::set<std::string> keysDesignX = {"emitx", "emitnx"};
163 G4int nSetDesignX = BDS::NBeamParametersSet(beam, keysDesignX);
164 BDS::ConflictingParametersSet(beam, keysDesignX, nSetDesignX);
165 if (BDS::IsFinite(beam.emitNX))
166 {
167 emittNormalisedX = G4double(beam.emitNX);
168 emittGeometricX = G4double(beam.emitNX) / beamParticle->Gamma();
169 }
170 else
171 {
172 emittGeometricX = G4double(beam.emitx);
173 emittNormalisedX = G4double(beam.emitx) * beamParticle->Gamma();
174 }
175
176 std::set<std::string> keysDesignY = {"emity", "emitny"};
177 G4int nSetDesignY = BDS::NBeamParametersSet(beam, keysDesignY);
178 BDS::ConflictingParametersSet(beam, keysDesignY, nSetDesignY);
179 if (BDS::IsFinite(beam.emitNY))
180 {
181 emittNormalisedY = G4double(beam.emitNY);
182 emittGeometricY = G4double(beam.emitNY) / beamParticle->Gamma();}
183 else
184 {
185 emittGeometricY = G4double(beam.emity);
186 emittNormalisedY = G4double(beam.emity) * beamParticle->Gamma();
187 }
188
189 G4cout << __METHOD_NAME__ << "Geometric (x): " << emittGeometricX
190 << ", Normalised (x): " << emittNormalisedX << G4endl;
191 G4cout << __METHOD_NAME__ << "Geometric (y): " << emittGeometricY
192 << ", Normalised (y): " << emittNormalisedY << G4endl;
193}
194
196{
197 if (sigmaE < 0)
198 {throw BDSException(__METHOD_NAME__, "sigmaE " + std::to_string(sigmaE) + " < 0!");}
199 if (sigmaT < 0)
200 {throw BDSException(__METHOD_NAME__, "sigmaT " + std::to_string(sigmaT) + " < 0!");}
201}
202
204{;}
205
207{
208 particleDefinitionHasBeenUpdated = false; // reset flag for this call
209 // use a separate flag to record whether the particle definitions has
210 // been updated as subsequent calls to GetNextParticle may reset it to
211 // false but it was updated in the first call
212 G4bool flag = false;
213
214 // continue generating particles until positive finite kinetic energy.
215 G4int n = 0;
217 while (n < maxTries) // prevent infinite loops
218 {
219 ++n;
220 coords = GetNextParticle();
222
223 // ensure total energy is greater than the rest mass
224 if ((coords.local.totalEnergy - particleDefinition->Mass()) > 0)
225 {break;}
226 }
227 if (n >= maxTries)
228 {throw BDSException(__METHOD_NAME__, "unable to generate coordinates above rest mass after 100 attempts.");}
229
231 return coords;
232}
233
235{
236 particleDefinitionHasBeenUpdated = false; // reset flag
238 if (finiteTilt)
239 {ApplyTilt(local);}
241 return all;
242}
243
245{
247 Xp0, Yp0, Zp0,
248 T0, S0, E0, /*weight=*/1.0);
249 return local;
250}
251
252void BDSBunch::BeginOfRunAction(G4int /*numberOfEvents*/)
253{;}
254
255void BDSBunch::SetGeneratePrimariesOnly(G4bool generatePrimariesOnlyIn)
256{generatePrimariesOnly = generatePrimariesOnlyIn;}
257
259{
260 if (useCurvilinear) // i.e. S0 is finite - use beam line
261 {return ApplyCurvilinearTransform(localIn);}
262 else // just general beam line offset
264}
265
267{
268 G4TwoVector xy(localIn.x, localIn.y);
269 G4TwoVector xpyp(localIn.xp, localIn.yp);
270 xy.rotate(tilt);
271 xpyp.rotate(tilt);
272 localIn.x = xy.x();
273 localIn.y = xy.y();
274 localIn.xp = xpyp.x();
275 localIn.yp = xpyp.y();
276}
277
279{
280 if (generatePrimariesOnly) // no beam line built so no possible transform
281 {return BDSParticleCoordsFullGlobal(localIn, (BDSParticleCoords)localIn);}
282
283 if (!beamline)
284 {// initialise cache of beam line pointer
285 beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
286 if (!beamline)
287 {throw BDSException(__METHOD_NAME__, "no beamline constructed!");}
288 }
289
290 // 'c' for curvilinear
291 G4int beamlineIndex = 0;
292 G4double S = S0 + localIn.z;
293 if (S < 0)
294 {throw BDSException(__METHOD_NAME__, "Negative S detected for particle.");}
295
296 G4Transform3D cTrans = beamline->GetGlobalEuclideanTransform(S,
297 localIn.x,
298 localIn.y,
299 &beamlineIndex);
300 // rotate the momentum vector
301 G4ThreeVector cMom = G4ThreeVector(localIn.xp, localIn.yp, localIn.zp).transform(cTrans.getRotation());
302 // translation contains displacement from origin already - including any local offset
303 G4ThreeVector cPos = cTrans.getTranslation();
304
305 BDSParticleCoords global = BDSParticleCoords(cPos.x(), cPos.y(), cPos.z(),
306 cMom.x(), cMom.y(), cMom.z(),
307 localIn.T);
308
310 result.beamlineIndex = beamlineIndex;
311#ifdef BDSDEBUG
312 G4cout << __METHOD_NAME__ << result << G4endl;
313#endif
314 return result;
315}
316
317G4double BDSBunch::CalculateZp(G4double xp, G4double yp, G4double Zp0In)
318{
319 G4double zp;
320 G4double transMom = std::pow(xp, 2) + std::pow(yp, 2);
321
322 if (transMom > (1 - std::numeric_limits<double>::epsilon()))
323 {throw BDSException(__METHOD_NAME__, "xp, yp too large, xp: " + std::to_string(xp) + " yp: " + std::to_string(yp));}
324 if (Zp0In < 0)
325 {zp = -std::sqrt(1.0 - transMom);}
326 else
327 {zp = std::sqrt(1.0 - transMom);}
328
329 return zp;
330}
331
333{
335 {return;}
336
337 G4IonTable* ionTable = G4ParticleTable::GetParticleTable()->GetIonTable();
339 G4ParticleDefinition* ionParticleDef = ionTable->GetIon(ionDefinition->Z(),
340 ionDefinition->A(),
341 ionDefinition->ExcitationEnergy());
343 // Note we don't need to take care of electrons here. These are automatically
344 // allocated by Geant4 when it converts the primary vertex to a dynamic particle
345 // (in the process of constructing a track from it) (done in G4PrimaryTransformer)
346 // this relies on the charge being set correctly - Geant4 detects this isn't the same
347 // as Z and adds electrons accordingly.
348#if G4VERSION_NUMBER > 1049
349 // in the case of ions the particle definition is only available now
350 // fix the looping thresholds now it's available
352#endif
353}
const BDSBeamline * BeamlineMain() const
Accessor.
G4Transform3D GetGlobalEuclideanTransform(G4double s, G4double x=0, G4double y=0, G4int *indexOfFoundElement=nullptr) const
Definition: BDSBeamline.cc:598
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:278
G4double sigmaE
Centre of distributions.
Definition: BDSBunch.hh:173
G4bool finiteTilt
Definition: BDSBunch.hh:187
const BDSBeamline * beamline
A reference to the fully constructed beamline that's lazily instantiated.
Definition: BDSBunch.hh:205
G4double Yp0
Centre of distributions.
Definition: BDSBunch.hh:166
G4double T0
Centre of distributions.
Definition: BDSBunch.hh:164
G4bool finiteSigmaE
Flags to ignore random number generator in case of no finite E or T.
Definition: BDSBunch.hh:189
G4double tilt
Centre of distributions.
Definition: BDSBunch.hh:170
BDSParticleCoordsFullGlobal ApplyTransform(const BDSParticleCoordsFull &localIn) const
Definition: BDSBunch.cc:258
void ApplyTilt(BDSParticleCoordsFull &localIn) const
Apply a rotation about unitZ for the local coordinates according to member variable tilt.
Definition: BDSBunch.cc:266
G4double S0
Centre of distributions.
Definition: BDSBunch.hh:163
G4double P0
central momentum
Definition: BDSBunch.hh:169
G4bool particleDefinitionHasBeenUpdated
Definition: BDSBunch.hh:185
G4bool generatePrimariesOnly
Definition: BDSBunch.hh:195
virtual void SetGeneratePrimariesOnly(G4bool generatePrimariesOnlyIn)
Definition: BDSBunch.cc:255
G4double sigmaP
Centre of distributions.
Definition: BDSBunch.hh:172
G4bool useCurvilinear
Whether to ignore z and use s and transform for curvilinear coordinates.
Definition: BDSBunch.hh:178
G4Transform3D beamlineTransform
Transform that beam line starts with that will also be applied to coordinates.
Definition: BDSBunch.hh:199
G4double sigmaT
Centre of distributions.
Definition: BDSBunch.hh:171
G4double Z0
Centre of distributions.
Definition: BDSBunch.hh:162
virtual BDSParticleCoordsFull GetNextParticleLocal()
Definition: BDSBunch.cc:244
G4double X0
Centre of distributions.
Definition: BDSBunch.hh:160
G4double Zp0
Centre of distributions.
Definition: BDSBunch.hh:167
virtual void Initialise()
Any initialisation - to be used after SetOptions, then CheckParameters.
Definition: BDSBunch.cc:203
G4double Xp0
Centre of distributions.
Definition: BDSBunch.hh:165
G4double E0
Centre of distributions.
Definition: BDSBunch.hh:168
G4double sigmaEk
Centre of distributions.
Definition: BDSBunch.hh:174
static void SetEmittances(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, G4double &emittGeometricX, G4double &emittGeometricY, G4double &emittNormalisedX, G4double &emittNormalisedY)
Definition: BDSBunch.cc:155
G4double Y0
Centre of distributions.
Definition: BDSBunch.hh:161
virtual void BeginOfRunAction(G4int numberOfEvents)
Definition: BDSBunch.cc:252
virtual void SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
Definition: BDSBunch.cc:78
virtual void UpdateIonDefinition()
Definition: BDSBunch.cc:332
G4double beamlineS
Beamline initial S position.
Definition: BDSBunch.hh:201
BDSParticleCoordsFullGlobal GetNextParticle()
Definition: BDSBunch.cc:234
G4bool finiteSigmaT
Flags to ignore random number generator in case of no finite E or T.
Definition: BDSBunch.hh:190
static G4double CalculateZp(G4double xp, G4double yp, G4double Zp0)
Calculate zp safely based on other components.
Definition: BDSBunch.cc:317
BDSParticleDefinition * particleDefinition
Particle definition for bunch - this class owns it.
Definition: BDSBunch.hh:181
virtual BDSParticleCoordsFullGlobal GetNextParticleValid(G4int maxTries=100)
Definition: BDSBunch.cc:206
virtual void CheckParameters()
Definition: BDSBunch.cc:195
General exception with possible name of object and message.
Definition: BDSException.hh:35
static BDSGlobalConstants * Instance()
Access method.
Class to parse an ion particle definition.
G4int A() const
Accessor.
G4double ExcitationEnergy() const
Accessor.
G4int Z() const
Accessor.
A set of particle coordinates in both local and global.
G4int beamlineIndex
Optional index for where transform was found.
A set of particle coordinates including energy and weight.
A set of particle coordinates.
BDSParticleCoords ApplyTransform(const G4Transform3D &transform) const
Apply a transform to the coordinates and return a copy of them transformed.
Wrapper for particle definition.
G4double Mass() const
Accessor.
G4double Momentum() const
Accessor.
G4double KineticEnergy() const
Accessor.
G4double Beta() const
Accessor.
G4double Gamma() const
Accessor.
BDSIonDefinition * IonDefinition() const
Accessor.
G4double TotalEnergy() const
Accessor.
G4bool IsAnIon() const
Accessor.
void UpdateG4ParticleDefinition(G4ParticleDefinition *particleIn)
Improve type-safety of native enum data type in C++.
double sigmaE
for the gaussian, elliptic shell, ring distributions
Definition: beamBase.h:77
double emity
initial twiss parameters
Definition: beamBase.h:82
double emitNY
initial twiss parameters
Definition: beamBase.h:83
double tilt
tilt of beam applied as rotation about unit local z
Definition: beamBase.h:71
double emitNX
initial twiss parameters
Definition: beamBase.h:83
double S0
initial beam centroid
Definition: beamBase.h:62
double emitx
initial twiss parameters
Definition: beamBase.h:82
double sigmaT
bunch length
Definition: beamBase.h:74
double Zp0
initial beam centroid
Definition: beamBase.h:63
double Z0
initial beam centroid
Definition: beamBase.h:62
double X0
initial beam centroid
Definition: beamBase.h:62
double Xp0
initial beam centroid
Definition: beamBase.h:63
double Yp0
initial beam centroid
Definition: beamBase.h:63
double Y0
initial beam centroid
Definition: beamBase.h:62
double T0
initial beam centroid
Definition: beamBase.h:64
Beam class.
Definition: beam.h:44
void ConflictingParametersSet(const GMAD::Beam &beamDefinition, const std::set< std::string > &keys, G4int nSet, G4bool warnZeroParamsSet=true, const G4String &unitString="")
Throw an exception if too few or too many parameters are set for the supplied keys.
void FixGeant105ThreshholdsForParticle(const G4ParticleDefinition *particleDefinition)
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())
G4int NBeamParametersSet(const GMAD::Beam &beamDefinition, const std::set< std::string > &keys)
Count how many out of the set of keys in a beam instance are set.