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-2019,2020,2021, 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/listed_forces.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/rf_util.h"
75 #include "gromacs/mdlib/wall.h"
76 #include "gromacs/mdlib/wholemoleculetransform.h"
77 #include "gromacs/mdtypes/commrec.h"
78 #include "gromacs/mdtypes/fcdata.h"
79 #include "gromacs/mdtypes/forcerec.h"
80 #include "gromacs/mdtypes/group.h"
81 #include "gromacs/mdtypes/iforceprovider.h"
82 #include "gromacs/mdtypes/inputrec.h"
83 #include "gromacs/mdtypes/interaction_const.h"
84 #include "gromacs/mdtypes/md_enums.h"
85 #include "gromacs/mdtypes/multipletimestepping.h"
86 #include "gromacs/mdtypes/nblist.h"
87 #include "gromacs/nbnxm/nbnxm.h"
88 #include "gromacs/pbcutil/ishift.h"
89 #include "gromacs/pbcutil/pbc.h"
90 #include "gromacs/tables/forcetable.h"
91 #include "gromacs/topology/mtop_util.h"
92 #include "gromacs/trajectory/trajectoryframe.h"
93 #include "gromacs/utility/cstringutil.h"
94 #include "gromacs/utility/exceptions.h"
95 #include "gromacs/utility/fatalerror.h"
96 #include "gromacs/utility/gmxassert.h"
97 #include "gromacs/utility/logger.h"
98 #include "gromacs/utility/physicalnodecommunicator.h"
99 #include "gromacs/utility/pleasecite.h"
100 #include "gromacs/utility/smalloc.h"
101 #include "gromacs/utility/strconvert.h"
103 #include "gpuforcereduction.h"
105 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
106 haveDirectVirialContributions_(haveDirectVirialContributions)
108 shiftForces_.resize(SHIFTS);
111 void ForceHelperBuffers::resize(int numAtoms)
113 if (haveDirectVirialContributions_)
115 forceBufferForDirectVirialContributions_.resize(numAtoms);
119 std::vector<real> makeNonBondedParameterLists(const gmx_ffparams_t& forceFieldParams, bool useBuckinghamPotential)
121 std::vector<real> nbfp;
124 atnr = forceFieldParams.atnr;
125 if (useBuckinghamPotential)
127 nbfp.resize(3 * atnr * atnr);
129 for (int i = 0; (i < atnr); i++)
131 for (int j = 0; (j < atnr); j++, k++)
133 BHAMA(nbfp, atnr, i, j) = forceFieldParams.iparams[k].bham.a;
134 BHAMB(nbfp, atnr, i, j) = forceFieldParams.iparams[k].bham.b;
135 /* nbfp now includes the 6.0 derivative prefactor */
136 BHAMC(nbfp, atnr, i, j) = forceFieldParams.iparams[k].bham.c * 6.0;
142 nbfp.resize(2 * atnr * atnr);
144 for (int i = 0; (i < atnr); i++)
146 for (int j = 0; (j < atnr); j++, k++)
148 /* nbfp now includes the 6.0/12.0 derivative prefactors */
149 C6(nbfp, atnr, i, j) = forceFieldParams.iparams[k].lj.c6 * 6.0;
150 C12(nbfp, atnr, i, j) = forceFieldParams.iparams[k].lj.c12 * 12.0;
158 std::vector<real> makeLJPmeC6GridCorrectionParameters(const gmx_ffparams_t& forceFieldParams,
159 const t_forcerec& forceRec)
162 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
164 /* For LJ-PME simulations, we correct the energies with the reciprocal space
165 * inside of the cut-off. To do this the non-bonded kernels needs to have
166 * access to the C6-values used on the reciprocal grid in pme.c
169 atnr = forceFieldParams.atnr;
170 std::vector<real> grid(2 * atnr * atnr, 0.0);
171 for (i = k = 0; (i < atnr); i++)
173 for (j = 0; (j < atnr); j++, k++)
175 c6i = forceFieldParams.iparams[i * (atnr + 1)].lj.c6;
176 c12i = forceFieldParams.iparams[i * (atnr + 1)].lj.c12;
177 c6j = forceFieldParams.iparams[j * (atnr + 1)].lj.c6;
178 c12j = forceFieldParams.iparams[j * (atnr + 1)].lj.c12;
179 c6 = std::sqrt(c6i * c6j);
180 if (forceRec.ljpme_combination_rule == LongRangeVdW::LB && !gmx_numzero(c6)
181 && !gmx_numzero(c12i) && !gmx_numzero(c12j))
183 sigmai = gmx::sixthroot(c12i / c6i);
184 sigmaj = gmx::sixthroot(c12j / c6j);
185 epsi = c6i * c6i / c12i;
186 epsj = c6j * c6j / c12j;
187 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
189 /* Store the elements at the same relative positions as C6 in nbfp in order
190 * to simplify access in the kernels
192 grid[2 * (atnr * i + j)] = c6 * 6.0;
205 static std::vector<cginfo_mb_t> init_cginfo_mb(const gmx_mtop_t& mtop, const t_forcerec* fr)
210 snew(type_VDW, fr->ntype);
211 for (int ai = 0; ai < fr->ntype; ai++)
213 type_VDW[ai] = FALSE;
214 for (int j = 0; j < fr->ntype; j++)
216 type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
217 || C12(fr->nbfp, fr->ntype, ai, j) != 0;
221 std::vector<cginfo_mb_t> cginfoPerMolblock;
223 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
225 const gmx_molblock_t& molb = mtop.molblock[mb];
226 const gmx_moltype_t& molt = mtop.moltype[molb.type];
227 const auto& excl = molt.excls;
229 /* Check if the cginfo is identical for all molecules in this block.
230 * If so, we only need an array of the size of one molecule.
231 * Otherwise we make an array of #mol times #cgs per molecule.
234 for (int m = 0; m < molb.nmol; m++)
236 const int am = m * molt.atoms.nr;
237 for (int a = 0; a < molt.atoms.nr; a++)
239 if (getGroupType(mtop.groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
240 != getGroupType(mtop.groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
244 if (!mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
246 if (mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
247 != mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
255 cginfo_mb_t cginfo_mb;
256 cginfo_mb.cg_start = a_offset;
257 cginfo_mb.cg_end = a_offset + molb.nmol * molt.atoms.nr;
258 cginfo_mb.cg_mod = (bId ? 1 : molb.nmol) * molt.atoms.nr;
259 cginfo_mb.cginfo.resize(cginfo_mb.cg_mod);
260 gmx::ArrayRef<int> cginfo = cginfo_mb.cginfo;
262 /* Set constraints flags for constrained atoms */
263 snew(a_con, molt.atoms.nr);
264 for (int ftype = 0; ftype < F_NRE; ftype++)
266 if (interaction_function[ftype].flags & IF_CONSTRAINT)
268 const int nral = NRAL(ftype);
269 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
273 for (a = 0; a < nral; a++)
275 a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
276 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
282 for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
284 const int molculeOffsetInBlock = m * molt.atoms.nr;
285 for (int a = 0; a < molt.atoms.nr; a++)
287 const t_atom& atom = molt.atoms.atom[a];
288 int& atomInfo = cginfo[molculeOffsetInBlock + a];
290 /* Store the energy group in cginfo */
291 int gid = getGroupType(mtop.groups,
292 SimulationAtomGroupType::EnergyOutput,
293 a_offset + molculeOffsetInBlock + a);
294 SET_CGINFO_GID(atomInfo, gid);
296 bool bHaveVDW = (type_VDW[atom.type] || type_VDW[atom.typeB]);
297 bool bHaveQ = (atom.q != 0 || atom.qB != 0);
299 bool haveExclusions = false;
300 /* Loop over all the exclusions of atom ai */
301 for (const int j : excl[a])
305 haveExclusions = true;
312 case acCONSTRAINT: SET_CGINFO_CONSTR(atomInfo); break;
313 case acSETTLE: SET_CGINFO_SETTLE(atomInfo); break;
318 SET_CGINFO_EXCL_INTER(atomInfo);
322 SET_CGINFO_HAS_VDW(atomInfo);
326 SET_CGINFO_HAS_Q(atomInfo);
328 if (fr->efep != FreeEnergyPerturbationType::No && PERTURBED(atom))
330 SET_CGINFO_FEP(atomInfo);
337 cginfoPerMolblock.push_back(cginfo_mb);
339 a_offset += molb.nmol * molt.atoms.nr;
343 return cginfoPerMolblock;
346 static std::vector<int> cginfo_expand(const int nmb, gmx::ArrayRef<const cginfo_mb_t> cgi_mb)
348 const int ncg = cgi_mb[nmb - 1].cg_end;
350 std::vector<int> cginfo(ncg);
353 for (int cg = 0; cg < ncg; cg++)
355 while (cg >= cgi_mb[mb].cg_end)
359 cginfo[cg] = cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
365 /* Sets the sum of charges (squared) and C6 in the system in fr.
366 * Returns whether the system has a net charge.
368 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t& mtop)
370 /*This now calculates sum for q and c6*/
371 double qsum, q2sum, q, c6sum, c6;
376 for (const gmx_molblock_t& molb : mtop.molblock)
378 int nmol = molb.nmol;
379 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
380 for (int i = 0; i < atoms->nr; i++)
382 q = atoms->atom[i].q;
384 q2sum += nmol * q * q;
385 c6 = mtop.ffparams.iparams[atoms->atom[i].type * (mtop.ffparams.atnr + 1)].lj.c6;
390 fr->q2sum[0] = q2sum;
391 fr->c6sum[0] = c6sum;
393 if (fr->efep != FreeEnergyPerturbationType::No)
398 for (const gmx_molblock_t& molb : mtop.molblock)
400 int nmol = molb.nmol;
401 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
402 for (int i = 0; i < atoms->nr; i++)
404 q = atoms->atom[i].qB;
406 q2sum += nmol * q * q;
407 c6 = mtop.ffparams.iparams[atoms->atom[i].typeB * (mtop.ffparams.atnr + 1)].lj.c6;
411 fr->q2sum[1] = q2sum;
412 fr->c6sum[1] = c6sum;
417 fr->qsum[1] = fr->qsum[0];
418 fr->q2sum[1] = fr->q2sum[0];
419 fr->c6sum[1] = fr->c6sum[0];
423 if (fr->efep == FreeEnergyPerturbationType::No)
425 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
429 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
433 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
434 return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
437 static real calcBuckinghamBMax(FILE* fplog, const gmx_mtop_t& mtop)
439 const t_atoms *at1, *at2;
440 int i, j, tpi, tpj, ntypes;
445 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
447 ntypes = mtop.ffparams.atnr;
451 for (size_t mt1 = 0; mt1 < mtop.moltype.size(); mt1++)
453 at1 = &mtop.moltype[mt1].atoms;
454 for (i = 0; (i < at1->nr); i++)
456 tpi = at1->atom[i].type;
459 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
462 for (size_t mt2 = mt1; mt2 < mtop.moltype.size(); mt2++)
464 at2 = &mtop.moltype[mt2].atoms;
465 for (j = 0; (j < at2->nr); j++)
467 tpj = at2->atom[j].type;
470 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
472 b = mtop.ffparams.iparams[tpi * ntypes + tpj].bham.b;
477 if ((b < bmin) || (bmin == -1))
487 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n", bmin, bham_b_max);
493 /*!\brief If there's bonded interactions of type \c ftype1 or \c
494 * ftype2 present in the topology, build an array of the number of
495 * interactions present for each bonded interaction index found in the
498 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
499 * valid type with that parameter.
501 * \c count will be reallocated as necessary to fit the largest bonded
502 * interaction index found, and its current size will be returned in
503 * \c ncount. It will contain zero for every bonded interaction index
504 * for which no interactions are present in the topology.
506 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t& mtop, int* ncount, int** count)
508 int ftype, i, j, tabnr;
510 // Loop over all moleculetypes
511 for (const gmx_moltype_t& molt : mtop.moltype)
513 // Loop over all interaction types
514 for (ftype = 0; ftype < F_NRE; ftype++)
516 // If the current interaction type is one of the types whose tables we're trying to count...
517 if (ftype == ftype1 || ftype == ftype2)
519 const InteractionList& il = molt.ilist[ftype];
520 const int stride = 1 + NRAL(ftype);
521 // ... and there are actually some interactions for this type
522 for (i = 0; i < il.size(); i += stride)
524 // Find out which table index the user wanted
525 tabnr = mtop.ffparams.iparams[il.iatoms[i]].tab.table;
528 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
530 // Make room for this index in the data structure
531 if (tabnr >= *ncount)
533 srenew(*count, tabnr + 1);
534 for (j = *ncount; j < tabnr + 1; j++)
540 // Record that this table index is used and must have a valid file
548 /*!\brief If there's bonded interactions of flavour \c tabext and type
549 * \c ftype1 or \c ftype2 present in the topology, seek them in the
550 * list of filenames passed to mdrun, and make bonded tables from
553 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
554 * valid type with that parameter.
556 * A fatal error occurs if no matching filename is found.
558 static std::vector<bondedtable_t> make_bonded_tables(FILE* fplog,
561 const gmx_mtop_t& mtop,
562 gmx::ArrayRef<const std::string> tabbfnm,
565 std::vector<bondedtable_t> tab;
568 int* count = nullptr;
569 count_tables(ftype1, ftype2, mtop, &ncount, &count);
571 // Are there any relevant tabulated bond interactions?
575 for (int i = 0; i < ncount; i++)
577 // Do any interactions exist that requires this table?
580 // This pattern enforces the current requirement that
581 // table filenames end in a characteristic sequence
582 // before the file type extension, and avoids table 13
583 // being recognized and used for table 1.
584 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
585 bool madeTable = false;
586 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
588 if (gmx::endsWith(tabbfnm[j], patternToFind))
590 // Finally read the table from the file found
591 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
597 bool isPlural = (ftype2 != -1);
599 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
600 "because no table file whose name matched '%s' was passed via the "
601 "gmx mdrun -tableb command-line option.",
602 interaction_function[ftype1].longname,
603 isPlural ? "' or '" : "",
604 isPlural ? interaction_function[ftype2].longname : "",
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 for (auto& forceHelperBuffers : fr->forceHelperBuffers)
623 forceHelperBuffers.resize(natoms_f_novirsum);
627 static real cutoff_inf(real cutoff)
631 cutoff = GMX_CUTOFF_INF;
637 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
638 static void initCoulombEwaldParameters(FILE* fp,
639 const t_inputrec& ir,
640 bool systemHasNetCharge,
641 interaction_const_t* ic)
643 if (!EEL_PME_EWALD(ir.coulombtype))
650 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
652 if (ir.coulombtype == CoulombInteractionType::P3mAD)
654 please_cite(fp, "Hockney1988");
655 please_cite(fp, "Ballenegger2012");
659 please_cite(fp, "Essmann95a");
662 if (ir.ewald_geometry == EwaldGeometry::ThreeDC)
667 "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
668 systemHasNetCharge ? " and net charge" : "");
670 please_cite(fp, "In-Chul99a");
671 if (systemHasNetCharge)
673 please_cite(fp, "Ballenegger2009");
678 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir.rcoulomb, ir.ewald_rtol);
681 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
684 if (ic->coulomb_modifier == InteractionModifiers::PotShift)
686 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
687 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
695 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
696 static void initVdwEwaldParameters(FILE* fp, const t_inputrec& ir, interaction_const_t* ic)
698 if (!EVDW_PME(ir.vdwtype))
705 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
706 please_cite(fp, "Essmann95a");
708 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir.rvdw, ir.ewald_rtol_lj);
711 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
714 if (ic->vdw_modifier == InteractionModifiers::PotShift)
716 real crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
717 ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
725 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
727 * Tables are generated for one or both, depending on if the pointers
728 * are non-null. The spacing for both table sets is the same and obeys
729 * both accuracy requirements, when relevant.
731 static void init_ewald_f_table(const interaction_const_t& ic,
734 EwaldCorrectionTables* coulombTables,
735 EwaldCorrectionTables* vdwTables)
737 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
738 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
740 /* Get the Ewald table spacing based on Coulomb and/or LJ
741 * Ewald coefficients and rtol.
743 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
745 const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
747 real tableLen = ic.rcoulomb;
748 if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
750 /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
751 * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
752 * The alternative is to look through all the exclusions and check if they come from
753 * couple-intramol == no. Meanwhile, always having larger tables should only affect
754 * memory consumption, not speed (barring cache issues).
756 tableLen = rlist + tabext;
758 const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
763 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
768 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
772 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
774 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
777 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
781 "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
782 1 / ic->coulombEwaldTables->scale,
783 ic->coulombEwaldTables->tableF.size());
788 static void clear_force_switch_constants(shift_consts_t* sc)
795 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
797 /* Here we determine the coefficient for shifting the force to zero
798 * between distance rsw and the cut-off rc.
799 * For a potential of r^-p, we have force p*r^-(p+1).
800 * But to save flops we absorb p in the coefficient.
802 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
803 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
805 sc->c2 = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
806 sc->c3 = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
807 sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
808 + p * sc->c3 / 4 * gmx::power4(rc - rsw);
811 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
813 /* The switch function is 1 at rsw and 0 at rc.
814 * The derivative and second derivate are zero at both ends.
815 * rsw = max(r - r_switch, 0)
816 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
817 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
818 * force = force*dsw - potential*sw
821 sc->c3 = -10 / gmx::power3(rc - rsw);
822 sc->c4 = 15 / gmx::power4(rc - rsw);
823 sc->c5 = -6 / gmx::power5(rc - rsw);
826 /*! \brief Construct interaction constants
828 * This data is used (particularly) by search and force code for
829 * short-range interactions. Many of these are constant for the whole
830 * simulation; some are constant only after PME tuning completes.
832 static interaction_const_t init_interaction_const(FILE* fp,
833 const t_inputrec& ir,
834 const gmx_mtop_t& mtop,
835 bool systemHasNetCharge)
837 interaction_const_t interactionConst;
839 interactionConst.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
840 interactionConst.vdwEwaldTables = std::make_unique<EwaldCorrectionTables>();
843 interactionConst.vdwtype = ir.vdwtype;
844 interactionConst.vdw_modifier = ir.vdw_modifier;
845 interactionConst.reppow = mtop.ffparams.reppow;
846 interactionConst.rvdw = cutoff_inf(ir.rvdw);
847 interactionConst.rvdw_switch = ir.rvdw_switch;
848 interactionConst.ljpme_comb_rule = ir.ljpme_combination_rule;
849 interactionConst.useBuckingham = (mtop.ffparams.functype[0] == F_BHAM);
850 if (interactionConst.useBuckingham)
852 interactionConst.buckinghamBMax = calcBuckinghamBMax(fp, mtop);
855 initVdwEwaldParameters(fp, ir, &interactionConst);
857 clear_force_switch_constants(&interactionConst.dispersion_shift);
858 clear_force_switch_constants(&interactionConst.repulsion_shift);
860 switch (interactionConst.vdw_modifier)
862 case InteractionModifiers::PotShift:
863 /* Only shift the potential, don't touch the force */
864 interactionConst.dispersion_shift.cpot = -1.0 / gmx::power6(interactionConst.rvdw);
865 interactionConst.repulsion_shift.cpot = -1.0 / gmx::power12(interactionConst.rvdw);
867 case InteractionModifiers::ForceSwitch:
868 /* Switch the force, switch and shift the potential */
869 force_switch_constants(
870 6.0, interactionConst.rvdw_switch, interactionConst.rvdw, &interactionConst.dispersion_shift);
871 force_switch_constants(
872 12.0, interactionConst.rvdw_switch, interactionConst.rvdw, &interactionConst.repulsion_shift);
874 case InteractionModifiers::PotSwitch:
875 /* Switch the potential and force */
876 potential_switch_constants(
877 interactionConst.rvdw_switch, interactionConst.rvdw, &interactionConst.vdw_switch);
879 case InteractionModifiers::None:
880 case InteractionModifiers::ExactCutoff:
881 /* Nothing to do here */
883 default: gmx_incons("unimplemented potential modifier");
887 interactionConst.eeltype = ir.coulombtype;
888 interactionConst.coulomb_modifier = ir.coulomb_modifier;
889 interactionConst.rcoulomb = cutoff_inf(ir.rcoulomb);
890 interactionConst.rcoulomb_switch = ir.rcoulomb_switch;
891 interactionConst.epsilon_r = ir.epsilon_r;
893 /* Set the Coulomb energy conversion factor */
894 if (interactionConst.epsilon_r != 0)
896 interactionConst.epsfac = gmx::c_one4PiEps0 / interactionConst.epsilon_r;
900 /* eps = 0 is infinite dieletric: no Coulomb interactions */
901 interactionConst.epsfac = 0;
905 if (EEL_RF(interactionConst.eeltype))
907 GMX_RELEASE_ASSERT(interactionConst.eeltype != CoulombInteractionType::GRFNotused,
908 "GRF is no longer supported");
909 interactionConst.reactionFieldPermitivity = ir.epsilon_rf;
911 interactionConst.epsilon_r,
912 interactionConst.reactionFieldPermitivity,
913 interactionConst.rcoulomb,
914 &interactionConst.reactionFieldCoefficient,
915 &interactionConst.reactionFieldShift);
919 /* For plain cut-off we might use the reaction-field kernels */
920 interactionConst.reactionFieldPermitivity = interactionConst.epsilon_r;
921 interactionConst.reactionFieldCoefficient = 0;
922 if (ir.coulomb_modifier == InteractionModifiers::PotShift)
924 interactionConst.reactionFieldShift = 1 / interactionConst.rcoulomb;
928 interactionConst.reactionFieldShift = 0;
932 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, &interactionConst);
936 real dispersion_shift;
938 dispersion_shift = interactionConst.dispersion_shift.cpot;
939 if (EVDW_PME(interactionConst.vdwtype))
941 dispersion_shift -= interactionConst.sh_lj_ewald;
944 "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
945 interactionConst.repulsion_shift.cpot,
948 if (interactionConst.eeltype == CoulombInteractionType::Cut)
950 fprintf(fp, ", Coulomb %.e", -interactionConst.reactionFieldShift);
952 else if (EEL_PME(interactionConst.eeltype))
954 fprintf(fp, ", Ewald %.3e", -interactionConst.sh_ewald);
959 if (ir.efep != FreeEnergyPerturbationType::No)
961 GMX_RELEASE_ASSERT(ir.fepvals, "ir.fepvals should be set wth free-energy");
962 interactionConst.softCoreParameters =
963 std::make_unique<interaction_const_t::SoftCoreParameters>(*ir.fepvals);
966 return interactionConst;
969 void init_forcerec(FILE* fplog,
970 const gmx::MDLogger& mdlog,
971 t_forcerec* forcerec,
972 const t_inputrec& inputrec,
973 const gmx_mtop_t& mtop,
974 const t_commrec* commrec,
978 gmx::ArrayRef<const std::string> tabbfnm,
981 /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
982 forcerec->use_simd_kernels = GMX_USE_SIMD_KERNELS;
984 if (check_box(inputrec.pbcType, box))
986 gmx_fatal(FARGS, "%s", check_box(inputrec.pbcType, box));
989 /* Test particle insertion ? */
990 if (EI_TPI(inputrec.eI))
992 /* Set to the size of the molecule to be inserted (the last one) */
993 gmx::RangePartitioning molecules = gmx_mtop_molecules(mtop);
994 forcerec->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
1001 if (inputrec.coulombtype == CoulombInteractionType::RFNecUnsupported
1002 || inputrec.coulombtype == CoulombInteractionType::GRFNotused)
1004 gmx_fatal(FARGS, "%s electrostatics is no longer supported", enumValueToString(inputrec.coulombtype));
1007 if (inputrec.bAdress)
1009 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1011 if (inputrec.useTwinRange)
1013 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1015 /* Copy the user determined parameters */
1016 forcerec->userint1 = inputrec.userint1;
1017 forcerec->userint2 = inputrec.userint2;
1018 forcerec->userint3 = inputrec.userint3;
1019 forcerec->userint4 = inputrec.userint4;
1020 forcerec->userreal1 = inputrec.userreal1;
1021 forcerec->userreal2 = inputrec.userreal2;
1022 forcerec->userreal3 = inputrec.userreal3;
1023 forcerec->userreal4 = inputrec.userreal4;
1026 forcerec->fc_stepsize = inputrec.fc_stepsize;
1029 forcerec->efep = inputrec.efep;
1031 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1033 forcerec->use_simd_kernels = FALSE;
1034 if (fplog != nullptr)
1037 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1038 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1039 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1043 forcerec->bBHAM = (mtop.ffparams.functype[0] == F_BHAM);
1045 /* Neighbour searching stuff */
1046 forcerec->pbcType = inputrec.pbcType;
1048 /* Determine if we will do PBC for distances in bonded interactions */
1049 if (forcerec->pbcType == PbcType::No)
1051 forcerec->bMolPBC = FALSE;
1055 const bool useEwaldSurfaceCorrection =
1056 (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0);
1057 const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
1058 if (!DOMAINDECOMP(commrec))
1060 forcerec->bMolPBC = true;
1062 if (useEwaldSurfaceCorrection || haveOrientationRestraints)
1064 forcerec->wholeMoleculeTransform =
1065 std::make_unique<gmx::WholeMoleculeTransform>(mtop, inputrec.pbcType);
1070 forcerec->bMolPBC = dd_bonded_molpbc(*commrec->dd, forcerec->pbcType);
1072 /* With Ewald surface correction it is useful to support e.g. running water
1073 * in parallel with update groups.
1074 * With orientation restraints there is no sensible use case supported with DD.
1076 if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*commrec->dd))
1077 || haveOrientationRestraints)
1080 "You requested Ewald surface correction or orientation restraints, "
1081 "but molecules are broken "
1082 "over periodic boundary conditions by the domain decomposition. "
1083 "Run without domain decomposition instead.");
1087 if (useEwaldSurfaceCorrection)
1089 GMX_RELEASE_ASSERT(!DOMAINDECOMP(commrec) || dd_moleculesAreAlwaysWhole(*commrec->dd),
1090 "Molecules can not be broken by PBC with epsilon_surface > 0");
1094 forcerec->rc_scaling = inputrec.refcoord_scaling;
1095 copy_rvec(inputrec.posres_com, forcerec->posres_com);
1096 copy_rvec(inputrec.posres_comB, forcerec->posres_comB);
1097 forcerec->rlist = cutoff_inf(inputrec.rlist);
1098 forcerec->ljpme_combination_rule = inputrec.ljpme_combination_rule;
1100 /* This now calculates sum for q and c6*/
1101 bool systemHasNetCharge = set_chargesum(fplog, forcerec, mtop);
1103 /* Make data structure used by kernels */
1104 forcerec->ic = std::make_unique<interaction_const_t>(
1105 init_interaction_const(fplog, inputrec, mtop, systemHasNetCharge));
1106 init_interaction_const_tables(fplog, forcerec->ic.get(), forcerec->rlist, inputrec.tabext);
1108 const interaction_const_t* interactionConst = forcerec->ic.get();
1110 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1111 if (inputrec.coulombtype == CoulombInteractionType::Ewald)
1113 forcerec->ewald_table = std::make_unique<gmx_ewald_tab_t>(inputrec, fplog);
1116 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1117 switch (interactionConst->eeltype)
1119 case CoulombInteractionType::Cut:
1120 forcerec->nbkernel_elec_interaction = NbkernelElecType::Coulomb;
1123 case CoulombInteractionType::RF:
1124 case CoulombInteractionType::RFZero:
1125 forcerec->nbkernel_elec_interaction = NbkernelElecType::ReactionField;
1128 case CoulombInteractionType::Switch:
1129 case CoulombInteractionType::Shift:
1130 case CoulombInteractionType::User:
1131 case CoulombInteractionType::PmeSwitch:
1132 case CoulombInteractionType::PmeUser:
1133 case CoulombInteractionType::PmeUserSwitch:
1134 forcerec->nbkernel_elec_interaction = NbkernelElecType::CubicSplineTable;
1137 case CoulombInteractionType::Pme:
1138 case CoulombInteractionType::P3mAD:
1139 case CoulombInteractionType::Ewald:
1140 forcerec->nbkernel_elec_interaction = NbkernelElecType::Ewald;
1145 "Unsupported electrostatic interaction: %s",
1146 enumValueToString(interactionConst->eeltype));
1148 forcerec->nbkernel_elec_modifier = interactionConst->coulomb_modifier;
1150 /* Vdw: Translate from mdp settings to kernel format */
1151 switch (interactionConst->vdwtype)
1153 case VanDerWaalsType::Cut:
1154 if (forcerec->bBHAM)
1156 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::Buckingham;
1160 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LennardJones;
1163 case VanDerWaalsType::Pme:
1164 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LJEwald;
1167 case VanDerWaalsType::Switch:
1168 case VanDerWaalsType::Shift:
1169 case VanDerWaalsType::User:
1170 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::CubicSplineTable;
1174 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", enumValueToString(interactionConst->vdwtype));
1176 forcerec->nbkernel_vdw_modifier = interactionConst->vdw_modifier;
1178 if (!gmx_within_tol(interactionConst->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1180 gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
1182 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1183 * while mdrun does not (and never did) support this.
1185 if (EEL_USER(forcerec->ic->eeltype))
1188 "Electrostatics type %s is currently not supported",
1189 enumValueToString(inputrec.coulombtype));
1192 forcerec->bvdwtab = FALSE;
1193 forcerec->bcoultab = FALSE;
1195 /* 1-4 interaction electrostatics */
1196 forcerec->fudgeQQ = mtop.ffparams.fudgeQQ;
1198 // Multiple time stepping
1199 forcerec->useMts = inputrec.useMts;
1201 if (forcerec->useMts)
1203 GMX_ASSERT(gmx::checkMtsRequirements(inputrec).empty(),
1204 "All MTS requirements should be met here");
1207 const bool haveDirectVirialContributionsFast =
1208 forcerec->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
1209 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || inputrec.nwall > 0 || inputrec.bPull
1210 || inputrec.bRot || inputrec.bIMD;
1211 const bool haveDirectVirialContributionsSlow =
1212 EEL_FULL(interactionConst->eeltype) || EVDW_PME(interactionConst->vdwtype);
1213 for (int i = 0; i < (forcerec->useMts ? 2 : 1); i++)
1215 bool haveDirectVirialContributions =
1216 (((!forcerec->useMts || i == 0) && haveDirectVirialContributionsFast)
1217 || ((!forcerec->useMts || i == 1) && haveDirectVirialContributionsSlow));
1218 forcerec->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
1221 if (forcerec->shift_vec == nullptr)
1223 snew(forcerec->shift_vec, SHIFTS);
1226 if (forcerec->nbfp.empty())
1228 forcerec->ntype = mtop.ffparams.atnr;
1229 forcerec->nbfp = makeNonBondedParameterLists(mtop.ffparams, forcerec->bBHAM);
1230 if (EVDW_PME(interactionConst->vdwtype))
1232 forcerec->ljpme_c6grid = makeLJPmeC6GridCorrectionParameters(mtop.ffparams, *forcerec);
1236 /* Copy the energy group exclusions */
1237 forcerec->egp_flags = inputrec.opts.egp_flags;
1239 /* Van der Waals stuff */
1240 if ((interactionConst->vdwtype != VanDerWaalsType::Cut)
1241 && (interactionConst->vdwtype != VanDerWaalsType::User) && !forcerec->bBHAM)
1243 if (interactionConst->rvdw_switch >= interactionConst->rvdw)
1246 "rvdw_switch (%f) must be < rvdw (%f)",
1247 interactionConst->rvdw_switch,
1248 interactionConst->rvdw);
1253 "Using %s Lennard-Jones, switch between %g and %g nm\n",
1254 (interactionConst->eeltype == CoulombInteractionType::Switch) ? "switched" : "shifted",
1255 interactionConst->rvdw_switch,
1256 interactionConst->rvdw);
1260 if (forcerec->bBHAM && EVDW_PME(interactionConst->vdwtype))
1262 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1266 && (interactionConst->vdwtype == VanDerWaalsType::Shift
1267 || interactionConst->vdwtype == VanDerWaalsType::Switch))
1269 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1272 if (forcerec->bBHAM)
1274 gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
1277 if (inputrec.implicit_solvent)
1279 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1283 /* This code automatically gives table length tabext without cut-off's,
1284 * in that case grompp should already have checked that we do not need
1285 * normal tables and we only generate tables for 1-4 interactions.
1287 real rtab = inputrec.rlist + inputrec.tabext;
1289 /* We want to use unmodified tables for 1-4 coulombic
1290 * interactions, so we must in general have an extra set of
1292 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1293 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1295 forcerec->pairsTable = make_tables(fplog, interactionConst, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1299 forcerec->nwall = inputrec.nwall;
1300 if (inputrec.nwall && inputrec.wall_type == WallType::Table)
1302 make_wall_tables(fplog, inputrec, tabfn, &mtop.groups, forcerec);
1305 forcerec->fcdata = std::make_unique<t_fcdata>();
1307 if (!tabbfnm.empty())
1309 t_fcdata& fcdata = *forcerec->fcdata;
1310 // Need to catch std::bad_alloc
1311 // TODO Don't need to catch this here, when merging with master branch
1314 // TODO move these tables into a separate struct and store reference in ListedForces
1315 fcdata.bondtab = make_bonded_tables(fplog, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1316 fcdata.angletab = make_bonded_tables(fplog, F_TABANGLES, -1, mtop, tabbfnm, "a");
1317 fcdata.dihtab = make_bonded_tables(fplog, F_TABDIHS, -1, mtop, tabbfnm, "d");
1319 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1326 "No fcdata or table file name passed, can not read table, can not do bonded "
1331 /* Initialize the thread working data for bonded interactions */
1332 if (forcerec->useMts)
1334 // Add one ListedForces object for each MTS level
1335 bool isFirstLevel = true;
1336 for (const auto& mtsLevel : inputrec.mtsLevels)
1338 ListedForces::InteractionSelection interactionSelection;
1339 const auto& forceGroups = mtsLevel.forceGroups;
1340 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1342 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1344 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1346 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1348 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1350 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1354 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1355 isFirstLevel = false;
1357 forcerec->listedForces.emplace_back(
1359 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1360 gmx_omp_nthreads_get(emntBonded),
1361 interactionSelection,
1367 // Add one ListedForces object with all listed interactions
1368 forcerec->listedForces.emplace_back(
1370 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1371 gmx_omp_nthreads_get(emntBonded),
1372 ListedForces::interactionSelectionAll(),
1376 // QM/MM initialization if requested
1379 gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1382 /* Set all the static charge group info */
1383 forcerec->cginfo_mb = init_cginfo_mb(mtop, forcerec);
1384 if (!DOMAINDECOMP(commrec))
1386 forcerec->cginfo = cginfo_expand(mtop.molblock.size(), forcerec->cginfo_mb);
1389 if (!DOMAINDECOMP(commrec))
1391 forcerec_set_ranges(forcerec, mtop.natoms, mtop.natoms, mtop.natoms);
1394 forcerec->print_force = print_force;
1396 forcerec->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1397 snew(forcerec->ewc_t, forcerec->nthread_ewc);
1399 if (inputrec.eDispCorr != DispersionCorrectionType::No)
1401 forcerec->dispersionCorrection = std::make_unique<DispersionCorrection>(
1402 mtop, inputrec, forcerec->bBHAM, forcerec->ntype, forcerec->nbfp, *forcerec->ic, tabfn);
1403 forcerec->dispersionCorrection->print(mdlog);
1406 if (fplog != nullptr)
1408 /* Here we switch from using mdlog, which prints the newline before
1409 * the paragraph, to our old fprintf logging, which prints the newline
1410 * after the paragraph, so we should add a newline here.
1412 fprintf(fplog, "\n");
1416 t_forcerec::t_forcerec() = default;
1418 t_forcerec::~t_forcerec()
1420 /* Note: This code will disappear when types are converted to C++ */