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/pleasecite.h"
99 #include "gromacs/utility/smalloc.h"
100 #include "gromacs/utility/strconvert.h"
102 #include "nbnxn_gpu_jit_support.h"
104 const char *egrp_nm[egNR+1] = {
105 "Coul-SR", "LJ-SR", "Buck-SR",
106 "Coul-14", "LJ-14", nullptr
109 t_forcerec *mk_forcerec(void)
119 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
123 for (i = 0; (i < atnr); i++)
125 for (j = 0; (j < atnr); j++)
127 fprintf(fp, "%2d - %2d", i, j);
130 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
131 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
135 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
136 C12(nbfp, atnr, i, j)/12.0);
143 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
151 snew(nbfp, 3*atnr*atnr);
152 for (i = k = 0; (i < atnr); i++)
154 for (j = 0; (j < atnr); j++, k++)
156 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
157 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
158 /* nbfp now includes the 6.0 derivative prefactor */
159 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
165 snew(nbfp, 2*atnr*atnr);
166 for (i = k = 0; (i < atnr); i++)
168 for (j = 0; (j < atnr); j++, k++)
170 /* nbfp now includes the 6.0/12.0 derivative prefactors */
171 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
172 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
180 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
183 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
186 /* For LJ-PME simulations, we correct the energies with the reciprocal space
187 * inside of the cut-off. To do this the non-bonded kernels needs to have
188 * access to the C6-values used on the reciprocal grid in pme.c
192 snew(grid, 2*atnr*atnr);
193 for (i = k = 0; (i < atnr); i++)
195 for (j = 0; (j < atnr); j++, k++)
197 c6i = idef->iparams[i*(atnr+1)].lj.c6;
198 c12i = idef->iparams[i*(atnr+1)].lj.c12;
199 c6j = idef->iparams[j*(atnr+1)].lj.c6;
200 c12j = idef->iparams[j*(atnr+1)].lj.c12;
201 c6 = std::sqrt(c6i * c6j);
202 if (fr->ljpme_combination_rule == eljpmeLB
203 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
205 sigmai = gmx::sixthroot(c12i / c6i);
206 sigmaj = gmx::sixthroot(c12j / c6j);
207 epsi = c6i * c6i / c12i;
208 epsj = c6j * c6j / c12j;
209 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
211 /* Store the elements at the same relative positions as C6 in nbfp in order
212 * to simplify access in the kernels
214 grid[2*(atnr*i+j)] = c6*6.0;
220 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
224 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
228 snew(nbfp, 2*atnr*atnr);
229 for (i = 0; i < atnr; ++i)
231 for (j = 0; j < atnr; ++j)
233 c6i = idef->iparams[i*(atnr+1)].lj.c6;
234 c12i = idef->iparams[i*(atnr+1)].lj.c12;
235 c6j = idef->iparams[j*(atnr+1)].lj.c6;
236 c12j = idef->iparams[j*(atnr+1)].lj.c12;
237 c6 = std::sqrt(c6i * c6j);
238 c12 = std::sqrt(c12i * c12j);
239 if (comb_rule == eCOMB_ARITHMETIC
240 && !gmx_numzero(c6) && !gmx_numzero(c12))
242 sigmai = gmx::sixthroot(c12i / c6i);
243 sigmaj = gmx::sixthroot(c12j / c6j);
244 epsi = c6i * c6i / c12i;
245 epsj = c6j * c6j / c12j;
246 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
247 c12 = std::sqrt(epsi * epsj) * gmx::power12(0.5*(sigmai+sigmaj));
249 C6(nbfp, atnr, i, j) = c6*6.0;
250 C12(nbfp, atnr, i, j) = c12*12.0;
256 /* This routine sets fr->solvent_opt to the most common solvent in the
257 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
258 * the fr->solvent_type array with the correct type (or esolNO).
260 * Charge groups that fulfill the conditions but are not identical to the
261 * most common one will be marked as esolNO in the solvent_type array.
263 * TIP3p is identical to SPC for these purposes, so we call it
264 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
266 * NOTE: QM particle should not
267 * become an optimized solvent. Not even if there is only one charge
277 } solvent_parameters_t;
280 check_solvent_cg(const gmx_moltype_t *molt,
283 const unsigned char *qm_grpnr,
284 const t_grps *qm_grps,
286 int *n_solvent_parameters,
287 solvent_parameters_t **solvent_parameters_p,
297 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
298 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
301 solvent_parameters_t *solvent_parameters;
303 /* We use a list with parameters for each solvent type.
304 * Every time we discover a new molecule that fulfills the basic
305 * conditions for a solvent we compare with the previous entries
306 * in these lists. If the parameters are the same we just increment
307 * the counter for that type, and otherwise we create a new type
308 * based on the current molecule.
310 * Once we've finished going through all molecules we check which
311 * solvent is most common, and mark all those molecules while we
312 * clear the flag on all others.
315 solvent_parameters = *solvent_parameters_p;
317 /* Mark the cg first as non optimized */
320 /* Check if this cg has no exclusions with atoms in other charge groups
321 * and all atoms inside the charge group excluded.
322 * We only have 3 or 4 atom solvent loops.
324 if (GET_CGINFO_EXCL_INTER(cginfo) ||
325 !GET_CGINFO_EXCL_INTRA(cginfo))
330 /* Get the indices of the first atom in this charge group */
331 j0 = molt->cgs.index[cg0];
332 j1 = molt->cgs.index[cg0+1];
334 /* Number of atoms in our molecule */
340 "Moltype '%s': there are %d atoms in this charge group\n",
344 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
347 if (nj < 3 || nj > 4)
352 /* Check if we are doing QM on this group */
354 if (qm_grpnr != nullptr)
356 for (j = j0; j < j1 && !qm; j++)
358 qm = (qm_grpnr[j] < qm_grps->nr - 1);
361 /* Cannot use solvent optimization with QM */
367 atom = molt->atoms.atom;
369 /* Still looks like a solvent, time to check parameters */
371 /* If it is perturbed (free energy) we can't use the solvent loops,
372 * so then we just skip to the next molecule.
376 for (j = j0; j < j1 && !perturbed; j++)
378 perturbed = PERTURBED(atom[j]);
386 /* Now it's only a question if the VdW and charge parameters
387 * are OK. Before doing the check we compare and see if they are
388 * identical to a possible previous solvent type.
389 * First we assign the current types and charges.
391 for (j = 0; j < nj; j++)
393 tmp_vdwtype[j] = atom[j0+j].type;
394 tmp_charge[j] = atom[j0+j].q;
397 /* Does it match any previous solvent type? */
398 for (k = 0; k < *n_solvent_parameters; k++)
403 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
404 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
405 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
410 /* Check that types & charges match for all atoms in molecule */
411 for (j = 0; j < nj && match == TRUE; j++)
413 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
417 if (tmp_charge[j] != solvent_parameters[k].charge[j])
424 /* Congratulations! We have a matched solvent.
425 * Flag it with this type for later processing.
428 solvent_parameters[k].count += nmol;
430 /* We are done with this charge group */
435 /* If we get here, we have a tentative new solvent type.
436 * Before we add it we must check that it fulfills the requirements
437 * of the solvent optimized loops. First determine which atoms have
440 for (j = 0; j < nj; j++)
443 tjA = tmp_vdwtype[j];
445 /* Go through all other tpes and see if any have non-zero
446 * VdW parameters when combined with this one.
448 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
450 /* We already checked that the atoms weren't perturbed,
451 * so we only need to check state A now.
455 has_vdw[j] = (has_vdw[j] ||
456 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
457 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
458 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
463 has_vdw[j] = (has_vdw[j] ||
464 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
465 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
470 /* Now we know all we need to make the final check and assignment. */
474 * For this we require thatn all atoms have charge,
475 * the charges on atom 2 & 3 should be the same, and only
476 * atom 1 might have VdW.
478 if (has_vdw[1] == FALSE &&
479 has_vdw[2] == FALSE &&
480 tmp_charge[0] != 0 &&
481 tmp_charge[1] != 0 &&
482 tmp_charge[2] == tmp_charge[1])
484 srenew(solvent_parameters, *n_solvent_parameters+1);
485 solvent_parameters[*n_solvent_parameters].model = esolSPC;
486 solvent_parameters[*n_solvent_parameters].count = nmol;
487 for (k = 0; k < 3; k++)
489 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
490 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
493 *cg_sp = *n_solvent_parameters;
494 (*n_solvent_parameters)++;
499 /* Or could it be a TIP4P?
500 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
501 * Only atom 1 mght have VdW.
503 if (has_vdw[1] == FALSE &&
504 has_vdw[2] == FALSE &&
505 has_vdw[3] == FALSE &&
506 tmp_charge[0] == 0 &&
507 tmp_charge[1] != 0 &&
508 tmp_charge[2] == tmp_charge[1] &&
511 srenew(solvent_parameters, *n_solvent_parameters+1);
512 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
513 solvent_parameters[*n_solvent_parameters].count = nmol;
514 for (k = 0; k < 4; k++)
516 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
517 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
520 *cg_sp = *n_solvent_parameters;
521 (*n_solvent_parameters)++;
525 *solvent_parameters_p = solvent_parameters;
529 check_solvent(FILE * fp,
530 const gmx_mtop_t * mtop,
532 cginfo_mb_t *cginfo_mb)
535 const gmx_moltype_t *molt;
536 int mb, mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
537 int n_solvent_parameters;
538 solvent_parameters_t *solvent_parameters;
544 fprintf(debug, "Going to determine what solvent types we have.\n");
547 n_solvent_parameters = 0;
548 solvent_parameters = nullptr;
549 /* Allocate temporary array for solvent type */
550 snew(cg_sp, mtop->nmolblock);
553 for (mb = 0; mb < mtop->nmolblock; mb++)
555 molt = &mtop->moltype[mtop->molblock[mb].type];
557 /* Here we have to loop over all individual molecules
558 * because we need to check for QMMM particles.
560 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
561 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
562 nmol = mtop->molblock[mb].nmol/nmol_ch;
563 for (mol = 0; mol < nmol_ch; mol++)
566 am = mol*cgs->index[cgs->nr];
567 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
569 check_solvent_cg(molt, cg_mol, nmol,
570 mtop->groups.grpnr[egcQMMM] ?
571 mtop->groups.grpnr[egcQMMM]+at_offset+am : nullptr,
572 &mtop->groups.grps[egcQMMM],
574 &n_solvent_parameters, &solvent_parameters,
575 cginfo_mb[mb].cginfo[cgm+cg_mol],
576 &cg_sp[mb][cgm+cg_mol]);
579 at_offset += cgs->index[cgs->nr];
582 /* Puh! We finished going through all charge groups.
583 * Now find the most common solvent model.
586 /* Most common solvent this far */
588 for (i = 0; i < n_solvent_parameters; i++)
591 solvent_parameters[i].count > solvent_parameters[bestsp].count)
599 bestsol = solvent_parameters[bestsp].model;
607 for (mb = 0; mb < mtop->nmolblock; mb++)
609 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
610 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
611 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
613 if (cg_sp[mb][i] == bestsp)
615 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
620 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
627 if (bestsol != esolNO && fp != nullptr)
629 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
631 solvent_parameters[bestsp].count);
634 sfree(solvent_parameters);
635 fr->solvent_opt = bestsol;
639 acNONE = 0, acCONSTRAINT, acSETTLE
642 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
643 t_forcerec *fr, gmx_bool bNoSolvOpt,
644 gmx_bool *bFEP_NonBonded,
645 gmx_bool *bExcl_IntraCGAll_InterCGNone)
648 const t_blocka *excl;
649 const gmx_moltype_t *molt;
650 const gmx_molblock_t *molb;
651 cginfo_mb_t *cginfo_mb;
654 int cg_offset, a_offset;
655 int mb, m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
659 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
661 snew(cginfo_mb, mtop->nmolblock);
663 snew(type_VDW, fr->ntype);
664 for (ai = 0; ai < fr->ntype; ai++)
666 type_VDW[ai] = FALSE;
667 for (j = 0; j < fr->ntype; j++)
669 type_VDW[ai] = type_VDW[ai] ||
671 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
672 C12(fr->nbfp, fr->ntype, ai, j) != 0;
676 *bFEP_NonBonded = FALSE;
677 *bExcl_IntraCGAll_InterCGNone = TRUE;
680 snew(bExcl, excl_nalloc);
683 for (mb = 0; mb < mtop->nmolblock; mb++)
685 molb = &mtop->molblock[mb];
686 molt = &mtop->moltype[molb->type];
690 /* Check if the cginfo is identical for all molecules in this block.
691 * If so, we only need an array of the size of one molecule.
692 * Otherwise we make an array of #mol times #cgs per molecule.
695 for (m = 0; m < molb->nmol; m++)
697 int am = m*cgs->index[cgs->nr];
698 for (cg = 0; cg < cgs->nr; cg++)
701 a1 = cgs->index[cg+1];
702 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
703 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
707 if (mtop->groups.grpnr[egcQMMM] != nullptr)
709 for (ai = a0; ai < a1; ai++)
711 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
712 mtop->groups.grpnr[egcQMMM][a_offset +ai])
721 cginfo_mb[mb].cg_start = cg_offset;
722 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
723 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
724 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
725 cginfo = cginfo_mb[mb].cginfo;
727 /* Set constraints flags for constrained atoms */
728 snew(a_con, molt->atoms.nr);
729 for (ftype = 0; ftype < F_NRE; ftype++)
731 if (interaction_function[ftype].flags & IF_CONSTRAINT)
736 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
740 for (a = 0; a < nral; a++)
742 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
743 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
749 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
752 int am = m*cgs->index[cgs->nr];
753 for (cg = 0; cg < cgs->nr; cg++)
756 a1 = cgs->index[cg+1];
758 /* Store the energy group in cginfo */
759 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
760 SET_CGINFO_GID(cginfo[cgm+cg], gid);
762 /* Check the intra/inter charge group exclusions */
763 if (a1-a0 > excl_nalloc)
765 excl_nalloc = a1 - a0;
766 srenew(bExcl, excl_nalloc);
768 /* bExclIntraAll: all intra cg interactions excluded
769 * bExclInter: any inter cg interactions excluded
771 bExclIntraAll = TRUE;
775 bHavePerturbedAtoms = FALSE;
776 for (ai = a0; ai < a1; ai++)
778 /* Check VDW and electrostatic interactions */
779 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
780 type_VDW[molt->atoms.atom[ai].typeB]);
781 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
782 molt->atoms.atom[ai].qB != 0);
784 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
786 /* Clear the exclusion list for atom ai */
787 for (aj = a0; aj < a1; aj++)
789 bExcl[aj-a0] = FALSE;
791 /* Loop over all the exclusions of atom ai */
792 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
795 if (aj < a0 || aj >= a1)
804 /* Check if ai excludes a0 to a1 */
805 for (aj = a0; aj < a1; aj++)
809 bExclIntraAll = FALSE;
816 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
819 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
827 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
831 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
833 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
835 /* The size in cginfo is currently only read with DD */
836 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
840 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
844 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
846 if (bHavePerturbedAtoms && fr->efep != efepNO)
848 SET_CGINFO_FEP(cginfo[cgm+cg]);
849 *bFEP_NonBonded = TRUE;
851 /* Store the charge group size */
852 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
854 if (!bExclIntraAll || bExclInter)
856 *bExcl_IntraCGAll_InterCGNone = FALSE;
863 cg_offset += molb->nmol*cgs->nr;
864 a_offset += molb->nmol*cgs->index[cgs->nr];
868 /* the solvent optimizer is called after the QM is initialized,
869 * because we don't want to have the QM subsystemto become an
873 check_solvent(fplog, mtop, fr, cginfo_mb);
875 if (getenv("GMX_NO_SOLV_OPT"))
879 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
880 "Disabling all solvent optimization\n");
882 fr->solvent_opt = esolNO;
886 fr->solvent_opt = esolNO;
888 if (!fr->solvent_opt)
890 for (mb = 0; mb < mtop->nmolblock; mb++)
892 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
894 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
902 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
907 ncg = cgi_mb[nmb-1].cg_end;
910 for (cg = 0; cg < ncg; cg++)
912 while (cg >= cgi_mb[mb].cg_end)
917 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
923 /* Sets the sum of charges (squared) and C6 in the system in fr.
924 * Returns whether the system has a net charge.
926 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
928 /*This now calculates sum for q and c6*/
929 double qsum, q2sum, q, c6sum, c6;
931 const t_atoms *atoms;
936 for (mb = 0; mb < mtop->nmolblock; mb++)
938 nmol = mtop->molblock[mb].nmol;
939 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
940 for (i = 0; i < atoms->nr; i++)
942 q = atoms->atom[i].q;
945 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
950 fr->q2sum[0] = q2sum;
951 fr->c6sum[0] = c6sum;
953 if (fr->efep != efepNO)
958 for (mb = 0; mb < mtop->nmolblock; mb++)
960 nmol = mtop->molblock[mb].nmol;
961 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
962 for (i = 0; i < atoms->nr; i++)
964 q = atoms->atom[i].qB;
967 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
971 fr->q2sum[1] = q2sum;
972 fr->c6sum[1] = c6sum;
977 fr->qsum[1] = fr->qsum[0];
978 fr->q2sum[1] = fr->q2sum[0];
979 fr->c6sum[1] = fr->c6sum[0];
983 if (fr->efep == efepNO)
985 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
989 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
990 fr->qsum[0], fr->qsum[1]);
994 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
995 return (std::abs(fr->qsum[0]) > 1e-4 ||
996 std::abs(fr->qsum[1]) > 1e-4);
999 void update_forcerec(t_forcerec *fr, matrix box)
1001 if (fr->ic->eeltype == eelGRF)
1003 calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
1004 fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
1005 &fr->ic->k_rf, &fr->ic->c_rf);
1009 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
1011 const t_atoms *atoms, *atoms_tpi;
1012 const t_blocka *excl;
1013 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
1014 gmx_int64_t npair, npair_ij, tmpi, tmpj;
1015 double csix, ctwelve;
1016 int ntp, *typecount;
1019 real *nbfp_comb = nullptr;
1025 /* For LJ-PME, we want to correct for the difference between the
1026 * actual C6 values and the C6 values used by the LJ-PME based on
1027 * combination rules. */
1029 if (EVDW_PME(fr->ic->vdwtype))
1031 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1032 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1033 for (tpi = 0; tpi < ntp; ++tpi)
1035 for (tpj = 0; tpj < ntp; ++tpj)
1037 C6(nbfp_comb, ntp, tpi, tpj) =
1038 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1039 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1044 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1052 /* Count the types so we avoid natoms^2 operations */
1053 snew(typecount, ntp);
1054 gmx_mtop_count_atomtypes(mtop, q, typecount);
1056 for (tpi = 0; tpi < ntp; tpi++)
1058 for (tpj = tpi; tpj < ntp; tpj++)
1060 tmpi = typecount[tpi];
1061 tmpj = typecount[tpj];
1064 npair_ij = tmpi*tmpj;
1068 npair_ij = tmpi*(tmpi - 1)/2;
1072 /* nbfp now includes the 6.0 derivative prefactor */
1073 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1077 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1078 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1079 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1085 /* Subtract the excluded pairs.
1086 * The main reason for substracting exclusions is that in some cases
1087 * some combinations might never occur and the parameters could have
1088 * any value. These unused values should not influence the dispersion
1091 for (mb = 0; mb < mtop->nmolblock; mb++)
1093 nmol = mtop->molblock[mb].nmol;
1094 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1095 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1096 for (i = 0; (i < atoms->nr); i++)
1100 tpi = atoms->atom[i].type;
1104 tpi = atoms->atom[i].typeB;
1106 j1 = excl->index[i];
1107 j2 = excl->index[i+1];
1108 for (j = j1; j < j2; j++)
1115 tpj = atoms->atom[k].type;
1119 tpj = atoms->atom[k].typeB;
1123 /* nbfp now includes the 6.0 derivative prefactor */
1124 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1128 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1129 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1130 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1140 /* Only correct for the interaction of the test particle
1141 * with the rest of the system.
1144 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1147 for (mb = 0; mb < mtop->nmolblock; mb++)
1149 nmol = mtop->molblock[mb].nmol;
1150 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1151 for (j = 0; j < atoms->nr; j++)
1154 /* Remove the interaction of the test charge group
1157 if (mb == mtop->nmolblock-1)
1161 if (mb == 0 && nmol == 1)
1163 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1168 tpj = atoms->atom[j].type;
1172 tpj = atoms->atom[j].typeB;
1174 for (i = 0; i < fr->n_tpi; i++)
1178 tpi = atoms_tpi->atom[i].type;
1182 tpi = atoms_tpi->atom[i].typeB;
1186 /* nbfp now includes the 6.0 derivative prefactor */
1187 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1191 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1192 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1193 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1200 if (npair - nexcl <= 0 && fplog)
1202 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1208 csix /= npair - nexcl;
1209 ctwelve /= npair - nexcl;
1213 fprintf(debug, "Counted %d exclusions\n", nexcl);
1214 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1215 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1217 fr->avcsix[q] = csix;
1218 fr->avctwelve[q] = ctwelve;
1221 if (EVDW_PME(fr->ic->vdwtype))
1226 if (fplog != nullptr)
1228 if (fr->eDispCorr == edispcAllEner ||
1229 fr->eDispCorr == edispcAllEnerPres)
1231 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1232 fr->avcsix[0], fr->avctwelve[0]);
1236 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1242 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
1244 const t_atoms *at1, *at2;
1245 int mt1, mt2, i, j, tpi, tpj, ntypes;
1250 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1252 ntypes = mtop->ffparams.atnr;
1255 real bham_b_max = 0;
1256 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1258 at1 = &mtop->moltype[mt1].atoms;
1259 for (i = 0; (i < at1->nr); i++)
1261 tpi = at1->atom[i].type;
1264 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1267 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1269 at2 = &mtop->moltype[mt2].atoms;
1270 for (j = 0; (j < at2->nr); j++)
1272 tpj = at2->atom[j].type;
1275 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1277 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
1282 if ((b < bmin) || (bmin == -1))
1292 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1299 static void make_nbf_tables(FILE *fp,
1300 const interaction_const_t *ic, real rtab,
1301 const char *tabfn, char *eg1, char *eg2,
1307 if (tabfn == nullptr)
1311 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1316 sprintf(buf, "%s", tabfn);
1319 /* Append the two energy group names */
1320 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1321 eg1, eg2, ftp2ext(efXVG));
1323 nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1324 /* Copy the contents of the table to separate coulomb and LJ tables too,
1325 * to improve cache performance.
1327 /* For performance reasons we want
1328 * the table data to be aligned to 16-byte. The pointers could be freed
1329 * but currently aren't.
1331 snew(nbl->table_elec, 1);
1332 nbl->table_elec->interaction = GMX_TABLE_INTERACTION_ELEC;
1333 nbl->table_elec->format = nbl->table_elec_vdw->format;
1334 nbl->table_elec->r = nbl->table_elec_vdw->r;
1335 nbl->table_elec->n = nbl->table_elec_vdw->n;
1336 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1337 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1338 nbl->table_elec->ninteractions = 1;
1339 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1340 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1342 snew(nbl->table_vdw, 1);
1343 nbl->table_vdw->interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1344 nbl->table_vdw->format = nbl->table_elec_vdw->format;
1345 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1346 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1347 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1348 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1349 nbl->table_vdw->ninteractions = 2;
1350 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1351 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1353 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1355 for (j = 0; j < 4; j++)
1357 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1359 for (j = 0; j < 8; j++)
1361 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1366 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1367 * ftype2 present in the topology, build an array of the number of
1368 * interactions present for each bonded interaction index found in the
1371 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1372 * valid type with that parameter.
1374 * \c count will be reallocated as necessary to fit the largest bonded
1375 * interaction index found, and its current size will be returned in
1376 * \c ncount. It will contain zero for every bonded interaction index
1377 * for which no interactions are present in the topology.
1379 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1380 int *ncount, int **count)
1382 const gmx_moltype_t *molt;
1384 int mt, ftype, stride, i, j, tabnr;
1386 // Loop over all moleculetypes
1387 for (mt = 0; mt < mtop->nmoltype; mt++)
1389 molt = &mtop->moltype[mt];
1390 // Loop over all interaction types
1391 for (ftype = 0; ftype < F_NRE; ftype++)
1393 // If the current interaction type is one of the types whose tables we're trying to count...
1394 if (ftype == ftype1 || ftype == ftype2)
1396 il = &molt->ilist[ftype];
1397 stride = 1 + NRAL(ftype);
1398 // ... and there are actually some interactions for this type
1399 for (i = 0; i < il->nr; i += stride)
1401 // Find out which table index the user wanted
1402 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1405 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1407 // Make room for this index in the data structure
1408 if (tabnr >= *ncount)
1410 srenew(*count, tabnr+1);
1411 for (j = *ncount; j < tabnr+1; j++)
1417 // Record that this table index is used and must have a valid file
1425 /*!\brief If there's bonded interactions of flavour \c tabext and type
1426 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1427 * list of filenames passed to mdrun, and make bonded tables from
1430 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1431 * valid type with that parameter.
1433 * A fatal error occurs if no matching filename is found.
1435 static bondedtable_t *make_bonded_tables(FILE *fplog,
1436 int ftype1, int ftype2,
1437 const gmx_mtop_t *mtop,
1438 const t_filenm *tabbfnm,
1448 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1450 // Are there any relevant tabulated bond interactions?
1454 for (int i = 0; i < ncount; i++)
1456 // Do any interactions exist that requires this table?
1459 // This pattern enforces the current requirement that
1460 // table filenames end in a characteristic sequence
1461 // before the file type extension, and avoids table 13
1462 // being recognized and used for table 1.
1463 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1464 bool madeTable = false;
1465 for (int j = 0; j < tabbfnm->nfiles && !madeTable; ++j)
1467 std::string filename(tabbfnm->fns[j]);
1468 if (gmx::endsWith(filename, patternToFind))
1470 // Finally read the table from the file found
1471 tab[i] = make_bonded_table(fplog, tabbfnm->fns[j], NRAL(ftype1)-2);
1477 bool isPlural = (ftype2 != -1);
1478 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.",
1479 interaction_function[ftype1].longname,
1480 isPlural ? "' or '" : "",
1481 isPlural ? interaction_function[ftype2].longname : "",
1483 patternToFind.c_str());
1493 void forcerec_set_ranges(t_forcerec *fr,
1494 int ncg_home, int ncg_force,
1496 int natoms_force_constr, int natoms_f_novirsum)
1501 /* fr->ncg_force is unused in the standard code,
1502 * but it can be useful for modified code dealing with charge groups.
1504 fr->ncg_force = ncg_force;
1505 fr->natoms_force = natoms_force;
1506 fr->natoms_force_constr = natoms_force_constr;
1508 if (fr->natoms_force_constr > fr->nalloc_force)
1510 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1513 if (fr->haveDirectVirialContributions)
1515 fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1519 static real cutoff_inf(real cutoff)
1523 cutoff = GMX_CUTOFF_INF;
1529 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1536 ir->rcoulomb == 0 &&
1538 ir->ePBC == epbcNONE &&
1539 ir->vdwtype == evdwCUT &&
1540 ir->coulombtype == eelCUT &&
1541 ir->efep == efepNO &&
1542 getenv("GMX_NO_ALLVSALL") == nullptr
1545 if (bAllvsAll && ir->opts.ngener > 1)
1547 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";
1553 fprintf(fp, "\n%s\n", note);
1559 if (bAllvsAll && fp && MASTER(cr))
1561 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1568 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
1569 const t_inputrec *ir)
1571 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1573 /* LJ PME with LB combination rule does 7 mesh operations.
1574 * This so slow that we don't compile SIMD non-bonded kernels
1576 GMX_LOG(mdlog.warning).asParagraph().appendText("LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels");
1584 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1587 const gmx_hw_info_t gmx_unused &hardwareInfo)
1589 *kernel_type = nbnxnk4x4_PlainC;
1590 *ewald_excl = ewaldexclTable;
1594 #ifdef GMX_NBNXN_SIMD_4XN
1595 *kernel_type = nbnxnk4xN_SIMD_4xN;
1597 #ifdef GMX_NBNXN_SIMD_2XNN
1598 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1601 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1602 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1603 * This is based on the SIMD acceleration choice and CPU information
1604 * detected at runtime.
1606 * 4xN calculates more (zero) interactions, but has less pair-search
1607 * work and much better kernel instruction scheduling.
1609 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1610 * which doesn't have FMA, both the analytical and tabulated Ewald
1611 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1612 * 2x(4+4) because it results in significantly fewer pairs.
1613 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1614 * 10% with HT, 50% without HT. As we currently don't detect the actual
1615 * use of HT, use 4x8 to avoid a potential performance hit.
1616 * On Intel Haswell 4x8 is always faster.
1618 *kernel_type = nbnxnk4xN_SIMD_4xN;
1620 #if !GMX_SIMD_HAVE_FMA
1621 if (EEL_PME_EWALD(ir->coulombtype) ||
1622 EVDW_PME(ir->vdwtype))
1624 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1625 * There are enough instructions to make 2x(4+4) efficient.
1627 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1630 if (hardwareInfo.haveAmdZenCpu)
1632 /* One 256-bit FMA per cycle makes 2xNN faster */
1633 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1635 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1638 if (getenv("GMX_NBNXN_SIMD_4XN") != nullptr)
1640 #ifdef GMX_NBNXN_SIMD_4XN
1641 *kernel_type = nbnxnk4xN_SIMD_4xN;
1643 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1646 if (getenv("GMX_NBNXN_SIMD_2XNN") != nullptr)
1648 #ifdef GMX_NBNXN_SIMD_2XNN
1649 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1651 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1655 /* Analytical Ewald exclusion correction is only an option in
1657 * Since table lookup's don't parallelize with SIMD, analytical
1658 * will probably always be faster for a SIMD width of 8 or more.
1659 * With FMA analytical is sometimes faster for a width if 4 as well.
1660 * In single precision, this is faster on Bulldozer.
1662 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1663 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE)
1664 /* On AMD Zen, tabulated Ewald kernels are faster on all 4 combinations
1665 * of single or double precision and 128 or 256-bit AVX2.
1667 if (!hardwareInfo.haveAmdZenCpu)
1669 *ewald_excl = ewaldexclAnalytical;
1672 if (getenv("GMX_NBNXN_EWALD_TABLE") != nullptr)
1674 *ewald_excl = ewaldexclTable;
1676 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != nullptr)
1678 *ewald_excl = ewaldexclAnalytical;
1686 const char *lookup_nbnxn_kernel_name(int kernel_type)
1688 const char *returnvalue = nullptr;
1689 switch (kernel_type)
1692 returnvalue = "not set";
1694 case nbnxnk4x4_PlainC:
1695 returnvalue = "plain C";
1697 case nbnxnk4xN_SIMD_4xN:
1698 case nbnxnk4xN_SIMD_2xNN:
1700 returnvalue = "SIMD";
1702 returnvalue = "not available";
1705 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1706 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1710 gmx_fatal(FARGS, "Illegal kernel type selected");
1711 returnvalue = nullptr;
1717 static void pick_nbnxn_kernel(const gmx::MDLogger &mdlog,
1718 gmx_bool use_simd_kernels,
1719 const gmx_hw_info_t &hardwareInfo,
1721 EmulateGpuNonbonded emulateGpu,
1722 const t_inputrec *ir,
1725 gmx_bool bDoNonbonded)
1727 assert(kernel_type);
1729 *kernel_type = nbnxnkNotSet;
1730 *ewald_excl = ewaldexclTable;
1732 if (emulateGpu == EmulateGpuNonbonded::Yes)
1734 *kernel_type = nbnxnk8x8x8_PlainC;
1738 GMX_LOG(mdlog.warning).asParagraph().appendText("Emulating a GPU run on the CPU (slow)");
1743 *kernel_type = nbnxnk8x8x8_GPU;
1746 if (*kernel_type == nbnxnkNotSet)
1748 if (use_simd_kernels &&
1749 nbnxn_simd_supported(mdlog, ir))
1751 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl, hardwareInfo);
1755 *kernel_type = nbnxnk4x4_PlainC;
1761 GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
1762 "Using %s %dx%d nonbonded short-range kernels",
1763 lookup_nbnxn_kernel_name(*kernel_type),
1764 nbnxn_kernel_to_cluster_i_size(*kernel_type),
1765 nbnxn_kernel_to_cluster_j_size(*kernel_type));
1767 if (nbnxnk4x4_PlainC == *kernel_type ||
1768 nbnxnk8x8x8_PlainC == *kernel_type)
1770 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1771 "WARNING: Using the slow %s kernels. This should\n"
1772 "not happen during routine usage on supported platforms.",
1773 lookup_nbnxn_kernel_name(*kernel_type));
1778 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1779 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1780 bool systemHasNetCharge,
1781 interaction_const_t *ic)
1783 if (!EEL_PME_EWALD(ir->coulombtype))
1790 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1792 if (ir->coulombtype == eelP3M_AD)
1794 please_cite(fp, "Hockney1988");
1795 please_cite(fp, "Ballenegger2012");
1799 please_cite(fp, "Essmann95a");
1802 if (ir->ewald_geometry == eewg3DC)
1806 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1807 systemHasNetCharge ? " and net charge" : "");
1809 please_cite(fp, "In-Chul99a");
1810 if (systemHasNetCharge)
1812 please_cite(fp, "Ballenegger2009");
1817 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1820 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1821 1/ic->ewaldcoeff_q);
1824 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1826 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1827 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1835 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1836 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1837 interaction_const_t *ic)
1839 if (!EVDW_PME(ir->vdwtype))
1846 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1847 please_cite(fp, "Essmann95a");
1849 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1852 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1853 1/ic->ewaldcoeff_lj);
1856 if (ic->vdw_modifier == eintmodPOTSHIFT)
1858 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1859 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1863 ic->sh_lj_ewald = 0;
1867 gmx_bool uses_simple_tables(int cutoff_scheme,
1868 nonbonded_verlet_t *nbv,
1871 gmx_bool bUsesSimpleTables = TRUE;
1874 switch (cutoff_scheme)
1877 bUsesSimpleTables = TRUE;
1880 assert(NULL != nbv && NULL != nbv->grp);
1881 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1882 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1885 gmx_incons("unimplemented");
1887 return bUsesSimpleTables;
1890 static void init_ewald_f_table(interaction_const_t *ic,
1895 /* Get the Ewald table spacing based on Coulomb and/or LJ
1896 * Ewald coefficients and rtol.
1898 ic->tabq_scale = ewald_spline3_table_scale(ic);
1900 if (ic->cutoff_scheme == ecutsVERLET)
1902 maxr = ic->rcoulomb;
1906 maxr = std::max(ic->rcoulomb, rtab);
1908 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1910 sfree_aligned(ic->tabq_coul_FDV0);
1911 sfree_aligned(ic->tabq_coul_F);
1912 sfree_aligned(ic->tabq_coul_V);
1914 sfree_aligned(ic->tabq_vdw_FDV0);
1915 sfree_aligned(ic->tabq_vdw_F);
1916 sfree_aligned(ic->tabq_vdw_V);
1918 if (EEL_PME_EWALD(ic->eeltype))
1920 /* Create the original table data in FDV0 */
1921 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1922 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1923 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1924 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1925 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1928 if (EVDW_PME(ic->vdwtype))
1930 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1931 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1932 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1933 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1934 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1938 void init_interaction_const_tables(FILE *fp,
1939 interaction_const_t *ic,
1942 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1944 init_ewald_f_table(ic, rtab);
1948 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1949 1/ic->tabq_scale, ic->tabq_size);
1954 static void clear_force_switch_constants(shift_consts_t *sc)
1961 static void force_switch_constants(real p,
1965 /* Here we determine the coefficient for shifting the force to zero
1966 * between distance rsw and the cut-off rc.
1967 * For a potential of r^-p, we have force p*r^-(p+1).
1968 * But to save flops we absorb p in the coefficient.
1970 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1971 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1973 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1974 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1975 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1978 static void potential_switch_constants(real rsw, real rc,
1979 switch_consts_t *sc)
1981 /* The switch function is 1 at rsw and 0 at rc.
1982 * The derivative and second derivate are zero at both ends.
1983 * rsw = max(r - r_switch, 0)
1984 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1985 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1986 * force = force*dsw - potential*sw
1989 sc->c3 = -10/gmx::power3(rc - rsw);
1990 sc->c4 = 15/gmx::power4(rc - rsw);
1991 sc->c5 = -6/gmx::power5(rc - rsw);
1994 /*! \brief Construct interaction constants
1996 * This data is used (particularly) by search and force code for
1997 * short-range interactions. Many of these are constant for the whole
1998 * simulation; some are constant only after PME tuning completes.
2001 init_interaction_const(FILE *fp,
2002 interaction_const_t **interaction_const,
2003 const t_inputrec *ir,
2004 const gmx_mtop_t *mtop,
2005 bool systemHasNetCharge)
2007 interaction_const_t *ic;
2011 ic->cutoff_scheme = ir->cutoff_scheme;
2013 /* Just allocate something so we can free it */
2014 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
2015 snew_aligned(ic->tabq_coul_F, 16, 32);
2016 snew_aligned(ic->tabq_coul_V, 16, 32);
2019 ic->vdwtype = ir->vdwtype;
2020 ic->vdw_modifier = ir->vdw_modifier;
2021 ic->reppow = mtop->ffparams.reppow;
2022 ic->rvdw = cutoff_inf(ir->rvdw);
2023 ic->rvdw_switch = ir->rvdw_switch;
2024 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
2025 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
2026 if (ic->useBuckingham)
2028 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
2031 initVdwEwaldParameters(fp, ir, ic);
2033 clear_force_switch_constants(&ic->dispersion_shift);
2034 clear_force_switch_constants(&ic->repulsion_shift);
2036 switch (ic->vdw_modifier)
2038 case eintmodPOTSHIFT:
2039 /* Only shift the potential, don't touch the force */
2040 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
2041 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
2043 case eintmodFORCESWITCH:
2044 /* Switch the force, switch and shift the potential */
2045 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2046 &ic->dispersion_shift);
2047 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2048 &ic->repulsion_shift);
2050 case eintmodPOTSWITCH:
2051 /* Switch the potential and force */
2052 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2056 case eintmodEXACTCUTOFF:
2057 /* Nothing to do here */
2060 gmx_incons("unimplemented potential modifier");
2063 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2065 /* Electrostatics */
2066 ic->eeltype = ir->coulombtype;
2067 ic->coulomb_modifier = ir->coulomb_modifier;
2068 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
2069 ic->rcoulomb_switch = ir->rcoulomb_switch;
2070 ic->epsilon_r = ir->epsilon_r;
2072 /* Set the Coulomb energy conversion factor */
2073 if (ic->epsilon_r != 0)
2075 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
2079 /* eps = 0 is infinite dieletric: no Coulomb interactions */
2083 /* Reaction-field */
2084 if (EEL_RF(ic->eeltype))
2086 ic->epsilon_rf = ir->epsilon_rf;
2087 /* Generalized reaction field parameters are updated every step */
2088 if (ic->eeltype != eelGRF)
2090 calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
2091 ic->rcoulomb, 0, 0, NULL,
2092 &ic->k_rf, &ic->c_rf);
2095 if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
2097 /* grompp should have done this, but this scheme is obsolete */
2098 ic->coulomb_modifier = eintmodEXACTCUTOFF;
2103 /* For plain cut-off we might use the reaction-field kernels */
2104 ic->epsilon_rf = ic->epsilon_r;
2106 if (ir->coulomb_modifier == eintmodPOTSHIFT)
2108 ic->c_rf = 1/ic->rcoulomb;
2116 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
2120 real dispersion_shift;
2122 dispersion_shift = ic->dispersion_shift.cpot;
2123 if (EVDW_PME(ic->vdwtype))
2125 dispersion_shift -= ic->sh_lj_ewald;
2127 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2128 ic->repulsion_shift.cpot, dispersion_shift);
2130 if (ic->eeltype == eelCUT)
2132 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2134 else if (EEL_PME(ic->eeltype))
2136 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2141 *interaction_const = ic;
2144 static void init_nb_verlet(const gmx::MDLogger &mdlog,
2145 nonbonded_verlet_t **nb_verlet,
2146 gmx_bool bFEP_NonBonded,
2147 const t_inputrec *ir,
2148 const t_forcerec *fr,
2149 const t_commrec *cr,
2150 const gmx_hw_info_t &hardwareInfo,
2151 const gmx_device_info_t *deviceInfo,
2152 const gmx_mtop_t *mtop,
2155 nonbonded_verlet_t *nbv;
2158 nbnxn_alloc_t *nb_alloc;
2159 nbnxn_free_t *nb_free;
2161 nbv = new nonbonded_verlet_t();
2163 nbv->emulateGpu = ((getenv("GMX_EMULATE_GPU") != nullptr) ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
2164 nbv->bUseGPU = deviceInfo != nullptr;
2166 GMX_RELEASE_ASSERT(!(nbv->emulateGpu == EmulateGpuNonbonded::Yes && nbv->bUseGPU), "When GPU emulation is active, there cannot be a GPU assignment");
2169 nbv->min_ci_balanced = 0;
2171 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2172 for (int i = 0; i < nbv->ngrp; i++)
2174 nbv->grp[i].nbl_lists.nnbl = 0;
2175 nbv->grp[i].kernel_type = nbnxnkNotSet;
2177 if (i == 0) /* local */
2179 pick_nbnxn_kernel(mdlog, fr->use_simd_kernels, hardwareInfo,
2180 nbv->bUseGPU, nbv->emulateGpu, ir,
2181 &nbv->grp[i].kernel_type,
2182 &nbv->grp[i].ewald_excl,
2185 else /* non-local */
2187 /* Use the same kernel for local and non-local interactions */
2188 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2189 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2193 nbv->listParams = std::unique_ptr<NbnxnListParameters>(new NbnxnListParameters(ir->rlist));
2194 setupDynamicPairlistPruning(mdlog, ir, mtop, box, nbv->grp[0].kernel_type, fr->ic,
2195 nbv->listParams.get());
2197 nbnxn_init_search(&nbv->nbs,
2198 DOMAINDECOMP(cr) ? &cr->dd->nc : nullptr,
2199 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
2201 gmx_omp_nthreads_get(emntPairsearch));
2203 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2204 &nb_alloc, &nb_free);
2206 for (int i = 0; i < nbv->ngrp; i++)
2208 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2209 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2210 /* 8x8x8 "non-simple" lists are ATM always combined */
2211 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2215 int enbnxninitcombrule;
2216 if (fr->ic->vdwtype == evdwCUT &&
2217 (fr->ic->vdw_modifier == eintmodNONE ||
2218 fr->ic->vdw_modifier == eintmodPOTSHIFT) &&
2219 getenv("GMX_NO_LJ_COMB_RULE") == nullptr)
2221 /* Plain LJ cut-off: we can optimize with combination rules */
2222 enbnxninitcombrule = enbnxninitcombruleDETECT;
2224 else if (fr->ic->vdwtype == evdwPME)
2226 /* LJ-PME: we need to use a combination rule for the grid */
2227 if (fr->ljpme_combination_rule == eljpmeGEOM)
2229 enbnxninitcombrule = enbnxninitcombruleGEOM;
2233 enbnxninitcombrule = enbnxninitcombruleLB;
2238 /* We use a full combination matrix: no rule required */
2239 enbnxninitcombrule = enbnxninitcombruleNONE;
2243 int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
2244 if (ir->opts.ngener - ir->nwall == 1)
2246 /* We have only one non-wall energy group, we do not need energy group
2247 * support in the non-bondeds kernels, since all non-bonded energy
2248 * contributions go to the first element of the energy group matrix.
2250 mimimumNumEnergyGroupNonbonded = 1;
2252 bool bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[0].kernel_type);
2253 nbnxn_atomdata_init(mdlog,
2255 nbv->grp[0].kernel_type,
2257 fr->ntype, fr->nbfp,
2258 mimimumNumEnergyGroupNonbonded,
2259 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2264 /* init the NxN GPU data; the last argument tells whether we'll have
2265 * both local and non-local NB calculation on GPU */
2266 nbnxn_gpu_init(&nbv->gpu_nbv,
2269 nbv->listParams.get(),
2274 // TODO: with the texture reference support removed, this barrier is
2275 // in principle not needed. Remove now or do it in a follow-up?
2276 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2277 * also sharing texture references. To keep the code simple, we don't
2278 * treat texture references as shared resources, but this means that
2279 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2280 * Hence, to ensure that the non-bonded kernels don't start before all
2281 * texture binding operations are finished, we need to wait for all ranks
2282 * to arrive here before continuing.
2284 * Note that we could omit this barrier if GPUs are not shared (or
2285 * texture objects are used), but as this is initialization code, there
2286 * is no point in complicating things.
2293 #endif /* GMX_THREAD_MPI */
2295 if ((env = getenv("GMX_NB_MIN_CI")) != nullptr)
2299 nbv->min_ci_balanced = strtol(env, &end, 10);
2300 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
2302 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
2307 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2308 nbv->min_ci_balanced);
2313 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2316 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2317 nbv->min_ci_balanced);
2326 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2328 return nbv != nullptr && nbv->bUseGPU;
2331 void init_forcerec(FILE *fp,
2332 const gmx::MDLogger &mdlog,
2335 const t_inputrec *ir,
2336 const gmx_mtop_t *mtop,
2337 const t_commrec *cr,
2341 const t_filenm *tabbfnm,
2342 const gmx_hw_info_t &hardwareInfo,
2343 const gmx_device_info_t *deviceInfo,
2344 gmx_bool bNoSolvOpt,
2347 int m, negp_pp, negptable, egi, egj;
2352 gmx_bool bGenericKernelOnly;
2353 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
2354 gmx_bool bFEP_NonBonded;
2355 int *nm_ind, egp_flags;
2357 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2358 fr->use_simd_kernels = TRUE;
2360 fr->bDomDec = DOMAINDECOMP(cr);
2362 if (check_box(ir->ePBC, box))
2364 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2367 /* Test particle insertion ? */
2370 /* Set to the size of the molecule to be inserted (the last one) */
2371 /* Because of old style topologies, we have to use the last cg
2372 * instead of the last molecule type.
2374 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2375 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2376 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2378 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2386 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
2388 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
2389 eel_names[ir->coulombtype]);
2394 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
2396 if (ir->useTwinRange)
2398 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
2400 /* Copy the user determined parameters */
2401 fr->userint1 = ir->userint1;
2402 fr->userint2 = ir->userint2;
2403 fr->userint3 = ir->userint3;
2404 fr->userint4 = ir->userint4;
2405 fr->userreal1 = ir->userreal1;
2406 fr->userreal2 = ir->userreal2;
2407 fr->userreal3 = ir->userreal3;
2408 fr->userreal4 = ir->userreal4;
2411 fr->fc_stepsize = ir->fc_stepsize;
2414 fr->efep = ir->efep;
2415 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2416 if (ir->fepvals->bScCoul)
2418 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2419 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
2423 fr->sc_alphacoul = 0;
2424 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2426 fr->sc_power = ir->fepvals->sc_power;
2427 fr->sc_r_power = ir->fepvals->sc_r_power;
2428 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
2430 env = getenv("GMX_SCSIGMA_MIN");
2434 sscanf(env, "%20lf", &dbl);
2435 fr->sc_sigma6_min = gmx::power6(dbl);
2438 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2442 fr->bNonbonded = TRUE;
2443 if (getenv("GMX_NO_NONBONDED") != nullptr)
2445 /* turn off non-bonded calculations */
2446 fr->bNonbonded = FALSE;
2447 GMX_LOG(mdlog.warning).asParagraph().appendText(
2448 "Found environment variable GMX_NO_NONBONDED.\n"
2449 "Disabling nonbonded calculations.");
2452 bGenericKernelOnly = FALSE;
2454 /* We now check in the NS code whether a particular combination of interactions
2455 * can be used with water optimization, and disable it if that is not the case.
2458 if (getenv("GMX_NB_GENERIC") != nullptr)
2463 "Found environment variable GMX_NB_GENERIC.\n"
2464 "Disabling all interaction-specific nonbonded kernels, will only\n"
2465 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2467 bGenericKernelOnly = TRUE;
2470 if (bGenericKernelOnly == TRUE)
2475 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
2477 fr->use_simd_kernels = FALSE;
2481 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2482 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2483 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2487 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2489 /* Check if we can/should do all-vs-all kernels */
2490 fr->bAllvsAll = can_use_allvsall(ir, FALSE, nullptr, nullptr);
2491 fr->AllvsAll_work = nullptr;
2493 /* All-vs-all kernels have not been implemented in 4.6 and later.
2494 * See Redmine #1249. */
2497 fr->bAllvsAll = FALSE;
2501 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2502 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2503 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2504 "or try cutoff-scheme = Verlet.\n\n");
2508 /* Neighbour searching stuff */
2509 fr->cutoff_scheme = ir->cutoff_scheme;
2510 fr->bGrid = (ir->ns_type == ensGRID);
2511 fr->ePBC = ir->ePBC;
2513 if (fr->cutoff_scheme == ecutsGROUP)
2515 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2516 "removed in a future release when 'verlet' supports all interaction forms.\n";
2520 fprintf(stderr, "\n%s\n", note);
2524 fprintf(fp, "\n%s\n", note);
2528 /* Determine if we will do PBC for distances in bonded interactions */
2529 if (fr->ePBC == epbcNONE)
2531 fr->bMolPBC = FALSE;
2535 if (!DOMAINDECOMP(cr))
2539 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2540 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2541 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2543 /* The group cut-off scheme and SHAKE assume charge groups
2544 * are whole, but not using molpbc is faster in most cases.
2545 * With intermolecular interactions we need PBC for calculating
2546 * distances between atoms in different molecules.
2548 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2549 !mtop->bIntermolecularInteractions)
2551 fr->bMolPBC = ir->bPeriodicMols;
2553 if (bSHAKE && fr->bMolPBC)
2555 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2560 /* Not making molecules whole is faster in most cases,
2561 * but With orientation restraints we need whole molecules.
2563 fr->bMolPBC = (fcd->orires.nr == 0);
2565 if (getenv("GMX_USE_GRAPH") != nullptr)
2567 fr->bMolPBC = FALSE;
2570 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
2573 if (mtop->bIntermolecularInteractions)
2575 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
2579 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
2581 if (bSHAKE && fr->bMolPBC)
2583 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");
2589 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2593 fr->rc_scaling = ir->refcoord_scaling;
2594 copy_rvec(ir->posres_com, fr->posres_com);
2595 copy_rvec(ir->posres_comB, fr->posres_comB);
2596 fr->rlist = cutoff_inf(ir->rlist);
2597 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2599 /* This now calculates sum for q and c6*/
2600 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
2602 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2603 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
2604 init_interaction_const_tables(fp, fr->ic, ir->rlist + ir->tabext);
2606 const interaction_const_t *ic = fr->ic;
2608 /* TODO: Replace this Ewald table or move it into interaction_const_t */
2609 if (ir->coulombtype == eelEWALD)
2611 init_ewald_tab(&(fr->ewald_table), ir, fp);
2614 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2615 switch (ic->eeltype)
2618 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
2623 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2627 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2628 GMX_RELEASE_ASSERT(ic->coulomb_modifier == eintmodEXACTCUTOFF, "With the group scheme RF-zero needs the exact cut-off modifier");
2637 case eelPMEUSERSWITCH:
2638 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2644 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2648 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
2651 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
2653 /* Vdw: Translate from mdp settings to kernel format */
2654 switch (ic->vdwtype)
2659 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2663 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2667 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2673 case evdwENCADSHIFT:
2674 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2678 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
2681 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
2683 if (ir->cutoff_scheme == ecutsGROUP)
2685 fr->bvdwtab = ((ic->vdwtype != evdwCUT || !gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2686 && !EVDW_PME(ic->vdwtype));
2687 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2688 fr->bcoultab = !(ic->eeltype == eelCUT ||
2689 ic->eeltype == eelEWALD ||
2690 ic->eeltype == eelPME ||
2691 ic->eeltype == eelP3M_AD ||
2692 ic->eeltype == eelRF ||
2693 ic->eeltype == eelRF_ZERO);
2695 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2696 * going to be faster to tabulate the interaction than calling the generic kernel.
2697 * However, if generic kernels have been requested we keep things analytically.
2699 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2700 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2701 bGenericKernelOnly == FALSE)
2703 if ((ic->rcoulomb_switch != ic->rvdw_switch) || (ic->rcoulomb != ic->rvdw))
2705 fr->bcoultab = TRUE;
2706 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2707 * which would otherwise need two tables.
2711 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2712 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2713 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2714 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2716 if ((ic->rcoulomb != ic->rvdw) && (bGenericKernelOnly == FALSE))
2718 fr->bcoultab = TRUE;
2722 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2724 fr->bcoultab = TRUE;
2726 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2731 if (getenv("GMX_REQUIRE_TABLES"))
2734 fr->bcoultab = TRUE;
2739 fprintf(fp, "Table routines are used for coulomb: %s\n",
2740 gmx::boolToString(fr->bcoultab));
2741 fprintf(fp, "Table routines are used for vdw: %s\n",
2742 gmx::boolToString(fr->bvdwtab));
2745 if (fr->bvdwtab == TRUE)
2747 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2748 fr->nbkernel_vdw_modifier = eintmodNONE;
2750 if (fr->bcoultab == TRUE)
2752 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2753 fr->nbkernel_elec_modifier = eintmodNONE;
2757 if (ir->cutoff_scheme == ecutsVERLET)
2759 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2761 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2763 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2764 * while mdrun does not (and never did) support this.
2766 if (EEL_USER(fr->ic->eeltype))
2768 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2769 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2772 fr->bvdwtab = FALSE;
2773 fr->bcoultab = FALSE;
2776 /* 1-4 interaction electrostatics */
2777 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2779 /* Parameters for generalized RF */
2783 if (ic->eeltype == eelGRF)
2785 init_generalized_rf(fp, mtop, ir, fr);
2788 fr->haveDirectVirialContributions =
2789 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
2790 fr->forceProviders->hasForceProvider() ||
2791 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2792 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2797 if (fr->haveDirectVirialContributions)
2799 fr->forceBufferForDirectVirialContributions = new std::vector<gmx::RVec>;
2802 if (fr->cutoff_scheme == ecutsGROUP &&
2803 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2805 /* Count the total number of charge groups */
2806 fr->cg_nalloc = ncg_mtop(mtop);
2807 srenew(fr->cg_cm, fr->cg_nalloc);
2809 if (fr->shift_vec == nullptr)
2811 snew(fr->shift_vec, SHIFTS);
2814 if (fr->fshift == nullptr)
2816 snew(fr->fshift, SHIFTS);
2819 if (fr->nbfp == nullptr)
2821 fr->ntype = mtop->ffparams.atnr;
2822 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2823 if (EVDW_PME(ic->vdwtype))
2825 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2829 /* Copy the energy group exclusions */
2830 fr->egp_flags = ir->opts.egp_flags;
2832 /* Van der Waals stuff */
2833 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
2835 if (ic->rvdw_switch >= ic->rvdw)
2837 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2838 ic->rvdw_switch, ic->rvdw);
2842 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2843 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
2844 ic->rvdw_switch, ic->rvdw);
2848 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
2850 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2853 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
2855 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2858 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2860 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2863 if (fp && fr->cutoff_scheme == ecutsGROUP)
2865 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2866 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
2869 fr->eDispCorr = ir->eDispCorr;
2870 fr->numAtomsForDispersionCorrection = mtop->natoms;
2871 if (ir->eDispCorr != edispcNO)
2873 set_avcsixtwelve(fp, fr, mtop);
2876 if (ir->implicit_solvent)
2878 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
2881 /* Construct tables for the group scheme. A little unnecessary to
2882 * make both vdw and coul tables sometimes, but what the
2883 * heck. Note that both cutoff schemes construct Ewald tables in
2884 * init_interaction_const_tables. */
2885 needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
2886 (fr->bcoultab || fr->bvdwtab));
2888 negp_pp = ir->opts.ngener - ir->nwall;
2890 if (!needGroupSchemeTables)
2892 bSomeNormalNbListsAreInUse = TRUE;
2897 bSomeNormalNbListsAreInUse = FALSE;
2898 for (egi = 0; egi < negp_pp; egi++)
2900 for (egj = egi; egj < negp_pp; egj++)
2902 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2903 if (!(egp_flags & EGP_EXCL))
2905 if (egp_flags & EGP_TABLE)
2911 bSomeNormalNbListsAreInUse = TRUE;
2916 if (bSomeNormalNbListsAreInUse)
2918 fr->nnblists = negptable + 1;
2922 fr->nnblists = negptable;
2924 if (fr->nnblists > 1)
2926 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
2930 snew(fr->nblists, fr->nnblists);
2932 /* This code automatically gives table length tabext without cut-off's,
2933 * in that case grompp should already have checked that we do not need
2934 * normal tables and we only generate tables for 1-4 interactions.
2936 rtab = ir->rlist + ir->tabext;
2938 if (needGroupSchemeTables)
2940 /* make tables for ordinary interactions */
2941 if (bSomeNormalNbListsAreInUse)
2943 make_nbf_tables(fp, ic, rtab, tabfn, nullptr, nullptr, &fr->nblists[0]);
2952 /* Read the special tables for certain energy group pairs */
2953 nm_ind = mtop->groups.grps[egcENER].nm_ind;
2954 for (egi = 0; egi < negp_pp; egi++)
2956 for (egj = egi; egj < negp_pp; egj++)
2958 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2959 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2961 if (fr->nnblists > 1)
2963 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
2965 /* Read the table file with the two energy groups names appended */
2966 make_nbf_tables(fp, ic, rtab, tabfn,
2967 *mtop->groups.grpname[nm_ind[egi]],
2968 *mtop->groups.grpname[nm_ind[egj]],
2972 else if (fr->nnblists > 1)
2974 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
2981 /* Tables might not be used for the potential modifier
2982 * interactions per se, but we still need them to evaluate
2983 * switch/shift dispersion corrections in this case. */
2984 if (fr->eDispCorr != edispcNO)
2986 fr->dispersionCorrectionTable = makeDispersionCorrectionTable(fp, ic, rtab, tabfn);
2989 /* We want to use unmodified tables for 1-4 coulombic
2990 * interactions, so we must in general have an extra set of
2992 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2993 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2994 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
2996 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
2997 GMX_MAKETABLES_14ONLY);
3001 fr->nwall = ir->nwall;
3002 if (ir->nwall && ir->wall_type == ewtTABLE)
3004 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
3009 // Need to catch std::bad_alloc
3010 // TODO Don't need to catch this here, when merging with master branch
3013 fcd->bondtab = make_bonded_tables(fp,
3014 F_TABBONDS, F_TABBONDSNC,
3015 mtop, tabbfnm, "b");
3016 fcd->angletab = make_bonded_tables(fp,
3018 mtop, tabbfnm, "a");
3019 fcd->dihtab = make_bonded_tables(fp,
3021 mtop, tabbfnm, "d");
3023 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
3029 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3033 /* QM/MM initialization if requested
3037 fprintf(stderr, "QM/MM calculation requested.\n");
3040 fr->bQMMM = ir->bQMMM;
3041 fr->qr = mk_QMMMrec();
3043 /* Set all the static charge group info */
3044 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3046 &fr->bExcl_IntraCGAll_InterCGNone);
3047 if (DOMAINDECOMP(cr))
3049 fr->cginfo = nullptr;
3053 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3056 if (!DOMAINDECOMP(cr))
3058 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3059 mtop->natoms, mtop->natoms, mtop->natoms);
3062 fr->print_force = print_force;
3065 /* coarse load balancing vars */
3070 /* Initialize neighbor search */
3072 init_ns(fp, cr, fr->ns, fr, mtop);
3074 if (thisRankHasDuty(cr, DUTY_PP))
3076 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3079 /* Initialize the thread working data for bonded interactions */
3080 init_bonded_threading(fp, mtop->groups.grps[egcENER].nr,
3081 &fr->bonded_threading);
3083 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
3084 snew(fr->ewc_t, fr->nthread_ewc);
3086 if (fr->cutoff_scheme == ecutsVERLET)
3088 // We checked the cut-offs in grompp, but double-check here.
3089 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
3090 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
3092 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
3096 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
3099 init_nb_verlet(mdlog, &fr->nbv, bFEP_NonBonded, ir, fr,
3100 cr, hardwareInfo, deviceInfo,
3106 /* Here we switch from using mdlog, which prints the newline before
3107 * the paragraph, to our old fprintf logging, which prints the newline
3108 * after the paragraph, so we should add a newline here.
3113 if (ir->eDispCorr != edispcNO)
3115 calc_enervirdiff(fp, ir->eDispCorr, fr);
3119 /* Frees GPU memory and sets a tMPI node barrier.
3121 * Note that this function needs to be called even if GPUs are not used
3122 * in this run because the PME ranks have no knowledge of whether GPUs
3123 * are used or not, but all ranks need to enter the barrier below.
3124 * \todo Remove physical node barrier from this function after making sure
3125 * that it's not needed anymore (with a shared GPU run).
3127 void free_gpu_resources(const t_forcerec *fr,
3128 const t_commrec *cr)
3130 bool isPPrankUsingGPU = fr && fr->nbv && fr->nbv->bUseGPU;
3132 /* stop the GPU profiler (only CUDA) */
3135 if (isPPrankUsingGPU)
3137 /* free nbnxn data in GPU memory */
3138 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3141 /* With tMPI we need to wait for all ranks to finish deallocation before
3142 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
3145 * This is not a concern in OpenCL where we use one context per rank which
3146 * is freed in nbnxn_gpu_free().
3148 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
3149 * but it is easier and more futureproof to call it on the whole node.
3151 if (GMX_THREAD_MPI && (PAR(cr) || isMultiSim(cr->ms)))
3153 gmx_barrier_physical_node(cr);