ImpactX
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/beamoptic.H"
16 #include "mixin/thick.H"
17 #include "mixin/nofinalize.H"
18 
19 #include <AMReX_Extension.H>
20 #include <AMReX_REAL.H>
21 #include <AMReX_Print.H> // for PrintToFile
22 
23 #include <cmath>
24 
25 
26 namespace impactx
27 {
28  struct ChrQuad
29  : public elements::BeamOptic<ChrQuad>,
30  public elements::Thick,
31  public elements::Alignment,
33  {
34  static constexpr auto name = "ChrQuad";
36 
57  amrex::ParticleReal ds,
58  amrex::ParticleReal k,
59  int unit,
60  amrex::ParticleReal dx = 0,
61  amrex::ParticleReal dy = 0,
62  amrex::ParticleReal rotation_degree = 0,
63  int nslice = 1
64  )
65  : Thick(ds, nslice),
66  Alignment(dx, dy, rotation_degree),
67  m_k(k), m_unit(unit)
68  {
69  }
70 
72  using BeamOptic::operator();
73 
84  PType& AMREX_RESTRICT p,
85  amrex::ParticleReal & AMREX_RESTRICT px,
86  amrex::ParticleReal & AMREX_RESTRICT py,
87  amrex::ParticleReal & AMREX_RESTRICT pt,
88  RefPart const & refpart) const {
89 
90  using namespace amrex::literals; // for _rt and _prt
91 
92  // shift due to alignment errors of the element
93  shift_in(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
94 
95  // access AoS data such as positions and cpu/id
96  amrex::ParticleReal const x = p.pos(RealAoS::x);
97  amrex::ParticleReal const y = p.pos(RealAoS::y);
98  amrex::ParticleReal const t = p.pos(RealAoS::t);
99 
100  // length of the current slice
101  amrex::ParticleReal const slice_ds = m_ds / nslice();
102 
103  // access reference particle values to find beta
104  amrex::ParticleReal const bet = refpart.beta();
105 
106  // normalize quad units to MAD-X convention if needed
107  amrex::ParticleReal g = m_k;
108  if (m_unit == 1) {
109  g = m_k / refpart.rigidity_Tm();
110  }
111 
112  // compute particle momentum deviation delta + 1
113  amrex::ParticleReal delta1;
114  delta1 = sqrt(1_prt - 2_prt*pt/bet + pow(pt,2));
115  amrex::ParticleReal const delta = delta1 - 1_prt;
116 
117  // compute phase advance per unit length in s (in rad/m)
118  // chromatic dependence on delta is included
119  amrex::ParticleReal const omega = sqrt(std::abs(g)/delta1);
120 
121  // intialize output values of momenta
122  amrex::ParticleReal pxout = px;
123  amrex::ParticleReal pyout = py;
124  amrex::ParticleReal const ptout = pt;
125 
126  // paceholder variables
127  amrex::ParticleReal q1 = x;
128  amrex::ParticleReal q2 = y;
129  amrex::ParticleReal p1 = px;
130  amrex::ParticleReal p2 = py;
131 
132  if(g > 0.0) {
133  // advance transverse position and momentum (focusing quad)
134  p.pos(RealAoS::x) = cos(omega*slice_ds)*x +
135  sin(omega*slice_ds)/(omega*delta1)*px;
136  pxout = -omega*delta1*sin(omega*slice_ds)*x + cos(omega*slice_ds)*px;
137 
138  p.pos(RealAoS::y) = cosh(omega*slice_ds)*y +
139  sinh(omega*slice_ds)/(omega*delta1)*py;
140  pyout = omega*delta1*sinh(omega*slice_ds)*y + cosh(omega*slice_ds)*py;
141 
142  } else {
143  // advance transverse position and momentum (defocusing quad)
144  p.pos(RealAoS::x) = cosh(omega*slice_ds)*x +
145  sinh(omega*slice_ds)/(omega*delta1)*px;
146  pxout = omega*delta1*sinh(omega*slice_ds)*x + cosh(omega*slice_ds)*px;
147 
148  p.pos(RealAoS::y) = cos(omega*slice_ds)*y +
149  sin(omega*slice_ds)/(omega*delta1)*py;
150  pyout = -omega*delta1*sin(omega*slice_ds)*y + cos(omega*slice_ds)*py;
151 
152  q1 = y;
153  q2 = x;
154  p1 = py;
155  p2 = px;
156 
157  }
158 
159  // advance longitudinal position and momentum
160 
161  // the corresponding symplectic update to t
162  amrex::ParticleReal const term = pt + delta/bet;
163  amrex::ParticleReal const t0 = t - term*slice_ds/delta1;
164 
165  amrex::ParticleReal const w = omega*delta1;
166  amrex::ParticleReal const term1 = -(pow(p2,2)+pow(q2,2)*pow(w,2))*sinh(2_prt*slice_ds*omega);
167  amrex::ParticleReal const term2 = -(pow(p1,2)-pow(q1,2)*pow(w,2))*sin(2_prt*slice_ds*omega);
168  amrex::ParticleReal const term3 = -2_prt*q2*p2*w*cosh(2_prt*slice_ds*omega);
169  amrex::ParticleReal const term4 = -2_prt*q1*p1*w*cos(2_prt*slice_ds*omega);
170  amrex::ParticleReal const term5 = 2_prt*omega*(q1*p1*delta1 + q2*p2*delta1
171  -(pow(p1,2)+pow(p2,2))*slice_ds - (pow(q1,2)-pow(q2,2))*pow(w,2)*slice_ds);
172  p.pos(RealAoS::t) = t0 + (-1_prt+bet*pt)/(8_prt*bet*pow(delta1,3)*omega)
173  *(term1+term2+term3+term4+term5);
174 
175  // ptout = pt;
176 
177  // assign updated momenta
178  px = pxout;
179  py = pyout;
180  pt = ptout;
181 
182  // undo shift due to alignment errors of the element
183  shift_out(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
184  }
185 
191  void operator() (RefPart & AMREX_RESTRICT refpart) const {
192 
193  using namespace amrex::literals; // for _rt and _prt
194 
195  // assign input reference particle values
196  amrex::ParticleReal const x = refpart.x;
197  amrex::ParticleReal const px = refpart.px;
198  amrex::ParticleReal const y = refpart.y;
199  amrex::ParticleReal const py = refpart.py;
200  amrex::ParticleReal const z = refpart.z;
201  amrex::ParticleReal const pz = refpart.pz;
202  amrex::ParticleReal const t = refpart.t;
203  amrex::ParticleReal const pt = refpart.pt;
204  amrex::ParticleReal const s = refpart.s;
205 
206  // length of the current slice
207  amrex::ParticleReal const slice_ds = m_ds / nslice();
208 
209  // assign intermediate parameter
210  amrex::ParticleReal const step = slice_ds / sqrt(pow(pt,2)-1.0_prt);
211 
212  // advance position and momentum (straight element)
213  refpart.x = x + step*px;
214  refpart.y = y + step*py;
215  refpart.z = z + step*pz;
216  refpart.t = t - step*pt;
217 
218  // advance integrated path length
219  refpart.s = s + slice_ds;
220  }
221 
222  private:
223  amrex::ParticleReal m_k;
224  int m_unit;
225  };
226 
227 } // namespace impactx
228 
229 #endif // IMPACTX_CHRQUAD_H
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > pow(const GpuComplex< T > &a_z, const T &a_y) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
Definition: ImpactX.cpp:36
@ t
fixed t as the independent variable
s
tuple w
Definition: ChrQuad.H:33
ImpactXParticleContainer::ParticleType PType
Definition: ChrQuad.H:35
static constexpr auto name
Definition: ChrQuad.H:34
ChrQuad(amrex::ParticleReal ds, amrex::ParticleReal k, int unit, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, int nslice=1)
Definition: ChrQuad.H:56
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, RefPart const &refpart) const
Definition: ChrQuad.H:83
amrex::ParticleReal m_k
Definition: ChrQuad.H:223
int m_unit
quadrupole strength in 1/m^2 (or T/m)
Definition: ChrQuad.H:224
@ x
position in x [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:48
@ y
position in y [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:49
@ t
c * time-of-flight [m] (at fixed s)
Definition: ImpactXParticleContainer.H:50
Definition: ReferenceParticle.H:30
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal rigidity_Tm() const
Definition: ReferenceParticle.H:169
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal beta() const
Definition: ReferenceParticle.H:64
Definition: alignment.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void shift_out(amrex::ParticleReal &AMREX_RESTRICT x, amrex::ParticleReal &AMREX_RESTRICT y, amrex::ParticleReal &AMREX_RESTRICT px, amrex::ParticleReal &AMREX_RESTRICT py) const
Definition: alignment.H:91
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dx() const
Definition: alignment.H:120
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void shift_in(amrex::ParticleReal &AMREX_RESTRICT x, amrex::ParticleReal &AMREX_RESTRICT y, amrex::ParticleReal &AMREX_RESTRICT px, amrex::ParticleReal &AMREX_RESTRICT py) const
Definition: alignment.H:61
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dy() const
Definition: alignment.H:130
Definition: beamoptic.H:135
Definition: nofinalize.H:22
Definition: thick.H:24
Thick(amrex::ParticleReal ds, int nslice)
Definition: thick.H:30
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
amrex::ParticleReal m_ds
Definition: thick.H:58