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-2020, 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/ewald/pme_pp_comm_gpu.h"
57 #include "gromacs/fileio/filetypes.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/hardware/hw_info.h"
62 #include "gromacs/listed_forces/gpubonded.h"
63 #include "gromacs/listed_forces/manage_threading.h"
64 #include "gromacs/listed_forces/pairs.h"
65 #include "gromacs/math/functions.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/utilities.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/mdlib/dispersioncorrection.h"
70 #include "gromacs/mdlib/force.h"
71 #include "gromacs/mdlib/forcerec_threading.h"
72 #include "gromacs/mdlib/gmx_omp_nthreads.h"
73 #include "gromacs/mdlib/md_support.h"
74 #include "gromacs/mdlib/qmmm.h"
75 #include "gromacs/mdlib/rf_util.h"
76 #include "gromacs/mdlib/wall.h"
77 #include "gromacs/mdtypes/commrec.h"
78 #include "gromacs/mdtypes/fcdata.h"
79 #include "gromacs/mdtypes/group.h"
80 #include "gromacs/mdtypes/iforceprovider.h"
81 #include "gromacs/mdtypes/inputrec.h"
82 #include "gromacs/mdtypes/md_enums.h"
83 #include "gromacs/nbnxm/gpu_data_mgmt.h"
84 #include "gromacs/nbnxm/nbnxm.h"
85 #include "gromacs/nbnxm/nbnxm_geometry.h"
86 #include "gromacs/pbcutil/ishift.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/tables/forcetable.h"
89 #include "gromacs/topology/mtop_util.h"
90 #include "gromacs/trajectory/trajectoryframe.h"
91 #include "gromacs/utility/cstringutil.h"
92 #include "gromacs/utility/exceptions.h"
93 #include "gromacs/utility/fatalerror.h"
94 #include "gromacs/utility/gmxassert.h"
95 #include "gromacs/utility/logger.h"
96 #include "gromacs/utility/physicalnodecommunicator.h"
97 #include "gromacs/utility/pleasecite.h"
98 #include "gromacs/utility/smalloc.h"
99 #include "gromacs/utility/strconvert.h"
101 /*! \brief environment variable to enable GPU P2P communication */
102 static const bool c_enableGpuPmePpComms =
103 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr) && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA);
105 static real* mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
113 snew(nbfp, 3 * atnr * atnr);
114 for (i = k = 0; (i < atnr); i++)
116 for (j = 0; (j < atnr); j++, k++)
118 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
119 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
120 /* nbfp now includes the 6.0 derivative prefactor */
121 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c * 6.0;
127 snew(nbfp, 2 * atnr * atnr);
128 for (i = k = 0; (i < atnr); i++)
130 for (j = 0; (j < atnr); j++, k++)
132 /* nbfp now includes the 6.0/12.0 derivative prefactors */
133 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6 * 6.0;
134 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12 * 12.0;
142 static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
145 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
148 /* For LJ-PME simulations, we correct the energies with the reciprocal space
149 * inside of the cut-off. To do this the non-bonded kernels needs to have
150 * access to the C6-values used on the reciprocal grid in pme.c
154 snew(grid, 2 * atnr * atnr);
155 for (i = k = 0; (i < atnr); i++)
157 for (j = 0; (j < atnr); j++, k++)
159 c6i = idef->iparams[i * (atnr + 1)].lj.c6;
160 c12i = idef->iparams[i * (atnr + 1)].lj.c12;
161 c6j = idef->iparams[j * (atnr + 1)].lj.c6;
162 c12j = idef->iparams[j * (atnr + 1)].lj.c12;
163 c6 = std::sqrt(c6i * c6j);
164 if (fr->ljpme_combination_rule == eljpmeLB && !gmx_numzero(c6) && !gmx_numzero(c12i)
165 && !gmx_numzero(c12j))
167 sigmai = gmx::sixthroot(c12i / c6i);
168 sigmaj = gmx::sixthroot(c12j / c6j);
169 epsi = c6i * c6i / c12i;
170 epsj = c6j * c6j / c12j;
171 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
173 /* Store the elements at the same relative positions as C6 in nbfp in order
174 * to simplify access in the kernels
176 grid[2 * (atnr * i + j)] = c6 * 6.0;
189 static cginfo_mb_t* init_cginfo_mb(const gmx_mtop_t* mtop, const t_forcerec* fr, gmx_bool* bFEP_NonBonded)
191 cginfo_mb_t* cginfo_mb;
196 snew(cginfo_mb, mtop->molblock.size());
198 snew(type_VDW, fr->ntype);
199 for (int ai = 0; ai < fr->ntype; ai++)
201 type_VDW[ai] = FALSE;
202 for (int j = 0; j < fr->ntype; j++)
204 type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
205 || C12(fr->nbfp, fr->ntype, ai, j) != 0;
209 *bFEP_NonBonded = FALSE;
212 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
214 const gmx_molblock_t& molb = mtop->molblock[mb];
215 const gmx_moltype_t& molt = mtop->moltype[molb.type];
216 const t_blocka& excl = molt.excls;
218 /* Check if the cginfo is identical for all molecules in this block.
219 * If so, we only need an array of the size of one molecule.
220 * Otherwise we make an array of #mol times #cgs per molecule.
223 for (int m = 0; m < molb.nmol; m++)
225 const int am = m * molt.atoms.nr;
226 for (int a = 0; a < molt.atoms.nr; a++)
228 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
229 != getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
233 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
235 if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
236 != mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
244 cginfo_mb[mb].cg_start = a_offset;
245 cginfo_mb[mb].cg_end = a_offset + molb.nmol * molt.atoms.nr;
246 cginfo_mb[mb].cg_mod = (bId ? 1 : molb.nmol) * molt.atoms.nr;
247 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
248 cginfo = cginfo_mb[mb].cginfo;
250 /* Set constraints flags for constrained atoms */
251 snew(a_con, molt.atoms.nr);
252 for (int ftype = 0; ftype < F_NRE; ftype++)
254 if (interaction_function[ftype].flags & IF_CONSTRAINT)
256 const int nral = NRAL(ftype);
257 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
261 for (a = 0; a < nral; a++)
263 a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
264 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
270 for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
272 const int molculeOffsetInBlock = m * molt.atoms.nr;
273 for (int a = 0; a < molt.atoms.nr; a++)
275 const t_atom& atom = molt.atoms.atom[a];
276 int& atomInfo = cginfo[molculeOffsetInBlock + a];
278 /* Store the energy group in cginfo */
279 int gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput,
280 a_offset + molculeOffsetInBlock + a);
281 SET_CGINFO_GID(atomInfo, gid);
283 bool bHaveVDW = (type_VDW[atom.type] || type_VDW[atom.typeB]);
284 bool bHaveQ = (atom.q != 0 || atom.qB != 0);
286 bool haveExclusions = false;
287 /* Loop over all the exclusions of atom ai */
288 for (int j = excl.index[a]; j < excl.index[a + 1]; j++)
292 haveExclusions = true;
299 case acCONSTRAINT: SET_CGINFO_CONSTR(atomInfo); break;
300 case acSETTLE: SET_CGINFO_SETTLE(atomInfo); break;
305 SET_CGINFO_EXCL_INTER(atomInfo);
309 SET_CGINFO_HAS_VDW(atomInfo);
313 SET_CGINFO_HAS_Q(atomInfo);
315 if (fr->efep != efepNO && PERTURBED(atom))
317 SET_CGINFO_FEP(atomInfo);
318 *bFEP_NonBonded = TRUE;
325 a_offset += molb.nmol * molt.atoms.nr;
332 static std::vector<int> cginfo_expand(const int nmb, const cginfo_mb_t* cgi_mb)
334 const int ncg = cgi_mb[nmb - 1].cg_end;
336 std::vector<int> cginfo(ncg);
339 for (int cg = 0; cg < ncg; cg++)
341 while (cg >= cgi_mb[mb].cg_end)
345 cginfo[cg] = cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
351 static void done_cginfo_mb(cginfo_mb_t* cginfo_mb, int numMolBlocks)
353 if (cginfo_mb == nullptr)
357 for (int mb = 0; mb < numMolBlocks; ++mb)
359 sfree(cginfo_mb[mb].cginfo);
364 /* Sets the sum of charges (squared) and C6 in the system in fr.
365 * Returns whether the system has a net charge.
367 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t* mtop)
369 /*This now calculates sum for q and c6*/
370 double qsum, q2sum, q, c6sum, c6;
375 for (const gmx_molblock_t& molb : mtop->molblock)
377 int nmol = molb.nmol;
378 const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
379 for (int i = 0; i < atoms->nr; i++)
381 q = atoms->atom[i].q;
383 q2sum += nmol * q * q;
384 c6 = mtop->ffparams.iparams[atoms->atom[i].type * (mtop->ffparams.atnr + 1)].lj.c6;
389 fr->q2sum[0] = q2sum;
390 fr->c6sum[0] = c6sum;
392 if (fr->efep != efepNO)
397 for (const gmx_molblock_t& molb : mtop->molblock)
399 int nmol = molb.nmol;
400 const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
401 for (int i = 0; i < atoms->nr; i++)
403 q = atoms->atom[i].qB;
405 q2sum += nmol * q * q;
406 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB * (mtop->ffparams.atnr + 1)].lj.c6;
410 fr->q2sum[1] = q2sum;
411 fr->c6sum[1] = c6sum;
416 fr->qsum[1] = fr->qsum[0];
417 fr->q2sum[1] = fr->q2sum[0];
418 fr->c6sum[1] = fr->c6sum[0];
422 if (fr->efep == efepNO)
424 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
428 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
432 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
433 return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
436 static real calcBuckinghamBMax(FILE* fplog, const gmx_mtop_t* mtop)
438 const t_atoms *at1, *at2;
439 int i, j, tpi, tpj, ntypes;
444 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
446 ntypes = mtop->ffparams.atnr;
450 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
452 at1 = &mtop->moltype[mt1].atoms;
453 for (i = 0; (i < at1->nr); i++)
455 tpi = at1->atom[i].type;
458 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
461 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
463 at2 = &mtop->moltype[mt2].atoms;
464 for (j = 0; (j < at2->nr); j++)
466 tpj = at2->atom[j].type;
469 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
471 b = mtop->ffparams.iparams[tpi * ntypes + tpj].bham.b;
476 if ((b < bmin) || (bmin == -1))
486 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n", bmin, bham_b_max);
492 /*!\brief If there's bonded interactions of type \c ftype1 or \c
493 * ftype2 present in the topology, build an array of the number of
494 * interactions present for each bonded interaction index found in the
497 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
498 * valid type with that parameter.
500 * \c count will be reallocated as necessary to fit the largest bonded
501 * interaction index found, and its current size will be returned in
502 * \c ncount. It will contain zero for every bonded interaction index
503 * for which no interactions are present in the topology.
505 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t* mtop, int* ncount, int** count)
507 int ftype, i, j, tabnr;
509 // Loop over all moleculetypes
510 for (const gmx_moltype_t& molt : mtop->moltype)
512 // Loop over all interaction types
513 for (ftype = 0; ftype < F_NRE; ftype++)
515 // If the current interaction type is one of the types whose tables we're trying to count...
516 if (ftype == ftype1 || ftype == ftype2)
518 const InteractionList& il = molt.ilist[ftype];
519 const int stride = 1 + NRAL(ftype);
520 // ... and there are actually some interactions for this type
521 for (i = 0; i < il.size(); i += stride)
523 // Find out which table index the user wanted
524 tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
527 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
529 // Make room for this index in the data structure
530 if (tabnr >= *ncount)
532 srenew(*count, tabnr + 1);
533 for (j = *ncount; j < tabnr + 1; j++)
539 // Record that this table index is used and must have a valid file
547 /*!\brief If there's bonded interactions of flavour \c tabext and type
548 * \c ftype1 or \c ftype2 present in the topology, seek them in the
549 * list of filenames passed to mdrun, and make bonded tables from
552 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
553 * valid type with that parameter.
555 * A fatal error occurs if no matching filename is found.
557 static bondedtable_t* make_bonded_tables(FILE* fplog,
560 const gmx_mtop_t* mtop,
561 gmx::ArrayRef<const std::string> tabbfnm,
571 count_tables(ftype1, ftype2, mtop, &ncount, &count);
573 // Are there any relevant tabulated bond interactions?
577 for (int i = 0; i < ncount; i++)
579 // Do any interactions exist that requires this table?
582 // This pattern enforces the current requirement that
583 // table filenames end in a characteristic sequence
584 // before the file type extension, and avoids table 13
585 // being recognized and used for table 1.
586 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
587 bool madeTable = false;
588 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
590 if (gmx::endsWith(tabbfnm[j], patternToFind))
592 // Finally read the table from the file found
593 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
599 bool isPlural = (ftype2 != -1);
601 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
602 "because no table file whose name matched '%s' was passed via the "
603 "gmx mdrun -tableb command-line option.",
604 interaction_function[ftype1].longname, isPlural ? "' or '" : "",
605 isPlural ? interaction_function[ftype2].longname : "", i,
606 patternToFind.c_str());
616 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
618 fr->natoms_force = natoms_force;
619 fr->natoms_force_constr = natoms_force_constr;
621 if (fr->natoms_force_constr > fr->nalloc_force)
623 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
626 if (fr->haveDirectVirialContributions)
628 fr->forceBufferForDirectVirialContributions.resize(natoms_f_novirsum);
632 static real cutoff_inf(real cutoff)
636 cutoff = GMX_CUTOFF_INF;
642 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
643 static void initCoulombEwaldParameters(FILE* fp,
644 const t_inputrec* ir,
645 bool systemHasNetCharge,
646 interaction_const_t* ic)
648 if (!EEL_PME_EWALD(ir->coulombtype))
655 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
657 if (ir->coulombtype == eelP3M_AD)
659 please_cite(fp, "Hockney1988");
660 please_cite(fp, "Ballenegger2012");
664 please_cite(fp, "Essmann95a");
667 if (ir->ewald_geometry == eewg3DC)
671 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
672 systemHasNetCharge ? " and net charge" : "");
674 please_cite(fp, "In-Chul99a");
675 if (systemHasNetCharge)
677 please_cite(fp, "Ballenegger2009");
682 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
685 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
688 if (ic->coulomb_modifier == eintmodPOTSHIFT)
690 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
691 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
699 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
700 static void initVdwEwaldParameters(FILE* fp, const t_inputrec* ir, interaction_const_t* ic)
702 if (!EVDW_PME(ir->vdwtype))
709 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
710 please_cite(fp, "Essmann95a");
712 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
715 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
718 if (ic->vdw_modifier == eintmodPOTSHIFT)
720 real crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
721 ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
729 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
731 * Tables are generated for one or both, depending on if the pointers
732 * are non-null. The spacing for both table sets is the same and obeys
733 * both accuracy requirements, when relevant.
735 static void init_ewald_f_table(const interaction_const_t& ic,
736 EwaldCorrectionTables* coulombTables,
737 EwaldCorrectionTables* vdwTables)
739 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
740 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
742 /* Get the Ewald table spacing based on Coulomb and/or LJ
743 * Ewald coefficients and rtol.
745 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
747 const int tableSize = static_cast<int>(ic.rcoulomb * tableScale) + 2;
752 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
757 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
761 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic)
763 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
765 init_ewald_f_table(*ic, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
768 fprintf(fp, "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
769 1 / ic->coulombEwaldTables->scale, ic->coulombEwaldTables->tableF.size());
774 static void clear_force_switch_constants(shift_consts_t* sc)
781 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
783 /* Here we determine the coefficient for shifting the force to zero
784 * between distance rsw and the cut-off rc.
785 * For a potential of r^-p, we have force p*r^-(p+1).
786 * But to save flops we absorb p in the coefficient.
788 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
789 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
791 sc->c2 = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
792 sc->c3 = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
793 sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
794 + p * sc->c3 / 4 * gmx::power4(rc - rsw);
797 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
799 /* The switch function is 1 at rsw and 0 at rc.
800 * The derivative and second derivate are zero at both ends.
801 * rsw = max(r - r_switch, 0)
802 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
803 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
804 * force = force*dsw - potential*sw
807 sc->c3 = -10 / gmx::power3(rc - rsw);
808 sc->c4 = 15 / gmx::power4(rc - rsw);
809 sc->c5 = -6 / gmx::power5(rc - rsw);
812 /*! \brief Construct interaction constants
814 * This data is used (particularly) by search and force code for
815 * short-range interactions. Many of these are constant for the whole
816 * simulation; some are constant only after PME tuning completes.
818 static void init_interaction_const(FILE* fp,
819 interaction_const_t** interaction_const,
820 const t_inputrec* ir,
821 const gmx_mtop_t* mtop,
822 bool systemHasNetCharge)
824 interaction_const_t* ic = new interaction_const_t;
826 ic->cutoff_scheme = ir->cutoff_scheme;
828 ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
829 ic->vdwEwaldTables = std::make_unique<EwaldCorrectionTables>();
832 ic->vdwtype = ir->vdwtype;
833 ic->vdw_modifier = ir->vdw_modifier;
834 ic->reppow = mtop->ffparams.reppow;
835 ic->rvdw = cutoff_inf(ir->rvdw);
836 ic->rvdw_switch = ir->rvdw_switch;
837 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
838 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
839 if (ic->useBuckingham)
841 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
844 initVdwEwaldParameters(fp, ir, ic);
846 clear_force_switch_constants(&ic->dispersion_shift);
847 clear_force_switch_constants(&ic->repulsion_shift);
849 switch (ic->vdw_modifier)
851 case eintmodPOTSHIFT:
852 /* Only shift the potential, don't touch the force */
853 ic->dispersion_shift.cpot = -1.0 / gmx::power6(ic->rvdw);
854 ic->repulsion_shift.cpot = -1.0 / gmx::power12(ic->rvdw);
856 case eintmodFORCESWITCH:
857 /* Switch the force, switch and shift the potential */
858 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw, &ic->dispersion_shift);
859 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw, &ic->repulsion_shift);
861 case eintmodPOTSWITCH:
862 /* Switch the potential and force */
863 potential_switch_constants(ic->rvdw_switch, ic->rvdw, &ic->vdw_switch);
866 case eintmodEXACTCUTOFF:
867 /* Nothing to do here */
869 default: gmx_incons("unimplemented potential modifier");
873 ic->eeltype = ir->coulombtype;
874 ic->coulomb_modifier = ir->coulomb_modifier;
875 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
876 ic->rcoulomb_switch = ir->rcoulomb_switch;
877 ic->epsilon_r = ir->epsilon_r;
879 /* Set the Coulomb energy conversion factor */
880 if (ic->epsilon_r != 0)
882 ic->epsfac = ONE_4PI_EPS0 / ic->epsilon_r;
886 /* eps = 0 is infinite dieletric: no Coulomb interactions */
891 if (EEL_RF(ic->eeltype))
893 GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
894 ic->epsilon_rf = ir->epsilon_rf;
896 calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf, ic->rcoulomb, &ic->k_rf, &ic->c_rf);
900 /* For plain cut-off we might use the reaction-field kernels */
901 ic->epsilon_rf = ic->epsilon_r;
903 if (ir->coulomb_modifier == eintmodPOTSHIFT)
905 ic->c_rf = 1 / ic->rcoulomb;
913 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
917 real dispersion_shift;
919 dispersion_shift = ic->dispersion_shift.cpot;
920 if (EVDW_PME(ic->vdwtype))
922 dispersion_shift -= ic->sh_lj_ewald;
924 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e", ic->repulsion_shift.cpot, dispersion_shift);
926 if (ic->eeltype == eelCUT)
928 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
930 else if (EEL_PME(ic->eeltype))
932 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
937 *interaction_const = ic;
940 bool areMoleculesDistributedOverPbc(const t_inputrec& ir, const gmx_mtop_t& mtop, const gmx::MDLogger& mdlog)
942 bool areMoleculesDistributedOverPbc = false;
943 const bool useEwaldSurfaceCorrection = (EEL_PME_EWALD(ir.coulombtype) && ir.epsilon_surface != 0);
946 (ir.eConstrAlg == econtSHAKE
947 && (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 || gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
949 /* The group cut-off scheme and SHAKE assume charge groups
950 * are whole, but not using molpbc is faster in most cases.
951 * With intermolecular interactions we need PBC for calculating
952 * distances between atoms in different molecules.
954 if (bSHAKE && !mtop.bIntermolecularInteractions)
956 areMoleculesDistributedOverPbc = ir.bPeriodicMols;
958 if (areMoleculesDistributedOverPbc)
960 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
965 /* Not making molecules whole is faster in most cases,
966 * but with orientation restraints or non-tinfoil boundary
967 * conditions we need whole molecules.
969 areMoleculesDistributedOverPbc =
970 (gmx_mtop_ftype_count(mtop, F_ORIRES) == 0 && !useEwaldSurfaceCorrection);
972 if (getenv("GMX_USE_GRAPH") != nullptr)
974 areMoleculesDistributedOverPbc = false;
976 GMX_LOG(mdlog.warning)
979 "GMX_USE_GRAPH is set, using the graph for bonded "
982 if (mtop.bIntermolecularInteractions)
984 GMX_LOG(mdlog.warning)
987 "WARNING: Molecules linked by intermolecular interactions "
988 "have to reside in the same periodic image, otherwise "
989 "artifacts will occur!");
993 GMX_RELEASE_ASSERT(areMoleculesDistributedOverPbc || !mtop.bIntermolecularInteractions,
994 "We need to use PBC within molecules with inter-molecular interactions");
996 if (bSHAKE && areMoleculesDistributedOverPbc)
999 "SHAKE is not properly supported with intermolecular interactions. "
1000 "For short simulations where linked molecules remain in the same "
1001 "periodic image, the environment variable GMX_USE_GRAPH can be used "
1002 "to override this check.\n");
1006 return areMoleculesDistributedOverPbc;
1009 void init_forcerec(FILE* fp,
1010 const gmx::MDLogger& mdlog,
1013 const t_inputrec* ir,
1014 const gmx_mtop_t* mtop,
1015 const t_commrec* cr,
1019 gmx::ArrayRef<const std::string> tabbfnm,
1020 const gmx_hw_info_t& hardwareInfo,
1021 const gmx_device_info_t* deviceInfo,
1022 const bool useGpuForBonded,
1023 const bool pmeOnlyRankUsesGpu,
1025 gmx_wallcycle* wcycle)
1030 gmx_bool bFEP_NonBonded;
1032 /* By default we turn SIMD kernels on, but it might be turned off further down... */
1033 fr->use_simd_kernels = TRUE;
1035 if (check_box(ir->ePBC, box))
1037 gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
1040 /* Test particle insertion ? */
1043 /* Set to the size of the molecule to be inserted (the last one) */
1044 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
1045 fr->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
1052 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED || ir->coulombtype == eelGRF_NOTUSED)
1054 gmx_fatal(FARGS, "%s electrostatics is no longer supported", eel_names[ir->coulombtype]);
1059 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1061 if (ir->useTwinRange)
1063 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1065 /* Copy the user determined parameters */
1066 fr->userint1 = ir->userint1;
1067 fr->userint2 = ir->userint2;
1068 fr->userint3 = ir->userint3;
1069 fr->userint4 = ir->userint4;
1070 fr->userreal1 = ir->userreal1;
1071 fr->userreal2 = ir->userreal2;
1072 fr->userreal3 = ir->userreal3;
1073 fr->userreal4 = ir->userreal4;
1076 fr->fc_stepsize = ir->fc_stepsize;
1079 fr->efep = ir->efep;
1080 fr->sc_alphavdw = ir->fepvals->sc_alpha;
1081 if (ir->fepvals->bScCoul)
1083 fr->sc_alphacoul = ir->fepvals->sc_alpha;
1084 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1088 fr->sc_alphacoul = 0;
1089 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1091 fr->sc_power = ir->fepvals->sc_power;
1092 fr->sc_r_power = ir->fepvals->sc_r_power;
1093 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1095 env = getenv("GMX_SCSIGMA_MIN");
1099 sscanf(env, "%20lf", &dbl);
1100 fr->sc_sigma6_min = gmx::power6(dbl);
1103 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1107 fr->bNonbonded = TRUE;
1108 if (getenv("GMX_NO_NONBONDED") != nullptr)
1110 /* turn off non-bonded calculations */
1111 fr->bNonbonded = FALSE;
1112 GMX_LOG(mdlog.warning)
1115 "Found environment variable GMX_NO_NONBONDED.\n"
1116 "Disabling nonbonded calculations.");
1119 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1121 fr->use_simd_kernels = FALSE;
1125 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1126 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1127 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1131 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1133 /* Neighbour searching stuff */
1134 fr->cutoff_scheme = ir->cutoff_scheme;
1135 fr->ePBC = ir->ePBC;
1137 /* Determine if we will do PBC for distances in bonded interactions */
1138 if (fr->ePBC == epbcNONE)
1140 fr->bMolPBC = FALSE;
1144 const bool useEwaldSurfaceCorrection =
1145 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1146 if (!DOMAINDECOMP(cr))
1148 fr->bMolPBC = areMoleculesDistributedOverPbc(*ir, *mtop, mdlog);
1149 // The assert below is equivalent to fcd->orires.nr != gmx_mtop_ftype_count(mtop, F_ORIRES)
1150 GMX_RELEASE_ASSERT(!fr->bMolPBC || fcd->orires.nr == 0,
1151 "Molecules broken over PBC exist in a simulation including "
1152 "orientation restraints. "
1153 "This likely means that the global topology and the force constant "
1154 "data have gotten out of sync.");
1155 if (useEwaldSurfaceCorrection)
1158 "In GROMACS 2020, Ewald dipole correction is disabled when not "
1159 "using domain decomposition. With domain decomposition, it only works "
1160 "when each molecule consists of a single update group (e.g. water). "
1161 "This will be fixed in GROMACS 2021.");
1166 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
1168 if (useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd))
1171 "You requested dipole correction (epsilon_surface > 0), but molecules "
1173 "over periodic boundary conditions by the domain decomposition. "
1174 "Run without domain decomposition instead.");
1178 if (useEwaldSurfaceCorrection)
1180 GMX_RELEASE_ASSERT((!DOMAINDECOMP(cr) && !fr->bMolPBC)
1181 || (DOMAINDECOMP(cr) && dd_moleculesAreAlwaysWhole(*cr->dd)),
1182 "Molecules can not be broken by PBC with epsilon_surface > 0");
1186 fr->rc_scaling = ir->refcoord_scaling;
1187 copy_rvec(ir->posres_com, fr->posres_com);
1188 copy_rvec(ir->posres_comB, fr->posres_comB);
1189 fr->rlist = cutoff_inf(ir->rlist);
1190 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1192 /* This now calculates sum for q and c6*/
1193 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1195 /* fr->ic is used both by verlet and group kernels (to some extent) now */
1196 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1197 init_interaction_const_tables(fp, fr->ic);
1199 const interaction_const_t* ic = fr->ic;
1201 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1202 if (ir->coulombtype == eelEWALD)
1204 init_ewald_tab(&(fr->ewald_table), ir, fp);
1207 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1208 switch (ic->eeltype)
1210 case eelCUT: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB; break;
1213 case eelRF_ZERO: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD; break;
1221 case eelPMEUSERSWITCH:
1222 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1227 case eelEWALD: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD; break;
1230 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1232 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1234 /* Vdw: Translate from mdp settings to kernel format */
1235 switch (ic->vdwtype)
1240 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1244 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1247 case evdwPME: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD; break;
1252 case evdwENCADSHIFT:
1253 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1256 default: gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1258 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1260 if (ir->cutoff_scheme == ecutsVERLET)
1262 if (!gmx_within_tol(ic->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1264 gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12",
1265 ecutscheme_names[ir->cutoff_scheme]);
1267 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1268 * while mdrun does not (and never did) support this.
1270 if (EEL_USER(fr->ic->eeltype))
1272 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
1273 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
1276 fr->bvdwtab = FALSE;
1277 fr->bcoultab = FALSE;
1280 /* 1-4 interaction electrostatics */
1281 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1283 fr->haveDirectVirialContributions =
1284 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) || fr->forceProviders->hasForceProvider()
1285 || gmx_mtop_ftype_count(mtop, F_POSRES) > 0 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0
1286 || ir->nwall > 0 || ir->bPull || ir->bRot || ir->bIMD);
1288 if (fr->shift_vec == nullptr)
1290 snew(fr->shift_vec, SHIFTS);
1293 fr->shiftForces.resize(SHIFTS);
1295 if (fr->nbfp == nullptr)
1297 fr->ntype = mtop->ffparams.atnr;
1298 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1299 if (EVDW_PME(ic->vdwtype))
1301 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1305 /* Copy the energy group exclusions */
1306 fr->egp_flags = ir->opts.egp_flags;
1308 /* Van der Waals stuff */
1309 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1311 if (ic->rvdw_switch >= ic->rvdw)
1313 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)", ic->rvdw_switch, ic->rvdw);
1317 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
1318 (ic->eeltype == eelSWITCH) ? "switched" : "shifted", ic->rvdw_switch, ic->rvdw);
1322 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1324 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1327 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1329 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1332 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
1334 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
1337 if (fp && fr->cutoff_scheme == ecutsGROUP)
1339 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n", fr->rlist, ic->rcoulomb,
1340 fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
1343 if (ir->implicit_solvent)
1345 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1349 /* This code automatically gives table length tabext without cut-off's,
1350 * in that case grompp should already have checked that we do not need
1351 * normal tables and we only generate tables for 1-4 interactions.
1353 rtab = ir->rlist + ir->tabext;
1355 /* We want to use unmodified tables for 1-4 coulombic
1356 * interactions, so we must in general have an extra set of
1358 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1359 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1361 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1365 fr->nwall = ir->nwall;
1366 if (ir->nwall && ir->wall_type == ewtTABLE)
1368 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1371 if (fcd && !tabbfnm.empty())
1373 // Need to catch std::bad_alloc
1374 // TODO Don't need to catch this here, when merging with master branch
1377 fcd->bondtab = make_bonded_tables(fp, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1378 fcd->angletab = make_bonded_tables(fp, F_TABANGLES, -1, mtop, tabbfnm, "a");
1379 fcd->dihtab = make_bonded_tables(fp, F_TABDIHS, -1, mtop, tabbfnm, "d");
1381 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1388 "No fcdata or table file name passed, can not read table, can not do bonded "
1393 // QM/MM initialization if requested
1394 fr->bQMMM = ir->bQMMM;
1397 // Initialize QM/MM if supported
1403 "Large parts of the QM/MM support is deprecated, and may be removed in "
1405 "version. Please get in touch with the developers if you find the "
1407 "as help is needed if the functionality is to continue to be "
1409 fr->qr = mk_QMMMrec();
1410 init_QMMMrec(cr, mtop, ir, fr);
1415 "QM/MM was requested, but is only available when GROMACS "
1416 "is configured with QM/MM support");
1420 /* Set all the static charge group info */
1421 fr->cginfo_mb = init_cginfo_mb(mtop, fr, &bFEP_NonBonded);
1422 if (!DOMAINDECOMP(cr))
1424 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1427 if (!DOMAINDECOMP(cr))
1429 forcerec_set_ranges(fr, mtop->natoms, mtop->natoms, mtop->natoms);
1432 fr->print_force = print_force;
1434 /* Initialize the thread working data for bonded interactions */
1435 fr->bondedThreading = init_bonded_threading(
1436 fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size());
1438 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1439 snew(fr->ewc_t, fr->nthread_ewc);
1441 if (fr->cutoff_scheme == ecutsVERLET)
1443 // We checked the cut-offs in grompp, but double-check here.
1444 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
1445 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
1447 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
1448 "With Verlet lists and PME we should have rcoulomb>=rvdw");
1453 ir->rcoulomb == ir->rvdw,
1454 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
1457 fr->nbv = Nbnxm::init_nb_verlet(mdlog, bFEP_NonBonded, ir, fr, cr, hardwareInfo, deviceInfo,
1460 if (useGpuForBonded)
1462 auto stream = havePPDomainDecomposition(cr)
1463 ? Nbnxm::gpu_get_command_stream(
1464 fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1465 : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1466 gmx::InteractionLocality::Local);
1467 // TODO the heap allocation is only needed while
1468 // t_forcerec lacks a constructor.
1469 fr->gpuBonded = new gmx::GpuBonded(mtop->ffparams, stream, wcycle);
1473 if (ir->eDispCorr != edispcNO)
1475 fr->dispersionCorrection = std::make_unique<DispersionCorrection>(
1476 *mtop, *ir, fr->bBHAM, fr->ntype,
1477 gmx::arrayRefFromArray(fr->nbfp, fr->ntype * fr->ntype * 2), *fr->ic, tabfn);
1478 fr->dispersionCorrection->print(mdlog);
1483 /* Here we switch from using mdlog, which prints the newline before
1484 * the paragraph, to our old fprintf logging, which prints the newline
1485 * after the paragraph, so we should add a newline here.
1490 if (pmeOnlyRankUsesGpu && c_enableGpuPmePpComms)
1492 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(cr->mpi_comm_mysim, cr->dd->pme_nodeid);
1496 t_forcerec::t_forcerec() = default;
1498 t_forcerec::~t_forcerec() = default;
1500 /* Frees GPU memory and sets a tMPI node barrier.
1502 * Note that this function needs to be called even if GPUs are not used
1503 * in this run because the PME ranks have no knowledge of whether GPUs
1504 * are used or not, but all ranks need to enter the barrier below.
1505 * \todo Remove physical node barrier from this function after making sure
1506 * that it's not needed anymore (with a shared GPU run).
1508 void free_gpu_resources(t_forcerec* fr,
1509 const gmx::PhysicalNodeCommunicator& physicalNodeCommunicator,
1510 const gmx_gpu_info_t& gpu_info)
1512 bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->useGpu();
1514 /* stop the GPU profiler (only CUDA) */
1515 if (gpu_info.n_dev > 0)
1520 if (isPPrankUsingGPU)
1522 /* Free data in GPU memory and pinned memory before destroying the GPU context */
1525 delete fr->gpuBonded;
1526 fr->gpuBonded = nullptr;
1529 /* With tMPI we need to wait for all ranks to finish deallocation before
1530 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1533 * This is not a concern in OpenCL where we use one context per rank which
1534 * is freed in nbnxn_gpu_free().
1536 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1537 * but it is easier and more futureproof to call it on the whole node.
1541 physicalNodeCommunicator.barrier();
1545 void done_forcerec(t_forcerec* fr, int numMolBlocks)
1549 // PME-only ranks don't have a forcerec
1552 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
1555 sfree(fr->shift_vec);
1557 tear_down_bonded_threading(fr->bondedThreading);
1558 GMX_RELEASE_ASSERT(fr->gpuBonded == nullptr, "Should have been deleted earlier, when used");
1559 fr->bondedThreading = nullptr;