2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file defines functions used in making the
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "gromacs/domdec/reversetopology.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/domdec/domdec_constraints.h"
56 #include "gromacs/domdec/domdec_internal.h"
57 #include "gromacs/domdec/domdec_vsite.h"
58 #include "gromacs/domdec/options.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/mdlib/vsite.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/topology/atoms.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/topsort.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/fatalerror.h"
69 using gmx::DDBondedChecking;
70 using gmx::ListOfLists;
73 /*! \brief Struct for the reverse topology: links bonded interactions to atoms */
74 struct gmx_reverse_top_t::Impl
76 //! Constructs a reverse topology from \p mtop
77 Impl(const gmx_mtop_t& mtop, bool useFreeEnergy, const ReverseTopOptions& reverseTopOptions);
79 //! @cond Doxygen_Suppress
80 //! Options for the setup of this reverse topology
81 const ReverseTopOptions options;
82 //! Are there interaction of type F_POSRES and/or F_FBPOSRES
83 bool hasPositionRestraints;
84 //! \brief Are there bondeds/exclusions between atoms?
85 bool bInterAtomicInteractions = false;
86 //! \brief Reverse ilist for all moltypes
87 std::vector<reverse_ilist_t> ril_mt;
88 //! \brief The size of ril_mt[?].index summed over all entries
89 int ril_mt_tot_size = 0;
90 //! \brief The sorting state of bondeds for free energy
91 int ilsort = ilsortUNKNOWN;
92 //! \brief molblock to global atom index for quick lookup of molblocks on atom index
93 std::vector<MolblockIndices> mbi;
95 //! \brief Do we have intermolecular interactions?
96 bool bIntermolecularInteractions = false;
97 //! \brief Intermolecular reverse ilist
98 reverse_ilist_t ril_intermol;
100 //! The number of bonded interactions computed from the full topology
101 int expectedNumGlobalBondedInteractions = 0;
103 /* Work data structures for multi-threading */
104 //! \brief Thread work array for local topology generation
105 std::vector<thread_work_t> th_work;
110 int nral_rt(int ftype)
112 int nral = NRAL(ftype);
113 if (interaction_function[ftype].flags & IF_VSITE)
115 /* With vsites the reverse topology contains an extra entry
116 * for storing if constructing atoms are vsites.
124 bool dd_check_ftype(const int ftype, const ReverseTopOptions& rtOptions)
126 return ((((interaction_function[ftype].flags & IF_BOND) != 0U)
127 && ((interaction_function[ftype].flags & IF_VSITE) == 0U)
128 && ((rtOptions.ddBondedChecking == DDBondedChecking::All)
129 || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U)))
130 || (rtOptions.includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
131 || (rtOptions.includeSettles && ftype == F_SETTLE));
134 void global_atomnr_to_moltype_ind(ArrayRef<const MolblockIndices> molblockIndices,
141 const MolblockIndices* mbi = molblockIndices.data();
143 int end = molblockIndices.size(); /* exclusive */
146 /* binary search for molblock_ind */
149 mid = (start + end) / 2;
150 if (i_gl >= mbi[mid].a_end)
154 else if (i_gl < mbi[mid].a_start)
168 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
169 *i_mol = (i_gl - mbi->a_start) - (*mol) * mbi->natoms_mol;
172 /*! \brief Returns the maximum number of exclusions per atom */
173 static int getMaxNumExclusionsPerAtom(const ListOfLists<int>& excls)
176 for (gmx::index at = 0; at < excls.ssize(); at++)
178 const auto list = excls[at];
179 const int numExcls = list.ssize();
181 GMX_RELEASE_ASSERT(numExcls != 1 || list[0] == at,
182 "With 1 exclusion we expect a self-exclusion");
184 maxNumExcls = std::max(maxNumExcls, numExcls);
190 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
191 static void low_make_reverse_ilist(const InteractionLists& il_mt,
194 const ReverseTopOptions& rtOptions,
195 gmx::ArrayRef<const int> r_index,
196 gmx::ArrayRef<int> r_il,
197 const AtomLinkRule atomLinkRule,
198 const bool assignReverseIlist)
200 const bool includeConstraints = rtOptions.includeConstraints;
201 const bool includeSettles = rtOptions.includeSettles;
203 for (int ftype = 0; ftype < F_NRE; ftype++)
205 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
206 || (includeConstraints && (ftype == F_CONSTR || ftype == F_CONSTRNC))
207 || (includeSettles && ftype == F_SETTLE))
209 const bool isVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
210 const int nral = NRAL(ftype);
211 const auto& il = il_mt[ftype];
212 for (int i = 0; i < il.size(); i += 1 + nral)
214 const int* ia = il.iatoms.data() + i;
215 // Virtual sites should not be linked for bonded interactions
216 const int nlink = (atomLinkRule == AtomLinkRule::FirstAtom) ? 1 : (isVSite ? 0 : nral);
217 for (int link = 0; link < nlink; link++)
219 const int a = ia[1 + link];
220 if (assignReverseIlist)
222 GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
223 GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
224 r_il[r_index[a] + count[a]] = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
225 r_il[r_index[a] + count[a] + 1] = ia[0];
226 for (int j = 1; j < 1 + nral; j++)
228 /* Store the molecular atom number */
229 r_il[r_index[a] + count[a] + 1 + j] = ia[j];
232 if (interaction_function[ftype].flags & IF_VSITE)
234 if (assignReverseIlist)
236 /* Add an entry to iatoms for storing
237 * which of the constructing atoms are
240 r_il[r_index[a] + count[a] + 2 + nral] = 0;
241 for (int j = 2; j < 1 + nral; j++)
243 if (atom[ia[j]].ptype == ParticleType::VSite)
245 r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
250 count[a] += 2 + nral_rt(ftype);
257 void make_reverse_ilist(const InteractionLists& ilist,
258 const t_atoms* atoms,
259 const ReverseTopOptions& rtOptions,
260 const AtomLinkRule atomLinkRule,
261 reverse_ilist_t* ril_mt)
263 /* Count the interactions */
264 const int nat_mt = atoms->nr;
265 std::vector<int> count(nat_mt);
266 low_make_reverse_ilist(ilist, atoms->atom, count.data(), rtOptions, {}, {}, atomLinkRule, FALSE);
268 ril_mt->index.push_back(0);
269 for (int i = 0; i < nat_mt; i++)
271 ril_mt->index.push_back(ril_mt->index[i] + count[i]);
274 ril_mt->il.resize(ril_mt->index[nat_mt]);
276 /* Store the interactions */
277 low_make_reverse_ilist(
278 ilist, atoms->atom, count.data(), rtOptions, ril_mt->index, ril_mt->il, atomLinkRule, TRUE);
280 ril_mt->numAtomsInMolecule = atoms->nr;
283 gmx_reverse_top_t::gmx_reverse_top_t(const gmx_mtop_t& mtop,
285 const ReverseTopOptions& reverseTopOptions) :
286 impl_(std::make_unique<Impl>(mtop, useFreeEnergy, reverseTopOptions))
290 gmx_reverse_top_t::~gmx_reverse_top_t() {}
292 const ReverseTopOptions& gmx_reverse_top_t::options() const
294 return impl_->options;
297 const reverse_ilist_t& gmx_reverse_top_t::interactionListForMoleculeType(int moleculeType) const
299 return impl_->ril_mt[moleculeType];
302 int gmx_reverse_top_t::expectedNumGlobalBondedInteractions() const
304 return impl_->expectedNumGlobalBondedInteractions;
307 ArrayRef<const MolblockIndices> gmx_reverse_top_t::molblockIndices() const
312 bool gmx_reverse_top_t::hasIntermolecularInteractions() const
314 return impl_->bIntermolecularInteractions;
317 const reverse_ilist_t& gmx_reverse_top_t::interactionListForIntermolecularInteractions() const
319 return impl_->ril_intermol;
322 bool gmx_reverse_top_t::hasInterAtomicInteractions() const
324 return impl_->bInterAtomicInteractions;
327 bool gmx_reverse_top_t::hasPositionRestraints() const
329 return impl_->hasPositionRestraints;
332 ArrayRef<thread_work_t> gmx_reverse_top_t::threadWorkObjects() const
334 return impl_->th_work;
337 bool gmx_reverse_top_t::doSorting() const
339 return impl_->ilsort != ilsortNO_FE;
342 /*! \brief Compute the total bonded interaction count
344 * \param[in] mtop The global system topology
345 * \param[in] includeConstraints Whether constraint interactions are included in the count
346 * \param[in] includeSettles Whether settle interactions are included in the count
348 * When using domain decomposition without update groups,
349 * constraint-type interactions can be split across domains, and so we
350 * do not consider them in this correctness check. Otherwise, we
353 static int computeExpectedNumGlobalBondedInteractions(const gmx_mtop_t& mtop,
354 const bool includeConstraints,
355 const bool includeSettles)
357 int expectedNumGlobalBondedInteractions = gmx_mtop_interaction_count(mtop, IF_BOND);
358 if (includeConstraints)
360 expectedNumGlobalBondedInteractions +=
361 gmx_mtop_ftype_count(mtop, F_CONSTR) + gmx_mtop_ftype_count(mtop, F_CONSTRNC);
365 expectedNumGlobalBondedInteractions += gmx_mtop_ftype_count(mtop, F_SETTLE);
367 return expectedNumGlobalBondedInteractions;
370 /*! \brief Generate the reverse topology */
371 gmx_reverse_top_t::Impl::Impl(const gmx_mtop_t& mtop,
372 const bool useFreeEnergy,
373 const ReverseTopOptions& reverseTopOptions) :
374 options(reverseTopOptions),
375 hasPositionRestraints(gmx_mtop_ftype_count(mtop, F_POSRES) + gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0),
376 bInterAtomicInteractions(mtop.bIntermolecularInteractions),
377 expectedNumGlobalBondedInteractions(
378 computeExpectedNumGlobalBondedInteractions(mtop,
379 reverseTopOptions.includeConstraints,
380 reverseTopOptions.includeSettles))
382 bInterAtomicInteractions = mtop.bIntermolecularInteractions;
383 ril_mt.resize(mtop.moltype.size());
385 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
387 const gmx_moltype_t& molt = mtop.moltype[mt];
388 if (molt.atoms.nr > 1)
390 bInterAtomicInteractions = true;
393 /* Make the atom to interaction list for this molecule type */
394 make_reverse_ilist(molt.ilist, &molt.atoms, options, AtomLinkRule::FirstAtom, &ril_mt[mt]);
396 ril_mt_tot_size += ril_mt[mt].index[molt.atoms.nr];
400 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", ril_mt_tot_size);
403 /* Make an intermolecular reverse top, if necessary */
404 bIntermolecularInteractions = mtop.bIntermolecularInteractions;
405 if (bIntermolecularInteractions)
407 t_atoms atoms_global;
409 atoms_global.nr = mtop.natoms;
410 atoms_global.atom = nullptr; /* Only used with virtual sites */
412 GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
413 "We should have an ilist when intermolecular interactions are on");
416 *mtop.intermolecular_ilist, &atoms_global, options, AtomLinkRule::FirstAtom, &ril_intermol);
419 if (useFreeEnergy && gmx_mtop_bondeds_free_energy(&mtop))
421 ilsort = ilsortFE_UNSORTED;
425 ilsort = ilsortNO_FE;
428 /* Make a molblock index for fast searching */
430 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
432 const gmx_molblock_t& molb = mtop.molblock[mb];
433 const int numAtomsPerMol = mtop.moltype[molb.type].atoms.nr;
434 MolblockIndices mbiMolblock;
435 mbiMolblock.a_start = i;
436 i += molb.nmol * numAtomsPerMol;
437 mbiMolblock.a_end = i;
438 mbiMolblock.natoms_mol = numAtomsPerMol;
439 mbiMolblock.type = molb.type;
440 mbi.push_back(mbiMolblock);
443 for (int th = 0; th < gmx_omp_nthreads_get(ModuleMultiThread::Domdec); th++)
445 th_work.emplace_back(mtop.ffparams);
449 void dd_make_reverse_top(FILE* fplog,
451 const gmx_mtop_t& mtop,
452 const gmx::VirtualSitesHandler* vsite,
453 const t_inputrec& inputrec,
454 const DDBondedChecking ddBondedChecking)
458 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
461 /* If normal and/or settle constraints act only within charge groups,
462 * we can store them in the reverse top and simply assign them to domains.
463 * Otherwise we need to assign them to multiple domains and set up
464 * the parallel version constraint algorithm(s).
466 GMX_RELEASE_ASSERT(ddBondedChecking == DDBondedChecking::ExcludeZeroLimit
467 || ddBondedChecking == DDBondedChecking::All,
468 "Invalid enum value for mdrun -ddcheck");
469 const ReverseTopOptions rtOptions(ddBondedChecking,
470 !dd->comm->systemInfo.mayHaveSplitConstraints,
471 !dd->comm->systemInfo.mayHaveSplitSettles);
473 dd->reverse_top = std::make_unique<gmx_reverse_top_t>(
474 mtop, inputrec.efep != FreeEnergyPerturbationType::No, rtOptions);
476 dd->haveExclusions = false;
477 for (const gmx_molblock_t& molb : mtop.molblock)
479 const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop.moltype[molb.type].excls);
480 // We checked above that max 1 exclusion means only self exclusions
481 if (maxNumExclusionsPerAtom > 1)
483 dd->haveExclusions = true;
487 const int numInterUpdategroupVirtualSites =
488 (vsite == nullptr ? 0 : vsite->numInterUpdategroupVirtualSites());
489 if (numInterUpdategroupVirtualSites > 0)
494 "There are %d inter update-group virtual sites,\n"
495 "will perform an extra communication step for selected coordinates and "
497 numInterUpdategroupVirtualSites);
499 init_domdec_vsites(dd, numInterUpdategroupVirtualSites);
502 if (dd->comm->systemInfo.mayHaveSplitConstraints || dd->comm->systemInfo.mayHaveSplitSettles)
504 init_domdec_constraints(dd, mtop);
508 fprintf(fplog, "\n");