ImpactX
Sbend.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_SBEND_H
11 #define IMPACTX_SBEND_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 Sbend
29  : public elements::BeamOptic<Sbend>,
30  public elements::Thick,
31  public elements::Alignment,
33  {
34  static constexpr auto name = "Sbend";
36 
47  amrex::ParticleReal ds,
48  amrex::ParticleReal rc,
49  amrex::ParticleReal dx = 0,
50  amrex::ParticleReal dy = 0,
51  amrex::ParticleReal rotation_degree = 0,
52  int nslice = 1
53  )
54  : Thick(ds, nslice),
55  Alignment(dx, dy, rotation_degree),
56  m_rc(rc)
57  {
58  }
59 
61  using BeamOptic::operator();
62 
73  PType& AMREX_RESTRICT p,
74  amrex::ParticleReal & AMREX_RESTRICT px,
75  amrex::ParticleReal & AMREX_RESTRICT py,
76  amrex::ParticleReal & AMREX_RESTRICT pt,
77  RefPart const & refpart) const {
78 
79  using namespace amrex::literals; // for _rt and _prt
80 
81  // shift due to alignment errors of the element
82  shift_in(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
83 
84  // access AoS data such as positions and cpu/id
85  amrex::ParticleReal const x = p.pos(RealAoS::x);
86  amrex::ParticleReal const y = p.pos(RealAoS::y);
87  amrex::ParticleReal const t = p.pos(RealAoS::t);
88 
89  // initialize output values of momenta
90  amrex::ParticleReal pxout = px;
91  amrex::ParticleReal const pyout = py;
92  amrex::ParticleReal const ptout = pt;
93 
94  // length of the current slice
95  amrex::ParticleReal const slice_ds = m_ds / nslice();
96 
97  // access reference particle values to find beta*gamma^2
98  amrex::ParticleReal const pt_ref = refpart.pt;
99  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
100  amrex::ParticleReal const bet = sqrt(betgam2/(1.0_prt + betgam2));
101 
102  // calculate expensive terms once
103  amrex::ParticleReal const theta = slice_ds/m_rc;
104  auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
105 
106  // advance position and momentum (sector bend)
107  p.pos(RealAoS::x) = cos_theta*x + m_rc*sin_theta*px
108  - (m_rc/bet)*(1.0_prt - cos_theta)*pt;
109 
110  pxout = -sin_theta/m_rc*x + cos_theta*px - sin_theta/bet*pt;
111 
112  p.pos(RealAoS::y) = y + m_rc*theta*py;
113 
114  // pyout = py;
115 
116  p.pos(RealAoS::t) = sin_theta/bet*x + m_rc/bet*(1.0_prt - cos_theta)*px + t
117  + m_rc*(-theta+sin_theta/(bet*bet))*pt;
118 
119  // ptout = pt;
120 
121  // assign updated momenta
122  px = pxout;
123  py = pyout;
124  pt = ptout;
125 
126  // undo shift due to alignment errors of the element
127  shift_out(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
128  }
129 
135  void operator() (RefPart & AMREX_RESTRICT refpart) const {
136 
137  using namespace amrex::literals; // for _rt and _prt
138 
139  // assign input reference particle values
140  amrex::ParticleReal const x = refpart.x;
141  amrex::ParticleReal const px = refpart.px;
142  amrex::ParticleReal const y = refpart.y;
143  amrex::ParticleReal const py = refpart.py;
144  amrex::ParticleReal const z = refpart.z;
145  amrex::ParticleReal const pz = refpart.pz;
146  amrex::ParticleReal const t = refpart.t;
147  amrex::ParticleReal const pt = refpart.pt;
148  amrex::ParticleReal const s = refpart.s;
149 
150  // length of the current slice
151  amrex::ParticleReal const slice_ds = m_ds / nslice();
152 
153  // assign intermediate parameter
154  amrex::ParticleReal const theta = slice_ds/m_rc;
155  amrex::ParticleReal const B = sqrt(pow(pt,2)-1.0_prt)/m_rc;
156 
157  // calculate expensive terms once
158  auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
159 
160  // advance position and momentum (bend)
161  refpart.px = px*cos_theta - pz*sin_theta;
162  refpart.py = py;
163  refpart.pz = pz*cos_theta + px*sin_theta;
164  refpart.pt = pt;
165 
166  refpart.x = x + (refpart.pz - pz)/B;
167  refpart.y = y + (theta/B)*py;
168  refpart.z = z - (refpart.px - px)/B;
169  refpart.t = t - (theta/B)*pt;
170 
171  // advance integrated path length
172  refpart.s = s + slice_ds;
173 
174  }
175 
176  private:
177  amrex::ParticleReal m_rc;
178  };
179 
180 } // namespace impactx
181 
182 #endif // IMPACTX_SBEND_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:36
@ t
fixed t as the independent variable
s
@ 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: Sbend.H:33
Sbend(amrex::ParticleReal ds, amrex::ParticleReal rc, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, int nslice=1)
Definition: Sbend.H:46
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: Sbend.H:72
static constexpr auto name
Definition: Sbend.H:34
amrex::ParticleReal m_rc
Definition: Sbend.H:177
ImpactXParticleContainer::ParticleType PType
Definition: Sbend.H:35
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