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/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/ewald/ewald.h"
55 #include "gromacs/ewald/ewald-utils.h"
56 #include "gromacs/fileio/filetypes.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/hardware/hw_info.h"
61 #include "gromacs/listed-forces/manage-threading.h"
62 #include "gromacs/listed-forces/pairs.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/units.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/force.h"
68 #include "gromacs/mdlib/forcerec-threading.h"
69 #include "gromacs/mdlib/gmx_omp_nthreads.h"
70 #include "gromacs/mdlib/md_support.h"
71 #include "gromacs/mdlib/nb_verlet.h"
72 #include "gromacs/mdlib/nbnxn_atomdata.h"
73 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
74 #include "gromacs/mdlib/nbnxn_search.h"
75 #include "gromacs/mdlib/nbnxn_simd.h"
76 #include "gromacs/mdlib/nbnxn_tuning.h"
77 #include "gromacs/mdlib/nbnxn_util.h"
78 #include "gromacs/mdlib/ns.h"
79 #include "gromacs/mdlib/qmmm.h"
80 #include "gromacs/mdlib/sim_util.h"
81 #include "gromacs/mdtypes/commrec.h"
82 #include "gromacs/mdtypes/fcdata.h"
83 #include "gromacs/mdtypes/group.h"
84 #include "gromacs/mdtypes/iforceprovider.h"
85 #include "gromacs/mdtypes/inputrec.h"
86 #include "gromacs/mdtypes/md_enums.h"
87 #include "gromacs/pbcutil/ishift.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/simd/simd.h"
90 #include "gromacs/tables/forcetable.h"
91 #include "gromacs/topology/mtop_util.h"
92 #include "gromacs/trajectory/trajectoryframe.h"
93 #include "gromacs/utility/cstringutil.h"
94 #include "gromacs/utility/exceptions.h"
95 #include "gromacs/utility/fatalerror.h"
96 #include "gromacs/utility/gmxassert.h"
97 #include "gromacs/utility/logger.h"
98 #include "gromacs/utility/physicalnodecommunicator.h"
99 #include "gromacs/utility/pleasecite.h"
100 #include "gromacs/utility/smalloc.h"
101 #include "gromacs/utility/strconvert.h"
103 #include "nbnxn_gpu_jit_support.h"
105 const char *egrp_nm[egNR+1] = {
106 "Coul-SR", "LJ-SR", "Buck-SR",
107 "Coul-14", "LJ-14", nullptr
110 t_forcerec *mk_forcerec(void)
120 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
124 for (i = 0; (i < atnr); i++)
126 for (j = 0; (j < atnr); j++)
128 fprintf(fp, "%2d - %2d", i, j);
131 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
132 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
136 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
137 C12(nbfp, atnr, i, j)/12.0);
144 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
152 snew(nbfp, 3*atnr*atnr);
153 for (i = k = 0; (i < atnr); i++)
155 for (j = 0; (j < atnr); j++, k++)
157 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
158 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
159 /* nbfp now includes the 6.0 derivative prefactor */
160 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
166 snew(nbfp, 2*atnr*atnr);
167 for (i = k = 0; (i < atnr); i++)
169 for (j = 0; (j < atnr); j++, k++)
171 /* nbfp now includes the 6.0/12.0 derivative prefactors */
172 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
173 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
181 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
184 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
187 /* For LJ-PME simulations, we correct the energies with the reciprocal space
188 * inside of the cut-off. To do this the non-bonded kernels needs to have
189 * access to the C6-values used on the reciprocal grid in pme.c
193 snew(grid, 2*atnr*atnr);
194 for (i = k = 0; (i < atnr); i++)
196 for (j = 0; (j < atnr); j++, k++)
198 c6i = idef->iparams[i*(atnr+1)].lj.c6;
199 c12i = idef->iparams[i*(atnr+1)].lj.c12;
200 c6j = idef->iparams[j*(atnr+1)].lj.c6;
201 c12j = idef->iparams[j*(atnr+1)].lj.c12;
202 c6 = std::sqrt(c6i * c6j);
203 if (fr->ljpme_combination_rule == eljpmeLB
204 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
206 sigmai = gmx::sixthroot(c12i / c6i);
207 sigmaj = gmx::sixthroot(c12j / c6j);
208 epsi = c6i * c6i / c12i;
209 epsj = c6j * c6j / c12j;
210 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
212 /* Store the elements at the same relative positions as C6 in nbfp in order
213 * to simplify access in the kernels
215 grid[2*(atnr*i+j)] = c6*6.0;
221 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
225 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
229 snew(nbfp, 2*atnr*atnr);
230 for (i = 0; i < atnr; ++i)
232 for (j = 0; j < atnr; ++j)
234 c6i = idef->iparams[i*(atnr+1)].lj.c6;
235 c12i = idef->iparams[i*(atnr+1)].lj.c12;
236 c6j = idef->iparams[j*(atnr+1)].lj.c6;
237 c12j = idef->iparams[j*(atnr+1)].lj.c12;
238 c6 = std::sqrt(c6i * c6j);
239 c12 = std::sqrt(c12i * c12j);
240 if (comb_rule == eCOMB_ARITHMETIC
241 && !gmx_numzero(c6) && !gmx_numzero(c12))
243 sigmai = gmx::sixthroot(c12i / c6i);
244 sigmaj = gmx::sixthroot(c12j / c6j);
245 epsi = c6i * c6i / c12i;
246 epsj = c6j * c6j / c12j;
247 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
248 c12 = std::sqrt(epsi * epsj) * gmx::power12(0.5*(sigmai+sigmaj));
250 C6(nbfp, atnr, i, j) = c6*6.0;
251 C12(nbfp, atnr, i, j) = c12*12.0;
257 /* This routine sets fr->solvent_opt to the most common solvent in the
258 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
259 * the fr->solvent_type array with the correct type (or esolNO).
261 * Charge groups that fulfill the conditions but are not identical to the
262 * most common one will be marked as esolNO in the solvent_type array.
264 * TIP3p is identical to SPC for these purposes, so we call it
265 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
267 * NOTE: QM particle should not
268 * become an optimized solvent. Not even if there is only one charge
278 } solvent_parameters_t;
281 check_solvent_cg(const gmx_moltype_t *molt,
284 const unsigned char *qm_grpnr,
285 const t_grps *qm_grps,
287 int *n_solvent_parameters,
288 solvent_parameters_t **solvent_parameters_p,
298 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
299 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
302 solvent_parameters_t *solvent_parameters;
304 /* We use a list with parameters for each solvent type.
305 * Every time we discover a new molecule that fulfills the basic
306 * conditions for a solvent we compare with the previous entries
307 * in these lists. If the parameters are the same we just increment
308 * the counter for that type, and otherwise we create a new type
309 * based on the current molecule.
311 * Once we've finished going through all molecules we check which
312 * solvent is most common, and mark all those molecules while we
313 * clear the flag on all others.
316 solvent_parameters = *solvent_parameters_p;
318 /* Mark the cg first as non optimized */
321 /* Check if this cg has no exclusions with atoms in other charge groups
322 * and all atoms inside the charge group excluded.
323 * We only have 3 or 4 atom solvent loops.
325 if (GET_CGINFO_EXCL_INTER(cginfo) ||
326 !GET_CGINFO_EXCL_INTRA(cginfo))
331 /* Get the indices of the first atom in this charge group */
332 j0 = molt->cgs.index[cg0];
333 j1 = molt->cgs.index[cg0+1];
335 /* Number of atoms in our molecule */
341 "Moltype '%s': there are %d atoms in this charge group\n",
345 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
348 if (nj < 3 || nj > 4)
353 /* Check if we are doing QM on this group */
355 if (qm_grpnr != nullptr)
357 for (j = j0; j < j1 && !qm; j++)
359 qm = (qm_grpnr[j] < qm_grps->nr - 1);
362 /* Cannot use solvent optimization with QM */
368 atom = molt->atoms.atom;
370 /* Still looks like a solvent, time to check parameters */
372 /* If it is perturbed (free energy) we can't use the solvent loops,
373 * so then we just skip to the next molecule.
377 for (j = j0; j < j1 && !perturbed; j++)
379 perturbed = PERTURBED(atom[j]);
387 /* Now it's only a question if the VdW and charge parameters
388 * are OK. Before doing the check we compare and see if they are
389 * identical to a possible previous solvent type.
390 * First we assign the current types and charges.
392 for (j = 0; j < nj; j++)
394 tmp_vdwtype[j] = atom[j0+j].type;
395 tmp_charge[j] = atom[j0+j].q;
398 /* Does it match any previous solvent type? */
399 for (k = 0; k < *n_solvent_parameters; k++)
404 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
405 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
406 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
411 /* Check that types & charges match for all atoms in molecule */
412 for (j = 0; j < nj && match == TRUE; j++)
414 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
418 if (tmp_charge[j] != solvent_parameters[k].charge[j])
425 /* Congratulations! We have a matched solvent.
426 * Flag it with this type for later processing.
429 solvent_parameters[k].count += nmol;
431 /* We are done with this charge group */
436 /* If we get here, we have a tentative new solvent type.
437 * Before we add it we must check that it fulfills the requirements
438 * of the solvent optimized loops. First determine which atoms have
441 for (j = 0; j < nj; j++)
444 tjA = tmp_vdwtype[j];
446 /* Go through all other tpes and see if any have non-zero
447 * VdW parameters when combined with this one.
449 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
451 /* We already checked that the atoms weren't perturbed,
452 * so we only need to check state A now.
456 has_vdw[j] = (has_vdw[j] ||
457 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
458 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
459 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
464 has_vdw[j] = (has_vdw[j] ||
465 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
466 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
471 /* Now we know all we need to make the final check and assignment. */
475 * For this we require thatn all atoms have charge,
476 * the charges on atom 2 & 3 should be the same, and only
477 * atom 1 might have VdW.
479 if (has_vdw[1] == FALSE &&
480 has_vdw[2] == FALSE &&
481 tmp_charge[0] != 0 &&
482 tmp_charge[1] != 0 &&
483 tmp_charge[2] == tmp_charge[1])
485 srenew(solvent_parameters, *n_solvent_parameters+1);
486 solvent_parameters[*n_solvent_parameters].model = esolSPC;
487 solvent_parameters[*n_solvent_parameters].count = nmol;
488 for (k = 0; k < 3; k++)
490 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
491 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
494 *cg_sp = *n_solvent_parameters;
495 (*n_solvent_parameters)++;
500 /* Or could it be a TIP4P?
501 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
502 * Only atom 1 mght have VdW.
504 if (has_vdw[1] == FALSE &&
505 has_vdw[2] == FALSE &&
506 has_vdw[3] == FALSE &&
507 tmp_charge[0] == 0 &&
508 tmp_charge[1] != 0 &&
509 tmp_charge[2] == tmp_charge[1] &&
512 srenew(solvent_parameters, *n_solvent_parameters+1);
513 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
514 solvent_parameters[*n_solvent_parameters].count = nmol;
515 for (k = 0; k < 4; k++)
517 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
518 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
521 *cg_sp = *n_solvent_parameters;
522 (*n_solvent_parameters)++;
526 *solvent_parameters_p = solvent_parameters;
530 check_solvent(FILE * fp,
531 const gmx_mtop_t * mtop,
533 cginfo_mb_t *cginfo_mb)
536 const gmx_moltype_t *molt;
537 int mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
538 int n_solvent_parameters;
539 solvent_parameters_t *solvent_parameters;
545 fprintf(debug, "Going to determine what solvent types we have.\n");
548 n_solvent_parameters = 0;
549 solvent_parameters = nullptr;
550 /* Allocate temporary array for solvent type */
551 snew(cg_sp, mtop->molblock.size());
554 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
556 molt = &mtop->moltype[mtop->molblock[mb].type];
558 /* Here we have to loop over all individual molecules
559 * because we need to check for QMMM particles.
561 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
562 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
563 nmol = mtop->molblock[mb].nmol/nmol_ch;
564 for (mol = 0; mol < nmol_ch; mol++)
567 am = mol*cgs->index[cgs->nr];
568 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
570 check_solvent_cg(molt, cg_mol, nmol,
571 mtop->groups.grpnr[egcQMMM] ?
572 mtop->groups.grpnr[egcQMMM]+at_offset+am : nullptr,
573 &mtop->groups.grps[egcQMMM],
575 &n_solvent_parameters, &solvent_parameters,
576 cginfo_mb[mb].cginfo[cgm+cg_mol],
577 &cg_sp[mb][cgm+cg_mol]);
580 at_offset += cgs->index[cgs->nr];
583 /* Puh! We finished going through all charge groups.
584 * Now find the most common solvent model.
587 /* Most common solvent this far */
589 for (i = 0; i < n_solvent_parameters; i++)
592 solvent_parameters[i].count > solvent_parameters[bestsp].count)
600 bestsol = solvent_parameters[bestsp].model;
608 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
610 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
611 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
612 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
614 if (cg_sp[mb][i] == bestsp)
616 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
621 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
628 if (bestsol != esolNO && fp != nullptr)
630 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
632 solvent_parameters[bestsp].count);
635 sfree(solvent_parameters);
636 fr->solvent_opt = bestsol;
640 acNONE = 0, acCONSTRAINT, acSETTLE
643 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
644 t_forcerec *fr, gmx_bool bNoSolvOpt,
645 gmx_bool *bFEP_NonBonded,
646 gmx_bool *bExcl_IntraCGAll_InterCGNone)
649 const t_blocka *excl;
650 const gmx_moltype_t *molt;
651 const gmx_molblock_t *molb;
652 cginfo_mb_t *cginfo_mb;
655 int cg_offset, a_offset;
656 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
660 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
662 snew(cginfo_mb, mtop->molblock.size());
664 snew(type_VDW, fr->ntype);
665 for (ai = 0; ai < fr->ntype; ai++)
667 type_VDW[ai] = FALSE;
668 for (j = 0; j < fr->ntype; j++)
670 type_VDW[ai] = type_VDW[ai] ||
672 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
673 C12(fr->nbfp, fr->ntype, ai, j) != 0;
677 *bFEP_NonBonded = FALSE;
678 *bExcl_IntraCGAll_InterCGNone = TRUE;
681 snew(bExcl, excl_nalloc);
684 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
686 molb = &mtop->molblock[mb];
687 molt = &mtop->moltype[molb->type];
691 /* Check if the cginfo is identical for all molecules in this block.
692 * If so, we only need an array of the size of one molecule.
693 * Otherwise we make an array of #mol times #cgs per molecule.
696 for (m = 0; m < molb->nmol; m++)
698 int am = m*cgs->index[cgs->nr];
699 for (cg = 0; cg < cgs->nr; cg++)
702 a1 = cgs->index[cg+1];
703 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
704 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
708 if (mtop->groups.grpnr[egcQMMM] != nullptr)
710 for (ai = a0; ai < a1; ai++)
712 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
713 mtop->groups.grpnr[egcQMMM][a_offset +ai])
722 cginfo_mb[mb].cg_start = cg_offset;
723 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
724 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
725 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
726 cginfo = cginfo_mb[mb].cginfo;
728 /* Set constraints flags for constrained atoms */
729 snew(a_con, molt->atoms.nr);
730 for (ftype = 0; ftype < F_NRE; ftype++)
732 if (interaction_function[ftype].flags & IF_CONSTRAINT)
737 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
741 for (a = 0; a < nral; a++)
743 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
744 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
750 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
753 int am = m*cgs->index[cgs->nr];
754 for (cg = 0; cg < cgs->nr; cg++)
757 a1 = cgs->index[cg+1];
759 /* Store the energy group in cginfo */
760 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
761 SET_CGINFO_GID(cginfo[cgm+cg], gid);
763 /* Check the intra/inter charge group exclusions */
764 if (a1-a0 > excl_nalloc)
766 excl_nalloc = a1 - a0;
767 srenew(bExcl, excl_nalloc);
769 /* bExclIntraAll: all intra cg interactions excluded
770 * bExclInter: any inter cg interactions excluded
772 bExclIntraAll = TRUE;
776 bHavePerturbedAtoms = FALSE;
777 for (ai = a0; ai < a1; ai++)
779 /* Check VDW and electrostatic interactions */
780 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
781 type_VDW[molt->atoms.atom[ai].typeB]);
782 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
783 molt->atoms.atom[ai].qB != 0);
785 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
787 /* Clear the exclusion list for atom ai */
788 for (aj = a0; aj < a1; aj++)
790 bExcl[aj-a0] = FALSE;
792 /* Loop over all the exclusions of atom ai */
793 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
796 if (aj < a0 || aj >= a1)
805 /* Check if ai excludes a0 to a1 */
806 for (aj = a0; aj < a1; aj++)
810 bExclIntraAll = FALSE;
817 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
820 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
828 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
832 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
834 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
836 /* The size in cginfo is currently only read with DD */
837 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
841 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
845 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
847 if (bHavePerturbedAtoms && fr->efep != efepNO)
849 SET_CGINFO_FEP(cginfo[cgm+cg]);
850 *bFEP_NonBonded = TRUE;
852 /* Store the charge group size */
853 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
855 if (!bExclIntraAll || bExclInter)
857 *bExcl_IntraCGAll_InterCGNone = FALSE;
864 cg_offset += molb->nmol*cgs->nr;
865 a_offset += molb->nmol*cgs->index[cgs->nr];
870 /* the solvent optimizer is called after the QM is initialized,
871 * because we don't want to have the QM subsystemto become an
875 check_solvent(fplog, mtop, fr, cginfo_mb);
877 if (getenv("GMX_NO_SOLV_OPT"))
881 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
882 "Disabling all solvent optimization\n");
884 fr->solvent_opt = esolNO;
888 fr->solvent_opt = esolNO;
890 if (!fr->solvent_opt)
892 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
894 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
896 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
904 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
909 ncg = cgi_mb[nmb-1].cg_end;
912 for (cg = 0; cg < ncg; cg++)
914 while (cg >= cgi_mb[mb].cg_end)
919 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
925 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
927 if (cginfo_mb == nullptr)
931 for (int mb = 0; mb < numMolBlocks; ++mb)
933 sfree(cginfo_mb[mb].cginfo);
938 /* Sets the sum of charges (squared) and C6 in the system in fr.
939 * Returns whether the system has a net charge.
941 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
943 /*This now calculates sum for q and c6*/
944 double qsum, q2sum, q, c6sum, c6;
949 for (const gmx_molblock_t &molb : mtop->molblock)
951 int nmol = molb.nmol;
952 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
953 for (int i = 0; i < atoms->nr; i++)
955 q = atoms->atom[i].q;
958 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
963 fr->q2sum[0] = q2sum;
964 fr->c6sum[0] = c6sum;
966 if (fr->efep != efepNO)
971 for (const gmx_molblock_t &molb : mtop->molblock)
973 int nmol = molb.nmol;
974 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
975 for (int i = 0; i < atoms->nr; i++)
977 q = atoms->atom[i].qB;
980 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
984 fr->q2sum[1] = q2sum;
985 fr->c6sum[1] = c6sum;
990 fr->qsum[1] = fr->qsum[0];
991 fr->q2sum[1] = fr->q2sum[0];
992 fr->c6sum[1] = fr->c6sum[0];
996 if (fr->efep == efepNO)
998 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
1002 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
1003 fr->qsum[0], fr->qsum[1]);
1007 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
1008 return (std::abs(fr->qsum[0]) > 1e-4 ||
1009 std::abs(fr->qsum[1]) > 1e-4);
1012 void update_forcerec(t_forcerec *fr, matrix box)
1014 if (fr->ic->eeltype == eelGRF)
1016 calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
1017 fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
1018 &fr->ic->k_rf, &fr->ic->c_rf);
1022 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
1024 const t_atoms *atoms, *atoms_tpi;
1025 const t_blocka *excl;
1026 int nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
1027 gmx_int64_t npair, npair_ij, tmpi, tmpj;
1028 double csix, ctwelve;
1029 int ntp, *typecount;
1032 real *nbfp_comb = nullptr;
1038 /* For LJ-PME, we want to correct for the difference between the
1039 * actual C6 values and the C6 values used by the LJ-PME based on
1040 * combination rules. */
1042 if (EVDW_PME(fr->ic->vdwtype))
1044 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1045 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1046 for (tpi = 0; tpi < ntp; ++tpi)
1048 for (tpj = 0; tpj < ntp; ++tpj)
1050 C6(nbfp_comb, ntp, tpi, tpj) =
1051 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1052 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1057 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1065 /* Count the types so we avoid natoms^2 operations */
1066 snew(typecount, ntp);
1067 gmx_mtop_count_atomtypes(mtop, q, typecount);
1069 for (tpi = 0; tpi < ntp; tpi++)
1071 for (tpj = tpi; tpj < ntp; tpj++)
1073 tmpi = typecount[tpi];
1074 tmpj = typecount[tpj];
1077 npair_ij = tmpi*tmpj;
1081 npair_ij = tmpi*(tmpi - 1)/2;
1085 /* nbfp now includes the 6.0 derivative prefactor */
1086 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1090 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1091 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1092 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1098 /* Subtract the excluded pairs.
1099 * The main reason for substracting exclusions is that in some cases
1100 * some combinations might never occur and the parameters could have
1101 * any value. These unused values should not influence the dispersion
1104 for (const gmx_molblock_t &molb : mtop->molblock)
1106 int nmol = molb.nmol;
1107 atoms = &mtop->moltype[molb.type].atoms;
1108 excl = &mtop->moltype[molb.type].excls;
1109 for (int i = 0; (i < atoms->nr); i++)
1113 tpi = atoms->atom[i].type;
1117 tpi = atoms->atom[i].typeB;
1119 j1 = excl->index[i];
1120 j2 = excl->index[i+1];
1121 for (j = j1; j < j2; j++)
1128 tpj = atoms->atom[k].type;
1132 tpj = atoms->atom[k].typeB;
1136 /* nbfp now includes the 6.0 derivative prefactor */
1137 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1141 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1142 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1143 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1153 /* Only correct for the interaction of the test particle
1154 * with the rest of the system.
1157 &mtop->moltype[mtop->molblock.back().type].atoms;
1160 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
1162 const gmx_molblock_t &molb = mtop->molblock[mb];
1163 atoms = &mtop->moltype[molb.type].atoms;
1164 for (j = 0; j < atoms->nr; j++)
1167 /* Remove the interaction of the test charge group
1170 if (mb == mtop->molblock.size() - 1)
1174 if (mb == 0 && molb.nmol == 1)
1176 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1181 tpj = atoms->atom[j].type;
1185 tpj = atoms->atom[j].typeB;
1187 for (i = 0; i < fr->n_tpi; i++)
1191 tpi = atoms_tpi->atom[i].type;
1195 tpi = atoms_tpi->atom[i].typeB;
1199 /* nbfp now includes the 6.0 derivative prefactor */
1200 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1204 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1205 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1206 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1213 if (npair - nexcl <= 0 && fplog)
1215 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1221 csix /= npair - nexcl;
1222 ctwelve /= npair - nexcl;
1226 fprintf(debug, "Counted %d exclusions\n", nexcl);
1227 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1228 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1230 fr->avcsix[q] = csix;
1231 fr->avctwelve[q] = ctwelve;
1234 if (EVDW_PME(fr->ic->vdwtype))
1239 if (fplog != nullptr)
1241 if (fr->eDispCorr == edispcAllEner ||
1242 fr->eDispCorr == edispcAllEnerPres)
1244 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1245 fr->avcsix[0], fr->avctwelve[0]);
1249 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1255 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
1257 const t_atoms *at1, *at2;
1258 int i, j, tpi, tpj, ntypes;
1263 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1265 ntypes = mtop->ffparams.atnr;
1268 real bham_b_max = 0;
1269 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
1271 at1 = &mtop->moltype[mt1].atoms;
1272 for (i = 0; (i < at1->nr); i++)
1274 tpi = at1->atom[i].type;
1277 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1280 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
1282 at2 = &mtop->moltype[mt2].atoms;
1283 for (j = 0; (j < at2->nr); j++)
1285 tpj = at2->atom[j].type;
1288 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1290 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
1295 if ((b < bmin) || (bmin == -1))
1305 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1312 static void make_nbf_tables(FILE *fp,
1313 const interaction_const_t *ic, real rtab,
1314 const char *tabfn, char *eg1, char *eg2,
1320 if (tabfn == nullptr)
1324 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1329 sprintf(buf, "%s", tabfn);
1332 /* Append the two energy group names */
1333 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1334 eg1, eg2, ftp2ext(efXVG));
1336 nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1337 /* Copy the contents of the table to separate coulomb and LJ tables too,
1338 * to improve cache performance.
1340 /* For performance reasons we want
1341 * the table data to be aligned to 16-byte. The pointers could be freed
1342 * but currently aren't.
1344 snew(nbl->table_elec, 1);
1345 nbl->table_elec->interaction = GMX_TABLE_INTERACTION_ELEC;
1346 nbl->table_elec->format = nbl->table_elec_vdw->format;
1347 nbl->table_elec->r = nbl->table_elec_vdw->r;
1348 nbl->table_elec->n = nbl->table_elec_vdw->n;
1349 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1350 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1351 nbl->table_elec->ninteractions = 1;
1352 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1353 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1355 snew(nbl->table_vdw, 1);
1356 nbl->table_vdw->interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1357 nbl->table_vdw->format = nbl->table_elec_vdw->format;
1358 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1359 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1360 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1361 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1362 nbl->table_vdw->ninteractions = 2;
1363 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1364 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1366 /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1367 * with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1369 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1371 for (j = 0; j < 4; j++)
1373 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1376 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1378 for (j = 0; j < 8; j++)
1380 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1385 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1386 * ftype2 present in the topology, build an array of the number of
1387 * interactions present for each bonded interaction index found in the
1390 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1391 * valid type with that parameter.
1393 * \c count will be reallocated as necessary to fit the largest bonded
1394 * interaction index found, and its current size will be returned in
1395 * \c ncount. It will contain zero for every bonded interaction index
1396 * for which no interactions are present in the topology.
1398 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1399 int *ncount, int **count)
1402 int ftype, stride, i, j, tabnr;
1404 // Loop over all moleculetypes
1405 for (const gmx_moltype_t &molt : mtop->moltype)
1407 // Loop over all interaction types
1408 for (ftype = 0; ftype < F_NRE; ftype++)
1410 // If the current interaction type is one of the types whose tables we're trying to count...
1411 if (ftype == ftype1 || ftype == ftype2)
1413 il = &molt.ilist[ftype];
1414 stride = 1 + NRAL(ftype);
1415 // ... and there are actually some interactions for this type
1416 for (i = 0; i < il->nr; i += stride)
1418 // Find out which table index the user wanted
1419 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1422 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1424 // Make room for this index in the data structure
1425 if (tabnr >= *ncount)
1427 srenew(*count, tabnr+1);
1428 for (j = *ncount; j < tabnr+1; j++)
1434 // Record that this table index is used and must have a valid file
1442 /*!\brief If there's bonded interactions of flavour \c tabext and type
1443 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1444 * list of filenames passed to mdrun, and make bonded tables from
1447 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1448 * valid type with that parameter.
1450 * A fatal error occurs if no matching filename is found.
1452 static bondedtable_t *make_bonded_tables(FILE *fplog,
1453 int ftype1, int ftype2,
1454 const gmx_mtop_t *mtop,
1455 gmx::ArrayRef<const std::string> tabbfnm,
1465 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1467 // Are there any relevant tabulated bond interactions?
1471 for (int i = 0; i < ncount; i++)
1473 // Do any interactions exist that requires this table?
1476 // This pattern enforces the current requirement that
1477 // table filenames end in a characteristic sequence
1478 // before the file type extension, and avoids table 13
1479 // being recognized and used for table 1.
1480 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1481 bool madeTable = false;
1482 for (size_t j = 0; j < tabbfnm.size() && !madeTable; ++j)
1484 if (gmx::endsWith(tabbfnm[j].c_str(), patternToFind))
1486 // Finally read the table from the file found
1487 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1493 bool isPlural = (ftype2 != -1);
1494 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.",
1495 interaction_function[ftype1].longname,
1496 isPlural ? "' or '" : "",
1497 isPlural ? interaction_function[ftype2].longname : "",
1499 patternToFind.c_str());
1509 void forcerec_set_ranges(t_forcerec *fr,
1510 int ncg_home, int ncg_force,
1512 int natoms_force_constr, int natoms_f_novirsum)
1517 /* fr->ncg_force is unused in the standard code,
1518 * but it can be useful for modified code dealing with charge groups.
1520 fr->ncg_force = ncg_force;
1521 fr->natoms_force = natoms_force;
1522 fr->natoms_force_constr = natoms_force_constr;
1524 if (fr->natoms_force_constr > fr->nalloc_force)
1526 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1529 if (fr->haveDirectVirialContributions)
1531 fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1535 static real cutoff_inf(real cutoff)
1539 cutoff = GMX_CUTOFF_INF;
1545 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, const t_commrec *cr, FILE *fp)
1552 ir->rcoulomb == 0 &&
1554 ir->ePBC == epbcNONE &&
1555 ir->vdwtype == evdwCUT &&
1556 ir->coulombtype == eelCUT &&
1557 ir->efep == efepNO &&
1558 getenv("GMX_NO_ALLVSALL") == nullptr
1561 if (bAllvsAll && ir->opts.ngener > 1)
1563 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";
1569 fprintf(fp, "\n%s\n", note);
1575 if (bAllvsAll && fp && MASTER(cr))
1577 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1584 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
1585 const t_inputrec *ir)
1587 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1589 /* LJ PME with LB combination rule does 7 mesh operations.
1590 * This so slow that we don't compile SIMD non-bonded kernels
1592 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
1600 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1603 const gmx_hw_info_t gmx_unused &hardwareInfo)
1605 *kernel_type = nbnxnk4x4_PlainC;
1606 *ewald_excl = ewaldexclTable;
1610 #ifdef GMX_NBNXN_SIMD_4XN
1611 *kernel_type = nbnxnk4xN_SIMD_4xN;
1613 #ifdef GMX_NBNXN_SIMD_2XNN
1614 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1617 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1618 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1619 * This is based on the SIMD acceleration choice and CPU information
1620 * detected at runtime.
1622 * 4xN calculates more (zero) interactions, but has less pair-search
1623 * work and much better kernel instruction scheduling.
1625 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1626 * which doesn't have FMA, both the analytical and tabulated Ewald
1627 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1628 * 2x(4+4) because it results in significantly fewer pairs.
1629 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1630 * 10% with HT, 50% without HT. As we currently don't detect the actual
1631 * use of HT, use 4x8 to avoid a potential performance hit.
1632 * On Intel Haswell 4x8 is always faster.
1634 *kernel_type = nbnxnk4xN_SIMD_4xN;
1636 #if !GMX_SIMD_HAVE_FMA
1637 if (EEL_PME_EWALD(ir->coulombtype) ||
1638 EVDW_PME(ir->vdwtype))
1640 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1641 * There are enough instructions to make 2x(4+4) efficient.
1643 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1646 if (hardwareInfo.haveAmdZenCpu)
1648 /* One 256-bit FMA per cycle makes 2xNN faster */
1649 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1651 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1654 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
1656 #ifdef GMX_NBNXN_SIMD_4XN
1657 *kernel_type = nbnxnk4xN_SIMD_4xN;
1659 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1662 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
1664 #ifdef GMX_NBNXN_SIMD_2XNN
1665 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1667 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1671 /* Analytical Ewald exclusion correction is only an option in
1673 * Since table lookup's don't parallelize with SIMD, analytical
1674 * will probably always be faster for a SIMD width of 8 or more.
1675 * With FMA analytical is sometimes faster for a width if 4 as well.
1676 * In single precision, this is faster on Bulldozer.
1678 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1679 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
1680 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
1681 * of single or double precision and 128 or 256-bit AVX2.
1683 if (!hardwareInfo.haveAmdZenCpu)
1685 *ewald_excl = ewaldexclAnalytical;
1688 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
1690 *ewald_excl = ewaldexclTable;
1692 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
1694 *ewald_excl = ewaldexclAnalytical;
1702 const char *lookup_nbnxn_kernel_name(int kernel_type)
1704 const char *returnvalue = nullptr;
1705 switch (kernel_type)
1708 returnvalue = "not set";
1710 case nbnxnk4x4_PlainC:
1711 returnvalue = "plain C";
1713 case nbnxnk4xN_SIMD_4xN:
1714 case nbnxnk4xN_SIMD_2xNN:
1716 returnvalue = "SIMD";
1718 returnvalue = "not available";
1721 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1722 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1726 gmx_fatal(FARGS, "Illegal kernel type selected");
1727 returnvalue = nullptr;
1733 static void pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
1734 gmx_bool use_simd_kernels,
1735 const gmx_hw_info_t &hardwareInfo,
1737 EmulateGpuNonbonded emulateGpu,
1738 const t_inputrec *ir,
1741 gmx_bool bDoNonbonded)
1743 assert(kernel_type);
1745 *kernel_type = nbnxnkNotSet;
1746 *ewald_excl = ewaldexclTable;
1748 if (emulateGpu == EmulateGpuNonbonded::Yes)
1750 *kernel_type = nbnxnk8x8x8_PlainC;
1754 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
1759 *kernel_type = nbnxnk8x8x8_GPU;
1762 if (*kernel_type == nbnxnkNotSet)
1764 if (use_simd_kernels &&
1765 nbnxn_simd_supported(mdlog, ir))
1767 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl, hardwareInfo);
1771 *kernel_type = nbnxnk4x4_PlainC;
1777 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
1778 "Using %s %dx%d nonbonded short-range kernels",
1779 lookup_nbnxn_kernel_name(*kernel_type),
1780 nbnxn_kernel_to_cluster_i_size(*kernel_type),
1781 nbnxn_kernel_to_cluster_j_size(*kernel_type));
1783 if (nbnxnk4x4_PlainC == *kernel_type ||
1784 nbnxnk8x8x8_PlainC == *kernel_type)
1786 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1787 "WARNING: Using the slow %s kernels. This should\n"
1788 "not happen during routine usage on supported platforms.",
1789 lookup_nbnxn_kernel_name(*kernel_type));
1794 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1795 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1796 bool systemHasNetCharge,
1797 interaction_const_t *ic)
1799 if (!EEL_PME_EWALD(ir->coulombtype))
1806 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1808 if (ir->coulombtype == eelP3M_AD)
1810 please_cite(fp, "Hockney1988");
1811 please_cite(fp, "Ballenegger2012");
1815 please_cite(fp, "Essmann95a");
1818 if (ir->ewald_geometry == eewg3DC)
1822 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1823 systemHasNetCharge ? " and net charge" : "");
1825 please_cite(fp, "In-Chul99a");
1826 if (systemHasNetCharge)
1828 please_cite(fp, "Ballenegger2009");
1833 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1836 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1837 1/ic->ewaldcoeff_q);
1840 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1842 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1843 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1851 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1852 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1853 interaction_const_t *ic)
1855 if (!EVDW_PME(ir->vdwtype))
1862 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1863 please_cite(fp, "Essmann95a");
1865 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1868 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1869 1/ic->ewaldcoeff_lj);
1872 if (ic->vdw_modifier == eintmodPOTSHIFT)
1874 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1875 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1879 ic->sh_lj_ewald = 0;
1883 gmx_bool uses_simple_tables(int cutoff_scheme,
1884 nonbonded_verlet_t *nbv,
1887 gmx_bool bUsesSimpleTables = TRUE;
1890 switch (cutoff_scheme)
1893 bUsesSimpleTables = TRUE;
1896 assert(NULL != nbv && NULL != nbv->grp);
1897 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1898 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1901 gmx_incons("unimplemented");
1903 return bUsesSimpleTables;
1906 static void init_ewald_f_table(interaction_const_t *ic,
1911 /* Get the Ewald table spacing based on Coulomb and/or LJ
1912 * Ewald coefficients and rtol.
1914 ic->tabq_scale = ewald_spline3_table_scale(ic);
1916 if (ic->cutoff_scheme == ecutsVERLET)
1918 maxr = ic->rcoulomb;
1922 maxr = std::max(ic->rcoulomb, rtab);
1924 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1926 sfree_aligned(ic->tabq_coul_FDV0);
1927 sfree_aligned(ic->tabq_coul_F);
1928 sfree_aligned(ic->tabq_coul_V);
1930 sfree_aligned(ic->tabq_vdw_FDV0);
1931 sfree_aligned(ic->tabq_vdw_F);
1932 sfree_aligned(ic->tabq_vdw_V);
1934 if (EEL_PME_EWALD(ic->eeltype))
1936 /* Create the original table data in FDV0 */
1937 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1938 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1939 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1940 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1941 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1944 if (EVDW_PME(ic->vdwtype))
1946 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1947 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1948 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1949 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1950 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1954 void init_interaction_const_tables(FILE *fp,
1955 interaction_const_t *ic,
1958 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1960 init_ewald_f_table(ic, rtab);
1964 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1965 1/ic->tabq_scale, ic->tabq_size);
1970 static void clear_force_switch_constants(shift_consts_t *sc)
1977 static void force_switch_constants(real p,
1981 /* Here we determine the coefficient for shifting the force to zero
1982 * between distance rsw and the cut-off rc.
1983 * For a potential of r^-p, we have force p*r^-(p+1).
1984 * But to save flops we absorb p in the coefficient.
1986 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1987 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1989 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1990 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1991 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1994 static void potential_switch_constants(real rsw, real rc,
1995 switch_consts_t *sc)
1997 /* The switch function is 1 at rsw and 0 at rc.
1998 * The derivative and second derivate are zero at both ends.
1999 * rsw = max(r - r_switch, 0)
2000 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
2001 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
2002 * force = force*dsw - potential*sw
2005 sc->c3 = -10/gmx::power3(rc - rsw);
2006 sc->c4 = 15/gmx::power4(rc - rsw);
2007 sc->c5 = -6/gmx::power5(rc - rsw);
2010 /*! \brief Construct interaction constants
2012 * This data is used (particularly) by search and force code for
2013 * short-range interactions. Many of these are constant for the whole
2014 * simulation; some are constant only after PME tuning completes.
2017 init_interaction_const(FILE *fp,
2018 interaction_const_t **interaction_const,
2019 const t_inputrec *ir,
2020 const gmx_mtop_t *mtop,
2021 bool systemHasNetCharge)
2023 interaction_const_t *ic;
2027 ic->cutoff_scheme = ir->cutoff_scheme;
2029 /* Just allocate something so we can free it */
2030 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
2031 snew_aligned(ic->tabq_coul_F, 16, 32);
2032 snew_aligned(ic->tabq_coul_V, 16, 32);
2035 ic->vdwtype = ir->vdwtype;
2036 ic->vdw_modifier = ir->vdw_modifier;
2037 ic->reppow = mtop->ffparams.reppow;
2038 ic->rvdw = cutoff_inf(ir->rvdw);
2039 ic->rvdw_switch = ir->rvdw_switch;
2040 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
2041 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
2042 if (ic->useBuckingham)
2044 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
2047 initVdwEwaldParameters(fp, ir, ic);
2049 clear_force_switch_constants(&ic->dispersion_shift);
2050 clear_force_switch_constants(&ic->repulsion_shift);
2052 switch (ic->vdw_modifier)
2054 case eintmodPOTSHIFT:
2055 /* Only shift the potential, don't touch the force */
2056 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
2057 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
2059 case eintmodFORCESWITCH:
2060 /* Switch the force, switch and shift the potential */
2061 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2062 &ic->dispersion_shift);
2063 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2064 &ic->repulsion_shift);
2066 case eintmodPOTSWITCH:
2067 /* Switch the potential and force */
2068 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2072 case eintmodEXACTCUTOFF:
2073 /* Nothing to do here */
2076 gmx_incons("unimplemented potential modifier");
2079 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2081 /* Electrostatics */
2082 ic->eeltype = ir->coulombtype;
2083 ic->coulomb_modifier = ir->coulomb_modifier;
2084 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
2085 ic->rcoulomb_switch = ir->rcoulomb_switch;
2086 ic->epsilon_r = ir->epsilon_r;
2088 /* Set the Coulomb energy conversion factor */
2089 if (ic->epsilon_r != 0)
2091 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
2095 /* eps = 0 is infinite dieletric: no Coulomb interactions */
2099 /* Reaction-field */
2100 if (EEL_RF(ic->eeltype))
2102 ic->epsilon_rf = ir->epsilon_rf;
2103 /* Generalized reaction field parameters are updated every step */
2104 if (ic->eeltype != eelGRF)
2106 calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
2107 ic->rcoulomb, 0, 0, NULL,
2108 &ic->k_rf, &ic->c_rf);
2111 if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
2113 /* grompp should have done this, but this scheme is obsolete */
2114 ic->coulomb_modifier = eintmodEXACTCUTOFF;
2119 /* For plain cut-off we might use the reaction-field kernels */
2120 ic->epsilon_rf = ic->epsilon_r;
2122 if (ir->coulomb_modifier == eintmodPOTSHIFT)
2124 ic->c_rf = 1/ic->rcoulomb;
2132 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
2136 real dispersion_shift;
2138 dispersion_shift = ic->dispersion_shift.cpot;
2139 if (EVDW_PME(ic->vdwtype))
2141 dispersion_shift -= ic->sh_lj_ewald;
2143 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2144 ic->repulsion_shift.cpot, dispersion_shift);
2146 if (ic->eeltype == eelCUT)
2148 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2150 else if (EEL_PME(ic->eeltype))
2152 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2157 *interaction_const = ic;
2161 done_interaction_const(interaction_const_t *interaction_const)
2163 sfree_aligned(interaction_const->tabq_coul_FDV0);
2164 sfree_aligned(interaction_const->tabq_coul_F);
2165 sfree_aligned(interaction_const->tabq_coul_V);
2166 sfree(interaction_const);
2169 static void init_nb_verlet(const gmx::MDLogger &mdlog,
2170 nonbonded_verlet_t **nb_verlet,
2171 gmx_bool bFEP_NonBonded,
2172 const t_inputrec *ir,
2173 const t_forcerec *fr,
2174 const t_commrec *cr,
2175 const gmx_hw_info_t &hardwareInfo,
2176 const gmx_device_info_t *deviceInfo,
2177 const gmx_mtop_t *mtop,
2180 nonbonded_verlet_t *nbv;
2183 nbnxn_alloc_t *nb_alloc;
2184 nbnxn_free_t *nb_free;
2186 nbv = new nonbonded_verlet_t();
2188 nbv->emulateGpu = ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
2189 nbv->bUseGPU = deviceInfo != nullptr;
2191 GMX_RELEASE_ASSERT(!(nbv->emulateGpu == EmulateGpuNonbonded::Yes && nbv->bUseGPU), "When GPU emulation is active, there cannot be a GPU assignment");
2194 nbv->min_ci_balanced = 0;
2196 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2197 for (int i = 0; i < nbv->ngrp; i++)
2199 nbv->grp[i].nbl_lists.nnbl = 0;
2200 nbv->grp[i].kernel_type = nbnxnkNotSet;
2202 if (i == 0) /* local */
2204 pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
2205 nbv->bUseGPU, nbv->emulateGpu, ir,
2206 &nbv->grp[i].kernel_type,
2207 &nbv->grp[i].ewald_excl,
2210 else /* non-local */
2212 /* Use the same kernel for local and non-local interactions */
2213 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2214 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2218 nbv->listParams = std::unique_ptr<NbnxnListParameters>(new NbnxnListParameters(ir->rlist));
2219 setupDynamicPairlistPruning(mdlog, ir, mtop, box, nbv->grp[0].kernel_type, fr->ic,
2220 nbv->listParams.get());
2222 nbnxn_init_search(&nbv->nbs,
2223 DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
2224 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
2226 gmx_omp_nthreads_get(emntPairsearch));
2228 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2229 &nb_alloc, &nb_free);
2231 for (int i = 0; i < nbv->ngrp; i++)
2233 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2234 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2235 /* 8x8x8 "non-simple" lists are ATM always combined */
2236 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2240 int enbnxninitcombrule;
2241 if (fr->ic->vdwtype == evdwCUT &&
2242 (fr->ic->vdw_modifier == eintmodNONE ||
2243 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
2244 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
2246 /* Plain LJ cut-off: we can optimize with combination rules */
2247 enbnxninitcombrule = enbnxninitcombruleDETECT;
2249 else if (fr->ic->vdwtype == evdwPME)
2251 /* LJ-PME: we need to use a combination rule for the grid */
2252 if (fr->ljpme_combination_rule == eljpmeGEOM)
2254 enbnxninitcombrule = enbnxninitcombruleGEOM;
2258 enbnxninitcombrule = enbnxninitcombruleLB;
2263 /* We use a full combination matrix: no rule required */
2264 enbnxninitcombrule = enbnxninitcombruleNONE;
2268 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
2269 if (ir->opts.ngener - ir->nwall == 1)
2271 /* We have only one non-wall energy group, we do not need energy group
2272 * support in the non-bondeds kernels, since all non-bonded energy
2273 * contributions go to the first element of the energy group matrix.
2275 mimimumNumEnergyGroupNonbonded = 1;
2277 bool bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[0].kernel_type);
2278 nbnxn_atomdata_init(mdlog,
2280 nbv->grp[0].kernel_type,
2282 fr->ntype, fr->nbfp,
2283 mimimumNumEnergyGroupNonbonded,
2284 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2289 /* init the NxN GPU data; the last argument tells whether we'll have
2290 * both local and non-local NB calculation on GPU */
2291 nbnxn_gpu_init(&nbv->gpu_nbv,
2294 nbv->listParams.get(),
2299 if ((env = getenv("GMX_NB_MIN_CI")) != nullptr)
2303 nbv->min_ci_balanced = strtol(env, &end, 10);
2304 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
2306 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
2311 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2312 nbv->min_ci_balanced);
2317 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2320 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2321 nbv->min_ci_balanced);
2330 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2332 return nbv != nullptr && nbv->bUseGPU;
2335 void init_forcerec(FILE *fp,
2336 const gmx::MDLogger &mdlog,
2339 const t_inputrec *ir,
2340 const gmx_mtop_t *mtop,
2341 const t_commrec *cr,
2345 gmx::ArrayRef<const std::string> tabbfnm,
2346 const gmx_hw_info_t &hardwareInfo,
2347 const gmx_device_info_t *deviceInfo,
2348 gmx_bool bNoSolvOpt,
2351 int m, negp_pp, negptable, egi, egj;
2356 gmx_bool bGenericKernelOnly;
2357 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
2358 gmx_bool bFEP_NonBonded;
2359 int *nm_ind, egp_flags;
2361 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2362 fr->use_simd_kernels = TRUE;
2364 fr->bDomDec = DOMAINDECOMP(cr);
2366 if (check_box(ir->ePBC, box))
2368 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2371 /* Test particle insertion ? */
2374 /* Set to the size of the molecule to be inserted (the last one) */
2375 /* Because of old style topologies, we have to use the last cg
2376 * instead of the last molecule type.
2378 cgs = &mtop->moltype[mtop->molblock.back().type].cgs;
2379 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2380 gmx::BlockRanges molecules = gmx_mtop_molecules(*mtop);
2381 if (fr->n_tpi != molecules.index[molecules.numBlocks()] - molecules.index[molecules.numBlocks() - 1])
2383 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2391 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
2393 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
2394 eel_names[ir->coulombtype]);
2399 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
2401 if (ir->useTwinRange)
2403 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
2405 /* Copy the user determined parameters */
2406 fr->userint1 = ir->userint1;
2407 fr->userint2 = ir->userint2;
2408 fr->userint3 = ir->userint3;
2409 fr->userint4 = ir->userint4;
2410 fr->userreal1 = ir->userreal1;
2411 fr->userreal2 = ir->userreal2;
2412 fr->userreal3 = ir->userreal3;
2413 fr->userreal4 = ir->userreal4;
2416 fr->fc_stepsize = ir->fc_stepsize;
2419 fr->efep = ir->efep;
2420 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2421 if (ir->fepvals->bScCoul)
2423 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2424 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
2428 fr->sc_alphacoul = 0;
2429 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2431 fr->sc_power = ir->fepvals->sc_power;
2432 fr->sc_r_power = ir->fepvals->sc_r_power;
2433 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
2435 env = getenv("GMX_SCSIGMA_MIN");
2439 sscanf(env, "%20lf", &dbl);
2440 fr->sc_sigma6_min = gmx::power6(dbl);
2443 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2447 fr->bNonbonded = TRUE;
2448 if (getenv("GMX_NO_NONBONDED") != nullptr)
2450 /* turn off non-bonded calculations */
2451 fr->bNonbonded = FALSE;
2452 GMX_LOG(mdlog.warning).asParagraph().appendText(
2453 "Found environment variable GMX_NO_NONBONDED.\n"
2454 "Disabling nonbonded calculations.");
2457 bGenericKernelOnly = FALSE;
2459 /* We now check in the NS code whether a particular combination of interactions
2460 * can be used with water optimization, and disable it if that is not the case.
2463 if (getenv("GMX_NB_GENERIC") != nullptr)
2468 "Found environment variable GMX_NB_GENERIC.\n"
2469 "Disabling all interaction-specific nonbonded kernels, will only\n"
2470 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2472 bGenericKernelOnly = TRUE;
2475 if (bGenericKernelOnly == TRUE)
2480 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
2482 fr->use_simd_kernels = FALSE;
2486 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2487 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2488 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2492 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2494 /* Check if we can/should do all-vs-all kernels */
2495 fr->bAllvsAll = can_use_allvsall(ir, FALSE, nullptr, nullptr);
2496 fr->AllvsAll_work = nullptr;
2498 /* All-vs-all kernels have not been implemented in 4.6 and later.
2499 * See Redmine #1249. */
2502 fr->bAllvsAll = FALSE;
2506 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2507 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2508 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2509 "or try cutoff-scheme = Verlet.\n\n");
2513 /* Neighbour searching stuff */
2514 fr->cutoff_scheme = ir->cutoff_scheme;
2515 fr->bGrid = (ir->ns_type == ensGRID);
2516 fr->ePBC = ir->ePBC;
2518 if (fr->cutoff_scheme == ecutsGROUP)
2520 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2521 "removed in a future release when 'verlet' supports all interaction forms.\n";
2525 fprintf(stderr, "\n%s\n", note);
2529 fprintf(fp, "\n%s\n", note);
2533 /* Determine if we will do PBC for distances in bonded interactions */
2534 if (fr->ePBC == epbcNONE)
2536 fr->bMolPBC = FALSE;
2540 if (!DOMAINDECOMP(cr))
2544 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2545 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2546 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2548 /* The group cut-off scheme and SHAKE assume charge groups
2549 * are whole, but not using molpbc is faster in most cases.
2550 * With intermolecular interactions we need PBC for calculating
2551 * distances between atoms in different molecules.
2553 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2554 !mtop->bIntermolecularInteractions)
2556 fr->bMolPBC = ir->bPeriodicMols;
2558 if (bSHAKE && fr->bMolPBC)
2560 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2565 /* Not making molecules whole is faster in most cases,
2566 * but With orientation restraints we need whole molecules.
2568 fr->bMolPBC = (fcd->orires.nr == 0);
2570 if (getenv("GMX_USE_GRAPH") != nullptr)
2572 fr->bMolPBC = FALSE;
2575 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
2578 if (mtop->bIntermolecularInteractions)
2580 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
2584 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
2586 if (bSHAKE && fr->bMolPBC)
2588 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");
2594 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2598 fr->rc_scaling = ir->refcoord_scaling;
2599 copy_rvec(ir->posres_com, fr->posres_com);
2600 copy_rvec(ir->posres_comB, fr->posres_comB);
2601 fr->rlist = cutoff_inf(ir->rlist);
2602 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2604 /* This now calculates sum for q and c6*/
2605 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
2607 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2608 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
2609 init_interaction_const_tables(fp, fr->ic, ir->rlist + ir->tabext);
2611 const interaction_const_t *ic = fr->ic;
2613 /* TODO: Replace this Ewald table or move it into interaction_const_t */
2614 if (ir->coulombtype == eelEWALD)
2616 init_ewald_tab(&(fr->ewald_table), ir, fp);
2619 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2620 switch (ic->eeltype)
2623 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
2628 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2632 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2633 GMX_RELEASE_ASSERT(ic->coulomb_modifier == eintmodEXACTCUTOFF, "With the group scheme RF-zero needs the exact cut-off modifier");
2642 case eelPMEUSERSWITCH:
2643 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2649 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2653 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
2656 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
2658 /* Vdw: Translate from mdp settings to kernel format */
2659 switch (ic->vdwtype)
2664 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2668 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2672 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2678 case evdwENCADSHIFT:
2679 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2683 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
2686 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
2688 if (ir->cutoff_scheme == ecutsGROUP)
2690 fr->bvdwtab = ((ic->vdwtype != evdwCUT || !gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2691 && !EVDW_PME(ic->vdwtype));
2692 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2693 fr->bcoultab = !(ic->eeltype == eelCUT ||
2694 ic->eeltype == eelEWALD ||
2695 ic->eeltype == eelPME ||
2696 ic->eeltype == eelP3M_AD ||
2697 ic->eeltype == eelRF ||
2698 ic->eeltype == eelRF_ZERO);
2700 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2701 * going to be faster to tabulate the interaction than calling the generic kernel.
2702 * However, if generic kernels have been requested we keep things analytically.
2704 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2705 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2706 bGenericKernelOnly == FALSE)
2708 if ((ic->rcoulomb_switch != ic->rvdw_switch) || (ic->rcoulomb != ic->rvdw))
2710 fr->bcoultab = TRUE;
2711 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2712 * which would otherwise need two tables.
2716 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2717 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2718 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2719 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2721 if ((ic->rcoulomb != ic->rvdw) && (bGenericKernelOnly == FALSE))
2723 fr->bcoultab = TRUE;
2727 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2729 fr->bcoultab = TRUE;
2731 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2736 if (getenv("GMX_REQUIRE_TABLES"))
2739 fr->bcoultab = TRUE;
2744 fprintf(fp, "Table routines are used for coulomb: %s\n",
2745 gmx::boolToString(fr->bcoultab));
2746 fprintf(fp, "Table routines are used for vdw: %s\n",
2747 gmx::boolToString(fr->bvdwtab));
2750 if (fr->bvdwtab == TRUE)
2752 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2753 fr->nbkernel_vdw_modifier = eintmodNONE;
2755 if (fr->bcoultab == TRUE)
2757 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2758 fr->nbkernel_elec_modifier = eintmodNONE;
2762 if (ir->cutoff_scheme == ecutsVERLET)
2764 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2766 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2768 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2769 * while mdrun does not (and never did) support this.
2771 if (EEL_USER(fr->ic->eeltype))
2773 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2774 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2777 fr->bvdwtab = FALSE;
2778 fr->bcoultab = FALSE;
2781 /* 1-4 interaction electrostatics */
2782 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2784 /* Parameters for generalized RF */
2788 if (ic->eeltype == eelGRF)
2790 init_generalized_rf(fp, mtop, ir, fr);
2793 fr->haveDirectVirialContributions =
2794 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
2795 fr->forceProviders->hasForceProvider() ||
2796 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2797 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2802 if (fr->haveDirectVirialContributions)
2804 fr->forceBufferForDirectVirialContributions = new std::vector<gmx::RVec>;
2807 if (fr->cutoff_scheme == ecutsGROUP &&
2808 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2810 /* Count the total number of charge groups */
2811 fr->cg_nalloc = ncg_mtop(mtop);
2812 srenew(fr->cg_cm, fr->cg_nalloc);
2814 if (fr->shift_vec == nullptr)
2816 snew(fr->shift_vec, SHIFTS);
2819 if (fr->fshift == nullptr)
2821 snew(fr->fshift, SHIFTS);
2824 if (fr->nbfp == nullptr)
2826 fr->ntype = mtop->ffparams.atnr;
2827 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2828 if (EVDW_PME(ic->vdwtype))
2830 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2834 /* Copy the energy group exclusions */
2835 fr->egp_flags = ir->opts.egp_flags;
2837 /* Van der Waals stuff */
2838 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
2840 if (ic->rvdw_switch >= ic->rvdw)
2842 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2843 ic->rvdw_switch, ic->rvdw);
2847 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2848 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
2849 ic->rvdw_switch, ic->rvdw);
2853 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
2855 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2858 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
2860 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2863 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2865 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2868 if (fp && fr->cutoff_scheme == ecutsGROUP)
2870 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2871 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
2874 fr->eDispCorr = ir->eDispCorr;
2875 fr->numAtomsForDispersionCorrection = mtop->natoms;
2876 if (ir->eDispCorr != edispcNO)
2878 set_avcsixtwelve(fp, fr, mtop);
2881 if (ir->implicit_solvent)
2883 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
2886 /* Construct tables for the group scheme. A little unnecessary to
2887 * make both vdw and coul tables sometimes, but what the
2888 * heck. Note that both cutoff schemes construct Ewald tables in
2889 * init_interaction_const_tables. */
2890 needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
2891 (fr->bcoultab || fr->bvdwtab));
2893 negp_pp = ir->opts.ngener - ir->nwall;
2895 if (!needGroupSchemeTables)
2897 bSomeNormalNbListsAreInUse = TRUE;
2902 bSomeNormalNbListsAreInUse = FALSE;
2903 for (egi = 0; egi < negp_pp; egi++)
2905 for (egj = egi; egj < negp_pp; egj++)
2907 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2908 if (!(egp_flags & EGP_EXCL))
2910 if (egp_flags & EGP_TABLE)
2916 bSomeNormalNbListsAreInUse = TRUE;
2921 if (bSomeNormalNbListsAreInUse)
2923 fr->nnblists = negptable + 1;
2927 fr->nnblists = negptable;
2929 if (fr->nnblists > 1)
2931 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
2935 snew(fr->nblists, fr->nnblists);
2937 /* This code automatically gives table length tabext without cut-off's,
2938 * in that case grompp should already have checked that we do not need
2939 * normal tables and we only generate tables for 1-4 interactions.
2941 rtab = ir->rlist + ir->tabext;
2943 if (needGroupSchemeTables)
2945 /* make tables for ordinary interactions */
2946 if (bSomeNormalNbListsAreInUse)
2948 make_nbf_tables(fp, ic, rtab, tabfn, nullptr, nullptr, &fr->nblists[0]);
2957 /* Read the special tables for certain energy group pairs */
2958 nm_ind = mtop->groups.grps[egcENER].nm_ind;
2959 for (egi = 0; egi < negp_pp; egi++)
2961 for (egj = egi; egj < negp_pp; egj++)
2963 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2964 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2966 if (fr->nnblists > 1)
2968 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
2970 /* Read the table file with the two energy groups names appended */
2971 make_nbf_tables(fp, ic, rtab, tabfn,
2972 *mtop->groups.grpname[nm_ind[egi]],
2973 *mtop->groups.grpname[nm_ind[egj]],
2977 else if (fr->nnblists > 1)
2979 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
2986 /* Tables might not be used for the potential modifier
2987 * interactions per se, but we still need them to evaluate
2988 * switch/shift dispersion corrections in this case. */
2989 if (fr->eDispCorr != edispcNO)
2991 fr->dispersionCorrectionTable = makeDispersionCorrectionTable(fp, ic, rtab, tabfn);
2994 /* We want to use unmodified tables for 1-4 coulombic
2995 * interactions, so we must in general have an extra set of
2997 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2998 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2999 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
3001 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
3002 GMX_MAKETABLES_14ONLY);
3006 fr->nwall = ir->nwall;
3007 if (ir->nwall && ir->wall_type == ewtTABLE)
3009 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
3012 if (fcd && !tabbfnm.empty())
3014 // Need to catch std::bad_alloc
3015 // TODO Don't need to catch this here, when merging with master branch
3018 fcd->bondtab = make_bonded_tables(fp,
3019 F_TABBONDS, F_TABBONDSNC,
3020 mtop, tabbfnm, "b");
3021 fcd->angletab = make_bonded_tables(fp,
3023 mtop, tabbfnm, "a");
3024 fcd->dihtab = make_bonded_tables(fp,
3026 mtop, tabbfnm, "d");
3028 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
3034 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3038 /* QM/MM initialization if requested
3042 fprintf(stderr, "QM/MM calculation requested.\n");
3045 fr->bQMMM = ir->bQMMM;
3046 fr->qr = mk_QMMMrec();
3048 /* Set all the static charge group info */
3049 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3051 &fr->bExcl_IntraCGAll_InterCGNone);
3052 if (DOMAINDECOMP(cr))
3054 fr->cginfo = nullptr;
3058 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
3061 if (!DOMAINDECOMP(cr))
3063 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3064 mtop->natoms, mtop->natoms, mtop->natoms);
3067 fr->print_force = print_force;
3070 /* coarse load balancing vars */
3075 /* Initialize neighbor search */
3077 init_ns(fp, cr, fr->ns, fr, mtop);
3079 if (thisRankHasDuty(cr, DUTY_PP))
3081 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3084 /* Initialize the thread working data for bonded interactions */
3085 init_bonded_threading(fp, mtop->groups.grps[egcENER].nr,
3086 &fr->bonded_threading);
3088 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
3089 snew(fr->ewc_t, fr->nthread_ewc);
3091 if (fr->cutoff_scheme == ecutsVERLET)
3093 // We checked the cut-offs in grompp, but double-check here.
3094 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
3095 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
3097 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
3101 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
3104 init_nb_verlet(mdlog, &fr->nbv, bFEP_NonBonded, ir, fr,
3105 cr, hardwareInfo, deviceInfo,
3111 /* Here we switch from using mdlog, which prints the newline before
3112 * the paragraph, to our old fprintf logging, which prints the newline
3113 * after the paragraph, so we should add a newline here.
3118 if (ir->eDispCorr != edispcNO)
3120 calc_enervirdiff(fp, ir->eDispCorr, fr);
3124 /* Frees GPU memory and sets a tMPI node barrier.
3126 * Note that this function needs to be called even if GPUs are not used
3127 * in this run because the PME ranks have no knowledge of whether GPUs
3128 * are used or not, but all ranks need to enter the barrier below.
3129 * \todo Remove physical node barrier from this function after making sure
3130 * that it's not needed anymore (with a shared GPU run).
3132 void free_gpu_resources(const t_forcerec *fr,
3133 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator)
3135 bool isPPrankUsingGPU = fr && fr->nbv && fr->nbv->bUseGPU;
3137 /* stop the GPU profiler (only CUDA) */
3140 if (isPPrankUsingGPU)
3142 /* free nbnxn data in GPU memory */
3143 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3146 /* With tMPI we need to wait for all ranks to finish deallocation before
3147 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
3150 * This is not a concern in OpenCL where we use one context per rank which
3151 * is freed in nbnxn_gpu_free().
3153 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
3154 * but it is easier and more futureproof to call it on the whole node.
3158 physicalNodeCommunicator.barrier();
3162 void done_forcerec(t_forcerec *fr, int numMolBlocks, int numEnergyGroups)
3166 // PME-only ranks don't have a forcerec
3169 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
3171 done_interaction_const(fr->ic);
3172 sfree(fr->shift_vec);
3175 done_ns(fr->ns, numEnergyGroups);