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 
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  // initialize output values of momenta
96  amrex::ParticleReal pxout = px;
97  amrex::ParticleReal pyout = py;
98  amrex::ParticleReal ptout = pt;
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*gamma^2
104  amrex::ParticleReal const pt_ref = refpart.pt;
105  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
106  amrex::ParticleReal const bet = sqrt(betgam2/(1.0_prt + betgam2));
107 
108  // update horizontal and longitudinal phase space variables
109  amrex::ParticleReal const gx = m_k + pow(m_rc,-2);
110  amrex::ParticleReal const omegax = sqrt(std::abs(gx));
111 
112  if(gx > 0.0) {
113  // calculate expensive terms once
114  auto const [sinx, cosx] = amrex::Math::sincos(omegax * slice_ds);
115  amrex::ParticleReal const r56 = slice_ds/betgam2
116  + (sinx - omegax*slice_ds)/(gx*omegax*pow(bet,2)*pow(m_rc,2));
117 
118  // advance position and momentum (focusing)
119  p.pos(RealAoS::x) = cosx*x + sinx/omegax*px - (1.0_prt - cosx)/(gx*bet*m_rc)*pt;
120  pxout = -omegax*sinx*x + cosx*px - sinx/(omegax*bet*m_rc)*pt;
121 
122  p.pos(RealAoS::t) = sinx/(omegax*bet*m_rc)*x + (1.0_prt - cosx)/(gx*bet*m_rc)*px
123  + t + r56*pt;
124  ptout = pt;
125  } else {
126  // calculate expensive terms once
127  amrex::ParticleReal const sinhx = sinh(omegax * slice_ds);
128  amrex::ParticleReal const coshx = cosh(omegax * slice_ds);
129  amrex::ParticleReal const r56 = slice_ds/betgam2
130  + (sinhx - omegax*slice_ds)/(gx*omegax*pow(bet,2)*pow(m_rc,2));
131 
132  // advance position and momentum (defocusing)
133  p.pos(RealAoS::x) = coshx*x + sinhx/omegax*px - (1.0_prt - coshx)/(gx*bet*m_rc)*pt;
134  pxout = omegax*sinhx*x + coshx*px - sinhx/(omegax*bet*m_rc)*pt;
135 
136  p.pos(RealAoS::t) = sinhx/(omegax*bet*m_rc)*x + (1.0_prt - coshx)/(gx*bet*m_rc)*px
137  + t + r56*pt;
138  ptout = pt;
139  }
140 
141  // update vertical phase space variables
142  amrex::ParticleReal const gy = -m_k;
143  amrex::ParticleReal const omegay = sqrt(std::abs(gy));
144 
145  if(gy > 0.0) {
146  // calculate expensive terms once
147  auto const [siny, cosy] = amrex::Math::sincos(omegay * slice_ds);
148 
149  // advance position and momentum (focusing)
150  p.pos(RealAoS::y) = cosy*y + siny/omegay*py;
151  pyout = -omegay*siny*y + cosy*py;
152 
153  } else {
154  // calculate expensive terms once
155  amrex::ParticleReal const sinhy = sinh(omegay * slice_ds);
156  amrex::ParticleReal const coshy = cosh(omegay * slice_ds);
157 
158  // advance position and momentum (defocusing)
159  p.pos(RealAoS::y) = coshy*y + sinhy/omegay*py;
160  pyout = omegay*sinhy*y + coshy*py;
161  }
162 
163  // assign updated momenta
164  px = pxout;
165  py = pyout;
166  pt = ptout;
167 
168  // undo shift due to alignment errors of the element
169  shift_out(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
170  }
171 
177  void operator() (RefPart & AMREX_RESTRICT refpart) const
178  {
179  using namespace amrex::literals; // for _rt and _prt
180 
181  // assign input reference particle values
182  amrex::ParticleReal const x = refpart.x;
183  amrex::ParticleReal const px = refpart.px;
184  amrex::ParticleReal const y = refpart.y;
185  amrex::ParticleReal const py = refpart.py;
186  amrex::ParticleReal const z = refpart.z;
187  amrex::ParticleReal const pz = refpart.pz;
188  amrex::ParticleReal const t = refpart.t;
189  amrex::ParticleReal const pt = refpart.pt;
190  amrex::ParticleReal const s = refpart.s;
191 
192  // length of the current slice
193  amrex::ParticleReal const slice_ds = m_ds / nslice();
194 
195  // assign intermediate parameter
196  amrex::ParticleReal const theta = slice_ds/m_rc;
197  amrex::ParticleReal const B = sqrt(pow(pt,2)-1.0_prt)/m_rc;
198 
199  // calculate expensive terms once
200  auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
201 
202  // advance position and momentum (bend)
203  refpart.px = px*cos_theta - pz*sin_theta;
204  refpart.py = py;
205  refpart.pz = pz*cos_theta + px*sin_theta;
206  refpart.pt = pt;
207 
208  refpart.x = x + (refpart.pz - pz)/B;
209  refpart.y = y + (theta/B)*py;
210  refpart.z = z - (refpart.px - px)/B;
211  refpart.t = t - (theta/B)*pt;
212 
213  // advance integrated path length
214  refpart.s = s + slice_ds;
215  }
216 
217  private:
218  amrex::ParticleReal m_rc;
219  amrex::ParticleReal m_k;
220  };
221 
222 } // namespace impactx
223 
224 #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:35
@ 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:218
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: CFbend.H:78
amrex::ParticleReal m_k
bend radius in m
Definition: CFbend.H:219
ImpactXParticleContainer::ParticleType PType
Definition: CFbend.H:35
static constexpr auto name
Definition: CFbend.H:34
@ x
position in x [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:49
@ y
position in y [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:50
@ t
c * time-of-flight [m] (at fixed s)
Definition: ImpactXParticleContainer.H:51
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