ImpactX
Loading...
Searching...
No Matches
Drift.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_DRIFT_H
11#define IMPACTX_DRIFT_H
12
14#include "mixin/alignment.H"
15#include "mixin/beamoptic.H"
17#include "mixin/spintransport.H"
18#include "mixin/named.H"
19#include "mixin/nofinalize.H"
20#include "mixin/pipeaperture.H"
21#include "mixin/thick.H"
22
23#include <AMReX_Extension.H>
24#include <AMReX_Math.H>
25#include <AMReX_REAL.H>
26#include <AMReX_SIMD.H>
27#include <AMReX_SmallMatrix.H>
28
29#include <cmath>
30#include <utility>
31
32
33namespace impactx::elements
34{
35 struct Drift
36 : public mixin::Named,
37 public mixin::BeamOptic<Drift>,
38 public mixin::LinearTransport<Drift>,
39 public mixin::Thick,
40 public mixin::Alignment,
43 public mixin::NoFinalize,
44 public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
45 {
46 static constexpr auto type = "Drift";
48
64 amrex::ParticleReal rotation_degree = 0,
67 int nslice = 1,
68 std::optional<std::string> name = std::nullopt
69 )
70 : Named(std::move(name)),
71 Thick(ds, nslice),
72 Alignment(dx, dy, rotation_degree),
74 {
75 }
76
78 void reverse () { Thick::reverse(); }
79
81 using BeamOptic::operator();
82
90 void compute_constants (RefPart const & refpart)
91 {
92 using namespace amrex::literals; // for _rt and _prt
94
95 Alignment::compute_constants(refpart);
96
97 // length of the current slice
99
100 // find beta*gamma^2
101 m_betgam2 = powi<2>(refpart.pt) - 1.0_prt;
102
104 }
105
119 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
122 T_Real & AMREX_RESTRICT x,
123 T_Real & AMREX_RESTRICT y,
124 T_Real & AMREX_RESTRICT t,
125 T_Real const & AMREX_RESTRICT px,
126 T_Real const & AMREX_RESTRICT py,
127 T_Real const & AMREX_RESTRICT pt,
128 [[maybe_unused]] T_IdCpu const & AMREX_RESTRICT idcpu,
129 [[maybe_unused]] RefPart const & AMREX_RESTRICT refpart
130 ) const
131 {
132 using namespace amrex::literals; // for _rt and _prt
133
134 // initialize output values
135 T_Real xout = x;
136 T_Real yout = y;
137 T_Real tout = t;
138 // T_Real pxout = px;
139 // T_Real pyout = py;
140 // T_Real ptout = pt;
141
142 // advance position and momentum (drift)
143 xout = x + m_slice_ds * px;
144 // pxout = px;
145 yout = y + m_slice_ds * py;
146 // pyout = py;
147 tout = t + m_slice_bg * pt;
148 // ptout = pt;
149
150 // assign updated values
151 x = xout;
152 y = yout;
153 t = tout;
154 // px = pxout;
155 // py = pyout;
156 // pt = ptout;
157 }
158
164 void operator() (RefPart & AMREX_RESTRICT refpart) const
165 {
166 using namespace amrex::literals; // for _rt and _prt
167 using amrex::Math::powi;
168
169 // assign input reference particle values
170 amrex::ParticleReal const x = refpart.x;
171 amrex::ParticleReal const px = refpart.px;
172 amrex::ParticleReal const y = refpart.y;
173 amrex::ParticleReal const py = refpart.py;
174 amrex::ParticleReal const z = refpart.z;
175 amrex::ParticleReal const pz = refpart.pz;
176 amrex::ParticleReal const t = refpart.t;
177 amrex::ParticleReal const pt = refpart.pt;
178 amrex::ParticleReal const s = refpart.s;
179
180 // length of the current slice
181 amrex::ParticleReal const slice_ds = m_ds / nslice();
182
183 // assign intermediate parameter
184 amrex::ParticleReal const step = slice_ds / std::sqrt(powi<2>(pt)-1.0_prt);
185
186 // advance position and momentum (drift)
187 refpart.x = x + step*px;
188 refpart.y = y + step*py;
189 refpart.z = z + step*pz;
190 refpart.t = t - step*pt;
191
192 // advance integrated path length
193 refpart.s = s + slice_ds;
194 }
195
196
211 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
214 T_Real & AMREX_RESTRICT x,
215 T_Real & AMREX_RESTRICT y,
216 T_Real & AMREX_RESTRICT t,
217 T_Real const & AMREX_RESTRICT px,
218 T_Real const & AMREX_RESTRICT py,
219 T_Real const & AMREX_RESTRICT pt,
220 [[maybe_unused]] T_Real const & AMREX_RESTRICT sx,
221 [[maybe_unused]] T_Real const & AMREX_RESTRICT sy,
222 [[maybe_unused]] T_Real const & AMREX_RESTRICT sz,
223 T_IdCpu const & AMREX_RESTRICT idcpu,
224 RefPart const & AMREX_RESTRICT refpart
225 ) const
226 {
227 // spin is unaffected
228
229 // phase space push
230 (*this)(x, y, t, px, py, pt, idcpu, refpart);
231
232 }
233
234
236 using LinearTransport::operator();
237
243 Map6x6
244 transport_map (RefPart const & AMREX_RESTRICT refpart) const
245 {
246 using namespace amrex::literals; // for _rt and _prt
247 using amrex::Math::powi;
248
249 // length of the current slice
250 amrex::ParticleReal const slice_ds = m_ds / nslice();
251
252 // access reference particle values to find beta*gamma^2
253 amrex::ParticleReal const pt_ref = refpart.pt;
254 amrex::ParticleReal const betgam2 = powi<2>(pt_ref) - 1.0_prt;
255
256 // assign linear map matrix elements
258 R(1,2) = slice_ds;
259 R(3,4) = slice_ds;
260 R(5,6) = slice_ds / betgam2;
261
262 // apply the transverse rotation (roll) alignment error
263 return rotate_aligned_map(R);
264 }
265
266 private:
267 // constants that are independent of the individually tracked particle,
268 // see: compute_constants() to refresh
272 };
273
274} // namespace impactx
275
277
278#endif // IMPACTX_DRIFT_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
#define AMREX_GPU_HOST
#define IMPACTX_PUSH_EXTERN_TEMPLATE(ElementType)
Definition PushAll.H:78
amrex_particle_real ParticleReal
constexpr T powi(T x) noexcept
Definition All.H:56
@ s
fixed s as the independent variable
Definition ImpactXParticleContainer.H:37
@ t
fixed t as the independent variable
Definition ImpactXParticleContainer.H:38
amrex::SmallMatrix< amrex::ParticleReal, 6, 6, amrex::Order::F, 1 > Map6x6
Definition CovarianceMatrix.H:20
static constexpr __host__ __device__ SmallMatrix< T, NRows, NCols, ORDER, StartIndex > Identity() noexcept
Definition ReferenceParticle.H:33
amrex::ParticleReal pt
energy, normalized by rest energy
Definition ReferenceParticle.H:42
Definition Drift.H:45
amrex::ParticleReal m_betgam2
m_ds / nslice();
Definition Drift.H:270
void reverse()
Definition Drift.H:78
ImpactXParticleContainer::ParticleType PType
Definition Drift.H:47
Drift(amrex::ParticleReal ds, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, amrex::ParticleReal aperture_x=0, amrex::ParticleReal aperture_y=0, int nslice=1, std::optional< std::string > name=std::nullopt)
Definition Drift.H:60
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT t, T_Real const &AMREX_RESTRICT px, T_Real const &AMREX_RESTRICT py, T_Real const &AMREX_RESTRICT pt, T_IdCpu const &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition Drift.H:121
amrex::ParticleReal m_slice_ds
Definition Drift.H:269
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void spin_and_phasespace_push(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT t, T_Real const &AMREX_RESTRICT px, T_Real const &AMREX_RESTRICT py, T_Real const &AMREX_RESTRICT pt, T_Real const &AMREX_RESTRICT sx, T_Real const &AMREX_RESTRICT sy, T_Real const &AMREX_RESTRICT sz, T_IdCpu const &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition Drift.H:213
amrex::ParticleReal m_slice_bg
beta*gamma^2
Definition Drift.H:271
void compute_constants(RefPart const &refpart)
Definition Drift.H:90
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition Drift.H:244
static constexpr auto type
Definition Drift.H:46
Definition alignment.H:29
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dy() const
Definition alignment.H:189
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dx() const
Definition alignment.H:179
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 rotate_aligned_map(Map6x6 const &R) const
Definition alignment.H:263
Alignment(amrex::ParticleReal dx, amrex::ParticleReal dy, amrex::ParticleReal rotation_degree)
Definition alignment.H:39
Definition beamoptic.H:529
Definition lineartransport.H:50
Definition named.H:29
AMREX_GPU_HOST Named(std::optional< std::string > name)
Definition named.H:57
AMREX_FORCE_INLINE std::string name() const
Definition named.H:122
Definition nofinalize.H:22
Definition pipeaperture.H:26
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal aperture_x() const
Definition pipeaperture.H:90
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal aperture_y() const
Definition pipeaperture.H:101
PipeAperture(amrex::ParticleReal aperture_x, amrex::ParticleReal aperture_y)
Definition pipeaperture.H:32
Definition spintransport.H:36
Definition thick.H:24
Thick(amrex::ParticleReal ds, int nslice)
Definition thick.H:30
amrex::ParticleReal m_ds
Definition thick.H:68
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