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/nbnxm/gpu_data_mgmt.h"
86 #include "gromacs/nbnxm/nbnxm.h"
87 #include "gromacs/nbnxm/nbnxm_geometry.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 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
104 haveDirectVirialContributions_(haveDirectVirialContributions)
106 shiftForces_.resize(SHIFTS);
109 void ForceHelperBuffers::resize(int numAtoms)
111 if (haveDirectVirialContributions_)
113 forceBufferForDirectVirialContributions_.resize(numAtoms);
117 static std::vector<real> mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
119 std::vector<real> nbfp;
125 nbfp.resize(3 * atnr * atnr);
127 for (int i = 0; (i < atnr); i++)
129 for (int j = 0; (j < atnr); j++, k++)
131 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
132 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
133 /* nbfp now includes the 6.0 derivative prefactor */
134 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c * 6.0;
140 nbfp.resize(2 * atnr * atnr);
142 for (int i = 0; (i < atnr); i++)
144 for (int j = 0; (j < atnr); j++, k++)
146 /* nbfp now includes the 6.0/12.0 derivative prefactors */
147 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6 * 6.0;
148 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12 * 12.0;
156 static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
159 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
162 /* For LJ-PME simulations, we correct the energies with the reciprocal space
163 * inside of the cut-off. To do this the non-bonded kernels needs to have
164 * access to the C6-values used on the reciprocal grid in pme.c
168 snew(grid, 2 * atnr * atnr);
169 for (i = k = 0; (i < atnr); i++)
171 for (j = 0; (j < atnr); j++, k++)
173 c6i = idef->iparams[i * (atnr + 1)].lj.c6;
174 c12i = idef->iparams[i * (atnr + 1)].lj.c12;
175 c6j = idef->iparams[j * (atnr + 1)].lj.c6;
176 c12j = idef->iparams[j * (atnr + 1)].lj.c12;
177 c6 = std::sqrt(c6i * c6j);
178 if (fr->ljpme_combination_rule == eljpmeLB && !gmx_numzero(c6) && !gmx_numzero(c12i)
179 && !gmx_numzero(c12j))
181 sigmai = gmx::sixthroot(c12i / c6i);
182 sigmaj = gmx::sixthroot(c12j / c6j);
183 epsi = c6i * c6i / c12i;
184 epsj = c6j * c6j / c12j;
185 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
187 /* Store the elements at the same relative positions as C6 in nbfp in order
188 * to simplify access in the kernels
190 grid[2 * (atnr * i + j)] = c6 * 6.0;
203 static std::vector<cginfo_mb_t> init_cginfo_mb(const gmx_mtop_t* mtop, const t_forcerec* fr)
208 snew(type_VDW, fr->ntype);
209 for (int ai = 0; ai < fr->ntype; ai++)
211 type_VDW[ai] = FALSE;
212 for (int j = 0; j < fr->ntype; j++)
214 type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
215 || C12(fr->nbfp, fr->ntype, ai, j) != 0;
219 std::vector<cginfo_mb_t> cginfoPerMolblock;
221 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
223 const gmx_molblock_t& molb = mtop->molblock[mb];
224 const gmx_moltype_t& molt = mtop->moltype[molb.type];
225 const auto& excl = molt.excls;
227 /* Check if the cginfo is identical for all molecules in this block.
228 * If so, we only need an array of the size of one molecule.
229 * Otherwise we make an array of #mol times #cgs per molecule.
232 for (int m = 0; m < molb.nmol; m++)
234 const int am = m * molt.atoms.nr;
235 for (int a = 0; a < molt.atoms.nr; a++)
237 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
238 != getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
242 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
244 if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
245 != mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
253 cginfo_mb_t cginfo_mb;
254 cginfo_mb.cg_start = a_offset;
255 cginfo_mb.cg_end = a_offset + molb.nmol * molt.atoms.nr;
256 cginfo_mb.cg_mod = (bId ? 1 : molb.nmol) * molt.atoms.nr;
257 cginfo_mb.cginfo.resize(cginfo_mb.cg_mod);
258 gmx::ArrayRef<int> cginfo = cginfo_mb.cginfo;
260 /* Set constraints flags for constrained atoms */
261 snew(a_con, molt.atoms.nr);
262 for (int ftype = 0; ftype < F_NRE; ftype++)
264 if (interaction_function[ftype].flags & IF_CONSTRAINT)
266 const int nral = NRAL(ftype);
267 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
271 for (a = 0; a < nral; a++)
273 a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
274 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
280 for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
282 const int molculeOffsetInBlock = m * molt.atoms.nr;
283 for (int a = 0; a < molt.atoms.nr; a++)
285 const t_atom& atom = molt.atoms.atom[a];
286 int& atomInfo = cginfo[molculeOffsetInBlock + a];
288 /* Store the energy group in cginfo */
289 int gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput,
290 a_offset + molculeOffsetInBlock + a);
291 SET_CGINFO_GID(atomInfo, gid);
293 bool bHaveVDW = (type_VDW[atom.type] || type_VDW[atom.typeB]);
294 bool bHaveQ = (atom.q != 0 || atom.qB != 0);
296 bool haveExclusions = false;
297 /* Loop over all the exclusions of atom ai */
298 for (const int j : excl[a])
302 haveExclusions = true;
309 case acCONSTRAINT: SET_CGINFO_CONSTR(atomInfo); break;
310 case acSETTLE: SET_CGINFO_SETTLE(atomInfo); break;
315 SET_CGINFO_EXCL_INTER(atomInfo);
319 SET_CGINFO_HAS_VDW(atomInfo);
323 SET_CGINFO_HAS_Q(atomInfo);
325 if (fr->efep != efepNO && PERTURBED(atom))
327 SET_CGINFO_FEP(atomInfo);
334 cginfoPerMolblock.push_back(cginfo_mb);
336 a_offset += molb.nmol * molt.atoms.nr;
340 return cginfoPerMolblock;
343 static std::vector<int> cginfo_expand(const int nmb, gmx::ArrayRef<const cginfo_mb_t> cgi_mb)
345 const int ncg = cgi_mb[nmb - 1].cg_end;
347 std::vector<int> cginfo(ncg);
350 for (int cg = 0; cg < ncg; cg++)
352 while (cg >= cgi_mb[mb].cg_end)
356 cginfo[cg] = cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
362 /* Sets the sum of charges (squared) and C6 in the system in fr.
363 * Returns whether the system has a net charge.
365 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t* mtop)
367 /*This now calculates sum for q and c6*/
368 double qsum, q2sum, q, c6sum, c6;
373 for (const gmx_molblock_t& molb : mtop->molblock)
375 int nmol = molb.nmol;
376 const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
377 for (int i = 0; i < atoms->nr; i++)
379 q = atoms->atom[i].q;
381 q2sum += nmol * q * q;
382 c6 = mtop->ffparams.iparams[atoms->atom[i].type * (mtop->ffparams.atnr + 1)].lj.c6;
387 fr->q2sum[0] = q2sum;
388 fr->c6sum[0] = c6sum;
390 if (fr->efep != efepNO)
395 for (const gmx_molblock_t& molb : mtop->molblock)
397 int nmol = molb.nmol;
398 const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
399 for (int i = 0; i < atoms->nr; i++)
401 q = atoms->atom[i].qB;
403 q2sum += nmol * q * q;
404 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB * (mtop->ffparams.atnr + 1)].lj.c6;
408 fr->q2sum[1] = q2sum;
409 fr->c6sum[1] = c6sum;
414 fr->qsum[1] = fr->qsum[0];
415 fr->q2sum[1] = fr->q2sum[0];
416 fr->c6sum[1] = fr->c6sum[0];
420 if (fr->efep == efepNO)
422 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
426 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
430 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
431 return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
434 static real calcBuckinghamBMax(FILE* fplog, const gmx_mtop_t* mtop)
436 const t_atoms *at1, *at2;
437 int i, j, tpi, tpj, ntypes;
442 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
444 ntypes = mtop->ffparams.atnr;
448 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
450 at1 = &mtop->moltype[mt1].atoms;
451 for (i = 0; (i < at1->nr); i++)
453 tpi = at1->atom[i].type;
456 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
459 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
461 at2 = &mtop->moltype[mt2].atoms;
462 for (j = 0; (j < at2->nr); j++)
464 tpj = at2->atom[j].type;
467 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
469 b = mtop->ffparams.iparams[tpi * ntypes + tpj].bham.b;
474 if ((b < bmin) || (bmin == -1))
484 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n", bmin, bham_b_max);
490 /*!\brief If there's bonded interactions of type \c ftype1 or \c
491 * ftype2 present in the topology, build an array of the number of
492 * interactions present for each bonded interaction index found in the
495 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
496 * valid type with that parameter.
498 * \c count will be reallocated as necessary to fit the largest bonded
499 * interaction index found, and its current size will be returned in
500 * \c ncount. It will contain zero for every bonded interaction index
501 * for which no interactions are present in the topology.
503 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t* mtop, int* ncount, int** count)
505 int ftype, i, j, tabnr;
507 // Loop over all moleculetypes
508 for (const gmx_moltype_t& molt : mtop->moltype)
510 // Loop over all interaction types
511 for (ftype = 0; ftype < F_NRE; ftype++)
513 // If the current interaction type is one of the types whose tables we're trying to count...
514 if (ftype == ftype1 || ftype == ftype2)
516 const InteractionList& il = molt.ilist[ftype];
517 const int stride = 1 + NRAL(ftype);
518 // ... and there are actually some interactions for this type
519 for (i = 0; i < il.size(); i += stride)
521 // Find out which table index the user wanted
522 tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
525 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
527 // Make room for this index in the data structure
528 if (tabnr >= *ncount)
530 srenew(*count, tabnr + 1);
531 for (j = *ncount; j < tabnr + 1; j++)
537 // Record that this table index is used and must have a valid file
545 /*!\brief If there's bonded interactions of flavour \c tabext and type
546 * \c ftype1 or \c ftype2 present in the topology, seek them in the
547 * list of filenames passed to mdrun, and make bonded tables from
550 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
551 * valid type with that parameter.
553 * A fatal error occurs if no matching filename is found.
555 static std::vector<bondedtable_t> make_bonded_tables(FILE* fplog,
558 const gmx_mtop_t* mtop,
559 gmx::ArrayRef<const std::string> tabbfnm,
562 std::vector<bondedtable_t> tab;
565 int* count = nullptr;
566 count_tables(ftype1, ftype2, mtop, &ncount, &count);
568 // Are there any relevant tabulated bond interactions?
572 for (int i = 0; i < ncount; i++)
574 // Do any interactions exist that requires this table?
577 // This pattern enforces the current requirement that
578 // table filenames end in a characteristic sequence
579 // before the file type extension, and avoids table 13
580 // being recognized and used for table 1.
581 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
582 bool madeTable = false;
583 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
585 if (gmx::endsWith(tabbfnm[j], patternToFind))
587 // Finally read the table from the file found
588 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
594 bool isPlural = (ftype2 != -1);
596 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
597 "because no table file whose name matched '%s' was passed via the "
598 "gmx mdrun -tableb command-line option.",
599 interaction_function[ftype1].longname, isPlural ? "' or '" : "",
600 isPlural ? interaction_function[ftype2].longname : "", i,
601 patternToFind.c_str());
611 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
613 fr->natoms_force = natoms_force;
614 fr->natoms_force_constr = natoms_force_constr;
616 fr->forceHelperBuffers->resize(natoms_f_novirsum);
619 static real cutoff_inf(real cutoff)
623 cutoff = GMX_CUTOFF_INF;
629 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
630 static void initCoulombEwaldParameters(FILE* fp,
631 const t_inputrec* ir,
632 bool systemHasNetCharge,
633 interaction_const_t* ic)
635 if (!EEL_PME_EWALD(ir->coulombtype))
642 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
644 if (ir->coulombtype == eelP3M_AD)
646 please_cite(fp, "Hockney1988");
647 please_cite(fp, "Ballenegger2012");
651 please_cite(fp, "Essmann95a");
654 if (ir->ewald_geometry == eewg3DC)
658 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
659 systemHasNetCharge ? " and net charge" : "");
661 please_cite(fp, "In-Chul99a");
662 if (systemHasNetCharge)
664 please_cite(fp, "Ballenegger2009");
669 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
672 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
675 if (ic->coulomb_modifier == eintmodPOTSHIFT)
677 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
678 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
686 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
687 static void initVdwEwaldParameters(FILE* fp, const t_inputrec* ir, interaction_const_t* ic)
689 if (!EVDW_PME(ir->vdwtype))
696 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
697 please_cite(fp, "Essmann95a");
699 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
702 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
705 if (ic->vdw_modifier == eintmodPOTSHIFT)
707 real crc2 = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
708 ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
716 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
718 * Tables are generated for one or both, depending on if the pointers
719 * are non-null. The spacing for both table sets is the same and obeys
720 * both accuracy requirements, when relevant.
722 static void init_ewald_f_table(const interaction_const_t& ic,
723 const real tableExtensionLength,
724 EwaldCorrectionTables* coulombTables,
725 EwaldCorrectionTables* vdwTables)
727 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
728 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
730 /* Get the Ewald table spacing based on Coulomb and/or LJ
731 * Ewald coefficients and rtol.
733 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
735 const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
737 real tableLen = ic.rcoulomb;
738 if (useCoulombTable && havePerturbedNonbondeds && tableExtensionLength > 0.0)
740 /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
741 * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
742 * The alternative is to look through all the exclusions and check if they come from
743 * couple-intramol == no. Meanwhile, always having larger tables should only affect
744 * memory consumption, not speed (barring cache issues).
746 tableLen = ic.rcoulomb + tableExtensionLength;
748 const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
753 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
758 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
762 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real tableExtensionLength)
764 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
766 init_ewald_f_table(*ic, tableExtensionLength, ic->coulombEwaldTables.get(),
767 ic->vdwEwaldTables.get());
770 fprintf(fp, "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
771 1 / ic->coulombEwaldTables->scale, ic->coulombEwaldTables->tableF.size());
776 static void clear_force_switch_constants(shift_consts_t* sc)
783 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
785 /* Here we determine the coefficient for shifting the force to zero
786 * between distance rsw and the cut-off rc.
787 * For a potential of r^-p, we have force p*r^-(p+1).
788 * But to save flops we absorb p in the coefficient.
790 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
791 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
793 sc->c2 = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
794 sc->c3 = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
795 sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
796 + p * sc->c3 / 4 * gmx::power4(rc - rsw);
799 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
801 /* The switch function is 1 at rsw and 0 at rc.
802 * The derivative and second derivate are zero at both ends.
803 * rsw = max(r - r_switch, 0)
804 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
805 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
806 * force = force*dsw - potential*sw
809 sc->c3 = -10 / gmx::power3(rc - rsw);
810 sc->c4 = 15 / gmx::power4(rc - rsw);
811 sc->c5 = -6 / gmx::power5(rc - rsw);
814 /*! \brief Construct interaction constants
816 * This data is used (particularly) by search and force code for
817 * short-range interactions. Many of these are constant for the whole
818 * simulation; some are constant only after PME tuning completes.
820 static void init_interaction_const(FILE* fp,
821 interaction_const_t** interaction_const,
822 const t_inputrec* ir,
823 const gmx_mtop_t* mtop,
824 bool systemHasNetCharge)
826 interaction_const_t* ic = new interaction_const_t;
828 ic->cutoff_scheme = ir->cutoff_scheme;
830 ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
831 ic->vdwEwaldTables = std::make_unique<EwaldCorrectionTables>();
834 ic->vdwtype = ir->vdwtype;
835 ic->vdw_modifier = ir->vdw_modifier;
836 ic->reppow = mtop->ffparams.reppow;
837 ic->rvdw = cutoff_inf(ir->rvdw);
838 ic->rvdw_switch = ir->rvdw_switch;
839 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
840 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
841 if (ic->useBuckingham)
843 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
846 initVdwEwaldParameters(fp, ir, ic);
848 clear_force_switch_constants(&ic->dispersion_shift);
849 clear_force_switch_constants(&ic->repulsion_shift);
851 switch (ic->vdw_modifier)
853 case eintmodPOTSHIFT:
854 /* Only shift the potential, don't touch the force */
855 ic->dispersion_shift.cpot = -1.0 / gmx::power6(ic->rvdw);
856 ic->repulsion_shift.cpot = -1.0 / gmx::power12(ic->rvdw);
858 case eintmodFORCESWITCH:
859 /* Switch the force, switch and shift the potential */
860 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw, &ic->dispersion_shift);
861 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw, &ic->repulsion_shift);
863 case eintmodPOTSWITCH:
864 /* Switch the potential and force */
865 potential_switch_constants(ic->rvdw_switch, ic->rvdw, &ic->vdw_switch);
868 case eintmodEXACTCUTOFF:
869 /* Nothing to do here */
871 default: gmx_incons("unimplemented potential modifier");
875 ic->eeltype = ir->coulombtype;
876 ic->coulomb_modifier = ir->coulomb_modifier;
877 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
878 ic->rcoulomb_switch = ir->rcoulomb_switch;
879 ic->epsilon_r = ir->epsilon_r;
881 /* Set the Coulomb energy conversion factor */
882 if (ic->epsilon_r != 0)
884 ic->epsfac = ONE_4PI_EPS0 / ic->epsilon_r;
888 /* eps = 0 is infinite dieletric: no Coulomb interactions */
893 if (EEL_RF(ic->eeltype))
895 GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
896 ic->epsilon_rf = ir->epsilon_rf;
898 calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf, ic->rcoulomb, &ic->k_rf, &ic->c_rf);
902 /* For plain cut-off we might use the reaction-field kernels */
903 ic->epsilon_rf = ic->epsilon_r;
905 if (ir->coulomb_modifier == eintmodPOTSHIFT)
907 ic->c_rf = 1 / ic->rcoulomb;
915 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
919 real dispersion_shift;
921 dispersion_shift = ic->dispersion_shift.cpot;
922 if (EVDW_PME(ic->vdwtype))
924 dispersion_shift -= ic->sh_lj_ewald;
926 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e", ic->repulsion_shift.cpot, dispersion_shift);
928 if (ic->eeltype == eelCUT)
930 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
932 else if (EEL_PME(ic->eeltype))
934 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
939 if (ir->efep != efepNO)
941 GMX_RELEASE_ASSERT(ir->fepvals, "ir->fepvals should be set wth free-energy");
942 ic->softCoreParameters = std::make_unique<interaction_const_t::SoftCoreParameters>(*ir->fepvals);
945 *interaction_const = ic;
948 void init_forcerec(FILE* fp,
949 const gmx::MDLogger& mdlog,
951 const t_inputrec* ir,
952 const gmx_mtop_t* mtop,
957 gmx::ArrayRef<const std::string> tabbfnm,
960 /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
961 fr->use_simd_kernels = GMX_USE_SIMD_KERNELS;
963 if (check_box(ir->pbcType, box))
965 gmx_fatal(FARGS, "%s", check_box(ir->pbcType, box));
968 /* Test particle insertion ? */
971 /* Set to the size of the molecule to be inserted (the last one) */
972 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
973 fr->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
980 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED || ir->coulombtype == eelGRF_NOTUSED)
982 gmx_fatal(FARGS, "%s electrostatics is no longer supported", eel_names[ir->coulombtype]);
987 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
989 if (ir->useTwinRange)
991 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
993 /* Copy the user determined parameters */
994 fr->userint1 = ir->userint1;
995 fr->userint2 = ir->userint2;
996 fr->userint3 = ir->userint3;
997 fr->userint4 = ir->userint4;
998 fr->userreal1 = ir->userreal1;
999 fr->userreal2 = ir->userreal2;
1000 fr->userreal3 = ir->userreal3;
1001 fr->userreal4 = ir->userreal4;
1004 fr->fc_stepsize = ir->fc_stepsize;
1007 fr->efep = ir->efep;
1009 fr->bNonbonded = TRUE;
1010 if (getenv("GMX_NO_NONBONDED") != nullptr)
1012 /* turn off non-bonded calculations */
1013 fr->bNonbonded = FALSE;
1014 GMX_LOG(mdlog.warning)
1017 "Found environment variable GMX_NO_NONBONDED.\n"
1018 "Disabling nonbonded calculations.");
1021 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1023 fr->use_simd_kernels = FALSE;
1027 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1028 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1029 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1033 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1035 /* Neighbour searching stuff */
1036 fr->cutoff_scheme = ir->cutoff_scheme;
1037 fr->pbcType = ir->pbcType;
1039 /* Determine if we will do PBC for distances in bonded interactions */
1040 if (fr->pbcType == PbcType::No)
1042 fr->bMolPBC = FALSE;
1046 const bool useEwaldSurfaceCorrection =
1047 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1048 const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
1049 if (!DOMAINDECOMP(cr))
1053 if (useEwaldSurfaceCorrection || haveOrientationRestraints)
1055 fr->wholeMoleculeTransform =
1056 std::make_unique<gmx::WholeMoleculeTransform>(*mtop, ir->pbcType);
1061 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->pbcType);
1063 /* With Ewald surface correction it is useful to support e.g. running water
1064 * in parallel with update groups.
1065 * With orientation restraints there is no sensible use case supported with DD.
1067 if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd)) || haveOrientationRestraints)
1070 "You requested Ewald surface correction or orientation restraints, "
1071 "but molecules are broken "
1072 "over periodic boundary conditions by the domain decomposition. "
1073 "Run without domain decomposition instead.");
1077 if (useEwaldSurfaceCorrection)
1079 GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || dd_moleculesAreAlwaysWhole(*cr->dd),
1080 "Molecules can not be broken by PBC with epsilon_surface > 0");
1084 fr->rc_scaling = ir->refcoord_scaling;
1085 copy_rvec(ir->posres_com, fr->posres_com);
1086 copy_rvec(ir->posres_comB, fr->posres_comB);
1087 fr->rlist = cutoff_inf(ir->rlist);
1088 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1090 /* This now calculates sum for q and c6*/
1091 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1093 /* fr->ic is used both by verlet and group kernels (to some extent) now */
1094 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1095 init_interaction_const_tables(fp, fr->ic, ir->tabext);
1097 const interaction_const_t* ic = fr->ic;
1099 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1100 if (ir->coulombtype == eelEWALD)
1102 init_ewald_tab(&(fr->ewald_table), ir, fp);
1105 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1106 switch (ic->eeltype)
1108 case eelCUT: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB; break;
1111 case eelRF_ZERO: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD; break;
1118 case eelPMEUSERSWITCH:
1119 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1124 case eelEWALD: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD; break;
1127 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1129 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1131 /* Vdw: Translate from mdp settings to kernel format */
1132 switch (ic->vdwtype)
1137 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1141 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1144 case evdwPME: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD; break;
1148 case evdwUSER: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE; break;
1150 default: gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1152 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1154 if (ir->cutoff_scheme == ecutsVERLET)
1156 if (!gmx_within_tol(ic->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1158 gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12",
1159 ecutscheme_names[ir->cutoff_scheme]);
1161 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1162 * while mdrun does not (and never did) support this.
1164 if (EEL_USER(fr->ic->eeltype))
1166 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
1167 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
1170 fr->bvdwtab = FALSE;
1171 fr->bcoultab = FALSE;
1174 /* 1-4 interaction electrostatics */
1175 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1177 const bool haveDirectVirialContributions =
1178 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) || fr->forceProviders->hasForceProvider()
1179 || gmx_mtop_ftype_count(mtop, F_POSRES) > 0 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0
1180 || ir->nwall > 0 || ir->bPull || ir->bRot || ir->bIMD);
1181 fr->forceHelperBuffers = std::make_unique<ForceHelperBuffers>(haveDirectVirialContributions);
1183 if (fr->shift_vec == nullptr)
1185 snew(fr->shift_vec, SHIFTS);
1188 if (fr->nbfp.empty())
1190 fr->ntype = mtop->ffparams.atnr;
1191 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1192 if (EVDW_PME(ic->vdwtype))
1194 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1198 /* Copy the energy group exclusions */
1199 fr->egp_flags = ir->opts.egp_flags;
1201 /* Van der Waals stuff */
1202 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1204 if (ic->rvdw_switch >= ic->rvdw)
1206 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)", ic->rvdw_switch, ic->rvdw);
1210 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
1211 (ic->eeltype == eelSWITCH) ? "switched" : "shifted", ic->rvdw_switch, ic->rvdw);
1215 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1217 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1220 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1222 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1225 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
1227 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
1230 if (fp && fr->cutoff_scheme == ecutsGROUP)
1232 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n", fr->rlist, ic->rcoulomb,
1233 fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
1236 if (ir->implicit_solvent)
1238 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1242 /* This code automatically gives table length tabext without cut-off's,
1243 * in that case grompp should already have checked that we do not need
1244 * normal tables and we only generate tables for 1-4 interactions.
1246 real rtab = ir->rlist + ir->tabext;
1248 /* We want to use unmodified tables for 1-4 coulombic
1249 * interactions, so we must in general have an extra set of
1251 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1252 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1254 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1258 fr->nwall = ir->nwall;
1259 if (ir->nwall && ir->wall_type == ewtTABLE)
1261 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1264 fr->fcdata = std::make_unique<t_fcdata>();
1266 if (!tabbfnm.empty())
1268 t_fcdata& fcdata = *fr->fcdata;
1269 // Need to catch std::bad_alloc
1270 // TODO Don't need to catch this here, when merging with master branch
1273 // TODO move these tables into a separate struct and store reference in ListedForces
1274 fcdata.bondtab = make_bonded_tables(fp, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1275 fcdata.angletab = make_bonded_tables(fp, F_TABANGLES, -1, mtop, tabbfnm, "a");
1276 fcdata.dihtab = make_bonded_tables(fp, F_TABDIHS, -1, mtop, tabbfnm, "d");
1278 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1285 "No fcdata or table file name passed, can not read table, can not do bonded "
1290 /* Initialize the thread working data for bonded interactions */
1291 fr->listedForces.emplace_back(
1292 mtop->ffparams, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1293 gmx_omp_nthreads_get(emntBonded), ListedForces::interactionSelectionAll(), fp);
1295 // QM/MM initialization if requested
1298 gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1301 /* Set all the static charge group info */
1302 fr->cginfo_mb = init_cginfo_mb(mtop, fr);
1303 if (!DOMAINDECOMP(cr))
1305 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1308 if (!DOMAINDECOMP(cr))
1310 forcerec_set_ranges(fr, mtop->natoms, mtop->natoms, mtop->natoms);
1313 fr->print_force = print_force;
1315 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1316 snew(fr->ewc_t, fr->nthread_ewc);
1318 if (ir->eDispCorr != edispcNO)
1320 fr->dispersionCorrection = std::make_unique<DispersionCorrection>(
1321 *mtop, *ir, fr->bBHAM, fr->ntype, fr->nbfp, *fr->ic, tabfn);
1322 fr->dispersionCorrection->print(mdlog);
1327 /* Here we switch from using mdlog, which prints the newline before
1328 * the paragraph, to our old fprintf logging, which prints the newline
1329 * after the paragraph, so we should add a newline here.
1335 t_forcerec::t_forcerec() = default;
1337 t_forcerec::~t_forcerec()
1339 /* Note: This code will disappear when types are converted to C++ */