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,2019, 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/gpubonded.h"
62 #include "gromacs/listed_forces/manage_threading.h"
63 #include "gromacs/listed_forces/pairs.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/dispersioncorrection.h"
69 #include "gromacs/mdlib/force.h"
70 #include "gromacs/mdlib/forcerec_threading.h"
71 #include "gromacs/mdlib/gmx_omp_nthreads.h"
72 #include "gromacs/mdlib/md_support.h"
73 #include "gromacs/mdlib/qmmm.h"
74 #include "gromacs/mdlib/rf_util.h"
75 #include "gromacs/mdlib/wall.h"
76 #include "gromacs/mdtypes/commrec.h"
77 #include "gromacs/mdtypes/fcdata.h"
78 #include "gromacs/mdtypes/group.h"
79 #include "gromacs/mdtypes/iforceprovider.h"
80 #include "gromacs/mdtypes/inputrec.h"
81 #include "gromacs/mdtypes/md_enums.h"
82 #include "gromacs/nbnxm/gpu_data_mgmt.h"
83 #include "gromacs/nbnxm/nbnxm.h"
84 #include "gromacs/nbnxm/nbnxm_geometry.h"
85 #include "gromacs/pbcutil/ishift.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/tables/forcetable.h"
88 #include "gromacs/topology/mtop_util.h"
89 #include "gromacs/trajectory/trajectoryframe.h"
90 #include "gromacs/utility/cstringutil.h"
91 #include "gromacs/utility/exceptions.h"
92 #include "gromacs/utility/fatalerror.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/logger.h"
95 #include "gromacs/utility/physicalnodecommunicator.h"
96 #include "gromacs/utility/pleasecite.h"
97 #include "gromacs/utility/smalloc.h"
98 #include "gromacs/utility/strconvert.h"
100 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
108 snew(nbfp, 3*atnr*atnr);
109 for (i = k = 0; (i < atnr); i++)
111 for (j = 0; (j < atnr); j++, k++)
113 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
114 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
115 /* nbfp now includes the 6.0 derivative prefactor */
116 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
122 snew(nbfp, 2*atnr*atnr);
123 for (i = k = 0; (i < atnr); i++)
125 for (j = 0; (j < atnr); j++, k++)
127 /* nbfp now includes the 6.0/12.0 derivative prefactors */
128 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
129 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
137 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
140 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
143 /* For LJ-PME simulations, we correct the energies with the reciprocal space
144 * inside of the cut-off. To do this the non-bonded kernels needs to have
145 * access to the C6-values used on the reciprocal grid in pme.c
149 snew(grid, 2*atnr*atnr);
150 for (i = k = 0; (i < atnr); i++)
152 for (j = 0; (j < atnr); j++, k++)
154 c6i = idef->iparams[i*(atnr+1)].lj.c6;
155 c12i = idef->iparams[i*(atnr+1)].lj.c12;
156 c6j = idef->iparams[j*(atnr+1)].lj.c6;
157 c12j = idef->iparams[j*(atnr+1)].lj.c12;
158 c6 = std::sqrt(c6i * c6j);
159 if (fr->ljpme_combination_rule == eljpmeLB
160 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
162 sigmai = gmx::sixthroot(c12i / c6i);
163 sigmaj = gmx::sixthroot(c12j / c6j);
164 epsi = c6i * c6i / c12i;
165 epsj = c6j * c6j / c12j;
166 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
168 /* Store the elements at the same relative positions as C6 in nbfp in order
169 * to simplify access in the kernels
171 grid[2*(atnr*i+j)] = c6*6.0;
177 /* This routine sets fr->solvent_opt to the most common solvent in the
178 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
179 * the fr->solvent_type array with the correct type (or esolNO).
181 * Charge groups that fulfill the conditions but are not identical to the
182 * most common one will be marked as esolNO in the solvent_type array.
184 * TIP3p is identical to SPC for these purposes, so we call it
185 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
187 * NOTE: QM particle should not
188 * become an optimized solvent. Not even if there is only one charge
198 } solvent_parameters_t;
201 check_solvent_cg(const gmx_moltype_t *molt,
204 const unsigned char *qm_grpnr,
205 const t_grps *qm_grps,
207 int *n_solvent_parameters,
208 solvent_parameters_t **solvent_parameters_p,
218 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc 7 happy */
219 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc 7 happy */
222 solvent_parameters_t *solvent_parameters;
224 /* We use a list with parameters for each solvent type.
225 * Every time we discover a new molecule that fulfills the basic
226 * conditions for a solvent we compare with the previous entries
227 * in these lists. If the parameters are the same we just increment
228 * the counter for that type, and otherwise we create a new type
229 * based on the current molecule.
231 * Once we've finished going through all molecules we check which
232 * solvent is most common, and mark all those molecules while we
233 * clear the flag on all others.
236 solvent_parameters = *solvent_parameters_p;
238 /* Mark the cg first as non optimized */
241 /* Check if this cg has no exclusions with atoms in other charge groups
242 * and all atoms inside the charge group excluded.
243 * We only have 3 or 4 atom solvent loops.
245 if (GET_CGINFO_EXCL_INTER(cginfo) ||
246 !GET_CGINFO_EXCL_INTRA(cginfo))
251 /* Get the indices of the first atom in this charge group */
252 j0 = molt->cgs.index[cg0];
253 j1 = molt->cgs.index[cg0+1];
255 /* Number of atoms in our molecule */
261 "Moltype '%s': there are %d atoms in this charge group\n",
265 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
268 if (nj < 3 || nj > 4)
273 /* Check if we are doing QM on this group */
275 if (qm_grpnr != nullptr)
277 for (j = j0; j < j1 && !qm; j++)
279 qm = (qm_grpnr[j] < qm_grps->nr - 1);
282 /* Cannot use solvent optimization with QM */
288 atom = molt->atoms.atom;
290 /* Still looks like a solvent, time to check parameters */
292 /* If it is perturbed (free energy) we can't use the solvent loops,
293 * so then we just skip to the next molecule.
297 for (j = j0; j < j1 && !perturbed; j++)
299 perturbed = PERTURBED(atom[j]);
307 /* Now it's only a question if the VdW and charge parameters
308 * are OK. Before doing the check we compare and see if they are
309 * identical to a possible previous solvent type.
310 * First we assign the current types and charges.
312 for (j = 0; j < nj; j++)
314 tmp_vdwtype[j] = atom[j0+j].type;
315 tmp_charge[j] = atom[j0+j].q;
318 /* Does it match any previous solvent type? */
319 for (k = 0; k < *n_solvent_parameters; k++)
324 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
325 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
326 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
331 /* Check that types & charges match for all atoms in molecule */
332 for (j = 0; j < nj && match; j++)
334 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
338 if (tmp_charge[j] != solvent_parameters[k].charge[j])
345 /* Congratulations! We have a matched solvent.
346 * Flag it with this type for later processing.
349 solvent_parameters[k].count += nmol;
351 /* We are done with this charge group */
356 /* If we get here, we have a tentative new solvent type.
357 * Before we add it we must check that it fulfills the requirements
358 * of the solvent optimized loops. First determine which atoms have
361 for (j = 0; j < nj; j++)
364 tjA = tmp_vdwtype[j];
366 /* Go through all other tpes and see if any have non-zero
367 * VdW parameters when combined with this one.
369 for (k = 0; k < fr->ntype && (!has_vdw[j]); k++)
371 /* We already checked that the atoms weren't perturbed,
372 * so we only need to check state A now.
376 has_vdw[j] = (has_vdw[j] ||
377 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
378 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
379 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
384 has_vdw[j] = (has_vdw[j] ||
385 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
386 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
391 /* Now we know all we need to make the final check and assignment. */
395 * For this we require thatn all atoms have charge,
396 * the charges on atom 2 & 3 should be the same, and only
397 * atom 1 might have VdW.
401 tmp_charge[0] != 0 &&
402 tmp_charge[1] != 0 &&
403 tmp_charge[2] == tmp_charge[1])
405 srenew(solvent_parameters, *n_solvent_parameters+1);
406 solvent_parameters[*n_solvent_parameters].model = esolSPC;
407 solvent_parameters[*n_solvent_parameters].count = nmol;
408 for (k = 0; k < 3; k++)
410 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
411 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
414 *cg_sp = *n_solvent_parameters;
415 (*n_solvent_parameters)++;
420 /* Or could it be a TIP4P?
421 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
422 * Only atom 1 mght have VdW.
427 tmp_charge[0] == 0 &&
428 tmp_charge[1] != 0 &&
429 tmp_charge[2] == tmp_charge[1] &&
432 srenew(solvent_parameters, *n_solvent_parameters+1);
433 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
434 solvent_parameters[*n_solvent_parameters].count = nmol;
435 for (k = 0; k < 4; k++)
437 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
438 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
441 *cg_sp = *n_solvent_parameters;
442 (*n_solvent_parameters)++;
446 *solvent_parameters_p = solvent_parameters;
450 check_solvent(FILE * fp,
451 const gmx_mtop_t * mtop,
453 cginfo_mb_t *cginfo_mb)
456 const gmx_moltype_t *molt;
457 int mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
458 int n_solvent_parameters;
459 solvent_parameters_t *solvent_parameters;
465 fprintf(debug, "Going to determine what solvent types we have.\n");
468 n_solvent_parameters = 0;
469 solvent_parameters = nullptr;
470 /* Allocate temporary array for solvent type */
471 snew(cg_sp, mtop->molblock.size());
474 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
476 molt = &mtop->moltype[mtop->molblock[mb].type];
478 /* Here we have to loop over all individual molecules
479 * because we need to check for QMMM particles.
481 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
482 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
483 nmol = mtop->molblock[mb].nmol/nmol_ch;
484 for (mol = 0; mol < nmol_ch; mol++)
487 am = mol*cgs->index[cgs->nr];
488 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
490 check_solvent_cg(molt, cg_mol, nmol,
491 mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty() ?
492 nullptr : mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].data()+at_offset+am,
493 &mtop->groups.groups[SimulationAtomGroupType::QuantumMechanics],
495 &n_solvent_parameters, &solvent_parameters,
496 cginfo_mb[mb].cginfo[cgm+cg_mol],
497 &cg_sp[mb][cgm+cg_mol]);
500 at_offset += cgs->index[cgs->nr];
503 /* Puh! We finished going through all charge groups.
504 * Now find the most common solvent model.
507 /* Most common solvent this far */
509 for (i = 0; i < n_solvent_parameters; i++)
512 solvent_parameters[i].count > solvent_parameters[bestsp].count)
520 bestsol = solvent_parameters[bestsp].model;
528 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
530 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
531 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
532 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
534 if (cg_sp[mb][i] == bestsp)
536 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
541 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
548 if (bestsol != esolNO && fp != nullptr)
550 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
552 solvent_parameters[bestsp].count);
555 sfree(solvent_parameters);
556 fr->solvent_opt = bestsol;
560 acNONE = 0, acCONSTRAINT, acSETTLE
563 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
564 t_forcerec *fr, gmx_bool bNoSolvOpt,
565 gmx_bool *bFEP_NonBonded,
566 gmx_bool *bExcl_IntraCGAll_InterCGNone)
569 const t_blocka *excl;
570 const gmx_moltype_t *molt;
571 const gmx_molblock_t *molb;
572 cginfo_mb_t *cginfo_mb;
575 int cg_offset, a_offset;
576 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
580 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
582 snew(cginfo_mb, mtop->molblock.size());
584 snew(type_VDW, fr->ntype);
585 for (ai = 0; ai < fr->ntype; ai++)
587 type_VDW[ai] = FALSE;
588 for (j = 0; j < fr->ntype; j++)
590 type_VDW[ai] = type_VDW[ai] ||
592 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
593 C12(fr->nbfp, fr->ntype, ai, j) != 0;
597 *bFEP_NonBonded = FALSE;
598 *bExcl_IntraCGAll_InterCGNone = TRUE;
601 snew(bExcl, excl_nalloc);
604 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
606 molb = &mtop->molblock[mb];
607 molt = &mtop->moltype[molb->type];
611 /* Check if the cginfo is identical for all molecules in this block.
612 * If so, we only need an array of the size of one molecule.
613 * Otherwise we make an array of #mol times #cgs per molecule.
616 for (m = 0; m < molb->nmol; m++)
618 int am = m*cgs->index[cgs->nr];
619 for (cg = 0; cg < cgs->nr; cg++)
622 a1 = cgs->index[cg+1];
623 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset+am+a0) !=
624 getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset +a0))
628 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
630 for (ai = a0; ai < a1; ai++)
632 if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset+am+ai] !=
633 mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset +ai])
642 cginfo_mb[mb].cg_start = cg_offset;
643 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
644 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
645 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
646 cginfo = cginfo_mb[mb].cginfo;
648 /* Set constraints flags for constrained atoms */
649 snew(a_con, molt->atoms.nr);
650 for (ftype = 0; ftype < F_NRE; ftype++)
652 if (interaction_function[ftype].flags & IF_CONSTRAINT)
657 for (ia = 0; ia < molt->ilist[ftype].size(); ia += 1+nral)
661 for (a = 0; a < nral; a++)
663 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
664 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
670 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
673 int am = m*cgs->index[cgs->nr];
674 for (cg = 0; cg < cgs->nr; cg++)
677 a1 = cgs->index[cg+1];
679 /* Store the energy group in cginfo */
680 gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput, a_offset+am+a0);
681 SET_CGINFO_GID(cginfo[cgm+cg], gid);
683 /* Check the intra/inter charge group exclusions */
684 if (a1-a0 > excl_nalloc)
686 excl_nalloc = a1 - a0;
687 srenew(bExcl, excl_nalloc);
689 /* bExclIntraAll: all intra cg interactions excluded
690 * bExclInter: any inter cg interactions excluded
692 bExclIntraAll = TRUE;
696 bHavePerturbedAtoms = FALSE;
697 for (ai = a0; ai < a1; ai++)
699 /* Check VDW and electrostatic interactions */
700 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
701 type_VDW[molt->atoms.atom[ai].typeB]);
702 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
703 molt->atoms.atom[ai].qB != 0);
705 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
707 /* Clear the exclusion list for atom ai */
708 for (aj = a0; aj < a1; aj++)
710 bExcl[aj-a0] = FALSE;
712 /* Loop over all the exclusions of atom ai */
713 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
716 if (aj < a0 || aj >= a1)
725 /* Check if ai excludes a0 to a1 */
726 for (aj = a0; aj < a1; aj++)
730 bExclIntraAll = FALSE;
737 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
740 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
748 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
752 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
754 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
756 /* The size in cginfo is currently only read with DD */
757 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
761 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
765 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
767 if (bHavePerturbedAtoms && fr->efep != efepNO)
769 SET_CGINFO_FEP(cginfo[cgm+cg]);
770 *bFEP_NonBonded = TRUE;
772 /* Store the charge group size */
773 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
775 if (!bExclIntraAll || bExclInter)
777 *bExcl_IntraCGAll_InterCGNone = FALSE;
784 cg_offset += molb->nmol*cgs->nr;
785 a_offset += molb->nmol*cgs->index[cgs->nr];
790 /* the solvent optimizer is called after the QM is initialized,
791 * because we don't want to have the QM subsystemto become an
795 check_solvent(fplog, mtop, fr, cginfo_mb);
797 if (getenv("GMX_NO_SOLV_OPT"))
801 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
802 "Disabling all solvent optimization\n");
804 fr->solvent_opt = esolNO;
808 fr->solvent_opt = esolNO;
810 if (!fr->solvent_opt)
812 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
814 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
816 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
824 static std::vector<int> cginfo_expand(const int nmb,
825 const cginfo_mb_t *cgi_mb)
827 const int ncg = cgi_mb[nmb - 1].cg_end;
829 std::vector<int> cginfo(ncg);
832 for (int cg = 0; cg < ncg; cg++)
834 while (cg >= cgi_mb[mb].cg_end)
839 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
845 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
847 if (cginfo_mb == nullptr)
851 for (int mb = 0; mb < numMolBlocks; ++mb)
853 sfree(cginfo_mb[mb].cginfo);
858 /* Sets the sum of charges (squared) and C6 in the system in fr.
859 * Returns whether the system has a net charge.
861 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
863 /*This now calculates sum for q and c6*/
864 double qsum, q2sum, q, c6sum, c6;
869 for (const gmx_molblock_t &molb : mtop->molblock)
871 int nmol = molb.nmol;
872 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
873 for (int i = 0; i < atoms->nr; i++)
875 q = atoms->atom[i].q;
878 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
883 fr->q2sum[0] = q2sum;
884 fr->c6sum[0] = c6sum;
886 if (fr->efep != efepNO)
891 for (const gmx_molblock_t &molb : mtop->molblock)
893 int nmol = molb.nmol;
894 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
895 for (int i = 0; i < atoms->nr; i++)
897 q = atoms->atom[i].qB;
900 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
904 fr->q2sum[1] = q2sum;
905 fr->c6sum[1] = c6sum;
910 fr->qsum[1] = fr->qsum[0];
911 fr->q2sum[1] = fr->q2sum[0];
912 fr->c6sum[1] = fr->c6sum[0];
916 if (fr->efep == efepNO)
918 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
922 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
923 fr->qsum[0], fr->qsum[1]);
927 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
928 return (std::abs(fr->qsum[0]) > 1e-4 ||
929 std::abs(fr->qsum[1]) > 1e-4);
932 void update_forcerec(t_forcerec *fr, matrix box)
934 if (fr->ic->eeltype == eelGRF)
936 calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
937 fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
938 &fr->ic->k_rf, &fr->ic->c_rf);
942 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
944 const t_atoms *at1, *at2;
945 int i, j, tpi, tpj, ntypes;
950 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
952 ntypes = mtop->ffparams.atnr;
956 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
958 at1 = &mtop->moltype[mt1].atoms;
959 for (i = 0; (i < at1->nr); i++)
961 tpi = at1->atom[i].type;
964 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
967 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
969 at2 = &mtop->moltype[mt2].atoms;
970 for (j = 0; (j < at2->nr); j++)
972 tpj = at2->atom[j].type;
975 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
977 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
982 if ((b < bmin) || (bmin == -1))
992 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
999 static void make_nbf_tables(FILE *fp,
1000 const interaction_const_t *ic, real rtab,
1001 const char *tabfn, char *eg1, char *eg2,
1007 if (tabfn == nullptr)
1011 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1016 sprintf(buf, "%s", tabfn);
1019 /* Append the two energy group names */
1020 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1021 eg1, eg2, ftp2ext(efXVG));
1023 nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1024 /* Copy the contents of the table to separate coulomb and LJ tables too,
1025 * to improve cache performance.
1027 /* For performance reasons we want
1028 * the table data to be aligned to 16-byte.
1031 new t_forcetable(GMX_TABLE_INTERACTION_ELEC,
1032 nbl->table_elec_vdw->format);
1033 nbl->table_elec->r = nbl->table_elec_vdw->r;
1034 nbl->table_elec->n = nbl->table_elec_vdw->n;
1035 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1036 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1037 nbl->table_elec->ninteractions = 1;
1038 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1039 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1042 new t_forcetable(GMX_TABLE_INTERACTION_VDWREP_VDWDISP,
1043 nbl->table_elec_vdw->format);
1044 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1045 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1046 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1047 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1048 nbl->table_vdw->ninteractions = 2;
1049 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1050 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1052 /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1053 * with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1055 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1057 for (j = 0; j < 4; j++)
1059 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1062 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1064 for (j = 0; j < 8; j++)
1066 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1071 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1072 * ftype2 present in the topology, build an array of the number of
1073 * interactions present for each bonded interaction index found in the
1076 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1077 * valid type with that parameter.
1079 * \c count will be reallocated as necessary to fit the largest bonded
1080 * interaction index found, and its current size will be returned in
1081 * \c ncount. It will contain zero for every bonded interaction index
1082 * for which no interactions are present in the topology.
1084 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1085 int *ncount, int **count)
1087 int ftype, i, j, tabnr;
1089 // Loop over all moleculetypes
1090 for (const gmx_moltype_t &molt : mtop->moltype)
1092 // Loop over all interaction types
1093 for (ftype = 0; ftype < F_NRE; ftype++)
1095 // If the current interaction type is one of the types whose tables we're trying to count...
1096 if (ftype == ftype1 || ftype == ftype2)
1098 const InteractionList &il = molt.ilist[ftype];
1099 const int stride = 1 + NRAL(ftype);
1100 // ... and there are actually some interactions for this type
1101 for (i = 0; i < il.size(); i += stride)
1103 // Find out which table index the user wanted
1104 tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
1107 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1109 // Make room for this index in the data structure
1110 if (tabnr >= *ncount)
1112 srenew(*count, tabnr+1);
1113 for (j = *ncount; j < tabnr+1; j++)
1119 // Record that this table index is used and must have a valid file
1127 /*!\brief If there's bonded interactions of flavour \c tabext and type
1128 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1129 * list of filenames passed to mdrun, and make bonded tables from
1132 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1133 * valid type with that parameter.
1135 * A fatal error occurs if no matching filename is found.
1137 static bondedtable_t *make_bonded_tables(FILE *fplog,
1138 int ftype1, int ftype2,
1139 const gmx_mtop_t *mtop,
1140 gmx::ArrayRef<const std::string> tabbfnm,
1150 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1152 // Are there any relevant tabulated bond interactions?
1156 for (int i = 0; i < ncount; i++)
1158 // Do any interactions exist that requires this table?
1161 // This pattern enforces the current requirement that
1162 // table filenames end in a characteristic sequence
1163 // before the file type extension, and avoids table 13
1164 // being recognized and used for table 1.
1165 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1166 bool madeTable = false;
1167 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
1169 if (gmx::endsWith(tabbfnm[j], patternToFind))
1171 // Finally read the table from the file found
1172 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1178 bool isPlural = (ftype2 != -1);
1179 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.",
1180 interaction_function[ftype1].longname,
1181 isPlural ? "' or '" : "",
1182 isPlural ? interaction_function[ftype2].longname : "",
1184 patternToFind.c_str());
1194 void forcerec_set_ranges(t_forcerec *fr,
1195 int ncg_home, int ncg_force,
1197 int natoms_force_constr, int natoms_f_novirsum)
1202 /* fr->ncg_force is unused in the standard code,
1203 * but it can be useful for modified code dealing with charge groups.
1205 fr->ncg_force = ncg_force;
1206 fr->natoms_force = natoms_force;
1207 fr->natoms_force_constr = natoms_force_constr;
1209 if (fr->natoms_force_constr > fr->nalloc_force)
1211 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1214 if (fr->haveDirectVirialContributions)
1216 fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1220 static real cutoff_inf(real cutoff)
1224 cutoff = GMX_CUTOFF_INF;
1230 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1231 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1232 bool systemHasNetCharge,
1233 interaction_const_t *ic)
1235 if (!EEL_PME_EWALD(ir->coulombtype))
1242 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1244 if (ir->coulombtype == eelP3M_AD)
1246 please_cite(fp, "Hockney1988");
1247 please_cite(fp, "Ballenegger2012");
1251 please_cite(fp, "Essmann95a");
1254 if (ir->ewald_geometry == eewg3DC)
1258 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1259 systemHasNetCharge ? " and net charge" : "");
1261 please_cite(fp, "In-Chul99a");
1262 if (systemHasNetCharge)
1264 please_cite(fp, "Ballenegger2009");
1269 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1272 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1273 1/ic->ewaldcoeff_q);
1276 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1278 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1279 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1287 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1288 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1289 interaction_const_t *ic)
1291 if (!EVDW_PME(ir->vdwtype))
1298 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1299 please_cite(fp, "Essmann95a");
1301 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1304 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1305 1/ic->ewaldcoeff_lj);
1308 if (ic->vdw_modifier == eintmodPOTSHIFT)
1310 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1311 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1315 ic->sh_lj_ewald = 0;
1319 gmx_bool uses_simple_tables(int cutoff_scheme,
1320 const nonbonded_verlet_t *nbv)
1322 gmx_bool bUsesSimpleTables = TRUE;
1324 switch (cutoff_scheme)
1327 bUsesSimpleTables = TRUE;
1330 GMX_RELEASE_ASSERT(nullptr != nbv, "A non-bonded verlet object is required with the Verlet cutoff-scheme");
1331 bUsesSimpleTables = nbv->pairlistIsSimple();
1334 gmx_incons("unimplemented");
1336 return bUsesSimpleTables;
1339 static void init_ewald_f_table(interaction_const_t *ic,
1344 /* Get the Ewald table spacing based on Coulomb and/or LJ
1345 * Ewald coefficients and rtol.
1347 ic->tabq_scale = ewald_spline3_table_scale(ic);
1349 if (ic->cutoff_scheme == ecutsVERLET)
1351 maxr = ic->rcoulomb;
1355 maxr = std::max(ic->rcoulomb, rtab);
1357 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1359 sfree_aligned(ic->tabq_coul_FDV0);
1360 sfree_aligned(ic->tabq_coul_F);
1361 sfree_aligned(ic->tabq_coul_V);
1363 sfree_aligned(ic->tabq_vdw_FDV0);
1364 sfree_aligned(ic->tabq_vdw_F);
1365 sfree_aligned(ic->tabq_vdw_V);
1367 if (EEL_PME_EWALD(ic->eeltype))
1369 /* Create the original table data in FDV0 */
1370 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1371 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1372 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1373 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1374 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1377 if (EVDW_PME(ic->vdwtype))
1379 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1380 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1381 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1382 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1383 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1387 void init_interaction_const_tables(FILE *fp,
1388 interaction_const_t *ic,
1391 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1393 init_ewald_f_table(ic, rtab);
1397 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1398 1/ic->tabq_scale, ic->tabq_size);
1403 static void clear_force_switch_constants(shift_consts_t *sc)
1410 static void force_switch_constants(real p,
1414 /* Here we determine the coefficient for shifting the force to zero
1415 * between distance rsw and the cut-off rc.
1416 * For a potential of r^-p, we have force p*r^-(p+1).
1417 * But to save flops we absorb p in the coefficient.
1419 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1420 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1422 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1423 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1424 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1427 static void potential_switch_constants(real rsw, real rc,
1428 switch_consts_t *sc)
1430 /* The switch function is 1 at rsw and 0 at rc.
1431 * The derivative and second derivate are zero at both ends.
1432 * rsw = max(r - r_switch, 0)
1433 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1434 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1435 * force = force*dsw - potential*sw
1438 sc->c3 = -10/gmx::power3(rc - rsw);
1439 sc->c4 = 15/gmx::power4(rc - rsw);
1440 sc->c5 = -6/gmx::power5(rc - rsw);
1443 /*! \brief Construct interaction constants
1445 * This data is used (particularly) by search and force code for
1446 * short-range interactions. Many of these are constant for the whole
1447 * simulation; some are constant only after PME tuning completes.
1450 init_interaction_const(FILE *fp,
1451 interaction_const_t **interaction_const,
1452 const t_inputrec *ir,
1453 const gmx_mtop_t *mtop,
1454 bool systemHasNetCharge)
1456 interaction_const_t *ic;
1460 ic->cutoff_scheme = ir->cutoff_scheme;
1462 /* Just allocate something so we can free it */
1463 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1464 snew_aligned(ic->tabq_coul_F, 16, 32);
1465 snew_aligned(ic->tabq_coul_V, 16, 32);
1468 ic->vdwtype = ir->vdwtype;
1469 ic->vdw_modifier = ir->vdw_modifier;
1470 ic->reppow = mtop->ffparams.reppow;
1471 ic->rvdw = cutoff_inf(ir->rvdw);
1472 ic->rvdw_switch = ir->rvdw_switch;
1473 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
1474 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
1475 if (ic->useBuckingham)
1477 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
1480 initVdwEwaldParameters(fp, ir, ic);
1482 clear_force_switch_constants(&ic->dispersion_shift);
1483 clear_force_switch_constants(&ic->repulsion_shift);
1485 switch (ic->vdw_modifier)
1487 case eintmodPOTSHIFT:
1488 /* Only shift the potential, don't touch the force */
1489 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
1490 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
1492 case eintmodFORCESWITCH:
1493 /* Switch the force, switch and shift the potential */
1494 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
1495 &ic->dispersion_shift);
1496 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
1497 &ic->repulsion_shift);
1499 case eintmodPOTSWITCH:
1500 /* Switch the potential and force */
1501 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
1505 case eintmodEXACTCUTOFF:
1506 /* Nothing to do here */
1509 gmx_incons("unimplemented potential modifier");
1512 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
1514 /* Electrostatics */
1515 ic->eeltype = ir->coulombtype;
1516 ic->coulomb_modifier = ir->coulomb_modifier;
1517 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
1518 ic->rcoulomb_switch = ir->rcoulomb_switch;
1519 ic->epsilon_r = ir->epsilon_r;
1521 /* Set the Coulomb energy conversion factor */
1522 if (ic->epsilon_r != 0)
1524 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
1528 /* eps = 0 is infinite dieletric: no Coulomb interactions */
1532 /* Reaction-field */
1533 if (EEL_RF(ic->eeltype))
1535 ic->epsilon_rf = ir->epsilon_rf;
1536 /* Generalized reaction field parameters are updated every step */
1537 if (ic->eeltype != eelGRF)
1539 calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
1540 ic->rcoulomb, 0, 0, nullptr,
1541 &ic->k_rf, &ic->c_rf);
1544 if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
1546 /* grompp should have done this, but this scheme is obsolete */
1547 ic->coulomb_modifier = eintmodEXACTCUTOFF;
1552 /* For plain cut-off we might use the reaction-field kernels */
1553 ic->epsilon_rf = ic->epsilon_r;
1555 if (ir->coulomb_modifier == eintmodPOTSHIFT)
1557 ic->c_rf = 1/ic->rcoulomb;
1565 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
1569 real dispersion_shift;
1571 dispersion_shift = ic->dispersion_shift.cpot;
1572 if (EVDW_PME(ic->vdwtype))
1574 dispersion_shift -= ic->sh_lj_ewald;
1576 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
1577 ic->repulsion_shift.cpot, dispersion_shift);
1579 if (ic->eeltype == eelCUT)
1581 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
1583 else if (EEL_PME(ic->eeltype))
1585 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
1590 *interaction_const = ic;
1594 done_interaction_const(interaction_const_t *interaction_const)
1596 sfree_aligned(interaction_const->tabq_coul_FDV0);
1597 sfree_aligned(interaction_const->tabq_coul_F);
1598 sfree_aligned(interaction_const->tabq_coul_V);
1599 sfree(interaction_const);
1602 void init_forcerec(FILE *fp,
1603 const gmx::MDLogger &mdlog,
1606 const t_inputrec *ir,
1607 const gmx_mtop_t *mtop,
1608 const t_commrec *cr,
1612 gmx::ArrayRef<const std::string> tabbfnm,
1613 const gmx_hw_info_t &hardwareInfo,
1614 const gmx_device_info_t *deviceInfo,
1615 const bool useGpuForBonded,
1616 gmx_bool bNoSolvOpt,
1619 int m, negp_pp, negptable, egi, egj;
1624 gmx_bool bGenericKernelOnly;
1625 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
1626 gmx_bool bFEP_NonBonded;
1627 int *nm_ind, egp_flags;
1630 /* By default we turn SIMD kernels on, but it might be turned off further down... */
1631 fr->use_simd_kernels = TRUE;
1633 if (check_box(ir->ePBC, box))
1635 gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
1638 /* Test particle insertion ? */
1641 /* Set to the size of the molecule to be inserted (the last one) */
1642 /* Because of old style topologies, we have to use the last cg
1643 * instead of the last molecule type.
1645 cgs = &mtop->moltype[mtop->molblock.back().type].cgs;
1646 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
1647 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
1648 if (fr->n_tpi != molecules.block(molecules.numBlocks() - 1).size())
1650 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
1658 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
1660 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
1661 eel_names[ir->coulombtype]);
1666 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1668 if (ir->useTwinRange)
1670 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1672 /* Copy the user determined parameters */
1673 fr->userint1 = ir->userint1;
1674 fr->userint2 = ir->userint2;
1675 fr->userint3 = ir->userint3;
1676 fr->userint4 = ir->userint4;
1677 fr->userreal1 = ir->userreal1;
1678 fr->userreal2 = ir->userreal2;
1679 fr->userreal3 = ir->userreal3;
1680 fr->userreal4 = ir->userreal4;
1683 fr->fc_stepsize = ir->fc_stepsize;
1686 fr->efep = ir->efep;
1687 fr->sc_alphavdw = ir->fepvals->sc_alpha;
1688 if (ir->fepvals->bScCoul)
1690 fr->sc_alphacoul = ir->fepvals->sc_alpha;
1691 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1695 fr->sc_alphacoul = 0;
1696 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1698 fr->sc_power = ir->fepvals->sc_power;
1699 fr->sc_r_power = ir->fepvals->sc_r_power;
1700 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1702 env = getenv("GMX_SCSIGMA_MIN");
1706 sscanf(env, "%20lf", &dbl);
1707 fr->sc_sigma6_min = gmx::power6(dbl);
1710 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1714 fr->bNonbonded = TRUE;
1715 if (getenv("GMX_NO_NONBONDED") != nullptr)
1717 /* turn off non-bonded calculations */
1718 fr->bNonbonded = FALSE;
1719 GMX_LOG(mdlog.warning).asParagraph().appendText(
1720 "Found environment variable GMX_NO_NONBONDED.\n"
1721 "Disabling nonbonded calculations.");
1724 bGenericKernelOnly = FALSE;
1726 /* We now check in the NS code whether a particular combination of interactions
1727 * can be used with water optimization, and disable it if that is not the case.
1730 if (getenv("GMX_NB_GENERIC") != nullptr)
1735 "Found environment variable GMX_NB_GENERIC.\n"
1736 "Disabling all interaction-specific nonbonded kernels, will only\n"
1737 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
1739 bGenericKernelOnly = TRUE;
1742 if (bGenericKernelOnly)
1747 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
1749 fr->use_simd_kernels = FALSE;
1753 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1754 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1755 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1759 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1761 /* Neighbour searching stuff */
1762 fr->cutoff_scheme = ir->cutoff_scheme;
1763 fr->bGrid = (ir->ns_type == ensGRID);
1764 fr->ePBC = ir->ePBC;
1766 if (fr->cutoff_scheme == ecutsGROUP)
1768 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
1769 "removed in a future release when 'verlet' supports all interaction forms.\n";
1773 fprintf(stderr, "\n%s\n", note);
1777 fprintf(fp, "\n%s\n", note);
1781 /* Determine if we will do PBC for distances in bonded interactions */
1782 if (fr->ePBC == epbcNONE)
1784 fr->bMolPBC = FALSE;
1788 if (!DOMAINDECOMP(cr))
1792 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
1793 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
1794 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
1796 /* The group cut-off scheme and SHAKE assume charge groups
1797 * are whole, but not using molpbc is faster in most cases.
1798 * With intermolecular interactions we need PBC for calculating
1799 * distances between atoms in different molecules.
1801 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
1802 !mtop->bIntermolecularInteractions)
1804 fr->bMolPBC = ir->bPeriodicMols;
1806 if (bSHAKE && fr->bMolPBC)
1808 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
1813 /* Not making molecules whole is faster in most cases,
1814 * but With orientation restraints we need whole molecules.
1816 fr->bMolPBC = (fcd->orires.nr == 0);
1818 if (getenv("GMX_USE_GRAPH") != nullptr)
1820 fr->bMolPBC = FALSE;
1823 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
1826 if (mtop->bIntermolecularInteractions)
1828 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
1832 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
1834 if (bSHAKE && fr->bMolPBC)
1836 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");
1842 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
1846 fr->rc_scaling = ir->refcoord_scaling;
1847 copy_rvec(ir->posres_com, fr->posres_com);
1848 copy_rvec(ir->posres_comB, fr->posres_comB);
1849 fr->rlist = cutoff_inf(ir->rlist);
1850 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1852 /* This now calculates sum for q and c6*/
1853 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1855 /* fr->ic is used both by verlet and group kernels (to some extent) now */
1856 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1857 init_interaction_const_tables(fp, fr->ic, ir->rlist + ir->tabext);
1859 const interaction_const_t *ic = fr->ic;
1861 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1862 if (ir->coulombtype == eelEWALD)
1864 init_ewald_tab(&(fr->ewald_table), ir, fp);
1867 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1868 switch (ic->eeltype)
1871 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
1876 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
1880 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
1881 GMX_RELEASE_ASSERT(ic->coulomb_modifier == eintmodEXACTCUTOFF, "With the group scheme RF-zero needs the exact cut-off modifier");
1890 case eelPMEUSERSWITCH:
1891 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1897 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
1901 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1903 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1905 /* Vdw: Translate from mdp settings to kernel format */
1906 switch (ic->vdwtype)
1911 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1915 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1919 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
1925 case evdwENCADSHIFT:
1926 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1930 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1932 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1934 if (ir->cutoff_scheme == ecutsGROUP)
1936 fr->bvdwtab = ((ic->vdwtype != evdwCUT || !gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
1937 && !EVDW_PME(ic->vdwtype));
1938 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
1939 fr->bcoultab = !(ic->eeltype == eelCUT ||
1940 ic->eeltype == eelEWALD ||
1941 ic->eeltype == eelPME ||
1942 ic->eeltype == eelP3M_AD ||
1943 ic->eeltype == eelRF ||
1944 ic->eeltype == eelRF_ZERO);
1946 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
1947 * going to be faster to tabulate the interaction than calling the generic kernel.
1948 * However, if generic kernels have been requested we keep things analytically.
1950 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
1951 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
1952 !bGenericKernelOnly)
1954 if ((ic->rcoulomb_switch != ic->rvdw_switch) || (ic->rcoulomb != ic->rvdw))
1956 fr->bcoultab = TRUE;
1957 /* Once we tabulate electrostatics, we can use the switch function for LJ,
1958 * which would otherwise need two tables.
1962 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
1963 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
1964 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
1965 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
1967 if ((ic->rcoulomb != ic->rvdw) && (!bGenericKernelOnly))
1969 fr->bcoultab = TRUE;
1973 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
1975 fr->bcoultab = TRUE;
1977 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
1982 if (getenv("GMX_REQUIRE_TABLES"))
1985 fr->bcoultab = TRUE;
1990 fprintf(fp, "Table routines are used for coulomb: %s\n",
1991 gmx::boolToString(fr->bcoultab));
1992 fprintf(fp, "Table routines are used for vdw: %s\n",
1993 gmx::boolToString(fr->bvdwtab));
1998 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1999 fr->nbkernel_vdw_modifier = eintmodNONE;
2003 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2004 fr->nbkernel_elec_modifier = eintmodNONE;
2008 if (ir->cutoff_scheme == ecutsVERLET)
2010 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2012 gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2014 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2015 * while mdrun does not (and never did) support this.
2017 if (EEL_USER(fr->ic->eeltype))
2019 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2020 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2023 fr->bvdwtab = FALSE;
2024 fr->bcoultab = FALSE;
2027 /* 1-4 interaction electrostatics */
2028 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2030 /* Parameters for generalized RF */
2034 if (ic->eeltype == eelGRF)
2036 init_generalized_rf(fp, mtop, ir, fr);
2039 fr->haveDirectVirialContributions =
2040 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
2041 fr->forceProviders->hasForceProvider() ||
2042 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2043 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2049 if (fr->haveDirectVirialContributions)
2051 fr->forceBufferForDirectVirialContributions = new std::vector<gmx::RVec>;
2054 if (fr->shift_vec == nullptr)
2056 snew(fr->shift_vec, SHIFTS);
2059 if (fr->fshift == nullptr)
2061 snew(fr->fshift, SHIFTS);
2064 if (fr->nbfp == nullptr)
2066 fr->ntype = mtop->ffparams.atnr;
2067 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2068 if (EVDW_PME(ic->vdwtype))
2070 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2074 /* Copy the energy group exclusions */
2075 fr->egp_flags = ir->opts.egp_flags;
2077 /* Van der Waals stuff */
2078 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
2080 if (ic->rvdw_switch >= ic->rvdw)
2082 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2083 ic->rvdw_switch, ic->rvdw);
2087 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2088 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
2089 ic->rvdw_switch, ic->rvdw);
2093 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
2095 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2098 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
2100 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2103 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2105 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2108 if (fp && fr->cutoff_scheme == ecutsGROUP)
2110 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2111 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
2114 if (ir->implicit_solvent)
2116 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
2119 /* Construct tables for the group scheme. A little unnecessary to
2120 * make both vdw and coul tables sometimes, but what the
2121 * heck. Note that both cutoff schemes construct Ewald tables in
2122 * init_interaction_const_tables. */
2123 needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
2124 (fr->bcoultab || fr->bvdwtab));
2126 negp_pp = ir->opts.ngener - ir->nwall;
2128 if (!needGroupSchemeTables)
2130 bSomeNormalNbListsAreInUse = TRUE;
2134 bSomeNormalNbListsAreInUse = FALSE;
2135 for (egi = 0; egi < negp_pp; egi++)
2137 for (egj = egi; egj < negp_pp; egj++)
2139 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2140 if (!(egp_flags & EGP_EXCL))
2142 if (egp_flags & EGP_TABLE)
2148 bSomeNormalNbListsAreInUse = TRUE;
2155 /* This code automatically gives table length tabext without cut-off's,
2156 * in that case grompp should already have checked that we do not need
2157 * normal tables and we only generate tables for 1-4 interactions.
2159 rtab = ir->rlist + ir->tabext;
2161 if (needGroupSchemeTables)
2163 /* make tables for ordinary interactions */
2164 if (bSomeNormalNbListsAreInUse)
2166 make_nbf_tables(fp, ic, rtab, tabfn, nullptr, nullptr, nullptr);
2175 /* Read the special tables for certain energy group pairs */
2176 nm_ind = mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].nm_ind;
2177 for (egi = 0; egi < negp_pp; egi++)
2179 for (egj = egi; egj < negp_pp; egj++)
2181 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2182 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2184 /* Read the table file with the two energy groups names appended */
2185 make_nbf_tables(fp, ic, rtab, tabfn,
2186 *mtop->groups.groupNames[nm_ind[egi]],
2187 *mtop->groups.groupNames[nm_ind[egj]],
2196 /* We want to use unmodified tables for 1-4 coulombic
2197 * interactions, so we must in general have an extra set of
2199 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2200 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2201 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
2203 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
2204 GMX_MAKETABLES_14ONLY);
2208 fr->nwall = ir->nwall;
2209 if (ir->nwall && ir->wall_type == ewtTABLE)
2211 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
2214 if (fcd && !tabbfnm.empty())
2216 // Need to catch std::bad_alloc
2217 // TODO Don't need to catch this here, when merging with master branch
2220 fcd->bondtab = make_bonded_tables(fp,
2221 F_TABBONDS, F_TABBONDSNC,
2222 mtop, tabbfnm, "b");
2223 fcd->angletab = make_bonded_tables(fp,
2225 mtop, tabbfnm, "a");
2226 fcd->dihtab = make_bonded_tables(fp,
2228 mtop, tabbfnm, "d");
2230 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2236 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
2240 // QM/MM initialization if requested
2241 fr->bQMMM = ir->bQMMM;
2244 // Initialize QM/MM if supported
2247 GMX_LOG(mdlog.info).asParagraph().
2248 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
2249 "version. Please get in touch with the developers if you find the support useful, "
2250 "as help is needed if the functionality is to continue to be available.");
2251 fr->qr = mk_QMMMrec();
2252 init_QMMMrec(cr, mtop, ir, fr);
2256 gmx_incons("QM/MM was requested, but is only available when GROMACS "
2257 "is configured with QM/MM support");
2261 /* Set all the static charge group info */
2262 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
2264 &fr->bExcl_IntraCGAll_InterCGNone);
2265 if (!DOMAINDECOMP(cr))
2267 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
2270 if (!DOMAINDECOMP(cr))
2272 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
2273 mtop->natoms, mtop->natoms, mtop->natoms);
2276 fr->print_force = print_force;
2278 /* Initialize the thread working data for bonded interactions */
2279 fr->bondedThreading =
2280 init_bonded_threading(fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].nr);
2282 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
2283 snew(fr->ewc_t, fr->nthread_ewc);
2285 if (fr->cutoff_scheme == ecutsVERLET)
2287 // We checked the cut-offs in grompp, but double-check here.
2288 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
2289 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
2291 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
2295 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
2298 fr->nbv = Nbnxm::init_nb_verlet(mdlog, bFEP_NonBonded, ir, fr,
2299 cr, hardwareInfo, deviceInfo,
2302 if (useGpuForBonded)
2304 auto stream = DOMAINDECOMP(cr) ?
2305 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::NonLocal) :
2306 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::Local);
2307 // TODO the heap allocation is only needed while
2308 // t_forcerec lacks a constructor.
2309 fr->gpuBonded = new gmx::GpuBonded(mtop->ffparams,
2314 if (ir->eDispCorr != edispcNO)
2316 fr->dispersionCorrection =
2317 std::make_unique<DispersionCorrection>(*mtop, *ir, fr->bBHAM,
2319 gmx::arrayRefFromArray(fr->nbfp, fr->ntype*fr->ntype*2),
2321 fr->dispersionCorrection->print(mdlog);
2326 /* Here we switch from using mdlog, which prints the newline before
2327 * the paragraph, to our old fprintf logging, which prints the newline
2328 * after the paragraph, so we should add a newline here.
2334 /* Frees GPU memory and sets a tMPI node barrier.
2336 * Note that this function needs to be called even if GPUs are not used
2337 * in this run because the PME ranks have no knowledge of whether GPUs
2338 * are used or not, but all ranks need to enter the barrier below.
2339 * \todo Remove physical node barrier from this function after making sure
2340 * that it's not needed anymore (with a shared GPU run).
2342 void free_gpu_resources(t_forcerec *fr,
2343 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator)
2345 bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->useGpu();
2347 /* stop the GPU profiler (only CUDA) */
2350 if (isPPrankUsingGPU)
2352 /* Free data in GPU memory and pinned memory before destroying the GPU context */
2355 delete fr->gpuBonded;
2356 fr->gpuBonded = nullptr;
2359 /* With tMPI we need to wait for all ranks to finish deallocation before
2360 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
2363 * This is not a concern in OpenCL where we use one context per rank which
2364 * is freed in nbnxn_gpu_free().
2366 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
2367 * but it is easier and more futureproof to call it on the whole node.
2371 physicalNodeCommunicator.barrier();
2375 void done_forcerec(t_forcerec *fr, int numMolBlocks)
2379 // PME-only ranks don't have a forcerec
2382 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
2384 done_interaction_const(fr->ic);
2385 sfree(fr->shift_vec);
2388 tear_down_bonded_threading(fr->bondedThreading);
2389 GMX_RELEASE_ASSERT(fr->gpuBonded == nullptr, "Should have been deleted earlier, when used");
2390 fr->bondedThreading = nullptr;