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.
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/ewald/ewald.h"
54 #include "gromacs/ewald/pme_pp_comm_gpu.h"
55 #include "gromacs/fileio/filetypes.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
58 #include "gromacs/gpu_utils/gpu_utils.h"
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/listed_forces/listed_forces_gpu.h"
61 #include "gromacs/listed_forces/listed_forces.h"
62 #include "gromacs/listed_forces/pairs.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/utilities.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/mdlib/dispersioncorrection.h"
67 #include "gromacs/mdlib/force.h"
68 #include "gromacs/mdlib/forcerec_threading.h"
69 #include "gromacs/mdlib/gmx_omp_nthreads.h"
70 #include "gromacs/mdlib/md_support.h"
71 #include "gromacs/mdlib/wall.h"
72 #include "gromacs/mdlib/wholemoleculetransform.h"
73 #include "gromacs/mdtypes/commrec.h"
74 #include "gromacs/mdtypes/fcdata.h"
75 #include "gromacs/mdtypes/forcerec.h"
76 #include "gromacs/mdtypes/group.h"
77 #include "gromacs/mdtypes/iforceprovider.h"
78 #include "gromacs/mdtypes/inputrec.h"
79 #include "gromacs/mdtypes/interaction_const.h"
80 #include "gromacs/mdtypes/md_enums.h"
81 #include "gromacs/mdtypes/multipletimestepping.h"
82 #include "gromacs/mdtypes/nblist.h"
83 #include "gromacs/mdtypes/simulation_workload.h"
84 #include "gromacs/nbnxm/nbnxm.h"
85 #include "gromacs/pbcutil/ishift.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/tables/forcetable.h"
88 #include "gromacs/topology/mtop_util.h"
89 #include "gromacs/topology/idef.h"
90 #include "gromacs/trajectory/trajectoryframe.h"
91 #include "gromacs/utility/cstringutil.h"
92 #include "gromacs/utility/exceptions.h"
93 #include "gromacs/utility/fatalerror.h"
94 #include "gromacs/utility/gmxassert.h"
95 #include "gromacs/utility/logger.h"
96 #include "gromacs/utility/physicalnodecommunicator.h"
97 #include "gromacs/utility/pleasecite.h"
98 #include "gromacs/utility/smalloc.h"
99 #include "gromacs/utility/strconvert.h"
101 #include "gpuforcereduction.h"
103 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
104 haveDirectVirialContributions_(haveDirectVirialContributions)
106 shiftForces_.resize(gmx::c_numShiftVectors);
109 void ForceHelperBuffers::resize(int numAtoms)
111 if (haveDirectVirialContributions_)
113 forceBufferForDirectVirialContributions_.resize(numAtoms);
117 std::vector<real> makeNonBondedParameterLists(const int numAtomTypes,
118 gmx::ArrayRef<const t_iparams> iparams,
119 bool useBuckinghamPotential)
121 std::vector<real> nbfp;
123 if (useBuckinghamPotential)
125 nbfp.resize(3 * numAtomTypes * numAtomTypes);
127 for (int i = 0; (i < numAtomTypes); i++)
129 for (int j = 0; (j < numAtomTypes); j++, k++)
131 BHAMA(nbfp, numAtomTypes, i, j) = iparams[k].bham.a;
132 BHAMB(nbfp, numAtomTypes, i, j) = iparams[k].bham.b;
133 /* nbfp now includes the 6.0 derivative prefactor */
134 BHAMC(nbfp, numAtomTypes, i, j) = iparams[k].bham.c * 6.0;
140 nbfp.resize(2 * numAtomTypes * numAtomTypes);
142 for (int i = 0; (i < numAtomTypes); i++)
144 for (int j = 0; (j < numAtomTypes); j++, k++)
146 /* nbfp now includes the 6.0/12.0 derivative prefactors */
147 C6(nbfp, numAtomTypes, i, j) = iparams[k].lj.c6 * 6.0;
148 C12(nbfp, numAtomTypes, i, j) = iparams[k].lj.c12 * 12.0;
156 std::vector<real> makeLJPmeC6GridCorrectionParameters(const int numAtomTypes,
157 gmx::ArrayRef<const t_iparams> iparams,
158 LongRangeVdW ljpme_combination_rule)
160 /* For LJ-PME simulations, we correct the energies with the reciprocal space
161 * inside of the cut-off. To do this the non-bonded kernels needs to have
162 * access to the C6-values used on the reciprocal grid in pme.c
165 std::vector<real> grid(2 * numAtomTypes * numAtomTypes, 0.0);
167 for (int i = 0; (i < numAtomTypes); i++)
169 for (int j = 0; (j < numAtomTypes); j++, k++)
171 real c6i = iparams[i * (numAtomTypes + 1)].lj.c6;
172 real c12i = iparams[i * (numAtomTypes + 1)].lj.c12;
173 real c6j = iparams[j * (numAtomTypes + 1)].lj.c6;
174 real c12j = iparams[j * (numAtomTypes + 1)].lj.c12;
175 real c6 = std::sqrt(c6i * c6j);
176 if (ljpme_combination_rule == LongRangeVdW::LB && !gmx_numzero(c6) && !gmx_numzero(c12i)
177 && !gmx_numzero(c12j))
179 real sigmai = gmx::sixthroot(c12i / c6i);
180 real sigmaj = gmx::sixthroot(c12j / c6j);
181 real epsi = c6i * c6i / c12i;
182 real epsj = c6j * c6j / c12j;
183 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
185 /* Store the elements at the same relative positions as C6 in nbfp in order
186 * to simplify access in the kernels
188 grid[2 * (numAtomTypes * i + j)] = c6 * 6.0;
194 //! What kind of constraint affects an atom
195 enum class ConstraintTypeForAtom : int
197 None, //!< No constraint active
198 Constraint, //!< F_CONSTR or F_CONSTRNC active
199 Settle, //! F_SETTLE active
202 static std::vector<gmx::AtomInfoWithinMoleculeBlock>
203 makeAtomInfoForEachMoleculeBlock(const gmx_mtop_t& mtop, const t_forcerec* fr)
205 std::vector<bool> atomUsesVdw(fr->ntype, false);
206 for (int ai = 0; ai < fr->ntype; ai++)
208 for (int j = 0; j < fr->ntype; j++)
210 atomUsesVdw[ai] = atomUsesVdw[ai] || fr->haveBuckingham || C6(fr->nbfp, fr->ntype, ai, j) != 0
211 || C12(fr->nbfp, fr->ntype, ai, j) != 0;
215 std::vector<gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock;
216 int indexOfFirstAtomInMoleculeBlock = 0;
217 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
219 const gmx_molblock_t& molb = mtop.molblock[mb];
220 const gmx_moltype_t& molt = mtop.moltype[molb.type];
221 const auto& excl = molt.excls;
223 /* Check if all molecules in this block have identical
224 * atominfo. (That's true unless some kind of group- or
225 * distance-based algorithm is involved, e.g. QM/MM groups
226 * affecting multiple molecules within a block differently.)
227 * If so, we only need an array of the size of one molecule.
228 * Otherwise we make an array of #mol times #atoms per
231 bool allMoleculesWithinBlockAreIdentical = true;
232 for (int m = 0; m < molb.nmol; m++)
234 const int numAtomsInAllMolecules = m * molt.atoms.nr;
235 for (int a = 0; a < molt.atoms.nr; a++)
237 if (getGroupType(mtop.groups,
238 SimulationAtomGroupType::QuantumMechanics,
239 indexOfFirstAtomInMoleculeBlock + numAtomsInAllMolecules + a)
240 != getGroupType(mtop.groups,
241 SimulationAtomGroupType::QuantumMechanics,
242 indexOfFirstAtomInMoleculeBlock + a))
244 allMoleculesWithinBlockAreIdentical = false;
246 if (!mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
248 if (mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics]
249 [indexOfFirstAtomInMoleculeBlock + numAtomsInAllMolecules + a]
250 != mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics]
251 [indexOfFirstAtomInMoleculeBlock + a])
253 allMoleculesWithinBlockAreIdentical = false;
259 gmx::AtomInfoWithinMoleculeBlock atomInfoOfMoleculeBlock;
260 atomInfoOfMoleculeBlock.indexOfFirstAtomInMoleculeBlock = indexOfFirstAtomInMoleculeBlock;
261 atomInfoOfMoleculeBlock.indexOfLastAtomInMoleculeBlock =
262 indexOfFirstAtomInMoleculeBlock + molb.nmol * molt.atoms.nr;
263 int atomInfoSize = (allMoleculesWithinBlockAreIdentical ? 1 : molb.nmol) * molt.atoms.nr;
264 atomInfoOfMoleculeBlock.atomInfo.resize(atomInfoSize);
266 /* Set constraints flags for constrained atoms */
267 std::vector<ConstraintTypeForAtom> constraintTypeOfAtom(molt.atoms.nr, ConstraintTypeForAtom::None);
268 for (int ftype = 0; ftype < F_NRE; ftype++)
270 if (interaction_function[ftype].flags & IF_CONSTRAINT)
272 const int nral = NRAL(ftype);
273 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
275 for (int a = 0; a < nral; a++)
277 constraintTypeOfAtom[molt.ilist[ftype].iatoms[ia + 1 + a]] =
278 (ftype == F_SETTLE ? ConstraintTypeForAtom::Settle
279 : ConstraintTypeForAtom::Constraint);
285 for (int m = 0; m < (allMoleculesWithinBlockAreIdentical ? 1 : molb.nmol); m++)
287 const int moleculeOffsetInBlock = m * molt.atoms.nr;
288 for (int a = 0; a < molt.atoms.nr; a++)
290 const t_atom& atom = molt.atoms.atom[a];
291 int64_t& atomInfo = atomInfoOfMoleculeBlock.atomInfo[moleculeOffsetInBlock + a];
293 /* Store the energy group in atomInfo */
294 int gid = getGroupType(mtop.groups,
295 SimulationAtomGroupType::EnergyOutput,
296 indexOfFirstAtomInMoleculeBlock + moleculeOffsetInBlock + a);
297 atomInfo = (atomInfo & ~gmx::sc_atomInfo_EnergyGroupIdMask) | gid;
299 bool bHaveVDW = (atomUsesVdw[atom.type] || atomUsesVdw[atom.typeB]);
300 bool bHaveQ = (atom.q != 0 || atom.qB != 0);
302 bool haveExclusions = false;
303 /* Loop over all the exclusions of atom ai */
304 for (const int j : excl[a])
308 haveExclusions = true;
313 switch (constraintTypeOfAtom[a])
315 case ConstraintTypeForAtom::Constraint:
316 atomInfo |= gmx::sc_atomInfo_Constraint;
318 case ConstraintTypeForAtom::Settle: atomInfo |= gmx::sc_atomInfo_Settle; break;
323 atomInfo |= gmx::sc_atomInfo_Exclusion;
327 atomInfo |= gmx::sc_atomInfo_HasVdw;
331 atomInfo |= gmx::sc_atomInfo_HasCharge;
333 if (fr->efep != FreeEnergyPerturbationType::No && PERTURBED(atom))
335 atomInfo |= gmx::sc_atomInfo_FreeEnergyPerturbation;
340 atomInfoForEachMoleculeBlock.push_back(atomInfoOfMoleculeBlock);
342 indexOfFirstAtomInMoleculeBlock += molb.nmol * molt.atoms.nr;
345 return atomInfoForEachMoleculeBlock;
348 static std::vector<int64_t> expandAtomInfo(const int nmb,
349 gmx::ArrayRef<const gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock)
351 const int numAtoms = atomInfoForEachMoleculeBlock[nmb - 1].indexOfLastAtomInMoleculeBlock;
353 std::vector<int64_t> atomInfo(numAtoms);
356 for (int a = 0; a < numAtoms; a++)
358 while (a >= atomInfoForEachMoleculeBlock[mb].indexOfLastAtomInMoleculeBlock)
362 atomInfo[a] = atomInfoForEachMoleculeBlock[mb]
363 .atomInfo[(a - atomInfoForEachMoleculeBlock[mb].indexOfFirstAtomInMoleculeBlock)
364 % atomInfoForEachMoleculeBlock[mb].atomInfo.size()];
370 /* Sets the sum of charges (squared) and C6 in the system in fr.
371 * Returns whether the system has a net charge.
373 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t& mtop)
375 /*This now calculates sum for q and c6*/
376 double qsum, q2sum, q, c6sum, c6;
381 for (const gmx_molblock_t& molb : mtop.molblock)
383 int nmol = molb.nmol;
384 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
385 for (int i = 0; i < atoms->nr; i++)
387 q = atoms->atom[i].q;
389 q2sum += nmol * q * q;
390 c6 = mtop.ffparams.iparams[atoms->atom[i].type * (mtop.ffparams.atnr + 1)].lj.c6;
395 fr->q2sum[0] = q2sum;
396 fr->c6sum[0] = c6sum;
398 if (fr->efep != FreeEnergyPerturbationType::No)
403 for (const gmx_molblock_t& molb : mtop.molblock)
405 int nmol = molb.nmol;
406 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
407 for (int i = 0; i < atoms->nr; i++)
409 q = atoms->atom[i].qB;
411 q2sum += nmol * q * q;
412 c6 = mtop.ffparams.iparams[atoms->atom[i].typeB * (mtop.ffparams.atnr + 1)].lj.c6;
416 fr->q2sum[1] = q2sum;
417 fr->c6sum[1] = c6sum;
422 fr->qsum[1] = fr->qsum[0];
423 fr->q2sum[1] = fr->q2sum[0];
424 fr->c6sum[1] = fr->c6sum[0];
428 if (fr->efep == FreeEnergyPerturbationType::No)
430 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
434 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
438 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
439 return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
442 /*!\brief If there's bonded interactions of type \c ftype1 or \c
443 * ftype2 present in the topology, build an array of the number of
444 * interactions present for each bonded interaction index found in the
447 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
448 * valid type with that parameter.
450 * \c count will be reallocated as necessary to fit the largest bonded
451 * interaction index found, and its current size will be returned in
452 * \c ncount. It will contain zero for every bonded interaction index
453 * for which no interactions are present in the topology.
455 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t& mtop, int* ncount, int** count)
457 int ftype, i, j, tabnr;
459 // Loop over all moleculetypes
460 for (const gmx_moltype_t& molt : mtop.moltype)
462 // Loop over all interaction types
463 for (ftype = 0; ftype < F_NRE; ftype++)
465 // If the current interaction type is one of the types whose tables we're trying to count...
466 if (ftype == ftype1 || ftype == ftype2)
468 const InteractionList& il = molt.ilist[ftype];
469 const int stride = 1 + NRAL(ftype);
470 // ... and there are actually some interactions for this type
471 for (i = 0; i < il.size(); i += stride)
473 // Find out which table index the user wanted
474 tabnr = mtop.ffparams.iparams[il.iatoms[i]].tab.table;
477 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
479 // Make room for this index in the data structure
480 if (tabnr >= *ncount)
482 srenew(*count, tabnr + 1);
483 for (j = *ncount; j < tabnr + 1; j++)
489 // Record that this table index is used and must have a valid file
497 /*!\brief If there's bonded interactions of flavour \c tabext and type
498 * \c ftype1 or \c ftype2 present in the topology, seek them in the
499 * list of filenames passed to mdrun, and make bonded tables from
502 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
503 * valid type with that parameter.
505 * A fatal error occurs if no matching filename is found.
507 static std::vector<bondedtable_t> make_bonded_tables(FILE* fplog,
510 const gmx_mtop_t& mtop,
511 gmx::ArrayRef<const std::string> tabbfnm,
514 std::vector<bondedtable_t> tab;
517 int* count = nullptr;
518 count_tables(ftype1, ftype2, mtop, &ncount, &count);
520 // Are there any relevant tabulated bond interactions?
524 for (int i = 0; i < ncount; i++)
526 // Do any interactions exist that requires this table?
529 // This pattern enforces the current requirement that
530 // table filenames end in a characteristic sequence
531 // before the file type extension, and avoids table 13
532 // being recognized and used for table 1.
533 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
534 bool madeTable = false;
535 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
537 if (gmx::endsWith(tabbfnm[j], patternToFind))
539 // Finally read the table from the file found
540 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
546 bool isPlural = (ftype2 != -1);
548 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
549 "because no table file whose name matched '%s' was passed via the "
550 "gmx mdrun -tableb command-line option.",
551 interaction_function[ftype1].longname,
552 isPlural ? "' or '" : "",
553 isPlural ? interaction_function[ftype2].longname : "",
555 patternToFind.c_str());
565 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
567 fr->natoms_force = natoms_force;
568 fr->natoms_force_constr = natoms_force_constr;
570 for (auto& forceHelperBuffers : fr->forceHelperBuffers)
572 forceHelperBuffers.resize(natoms_f_novirsum);
576 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
578 * Tables are generated for one or both, depending on if the pointers
579 * are non-null. The spacing for both table sets is the same and obeys
580 * both accuracy requirements, when relevant.
582 static void init_ewald_f_table(const interaction_const_t& ic,
585 EwaldCorrectionTables* coulombTables,
586 EwaldCorrectionTables* vdwTables)
588 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
589 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
591 /* Get the Ewald table spacing based on Coulomb and/or LJ
592 * Ewald coefficients and rtol.
594 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
596 const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
598 real tableLen = ic.rcoulomb;
599 if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
601 /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
602 * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
603 * The alternative is to look through all the exclusions and check if they come from
604 * couple-intramol == no. Meanwhile, always having larger tables should only affect
605 * memory consumption, not speed (barring cache issues).
607 tableLen = rlist + tabext;
609 const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
614 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
619 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
623 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
625 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
628 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
632 "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
633 1 / ic->coulombEwaldTables->scale,
634 ic->coulombEwaldTables->tableF.size());
639 real cutoff_inf(real cutoff)
643 cutoff = GMX_CUTOFF_INF;
649 void init_forcerec(FILE* fplog,
650 const gmx::MDLogger& mdlog,
651 const gmx::SimulationWorkload& simulationWork,
652 t_forcerec* forcerec,
653 const t_inputrec& inputrec,
654 const gmx_mtop_t& mtop,
655 const t_commrec* commrec,
659 gmx::ArrayRef<const std::string> tabbfnm,
662 /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
663 forcerec->use_simd_kernels = GMX_USE_SIMD_KERNELS;
665 if (check_box(inputrec.pbcType, box))
667 gmx_fatal(FARGS, "%s", check_box(inputrec.pbcType, box));
670 /* Test particle insertion ? */
671 if (EI_TPI(inputrec.eI))
673 /* Set to the size of the molecule to be inserted (the last one) */
674 gmx::RangePartitioning molecules = gmx_mtop_molecules(mtop);
675 forcerec->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
682 if (inputrec.coulombtype == CoulombInteractionType::RFNecUnsupported
683 || inputrec.coulombtype == CoulombInteractionType::GRFNotused)
685 gmx_fatal(FARGS, "%s electrostatics is no longer supported", enumValueToString(inputrec.coulombtype));
688 if (inputrec.bAdress)
690 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
692 if (inputrec.useTwinRange)
694 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
696 /* Copy the user determined parameters */
697 forcerec->userint1 = inputrec.userint1;
698 forcerec->userint2 = inputrec.userint2;
699 forcerec->userint3 = inputrec.userint3;
700 forcerec->userint4 = inputrec.userint4;
701 forcerec->userreal1 = inputrec.userreal1;
702 forcerec->userreal2 = inputrec.userreal2;
703 forcerec->userreal3 = inputrec.userreal3;
704 forcerec->userreal4 = inputrec.userreal4;
707 forcerec->fc_stepsize = inputrec.fc_stepsize;
710 forcerec->efep = inputrec.efep;
712 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
714 forcerec->use_simd_kernels = FALSE;
715 if (fplog != nullptr)
718 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
719 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
720 "(e.g. SSE2/SSE4.1/AVX).\n\n");
724 forcerec->haveBuckingham = (mtop.ffparams.functype[0] == F_BHAM);
726 /* Neighbour searching stuff */
727 forcerec->pbcType = inputrec.pbcType;
729 /* Determine if we will do PBC for distances in bonded interactions */
730 if (forcerec->pbcType == PbcType::No)
732 forcerec->bMolPBC = FALSE;
736 const bool useEwaldSurfaceCorrection =
737 (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0);
738 const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
739 if (!DOMAINDECOMP(commrec))
741 forcerec->bMolPBC = true;
743 if (useEwaldSurfaceCorrection || haveOrientationRestraints)
745 forcerec->wholeMoleculeTransform =
746 std::make_unique<gmx::WholeMoleculeTransform>(mtop, inputrec.pbcType);
751 forcerec->bMolPBC = dd_bonded_molpbc(*commrec->dd, forcerec->pbcType);
753 /* With Ewald surface correction it is useful to support e.g. running water
754 * in parallel with update groups.
755 * With orientation restraints there is no sensible use case supported with DD.
757 if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*commrec->dd))
758 || haveOrientationRestraints)
761 "You requested Ewald surface correction or orientation restraints, "
762 "but molecules are broken "
763 "over periodic boundary conditions by the domain decomposition. "
764 "Run without domain decomposition instead.");
768 if (useEwaldSurfaceCorrection)
770 GMX_RELEASE_ASSERT(!DOMAINDECOMP(commrec) || dd_moleculesAreAlwaysWhole(*commrec->dd),
771 "Molecules can not be broken by PBC with epsilon_surface > 0");
775 forcerec->rc_scaling = inputrec.refcoord_scaling;
776 copy_rvec(inputrec.posres_com, forcerec->posres_com);
777 copy_rvec(inputrec.posres_comB, forcerec->posres_comB);
778 forcerec->rlist = cutoff_inf(inputrec.rlist);
779 forcerec->ljpme_combination_rule = inputrec.ljpme_combination_rule;
781 /* This now calculates sum for q and c6*/
782 bool systemHasNetCharge = set_chargesum(fplog, forcerec, mtop);
784 /* Make data structure used by kernels */
785 forcerec->ic = std::make_unique<interaction_const_t>(
786 init_interaction_const(fplog, inputrec, mtop, systemHasNetCharge));
787 init_interaction_const_tables(fplog, forcerec->ic.get(), forcerec->rlist, inputrec.tabext);
789 const interaction_const_t* interactionConst = forcerec->ic.get();
791 /* TODO: Replace this Ewald table or move it into interaction_const_t */
792 if (inputrec.coulombtype == CoulombInteractionType::Ewald)
794 forcerec->ewald_table = std::make_unique<gmx_ewald_tab_t>(inputrec, fplog);
797 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
798 switch (interactionConst->eeltype)
800 case CoulombInteractionType::Cut:
801 forcerec->nbkernel_elec_interaction = NbkernelElecType::Coulomb;
804 case CoulombInteractionType::RF:
805 case CoulombInteractionType::RFZero:
806 forcerec->nbkernel_elec_interaction = NbkernelElecType::ReactionField;
809 case CoulombInteractionType::Switch:
810 case CoulombInteractionType::Shift:
811 case CoulombInteractionType::User:
812 case CoulombInteractionType::PmeSwitch:
813 case CoulombInteractionType::PmeUser:
814 case CoulombInteractionType::PmeUserSwitch:
815 forcerec->nbkernel_elec_interaction = NbkernelElecType::CubicSplineTable;
818 case CoulombInteractionType::Pme:
819 case CoulombInteractionType::P3mAD:
820 case CoulombInteractionType::Ewald:
821 forcerec->nbkernel_elec_interaction = NbkernelElecType::Ewald;
826 "Unsupported electrostatic interaction: %s",
827 enumValueToString(interactionConst->eeltype));
829 forcerec->nbkernel_elec_modifier = interactionConst->coulomb_modifier;
831 /* Vdw: Translate from mdp settings to kernel format */
832 switch (interactionConst->vdwtype)
834 case VanDerWaalsType::Cut:
835 if (forcerec->haveBuckingham)
837 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::Buckingham;
841 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LennardJones;
844 case VanDerWaalsType::Pme:
845 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LJEwald;
848 case VanDerWaalsType::Switch:
849 case VanDerWaalsType::Shift:
850 case VanDerWaalsType::User:
851 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::CubicSplineTable;
855 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", enumValueToString(interactionConst->vdwtype));
857 forcerec->nbkernel_vdw_modifier = interactionConst->vdw_modifier;
859 if (!gmx_within_tol(interactionConst->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
861 gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
863 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
864 * while mdrun does not (and never did) support this.
866 if (EEL_USER(forcerec->ic->eeltype))
869 "Electrostatics type %s is currently not supported",
870 enumValueToString(inputrec.coulombtype));
873 /* 1-4 interaction electrostatics */
874 forcerec->fudgeQQ = mtop.ffparams.fudgeQQ;
876 if (simulationWork.useMts)
878 GMX_ASSERT(gmx::checkMtsRequirements(inputrec).empty(),
879 "All MTS requirements should be met here");
882 const bool haveDirectVirialContributionsFast =
883 forcerec->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
884 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || inputrec.nwall > 0 || inputrec.bPull
885 || inputrec.bRot || inputrec.bIMD;
886 const bool haveDirectVirialContributionsSlow =
887 EEL_FULL(interactionConst->eeltype) || EVDW_PME(interactionConst->vdwtype);
888 for (int i = 0; i < (simulationWork.useMts ? 2 : 1); i++)
890 bool haveDirectVirialContributions =
891 (((!simulationWork.useMts || i == 0) && haveDirectVirialContributionsFast)
892 || ((!simulationWork.useMts || i == 1) && haveDirectVirialContributionsSlow));
893 forcerec->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
896 if (forcerec->shift_vec.empty())
898 forcerec->shift_vec.resize(gmx::c_numShiftVectors);
901 GMX_ASSERT(forcerec->nbfp.empty(), "The nonbonded force parameters should not be set up yet.");
902 forcerec->ntype = mtop.ffparams.atnr;
903 forcerec->nbfp = makeNonBondedParameterLists(
904 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->haveBuckingham);
905 if (EVDW_PME(interactionConst->vdwtype))
907 forcerec->ljpme_c6grid = makeLJPmeC6GridCorrectionParameters(
908 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->ljpme_combination_rule);
911 /* Copy the energy group exclusions */
912 forcerec->egp_flags = inputrec.opts.egp_flags;
914 /* Van der Waals stuff */
915 if ((interactionConst->vdwtype != VanDerWaalsType::Cut)
916 && (interactionConst->vdwtype != VanDerWaalsType::User) && !forcerec->haveBuckingham)
918 if (interactionConst->rvdw_switch >= interactionConst->rvdw)
921 "rvdw_switch (%f) must be < rvdw (%f)",
922 interactionConst->rvdw_switch,
923 interactionConst->rvdw);
928 "Using %s Lennard-Jones, switch between %g and %g nm\n",
929 (interactionConst->eeltype == CoulombInteractionType::Switch) ? "switched" : "shifted",
930 interactionConst->rvdw_switch,
931 interactionConst->rvdw);
935 if (forcerec->haveBuckingham && EVDW_PME(interactionConst->vdwtype))
937 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
940 if (forcerec->haveBuckingham
941 && (interactionConst->vdwtype == VanDerWaalsType::Shift
942 || interactionConst->vdwtype == VanDerWaalsType::Switch))
944 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
947 if (forcerec->haveBuckingham)
949 gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
952 if (inputrec.implicit_solvent)
954 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
958 /* This code automatically gives table length tabext without cut-off's,
959 * in that case grompp should already have checked that we do not need
960 * normal tables and we only generate tables for 1-4 interactions.
962 real rtab = inputrec.rlist + inputrec.tabext;
964 /* We want to use unmodified tables for 1-4 coulombic
965 * interactions, so we must in general have an extra set of
967 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
968 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
970 forcerec->pairsTable = make_tables(fplog, interactionConst, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
974 forcerec->nwall = inputrec.nwall;
975 if (inputrec.nwall && inputrec.wall_type == WallType::Table)
977 make_wall_tables(fplog, inputrec, tabfn, &mtop.groups, forcerec);
980 forcerec->fcdata = std::make_unique<t_fcdata>();
982 if (!tabbfnm.empty())
984 t_fcdata& fcdata = *forcerec->fcdata;
985 // Need to catch std::bad_alloc
986 // TODO Don't need to catch this here, when merging with master branch
989 // TODO move these tables into a separate struct and store reference in ListedForces
990 fcdata.bondtab = make_bonded_tables(fplog, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
991 fcdata.angletab = make_bonded_tables(fplog, F_TABANGLES, -1, mtop, tabbfnm, "a");
992 fcdata.dihtab = make_bonded_tables(fplog, F_TABDIHS, -1, mtop, tabbfnm, "d");
994 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1001 "No fcdata or table file name passed, can not read table, can not do bonded "
1006 /* Initialize the thread working data for bonded interactions */
1007 if (simulationWork.useMts)
1009 // Add one ListedForces object for each MTS level
1010 bool isFirstLevel = true;
1011 for (const auto& mtsLevel : inputrec.mtsLevels)
1013 ListedForces::InteractionSelection interactionSelection;
1014 const auto& forceGroups = mtsLevel.forceGroups;
1015 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1017 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1019 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1021 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1023 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1025 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1029 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1030 isFirstLevel = false;
1032 forcerec->listedForces.emplace_back(
1034 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1035 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1036 interactionSelection,
1042 // Add one ListedForces object with all listed interactions
1043 forcerec->listedForces.emplace_back(
1045 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1046 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1047 ListedForces::interactionSelectionAll(),
1051 // QM/MM initialization if requested
1054 gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1057 /* Set all the static charge group info */
1058 forcerec->atomInfoForEachMoleculeBlock = makeAtomInfoForEachMoleculeBlock(mtop, forcerec);
1059 if (!DOMAINDECOMP(commrec))
1061 forcerec->atomInfo = expandAtomInfo(mtop.molblock.size(), forcerec->atomInfoForEachMoleculeBlock);
1064 if (!DOMAINDECOMP(commrec))
1066 forcerec_set_ranges(forcerec, mtop.natoms, mtop.natoms, mtop.natoms);
1069 forcerec->print_force = print_force;
1071 forcerec->nthread_ewc = gmx_omp_nthreads_get(ModuleMultiThread::Bonded);
1072 forcerec->ewc_t.resize(forcerec->nthread_ewc);
1074 if (inputrec.eDispCorr != DispersionCorrectionType::No)
1076 forcerec->dispersionCorrection = std::make_unique<DispersionCorrection>(
1077 mtop, inputrec, forcerec->haveBuckingham, forcerec->ntype, forcerec->nbfp, *forcerec->ic, tabfn);
1078 forcerec->dispersionCorrection->print(mdlog);
1081 if (fplog != nullptr)
1083 /* Here we switch from using mdlog, which prints the newline before
1084 * the paragraph, to our old fprintf logging, which prints the newline
1085 * after the paragraph, so we should add a newline here.
1087 fprintf(fplog, "\n");
1091 t_forcerec::t_forcerec() = default;
1093 t_forcerec::~t_forcerec() = default;