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 
15 #include <AMReX_Extension.H>
16 #include <AMReX_REAL.H>
17 
18 #include <cmath>
19 
20 
21 namespace impactx
22 {
23  struct ConstF
24  {
25  static constexpr auto name = "ConstF";
27 
36  ConstF( amrex::ParticleReal const ds, amrex::ParticleReal const kx,
37  amrex::ParticleReal const ky, amrex::ParticleReal const kt,
38  int const nslice )
39  : m_ds(ds), m_kx(kx), m_ky(ky), m_kt(kt), m_nslice(nslice)
40  {
41  }
42 
54  PType& AMREX_RESTRICT p,
55  amrex::ParticleReal & AMREX_RESTRICT px,
56  amrex::ParticleReal & AMREX_RESTRICT py,
57  amrex::ParticleReal & AMREX_RESTRICT pt,
58  RefPart const & refpart) const {
59 
60  using namespace amrex::literals; // for _rt and _prt
61 
62  // access AoS data such as positions and cpu/id
63  amrex::ParticleReal const x = p.pos(0);
64  amrex::ParticleReal const y = p.pos(1);
65  amrex::ParticleReal const t = p.pos(2);
66 
67  // access reference particle values to find beta*gamma^2
68  amrex::ParticleReal const pt_ref = refpart.pt;
69  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
70 
71  // intialize output values of momenta
72  amrex::ParticleReal pxout = px;
73  amrex::ParticleReal pyout = py;
74  amrex::ParticleReal ptout = pt;
75 
76  // length of the current slice
77  amrex::ParticleReal const slice_ds = m_ds / nslice();
78 
79  // advance position and momentum
80  p.pos(0) = cos(m_kx*slice_ds)*x + sin(m_kx*slice_ds)/m_kx*px;
81  pxout = -m_kx*sin(m_kx*slice_ds)*x + cos(m_kx*slice_ds)*px;
82 
83  p.pos(1) = cos(m_ky*slice_ds)*y + sin(m_ky*slice_ds)/m_ky*py;
84  pyout = -m_ky*sin(m_ky*slice_ds)*y + cos(m_ky*slice_ds)*py;
85 
86  p.pos(2) = cos(m_kt*slice_ds)*t + sin(m_kt*slice_ds)/(betgam2*m_kt)*pt;
87  ptout = -(m_kt*betgam2)*sin(m_kt*slice_ds)*t + cos(m_kt*slice_ds)*pt;
88 
89  // assign updated momenta
90  px = pxout;
91  py = pyout;
92  pt = ptout;
93 
94  }
95 
101  void operator() (RefPart & AMREX_RESTRICT refpart) const {
102 
103  using namespace amrex::literals; // for _rt and _prt
104 
105  // assign input reference particle values
106  amrex::ParticleReal const x = refpart.x;
107  amrex::ParticleReal const px = refpart.px;
108  amrex::ParticleReal const y = refpart.y;
109  amrex::ParticleReal const py = refpart.py;
110  amrex::ParticleReal const z = refpart.z;
111  amrex::ParticleReal const pz = refpart.pz;
112  amrex::ParticleReal const t = refpart.t;
113  amrex::ParticleReal const pt = refpart.pt;
114  amrex::ParticleReal const s = refpart.s;
115 
116  // length of the current slice
117  amrex::ParticleReal const slice_ds = m_ds / nslice();
118 
119  // assign intermediate parameter
120  amrex::ParticleReal const step = slice_ds / sqrt(pow(pt, 2)-1.0_prt);
121 
122  // advance position and momentum (straight element)
123  refpart.x = x + step*px;
124  refpart.y = y + step*py;
125  refpart.z = z + step*pz;
126  refpart.t = t - step*pt;
127 
128  // advance integrated path length
129  refpart.s = s + slice_ds;
130 
131  }
132 
138  int nslice () const
139  {
140  return m_nslice;
141  }
142 
148  amrex::ParticleReal ds () const
149  {
150  return m_ds;
151  }
152 
153  private:
154  amrex::ParticleReal m_ds;
155  amrex::ParticleReal m_kx;
156  amrex::ParticleReal m_ky;
157  amrex::ParticleReal m_kt;
158  int m_nslice;
159  };
160 
161 } // namespace impactx
162 
163 #endif // IMPACTX_CONSTF_H
amrex::ParticleReal m_ky
focusing x strength in 1/m
Definition: ConstF.H:156
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:53
Definition: ImpactX.cpp:31
static constexpr auto name
Definition: ConstF.H:25
Definition: ConstF.H:23
int m_nslice
focusing t strength in 1/m
Definition: ConstF.H:158
amrex::ParticleReal pt
energy deviation, normalized by rest energy
Definition: ReferenceParticle.H:39
#define AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal ds() const
Definition: ConstF.H:148
#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:36
amrex::ParticleReal m_kt
focusing y strength in 1/m
Definition: ConstF.H:157
amrex::ParticleReal m_kx
segment length in m
Definition: ConstF.H:155
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 int nslice() const
Definition: ConstF.H:138
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
amrex::ParticleReal m_ds
Definition: ConstF.H:154