2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file defines functions for "pair" interactions
39 * (i.e. listed non-bonded interactions, e.g. 1-4 interactions)
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
43 * \ingroup module_listed_forces
51 #include "gromacs/listed_forces/bonded.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdtypes/forcerec.h"
55 #include "gromacs/mdtypes/group.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/mdtypes/mdatom.h"
59 #include "gromacs/mdtypes/nblist.h"
60 #include "gromacs/mdtypes/simulation_workload.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/pbcutil/pbc_simd.h"
64 #include "gromacs/simd/simd.h"
65 #include "gromacs/simd/simd_math.h"
66 #include "gromacs/simd/vector_operations.h"
67 #include "gromacs/tables/forcetable.h"
68 #include "gromacs/utility/basedefinitions.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/gmxassert.h"
72 #include "listed_internal.h"
74 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
76 /*! \brief Issue a warning if a listed interaction is beyond a table limit */
77 static void warning_rlimit(const rvec* x, int ai, int aj, int* global_atom_index, real r, real rlimit)
80 "Listed nonbonded interaction between particles %d and %d\n"
81 "at distance %.3f which is larger than the table limit %.3f nm.\n\n"
82 "This is likely either a 1,4 interaction, or a listed interaction inside\n"
83 "a smaller molecule you are decoupling during a free energy calculation.\n"
84 "Since interactions at distances beyond the table cannot be computed,\n"
85 "they are skipped until they are inside the table limit again. You will\n"
86 "only see this message once, even if it occurs for several interactions.\n\n"
87 "IMPORTANT: This should not happen in a stable simulation, so there is\n"
88 "probably something wrong with your system. Only change the table-extension\n"
89 "distance in the mdp file if you are really sure that is the reason.\n",
90 glatnr(global_atom_index, ai), glatnr(global_atom_index, aj), r, rlimit);
95 "%8f %8f %8f\n%8f %8f %8f\n1-4 (%d,%d) interaction not within cut-off! r=%g. "
97 x[ai][XX], x[ai][YY], x[ai][ZZ], x[aj][XX], x[aj][YY], x[aj][ZZ],
98 glatnr(global_atom_index, ai), glatnr(global_atom_index, aj), r);
102 /*! \brief Compute the energy and force for a single pair interaction */
103 static real evaluate_single(real r2,
113 real rinv, r, rtab, eps, eps2, Y, F, Geps, Heps2, Fp, VVe, FFe, VVd, FFd, VVr, FFr, fscal;
116 /* Do the tabulated interactions - first table lookup */
117 rinv = gmx::invsqrt(r2);
120 ntab = static_cast<int>(rtab);
123 ntab = static_cast<int>(tableStride * ntab);
127 Geps = eps * vftab[ntab + 2];
128 Heps2 = eps2 * vftab[ntab + 3];
129 Fp = F + Geps + Heps2;
131 FFe = Fp + Geps + 2.0 * Heps2;
135 Geps = eps * vftab[ntab + 6];
136 Heps2 = eps2 * vftab[ntab + 7];
137 Fp = F + Geps + Heps2;
139 FFd = Fp + Geps + 2.0 * Heps2;
143 Geps = eps * vftab[ntab + 10];
144 Heps2 = eps2 * vftab[ntab + 11];
145 Fp = F + Geps + Heps2;
147 FFr = Fp + Geps + 2.0 * Heps2;
150 *vvdw = c6 * VVd + c12 * VVr;
152 fscal = -(qq * FFe + c6 * FFd + c12 * FFr) * tabscale * rinv;
157 /*! \brief Compute the energy and force for a single pair interaction under FEP */
158 static real free_energy_evaluate_single(real r2,
174 const real lfac_coul[2],
175 const real lfac_vdw[2],
176 const real dlfac_coul[2],
177 const real dlfac_vdw[2],
186 real rp, rpm2, rtab, eps, eps2, Y, F, Geps, Heps2, Fp, VV, FF, fscal;
187 real qq[2], c6[2], c12[2], sigma6[2], sigma2[2], sigma_pow[2];
188 real alpha_coul_eff, alpha_vdw_eff, dvdl_coul, dvdl_vdw;
189 real rpinv, r_coul, r_vdw, velecsum, vvdwsum;
190 real fscal_vdw[2], fscal_elec[2];
191 real velec[2], vvdw[2];
193 const real half = 0.5;
194 const real minusOne = -1.0;
195 const real one = 1.0;
196 const real two = 2.0;
197 const real six = 6.0;
206 if (sc_r_power == six)
208 rpm2 = r2 * r2; /* r4 */
209 rp = rpm2 * r2; /* r6 */
213 rp = std::pow(r2, half * sc_r_power); /* not currently supported as input, but can handle it */
217 /* Loop over state A(0) and B(1) */
218 for (i = 0; i < 2; i++)
220 if ((c6[i] > 0) && (c12[i] > 0))
222 /* The c6 & c12 coefficients now contain the constants 6.0 and 12.0, respectively.
223 * Correct for this by multiplying with (1/12.0)/(1/6.0)=6.0/12.0=0.5.
225 sigma6[i] = half * c12[i] / c6[i];
226 sigma2[i] = std::cbrt(half * c12[i] / c6[i]);
227 /* should be able to get rid of this ^^^ internal pow call eventually. Will require agreement on
228 what data to store externally. Can't be fixed without larger scale changes, so not 5.0 */
229 if (sigma6[i] < sigma6_min) /* for disappearing coul and vdw with soft core at the same time */
231 sigma6[i] = sigma6_min;
232 sigma2[i] = sigma2_min;
237 sigma6[i] = sigma6_def;
238 sigma2[i] = sigma2_def;
240 if (sc_r_power == six)
242 sigma_pow[i] = sigma6[i];
245 { /* not really supported as input, but in here for testing the general case*/
246 sigma_pow[i] = std::pow(sigma2[i], sc_r_power / 2);
250 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
251 if ((c12[0] > 0) && (c12[1] > 0))
258 alpha_vdw_eff = alpha_vdw;
259 alpha_coul_eff = alpha_coul;
262 /* Loop over A and B states again */
263 for (i = 0; i < 2; i++)
270 /* Only spend time on A or B state if it is non-zero */
271 if ((qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0))
274 rpinv = one / (alpha_coul_eff * lfac_coul[i] * sigma_pow[i] + rp);
275 r_coul = std::pow(rpinv, minusOne / sc_r_power);
277 /* Electrostatics table lookup data */
278 rtab = r_coul * tabscale;
279 ntab = static_cast<int>(rtab);
282 ntab = static_cast<int>(tableStride * ntab);
286 Geps = eps * vftab[ntab + 2];
287 Heps2 = eps2 * vftab[ntab + 3];
288 Fp = F + Geps + Heps2;
290 FF = Fp + Geps + two * Heps2;
291 velec[i] = qq[i] * VV;
292 fscal_elec[i] = -qq[i] * FF * r_coul * rpinv * tabscale;
295 rpinv = one / (alpha_vdw_eff * lfac_vdw[i] * sigma_pow[i] + rp);
296 r_vdw = std::pow(rpinv, minusOne / sc_r_power);
297 /* Vdw table lookup data */
298 rtab = r_vdw * tabscale;
299 ntab = static_cast<int>(rtab);
306 Geps = eps * vftab[ntab + 6];
307 Heps2 = eps2 * vftab[ntab + 7];
308 Fp = F + Geps + Heps2;
310 FF = Fp + Geps + two * Heps2;
311 vvdw[i] = c6[i] * VV;
312 fscal_vdw[i] = -c6[i] * FF;
317 Geps = eps * vftab[ntab + 10];
318 Heps2 = eps2 * vftab[ntab + 11];
319 Fp = F + Geps + Heps2;
321 FF = Fp + Geps + two * Heps2;
322 vvdw[i] += c12[i] * VV;
323 fscal_vdw[i] -= c12[i] * FF;
324 fscal_vdw[i] *= r_vdw * rpinv * tabscale;
327 /* Now we have velec[i], vvdw[i], and fscal[i] for both states */
328 /* Assemble A and B states */
334 for (i = 0; i < 2; i++)
336 velecsum += LFC[i] * velec[i];
337 vvdwsum += LFV[i] * vvdw[i];
339 fscal += (LFC[i] * fscal_elec[i] + LFV[i] * fscal_vdw[i]) * rpm2;
341 dvdl_coul += velec[i] * DLF[i]
342 + LFC[i] * alpha_coul_eff * dlfac_coul[i] * fscal_elec[i] * sigma_pow[i];
343 dvdl_vdw += vvdw[i] * DLF[i]
344 + LFV[i] * alpha_vdw_eff * dlfac_vdw[i] * fscal_vdw[i] * sigma_pow[i];
347 dvdl[efptCOUL] += dvdl_coul;
348 dvdl[efptVDW] += dvdl_vdw;
350 *velectot = velecsum;
356 /*! \brief Calculate pair interactions, supports all types and conditions. */
357 template<BondedKernelFlavor flavor>
358 static real do_pairs_general(int ftype,
360 const t_iatom iatoms[],
361 const t_iparams iparams[],
365 const struct t_pbc* pbc,
369 const t_forcerec* fr,
370 gmx_grppairener_t* grppener,
371 int* global_atom_index)
375 int i, itype, ai, aj, gid;
378 real fscal, velec, vvdw;
379 real* energygrp_elec;
381 static gmx_bool warned_rlimit = FALSE;
382 /* Free energy stuff */
383 gmx_bool bFreeEnergy;
384 real LFC[2], LFV[2], DLF[2], lfac_coul[2], lfac_vdw[2], dlfac_coul[2], dlfac_vdw[2];
385 real qqB, c6B, c12B, sigma2_def, sigma2_min;
391 energygrp_elec = grppener->ener[egCOUL14].data();
392 energygrp_vdw = grppener->ener[egLJ14].data();
395 energygrp_elec = grppener->ener[egCOULSR].data();
396 energygrp_vdw = grppener->ener[egLJSR].data();
399 energygrp_elec = nullptr; /* Keep compiler happy */
400 energygrp_vdw = nullptr; /* Keep compiler happy */
401 gmx_fatal(FARGS, "Unknown function type %d in do_nonbonded14", ftype);
404 if (fr->efep != efepNO)
406 /* Lambda factor for state A=1-lambda and B=lambda */
407 LFC[0] = 1.0 - lambda[efptCOUL];
408 LFV[0] = 1.0 - lambda[efptVDW];
409 LFC[1] = lambda[efptCOUL];
410 LFV[1] = lambda[efptVDW];
412 /*derivative of the lambda factor for state A and B */
417 sigma2_def = std::cbrt(fr->sc_sigma6_def);
418 sigma2_min = std::cbrt(fr->sc_sigma6_min);
420 for (i = 0; i < 2; i++)
422 lfac_coul[i] = (fr->sc_power == 2 ? (1 - LFC[i]) * (1 - LFC[i]) : (1 - LFC[i]));
424 DLF[i] * fr->sc_power / fr->sc_r_power * (fr->sc_power == 2 ? (1 - LFC[i]) : 1);
425 lfac_vdw[i] = (fr->sc_power == 2 ? (1 - LFV[i]) * (1 - LFV[i]) : (1 - LFV[i]));
427 DLF[i] * fr->sc_power / fr->sc_r_power * (fr->sc_power == 2 ? (1 - LFV[i]) : 1);
432 sigma2_min = sigma2_def = 0;
435 /* TODO This code depends on the logic in tables.c that constructs
436 the table layout, which should be made explicit in future
440 "Pair-interaction code that uses GROMACS interaction tables supports exactly 3 tables");
442 fr->pairsTable->interaction == GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP,
443 "Pair interaction kernels need a table with Coulomb, repulsion and dispersion entries");
445 const real epsfac = fr->ic->epsfac;
448 for (i = 0; (i < nbonds);)
453 gid = GID(md->cENER[ai], md->cENER[aj], md->nenergrp);
459 bFreeEnergy = (fr->efep != efepNO
460 && (((md->nPerturbed != 0) && (md->bPerturbed[ai] || md->bPerturbed[aj]))
461 || iparams[itype].lj14.c6A != iparams[itype].lj14.c6B
462 || iparams[itype].lj14.c12A != iparams[itype].lj14.c12B));
463 qq = md->chargeA[ai] * md->chargeA[aj] * epsfac * fr->fudgeQQ;
464 c6 = iparams[itype].lj14.c6A;
465 c12 = iparams[itype].lj14.c12A;
468 qq = iparams[itype].ljc14.qi * iparams[itype].ljc14.qj * epsfac
469 * iparams[itype].ljc14.fqq;
470 c6 = iparams[itype].ljc14.c6;
471 c12 = iparams[itype].ljc14.c12;
474 qq = iparams[itype].ljcnb.qi * iparams[itype].ljcnb.qj * epsfac;
475 c6 = iparams[itype].ljcnb.c6;
476 c12 = iparams[itype].ljcnb.c12;
479 /* Cannot happen since we called gmx_fatal() above in this case */
480 qq = c6 = c12 = 0; /* Keep compiler happy */
484 /* To save flops in the optimized kernels, c6/c12 have 6.0/12.0 derivative prefactors
485 * included in the general nfbp array now. This means the tables are scaled down by the
486 * same factor, so when we use the original c6/c12 parameters from iparams[] they must
492 /* Do we need to apply full periodic boundary conditions? */
495 fshift_index = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
499 fshift_index = CENTRAL;
500 rvec_sub(x[ai], x[aj], dx);
504 if (r2 >= fr->pairsTable->r * fr->pairsTable->r)
506 /* This check isn't race free. But it doesn't matter because if a race occurs the only
507 * disadvantage is that the warning is printed twice */
510 warning_rlimit(x, ai, aj, global_atom_index, sqrt(r2), fr->pairsTable->r);
511 warned_rlimit = TRUE;
518 /* Currently free energy is only supported for F_LJ14, so no need to check for that if we got here */
519 qqB = md->chargeB[ai] * md->chargeB[aj] * epsfac * fr->fudgeQQ;
520 c6B = iparams[itype].lj14.c6B * 6.0;
521 c12B = iparams[itype].lj14.c12B * 12.0;
523 fscal = free_energy_evaluate_single(
524 r2, fr->sc_r_power, fr->sc_alphacoul, fr->sc_alphavdw, fr->pairsTable->scale,
525 fr->pairsTable->data, fr->pairsTable->stride, qq, c6, c12, qqB, c6B, c12B, LFC,
526 LFV, DLF, lfac_coul, lfac_vdw, dlfac_coul, dlfac_vdw, fr->sc_sigma6_def,
527 fr->sc_sigma6_min, sigma2_def, sigma2_min, &velec, &vvdw, dvdl);
531 /* Evaluate tabulated interaction without free energy */
532 fscal = evaluate_single(r2, fr->pairsTable->scale, fr->pairsTable->data,
533 fr->pairsTable->stride, qq, c6, c12, &velec, &vvdw);
536 energygrp_elec[gid] += velec;
537 energygrp_vdw[gid] += vvdw;
538 svmul(fscal, dx, dx);
544 if (computeVirial(flavor))
546 if (fshift_index != CENTRAL)
548 rvec_inc(fshift[fshift_index], dx);
549 rvec_dec(fshift[CENTRAL], dx);
556 /*! \brief Calculate pairs, only for plain-LJ + plain Coulomb normal type.
558 * This function is templated for real/SimdReal and for optimization.
560 template<typename T, int pack_size, typename pbc_type>
561 static void do_pairs_simple(int nbonds,
562 const t_iatom iatoms[],
563 const t_iparams iparams[],
568 const real scale_factor)
570 const int nfa1 = 1 + 2;
576 #if GMX_SIMD_HAVE_REAL
577 alignas(GMX_SIMD_ALIGNMENT) std::int32_t ai[pack_size];
578 alignas(GMX_SIMD_ALIGNMENT) std::int32_t aj[pack_size];
579 alignas(GMX_SIMD_ALIGNMENT) real coeff[3 * pack_size];
581 std::int32_t ai[pack_size];
582 std::int32_t aj[pack_size];
583 real coeff[3 * pack_size];
586 /* nbonds is #pairs*nfa1, here we step pack_size pairs */
587 for (int i = 0; i < nbonds; i += pack_size * nfa1)
589 /* Collect atoms for pack_size pairs.
590 * iu indexes into iatoms, we should not let iu go beyond nbonds.
593 for (int s = 0; s < pack_size; s++)
595 int itype = iatoms[iu];
596 ai[s] = iatoms[iu + 1];
597 aj[s] = iatoms[iu + 2];
599 if (i + s * nfa1 < nbonds)
601 coeff[0 * pack_size + s] = iparams[itype].lj14.c6A;
602 coeff[1 * pack_size + s] = iparams[itype].lj14.c12A;
603 coeff[2 * pack_size + s] = md->chargeA[ai[s]] * md->chargeA[aj[s]];
605 /* Avoid indexing the iatoms array out of bounds.
606 * We pad the coordinate indices with the last atom pair.
608 if (iu + nfa1 < nbonds)
615 /* Pad the coefficient arrays with zeros to get zero forces */
616 coeff[0 * pack_size + s] = 0;
617 coeff[1 * pack_size + s] = 0;
618 coeff[2 * pack_size + s] = 0;
622 /* Load the coordinates */
624 gatherLoadUTranspose<3>(reinterpret_cast<const real*>(x), ai, &xi[XX], &xi[YY], &xi[ZZ]);
625 gatherLoadUTranspose<3>(reinterpret_cast<const real*>(x), aj, &xj[XX], &xj[YY], &xj[ZZ]);
627 T c6 = load<T>(coeff + 0 * pack_size);
628 T c12 = load<T>(coeff + 1 * pack_size);
629 T qq = load<T>(coeff + 2 * pack_size);
631 /* We could save these operations by storing 6*C6,12*C12 */
636 pbc_dx_aiuc(pbc, xi, xj, dr);
638 T rsq = dr[XX] * dr[XX] + dr[YY] * dr[YY] + dr[ZZ] * dr[ZZ];
639 T rinv = gmx::invsqrt(rsq);
640 T rinv2 = rinv * rinv;
641 T rinv6 = rinv2 * rinv2 * rinv2;
643 /* Calculate the Coulomb force * r */
644 T cfr = ef * qq * rinv;
646 /* Calculate the LJ force * r and add it to the Coulomb part */
647 T fr = gmx::fma(fms(c12, rinv6, c6), rinv6, cfr);
649 T finvr = fr * rinv2;
650 T fx = finvr * dr[XX];
651 T fy = finvr * dr[YY];
652 T fz = finvr * dr[ZZ];
654 /* Add the pair forces to the force array.
655 * Note that here we might add multiple force components for some atoms
656 * due to the SIMD padding. But the extra force components are zero.
658 transposeScatterIncrU<4>(reinterpret_cast<real*>(f), ai, fx, fy, fz);
659 transposeScatterDecrU<4>(reinterpret_cast<real*>(f), aj, fx, fy, fz);
663 /*! \brief Calculate all listed pair interactions */
664 void do_pairs(int ftype,
666 const t_iatom iatoms[],
667 const t_iparams iparams[],
671 const struct t_pbc* pbc,
675 const t_forcerec* fr,
676 const bool havePerturbedInteractions,
677 const gmx::StepWorkload& stepWork,
678 gmx_grppairener_t* grppener,
679 int* global_atom_index)
681 if (ftype == F_LJ14 && fr->ic->vdwtype != evdwUSER && !EEL_USER(fr->ic->eeltype)
682 && !havePerturbedInteractions && (!stepWork.computeVirial && !stepWork.computeEnergy))
684 /* We use a fast code-path for plain LJ 1-4 without FEP.
686 * TODO: Add support for energies (straightforward) and virial
687 * in the SIMD template. For the virial it's inconvenient to store
688 * the force sums for the shifts and we should directly calculate
689 * and sum the virial for the shifts. But we should do this
690 * at once for the angles and dihedrals as well.
692 #if GMX_SIMD_HAVE_REAL
693 if (fr->use_simd_kernels)
695 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
696 set_pbc_simd(pbc, pbc_simd);
698 do_pairs_simple<SimdReal, GMX_SIMD_REAL_WIDTH, const real*>(
699 nbonds, iatoms, iparams, x, f, pbc_simd, md, fr->ic->epsfac * fr->fudgeQQ);
704 /* This construct is needed because pbc_dx_aiuc doesn't accept pbc=NULL */
706 const t_pbc* pbc_nonnull;
714 set_pbc(&pbc_no, PbcType::No, nullptr);
715 pbc_nonnull = &pbc_no;
718 do_pairs_simple<real, 1, const t_pbc*>(nbonds, iatoms, iparams, x, f, pbc_nonnull, md,
719 fr->ic->epsfac * fr->fudgeQQ);
722 else if (stepWork.computeVirial)
724 do_pairs_general<BondedKernelFlavor::ForcesAndVirialAndEnergy>(
725 ftype, nbonds, iatoms, iparams, x, f, fshift, pbc, lambda, dvdl, md, fr, grppener,
730 do_pairs_general<BondedKernelFlavor::ForcesAndEnergy>(ftype, nbonds, iatoms, iparams, x, f,
731 fshift, pbc, lambda, dvdl, md, fr,
732 grppener, global_atom_index);