Tinker9 70bd052 (Thu Nov 9 12:11:35 2023 -0800)
Loading...
Searching...
No Matches
emselfamoeba.h
1#pragma once
2#include "ff/amoeba/mpole.h"
3#include "ff/energybuffer.h"
4#include "seq/add.h"
5
6namespace tinker {
7template <bool do_a>
8__global__
10 const real (*restrict rpole)[10], int n, real f, real aewald)
11{
12 real aewald_sq_2 = 2 * aewald * aewald;
13 real fterm = -f * aewald * 0.5f * (real)(M_2_SQRTPI);
14
15 for (int i = threadIdx.x + blockIdx.x * blockDim.x; i < n;
16 i += blockDim.x * gridDim.x) {
17 real ci = rpole[i][MPL_PME_0];
18 real dix = rpole[i][MPL_PME_X];
19 real diy = rpole[i][MPL_PME_Y];
20 real diz = rpole[i][MPL_PME_Z];
21 real qixx = rpole[i][MPL_PME_XX];
22 real qixy = rpole[i][MPL_PME_XY];
23 real qixz = rpole[i][MPL_PME_XZ];
24 real qiyy = rpole[i][MPL_PME_YY];
25 real qiyz = rpole[i][MPL_PME_YZ];
26 real qizz = rpole[i][MPL_PME_ZZ];
27
28 real cii = ci * ci;
29 real dii = dix * dix + diy * diy + diz * diz;
30 real qii = 2 * (qixy * qixy + qixz * qixz + qiyz * qiyz) + qixx * qixx
31 + qiyy * qiyy + qizz * qizz;
32
33 int offset = threadIdx.x + blockIdx.x * blockDim.x;
34 real e = fterm
35 * (cii + aewald_sq_2 * (dii / 3 + 2 * aewald_sq_2 * qii * (real)0.2));
36 atomic_add(e, em, offset);
37 if CONSTEXPR (do_a) atomic_add(1, nem, offset);
38 }
39}
40}
void atomic_add(T value, T *buffer, size_t offset=0)
Definition: acc/adddef.h:14
#define restrict
Definition: macro.h:51
#define CONSTEXPR
Definition: macro.h:61
CountBufferTraits::type * CountBuffer
Definition: energybuffer.h:92
EnergyBufferTraits::type * EnergyBuffer
Definition: energybuffer.h:93
int n
Number of atoms padded by WARP_SIZE.
float real
Definition: precision.h:80
Definition: testrt.h:9
@ MPL_PME_XY
Definition: mpole.h:13
@ MPL_PME_Y
Definition: mpole.h:8
@ MPL_PME_0
Definition: mpole.h:6
@ MPL_PME_ZZ
Definition: mpole.h:12
@ MPL_PME_YZ
Definition: mpole.h:15
@ MPL_PME_XX
Definition: mpole.h:10
@ MPL_PME_XZ
Definition: mpole.h:14
@ MPL_PME_X
Definition: mpole.h:7
@ MPL_PME_YY
Definition: mpole.h:11
@ MPL_PME_Z
Definition: mpole.h:9
CountBuffer nem
__global__ void empoleSelf_cu(CountBuffer __restrict__ nem, EnergyBuffer __restrict__ em, const real(*__restrict__ rpole)[10], int n, real f, real aewald)
Definition: emselfamoeba.h:9
real(* rpole)[MPL_TOTAL]
EnergyBuffer em