ImpactX
Loading...
Searching...
No Matches
LinearMap.H
Go to the documentation of this file.
1/* Copyright 2022-2024 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_ELEMENT_LINEAR_MAP_H
11#define IMPACTX_ELEMENT_LINEAR_MAP_H
12
14#include "mixin/alignment.H"
15#include "mixin/beamoptic.H"
17#include "mixin/named.H"
18#include "mixin/nofinalize.H"
19
20#include <AMReX_Extension.H>
21#include <AMReX_Math.H>
22#include <AMReX_REAL.H>
23#include <AMReX_SIMD.H>
24
25#include <cmath>
26#include <stdexcept>
27
28
29namespace impactx::elements
30{
31 struct LinearMap
32 : public mixin::Named,
33 public mixin::BeamOptic<LinearMap>,
34 public mixin::LinearTransport<LinearMap>,
35 public mixin::Alignment,
36 public mixin::NoFinalize,
37 public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
38 {
39 static constexpr auto type = "LinearMap";
41
55 Map6x6 const & R,
56 amrex::ParticleReal ds = 0,
57 amrex::ParticleReal dx = 0,
58 amrex::ParticleReal dy = 0,
59 amrex::ParticleReal rotation_degree = 0,
60 std::optional<std::string> name = std::nullopt
61 )
62 : Named(std::move(name)),
63 Alignment(dx, dy, rotation_degree)
64 {
66 m_ds = ds;
67 }
68
70 using BeamOptic::operator();
71
84 template<typename T_Real=amrex::ParticleReal, typename T_IdCpu=uint64_t>
87 T_Real & AMREX_RESTRICT x,
88 T_Real & AMREX_RESTRICT y,
89 T_Real & AMREX_RESTRICT t,
90 T_Real & AMREX_RESTRICT px,
91 T_Real & AMREX_RESTRICT py,
92 T_Real & AMREX_RESTRICT pt,
93 [[maybe_unused]] T_IdCpu & AMREX_RESTRICT idcpu,
94 [[maybe_unused]] RefPart const & AMREX_RESTRICT refpart
95 ) const
96 {
97 using namespace amrex::literals; // for _rt and _prt
98
99 // shift due to alignment errors of the element
100 shift_in(x, y, px, py);
101
102 // input and output phase space vectors
104 x, px, y, py, t, pt
105 };
106
107 // TODO amrex::SmallMatrix::operator* overloads on T for SIMD*scalar / scalar*SIMD
108 amrex::SmallVector<T_Real, 6, 1> const vectorout = m_transport_map * vectorin;
109
110 // assign updated values
111 x = vectorout(1);
112 px = vectorout(2);
113 y = vectorout(3);
114 py = vectorout(4);
115 t = vectorout(5);
116 pt = vectorout(6);
117
118 // undo shift due to alignment errors of the element
119 shift_out(x, y, px, py);
120 }
121
127 void operator() (RefPart & AMREX_RESTRICT refpart) const
128 {
129 if (m_ds > 0) // Drift
130 {
131 using namespace amrex::literals; // for _rt and _prt
132 using amrex::Math::powi;
133
134 // assign input reference particle values
135 amrex::ParticleReal const x = refpart.x;
136 amrex::ParticleReal const px = refpart.px;
137 amrex::ParticleReal const y = refpart.y;
138 amrex::ParticleReal const py = refpart.py;
139 amrex::ParticleReal const z = refpart.z;
140 amrex::ParticleReal const pz = refpart.pz;
141 amrex::ParticleReal const t = refpart.t;
142 amrex::ParticleReal const pt = refpart.pt;
143 amrex::ParticleReal const s = refpart.s;
144
145 // length of the current slice
146 amrex::ParticleReal const slice_ds = m_ds / nslice();
147
148 // assign intermediate parameter
149 amrex::ParticleReal const step = slice_ds / std::sqrt(powi<2>(pt)-1.0_prt);
150
151 // advance position and momentum (drift)
152 refpart.x = x + step*px;
153 refpart.y = y + step*py;
154 refpart.z = z + step*pz;
155 refpart.t = t - step*pt;
156
157 // advance integrated path length
158 refpart.s = s + slice_ds;
159 }
160 // else nothing to do for a zero-length element
161 }
162
168 int nslice () const
169 {
170 return 1;
171 }
172
178 amrex::ParticleReal ds () const
179 {
180 return m_ds;
181 }
182
184 using LinearTransport::operator();
185
192 Map6x6
193 transport_map ([[maybe_unused]] RefPart const & AMREX_RESTRICT refpart) const
194 {
196 R = m_transport_map;
197
198 return R;
199 }
200
201 Map6x6 m_transport_map; // 6x6 transport map
202 amrex::ParticleReal m_ds; // finite ds allowed for bookkeeping, but we do not allow slicing
203 };
204
205} // namespace impactx
206
207#endif // IMPACTX_ELEMENT_LINEAR_MAP_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
#define AMREX_GPU_HOST
constexpr T powi(T x) noexcept
SmallMatrix< T, N, 1, Order::F, StartIndex > SmallVector
Definition All.H:54
@ 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:31
Map6x6 m_transport_map
Definition LinearMap.H:201
AMREX_GPU_HOST AMREX_FORCE_INLINE Map6x6 transport_map(RefPart const &AMREX_RESTRICT refpart) const
Definition LinearMap.H:193
amrex::ParticleReal m_ds
Definition LinearMap.H:202
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal ds() const
Definition LinearMap.H:178
static constexpr auto type
Definition LinearMap.H:39
ImpactXParticleContainer::ParticleType PType
Definition LinearMap.H:40
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 &AMREX_RESTRICT px, T_Real &AMREX_RESTRICT py, T_Real &AMREX_RESTRICT pt, T_IdCpu &AMREX_RESTRICT idcpu, RefPart const &AMREX_RESTRICT refpart) const
Definition LinearMap.H:86
LinearMap(Map6x6 const &R, amrex::ParticleReal ds=0, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, std::optional< std::string > name=std::nullopt)
Definition LinearMap.H:54
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int nslice() const
Definition LinearMap.H:168
Definition alignment.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void shift_out(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT px, T_Real &AMREX_RESTRICT py) const
Definition alignment.H:109
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dy() const
Definition alignment.H:146
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal dx() const
Definition alignment.H:136
Alignment(amrex::ParticleReal dx, amrex::ParticleReal dy, amrex::ParticleReal rotation_degree)
Definition alignment.H:36
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void shift_in(T_Real &AMREX_RESTRICT x, T_Real &AMREX_RESTRICT y, T_Real &AMREX_RESTRICT px, T_Real &AMREX_RESTRICT py) const
Definition alignment.H:78
Definition beamoptic.H:219
Definition lineartransport.H:29
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