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/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 ConstF
26  : public elements::BeamOptic<ConstF>,
27  public elements::Thick
28  {
29  static constexpr auto name = "ConstF";
31 
40  ConstF( amrex::ParticleReal const ds, amrex::ParticleReal const kx,
41  amrex::ParticleReal const ky, amrex::ParticleReal const kt,
42  int const nslice )
43  : Thick(ds, nslice), m_kx(kx), m_ky(ky), m_kt(kt)
44  {
45  }
46 
48  using BeamOptic::operator();
49 
60  PType& AMREX_RESTRICT p,
61  amrex::ParticleReal & AMREX_RESTRICT px,
62  amrex::ParticleReal & AMREX_RESTRICT py,
63  amrex::ParticleReal & AMREX_RESTRICT pt,
64  RefPart const & refpart) const {
65 
66  using namespace amrex::literals; // for _rt and _prt
67 
68  // access AoS data such as positions and cpu/id
69  amrex::ParticleReal const x = p.pos(0);
70  amrex::ParticleReal const y = p.pos(1);
71  amrex::ParticleReal const t = p.pos(2);
72 
73  // access reference particle values to find beta*gamma^2
74  amrex::ParticleReal const pt_ref = refpart.pt;
75  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
76 
77  // intialize output values of momenta
78  amrex::ParticleReal pxout = px;
79  amrex::ParticleReal pyout = py;
80  amrex::ParticleReal ptout = pt;
81 
82  // length of the current slice
83  amrex::ParticleReal const slice_ds = m_ds / nslice();
84 
85  // advance position and momentum
86  p.pos(0) = cos(m_kx*slice_ds)*x + sin(m_kx*slice_ds)/m_kx*px;
87  pxout = -m_kx*sin(m_kx*slice_ds)*x + cos(m_kx*slice_ds)*px;
88 
89  p.pos(1) = cos(m_ky*slice_ds)*y + sin(m_ky*slice_ds)/m_ky*py;
90  pyout = -m_ky*sin(m_ky*slice_ds)*y + cos(m_ky*slice_ds)*py;
91 
92  p.pos(2) = cos(m_kt*slice_ds)*t + sin(m_kt*slice_ds)/(betgam2*m_kt)*pt;
93  ptout = -(m_kt*betgam2)*sin(m_kt*slice_ds)*t + cos(m_kt*slice_ds)*pt;
94 
95  // assign updated momenta
96  px = pxout;
97  py = pyout;
98  pt = ptout;
99 
100  }
101 
107  void operator() (RefPart & AMREX_RESTRICT refpart) const {
108 
109  using namespace amrex::literals; // for _rt and _prt
110 
111  // assign input reference particle values
112  amrex::ParticleReal const x = refpart.x;
113  amrex::ParticleReal const px = refpart.px;
114  amrex::ParticleReal const y = refpart.y;
115  amrex::ParticleReal const py = refpart.py;
116  amrex::ParticleReal const z = refpart.z;
117  amrex::ParticleReal const pz = refpart.pz;
118  amrex::ParticleReal const t = refpart.t;
119  amrex::ParticleReal const pt = refpart.pt;
120  amrex::ParticleReal const s = refpart.s;
121 
122  // length of the current slice
123  amrex::ParticleReal const slice_ds = m_ds / nslice();
124 
125  // assign intermediate parameter
126  amrex::ParticleReal const step = slice_ds / sqrt(pow(pt, 2)-1.0_prt);
127 
128  // advance position and momentum (straight element)
129  refpart.x = x + step*px;
130  refpart.y = y + step*py;
131  refpart.z = z + step*pz;
132  refpart.t = t - step*pt;
133 
134  // advance integrated path length
135  refpart.s = s + slice_ds;
136 
137  }
138 
139  private:
140  amrex::ParticleReal m_kx;
141  amrex::ParticleReal m_ky;
142  amrex::ParticleReal m_kt;
143  };
144 
145 } // namespace impactx
146 
147 #endif // IMPACTX_CONSTF_H
amrex::ParticleReal m_ky
focusing x strength in 1/m
Definition: ConstF.H:141
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:59
Definition: ImpactX.cpp:31
static constexpr auto name
Definition: ConstF.H:29
Definition: ConstF.H:25
Definition: beamoptic.H:134
amrex::ParticleReal pt
energy deviation, normalized by rest energy
Definition: ReferenceParticle.H:39
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
ConstF(amrex::ParticleReal const ds, amrex::ParticleReal const kx, amrex::ParticleReal const ky, amrex::ParticleReal const kt, int const nslice)
Definition: ConstF.H:40
Definition: thick.H:23
amrex::ParticleReal m_kt
focusing y strength in 1/m
Definition: ConstF.H:142
amrex::ParticleReal m_kx
Definition: ConstF.H:140
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
amrex::ParticleReal m_ds
Definition: thick.H:56