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.
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/compat/make_unique.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/ewald/ewald.h"
56 #include "gromacs/ewald/ewald-utils.h"
57 #include "gromacs/fileio/filetypes.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/hardware/hw_info.h"
62 #include "gromacs/listed-forces/manage-threading.h"
63 #include "gromacs/listed-forces/pairs.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/force.h"
69 #include "gromacs/mdlib/forcerec-threading.h"
70 #include "gromacs/mdlib/gmx_omp_nthreads.h"
71 #include "gromacs/mdlib/md_support.h"
72 #include "gromacs/mdlib/nb_verlet.h"
73 #include "gromacs/mdlib/nbnxn_atomdata.h"
74 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
75 #include "gromacs/mdlib/nbnxn_internal.h"
76 #include "gromacs/mdlib/nbnxn_search.h"
77 #include "gromacs/mdlib/nbnxn_simd.h"
78 #include "gromacs/mdlib/nbnxn_tuning.h"
79 #include "gromacs/mdlib/nbnxn_util.h"
80 #include "gromacs/mdlib/ns.h"
81 #include "gromacs/mdlib/qmmm.h"
82 #include "gromacs/mdlib/rf_util.h"
83 #include "gromacs/mdlib/sim_util.h"
84 #include "gromacs/mdlib/wall.h"
85 #include "gromacs/mdtypes/commrec.h"
86 #include "gromacs/mdtypes/fcdata.h"
87 #include "gromacs/mdtypes/group.h"
88 #include "gromacs/mdtypes/iforceprovider.h"
89 #include "gromacs/mdtypes/inputrec.h"
90 #include "gromacs/mdtypes/md_enums.h"
91 #include "gromacs/pbcutil/ishift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/simd/simd.h"
94 #include "gromacs/tables/forcetable.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/trajectory/trajectoryframe.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/exceptions.h"
99 #include "gromacs/utility/fatalerror.h"
100 #include "gromacs/utility/gmxassert.h"
101 #include "gromacs/utility/logger.h"
102 #include "gromacs/utility/physicalnodecommunicator.h"
103 #include "gromacs/utility/pleasecite.h"
104 #include "gromacs/utility/smalloc.h"
105 #include "gromacs/utility/strconvert.h"
107 #include "nbnxn_gpu_jit_support.h"
109 const char *egrp_nm[egNR+1] = {
110 "Coul-SR", "LJ-SR", "Buck-SR",
111 "Coul-14", "LJ-14", nullptr
114 t_forcerec *mk_forcerec(void)
123 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
131 snew(nbfp, 3*atnr*atnr);
132 for (i = k = 0; (i < atnr); i++)
134 for (j = 0; (j < atnr); j++, k++)
136 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
137 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
138 /* nbfp now includes the 6.0 derivative prefactor */
139 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
145 snew(nbfp, 2*atnr*atnr);
146 for (i = k = 0; (i < atnr); i++)
148 for (j = 0; (j < atnr); j++, k++)
150 /* nbfp now includes the 6.0/12.0 derivative prefactors */
151 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
152 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
160 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
163 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
166 /* For LJ-PME simulations, we correct the energies with the reciprocal space
167 * inside of the cut-off. To do this the non-bonded kernels needs to have
168 * access to the C6-values used on the reciprocal grid in pme.c
172 snew(grid, 2*atnr*atnr);
173 for (i = k = 0; (i < atnr); i++)
175 for (j = 0; (j < atnr); j++, k++)
177 c6i = idef->iparams[i*(atnr+1)].lj.c6;
178 c12i = idef->iparams[i*(atnr+1)].lj.c12;
179 c6j = idef->iparams[j*(atnr+1)].lj.c6;
180 c12j = idef->iparams[j*(atnr+1)].lj.c12;
181 c6 = std::sqrt(c6i * c6j);
182 if (fr->ljpme_combination_rule == eljpmeLB
183 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
185 sigmai = gmx::sixthroot(c12i / c6i);
186 sigmaj = gmx::sixthroot(c12j / c6j);
187 epsi = c6i * c6i / c12i;
188 epsj = c6j * c6j / c12j;
189 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
191 /* Store the elements at the same relative positions as C6 in nbfp in order
192 * to simplify access in the kernels
194 grid[2*(atnr*i+j)] = c6*6.0;
200 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
204 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
208 snew(nbfp, 2*atnr*atnr);
209 for (i = 0; i < atnr; ++i)
211 for (j = 0; j < atnr; ++j)
213 c6i = idef->iparams[i*(atnr+1)].lj.c6;
214 c12i = idef->iparams[i*(atnr+1)].lj.c12;
215 c6j = idef->iparams[j*(atnr+1)].lj.c6;
216 c12j = idef->iparams[j*(atnr+1)].lj.c12;
217 c6 = std::sqrt(c6i * c6j);
218 c12 = std::sqrt(c12i * c12j);
219 if (comb_rule == eCOMB_ARITHMETIC
220 && !gmx_numzero(c6) && !gmx_numzero(c12))
222 sigmai = gmx::sixthroot(c12i / c6i);
223 sigmaj = gmx::sixthroot(c12j / c6j);
224 epsi = c6i * c6i / c12i;
225 epsj = c6j * c6j / c12j;
226 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
227 c12 = std::sqrt(epsi * epsj) * gmx::power12(0.5*(sigmai+sigmaj));
229 C6(nbfp, atnr, i, j) = c6*6.0;
230 C12(nbfp, atnr, i, j) = c12*12.0;
236 /* This routine sets fr->solvent_opt to the most common solvent in the
237 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
238 * the fr->solvent_type array with the correct type (or esolNO).
240 * Charge groups that fulfill the conditions but are not identical to the
241 * most common one will be marked as esolNO in the solvent_type array.
243 * TIP3p is identical to SPC for these purposes, so we call it
244 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
246 * NOTE: QM particle should not
247 * become an optimized solvent. Not even if there is only one charge
257 } solvent_parameters_t;
260 check_solvent_cg(const gmx_moltype_t *molt,
263 const unsigned char *qm_grpnr,
264 const t_grps *qm_grps,
266 int *n_solvent_parameters,
267 solvent_parameters_t **solvent_parameters_p,
277 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
278 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
281 solvent_parameters_t *solvent_parameters;
283 /* We use a list with parameters for each solvent type.
284 * Every time we discover a new molecule that fulfills the basic
285 * conditions for a solvent we compare with the previous entries
286 * in these lists. If the parameters are the same we just increment
287 * the counter for that type, and otherwise we create a new type
288 * based on the current molecule.
290 * Once we've finished going through all molecules we check which
291 * solvent is most common, and mark all those molecules while we
292 * clear the flag on all others.
295 solvent_parameters = *solvent_parameters_p;
297 /* Mark the cg first as non optimized */
300 /* Check if this cg has no exclusions with atoms in other charge groups
301 * and all atoms inside the charge group excluded.
302 * We only have 3 or 4 atom solvent loops.
304 if (GET_CGINFO_EXCL_INTER(cginfo) ||
305 !GET_CGINFO_EXCL_INTRA(cginfo))
310 /* Get the indices of the first atom in this charge group */
311 j0 = molt->cgs.index[cg0];
312 j1 = molt->cgs.index[cg0+1];
314 /* Number of atoms in our molecule */
320 "Moltype '%s': there are %d atoms in this charge group\n",
324 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
327 if (nj < 3 || nj > 4)
332 /* Check if we are doing QM on this group */
334 if (qm_grpnr != nullptr)
336 for (j = j0; j < j1 && !qm; j++)
338 qm = (qm_grpnr[j] < qm_grps->nr - 1);
341 /* Cannot use solvent optimization with QM */
347 atom = molt->atoms.atom;
349 /* Still looks like a solvent, time to check parameters */
351 /* If it is perturbed (free energy) we can't use the solvent loops,
352 * so then we just skip to the next molecule.
356 for (j = j0; j < j1 && !perturbed; j++)
358 perturbed = PERTURBED(atom[j]);
366 /* Now it's only a question if the VdW and charge parameters
367 * are OK. Before doing the check we compare and see if they are
368 * identical to a possible previous solvent type.
369 * First we assign the current types and charges.
371 for (j = 0; j < nj; j++)
373 tmp_vdwtype[j] = atom[j0+j].type;
374 tmp_charge[j] = atom[j0+j].q;
377 /* Does it match any previous solvent type? */
378 for (k = 0; k < *n_solvent_parameters; k++)
383 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
384 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
385 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
390 /* Check that types & charges match for all atoms in molecule */
391 for (j = 0; j < nj && match == TRUE; j++)
393 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
397 if (tmp_charge[j] != solvent_parameters[k].charge[j])
404 /* Congratulations! We have a matched solvent.
405 * Flag it with this type for later processing.
408 solvent_parameters[k].count += nmol;
410 /* We are done with this charge group */
415 /* If we get here, we have a tentative new solvent type.
416 * Before we add it we must check that it fulfills the requirements
417 * of the solvent optimized loops. First determine which atoms have
420 for (j = 0; j < nj; j++)
423 tjA = tmp_vdwtype[j];
425 /* Go through all other tpes and see if any have non-zero
426 * VdW parameters when combined with this one.
428 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
430 /* We already checked that the atoms weren't perturbed,
431 * so we only need to check state A now.
435 has_vdw[j] = (has_vdw[j] ||
436 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
437 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
438 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
443 has_vdw[j] = (has_vdw[j] ||
444 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
445 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
450 /* Now we know all we need to make the final check and assignment. */
454 * For this we require thatn all atoms have charge,
455 * the charges on atom 2 & 3 should be the same, and only
456 * atom 1 might have VdW.
458 if (has_vdw[1] == FALSE &&
459 has_vdw[2] == FALSE &&
460 tmp_charge[0] != 0 &&
461 tmp_charge[1] != 0 &&
462 tmp_charge[2] == tmp_charge[1])
464 srenew(solvent_parameters, *n_solvent_parameters+1);
465 solvent_parameters[*n_solvent_parameters].model = esolSPC;
466 solvent_parameters[*n_solvent_parameters].count = nmol;
467 for (k = 0; k < 3; k++)
469 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
470 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
473 *cg_sp = *n_solvent_parameters;
474 (*n_solvent_parameters)++;
479 /* Or could it be a TIP4P?
480 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
481 * Only atom 1 mght have VdW.
483 if (has_vdw[1] == FALSE &&
484 has_vdw[2] == FALSE &&
485 has_vdw[3] == FALSE &&
486 tmp_charge[0] == 0 &&
487 tmp_charge[1] != 0 &&
488 tmp_charge[2] == tmp_charge[1] &&
491 srenew(solvent_parameters, *n_solvent_parameters+1);
492 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
493 solvent_parameters[*n_solvent_parameters].count = nmol;
494 for (k = 0; k < 4; k++)
496 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
497 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
500 *cg_sp = *n_solvent_parameters;
501 (*n_solvent_parameters)++;
505 *solvent_parameters_p = solvent_parameters;
509 check_solvent(FILE * fp,
510 const gmx_mtop_t * mtop,
512 cginfo_mb_t *cginfo_mb)
515 const gmx_moltype_t *molt;
516 int mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
517 int n_solvent_parameters;
518 solvent_parameters_t *solvent_parameters;
524 fprintf(debug, "Going to determine what solvent types we have.\n");
527 n_solvent_parameters = 0;
528 solvent_parameters = nullptr;
529 /* Allocate temporary array for solvent type */
530 snew(cg_sp, mtop->molblock.size());
533 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
535 molt = &mtop->moltype[mtop->molblock[mb].type];
537 /* Here we have to loop over all individual molecules
538 * because we need to check for QMMM particles.
540 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
541 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
542 nmol = mtop->molblock[mb].nmol/nmol_ch;
543 for (mol = 0; mol < nmol_ch; mol++)
546 am = mol*cgs->index[cgs->nr];
547 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
549 check_solvent_cg(molt, cg_mol, nmol,
550 mtop->groups.grpnr[egcQMMM] ?
551 mtop->groups.grpnr[egcQMMM]+at_offset+am : nullptr,
552 &mtop->groups.grps[egcQMMM],
554 &n_solvent_parameters, &solvent_parameters,
555 cginfo_mb[mb].cginfo[cgm+cg_mol],
556 &cg_sp[mb][cgm+cg_mol]);
559 at_offset += cgs->index[cgs->nr];
562 /* Puh! We finished going through all charge groups.
563 * Now find the most common solvent model.
566 /* Most common solvent this far */
568 for (i = 0; i < n_solvent_parameters; i++)
571 solvent_parameters[i].count > solvent_parameters[bestsp].count)
579 bestsol = solvent_parameters[bestsp].model;
587 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
589 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
590 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
591 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
593 if (cg_sp[mb][i] == bestsp)
595 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
600 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
607 if (bestsol != esolNO && fp != nullptr)
609 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
611 solvent_parameters[bestsp].count);
614 sfree(solvent_parameters);
615 fr->solvent_opt = bestsol;
619 acNONE = 0, acCONSTRAINT, acSETTLE
622 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
623 t_forcerec *fr, gmx_bool bNoSolvOpt,
624 gmx_bool *bFEP_NonBonded,
625 gmx_bool *bExcl_IntraCGAll_InterCGNone)
628 const t_blocka *excl;
629 const gmx_moltype_t *molt;
630 const gmx_molblock_t *molb;
631 cginfo_mb_t *cginfo_mb;
634 int cg_offset, a_offset;
635 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
639 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
641 snew(cginfo_mb, mtop->molblock.size());
643 snew(type_VDW, fr->ntype);
644 for (ai = 0; ai < fr->ntype; ai++)
646 type_VDW[ai] = FALSE;
647 for (j = 0; j < fr->ntype; j++)
649 type_VDW[ai] = type_VDW[ai] ||
651 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
652 C12(fr->nbfp, fr->ntype, ai, j) != 0;
656 *bFEP_NonBonded = FALSE;
657 *bExcl_IntraCGAll_InterCGNone = TRUE;
660 snew(bExcl, excl_nalloc);
663 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
665 molb = &mtop->molblock[mb];
666 molt = &mtop->moltype[molb->type];
670 /* Check if the cginfo is identical for all molecules in this block.
671 * If so, we only need an array of the size of one molecule.
672 * Otherwise we make an array of #mol times #cgs per molecule.
675 for (m = 0; m < molb->nmol; m++)
677 int am = m*cgs->index[cgs->nr];
678 for (cg = 0; cg < cgs->nr; cg++)
681 a1 = cgs->index[cg+1];
682 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
683 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
687 if (mtop->groups.grpnr[egcQMMM] != nullptr)
689 for (ai = a0; ai < a1; ai++)
691 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
692 mtop->groups.grpnr[egcQMMM][a_offset +ai])
701 cginfo_mb[mb].cg_start = cg_offset;
702 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
703 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
704 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
705 cginfo = cginfo_mb[mb].cginfo;
707 /* Set constraints flags for constrained atoms */
708 snew(a_con, molt->atoms.nr);
709 for (ftype = 0; ftype < F_NRE; ftype++)
711 if (interaction_function[ftype].flags & IF_CONSTRAINT)
716 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
720 for (a = 0; a < nral; a++)
722 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
723 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
729 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
732 int am = m*cgs->index[cgs->nr];
733 for (cg = 0; cg < cgs->nr; cg++)
736 a1 = cgs->index[cg+1];
738 /* Store the energy group in cginfo */
739 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
740 SET_CGINFO_GID(cginfo[cgm+cg], gid);
742 /* Check the intra/inter charge group exclusions */
743 if (a1-a0 > excl_nalloc)
745 excl_nalloc = a1 - a0;
746 srenew(bExcl, excl_nalloc);
748 /* bExclIntraAll: all intra cg interactions excluded
749 * bExclInter: any inter cg interactions excluded
751 bExclIntraAll = TRUE;
755 bHavePerturbedAtoms = FALSE;
756 for (ai = a0; ai < a1; ai++)
758 /* Check VDW and electrostatic interactions */
759 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
760 type_VDW[molt->atoms.atom[ai].typeB]);
761 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
762 molt->atoms.atom[ai].qB != 0);
764 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
766 /* Clear the exclusion list for atom ai */
767 for (aj = a0; aj < a1; aj++)
769 bExcl[aj-a0] = FALSE;
771 /* Loop over all the exclusions of atom ai */
772 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
775 if (aj < a0 || aj >= a1)
784 /* Check if ai excludes a0 to a1 */
785 for (aj = a0; aj < a1; aj++)
789 bExclIntraAll = FALSE;
796 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
799 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
807 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
811 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
813 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
815 /* The size in cginfo is currently only read with DD */
816 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
820 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
824 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
826 if (bHavePerturbedAtoms && fr->efep != efepNO)
828 SET_CGINFO_FEP(cginfo[cgm+cg]);
829 *bFEP_NonBonded = TRUE;
831 /* Store the charge group size */
832 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
834 if (!bExclIntraAll || bExclInter)
836 *bExcl_IntraCGAll_InterCGNone = FALSE;
843 cg_offset += molb->nmol*cgs->nr;
844 a_offset += molb->nmol*cgs->index[cgs->nr];
849 /* the solvent optimizer is called after the QM is initialized,
850 * because we don't want to have the QM subsystemto become an
854 check_solvent(fplog, mtop, fr, cginfo_mb);
856 if (getenv("GMX_NO_SOLV_OPT"))
860 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
861 "Disabling all solvent optimization\n");
863 fr->solvent_opt = esolNO;
867 fr->solvent_opt = esolNO;
869 if (!fr->solvent_opt)
871 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
873 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
875 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
883 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
888 ncg = cgi_mb[nmb-1].cg_end;
891 for (cg = 0; cg < ncg; cg++)
893 while (cg >= cgi_mb[mb].cg_end)
898 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
904 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
906 if (cginfo_mb == nullptr)
910 for (int mb = 0; mb < numMolBlocks; ++mb)
912 sfree(cginfo_mb[mb].cginfo);
917 /* Sets the sum of charges (squared) and C6 in the system in fr.
918 * Returns whether the system has a net charge.
920 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
922 /*This now calculates sum for q and c6*/
923 double qsum, q2sum, q, c6sum, c6;
928 for (const gmx_molblock_t &molb : mtop->molblock)
930 int nmol = molb.nmol;
931 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
932 for (int i = 0; i < atoms->nr; i++)
934 q = atoms->atom[i].q;
937 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
942 fr->q2sum[0] = q2sum;
943 fr->c6sum[0] = c6sum;
945 if (fr->efep != efepNO)
950 for (const gmx_molblock_t &molb : mtop->molblock)
952 int nmol = molb.nmol;
953 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
954 for (int i = 0; i < atoms->nr; i++)
956 q = atoms->atom[i].qB;
959 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
963 fr->q2sum[1] = q2sum;
964 fr->c6sum[1] = c6sum;
969 fr->qsum[1] = fr->qsum[0];
970 fr->q2sum[1] = fr->q2sum[0];
971 fr->c6sum[1] = fr->c6sum[0];
975 if (fr->efep == efepNO)
977 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
981 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
982 fr->qsum[0], fr->qsum[1]);
986 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
987 return (std::abs(fr->qsum[0]) > 1e-4 ||
988 std::abs(fr->qsum[1]) > 1e-4);
991 void update_forcerec(t_forcerec *fr, matrix box)
993 if (fr->ic->eeltype == eelGRF)
995 calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
996 fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
997 &fr->ic->k_rf, &fr->ic->c_rf);
1001 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
1003 const t_atoms *atoms, *atoms_tpi;
1004 const t_blocka *excl;
1005 int nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
1006 gmx_int64_t npair, npair_ij, tmpi, tmpj;
1007 double csix, ctwelve;
1008 int ntp, *typecount;
1011 real *nbfp_comb = nullptr;
1017 /* For LJ-PME, we want to correct for the difference between the
1018 * actual C6 values and the C6 values used by the LJ-PME based on
1019 * combination rules. */
1021 if (EVDW_PME(fr->ic->vdwtype))
1023 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1024 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1025 for (tpi = 0; tpi < ntp; ++tpi)
1027 for (tpj = 0; tpj < ntp; ++tpj)
1029 C6(nbfp_comb, ntp, tpi, tpj) =
1030 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1031 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1036 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1044 /* Count the types so we avoid natoms^2 operations */
1045 snew(typecount, ntp);
1046 gmx_mtop_count_atomtypes(mtop, q, typecount);
1048 for (tpi = 0; tpi < ntp; tpi++)
1050 for (tpj = tpi; tpj < ntp; tpj++)
1052 tmpi = typecount[tpi];
1053 tmpj = typecount[tpj];
1056 npair_ij = tmpi*tmpj;
1060 npair_ij = tmpi*(tmpi - 1)/2;
1064 /* nbfp now includes the 6.0 derivative prefactor */
1065 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1069 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1070 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1071 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1077 /* Subtract the excluded pairs.
1078 * The main reason for substracting exclusions is that in some cases
1079 * some combinations might never occur and the parameters could have
1080 * any value. These unused values should not influence the dispersion
1083 for (const gmx_molblock_t &molb : mtop->molblock)
1085 int nmol = molb.nmol;
1086 atoms = &mtop->moltype[molb.type].atoms;
1087 excl = &mtop->moltype[molb.type].excls;
1088 for (int i = 0; (i < atoms->nr); i++)
1092 tpi = atoms->atom[i].type;
1096 tpi = atoms->atom[i].typeB;
1098 j1 = excl->index[i];
1099 j2 = excl->index[i+1];
1100 for (j = j1; j < j2; j++)
1107 tpj = atoms->atom[k].type;
1111 tpj = atoms->atom[k].typeB;
1115 /* nbfp now includes the 6.0 derivative prefactor */
1116 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1120 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1121 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1122 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1132 /* Only correct for the interaction of the test particle
1133 * with the rest of the system.
1136 &mtop->moltype[mtop->molblock.back().type].atoms;
1139 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
1141 const gmx_molblock_t &molb = mtop->molblock[mb];
1142 atoms = &mtop->moltype[molb.type].atoms;
1143 for (j = 0; j < atoms->nr; j++)
1146 /* Remove the interaction of the test charge group
1149 if (mb == mtop->molblock.size() - 1)
1153 if (mb == 0 && molb.nmol == 1)
1155 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1160 tpj = atoms->atom[j].type;
1164 tpj = atoms->atom[j].typeB;
1166 for (i = 0; i < fr->n_tpi; i++)
1170 tpi = atoms_tpi->atom[i].type;
1174 tpi = atoms_tpi->atom[i].typeB;
1178 /* nbfp now includes the 6.0 derivative prefactor */
1179 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1183 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1184 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1185 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1192 if (npair - nexcl <= 0 && fplog)
1194 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1200 csix /= npair - nexcl;
1201 ctwelve /= npair - nexcl;
1205 fprintf(debug, "Counted %d exclusions\n", nexcl);
1206 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1207 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1209 fr->avcsix[q] = csix;
1210 fr->avctwelve[q] = ctwelve;
1213 if (EVDW_PME(fr->ic->vdwtype))
1218 if (fplog != nullptr)
1220 if (fr->eDispCorr == edispcAllEner ||
1221 fr->eDispCorr == edispcAllEnerPres)
1223 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1224 fr->avcsix[0], fr->avctwelve[0]);
1228 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1234 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
1236 const t_atoms *at1, *at2;
1237 int i, j, tpi, tpj, ntypes;
1242 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1244 ntypes = mtop->ffparams.atnr;
1247 real bham_b_max = 0;
1248 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
1250 at1 = &mtop->moltype[mt1].atoms;
1251 for (i = 0; (i < at1->nr); i++)
1253 tpi = at1->atom[i].type;
1256 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1259 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
1261 at2 = &mtop->moltype[mt2].atoms;
1262 for (j = 0; (j < at2->nr); j++)
1264 tpj = at2->atom[j].type;
1267 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1269 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
1274 if ((b < bmin) || (bmin == -1))
1284 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1291 static void make_nbf_tables(FILE *fp,
1292 const interaction_const_t *ic, real rtab,
1293 const char *tabfn, char *eg1, char *eg2,
1299 if (tabfn == nullptr)
1303 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1308 sprintf(buf, "%s", tabfn);
1311 /* Append the two energy group names */
1312 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1313 eg1, eg2, ftp2ext(efXVG));
1315 nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1316 /* Copy the contents of the table to separate coulomb and LJ tables too,
1317 * to improve cache performance.
1319 /* For performance reasons we want
1320 * the table data to be aligned to 16-byte. The pointers could be freed
1321 * but currently aren't.
1323 snew(nbl->table_elec, 1);
1324 nbl->table_elec->interaction = GMX_TABLE_INTERACTION_ELEC;
1325 nbl->table_elec->format = nbl->table_elec_vdw->format;
1326 nbl->table_elec->r = nbl->table_elec_vdw->r;
1327 nbl->table_elec->n = nbl->table_elec_vdw->n;
1328 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1329 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1330 nbl->table_elec->ninteractions = 1;
1331 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1332 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1334 snew(nbl->table_vdw, 1);
1335 nbl->table_vdw->interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1336 nbl->table_vdw->format = nbl->table_elec_vdw->format;
1337 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1338 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1339 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1340 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1341 nbl->table_vdw->ninteractions = 2;
1342 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1343 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1345 /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1346 * with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1348 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1350 for (j = 0; j < 4; j++)
1352 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1355 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1357 for (j = 0; j < 8; j++)
1359 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1364 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1365 * ftype2 present in the topology, build an array of the number of
1366 * interactions present for each bonded interaction index found in the
1369 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1370 * valid type with that parameter.
1372 * \c count will be reallocated as necessary to fit the largest bonded
1373 * interaction index found, and its current size will be returned in
1374 * \c ncount. It will contain zero for every bonded interaction index
1375 * for which no interactions are present in the topology.
1377 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1378 int *ncount, int **count)
1381 int ftype, stride, i, j, tabnr;
1383 // Loop over all moleculetypes
1384 for (const gmx_moltype_t &molt : mtop->moltype)
1386 // Loop over all interaction types
1387 for (ftype = 0; ftype < F_NRE; ftype++)
1389 // If the current interaction type is one of the types whose tables we're trying to count...
1390 if (ftype == ftype1 || ftype == ftype2)
1392 il = &molt.ilist[ftype];
1393 stride = 1 + NRAL(ftype);
1394 // ... and there are actually some interactions for this type
1395 for (i = 0; i < il->nr; i += stride)
1397 // Find out which table index the user wanted
1398 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1401 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1403 // Make room for this index in the data structure
1404 if (tabnr >= *ncount)
1406 srenew(*count, tabnr+1);
1407 for (j = *ncount; j < tabnr+1; j++)
1413 // Record that this table index is used and must have a valid file
1421 /*!\brief If there's bonded interactions of flavour \c tabext and type
1422 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1423 * list of filenames passed to mdrun, and make bonded tables from
1426 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1427 * valid type with that parameter.
1429 * A fatal error occurs if no matching filename is found.
1431 static bondedtable_t *make_bonded_tables(FILE *fplog,
1432 int ftype1, int ftype2,
1433 const gmx_mtop_t *mtop,
1434 gmx::ArrayRef<const std::string> tabbfnm,
1444 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1446 // Are there any relevant tabulated bond interactions?
1450 for (int i = 0; i < ncount; i++)
1452 // Do any interactions exist that requires this table?
1455 // This pattern enforces the current requirement that
1456 // table filenames end in a characteristic sequence
1457 // before the file type extension, and avoids table 13
1458 // being recognized and used for table 1.
1459 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1460 bool madeTable = false;
1461 for (size_t j = 0; j < tabbfnm.size() && !madeTable; ++j)
1463 if (gmx::endsWith(tabbfnm[j], patternToFind))
1465 // Finally read the table from the file found
1466 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1472 bool isPlural = (ftype2 != -1);
1473 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.",
1474 interaction_function[ftype1].longname,
1475 isPlural ? "' or '" : "",
1476 isPlural ? interaction_function[ftype2].longname : "",
1478 patternToFind.c_str());
1488 void forcerec_set_ranges(t_forcerec *fr,
1489 int ncg_home, int ncg_force,
1491 int natoms_force_constr, int natoms_f_novirsum)
1496 /* fr->ncg_force is unused in the standard code,
1497 * but it can be useful for modified code dealing with charge groups.
1499 fr->ncg_force = ncg_force;
1500 fr->natoms_force = natoms_force;
1501 fr->natoms_force_constr = natoms_force_constr;
1503 if (fr->natoms_force_constr > fr->nalloc_force)
1505 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1508 if (fr->haveDirectVirialContributions)
1510 fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1514 static real cutoff_inf(real cutoff)
1518 cutoff = GMX_CUTOFF_INF;
1524 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, const t_commrec *cr, FILE *fp)
1531 ir->rcoulomb == 0 &&
1533 ir->ePBC == epbcNONE &&
1534 ir->vdwtype == evdwCUT &&
1535 ir->coulombtype == eelCUT &&
1536 ir->efep == efepNO &&
1537 getenv("GMX_NO_ALLVSALL") == nullptr
1540 if (bAllvsAll && ir->opts.ngener > 1)
1542 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";
1548 fprintf(fp, "\n%s\n", note);
1554 if (bAllvsAll && fp && MASTER(cr))
1556 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1563 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
1564 const t_inputrec *ir)
1566 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1568 /* LJ PME with LB combination rule does 7 mesh operations.
1569 * This so slow that we don't compile SIMD non-bonded kernels
1571 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
1579 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1582 const gmx_hw_info_t gmx_unused &hardwareInfo)
1584 *kernel_type = nbnxnk4x4_PlainC;
1585 *ewald_excl = ewaldexclTable;
1589 #ifdef GMX_NBNXN_SIMD_4XN
1590 *kernel_type = nbnxnk4xN_SIMD_4xN;
1592 #ifdef GMX_NBNXN_SIMD_2XNN
1593 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1596 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1597 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1598 * This is based on the SIMD acceleration choice and CPU information
1599 * detected at runtime.
1601 * 4xN calculates more (zero) interactions, but has less pair-search
1602 * work and much better kernel instruction scheduling.
1604 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1605 * which doesn't have FMA, both the analytical and tabulated Ewald
1606 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1607 * 2x(4+4) because it results in significantly fewer pairs.
1608 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1609 * 10% with HT, 50% without HT. As we currently don't detect the actual
1610 * use of HT, use 4x8 to avoid a potential performance hit.
1611 * On Intel Haswell 4x8 is always faster.
1613 *kernel_type = nbnxnk4xN_SIMD_4xN;
1615 #if !GMX_SIMD_HAVE_FMA
1616 if (EEL_PME_EWALD(ir->coulombtype) ||
1617 EVDW_PME(ir->vdwtype))
1619 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1620 * There are enough instructions to make 2x(4+4) efficient.
1622 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1625 if (hardwareInfo.haveAmdZenCpu)
1627 /* One 256-bit FMA per cycle makes 2xNN faster */
1628 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1630 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1633 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
1635 #ifdef GMX_NBNXN_SIMD_4XN
1636 *kernel_type = nbnxnk4xN_SIMD_4xN;
1638 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1641 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
1643 #ifdef GMX_NBNXN_SIMD_2XNN
1644 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1646 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1650 /* Analytical Ewald exclusion correction is only an option in
1652 * Since table lookup's don't parallelize with SIMD, analytical
1653 * will probably always be faster for a SIMD width of 8 or more.
1654 * With FMA analytical is sometimes faster for a width if 4 as well.
1655 * In single precision, this is faster on Bulldozer.
1657 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1658 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
1659 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
1660 * of single or double precision and 128 or 256-bit AVX2.
1662 if (!hardwareInfo.haveAmdZenCpu)
1664 *ewald_excl = ewaldexclAnalytical;
1667 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
1669 *ewald_excl = ewaldexclTable;
1671 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
1673 *ewald_excl = ewaldexclAnalytical;
1681 const char *lookup_nbnxn_kernel_name(int kernel_type)
1683 const char *returnvalue = nullptr;
1684 switch (kernel_type)
1687 returnvalue = "not set";
1689 case nbnxnk4x4_PlainC:
1690 returnvalue = "plain C";
1692 case nbnxnk4xN_SIMD_4xN:
1693 case nbnxnk4xN_SIMD_2xNN:
1695 returnvalue = "SIMD";
1697 returnvalue = "not available";
1700 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1701 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1705 gmx_fatal(FARGS, "Illegal kernel type selected");
1706 returnvalue = nullptr;
1712 static void pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
1713 gmx_bool use_simd_kernels,
1714 const gmx_hw_info_t &hardwareInfo,
1716 EmulateGpuNonbonded emulateGpu,
1717 const t_inputrec *ir,
1720 gmx_bool bDoNonbonded)
1722 assert(kernel_type);
1724 *kernel_type = nbnxnkNotSet;
1725 *ewald_excl = ewaldexclTable;
1727 if (emulateGpu == EmulateGpuNonbonded::Yes)
1729 *kernel_type = nbnxnk8x8x8_PlainC;
1733 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
1738 *kernel_type = nbnxnk8x8x8_GPU;
1741 if (*kernel_type == nbnxnkNotSet)
1743 if (use_simd_kernels &&
1744 nbnxn_simd_supported(mdlog, ir))
1746 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl, hardwareInfo);
1750 *kernel_type = nbnxnk4x4_PlainC;
1756 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
1757 "Using %s %dx%d nonbonded short-range kernels",
1758 lookup_nbnxn_kernel_name(*kernel_type),
1759 nbnxn_kernel_to_cluster_i_size(*kernel_type),
1760 nbnxn_kernel_to_cluster_j_size(*kernel_type));
1762 if (nbnxnk4x4_PlainC == *kernel_type ||
1763 nbnxnk8x8x8_PlainC == *kernel_type)
1765 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1766 "WARNING: Using the slow %s kernels. This should\n"
1767 "not happen during routine usage on supported platforms.",
1768 lookup_nbnxn_kernel_name(*kernel_type));
1773 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1774 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1775 bool systemHasNetCharge,
1776 interaction_const_t *ic)
1778 if (!EEL_PME_EWALD(ir->coulombtype))
1785 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1787 if (ir->coulombtype == eelP3M_AD)
1789 please_cite(fp, "Hockney1988");
1790 please_cite(fp, "Ballenegger2012");
1794 please_cite(fp, "Essmann95a");
1797 if (ir->ewald_geometry == eewg3DC)
1801 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1802 systemHasNetCharge ? " and net charge" : "");
1804 please_cite(fp, "In-Chul99a");
1805 if (systemHasNetCharge)
1807 please_cite(fp, "Ballenegger2009");
1812 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1815 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1816 1/ic->ewaldcoeff_q);
1819 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1821 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1822 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1830 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1831 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1832 interaction_const_t *ic)
1834 if (!EVDW_PME(ir->vdwtype))
1841 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1842 please_cite(fp, "Essmann95a");
1844 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1847 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1848 1/ic->ewaldcoeff_lj);
1851 if (ic->vdw_modifier == eintmodPOTSHIFT)
1853 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1854 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1858 ic->sh_lj_ewald = 0;
1862 gmx_bool uses_simple_tables(int cutoff_scheme,
1863 nonbonded_verlet_t *nbv,
1866 gmx_bool bUsesSimpleTables = TRUE;
1869 switch (cutoff_scheme)
1872 bUsesSimpleTables = TRUE;
1875 assert(nullptr != nbv);
1876 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1877 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1880 gmx_incons("unimplemented");
1882 return bUsesSimpleTables;
1885 static void init_ewald_f_table(interaction_const_t *ic,
1890 /* Get the Ewald table spacing based on Coulomb and/or LJ
1891 * Ewald coefficients and rtol.
1893 ic->tabq_scale = ewald_spline3_table_scale(ic);
1895 if (ic->cutoff_scheme == ecutsVERLET)
1897 maxr = ic->rcoulomb;
1901 maxr = std::max(ic->rcoulomb, rtab);
1903 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1905 sfree_aligned(ic->tabq_coul_FDV0);
1906 sfree_aligned(ic->tabq_coul_F);
1907 sfree_aligned(ic->tabq_coul_V);
1909 sfree_aligned(ic->tabq_vdw_FDV0);
1910 sfree_aligned(ic->tabq_vdw_F);
1911 sfree_aligned(ic->tabq_vdw_V);
1913 if (EEL_PME_EWALD(ic->eeltype))
1915 /* Create the original table data in FDV0 */
1916 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1917 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1918 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1919 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1920 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1923 if (EVDW_PME(ic->vdwtype))
1925 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1926 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1927 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1928 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1929 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1933 void init_interaction_const_tables(FILE *fp,
1934 interaction_const_t *ic,
1937 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1939 init_ewald_f_table(ic, rtab);
1943 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1944 1/ic->tabq_scale, ic->tabq_size);
1949 static void clear_force_switch_constants(shift_consts_t *sc)
1956 static void force_switch_constants(real p,
1960 /* Here we determine the coefficient for shifting the force to zero
1961 * between distance rsw and the cut-off rc.
1962 * For a potential of r^-p, we have force p*r^-(p+1).
1963 * But to save flops we absorb p in the coefficient.
1965 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1966 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1968 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1969 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1970 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1973 static void potential_switch_constants(real rsw, real rc,
1974 switch_consts_t *sc)
1976 /* The switch function is 1 at rsw and 0 at rc.
1977 * The derivative and second derivate are zero at both ends.
1978 * rsw = max(r - r_switch, 0)
1979 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1980 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1981 * force = force*dsw - potential*sw
1984 sc->c3 = -10/gmx::power3(rc - rsw);
1985 sc->c4 = 15/gmx::power4(rc - rsw);
1986 sc->c5 = -6/gmx::power5(rc - rsw);
1989 /*! \brief Construct interaction constants
1991 * This data is used (particularly) by search and force code for
1992 * short-range interactions. Many of these are constant for the whole
1993 * simulation; some are constant only after PME tuning completes.
1996 init_interaction_const(FILE *fp,
1997 interaction_const_t **interaction_const,
1998 const t_inputrec *ir,
1999 const gmx_mtop_t *mtop,
2000 bool systemHasNetCharge)
2002 interaction_const_t *ic;
2006 ic->cutoff_scheme = ir->cutoff_scheme;
2008 /* Just allocate something so we can free it */
2009 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
2010 snew_aligned(ic->tabq_coul_F, 16, 32);
2011 snew_aligned(ic->tabq_coul_V, 16, 32);
2014 ic->vdwtype = ir->vdwtype;
2015 ic->vdw_modifier = ir->vdw_modifier;
2016 ic->reppow = mtop->ffparams.reppow;
2017 ic->rvdw = cutoff_inf(ir->rvdw);
2018 ic->rvdw_switch = ir->rvdw_switch;
2019 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
2020 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
2021 if (ic->useBuckingham)
2023 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
2026 initVdwEwaldParameters(fp, ir, ic);
2028 clear_force_switch_constants(&ic->dispersion_shift);
2029 clear_force_switch_constants(&ic->repulsion_shift);
2031 switch (ic->vdw_modifier)
2033 case eintmodPOTSHIFT:
2034 /* Only shift the potential, don't touch the force */
2035 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
2036 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
2038 case eintmodFORCESWITCH:
2039 /* Switch the force, switch and shift the potential */
2040 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2041 &ic->dispersion_shift);
2042 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2043 &ic->repulsion_shift);
2045 case eintmodPOTSWITCH:
2046 /* Switch the potential and force */
2047 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2051 case eintmodEXACTCUTOFF:
2052 /* Nothing to do here */
2055 gmx_incons("unimplemented potential modifier");
2058 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2060 /* Electrostatics */
2061 ic->eeltype = ir->coulombtype;
2062 ic->coulomb_modifier = ir->coulomb_modifier;
2063 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
2064 ic->rcoulomb_switch = ir->rcoulomb_switch;
2065 ic->epsilon_r = ir->epsilon_r;
2067 /* Set the Coulomb energy conversion factor */
2068 if (ic->epsilon_r != 0)
2070 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
2074 /* eps = 0 is infinite dieletric: no Coulomb interactions */
2078 /* Reaction-field */
2079 if (EEL_RF(ic->eeltype))
2081 ic->epsilon_rf = ir->epsilon_rf;
2082 /* Generalized reaction field parameters are updated every step */
2083 if (ic->eeltype != eelGRF)
2085 calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
2086 ic->rcoulomb, 0, 0, nullptr,
2087 &ic->k_rf, &ic->c_rf);
2090 if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
2092 /* grompp should have done this, but this scheme is obsolete */
2093 ic->coulomb_modifier = eintmodEXACTCUTOFF;
2098 /* For plain cut-off we might use the reaction-field kernels */
2099 ic->epsilon_rf = ic->epsilon_r;
2101 if (ir->coulomb_modifier == eintmodPOTSHIFT)
2103 ic->c_rf = 1/ic->rcoulomb;
2111 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
2115 real dispersion_shift;
2117 dispersion_shift = ic->dispersion_shift.cpot;
2118 if (EVDW_PME(ic->vdwtype))
2120 dispersion_shift -= ic->sh_lj_ewald;
2122 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2123 ic->repulsion_shift.cpot, dispersion_shift);
2125 if (ic->eeltype == eelCUT)
2127 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2129 else if (EEL_PME(ic->eeltype))
2131 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2136 *interaction_const = ic;
2140 done_interaction_const(interaction_const_t *interaction_const)
2142 sfree_aligned(interaction_const->tabq_coul_FDV0);
2143 sfree_aligned(interaction_const->tabq_coul_F);
2144 sfree_aligned(interaction_const->tabq_coul_V);
2145 sfree(interaction_const);
2148 static void init_nb_verlet(const gmx::MDLogger &mdlog,
2149 nonbonded_verlet_t **nb_verlet,
2150 gmx_bool bFEP_NonBonded,
2151 const t_inputrec *ir,
2152 const t_forcerec *fr,
2153 const t_commrec *cr,
2154 const gmx_hw_info_t &hardwareInfo,
2155 const gmx_device_info_t *deviceInfo,
2156 const gmx_mtop_t *mtop,
2159 nonbonded_verlet_t *nbv;
2162 nbnxn_alloc_t *nb_alloc;
2163 nbnxn_free_t *nb_free;
2165 nbv = new nonbonded_verlet_t();
2167 nbv->emulateGpu = ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
2168 nbv->bUseGPU = deviceInfo != nullptr;
2170 GMX_RELEASE_ASSERT(!(nbv->emulateGpu == EmulateGpuNonbonded::Yes && nbv->bUseGPU), "When GPU emulation is active, there cannot be a GPU assignment");
2173 nbv->min_ci_balanced = 0;
2175 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2176 for (int i = 0; i < nbv->ngrp; i++)
2178 nbv->grp[i].nbl_lists.nnbl = 0;
2179 nbv->grp[i].kernel_type = nbnxnkNotSet;
2181 if (i == 0) /* local */
2183 pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
2184 nbv->bUseGPU, nbv->emulateGpu, ir,
2185 &nbv->grp[i].kernel_type,
2186 &nbv->grp[i].ewald_excl,
2189 else /* non-local */
2191 /* Use the same kernel for local and non-local interactions */
2192 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2193 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2197 nbv->listParams = std::unique_ptr<NbnxnListParameters>(new NbnxnListParameters(ir->rlist));
2198 setupDynamicPairlistPruning(mdlog, ir, mtop, box, nbv->grp[0].kernel_type, fr->ic,
2199 nbv->listParams.get());
2201 nbv->nbs = gmx::compat::make_unique<nbnxn_search>(DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
2202 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
2204 gmx_omp_nthreads_get(emntPairsearch));
2206 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2207 &nb_alloc, &nb_free);
2209 for (int i = 0; i < nbv->ngrp; i++)
2211 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2212 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2213 /* 8x8x8 "non-simple" lists are ATM always combined */
2214 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2218 int enbnxninitcombrule;
2219 if (fr->ic->vdwtype == evdwCUT &&
2220 (fr->ic->vdw_modifier == eintmodNONE ||
2221 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
2222 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
2224 /* Plain LJ cut-off: we can optimize with combination rules */
2225 enbnxninitcombrule = enbnxninitcombruleDETECT;
2227 else if (fr->ic->vdwtype == evdwPME)
2229 /* LJ-PME: we need to use a combination rule for the grid */
2230 if (fr->ljpme_combination_rule == eljpmeGEOM)
2232 enbnxninitcombrule = enbnxninitcombruleGEOM;
2236 enbnxninitcombrule = enbnxninitcombruleLB;
2241 /* We use a full combination matrix: no rule required */
2242 enbnxninitcombrule = enbnxninitcombruleNONE;
2246 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
2247 if (ir->opts.ngener - ir->nwall == 1)
2249 /* We have only one non-wall energy group, we do not need energy group
2250 * support in the non-bondeds kernels, since all non-bonded energy
2251 * contributions go to the first element of the energy group matrix.
2253 mimimumNumEnergyGroupNonbonded = 1;
2255 bool bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[0].kernel_type);
2256 nbnxn_atomdata_init(mdlog,
2258 nbv->grp[0].kernel_type,
2260 fr->ntype, fr->nbfp,
2261 mimimumNumEnergyGroupNonbonded,
2262 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2267 /* init the NxN GPU data; the last argument tells whether we'll have
2268 * both local and non-local NB calculation on GPU */
2269 nbnxn_gpu_init(&nbv->gpu_nbv,
2272 nbv->listParams.get(),
2277 if ((env = getenv("GMX_NB_MIN_CI")) != nullptr)
2281 nbv->min_ci_balanced = strtol(env, &end, 10);
2282 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
2284 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
2289 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2290 nbv->min_ci_balanced);
2295 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2298 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2299 nbv->min_ci_balanced);
2308 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2310 return nbv != nullptr && nbv->bUseGPU;
2313 void init_forcerec(FILE *fp,
2314 const gmx::MDLogger &mdlog,
2317 const t_inputrec *ir,
2318 const gmx_mtop_t *mtop,
2319 const t_commrec *cr,
2323 gmx::ArrayRef<const std::string> tabbfnm,
2324 const gmx_hw_info_t &hardwareInfo,
2325 const gmx_device_info_t *deviceInfo,
2326 gmx_bool bNoSolvOpt,
2329 int m, negp_pp, negptable, egi, egj;
2334 gmx_bool bGenericKernelOnly;
2335 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
2336 gmx_bool bFEP_NonBonded;
2337 int *nm_ind, egp_flags;
2339 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2340 fr->use_simd_kernels = TRUE;
2342 if (check_box(ir->ePBC, box))
2344 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2347 /* Test particle insertion ? */
2350 /* Set to the size of the molecule to be inserted (the last one) */
2351 /* Because of old style topologies, we have to use the last cg
2352 * instead of the last molecule type.
2354 cgs = &mtop->moltype[mtop->molblock.back().type].cgs;
2355 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2356 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
2357 if (fr->n_tpi != molecules.block(molecules.numBlocks() - 1).size())
2359 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2367 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
2369 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
2370 eel_names[ir->coulombtype]);
2375 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
2377 if (ir->useTwinRange)
2379 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
2381 /* Copy the user determined parameters */
2382 fr->userint1 = ir->userint1;
2383 fr->userint2 = ir->userint2;
2384 fr->userint3 = ir->userint3;
2385 fr->userint4 = ir->userint4;
2386 fr->userreal1 = ir->userreal1;
2387 fr->userreal2 = ir->userreal2;
2388 fr->userreal3 = ir->userreal3;
2389 fr->userreal4 = ir->userreal4;
2392 fr->fc_stepsize = ir->fc_stepsize;
2395 fr->efep = ir->efep;
2396 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2397 if (ir->fepvals->bScCoul)
2399 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2400 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
2404 fr->sc_alphacoul = 0;
2405 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2407 fr->sc_power = ir->fepvals->sc_power;
2408 fr->sc_r_power = ir->fepvals->sc_r_power;
2409 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
2411 env = getenv("GMX_SCSIGMA_MIN");
2415 sscanf(env, "%20lf", &dbl);
2416 fr->sc_sigma6_min = gmx::power6(dbl);
2419 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2423 fr->bNonbonded = TRUE;
2424 if (getenv("GMX_NO_NONBONDED") != nullptr)
2426 /* turn off non-bonded calculations */
2427 fr->bNonbonded = FALSE;
2428 GMX_LOG(mdlog.warning).asParagraph().appendText(
2429 "Found environment variable GMX_NO_NONBONDED.\n"
2430 "Disabling nonbonded calculations.");
2433 bGenericKernelOnly = FALSE;
2435 /* We now check in the NS code whether a particular combination of interactions
2436 * can be used with water optimization, and disable it if that is not the case.
2439 if (getenv("GMX_NB_GENERIC") != nullptr)
2444 "Found environment variable GMX_NB_GENERIC.\n"
2445 "Disabling all interaction-specific nonbonded kernels, will only\n"
2446 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2448 bGenericKernelOnly = TRUE;
2451 if (bGenericKernelOnly == TRUE)
2456 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
2458 fr->use_simd_kernels = FALSE;
2462 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2463 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2464 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2468 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2470 /* Check if we can/should do all-vs-all kernels */
2471 fr->bAllvsAll = can_use_allvsall(ir, FALSE, nullptr, nullptr);
2472 fr->AllvsAll_work = nullptr;
2474 /* All-vs-all kernels have not been implemented in 4.6 and later.
2475 * See Redmine #1249. */
2478 fr->bAllvsAll = FALSE;
2482 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2483 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2484 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2485 "or try cutoff-scheme = Verlet.\n\n");
2489 /* Neighbour searching stuff */
2490 fr->cutoff_scheme = ir->cutoff_scheme;
2491 fr->bGrid = (ir->ns_type == ensGRID);
2492 fr->ePBC = ir->ePBC;
2494 if (fr->cutoff_scheme == ecutsGROUP)
2496 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2497 "removed in a future release when 'verlet' supports all interaction forms.\n";
2501 fprintf(stderr, "\n%s\n", note);
2505 fprintf(fp, "\n%s\n", note);
2509 /* Determine if we will do PBC for distances in bonded interactions */
2510 if (fr->ePBC == epbcNONE)
2512 fr->bMolPBC = FALSE;
2516 if (!DOMAINDECOMP(cr))
2520 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2521 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2522 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2524 /* The group cut-off scheme and SHAKE assume charge groups
2525 * are whole, but not using molpbc is faster in most cases.
2526 * With intermolecular interactions we need PBC for calculating
2527 * distances between atoms in different molecules.
2529 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2530 !mtop->bIntermolecularInteractions)
2532 fr->bMolPBC = ir->bPeriodicMols;
2534 if (bSHAKE && fr->bMolPBC)
2536 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2541 /* Not making molecules whole is faster in most cases,
2542 * but With orientation restraints we need whole molecules.
2544 fr->bMolPBC = (fcd->orires.nr == 0);
2546 if (getenv("GMX_USE_GRAPH") != nullptr)
2548 fr->bMolPBC = FALSE;
2551 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
2554 if (mtop->bIntermolecularInteractions)
2556 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
2560 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
2562 if (bSHAKE && fr->bMolPBC)
2564 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");
2570 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2574 fr->rc_scaling = ir->refcoord_scaling;
2575 copy_rvec(ir->posres_com, fr->posres_com);
2576 copy_rvec(ir->posres_comB, fr->posres_comB);
2577 fr->rlist = cutoff_inf(ir->rlist);
2578 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2580 /* This now calculates sum for q and c6*/
2581 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
2583 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2584 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
2585 init_interaction_const_tables(fp, fr->ic, ir->rlist + ir->tabext);
2587 const interaction_const_t *ic = fr->ic;
2589 /* TODO: Replace this Ewald table or move it into interaction_const_t */
2590 if (ir->coulombtype == eelEWALD)
2592 init_ewald_tab(&(fr->ewald_table), ir, fp);
2595 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2596 switch (ic->eeltype)
2599 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
2604 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2608 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2609 GMX_RELEASE_ASSERT(ic->coulomb_modifier == eintmodEXACTCUTOFF, "With the group scheme RF-zero needs the exact cut-off modifier");
2618 case eelPMEUSERSWITCH:
2619 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2625 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2629 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
2632 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
2634 /* Vdw: Translate from mdp settings to kernel format */
2635 switch (ic->vdwtype)
2640 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2644 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2648 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2654 case evdwENCADSHIFT:
2655 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2659 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
2662 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
2664 if (ir->cutoff_scheme == ecutsGROUP)
2666 fr->bvdwtab = ((ic->vdwtype != evdwCUT || !gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2667 && !EVDW_PME(ic->vdwtype));
2668 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2669 fr->bcoultab = !(ic->eeltype == eelCUT ||
2670 ic->eeltype == eelEWALD ||
2671 ic->eeltype == eelPME ||
2672 ic->eeltype == eelP3M_AD ||
2673 ic->eeltype == eelRF ||
2674 ic->eeltype == eelRF_ZERO);
2676 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2677 * going to be faster to tabulate the interaction than calling the generic kernel.
2678 * However, if generic kernels have been requested we keep things analytically.
2680 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2681 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2682 bGenericKernelOnly == FALSE)
2684 if ((ic->rcoulomb_switch != ic->rvdw_switch) || (ic->rcoulomb != ic->rvdw))
2686 fr->bcoultab = TRUE;
2687 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2688 * which would otherwise need two tables.
2692 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2693 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2694 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2695 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2697 if ((ic->rcoulomb != ic->rvdw) && (bGenericKernelOnly == FALSE))
2699 fr->bcoultab = TRUE;
2703 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2705 fr->bcoultab = TRUE;
2707 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2712 if (getenv("GMX_REQUIRE_TABLES"))
2715 fr->bcoultab = TRUE;
2720 fprintf(fp, "Table routines are used for coulomb: %s\n",
2721 gmx::boolToString(fr->bcoultab));
2722 fprintf(fp, "Table routines are used for vdw: %s\n",
2723 gmx::boolToString(fr->bvdwtab));
2726 if (fr->bvdwtab == TRUE)
2728 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2729 fr->nbkernel_vdw_modifier = eintmodNONE;
2731 if (fr->bcoultab == TRUE)
2733 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2734 fr->nbkernel_elec_modifier = eintmodNONE;
2738 if (ir->cutoff_scheme == ecutsVERLET)
2740 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2742 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2744 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2745 * while mdrun does not (and never did) support this.
2747 if (EEL_USER(fr->ic->eeltype))
2749 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2750 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2753 fr->bvdwtab = FALSE;
2754 fr->bcoultab = FALSE;
2757 /* 1-4 interaction electrostatics */
2758 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2760 /* Parameters for generalized RF */
2764 if (ic->eeltype == eelGRF)
2766 init_generalized_rf(fp, mtop, ir, fr);
2769 fr->haveDirectVirialContributions =
2770 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
2771 fr->forceProviders->hasForceProvider() ||
2772 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2773 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2778 if (fr->haveDirectVirialContributions)
2780 fr->forceBufferForDirectVirialContributions = new std::vector<gmx::RVec>;
2783 if (fr->cutoff_scheme == ecutsGROUP &&
2784 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2786 /* Count the total number of charge groups */
2787 fr->cg_nalloc = ncg_mtop(mtop);
2788 srenew(fr->cg_cm, fr->cg_nalloc);
2790 if (fr->shift_vec == nullptr)
2792 snew(fr->shift_vec, SHIFTS);
2795 if (fr->fshift == nullptr)
2797 snew(fr->fshift, SHIFTS);
2800 if (fr->nbfp == nullptr)
2802 fr->ntype = mtop->ffparams.atnr;
2803 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2804 if (EVDW_PME(ic->vdwtype))
2806 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2810 /* Copy the energy group exclusions */
2811 fr->egp_flags = ir->opts.egp_flags;
2813 /* Van der Waals stuff */
2814 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
2816 if (ic->rvdw_switch >= ic->rvdw)
2818 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2819 ic->rvdw_switch, ic->rvdw);
2823 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2824 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
2825 ic->rvdw_switch, ic->rvdw);
2829 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
2831 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2834 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
2836 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2839 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2841 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2844 if (fp && fr->cutoff_scheme == ecutsGROUP)
2846 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2847 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
2850 fr->eDispCorr = ir->eDispCorr;
2851 fr->numAtomsForDispersionCorrection = mtop->natoms;
2852 if (ir->eDispCorr != edispcNO)
2854 set_avcsixtwelve(fp, fr, mtop);
2857 if (ir->implicit_solvent)
2859 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
2862 /* Construct tables for the group scheme. A little unnecessary to
2863 * make both vdw and coul tables sometimes, but what the
2864 * heck. Note that both cutoff schemes construct Ewald tables in
2865 * init_interaction_const_tables. */
2866 needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
2867 (fr->bcoultab || fr->bvdwtab));
2869 negp_pp = ir->opts.ngener - ir->nwall;
2871 if (!needGroupSchemeTables)
2873 bSomeNormalNbListsAreInUse = TRUE;
2878 bSomeNormalNbListsAreInUse = FALSE;
2879 for (egi = 0; egi < negp_pp; egi++)
2881 for (egj = egi; egj < negp_pp; egj++)
2883 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2884 if (!(egp_flags & EGP_EXCL))
2886 if (egp_flags & EGP_TABLE)
2892 bSomeNormalNbListsAreInUse = TRUE;
2897 if (bSomeNormalNbListsAreInUse)
2899 fr->nnblists = negptable + 1;
2903 fr->nnblists = negptable;
2905 if (fr->nnblists > 1)
2907 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
2911 snew(fr->nblists, fr->nnblists);
2913 /* This code automatically gives table length tabext without cut-off's,
2914 * in that case grompp should already have checked that we do not need
2915 * normal tables and we only generate tables for 1-4 interactions.
2917 rtab = ir->rlist + ir->tabext;
2919 if (needGroupSchemeTables)
2921 /* make tables for ordinary interactions */
2922 if (bSomeNormalNbListsAreInUse)
2924 make_nbf_tables(fp, ic, rtab, tabfn, nullptr, nullptr, &fr->nblists[0]);
2933 /* Read the special tables for certain energy group pairs */
2934 nm_ind = mtop->groups.grps[egcENER].nm_ind;
2935 for (egi = 0; egi < negp_pp; egi++)
2937 for (egj = egi; egj < negp_pp; egj++)
2939 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2940 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2942 if (fr->nnblists > 1)
2944 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
2946 /* Read the table file with the two energy groups names appended */
2947 make_nbf_tables(fp, ic, rtab, tabfn,
2948 *mtop->groups.grpname[nm_ind[egi]],
2949 *mtop->groups.grpname[nm_ind[egj]],
2953 else if (fr->nnblists > 1)
2955 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
2962 /* Tables might not be used for the potential modifier
2963 * interactions per se, but we still need them to evaluate
2964 * switch/shift dispersion corrections in this case. */
2965 if (fr->eDispCorr != edispcNO)
2967 fr->dispersionCorrectionTable = makeDispersionCorrectionTable(fp, ic, rtab, tabfn);
2970 /* We want to use unmodified tables for 1-4 coulombic
2971 * interactions, so we must in general have an extra set of
2973 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2974 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2975 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
2977 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
2978 GMX_MAKETABLES_14ONLY);
2982 fr->nwall = ir->nwall;
2983 if (ir->nwall && ir->wall_type == ewtTABLE)
2985 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
2988 if (fcd && !tabbfnm.empty())
2990 // Need to catch std::bad_alloc
2991 // TODO Don't need to catch this here, when merging with master branch
2994 fcd->bondtab = make_bonded_tables(fp,
2995 F_TABBONDS, F_TABBONDSNC,
2996 mtop, tabbfnm, "b");
2997 fcd->angletab = make_bonded_tables(fp,
2999 mtop, tabbfnm, "a");
3000 fcd->dihtab = make_bonded_tables(fp,
3002 mtop, tabbfnm, "d");
3004 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
3010 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3014 /* QM/MM initialization if requested
3018 fprintf(stderr, "QM/MM calculation requested.\n");
3021 fr->bQMMM = ir->bQMMM;
3022 fr->qr = mk_QMMMrec();
3024 /* Set all the static charge group info */
3025 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3027 &fr->bExcl_IntraCGAll_InterCGNone);
3028 if (DOMAINDECOMP(cr))
3030 fr->cginfo = nullptr;
3034 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
3037 if (!DOMAINDECOMP(cr))
3039 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3040 mtop->natoms, mtop->natoms, mtop->natoms);
3043 fr->print_force = print_force;
3046 /* coarse load balancing vars */
3051 /* Initialize neighbor search */
3053 init_ns(fp, cr, fr->ns, fr, mtop);
3055 if (thisRankHasDuty(cr, DUTY_PP))
3057 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3060 /* Initialize the thread working data for bonded interactions */
3061 init_bonded_threading(fp, mtop->groups.grps[egcENER].nr,
3062 &fr->bondedThreading);
3064 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
3065 snew(fr->ewc_t, fr->nthread_ewc);
3067 if (fr->cutoff_scheme == ecutsVERLET)
3069 // We checked the cut-offs in grompp, but double-check here.
3070 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
3071 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
3073 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
3077 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
3080 init_nb_verlet(mdlog, &fr->nbv, bFEP_NonBonded, ir, fr,
3081 cr, hardwareInfo, deviceInfo,
3087 /* Here we switch from using mdlog, which prints the newline before
3088 * the paragraph, to our old fprintf logging, which prints the newline
3089 * after the paragraph, so we should add a newline here.
3094 if (ir->eDispCorr != edispcNO)
3096 calc_enervirdiff(fp, ir->eDispCorr, fr);
3100 /* Frees GPU memory and sets a tMPI node barrier.
3102 * Note that this function needs to be called even if GPUs are not used
3103 * in this run because the PME ranks have no knowledge of whether GPUs
3104 * are used or not, but all ranks need to enter the barrier below.
3105 * \todo Remove physical node barrier from this function after making sure
3106 * that it's not needed anymore (with a shared GPU run).
3108 void free_gpu_resources(const t_forcerec *fr,
3109 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator)
3111 bool isPPrankUsingGPU = fr && fr->nbv && fr->nbv->bUseGPU;
3113 /* stop the GPU profiler (only CUDA) */
3116 if (isPPrankUsingGPU)
3118 /* free nbnxn data in GPU memory */
3119 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3122 /* With tMPI we need to wait for all ranks to finish deallocation before
3123 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
3126 * This is not a concern in OpenCL where we use one context per rank which
3127 * is freed in nbnxn_gpu_free().
3129 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
3130 * but it is easier and more futureproof to call it on the whole node.
3134 physicalNodeCommunicator.barrier();
3138 void done_forcerec(t_forcerec *fr, int numMolBlocks, int numEnergyGroups)
3142 // PME-only ranks don't have a forcerec
3145 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
3147 done_interaction_const(fr->ic);
3148 sfree(fr->shift_vec);
3151 done_ns(fr->ns, numEnergyGroups);
3153 tear_down_bonded_threading(fr->bondedThreading);
3154 fr->bondedThreading = nullptr;