ImpactX
Loading...
Searching...
No Matches
ExactCFbend.H
Go to the documentation of this file.
1/* Copyright 2022-2026 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_EXACTCFBEND_H
11#define IMPACTX_EXACTCFBEND_H
12
15#include "mixin/alignment.H"
16#include "mixin/beamoptic.H"
17#include "mixin/dynamicdata.H"
19#include "mixin/named.H"
20#include "mixin/nofinalize.H"
21#include "mixin/pipeaperture.H"
22#include "mixin/thick.H"
23#include "mixin/TrackedVector.H"
24
25#include <AMReX_Extension.H>
26#include <AMReX_REAL.H>
27#include <AMReX_Print.H>
28#include <AMReX_GpuComplex.H>
29#include <AMReX_Math.H>
30
31#include <cmath>
32#include <memory>
33#include <stdexcept>
34#include <vector>
35
36namespace impactx::elements
37{
38
49
51 : public mixin::Named,
52 public mixin::BeamOptic<ExactCFbend>,
53 public mixin::LinearTransport<ExactCFbend>,
54 public mixin::Thick,
55 public mixin::Alignment,
56 public mixin::NoFinalize,
58 {
59 static constexpr auto type = "ExactCFbend";
61
63
91 std::vector<amrex::ParticleReal> k_normal,
92 std::vector<amrex::ParticleReal> k_skew,
93 int unit,
96 amrex::ParticleReal rotation_degree = 0,
99 int int_order = 2,
100 int mapsteps = 1,
101 int nslice = 1,
102 std::optional<std::string> name = std::nullopt
103 )
104 : Named(std::move(name)),
105 Thick(ds, nslice),
106 Alignment(dx, dy, rotation_degree),
108 m_unit(unit), m_int_order(int_order), m_mapsteps(mapsteps),
109 m_id(DynamicData::allocate_id())
110 {
111 if (m_int_order != 2 && m_int_order != 4 && m_int_order != 6)
112 throw std::runtime_error("ExactCFbend: The order used for symplectic integration must be 2, 4 or 6.");
113
114 m_ncoef = int(k_normal.size());
115 if (m_ncoef != int(k_skew.size()))
116 throw std::runtime_error("ExactCFbend: normal and skew coefficient arrays must have same length!");
117
118 auto& coef = DynamicData::emplace(
119 m_id,
120 std::move(k_normal),
121 std::move(k_skew)
122 );
123 m_k_normal_h_data = coef.k_normal.host_const().data();
124 m_k_skew_h_data = coef.k_skew.host_const().data();
125 }
126
128 void reverse () { Thick::reverse(); }
129
131 using BeamOptic::operator();
132
140 void compute_constants (RefPart const & refpart)
141 {
142 using namespace amrex::literals; // for _rt and _prt
143
144 Alignment::compute_constants(refpart);
145
146 // Ensure dynamic coefficient data is on GPU and we have fresh pointers to it
147 auto const & coef = *DynamicData::get(m_id);
148 m_k_normal_d_data = coef.k_normal.device_const().data();
149 m_k_skew_d_data = coef.k_skew.device_const().data();
150
151 // access reference particle values to find relativistic factors
152 m_beta = refpart.beta();
153 m_ibeta = 1_prt / m_beta;
154 m_ibetgam2 = 1_prt / amrex::Math::powi<2>(refpart.beta_gamma());
155 m_brho = refpart.rigidity_Tm();
156
157 amrex::ParticleReal const * k_normal = m_k_normal_h_data;
158
159 // find magnetic field and curvature; normalize units to MAD-X convention if needed
160 amrex::ParticleReal curvature = m_unit == 1 ? k_normal[0] / m_brho : k_normal[0];
161 m_rc = curvature == 0.0 ? 0.0_prt : 1.0_prt/curvature;
162
163 // arc length and angle of the current slice
164 m_slice_ds = m_ds / nslice();
165
166 }
167
189 uint64_t & AMREX_RESTRICT idcpu,
190 [[maybe_unused]] RefPart const & AMREX_RESTRICT refpart
191 ) const
192 {
193 using namespace amrex::literals; // for _rt and _prt
194
195 // shift due to alignment errors of the element
196 shift_in(x, y, px, py);
197
198 // numerical integration parameters
199 amrex::ParticleReal const zin = 0_prt;
200 amrex::ParticleReal const zout = m_slice_ds;
201 int const nsteps = m_mapsteps;
202
203 // initialize phase space 6-vector
205 x, px, y, py, t, pt
206 };
207
208 // call integrator to advance through slice (int_order = 2 or 4, otherwise default 2)
209 if (m_int_order == 2) {
210 integrators::symp2_integrate_particle(particle,zin,zout,nsteps,*this);
211 } else if (m_int_order == 4) {
212 integrators::symp4_integrate_particle(particle,zin,zout,nsteps,*this);
213 } else if (m_int_order == 6) {
214 integrators::symp6_integrate_particle(particle,zin,zout,nsteps,*this);
215 }
216
217 // assign updated values
218 x = particle(1);
219 px = particle(2);
220 y = particle(3);
221 py = particle(4);
222 t = particle(5);
223 pt = particle(6);
224
225 // apply transverse aperture
226 apply_aperture(x, y, idcpu);
227
228 // undo shift due to alignment errors of the element
229 shift_out(x, y, px, py);
230 }
231
241 void map1 (amrex::ParticleReal const tau,
243 amrex::ParticleReal & zeval) const
244 {
245 using namespace amrex::literals; // for _rt and _prt
246 using namespace amrex::Math; // for powi
247
248 amrex::ParticleReal const x = particle(1);
249 amrex::ParticleReal const px = particle(2);
250 amrex::ParticleReal const y = particle(3);
251 amrex::ParticleReal const py = particle(4);
252 amrex::ParticleReal const t = particle(5);
253 amrex::ParticleReal const pt = particle(6);
254
255 amrex::ParticleReal xout = x;
256 amrex::ParticleReal pxout = px;
257 amrex::ParticleReal yout = y;
258 amrex::ParticleReal pyout = py;
259 amrex::ParticleReal tout = t;
260 amrex::ParticleReal ptout = pt;
261
262 // treat the special case of zero field
263 if (m_rc==0_prt) {
264 // compute the radical in the denominator (= pz):
265 amrex::ParticleReal const inv_pzden = 1_prt / std::sqrt(
266 powi<2>(pt - m_ibeta) -
267 m_ibetgam2 -
268 powi<2>(px) -
269 powi<2>(py)
270 );
271
272 // advance position and momentum (exact drift)
273 xout = x + tau * px * inv_pzden;
274 yout = y + tau * py * inv_pzden;
275 tout = t - tau * (m_ibeta + (pt - m_ibeta) * inv_pzden);
276
277 // assign updated momenta
278 pxout = px;
279 pyout = py;
280 ptout = pt;
281
282 } else {
283
284 // assign intermediate quantities
285 amrex::ParticleReal const pperp2 = powi<2>(pt)-2.0_prt * m_ibeta * pt - powi<2>(py)+1.0_prt;
286 amrex::ParticleReal const px2 = powi<2>(px);
287
288 // trigonometric evaluations
289 amrex::ParticleReal const slice_phi = tau / m_rc;
290 auto const [sin_phi, cos_phi] = amrex::Math::sincos(slice_phi);
291
292 // determine if particle lies within the domain of map definition
293 if (pperp2 > px2)
294 {
295 amrex::ParticleReal const pperp = std::sqrt(pperp2);
296 amrex::ParticleReal const pzi = std::sqrt(powi<2>(pperp) - powi<2>(px));
297 amrex::ParticleReal const rho = m_rc + x;
298
299 // update momenta
300 pxout = px * cos_phi + (pzi - rho / m_rc) * sin_phi;
301 pyout = py;
302 ptout = pt;
303
304 // angle of momentum rotation
305 amrex::ParticleReal const px2f = powi<2>(pxout);
306 // determine if particle lies within domain of map definition
307 if (pperp2 > px2f)
308 {
309 amrex::ParticleReal const pzf = std::sqrt(powi<2>(pperp)-powi<2>(pxout));
310 amrex::ParticleReal const theta = slice_phi + std::asin(px/pperp) - std::asin(pxout/pperp);
311
312 // update position coordinates
313 xout = -m_rc + rho*cos_phi + m_rc*(pzf + px*sin_phi - pzi*cos_phi);
314 yout = y + theta * m_rc * py;
315 tout = t - theta * m_rc * (pt - 1.0_prt * m_ibeta) - slice_phi * m_rc * m_ibeta;
316
317 } // else { amrex::Print() << "Warning: outside map domain." << "\n"; }
318 } // else { amrex::Print() << "Warning: outside map domain." << "\n"; }
319 }
320
321 // push particle coordinates
322 particle(1) = xout;
323 particle(2) = pxout;
324 particle(3) = yout;
325 particle(4) = pyout;
326 particle(5) = tout;
327 particle(6) = ptout;
328
329 zeval += 0_prt;
330 }
331
342 void map2 (amrex::ParticleReal const tau,
344 amrex::ParticleReal & zeval) const
345 {
346 using namespace amrex::literals; // for _rt and _prt
347 using namespace amrex::Math; // for powi
348
349 amrex::ParticleReal const x = particle(1);
350 amrex::ParticleReal const px = particle(2);
351 amrex::ParticleReal const y = particle(3);
352 amrex::ParticleReal const py = particle(4);
353 amrex::ParticleReal const t = particle(5);
354 amrex::ParticleReal const pt = particle(6);
355
356 amrex::ParticleReal xout = x;
357 amrex::ParticleReal pxout = px;
358 amrex::ParticleReal yout = y;
359 amrex::ParticleReal pyout = py;
360 amrex::ParticleReal tout = t;
361 amrex::ParticleReal ptout = pt;
362
363 amrex::ParticleReal const * k_normal = m_k_normal_d_data;
364 amrex::ParticleReal const * k_skew = m_k_skew_d_data;
365
366 // defined dimensionless variables scaled by radius of curvature
367 amrex::ParticleReal rho = 1_prt + x / m_rc;
368 amrex::ParticleReal y_scl = y / m_rc;
369 amrex::ParticleReal ln_rho = std::log(rho);
370
371 // initialize the momentum kick
372 amrex::ParticleReal dpx = 0_prt;
373 amrex::ParticleReal dpy = 0_prt;
374
375 // loop over array of multipole coefficients
376 for (int m = 1; m < m_ncoef; m++) {
377
378 // normalize units to MAD-X convention if needed
379 amrex::ParticleReal kn = m_unit == 1 ? k_normal[m] / m_brho : k_normal[m];
380 amrex::ParticleReal ks = m_unit == 1 ? k_skew[m] / m_brho : k_skew[m];
381
382 int n = m+1; // multipole order
383
384 //amrex::Print() << "n, kn, ks before scaling: " << n << ", " << kn << ", " << ks << "\n";
385
386 // define dimensionless multipole coefficients (scaled by radius of curvature)
387 kn = std::pow(m_rc, m) * kn;
388 ks = std::pow(m_rc, m) * ks;
389
390 switch (n) {
391
392 case 1: // dipole-like contribution (not used - contributing to map 1)
393 dpx += -rho * kn;
394 dpy += ks;
395 break;
396 case 2: // quadrupole-like contribution
397 dpx += -rho*ln_rho*kn + y_scl*rho*ks;
398 dpy += y_scl*kn + 0.5_prt*(rho*rho - 1_prt)*ks;
399 break;
400 case 3: // sextupole-like contribution
401 dpx += 0.25_prt*rho*(1_prt + 2_prt*y_scl*y_scl - rho*rho + 2_prt*ln_rho)*kn + y_scl*rho*ln_rho*ks;
402 dpy += 0.5_prt*y_scl*(rho*rho - 1_prt)*kn + 0.25_prt*(1_prt - 2_prt*y_scl*y_scl - rho*rho + 2_prt*rho*rho*ln_rho)*ks;
403 break;
404 case 4: // octupole-like contribution
405 dpx += -0.25_prt*(rho - powi<3>(rho) + (rho - 2_prt*y_scl*y_scl*rho + powi<3>(rho))*ln_rho)*kn
406 -1_prt/12_prt*y_scl*rho*(3_prt + 2_prt*y_scl*y_scl - 3_prt*rho*rho + 6_prt*ln_rho)*ks;
407 dpy += -1_prt/12_prt*y_scl*(-3_prt + 2_prt*y_scl*y_scl + 3_prt*rho*rho - 6_prt*rho*rho*ln_rho)*kn
408 -1_prt/16_prt*((-1_prt + 4_prt*y_scl*y_scl - rho*rho)*(-1_prt + rho*rho) + 4_prt*rho*rho*ln_rho)*ks;
409 //amrex::Print() << "x, y, dpx, dpy = " << x << ", " << y << ", " << dpx << ", " << dpy << "\n";
410 break;
411 case 5: // decapole-like contribution
412 dpx += -1_prt/192_prt*rho*(8_prt*powi<4>(y_scl) - 24_prt*powi<2>(y_scl)*(-1_prt + rho*rho)
413 + 3_prt*(-5_prt + 4_prt*powi<2>(rho) + powi<4>(rho))
414 +12_prt*(-1_prt + 4_prt*y_scl*y_scl - 2_prt*rho*rho)*ln_rho)*kn
415 +1_prt/12_prt*y_scl*rho*(3_prt*(-1_prt + rho*rho) + (2_prt*y_scl*y_scl - 3_prt*(1_prt + rho*rho))*ln_rho)*ks;
416 dpy += -1_prt/48_prt*y_scl*((-1_prt + rho*rho)*(4_prt*y_scl*y_scl - 3_prt*(1_prt + rho*rho)) + 12_prt*rho*rho*ln_rho)*kn
417 + 1_prt/192_prt*(3_prt + 8_prt*powi<4>(y_scl) + 12_prt*rho*rho
418 - 15_prt*powi<4>(rho) + 24_prt*y_scl*y_scl*(-1_prt + rho*rho)
419 +12_prt*rho*rho*(2_prt - 4_prt*y_scl*y_scl + rho*rho)*ln_rho)*ks;
420 break;
421 // default:
422 // TODO: print a warning here that order > 5 is not currently supported
423 }
424
425 }
426
427 // advance momentum by a step of length tau
428 pxout = pxout + tau * dpx;
429 pyout = pyout + tau * dpy;
430
431 //update the particle coordinates
432 particle(1) = xout;
433 particle(2) = pxout;
434 particle(3) = yout;
435 particle(4) = pyout;
436 particle(5) = tout;
437 particle(6) = ptout;
438
439 zeval += tau;
440 }
441
442
448 void operator() (RefPart & AMREX_RESTRICT refpart) const
449 {
450 using namespace amrex::literals; // for _rt and _prt
451
452 // assign input reference particle values
453 amrex::ParticleReal const x = refpart.x;
454 amrex::ParticleReal const px = refpart.px;
455 amrex::ParticleReal const y = refpart.y;
456 amrex::ParticleReal const py = refpart.py;
457 amrex::ParticleReal const z = refpart.z;
458 amrex::ParticleReal const pz = refpart.pz;
459 amrex::ParticleReal const t = refpart.t;
460 amrex::ParticleReal const pt = refpart.pt;
461 amrex::ParticleReal const s = refpart.s;
462 amrex::ParticleReal const brho = refpart.rigidity_Tm();
463
464#if AMREX_DEVICE_COMPILE
465 amrex::ParticleReal const * k_normal = m_k_normal_d_data;
466#else
467 amrex::ParticleReal const * k_normal = m_k_normal_h_data;
468#endif
469 // find magnetic field and curvature; normalize units to MAD-X convention if needed
470 amrex::ParticleReal const curvature = m_unit == 1 ? k_normal[0] / brho : k_normal[0];
471 amrex::ParticleReal const rc = curvature == 0.0 ? 0.0_prt : 1.0_prt/curvature;
472
473 // length of the current slice
474 amrex::ParticleReal const slice_ds = m_ds / nslice();
475
476 // special case of zero field = an exact drift
477 if (rc==0_prt) {
478 // advance position and momentum (drift)
479 amrex::ParticleReal const step = slice_ds / std::sqrt(amrex::Math::powi<2>(pt) - 1.0_prt);
480 refpart.x = x + step*px;
481 refpart.y = y + step*py;
482 refpart.z = z + step*pz;
483 refpart.t = t - step*pt;
484
485 } else {
486
487 // assign intermediate parameters
488 amrex::ParticleReal const B = refpart.beta_gamma() /rc;
489 amrex::ParticleReal const theta = slice_ds / rc;
490
491 // calculate expensive terms once
492 auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
493
494 // advance position and momentum (bend)
495 refpart.px = px*cos_theta - pz*sin_theta;
496 refpart.py = py;
497 refpart.pz = pz*cos_theta + px*sin_theta;
498 refpart.pt = pt;
499
500 refpart.x = x + (refpart.pz - pz)/B;
501 refpart.y = y + (theta/B)*py;
502 refpart.z = z - (refpart.px - px)/B;
503 refpart.t = t - (theta/B)*pt;
504
505 }
506
507 // advance integrated path length
508 refpart.s = s + slice_ds;
509 }
510
511
512
518 Map6x6
519 transport_map (RefPart const & AMREX_RESTRICT refpart) const // TODO: update as well, but needs more careful placement of calc_constants
520 {
521 using namespace amrex::literals; // for _rt and _prt
522
523 // length of the current slice
524 amrex::ParticleReal const slice_ds = m_ds / nslice();
525
526 // find beta*gamma^2, beta
527 amrex::ParticleReal const betgam2 = amrex::Math::powi<2>(refpart.pt) - 1_prt;
528 amrex::ParticleReal const bet = refpart.beta();
529 amrex::ParticleReal const ibetgam2 = 1_prt / betgam2;
530
531 amrex::ParticleReal const * k_normal = m_k_normal_h_data;
532
533 // find magnetic field and curvature; normalize units to MAD-X convention if needed
534 amrex::ParticleReal const curvature = m_unit == 1 ? k_normal[0] / m_brho : k_normal[0];
535 amrex::ParticleReal const rc = curvature == 0.0 ? 0.0_prt : 1.0_prt/curvature;
536 amrex::ParticleReal const kn = m_unit == 1 ? k_normal[1] / m_brho : k_normal[1];
538
539 // update horizontal and longitudinal phase space variables
540 amrex::ParticleReal const gx = kn + amrex::Math::powi<-2>(rc);
541 amrex::ParticleReal const omega_x = std::sqrt(std::abs(gx));
542
543 // update vertical phase space variables
544 amrex::ParticleReal const gy = -kn;
545 amrex::ParticleReal const omega_y = std::sqrt(std::abs(gy));
546
547 // trigonometry
548 auto const [sinx, cosx] = amrex::Math::sincos(omega_x * slice_ds);
549 amrex::ParticleReal const sinhx = std::sinh(omega_x * slice_ds);
550 amrex::ParticleReal const coshx = std::cosh(omega_x * slice_ds);
551 auto const [siny, cosy] = amrex::Math::sincos(omega_y * slice_ds);
552 amrex::ParticleReal const sinhy = std::sinh(omega_y * slice_ds);
553 amrex::ParticleReal const coshy = std::cosh(omega_y * slice_ds);
554
555 amrex::ParticleReal const igbrc = 1_prt / (gx * bet * rc);
556 amrex::ParticleReal const iobrc = 1_prt / (omega_x * bet * rc);
557 amrex::ParticleReal const igobr = 1_prt / (gx * omega_x * b2rc2);
558
559 // initialize linear map matrix elements
561
562 R(1,1) = gx > 0_prt ? cosx : coshx;
563 R(1,2) = gx > 0_prt ? sinx / omega_x : sinhx / omega_x;
564 R(1,6) = gx > 0_prt ? -(1_prt - cosx) * igbrc : -(1_prt - coshx) * igbrc;
565 R(2,1) = gx > 0_prt ? -omega_x * sinx : omega_x * sinhx;
566 R(2,2) = gx > 0_prt ? cosx : coshx;
567 R(2,6) = gx > 0_prt ? -sinx * iobrc : -sinhx * iobrc;
568 R(3,3) = gy > 0_prt ? cosy : coshy;
569 R(3,4) = gy > 0_prt ? siny / omega_y : sinhy / omega_y;
570 R(4,3) = gy > 0_prt ? -omega_y * siny : omega_y * sinhy;
571 R(4,4) = gy > 0_prt ? cosy : coshy;
572 R(5,1) = gx > 0_prt ? sinx * iobrc : sinhx * iobrc;
573 R(5,2) = gx > 0_prt ? (1_prt - cosx) * igbrc : (1_prt - coshx) * igbrc;
574 R(5,6) = gx > 0_prt ?
575 slice_ds * ibetgam2 + (sinx - omega_x * slice_ds) * igobr :
576 slice_ds * ibetgam2 + (sinhx - omega_x * slice_ds) * igobr;
577
578 return R;
579 }
580
582 using LinearTransport::operator();
583
584
585 int m_unit;
588 int m_id;
589
590 int m_ncoef = 0;
595
596 private:
597 // constants that are independent of the individually tracked particle,
598 // see: compute_constants() to refresh
605 };
606
607} // namespace impactx
608
611
612#endif // IMPACTX_EXACTCFBEND_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
#define IMPACTX_GPUDATA_EXTERN(ElementType)
Definition dynamicdata.H:169
amrex_particle_real ParticleReal
__host__ __device__ std::pair< double, double > sincos(double x)
constexpr T powi(T x) noexcept
SmallMatrix< T, N, 1, Order::F, StartIndex > SmallVector
Definition All.H:56
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:33
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal beta_gamma() const
Definition ReferenceParticle.H:168
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal rigidity_Tm() const
Definition ReferenceParticle.H:261
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal beta() const
Definition ReferenceParticle.H:152
Definition ExactCFbend.H:45
mixin::TrackedVector< amrex::ParticleReal > k_normal
Definition ExactCFbend.H:46
mixin::TrackedVector< amrex::ParticleReal > k_skew
Definition ExactCFbend.H:47
Definition ExactCFbend.H:58
ImpactXParticleContainer::ParticleType PType
Definition ExactCFbend.H:60
int m_int_order
unit specification for Multipole strength
Definition ExactCFbend.H:586
amrex::ParticleReal m_ibeta
beta
Definition ExactCFbend.H:602
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void map1(amrex::ParticleReal const tau, amrex::SmallVector< amrex::ParticleReal, 6, 1 > &particle, amrex::ParticleReal &zeval) const
Definition ExactCFbend.H:241
amrex::ParticleReal const * m_k_skew_d_data
non-owning pointer to device normal coefficients
Definition ExactCFbend.H:594
int m_id
number of integration steps per slice
Definition ExactCFbend.H:588
amrex::ParticleReal const * m_k_skew_h_data
non-owning pointer to host normal coefficients
Definition ExactCFbend.H:592
amrex::ParticleReal m_rc
magnetic ridigity in T-m
Definition ExactCFbend.H:604
int m_mapsteps
order used for the symplectic integration (2 or 4)
Definition ExactCFbend.H:587
mixin::GPUDataRegistry< CFbendCoefficients > DynamicData
Definition ExactCFbend.H:62
amrex::ParticleReal m_beta
1 / m_betgam2
Definition ExactCFbend.H:601
amrex::ParticleReal m_slice_ds
non-owning pointer to device skew coefficients
Definition ExactCFbend.H:599
int m_ncoef
unique ExactMultipole id used for data lookup map
Definition ExactCFbend.H:590
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void map2(amrex::ParticleReal const tau, amrex::SmallVector< amrex::ParticleReal, 6, 1 > &particle, amrex::ParticleReal &zeval) const
Definition ExactCFbend.H:342
amrex::ParticleReal m_ibetgam2
m_ds / nslice();
Definition ExactCFbend.H:600
void compute_constants(RefPart const &refpart)
Definition ExactCFbend.H:140
amrex::ParticleReal const * m_k_normal_d_data
non-owning pointer to host skew coefficients
Definition ExactCFbend.H:593
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition ExactCFbend.H:519
amrex::ParticleReal m_brho
1 / m_beta
Definition ExactCFbend.H:603
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(amrex::ParticleReal &AMREX_RESTRICT x, amrex::ParticleReal &AMREX_RESTRICT y, amrex::ParticleReal &AMREX_RESTRICT t, amrex::ParticleReal &AMREX_RESTRICT px, amrex::ParticleReal &AMREX_RESTRICT py, amrex::ParticleReal &AMREX_RESTRICT pt, uint64_t &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition ExactCFbend.H:182
amrex::ParticleReal const * m_k_normal_h_data
number of Fourier coefficients
Definition ExactCFbend.H:591
ExactCFbend(amrex::ParticleReal ds, std::vector< amrex::ParticleReal > k_normal, std::vector< amrex::ParticleReal > k_skew, 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=1, int nslice=1, std::optional< std::string > name=std::nullopt)
Definition ExactCFbend.H:89
void reverse()
Definition ExactCFbend.H:128
int m_unit
Definition ExactCFbend.H:585
static constexpr auto type
Definition ExactCFbend.H:59
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:436
static std::shared_ptr< CFbendCoefficients > const & get(int id)
Definition dynamicdata.H:98
static CFbendCoefficients & emplace(int id, Args &&... args)
Definition dynamicdata.H:125
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 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: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
Definition TrackedVector.H:49