32#define ESPRESSO_P3M_GPU_FLOAT
35#ifdef ESPRESSO_P3M_GPU_FLOAT
36#define REAL_TYPE float
37#define FFT_TYPE_COMPLEX cufftComplex
38#define FFT_FORW_FFT cufftExecR2C
39#define FFT_BACK_FFT cufftExecC2R
40#define FFT_PLAN_FORW_FLAG CUFFT_R2C
41#define FFT_PLAN_BACK_FLAG CUFFT_C2R
44#ifdef ESPRESSO_P3M_GPU_REAL_DOUBLE
45#define REAL_TYPE double
46#define FFT_TYPE_COMPLEX cufftDoubleComplex
47#define FFT_FORW_FFT cufftExecD2Z
48#define FFT_BACK_FFT cufftExecZ2D
49#define FFT_PLAN_FORW_FLAG CUFFT_D2Z
50#define FFT_PLAN_BACK_FLAG CUFFT_Z2D
57#include "system/System.hpp"
75#if defined(OMPI_MPI_H) || defined(_MPI_H)
76#error CU-file includes mpi.h! This should not happen!
126 auto const free_device_pointer = [](
auto *&ptr) {
127 if (ptr !=
nullptr) {
144 auto const cao3 = Utils::int_pow<3>(cao);
145 auto parts_per_block = 1u;
146 while ((parts_per_block + 1u) * cao3 <= 1024u) {
149 assert((n_part / parts_per_block) <= std::numeric_limits<unsigned>::max());
150 auto n =
static_cast<unsigned int>(n_part / parts_per_block);
151 auto n_blocks = ((n_part % parts_per_block) == 0u) ? std::max(1u, n) : n + 1u;
152 assert(n_blocks <= std::numeric_limits<unsigned>::max());
153 return std::make_pair(parts_per_block,
static_cast<unsigned>(n_blocks));
157 dim3 grid(n_blocks, 1u, 1u);
158 while (grid.x > 65536u) {
160 if ((n_blocks % grid.y) == 0u)
161 grid.x = std::max(1u, n_blocks / grid.y);
163 grid.x = n_blocks / grid.y + 1u;
180 for (
int i = 0; i < 3; ++i) {
181 Leni[i] = 1.0f / p.
box[i];
185 Zaehler[0] = Zaehler[1] = Zaehler[2] = *Nenner = 0.0;
190 S1 = int_pow<2 * cao>(
math::sinc(Meshi[0] * NMX));
193 ((NY > p.
mesh[1] / 2) ? NY - p.
mesh[1] : NY) + p.
mesh[1] * MY);
194 S2 = S1 * int_pow<2 * cao>(
math::sinc(Meshi[1] * NMY));
197 ((NZ > p.
mesh[2] / 2) ? NZ - p.
mesh[2] : NZ) + p.
mesh[2] * MZ);
198 S3 = S2 * int_pow<2 * cao>(
math::sinc(Meshi[2] * NMZ));
200 NM2 = sqr(NMX * Leni[0]) + sqr(NMY * Leni[1]) + sqr(NMZ * Leni[2]);
203 TE = exp(-sqr(std::numbers::pi_v<REAL_TYPE> / (p.
alpha)) * NM2);
205 Zaehler[0] += NMX * zwi * Leni[0];
206 Zaehler[1] += NMY * zwi * Leni[1];
207 Zaehler[2] += NMZ * zwi * Leni[2];
217 const auto NX =
static_cast<int>(blockDim.x * blockIdx.x + threadIdx.x);
218 const auto NY =
static_cast<int>(blockDim.y * blockIdx.y + threadIdx.y);
219 const auto NZ =
static_cast<int>(blockDim.z * blockIdx.z + threadIdx.z);
221 REAL_TYPE Zaehler[3] = {0.0, 0.0, 0.0}, Nenner = 0.0;
225 for (
int i = 0; i < 3; ++i) {
229 if ((NX >= p.
mesh[0]) || (NY >= p.
mesh[1]) || (NZ >= (p.
mesh[2] / 2 + 1)))
232 index = NX * p.
mesh[1] * (p.
mesh[2] / 2 + 1) + NY * (p.
mesh[2] / 2 + 1) + NZ;
234 if (((NX == 0) && (NY == 0) && (NZ == 0)) ||
235 ((NX % (p.
mesh[0] / 2) == 0) && (NY % (p.
mesh[1] / 2) == 0) &&
236 (NZ % (p.
mesh[2] / 2) == 0))) {
239 Aliasing_sums_ik<cao>(p, NX, NY, NZ, Zaehler, &Nenner);
245 zwi = Dnx * Zaehler[0] * Leni[0] + Dny * Zaehler[1] * Leni[1] +
246 Dnz * Zaehler[2] * Leni[2];
247 zwi /= ((sqr(Dnx * Leni[0]) + sqr(Dny * Leni[1]) + sqr(Dnz * Leni[2])) *
249 p.
G_hat[index] =
REAL_TYPE{2} * zwi / std::numbers::pi_v<REAL_TYPE>;
262 return static_cast<unsigned int>(p.
mesh[1] * (p.
mesh[2] / 2 + 1) * i +
263 (p.
mesh[2] / 2 + 1) * j + k);
268 auto const linear_index =
linear_index_k(p,
static_cast<int>(blockIdx.x),
269 static_cast<int>(blockIdx.y),
270 static_cast<int>(threadIdx.x));
272 auto const bidx =
static_cast<int>(blockIdx.x);
273 auto const bidy =
static_cast<int>(blockIdx.y);
274 auto const nx = (bidx > p.
mesh[0] / 2) ? bidx - p.
mesh[0] : bidx;
275 auto const ny = (bidy > p.
mesh[1] / 2) ? bidy - p.
mesh[1] : bidy;
276 auto const nz =
static_cast<int>(threadIdx.x);
280 buf.x =
REAL_TYPE(-2) * std::numbers::pi_v<REAL_TYPE> * meshw.y;
281 buf.y =
REAL_TYPE(+2) * std::numbers::pi_v<REAL_TYPE> * meshw.x;
284 static_cast<decltype(FFT_TYPE_COMPLEX::x)
>(nx) * buf.x / p.
box[0];
286 static_cast<decltype(FFT_TYPE_COMPLEX::x)
>(nx) * buf.y / p.
box[0];
289 static_cast<decltype(FFT_TYPE_COMPLEX::x)
>(ny) * buf.x / p.
box[1];
291 static_cast<decltype(FFT_TYPE_COMPLEX::x)
>(ny) * buf.y / p.
box[1];
294 static_cast<decltype(FFT_TYPE_COMPLEX::x)
>(nz) * buf.x / p.
box[2];
296 static_cast<decltype(FFT_TYPE_COMPLEX::x)
>(nz) * buf.y / p.
box[2];
299__device__
inline int wrap_index(
const int ind,
const int mesh) {
308 auto const linear_index =
linear_index_k(p,
static_cast<int>(blockIdx.x),
309 static_cast<int>(blockIdx.y),
310 static_cast<int>(threadIdx.x));
316template <
int cao,
bool shared>
318 float const *
const __restrict__ part_pos,
319 float const *
const __restrict__ part_q,
320 unsigned int const parts_per_block) {
321 auto const part_in_block = threadIdx.x /
static_cast<unsigned int>(cao);
322 auto const cao_id_x =
323 threadIdx.x - part_in_block *
static_cast<unsigned int>(cao);
326 parts_per_block * (blockIdx.x * gridDim.y + blockIdx.y) + part_in_block;
332 int nmp_x, nmp_y, nmp_z;
336 m_pos[0] = part_pos[3 *
id + 0] * params.
hi[0] - params.
pos_shift;
337 m_pos[1] = part_pos[3 *
id + 1] * params.
hi[1] - params.
pos_shift;
338 m_pos[2] = part_pos[3 *
id + 2] * params.
hi[2] - params.
pos_shift;
340 nmp_x =
static_cast<int>(floorf(m_pos[0] + 0.5f));
341 nmp_y =
static_cast<int>(floorf(m_pos[1] + 0.5f));
342 nmp_z =
static_cast<int>(floorf(m_pos[2] + 0.5f));
344 m_pos[0] -=
static_cast<REAL_TYPE>(nmp_x);
345 m_pos[1] -=
static_cast<REAL_TYPE>(nmp_y);
346 m_pos[2] -=
static_cast<REAL_TYPE>(nmp_z);
348 nmp_x =
wrap_index(nmp_x +
static_cast<int>(cao_id_x), params.
mesh[0]);
349 nmp_y =
wrap_index(nmp_y +
static_cast<int>(threadIdx.y), params.
mesh[1]);
350 nmp_z =
wrap_index(nmp_z +
static_cast<int>(threadIdx.z), params.
mesh[2]);
354 extern __shared__
float weights[];
357 auto const offset =
static_cast<unsigned int>(cao) * part_in_block;
358 if ((threadIdx.y < 3u) && (threadIdx.z == 0u)) {
359 weights[3u * offset + 3u * cao_id_x + threadIdx.y] =
360 Utils::bspline<cao>(
static_cast<int>(cao_id_x), m_pos[threadIdx.y]);
365 auto const c = weights[3u * offset + 3u * cao_id_x] *
366 weights[3u * offset + 3u * threadIdx.y + 1u] *
367 weights[3u * offset + 3u * threadIdx.z + 2u] * part_q[id];
368 atomicAdd(&(charge_mesh[index]),
REAL_TYPE(c));
372 Utils::bspline<cao>(
static_cast<int>(cao_id_x), m_pos[0]) * part_q[id] *
373 Utils::bspline<cao>(
static_cast<int>(threadIdx.y), m_pos[1]) *
374 Utils::bspline<cao>(
static_cast<int>(threadIdx.z), m_pos[2]);
375 atomicAdd(&(charge_mesh[index]),
REAL_TYPE(c));
380 float const *
const __restrict__ part_pos,
381 float const *
const __restrict__ part_q) {
382 auto const cao =
static_cast<unsigned int>(params.
cao);
384 dim3
block(parts_per_block * cao, cao, cao);
387 auto const data_length = std::size_t(3u) *
388 static_cast<std::size_t
>(parts_per_block) *
389 static_cast<std::size_t
>(cao) *
sizeof(
REAL_TYPE);
390 switch (params.
cao) {
392 (assign_charge_kernel<1, false>)<<<grid,
block, std::size_t(0u),
nullptr>>>(
393 params, part_pos, part_q, parts_per_block);
396 (assign_charge_kernel<2, false>)<<<grid,
block, std::size_t(0u),
nullptr>>>(
397 params, part_pos, part_q, parts_per_block);
400 (assign_charge_kernel<3, true>)<<<grid,
block, data_length,
nullptr>>>(
401 params, part_pos, part_q, parts_per_block);
404 (assign_charge_kernel<4, true>)<<<grid,
block, data_length,
nullptr>>>(
405 params, part_pos, part_q, parts_per_block);
408 (assign_charge_kernel<5, true>)<<<grid,
block, data_length,
nullptr>>>(
409 params, part_pos, part_q, parts_per_block);
412 (assign_charge_kernel<6, true>)<<<grid,
block, data_length,
nullptr>>>(
413 params, part_pos, part_q, parts_per_block);
416 (assign_charge_kernel<7, true>)<<<grid,
block, data_length,
nullptr>>>(
417 params, part_pos, part_q, parts_per_block);
425template <
int cao,
bool shared>
427 float const *
const __restrict__ part_pos,
428 float const *
const __restrict__ part_q,
429 float *
const __restrict__ part_f,
431 unsigned int const parts_per_block) {
432 auto const part_in_block = threadIdx.x /
static_cast<unsigned int>(cao);
433 auto const cao_id_x =
434 threadIdx.x - part_in_block *
static_cast<unsigned int>(cao);
437 parts_per_block * (blockIdx.x * gridDim.y + blockIdx.y) + part_in_block;
438 if (
id >=
static_cast<unsigned>(params.
n_part))
443 int nmp_x, nmp_y, nmp_z;
445 m_pos[0] = part_pos[3u *
id + 0u] * params.
hi[0] - params.
pos_shift;
446 m_pos[1] = part_pos[3u *
id + 1u] * params.
hi[1] - params.
pos_shift;
447 m_pos[2] = part_pos[3u *
id + 2u] * params.
hi[2] - params.
pos_shift;
449 nmp_x =
static_cast<int>(floorf(m_pos[0] +
REAL_TYPE{0.5}));
450 nmp_y =
static_cast<int>(floorf(m_pos[1] +
REAL_TYPE{0.5}));
451 nmp_z =
static_cast<int>(floorf(m_pos[2] +
REAL_TYPE{0.5}));
453 m_pos[0] -=
static_cast<REAL_TYPE>(nmp_x);
454 m_pos[1] -=
static_cast<REAL_TYPE>(nmp_y);
455 m_pos[2] -=
static_cast<REAL_TYPE>(nmp_z);
457 nmp_x =
wrap_index(nmp_x +
static_cast<int>(cao_id_x), params.
mesh[0]);
458 nmp_y =
wrap_index(nmp_y +
static_cast<int>(threadIdx.y), params.
mesh[1]);
459 nmp_z =
wrap_index(nmp_z +
static_cast<int>(threadIdx.z), params.
mesh[2]);
463 extern __shared__
float weights[];
467 auto const offset =
static_cast<unsigned int>(cao) * part_in_block;
468 if ((threadIdx.y < 3u) && (threadIdx.z == 0u)) {
469 weights[3u * offset + 3u * cao_id_x + threadIdx.y] =
470 Utils::bspline<cao>(
static_cast<int>(cao_id_x), m_pos[threadIdx.y]);
475 c *=
REAL_TYPE(weights[3u * offset + 3u * cao_id_x] *
476 weights[3u * offset + 3u * threadIdx.y + 1u] *
477 weights[3u * offset + 3u * threadIdx.z + 2u]);
480 REAL_TYPE(Utils::bspline<cao>(
static_cast<int>(cao_id_x), m_pos[0]) *
481 Utils::bspline<cao>(
static_cast<int>(threadIdx.y), m_pos[1]) *
482 Utils::bspline<cao>(
static_cast<int>(threadIdx.z), m_pos[2]));
489 atomicAdd(&(part_f[3u *
id + 0u]),
float(c * force_mesh_x[index]));
490 atomicAdd(&(part_f[3u *
id + 1u]),
float(c * force_mesh_y[index]));
491 atomicAdd(&(part_f[3u *
id + 2u]),
float(c * force_mesh_z[index]));
495 float const *
const __restrict__ part_pos,
496 float const *
const __restrict__ part_q,
497 float *
const __restrict__ part_f,
499 auto const cao =
static_cast<unsigned int>(params.
cao);
501 dim3
block(parts_per_block * cao, cao, cao);
506 auto const data_length = std::size_t(3u) *
507 static_cast<std::size_t
>(parts_per_block) *
508 static_cast<std::size_t
>(cao) *
sizeof(float);
509 switch (params.
cao) {
511 (assign_forces_kernel<1, false>)<<<grid,
block, std::size_t(0u),
nullptr>>>(
512 params, part_pos, part_q, part_f, prefactor, parts_per_block);
515 (assign_forces_kernel<2, false>)<<<grid,
block, std::size_t(0u),
nullptr>>>(
516 params, part_pos, part_q, part_f, prefactor, parts_per_block);
519 (assign_forces_kernel<3, true>)<<<grid,
block, data_length,
nullptr>>>(
520 params, part_pos, part_q, part_f, prefactor, parts_per_block);
523 (assign_forces_kernel<4, true>)<<<grid,
block, data_length,
nullptr>>>(
524 params, part_pos, part_q, part_f, prefactor, parts_per_block);
527 (assign_forces_kernel<5, true>)<<<grid,
block, data_length,
nullptr>>>(
528 params, part_pos, part_q, part_f, prefactor, parts_per_block);
531 (assign_forces_kernel<6, true>)<<<grid,
block, data_length,
nullptr>>>(
532 params, part_pos, part_q, part_f, prefactor, parts_per_block);
535 (assign_forces_kernel<7, true>)<<<grid,
block, data_length,
nullptr>>>(
536 params, part_pos, part_q, part_f, prefactor, parts_per_block);
557 data = std::make_shared<P3MGpuParams>();
560 auto &p3m_gpu_data = data->p3m_gpu_data;
561 bool do_reinit =
false, mesh_changed =
false;
562 assert(n_part <= std::numeric_limits<unsigned int>::max());
563 p3m_gpu_data.n_part =
static_cast<unsigned>(n_part);
565 if (not data->is_initialized or p3m_gpu_data.alpha != alpha) {
566 p3m_gpu_data.alpha =
static_cast<REAL_TYPE>(alpha);
570 if (not data->is_initialized or p3m_gpu_data.cao != cao) {
571 p3m_gpu_data.cao = cao;
573 p3m_gpu_data.pos_shift =
static_cast<REAL_TYPE>((p3m_gpu_data.cao - 1) / 2);
577 if (not data->is_initialized or mesh !=
Utils::Vector3i(p3m_gpu_data.mesh)) {
578 std::ranges::copy(mesh, p3m_gpu_data.mesh);
583 if (
auto constexpr eps =
584 static_cast<double>(std::numeric_limits<float>::epsilon());
585 not data->is_initialized or
587 std::ranges::copy(box_l, p3m_gpu_data.box);
591 p3m_gpu_data.mesh_z_padded = (mesh[2] / 2 + 1) * 2;
592 p3m_gpu_data.mesh_size = mesh[0] * mesh[1] * p3m_gpu_data.mesh_z_padded;
594 for (
auto i = 0u; i < 3u; ++i) {
596 static_cast<REAL_TYPE>(p3m_gpu_data.mesh[i]) / p3m_gpu_data.box[i];
599 if (data->is_initialized and mesh_changed) {
600 data->free_device_memory();
601 data->is_initialized =
false;
604 if (not data->is_initialized and p3m_gpu_data.mesh_size > 0) {
606 auto const cmesh_size =
607 static_cast<std::size_t
>(p3m_gpu_data.mesh[0]) *
608 static_cast<std::size_t
>(p3m_gpu_data.mesh[1]) *
609 static_cast<std::size_t
>(p3m_gpu_data.mesh[2] / 2 + 1);
611 cuda_safe_mem(cudaMalloc((
void **)&(p3m_gpu_data.charge_mesh), mesh_len));
612 cuda_safe_mem(cudaMalloc((
void **)&(p3m_gpu_data.force_mesh_x), mesh_len));
613 cuda_safe_mem(cudaMalloc((
void **)&(p3m_gpu_data.force_mesh_y), mesh_len));
614 cuda_safe_mem(cudaMalloc((
void **)&(p3m_gpu_data.force_mesh_z), mesh_len));
626 if (cufftPlan3d(&(data->p3m_fft.forw_plan), mesh[0], mesh[1], mesh[2],
628 cufftPlan3d(&(data->p3m_fft.back_plan), mesh[0], mesh[1], mesh[2],
630 throw std::runtime_error(
"Unable to create fft plan");
635 if ((do_reinit or not data->is_initialized) and p3m_gpu_data.mesh_size > 0) {
638 block.x =
static_cast<unsigned>(512 / mesh[0] + 1);
639 block.y =
static_cast<unsigned>(mesh[1]);
641 grid.x =
static_cast<unsigned>(mesh[0]) /
block.x + 1;
642 grid.z =
static_cast<unsigned>(mesh[2]) / 2 + 1;
644 switch (p3m_gpu_data.cao) {
675 if (p3m_gpu_data.mesh_size > 0)
676 data->is_initialized =
true;
683 double prefactor, std::size_t n_part) {
685 p3m_gpu_data.
n_part =
static_cast<unsigned>(n_part);
694 dim3 gridConv(
static_cast<unsigned>(p3m_gpu_data.mesh[0]),
695 static_cast<unsigned>(p3m_gpu_data.mesh[1]), 1u);
696 dim3 threadsConv(
static_cast<unsigned>(p3m_gpu_data.mesh[2] / 2 + 1), 1u, 1u);
703 static_cast<std::size_t
>(p3m_gpu_data.mesh_size) *
712 p3m_gpu_data.charge_mesh) != CUFFT_SUCCESS) {
713 throw std::runtime_error(
"CUFFT error: Forward FFT failed");
731 assign_forces(p3m_gpu_data, positions_device, charges_device, forces_device,
Particle data communication manager for the GPU.
float * get_particle_charges_device() const
float * get_particle_forces_device() const
float * get_particle_positions_device() const
static DEVICE_QUALIFIER constexpr Vector< T, N > broadcast(typename Base::value_type const &value) noexcept
Create a vector that has all entries set to the same value.
static std::shared_ptr< scoped_pause > make_shared_pause_scoped()
Generate a shared handle to temporarily disable any currently active exception trap for the lifetime ...
void cuda_check_errors_exit(const dim3 &block, const dim3 &grid, const char *function, const char *file, unsigned int line)
In case of error during a CUDA operation, print the error message and exit.
static double * block(double *p, std::size_t index, std::size_t size)
#define P3M_BRILLOUIN
P3M: Number of Brillouin zones taken into account in the calculation of the optimal influence functio...
T product(Vector< T, N > const &v)
DEVICE_QUALIFIER constexpr T sqr(T x)
Calculates the SQuaRe of x.
DEVICE_QUALIFIER constexpr T int_pow(T x)
Calculate integer powers.
__device__ auto linear_index_k(P3MGpuData const &p, int i, int j, int k)
__device__ auto linear_index_r(P3MGpuData const &p, int i, int j, int k)
DEVICE_QUALIFIER auto sinc(T x)
Calculate the function .
__device__ static void Aliasing_sums_ik(const P3MGpuData p, int NX, int NY, int NZ, REAL_TYPE *Zaehler, REAL_TYPE *Nenner)
__global__ void assign_charge_kernel(P3MGpuData const params, float const *const __restrict__ part_pos, float const *const __restrict__ part_q, unsigned int const parts_per_block)
void p3m_gpu_add_farfield_force(P3MGpuParams &data, GpuParticleData &gpu, double prefactor, std::size_t n_part)
The long-range part of the P3M algorithm.
__global__ void apply_influence_function(const P3MGpuData p)
void assign_charges(P3MGpuData const ¶ms, float const *const __restrict__ part_pos, float const *const __restrict__ part_q)
#define FFT_PLAN_FORW_FLAG
__global__ void apply_diff_op(const P3MGpuData p)
void p3m_gpu_init(std::shared_ptr< P3MGpuParams > &data, int cao, Utils::Vector3i const &mesh, double alpha, Utils::Vector3d const &box_l, std::size_t n_part)
Initialize the internal data structure of the P3M GPU.
static auto p3m_calc_blocks(unsigned int cao, std::size_t n_part)
void assign_forces(P3MGpuData const ¶ms, float const *const __restrict__ part_pos, float const *const __restrict__ part_q, float *const __restrict__ part_f, REAL_TYPE const prefactor)
__global__ void calculate_influence_function_device(const P3MGpuData p)
__device__ int wrap_index(const int ind, const int mesh)
#define FFT_PLAN_BACK_FLAG
dim3 p3m_make_grid(unsigned int n_blocks)
__global__ void assign_forces_kernel(P3MGpuData const params, float const *const __restrict__ part_pos, float const *const __restrict__ part_q, float *const __restrict__ part_f, REAL_TYPE prefactor, unsigned int const parts_per_block)
int mesh[3]
Mesh dimensions.
REAL_TYPE pos_shift
Position shift.
REAL_TYPE * G_hat
Influence Function.
FFT_TYPE_COMPLEX * force_mesh_x
Force meshes.
int mesh_z_padded
Padded size.
int mesh_size
Total number of mesh points (including padding)
int cao
Charge assignment order.
unsigned int n_part
Number of particles.
FFT_TYPE_COMPLEX * force_mesh_z
REAL_TYPE hi[3]
Inverse mesh spacing.
FFT_TYPE_COMPLEX * charge_mesh
Charge mesh.
REAL_TYPE box[3]
Box size.
FFT_TYPE_COMPLEX * force_mesh_y
REAL_TYPE alpha
Ewald parameter.
cufftHandle forw_plan
Forward FFT plan.
cufftHandle back_plan
Backward FFT plan.
void free_device_memory()
#define KERNELCALL(_function, _grid, _block,...)