ImpactX
Sol.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_SOL_H
11 #define IMPACTX_SOL_H
12 
14 #include "mixin/alignment.H"
15 #include "mixin/beamoptic.H"
16 #include "mixin/thick.H"
17 #include "mixin/nofinalize.H"
18 
19 #include <AMReX_Extension.H>
20 #include <AMReX_Math.H>
21 #include <AMReX_REAL.H>
22 
23 #include <cmath>
24 
25 
26 namespace impactx
27 {
28  struct Sol
29  : public elements::BeamOptic<Sol>,
30  public elements::Thick,
31  public elements::Alignment,
33  {
34  static constexpr auto name = "Sol";
36 
47  Sol (
48  amrex::ParticleReal ds,
49  amrex::ParticleReal ks,
50  amrex::ParticleReal dx = 0,
51  amrex::ParticleReal dy = 0,
52  amrex::ParticleReal rotation_degree = 0,
53  int nslice = 1
54  )
55  : Thick(ds, nslice),
56  Alignment(dx, dy, rotation_degree),
57  m_ks(ks)
58  {
59  }
60 
62  using BeamOptic::operator();
63 
74  PType & AMREX_RESTRICT p,
75  amrex::ParticleReal & AMREX_RESTRICT px,
76  amrex::ParticleReal & AMREX_RESTRICT py,
77  amrex::ParticleReal & AMREX_RESTRICT pt,
78  RefPart const & refpart
79  ) const
80  {
81  using namespace amrex::literals; // for _rt and _prt
82 
83  // shift due to alignment errors of the element
84  shift_in(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
85 
86  // access AoS data such as positions and cpu/id
87  amrex::ParticleReal const x = p.pos(RealAoS::x);
88  amrex::ParticleReal const y = p.pos(RealAoS::y);
89  amrex::ParticleReal const t = p.pos(RealAoS::t);
90 
91  // length of the current slice
92  amrex::ParticleReal const slice_ds = m_ds / nslice();
93 
94  // access reference particle values to find beta*gamma^2
95  amrex::ParticleReal const pt_ref = refpart.pt;
96  amrex::ParticleReal const betgam2 = pow(pt_ref, 2) - 1.0_prt;
97 
98  // compute phase advance per unit length (in rad/m) and
99  // rotation angle (in rad)
100  amrex::ParticleReal const alpha = m_ks/2.0_prt;
101  amrex::ParticleReal const theta = alpha*slice_ds;
102 
103  // intialize output values
104  amrex::ParticleReal xout = x;
105  amrex::ParticleReal yout = y;
106  amrex::ParticleReal tout = t;
107  amrex::ParticleReal pxout = px;
108  amrex::ParticleReal pyout = py;
109  amrex::ParticleReal ptout = pt;
110 
111  // advance positions and momenta using map for focusing
112  auto const [sin_theta, cos_theta] = amrex::Math::sincos(theta);
113  xout = cos_theta*x + sin_theta/alpha*px;
114  pxout = -alpha*sin_theta*x + cos_theta*px;
115 
116  yout = cos_theta*y + sin_theta/alpha*py;
117  pyout = -alpha*sin_theta*y + cos_theta*py;
118 
119  tout = t + (slice_ds/betgam2)*pt;
120  ptout = pt;
121 
122  // assign intermediate momenta
123  px = pxout;
124  py = pyout;
125  pt = ptout;
126 
127  // advance positions and momenta using map for rotation
128  p.pos(RealAoS::x) = cos_theta*xout + sin_theta*yout;
129  pxout = cos_theta*px + sin_theta*py;
130 
131  p.pos(RealAoS::y) = -sin_theta*xout + cos_theta*yout;
132  pyout = -sin_theta*px + cos_theta*py;
133 
134  p.pos(RealAoS::t) = tout;
135  ptout = pt;
136 
137  // assign updated momenta
138  px = pxout;
139  py = pyout;
140  pt = ptout;
141 
142  // undo shift due to alignment errors of the element
143  shift_out(p.pos(RealAoS::x), p.pos(RealAoS::y), px, py);
144  }
145 
151  void operator() (RefPart & AMREX_RESTRICT refpart) const
152  {
153  using namespace amrex::literals; // for _rt and _prt
154 
155  // assign input reference particle values
156  amrex::ParticleReal const x = refpart.x;
157  amrex::ParticleReal const px = refpart.px;
158  amrex::ParticleReal const y = refpart.y;
159  amrex::ParticleReal const py = refpart.py;
160  amrex::ParticleReal const z = refpart.z;
161  amrex::ParticleReal const pz = refpart.pz;
162  amrex::ParticleReal const t = refpart.t;
163  amrex::ParticleReal const pt = refpart.pt;
164  amrex::ParticleReal const s = refpart.s;
165 
166  // length of the current slice
167  amrex::ParticleReal const slice_ds = m_ds / nslice();
168 
169  // assign intermediate parameter
170  amrex::ParticleReal const step = slice_ds / sqrt(pow(pt,2)-1.0_prt);
171 
172  // advance position and momentum (straight element)
173  refpart.x = x + step*px;
174  refpart.y = y + step*py;
175  refpart.z = z + step*pz;
176  refpart.t = t - step*pt;
177 
178  // advance integrated path length
179  refpart.s = s + slice_ds;
180  }
181 
182  private:
183  amrex::ParticleReal m_ks;
184  };
185 
186 } // namespace impactx
187 
188 #endif // IMPACTX_SOL_H
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE std::pair< double, double > sincos(double x)
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:36
@ t
fixed t as the independent variable
s
@ x
position in x [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:48
@ y
position in y [m] (at fixed s OR fixed t)
Definition: ImpactXParticleContainer.H:49
@ t
c * time-of-flight [m] (at fixed s)
Definition: ImpactXParticleContainer.H:50
Definition: ReferenceParticle.H:30
amrex::ParticleReal pt
energy, normalized by rest energy
Definition: ReferenceParticle.H:39
Definition: Sol.H:33
static constexpr auto name
Definition: Sol.H:34
ImpactXParticleContainer::ParticleType PType
Definition: Sol.H:35
Sol(amrex::ParticleReal ds, amrex::ParticleReal ks, amrex::ParticleReal dx=0, amrex::ParticleReal dy=0, amrex::ParticleReal rotation_degree=0, int nslice=1)
Definition: Sol.H:47
amrex::ParticleReal m_ks
Definition: Sol.H:183
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: Sol.H:73
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:135
Definition: nofinalize.H:22
Definition: thick.H:24
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