#include <ConstF.H>
|
| 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) |
|
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 |
|
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void | operator() (RefPart &AMREX_RESTRICT refpart) const |
|
void | operator() (ImpactXParticleContainer &pc, int step) |
|
void | operator() (ImpactXParticleContainer::iterator &pti, RefPart &AMREX_RESTRICT ref_part) |
|
| Thick (amrex::ParticleReal ds, int nslice) |
|
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int | nslice () const |
|
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal | ds () const |
|
| Alignment (amrex::ParticleReal dx, amrex::ParticleReal dy, amrex::ParticleReal rotation_degree) |
|
| Alignment ()=default |
|
| Alignment (Alignment const &)=default |
|
Alignment & | operator= (Alignment const &)=default |
|
| Alignment (Alignment &&)=default |
|
Alignment & | operator= (Alignment &&rhs)=default |
|
| ~Alignment ()=default |
|
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 |
|
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 |
|
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal | dx () const |
|
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal | dy () const |
|
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::ParticleReal | rotation () const |
|
void | finalize () |
|
◆ PType
◆ ConstF()
impactx::ConstF::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 |
|
) |
| |
|
inline |
A linear Constant Focusing element
- Parameters
-
ds | Segment length in m. |
kx | Focusing strength for x in 1/m. |
ky | Focusing strength for y in 1/m. |
kt | Focusing strength for t in 1/m. |
dx | horizontal translation error in m |
dy | vertical translation error in m |
rotation_degree | rotation error in the transverse plane [degrees] |
nslice | number of slices used for the application of space charge |
◆ operator()() [1/2]
This pushes a single particle, relative to the reference particle
- Parameters
-
x | particle position in x |
y | particle position in y |
t | particle position in t |
px | particle momentum in x |
py | particle momentum in y |
pt | particle momentum in t |
idcpu | particle global index (unused) |
refpart | reference particle |
◆ operator()() [2/2]
This pushes the reference particle.
- Parameters
-
[in,out] | refpart | reference particle |
◆ m_kt
amrex::ParticleReal impactx::ConstF::m_kt |
focusing y strength in 1/m
◆ m_kx
amrex::ParticleReal impactx::ConstF::m_kx |
◆ m_ky
amrex::ParticleReal impactx::ConstF::m_ky |
focusing x strength in 1/m
◆ name
constexpr auto impactx::ConstF::name = "ConstF" |
|
staticconstexpr |
The documentation for this struct was generated from the following file:
- /home/docs/checkouts/readthedocs.org/user_builds/impactx/checkouts/latest/src/particles/elements/ConstF.H