ImpactX
CFbend.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_CFBEND_H
11 #define IMPACTX_CFBEND_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_Math.H>
21 #include <AMReX_REAL.H>
22 
23 #include <cmath>
24 
25 
26 namespace impactx
27 {
28  struct CFbend
29  : public elements::BeamOptic<CFbend>,
30  public elements::Thick,
31  public elements::Alignment,
33  {
34  static constexpr auto name = "CFbend";
36 
52  amrex::ParticleReal ds,
53  amrex::ParticleReal rc,
54  amrex::ParticleReal k,
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_rc(rc), m_k(k)
63  {
64  }
65 
67  using BeamOptic::operator();
68 
82  amrex::ParticleReal & AMREX_RESTRICT x,
83  amrex::ParticleReal & AMREX_RESTRICT y,
84  amrex::ParticleReal & AMREX_RESTRICT t,
85  amrex::ParticleReal & AMREX_RESTRICT px,
86  amrex::ParticleReal & AMREX_RESTRICT py,
87  amrex::ParticleReal & AMREX_RESTRICT pt,
88  [[maybe_unused]] uint64_t & AMREX_RESTRICT idcpu,
89  RefPart const & refpart
90  ) const
91  {
92  using namespace amrex::literals; // for _rt and _prt
93 
94  // shift due to alignment errors of the element
95  shift_in(x, y, px, py);
96 
97  // initialize output values
98  amrex::ParticleReal xout = x;
99  amrex::ParticleReal yout = y;
100  amrex::ParticleReal tout = t;
101 
102  // initialize output values of momenta
103  amrex::ParticleReal pxout = px;
104  amrex::ParticleReal pyout = py;
105  amrex::ParticleReal ptout = pt;
106 
107  // length of the current slice
108  amrex::ParticleReal const slice_ds = m_ds / nslice();
109 
110  // access reference particle values to find beta*gamma^2
111  amrex::ParticleReal const pt_ref = refpart.pt;
112  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
113  amrex::ParticleReal const bet = sqrt(betgam2/(1.0_prt + betgam2));
114 
115  // update horizontal and longitudinal phase space variables
116  amrex::ParticleReal const gx = m_k + pow(m_rc,-2);
117  amrex::ParticleReal const omegax = sqrt(std::abs(gx));
118 
119  if(gx > 0.0) {
120  // calculate expensive terms once
121  auto const [sinx, cosx] = amrex::Math::sincos(omegax * slice_ds);
122  amrex::ParticleReal const r56 = slice_ds/betgam2
123  + (sinx - omegax*slice_ds)/(gx*omegax*pow(bet,2)*pow(m_rc,2));
124 
125  // advance position and momentum (focusing)
126  x = cosx*xout + sinx/omegax*px - (1.0_prt - cosx)/(gx*bet*m_rc)*pt;
127  pxout = -omegax*sinx*xout + cosx*px - sinx/(omegax*bet*m_rc)*pt;
128 
129  y = sinx/(omegax*bet*m_rc)*xout + (1.0_prt - cosx)/(gx*bet*m_rc)*px
130  + tout + r56*pt;
131  ptout = pt;
132  } else {
133  // calculate expensive terms once
134  amrex::ParticleReal const sinhx = sinh(omegax * slice_ds);
135  amrex::ParticleReal const coshx = cosh(omegax * slice_ds);
136  amrex::ParticleReal const r56 = slice_ds/betgam2
137  + (sinhx - omegax*slice_ds)/(gx*omegax*pow(bet,2)*pow(m_rc,2));
138 
139  // advance position and momentum (defocusing)
140  x = coshx*xout + sinhx/omegax*px - (1.0_prt - coshx)/(gx*bet*m_rc)*pt;
141  pxout = omegax*sinhx*xout + coshx*px - sinhx/(omegax*bet*m_rc)*pt;
142 
143  t = sinhx/(omegax*bet*m_rc)*xout + (1.0_prt - coshx)/(gx*bet*m_rc)*px
144  + tout + r56*pt;
145  ptout = pt;
146  }
147 
148  // update vertical phase space variables
149  amrex::ParticleReal const gy = -m_k;
150  amrex::ParticleReal const omegay = sqrt(std::abs(gy));
151 
152  if(gy > 0.0) {
153  // calculate expensive terms once
154  auto const [siny, cosy] = amrex::Math::sincos(omegay * slice_ds);
155 
156  // advance position and momentum (focusing)
157  y = cosy*yout + siny/omegay*py;
158  pyout = -omegay*siny*yout + cosy*py;
159 
160  } else {
161  // calculate expensive terms once
162  amrex::ParticleReal const sinhy = sinh(omegay * slice_ds);
163  amrex::ParticleReal const coshy = cosh(omegay * slice_ds);
164 
165  // advance position and momentum (defocusing)
166  y = coshy*yout + sinhy/omegay*py;
167  pyout = omegay*sinhy*yout + coshy*py;
168  }
169 
170  // assign updated momenta
171  px = pxout;
172  py = pyout;
173  pt = ptout;
174 
175  // undo shift due to alignment errors of the element
176  shift_out(x, y, px, py);
177  }
178 
184  void operator() (RefPart & AMREX_RESTRICT refpart) const
185  {
186  using namespace amrex::literals; // for _rt and _prt
187 
188  // assign input reference particle values
189  amrex::ParticleReal const x = refpart.x;
190  amrex::ParticleReal const px = refpart.px;
191  amrex::ParticleReal const y = refpart.y;
192  amrex::ParticleReal const py = refpart.py;
193  amrex::ParticleReal const z = refpart.z;
194  amrex::ParticleReal const pz = refpart.pz;
195  amrex::ParticleReal const t = refpart.t;
196  amrex::ParticleReal const pt = refpart.pt;
197  amrex::ParticleReal const s = refpart.s;
198 
199  // length of the current slice
200  amrex::ParticleReal const slice_ds = m_ds / nslice();
201 
202  // assign intermediate parameter
203  amrex::ParticleReal const theta = slice_ds/m_rc;
204  amrex::ParticleReal const B = sqrt(pow(pt,2)-1.0_prt)/m_rc;
205 
206  // calculate expensive terms once
207  auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
208 
209  // advance position and momentum (bend)
210  refpart.px = px*cos_theta - pz*sin_theta;
211  refpart.py = py;
212  refpart.pz = pz*cos_theta + px*sin_theta;
213  refpart.pt = pt;
214 
215  refpart.x = x + (refpart.pz - pz)/B;
216  refpart.y = y + (theta/B)*py;
217  refpart.z = z - (refpart.px - px)/B;
218  refpart.t = t - (theta/B)*pt;
219 
220  // advance integrated path length
221  refpart.s = s + slice_ds;
222  }
223 
224  amrex::ParticleReal m_rc;
225  amrex::ParticleReal m_k;
226  };
227 
228 } // namespace impactx
229 
230 #endif // IMPACTX_CFBEND_H
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE std::pair< double, double > sincos(double x)
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:33
@ t
fixed t as the independent variable
s
Definition: CFbend.H:33
CFbend(amrex::ParticleReal ds, amrex::ParticleReal rc, amrex::ParticleReal k, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, int nslice=1)
Definition: CFbend.H:51
amrex::ParticleReal m_rc
Definition: CFbend.H:224
amrex::ParticleReal m_k
bend radius in m
Definition: CFbend.H:225
ImpactXParticleContainer::ParticleType PType
Definition: CFbend.H:35
static constexpr auto name
Definition: CFbend.H:34
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, [[maybe_unused]] uint64_t &AMREX_RESTRICT idcpu, RefPart const &refpart) const
Definition: CFbend.H:81
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:149
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