ImpactX
Loading...
Searching...
No Matches
Quad.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, Kyrre Ness Sjobak
8 * License: BSD-3-Clause-LBNL
9 */
10#ifndef IMPACTX_QUAD_H
11#define IMPACTX_QUAD_H
12
14#include "mixin/alignment.H"
15#include "mixin/pipeaperture.H"
16#include "mixin/beamoptic.H"
17#include "mixin/thick.H"
18#include "mixin/named.H"
19#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
29
30namespace impactx::elements
31{
32 struct Quad
33 : public mixin::Named,
34 public mixin::BeamOptic<Quad>,
35 public mixin::LinearTransport<Quad>,
36 public mixin::Thick,
37 public mixin::Alignment,
40 // At least on Intel AVX512, there is a small overhead to vectorize this element for DP and a gain for SP, see
41 // https://github.com/BLAST-ImpactX/impactx/pull/1002
42#ifdef AMREX_SINGLE_PRECISION_PARTICLES
43 , public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
44#endif
45 {
46 static constexpr auto type = "Quad";
48
65 amrex::ParticleReal ds,
66 amrex::ParticleReal k,
67 amrex::ParticleReal dx = 0,
68 amrex::ParticleReal dy = 0,
69 amrex::ParticleReal rotation_degree = 0,
70 amrex::ParticleReal aperture_x = 0,
71 amrex::ParticleReal aperture_y = 0,
72 int nslice = 1,
73 std::optional<std::string> name = std::nullopt
74 )
75 : Named(std::move(name)),
76 Thick(ds, nslice),
77 Alignment(dx, dy, rotation_degree),
79 m_k(k)
80 {
81 }
82
84 using BeamOptic::operator();
85
93 void compute_constants (RefPart const & refpart)
94 {
95 using namespace amrex::literals; // for _rt and _prt
97
98 Alignment::compute_constants(refpart);
99
100 // length of the current slice
101 m_slice_ds = m_ds / nslice();
102
103 amrex::ParticleReal const pt_ref = refpart.pt;
104 // find beta*gamma^2
105 m_betgam2 = powi<2>(pt_ref) - 1.0_prt;
107
108 // compute phase advance per unit length in s (in rad/m)
109 m_omega = std::sqrt(std::abs(m_k));
114 }
115
129 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
132 T_Real & AMREX_RESTRICT x,
133 T_Real & AMREX_RESTRICT y,
134 T_Real & AMREX_RESTRICT t,
135 T_Real & AMREX_RESTRICT px,
136 T_Real & AMREX_RESTRICT py,
137 T_Real & AMREX_RESTRICT pt,
138 T_IdCpu & AMREX_RESTRICT idcpu,
139 [[maybe_unused]] RefPart const & AMREX_RESTRICT refpart
140 ) const
141 {
142 using namespace amrex::literals; // for _rt and _prt
143
144 // shift due to alignment errors of the element
145 shift_in(x, y, px, py);
146
147 // initialize output values
148 T_Real xout = x;
149 T_Real yout = y;
150 T_Real tout = t;
151 T_Real pxout = px;
152 T_Real pyout = py;
153 T_Real const ptout = pt;
154
155 if (m_k > 0.0_prt)
156 {
157 // advance position and momentum (focusing quad)
159 pxout = -m_omega*m_sin_omega_ds*x + m_cos_omega_ds*px;
160
163
164 tout = t + (m_slice_bg)*pt;
165 // ptout = pt;
166 } else if (m_k < 0.0_prt)
167 {
168 // advance position and momentum (defocusing quad)
171
173 pyout = -m_omega*m_sin_omega_ds*y + m_cos_omega_ds*py;
174
175 tout = t + m_slice_bg*pt;
176 // ptout = pt;
177 } else {
178 // advance position and momentum (zero strength = drift)
179 xout = x + m_slice_ds * px;
180 // pxout = px;
181 yout = y + m_slice_ds * py;
182 // pyout = py;
183 tout = t + m_slice_bg * pt;
184 // ptout = pt;
185 }
186
187 // assign updated values
188 x = xout;
189 y = yout;
190 t = tout;
191 px = pxout;
192 py = pyout;
193 pt = ptout;
194
195 // apply transverse aperture
196 apply_aperture(x, y, idcpu);
197
198 // undo shift due to alignment errors of the element
199 shift_out(x, y, px, py);
200 }
201
207 void operator() (RefPart & AMREX_RESTRICT refpart) const { // TODO: update as well, but needs more careful placement of calc_constants
208
209 using namespace amrex::literals; // for _rt and _prt
210 using amrex::Math::powi;
211
212 // assign input reference particle values
213 amrex::ParticleReal const x = refpart.x;
214 amrex::ParticleReal const px = refpart.px;
215 amrex::ParticleReal const y = refpart.y;
216 amrex::ParticleReal const py = refpart.py;
217 amrex::ParticleReal const z = refpart.z;
218 amrex::ParticleReal const pz = refpart.pz;
219 amrex::ParticleReal const t = refpart.t;
220 amrex::ParticleReal const pt = refpart.pt;
221 amrex::ParticleReal const s = refpart.s;
222
223 // length of the current slice
224 amrex::ParticleReal const slice_ds = m_ds / nslice();
225
226 // assign intermediate parameter
227 amrex::ParticleReal const step = slice_ds / std::sqrt(powi<2>(pt)-1.0_prt);
228
229 // advance position and momentum (straight element)
230 refpart.x = x + step*px;
231 refpart.y = y + step*py;
232 refpart.z = z + step*pz;
233 refpart.t = t - step*pt;
234
235 // advance integrated path length
236 refpart.s = s + slice_ds;
237 }
238
244 Map6x6
245 transport_map (RefPart const & AMREX_RESTRICT refpart) const // TODO: update as well, but needs more careful placement of calc_constants
246 {
247 using namespace amrex::literals; // for _rt and _prt
248 using amrex::Math::powi;
249
250 // length of the current slice
251 amrex::ParticleReal const slice_ds = m_ds / nslice();
252
253 // access reference particle values to find beta*gamma^2
254 amrex::ParticleReal const pt_ref = refpart.pt;
255 amrex::ParticleReal const betgam2 = powi<2>(pt_ref) - 1.0_prt;
256
257 // compute phase advance per unit length in s (in rad/m)
258 amrex::ParticleReal const omega = std::sqrt(std::abs(m_k));
259
260 // initialize linear map matrix elements
262
263 if (m_k > 0.0) {
264 R(1,1) = std::cos(omega*slice_ds);
265 R(1,2) = std::sin(omega*slice_ds)/omega;
266 R(2,1) = -omega*std::sin(omega*slice_ds);
267 R(2,2) = std::cos(omega*slice_ds);
268 R(3,3) = std::cosh(omega*slice_ds);
269 R(3,4) = std::sinh(omega*slice_ds)/omega;
270 R(4,3) = omega*std::sinh(omega*slice_ds);
271 R(4,4) = std::cosh(omega*slice_ds);
272 R(5,6) = slice_ds/betgam2;
273 } else if (m_k < 0.0) {
274 R(1,1) = std::cosh(omega*slice_ds);
275 R(1,2) = std::sinh(omega*slice_ds)/omega;
276 R(2,1) = omega*std::sinh(omega*slice_ds);
277 R(2,2) = std::cosh(omega*slice_ds);
278 R(3,3) = std::cos(omega*slice_ds);
279 R(3,4) = std::sin(omega*slice_ds)/omega;
280 R(4,3) = -omega*std::sin(omega*slice_ds);
281 R(4,4) = std::cos(omega*slice_ds);
282 R(5,6) = slice_ds/betgam2;
283 } else {
284 R(1,2) = slice_ds;
285 R(3,4) = slice_ds;
286 R(5,6) = slice_ds / betgam2;
287 }
288
289 return R;
290 }
291
293 using LinearTransport::operator();
294
295 amrex::ParticleReal m_k;
296
297 private:
298 // constants that are independent of the individually tracked particle,
299 // see: compute_constants() to refresh
300 amrex::ParticleReal m_slice_ds;
301 amrex::ParticleReal m_betgam2;
302 amrex::ParticleReal m_slice_bg;
303 amrex::ParticleReal m_omega;
304 amrex::ParticleReal m_sin_omega_ds;
305 amrex::ParticleReal m_cos_omega_ds;
306 amrex::ParticleReal m_sinh_omega_ds;
307 amrex::ParticleReal m_cosh_omega_ds;
308 };
309
310} // namespace impactx
311
312#endif // IMPACTX_QUAD_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
#define AMREX_GPU_HOST
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
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::ParticleReal m_cosh_omega_ds
std::sinh(omega*slice_ds)
Definition Quad.H:307
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 Quad.H:131
amrex::ParticleReal m_k
Definition Quad.H:295
static constexpr auto type
Definition Quad.H:46
amrex::ParticleReal m_slice_ds
quadrupole strength in 1/m
Definition Quad.H:300
amrex::ParticleReal m_cos_omega_ds
std::sin(omega*slice_ds)
Definition Quad.H:305
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition Quad.H:245
void compute_constants(RefPart const &refpart)
Definition Quad.H:93
amrex::ParticleReal m_slice_bg
beta*gamma^2
Definition Quad.H:302
ImpactXParticleContainer::ParticleType PType
Definition Quad.H:47
amrex::ParticleReal m_sin_omega_ds
std::sqrt(std::abs(m_k)) compute phase advance per unit length in s (in rad/m)
Definition Quad.H:304
Quad(amrex::ParticleReal ds, amrex::ParticleReal k, 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 Quad.H:64
amrex::ParticleReal m_sinh_omega_ds
std::cos(omega*slice_ds)
Definition Quad.H:306
amrex::ParticleReal m_betgam2
m_ds / nslice();
Definition Quad.H:301
amrex::ParticleReal m_omega
m_slice_ds / m_betgam2
Definition Quad.H:303
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