ImpactX
Loading...
Searching...
No Matches
ExactQuad.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_EXACTQUAD_H
11#define IMPACTX_EXACTQUAD_H
12
15#include "mixin/alignment.H"
16#include "mixin/pipeaperture.H"
17#include "mixin/beamoptic.H"
18#include "mixin/thick.H"
19#include "mixin/named.H"
20#include "mixin/nofinalize.H"
22
23#include <AMReX_Extension.H>
24#include <AMReX_Math.H>
25#include <AMReX_REAL.H>
26#include <AMReX_SIMD.H>
27
28#include <cmath>
29#include <stdexcept>
30
31
32namespace impactx::elements
33{
34 struct ExactQuad
35 : public mixin::Named,
36 public mixin::BeamOptic<ExactQuad>,
37 public mixin::LinearTransport<ExactQuad>,
38 public mixin::Thick,
39 public mixin::Alignment,
41 public mixin::NoFinalize,
42 public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
43 {
44 static constexpr auto type = "ExactQuad";
46
74 amrex::ParticleReal ds,
75 amrex::ParticleReal k,
76 int unit,
77 amrex::ParticleReal dx = 0,
78 amrex::ParticleReal dy = 0,
79 amrex::ParticleReal rotation_degree = 0,
80 amrex::ParticleReal aperture_x = 0,
81 amrex::ParticleReal aperture_y = 0,
82 int int_order = 2,
83 int mapsteps = 5,
84 int nslice = 1,
85 std::optional<std::string> name = std::nullopt
86 )
87 : Named(std::move(name)),
88 Thick(ds, nslice),
89 Alignment(dx, dy, rotation_degree),
91 m_k(k),m_unit(unit),m_int_order(int_order),m_mapsteps(mapsteps)
92 {
93 }
94
96 using BeamOptic::operator();
97
105 void compute_constants (RefPart const & refpart)
106 {
107 using namespace amrex::literals; // for _rt and _prt
108 using amrex::Math::powi;
109
110 Alignment::compute_constants(refpart);
111
112 // length of the current slice
113 m_slice_ds = m_ds / nslice();
114
115 // find beta*gamma^2
116 amrex::ParticleReal const pt_ref = refpart.pt;
117 m_betgam2 = powi<2>(pt_ref) - 1.0_prt;
118 m_ibetgam2 = 1_prt / m_betgam2;
119 m_beta = refpart.beta();
120 m_ibeta = 1_prt / m_beta;
121
122 // normalize quad units to MAD-X convention if needed
123 m_g = m_unit == 1 ? m_k / refpart.rigidity_Tm() : m_k;
124
125 // compute phase advance per unit length in s (in rad/m)
126 m_omega = std::sqrt(std::abs(m_g));
127
128 }
129
143 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
146 T_Real & AMREX_RESTRICT x,
147 T_Real & AMREX_RESTRICT y,
148 T_Real & AMREX_RESTRICT t,
149 T_Real & AMREX_RESTRICT px,
150 T_Real & AMREX_RESTRICT py,
151 T_Real & AMREX_RESTRICT pt,
152 T_IdCpu & AMREX_RESTRICT idcpu,
153 [[maybe_unused]] RefPart const & AMREX_RESTRICT refpart
154 ) const
155 {
156 using namespace amrex::literals; // for _rt and _prt
157
158 // shift due to alignment errors of the element
159 shift_in(x, y, px, py);
160
161 // coordinates within the quad (not needed)
162 //amrex::ParticleReal const s = refpart.s;
163 //amrex::ParticleReal const sedge = refpart.sedge;
164
165 // numerical integration parameters
166 amrex::ParticleReal const zin = 0_prt;
167 amrex::ParticleReal const zout = m_slice_ds;
168 int const nsteps = m_mapsteps;
169
170 // initialize phase space 6-vector
172 x, px, y, py, t, pt
173 };
174
175 // call integrator to advance through slice (int_order = 2 or 4, otherwise default 2)
176 if (m_int_order == 2) {
177 integrators::symp2_integrate_particle(particle,zin,zout,nsteps,*this);
178 } else if (m_int_order == 4) {
179 integrators::symp4_integrate_particle(particle,zin,zout,nsteps,*this);
180 } else if (m_int_order == 6) {
181 integrators::symp6_integrate_particle(particle,zin,zout,nsteps,*this);
182 } else {
183 integrators::symp2_integrate_particle(particle,zin,zout,nsteps,*this);
184 }
185
186 // assign updated values
187 x = particle(1);
188 px = particle(2);
189 y = particle(3);
190 py = particle(4);
191 t = particle(5);
192 pt = particle(6);
193
194 // apply transverse aperture
195 apply_aperture(x, y, idcpu);
196
197 // undo shift due to alignment errors of the element
198 shift_out(x, y, px, py);
199 }
200
210 template<typename T_Real>
212 void map1 (amrex::ParticleReal const tau,
214 amrex::ParticleReal & zeval) const
215 {
216 using namespace amrex::literals; // for _rt and _prt
217
218 T_Real const x = particle(1);
219 T_Real const px = particle(2);
220 T_Real const y = particle(3);
221 T_Real const py = particle(4);
222 T_Real const t = particle(5);
223 T_Real const pt = particle(6);
224
225 T_Real xout = x;
226 T_Real pxout = px;
227 T_Real yout = y;
228 T_Real pyout = py;
229 T_Real tout = t;
230 T_Real ptout = pt;
231
232 amrex::ParticleReal const sin_omega_ds = sin(m_omega*tau);
233 amrex::ParticleReal const cos_omega_ds = cos(m_omega*tau);
234 amrex::ParticleReal const sinh_omega_ds = sinh(m_omega*tau);
235 amrex::ParticleReal const cosh_omega_ds = cosh(m_omega*tau);
236 amrex::ParticleReal const slice_bg = tau / m_betgam2;
237
238 if (m_g > 0.0_prt)
239 {
240 // advance position and momentum (focusing quad)
241 xout = cos_omega_ds*x + sin_omega_ds/m_omega*px;
242 pxout = -m_omega*sin_omega_ds*x + cos_omega_ds*px;
243
244 yout = cosh_omega_ds*y + sinh_omega_ds/m_omega*py;
245 pyout = m_omega*sinh_omega_ds*y + cosh_omega_ds*py;
246
247 tout = t + (slice_bg)*pt;
248 // ptout = pt;
249 } else if (m_g < 0.0_prt)
250 {
251 // advance position and momentum (defocusing quad)
252 xout = cosh_omega_ds*x + sinh_omega_ds/m_omega*px;
253 pxout = m_omega*sinh_omega_ds*x + cosh_omega_ds*px;
254
255 yout = cos_omega_ds*y + sin_omega_ds/m_omega*py;
256 pyout = -m_omega*sin_omega_ds*y + cos_omega_ds*py;
257
258 tout = t + slice_bg*pt;
259 // ptout = pt;
260 } else {
261 // advance position and momentum (zero strength = drift)
262 xout = x + tau * px;
263 // pxout = px;
264 yout = y + tau * py;
265 // pyout = py;
266 tout = t + slice_bg * pt;
267 // ptout = pt;
268 }
269
270 // push particle coordinates
271 particle(1) = xout;
272 particle(2) = pxout;
273 particle(3) = yout;
274 particle(4) = pyout;
275 particle(5) = tout;
276 particle(6) = ptout;
277
278 zeval += 0_prt;
279 }
280
289 template<typename T_Real>
291 void map2 (amrex::ParticleReal const tau,
293 amrex::ParticleReal & zeval) const
294 {
295 using namespace amrex::literals; // for _rt and _prt
296 using namespace std; // for cmath(float)
297 using amrex::Math::powi;
298
299 T_Real const x = particle(1);
300 T_Real const px = particle(2);
301 T_Real const y = particle(3);
302 T_Real const py = particle(4);
303 T_Real const t = particle(5);
304 T_Real const pt = particle(6);
305
306 // compute the radical in the denominator (= pz):
307 T_Real const inv_pzden = 1_prt / sqrt(
308 powi<2>(pt - m_ibeta) -
309 m_ibetgam2 -
310 powi<2>(px) -
311 powi<2>(py)
312 );
313
314 // The momenta are not affected.
315 particle(1) = x + tau * px * (-1_prt + inv_pzden);
316 particle(2) = px;
317 particle(3) = y + tau * py * (-1_prt + inv_pzden);
318 particle(4) = py;
319 particle(5) = t + tau * (-m_ibeta - m_ibetgam2*pt + inv_pzden*(1_prt - pt*m_beta));
320 particle(6) = pt;
321
322 zeval += tau;
323 }
324
325
331 void operator() (RefPart & AMREX_RESTRICT refpart) const { // TODO: update as well, but needs more careful placement of calc_constants
332
333 using namespace amrex::literals; // for _rt and _prt
334 using amrex::Math::powi;
335
336 // assign input reference particle values
337 amrex::ParticleReal const x = refpart.x;
338 amrex::ParticleReal const px = refpart.px;
339 amrex::ParticleReal const y = refpart.y;
340 amrex::ParticleReal const py = refpart.py;
341 amrex::ParticleReal const z = refpart.z;
342 amrex::ParticleReal const pz = refpart.pz;
343 amrex::ParticleReal const t = refpart.t;
344 amrex::ParticleReal const pt = refpart.pt;
345 amrex::ParticleReal const s = refpart.s;
346
347 // length of the current slice
348 amrex::ParticleReal const slice_ds = m_ds / nslice();
349
350 // assign intermediate parameter
351 amrex::ParticleReal const step = slice_ds / std::sqrt(powi<2>(pt)-1.0_prt);
352
353 // advance position and momentum (straight element)
354 refpart.x = x + step*px;
355 refpart.y = y + step*py;
356 refpart.z = z + step*pz;
357 refpart.t = t - step*pt;
358
359 // advance integrated path length
360 refpart.s = s + slice_ds;
361 }
362
368 Map6x6
369 transport_map (RefPart const & AMREX_RESTRICT refpart) const // TODO: update as well, but needs more careful placement of calc_constants
370 {
371 using namespace amrex::literals; // for _rt and _prt
372 using amrex::Math::powi;
373
374 // length of the current slice
375 amrex::ParticleReal const slice_ds = m_ds / nslice();
376
377 // access reference particle values to find beta*gamma^2
378 amrex::ParticleReal const pt_ref = refpart.pt;
379 amrex::ParticleReal const betgam2 = powi<2>(pt_ref) - 1.0_prt;
380
381 // compute phase advance per unit length in s (in rad/m)
382 amrex::ParticleReal const omega = std::sqrt(std::abs(m_k));
383
384 // initialize linear map matrix elements
386
387 if (m_k > 0.0) {
388 R(1,1) = std::cos(omega*slice_ds);
389 R(1,2) = std::sin(omega*slice_ds)/omega;
390 R(2,1) = -omega*std::sin(omega*slice_ds);
391 R(2,2) = std::cos(omega*slice_ds);
392 R(3,3) = std::cosh(omega*slice_ds);
393 R(3,4) = std::sinh(omega*slice_ds)/omega;
394 R(4,3) = omega*std::sinh(omega*slice_ds);
395 R(4,4) = std::cosh(omega*slice_ds);
396 R(5,6) = slice_ds/betgam2;
397 } else if (m_k < 0.0) {
398 R(1,1) = std::cosh(omega*slice_ds);
399 R(1,2) = std::sinh(omega*slice_ds)/omega;
400 R(2,1) = omega*std::sinh(omega*slice_ds);
401 R(2,2) = std::cosh(omega*slice_ds);
402 R(3,3) = std::cos(omega*slice_ds);
403 R(3,4) = std::sin(omega*slice_ds)/omega;
404 R(4,3) = -omega*std::sin(omega*slice_ds);
405 R(4,4) = std::cos(omega*slice_ds);
406 R(5,6) = slice_ds/betgam2;
407 } else {
408 R(1,2) = m_slice_ds;
409 R(3,4) = m_slice_ds;
410 R(5,6) = m_slice_ds / betgam2;
411 }
412
413 return R;
414 }
415
417 using LinearTransport::operator();
418
419 amrex::ParticleReal m_k;
420 int m_unit;
423
424 private:
425 // constants that are independent of the individually tracked particle,
426 // see: compute_constants() to refresh
427 amrex::ParticleReal m_slice_ds;
428 amrex::ParticleReal m_betgam2;
429 amrex::ParticleReal m_ibetgam2;
430 amrex::ParticleReal m_beta;
431 amrex::ParticleReal m_ibeta;
432 amrex::ParticleReal m_g;
433 amrex::ParticleReal m_omega;
434 };
435
436} // namespace impactx
437
438#endif // IMPACTX_EXACTQUAD_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
#define AMREX_GPU_HOST
constexpr T powi(T x) noexcept
SmallMatrix< T, N, 1, Order::F, StartIndex > SmallVector
Definition All.H:54
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void symp4_integrate_particle(amrex::SmallVector< T_Real, 6, 1 > &particle, amrex::ParticleReal const zin, amrex::ParticleReal const zout, int const nsteps, T_Element const &element)
Definition Integrators.H:225
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void symp6_integrate_particle(amrex::SmallVector< T_Real, 6, 1 > &particle, amrex::ParticleReal const zin, amrex::ParticleReal const zout, int const nsteps, T_Element const &element)
Definition Integrators.H:284
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void symp2_integrate_particle(amrex::SmallVector< T_Real, 6, 1 > &particle, amrex::ParticleReal const zin, amrex::ParticleReal const zout, int const nsteps, T_Element const &element)
Definition Integrators.H:178
@ 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
static constexpr __host__ __device__ SmallMatrix< T, NRows, NCols, ORDER, StartIndex > Identity() noexcept
Definition ReferenceParticle.H:31
amrex::ParticleReal pt
energy, normalized by rest energy
Definition ReferenceParticle.H:40
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
void compute_constants(RefPart const &refpart)
Definition ExactQuad.H:105
ExactQuad(amrex::ParticleReal ds, amrex::ParticleReal k, int unit, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, amrex::ParticleReal aperture_x=0, amrex::ParticleReal aperture_y=0, int int_order=2, int mapsteps=5, int nslice=1, std::optional< std::string > name=std::nullopt)
Definition ExactQuad.H:73
amrex::ParticleReal m_ibetgam2
beta*gamma^2
Definition ExactQuad.H:429
amrex::ParticleReal m_slice_ds
number of integration steps per slice
Definition ExactQuad.H:427
int m_int_order
unit specification for quad strength
Definition ExactQuad.H:421
int m_mapsteps
order used for the symplectic integration (2 or 4)
Definition ExactQuad.H:422
amrex::ParticleReal m_omega
quadrupole strength in 1/m^2
Definition ExactQuad.H:433
amrex::ParticleReal m_ibeta
beta
Definition ExactQuad.H:431
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void map1(amrex::ParticleReal const tau, amrex::SmallVector< T_Real, 6, 1 > &particle, amrex::ParticleReal &zeval) const
Definition ExactQuad.H:212
amrex::ParticleReal m_beta
1 / m_betgam2
Definition ExactQuad.H:430
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void map2(amrex::ParticleReal const tau, amrex::SmallVector< T_Real, 6, 1 > &particle, amrex::ParticleReal &zeval) const
Definition ExactQuad.H:291
amrex::ParticleReal m_g
1 / m_beta
Definition ExactQuad.H:432
ImpactXParticleContainer::ParticleType PType
Definition ExactQuad.H:45
int m_unit
quadrupole strength in 1/m^2 (or T/m)
Definition ExactQuad.H:420
amrex::ParticleReal m_k
Definition ExactQuad.H:419
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition ExactQuad.H:369
amrex::ParticleReal m_betgam2
m_ds / nslice();
Definition ExactQuad.H:428
static constexpr auto type
Definition ExactQuad.H:44
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 ExactQuad.H:145
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