10 #ifndef IMPACTX_CONSTF_H
11 #define IMPACTX_CONSTF_H
33 static constexpr
auto name =
"ConstF";
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,
64 using BeamOptic::operator();
79 amrex::ParticleReal & AMREX_RESTRICT x,
80 amrex::ParticleReal & AMREX_RESTRICT y,
81 amrex::ParticleReal & AMREX_RESTRICT t,
82 amrex::ParticleReal & AMREX_RESTRICT px,
83 amrex::ParticleReal & AMREX_RESTRICT py,
84 amrex::ParticleReal & AMREX_RESTRICT pt,
85 [[maybe_unused]] uint64_t & AMREX_RESTRICT idcpu,
86 RefPart const & refpart)
const {
88 using namespace amrex::literals;
94 amrex::ParticleReal
const pt_ref = refpart.
pt;
95 amrex::ParticleReal
const betgam2 =
pow(pt_ref, 2) - 1.0_prt;
98 amrex::ParticleReal xout =
x;
99 amrex::ParticleReal yout = y;
100 amrex::ParticleReal tout =
t;
101 amrex::ParticleReal pxout = px;
102 amrex::ParticleReal pyout = py;
103 amrex::ParticleReal ptout = pt;
106 amrex::ParticleReal
const slice_ds =
m_ds /
nslice();
112 yout = cos(
m_ky*slice_ds)*y + sin(
m_ky*slice_ds)/
m_ky*py;
113 pyout = -
m_ky*sin(
m_ky*slice_ds)*y + cos(
m_ky*slice_ds)*py;
115 tout = cos(
m_kt*slice_ds)*
t + sin(
m_kt*slice_ds)/(betgam2*
m_kt)*pt;
116 ptout = -(
m_kt*betgam2)*sin(
m_kt*slice_ds)*
t + cos(
m_kt*slice_ds)*pt;
137 using namespace amrex::literals;
140 amrex::ParticleReal
const x = refpart.x;
141 amrex::ParticleReal
const px = refpart.px;
142 amrex::ParticleReal
const y = refpart.y;
143 amrex::ParticleReal
const py = refpart.py;
144 amrex::ParticleReal
const z = refpart.z;
145 amrex::ParticleReal
const pz = refpart.pz;
146 amrex::ParticleReal
const t = refpart.t;
147 amrex::ParticleReal
const pt = refpart.pt;
148 amrex::ParticleReal
const s = refpart.s;
151 amrex::ParticleReal
const slice_ds =
m_ds /
nslice();
154 amrex::ParticleReal
const step = slice_ds /
sqrt(
pow(pt, 2)-1.0_prt);
157 refpart.x =
x + step*px;
158 refpart.y = y + step*py;
159 refpart.z = z + step*pz;
160 refpart.t =
t - step*pt;
163 refpart.s =
s + slice_ds;
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
T_ParticleType ParticleType
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
@ t
fixed t as the independent variable
ImpactXParticleContainer::ParticleType PType
Definition: ConstF.H:34
amrex::ParticleReal m_ky
focusing x strength in 1/m
Definition: ConstF.H:168
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()(amrex::ParticleReal &AMREX_RESTRICT x, amrex::ParticleReal &AMREX_RESTRICT y, amrex::ParticleReal &AMREX_RESTRICT t, amrex::ParticleReal &AMREX_RESTRICT px, amrex::ParticleReal &AMREX_RESTRICT py, amrex::ParticleReal &AMREX_RESTRICT pt, [[maybe_unused]] uint64_t &AMREX_RESTRICT idcpu, RefPart const &refpart) const
Definition: ConstF.H:78
static constexpr auto name
Definition: ConstF.H:33
amrex::ParticleReal m_kt
focusing y strength in 1/m
Definition: ConstF.H:169
amrex::ParticleReal m_kx
Definition: ConstF.H:167
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:149
Definition: nofinalize.H:22
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