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/beamoptic.H"
15 #include "mixin/thick.H"
16 
17 #include <AMReX_Extension.H>
18 #include <AMReX_REAL.H>
19 
20 #include <cmath>
21 
22 
23 namespace impactx
24 {
25  struct Sbend
26  : public elements::BeamOptic<Sbend>,
27  public elements::Thick
28  {
29  static constexpr auto name = "Sbend";
31 
38  Sbend( amrex::ParticleReal const ds, amrex::ParticleReal const rc,
39  int const nslice)
40  : Thick(ds, nslice), m_rc(rc)
41  {
42  }
43 
45  using BeamOptic::operator();
46 
57  PType& AMREX_RESTRICT p,
58  amrex::ParticleReal & AMREX_RESTRICT px,
59  amrex::ParticleReal & AMREX_RESTRICT py,
60  amrex::ParticleReal & AMREX_RESTRICT pt,
61  RefPart const & refpart) const {
62 
63  using namespace amrex::literals; // for _rt and _prt
64 
65  // access AoS data such as positions and cpu/id
66  amrex::ParticleReal const x = p.pos(0);
67  amrex::ParticleReal const y = p.pos(1);
68  amrex::ParticleReal const t = p.pos(2);
69 
70  // initialize output values of momenta
71  amrex::ParticleReal pxout = px;
72  amrex::ParticleReal pyout = py;
73  amrex::ParticleReal ptout = pt;
74 
75  // length of the current slice
76  amrex::ParticleReal const slice_ds = m_ds / nslice();
77 
78  // access reference particle values to find beta*gamma^2
79  amrex::ParticleReal const pt_ref = refpart.pt;
80  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
81  amrex::ParticleReal const bet = sqrt(betgam2/(1.0_prt + betgam2));
82 
83  // calculate expensive terms once
84  // TODO: use sincos function once wrapped in AMReX
85  amrex::ParticleReal const theta = slice_ds/m_rc;
86  amrex::ParticleReal const sin_theta = sin(theta);
87  amrex::ParticleReal const cos_theta = cos(theta);
88 
89  // advance position and momentum (sector bend)
90  p.pos(0) = cos_theta*x + m_rc*sin_theta*px
91  - (m_rc/bet)*(1.0_prt - cos_theta)*pt;
92 
93  pxout = -sin_theta/m_rc*x + cos_theta*px - sin_theta/bet*pt;
94 
95  p.pos(1) = y + m_rc*theta*py;
96 
97  // pyout = py;
98 
99  p.pos(2) = sin_theta/bet*x + m_rc/bet*(1.0_prt - cos_theta)*px + t
100  + m_rc*(-theta+sin_theta/(bet*bet))*pt;
101 
102  // ptout = pt;
103 
104  // assign updated momenta
105  px = pxout;
106  py = pyout;
107  pt = ptout;
108 
109  }
110 
116  void operator() (RefPart & AMREX_RESTRICT refpart) const {
117 
118  using namespace amrex::literals; // for _rt and _prt
119 
120  // assign input reference particle values
121  amrex::ParticleReal const x = refpart.x;
122  amrex::ParticleReal const px = refpart.px;
123  amrex::ParticleReal const y = refpart.y;
124  amrex::ParticleReal const py = refpart.py;
125  amrex::ParticleReal const z = refpart.z;
126  amrex::ParticleReal const pz = refpart.pz;
127  amrex::ParticleReal const t = refpart.t;
128  amrex::ParticleReal const pt = refpart.pt;
129  amrex::ParticleReal const s = refpart.s;
130 
131  // length of the current slice
132  amrex::ParticleReal const slice_ds = m_ds / nslice();
133 
134  // assign intermediate parameter
135  amrex::ParticleReal const theta = slice_ds/m_rc;
136  amrex::ParticleReal const B = sqrt(pow(pt,2)-1.0_prt)/m_rc;
137 
138  // calculate expensive terms once
139  // TODO: use sincos function once wrapped in AMReX
140  amrex::ParticleReal const sin_theta = sin(theta);
141  amrex::ParticleReal const cos_theta = cos(theta);
142 
143  // advance position and momentum (bend)
144  refpart.px = px*cos_theta - pz*sin_theta;
145  refpart.py = py;
146  refpart.pz = pz*cos_theta + px*sin_theta;
147  refpart.pt = pt;
148 
149  refpart.x = x + (refpart.pz - pz)/B;
150  refpart.y = y + (theta/B)*py;
151  refpart.z = z - (refpart.px - px)/B;
152  refpart.t = t - (theta/B)*pt;
153 
154  // advance integrated path length
155  refpart.s = s + slice_ds;
156 
157  }
158 
159  private:
160  amrex::ParticleReal m_rc;
161  };
162 
163 } // namespace impactx
164 
165 #endif // IMPACTX_SBEND_H
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:56
Definition: ImpactX.cpp:31
Definition: beamoptic.H:134
amrex::ParticleReal pt
energy deviation, normalized by rest energy
Definition: ReferenceParticle.H:39
#define AMREX_FORCE_INLINE
static constexpr auto name
Definition: Sbend.H:29
Definition: Sbend.H:25
#define AMREX_GPU_HOST_DEVICE
amrex::ParticleReal m_rc
Definition: Sbend.H:160
Definition: thick.H:23
Thick(amrex::ParticleReal const ds, int const nslice)
Definition: thick.H:30
Definition: ReferenceParticle.H:29
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 amrex::ParticleReal ds() const
Definition: thick.H:50
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int nslice() const
Definition: thick.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
Sbend(amrex::ParticleReal const ds, amrex::ParticleReal const rc, int const nslice)
Definition: Sbend.H:38
amrex::ParticleReal m_ds
Definition: thick.H:56