ImpactX
Loading...
Searching...
No Matches
ChrQuad.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_CHRQUAD_H
11#define IMPACTX_CHRQUAD_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/spintransport.H"
20#include "mixin/named.H"
21#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 ChrQuad
35 : public mixin::Named,
36 public mixin::BeamOptic<ChrQuad>,
37 public mixin::LinearTransport<ChrQuad>,
38 public mixin::Thick,
39 public mixin::Alignment,
42 public mixin::NoFinalize,
43 public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
44 {
45 static constexpr auto type = "ChrQuad";
47
73 int unit,
76 amrex::ParticleReal rotation_degree = 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_k(k), m_unit(unit)
87 {
88 }
89
91 void reverse () { Thick::reverse(); }
92
94 using BeamOptic::operator();
95
96
104 void compute_constants (RefPart const & refpart)
105 {
106 using namespace amrex::literals; // for _rt and _prt
107 using amrex::Math::powi;
108
109 Alignment::compute_constants(refpart);
110
111 // length of the current slice
112 m_slice_ds = m_ds / nslice();
113
114 // access reference particle values to find beta and gamma
115 m_beta = refpart.beta();
116 amrex::ParticleReal const gamma = refpart.gamma();
117
118 // normalize quad units to MAD-X convention if needed
119 m_g = m_unit == 1 ? m_k / refpart.rigidity_Tm() : m_k;
120
121 m_const1 = 1_prt / (2_prt * powi<3>(m_beta) * powi<2>(gamma));
122 }
123
135 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
138 T_Real & AMREX_RESTRICT x,
139 T_Real & AMREX_RESTRICT y,
140 T_Real & AMREX_RESTRICT t,
141 T_Real & AMREX_RESTRICT px,
142 T_Real & AMREX_RESTRICT py,
143 T_Real const & AMREX_RESTRICT pt,
144 [[maybe_unused]] T_IdCpu const & AMREX_RESTRICT idcpu,
145 RefPart const & AMREX_RESTRICT refpart
146 ) const
147 {
148 using namespace amrex::literals; // for _rt and _prt
149 using amrex::Math::powi;
150 using namespace std; // for cmath(float)
151
152 // access position data
153 T_Real const xout = x;
154 T_Real const yout = y;
155 T_Real const tout = t;
156
157 // compute particle momentum deviation delta + 1
158 T_Real const delta1 = sqrt(1_prt - 2_prt*pt/m_beta + powi<2>(pt));
159 T_Real const delta = delta1 - 1_prt;
160
161 // compute phase advance per unit length in s (in rad/m)
162 // chromatic dependence on delta is included
163 T_Real const omega = sqrt(abs(m_g)/delta1);
164
165 // initialize output values of momenta
166 T_Real pxout = px;
167 T_Real pyout = py;
168 // T_Real const ptout = pt;
169
170 // paceholder variables
171 T_Real q1 = xout;
172 T_Real q2 = yout;
173 T_Real p1 = px;
174 T_Real p2 = py;
175
176 if (m_g > 0_prt)
177 {
178 // advance transverse position and momentum (focusing quad)
179 x = cos(omega*m_slice_ds) * xout +
180 sin(omega*m_slice_ds) / (omega*delta1)*px;
181 pxout = -omega * delta1 * sin(omega*m_slice_ds) * xout + cos(omega * m_slice_ds) * px;
182
183 y = cosh(omega*m_slice_ds) * yout +
184 sinh(omega*m_slice_ds) / (omega*delta1)*py;
185 pyout = omega * delta1 * sinh(omega*m_slice_ds) * yout + cosh(omega * m_slice_ds) * py;
186
187 } else if (m_g < 0_prt)
188 {
189 // advance transverse position and momentum (defocusing quad)
190 x = cosh(omega*m_slice_ds) * xout +
191 sinh(omega*m_slice_ds) / (omega*delta1)*px;
192 pxout = omega * delta1 * sinh(omega*m_slice_ds) * xout + cosh(omega * m_slice_ds) * px;
193
194 y = cos(omega*m_slice_ds) * yout +
195 sin(omega*m_slice_ds) / (omega*delta1) * py;
196 pyout = -omega * delta1 * sin(omega*m_slice_ds) * yout + cos(omega * m_slice_ds) * py;
197
198 q1 = yout;
199 q2 = xout;
200 p1 = py;
201 p2 = px;
202 } else
203 {
204 // advance transverse position and momentum (zero focusing strength = drift)
205 x = xout + m_slice_ds * px / delta1;
206 // pxout = px;
207 y = yout + m_slice_ds * py / delta1;
208 // pyout = py;
209 }
210
211
212 // advance longitudinal position and momentum
213
214 if (m_g == 0.0)
215 {
216 // the corresponding symplectic update to t (zero strength = drift)
217 T_Real term = 2_prt * powi<2>(pt) + powi<2>(px) + powi<2>(py);
218 term = 2_prt - 4_prt * m_beta * pt + powi<2>(m_beta) * term;
219 term = -2_prt + powi<2>(refpart.gamma())*term;
220 term = (-1_prt + m_beta * pt) * term;
221 term = term * m_const1;
222 t = tout - m_slice_ds * (1_prt / m_beta + term / powi<3>(delta1));
223 // ptout = pt;
224
225 } else
226 {
227 // the corresponding symplectic update to t (nonzero strength)
228 T_Real const term = pt + delta / m_beta;
229 T_Real const t0 = tout - term * m_slice_ds / delta1;
230
231 T_Real const w = omega * delta1;
232 T_Real const term1 = -(powi<2>(p2) + powi<2>(q2) * powi<2>(w)) * sinh(2_prt*m_slice_ds*omega);
233 T_Real const term2 = -(powi<2>(p1) - powi<2>(q1) * powi<2>(w)) * sin(2_prt*m_slice_ds*omega);
234 T_Real const term3 = -2_prt * q2 * p2 * w * cosh(2_prt*m_slice_ds*omega);
235 T_Real const term4 = -2_prt * q1 * p1 * w * cos(2_prt*m_slice_ds*omega);
236 T_Real const term5 = 2_prt * omega * (q1*p1*delta1 + q2*p2*delta1
237 -(powi<2>(p1) + powi<2>(p2))*m_slice_ds - (powi<2>(q1)-powi<2>(q2)) * powi<2>(w)*m_slice_ds);
238 t = t0 + (-1_prt+m_beta*pt)
239 / (8_prt*m_beta * powi<3>(delta1)*omega)
240 * (term1+term2+term3+term4+term5);
241 // ptout = pt;
242 }
243
244 // assign updated momenta
245 px = pxout;
246 py = pyout;
247 // pt = ptout;
248 }
249
255 void operator() (RefPart & AMREX_RESTRICT refpart) const
256 {
257 using namespace amrex::literals; // for _rt and _prt
258 using amrex::Math::powi;
259
260 // assign input reference particle values
261 amrex::ParticleReal const x = refpart.x;
262 amrex::ParticleReal const px = refpart.px;
263 amrex::ParticleReal const y = refpart.y;
264 amrex::ParticleReal const py = refpart.py;
265 amrex::ParticleReal const z = refpart.z;
266 amrex::ParticleReal const pz = refpart.pz;
267 amrex::ParticleReal const t = refpart.t;
268 amrex::ParticleReal const pt = refpart.pt;
269 amrex::ParticleReal const s = refpart.s;
270
271 // length of the current slice
272 amrex::ParticleReal const slice_ds = m_ds / nslice();
273
274 // assign intermediate parameter
275 amrex::ParticleReal const step = slice_ds / std::sqrt(powi<2>(pt)-1.0_prt);
276
277 // advance position and momentum (straight element)
278 refpart.x = x + step*px;
279 refpart.y = y + step*py;
280 refpart.z = z + step*pz;
281 refpart.t = t - step*pt;
282
283 // advance integrated path length
284 refpart.s = s + slice_ds;
285 }
286
287
288
303 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
306 T_Real & AMREX_RESTRICT x,
307 T_Real & AMREX_RESTRICT y,
308 T_Real & AMREX_RESTRICT t,
309 T_Real & AMREX_RESTRICT px,
310 T_Real & AMREX_RESTRICT py,
311 T_Real const & AMREX_RESTRICT pt,
312 T_Real & AMREX_RESTRICT sx,
313 T_Real & AMREX_RESTRICT sy,
314 T_Real & AMREX_RESTRICT sz,
315 T_IdCpu const & AMREX_RESTRICT idcpu,
316 RefPart const & AMREX_RESTRICT refpart
317 ) const
318 {
319 using namespace amrex::literals; // for _rt and _prt
320 using amrex::Math::powi;
321 using namespace std; // for cmath(float)
322
323 // initialize the three components of the axis-angle vector
324 T_Real lambdax = 0_prt;
325 T_Real lambday = 0_prt;
326 T_Real lambdaz = 0_prt;
327
328 // compute particle momentum deviation delta + 1
329 T_Real const delta1 = sqrt(1_prt - 2_prt*pt/m_beta + powi<2>(pt));
330 T_Real const inv_delta1 = 1_prt / delta1;
331
332 // compute particle relativistic gamma
333 T_Real const gamma = refpart.gamma() * (1_prt - m_beta*pt);
334 T_Real const gyro_const = 1_prt + refpart.gyromagnetic_anomaly * gamma;
335
336 // compute phase advance per unit length in s (in rad/m)
337 // chromatic dependence on delta is included
338 T_Real const omega = sqrt(abs(m_g)*inv_delta1);
339
340 // compute trigonometric quantities
341 auto const [sin_omega_ds, cos_omega_ds] = amrex::Math::sincos(omega*m_slice_ds);
342 T_Real const sinh_omega_ds = sinh(omega*m_slice_ds);
343 T_Real const cosh_omega_ds = cosh(omega*m_slice_ds);
344
345 if (m_g > 0.0_prt)
346 {
347
348 // horizontally focusing quad case
349 lambdax = -gyro_const * (
350 py * inv_delta1 * (cosh_omega_ds - 1_prt) +
351 y * omega * sinh_omega_ds
352 );
353 lambday = -gyro_const * (
354 px * inv_delta1 * (1_prt - cos_omega_ds) +
355 x * omega * sin_omega_ds
356 );
357
358 } else if (m_g < 0.0_prt)
359 {
360
361 // horizontally defocusing quad case
362 lambdax = gyro_const * (
363 py * inv_delta1 * (1_prt - cos_omega_ds) +
364 y * omega * sin_omega_ds
365 );
366 lambday = gyro_const * (
367 px * inv_delta1 * (cosh_omega_ds - 1_prt) +
368 x * omega * sinh_omega_ds
369 );
370
371 } else {
372 // treat as a drift
373 }
374
375 // push the spin vector using the generator just determined
376 rotate_spin(lambdax,lambday,lambdaz,sx,sy,sz);
377
378 // phase space push
379 (*this)(x, y, t, px, py, pt, idcpu, refpart);
380
381 }
382
383
385 using LinearTransport::operator();
386
393 Map6x6
394 transport_map (RefPart const & AMREX_RESTRICT refpart) const
395 {
396 using namespace amrex::literals; // for _rt and _prt
397 using amrex::Math::powi;
398
399 // length of the current slice
400 amrex::ParticleReal const slice_ds = m_ds / nslice();
401
402 // access reference particle values to find beta*gamma^2
403 amrex::ParticleReal const pt_ref = refpart.pt;
404 amrex::ParticleReal const betgam2 = powi<2>(pt_ref) - 1_prt;
405
406 // normalize quad units to MAD-X convention if needed
407 amrex::ParticleReal const g = m_unit == 1 ? m_k / refpart.rigidity_Tm() : m_k;
408
409 // compute phase advance per unit length in s (in rad/m)
410 amrex::ParticleReal const omega = std::sqrt(std::abs(g));
411
412 // initialize linear map matrix elements
414
415 if (g > 0.0) {
416 R(1,1) = std::cos(omega*slice_ds);
417 R(1,2) = std::sin(omega*slice_ds)/omega;
418 R(2,1) = -omega*std::sin(omega*slice_ds);
419 R(2,2) = std::cos(omega*slice_ds);
420 R(3,3) = std::cosh(omega*slice_ds);
421 R(3,4) = std::sinh(omega*slice_ds)/omega;
422 R(4,3) = omega*std::sinh(omega*slice_ds);
423 R(4,4) = std::cosh(omega*slice_ds);
424 R(5,6) = slice_ds/betgam2;
425 } else if (g < 0.0) {
426 R(1,1) = std::cosh(omega*slice_ds);
427 R(1,2) = std::sinh(omega*slice_ds)/omega;
428 R(2,1) = omega*std::sinh(omega*slice_ds);
429 R(2,2) = std::cosh(omega*slice_ds);
430 R(3,3) = std::cos(omega*slice_ds);
431 R(3,4) = std::sin(omega*slice_ds)/omega;
432 R(4,3) = -omega*std::sin(omega*slice_ds);
433 R(4,4) = std::cos(omega*slice_ds);
434 R(5,6) = slice_ds/betgam2;
435 } else {
436 R(1,2) = m_slice_ds;
437 R(3,4) = m_slice_ds;
438 R(5,6) = m_slice_ds / betgam2;
439 }
440
441 // apply the transverse rotation (roll) alignment error
442 return rotate_aligned_map(R);
443 }
444
446 int m_unit;
447
448 private:
449 // constants that are independent of the individually tracked particle,
450 // see: compute_constants() to refresh
455 };
456
457} // namespace impactx
458
460
461#endif // IMPACTX_CHRQUAD_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
#define AMREX_GPU_HOST
#define IMPACTX_PUSH_EXTERN_TEMPLATE(ElementType)
Definition PushAll.H:78
amrex_particle_real ParticleReal
constexpr T powi(T x) noexcept
__host__ __device__ std::pair< double, double > sincos(double x)
__host__ __device__ T abs(const GpuComplex< T > &a_z) noexcept
__host__ __device__ GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
Definition All.H:56
@ 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:33
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal rigidity_Tm() const
Definition ReferenceParticle.H:260
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal beta() const
Definition ReferenceParticle.H:151
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal gamma() const
Definition ReferenceParticle.H:139
Definition ChrQuad.H:44
int m_unit
quadrupole strength in 1/m^2 (or T/m)
Definition ChrQuad.H:446
amrex::ParticleReal m_const1
Definition ChrQuad.H:454
amrex::ParticleReal m_slice_ds
unit specification for quad strength
Definition ChrQuad.H:451
ImpactXParticleContainer::ParticleType PType
Definition ChrQuad.H:46
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition ChrQuad.H:394
void compute_constants(RefPart const &refpart)
Definition ChrQuad.H:104
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 const &AMREX_RESTRICT pt, T_IdCpu const &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition ChrQuad.H:137
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void spin_and_phasespace_push(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 const &AMREX_RESTRICT pt, T_Real &AMREX_RESTRICT sx, T_Real &AMREX_RESTRICT sy, T_Real &AMREX_RESTRICT sz, T_IdCpu const &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition ChrQuad.H:305
ChrQuad(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 nslice=1, std::optional< std::string > name=std::nullopt)
Definition ChrQuad.H:70
amrex::ParticleReal m_k
Definition ChrQuad.H:445
static constexpr auto type
Definition ChrQuad.H:45
void reverse()
Definition ChrQuad.H:91
amrex::ParticleReal m_beta
m_ds / nslice();
Definition ChrQuad.H:452
amrex::ParticleReal m_g
Definition ChrQuad.H:453
Definition alignment.H:29
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dy() const
Definition alignment.H:189
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dx() const
Definition alignment.H:179
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 rotate_aligned_map(Map6x6 const &R) const
Definition alignment.H:263
Alignment(amrex::ParticleReal dx, amrex::ParticleReal dy, amrex::ParticleReal rotation_degree)
Definition alignment.H:39
Definition beamoptic.H:529
Definition lineartransport.H:50
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 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 spintransport.H:36
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void rotate_spin(T_Real const &AMREX_RESTRICT lambdax, T_Real const &AMREX_RESTRICT lambday, T_Real const &AMREX_RESTRICT lambdaz, T_Real &AMREX_RESTRICT sx, T_Real &AMREX_RESTRICT sy, T_Real &AMREX_RESTRICT sz) const
Definition spintransport.H:48
Definition thick.H:24
Thick(amrex::ParticleReal ds, int nslice)
Definition thick.H:30
amrex::ParticleReal m_ds
Definition thick.H:68
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