ImpactX
Sbend.H
Go to the documentation of this file.
1 /* Copyright 2022 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 
15 #include <AMReX_Extension.H>
16 #include <AMReX_REAL.H>
17 
18 #include <cmath>
19 
20 
21 namespace impactx
22 {
23  struct Sbend
24  {
25  static constexpr auto name = "Sbend";
27 
34  Sbend( amrex::ParticleReal const ds, amrex::ParticleReal const rc,
35  int const nslice)
36  : m_ds(ds), m_rc(rc), m_nslice(nslice)
37  {
38  }
39 
50  PType& AMREX_RESTRICT p,
51  amrex::ParticleReal & AMREX_RESTRICT px,
52  amrex::ParticleReal & AMREX_RESTRICT py,
53  amrex::ParticleReal & AMREX_RESTRICT pt,
54  RefPart const refpart) const {
55 
56  using namespace amrex::literals; // for _rt and _prt
57 
58  // access AoS data such as positions and cpu/id
59  amrex::ParticleReal const x = p.pos(0);
60  amrex::ParticleReal const y = p.pos(1);
61  amrex::ParticleReal const t = p.pos(2);
62 
63  // initialize output values of momenta
64  amrex::ParticleReal pxout = px;
65  amrex::ParticleReal pyout = py;
66  amrex::ParticleReal ptout = pt;
67 
68  // length of the current slice
69  amrex::ParticleReal const slice_ds = m_ds / nslice();
70 
71  // access reference particle values to find beta*gamma^2
72  amrex::ParticleReal const pt_ref = refpart.pt;
73  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
74  amrex::ParticleReal const bet = sqrt(betgam2/(1.0_prt + betgam2));
75 
76  // calculate expensive terms once
77  // TODO: use sincos function once wrapped in AMReX
78  amrex::ParticleReal const theta = slice_ds/m_rc;
79  amrex::ParticleReal const sin_theta = sin(theta);
80  amrex::ParticleReal const cos_theta = cos(theta);
81 
82  // advance position and momentum (sector bend)
83  p.pos(0) = cos_theta*x + m_rc*sin_theta*px
84  - (m_rc/bet)*(1.0_prt - cos_theta)*pt;
85 
86  pxout = -sin_theta/m_rc*x + cos_theta*px - sin_theta/bet*pt;
87 
88  p.pos(1) = y + m_rc*theta*py;
89 
90  // pyout = py;
91 
92  p.pos(2) = sin_theta/bet*x + m_rc/bet*(1.0_prt - cos_theta)*px + t
93  + m_rc*(-theta+sin_theta/(bet*bet))*pt;
94 
95  // ptout = pt;
96 
97  // assign updated momenta
98  px = pxout;
99  py = pyout;
100  pt = ptout;
101 
102  }
103 
110  RefPart & AMREX_RESTRICT refpart) const {
111 
112  using namespace amrex::literals; // for _rt and _prt
113 
114  // assign input reference particle values
115  amrex::ParticleReal const x = refpart.x;
116  amrex::ParticleReal const px = refpart.px;
117  amrex::ParticleReal const y = refpart.y;
118  amrex::ParticleReal const py = refpart.py;
119  amrex::ParticleReal const z = refpart.z;
120  amrex::ParticleReal const pz = refpart.pz;
121  amrex::ParticleReal const t = refpart.t;
122  amrex::ParticleReal const pt = refpart.pt;
123  amrex::ParticleReal const s = refpart.s;
124 
125  // length of the current slice
126  amrex::ParticleReal const slice_ds = m_ds / nslice();
127 
128  // assign intermediate parameter
129  amrex::ParticleReal const theta = slice_ds/m_rc;
130  amrex::ParticleReal const B = sqrt(pow(pt,2)-1.0_prt)/m_rc;
131 
132  // calculate expensive terms once
133  // TODO: use sincos function once wrapped in AMReX
134  amrex::ParticleReal const sin_theta = sin(theta);
135  amrex::ParticleReal const cos_theta = cos(theta);
136 
137  // advance position and momentum (bend)
138  refpart.px = px*cos_theta - pz*sin_theta;
139  refpart.py = py;
140  refpart.pz = pz*cos_theta + px*sin_theta;
141  refpart.pt = pt;
142 
143  refpart.x = x + (refpart.pz - pz)/B;
144  refpart.y = y + (theta/B)*py;
145  refpart.z = z - (refpart.px - px)/B;
146  refpart.t = t - (theta/B)*pt;
147 
148  // advance integrated path length
149  refpart.s = s + slice_ds;
150 
151  }
152 
158  int nslice () const
159  {
160  return m_nslice;
161  }
162 
168  amrex::ParticleReal ds () const
169  {
170  return m_ds;
171  }
172 
173  private:
174  amrex::ParticleReal m_ds;
175  amrex::ParticleReal m_rc;
176  int m_nslice;
177  };
178 
179 } // namespace impactx
180 
181 #endif // IMPACTX_SBEND_H
def x
Definition: ImpactX.cpp:31
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:49
amrex::ParticleReal pt
energy deviation, normalized by rest energy
Definition: ReferenceParticle.H:38
def z
int m_nslice
bend radius in m
Definition: Sbend.H:176
#define AMREX_FORCE_INLINE
static constexpr auto name
Definition: Sbend.H:25
Definition: Sbend.H:23
#define AMREX_GPU_HOST_DEVICE
amrex::ParticleReal m_rc
segment length in m
Definition: Sbend.H:175
amrex::ParticleReal m_ds
Definition: Sbend.H:174
Definition: ReferenceParticle.H:28
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > pow(const GpuComplex< T > &a_z, const T &a_y) noexcept
s
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int nslice() const
Definition: Sbend.H:158
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal ds() const
Definition: Sbend.H:168
Sbend(amrex::ParticleReal const ds, amrex::ParticleReal const rc, int const nslice)
Definition: Sbend.H:34