ImpactX
Loading...
Searching...
No Matches
InitDistribution.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: Axel Huebl, Andrew Myers, Marco Garten
8 * License: BSD-3-Clause-LBNL
9 */
10#ifndef IMPACTX_INITIALIZATION_INITDISTRIBUTION_H
11#define IMPACTX_INITIALIZATION_INITDISTRIBUTION_H
12
16
17#include <AMReX_Extension.H> // for AMREX_RESTRICT
18#include <AMReX_REAL.H>
19
20#include <optional> // for std::optional
21#include <utility> // for std::move
22
23
25{
31 RefPart
32 read_reference_particle (amrex::ParmParse const & pp_dist);
33
40 read_distribution (amrex::ParmParse const & pp_dist);
41
47 Envelope
50 std::optional<amrex::ParticleReal> intensity = std::nullopt
51 );
52
66 template <typename T_Distribution>
68 {
80 T_Distribution distribution,
81 amrex::ParticleReal* AMREX_RESTRICT part_x,
82 amrex::ParticleReal* AMREX_RESTRICT part_y,
83 amrex::ParticleReal* AMREX_RESTRICT part_t,
84 amrex::ParticleReal* AMREX_RESTRICT part_px,
85 amrex::ParticleReal* AMREX_RESTRICT part_py,
86 amrex::ParticleReal* AMREX_RESTRICT part_pt
87 )
88 : m_distribution(std::move(distribution)),
89 m_part_x(part_x), m_part_y(part_y), m_part_t(part_t),
90 m_part_px(part_px), m_part_py(part_py), m_part_pt(part_pt)
91 {
92 }
93
98
105 void
107 amrex::Long i,
108 amrex::RandomEngine const & engine
109 ) const
110 {
112 m_part_x[i],
113 m_part_y[i],
114 m_part_t[i],
115 m_part_px[i],
116 m_part_py[i],
117 m_part_pt[i],
118 engine
119 );
120 }
121
122 private:
123 T_Distribution const m_distribution;
124 amrex::ParticleReal* const AMREX_RESTRICT m_part_x;
125 amrex::ParticleReal* const AMREX_RESTRICT m_part_y;
126 amrex::ParticleReal* const AMREX_RESTRICT m_part_t;
127 amrex::ParticleReal* const AMREX_RESTRICT m_part_px;
128 amrex::ParticleReal* const AMREX_RESTRICT m_part_py;
129 amrex::ParticleReal* const AMREX_RESTRICT m_part_pt;
130 };
131
132
160 void
162 amrex::ParmParse const & pp_dist,
163 amrex::ParticleReal& lambdax, amrex::ParticleReal& lambday, amrex::ParticleReal& lambdat,
164 amrex::ParticleReal& lambdapx,amrex::ParticleReal& lambdapy, amrex::ParticleReal& lambdapt,
165 amrex::ParticleReal& muxpx, amrex::ParticleReal& muypy, amrex::ParticleReal& mutpt,
166 amrex::ParticleReal& meanx, amrex::ParticleReal& meany, amrex::ParticleReal& meant,
167 amrex::ParticleReal& meanpx, amrex::ParticleReal& meanpy, amrex::ParticleReal& meanpt,
168 amrex::ParticleReal& dispx, amrex::ParticleReal& disppx, amrex::ParticleReal& dispy, amrex::ParticleReal& disppy
169 );
170
171
198 void
200 amrex::ParmParse const & pp_dist,
201 amrex::ParticleReal& lambdax, amrex::ParticleReal& lambday, amrex::ParticleReal& lambdat,
202 amrex::ParticleReal& lambdapx,amrex::ParticleReal& lambdapy, amrex::ParticleReal& lambdapt,
203 amrex::ParticleReal& muxpx, amrex::ParticleReal& muypy, amrex::ParticleReal& mutpt,
204 amrex::ParticleReal& meanx, amrex::ParticleReal& meany, amrex::ParticleReal& meant,
205 amrex::ParticleReal& meanpx, amrex::ParticleReal& meanpy, amrex::ParticleReal& meanpt,
206 amrex::ParticleReal& dispx, amrex::ParticleReal& disppx, amrex::ParticleReal& dispy, amrex::ParticleReal& disppy
207 );
208
209} // namespace impactx::initialization
210
211#endif // IMPACTX_INITIALIZATION_INITDISTRIBUTION_H
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_DEVICE
Definition All.H:27
std::variant< Empty, Gaussian, Kurth4D, Kurth6D, KVdist, Thermal, Triangle, Semigaussian, Waterbag > KnownDistributions
Definition All.H:28
Definition AmrCoreData.cpp:18
void set_distribution_parameters_from_phase_space_inputs(amrex::ParmParse const &pp_dist, amrex::ParticleReal &lambdax, amrex::ParticleReal &lambday, amrex::ParticleReal &lambdat, amrex::ParticleReal &lambdapx, amrex::ParticleReal &lambdapy, amrex::ParticleReal &lambdapt, amrex::ParticleReal &muxpx, amrex::ParticleReal &muypy, amrex::ParticleReal &mutpt, amrex::ParticleReal &meanx, amrex::ParticleReal &meany, amrex::ParticleReal &meant, amrex::ParticleReal &meanpx, amrex::ParticleReal &meanpy, amrex::ParticleReal &meanpt, amrex::ParticleReal &dispx, amrex::ParticleReal &disppx, amrex::ParticleReal &dispy, amrex::ParticleReal &disppy)
Definition InitDistribution.cpp:491
RefPart read_reference_particle(amrex::ParmParse const &pp_dist)
Definition InitDistribution.cpp:39
distribution::KnownDistributions read_distribution(amrex::ParmParse const &pp_dist)
Definition InitDistribution.cpp:79
Envelope create_envelope(distribution::KnownDistributions const &distr, std::optional< amrex::ParticleReal > intensity=std::nullopt)
Definition InitDistribution.cpp:222
void set_distribution_parameters_from_twiss_inputs(amrex::ParmParse const &pp_dist, amrex::ParticleReal &lambdax, amrex::ParticleReal &lambday, amrex::ParticleReal &lambdat, amrex::ParticleReal &lambdapx, amrex::ParticleReal &lambdapy, amrex::ParticleReal &lambdapt, amrex::ParticleReal &muxpx, amrex::ParticleReal &muypy, amrex::ParticleReal &mutpt, amrex::ParticleReal &meanx, amrex::ParticleReal &meany, amrex::ParticleReal &meant, amrex::ParticleReal &meanpx, amrex::ParticleReal &meanpy, amrex::ParticleReal &meanpt, amrex::ParticleReal &dispx, amrex::ParticleReal &disppx, amrex::ParticleReal &dispy, amrex::ParticleReal &disppy)
Definition InitDistribution.cpp:411
amrex::ParticleReal *const AMREX_RESTRICT m_part_px
Definition InitDistribution.H:127
amrex::ParticleReal *const AMREX_RESTRICT m_part_pt
Definition InitDistribution.H:129
amrex::ParticleReal *const AMREX_RESTRICT m_part_t
Definition InitDistribution.H:126
amrex::ParticleReal *const AMREX_RESTRICT m_part_py
Definition InitDistribution.H:128
amrex::ParticleReal *const AMREX_RESTRICT m_part_y
Definition InitDistribution.H:125
amrex::ParticleReal *const AMREX_RESTRICT m_part_x
Definition InitDistribution.H:124
InitSingleParticleData(InitSingleParticleData const &)=default
InitSingleParticleData(InitSingleParticleData &&)=default
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(amrex::Long i, amrex::RandomEngine const &engine) const
Definition InitDistribution.H:106
InitSingleParticleData(T_Distribution distribution, amrex::ParticleReal *AMREX_RESTRICT part_x, amrex::ParticleReal *AMREX_RESTRICT part_y, amrex::ParticleReal *AMREX_RESTRICT part_t, amrex::ParticleReal *AMREX_RESTRICT part_px, amrex::ParticleReal *AMREX_RESTRICT part_py, amrex::ParticleReal *AMREX_RESTRICT part_pt)
Definition InitDistribution.H:79
T_Distribution const m_distribution
Definition InitDistribution.H:123