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/gmx_omp_nthreads.h"
69 #include "gromacs/mdlib/md_support.h"
70 #include "gromacs/mdlib/wall.h"
71 #include "gromacs/mdlib/wholemoleculetransform.h"
72 #include "gromacs/mdtypes/commrec.h"
73 #include "gromacs/mdtypes/fcdata.h"
74 #include "gromacs/mdtypes/forcerec.h"
75 #include "gromacs/mdtypes/group.h"
76 #include "gromacs/mdtypes/iforceprovider.h"
77 #include "gromacs/mdtypes/inputrec.h"
78 #include "gromacs/mdtypes/interaction_const.h"
79 #include "gromacs/mdtypes/md_enums.h"
80 #include "gromacs/mdtypes/multipletimestepping.h"
81 #include "gromacs/mdtypes/nblist.h"
82 #include "gromacs/mdtypes/simulation_workload.h"
83 #include "gromacs/nbnxm/nbnxm.h"
84 #include "gromacs/pbcutil/ishift.h"
85 #include "gromacs/pbcutil/pbc.h"
86 #include "gromacs/tables/forcetable.h"
87 #include "gromacs/topology/mtop_util.h"
88 #include "gromacs/topology/idef.h"
89 #include "gromacs/trajectory/trajectoryframe.h"
90 #include "gromacs/utility/cstringutil.h"
91 #include "gromacs/utility/exceptions.h"
92 #include "gromacs/utility/fatalerror.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/logger.h"
95 #include "gromacs/utility/physicalnodecommunicator.h"
96 #include "gromacs/utility/pleasecite.h"
97 #include "gromacs/utility/smalloc.h"
98 #include "gromacs/utility/strconvert.h"
100 #include "gpuforcereduction.h"
102 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
103 haveDirectVirialContributions_(haveDirectVirialContributions)
105 shiftForces_.resize(gmx::c_numShiftVectors);
108 void ForceHelperBuffers::resize(int numAtoms)
110 if (haveDirectVirialContributions_)
112 forceBufferForDirectVirialContributions_.resize(numAtoms);
116 std::vector<real> makeNonBondedParameterLists(const int numAtomTypes,
117 gmx::ArrayRef<const t_iparams> iparams,
118 bool useBuckinghamPotential)
120 std::vector<real> nbfp;
122 if (useBuckinghamPotential)
124 nbfp.resize(3 * numAtomTypes * numAtomTypes);
126 for (int i = 0; (i < numAtomTypes); i++)
128 for (int j = 0; (j < numAtomTypes); j++, k++)
130 BHAMA(nbfp, numAtomTypes, i, j) = iparams[k].bham.a;
131 BHAMB(nbfp, numAtomTypes, i, j) = iparams[k].bham.b;
132 /* nbfp now includes the 6.0 derivative prefactor */
133 BHAMC(nbfp, numAtomTypes, i, j) = iparams[k].bham.c * 6.0;
139 nbfp.resize(2 * numAtomTypes * numAtomTypes);
141 for (int i = 0; (i < numAtomTypes); i++)
143 for (int j = 0; (j < numAtomTypes); j++, k++)
145 /* nbfp now includes the 6.0/12.0 derivative prefactors */
146 C6(nbfp, numAtomTypes, i, j) = iparams[k].lj.c6 * 6.0;
147 C12(nbfp, numAtomTypes, i, j) = iparams[k].lj.c12 * 12.0;
155 std::vector<real> makeLJPmeC6GridCorrectionParameters(const int numAtomTypes,
156 gmx::ArrayRef<const t_iparams> iparams,
157 LongRangeVdW ljpme_combination_rule)
159 /* For LJ-PME simulations, we correct the energies with the reciprocal space
160 * inside of the cut-off. To do this the non-bonded kernels needs to have
161 * access to the C6-values used on the reciprocal grid in pme.c
164 std::vector<real> grid(2 * numAtomTypes * numAtomTypes, 0.0);
166 for (int i = 0; (i < numAtomTypes); i++)
168 for (int j = 0; (j < numAtomTypes); j++, k++)
170 real c6i = iparams[i * (numAtomTypes + 1)].lj.c6;
171 real c12i = iparams[i * (numAtomTypes + 1)].lj.c12;
172 real c6j = iparams[j * (numAtomTypes + 1)].lj.c6;
173 real c12j = iparams[j * (numAtomTypes + 1)].lj.c12;
174 real c6 = std::sqrt(c6i * c6j);
175 if (ljpme_combination_rule == LongRangeVdW::LB && !gmx_numzero(c6) && !gmx_numzero(c12i)
176 && !gmx_numzero(c12j))
178 real sigmai = gmx::sixthroot(c12i / c6i);
179 real sigmaj = gmx::sixthroot(c12j / c6j);
180 real epsi = c6i * c6i / c12i;
181 real epsj = c6j * c6j / c12j;
182 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
184 /* Store the elements at the same relative positions as C6 in nbfp in order
185 * to simplify access in the kernels
187 grid[2 * (numAtomTypes * i + j)] = c6 * 6.0;
193 //! What kind of constraint affects an atom
194 enum class ConstraintTypeForAtom : int
196 None, //!< No constraint active
197 Constraint, //!< F_CONSTR or F_CONSTRNC active
198 Settle, //! F_SETTLE active
201 static std::vector<gmx::AtomInfoWithinMoleculeBlock>
202 makeAtomInfoForEachMoleculeBlock(const gmx_mtop_t& mtop, const t_forcerec* fr)
204 std::vector<bool> atomUsesVdw(fr->ntype, false);
205 for (int ai = 0; ai < fr->ntype; ai++)
207 for (int j = 0; j < fr->ntype; j++)
209 atomUsesVdw[ai] = atomUsesVdw[ai] || fr->haveBuckingham || C6(fr->nbfp, fr->ntype, ai, j) != 0
210 || C12(fr->nbfp, fr->ntype, ai, j) != 0;
214 std::vector<gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock;
215 int indexOfFirstAtomInMoleculeBlock = 0;
216 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
218 const gmx_molblock_t& molb = mtop.molblock[mb];
219 const gmx_moltype_t& molt = mtop.moltype[molb.type];
220 const auto& excl = molt.excls;
222 /* Check if all molecules in this block have identical
223 * atominfo. (That's true unless some kind of group- or
224 * distance-based algorithm is involved, e.g. QM/MM groups
225 * affecting multiple molecules within a block differently.)
226 * If so, we only need an array of the size of one molecule.
227 * Otherwise we make an array of #mol times #atoms per
230 bool allMoleculesWithinBlockAreIdentical = true;
231 for (int m = 0; m < molb.nmol; m++)
233 const int numAtomsInAllMolecules = m * molt.atoms.nr;
234 for (int a = 0; a < molt.atoms.nr; a++)
236 if (getGroupType(mtop.groups,
237 SimulationAtomGroupType::QuantumMechanics,
238 indexOfFirstAtomInMoleculeBlock + numAtomsInAllMolecules + a)
239 != getGroupType(mtop.groups,
240 SimulationAtomGroupType::QuantumMechanics,
241 indexOfFirstAtomInMoleculeBlock + a))
243 allMoleculesWithinBlockAreIdentical = false;
245 if (!mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
247 if (mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics]
248 [indexOfFirstAtomInMoleculeBlock + numAtomsInAllMolecules + a]
249 != mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics]
250 [indexOfFirstAtomInMoleculeBlock + a])
252 allMoleculesWithinBlockAreIdentical = false;
258 gmx::AtomInfoWithinMoleculeBlock atomInfoOfMoleculeBlock;
259 atomInfoOfMoleculeBlock.indexOfFirstAtomInMoleculeBlock = indexOfFirstAtomInMoleculeBlock;
260 atomInfoOfMoleculeBlock.indexOfLastAtomInMoleculeBlock =
261 indexOfFirstAtomInMoleculeBlock + molb.nmol * molt.atoms.nr;
262 int atomInfoSize = (allMoleculesWithinBlockAreIdentical ? 1 : molb.nmol) * molt.atoms.nr;
263 atomInfoOfMoleculeBlock.atomInfo.resize(atomInfoSize);
265 /* Set constraints flags for constrained atoms */
266 std::vector<ConstraintTypeForAtom> constraintTypeOfAtom(molt.atoms.nr, ConstraintTypeForAtom::None);
267 for (int ftype = 0; ftype < F_NRE; ftype++)
269 if (interaction_function[ftype].flags & IF_CONSTRAINT)
271 const int nral = NRAL(ftype);
272 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
274 for (int a = 0; a < nral; a++)
276 constraintTypeOfAtom[molt.ilist[ftype].iatoms[ia + 1 + a]] =
277 (ftype == F_SETTLE ? ConstraintTypeForAtom::Settle
278 : ConstraintTypeForAtom::Constraint);
284 for (int m = 0; m < (allMoleculesWithinBlockAreIdentical ? 1 : molb.nmol); m++)
286 const int moleculeOffsetInBlock = m * molt.atoms.nr;
287 for (int a = 0; a < molt.atoms.nr; a++)
289 const t_atom& atom = molt.atoms.atom[a];
290 int64_t& atomInfo = atomInfoOfMoleculeBlock.atomInfo[moleculeOffsetInBlock + a];
292 /* Store the energy group in atomInfo */
293 int gid = getGroupType(mtop.groups,
294 SimulationAtomGroupType::EnergyOutput,
295 indexOfFirstAtomInMoleculeBlock + moleculeOffsetInBlock + a);
296 atomInfo = (atomInfo & ~gmx::sc_atomInfo_EnergyGroupIdMask) | gid;
298 bool bHaveVDW = (atomUsesVdw[atom.type] || atomUsesVdw[atom.typeB]);
299 bool bHaveQ = (atom.q != 0 || atom.qB != 0);
301 bool haveExclusions = false;
302 /* Loop over all the exclusions of atom ai */
303 for (const int j : excl[a])
307 haveExclusions = true;
312 switch (constraintTypeOfAtom[a])
314 case ConstraintTypeForAtom::Constraint:
315 atomInfo |= gmx::sc_atomInfo_Constraint;
317 case ConstraintTypeForAtom::Settle: atomInfo |= gmx::sc_atomInfo_Settle; break;
322 atomInfo |= gmx::sc_atomInfo_Exclusion;
326 atomInfo |= gmx::sc_atomInfo_HasVdw;
330 atomInfo |= gmx::sc_atomInfo_HasCharge;
332 if (fr->efep != FreeEnergyPerturbationType::No)
336 atomInfo |= gmx::sc_atomInfo_FreeEnergyPerturbation;
338 if (atomHasPerturbedChargeIn14Interaction(a, molt))
340 atomInfo |= gmx::sc_atomInfo_HasPerturbedChargeIn14Interaction;
346 atomInfoForEachMoleculeBlock.push_back(atomInfoOfMoleculeBlock);
348 indexOfFirstAtomInMoleculeBlock += molb.nmol * molt.atoms.nr;
351 return atomInfoForEachMoleculeBlock;
354 static std::vector<int64_t> expandAtomInfo(const int nmb,
355 gmx::ArrayRef<const gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock)
357 const int numAtoms = atomInfoForEachMoleculeBlock[nmb - 1].indexOfLastAtomInMoleculeBlock;
359 std::vector<int64_t> atomInfo(numAtoms);
362 for (int a = 0; a < numAtoms; a++)
364 while (a >= atomInfoForEachMoleculeBlock[mb].indexOfLastAtomInMoleculeBlock)
368 atomInfo[a] = atomInfoForEachMoleculeBlock[mb]
369 .atomInfo[(a - atomInfoForEachMoleculeBlock[mb].indexOfFirstAtomInMoleculeBlock)
370 % atomInfoForEachMoleculeBlock[mb].atomInfo.size()];
376 /* Sets the sum of charges (squared) and C6 in the system in fr.
377 * Returns whether the system has a net charge.
379 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t& mtop)
381 /*This now calculates sum for q and c6*/
382 double qsum, q2sum, q, c6sum, c6;
387 for (const gmx_molblock_t& molb : mtop.molblock)
389 int nmol = molb.nmol;
390 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
391 for (int i = 0; i < atoms->nr; i++)
393 q = atoms->atom[i].q;
395 q2sum += nmol * q * q;
396 c6 = mtop.ffparams.iparams[atoms->atom[i].type * (mtop.ffparams.atnr + 1)].lj.c6;
401 fr->q2sum[0] = q2sum;
402 fr->c6sum[0] = c6sum;
404 if (fr->efep != FreeEnergyPerturbationType::No)
409 for (const gmx_molblock_t& molb : mtop.molblock)
411 int nmol = molb.nmol;
412 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
413 for (int i = 0; i < atoms->nr; i++)
415 q = atoms->atom[i].qB;
417 q2sum += nmol * q * q;
418 c6 = mtop.ffparams.iparams[atoms->atom[i].typeB * (mtop.ffparams.atnr + 1)].lj.c6;
422 fr->q2sum[1] = q2sum;
423 fr->c6sum[1] = c6sum;
428 fr->qsum[1] = fr->qsum[0];
429 fr->q2sum[1] = fr->q2sum[0];
430 fr->c6sum[1] = fr->c6sum[0];
434 if (fr->efep == FreeEnergyPerturbationType::No)
436 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
440 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
444 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
445 return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
448 /*!\brief If there's bonded interactions of type \c ftype1 or \c
449 * ftype2 present in the topology, build an array of the number of
450 * interactions present for each bonded interaction index found in the
453 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
454 * valid type with that parameter.
456 * \c count will be reallocated as necessary to fit the largest bonded
457 * interaction index found, and its current size will be returned in
458 * \c ncount. It will contain zero for every bonded interaction index
459 * for which no interactions are present in the topology.
461 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t& mtop, int* ncount, int** count)
463 int ftype, i, j, tabnr;
465 // Loop over all moleculetypes
466 for (const gmx_moltype_t& molt : mtop.moltype)
468 // Loop over all interaction types
469 for (ftype = 0; ftype < F_NRE; ftype++)
471 // If the current interaction type is one of the types whose tables we're trying to count...
472 if (ftype == ftype1 || ftype == ftype2)
474 const InteractionList& il = molt.ilist[ftype];
475 const int stride = 1 + NRAL(ftype);
476 // ... and there are actually some interactions for this type
477 for (i = 0; i < il.size(); i += stride)
479 // Find out which table index the user wanted
480 tabnr = mtop.ffparams.iparams[il.iatoms[i]].tab.table;
483 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
485 // Make room for this index in the data structure
486 if (tabnr >= *ncount)
488 srenew(*count, tabnr + 1);
489 for (j = *ncount; j < tabnr + 1; j++)
495 // Record that this table index is used and must have a valid file
503 /*!\brief If there's bonded interactions of flavour \c tabext and type
504 * \c ftype1 or \c ftype2 present in the topology, seek them in the
505 * list of filenames passed to mdrun, and make bonded tables from
508 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
509 * valid type with that parameter.
511 * A fatal error occurs if no matching filename is found.
513 static std::vector<bondedtable_t> make_bonded_tables(FILE* fplog,
516 const gmx_mtop_t& mtop,
517 gmx::ArrayRef<const std::string> tabbfnm,
520 std::vector<bondedtable_t> tab;
523 int* count = nullptr;
524 count_tables(ftype1, ftype2, mtop, &ncount, &count);
526 // Are there any relevant tabulated bond interactions?
530 for (int i = 0; i < ncount; i++)
532 // Do any interactions exist that requires this table?
535 // This pattern enforces the current requirement that
536 // table filenames end in a characteristic sequence
537 // before the file type extension, and avoids table 13
538 // being recognized and used for table 1.
539 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
540 bool madeTable = false;
541 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
543 if (gmx::endsWith(tabbfnm[j], patternToFind))
545 // Finally read the table from the file found
546 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
552 bool isPlural = (ftype2 != -1);
554 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
555 "because no table file whose name matched '%s' was passed via the "
556 "gmx mdrun -tableb command-line option.",
557 interaction_function[ftype1].longname,
558 isPlural ? "' or '" : "",
559 isPlural ? interaction_function[ftype2].longname : "",
561 patternToFind.c_str());
571 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
573 fr->natoms_force = natoms_force;
574 fr->natoms_force_constr = natoms_force_constr;
576 for (auto& forceHelperBuffers : fr->forceHelperBuffers)
578 forceHelperBuffers.resize(natoms_f_novirsum);
582 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
584 * Tables are generated for one or both, depending on if the pointers
585 * are non-null. The spacing for both table sets is the same and obeys
586 * both accuracy requirements, when relevant.
588 static void init_ewald_f_table(const interaction_const_t& ic,
591 EwaldCorrectionTables* coulombTables,
592 EwaldCorrectionTables* vdwTables)
594 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
595 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
597 /* Get the Ewald table spacing based on Coulomb and/or LJ
598 * Ewald coefficients and rtol.
600 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
602 const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
604 real tableLen = ic.rcoulomb;
605 if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
607 /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
608 * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
609 * The alternative is to look through all the exclusions and check if they come from
610 * couple-intramol == no. Meanwhile, always having larger tables should only affect
611 * memory consumption, not speed (barring cache issues).
613 tableLen = rlist + tabext;
615 const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
620 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
625 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
629 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
631 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
634 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
638 "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
639 1 / ic->coulombEwaldTables->scale,
640 ic->coulombEwaldTables->tableF.size());
645 real cutoff_inf(real cutoff)
649 cutoff = GMX_CUTOFF_INF;
655 void init_forcerec(FILE* fplog,
656 const gmx::MDLogger& mdlog,
657 const gmx::SimulationWorkload& simulationWork,
658 t_forcerec* forcerec,
659 const t_inputrec& inputrec,
660 const gmx_mtop_t& mtop,
661 const t_commrec* commrec,
665 gmx::ArrayRef<const std::string> tabbfnm,
668 /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
669 forcerec->use_simd_kernels = GMX_USE_SIMD_KERNELS;
671 if (check_box(inputrec.pbcType, box))
673 gmx_fatal(FARGS, "%s", check_box(inputrec.pbcType, box));
676 /* Test particle insertion ? */
677 if (EI_TPI(inputrec.eI))
679 /* Set to the size of the molecule to be inserted (the last one) */
680 gmx::RangePartitioning molecules = gmx_mtop_molecules(mtop);
681 forcerec->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
688 if (inputrec.coulombtype == CoulombInteractionType::RFNecUnsupported
689 || inputrec.coulombtype == CoulombInteractionType::GRFNotused)
691 gmx_fatal(FARGS, "%s electrostatics is no longer supported", enumValueToString(inputrec.coulombtype));
694 if (inputrec.bAdress)
696 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
698 if (inputrec.useTwinRange)
700 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
702 /* Copy the user determined parameters */
703 forcerec->userint1 = inputrec.userint1;
704 forcerec->userint2 = inputrec.userint2;
705 forcerec->userint3 = inputrec.userint3;
706 forcerec->userint4 = inputrec.userint4;
707 forcerec->userreal1 = inputrec.userreal1;
708 forcerec->userreal2 = inputrec.userreal2;
709 forcerec->userreal3 = inputrec.userreal3;
710 forcerec->userreal4 = inputrec.userreal4;
713 forcerec->fc_stepsize = inputrec.fc_stepsize;
716 forcerec->efep = inputrec.efep;
718 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
720 forcerec->use_simd_kernels = FALSE;
721 if (fplog != nullptr)
724 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
725 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
726 "(e.g. SSE2/SSE4.1/AVX).\n\n");
730 forcerec->haveBuckingham = (mtop.ffparams.functype[0] == F_BHAM);
732 /* Neighbour searching stuff */
733 forcerec->pbcType = inputrec.pbcType;
735 /* Determine if we will do PBC for distances in bonded interactions */
736 if (forcerec->pbcType == PbcType::No)
738 forcerec->bMolPBC = FALSE;
742 forcerec->bMolPBC = (!DOMAINDECOMP(commrec) || dd_bonded_molpbc(*commrec->dd, forcerec->pbcType));
744 // Check and set up PBC for Ewald surface corrections or orientation restraints
745 const bool useEwaldSurfaceCorrection =
746 (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0);
747 const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
748 const bool moleculesAreAlwaysWhole =
749 (DOMAINDECOMP(commrec) && dd_moleculesAreAlwaysWhole(*commrec->dd));
750 // WholeMoleculeTransform is only supported with a single PP rank
751 if (!moleculesAreAlwaysWhole && !havePPDomainDecomposition(commrec)
752 && (useEwaldSurfaceCorrection || haveOrientationRestraints))
754 if (havePPDomainDecomposition(commrec))
757 "You requested Ewald surface correction or orientation restraints, "
758 "but molecules are broken "
759 "over periodic boundary conditions by the domain decomposition. "
760 "Run without domain decomposition instead.");
763 forcerec->wholeMoleculeTransform = std::make_unique<gmx::WholeMoleculeTransform>(
764 mtop, inputrec.pbcType, DOMAINDECOMP(commrec));
767 forcerec->bMolPBC = !DOMAINDECOMP(commrec) || dd_bonded_molpbc(*commrec->dd, forcerec->pbcType);
769 if (useEwaldSurfaceCorrection)
771 GMX_RELEASE_ASSERT(moleculesAreAlwaysWhole || forcerec->wholeMoleculeTransform,
772 "Molecules can not be broken by PBC with epsilon_surface > 0");
776 forcerec->rc_scaling = inputrec.refcoord_scaling;
777 copy_rvec(inputrec.posres_com, forcerec->posres_com);
778 copy_rvec(inputrec.posres_comB, forcerec->posres_comB);
779 forcerec->rlist = cutoff_inf(inputrec.rlist);
780 forcerec->ljpme_combination_rule = inputrec.ljpme_combination_rule;
782 /* This now calculates sum for q and c6*/
783 bool systemHasNetCharge = set_chargesum(fplog, forcerec, mtop);
785 /* Make data structure used by kernels */
786 forcerec->ic = std::make_unique<interaction_const_t>(
787 init_interaction_const(fplog, inputrec, mtop, systemHasNetCharge));
788 init_interaction_const_tables(fplog, forcerec->ic.get(), forcerec->rlist, inputrec.tabext);
790 const interaction_const_t* interactionConst = forcerec->ic.get();
792 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
793 switch (interactionConst->eeltype)
795 case CoulombInteractionType::Cut:
796 forcerec->nbkernel_elec_interaction = NbkernelElecType::Coulomb;
799 case CoulombInteractionType::RF:
800 case CoulombInteractionType::RFZero:
801 forcerec->nbkernel_elec_interaction = NbkernelElecType::ReactionField;
804 case CoulombInteractionType::Switch:
805 case CoulombInteractionType::Shift:
806 case CoulombInteractionType::User:
807 case CoulombInteractionType::PmeSwitch:
808 case CoulombInteractionType::PmeUser:
809 case CoulombInteractionType::PmeUserSwitch:
810 forcerec->nbkernel_elec_interaction = NbkernelElecType::CubicSplineTable;
813 case CoulombInteractionType::Pme:
814 case CoulombInteractionType::P3mAD:
815 case CoulombInteractionType::Ewald:
816 forcerec->nbkernel_elec_interaction = NbkernelElecType::Ewald;
821 "Unsupported electrostatic interaction: %s",
822 enumValueToString(interactionConst->eeltype));
824 forcerec->nbkernel_elec_modifier = interactionConst->coulomb_modifier;
826 /* Vdw: Translate from mdp settings to kernel format */
827 switch (interactionConst->vdwtype)
829 case VanDerWaalsType::Cut:
830 if (forcerec->haveBuckingham)
832 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::Buckingham;
836 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LennardJones;
839 case VanDerWaalsType::Pme:
840 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LJEwald;
843 case VanDerWaalsType::Switch:
844 case VanDerWaalsType::Shift:
845 case VanDerWaalsType::User:
846 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::CubicSplineTable;
850 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", enumValueToString(interactionConst->vdwtype));
852 forcerec->nbkernel_vdw_modifier = interactionConst->vdw_modifier;
854 if (!gmx_within_tol(interactionConst->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
856 gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
858 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
859 * while mdrun does not (and never did) support this.
861 if (EEL_USER(forcerec->ic->eeltype))
864 "Electrostatics type %s is currently not supported",
865 enumValueToString(inputrec.coulombtype));
868 /* 1-4 interaction electrostatics */
869 forcerec->fudgeQQ = mtop.ffparams.fudgeQQ;
871 if (simulationWork.useMts)
873 GMX_ASSERT(gmx::checkMtsRequirements(inputrec).empty(),
874 "All MTS requirements should be met here");
877 const bool haveDirectVirialContributionsFast =
878 forcerec->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
879 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || inputrec.nwall > 0 || inputrec.bPull
880 || inputrec.bRot || inputrec.bIMD;
881 const bool haveDirectVirialContributionsSlow =
882 EEL_FULL(interactionConst->eeltype) || EVDW_PME(interactionConst->vdwtype);
883 for (int i = 0; i < (simulationWork.useMts ? 2 : 1); i++)
885 bool haveDirectVirialContributions =
886 (((!simulationWork.useMts || i == 0) && haveDirectVirialContributionsFast)
887 || ((!simulationWork.useMts || i == 1) && haveDirectVirialContributionsSlow));
888 forcerec->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
891 if (forcerec->shift_vec.empty())
893 forcerec->shift_vec.resize(gmx::c_numShiftVectors);
896 GMX_ASSERT(forcerec->nbfp.empty(), "The nonbonded force parameters should not be set up yet.");
897 forcerec->ntype = mtop.ffparams.atnr;
898 forcerec->nbfp = makeNonBondedParameterLists(
899 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->haveBuckingham);
900 if (EVDW_PME(interactionConst->vdwtype))
902 forcerec->ljpme_c6grid = makeLJPmeC6GridCorrectionParameters(
903 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->ljpme_combination_rule);
906 /* Copy the energy group exclusions */
907 forcerec->egp_flags = inputrec.opts.egp_flags;
909 /* Van der Waals stuff */
910 if ((interactionConst->vdwtype != VanDerWaalsType::Cut)
911 && (interactionConst->vdwtype != VanDerWaalsType::User) && !forcerec->haveBuckingham)
913 if (interactionConst->rvdw_switch >= interactionConst->rvdw)
916 "rvdw_switch (%f) must be < rvdw (%f)",
917 interactionConst->rvdw_switch,
918 interactionConst->rvdw);
923 "Using %s Lennard-Jones, switch between %g and %g nm\n",
924 (interactionConst->eeltype == CoulombInteractionType::Switch) ? "switched" : "shifted",
925 interactionConst->rvdw_switch,
926 interactionConst->rvdw);
930 if (forcerec->haveBuckingham && EVDW_PME(interactionConst->vdwtype))
932 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
935 if (forcerec->haveBuckingham
936 && (interactionConst->vdwtype == VanDerWaalsType::Shift
937 || interactionConst->vdwtype == VanDerWaalsType::Switch))
939 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
942 if (forcerec->haveBuckingham)
944 gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
947 if (inputrec.implicit_solvent)
949 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
953 /* This code automatically gives table length tabext without cut-off's,
954 * in that case grompp should already have checked that we do not need
955 * normal tables and we only generate tables for 1-4 interactions.
957 real rtab = inputrec.rlist + inputrec.tabext;
959 /* We want to use unmodified tables for 1-4 coulombic
960 * interactions, so we must in general have an extra set of
962 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
963 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
965 forcerec->pairsTable = make_tables(fplog, interactionConst, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
969 forcerec->nwall = inputrec.nwall;
970 if (inputrec.nwall && inputrec.wall_type == WallType::Table)
972 make_wall_tables(fplog, inputrec, tabfn, &mtop.groups, forcerec);
975 forcerec->fcdata = std::make_unique<t_fcdata>();
977 if (!tabbfnm.empty())
979 t_fcdata& fcdata = *forcerec->fcdata;
980 // Need to catch std::bad_alloc
981 // TODO Don't need to catch this here, when merging with master branch
984 // TODO move these tables into a separate struct and store reference in ListedForces
985 fcdata.bondtab = make_bonded_tables(fplog, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
986 fcdata.angletab = make_bonded_tables(fplog, F_TABANGLES, -1, mtop, tabbfnm, "a");
987 fcdata.dihtab = make_bonded_tables(fplog, F_TABDIHS, -1, mtop, tabbfnm, "d");
989 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
996 "No fcdata or table file name passed, can not read table, can not do bonded "
1001 /* Initialize the thread working data for bonded interactions */
1002 if (simulationWork.useMts)
1004 // Add one ListedForces object for each MTS level
1005 bool isFirstLevel = true;
1006 for (const auto& mtsLevel : inputrec.mtsLevels)
1008 ListedForces::InteractionSelection interactionSelection;
1009 const auto& forceGroups = mtsLevel.forceGroups;
1010 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1012 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1014 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1016 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1018 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1020 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1024 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1025 isFirstLevel = false;
1027 forcerec->listedForces.emplace_back(
1029 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1030 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1031 interactionSelection,
1037 // Add one ListedForces object with all listed interactions
1038 forcerec->listedForces.emplace_back(
1040 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1041 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1042 ListedForces::interactionSelectionAll(),
1046 // QM/MM initialization if requested
1049 gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1052 /* Set all the static charge group info */
1053 forcerec->atomInfoForEachMoleculeBlock = makeAtomInfoForEachMoleculeBlock(mtop, forcerec);
1054 if (!DOMAINDECOMP(commrec))
1056 forcerec->atomInfo = expandAtomInfo(mtop.molblock.size(), forcerec->atomInfoForEachMoleculeBlock);
1059 if (!DOMAINDECOMP(commrec))
1061 forcerec_set_ranges(forcerec, mtop.natoms, mtop.natoms, mtop.natoms);
1064 forcerec->print_force = print_force;
1066 if (inputrec.eDispCorr != DispersionCorrectionType::No)
1068 forcerec->dispersionCorrection = std::make_unique<DispersionCorrection>(
1069 mtop, inputrec, forcerec->haveBuckingham, forcerec->ntype, forcerec->nbfp, *forcerec->ic, tabfn);
1070 forcerec->dispersionCorrection->print(mdlog);
1073 if (fplog != nullptr)
1075 /* Here we switch from using mdlog, which prints the newline before
1076 * the paragraph, to our old fprintf logging, which prints the newline
1077 * after the paragraph, so we should add a newline here.
1079 fprintf(fplog, "\n");
1083 t_forcerec::t_forcerec() = default;
1085 t_forcerec::~t_forcerec() = default;