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, 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/nbnxm/nbnxm.h"
87 #include "gromacs/pbcutil/ishift.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/tables/forcetable.h"
90 #include "gromacs/topology/mtop_util.h"
91 #include "gromacs/trajectory/trajectoryframe.h"
92 #include "gromacs/utility/cstringutil.h"
93 #include "gromacs/utility/exceptions.h"
94 #include "gromacs/utility/fatalerror.h"
95 #include "gromacs/utility/gmxassert.h"
96 #include "gromacs/utility/logger.h"
97 #include "gromacs/utility/physicalnodecommunicator.h"
98 #include "gromacs/utility/pleasecite.h"
99 #include "gromacs/utility/smalloc.h"
100 #include "gromacs/utility/strconvert.h"
102 #include "gpuforcereduction.h"
104 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
105 haveDirectVirialContributions_(haveDirectVirialContributions)
107 shiftForces_.resize(SHIFTS);
110 void ForceHelperBuffers::resize(int numAtoms)
112 if (haveDirectVirialContributions_)
114 forceBufferForDirectVirialContributions_.resize(numAtoms);
118 static std::vector<real> mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
120 std::vector<real> nbfp;
126 nbfp.resize(3 * atnr * atnr);
128 for (int i = 0; (i < atnr); i++)
130 for (int j = 0; (j < atnr); j++, k++)
132 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
133 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
134 /* nbfp now includes the 6.0 derivative prefactor */
135 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c * 6.0;
141 nbfp.resize(2 * atnr * atnr);
143 for (int i = 0; (i < atnr); i++)
145 for (int j = 0; (j < atnr); j++, k++)
147 /* nbfp now includes the 6.0/12.0 derivative prefactors */
148 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6 * 6.0;
149 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12 * 12.0;
157 static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
160 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
163 /* For LJ-PME simulations, we correct the energies with the reciprocal space
164 * inside of the cut-off. To do this the non-bonded kernels needs to have
165 * access to the C6-values used on the reciprocal grid in pme.c
169 snew(grid, 2 * atnr * atnr);
170 for (i = k = 0; (i < atnr); i++)
172 for (j = 0; (j < atnr); j++, k++)
174 c6i = idef->iparams[i * (atnr + 1)].lj.c6;
175 c12i = idef->iparams[i * (atnr + 1)].lj.c12;
176 c6j = idef->iparams[j * (atnr + 1)].lj.c6;
177 c12j = idef->iparams[j * (atnr + 1)].lj.c12;
178 c6 = std::sqrt(c6i * c6j);
179 if (fr->ljpme_combination_rule == eljpmeLB && !gmx_numzero(c6) && !gmx_numzero(c12i)
180 && !gmx_numzero(c12j))
182 sigmai = gmx::sixthroot(c12i / c6i);
183 sigmaj = gmx::sixthroot(c12j / c6j);
184 epsi = c6i * c6i / c12i;
185 epsj = c6j * c6j / c12j;
186 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
188 /* Store the elements at the same relative positions as C6 in nbfp in order
189 * to simplify access in the kernels
191 grid[2 * (atnr * i + j)] = c6 * 6.0;
204 static std::vector<cginfo_mb_t> init_cginfo_mb(const gmx_mtop_t* mtop, const t_forcerec* fr)
209 snew(type_VDW, fr->ntype);
210 for (int ai = 0; ai < fr->ntype; ai++)
212 type_VDW[ai] = FALSE;
213 for (int j = 0; j < fr->ntype; j++)
215 type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
216 || C12(fr->nbfp, fr->ntype, ai, j) != 0;
220 std::vector<cginfo_mb_t> cginfoPerMolblock;
222 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
224 const gmx_molblock_t& molb = mtop->molblock[mb];
225 const gmx_moltype_t& molt = mtop->moltype[molb.type];
226 const auto& excl = molt.excls;
228 /* Check if the cginfo is identical for all molecules in this block.
229 * If so, we only need an array of the size of one molecule.
230 * Otherwise we make an array of #mol times #cgs per molecule.
233 for (int m = 0; m < molb.nmol; m++)
235 const int am = m * molt.atoms.nr;
236 for (int a = 0; a < molt.atoms.nr; a++)
238 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
239 != getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
243 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
245 if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
246 != mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
254 cginfo_mb_t cginfo_mb;
255 cginfo_mb.cg_start = a_offset;
256 cginfo_mb.cg_end = a_offset + molb.nmol * molt.atoms.nr;
257 cginfo_mb.cg_mod = (bId ? 1 : molb.nmol) * molt.atoms.nr;
258 cginfo_mb.cginfo.resize(cginfo_mb.cg_mod);
259 gmx::ArrayRef<int> cginfo = cginfo_mb.cginfo;
261 /* Set constraints flags for constrained atoms */
262 snew(a_con, molt.atoms.nr);
263 for (int ftype = 0; ftype < F_NRE; ftype++)
265 if (interaction_function[ftype].flags & IF_CONSTRAINT)
267 const int nral = NRAL(ftype);
268 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
272 for (a = 0; a < nral; a++)
274 a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
275 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
281 for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
283 const int molculeOffsetInBlock = m * molt.atoms.nr;
284 for (int a = 0; a < molt.atoms.nr; a++)
286 const t_atom& atom = molt.atoms.atom[a];
287 int& atomInfo = cginfo[molculeOffsetInBlock + a];
289 /* Store the energy group in cginfo */
290 int gid = getGroupType(mtop->groups,
291 SimulationAtomGroupType::EnergyOutput,
292 a_offset + molculeOffsetInBlock + a);
293 SET_CGINFO_GID(atomInfo, gid);
295 bool bHaveVDW = (type_VDW[atom.type] || type_VDW[atom.typeB]);
296 bool bHaveQ = (atom.q != 0 || atom.qB != 0);
298 bool haveExclusions = false;
299 /* Loop over all the exclusions of atom ai */
300 for (const int j : excl[a])
304 haveExclusions = true;
311 case acCONSTRAINT: SET_CGINFO_CONSTR(atomInfo); break;
312 case acSETTLE: SET_CGINFO_SETTLE(atomInfo); break;
317 SET_CGINFO_EXCL_INTER(atomInfo);
321 SET_CGINFO_HAS_VDW(atomInfo);
325 SET_CGINFO_HAS_Q(atomInfo);
327 if (fr->efep != efepNO && PERTURBED(atom))
329 SET_CGINFO_FEP(atomInfo);
336 cginfoPerMolblock.push_back(cginfo_mb);
338 a_offset += molb.nmol * molt.atoms.nr;
342 return cginfoPerMolblock;
345 static std::vector<int> cginfo_expand(const int nmb, gmx::ArrayRef<const cginfo_mb_t> cgi_mb)
347 const int ncg = cgi_mb[nmb - 1].cg_end;
349 std::vector<int> cginfo(ncg);
352 for (int cg = 0; cg < ncg; cg++)
354 while (cg >= cgi_mb[mb].cg_end)
358 cginfo[cg] = cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
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 std::vector<bondedtable_t> make_bonded_tables(FILE* fplog,
560 const gmx_mtop_t* mtop,
561 gmx::ArrayRef<const std::string> tabbfnm,
564 std::vector<bondedtable_t> tab;
567 int* count = nullptr;
568 count_tables(ftype1, ftype2, mtop, &ncount, &count);
570 // Are there any relevant tabulated bond interactions?
574 for (int i = 0; i < ncount; i++)
576 // Do any interactions exist that requires this table?
579 // This pattern enforces the current requirement that
580 // table filenames end in a characteristic sequence
581 // before the file type extension, and avoids table 13
582 // being recognized and used for table 1.
583 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
584 bool madeTable = false;
585 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
587 if (gmx::endsWith(tabbfnm[j], patternToFind))
589 // Finally read the table from the file found
590 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
596 bool isPlural = (ftype2 != -1);
598 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
599 "because no table file whose name matched '%s' was passed via the "
600 "gmx mdrun -tableb command-line option.",
601 interaction_function[ftype1].longname,
602 isPlural ? "' or '" : "",
603 isPlural ? interaction_function[ftype2].longname : "",
605 patternToFind.c_str());
615 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
617 fr->natoms_force = natoms_force;
618 fr->natoms_force_constr = natoms_force_constr;
620 for (auto& forceHelperBuffers : fr->forceHelperBuffers)
622 forceHelperBuffers.resize(natoms_f_novirsum);
626 static real cutoff_inf(real cutoff)
630 cutoff = GMX_CUTOFF_INF;
636 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
637 static void initCoulombEwaldParameters(FILE* fp,
638 const t_inputrec* ir,
639 bool systemHasNetCharge,
640 interaction_const_t* ic)
642 if (!EEL_PME_EWALD(ir->coulombtype))
649 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
651 if (ir->coulombtype == eelP3M_AD)
653 please_cite(fp, "Hockney1988");
654 please_cite(fp, "Ballenegger2012");
658 please_cite(fp, "Essmann95a");
661 if (ir->ewald_geometry == eewg3DC)
666 "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
667 systemHasNetCharge ? " and net charge" : "");
669 please_cite(fp, "In-Chul99a");
670 if (systemHasNetCharge)
672 please_cite(fp, "Ballenegger2009");
677 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
680 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
683 if (ic->coulomb_modifier == eintmodPOTSHIFT)
685 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
686 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
694 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
695 static void initVdwEwaldParameters(FILE* fp, const t_inputrec* ir, interaction_const_t* ic)
697 if (!EVDW_PME(ir->vdwtype))
704 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
705 please_cite(fp, "Essmann95a");
707 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
710 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
713 if (ic->vdw_modifier == eintmodPOTSHIFT)
715 real crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
716 ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
724 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
726 * Tables are generated for one or both, depending on if the pointers
727 * are non-null. The spacing for both table sets is the same and obeys
728 * both accuracy requirements, when relevant.
730 static void init_ewald_f_table(const interaction_const_t& ic,
733 EwaldCorrectionTables* coulombTables,
734 EwaldCorrectionTables* vdwTables)
736 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
737 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
739 /* Get the Ewald table spacing based on Coulomb and/or LJ
740 * Ewald coefficients and rtol.
742 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
744 const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
746 real tableLen = ic.rcoulomb;
747 if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
749 /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
750 * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
751 * The alternative is to look through all the exclusions and check if they come from
752 * couple-intramol == no. Meanwhile, always having larger tables should only affect
753 * memory consumption, not speed (barring cache issues).
755 tableLen = rlist + tabext;
757 const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
762 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
767 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
771 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
773 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
776 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
780 "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
781 1 / ic->coulombEwaldTables->scale,
782 ic->coulombEwaldTables->tableF.size());
787 static void clear_force_switch_constants(shift_consts_t* sc)
794 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
796 /* Here we determine the coefficient for shifting the force to zero
797 * between distance rsw and the cut-off rc.
798 * For a potential of r^-p, we have force p*r^-(p+1).
799 * But to save flops we absorb p in the coefficient.
801 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
802 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
804 sc->c2 = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
805 sc->c3 = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
806 sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
807 + p * sc->c3 / 4 * gmx::power4(rc - rsw);
810 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
812 /* The switch function is 1 at rsw and 0 at rc.
813 * The derivative and second derivate are zero at both ends.
814 * rsw = max(r - r_switch, 0)
815 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
816 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
817 * force = force*dsw - potential*sw
820 sc->c3 = -10 / gmx::power3(rc - rsw);
821 sc->c4 = 15 / gmx::power4(rc - rsw);
822 sc->c5 = -6 / gmx::power5(rc - rsw);
825 /*! \brief Construct interaction constants
827 * This data is used (particularly) by search and force code for
828 * short-range interactions. Many of these are constant for the whole
829 * simulation; some are constant only after PME tuning completes.
831 static void init_interaction_const(FILE* fp,
832 interaction_const_t** interaction_const,
833 const t_inputrec* ir,
834 const gmx_mtop_t* mtop,
835 bool systemHasNetCharge)
837 interaction_const_t* ic = new interaction_const_t;
839 ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
840 ic->vdwEwaldTables = std::make_unique<EwaldCorrectionTables>();
843 ic->vdwtype = ir->vdwtype;
844 ic->vdw_modifier = ir->vdw_modifier;
845 ic->reppow = mtop->ffparams.reppow;
846 ic->rvdw = cutoff_inf(ir->rvdw);
847 ic->rvdw_switch = ir->rvdw_switch;
848 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
849 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
850 if (ic->useBuckingham)
852 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
855 initVdwEwaldParameters(fp, ir, ic);
857 clear_force_switch_constants(&ic->dispersion_shift);
858 clear_force_switch_constants(&ic->repulsion_shift);
860 switch (ic->vdw_modifier)
862 case eintmodPOTSHIFT:
863 /* Only shift the potential, don't touch the force */
864 ic->dispersion_shift.cpot = -1.0 / gmx::power6(ic->rvdw);
865 ic->repulsion_shift.cpot = -1.0 / gmx::power12(ic->rvdw);
867 case eintmodFORCESWITCH:
868 /* Switch the force, switch and shift the potential */
869 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw, &ic->dispersion_shift);
870 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw, &ic->repulsion_shift);
872 case eintmodPOTSWITCH:
873 /* Switch the potential and force */
874 potential_switch_constants(ic->rvdw_switch, ic->rvdw, &ic->vdw_switch);
877 case eintmodEXACTCUTOFF:
878 /* Nothing to do here */
880 default: gmx_incons("unimplemented potential modifier");
884 ic->eeltype = ir->coulombtype;
885 ic->coulomb_modifier = ir->coulomb_modifier;
886 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
887 ic->rcoulomb_switch = ir->rcoulomb_switch;
888 ic->epsilon_r = ir->epsilon_r;
890 /* Set the Coulomb energy conversion factor */
891 if (ic->epsilon_r != 0)
893 ic->epsfac = ONE_4PI_EPS0 / ic->epsilon_r;
897 /* eps = 0 is infinite dieletric: no Coulomb interactions */
902 if (EEL_RF(ic->eeltype))
904 GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
905 ic->epsilon_rf = ir->epsilon_rf;
907 calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf, ic->rcoulomb, &ic->k_rf, &ic->c_rf);
911 /* For plain cut-off we might use the reaction-field kernels */
912 ic->epsilon_rf = ic->epsilon_r;
914 if (ir->coulomb_modifier == eintmodPOTSHIFT)
916 ic->c_rf = 1 / ic->rcoulomb;
924 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
928 real dispersion_shift;
930 dispersion_shift = ic->dispersion_shift.cpot;
931 if (EVDW_PME(ic->vdwtype))
933 dispersion_shift -= ic->sh_lj_ewald;
935 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e", ic->repulsion_shift.cpot, dispersion_shift);
937 if (ic->eeltype == eelCUT)
939 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
941 else if (EEL_PME(ic->eeltype))
943 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
948 if (ir->efep != efepNO)
950 GMX_RELEASE_ASSERT(ir->fepvals, "ir->fepvals should be set wth free-energy");
951 ic->softCoreParameters = std::make_unique<interaction_const_t::SoftCoreParameters>(*ir->fepvals);
954 *interaction_const = ic;
957 void init_forcerec(FILE* fp,
958 const gmx::MDLogger& mdlog,
960 const t_inputrec* ir,
961 const gmx_mtop_t* mtop,
966 gmx::ArrayRef<const std::string> tabbfnm,
969 /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
970 fr->use_simd_kernels = GMX_USE_SIMD_KERNELS;
972 if (check_box(ir->pbcType, box))
974 gmx_fatal(FARGS, "%s", check_box(ir->pbcType, box));
977 /* Test particle insertion ? */
980 /* Set to the size of the molecule to be inserted (the last one) */
981 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
982 fr->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
989 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED || ir->coulombtype == eelGRF_NOTUSED)
991 gmx_fatal(FARGS, "%s electrostatics is no longer supported", eel_names[ir->coulombtype]);
996 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
998 if (ir->useTwinRange)
1000 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1002 /* Copy the user determined parameters */
1003 fr->userint1 = ir->userint1;
1004 fr->userint2 = ir->userint2;
1005 fr->userint3 = ir->userint3;
1006 fr->userint4 = ir->userint4;
1007 fr->userreal1 = ir->userreal1;
1008 fr->userreal2 = ir->userreal2;
1009 fr->userreal3 = ir->userreal3;
1010 fr->userreal4 = ir->userreal4;
1013 fr->fc_stepsize = ir->fc_stepsize;
1016 fr->efep = ir->efep;
1018 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1020 fr->use_simd_kernels = FALSE;
1024 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1025 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1026 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1030 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1032 /* Neighbour searching stuff */
1033 fr->pbcType = ir->pbcType;
1035 /* Determine if we will do PBC for distances in bonded interactions */
1036 if (fr->pbcType == PbcType::No)
1038 fr->bMolPBC = FALSE;
1042 const bool useEwaldSurfaceCorrection =
1043 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1044 const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
1045 if (!DOMAINDECOMP(cr))
1049 if (useEwaldSurfaceCorrection || haveOrientationRestraints)
1051 fr->wholeMoleculeTransform =
1052 std::make_unique<gmx::WholeMoleculeTransform>(*mtop, ir->pbcType);
1057 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->pbcType);
1059 /* With Ewald surface correction it is useful to support e.g. running water
1060 * in parallel with update groups.
1061 * With orientation restraints there is no sensible use case supported with DD.
1063 if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd)) || haveOrientationRestraints)
1066 "You requested Ewald surface correction or orientation restraints, "
1067 "but molecules are broken "
1068 "over periodic boundary conditions by the domain decomposition. "
1069 "Run without domain decomposition instead.");
1073 if (useEwaldSurfaceCorrection)
1075 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || dd_moleculesAreAlwaysWhole(*cr->dd),
1076 "Molecules can not be broken by PBC with epsilon_surface > 0");
1080 fr->rc_scaling = ir->refcoord_scaling;
1081 copy_rvec(ir->posres_com, fr->posres_com);
1082 copy_rvec(ir->posres_comB, fr->posres_comB);
1083 fr->rlist = cutoff_inf(ir->rlist);
1084 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1086 /* This now calculates sum for q and c6*/
1087 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1089 /* fr->ic is used both by verlet and group kernels (to some extent) now */
1090 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1091 init_interaction_const_tables(fp, fr->ic, fr->rlist, ir->tabext);
1093 const interaction_const_t* ic = fr->ic;
1095 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1096 if (ir->coulombtype == eelEWALD)
1098 init_ewald_tab(&(fr->ewald_table), ir, fp);
1101 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1102 switch (ic->eeltype)
1104 case eelCUT: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB; break;
1107 case eelRF_ZERO: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD; break;
1114 case eelPMEUSERSWITCH:
1115 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1120 case eelEWALD: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD; break;
1123 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1125 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1127 /* Vdw: Translate from mdp settings to kernel format */
1128 switch (ic->vdwtype)
1133 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1137 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1140 case evdwPME: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD; break;
1144 case evdwUSER: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE; break;
1146 default: gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1148 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1150 if (!gmx_within_tol(ic->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1152 gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
1154 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1155 * while mdrun does not (and never did) support this.
1157 if (EEL_USER(fr->ic->eeltype))
1159 gmx_fatal(FARGS, "Electrostatics type %s is currently not supported", eel_names[ir->coulombtype]);
1162 fr->bvdwtab = FALSE;
1163 fr->bcoultab = FALSE;
1165 /* 1-4 interaction electrostatics */
1166 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1168 // Multiple time stepping
1169 fr->useMts = ir->useMts;
1173 GMX_ASSERT(gmx::checkMtsRequirements(*ir).empty(),
1174 "All MTS requirements should be met here");
1177 const bool haveDirectVirialContributionsFast =
1178 fr->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
1179 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || ir->nwall > 0 || ir->bPull || ir->bRot
1181 const bool haveDirectVirialContributionsSlow = EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype);
1182 for (int i = 0; i < (fr->useMts ? 2 : 1); i++)
1184 bool haveDirectVirialContributions =
1185 (((!fr->useMts || i == 0) && haveDirectVirialContributionsFast)
1186 || ((!fr->useMts || i == 1) && haveDirectVirialContributionsSlow));
1187 fr->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
1190 if (fr->shift_vec == nullptr)
1192 snew(fr->shift_vec, SHIFTS);
1195 if (fr->nbfp.empty())
1197 fr->ntype = mtop->ffparams.atnr;
1198 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1199 if (EVDW_PME(ic->vdwtype))
1201 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1205 /* Copy the energy group exclusions */
1206 fr->egp_flags = ir->opts.egp_flags;
1208 /* Van der Waals stuff */
1209 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1211 if (ic->rvdw_switch >= ic->rvdw)
1213 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)", ic->rvdw_switch, ic->rvdw);
1218 "Using %s Lennard-Jones, switch between %g and %g nm\n",
1219 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
1225 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1227 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1230 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1232 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1237 gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
1240 if (ir->implicit_solvent)
1242 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1246 /* This code automatically gives table length tabext without cut-off's,
1247 * in that case grompp should already have checked that we do not need
1248 * normal tables and we only generate tables for 1-4 interactions.
1250 real rtab = ir->rlist + ir->tabext;
1252 /* We want to use unmodified tables for 1-4 coulombic
1253 * interactions, so we must in general have an extra set of
1255 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1256 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1258 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1262 fr->nwall = ir->nwall;
1263 if (ir->nwall && ir->wall_type == ewtTABLE)
1265 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1268 fr->fcdata = std::make_unique<t_fcdata>();
1270 if (!tabbfnm.empty())
1272 t_fcdata& fcdata = *fr->fcdata;
1273 // Need to catch std::bad_alloc
1274 // TODO Don't need to catch this here, when merging with master branch
1277 // TODO move these tables into a separate struct and store reference in ListedForces
1278 fcdata.bondtab = make_bonded_tables(fp, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1279 fcdata.angletab = make_bonded_tables(fp, F_TABANGLES, -1, mtop, tabbfnm, "a");
1280 fcdata.dihtab = make_bonded_tables(fp, F_TABDIHS, -1, mtop, tabbfnm, "d");
1282 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1289 "No fcdata or table file name passed, can not read table, can not do bonded "
1294 /* Initialize the thread working data for bonded interactions */
1297 // Add one ListedForces object for each MTS level
1298 bool isFirstLevel = true;
1299 for (const auto& mtsLevel : ir->mtsLevels)
1301 ListedForces::InteractionSelection interactionSelection;
1302 const auto& forceGroups = mtsLevel.forceGroups;
1303 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1305 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1307 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1309 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1311 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1313 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1317 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1318 isFirstLevel = false;
1320 fr->listedForces.emplace_back(mtop->ffparams,
1321 mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1322 gmx_omp_nthreads_get(emntBonded),
1323 interactionSelection,
1329 // Add one ListedForces object with all listed interactions
1330 fr->listedForces.emplace_back(mtop->ffparams,
1331 mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1332 gmx_omp_nthreads_get(emntBonded),
1333 ListedForces::interactionSelectionAll(),
1337 // QM/MM initialization if requested
1340 gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1343 /* Set all the static charge group info */
1344 fr->cginfo_mb = init_cginfo_mb(mtop, fr);
1345 if (!DOMAINDECOMP(cr))
1347 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1350 if (!DOMAINDECOMP(cr))
1352 forcerec_set_ranges(fr, mtop->natoms, mtop->natoms, mtop->natoms);
1355 fr->print_force = print_force;
1357 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1358 snew(fr->ewc_t, fr->nthread_ewc);
1360 if (ir->eDispCorr != edispcNO)
1362 fr->dispersionCorrection = std::make_unique<DispersionCorrection>(
1363 *mtop, *ir, fr->bBHAM, fr->ntype, fr->nbfp, *fr->ic, tabfn);
1364 fr->dispersionCorrection->print(mdlog);
1369 /* Here we switch from using mdlog, which prints the newline before
1370 * the paragraph, to our old fprintf logging, which prints the newline
1371 * after the paragraph, so we should add a newline here.
1377 t_forcerec::t_forcerec() = default;
1379 t_forcerec::~t_forcerec()
1381 /* Note: This code will disappear when types are converted to C++ */