2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/compat/make_unique.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/ewald/ewald.h"
55 #include "gromacs/ewald/ewald-utils.h"
56 #include "gromacs/fileio/filetypes.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/hardware/hw_info.h"
61 #include "gromacs/listed-forces/manage-threading.h"
62 #include "gromacs/listed-forces/pairs.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/units.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/force.h"
68 #include "gromacs/mdlib/forcerec-threading.h"
69 #include "gromacs/mdlib/gmx_omp_nthreads.h"
70 #include "gromacs/mdlib/md_support.h"
71 #include "gromacs/mdlib/nb_verlet.h"
72 #include "gromacs/mdlib/nbnxn_atomdata.h"
73 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
74 #include "gromacs/mdlib/nbnxn_internal.h"
75 #include "gromacs/mdlib/nbnxn_search.h"
76 #include "gromacs/mdlib/nbnxn_simd.h"
77 #include "gromacs/mdlib/nbnxn_tuning.h"
78 #include "gromacs/mdlib/nbnxn_util.h"
79 #include "gromacs/mdlib/ns.h"
80 #include "gromacs/mdlib/qmmm.h"
81 #include "gromacs/mdlib/rf_util.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/wall.h"
84 #include "gromacs/mdtypes/commrec.h"
85 #include "gromacs/mdtypes/fcdata.h"
86 #include "gromacs/mdtypes/group.h"
87 #include "gromacs/mdtypes/iforceprovider.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/pbcutil/ishift.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/simd/simd.h"
93 #include "gromacs/tables/forcetable.h"
94 #include "gromacs/topology/mtop_util.h"
95 #include "gromacs/trajectory/trajectoryframe.h"
96 #include "gromacs/utility/cstringutil.h"
97 #include "gromacs/utility/exceptions.h"
98 #include "gromacs/utility/fatalerror.h"
99 #include "gromacs/utility/gmxassert.h"
100 #include "gromacs/utility/logger.h"
101 #include "gromacs/utility/physicalnodecommunicator.h"
102 #include "gromacs/utility/pleasecite.h"
103 #include "gromacs/utility/smalloc.h"
104 #include "gromacs/utility/strconvert.h"
106 #include "nbnxn_gpu_jit_support.h"
108 t_forcerec *mk_forcerec()
117 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
125 snew(nbfp, 3*atnr*atnr);
126 for (i = k = 0; (i < atnr); i++)
128 for (j = 0; (j < atnr); j++, k++)
130 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
131 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
132 /* nbfp now includes the 6.0 derivative prefactor */
133 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
139 snew(nbfp, 2*atnr*atnr);
140 for (i = k = 0; (i < atnr); i++)
142 for (j = 0; (j < atnr); j++, k++)
144 /* nbfp now includes the 6.0/12.0 derivative prefactors */
145 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
146 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
154 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
157 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
160 /* For LJ-PME simulations, we correct the energies with the reciprocal space
161 * inside of the cut-off. To do this the non-bonded kernels needs to have
162 * access to the C6-values used on the reciprocal grid in pme.c
166 snew(grid, 2*atnr*atnr);
167 for (i = k = 0; (i < atnr); i++)
169 for (j = 0; (j < atnr); j++, k++)
171 c6i = idef->iparams[i*(atnr+1)].lj.c6;
172 c12i = idef->iparams[i*(atnr+1)].lj.c12;
173 c6j = idef->iparams[j*(atnr+1)].lj.c6;
174 c12j = idef->iparams[j*(atnr+1)].lj.c12;
175 c6 = std::sqrt(c6i * c6j);
176 if (fr->ljpme_combination_rule == eljpmeLB
177 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
179 sigmai = gmx::sixthroot(c12i / c6i);
180 sigmaj = gmx::sixthroot(c12j / c6j);
181 epsi = c6i * c6i / c12i;
182 epsj = c6j * c6j / c12j;
183 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
185 /* Store the elements at the same relative positions as C6 in nbfp in order
186 * to simplify access in the kernels
188 grid[2*(atnr*i+j)] = c6*6.0;
194 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
198 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
202 snew(nbfp, 2*atnr*atnr);
203 for (i = 0; i < atnr; ++i)
205 for (j = 0; j < atnr; ++j)
207 c6i = idef->iparams[i*(atnr+1)].lj.c6;
208 c12i = idef->iparams[i*(atnr+1)].lj.c12;
209 c6j = idef->iparams[j*(atnr+1)].lj.c6;
210 c12j = idef->iparams[j*(atnr+1)].lj.c12;
211 c6 = std::sqrt(c6i * c6j);
212 c12 = std::sqrt(c12i * c12j);
213 if (comb_rule == eCOMB_ARITHMETIC
214 && !gmx_numzero(c6) && !gmx_numzero(c12))
216 sigmai = gmx::sixthroot(c12i / c6i);
217 sigmaj = gmx::sixthroot(c12j / c6j);
218 epsi = c6i * c6i / c12i;
219 epsj = c6j * c6j / c12j;
220 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
221 c12 = std::sqrt(epsi * epsj) * gmx::power12(0.5*(sigmai+sigmaj));
223 C6(nbfp, atnr, i, j) = c6*6.0;
224 C12(nbfp, atnr, i, j) = c12*12.0;
230 /* This routine sets fr->solvent_opt to the most common solvent in the
231 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
232 * the fr->solvent_type array with the correct type (or esolNO).
234 * Charge groups that fulfill the conditions but are not identical to the
235 * most common one will be marked as esolNO in the solvent_type array.
237 * TIP3p is identical to SPC for these purposes, so we call it
238 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
240 * NOTE: QM particle should not
241 * become an optimized solvent. Not even if there is only one charge
251 } solvent_parameters_t;
254 check_solvent_cg(const gmx_moltype_t *molt,
257 const unsigned char *qm_grpnr,
258 const t_grps *qm_grps,
260 int *n_solvent_parameters,
261 solvent_parameters_t **solvent_parameters_p,
271 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
272 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
275 solvent_parameters_t *solvent_parameters;
277 /* We use a list with parameters for each solvent type.
278 * Every time we discover a new molecule that fulfills the basic
279 * conditions for a solvent we compare with the previous entries
280 * in these lists. If the parameters are the same we just increment
281 * the counter for that type, and otherwise we create a new type
282 * based on the current molecule.
284 * Once we've finished going through all molecules we check which
285 * solvent is most common, and mark all those molecules while we
286 * clear the flag on all others.
289 solvent_parameters = *solvent_parameters_p;
291 /* Mark the cg first as non optimized */
294 /* Check if this cg has no exclusions with atoms in other charge groups
295 * and all atoms inside the charge group excluded.
296 * We only have 3 or 4 atom solvent loops.
298 if (GET_CGINFO_EXCL_INTER(cginfo) ||
299 !GET_CGINFO_EXCL_INTRA(cginfo))
304 /* Get the indices of the first atom in this charge group */
305 j0 = molt->cgs.index[cg0];
306 j1 = molt->cgs.index[cg0+1];
308 /* Number of atoms in our molecule */
314 "Moltype '%s': there are %d atoms in this charge group\n",
318 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
321 if (nj < 3 || nj > 4)
326 /* Check if we are doing QM on this group */
328 if (qm_grpnr != nullptr)
330 for (j = j0; j < j1 && !qm; j++)
332 qm = (qm_grpnr[j] < qm_grps->nr - 1);
335 /* Cannot use solvent optimization with QM */
341 atom = molt->atoms.atom;
343 /* Still looks like a solvent, time to check parameters */
345 /* If it is perturbed (free energy) we can't use the solvent loops,
346 * so then we just skip to the next molecule.
350 for (j = j0; j < j1 && !perturbed; j++)
352 perturbed = PERTURBED(atom[j]);
360 /* Now it's only a question if the VdW and charge parameters
361 * are OK. Before doing the check we compare and see if they are
362 * identical to a possible previous solvent type.
363 * First we assign the current types and charges.
365 for (j = 0; j < nj; j++)
367 tmp_vdwtype[j] = atom[j0+j].type;
368 tmp_charge[j] = atom[j0+j].q;
371 /* Does it match any previous solvent type? */
372 for (k = 0; k < *n_solvent_parameters; k++)
377 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
378 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
379 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
384 /* Check that types & charges match for all atoms in molecule */
385 for (j = 0; j < nj && match; j++)
387 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
391 if (tmp_charge[j] != solvent_parameters[k].charge[j])
398 /* Congratulations! We have a matched solvent.
399 * Flag it with this type for later processing.
402 solvent_parameters[k].count += nmol;
404 /* We are done with this charge group */
409 /* If we get here, we have a tentative new solvent type.
410 * Before we add it we must check that it fulfills the requirements
411 * of the solvent optimized loops. First determine which atoms have
414 for (j = 0; j < nj; j++)
417 tjA = tmp_vdwtype[j];
419 /* Go through all other tpes and see if any have non-zero
420 * VdW parameters when combined with this one.
422 for (k = 0; k < fr->ntype && (!has_vdw[j]); k++)
424 /* We already checked that the atoms weren't perturbed,
425 * so we only need to check state A now.
429 has_vdw[j] = (has_vdw[j] ||
430 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
431 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
432 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
437 has_vdw[j] = (has_vdw[j] ||
438 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
439 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
444 /* Now we know all we need to make the final check and assignment. */
448 * For this we require thatn all atoms have charge,
449 * the charges on atom 2 & 3 should be the same, and only
450 * atom 1 might have VdW.
454 tmp_charge[0] != 0 &&
455 tmp_charge[1] != 0 &&
456 tmp_charge[2] == tmp_charge[1])
458 srenew(solvent_parameters, *n_solvent_parameters+1);
459 solvent_parameters[*n_solvent_parameters].model = esolSPC;
460 solvent_parameters[*n_solvent_parameters].count = nmol;
461 for (k = 0; k < 3; k++)
463 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
464 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
467 *cg_sp = *n_solvent_parameters;
468 (*n_solvent_parameters)++;
473 /* Or could it be a TIP4P?
474 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
475 * Only atom 1 mght have VdW.
480 tmp_charge[0] == 0 &&
481 tmp_charge[1] != 0 &&
482 tmp_charge[2] == tmp_charge[1] &&
485 srenew(solvent_parameters, *n_solvent_parameters+1);
486 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
487 solvent_parameters[*n_solvent_parameters].count = nmol;
488 for (k = 0; k < 4; k++)
490 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
491 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
494 *cg_sp = *n_solvent_parameters;
495 (*n_solvent_parameters)++;
499 *solvent_parameters_p = solvent_parameters;
503 check_solvent(FILE * fp,
504 const gmx_mtop_t * mtop,
506 cginfo_mb_t *cginfo_mb)
509 const gmx_moltype_t *molt;
510 int mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
511 int n_solvent_parameters;
512 solvent_parameters_t *solvent_parameters;
518 fprintf(debug, "Going to determine what solvent types we have.\n");
521 n_solvent_parameters = 0;
522 solvent_parameters = nullptr;
523 /* Allocate temporary array for solvent type */
524 snew(cg_sp, mtop->molblock.size());
527 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
529 molt = &mtop->moltype[mtop->molblock[mb].type];
531 /* Here we have to loop over all individual molecules
532 * because we need to check for QMMM particles.
534 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
535 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
536 nmol = mtop->molblock[mb].nmol/nmol_ch;
537 for (mol = 0; mol < nmol_ch; mol++)
540 am = mol*cgs->index[cgs->nr];
541 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
543 check_solvent_cg(molt, cg_mol, nmol,
544 mtop->groups.grpnr[egcQMMM] ?
545 mtop->groups.grpnr[egcQMMM]+at_offset+am : nullptr,
546 &mtop->groups.grps[egcQMMM],
548 &n_solvent_parameters, &solvent_parameters,
549 cginfo_mb[mb].cginfo[cgm+cg_mol],
550 &cg_sp[mb][cgm+cg_mol]);
553 at_offset += cgs->index[cgs->nr];
556 /* Puh! We finished going through all charge groups.
557 * Now find the most common solvent model.
560 /* Most common solvent this far */
562 for (i = 0; i < n_solvent_parameters; i++)
565 solvent_parameters[i].count > solvent_parameters[bestsp].count)
573 bestsol = solvent_parameters[bestsp].model;
581 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
583 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
584 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
585 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
587 if (cg_sp[mb][i] == bestsp)
589 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
594 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
601 if (bestsol != esolNO && fp != nullptr)
603 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
605 solvent_parameters[bestsp].count);
608 sfree(solvent_parameters);
609 fr->solvent_opt = bestsol;
613 acNONE = 0, acCONSTRAINT, acSETTLE
616 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
617 t_forcerec *fr, gmx_bool bNoSolvOpt,
618 gmx_bool *bFEP_NonBonded,
619 gmx_bool *bExcl_IntraCGAll_InterCGNone)
622 const t_blocka *excl;
623 const gmx_moltype_t *molt;
624 const gmx_molblock_t *molb;
625 cginfo_mb_t *cginfo_mb;
628 int cg_offset, a_offset;
629 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
633 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
635 snew(cginfo_mb, mtop->molblock.size());
637 snew(type_VDW, fr->ntype);
638 for (ai = 0; ai < fr->ntype; ai++)
640 type_VDW[ai] = FALSE;
641 for (j = 0; j < fr->ntype; j++)
643 type_VDW[ai] = type_VDW[ai] ||
645 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
646 C12(fr->nbfp, fr->ntype, ai, j) != 0;
650 *bFEP_NonBonded = FALSE;
651 *bExcl_IntraCGAll_InterCGNone = TRUE;
654 snew(bExcl, excl_nalloc);
657 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
659 molb = &mtop->molblock[mb];
660 molt = &mtop->moltype[molb->type];
664 /* Check if the cginfo is identical for all molecules in this block.
665 * If so, we only need an array of the size of one molecule.
666 * Otherwise we make an array of #mol times #cgs per molecule.
669 for (m = 0; m < molb->nmol; m++)
671 int am = m*cgs->index[cgs->nr];
672 for (cg = 0; cg < cgs->nr; cg++)
675 a1 = cgs->index[cg+1];
676 if (getGroupType(&mtop->groups, egcENER, a_offset+am+a0) !=
677 getGroupType(&mtop->groups, egcENER, a_offset +a0))
681 if (mtop->groups.grpnr[egcQMMM] != nullptr)
683 for (ai = a0; ai < a1; ai++)
685 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
686 mtop->groups.grpnr[egcQMMM][a_offset +ai])
695 cginfo_mb[mb].cg_start = cg_offset;
696 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
697 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
698 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
699 cginfo = cginfo_mb[mb].cginfo;
701 /* Set constraints flags for constrained atoms */
702 snew(a_con, molt->atoms.nr);
703 for (ftype = 0; ftype < F_NRE; ftype++)
705 if (interaction_function[ftype].flags & IF_CONSTRAINT)
710 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
714 for (a = 0; a < nral; a++)
716 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
717 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
723 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
726 int am = m*cgs->index[cgs->nr];
727 for (cg = 0; cg < cgs->nr; cg++)
730 a1 = cgs->index[cg+1];
732 /* Store the energy group in cginfo */
733 gid = getGroupType(&mtop->groups, egcENER, a_offset+am+a0);
734 SET_CGINFO_GID(cginfo[cgm+cg], gid);
736 /* Check the intra/inter charge group exclusions */
737 if (a1-a0 > excl_nalloc)
739 excl_nalloc = a1 - a0;
740 srenew(bExcl, excl_nalloc);
742 /* bExclIntraAll: all intra cg interactions excluded
743 * bExclInter: any inter cg interactions excluded
745 bExclIntraAll = TRUE;
749 bHavePerturbedAtoms = FALSE;
750 for (ai = a0; ai < a1; ai++)
752 /* Check VDW and electrostatic interactions */
753 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
754 type_VDW[molt->atoms.atom[ai].typeB]);
755 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
756 molt->atoms.atom[ai].qB != 0);
758 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
760 /* Clear the exclusion list for atom ai */
761 for (aj = a0; aj < a1; aj++)
763 bExcl[aj-a0] = FALSE;
765 /* Loop over all the exclusions of atom ai */
766 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
769 if (aj < a0 || aj >= a1)
778 /* Check if ai excludes a0 to a1 */
779 for (aj = a0; aj < a1; aj++)
783 bExclIntraAll = FALSE;
790 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
793 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
801 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
805 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
807 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
809 /* The size in cginfo is currently only read with DD */
810 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
814 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
818 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
820 if (bHavePerturbedAtoms && fr->efep != efepNO)
822 SET_CGINFO_FEP(cginfo[cgm+cg]);
823 *bFEP_NonBonded = TRUE;
825 /* Store the charge group size */
826 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
828 if (!bExclIntraAll || bExclInter)
830 *bExcl_IntraCGAll_InterCGNone = FALSE;
837 cg_offset += molb->nmol*cgs->nr;
838 a_offset += molb->nmol*cgs->index[cgs->nr];
843 /* the solvent optimizer is called after the QM is initialized,
844 * because we don't want to have the QM subsystemto become an
848 check_solvent(fplog, mtop, fr, cginfo_mb);
850 if (getenv("GMX_NO_SOLV_OPT"))
854 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
855 "Disabling all solvent optimization\n");
857 fr->solvent_opt = esolNO;
861 fr->solvent_opt = esolNO;
863 if (!fr->solvent_opt)
865 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
867 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
869 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
877 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
882 ncg = cgi_mb[nmb-1].cg_end;
885 for (cg = 0; cg < ncg; cg++)
887 while (cg >= cgi_mb[mb].cg_end)
892 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
898 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
900 if (cginfo_mb == nullptr)
904 for (int mb = 0; mb < numMolBlocks; ++mb)
906 sfree(cginfo_mb[mb].cginfo);
911 /* Sets the sum of charges (squared) and C6 in the system in fr.
912 * Returns whether the system has a net charge.
914 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
916 /*This now calculates sum for q and c6*/
917 double qsum, q2sum, q, c6sum, c6;
922 for (const gmx_molblock_t &molb : mtop->molblock)
924 int nmol = molb.nmol;
925 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
926 for (int i = 0; i < atoms->nr; i++)
928 q = atoms->atom[i].q;
931 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
936 fr->q2sum[0] = q2sum;
937 fr->c6sum[0] = c6sum;
939 if (fr->efep != efepNO)
944 for (const gmx_molblock_t &molb : mtop->molblock)
946 int nmol = molb.nmol;
947 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
948 for (int i = 0; i < atoms->nr; i++)
950 q = atoms->atom[i].qB;
953 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
957 fr->q2sum[1] = q2sum;
958 fr->c6sum[1] = c6sum;
963 fr->qsum[1] = fr->qsum[0];
964 fr->q2sum[1] = fr->q2sum[0];
965 fr->c6sum[1] = fr->c6sum[0];
969 if (fr->efep == efepNO)
971 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
975 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
976 fr->qsum[0], fr->qsum[1]);
980 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
981 return (std::abs(fr->qsum[0]) > 1e-4 ||
982 std::abs(fr->qsum[1]) > 1e-4);
985 void update_forcerec(t_forcerec *fr, matrix box)
987 if (fr->ic->eeltype == eelGRF)
989 calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
990 fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
991 &fr->ic->k_rf, &fr->ic->c_rf);
995 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
997 const t_atoms *atoms, *atoms_tpi;
998 const t_blocka *excl;
999 int nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
1000 int64_t npair, npair_ij, tmpi, tmpj;
1001 double csix, ctwelve;
1002 int ntp, *typecount;
1005 real *nbfp_comb = nullptr;
1011 /* For LJ-PME, we want to correct for the difference between the
1012 * actual C6 values and the C6 values used by the LJ-PME based on
1013 * combination rules. */
1015 if (EVDW_PME(fr->ic->vdwtype))
1017 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1018 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1019 for (tpi = 0; tpi < ntp; ++tpi)
1021 for (tpj = 0; tpj < ntp; ++tpj)
1023 C6(nbfp_comb, ntp, tpi, tpj) =
1024 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1025 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1030 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1038 /* Count the types so we avoid natoms^2 operations */
1039 snew(typecount, ntp);
1040 gmx_mtop_count_atomtypes(mtop, q, typecount);
1042 for (tpi = 0; tpi < ntp; tpi++)
1044 for (tpj = tpi; tpj < ntp; tpj++)
1046 tmpi = typecount[tpi];
1047 tmpj = typecount[tpj];
1050 npair_ij = tmpi*tmpj;
1054 npair_ij = tmpi*(tmpi - 1)/2;
1058 /* nbfp now includes the 6.0 derivative prefactor */
1059 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1063 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1064 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1065 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1071 /* Subtract the excluded pairs.
1072 * The main reason for substracting exclusions is that in some cases
1073 * some combinations might never occur and the parameters could have
1074 * any value. These unused values should not influence the dispersion
1077 for (const gmx_molblock_t &molb : mtop->molblock)
1079 int nmol = molb.nmol;
1080 atoms = &mtop->moltype[molb.type].atoms;
1081 excl = &mtop->moltype[molb.type].excls;
1082 for (int i = 0; (i < atoms->nr); i++)
1086 tpi = atoms->atom[i].type;
1090 tpi = atoms->atom[i].typeB;
1092 j1 = excl->index[i];
1093 j2 = excl->index[i+1];
1094 for (j = j1; j < j2; j++)
1101 tpj = atoms->atom[k].type;
1105 tpj = atoms->atom[k].typeB;
1109 /* nbfp now includes the 6.0 derivative prefactor */
1110 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1114 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1115 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1116 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1126 /* Only correct for the interaction of the test particle
1127 * with the rest of the system.
1130 &mtop->moltype[mtop->molblock.back().type].atoms;
1133 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
1135 const gmx_molblock_t &molb = mtop->molblock[mb];
1136 atoms = &mtop->moltype[molb.type].atoms;
1137 for (j = 0; j < atoms->nr; j++)
1140 /* Remove the interaction of the test charge group
1143 if (mb == mtop->molblock.size() - 1)
1147 if (mb == 0 && molb.nmol == 1)
1149 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1154 tpj = atoms->atom[j].type;
1158 tpj = atoms->atom[j].typeB;
1160 for (i = 0; i < fr->n_tpi; i++)
1164 tpi = atoms_tpi->atom[i].type;
1168 tpi = atoms_tpi->atom[i].typeB;
1172 /* nbfp now includes the 6.0 derivative prefactor */
1173 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1177 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1178 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1179 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1186 if (npair - nexcl <= 0 && fplog)
1188 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1194 csix /= npair - nexcl;
1195 ctwelve /= npair - nexcl;
1199 fprintf(debug, "Counted %d exclusions\n", nexcl);
1200 fprintf(debug, "Average C6 parameter is: %10g\n", csix);
1201 fprintf(debug, "Average C12 parameter is: %10g\n", ctwelve);
1203 fr->avcsix[q] = csix;
1204 fr->avctwelve[q] = ctwelve;
1207 if (EVDW_PME(fr->ic->vdwtype))
1212 if (fplog != nullptr)
1214 if (fr->eDispCorr == edispcAllEner ||
1215 fr->eDispCorr == edispcAllEnerPres)
1217 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1218 fr->avcsix[0], fr->avctwelve[0]);
1222 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1228 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
1230 const t_atoms *at1, *at2;
1231 int i, j, tpi, tpj, ntypes;
1236 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1238 ntypes = mtop->ffparams.atnr;
1241 real bham_b_max = 0;
1242 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
1244 at1 = &mtop->moltype[mt1].atoms;
1245 for (i = 0; (i < at1->nr); i++)
1247 tpi = at1->atom[i].type;
1250 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1253 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
1255 at2 = &mtop->moltype[mt2].atoms;
1256 for (j = 0; (j < at2->nr); j++)
1258 tpj = at2->atom[j].type;
1261 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1263 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
1268 if ((b < bmin) || (bmin == -1))
1278 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1285 static void make_nbf_tables(FILE *fp,
1286 const interaction_const_t *ic, real rtab,
1287 const char *tabfn, char *eg1, char *eg2,
1293 if (tabfn == nullptr)
1297 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1302 sprintf(buf, "%s", tabfn);
1305 /* Append the two energy group names */
1306 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1307 eg1, eg2, ftp2ext(efXVG));
1309 nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1310 /* Copy the contents of the table to separate coulomb and LJ tables too,
1311 * to improve cache performance.
1313 /* For performance reasons we want
1314 * the table data to be aligned to 16-byte. The pointers could be freed
1315 * but currently aren't.
1317 snew(nbl->table_elec, 1);
1318 nbl->table_elec->interaction = GMX_TABLE_INTERACTION_ELEC;
1319 nbl->table_elec->format = nbl->table_elec_vdw->format;
1320 nbl->table_elec->r = nbl->table_elec_vdw->r;
1321 nbl->table_elec->n = nbl->table_elec_vdw->n;
1322 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1323 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1324 nbl->table_elec->ninteractions = 1;
1325 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1326 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1328 snew(nbl->table_vdw, 1);
1329 nbl->table_vdw->interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1330 nbl->table_vdw->format = nbl->table_elec_vdw->format;
1331 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1332 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1333 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1334 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1335 nbl->table_vdw->ninteractions = 2;
1336 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1337 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1339 /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1340 * with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1342 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1344 for (j = 0; j < 4; j++)
1346 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1349 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1351 for (j = 0; j < 8; j++)
1353 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1358 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1359 * ftype2 present in the topology, build an array of the number of
1360 * interactions present for each bonded interaction index found in the
1363 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1364 * valid type with that parameter.
1366 * \c count will be reallocated as necessary to fit the largest bonded
1367 * interaction index found, and its current size will be returned in
1368 * \c ncount. It will contain zero for every bonded interaction index
1369 * for which no interactions are present in the topology.
1371 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1372 int *ncount, int **count)
1375 int ftype, stride, i, j, tabnr;
1377 // Loop over all moleculetypes
1378 for (const gmx_moltype_t &molt : mtop->moltype)
1380 // Loop over all interaction types
1381 for (ftype = 0; ftype < F_NRE; ftype++)
1383 // If the current interaction type is one of the types whose tables we're trying to count...
1384 if (ftype == ftype1 || ftype == ftype2)
1386 il = &molt.ilist[ftype];
1387 stride = 1 + NRAL(ftype);
1388 // ... and there are actually some interactions for this type
1389 for (i = 0; i < il->nr; i += stride)
1391 // Find out which table index the user wanted
1392 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1395 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1397 // Make room for this index in the data structure
1398 if (tabnr >= *ncount)
1400 srenew(*count, tabnr+1);
1401 for (j = *ncount; j < tabnr+1; j++)
1407 // Record that this table index is used and must have a valid file
1415 /*!\brief If there's bonded interactions of flavour \c tabext and type
1416 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1417 * list of filenames passed to mdrun, and make bonded tables from
1420 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1421 * valid type with that parameter.
1423 * A fatal error occurs if no matching filename is found.
1425 static bondedtable_t *make_bonded_tables(FILE *fplog,
1426 int ftype1, int ftype2,
1427 const gmx_mtop_t *mtop,
1428 gmx::ArrayRef<const std::string> tabbfnm,
1438 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1440 // Are there any relevant tabulated bond interactions?
1444 for (int i = 0; i < ncount; i++)
1446 // Do any interactions exist that requires this table?
1449 // This pattern enforces the current requirement that
1450 // table filenames end in a characteristic sequence
1451 // before the file type extension, and avoids table 13
1452 // being recognized and used for table 1.
1453 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1454 bool madeTable = false;
1455 for (gmx::index j = 0; j < tabbfnm.size() && !madeTable; ++j)
1457 if (gmx::endsWith(tabbfnm[j], patternToFind))
1459 // Finally read the table from the file found
1460 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1466 bool isPlural = (ftype2 != -1);
1467 gmx_fatal(FARGS, "Tabulated interaction of type '%s%s%s' with index %d cannot be used because no table file whose name matched '%s' was passed via the gmx mdrun -tableb command-line option.",
1468 interaction_function[ftype1].longname,
1469 isPlural ? "' or '" : "",
1470 isPlural ? interaction_function[ftype2].longname : "",
1472 patternToFind.c_str());
1482 void forcerec_set_ranges(t_forcerec *fr,
1483 int ncg_home, int ncg_force,
1485 int natoms_force_constr, int natoms_f_novirsum)
1490 /* fr->ncg_force is unused in the standard code,
1491 * but it can be useful for modified code dealing with charge groups.
1493 fr->ncg_force = ncg_force;
1494 fr->natoms_force = natoms_force;
1495 fr->natoms_force_constr = natoms_force_constr;
1497 if (fr->natoms_force_constr > fr->nalloc_force)
1499 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1502 if (fr->haveDirectVirialContributions)
1504 fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1508 static real cutoff_inf(real cutoff)
1512 cutoff = GMX_CUTOFF_INF;
1518 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, const t_commrec *cr, FILE *fp)
1525 ir->rcoulomb == 0 &&
1527 ir->ePBC == epbcNONE &&
1528 ir->vdwtype == evdwCUT &&
1529 ir->coulombtype == eelCUT &&
1530 ir->efep == efepNO &&
1531 getenv("GMX_NO_ALLVSALL") == nullptr
1534 if (bAllvsAll && ir->opts.ngener > 1)
1536 const char *note = "NOTE: Can not use all-vs-all force loops, because there are multiple energy monitor groups; you might get significantly higher performance when using only a single energy monitor group.\n";
1542 fprintf(fp, "\n%s\n", note);
1548 if (bAllvsAll && fp && MASTER(cr))
1550 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1557 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
1558 const t_inputrec *ir)
1560 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1562 /* LJ PME with LB combination rule does 7 mesh operations.
1563 * This so slow that we don't compile SIMD non-bonded kernels
1565 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
1573 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1576 const gmx_hw_info_t gmx_unused &hardwareInfo)
1578 *kernel_type = nbnxnk4x4_PlainC;
1579 *ewald_excl = ewaldexclTable;
1583 #ifdef GMX_NBNXN_SIMD_4XN
1584 *kernel_type = nbnxnk4xN_SIMD_4xN;
1586 #ifdef GMX_NBNXN_SIMD_2XNN
1587 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1590 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1591 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1592 * This is based on the SIMD acceleration choice and CPU information
1593 * detected at runtime.
1595 * 4xN calculates more (zero) interactions, but has less pair-search
1596 * work and much better kernel instruction scheduling.
1598 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1599 * which doesn't have FMA, both the analytical and tabulated Ewald
1600 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1601 * 2x(4+4) because it results in significantly fewer pairs.
1602 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1603 * 10% with HT, 50% without HT. As we currently don't detect the actual
1604 * use of HT, use 4x8 to avoid a potential performance hit.
1605 * On Intel Haswell 4x8 is always faster.
1607 *kernel_type = nbnxnk4xN_SIMD_4xN;
1609 #if !GMX_SIMD_HAVE_FMA
1610 if (EEL_PME_EWALD(ir->coulombtype) ||
1611 EVDW_PME(ir->vdwtype))
1613 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1614 * There are enough instructions to make 2x(4+4) efficient.
1616 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1619 if (hardwareInfo.haveAmdZenCpu)
1621 /* One 256-bit FMA per cycle makes 2xNN faster */
1622 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1624 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1627 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
1629 #ifdef GMX_NBNXN_SIMD_4XN
1630 *kernel_type = nbnxnk4xN_SIMD_4xN;
1632 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1635 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
1637 #ifdef GMX_NBNXN_SIMD_2XNN
1638 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1640 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1644 /* Analytical Ewald exclusion correction is only an option in
1646 * Since table lookup's don't parallelize with SIMD, analytical
1647 * will probably always be faster for a SIMD width of 8 or more.
1648 * With FMA analytical is sometimes faster for a width if 4 as well.
1649 * In single precision, this is faster on Bulldozer.
1651 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1652 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
1653 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
1654 * of single or double precision and 128 or 256-bit AVX2.
1656 if (!hardwareInfo.haveAmdZenCpu)
1658 *ewald_excl = ewaldexclAnalytical;
1661 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
1663 *ewald_excl = ewaldexclTable;
1665 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
1667 *ewald_excl = ewaldexclAnalytical;
1675 const char *lookup_nbnxn_kernel_name(int kernel_type)
1677 const char *returnvalue = nullptr;
1678 switch (kernel_type)
1681 returnvalue = "not set";
1683 case nbnxnk4x4_PlainC:
1684 returnvalue = "plain C";
1686 case nbnxnk4xN_SIMD_4xN:
1687 case nbnxnk4xN_SIMD_2xNN:
1689 returnvalue = "SIMD";
1691 returnvalue = "not available";
1694 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1695 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1699 gmx_fatal(FARGS, "Illegal kernel type selected");
1704 static void pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
1705 gmx_bool use_simd_kernels,
1706 const gmx_hw_info_t &hardwareInfo,
1708 EmulateGpuNonbonded emulateGpu,
1709 const t_inputrec *ir,
1712 gmx_bool bDoNonbonded)
1714 assert(kernel_type);
1716 *kernel_type = nbnxnkNotSet;
1717 *ewald_excl = ewaldexclTable;
1719 if (emulateGpu == EmulateGpuNonbonded::Yes)
1721 *kernel_type = nbnxnk8x8x8_PlainC;
1725 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
1730 *kernel_type = nbnxnk8x8x8_GPU;
1733 if (*kernel_type == nbnxnkNotSet)
1735 if (use_simd_kernels &&
1736 nbnxn_simd_supported(mdlog, ir))
1738 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl, hardwareInfo);
1742 *kernel_type = nbnxnk4x4_PlainC;
1748 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
1749 "Using %s %dx%d nonbonded short-range kernels",
1750 lookup_nbnxn_kernel_name(*kernel_type),
1751 nbnxn_kernel_to_cluster_i_size(*kernel_type),
1752 nbnxn_kernel_to_cluster_j_size(*kernel_type));
1754 if (nbnxnk4x4_PlainC == *kernel_type ||
1755 nbnxnk8x8x8_PlainC == *kernel_type)
1757 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1758 "WARNING: Using the slow %s kernels. This should\n"
1759 "not happen during routine usage on supported platforms.",
1760 lookup_nbnxn_kernel_name(*kernel_type));
1765 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1766 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1767 bool systemHasNetCharge,
1768 interaction_const_t *ic)
1770 if (!EEL_PME_EWALD(ir->coulombtype))
1777 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1779 if (ir->coulombtype == eelP3M_AD)
1781 please_cite(fp, "Hockney1988");
1782 please_cite(fp, "Ballenegger2012");
1786 please_cite(fp, "Essmann95a");
1789 if (ir->ewald_geometry == eewg3DC)
1793 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1794 systemHasNetCharge ? " and net charge" : "");
1796 please_cite(fp, "In-Chul99a");
1797 if (systemHasNetCharge)
1799 please_cite(fp, "Ballenegger2009");
1804 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1807 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1808 1/ic->ewaldcoeff_q);
1811 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1813 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1814 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1822 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1823 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1824 interaction_const_t *ic)
1826 if (!EVDW_PME(ir->vdwtype))
1833 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1834 please_cite(fp, "Essmann95a");
1836 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1839 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1840 1/ic->ewaldcoeff_lj);
1843 if (ic->vdw_modifier == eintmodPOTSHIFT)
1845 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1846 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1850 ic->sh_lj_ewald = 0;
1854 gmx_bool uses_simple_tables(int cutoff_scheme,
1855 nonbonded_verlet_t *nbv,
1858 gmx_bool bUsesSimpleTables = TRUE;
1861 switch (cutoff_scheme)
1864 bUsesSimpleTables = TRUE;
1867 assert(nullptr != nbv);
1868 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1869 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1872 gmx_incons("unimplemented");
1874 return bUsesSimpleTables;
1877 static void init_ewald_f_table(interaction_const_t *ic,
1882 /* Get the Ewald table spacing based on Coulomb and/or LJ
1883 * Ewald coefficients and rtol.
1885 ic->tabq_scale = ewald_spline3_table_scale(ic);
1887 if (ic->cutoff_scheme == ecutsVERLET)
1889 maxr = ic->rcoulomb;
1893 maxr = std::max(ic->rcoulomb, rtab);
1895 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1897 sfree_aligned(ic->tabq_coul_FDV0);
1898 sfree_aligned(ic->tabq_coul_F);
1899 sfree_aligned(ic->tabq_coul_V);
1901 sfree_aligned(ic->tabq_vdw_FDV0);
1902 sfree_aligned(ic->tabq_vdw_F);
1903 sfree_aligned(ic->tabq_vdw_V);
1905 if (EEL_PME_EWALD(ic->eeltype))
1907 /* Create the original table data in FDV0 */
1908 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1909 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1910 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1911 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1912 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1915 if (EVDW_PME(ic->vdwtype))
1917 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1918 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1919 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1920 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1921 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1925 void init_interaction_const_tables(FILE *fp,
1926 interaction_const_t *ic,
1929 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1931 init_ewald_f_table(ic, rtab);
1935 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1936 1/ic->tabq_scale, ic->tabq_size);
1941 static void clear_force_switch_constants(shift_consts_t *sc)
1948 static void force_switch_constants(real p,
1952 /* Here we determine the coefficient for shifting the force to zero
1953 * between distance rsw and the cut-off rc.
1954 * For a potential of r^-p, we have force p*r^-(p+1).
1955 * But to save flops we absorb p in the coefficient.
1957 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1958 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1960 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1961 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1962 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1965 static void potential_switch_constants(real rsw, real rc,
1966 switch_consts_t *sc)
1968 /* The switch function is 1 at rsw and 0 at rc.
1969 * The derivative and second derivate are zero at both ends.
1970 * rsw = max(r - r_switch, 0)
1971 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1972 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1973 * force = force*dsw - potential*sw
1976 sc->c3 = -10/gmx::power3(rc - rsw);
1977 sc->c4 = 15/gmx::power4(rc - rsw);
1978 sc->c5 = -6/gmx::power5(rc - rsw);
1981 /*! \brief Construct interaction constants
1983 * This data is used (particularly) by search and force code for
1984 * short-range interactions. Many of these are constant for the whole
1985 * simulation; some are constant only after PME tuning completes.
1988 init_interaction_const(FILE *fp,
1989 interaction_const_t **interaction_const,
1990 const t_inputrec *ir,
1991 const gmx_mtop_t *mtop,
1992 bool systemHasNetCharge)
1994 interaction_const_t *ic;
1998 ic->cutoff_scheme = ir->cutoff_scheme;
2000 /* Just allocate something so we can free it */
2001 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
2002 snew_aligned(ic->tabq_coul_F, 16, 32);
2003 snew_aligned(ic->tabq_coul_V, 16, 32);
2006 ic->vdwtype = ir->vdwtype;
2007 ic->vdw_modifier = ir->vdw_modifier;
2008 ic->reppow = mtop->ffparams.reppow;
2009 ic->rvdw = cutoff_inf(ir->rvdw);
2010 ic->rvdw_switch = ir->rvdw_switch;
2011 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
2012 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
2013 if (ic->useBuckingham)
2015 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
2018 initVdwEwaldParameters(fp, ir, ic);
2020 clear_force_switch_constants(&ic->dispersion_shift);
2021 clear_force_switch_constants(&ic->repulsion_shift);
2023 switch (ic->vdw_modifier)
2025 case eintmodPOTSHIFT:
2026 /* Only shift the potential, don't touch the force */
2027 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
2028 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
2030 case eintmodFORCESWITCH:
2031 /* Switch the force, switch and shift the potential */
2032 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2033 &ic->dispersion_shift);
2034 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2035 &ic->repulsion_shift);
2037 case eintmodPOTSWITCH:
2038 /* Switch the potential and force */
2039 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2043 case eintmodEXACTCUTOFF:
2044 /* Nothing to do here */
2047 gmx_incons("unimplemented potential modifier");
2050 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2052 /* Electrostatics */
2053 ic->eeltype = ir->coulombtype;
2054 ic->coulomb_modifier = ir->coulomb_modifier;
2055 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
2056 ic->rcoulomb_switch = ir->rcoulomb_switch;
2057 ic->epsilon_r = ir->epsilon_r;
2059 /* Set the Coulomb energy conversion factor */
2060 if (ic->epsilon_r != 0)
2062 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
2066 /* eps = 0 is infinite dieletric: no Coulomb interactions */
2070 /* Reaction-field */
2071 if (EEL_RF(ic->eeltype))
2073 ic->epsilon_rf = ir->epsilon_rf;
2074 /* Generalized reaction field parameters are updated every step */
2075 if (ic->eeltype != eelGRF)
2077 calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
2078 ic->rcoulomb, 0, 0, nullptr,
2079 &ic->k_rf, &ic->c_rf);
2082 if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
2084 /* grompp should have done this, but this scheme is obsolete */
2085 ic->coulomb_modifier = eintmodEXACTCUTOFF;
2090 /* For plain cut-off we might use the reaction-field kernels */
2091 ic->epsilon_rf = ic->epsilon_r;
2093 if (ir->coulomb_modifier == eintmodPOTSHIFT)
2095 ic->c_rf = 1/ic->rcoulomb;
2103 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
2107 real dispersion_shift;
2109 dispersion_shift = ic->dispersion_shift.cpot;
2110 if (EVDW_PME(ic->vdwtype))
2112 dispersion_shift -= ic->sh_lj_ewald;
2114 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2115 ic->repulsion_shift.cpot, dispersion_shift);
2117 if (ic->eeltype == eelCUT)
2119 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2121 else if (EEL_PME(ic->eeltype))
2123 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2128 *interaction_const = ic;
2132 done_interaction_const(interaction_const_t *interaction_const)
2134 sfree_aligned(interaction_const->tabq_coul_FDV0);
2135 sfree_aligned(interaction_const->tabq_coul_F);
2136 sfree_aligned(interaction_const->tabq_coul_V);
2137 sfree(interaction_const);
2140 static void init_nb_verlet(const gmx::MDLogger &mdlog,
2141 nonbonded_verlet_t **nb_verlet,
2142 gmx_bool bFEP_NonBonded,
2143 const t_inputrec *ir,
2144 const t_forcerec *fr,
2145 const t_commrec *cr,
2146 const gmx_hw_info_t &hardwareInfo,
2147 const gmx_device_info_t *deviceInfo,
2148 const gmx_mtop_t *mtop,
2151 nonbonded_verlet_t *nbv;
2154 nbnxn_alloc_t *nb_alloc;
2155 nbnxn_free_t *nb_free;
2157 nbv = new nonbonded_verlet_t();
2159 nbv->emulateGpu = ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
2160 nbv->bUseGPU = deviceInfo != nullptr;
2162 GMX_RELEASE_ASSERT(!(nbv->emulateGpu == EmulateGpuNonbonded::Yes && nbv->bUseGPU), "When GPU emulation is active, there cannot be a GPU assignment");
2165 nbv->min_ci_balanced = 0;
2167 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2168 for (int i = 0; i < nbv->ngrp; i++)
2170 nbv->grp[i].nbl_lists.nnbl = 0;
2171 nbv->grp[i].kernel_type = nbnxnkNotSet;
2173 if (i == 0) /* local */
2175 pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
2176 nbv->bUseGPU, nbv->emulateGpu, ir,
2177 &nbv->grp[i].kernel_type,
2178 &nbv->grp[i].ewald_excl,
2181 else /* non-local */
2183 /* Use the same kernel for local and non-local interactions */
2184 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2185 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2189 nbv->listParams = gmx::compat::make_unique<NbnxnListParameters>(ir->rlist);
2190 setupDynamicPairlistPruning(mdlog, ir, mtop, box, nbv->grp[0].kernel_type, fr->ic,
2191 nbv->listParams.get());
2193 nbv->nbs = gmx::compat::make_unique<nbnxn_search>(DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
2194 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
2196 gmx_omp_nthreads_get(emntPairsearch));
2198 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2199 &nb_alloc, &nb_free);
2201 for (int i = 0; i < nbv->ngrp; i++)
2203 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2204 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2205 /* 8x8x8 "non-simple" lists are ATM always combined */
2206 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2210 int enbnxninitcombrule;
2211 if (fr->ic->vdwtype == evdwCUT &&
2212 (fr->ic->vdw_modifier == eintmodNONE ||
2213 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
2214 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
2216 /* Plain LJ cut-off: we can optimize with combination rules */
2217 enbnxninitcombrule = enbnxninitcombruleDETECT;
2219 else if (fr->ic->vdwtype == evdwPME)
2221 /* LJ-PME: we need to use a combination rule for the grid */
2222 if (fr->ljpme_combination_rule == eljpmeGEOM)
2224 enbnxninitcombrule = enbnxninitcombruleGEOM;
2228 enbnxninitcombrule = enbnxninitcombruleLB;
2233 /* We use a full combination matrix: no rule required */
2234 enbnxninitcombrule = enbnxninitcombruleNONE;
2238 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
2239 if (ir->opts.ngener - ir->nwall == 1)
2241 /* We have only one non-wall energy group, we do not need energy group
2242 * support in the non-bondeds kernels, since all non-bonded energy
2243 * contributions go to the first element of the energy group matrix.
2245 mimimumNumEnergyGroupNonbonded = 1;
2247 bool bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[0].kernel_type);
2248 nbnxn_atomdata_init(mdlog,
2250 nbv->grp[0].kernel_type,
2252 fr->ntype, fr->nbfp,
2253 mimimumNumEnergyGroupNonbonded,
2254 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2259 /* init the NxN GPU data; the last argument tells whether we'll have
2260 * both local and non-local NB calculation on GPU */
2261 nbnxn_gpu_init(&nbv->gpu_nbv,
2264 nbv->listParams.get(),
2269 if ((env = getenv("GMX_NB_MIN_CI")) != nullptr)
2273 nbv->min_ci_balanced = strtol(env, &end, 10);
2274 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
2276 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
2281 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2282 nbv->min_ci_balanced);
2287 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2290 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2291 nbv->min_ci_balanced);
2300 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2302 return nbv != nullptr && nbv->bUseGPU;
2305 void init_forcerec(FILE *fp,
2306 const gmx::MDLogger &mdlog,
2309 const t_inputrec *ir,
2310 const gmx_mtop_t *mtop,
2311 const t_commrec *cr,
2315 gmx::ArrayRef<const std::string> tabbfnm,
2316 const gmx_hw_info_t &hardwareInfo,
2317 const gmx_device_info_t *deviceInfo,
2318 gmx_bool bNoSolvOpt,
2321 int m, negp_pp, negptable, egi, egj;
2326 gmx_bool bGenericKernelOnly;
2327 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
2328 gmx_bool bFEP_NonBonded;
2329 int *nm_ind, egp_flags;
2331 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2332 fr->use_simd_kernels = TRUE;
2334 if (check_box(ir->ePBC, box))
2336 gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
2339 /* Test particle insertion ? */
2342 /* Set to the size of the molecule to be inserted (the last one) */
2343 /* Because of old style topologies, we have to use the last cg
2344 * instead of the last molecule type.
2346 cgs = &mtop->moltype[mtop->molblock.back().type].cgs;
2347 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2348 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
2349 if (fr->n_tpi != molecules.block(molecules.numBlocks() - 1).size())
2351 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2359 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
2361 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
2362 eel_names[ir->coulombtype]);
2367 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
2369 if (ir->useTwinRange)
2371 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
2373 /* Copy the user determined parameters */
2374 fr->userint1 = ir->userint1;
2375 fr->userint2 = ir->userint2;
2376 fr->userint3 = ir->userint3;
2377 fr->userint4 = ir->userint4;
2378 fr->userreal1 = ir->userreal1;
2379 fr->userreal2 = ir->userreal2;
2380 fr->userreal3 = ir->userreal3;
2381 fr->userreal4 = ir->userreal4;
2384 fr->fc_stepsize = ir->fc_stepsize;
2387 fr->efep = ir->efep;
2388 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2389 if (ir->fepvals->bScCoul)
2391 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2392 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
2396 fr->sc_alphacoul = 0;
2397 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2399 fr->sc_power = ir->fepvals->sc_power;
2400 fr->sc_r_power = ir->fepvals->sc_r_power;
2401 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
2403 env = getenv("GMX_SCSIGMA_MIN");
2407 sscanf(env, "%20lf", &dbl);
2408 fr->sc_sigma6_min = gmx::power6(dbl);
2411 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2415 fr->bNonbonded = TRUE;
2416 if (getenv("GMX_NO_NONBONDED") != nullptr)
2418 /* turn off non-bonded calculations */
2419 fr->bNonbonded = FALSE;
2420 GMX_LOG(mdlog.warning).asParagraph().appendText(
2421 "Found environment variable GMX_NO_NONBONDED.\n"
2422 "Disabling nonbonded calculations.");
2425 bGenericKernelOnly = FALSE;
2427 /* We now check in the NS code whether a particular combination of interactions
2428 * can be used with water optimization, and disable it if that is not the case.
2431 if (getenv("GMX_NB_GENERIC") != nullptr)
2436 "Found environment variable GMX_NB_GENERIC.\n"
2437 "Disabling all interaction-specific nonbonded kernels, will only\n"
2438 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2440 bGenericKernelOnly = TRUE;
2443 if (bGenericKernelOnly)
2448 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
2450 fr->use_simd_kernels = FALSE;
2454 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2455 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2456 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2460 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2462 /* Check if we can/should do all-vs-all kernels */
2463 fr->bAllvsAll = can_use_allvsall(ir, FALSE, nullptr, nullptr);
2464 fr->AllvsAll_work = nullptr;
2466 /* All-vs-all kernels have not been implemented in 4.6 and later.
2467 * See Redmine #1249. */
2470 fr->bAllvsAll = FALSE;
2474 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2475 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2476 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2477 "or try cutoff-scheme = Verlet.\n\n");
2481 /* Neighbour searching stuff */
2482 fr->cutoff_scheme = ir->cutoff_scheme;
2483 fr->bGrid = (ir->ns_type == ensGRID);
2484 fr->ePBC = ir->ePBC;
2486 if (fr->cutoff_scheme == ecutsGROUP)
2488 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2489 "removed in a future release when 'verlet' supports all interaction forms.\n";
2493 fprintf(stderr, "\n%s\n", note);
2497 fprintf(fp, "\n%s\n", note);
2501 /* Determine if we will do PBC for distances in bonded interactions */
2502 if (fr->ePBC == epbcNONE)
2504 fr->bMolPBC = FALSE;
2508 if (!DOMAINDECOMP(cr))
2512 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2513 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2514 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2516 /* The group cut-off scheme and SHAKE assume charge groups
2517 * are whole, but not using molpbc is faster in most cases.
2518 * With intermolecular interactions we need PBC for calculating
2519 * distances between atoms in different molecules.
2521 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2522 !mtop->bIntermolecularInteractions)
2524 fr->bMolPBC = ir->bPeriodicMols;
2526 if (bSHAKE && fr->bMolPBC)
2528 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2533 /* Not making molecules whole is faster in most cases,
2534 * but With orientation restraints we need whole molecules.
2536 fr->bMolPBC = (fcd->orires.nr == 0);
2538 if (getenv("GMX_USE_GRAPH") != nullptr)
2540 fr->bMolPBC = FALSE;
2543 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
2546 if (mtop->bIntermolecularInteractions)
2548 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
2552 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
2554 if (bSHAKE && fr->bMolPBC)
2556 gmx_fatal(FARGS, "SHAKE is not properly supported with intermolecular interactions. For short simulations where linked molecules remain in the same periodic image, the environment variable GMX_USE_GRAPH can be used to override this check.\n");
2562 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2566 fr->rc_scaling = ir->refcoord_scaling;
2567 copy_rvec(ir->posres_com, fr->posres_com);
2568 copy_rvec(ir->posres_comB, fr->posres_comB);
2569 fr->rlist = cutoff_inf(ir->rlist);
2570 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2572 /* This now calculates sum for q and c6*/
2573 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
2575 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2576 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
2577 init_interaction_const_tables(fp, fr->ic, ir->rlist + ir->tabext);
2579 const interaction_const_t *ic = fr->ic;
2581 /* TODO: Replace this Ewald table or move it into interaction_const_t */
2582 if (ir->coulombtype == eelEWALD)
2584 init_ewald_tab(&(fr->ewald_table), ir, fp);
2587 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2588 switch (ic->eeltype)
2591 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
2596 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2600 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2601 GMX_RELEASE_ASSERT(ic->coulomb_modifier == eintmodEXACTCUTOFF, "With the group scheme RF-zero needs the exact cut-off modifier");
2610 case eelPMEUSERSWITCH:
2611 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2617 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2621 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
2623 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
2625 /* Vdw: Translate from mdp settings to kernel format */
2626 switch (ic->vdwtype)
2631 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2635 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2639 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2645 case evdwENCADSHIFT:
2646 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2650 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
2652 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
2654 if (ir->cutoff_scheme == ecutsGROUP)
2656 fr->bvdwtab = ((ic->vdwtype != evdwCUT || !gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2657 && !EVDW_PME(ic->vdwtype));
2658 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2659 fr->bcoultab = !(ic->eeltype == eelCUT ||
2660 ic->eeltype == eelEWALD ||
2661 ic->eeltype == eelPME ||
2662 ic->eeltype == eelP3M_AD ||
2663 ic->eeltype == eelRF ||
2664 ic->eeltype == eelRF_ZERO);
2666 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2667 * going to be faster to tabulate the interaction than calling the generic kernel.
2668 * However, if generic kernels have been requested we keep things analytically.
2670 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2671 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2672 !bGenericKernelOnly)
2674 if ((ic->rcoulomb_switch != ic->rvdw_switch) || (ic->rcoulomb != ic->rvdw))
2676 fr->bcoultab = TRUE;
2677 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2678 * which would otherwise need two tables.
2682 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2683 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2684 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2685 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2687 if ((ic->rcoulomb != ic->rvdw) && (!bGenericKernelOnly))
2689 fr->bcoultab = TRUE;
2693 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2695 fr->bcoultab = TRUE;
2697 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2702 if (getenv("GMX_REQUIRE_TABLES"))
2705 fr->bcoultab = TRUE;
2710 fprintf(fp, "Table routines are used for coulomb: %s\n",
2711 gmx::boolToString(fr->bcoultab));
2712 fprintf(fp, "Table routines are used for vdw: %s\n",
2713 gmx::boolToString(fr->bvdwtab));
2718 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2719 fr->nbkernel_vdw_modifier = eintmodNONE;
2723 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2724 fr->nbkernel_elec_modifier = eintmodNONE;
2728 if (ir->cutoff_scheme == ecutsVERLET)
2730 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2732 gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2734 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2735 * while mdrun does not (and never did) support this.
2737 if (EEL_USER(fr->ic->eeltype))
2739 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2740 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2743 fr->bvdwtab = FALSE;
2744 fr->bcoultab = FALSE;
2747 /* 1-4 interaction electrostatics */
2748 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2750 /* Parameters for generalized RF */
2754 if (ic->eeltype == eelGRF)
2756 init_generalized_rf(fp, mtop, ir, fr);
2759 fr->haveDirectVirialContributions =
2760 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
2761 fr->forceProviders->hasForceProvider() ||
2762 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2763 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2769 if (fr->haveDirectVirialContributions)
2771 fr->forceBufferForDirectVirialContributions = new std::vector<gmx::RVec>;
2774 if (fr->cutoff_scheme == ecutsGROUP &&
2775 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2777 /* Count the total number of charge groups */
2778 fr->cg_nalloc = ncg_mtop(mtop);
2779 srenew(fr->cg_cm, fr->cg_nalloc);
2781 if (fr->shift_vec == nullptr)
2783 snew(fr->shift_vec, SHIFTS);
2786 if (fr->fshift == nullptr)
2788 snew(fr->fshift, SHIFTS);
2791 if (fr->nbfp == nullptr)
2793 fr->ntype = mtop->ffparams.atnr;
2794 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2795 if (EVDW_PME(ic->vdwtype))
2797 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2801 /* Copy the energy group exclusions */
2802 fr->egp_flags = ir->opts.egp_flags;
2804 /* Van der Waals stuff */
2805 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
2807 if (ic->rvdw_switch >= ic->rvdw)
2809 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2810 ic->rvdw_switch, ic->rvdw);
2814 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2815 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
2816 ic->rvdw_switch, ic->rvdw);
2820 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
2822 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2825 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
2827 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2830 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2832 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2835 if (fp && fr->cutoff_scheme == ecutsGROUP)
2837 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2838 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
2841 fr->eDispCorr = ir->eDispCorr;
2842 fr->numAtomsForDispersionCorrection = mtop->natoms;
2843 if (ir->eDispCorr != edispcNO)
2845 set_avcsixtwelve(fp, fr, mtop);
2848 if (ir->implicit_solvent)
2850 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
2853 /* Construct tables for the group scheme. A little unnecessary to
2854 * make both vdw and coul tables sometimes, but what the
2855 * heck. Note that both cutoff schemes construct Ewald tables in
2856 * init_interaction_const_tables. */
2857 needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
2858 (fr->bcoultab || fr->bvdwtab));
2860 negp_pp = ir->opts.ngener - ir->nwall;
2862 if (!needGroupSchemeTables)
2864 bSomeNormalNbListsAreInUse = TRUE;
2869 bSomeNormalNbListsAreInUse = FALSE;
2870 for (egi = 0; egi < negp_pp; egi++)
2872 for (egj = egi; egj < negp_pp; egj++)
2874 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2875 if (!(egp_flags & EGP_EXCL))
2877 if (egp_flags & EGP_TABLE)
2883 bSomeNormalNbListsAreInUse = TRUE;
2888 if (bSomeNormalNbListsAreInUse)
2890 fr->nnblists = negptable + 1;
2894 fr->nnblists = negptable;
2896 if (fr->nnblists > 1)
2898 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
2902 snew(fr->nblists, fr->nnblists);
2904 /* This code automatically gives table length tabext without cut-off's,
2905 * in that case grompp should already have checked that we do not need
2906 * normal tables and we only generate tables for 1-4 interactions.
2908 rtab = ir->rlist + ir->tabext;
2910 if (needGroupSchemeTables)
2912 /* make tables for ordinary interactions */
2913 if (bSomeNormalNbListsAreInUse)
2915 make_nbf_tables(fp, ic, rtab, tabfn, nullptr, nullptr, &fr->nblists[0]);
2924 /* Read the special tables for certain energy group pairs */
2925 nm_ind = mtop->groups.grps[egcENER].nm_ind;
2926 for (egi = 0; egi < negp_pp; egi++)
2928 for (egj = egi; egj < negp_pp; egj++)
2930 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2931 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2933 if (fr->nnblists > 1)
2935 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
2937 /* Read the table file with the two energy groups names appended */
2938 make_nbf_tables(fp, ic, rtab, tabfn,
2939 *mtop->groups.grpname[nm_ind[egi]],
2940 *mtop->groups.grpname[nm_ind[egj]],
2944 else if (fr->nnblists > 1)
2946 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
2953 /* Tables might not be used for the potential modifier
2954 * interactions per se, but we still need them to evaluate
2955 * switch/shift dispersion corrections in this case. */
2956 if (fr->eDispCorr != edispcNO)
2958 fr->dispersionCorrectionTable = makeDispersionCorrectionTable(fp, ic, rtab, tabfn);
2961 /* We want to use unmodified tables for 1-4 coulombic
2962 * interactions, so we must in general have an extra set of
2964 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2965 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2966 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
2968 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
2969 GMX_MAKETABLES_14ONLY);
2973 fr->nwall = ir->nwall;
2974 if (ir->nwall && ir->wall_type == ewtTABLE)
2976 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
2979 if (fcd && !tabbfnm.empty())
2981 // Need to catch std::bad_alloc
2982 // TODO Don't need to catch this here, when merging with master branch
2985 fcd->bondtab = make_bonded_tables(fp,
2986 F_TABBONDS, F_TABBONDSNC,
2987 mtop, tabbfnm, "b");
2988 fcd->angletab = make_bonded_tables(fp,
2990 mtop, tabbfnm, "a");
2991 fcd->dihtab = make_bonded_tables(fp,
2993 mtop, tabbfnm, "d");
2995 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
3001 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3005 /* QM/MM initialization if requested
3009 fprintf(stderr, "QM/MM calculation requested.\n");
3012 fr->bQMMM = ir->bQMMM;
3013 fr->qr = mk_QMMMrec();
3015 /* Set all the static charge group info */
3016 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3018 &fr->bExcl_IntraCGAll_InterCGNone);
3019 if (DOMAINDECOMP(cr))
3021 fr->cginfo = nullptr;
3025 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
3028 if (!DOMAINDECOMP(cr))
3030 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3031 mtop->natoms, mtop->natoms, mtop->natoms);
3034 fr->print_force = print_force;
3037 /* coarse load balancing vars */
3042 /* Initialize neighbor search */
3044 init_ns(fp, cr, fr->ns, fr, mtop);
3046 if (thisRankHasDuty(cr, DUTY_PP))
3048 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3051 /* Initialize the thread working data for bonded interactions */
3052 init_bonded_threading(fp, mtop->groups.grps[egcENER].nr,
3053 &fr->bondedThreading);
3055 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
3056 snew(fr->ewc_t, fr->nthread_ewc);
3058 if (fr->cutoff_scheme == ecutsVERLET)
3060 // We checked the cut-offs in grompp, but double-check here.
3061 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
3062 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
3064 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
3068 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
3071 init_nb_verlet(mdlog, &fr->nbv, bFEP_NonBonded, ir, fr,
3072 cr, hardwareInfo, deviceInfo,
3078 /* Here we switch from using mdlog, which prints the newline before
3079 * the paragraph, to our old fprintf logging, which prints the newline
3080 * after the paragraph, so we should add a newline here.
3085 if (ir->eDispCorr != edispcNO)
3087 calc_enervirdiff(fp, ir->eDispCorr, fr);
3091 /* Frees GPU memory and sets a tMPI node barrier.
3093 * Note that this function needs to be called even if GPUs are not used
3094 * in this run because the PME ranks have no knowledge of whether GPUs
3095 * are used or not, but all ranks need to enter the barrier below.
3096 * \todo Remove physical node barrier from this function after making sure
3097 * that it's not needed anymore (with a shared GPU run).
3099 void free_gpu_resources(const t_forcerec *fr,
3100 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator)
3102 bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->bUseGPU;
3104 /* stop the GPU profiler (only CUDA) */
3107 if (isPPrankUsingGPU)
3109 /* free nbnxn data in GPU memory */
3110 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3113 /* With tMPI we need to wait for all ranks to finish deallocation before
3114 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
3117 * This is not a concern in OpenCL where we use one context per rank which
3118 * is freed in nbnxn_gpu_free().
3120 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
3121 * but it is easier and more futureproof to call it on the whole node.
3125 physicalNodeCommunicator.barrier();
3129 void done_forcerec(t_forcerec *fr, int numMolBlocks, int numEnergyGroups)
3133 // PME-only ranks don't have a forcerec
3136 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
3138 done_interaction_const(fr->ic);
3139 sfree(fr->shift_vec);
3142 done_ns(fr->ns, numEnergyGroups);
3144 tear_down_bonded_threading(fr->bondedThreading);
3145 fr->bondedThreading = nullptr;