ImpactX
ConstF.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_CONSTF_H
11 #define IMPACTX_CONSTF_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_REAL.H>
21 
22 #include <cmath>
23 
24 
25 namespace impactx
26 {
27  struct ConstF
28  : public elements::BeamOptic<ConstF>,
29  public elements::Thick,
30  public elements::Alignment,
32  {
33  static constexpr auto name = "ConstF";
35 
48  amrex::ParticleReal ds,
49  amrex::ParticleReal kx,
50  amrex::ParticleReal ky,
51  amrex::ParticleReal kt,
52  amrex::ParticleReal dx = 0,
53  amrex::ParticleReal dy = 0,
54  amrex::ParticleReal rotation_degree = 0,
55  int nslice = 1
56  )
57  : Thick(ds, nslice),
58  Alignment(dx, dy, rotation_degree),
59  m_kx(kx), m_ky(ky), m_kt(kt)
60  {
61  }
62 
64  using BeamOptic::operator();
65 
76  PType& AMREX_RESTRICT p,
77  amrex::ParticleReal & AMREX_RESTRICT px,
78  amrex::ParticleReal & AMREX_RESTRICT py,
79  amrex::ParticleReal & AMREX_RESTRICT pt,
80  RefPart const & refpart) const {
81 
82  using namespace amrex::literals; // for _rt and _prt
83 
84  // shift due to alignment errors of the element
85  shift_in(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
86 
87  // access AoS data such as positions and cpu/id
88  amrex::ParticleReal const x = p.pos(RealAoS::x);
89  amrex::ParticleReal const y = p.pos(RealAoS::y);
90  amrex::ParticleReal const t = p.pos(RealAoS::t);
91 
92  // access reference particle values to find beta*gamma^2
93  amrex::ParticleReal const pt_ref = refpart.pt;
94  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
95 
96  // intialize output values of momenta
97  amrex::ParticleReal pxout = px;
98  amrex::ParticleReal pyout = py;
99  amrex::ParticleReal ptout = pt;
100 
101  // length of the current slice
102  amrex::ParticleReal const slice_ds = m_ds / nslice();
103 
104  // advance position and momentum
105  p.pos(RealAoS::x) = cos(m_kx*slice_ds)*x + sin(m_kx*slice_ds)/m_kx*px;
106  pxout = -m_kx*sin(m_kx*slice_ds)*x + cos(m_kx*slice_ds)*px;
107 
108  p.pos(RealAoS::y) = cos(m_ky*slice_ds)*y + sin(m_ky*slice_ds)/m_ky*py;
109  pyout = -m_ky*sin(m_ky*slice_ds)*y + cos(m_ky*slice_ds)*py;
110 
111  p.pos(RealAoS::t) = cos(m_kt*slice_ds)*t + sin(m_kt*slice_ds)/(betgam2*m_kt)*pt;
112  ptout = -(m_kt*betgam2)*sin(m_kt*slice_ds)*t + cos(m_kt*slice_ds)*pt;
113 
114  // assign updated momenta
115  px = pxout;
116  py = pyout;
117  pt = ptout;
118 
119  // undo shift due to alignment errors of the element
120  shift_out(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
121  }
122 
128  void operator() (RefPart & AMREX_RESTRICT refpart) const {
129 
130  using namespace amrex::literals; // for _rt and _prt
131 
132  // assign input reference particle values
133  amrex::ParticleReal const x = refpart.x;
134  amrex::ParticleReal const px = refpart.px;
135  amrex::ParticleReal const y = refpart.y;
136  amrex::ParticleReal const py = refpart.py;
137  amrex::ParticleReal const z = refpart.z;
138  amrex::ParticleReal const pz = refpart.pz;
139  amrex::ParticleReal const t = refpart.t;
140  amrex::ParticleReal const pt = refpart.pt;
141  amrex::ParticleReal const s = refpart.s;
142 
143  // length of the current slice
144  amrex::ParticleReal const slice_ds = m_ds / nslice();
145 
146  // assign intermediate parameter
147  amrex::ParticleReal const step = slice_ds / sqrt(pow(pt, 2)-1.0_prt);
148 
149  // advance position and momentum (straight element)
150  refpart.x = x + step*px;
151  refpart.y = y + step*py;
152  refpart.z = z + step*pz;
153  refpart.t = t - step*pt;
154 
155  // advance integrated path length
156  refpart.s = s + slice_ds;
157 
158  }
159 
160  amrex::ParticleReal m_kx;
161  amrex::ParticleReal m_ky;
162  amrex::ParticleReal m_kt;
163  };
164 
165 } // namespace impactx
166 
167 #endif // IMPACTX_CONSTF_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:36
@ t
fixed t as the independent variable
s
Definition: ConstF.H:32
ImpactXParticleContainer::ParticleType PType
Definition: ConstF.H:34
amrex::ParticleReal m_ky
focusing x strength in 1/m
Definition: ConstF.H:161
ConstF(amrex::ParticleReal ds, amrex::ParticleReal kx, amrex::ParticleReal ky, amrex::ParticleReal kt, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, int nslice=1)
Definition: ConstF.H:47
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: ConstF.H:75
static constexpr auto name
Definition: ConstF.H:33
amrex::ParticleReal m_kt
focusing y strength in 1/m
Definition: ConstF.H:162
amrex::ParticleReal m_kx
Definition: ConstF.H:160
@ 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: 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