ImpactX
Kicker.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_KICKER_H
11 #define IMPACTX_KICKER_H
12 
14 #include "mixin/beamoptic.H"
15 #include "mixin/thin.H"
16 #include "mixin/nofinalize.H"
17 
18 #include <AMReX_Extension.H>
19 #include <AMReX_REAL.H>
20 #include <AMReX_GpuComplex.H>
21 
22 #include <cmath>
23 
24 namespace impactx
25 {
26  struct Kicker
27  : public elements::BeamOptic<Kicker>,
28  public elements::Thin,
30  {
31  static constexpr auto name = "Kicker";
33 
35  {
37  Tm = 1
38  };
39 
47  Kicker (amrex::ParticleReal xkick,
48  amrex::ParticleReal ykick,
49  UnitSystem unit)
50  : m_xkick(xkick), m_ykick(ykick), m_unit(unit)
51  {
52  }
53 
55  using BeamOptic::operator();
56 
68  PType& AMREX_RESTRICT p,
69  amrex::ParticleReal & AMREX_RESTRICT px,
70  amrex::ParticleReal & AMREX_RESTRICT py,
71  amrex::ParticleReal & AMREX_RESTRICT pt,
72  [[maybe_unused]] RefPart const & refpart) const
73  {
74  using namespace amrex::literals; // for _rt and _prt
75 
76  // access AoS data such as positions and cpu/id
77  amrex::ParticleReal const x = p.pos(RealAoS::x);
78  amrex::ParticleReal const y = p.pos(RealAoS::y);
79  amrex::ParticleReal const t = p.pos(RealAoS::t);
80 
81  // normalize quad units to MAD-X convention if needed
82  amrex::ParticleReal dpx = m_xkick;
83  amrex::ParticleReal dpy = m_ykick;
84  if (m_unit == UnitSystem::Tm) {
85  dpx /= refpart.rigidity_Tm();
86  dpy /= refpart.rigidity_Tm();
87  }
88 
89  // intialize output values of momenta
90  amrex::ParticleReal pxout = px;
91  amrex::ParticleReal pyout = py;
92  amrex::ParticleReal ptout = pt;
93 
94  // advance position and momentum
95  p.pos(RealAoS::x) = x;
96  pxout = px + dpx;
97 
98  p.pos(RealAoS::y) = y;
99  pyout = py + dpy;
100 
101  p.pos(RealAoS::t) = t;
102  ptout = pt;
103 
104  // assign updated momenta
105  px = pxout;
106  py = pyout;
107  pt = ptout;
108  }
109 
111  using Thin::operator();
112 
113  private:
114  amrex::ParticleReal m_xkick;
115  amrex::ParticleReal m_ykick;
117  };
118 
119 } // namespace impactx
120 
121 #endif // IMPACTX_KICKER_H
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
Definition: ImpactX.cpp:32
Definition: Kicker.H:30
amrex::ParticleReal m_ykick
horizontal kick strength
Definition: Kicker.H:115
UnitSystem
Definition: Kicker.H:35
@ Tm
in units of the magnetic rigidity of the reference particle
Definition: Kicker.H:37
@ dimensionless
Definition: Kicker.H:36
static constexpr auto name
Definition: Kicker.H:31
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(PType &AMREX_RESTRICT p, amrex::ParticleReal &AMREX_RESTRICT px, amrex::ParticleReal &AMREX_RESTRICT py, amrex::ParticleReal &AMREX_RESTRICT pt, [[maybe_unused]] RefPart const &refpart) const
Definition: Kicker.H:67
UnitSystem m_unit
vertical kick strength
Definition: Kicker.H:116
amrex::ParticleReal m_xkick
Definition: Kicker.H:114
ImpactXParticleContainer::ParticleType PType
Definition: Kicker.H:32
Kicker(amrex::ParticleReal xkick, amrex::ParticleReal ykick, UnitSystem unit)
Definition: Kicker.H:47
@ x
position in x [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:42
@ y
position in y [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:43
@ t
c * time-of-flight [m] (at fixed s)
Definition: ImpactXParticleContainer.H:44
Definition: ReferenceParticle.H:30
Definition: beamoptic.H:135
Definition: nofinalize.H:22
Definition: thin.H:24