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