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/sim_util.h"
83 #include "gromacs/mdtypes/commrec.h"
84 #include "gromacs/mdtypes/fcdata.h"
85 #include "gromacs/mdtypes/group.h"
86 #include "gromacs/mdtypes/iforceprovider.h"
87 #include "gromacs/mdtypes/inputrec.h"
88 #include "gromacs/mdtypes/md_enums.h"
89 #include "gromacs/pbcutil/ishift.h"
90 #include "gromacs/pbcutil/pbc.h"
91 #include "gromacs/simd/simd.h"
92 #include "gromacs/tables/forcetable.h"
93 #include "gromacs/topology/mtop_util.h"
94 #include "gromacs/trajectory/trajectoryframe.h"
95 #include "gromacs/utility/cstringutil.h"
96 #include "gromacs/utility/exceptions.h"
97 #include "gromacs/utility/fatalerror.h"
98 #include "gromacs/utility/gmxassert.h"
99 #include "gromacs/utility/logger.h"
100 #include "gromacs/utility/physicalnodecommunicator.h"
101 #include "gromacs/utility/pleasecite.h"
102 #include "gromacs/utility/smalloc.h"
103 #include "gromacs/utility/strconvert.h"
105 #include "nbnxn_gpu_jit_support.h"
107 const char *egrp_nm[egNR+1] = {
108 "Coul-SR", "LJ-SR", "Buck-SR",
109 "Coul-14", "LJ-14", nullptr
112 t_forcerec *mk_forcerec(void)
121 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
129 snew(nbfp, 3*atnr*atnr);
130 for (i = k = 0; (i < atnr); i++)
132 for (j = 0; (j < atnr); j++, k++)
134 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
135 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
136 /* nbfp now includes the 6.0 derivative prefactor */
137 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
143 snew(nbfp, 2*atnr*atnr);
144 for (i = k = 0; (i < atnr); i++)
146 for (j = 0; (j < atnr); j++, k++)
148 /* nbfp now includes the 6.0/12.0 derivative prefactors */
149 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
150 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
158 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
161 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
164 /* For LJ-PME simulations, we correct the energies with the reciprocal space
165 * inside of the cut-off. To do this the non-bonded kernels needs to have
166 * access to the C6-values used on the reciprocal grid in pme.c
170 snew(grid, 2*atnr*atnr);
171 for (i = k = 0; (i < atnr); i++)
173 for (j = 0; (j < atnr); j++, k++)
175 c6i = idef->iparams[i*(atnr+1)].lj.c6;
176 c12i = idef->iparams[i*(atnr+1)].lj.c12;
177 c6j = idef->iparams[j*(atnr+1)].lj.c6;
178 c12j = idef->iparams[j*(atnr+1)].lj.c12;
179 c6 = std::sqrt(c6i * c6j);
180 if (fr->ljpme_combination_rule == eljpmeLB
181 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
183 sigmai = gmx::sixthroot(c12i / c6i);
184 sigmaj = gmx::sixthroot(c12j / c6j);
185 epsi = c6i * c6i / c12i;
186 epsj = c6j * c6j / c12j;
187 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
189 /* Store the elements at the same relative positions as C6 in nbfp in order
190 * to simplify access in the kernels
192 grid[2*(atnr*i+j)] = c6*6.0;
198 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
202 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
206 snew(nbfp, 2*atnr*atnr);
207 for (i = 0; i < atnr; ++i)
209 for (j = 0; j < atnr; ++j)
211 c6i = idef->iparams[i*(atnr+1)].lj.c6;
212 c12i = idef->iparams[i*(atnr+1)].lj.c12;
213 c6j = idef->iparams[j*(atnr+1)].lj.c6;
214 c12j = idef->iparams[j*(atnr+1)].lj.c12;
215 c6 = std::sqrt(c6i * c6j);
216 c12 = std::sqrt(c12i * c12j);
217 if (comb_rule == eCOMB_ARITHMETIC
218 && !gmx_numzero(c6) && !gmx_numzero(c12))
220 sigmai = gmx::sixthroot(c12i / c6i);
221 sigmaj = gmx::sixthroot(c12j / c6j);
222 epsi = c6i * c6i / c12i;
223 epsj = c6j * c6j / c12j;
224 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
225 c12 = std::sqrt(epsi * epsj) * gmx::power12(0.5*(sigmai+sigmaj));
227 C6(nbfp, atnr, i, j) = c6*6.0;
228 C12(nbfp, atnr, i, j) = c12*12.0;
234 /* This routine sets fr->solvent_opt to the most common solvent in the
235 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
236 * the fr->solvent_type array with the correct type (or esolNO).
238 * Charge groups that fulfill the conditions but are not identical to the
239 * most common one will be marked as esolNO in the solvent_type array.
241 * TIP3p is identical to SPC for these purposes, so we call it
242 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
244 * NOTE: QM particle should not
245 * become an optimized solvent. Not even if there is only one charge
255 } solvent_parameters_t;
258 check_solvent_cg(const gmx_moltype_t *molt,
261 const unsigned char *qm_grpnr,
262 const t_grps *qm_grps,
264 int *n_solvent_parameters,
265 solvent_parameters_t **solvent_parameters_p,
275 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
276 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
279 solvent_parameters_t *solvent_parameters;
281 /* We use a list with parameters for each solvent type.
282 * Every time we discover a new molecule that fulfills the basic
283 * conditions for a solvent we compare with the previous entries
284 * in these lists. If the parameters are the same we just increment
285 * the counter for that type, and otherwise we create a new type
286 * based on the current molecule.
288 * Once we've finished going through all molecules we check which
289 * solvent is most common, and mark all those molecules while we
290 * clear the flag on all others.
293 solvent_parameters = *solvent_parameters_p;
295 /* Mark the cg first as non optimized */
298 /* Check if this cg has no exclusions with atoms in other charge groups
299 * and all atoms inside the charge group excluded.
300 * We only have 3 or 4 atom solvent loops.
302 if (GET_CGINFO_EXCL_INTER(cginfo) ||
303 !GET_CGINFO_EXCL_INTRA(cginfo))
308 /* Get the indices of the first atom in this charge group */
309 j0 = molt->cgs.index[cg0];
310 j1 = molt->cgs.index[cg0+1];
312 /* Number of atoms in our molecule */
318 "Moltype '%s': there are %d atoms in this charge group\n",
322 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
325 if (nj < 3 || nj > 4)
330 /* Check if we are doing QM on this group */
332 if (qm_grpnr != nullptr)
334 for (j = j0; j < j1 && !qm; j++)
336 qm = (qm_grpnr[j] < qm_grps->nr - 1);
339 /* Cannot use solvent optimization with QM */
345 atom = molt->atoms.atom;
347 /* Still looks like a solvent, time to check parameters */
349 /* If it is perturbed (free energy) we can't use the solvent loops,
350 * so then we just skip to the next molecule.
354 for (j = j0; j < j1 && !perturbed; j++)
356 perturbed = PERTURBED(atom[j]);
364 /* Now it's only a question if the VdW and charge parameters
365 * are OK. Before doing the check we compare and see if they are
366 * identical to a possible previous solvent type.
367 * First we assign the current types and charges.
369 for (j = 0; j < nj; j++)
371 tmp_vdwtype[j] = atom[j0+j].type;
372 tmp_charge[j] = atom[j0+j].q;
375 /* Does it match any previous solvent type? */
376 for (k = 0; k < *n_solvent_parameters; k++)
381 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
382 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
383 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
388 /* Check that types & charges match for all atoms in molecule */
389 for (j = 0; j < nj && match == TRUE; j++)
391 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
395 if (tmp_charge[j] != solvent_parameters[k].charge[j])
402 /* Congratulations! We have a matched solvent.
403 * Flag it with this type for later processing.
406 solvent_parameters[k].count += nmol;
408 /* We are done with this charge group */
413 /* If we get here, we have a tentative new solvent type.
414 * Before we add it we must check that it fulfills the requirements
415 * of the solvent optimized loops. First determine which atoms have
418 for (j = 0; j < nj; j++)
421 tjA = tmp_vdwtype[j];
423 /* Go through all other tpes and see if any have non-zero
424 * VdW parameters when combined with this one.
426 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
428 /* We already checked that the atoms weren't perturbed,
429 * so we only need to check state A now.
433 has_vdw[j] = (has_vdw[j] ||
434 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
435 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
436 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
441 has_vdw[j] = (has_vdw[j] ||
442 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
443 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
448 /* Now we know all we need to make the final check and assignment. */
452 * For this we require thatn all atoms have charge,
453 * the charges on atom 2 & 3 should be the same, and only
454 * atom 1 might have VdW.
456 if (has_vdw[1] == FALSE &&
457 has_vdw[2] == FALSE &&
458 tmp_charge[0] != 0 &&
459 tmp_charge[1] != 0 &&
460 tmp_charge[2] == tmp_charge[1])
462 srenew(solvent_parameters, *n_solvent_parameters+1);
463 solvent_parameters[*n_solvent_parameters].model = esolSPC;
464 solvent_parameters[*n_solvent_parameters].count = nmol;
465 for (k = 0; k < 3; k++)
467 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
468 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
471 *cg_sp = *n_solvent_parameters;
472 (*n_solvent_parameters)++;
477 /* Or could it be a TIP4P?
478 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
479 * Only atom 1 mght have VdW.
481 if (has_vdw[1] == FALSE &&
482 has_vdw[2] == FALSE &&
483 has_vdw[3] == FALSE &&
484 tmp_charge[0] == 0 &&
485 tmp_charge[1] != 0 &&
486 tmp_charge[2] == tmp_charge[1] &&
489 srenew(solvent_parameters, *n_solvent_parameters+1);
490 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
491 solvent_parameters[*n_solvent_parameters].count = nmol;
492 for (k = 0; k < 4; k++)
494 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
495 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
498 *cg_sp = *n_solvent_parameters;
499 (*n_solvent_parameters)++;
503 *solvent_parameters_p = solvent_parameters;
507 check_solvent(FILE * fp,
508 const gmx_mtop_t * mtop,
510 cginfo_mb_t *cginfo_mb)
513 const gmx_moltype_t *molt;
514 int mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
515 int n_solvent_parameters;
516 solvent_parameters_t *solvent_parameters;
522 fprintf(debug, "Going to determine what solvent types we have.\n");
525 n_solvent_parameters = 0;
526 solvent_parameters = nullptr;
527 /* Allocate temporary array for solvent type */
528 snew(cg_sp, mtop->molblock.size());
531 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
533 molt = &mtop->moltype[mtop->molblock[mb].type];
535 /* Here we have to loop over all individual molecules
536 * because we need to check for QMMM particles.
538 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
539 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
540 nmol = mtop->molblock[mb].nmol/nmol_ch;
541 for (mol = 0; mol < nmol_ch; mol++)
544 am = mol*cgs->index[cgs->nr];
545 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
547 check_solvent_cg(molt, cg_mol, nmol,
548 mtop->groups.grpnr[egcQMMM] ?
549 mtop->groups.grpnr[egcQMMM]+at_offset+am : nullptr,
550 &mtop->groups.grps[egcQMMM],
552 &n_solvent_parameters, &solvent_parameters,
553 cginfo_mb[mb].cginfo[cgm+cg_mol],
554 &cg_sp[mb][cgm+cg_mol]);
557 at_offset += cgs->index[cgs->nr];
560 /* Puh! We finished going through all charge groups.
561 * Now find the most common solvent model.
564 /* Most common solvent this far */
566 for (i = 0; i < n_solvent_parameters; i++)
569 solvent_parameters[i].count > solvent_parameters[bestsp].count)
577 bestsol = solvent_parameters[bestsp].model;
585 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
587 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
588 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
589 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
591 if (cg_sp[mb][i] == bestsp)
593 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
598 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
605 if (bestsol != esolNO && fp != nullptr)
607 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
609 solvent_parameters[bestsp].count);
612 sfree(solvent_parameters);
613 fr->solvent_opt = bestsol;
617 acNONE = 0, acCONSTRAINT, acSETTLE
620 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
621 t_forcerec *fr, gmx_bool bNoSolvOpt,
622 gmx_bool *bFEP_NonBonded,
623 gmx_bool *bExcl_IntraCGAll_InterCGNone)
626 const t_blocka *excl;
627 const gmx_moltype_t *molt;
628 const gmx_molblock_t *molb;
629 cginfo_mb_t *cginfo_mb;
632 int cg_offset, a_offset;
633 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
637 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
639 snew(cginfo_mb, mtop->molblock.size());
641 snew(type_VDW, fr->ntype);
642 for (ai = 0; ai < fr->ntype; ai++)
644 type_VDW[ai] = FALSE;
645 for (j = 0; j < fr->ntype; j++)
647 type_VDW[ai] = type_VDW[ai] ||
649 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
650 C12(fr->nbfp, fr->ntype, ai, j) != 0;
654 *bFEP_NonBonded = FALSE;
655 *bExcl_IntraCGAll_InterCGNone = TRUE;
658 snew(bExcl, excl_nalloc);
661 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
663 molb = &mtop->molblock[mb];
664 molt = &mtop->moltype[molb->type];
668 /* Check if the cginfo is identical for all molecules in this block.
669 * If so, we only need an array of the size of one molecule.
670 * Otherwise we make an array of #mol times #cgs per molecule.
673 for (m = 0; m < molb->nmol; m++)
675 int am = m*cgs->index[cgs->nr];
676 for (cg = 0; cg < cgs->nr; cg++)
679 a1 = cgs->index[cg+1];
680 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
681 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
685 if (mtop->groups.grpnr[egcQMMM] != nullptr)
687 for (ai = a0; ai < a1; ai++)
689 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
690 mtop->groups.grpnr[egcQMMM][a_offset +ai])
699 cginfo_mb[mb].cg_start = cg_offset;
700 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
701 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
702 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
703 cginfo = cginfo_mb[mb].cginfo;
705 /* Set constraints flags for constrained atoms */
706 snew(a_con, molt->atoms.nr);
707 for (ftype = 0; ftype < F_NRE; ftype++)
709 if (interaction_function[ftype].flags & IF_CONSTRAINT)
714 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
718 for (a = 0; a < nral; a++)
720 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
721 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
727 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
730 int am = m*cgs->index[cgs->nr];
731 for (cg = 0; cg < cgs->nr; cg++)
734 a1 = cgs->index[cg+1];
736 /* Store the energy group in cginfo */
737 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
738 SET_CGINFO_GID(cginfo[cgm+cg], gid);
740 /* Check the intra/inter charge group exclusions */
741 if (a1-a0 > excl_nalloc)
743 excl_nalloc = a1 - a0;
744 srenew(bExcl, excl_nalloc);
746 /* bExclIntraAll: all intra cg interactions excluded
747 * bExclInter: any inter cg interactions excluded
749 bExclIntraAll = TRUE;
753 bHavePerturbedAtoms = FALSE;
754 for (ai = a0; ai < a1; ai++)
756 /* Check VDW and electrostatic interactions */
757 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
758 type_VDW[molt->atoms.atom[ai].typeB]);
759 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
760 molt->atoms.atom[ai].qB != 0);
762 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
764 /* Clear the exclusion list for atom ai */
765 for (aj = a0; aj < a1; aj++)
767 bExcl[aj-a0] = FALSE;
769 /* Loop over all the exclusions of atom ai */
770 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
773 if (aj < a0 || aj >= a1)
782 /* Check if ai excludes a0 to a1 */
783 for (aj = a0; aj < a1; aj++)
787 bExclIntraAll = FALSE;
794 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
797 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
805 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
809 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
811 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
813 /* The size in cginfo is currently only read with DD */
814 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
818 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
822 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
824 if (bHavePerturbedAtoms && fr->efep != efepNO)
826 SET_CGINFO_FEP(cginfo[cgm+cg]);
827 *bFEP_NonBonded = TRUE;
829 /* Store the charge group size */
830 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
832 if (!bExclIntraAll || bExclInter)
834 *bExcl_IntraCGAll_InterCGNone = FALSE;
841 cg_offset += molb->nmol*cgs->nr;
842 a_offset += molb->nmol*cgs->index[cgs->nr];
847 /* the solvent optimizer is called after the QM is initialized,
848 * because we don't want to have the QM subsystemto become an
852 check_solvent(fplog, mtop, fr, cginfo_mb);
854 if (getenv("GMX_NO_SOLV_OPT"))
858 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
859 "Disabling all solvent optimization\n");
861 fr->solvent_opt = esolNO;
865 fr->solvent_opt = esolNO;
867 if (!fr->solvent_opt)
869 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
871 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
873 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
881 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
886 ncg = cgi_mb[nmb-1].cg_end;
889 for (cg = 0; cg < ncg; cg++)
891 while (cg >= cgi_mb[mb].cg_end)
896 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
902 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
904 if (cginfo_mb == nullptr)
908 for (int mb = 0; mb < numMolBlocks; ++mb)
910 sfree(cginfo_mb[mb].cginfo);
915 /* Sets the sum of charges (squared) and C6 in the system in fr.
916 * Returns whether the system has a net charge.
918 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
920 /*This now calculates sum for q and c6*/
921 double qsum, q2sum, q, c6sum, c6;
926 for (const gmx_molblock_t &molb : mtop->molblock)
928 int nmol = molb.nmol;
929 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
930 for (int i = 0; i < atoms->nr; i++)
932 q = atoms->atom[i].q;
935 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
940 fr->q2sum[0] = q2sum;
941 fr->c6sum[0] = c6sum;
943 if (fr->efep != efepNO)
948 for (const gmx_molblock_t &molb : mtop->molblock)
950 int nmol = molb.nmol;
951 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
952 for (int i = 0; i < atoms->nr; i++)
954 q = atoms->atom[i].qB;
957 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
961 fr->q2sum[1] = q2sum;
962 fr->c6sum[1] = c6sum;
967 fr->qsum[1] = fr->qsum[0];
968 fr->q2sum[1] = fr->q2sum[0];
969 fr->c6sum[1] = fr->c6sum[0];
973 if (fr->efep == efepNO)
975 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
979 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
980 fr->qsum[0], fr->qsum[1]);
984 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
985 return (std::abs(fr->qsum[0]) > 1e-4 ||
986 std::abs(fr->qsum[1]) > 1e-4);
989 void update_forcerec(t_forcerec *fr, matrix box)
991 if (fr->ic->eeltype == eelGRF)
993 calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
994 fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
995 &fr->ic->k_rf, &fr->ic->c_rf);
999 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
1001 const t_atoms *atoms, *atoms_tpi;
1002 const t_blocka *excl;
1003 int nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
1004 gmx_int64_t npair, npair_ij, tmpi, tmpj;
1005 double csix, ctwelve;
1006 int ntp, *typecount;
1009 real *nbfp_comb = nullptr;
1015 /* For LJ-PME, we want to correct for the difference between the
1016 * actual C6 values and the C6 values used by the LJ-PME based on
1017 * combination rules. */
1019 if (EVDW_PME(fr->ic->vdwtype))
1021 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1022 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1023 for (tpi = 0; tpi < ntp; ++tpi)
1025 for (tpj = 0; tpj < ntp; ++tpj)
1027 C6(nbfp_comb, ntp, tpi, tpj) =
1028 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1029 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1034 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1042 /* Count the types so we avoid natoms^2 operations */
1043 snew(typecount, ntp);
1044 gmx_mtop_count_atomtypes(mtop, q, typecount);
1046 for (tpi = 0; tpi < ntp; tpi++)
1048 for (tpj = tpi; tpj < ntp; tpj++)
1050 tmpi = typecount[tpi];
1051 tmpj = typecount[tpj];
1054 npair_ij = tmpi*tmpj;
1058 npair_ij = tmpi*(tmpi - 1)/2;
1062 /* nbfp now includes the 6.0 derivative prefactor */
1063 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1067 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1068 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1069 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1075 /* Subtract the excluded pairs.
1076 * The main reason for substracting exclusions is that in some cases
1077 * some combinations might never occur and the parameters could have
1078 * any value. These unused values should not influence the dispersion
1081 for (const gmx_molblock_t &molb : mtop->molblock)
1083 int nmol = molb.nmol;
1084 atoms = &mtop->moltype[molb.type].atoms;
1085 excl = &mtop->moltype[molb.type].excls;
1086 for (int i = 0; (i < atoms->nr); i++)
1090 tpi = atoms->atom[i].type;
1094 tpi = atoms->atom[i].typeB;
1096 j1 = excl->index[i];
1097 j2 = excl->index[i+1];
1098 for (j = j1; j < j2; j++)
1105 tpj = atoms->atom[k].type;
1109 tpj = atoms->atom[k].typeB;
1113 /* nbfp now includes the 6.0 derivative prefactor */
1114 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1118 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1119 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1120 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1130 /* Only correct for the interaction of the test particle
1131 * with the rest of the system.
1134 &mtop->moltype[mtop->molblock.back().type].atoms;
1137 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
1139 const gmx_molblock_t &molb = mtop->molblock[mb];
1140 atoms = &mtop->moltype[molb.type].atoms;
1141 for (j = 0; j < atoms->nr; j++)
1144 /* Remove the interaction of the test charge group
1147 if (mb == mtop->molblock.size() - 1)
1151 if (mb == 0 && molb.nmol == 1)
1153 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1158 tpj = atoms->atom[j].type;
1162 tpj = atoms->atom[j].typeB;
1164 for (i = 0; i < fr->n_tpi; i++)
1168 tpi = atoms_tpi->atom[i].type;
1172 tpi = atoms_tpi->atom[i].typeB;
1176 /* nbfp now includes the 6.0 derivative prefactor */
1177 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1181 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1182 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1183 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1190 if (npair - nexcl <= 0 && fplog)
1192 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1198 csix /= npair - nexcl;
1199 ctwelve /= npair - nexcl;
1203 fprintf(debug, "Counted %d exclusions\n", nexcl);
1204 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1205 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1207 fr->avcsix[q] = csix;
1208 fr->avctwelve[q] = ctwelve;
1211 if (EVDW_PME(fr->ic->vdwtype))
1216 if (fplog != nullptr)
1218 if (fr->eDispCorr == edispcAllEner ||
1219 fr->eDispCorr == edispcAllEnerPres)
1221 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1222 fr->avcsix[0], fr->avctwelve[0]);
1226 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1232 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
1234 const t_atoms *at1, *at2;
1235 int i, j, tpi, tpj, ntypes;
1240 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1242 ntypes = mtop->ffparams.atnr;
1245 real bham_b_max = 0;
1246 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
1248 at1 = &mtop->moltype[mt1].atoms;
1249 for (i = 0; (i < at1->nr); i++)
1251 tpi = at1->atom[i].type;
1254 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1257 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
1259 at2 = &mtop->moltype[mt2].atoms;
1260 for (j = 0; (j < at2->nr); j++)
1262 tpj = at2->atom[j].type;
1265 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1267 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
1272 if ((b < bmin) || (bmin == -1))
1282 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1289 static void make_nbf_tables(FILE *fp,
1290 const interaction_const_t *ic, real rtab,
1291 const char *tabfn, char *eg1, char *eg2,
1297 if (tabfn == nullptr)
1301 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1306 sprintf(buf, "%s", tabfn);
1309 /* Append the two energy group names */
1310 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1311 eg1, eg2, ftp2ext(efXVG));
1313 nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1314 /* Copy the contents of the table to separate coulomb and LJ tables too,
1315 * to improve cache performance.
1317 /* For performance reasons we want
1318 * the table data to be aligned to 16-byte. The pointers could be freed
1319 * but currently aren't.
1321 snew(nbl->table_elec, 1);
1322 nbl->table_elec->interaction = GMX_TABLE_INTERACTION_ELEC;
1323 nbl->table_elec->format = nbl->table_elec_vdw->format;
1324 nbl->table_elec->r = nbl->table_elec_vdw->r;
1325 nbl->table_elec->n = nbl->table_elec_vdw->n;
1326 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1327 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1328 nbl->table_elec->ninteractions = 1;
1329 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1330 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1332 snew(nbl->table_vdw, 1);
1333 nbl->table_vdw->interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1334 nbl->table_vdw->format = nbl->table_elec_vdw->format;
1335 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1336 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1337 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1338 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1339 nbl->table_vdw->ninteractions = 2;
1340 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1341 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1343 /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1344 * with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1346 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1348 for (j = 0; j < 4; j++)
1350 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1353 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1355 for (j = 0; j < 8; j++)
1357 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1362 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1363 * ftype2 present in the topology, build an array of the number of
1364 * interactions present for each bonded interaction index found in the
1367 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1368 * valid type with that parameter.
1370 * \c count will be reallocated as necessary to fit the largest bonded
1371 * interaction index found, and its current size will be returned in
1372 * \c ncount. It will contain zero for every bonded interaction index
1373 * for which no interactions are present in the topology.
1375 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1376 int *ncount, int **count)
1379 int ftype, stride, i, j, tabnr;
1381 // Loop over all moleculetypes
1382 for (const gmx_moltype_t &molt : mtop->moltype)
1384 // Loop over all interaction types
1385 for (ftype = 0; ftype < F_NRE; ftype++)
1387 // If the current interaction type is one of the types whose tables we're trying to count...
1388 if (ftype == ftype1 || ftype == ftype2)
1390 il = &molt.ilist[ftype];
1391 stride = 1 + NRAL(ftype);
1392 // ... and there are actually some interactions for this type
1393 for (i = 0; i < il->nr; i += stride)
1395 // Find out which table index the user wanted
1396 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1399 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1401 // Make room for this index in the data structure
1402 if (tabnr >= *ncount)
1404 srenew(*count, tabnr+1);
1405 for (j = *ncount; j < tabnr+1; j++)
1411 // Record that this table index is used and must have a valid file
1419 /*!\brief If there's bonded interactions of flavour \c tabext and type
1420 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1421 * list of filenames passed to mdrun, and make bonded tables from
1424 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1425 * valid type with that parameter.
1427 * A fatal error occurs if no matching filename is found.
1429 static bondedtable_t *make_bonded_tables(FILE *fplog,
1430 int ftype1, int ftype2,
1431 const gmx_mtop_t *mtop,
1432 gmx::ArrayRef<const std::string> tabbfnm,
1442 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1444 // Are there any relevant tabulated bond interactions?
1448 for (int i = 0; i < ncount; i++)
1450 // Do any interactions exist that requires this table?
1453 // This pattern enforces the current requirement that
1454 // table filenames end in a characteristic sequence
1455 // before the file type extension, and avoids table 13
1456 // being recognized and used for table 1.
1457 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1458 bool madeTable = false;
1459 for (size_t j = 0; j < tabbfnm.size() && !madeTable; ++j)
1461 if (gmx::endsWith(tabbfnm[j].c_str(), patternToFind))
1463 // Finally read the table from the file found
1464 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1470 bool isPlural = (ftype2 != -1);
1471 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.",
1472 interaction_function[ftype1].longname,
1473 isPlural ? "' or '" : "",
1474 isPlural ? interaction_function[ftype2].longname : "",
1476 patternToFind.c_str());
1486 void forcerec_set_ranges(t_forcerec *fr,
1487 int ncg_home, int ncg_force,
1489 int natoms_force_constr, int natoms_f_novirsum)
1494 /* fr->ncg_force is unused in the standard code,
1495 * but it can be useful for modified code dealing with charge groups.
1497 fr->ncg_force = ncg_force;
1498 fr->natoms_force = natoms_force;
1499 fr->natoms_force_constr = natoms_force_constr;
1501 if (fr->natoms_force_constr > fr->nalloc_force)
1503 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1506 if (fr->haveDirectVirialContributions)
1508 fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1512 static real cutoff_inf(real cutoff)
1516 cutoff = GMX_CUTOFF_INF;
1522 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, const t_commrec *cr, FILE *fp)
1529 ir->rcoulomb == 0 &&
1531 ir->ePBC == epbcNONE &&
1532 ir->vdwtype == evdwCUT &&
1533 ir->coulombtype == eelCUT &&
1534 ir->efep == efepNO &&
1535 getenv("GMX_NO_ALLVSALL") == nullptr
1538 if (bAllvsAll && ir->opts.ngener > 1)
1540 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";
1546 fprintf(fp, "\n%s\n", note);
1552 if (bAllvsAll && fp && MASTER(cr))
1554 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1561 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
1562 const t_inputrec *ir)
1564 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1566 /* LJ PME with LB combination rule does 7 mesh operations.
1567 * This so slow that we don't compile SIMD non-bonded kernels
1569 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
1577 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1580 const gmx_hw_info_t gmx_unused &hardwareInfo)
1582 *kernel_type = nbnxnk4x4_PlainC;
1583 *ewald_excl = ewaldexclTable;
1587 #ifdef GMX_NBNXN_SIMD_4XN
1588 *kernel_type = nbnxnk4xN_SIMD_4xN;
1590 #ifdef GMX_NBNXN_SIMD_2XNN
1591 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1594 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1595 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1596 * This is based on the SIMD acceleration choice and CPU information
1597 * detected at runtime.
1599 * 4xN calculates more (zero) interactions, but has less pair-search
1600 * work and much better kernel instruction scheduling.
1602 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1603 * which doesn't have FMA, both the analytical and tabulated Ewald
1604 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1605 * 2x(4+4) because it results in significantly fewer pairs.
1606 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1607 * 10% with HT, 50% without HT. As we currently don't detect the actual
1608 * use of HT, use 4x8 to avoid a potential performance hit.
1609 * On Intel Haswell 4x8 is always faster.
1611 *kernel_type = nbnxnk4xN_SIMD_4xN;
1613 #if !GMX_SIMD_HAVE_FMA
1614 if (EEL_PME_EWALD(ir->coulombtype) ||
1615 EVDW_PME(ir->vdwtype))
1617 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1618 * There are enough instructions to make 2x(4+4) efficient.
1620 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1623 if (hardwareInfo.haveAmdZenCpu)
1625 /* One 256-bit FMA per cycle makes 2xNN faster */
1626 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1628 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1631 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
1633 #ifdef GMX_NBNXN_SIMD_4XN
1634 *kernel_type = nbnxnk4xN_SIMD_4xN;
1636 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1639 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
1641 #ifdef GMX_NBNXN_SIMD_2XNN
1642 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1644 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1648 /* Analytical Ewald exclusion correction is only an option in
1650 * Since table lookup's don't parallelize with SIMD, analytical
1651 * will probably always be faster for a SIMD width of 8 or more.
1652 * With FMA analytical is sometimes faster for a width if 4 as well.
1653 * In single precision, this is faster on Bulldozer.
1655 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1656 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
1657 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
1658 * of single or double precision and 128 or 256-bit AVX2.
1660 if (!hardwareInfo.haveAmdZenCpu)
1662 *ewald_excl = ewaldexclAnalytical;
1665 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
1667 *ewald_excl = ewaldexclTable;
1669 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
1671 *ewald_excl = ewaldexclAnalytical;
1679 const char *lookup_nbnxn_kernel_name(int kernel_type)
1681 const char *returnvalue = nullptr;
1682 switch (kernel_type)
1685 returnvalue = "not set";
1687 case nbnxnk4x4_PlainC:
1688 returnvalue = "plain C";
1690 case nbnxnk4xN_SIMD_4xN:
1691 case nbnxnk4xN_SIMD_2xNN:
1693 returnvalue = "SIMD";
1695 returnvalue = "not available";
1698 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1699 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1703 gmx_fatal(FARGS, "Illegal kernel type selected");
1704 returnvalue = nullptr;
1710 static void pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
1711 gmx_bool use_simd_kernels,
1712 const gmx_hw_info_t &hardwareInfo,
1714 EmulateGpuNonbonded emulateGpu,
1715 const t_inputrec *ir,
1718 gmx_bool bDoNonbonded)
1720 assert(kernel_type);
1722 *kernel_type = nbnxnkNotSet;
1723 *ewald_excl = ewaldexclTable;
1725 if (emulateGpu == EmulateGpuNonbonded::Yes)
1727 *kernel_type = nbnxnk8x8x8_PlainC;
1731 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
1736 *kernel_type = nbnxnk8x8x8_GPU;
1739 if (*kernel_type == nbnxnkNotSet)
1741 if (use_simd_kernels &&
1742 nbnxn_simd_supported(mdlog, ir))
1744 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl, hardwareInfo);
1748 *kernel_type = nbnxnk4x4_PlainC;
1754 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
1755 "Using %s %dx%d nonbonded short-range kernels",
1756 lookup_nbnxn_kernel_name(*kernel_type),
1757 nbnxn_kernel_to_cluster_i_size(*kernel_type),
1758 nbnxn_kernel_to_cluster_j_size(*kernel_type));
1760 if (nbnxnk4x4_PlainC == *kernel_type ||
1761 nbnxnk8x8x8_PlainC == *kernel_type)
1763 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1764 "WARNING: Using the slow %s kernels. This should\n"
1765 "not happen during routine usage on supported platforms.",
1766 lookup_nbnxn_kernel_name(*kernel_type));
1771 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1772 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1773 bool systemHasNetCharge,
1774 interaction_const_t *ic)
1776 if (!EEL_PME_EWALD(ir->coulombtype))
1783 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1785 if (ir->coulombtype == eelP3M_AD)
1787 please_cite(fp, "Hockney1988");
1788 please_cite(fp, "Ballenegger2012");
1792 please_cite(fp, "Essmann95a");
1795 if (ir->ewald_geometry == eewg3DC)
1799 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1800 systemHasNetCharge ? " and net charge" : "");
1802 please_cite(fp, "In-Chul99a");
1803 if (systemHasNetCharge)
1805 please_cite(fp, "Ballenegger2009");
1810 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1813 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1814 1/ic->ewaldcoeff_q);
1817 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1819 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1820 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1828 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1829 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1830 interaction_const_t *ic)
1832 if (!EVDW_PME(ir->vdwtype))
1839 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1840 please_cite(fp, "Essmann95a");
1842 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1845 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1846 1/ic->ewaldcoeff_lj);
1849 if (ic->vdw_modifier == eintmodPOTSHIFT)
1851 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1852 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1856 ic->sh_lj_ewald = 0;
1860 gmx_bool uses_simple_tables(int cutoff_scheme,
1861 nonbonded_verlet_t *nbv,
1864 gmx_bool bUsesSimpleTables = TRUE;
1867 switch (cutoff_scheme)
1870 bUsesSimpleTables = TRUE;
1873 assert(NULL != nbv && NULL != nbv->grp);
1874 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1875 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1878 gmx_incons("unimplemented");
1880 return bUsesSimpleTables;
1883 static void init_ewald_f_table(interaction_const_t *ic,
1888 /* Get the Ewald table spacing based on Coulomb and/or LJ
1889 * Ewald coefficients and rtol.
1891 ic->tabq_scale = ewald_spline3_table_scale(ic);
1893 if (ic->cutoff_scheme == ecutsVERLET)
1895 maxr = ic->rcoulomb;
1899 maxr = std::max(ic->rcoulomb, rtab);
1901 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1903 sfree_aligned(ic->tabq_coul_FDV0);
1904 sfree_aligned(ic->tabq_coul_F);
1905 sfree_aligned(ic->tabq_coul_V);
1907 sfree_aligned(ic->tabq_vdw_FDV0);
1908 sfree_aligned(ic->tabq_vdw_F);
1909 sfree_aligned(ic->tabq_vdw_V);
1911 if (EEL_PME_EWALD(ic->eeltype))
1913 /* Create the original table data in FDV0 */
1914 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1915 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1916 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1917 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1918 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1921 if (EVDW_PME(ic->vdwtype))
1923 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1924 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1925 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1926 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1927 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1931 void init_interaction_const_tables(FILE *fp,
1932 interaction_const_t *ic,
1935 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1937 init_ewald_f_table(ic, rtab);
1941 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1942 1/ic->tabq_scale, ic->tabq_size);
1947 static void clear_force_switch_constants(shift_consts_t *sc)
1954 static void force_switch_constants(real p,
1958 /* Here we determine the coefficient for shifting the force to zero
1959 * between distance rsw and the cut-off rc.
1960 * For a potential of r^-p, we have force p*r^-(p+1).
1961 * But to save flops we absorb p in the coefficient.
1963 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1964 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1966 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1967 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1968 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1971 static void potential_switch_constants(real rsw, real rc,
1972 switch_consts_t *sc)
1974 /* The switch function is 1 at rsw and 0 at rc.
1975 * The derivative and second derivate are zero at both ends.
1976 * rsw = max(r - r_switch, 0)
1977 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1978 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1979 * force = force*dsw - potential*sw
1982 sc->c3 = -10/gmx::power3(rc - rsw);
1983 sc->c4 = 15/gmx::power4(rc - rsw);
1984 sc->c5 = -6/gmx::power5(rc - rsw);
1987 /*! \brief Construct interaction constants
1989 * This data is used (particularly) by search and force code for
1990 * short-range interactions. Many of these are constant for the whole
1991 * simulation; some are constant only after PME tuning completes.
1994 init_interaction_const(FILE *fp,
1995 interaction_const_t **interaction_const,
1996 const t_inputrec *ir,
1997 const gmx_mtop_t *mtop,
1998 bool systemHasNetCharge)
2000 interaction_const_t *ic;
2004 ic->cutoff_scheme = ir->cutoff_scheme;
2006 /* Just allocate something so we can free it */
2007 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
2008 snew_aligned(ic->tabq_coul_F, 16, 32);
2009 snew_aligned(ic->tabq_coul_V, 16, 32);
2012 ic->vdwtype = ir->vdwtype;
2013 ic->vdw_modifier = ir->vdw_modifier;
2014 ic->reppow = mtop->ffparams.reppow;
2015 ic->rvdw = cutoff_inf(ir->rvdw);
2016 ic->rvdw_switch = ir->rvdw_switch;
2017 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
2018 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
2019 if (ic->useBuckingham)
2021 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
2024 initVdwEwaldParameters(fp, ir, ic);
2026 clear_force_switch_constants(&ic->dispersion_shift);
2027 clear_force_switch_constants(&ic->repulsion_shift);
2029 switch (ic->vdw_modifier)
2031 case eintmodPOTSHIFT:
2032 /* Only shift the potential, don't touch the force */
2033 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
2034 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
2036 case eintmodFORCESWITCH:
2037 /* Switch the force, switch and shift the potential */
2038 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2039 &ic->dispersion_shift);
2040 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2041 &ic->repulsion_shift);
2043 case eintmodPOTSWITCH:
2044 /* Switch the potential and force */
2045 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2049 case eintmodEXACTCUTOFF:
2050 /* Nothing to do here */
2053 gmx_incons("unimplemented potential modifier");
2056 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2058 /* Electrostatics */
2059 ic->eeltype = ir->coulombtype;
2060 ic->coulomb_modifier = ir->coulomb_modifier;
2061 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
2062 ic->rcoulomb_switch = ir->rcoulomb_switch;
2063 ic->epsilon_r = ir->epsilon_r;
2065 /* Set the Coulomb energy conversion factor */
2066 if (ic->epsilon_r != 0)
2068 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
2072 /* eps = 0 is infinite dieletric: no Coulomb interactions */
2076 /* Reaction-field */
2077 if (EEL_RF(ic->eeltype))
2079 ic->epsilon_rf = ir->epsilon_rf;
2080 /* Generalized reaction field parameters are updated every step */
2081 if (ic->eeltype != eelGRF)
2083 calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
2084 ic->rcoulomb, 0, 0, NULL,
2085 &ic->k_rf, &ic->c_rf);
2088 if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
2090 /* grompp should have done this, but this scheme is obsolete */
2091 ic->coulomb_modifier = eintmodEXACTCUTOFF;
2096 /* For plain cut-off we might use the reaction-field kernels */
2097 ic->epsilon_rf = ic->epsilon_r;
2099 if (ir->coulomb_modifier == eintmodPOTSHIFT)
2101 ic->c_rf = 1/ic->rcoulomb;
2109 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
2113 real dispersion_shift;
2115 dispersion_shift = ic->dispersion_shift.cpot;
2116 if (EVDW_PME(ic->vdwtype))
2118 dispersion_shift -= ic->sh_lj_ewald;
2120 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2121 ic->repulsion_shift.cpot, dispersion_shift);
2123 if (ic->eeltype == eelCUT)
2125 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2127 else if (EEL_PME(ic->eeltype))
2129 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2134 *interaction_const = ic;
2138 done_interaction_const(interaction_const_t *interaction_const)
2140 sfree_aligned(interaction_const->tabq_coul_FDV0);
2141 sfree_aligned(interaction_const->tabq_coul_F);
2142 sfree_aligned(interaction_const->tabq_coul_V);
2143 sfree(interaction_const);
2146 static void init_nb_verlet(const gmx::MDLogger &mdlog,
2147 nonbonded_verlet_t **nb_verlet,
2148 gmx_bool bFEP_NonBonded,
2149 const t_inputrec *ir,
2150 const t_forcerec *fr,
2151 const t_commrec *cr,
2152 const gmx_hw_info_t &hardwareInfo,
2153 const gmx_device_info_t *deviceInfo,
2154 const gmx_mtop_t *mtop,
2157 nonbonded_verlet_t *nbv;
2160 nbnxn_alloc_t *nb_alloc;
2161 nbnxn_free_t *nb_free;
2163 nbv = new nonbonded_verlet_t();
2165 nbv->emulateGpu = ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
2166 nbv->bUseGPU = deviceInfo != nullptr;
2168 GMX_RELEASE_ASSERT(!(nbv->emulateGpu == EmulateGpuNonbonded::Yes && nbv->bUseGPU), "When GPU emulation is active, there cannot be a GPU assignment");
2171 nbv->min_ci_balanced = 0;
2173 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2174 for (int i = 0; i < nbv->ngrp; i++)
2176 nbv->grp[i].nbl_lists.nnbl = 0;
2177 nbv->grp[i].kernel_type = nbnxnkNotSet;
2179 if (i == 0) /* local */
2181 pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
2182 nbv->bUseGPU, nbv->emulateGpu, ir,
2183 &nbv->grp[i].kernel_type,
2184 &nbv->grp[i].ewald_excl,
2187 else /* non-local */
2189 /* Use the same kernel for local and non-local interactions */
2190 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2191 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2195 nbv->listParams = std::unique_ptr<NbnxnListParameters>(new NbnxnListParameters(ir->rlist));
2196 setupDynamicPairlistPruning(mdlog, ir, mtop, box, nbv->grp[0].kernel_type, fr->ic,
2197 nbv->listParams.get());
2199 nbv->nbs = gmx::compat::make_unique<nbnxn_search>(DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
2200 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
2202 gmx_omp_nthreads_get(emntPairsearch));
2204 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2205 &nb_alloc, &nb_free);
2207 for (int i = 0; i < nbv->ngrp; i++)
2209 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2210 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2211 /* 8x8x8 "non-simple" lists are ATM always combined */
2212 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2216 int enbnxninitcombrule;
2217 if (fr->ic->vdwtype == evdwCUT &&
2218 (fr->ic->vdw_modifier == eintmodNONE ||
2219 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
2220 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
2222 /* Plain LJ cut-off: we can optimize with combination rules */
2223 enbnxninitcombrule = enbnxninitcombruleDETECT;
2225 else if (fr->ic->vdwtype == evdwPME)
2227 /* LJ-PME: we need to use a combination rule for the grid */
2228 if (fr->ljpme_combination_rule == eljpmeGEOM)
2230 enbnxninitcombrule = enbnxninitcombruleGEOM;
2234 enbnxninitcombrule = enbnxninitcombruleLB;
2239 /* We use a full combination matrix: no rule required */
2240 enbnxninitcombrule = enbnxninitcombruleNONE;
2244 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
2245 if (ir->opts.ngener - ir->nwall == 1)
2247 /* We have only one non-wall energy group, we do not need energy group
2248 * support in the non-bondeds kernels, since all non-bonded energy
2249 * contributions go to the first element of the energy group matrix.
2251 mimimumNumEnergyGroupNonbonded = 1;
2253 bool bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[0].kernel_type);
2254 nbnxn_atomdata_init(mdlog,
2256 nbv->grp[0].kernel_type,
2258 fr->ntype, fr->nbfp,
2259 mimimumNumEnergyGroupNonbonded,
2260 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2265 /* init the NxN GPU data; the last argument tells whether we'll have
2266 * both local and non-local NB calculation on GPU */
2267 nbnxn_gpu_init(&nbv->gpu_nbv,
2270 nbv->listParams.get(),
2275 if ((env = getenv("GMX_NB_MIN_CI")) != nullptr)
2279 nbv->min_ci_balanced = strtol(env, &end, 10);
2280 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
2282 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
2287 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2288 nbv->min_ci_balanced);
2293 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2296 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2297 nbv->min_ci_balanced);
2306 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2308 return nbv != nullptr && nbv->bUseGPU;
2311 void init_forcerec(FILE *fp,
2312 const gmx::MDLogger &mdlog,
2315 const t_inputrec *ir,
2316 const gmx_mtop_t *mtop,
2317 const t_commrec *cr,
2321 gmx::ArrayRef<const std::string> tabbfnm,
2322 const gmx_hw_info_t &hardwareInfo,
2323 const gmx_device_info_t *deviceInfo,
2324 gmx_bool bNoSolvOpt,
2327 int m, negp_pp, negptable, egi, egj;
2332 gmx_bool bGenericKernelOnly;
2333 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
2334 gmx_bool bFEP_NonBonded;
2335 int *nm_ind, egp_flags;
2337 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2338 fr->use_simd_kernels = TRUE;
2340 fr->bDomDec = DOMAINDECOMP(cr);
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::BlockRanges molecules = gmx_mtop_molecules(*mtop);
2357 if (fr->n_tpi != molecules.index[molecules.numBlocks()] - molecules.index[molecules.numBlocks() - 1])
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->bonded_threading);
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);