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