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/beamoptic.H"
15 #include "mixin/thick.H"
16 #include "mixin/nofinalize.H"
17 
18 #include <AMReX_Extension.H>
19 #include <AMReX_REAL.H>
20 
21 #include <cmath>
22 
23 
24 namespace impactx
25 {
26  struct CFbend
27  : public elements::BeamOptic<CFbend>,
28  public elements::Thick,
30  {
31  static constexpr auto name = "CFbend";
33 
46  amrex::ParticleReal const ds,
47  amrex::ParticleReal const rc,
48  amrex::ParticleReal const k,
49  int const nslice)
50  : Thick(ds, nslice), m_rc(rc), m_k(k)
51  {
52  }
53 
55  using BeamOptic::operator();
56 
67  PType& AMREX_RESTRICT p,
68  amrex::ParticleReal & AMREX_RESTRICT px,
69  amrex::ParticleReal & AMREX_RESTRICT py,
70  amrex::ParticleReal & AMREX_RESTRICT pt,
71  RefPart const & refpart) const
72  {
73  using namespace amrex::literals; // for _rt and _prt
74 
75  // access AoS data such as positions and cpu/id
76  amrex::ParticleReal const x = p.pos(RealAoS::x);
77  amrex::ParticleReal const y = p.pos(RealAoS::y);
78  amrex::ParticleReal const t = p.pos(RealAoS::t);
79 
80  // initialize output values of momenta
81  amrex::ParticleReal pxout = px;
82  amrex::ParticleReal pyout = py;
83  amrex::ParticleReal ptout = pt;
84 
85  // length of the current slice
86  amrex::ParticleReal const slice_ds = m_ds / nslice();
87 
88  // access reference particle values to find beta*gamma^2
89  amrex::ParticleReal const pt_ref = refpart.pt;
90  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
91  amrex::ParticleReal const bet = sqrt(betgam2/(1.0_prt + betgam2));
92 
93  // update horizontal and longitudinal phase space variables
94  amrex::ParticleReal const gx = m_k + pow(m_rc,-2);
95  amrex::ParticleReal const omegax = sqrt(std::abs(gx));
96 
97  if(gx > 0.0) {
98  // calculate expensive terms once
99  amrex::ParticleReal const sinx = sin(omegax * slice_ds);
100  amrex::ParticleReal const cosx = cos(omegax * slice_ds);
101  amrex::ParticleReal const r56 = slice_ds/betgam2
102  + (sinx - omegax*slice_ds)/(gx*omegax*pow(bet,2)*pow(m_rc,2));
103 
104  // advance position and momentum (focusing)
105  p.pos(RealAoS::x) = cosx*x + sinx/omegax*px - (1.0_prt - cosx)/(gx*bet*m_rc)*pt;
106  pxout = -omegax*sinx*x + cosx*px - sinx/(omegax*bet*m_rc)*pt;
107 
108  p.pos(RealAoS::t) = sinx/(omegax*bet*m_rc)*x + (1.0_prt - cosx)/(gx*bet*m_rc)*px
109  + t + r56*pt;
110  ptout = pt;
111  } else {
112  // calculate expensive terms once
113  amrex::ParticleReal const sinhx = sinh(omegax * slice_ds);
114  amrex::ParticleReal const coshx = cosh(omegax * slice_ds);
115  amrex::ParticleReal const r56 = slice_ds/betgam2
116  + (sinhx - omegax*slice_ds)/(gx*omegax*pow(bet,2)*pow(m_rc,2));
117 
118  // advance position and momentum (defocusing)
119  p.pos(RealAoS::x) = coshx*x + sinhx/omegax*px - (1.0_prt - coshx)/(gx*bet*m_rc)*pt;
120  pxout = omegax*sinhx*x + coshx*px - sinhx/(omegax*bet*m_rc)*pt;
121 
122  p.pos(RealAoS::t) = sinhx/(omegax*bet*m_rc)*x + (1.0_prt - coshx)/(gx*bet*m_rc)*px
123  + t + r56*pt;
124  ptout = pt;
125  }
126 
127  // update vertical phase space variables
128  amrex::ParticleReal const gy = -m_k;
129  amrex::ParticleReal const omegay = sqrt(std::abs(gy));
130 
131  if(gy > 0.0) {
132  // calculate expensive terms once
133  amrex::ParticleReal const siny = sin(omegay * slice_ds);
134  amrex::ParticleReal const cosy = cos(omegay * slice_ds);
135 
136  // advance position and momentum (focusing)
137  p.pos(RealAoS::y) = cosy*y + siny/omegay*py;
138  pyout = -omegay*siny*y + cosy*py;
139 
140  } else {
141  // calculate expensive terms once
142  amrex::ParticleReal const sinhy = sinh(omegay * slice_ds);
143  amrex::ParticleReal const coshy = cosh(omegay * slice_ds);
144 
145  // advance position and momentum (defocusing)
146  p.pos(RealAoS::y) = coshy*y + sinhy/omegay*py;
147  pyout = omegay*sinhy*y + coshy*py;
148  }
149 
150  // assign updated momenta
151  px = pxout;
152  py = pyout;
153  pt = ptout;
154  }
155 
161  void operator() (RefPart & AMREX_RESTRICT refpart) const
162  {
163  using namespace amrex::literals; // for _rt and _prt
164 
165  // assign input reference particle values
166  amrex::ParticleReal const x = refpart.x;
167  amrex::ParticleReal const px = refpart.px;
168  amrex::ParticleReal const y = refpart.y;
169  amrex::ParticleReal const py = refpart.py;
170  amrex::ParticleReal const z = refpart.z;
171  amrex::ParticleReal const pz = refpart.pz;
172  amrex::ParticleReal const t = refpart.t;
173  amrex::ParticleReal const pt = refpart.pt;
174  amrex::ParticleReal const s = refpart.s;
175 
176  // length of the current slice
177  amrex::ParticleReal const slice_ds = m_ds / nslice();
178 
179  // assign intermediate parameter
180  amrex::ParticleReal const theta = slice_ds/m_rc;
181  amrex::ParticleReal const B = sqrt(pow(pt,2)-1.0_prt)/m_rc;
182 
183  // calculate expensive terms once
184  // TODO: use sincos function once wrapped in AMReX
185  amrex::ParticleReal const sin_theta = sin(theta);
186  amrex::ParticleReal const cos_theta = cos(theta);
187 
188  // advance position and momentum (bend)
189  refpart.px = px*cos_theta - pz*sin_theta;
190  refpart.py = py;
191  refpart.pz = pz*cos_theta + px*sin_theta;
192  refpart.pt = pt;
193 
194  refpart.x = x + (refpart.pz - pz)/B;
195  refpart.y = y + (theta/B)*py;
196  refpart.z = z - (refpart.px - px)/B;
197  refpart.t = t - (theta/B)*pt;
198 
199  // advance integrated path length
200  refpart.s = s + slice_ds;
201  }
202 
203  private:
204  amrex::ParticleReal m_rc;
205  amrex::ParticleReal m_k;
206  };
207 
208 } // namespace impactx
209 
210 #endif // IMPACTX_CFBEND_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:33
s
Definition: CFbend.H:30
CFbend(amrex::ParticleReal const ds, amrex::ParticleReal const rc, amrex::ParticleReal const k, int const nslice)
Definition: CFbend.H:45
amrex::ParticleReal m_rc
Definition: CFbend.H:204
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:66
amrex::ParticleReal m_k
bend radius in m
Definition: CFbend.H:205
ImpactXParticleContainer::ParticleType PType
Definition: CFbend.H:32
static constexpr auto name
Definition: CFbend.H:31
@ 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::ParticleReal pt
energy, normalized by rest energy
Definition: ReferenceParticle.H:39
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