10#ifndef HANDLE_WAKEFIELD_H
11#define HANDLE_WAKEFIELD_H
34 template <
typename T_Element>
37 T_Element
const& element_variant,
39 bool print_wakefield =
false
42 BL_PROFILE(
"impactx::particles::wakefields::HandleWakefield")
52 if (csr && element_has_csr)
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.");
62 [[maybe_unused]]
auto const [x_min, y_min, t_min, x_max, y_max, t_max] =
69 std::string error_msg =
"Coherent Synchrotron Radiation (CSR) "
70 "modeling was enabled, but ";
72 error_msg.append(
"there is only one particle in the beam. ");
74 error_msg.append(
"the beam is flat in t. ");
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."
81 throw std::runtime_error(error_msg);
85 bool const is_unity_particle_weight =
false;
86 bool const GetNumberDensity =
true;
88 int const num_bins = csr_bins;
91 amrex::Real const bin_size = (bin_max - bin_min) / (num_bins - 1);
104 charge_distribution.
data(),
105 charge_distribution.
size(),
114 bin_size, is_unity_particle_weight);
130 else if (i > num_bins) {
142 convolved_wakefield.
data(),
143 convolved_wakefield.
size(),
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;
156 std::cout << std::endl;
size_type size() const 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
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