6#include "seq/dampaplus.h"
20 real r = REAL_SQRT(r2);
21 real invr1 = REAL_RECIP(r);
22 real rr2 = invr1 * invr1;
24 real scale3, scale5, scale7;
25 damp_aplus3(r, pdi, ddi, pdk, ddk, scale3, scale5, scale7);
28 if CONSTEXPR (eq<ETYP, EWALD>()) damp_ewald<4>(bn, r, invr1, rr2, aewald);
31 real rr5 = 3 * rr1 * rr2 * rr2;
32 real rr7 = 15 * rr1 * rr2 * rr2 * rr2;
34 real dir = dix * xr + diy * yr + diz * zr;
35 real qix = qixx * xr + qixy * yr + qixz * zr;
36 real qiy = qixy * xr + qiyy * yr + qiyz * zr;
37 real qiz = qixz * xr + qiyz * yr + qizz * zr;
38 real qir = qix * xr + qiy * yr + qiz * zr;
39 real dkr = dkx * xr + dky * yr + dkz * zr;
40 real qkx = qkxx * xr + qkxy * yr + qkxz * zr;
41 real qky = qkxy * xr + qkyy * yr + qkyz * zr;
42 real qkz = qkxz * xr + qkyz * yr + qkzz * zr;
43 real qkr = qkx * xr + qky * yr + qkz * zr;
45 real3 dixyz = make_real3(dix, diy, diz);
46 real3 dkxyz = make_real3(dkx, dky, dkz);
47 real3 qixyz = make_real3(qix, qiy, qiz);
48 real3 qkxyz = make_real3(qkx, qky, qkz);
49 real3 dr = make_real3(xr, yr, zr);
56 bn[1] -= (1 - scale3) * rr3;
57 bn[2] -= (1 - scale5) * rr5;
58 bn[3] -= (1 - scale7) * rr7;
59 }
else if CONSTEXPR (eq<ETYP, NON_EWALD>()) {
60 bn[1] = dscale * scale3 * rr3;
61 bn[2] = dscale * scale5 * rr5;
62 bn[3] = dscale * scale7 * rr7;
65 c1 = -(bn[1] * ck - bn[2] * dkr + bn[3] * qkr);
66 inci = c1 * dr - bn[1] * dkxyz + 2 * bn[2] * qkxyz;
71 c1 = (bn[1] * ci + bn[2] * dir + bn[3] * qir);
72 inck = c1 * dr - bn[1] * dixyz - 2 * bn[2] * qixyz;
78#pragma acc routine seq
89 real r = REAL_SQRT(r2);
90 real invr1 = REAL_RECIP(r);
91 real rr2 = invr1 * invr1;
93 real scale3, scale5, scale7;
94 damp_aplus3(r, pdi, ddi, pdk, ddk, scale3, scale5, scale7);
97 if CONSTEXPR (eq<ETYP, EWALD>()) damp_ewald<4>(bn, r, invr1, rr2, aewald);
100 real rr5 = 3 * rr1 * rr2 * rr2;
101 real rr7 = 15 * rr1 * rr2 * rr2 * rr2;
103 real dir = dix * xr + diy * yr + diz * zr;
104 real qix = qixx * xr + qixy * yr + qixz * zr;
105 real qiy = qixy * xr + qiyy * yr + qiyz * zr;
106 real qiz = qixz * xr + qiyz * yr + qizz * zr;
107 real qir = qix * xr + qiy * yr + qiz * zr;
108 real dkr = dkx * xr + dky * yr + dkz * zr;
109 real qkx = qkxx * xr + qkxy * yr + qkxz * zr;
110 real qky = qkxy * xr + qkyy * yr + qkyz * zr;
111 real qkz = qkxz * xr + qkyz * yr + qkzz * zr;
112 real qkr = qkx * xr + qky * yr + qkz * zr;
114 real3 dixyz = make_real3(dix, diy, diz);
115 real3 dkxyz = make_real3(dkx, dky, dkz);
116 real3 qixyz = make_real3(qix, qiy, qiz);
117 real3 qkxyz = make_real3(qkx, qky, qkz);
118 real3 dr = make_real3(xr, yr, zr);
124 bn[1] -= (1 - dscale * scale3) * rr3;
125 bn[2] -= (1 - dscale * scale5) * rr5;
126 bn[3] -= (1 - dscale * scale7) * rr7;
127 }
else if CONSTEXPR (eq<ETYP, NON_EWALD>()) {
128 bn[1] = dscale * scale3 * rr3;
129 bn[2] = dscale * scale5 * rr5;
130 bn[3] = dscale * scale7 * rr7;
133 c1 = -(bn[1] * ck - bn[2] * dkr + bn[3] * qkr);
134 inci = c1 * dr - bn[1] * dkxyz + 2 * bn[2] * qkxyz;
139 c1 = (bn[1] * ci + bn[2] * dir + bn[3] * qir);
140 inck = c1 * dr - bn[1] * dixyz - 2 * bn[2] * qixyz;
146#pragma acc routine seq
157 real r = REAL_SQRT(r2);
158 real invr1 = REAL_RECIP(r);
159 real rr2 = invr1 * invr1;
162 damp_thole2(r, pdi, pti, pdk, ptk, scale3, scale5);
165 if CONSTEXPR (eq<ETYP, EWALD>()) damp_ewald<3>(bn, r, invr1, rr2, aewald);
167 real rr3 = rr1 * rr2;
168 real rr5 = 3 * rr1 * rr2 * rr2;
171 bn[1] -= (1 - scale3) * rr3;
172 bn[2] -= (1 - scale5) * rr5;
173 }
else if CONSTEXPR (eq<ETYP, NON_EWALD>()) {
174 bn[1] = uscale * scale3 * rr3;
175 bn[2] = uscale * scale5 * rr5;
180 real3 dr = make_real3(xr, yr, zr);
181 real3 uid = make_real3(uindi0, uindi1, uindi2);
182 real3 ukd = make_real3(uindk0, uindk1, uindk2);
184 coef = bn[2] *
dot3(dr, ukd);
185 inci = coef * dr - bn[1] * ukd;
190 coef = bn[2] *
dot3(dr, uid);
191 inck = coef * dr - bn[1] * uid;
197#pragma acc routine seq
208 real r = REAL_SQRT(r2);
209 real invr1 = REAL_RECIP(r);
210 real rr2 = invr1 * invr1;
213 damp_thole2(r, pdi, pti, pdk, ptk, scale3, scale5);
216 if CONSTEXPR (eq<ETYP, EWALD>()) damp_ewald<3>(bn, r, invr1, rr2, aewald);
218 real rr3 = rr1 * rr2;
219 real rr5 = 3 * rr1 * rr2 * rr2;
222 bn[1] -= (1 - uscale * scale3) * rr3;
223 bn[2] -= (1 - uscale * scale5) * rr5;
224 }
else if CONSTEXPR (eq<ETYP, NON_EWALD>()) {
225 bn[1] = uscale * scale3 * rr3;
226 bn[2] = uscale * scale5 * rr5;
231 real3 dr = make_real3(xr, yr, zr);
232 real3 uid = make_real3(uindi0, uindi1, uindi2);
233 real3 ukd = make_real3(uindk0, uindk1, uindk2);
235 coef = bn[2] *
dot3(dr, ukd);
236 inci = coef * dr - bn[1] * ukd;
241 coef = bn[2] *
dot3(dr, uid);
242 inck = coef * dr - bn[1] * uid;
float x
Definition: acc/realndef.h:49
float y
Definition: acc/realndef.h:49
float z
Definition: acc/realndef.h:49
#define SEQ_CUDA
Definition: acc/seqdef.h:12
Definition: acc/realndef.h:48
#define restrict
Definition: macro.h:51
#define CONSTEXPR
Definition: macro.h:61
__device__ real dot3(real3 a, real3 b)
Definition: realn.h:84
float real
Definition: precision.h:80
__device__ void pair_dfield_aplus_v2(real r2, real xr, real yr, real zr, real dscale, real ci, real dix, real diy, real diz, real pdi, real ddi, real qixx, real qixy, real qixz, real qiyy, real qiyz, real qizz, real ck, real dkx, real dky, real dkz, real pdk, real ddk, real qkxx, real qkxy, real qkxz, real qkyy, real qkyz, real qkzz, real aewald, real &__restrict__ fidx, real &__restrict__ fidy, real &__restrict__ fidz, real &__restrict__ fkdx, real &__restrict__ fkdy, real &__restrict__ fkdz)
Definition: pairfieldaplus.h:81
__device__ void pair_ufield_aplus(real r2, real xr, real yr, real zr, real uscale, real uindi0, real uindi1, real uindi2, real pdi, real pti, real uindk0, real uindk1, real uindk2, real pdk, real ptk, real aewald, real &__restrict__ fidx, real &__restrict__ fidy, real &__restrict__ fidz, real &__restrict__ fkdx, real &__restrict__ fkdy, real &__restrict__ fkdz)
Definition: pairfieldaplus.h:149
__device__ void pair_ufield_aplus_v2(real r2, real xr, real yr, real zr, real uscale, real uindi0, real uindi1, real uindi2, real pdi, real pti, real uindk0, real uindk1, real uindk2, real pdk, real ptk, real aewald, real &__restrict__ fidx, real &__restrict__ fidy, real &__restrict__ fidz, real &__restrict__ fkdx, real &__restrict__ fkdy, real &__restrict__ fkdz)
Definition: pairfieldaplus.h:200
__device__ void pair_dfield_aplus(real r2, real xr, real yr, real zr, real dscale, real ci, real dix, real diy, real diz, real pdi, real ddi, real qixx, real qixy, real qixz, real qiyy, real qiyz, real qizz, real ck, real dkx, real dky, real dkz, real pdk, real ddk, real qkxx, real qkxy, real qkxz, real qkyy, real qkyz, real qkzz, real aewald, real &__restrict__ fidx, real &__restrict__ fidy, real &__restrict__ fidz, real &__restrict__ fkdx, real &__restrict__ fkdy, real &__restrict__ fkdz)
Definition: pairfieldaplus.h:12
__device__ void damp_aplus3(real r, real pdi, real ddi, real pdk, real ddk, real &__restrict__ scale3, real &__restrict__ scale5, real &__restrict__ scale7)
Definition: dampaplus.h:131
__device__ void damp_thole2(real r, real pdi, real pti, real pdk, real ptk, real &__restrict__ scale3, real &__restrict__ scale5)
Definition: damp.h:21