ImpactX
Loading...
Searching...
No Matches
HandleWakefield.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: Alex Bojanich, Chad Mitchell, Axel Huebl
8 * License: BSD-3-Clause-LBNL
9 */
10#ifndef HANDLE_WAKEFIELD_H
11#define HANDLE_WAKEFIELD_H
12
14#include "ChargeBinning.H"
15#include "CSRBendElement.H"
16#include "WakeConvolution.H"
17#include "WakePush.H"
18
19#include <cmath>
20#include <stdexcept>
21#include <vector>
22
23
25{
34 template <typename T_Element>
36 impactx::ImpactXParticleContainer& particle_container,
37 T_Element const& element_variant,
38 amrex::Real slice_ds,
39 bool print_wakefield = false
40 )
41 {
42 BL_PROFILE("impactx::particles::wakefields::HandleWakefield")
43
44 amrex::ParmParse pp_algo("algo");
45 bool csr = false;
46 pp_algo.queryAdd("csr", csr);
47
48 // Call the CSR bend function
49 auto const [element_has_csr, R] = impactx::particles::wakefields::CSRBendElement(element_variant, particle_container.GetRefParticle());
50
51 // Enter loop if lattice has bend element
52 if (csr && element_has_csr)
53 {
54#ifndef ImpactX_USE_FFT
55 throw std::runtime_error("algo.csr was requested but ImpactX was not compiled with FFT support. Recompile with ImpactX_FFT=ON.");
56#endif
57
58 int csr_bins = 150;
59 pp_algo.queryAddWithParser("csr_bins", csr_bins);
60
61 // Measure beam size, extract the min, max of particle positions
62 [[maybe_unused]] auto const [x_min, y_min, t_min, x_max, y_max, t_max] =
63 particle_container.MinAndMaxPositions();
64
65 // Ensure we have a beam with a finite extent in t
66 if (t_min == t_max)
67 {
68 auto const npart = particle_container.TotalNumberOfParticles();
69 std::string error_msg = "Coherent Synchrotron Radiation (CSR) "
70 "modeling was enabled, but ";
71 if (npart == 1) {
72 error_msg.append("there is only one particle in the beam. ");
73 } else {
74 error_msg.append("the beam is flat in t. ");
75 }
76 error_msg.append(
77 "CSR modeling is a collective effect and requires "
78 "a particle beam with a finite extent in t. "
79 "Please set the CSR option to false."
80 );
81 throw std::runtime_error(error_msg);
82 }
83
84 // Set parameters for charge deposition
85 bool const is_unity_particle_weight = false; // Only true if w = 1
86 bool const GetNumberDensity = true;
87
88 int const num_bins = csr_bins; // Set resolution
89 amrex::Real const bin_min = t_min;
90 amrex::Real const bin_max = t_max;
91 amrex::Real const bin_size = (bin_max - bin_min) / (num_bins - 1); // number of evaluation points
92
93 // Allocate memory for the charge profile
94 amrex::Gpu::DeviceVector<amrex::Real> charge_distribution(num_bins + 1, 0.0);
95 amrex::Gpu::DeviceVector<amrex::Real> mean_x(num_bins, 0.0);
96 amrex::Gpu::DeviceVector<amrex::Real> mean_y(num_bins, 0.0);
97
98 // Call charge deposition function
99 impactx::particles::wakefields::DepositCharge1D(particle_container, charge_distribution, bin_min, bin_size, is_unity_particle_weight);
100
101 // Sum up all partial charge histograms to one MPI process to calculate the global wakefield.
102 // Once calculated, we will distribute convolved_wakefield back to every MPI process.
104 charge_distribution.data(),
105 charge_distribution.size(),
108 );
109
110 amrex::Gpu::DeviceVector<amrex::Real> convolved_wakefield;
112 // Call the mean transverse position function
113 impactx::particles::wakefields::MeanTransversePosition(particle_container, mean_x, mean_y, bin_min,
114 bin_size, is_unity_particle_weight);
115
116 // Call charge density derivative function
117 amrex::Gpu::DeviceVector<amrex::Real> slopes(charge_distribution.size() - 1, 0.0);
118 impactx::particles::wakefields::DerivativeCharge1D(charge_distribution, slopes, bin_size,
119 GetNumberDensity); // Use number derivatives for convolution with CSR
120
121 // Construct CSR wake function on 2N support
122 amrex::Gpu::DeviceVector<amrex::Real> wake_function(num_bins*2, 0.0);
123 amrex::Real *const dptr_wake_function = wake_function.data();
124 auto const dR = R; // for NVCC capture
125 amrex::ParallelFor(num_bins*2, [=] AMREX_GPU_DEVICE(int i) {
126 if (i < num_bins) {
127 amrex::Real const s = static_cast<amrex::Real>(i) * bin_size;
128 dptr_wake_function[i] = impactx::particles::wakefields::w_l_csr(s, dR, bin_size);
129 }
130 else if (i > num_bins) {
131 amrex::Real const s = static_cast<amrex::Real>(i - 2*num_bins) * bin_size;
132 dptr_wake_function[i] = impactx::particles::wakefields::w_l_csr(s, dR, bin_size);
133 }
134 });
135
136 // Call convolution function
137 convolved_wakefield = impactx::particles::wakefields::convolve_fft(slopes, wake_function, bin_size);
138 }
139
140 // Broadcast the global wakefield to every MPI rank
142 convolved_wakefield.data(),
143 convolved_wakefield.size(),
145 );
146
147 // Check convolution output
148 if (print_wakefield && amrex::ParallelDescriptor::IOProcessor())
149 {
150 std::cout << "Convolved wakefield: ";
151 std::ofstream outfile("convolved_wakefield.txt");
152 for (double const i : convolved_wakefield) {
153 std::cout << i << " ";
154 outfile << i << std::endl;
155 }
156 std::cout << std::endl;
157 outfile.close();
158 }
159
160 // Call function to kick particles with wake
161 impactx::particles::wakefields::WakePush(particle_container, convolved_wakefield, slice_ds, bin_size, t_min);
162 }
163 }
164
165} // namespace impactx::particles::wakefields
166
167#endif // HANDLE_WAKEFIELD_H
#define BL_PROFILE(a)
#define AMREX_GPU_DEVICE
size_type size() const noexcept
T * data() noexcept
int queryAdd(std::string_view name, T &ref)
int queryAddWithParser(std::string_view name, T &ref)
Long TotalNumberOfParticles(bool only_valid=true, bool only_local=false) const
Definition ImpactXParticleContainer.H:136
std::tuple< amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal > MinAndMaxPositions()
Definition ImpactXParticleContainer.cpp:440
RefPart & GetRefParticle()
Definition ImpactXParticleContainer.cpp:419
amrex_real Real
PODVector< T, ArenaAllocator< T > > DeviceVector
int IOProcessorNumber() noexcept
void Sum(T &v, int root, MPI_Comm comm)
bool IOProcessor() noexcept
MPI_Comm Communicator() noexcept
void Bcast(void *, int, MPI_Datatype, int, MPI_Comm)
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition ChargeBinning.cpp:17
void DerivativeCharge1D(amrex::Gpu::DeviceVector< amrex::Real > const &charge_distribution, amrex::Gpu::DeviceVector< amrex::Real > &slopes, amrex::Real bin_size, bool GetNumberDensity)
Definition ChargeBinning.cpp:87
void HandleWakefield(impactx::ImpactXParticleContainer &particle_container, T_Element const &element_variant, amrex::Real slice_ds, bool print_wakefield=false)
Definition HandleWakefield.H:35
void DepositCharge1D(impactx::ImpactXParticleContainer &myspc, amrex::Gpu::DeviceVector< amrex::Real > &charge_distribution, amrex::Real bin_min, amrex::Real bin_size, bool is_unity_particle_weight)
Definition ChargeBinning.cpp:18
void MeanTransversePosition(impactx::ImpactXParticleContainer &myspc, amrex::Gpu::DeviceVector< amrex::Real > &mean_x, amrex::Gpu::DeviceVector< amrex::Real > &mean_y, amrex::Real bin_min, amrex::Real bin_size, bool is_unity_particle_weight)
Definition ChargeBinning.cpp:115
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real w_l_csr(amrex::Real s, amrex::Real R, amrex::Real const bin_size)
Definition WakeConvolution.H:106
std::tuple< bool, amrex::Real > CSRBendElement(VariantType const &element_variant, RefPart const &refpart)
Definition CSRBendElement.H:35
amrex::Gpu::DeviceVector< amrex::Real > convolve_fft(amrex::Gpu::DeviceVector< amrex::Real > const &beam_profile_slope, amrex::Gpu::DeviceVector< amrex::Real > const &wake_func, amrex::Real delta_t)
Definition WakeConvolution.cpp:66
void WakePush(ImpactXParticleContainer &pc, amrex::Gpu::DeviceVector< amrex::Real > const &convolved_wakefield, amrex::ParticleReal slice_ds, amrex::Real bin_size, amrex::Real bin_min)
Definition WakePush.cpp:21
@ s
fixed s as the independent variable
Definition ImpactXParticleContainer.H:37