ImpactX
ChrUniformAcc.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_CHRACC_H
11 #define IMPACTX_CHRACC_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 
22 #include <cmath>
23 
24 
25 namespace impactx
26 {
27  struct ChrAcc
28  : public elements::BeamOptic<ChrAcc>,
29  public elements::Thick,
30  public elements::Alignment,
32  {
33  static constexpr auto name = "ChrAcc";
35 
52  amrex::ParticleReal ds,
53  amrex::ParticleReal ez,
54  amrex::ParticleReal bz,
55  amrex::ParticleReal dx = 0,
56  amrex::ParticleReal dy = 0,
57  amrex::ParticleReal rotation_degree = 0,
58  int nslice = 1
59  )
60  : Thick(ds, nslice),
61  Alignment(dx, dy, rotation_degree),
62  m_ez(ez), m_bz(bz)
63  {
64  }
65 
67  using BeamOptic::operator();
68 
79  PType& AMREX_RESTRICT p,
80  amrex::ParticleReal & AMREX_RESTRICT px,
81  amrex::ParticleReal & AMREX_RESTRICT py,
82  amrex::ParticleReal & AMREX_RESTRICT pt,
83  RefPart const & refpart) const {
84 
85  using namespace amrex::literals; // for _rt and _prt
86 
87  // shift due to alignment errors of the element
88  shift_in(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
89 
90  // access AoS data such as positions and cpu/id
91  amrex::ParticleReal const x = p.pos(RealAoS::x);
92  amrex::ParticleReal const y = p.pos(RealAoS::y);
93  amrex::ParticleReal const t = p.pos(RealAoS::t);
94 
95  // length of the current slice
96  amrex::ParticleReal const slice_ds = m_ds / nslice();
97 
98  // access reference particle values (final, initial):
99  amrex::ParticleReal const ptf_ref = refpart.pt;
100  amrex::ParticleReal const pti_ref = ptf_ref + m_ez*slice_ds;
101  amrex::ParticleReal const bgf = sqrt(pow(ptf_ref, 2) - 1.0_prt);
102  amrex::ParticleReal const bgi = sqrt(pow(pti_ref, 2) - 1.0_prt);
103 
104  // initial conversion from static to dynamic units:
105  px = px*bgi;
106  py = py*bgi;
107  pt = pt*bgi;
108 
109  // compute intermediate quantities related to acceleration
110  amrex::ParticleReal const pti_tot = pti_ref + pt;
111  amrex::ParticleReal const ptf_tot = ptf_ref + pt;
112  amrex::ParticleReal const pzi_tot = sqrt(pow(pti_tot,2)-1_prt);
113  amrex::ParticleReal const pzf_tot = sqrt(pow(ptf_tot,2)-1_prt);
114  amrex::ParticleReal const pzi_ref = sqrt(pow(pti_ref,2)-1_prt);
115  amrex::ParticleReal const pzf_ref = sqrt(pow(ptf_ref,2)-1_prt);
116 
117  amrex::ParticleReal const numer = -ptf_tot + pzf_tot;
118  amrex::ParticleReal const denom = -pti_tot + pzi_tot;
119 
120  // compute focusing constant (1/m) and rotation angle (in rad)
121  amrex::ParticleReal const alpha = m_bz/2.0_prt;
122  amrex::ParticleReal const theta = alpha/m_ez*log(numer/denom);
123 
124  // intialize output values
125  amrex::ParticleReal xout = x;
126  amrex::ParticleReal yout = y;
127  amrex::ParticleReal tout = t;
128  amrex::ParticleReal pxout = px;
129  amrex::ParticleReal pyout = py;
130  amrex::ParticleReal ptout = pt;
131 
132  // advance positions and momenta using map for focusing
133  xout = cos(theta)*x + sin(theta)/alpha*px;
134  pxout = -alpha*sin(theta)*x + cos(theta)*px;
135 
136  yout = cos(theta)*y + sin(theta)/alpha*py;
137  pyout = -alpha*sin(theta)*y + cos(theta)*py;
138 
139  // the correct symplectic update for t
140  tout = t + (pzf_tot - pzf_ref - pzi_tot + pzi_ref)/m_ez;
141  tout = tout + (1_prt/pzi_tot - 1_prt/pzf_tot)*(pow(py-alpha*x,2)+pow(px+alpha*y,2))/(2_prt*m_ez);
142  ptout = pt;
143 
144  // assign intermediate momenta
145  px = pxout;
146  py = pyout;
147  pt = ptout;
148 
149  // advance positions and momenta using map for rotation
150  p.pos(RealAoS::x) = cos(theta)*xout + sin(theta)*yout;
151  pxout = cos(theta)*px + sin(theta)*py;
152 
153  p.pos(RealAoS::y) = -sin(theta)*xout + cos(theta)*yout;
154  pyout = -sin(theta)*px + cos(theta)*py;
155 
156  p.pos(RealAoS::t) = tout;
157  ptout = pt;
158 
159  // assign updated momenta
160  px = pxout;
161  py = pyout;
162  pt = ptout;
163 
164  // final conversion from dynamic to static units:
165  px = px/bgf;
166  py = py/bgf;
167  pt = pt/bgf;
168 
169  // undo shift due to alignment errors of the element
170  shift_out(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
171  }
172 
178  void operator() (RefPart & AMREX_RESTRICT refpart) const {
179 
180  using namespace amrex::literals; // for _rt and _prt
181 
182  // assign input reference particle values
183  amrex::ParticleReal const x = refpart.x;
184  amrex::ParticleReal const px = refpart.px;
185  amrex::ParticleReal const y = refpart.y;
186  amrex::ParticleReal const py = refpart.py;
187  amrex::ParticleReal const z = refpart.z;
188  amrex::ParticleReal const pz = refpart.pz;
189  amrex::ParticleReal const t = refpart.t;
190  amrex::ParticleReal const pt = refpart.pt;
191  amrex::ParticleReal const s = refpart.s;
192 
193  // length of the current slice
194  amrex::ParticleReal const slice_ds = m_ds / nslice();
195 
196  // compute intial value of beta*gamma
197  amrex::ParticleReal const bgi = sqrt(pow(pt, 2) - 1.0_prt);
198 
199  // advance pt (uniform acceleration)
200  refpart.pt = pt - m_ez*slice_ds;
201 
202  // compute final value of beta*gamma
203  amrex::ParticleReal const ptf = refpart.pt;
204  amrex::ParticleReal const bgf = sqrt(pow(ptf, 2) - 1.0_prt);
205 
206  // update t
207  refpart.t = t + (bgf - bgi)/m_ez;
208 
209  // advance position (x,y,z)
210  refpart.x = x + slice_ds*px/bgi;
211  refpart.y = y + slice_ds*py/bgi;
212  refpart.z = z + slice_ds*pz/bgi;
213 
214  // advance momentum (px,py,pz)
215  refpart.px = px*bgf/bgi;
216  refpart.py = py*bgf/bgi;
217  refpart.pz = pz*bgf/bgi;
218 
219  // advance integrated path length
220  refpart.s = s + slice_ds;
221  }
222 
223  private:
224  amrex::ParticleReal m_ez;
225  amrex::ParticleReal m_bz;
226  };
227 
228 } // namespace impactx
229 
230 #endif // IMPACTX_CHRACC_H
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > log(const GpuComplex< T > &a_z) noexcept
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
Definition: ChrUniformAcc.H:32
amrex::ParticleReal m_bz
electric field strength in 1/m
Definition: ChrUniformAcc.H:225
ImpactXParticleContainer::ParticleType PType
Definition: ChrUniformAcc.H:34
amrex::ParticleReal m_ez
Definition: ChrUniformAcc.H:224
ChrAcc(amrex::ParticleReal ds, amrex::ParticleReal ez, amrex::ParticleReal bz, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, int nslice=1)
Definition: ChrUniformAcc.H:51
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: ChrUniformAcc.H:78
static constexpr auto name
Definition: ChrUniformAcc.H:33
@ 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::ParticleReal pt
energy, normalized by rest energy
Definition: ReferenceParticle.H:39
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