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/beamoptic.H"
15 #include "mixin/thick.H"
16 #include "mixin/nofinalize.H"
17 
18 #include <AMReX_Extension.H>
19 #include <AMReX_REAL.H>
20 #include <AMReX_Print.H> // for PrintToFile
21 
22 #include <cmath>
23 
24 
25 namespace impactx
26 {
27  struct ChrQuad
28  : public elements::BeamOptic<ChrQuad>,
29  public elements::Thick,
31  {
32  static constexpr auto name = "ChrQuad";
34 
51  ChrQuad( amrex::ParticleReal const ds, amrex::ParticleReal const k,
52  int const unit, int const nslice )
53  : Thick(ds, nslice), m_k(k), m_unit(unit)
54  {
55  }
56 
58  using BeamOptic::operator();
59 
70  PType& AMREX_RESTRICT p,
71  amrex::ParticleReal & AMREX_RESTRICT px,
72  amrex::ParticleReal & AMREX_RESTRICT py,
73  amrex::ParticleReal & AMREX_RESTRICT pt,
74  RefPart const & refpart) const {
75 
76  using namespace amrex::literals; // for _rt and _prt
77 
78  // access AoS data such as positions and cpu/id
79  amrex::ParticleReal const x = p.pos(RealAoS::x);
80  amrex::ParticleReal const y = p.pos(RealAoS::y);
81  amrex::ParticleReal const t = p.pos(RealAoS::t);
82 
83  // length of the current slice
84  amrex::ParticleReal const slice_ds = m_ds / nslice();
85 
86  // access reference particle values to find beta
87  amrex::ParticleReal const bet = refpart.beta();
88 
89  // normalize quad units to MAD-X convention if needed
90  amrex::ParticleReal g = m_k;
91  if (m_unit == 1) {
92  g = m_k / refpart.rigidity_Tm();
93  }
94 
95  // compute particle momentum deviation delta + 1
96  amrex::ParticleReal delta1;
97  delta1 = sqrt(1_prt - 2_prt*pt/bet + pow(pt,2));
98  amrex::ParticleReal const delta = delta1 - 1_prt;
99 
100  // compute phase advance per unit length in s (in rad/m)
101  // chromatic dependence on delta is included
102  amrex::ParticleReal const omega = sqrt(std::abs(g)/delta1);
103 
104  // intialize output values of momenta
105  amrex::ParticleReal pxout = px;
106  amrex::ParticleReal pyout = py;
107  amrex::ParticleReal const ptout = pt;
108 
109  // paceholder variables
110  amrex::ParticleReal q1 = x;
111  amrex::ParticleReal q2 = y;
112  amrex::ParticleReal p1 = px;
113  amrex::ParticleReal p2 = py;
114 
115  if(g > 0.0) {
116  // advance transverse position and momentum (focusing quad)
117  p.pos(RealAoS::x) = cos(omega*slice_ds)*x +
118  sin(omega*slice_ds)/(omega*delta1)*px;
119  pxout = -omega*delta1*sin(omega*slice_ds)*x + cos(omega*slice_ds)*px;
120 
121  p.pos(RealAoS::y) = cosh(omega*slice_ds)*y +
122  sinh(omega*slice_ds)/(omega*delta1)*py;
123  pyout = omega*delta1*sinh(omega*slice_ds)*y + cosh(omega*slice_ds)*py;
124 
125  } else {
126  // advance transverse position and momentum (defocusing quad)
127  p.pos(RealAoS::x) = cosh(omega*slice_ds)*x +
128  sinh(omega*slice_ds)/(omega*delta1)*px;
129  pxout = omega*delta1*sinh(omega*slice_ds)*x + cosh(omega*slice_ds)*px;
130 
131  p.pos(RealAoS::y) = cos(omega*slice_ds)*y +
132  sin(omega*slice_ds)/(omega*delta1)*py;
133  pyout = -omega*delta1*sin(omega*slice_ds)*y + cos(omega*slice_ds)*py;
134 
135  q1 = y;
136  q2 = x;
137  p1 = py;
138  p2 = px;
139 
140  }
141 
142  // advance longitudinal position and momentum
143 
144  // the corresponding symplectic update to t
145  amrex::ParticleReal const term = pt + delta/bet;
146  amrex::ParticleReal const t0 = t - term*slice_ds/delta1;
147 
148  amrex::ParticleReal const w = omega*delta1;
149  amrex::ParticleReal const term1 = -(pow(p2,2)+pow(q2,2)*pow(w,2))*sinh(2_prt*slice_ds*omega);
150  amrex::ParticleReal const term2 = -(pow(p1,2)-pow(q1,2)*pow(w,2))*sin(2_prt*slice_ds*omega);
151  amrex::ParticleReal const term3 = -2_prt*q2*p2*w*cosh(2_prt*slice_ds*omega);
152  amrex::ParticleReal const term4 = -2_prt*q1*p1*w*cos(2_prt*slice_ds*omega);
153  amrex::ParticleReal const term5 = 2_prt*omega*(q1*p1*delta1 + q2*p2*delta1
154  -(pow(p1,2)+pow(p2,2))*slice_ds - (pow(q1,2)-pow(q2,2))*pow(w,2)*slice_ds);
155  p.pos(RealAoS::t) = t0 + (-1_prt+bet*pt)/(8_prt*bet*pow(delta1,3)*omega)
156  *(term1+term2+term3+term4+term5);
157 
158  // ptout = pt;
159 
160  // assign updated momenta
161  px = pxout;
162  py = pyout;
163  pt = ptout;
164 
165  }
166 
172  void operator() (RefPart & AMREX_RESTRICT refpart) const {
173 
174  using namespace amrex::literals; // for _rt and _prt
175 
176  // assign input reference particle values
177  amrex::ParticleReal const x = refpart.x;
178  amrex::ParticleReal const px = refpart.px;
179  amrex::ParticleReal const y = refpart.y;
180  amrex::ParticleReal const py = refpart.py;
181  amrex::ParticleReal const z = refpart.z;
182  amrex::ParticleReal const pz = refpart.pz;
183  amrex::ParticleReal const t = refpart.t;
184  amrex::ParticleReal const pt = refpart.pt;
185  amrex::ParticleReal const s = refpart.s;
186 
187  // length of the current slice
188  amrex::ParticleReal const slice_ds = m_ds / nslice();
189 
190  // assign intermediate parameter
191  amrex::ParticleReal const step = slice_ds / sqrt(pow(pt,2)-1.0_prt);
192 
193  // advance position and momentum (straight element)
194  refpart.x = x + step*px;
195  refpart.y = y + step*py;
196  refpart.z = z + step*pz;
197  refpart.t = t - step*pt;
198 
199  // advance integrated path length
200  refpart.s = s + slice_ds;
201  }
202 
203  private:
204  amrex::ParticleReal m_k;
205  int m_unit;
206  };
207 
208 } // namespace impactx
209 
210 #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:34
s
tuple w
Definition: ChrQuad.H:31
ImpactXParticleContainer::ParticleType PType
Definition: ChrQuad.H:33
static constexpr auto name
Definition: ChrQuad.H:32
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:69
amrex::ParticleReal m_k
Definition: ChrQuad.H:204
ChrQuad(amrex::ParticleReal const ds, amrex::ParticleReal const k, int const unit, int const nslice)
Definition: ChrQuad.H:51
int m_unit
quadrupole strength in 1/m^2 (or T/m)
Definition: ChrQuad.H:205
@ 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
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: beamoptic.H:135
Definition: nofinalize.H:22
Definition: thick.H:24
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal ds() const
Definition: thick.H:50
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int nslice() const
Definition: thick.H:40
amrex::ParticleReal m_ds
Definition: thick.H:56
Thick(amrex::ParticleReal const ds, int const nslice)
Definition: thick.H:30