ImpactX
Loading...
Searching...
No Matches
KVdist.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_DISTRIBUTION_KVDIST
11#define IMPACTX_DISTRIBUTION_KVDIST
12
14
15#include <ablastr/constant.H>
16
17#include <AMReX_Math.H>
18#include <AMReX_Random.H>
19#include <AMReX_REAL.H>
20
21#include <cmath>
22
23
25{
26 struct KVdist
27 {
42 amrex::ParticleReal lambdax,
43 amrex::ParticleReal lambday,
44 amrex::ParticleReal lambdat,
45 amrex::ParticleReal lambdapx,
46 amrex::ParticleReal lambdapy,
47 amrex::ParticleReal lambdapt,
48 amrex::ParticleReal muxpx=0.0,
49 amrex::ParticleReal muypy=0.0,
50 amrex::ParticleReal mutpt=0.0,
51 amrex::ParticleReal meanx=0.0,
52 amrex::ParticleReal meany=0.0,
53 amrex::ParticleReal meant=0.0,
54 amrex::ParticleReal meanpx=0.0,
55 amrex::ParticleReal meanpy=0.0,
56 amrex::ParticleReal meanpt=0.0,
57 amrex::ParticleReal dispx=0.0,
58 amrex::ParticleReal disppx=0.0,
59 amrex::ParticleReal dispy=0.0,
60 amrex::ParticleReal disppy=0.0
61 )
62 : m_lambdaX(lambdax), m_lambdaY(lambday), m_lambdaT(lambdat), m_lambdaPx(lambdapx), m_lambdaPy(lambdapy),
63 m_lambdaPt(lambdapt), m_muxpx(muxpx), m_muypy(muypy), m_mutpt(mutpt), m_meanx(meanx), m_meany(meany),
64 m_meant(meant), m_meanpx(meanpx), m_meanpy(meanpy), m_meanpt(meanpt), m_dispx(dispx), m_disppx(disppx),
65 m_dispy(dispy), m_disppy(disppy)
66 {
67 }
68
76 void initialize ([[maybe_unused]] amrex::ParticleReal bunch_charge, [[maybe_unused]] RefPart const & ref)
77 {
78 }
79
84 void
86 {
87 }
88
101 amrex::ParticleReal & AMREX_RESTRICT x,
102 amrex::ParticleReal & AMREX_RESTRICT y,
103 amrex::ParticleReal & AMREX_RESTRICT t,
104 amrex::ParticleReal & AMREX_RESTRICT px,
105 amrex::ParticleReal & AMREX_RESTRICT py,
106 amrex::ParticleReal & AMREX_RESTRICT pt,
107 amrex::RandomEngine const & engine
108 ) const
109 {
110 using namespace amrex::literals;
111 using amrex::Math::powi;
113
114 amrex::ParticleReal v,phi,r,beta,p,u1,u2,ln1;
115 amrex::ParticleReal root,a1,a2;
116
117 // Sample and transform to define (x,y):
118 v = amrex::Random(engine);
119 phi = amrex::Random(engine);
120 phi = 2_prt*pi*phi;
121 r = std::sqrt(v);
122 x = r * std::cos(phi);
123 y = r * std::sin(phi);
124
125 // Sample and transform to define (px,py):
126 beta = amrex::Random(engine);
127 beta = 2_prt*pi*beta;
128 p = std::sqrt(1_prt-powi<2>(r));
129 px = p * std::cos(beta);
130 py = p * std::sin(beta);
131
132 // Sample and transform to define (t,pt):
133 t = amrex::Random(engine);
134 t = 2.0_prt*(t-0.5_prt);
135 u1 = amrex::Random(engine);
136 u2 = amrex::Random(engine);
137 ln1 = std::sqrt(-2_prt*std::log(u1));
138 pt = ln1 * std::cos(2_prt*pi*u2);
139
140 // Scale to produce the identity covariance matrix:
141 amrex::ParticleReal const c = std::sqrt(3.0_prt);
142 x = 2_prt*x;
143 y = 2_prt*y;
144 t = c*t;
145 px = 2_prt*px;
146 py = 2_prt*py;
147 // pt = pt;
148
149 // Transform to produce the desired second moments/correlations:
150 root = std::sqrt(1.0_prt-m_muxpx*m_muxpx);
151 a1 = m_lambdaX * x / root;
152 a2 = m_lambdaPx * (-m_muxpx * x / root + px);
153 x = a1;
154 px = a2;
155 root = std::sqrt(1.0_prt-m_muypy*m_muypy);
156 a1 = m_lambdaY * y / root;
157 a2 = m_lambdaPy * (-m_muypy * y / root + py);
158 y = a1;
159 py = a2;
160 root = std::sqrt(1.0_prt-m_mutpt*m_mutpt);
161 a1 = m_lambdaT * t / root;
162 a2 = m_lambdaPt * (-m_mutpt * t / root + pt);
163 t = a1;
164 pt = a2;
165
166 // Apply dispersive correlations
167 x = x - m_dispx * pt;
168 px = px - m_disppx * pt;
169 y = y - m_dispy * pt;
170 py = py - m_disppy * pt;
171
172 // Apply centroid translation
173 x = x + m_meanx;
174 px = px + m_meanpx;
175 y = y + m_meany;
176 py = py + m_meanpy;
177 t = t + m_meant;
178 pt = pt + m_meanpt;
179 }
180
181 amrex::ParticleReal m_lambdaX, m_lambdaY, m_lambdaT;
182 amrex::ParticleReal m_lambdaPx, m_lambdaPy, m_lambdaPt;
183 amrex::ParticleReal m_muxpx, m_muypy, m_mutpt;
184 amrex::ParticleReal m_meanx, m_meany, m_meant;
185 amrex::ParticleReal m_meanpx, m_meanpy, m_meanpt;
186 amrex::ParticleReal m_dispx, m_disppx, m_dispy, m_disppy;
187 };
188
189} // namespace impactx::distribution
190
191#endif // IMPACTX_DISTRIBUTION_KVDIST
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
static constexpr amrex::Real pi
constexpr T powi(T x) noexcept
Real Random()
Definition All.H:27
@ t
fixed t as the independent variable
Definition ImpactXParticleContainer.H:38
Definition ReferenceParticle.H:31
void initialize(amrex::ParticleReal bunch_charge, RefPart const &ref)
Definition KVdist.H:76
amrex::ParticleReal m_disppx
Definition KVdist.H:186
amrex::ParticleReal m_meant
Definition KVdist.H:184
amrex::ParticleReal m_lambdaPt
Definition KVdist.H:182
amrex::ParticleReal m_dispx
momentum coordinates of centroid offset
Definition KVdist.H:186
amrex::ParticleReal m_lambdaY
Definition KVdist.H:181
amrex::ParticleReal m_lambdaT
Definition KVdist.H:181
amrex::ParticleReal m_mutpt
Definition KVdist.H:183
amrex::ParticleReal m_lambdaX
Definition KVdist.H:181
void finalize()
Definition KVdist.H:85
amrex::ParticleReal m_meanpt
Definition KVdist.H:185
amrex::ParticleReal m_disppy
Definition KVdist.H:186
amrex::ParticleReal m_meany
Definition KVdist.H:184
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, amrex::RandomEngine const &engine) const
Definition KVdist.H:100
amrex::ParticleReal m_muxpx
related momentum axis intercepts of the phase space ellipse
Definition KVdist.H:183
amrex::ParticleReal m_meanpx
spatial coordinates of centroid offset
Definition KVdist.H:185
KVdist(amrex::ParticleReal lambdax, amrex::ParticleReal lambday, amrex::ParticleReal lambdat, amrex::ParticleReal lambdapx, amrex::ParticleReal lambdapy, amrex::ParticleReal lambdapt, amrex::ParticleReal muxpx=0.0, amrex::ParticleReal muypy=0.0, amrex::ParticleReal mutpt=0.0, amrex::ParticleReal meanx=0.0, amrex::ParticleReal meany=0.0, amrex::ParticleReal meant=0.0, amrex::ParticleReal meanpx=0.0, amrex::ParticleReal meanpy=0.0, amrex::ParticleReal meanpt=0.0, amrex::ParticleReal dispx=0.0, amrex::ParticleReal disppx=0.0, amrex::ParticleReal dispy=0.0, amrex::ParticleReal disppy=0.0)
Definition KVdist.H:41
amrex::ParticleReal m_lambdaPx
related position axis intercepts (length) of the phase space ellipse
Definition KVdist.H:182
amrex::ParticleReal m_muypy
Definition KVdist.H:183
amrex::ParticleReal m_lambdaPy
Definition KVdist.H:182
amrex::ParticleReal m_dispy
Definition KVdist.H:186
amrex::ParticleReal m_meanpy
Definition KVdist.H:185
amrex::ParticleReal m_meanx
correlation length-momentum
Definition KVdist.H:184