BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBunchHaloFlatSigma.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 "BDSBunchHaloFlatSigma.hh"
20#include "BDSDebug.hh"
21#include "BDSException.hh"
22#include "BDSParticleCoordsFull.hh"
23
24#include "parser/beam.h"
25
26#include "globals.hh" // geant4 types / globals
27
28#include "CLHEP/Random/RandFlat.h"
29#include "CLHEP/Units/SystemOfUnits.h"
30#include "Randomize.hh"
31
32#include <algorithm>
33#include <cmath>
34#include <iterator>
35#include <utility>
36#include <vector>
37
38
39namespace BDS
40{
41 PhaseSpaceCoord PhaseSpaceCoordFromActionAngle(ActionAngleCoord aa, const TwissPair& tp)
42 {
43 G4double action = aa.action;
44 G4double angle = aa.angle;
45 G4double alpha = tp.alpha;
46 G4double beta = tp.beta;
47 G4double x = std::sqrt(2 * action * beta) * std::cos(angle);
48 G4double xp = (-std::sqrt(2 * action / beta) * (std::sin(angle) + alpha * std::cos(angle)));
49 PhaseSpaceCoord result = {x, xp};
50 return result;
51 }
52
53 EllipsePointGenerator::EllipsePointGenerator(G4double actionIn,
54 const TwissPair& tp,
55 CLHEP::RandFlat* flatRandomGeneratorIn):
56 action(actionIn),
57 twisspair(tp),
58 flatRandomGenerator(flatRandomGeneratorIn)
59 {
60 // Here: populate angles vector from 0 to 2pi and the corresponding arc
61 // lengths from angle=0 to each angle.
62 G4int npoints = 100;
63 std::vector<PhaseSpaceCoord> points;
64 points.reserve(npoints);
65 // Sample points on the ellipse equidistantly.
66 // leq (<=) so that can can interpolate up to 2pi.
67 for (G4int i = 0; i <= npoints; ++i)
68 {
69 G4double angle = i * CLHEP::twopi / npoints;
70 angles.push_back(angle);
71 ActionAngleCoord aa = {action, angle};
72 points.push_back(PhaseSpaceCoordFromActionAngle(aa, twisspair));
73 }
74 // Get the cumulative distances from angle=0 to each angle. The distance
75 // between the first and second point is of course non-zero. So we have to
76 // consider the case where we roll a pathLength between 0 and the first
77 // distance. To do this we insert 0 at the start.
78 std::vector<G4double> cumulativeDistances = CumulativeDistances(points);
79 cumulativeDistances.insert(std::begin(cumulativeDistances), 0);
80 pathLengths = std::move(cumulativeDistances);
81 }
82
83 PhaseSpaceCoord EllipsePointGenerator::GetRandomPointOnEllipse() const
84 {
85 // Select random point on the perimeter
86 G4double pathLength = EllipsePerimeter() * flatRandomGenerator->shoot();
87 // Invert to get the angle and then use that to get the (x, xp) pair.
88 G4double angle = PathLengthToAngle(pathLength);
89 ActionAngleCoord aa = {action, angle};
90 return PhaseSpaceCoordFromActionAngle(aa, twisspair);
91 }
92
93 G4double EllipsePointGenerator::PathLengthToAngle(G4double pathLength) const
94 {
95 // Find position of pathLength in pathLengths
96 auto it = std::lower_bound(std::begin(pathLengths), std::end(pathLengths), pathLength);
97 auto n = std::distance(std::begin(pathLengths), it);
98
99 // Get path lengths and angles either side of the point we are trying to
100 // evaluate.
101 G4double s1 = *it;
102 G4double s0 = *(--it);
103 G4double angle0 = *std::next(std::begin(angles), n - 1);
104 G4double angle1 = *std::next(std::begin(angles), n);
105
106 // Do the linear interpolation
107 G4double angle = angle0 + (pathLength - s0) * ((angle1 - angle0) / (s1 - s0));
108 return angle;
109 }
110}
111
112BDSBunchHaloFlatSigma::BDSBunchHaloFlatSigma():
113 BDSBunch("halo-nsigma-flat"),
114 alphaX(0.0), alphaY(0.0),
115 betaX(0.0), betaY(0.0),
116 emitX(0.0), emitY(0.0),
117 sigmaX(0.0), sigmaY(0.0),
118 haloNSigmaXInner(0.0), haloNSigmaXOuter(0.0),
119 haloNSigmaYInner(0.0), haloNSigmaYOuter(0.0)
120{
121 flatGen = new CLHEP::RandFlat(*CLHEP::HepRandom::getTheEngine());
122}
123
125 const GMAD::Beam& beam,
126 const BDSBunchType& distrType,
127 G4Transform3D beamlineTransformIn,
128 const G4double beamlineSIn)
129{
130 BDSBunch::SetOptions(beamParticle, beam, distrType, beamlineTransformIn, beamlineSIn);
131 alphaX = G4double(beam.alfx);
132 alphaY = G4double(beam.alfy);
133 betaX = G4double(beam.betx);
134 betaY = G4double(beam.bety);
135 haloNSigmaXInner = G4double(beam.haloNSigmaXInner);
136 haloNSigmaXOuter = G4double(beam.haloNSigmaXOuter);
137 haloNSigmaYInner = G4double(beam.haloNSigmaYInner);
138 haloNSigmaYOuter = G4double(beam.haloNSigmaYOuter);
139
140 G4double ex, ey; // dummy variables we don't need
141 SetEmittances(beamParticle, beam, emitX, emitY, ex, ey);
142
143 sigmaX = std::sqrt(emitX * betaX);
144 sigmaY = std::sqrt(emitY * betaY);
145
147}
148
150{
151 // Sampler flat in nsigma between lower and upper.
152 G4double xsigrange = (haloNSigmaXOuter - haloNSigmaXInner);
153 G4double ysigrange = (haloNSigmaYOuter - haloNSigmaYInner);
154 G4double xnsig = haloNSigmaXInner + flatGen->shoot() * xsigrange;
155 G4double ynsig = haloNSigmaYInner + flatGen->shoot() * ysigrange;
156
157 // Courant-snyder invariants for the particle, which are actually twice the
158 // particle actions...
159 G4double csix = std::pow(xnsig * sigmaX, 2) / betaX;
160 G4double csiy = std::pow(ynsig * sigmaY, 2) / betaY;
161 // The particle actions, which we use here.
162 G4double jx = csix / 2.0;
163 G4double jy = csiy / 2.0;
164
167
170 auto xps = epgx.GetRandomPointOnEllipse();
171 auto yps = epgy.GetRandomPointOnEllipse();
172
173 G4double x = (xps.position + X0) * CLHEP::m;
174 G4double xp = (xps.momentum + Xp0) * CLHEP::rad;
175 G4double y = (yps.position + Y0) * CLHEP::m;
176 G4double yp = (yps.momentum + Yp0) * CLHEP::rad;
177 G4double z = 0.;
178 G4double zp = CalculateZp(xp, yp, Zp0);
179 G4double weight = 1.0;
180 BDSParticleCoordsFull result(x, y, z, xp, yp, zp, T0, S0 + z, E0, weight);
181 return result;
182}
183
185{
187
188 if (emitX <= 0)
189 {throw BDSException(__METHOD_NAME__, "emitx must be > 0!");}
190 if (emitY <= 0)
191 {throw BDSException(__METHOD_NAME__, "emity must be > 0!");}
192
193 if (haloNSigmaXInner <= 0)
194 {throw BDSException(__METHOD_NAME__, "haloNSigmaXInner <= 0");}
195
196 if (haloNSigmaYInner <= 0)
197 {throw BDSException(__METHOD_NAME__, "haloNSigmaYInner <= 0");}
198
199 if (haloNSigmaXInner > haloNSigmaXOuter)
200 {throw BDSException(__METHOD_NAME__, "haloNSigmaXInner cannot be less than haloNSigmaXOuter");}
201
202 if (haloNSigmaYInner > haloNSigmaYOuter)
203 {throw BDSException(__METHOD_NAME__, "haloNSigmaYInner cannot be less than haloNSigmaYOuter");}
204}
G4double alphaX
Twiss parameter.
G4double sigmaY
Twiss parameter.
G4double betaY
Twiss parameter.
virtual void SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
G4double sigmaX
Twiss parameter.
G4double emitX
Twiss parameter.
G4double betaX
Twiss parameter.
G4double emitY
Twiss parameter.
virtual BDSParticleCoordsFull GetNextParticleLocal()
G4double alphaY
Twiss parameter.
The base class for bunch distribution generators.
Definition: BDSBunch.hh:47
G4double Yp0
Centre of distributions.
Definition: BDSBunch.hh:166
G4double T0
Centre of distributions.
Definition: BDSBunch.hh:164
G4double S0
Centre of distributions.
Definition: BDSBunch.hh:163
G4double X0
Centre of distributions.
Definition: BDSBunch.hh:160
G4double Zp0
Centre of distributions.
Definition: BDSBunch.hh:167
G4double Xp0
Centre of distributions.
Definition: BDSBunch.hh:165
G4double E0
Centre of distributions.
Definition: BDSBunch.hh:168
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 SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
Definition: BDSBunch.cc:78
static G4double CalculateZp(G4double xp, G4double yp, G4double Zp0)
Calculate zp safely based on other components.
Definition: BDSBunch.cc:317
virtual void CheckParameters()
Definition: BDSBunch.cc:195
General exception with possible name of object and message.
Definition: BDSException.hh:35
A set of particle coordinates including energy and weight.
Wrapper for particle definition.
Improve type-safety of native enum data type in C++.
Class for generating points uniformly on ellipse perimeters via interpolation. Part of implementation...
double alfy
initial twiss parameters
Definition: beamBase.h:82
double haloNSigmaYInner
for the halo distribution
Definition: beamBase.h:117
double betx
initial twiss parameters
Definition: beamBase.h:82
double haloNSigmaXInner
for the halo distribution
Definition: beamBase.h:115
double haloNSigmaYOuter
for the halo distribution
Definition: beamBase.h:118
double haloNSigmaXOuter
for the halo distribution
Definition: beamBase.h:116
double alfx
initial twiss parameters
Definition: beamBase.h:82
double bety
initial twiss parameters
Definition: beamBase.h:82
Beam class.
Definition: beam.h:44
Return either G4Tubs or G4CutTubs depending on flat face.
Simple struct for storing action/angle pairs to aid readability. Implementation detail.
Simple struct for storing position/momentum pairs to aid readability. Implementation detail.
Simple struct for storing Twiss alpha/beta pairs to aid readability. Implementation detail.