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)
337 atomInfo |= gmx::sc_atomInfo_FreeEnergyPerturbation;
339 if (atomHasPerturbedChargeIn14Interaction(a, molt))
341 atomInfo |= gmx::sc_atomInfo_HasPerturbedChargeIn14Interaction;
347 atomInfoForEachMoleculeBlock.push_back(atomInfoOfMoleculeBlock);
349 indexOfFirstAtomInMoleculeBlock += molb.nmol * molt.atoms.nr;
352 return atomInfoForEachMoleculeBlock;
355 static std::vector<int64_t> expandAtomInfo(const int nmb,
356 gmx::ArrayRef<const gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock)
358 const int numAtoms = atomInfoForEachMoleculeBlock[nmb - 1].indexOfLastAtomInMoleculeBlock;
360 std::vector<int64_t> atomInfo(numAtoms);
363 for (int a = 0; a < numAtoms; a++)
365 while (a >= atomInfoForEachMoleculeBlock[mb].indexOfLastAtomInMoleculeBlock)
369 atomInfo[a] = atomInfoForEachMoleculeBlock[mb]
370 .atomInfo[(a - atomInfoForEachMoleculeBlock[mb].indexOfFirstAtomInMoleculeBlock)
371 % atomInfoForEachMoleculeBlock[mb].atomInfo.size()];
377 /* Sets the sum of charges (squared) and C6 in the system in fr.
378 * Returns whether the system has a net charge.
380 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t& mtop)
382 /*This now calculates sum for q and c6*/
383 double qsum, q2sum, q, c6sum, c6;
388 for (const gmx_molblock_t& molb : mtop.molblock)
390 int nmol = molb.nmol;
391 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
392 for (int i = 0; i < atoms->nr; i++)
394 q = atoms->atom[i].q;
396 q2sum += nmol * q * q;
397 c6 = mtop.ffparams.iparams[atoms->atom[i].type * (mtop.ffparams.atnr + 1)].lj.c6;
402 fr->q2sum[0] = q2sum;
403 fr->c6sum[0] = c6sum;
405 if (fr->efep != FreeEnergyPerturbationType::No)
410 for (const gmx_molblock_t& molb : mtop.molblock)
412 int nmol = molb.nmol;
413 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
414 for (int i = 0; i < atoms->nr; i++)
416 q = atoms->atom[i].qB;
418 q2sum += nmol * q * q;
419 c6 = mtop.ffparams.iparams[atoms->atom[i].typeB * (mtop.ffparams.atnr + 1)].lj.c6;
423 fr->q2sum[1] = q2sum;
424 fr->c6sum[1] = c6sum;
429 fr->qsum[1] = fr->qsum[0];
430 fr->q2sum[1] = fr->q2sum[0];
431 fr->c6sum[1] = fr->c6sum[0];
435 if (fr->efep == FreeEnergyPerturbationType::No)
437 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
441 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
445 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
446 return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
449 /*!\brief If there's bonded interactions of type \c ftype1 or \c
450 * ftype2 present in the topology, build an array of the number of
451 * interactions present for each bonded interaction index found in the
454 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
455 * valid type with that parameter.
457 * \c count will be reallocated as necessary to fit the largest bonded
458 * interaction index found, and its current size will be returned in
459 * \c ncount. It will contain zero for every bonded interaction index
460 * for which no interactions are present in the topology.
462 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t& mtop, int* ncount, int** count)
464 int ftype, i, j, tabnr;
466 // Loop over all moleculetypes
467 for (const gmx_moltype_t& molt : mtop.moltype)
469 // Loop over all interaction types
470 for (ftype = 0; ftype < F_NRE; ftype++)
472 // If the current interaction type is one of the types whose tables we're trying to count...
473 if (ftype == ftype1 || ftype == ftype2)
475 const InteractionList& il = molt.ilist[ftype];
476 const int stride = 1 + NRAL(ftype);
477 // ... and there are actually some interactions for this type
478 for (i = 0; i < il.size(); i += stride)
480 // Find out which table index the user wanted
481 tabnr = mtop.ffparams.iparams[il.iatoms[i]].tab.table;
484 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
486 // Make room for this index in the data structure
487 if (tabnr >= *ncount)
489 srenew(*count, tabnr + 1);
490 for (j = *ncount; j < tabnr + 1; j++)
496 // Record that this table index is used and must have a valid file
504 /*!\brief If there's bonded interactions of flavour \c tabext and type
505 * \c ftype1 or \c ftype2 present in the topology, seek them in the
506 * list of filenames passed to mdrun, and make bonded tables from
509 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
510 * valid type with that parameter.
512 * A fatal error occurs if no matching filename is found.
514 static std::vector<bondedtable_t> make_bonded_tables(FILE* fplog,
517 const gmx_mtop_t& mtop,
518 gmx::ArrayRef<const std::string> tabbfnm,
521 std::vector<bondedtable_t> tab;
524 int* count = nullptr;
525 count_tables(ftype1, ftype2, mtop, &ncount, &count);
527 // Are there any relevant tabulated bond interactions?
531 for (int i = 0; i < ncount; i++)
533 // Do any interactions exist that requires this table?
536 // This pattern enforces the current requirement that
537 // table filenames end in a characteristic sequence
538 // before the file type extension, and avoids table 13
539 // being recognized and used for table 1.
540 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
541 bool madeTable = false;
542 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
544 if (gmx::endsWith(tabbfnm[j], patternToFind))
546 // Finally read the table from the file found
547 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
553 bool isPlural = (ftype2 != -1);
555 "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
556 "because no table file whose name matched '%s' was passed via the "
557 "gmx mdrun -tableb command-line option.",
558 interaction_function[ftype1].longname,
559 isPlural ? "' or '" : "",
560 isPlural ? interaction_function[ftype2].longname : "",
562 patternToFind.c_str());
572 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
574 fr->natoms_force = natoms_force;
575 fr->natoms_force_constr = natoms_force_constr;
577 for (auto& forceHelperBuffers : fr->forceHelperBuffers)
579 forceHelperBuffers.resize(natoms_f_novirsum);
583 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
585 * Tables are generated for one or both, depending on if the pointers
586 * are non-null. The spacing for both table sets is the same and obeys
587 * both accuracy requirements, when relevant.
589 static void init_ewald_f_table(const interaction_const_t& ic,
592 EwaldCorrectionTables* coulombTables,
593 EwaldCorrectionTables* vdwTables)
595 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
596 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
598 /* Get the Ewald table spacing based on Coulomb and/or LJ
599 * Ewald coefficients and rtol.
601 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
603 const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
605 real tableLen = ic.rcoulomb;
606 if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
608 /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
609 * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
610 * The alternative is to look through all the exclusions and check if they come from
611 * couple-intramol == no. Meanwhile, always having larger tables should only affect
612 * memory consumption, not speed (barring cache issues).
614 tableLen = rlist + tabext;
616 const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
621 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
626 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
630 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
632 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
635 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
639 "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
640 1 / ic->coulombEwaldTables->scale,
641 ic->coulombEwaldTables->tableF.size());
646 real cutoff_inf(real cutoff)
650 cutoff = GMX_CUTOFF_INF;
656 void init_forcerec(FILE* fplog,
657 const gmx::MDLogger& mdlog,
658 const gmx::SimulationWorkload& simulationWork,
659 t_forcerec* forcerec,
660 const t_inputrec& inputrec,
661 const gmx_mtop_t& mtop,
662 const t_commrec* commrec,
666 gmx::ArrayRef<const std::string> tabbfnm,
669 /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
670 forcerec->use_simd_kernels = GMX_USE_SIMD_KERNELS;
672 if (check_box(inputrec.pbcType, box))
674 gmx_fatal(FARGS, "%s", check_box(inputrec.pbcType, box));
677 /* Test particle insertion ? */
678 if (EI_TPI(inputrec.eI))
680 /* Set to the size of the molecule to be inserted (the last one) */
681 gmx::RangePartitioning molecules = gmx_mtop_molecules(mtop);
682 forcerec->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
689 if (inputrec.coulombtype == CoulombInteractionType::RFNecUnsupported
690 || inputrec.coulombtype == CoulombInteractionType::GRFNotused)
692 gmx_fatal(FARGS, "%s electrostatics is no longer supported", enumValueToString(inputrec.coulombtype));
695 if (inputrec.bAdress)
697 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
699 if (inputrec.useTwinRange)
701 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
703 /* Copy the user determined parameters */
704 forcerec->userint1 = inputrec.userint1;
705 forcerec->userint2 = inputrec.userint2;
706 forcerec->userint3 = inputrec.userint3;
707 forcerec->userint4 = inputrec.userint4;
708 forcerec->userreal1 = inputrec.userreal1;
709 forcerec->userreal2 = inputrec.userreal2;
710 forcerec->userreal3 = inputrec.userreal3;
711 forcerec->userreal4 = inputrec.userreal4;
714 forcerec->fc_stepsize = inputrec.fc_stepsize;
717 forcerec->efep = inputrec.efep;
719 if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
721 forcerec->use_simd_kernels = FALSE;
722 if (fplog != nullptr)
725 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
726 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
727 "(e.g. SSE2/SSE4.1/AVX).\n\n");
731 forcerec->haveBuckingham = (mtop.ffparams.functype[0] == F_BHAM);
733 /* Neighbour searching stuff */
734 forcerec->pbcType = inputrec.pbcType;
736 /* Determine if we will do PBC for distances in bonded interactions */
737 if (forcerec->pbcType == PbcType::No)
739 forcerec->bMolPBC = FALSE;
743 const bool useEwaldSurfaceCorrection =
744 (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0);
745 const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
746 if (!DOMAINDECOMP(commrec))
748 forcerec->bMolPBC = true;
750 if (useEwaldSurfaceCorrection || haveOrientationRestraints)
752 forcerec->wholeMoleculeTransform =
753 std::make_unique<gmx::WholeMoleculeTransform>(mtop, inputrec.pbcType);
758 forcerec->bMolPBC = dd_bonded_molpbc(*commrec->dd, forcerec->pbcType);
760 /* With Ewald surface correction it is useful to support e.g. running water
761 * in parallel with update groups.
762 * With orientation restraints there is no sensible use case supported with DD.
764 if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*commrec->dd))
765 || haveOrientationRestraints)
768 "You requested Ewald surface correction or orientation restraints, "
769 "but molecules are broken "
770 "over periodic boundary conditions by the domain decomposition. "
771 "Run without domain decomposition instead.");
775 if (useEwaldSurfaceCorrection)
777 GMX_RELEASE_ASSERT(!DOMAINDECOMP(commrec) || dd_moleculesAreAlwaysWhole(*commrec->dd),
778 "Molecules can not be broken by PBC with epsilon_surface > 0");
782 forcerec->rc_scaling = inputrec.refcoord_scaling;
783 copy_rvec(inputrec.posres_com, forcerec->posres_com);
784 copy_rvec(inputrec.posres_comB, forcerec->posres_comB);
785 forcerec->rlist = cutoff_inf(inputrec.rlist);
786 forcerec->ljpme_combination_rule = inputrec.ljpme_combination_rule;
788 /* This now calculates sum for q and c6*/
789 bool systemHasNetCharge = set_chargesum(fplog, forcerec, mtop);
791 /* Make data structure used by kernels */
792 forcerec->ic = std::make_unique<interaction_const_t>(
793 init_interaction_const(fplog, inputrec, mtop, systemHasNetCharge));
794 init_interaction_const_tables(fplog, forcerec->ic.get(), forcerec->rlist, inputrec.tabext);
796 const interaction_const_t* interactionConst = forcerec->ic.get();
798 /* TODO: Replace this Ewald table or move it into interaction_const_t */
799 if (inputrec.coulombtype == CoulombInteractionType::Ewald)
801 forcerec->ewald_table = std::make_unique<gmx_ewald_tab_t>(inputrec, fplog);
804 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
805 switch (interactionConst->eeltype)
807 case CoulombInteractionType::Cut:
808 forcerec->nbkernel_elec_interaction = NbkernelElecType::Coulomb;
811 case CoulombInteractionType::RF:
812 case CoulombInteractionType::RFZero:
813 forcerec->nbkernel_elec_interaction = NbkernelElecType::ReactionField;
816 case CoulombInteractionType::Switch:
817 case CoulombInteractionType::Shift:
818 case CoulombInteractionType::User:
819 case CoulombInteractionType::PmeSwitch:
820 case CoulombInteractionType::PmeUser:
821 case CoulombInteractionType::PmeUserSwitch:
822 forcerec->nbkernel_elec_interaction = NbkernelElecType::CubicSplineTable;
825 case CoulombInteractionType::Pme:
826 case CoulombInteractionType::P3mAD:
827 case CoulombInteractionType::Ewald:
828 forcerec->nbkernel_elec_interaction = NbkernelElecType::Ewald;
833 "Unsupported electrostatic interaction: %s",
834 enumValueToString(interactionConst->eeltype));
836 forcerec->nbkernel_elec_modifier = interactionConst->coulomb_modifier;
838 /* Vdw: Translate from mdp settings to kernel format */
839 switch (interactionConst->vdwtype)
841 case VanDerWaalsType::Cut:
842 if (forcerec->haveBuckingham)
844 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::Buckingham;
848 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LennardJones;
851 case VanDerWaalsType::Pme:
852 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LJEwald;
855 case VanDerWaalsType::Switch:
856 case VanDerWaalsType::Shift:
857 case VanDerWaalsType::User:
858 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::CubicSplineTable;
862 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", enumValueToString(interactionConst->vdwtype));
864 forcerec->nbkernel_vdw_modifier = interactionConst->vdw_modifier;
866 if (!gmx_within_tol(interactionConst->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
868 gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
870 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
871 * while mdrun does not (and never did) support this.
873 if (EEL_USER(forcerec->ic->eeltype))
876 "Electrostatics type %s is currently not supported",
877 enumValueToString(inputrec.coulombtype));
880 /* 1-4 interaction electrostatics */
881 forcerec->fudgeQQ = mtop.ffparams.fudgeQQ;
883 if (simulationWork.useMts)
885 GMX_ASSERT(gmx::checkMtsRequirements(inputrec).empty(),
886 "All MTS requirements should be met here");
889 const bool haveDirectVirialContributionsFast =
890 forcerec->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
891 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || inputrec.nwall > 0 || inputrec.bPull
892 || inputrec.bRot || inputrec.bIMD;
893 const bool haveDirectVirialContributionsSlow =
894 EEL_FULL(interactionConst->eeltype) || EVDW_PME(interactionConst->vdwtype);
895 for (int i = 0; i < (simulationWork.useMts ? 2 : 1); i++)
897 bool haveDirectVirialContributions =
898 (((!simulationWork.useMts || i == 0) && haveDirectVirialContributionsFast)
899 || ((!simulationWork.useMts || i == 1) && haveDirectVirialContributionsSlow));
900 forcerec->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
903 if (forcerec->shift_vec.empty())
905 forcerec->shift_vec.resize(gmx::c_numShiftVectors);
908 GMX_ASSERT(forcerec->nbfp.empty(), "The nonbonded force parameters should not be set up yet.");
909 forcerec->ntype = mtop.ffparams.atnr;
910 forcerec->nbfp = makeNonBondedParameterLists(
911 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->haveBuckingham);
912 if (EVDW_PME(interactionConst->vdwtype))
914 forcerec->ljpme_c6grid = makeLJPmeC6GridCorrectionParameters(
915 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->ljpme_combination_rule);
918 /* Copy the energy group exclusions */
919 forcerec->egp_flags = inputrec.opts.egp_flags;
921 /* Van der Waals stuff */
922 if ((interactionConst->vdwtype != VanDerWaalsType::Cut)
923 && (interactionConst->vdwtype != VanDerWaalsType::User) && !forcerec->haveBuckingham)
925 if (interactionConst->rvdw_switch >= interactionConst->rvdw)
928 "rvdw_switch (%f) must be < rvdw (%f)",
929 interactionConst->rvdw_switch,
930 interactionConst->rvdw);
935 "Using %s Lennard-Jones, switch between %g and %g nm\n",
936 (interactionConst->eeltype == CoulombInteractionType::Switch) ? "switched" : "shifted",
937 interactionConst->rvdw_switch,
938 interactionConst->rvdw);
942 if (forcerec->haveBuckingham && EVDW_PME(interactionConst->vdwtype))
944 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
947 if (forcerec->haveBuckingham
948 && (interactionConst->vdwtype == VanDerWaalsType::Shift
949 || interactionConst->vdwtype == VanDerWaalsType::Switch))
951 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
954 if (forcerec->haveBuckingham)
956 gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
959 if (inputrec.implicit_solvent)
961 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
965 /* This code automatically gives table length tabext without cut-off's,
966 * in that case grompp should already have checked that we do not need
967 * normal tables and we only generate tables for 1-4 interactions.
969 real rtab = inputrec.rlist + inputrec.tabext;
971 /* We want to use unmodified tables for 1-4 coulombic
972 * interactions, so we must in general have an extra set of
974 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
975 || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
977 forcerec->pairsTable = make_tables(fplog, interactionConst, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
981 forcerec->nwall = inputrec.nwall;
982 if (inputrec.nwall && inputrec.wall_type == WallType::Table)
984 make_wall_tables(fplog, inputrec, tabfn, &mtop.groups, forcerec);
987 forcerec->fcdata = std::make_unique<t_fcdata>();
989 if (!tabbfnm.empty())
991 t_fcdata& fcdata = *forcerec->fcdata;
992 // Need to catch std::bad_alloc
993 // TODO Don't need to catch this here, when merging with master branch
996 // TODO move these tables into a separate struct and store reference in ListedForces
997 fcdata.bondtab = make_bonded_tables(fplog, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
998 fcdata.angletab = make_bonded_tables(fplog, F_TABANGLES, -1, mtop, tabbfnm, "a");
999 fcdata.dihtab = make_bonded_tables(fplog, F_TABDIHS, -1, mtop, tabbfnm, "d");
1001 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1008 "No fcdata or table file name passed, can not read table, can not do bonded "
1013 /* Initialize the thread working data for bonded interactions */
1014 if (simulationWork.useMts)
1016 // Add one ListedForces object for each MTS level
1017 bool isFirstLevel = true;
1018 for (const auto& mtsLevel : inputrec.mtsLevels)
1020 ListedForces::InteractionSelection interactionSelection;
1021 const auto& forceGroups = mtsLevel.forceGroups;
1022 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1024 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1026 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1028 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1030 if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1032 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1036 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1037 isFirstLevel = false;
1039 forcerec->listedForces.emplace_back(
1041 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1042 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1043 interactionSelection,
1049 // Add one ListedForces object with all listed interactions
1050 forcerec->listedForces.emplace_back(
1052 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1053 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1054 ListedForces::interactionSelectionAll(),
1058 // QM/MM initialization if requested
1061 gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1064 /* Set all the static charge group info */
1065 forcerec->atomInfoForEachMoleculeBlock = makeAtomInfoForEachMoleculeBlock(mtop, forcerec);
1066 if (!DOMAINDECOMP(commrec))
1068 forcerec->atomInfo = expandAtomInfo(mtop.molblock.size(), forcerec->atomInfoForEachMoleculeBlock);
1071 if (!DOMAINDECOMP(commrec))
1073 forcerec_set_ranges(forcerec, mtop.natoms, mtop.natoms, mtop.natoms);
1076 forcerec->print_force = print_force;
1078 forcerec->nthread_ewc = gmx_omp_nthreads_get(ModuleMultiThread::Bonded);
1079 forcerec->ewc_t.resize(forcerec->nthread_ewc);
1081 if (inputrec.eDispCorr != DispersionCorrectionType::No)
1083 forcerec->dispersionCorrection = std::make_unique<DispersionCorrection>(
1084 mtop, inputrec, forcerec->haveBuckingham, forcerec->ntype, forcerec->nbfp, *forcerec->ic, tabfn);
1085 forcerec->dispersionCorrection->print(mdlog);
1088 if (fplog != nullptr)
1090 /* Here we switch from using mdlog, which prints the newline before
1091 * the paragraph, to our old fprintf logging, which prints the newline
1092 * after the paragraph, so we should add a newline here.
1094 fprintf(fplog, "\n");
1098 t_forcerec::t_forcerec() = default;
1100 t_forcerec::~t_forcerec() = default;