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#ifndef ImpactX_USE_FFT
21#include <stdexcept>
22#endif
23#include <vector>
24
25
27{
36 template <typename T_Element>
38 impactx::ImpactXParticleContainer& particle_container,
39 T_Element const& element_variant,
40 amrex::Real slice_ds,
41 bool print_wakefield = false
42 )
43 {
44 BL_PROFILE("impactx::particles::wakefields::HandleWakefield")
45
46 amrex::ParmParse pp_algo("algo");
47 bool csr = false;
48 pp_algo.queryAdd("csr", csr);
49
50 // Call the CSR bend function
51 auto const [element_has_csr, R] = impactx::particles::wakefields::CSRBendElement(element_variant, particle_container.GetRefParticle());
52
53 // Enter loop if lattice has bend element
54 if (csr && element_has_csr)
55 {
56#ifndef ImpactX_USE_FFT
57 throw std::runtime_error("algo.csr was requested but ImpactX was not compiled with FFT support. Recompile with ImpactX_FFT=ON.");
58#endif
59
60 int csr_bins = 150;
61 pp_algo.queryAddWithParser("csr_bins", csr_bins);
62
63 // Measure beam size, extract the min, max of particle positions
64 [[maybe_unused]] auto const [x_min, y_min, t_min, x_max, y_max, t_max] =
65 particle_container.MinAndMaxPositions();
66
67 // Set parameters for charge deposition
68 bool const is_unity_particle_weight = false; // Only true if w = 1
69 bool const GetNumberDensity = true;
70
71 int const num_bins = csr_bins; // Set resolution
72 amrex::Real const bin_min = t_min;
73 amrex::Real const bin_max = t_max;
74 amrex::Real const bin_size = (bin_max - bin_min) / (num_bins - 1); // number of evaluation points
75
76 // Allocate memory for the charge profile
77 amrex::Gpu::DeviceVector<amrex::Real> charge_distribution(num_bins + 1, 0.0);
78 amrex::Gpu::DeviceVector<amrex::Real> mean_x(num_bins, 0.0);
79 amrex::Gpu::DeviceVector<amrex::Real> mean_y(num_bins, 0.0);
80
81 // Call charge deposition function
82 impactx::particles::wakefields::DepositCharge1D(particle_container, charge_distribution, bin_min, bin_size, is_unity_particle_weight);
83
84 // Sum up all partial charge histograms to one MPI process to calculate the global wakefield.
85 // Once calculated, we will distribute convolved_wakefield back to every MPI process.
87 charge_distribution.data(),
88 charge_distribution.size(),
91 );
92
93 amrex::Gpu::DeviceVector<amrex::Real> convolved_wakefield;
95 // Call the mean transverse position function
96 impactx::particles::wakefields::MeanTransversePosition(particle_container, mean_x, mean_y, bin_min,
97 bin_size, is_unity_particle_weight);
98
99 // Call charge density derivative function
100 amrex::Gpu::DeviceVector<amrex::Real> slopes(charge_distribution.size() - 1, 0.0);
101 impactx::particles::wakefields::DerivativeCharge1D(charge_distribution, slopes, bin_size,
102 GetNumberDensity); // Use number derivatives for convolution with CSR
103
104 // Construct CSR wake function on 2N support
105 amrex::Gpu::DeviceVector<amrex::Real> wake_function(num_bins*2, 0.0);
106 amrex::Real *const dptr_wake_function = wake_function.data();
107 auto const dR = R; // for NVCC capture
108 amrex::ParallelFor(num_bins*2, [=] AMREX_GPU_DEVICE(int i) {
109 if (i < num_bins) {
110 amrex::Real const s = static_cast<amrex::Real>(i) * bin_size;
111 dptr_wake_function[i] = impactx::particles::wakefields::w_l_csr(s, dR, bin_size);
112 }
113 else if (i > num_bins) {
114 amrex::Real const s = static_cast<amrex::Real>(i - 2*num_bins) * bin_size;
115 dptr_wake_function[i] = impactx::particles::wakefields::w_l_csr(s, dR, bin_size);
116 }
117 });
118
119 // Call convolution function
120 convolved_wakefield = impactx::particles::wakefields::convolve_fft(slopes, wake_function, bin_size);
121 }
122
123 // Broadcast the global wakefield to every MPI rank
125 convolved_wakefield.data(),
126 convolved_wakefield.size(),
128 );
129
130 // Check convolution output
131 if (print_wakefield && amrex::ParallelDescriptor::IOProcessor())
132 {
133 std::cout << "Convolved wakefield: ";
134 std::ofstream outfile("convolved_wakefield.txt");
135 for (double const i : convolved_wakefield) {
136 std::cout << i << " ";
137 outfile << i << std::endl;
138 }
139 std::cout << std::endl;
140 outfile.close();
141 }
142
143 // Call function to kick particles with wake
144 impactx::particles::wakefields::WakePush(particle_container, convolved_wakefield, slice_ds, bin_size, t_min);
145 }
146 }
147
148} // namespace impactx::particles::wakefields
149
150#endif // HANDLE_WAKEFIELD_H
#define BL_PROFILE(a)
#define AMREX_GPU_DEVICE
size_type size() const noexcept
T * data() noexcept
int queryAddWithParser(const char *name, T &ref)
int queryAdd(const char *name, T &ref)
Definition ImpactXParticleContainer.H:133
std::tuple< amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal, amrex::ParticleReal > MinAndMaxPositions()
Definition ImpactXParticleContainer.cpp:399
RefPart & GetRefParticle()
Definition ImpactXParticleContainer.cpp:378
PODVector< T, ArenaAllocator< T > > DeviceVector
MPI_Comm Communicator() noexcept
void Bcast(void *, int, MPI_Datatype, int, MPI_Comm)
int IOProcessorNumber() noexcept
bool IOProcessor() noexcept
void Sum(T &v, int root, MPI_Comm 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:37
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