#include <Triangle.H>
|
| | Triangle (amrex::ParticleReal const lambdax, amrex::ParticleReal const lambday, amrex::ParticleReal const lambdat, amrex::ParticleReal const lambdapx, amrex::ParticleReal const lambdapy, amrex::ParticleReal const lambdapt, amrex::ParticleReal const muxpx=0.0, amrex::ParticleReal const muypy=0.0, amrex::ParticleReal const mutpt=0.0, amrex::ParticleReal const meanx=0.0, amrex::ParticleReal const meany=0.0, amrex::ParticleReal const meant=0.0, amrex::ParticleReal const meanpx=0.0, amrex::ParticleReal const meanpy=0.0, amrex::ParticleReal const meanpt=0.0, amrex::ParticleReal const dispx=0.0, amrex::ParticleReal const disppx=0.0, amrex::ParticleReal const dispy=0.0, amrex::ParticleReal const disppy=0.0) |
| |
| void | initialize (amrex::ParticleReal bunch_charge, RefPart const &ref) |
| |
| void | finalize () |
| |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void | operator() (amrex::ParticleReal &x, amrex::ParticleReal &y, amrex::ParticleReal &t, amrex::ParticleReal &px, amrex::ParticleReal &py, amrex::ParticleReal &pt, amrex::RandomEngine const &engine) const |
| |
|
| amrex::ParticleReal | m_lambdaX |
| |
| amrex::ParticleReal | m_lambdaY |
| |
| amrex::ParticleReal | m_lambdaT |
| |
| amrex::ParticleReal | m_lambdaPx |
| | related position axis intercepts (length) of the phase space ellipse
|
| |
| amrex::ParticleReal | m_lambdaPy |
| |
| amrex::ParticleReal | m_lambdaPt |
| |
| amrex::ParticleReal | m_muxpx |
| | related momentum axis intercepts of the phase space ellipse
|
| |
| amrex::ParticleReal | m_muypy |
| |
| amrex::ParticleReal | m_mutpt |
| |
| amrex::ParticleReal | m_meanx |
| | correlation length-momentum
|
| |
| amrex::ParticleReal | m_meany |
| |
| amrex::ParticleReal | m_meant |
| |
| amrex::ParticleReal | m_meanpx |
| | spatial coordinates of centroid offset
|
| |
| amrex::ParticleReal | m_meanpy |
| |
| amrex::ParticleReal | m_meanpt |
| |
| amrex::ParticleReal | m_dispx |
| | momentum coordinates of centroid offset
|
| |
| amrex::ParticleReal | m_disppx |
| |
| amrex::ParticleReal | m_dispy |
| |
| amrex::ParticleReal | m_disppy |
| |
◆ Triangle()
| impactx::distribution::Triangle::Triangle |
( |
amrex::ParticleReal const | lambdax, |
|
|
amrex::ParticleReal const | lambday, |
|
|
amrex::ParticleReal const | lambdat, |
|
|
amrex::ParticleReal const | lambdapx, |
|
|
amrex::ParticleReal const | lambdapy, |
|
|
amrex::ParticleReal const | lambdapt, |
|
|
amrex::ParticleReal const | muxpx = 0.0, |
|
|
amrex::ParticleReal const | muypy = 0.0, |
|
|
amrex::ParticleReal const | mutpt = 0.0, |
|
|
amrex::ParticleReal const | meanx = 0.0, |
|
|
amrex::ParticleReal const | meany = 0.0, |
|
|
amrex::ParticleReal const | meant = 0.0, |
|
|
amrex::ParticleReal const | meanpx = 0.0, |
|
|
amrex::ParticleReal const | meanpy = 0.0, |
|
|
amrex::ParticleReal const | meanpt = 0.0, |
|
|
amrex::ParticleReal const | dispx = 0.0, |
|
|
amrex::ParticleReal const | disppx = 0.0, |
|
|
amrex::ParticleReal const | dispy = 0.0, |
|
|
amrex::ParticleReal const | disppy = 0.0 ) |
|
inline |
A Triangle distribution for LPA applications.
Return sampling from a ramped, triangular current profile with a Gaussian energy spread (possibly correlated). The transverse distribution is a 4D waterbag.
- Parameters
-
| lambdax,lambday,lambdat | for zero correlation, these are the related RMS sizes (in meters) |
| lambdapx,lambdapy,lambdapt | RMS momentum |
| muxpx,muypy,mutpt | correlation length-momentum |
| meanx,meany,meant | offsets of the mean (centroid) positions from those of the reference particle |
| meanpx,meanpy,meanpt | offsets of the mean (centroid) momenta from those of the reference particle |
| dispx,disppx,dispy,disppy | dispersion and its derivative in horizontal and vertical directions |
◆ finalize()
| void impactx::distribution::Triangle::finalize |
( |
| ) |
|
|
inline |
Close and deallocate all data and handles.
Nothing to do here.
◆ initialize()
| void impactx::distribution::Triangle::initialize |
( |
amrex::ParticleReal | bunch_charge, |
|
|
RefPart const & | ref ) |
|
inline |
Initialize the distribution.
Nothing to do here.
- Parameters
-
| bunch_charge | charge of the beam in C |
| ref | the reference particle |
◆ operator()()
Return 1 6D particle coordinate
- 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 |
| engine | a random number engine (with associated state) |
◆ m_disppx
| amrex::ParticleReal impactx::distribution::Triangle::m_disppx |
◆ m_disppy
| amrex::ParticleReal impactx::distribution::Triangle::m_disppy |
◆ m_dispx
| amrex::ParticleReal impactx::distribution::Triangle::m_dispx |
momentum coordinates of centroid offset
◆ m_dispy
| amrex::ParticleReal impactx::distribution::Triangle::m_dispy |
◆ m_lambdaPt
| amrex::ParticleReal impactx::distribution::Triangle::m_lambdaPt |
◆ m_lambdaPx
| amrex::ParticleReal impactx::distribution::Triangle::m_lambdaPx |
related position axis intercepts (length) of the phase space ellipse
◆ m_lambdaPy
| amrex::ParticleReal impactx::distribution::Triangle::m_lambdaPy |
◆ m_lambdaT
| amrex::ParticleReal impactx::distribution::Triangle::m_lambdaT |
◆ m_lambdaX
| amrex::ParticleReal impactx::distribution::Triangle::m_lambdaX |
◆ m_lambdaY
| amrex::ParticleReal impactx::distribution::Triangle::m_lambdaY |
◆ m_meanpt
| amrex::ParticleReal impactx::distribution::Triangle::m_meanpt |
◆ m_meanpx
| amrex::ParticleReal impactx::distribution::Triangle::m_meanpx |
spatial coordinates of centroid offset
◆ m_meanpy
| amrex::ParticleReal impactx::distribution::Triangle::m_meanpy |
◆ m_meant
| amrex::ParticleReal impactx::distribution::Triangle::m_meant |
◆ m_meanx
| amrex::ParticleReal impactx::distribution::Triangle::m_meanx |
correlation length-momentum
◆ m_meany
| amrex::ParticleReal impactx::distribution::Triangle::m_meany |
◆ m_mutpt
| amrex::ParticleReal impactx::distribution::Triangle::m_mutpt |
◆ m_muxpx
| amrex::ParticleReal impactx::distribution::Triangle::m_muxpx |
related momentum axis intercepts of the phase space ellipse
◆ m_muypy
| amrex::ParticleReal impactx::distribution::Triangle::m_muypy |
The documentation for this struct was generated from the following file:
- /home/docs/checkouts/readthedocs.org/user_builds/impactx/checkouts/1182/src/particles/distribution/Triangle.H