ImpactX
Loading...
Searching...
No Matches
ExactSbend.H
Go to the documentation of this file.
1/* Copyright 2022-2023 The Regents of the University of California, through Lawrence
2 * Berkeley National Laboratory (subject to receipt of any required
3 * approvals from the U.S. Dept. of Energy). All rights reserved.
4 *
5 * This file is part of ImpactX.
6 *
7 * Authors: Chad Mitchell, Axel Huebl
8 * License: BSD-3-Clause-LBNL
9 */
10#ifndef IMPACTX_EXACTSBEND_H
11#define IMPACTX_EXACTSBEND_H
12
14#include "mixin/alignment.H"
15#include "mixin/pipeaperture.H"
16#include "mixin/beamoptic.H"
17#include "mixin/thick.H"
19#include "mixin/named.H"
20#include "mixin/nofinalize.H"
21
22#include <AMReX_Extension.H>
23#include <AMReX_Math.H>
24#include <AMReX_REAL.H>
25#include <AMReX_SIMD.H>
26
27#include <cmath>
28#include <stdexcept>
29
30
31namespace impactx::elements
32{
34 : public mixin::Named,
35 public mixin::BeamOptic<ExactSbend>,
36 public mixin::LinearTransport<ExactSbend>,
37 public mixin::Thick,
38 public mixin::Alignment,
40 public mixin::NoFinalize,
41 public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
42 {
43 static constexpr auto type = "ExactSbend";
44 static constexpr amrex::ParticleReal degree2rad = ablastr::constant::math::pi / 180.0;
46
71 amrex::ParticleReal ds,
72 amrex::ParticleReal phi,
73 amrex::ParticleReal B,
74 amrex::ParticleReal dx = 0,
75 amrex::ParticleReal dy = 0,
76 amrex::ParticleReal rotation_degree = 0,
77 amrex::ParticleReal aperture_x = 0,
78 amrex::ParticleReal aperture_y = 0,
79 int nslice = 1,
80 std::optional<std::string> name = std::nullopt
81 )
82 : Named(std::move(name)),
83 Thick(ds, nslice),
84 Alignment(dx, dy, rotation_degree),
86 m_phi(phi * degree2rad), m_B(B)
87 {
88 }
89
91 amrex::ParticleReal
92 rc (RefPart const & refpart) const
93 {
94 using namespace amrex::literals; // for _rt and _prt
95
96 return m_B != 0_prt ? refpart.rigidity_Tm() / m_B : m_ds / m_phi;
97 }
98
100 using BeamOptic::operator();
101
102
110 void compute_constants (RefPart const & refpart)
111 {
112 using namespace amrex::literals; // for _rt and _prt
113 using amrex::Math::powi;
114
115 Alignment::compute_constants(refpart);
116
117 // reference particle's orbital radius
118 m_rc = this->rc(refpart);
119
120 // arc length and angle of the current slice
121 m_slice_ds = m_ds / nslice();
122 amrex::ParticleReal const bend_sign = m_ds * m_rc / std::abs(m_ds * m_rc);
123 m_slice_phi = bend_sign * std::abs(m_phi) / nslice(); // the sign is determined by rc/ds
124
125 // trigonometric evaluations
126 auto const [sin_phi, cos_phi] = amrex::Math::sincos(m_slice_phi);
127 m_sin_phi = sin_phi;
128 m_cos_phi = cos_phi;
129
130 // access reference particle values to find relativistic factors
131 m_ibeta = 1_prt / refpart.beta();
132 m_ibetagam2 = 1_prt / powi<2>(refpart.beta_gamma());
133
134 }
135
147 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
150 T_Real & AMREX_RESTRICT x,
151 T_Real & AMREX_RESTRICT y,
152 T_Real & AMREX_RESTRICT t,
153 T_Real & AMREX_RESTRICT px,
154 T_Real & AMREX_RESTRICT py,
155 T_Real & AMREX_RESTRICT pt,
156 T_IdCpu & AMREX_RESTRICT idcpu,
157 [[maybe_unused]] RefPart const & AMREX_RESTRICT refpart
158 ) const
159 {
160 using namespace amrex::literals; // for _rt and _prt
161 using amrex::Math::powi;
162 using namespace std; // for cmath(float)
163 namespace stdx = amrex::simd::stdx;
164
165 // shift due to alignment errors of the element
166 shift_in(x, y, px, py);
167
168 // access position data
169 T_Real const xout = x;
170 T_Real const yout = y;
171 T_Real const tout = t;
172
173 // initialize output values of momenta
174 T_Real pxout = px;
175 T_Real pyout = py;
176 T_Real ptout = pt;
177
178 // treat the special case of zero field
179 if (m_phi==0_prt && m_B==0_prt) {
180 // compute the radical in the denominator (= pz):
181 T_Real const inv_pzden = 1_prt / sqrt(
182 powi<2>(pt - m_ibeta) -
184 powi<2>(px) -
185 powi<2>(py)
186 );
187
188 // advance position and momentum (exact drift)
189 x = xout + m_slice_ds * px * inv_pzden;
190 // pxout = px;
191 y = yout + m_slice_ds * py * inv_pzden;
192 // pyout = py;
193 t = tout - m_slice_ds * (m_ibeta + (pt - m_ibeta) * inv_pzden);
194 // ptout = pt;
195
196 // assign updated momenta
197 px = pxout;
198 py = pyout;
199 pt = ptout;
200
201 } else {
202
203 // assign intermediate quantities
204 T_Real const pperp2 = powi<2>(pt)-2.0_prt * m_ibeta * pt - powi<2>(py)+1.0_prt;
205 T_Real const px2 = powi<2>(px);
206 // determine if particle lies within the domain of map definition
207 auto const mask = pperp2 <= px2;
209 {
210 T_Real const pperp = sqrt(pperp2);
211 T_Real const pzi = sqrt(powi<2>(pperp) - powi<2>(px));
212 T_Real const rho = m_rc + xout;
213
214 // update momenta
215 pxout = px * m_cos_phi + (pzi - rho / m_rc) * m_sin_phi;
216 pyout = py;
217 ptout = pt;
218
219 // angle of momentum rotation
220 T_Real const px2f = powi<2>(pxout);
221 // determine if particle lies within domain of map definition
222 auto const mask2 = pperp2 <= px2f;
224 {
225 T_Real const pzf = sqrt(powi<2>(pperp)-powi<2>(pxout));
226 T_Real const theta = m_slice_phi + asin(px/pperp) - asin(pxout/pperp);
227
228 // update position coordinates
229 x = -m_rc + rho*m_cos_phi + m_rc*(pzf + px*m_sin_phi - pzi*m_cos_phi);
230 y = yout + theta * m_rc * py;
231 t = tout - theta * m_rc * (pt - 1.0_prt * m_ibeta) - m_slice_phi * m_rc * m_ibeta;
232
233 // assign updated momenta
234 px = pxout;
235 py = pyout;
236 pt = ptout;
237 }
238 }
239 }
240
241 // apply transverse aperture
242 apply_aperture(x, y, idcpu);
243
244 // undo shift due to alignment errors of the element
245 shift_out(x, y, px, py);
246 }
247
253 void operator() (RefPart & AMREX_RESTRICT refpart) const
254 {
255 using namespace amrex::literals; // for _rt and _prt
256 using amrex::Math::powi;
257
258 // assign input reference particle values
259 amrex::ParticleReal const x = refpart.x;
260 amrex::ParticleReal const px = refpart.px;
261 amrex::ParticleReal const y = refpart.y;
262 amrex::ParticleReal const py = refpart.py;
263 amrex::ParticleReal const z = refpart.z;
264 amrex::ParticleReal const pz = refpart.pz;
265 amrex::ParticleReal const t = refpart.t;
266 amrex::ParticleReal const pt = refpart.pt;
267 amrex::ParticleReal const s = refpart.s;
268
269 // length of the current slice
270 amrex::ParticleReal const slice_ds = m_ds / nslice();
271
272 // special case of zero field = an exact drift
273 if (m_phi==0_prt && m_B==0_prt) {
274 // advance position and momentum (drift)
275 amrex::ParticleReal const step = slice_ds /std::sqrt(powi<2>(pt)-1.0_prt);
276 refpart.x = x + step*px;
277 refpart.y = y + step*py;
278 refpart.z = z + step*pz;
279 refpart.t = t - step*pt;
280
281 } else {
282
283 // assign intermediate parameters
284 amrex::ParticleReal const rc = (m_B != 0_prt) ? refpart.rigidity_Tm() / m_B : m_ds / m_phi;
285 amrex::ParticleReal const B = refpart.beta_gamma() /rc;
286 amrex::ParticleReal const bend_sign = m_ds * rc / std::abs(m_ds * rc);
287 amrex::ParticleReal const theta = bend_sign * std::abs(m_phi) / nslice(); // sign is determined by rc*ds
288
289 // calculate expensive terms once
290 auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
291
292 // advance position and momentum (bend)
293 refpart.px = px*cos_theta - pz*sin_theta;
294 refpart.py = py;
295 refpart.pz = pz*cos_theta + px*sin_theta;
296 refpart.pt = pt;
297
298 refpart.x = x + (refpart.pz - pz)/B;
299 refpart.y = y + (theta/B)*py;
300 refpart.z = z - (refpart.px - px)/B;
301 refpart.t = t - (theta/B)*pt;
302
303 }
304
305 // advance integrated path length
306 refpart.s = s + slice_ds;
307 }
308
310 using LinearTransport::operator();
311
318 Map6x6
319 transport_map ([[maybe_unused]] RefPart const & AMREX_RESTRICT refpart) const
320 {
321 throw std::runtime_error(std::string(type) + ": Envelope tracking is not yet implemented!");
322 return Map6x6::Identity();
323 }
324
325 amrex::ParticleReal m_phi;
326 amrex::ParticleReal m_B;
327
328 private:
329 // constants that are independent of the individually tracked particle,
330 // see: compute_constants() to refresh
332 };
333
334} // namespace impactx
335
336#endif // IMPACTX_EXACTSBEND_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
#define AMREX_GPU_HOST
Array4< int const > mask
static constexpr amrex::Real pi
__host__ __device__ std::pair< double, double > sincos(double x)
constexpr T powi(T x) noexcept
Definition All.H:54
@ s
fixed s as the independent variable
Definition ImpactXParticleContainer.H:37
@ t
fixed t as the independent variable
Definition ImpactXParticleContainer.H:38
amrex::SmallMatrix< amrex::ParticleReal, 6, 6, amrex::Order::F, 1 > Map6x6
Definition CovarianceMatrix.H:20
__host__ __device__ void make_invalid() const noexcept
static constexpr __host__ __device__ SmallMatrix< T, NRows, NCols, ORDER, StartIndex > Identity() noexcept
Definition ReferenceParticle.H:31
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal beta_gamma() const
Definition ReferenceParticle.H:110
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal rigidity_Tm() const
Definition ReferenceParticle.H:203
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal beta() const
Definition ReferenceParticle.H:94
amrex::ParticleReal m_rc
Definition ExactSbend.H:331
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition ExactSbend.H:319
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal rc(RefPart const &refpart) const
Definition ExactSbend.H:92
amrex::ParticleReal m_cos_phi
Definition ExactSbend.H:331
amrex::ParticleReal m_slice_ds
magnetic field in T
Definition ExactSbend.H:331
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT t, T_Real &AMREX_RESTRICT px, T_Real &AMREX_RESTRICT py, T_Real &AMREX_RESTRICT pt, T_IdCpu &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition ExactSbend.H:149
amrex::ParticleReal m_sin_phi
Definition ExactSbend.H:331
ImpactXParticleContainer::ParticleType PType
Definition ExactSbend.H:45
ExactSbend(amrex::ParticleReal ds, amrex::ParticleReal phi, amrex::ParticleReal B, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, amrex::ParticleReal aperture_x=0, amrex::ParticleReal aperture_y=0, int nslice=1, std::optional< std::string > name=std::nullopt)
Definition ExactSbend.H:70
static constexpr amrex::ParticleReal degree2rad
Definition ExactSbend.H:44
amrex::ParticleReal m_ibetagam2
Definition ExactSbend.H:331
amrex::ParticleReal m_ibeta
Definition ExactSbend.H:331
void compute_constants(RefPart const &refpart)
Definition ExactSbend.H:110
amrex::ParticleReal m_slice_phi
Definition ExactSbend.H:331
static constexpr auto type
Definition ExactSbend.H:43
amrex::ParticleReal m_B
bend angle in radians
Definition ExactSbend.H:326
amrex::ParticleReal m_phi
Definition ExactSbend.H:325
Definition alignment.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void shift_out(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT px, T_Real &AMREX_RESTRICT py) const
Definition alignment.H:109
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dy() const
Definition alignment.H:146
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dx() const
Definition alignment.H:136
Alignment(amrex::ParticleReal dx, amrex::ParticleReal dy, amrex::ParticleReal rotation_degree)
Definition alignment.H:36
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void shift_in(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT px, T_Real &AMREX_RESTRICT py) const
Definition alignment.H:78
Definition beamoptic.H:219
Definition lineartransport.H:29
Definition named.H:29
AMREX_GPU_HOST Named(std::optional< std::string > name)
Definition named.H:57
AMREX_FORCE_INLINE std::string name() const
Definition named.H:122
Definition nofinalize.H:22
Definition pipeaperture.H:26
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void apply_aperture(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_IdCpu &AMREX_RESTRICT idcpu) const
Definition pipeaperture.H:59
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal aperture_x() const
Definition pipeaperture.H:90
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal aperture_y() const
Definition pipeaperture.H:101
PipeAperture(amrex::ParticleReal aperture_x, amrex::ParticleReal aperture_y)
Definition pipeaperture.H:32
Definition thick.H:24
Thick(amrex::ParticleReal ds, int nslice)
Definition thick.H:30
amrex::ParticleReal m_ds
Definition thick.H:58
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal ds() const
Definition thick.H:53
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int nslice() const
Definition thick.H:43