2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /* When calculating RF or Ewald interactions we calculate the electrostatic
37 * forces and energies on excluded atom pairs here in the non-bonded loops.
39 #if defined CHECK_EXCLS && (defined CALC_COULOMB || defined LJ_EWALD)
53 egp_cj = nbatParams.energrp[cj];
55 for (i = 0; i < UNROLLI; i++)
63 type_i_off = type[ai]*ntype2;
65 for (j = 0; j < UNROLLJ; j++)
72 real FrLJ6 = 0, FrLJ12 = 0, frLJ = 0;
74 #if defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
93 /* A multiply mask used to zero an interaction
94 * when either the distance cutoff is exceeded, or
95 * (if appropriate) the i and j indices are
96 * unsuitable for this kind of inner loop. */
100 /* A multiply mask used to zero an interaction
101 * when that interaction should be excluded
102 * (e.g. because of bonding). */
103 real interact = static_cast<real>((l_cj[cjind].excl>>(i*UNROLLI + j)) & 1);
107 skipmask = (cj == ci_sh && j <= i) ? 0.0 : 1.0;
110 constexpr real interact = 1.0;
118 dx = xi[i*XI_STRIDE+XX] - x[aj*X_STRIDE+XX];
119 dy = xi[i*XI_STRIDE+YY] - x[aj*X_STRIDE+YY];
120 dz = xi[i*XI_STRIDE+ZZ] - x[aj*X_STRIDE+ZZ];
122 rsq = dx*dx + dy*dy + dz*dz;
124 /* Prepare to enforce the cut-off. */
125 skipmask = (rsq >= rcut2) ? 0 : skipmask;
126 /* 9 flops for r^2 + cut-off check */
128 // Ensure the distances do not fall below the limit where r^-12 overflows.
129 // This should never happen for normal interactions.
130 rsq = std::max(rsq, NBNXN_MIN_RSQ);
136 rinv = gmx::invsqrt(rsq);
137 /* 5 flops for invsqrt */
139 /* Partially enforce the cut-off (and perhaps
140 * exclusions) to avoid possible overflow of
141 * rinvsix when computing LJ, and/or overflowing
142 * the Coulomb table during lookup. */
143 rinv = rinv * skipmask;
151 c6 = nbfp[type_i_off+type[aj]*2 ];
152 c12 = nbfp[type_i_off+type[aj]*2+1];
154 #if defined LJ_CUT || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
155 rinvsix = interact*rinvsq*rinvsq*rinvsq;
157 FrLJ12 = c12*rinvsix*rinvsix;
158 frLJ = FrLJ12 - FrLJ6;
159 /* 7 flops for r^-2 + LJ force */
160 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
161 VLJ = (FrLJ12 + c12*ic->repulsion_shift.cpot)/12 -
162 (FrLJ6 + c6*ic->dispersion_shift.cpot)/6;
163 /* 7 flops for LJ energy */
167 #if defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
168 /* Force or potential switching from ic->rvdw_switch */
170 rsw = r - ic->rvdw_switch;
171 rsw = (rsw >= 0.0 ? rsw : 0.0);
173 #ifdef LJ_FORCE_SWITCH
175 -c6*(ic->dispersion_shift.c2 + ic->dispersion_shift.c3*rsw)*rsw*rsw*r
176 + c12*(ic->repulsion_shift.c2 + ic->repulsion_shift.c3*rsw)*rsw*rsw*r;
177 #if defined CALC_ENERGIES
179 -c6*(-ic->dispersion_shift.c2/3 - ic->dispersion_shift.c3/4*rsw)*rsw*rsw*rsw
180 + c12*(-ic->repulsion_shift.c2/3 - ic->repulsion_shift.c3/4*rsw)*rsw*rsw*rsw;
184 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
185 /* Masking should be done after force switching,
186 * but before potential switching.
188 /* Need to zero the interaction if there should be exclusion. */
189 VLJ = VLJ * interact;
196 sw = 1.0 + (swV3 + (swV4+ swV5*rsw)*rsw)*rsw*rsw*rsw;
197 dsw = (swF2 + (swF3 + swF4*rsw)*rsw)*rsw*rsw;
199 frLJ = frLJ*sw - r*VLJ*dsw;
206 real c6grid, rinvsix_nm, cr2, expmcr2, poly;
211 #ifdef LJ_EWALD_COMB_GEOM
212 c6grid = ljc[type[ai]*2]*ljc[type[aj]*2];
213 #elif defined LJ_EWALD_COMB_LB
215 real sigma, sigma2, epsilon;
217 /* These sigma and epsilon are scaled to give 6*C6 */
218 sigma = ljc[type[ai]*2] + ljc[type[aj]*2];
219 epsilon = ljc[type[ai]*2+1]*ljc[type[aj]*2+1];
221 sigma2 = sigma*sigma;
222 c6grid = epsilon*sigma2*sigma2*sigma2;
225 #error "No LJ Ewald combination rule defined"
229 /* Recalculate rinvsix without exclusion mask */
230 rinvsix_nm = rinvsq*rinvsq*rinvsq;
232 rinvsix_nm = rinvsix;
234 cr2 = lje_coeff2*rsq;
238 expmcr2 = expf(-cr2);
240 poly = 1 + cr2 + 0.5*cr2*cr2;
242 /* Subtract the grid force from the total LJ force */
243 frLJ += c6grid*(rinvsix_nm - expmcr2*(rinvsix_nm*poly + lje_coeff6_6));
245 /* Shift should only be applied to real LJ pairs */
246 sh_mask = lje_vc*interact;
248 VLJ += c6grid/6*(rinvsix_nm*(1 - expmcr2*poly) + sh_mask);
251 #endif /* LJ_EWALD */
253 #ifdef VDW_CUTOFF_CHECK
254 /* Mask for VdW cut-off shorter than Coulomb cut-off */
258 skipmask_rvdw = (rsq < rvdw2) ? 1.0 : 0.0;
259 frLJ *= skipmask_rvdw;
261 VLJ *= skipmask_rvdw;
265 #if defined CALC_ENERGIES
266 /* Need to zero the interaction if r >= rcut */
267 VLJ = VLJ * skipmask;
268 /* 1 more flop for LJ energy */
270 #endif /* VDW_CUTOFF_CHECK */
275 Vvdw[egp_sh_i[i] + ((egp_cj >> (nbatParams.neg_2log*j)) & egp_mask)] += VLJ;
278 /* 1 flop for LJ energy addition */
284 /* Enforce the cut-off and perhaps exclusions. In
285 * those cases, rinv is zero because of skipmask,
286 * but fcoul and vcoul will later be non-zero (in
287 * both RF and table cases) because of the
288 * contributions that do not depend on rinv. These
289 * contributions cannot be allowed to accumulate
290 * to the force and potential, and the easiest way
291 * to do this is to zero the charges in
293 qq = skipmask * qi[i] * q[aj];
296 fcoul = qq*(interact*rinv*rinvsq - k_rf2);
297 /* 4 flops for RF force */
299 vcoul = qq*(interact*rinv + k_rf*rsq - c_rf);
300 /* 4 flops for RF energy */
305 rs = rsq*rinv*tab_coul_scale;
307 frac = rs - static_cast<real>(ri);
309 /* fexcl = F_i + frac * (F_(i+1)-F_i) */
310 fexcl = tab_coul_FDV0[ri*4] + frac*tab_coul_FDV0[ri*4+1];
312 /* fexcl = (1-frac) * F_i + frac * F_(i+1) */
313 fexcl = (1 - frac)*tab_coul_F[ri] + frac*tab_coul_F[ri+1];
315 fcoul = interact*rinvsq - fexcl;
316 /* 7 flops for float 1/r-table force */
319 vcoul = qq*(interact*(rinv - ic->sh_ewald)
320 -(tab_coul_FDV0[ri*4+2]
321 -halfsp*frac*(tab_coul_FDV0[ri*4] + fexcl)));
322 /* 7 flops for float 1/r-table energy (8 with excls) */
324 vcoul = qq*(interact*(rinv - ic->sh_ewald)
326 -halfsp*frac*(tab_coul_F[ri] + fexcl)));
334 Vc[egp_sh_i[i] + ((egp_cj >> (nbatParams.neg_2log*j)) & egp_mask)] += vcoul;
337 /* 1 flop for Coulomb energy addition */
347 fscal = frLJ*rinvsq + fcoul;
348 /* 2 flops for scalar LJ+Coulomb force */
363 /* Increment i-atom force */
364 fi[i*FI_STRIDE+XX] += fx;
365 fi[i*FI_STRIDE+YY] += fy;
366 fi[i*FI_STRIDE+ZZ] += fz;
367 /* Decrement j-atom force */
368 f[aj*F_STRIDE+XX] -= fx;
369 f[aj*F_STRIDE+YY] -= fy;
370 f[aj*F_STRIDE+ZZ] -= fz;
371 /* 9 flops for force addition */