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