2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006 - 2014, The GROMACS development team.
5 * Copyright (c) 2015,2016,2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines functions used by the domdec module
40 * while managing the construction, use and error checking for
41 * topologies local to a DD rank.
43 * \author Berk Hess <hess@kth.se>
44 * \ingroup module_domdec
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/forcerec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/vsite.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/ifunc.h"
74 #include "gromacs/topology/mtop_util.h"
75 #include "gromacs/topology/topsort.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/exceptions.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/gmxassert.h"
80 #include "gromacs/utility/listoflists.h"
81 #include "gromacs/utility/logger.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/utility/strconvert.h"
84 #include "gromacs/utility/stringstream.h"
85 #include "gromacs/utility/stringutil.h"
86 #include "gromacs/utility/textwriter.h"
88 #include "domdec_constraints.h"
89 #include "domdec_internal.h"
90 #include "domdec_vsite.h"
93 using gmx::ListOfLists;
95 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
96 #define NITEM_DD_INIT_LOCAL_STATE 5
98 struct reverse_ilist_t
100 std::vector<int> index; /* Index for each atom into il */
101 std::vector<int> il; /* ftype|type|a0|...|an|ftype|... */
102 int numAtomsInMolecule; /* The number of atoms in this molecule */
105 struct MolblockIndices
113 /*! \brief Struct for thread local work data for local topology generation */
116 t_idef idef; /**< Partial local topology */
117 std::unique_ptr<VsitePbc> vsitePbc; /**< vsite PBC structure */
118 int nbonded; /**< The number of bondeds in this struct */
119 ListOfLists<int> excl; /**< List of exclusions */
120 int excl_count; /**< The total exclusion count for \p excl */
123 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
124 struct gmx_reverse_top_t
126 //! @cond Doxygen_Suppress
127 //! \brief Are there constraints in this revserse top?
128 bool bConstr = false;
129 //! \brief Are there settles in this revserse top?
130 bool bSettle = false;
131 //! \brief All bonded interactions have to be assigned?
132 bool bBCheck = false;
133 //! \brief Are there bondeds/exclusions between atoms?
134 bool bInterAtomicInteractions = false;
135 //! \brief Reverse ilist for all moltypes
136 std::vector<reverse_ilist_t> ril_mt;
137 //! \brief The size of ril_mt[?].index summed over all entries
138 int ril_mt_tot_size = 0;
139 //! \brief The sorting state of bondeds for free energy
140 int ilsort = ilsortUNKNOWN;
141 //! \brief molblock to global atom index for quick lookup of molblocks on atom index
142 std::vector<MolblockIndices> mbi;
144 //! \brief Do we have intermolecular interactions?
145 bool bIntermolecularInteractions = false;
146 //! \brief Intermolecular reverse ilist
147 reverse_ilist_t ril_intermol;
149 /* Work data structures for multi-threading */
150 //! \brief Thread work array for local topology generation
151 std::vector<thread_work_t> th_work;
155 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
156 static int nral_rt(int ftype)
161 if (interaction_function[ftype].flags & IF_VSITE)
163 /* With vsites the reverse topology contains an extra entry
164 * for storing if constructing atoms are vsites.
172 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
173 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck, gmx_bool bConstr, gmx_bool bSettle)
175 return ((((interaction_function[ftype].flags & IF_BOND) != 0U)
176 && ((interaction_function[ftype].flags & IF_VSITE) == 0U)
177 && (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U)))
178 || (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) || (bSettle && ftype == F_SETTLE));
181 /*! \brief Help print error output when interactions are missing */
182 static std::string print_missing_interactions_mb(t_commrec* cr,
183 const gmx_reverse_top_t* rt,
184 const char* moltypename,
185 const reverse_ilist_t* ril,
193 int nril_mol = ril->index[nat_mol];
194 snew(assigned, nmol * nril_mol);
195 gmx::StringOutputStream stream;
196 gmx::TextWriter log(&stream);
198 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
199 for (int ftype = 0; ftype < F_NRE; ftype++)
201 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
203 int nral = NRAL(ftype);
204 const t_ilist* il = &idef->il[ftype];
205 const t_iatom* ia = il->iatoms;
206 for (int i = 0; i < il->nr; i += 1 + nral)
208 int a0 = gatindex[ia[1]];
209 /* Check if this interaction is in
210 * the currently checked molblock.
212 if (a0 >= a_start && a0 < a_end)
214 int mol = (a0 - a_start) / nat_mol;
215 int a0_mol = (a0 - a_start) - mol * nat_mol;
216 int j_mol = ril->index[a0_mol];
218 while (j_mol < ril->index[a0_mol + 1] && !found)
220 int j = mol * nril_mol + j_mol;
221 int ftype_j = ril->il[j_mol];
222 /* Here we need to check if this interaction has
223 * not already been assigned, since we could have
224 * multiply defined interactions.
226 if (ftype == ftype_j && ia[0] == ril->il[j_mol + 1] && assigned[j] == 0)
228 /* Check the atoms */
230 for (int a = 0; a < nral; a++)
232 if (gatindex[ia[1 + a]] != a_start + mol * nat_mol + ril->il[j_mol + 2 + a])
242 j_mol += 2 + nral_rt(ftype_j);
246 gmx_incons("Some interactions seem to be assigned multiple times");
254 gmx_sumi(nmol * nril_mol, assigned, cr);
258 for (int mol = 0; mol < nmol; mol++)
261 while (j_mol < nril_mol)
263 int ftype = ril->il[j_mol];
264 int nral = NRAL(ftype);
265 int j = mol * nril_mol + j_mol;
266 if (assigned[j] == 0 && !(interaction_function[ftype].flags & IF_VSITE))
268 if (DDMASTER(cr->dd))
272 log.writeLineFormatted("Molecule type '%s'", moltypename);
273 log.writeLineFormatted(
274 "the first %d missing interactions, except for exclusions:", nprint);
276 log.writeStringFormatted("%20s atoms", interaction_function[ftype].longname);
278 for (a = 0; a < nral; a++)
280 log.writeStringFormatted("%5d", ril->il[j_mol + 2 + a] + 1);
284 log.writeString(" ");
287 log.writeString(" global");
288 for (a = 0; a < nral; a++)
290 log.writeStringFormatted(
291 "%6d", a_start + mol * nat_mol + ril->il[j_mol + 2 + a] + 1);
293 log.ensureLineBreak();
301 j_mol += 2 + nral_rt(ftype);
306 return stream.toString();
309 /*! \brief Help print error output when interactions are missing */
310 static void print_missing_interactions_atoms(const gmx::MDLogger& mdlog,
312 const gmx_mtop_t* mtop,
315 const gmx_reverse_top_t* rt = cr->dd->reverse_top;
317 /* Print the atoms in the missing interactions per molblock */
319 for (const gmx_molblock_t& molb : mtop->molblock)
321 const gmx_moltype_t& moltype = mtop->moltype[molb.type];
323 a_end = a_start + molb.nmol * moltype.atoms.nr;
325 GMX_LOG(mdlog.warning)
326 .appendText(print_missing_interactions_mb(cr, rt, *(moltype.name),
327 &rt->ril_mt[molb.type], a_start, a_end,
328 moltype.atoms.nr, molb.nmol, idef));
332 void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
335 const gmx_mtop_t* top_global,
336 const gmx_localtop_t* top_local,
337 gmx::ArrayRef<const gmx::RVec> x,
340 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
346 GMX_LOG(mdlog.warning)
348 "Not all bonded interactions have been properly assigned to the domain "
349 "decomposition cells");
351 ndiff_tot = local_count - dd->nbonded_global;
353 for (ftype = 0; ftype < F_NRE; ftype++)
356 cl[ftype] = top_local->idef.il[ftype].nr / (1 + nral);
359 gmx_sumi(F_NRE, cl, cr);
363 GMX_LOG(mdlog.warning).appendText("A list of missing interactions:");
364 rest_global = dd->nbonded_global;
365 rest_local = local_count;
366 for (ftype = 0; ftype < F_NRE; ftype++)
368 /* In the reverse and local top all constraints are merged
369 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
370 * and add these constraints when doing F_CONSTR.
372 if (((interaction_function[ftype].flags & IF_BOND)
373 && (dd->reverse_top->bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO)))
374 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
375 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
377 n = gmx_mtop_ftype_count(top_global, ftype);
378 if (ftype == F_CONSTR)
380 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
382 ndiff = cl[ftype] - n;
385 GMX_LOG(mdlog.warning)
386 .appendTextFormatted("%20s of %6d missing %6d",
387 interaction_function[ftype].longname, n, -ndiff);
390 rest_local -= cl[ftype];
394 ndiff = rest_local - rest_global;
397 GMX_LOG(mdlog.warning).appendTextFormatted("%20s of %6d missing %6d", "exclusions", rest_global, -ndiff);
401 print_missing_interactions_atoms(mdlog, cr, top_global, &top_local->idef);
402 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr, -1, as_rvec_array(x.data()), box);
404 std::string errorMessage;
409 "One or more interactions were assigned to multiple domains of the domain "
410 "decompostion. Please report this bug.";
414 errorMessage = gmx::formatString(
415 "%d of the %d bonded interactions could not be calculated because some atoms "
416 "involved moved further apart than the multi-body cut-off distance (%g nm) or the "
417 "two-body cut-off distance (%g nm), see option -rdd, for pairs and tabulated bonds "
418 "also see option -ddcheck",
419 -ndiff_tot, cr->dd->nbonded_global, dd_cutoff_multibody(dd), dd_cutoff_twobody(dd));
421 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
424 /*! \brief Return global topology molecule information for global atom index \p i_gl */
425 static void global_atomnr_to_moltype_ind(const gmx_reverse_top_t* rt, int i_gl, int* mb, int* mt, int* mol, int* i_mol)
427 const MolblockIndices* mbi = rt->mbi.data();
429 int end = rt->mbi.size(); /* exclusive */
432 /* binary search for molblock_ind */
435 mid = (start + end) / 2;
436 if (i_gl >= mbi[mid].a_end)
440 else if (i_gl < mbi[mid].a_start)
454 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
455 *i_mol = (i_gl - mbi->a_start) - (*mol) * mbi->natoms_mol;
458 /*! \brief Returns the maximum number of exclusions per atom */
459 static int getMaxNumExclusionsPerAtom(const ListOfLists<int>& excls)
462 for (gmx::index at = 0; at < excls.ssize(); at++)
464 const auto list = excls[at];
465 const int numExcls = list.ssize();
467 GMX_RELEASE_ASSERT(numExcls != 1 || list[0] == at,
468 "With 1 exclusion we expect a self-exclusion");
470 maxNumExcls = std::max(maxNumExcls, numExcls);
476 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
477 static int low_make_reverse_ilist(const InteractionLists& il_mt,
483 gmx::ArrayRef<const int> r_index,
484 gmx::ArrayRef<int> r_il,
485 gmx_bool bLinkToAllAtoms,
488 int ftype, j, nlink, link;
493 for (ftype = 0; ftype < F_NRE; ftype++)
495 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
496 || (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) || (bSettle && ftype == F_SETTLE))
498 const bool bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
499 const int nral = NRAL(ftype);
500 const auto& il = il_mt[ftype];
501 for (int i = 0; i < il.size(); i += 1 + nral)
503 const int* ia = il.iatoms.data() + i;
508 /* We don't need the virtual sites for the cg-links */
518 /* Couple to the first atom in the interaction */
521 for (link = 0; link < nlink; link++)
526 GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
527 GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
528 r_il[r_index[a] + count[a]] = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
529 r_il[r_index[a] + count[a] + 1] = ia[0];
530 for (j = 1; j < 1 + nral; j++)
532 /* Store the molecular atom number */
533 r_il[r_index[a] + count[a] + 1 + j] = ia[j];
536 if (interaction_function[ftype].flags & IF_VSITE)
540 /* Add an entry to iatoms for storing
541 * which of the constructing atoms are
544 r_il[r_index[a] + count[a] + 2 + nral] = 0;
545 for (j = 2; j < 1 + nral; j++)
547 if (atom[ia[j]].ptype == eptVSite)
549 r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
556 /* We do not count vsites since they are always
557 * uniquely assigned and can be assigned
558 * to multiple nodes with recursive vsites.
560 if (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))
565 count[a] += 2 + nral_rt(ftype);
574 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
575 static int make_reverse_ilist(const InteractionLists& ilist,
576 const t_atoms* atoms,
580 gmx_bool bLinkToAllAtoms,
581 reverse_ilist_t* ril_mt)
583 int nat_mt, *count, i, nint_mt;
585 /* Count the interactions */
588 low_make_reverse_ilist(ilist, atoms->atom, count, bConstr, bSettle, bBCheck, {}, {},
589 bLinkToAllAtoms, FALSE);
591 ril_mt->index.push_back(0);
592 for (i = 0; i < nat_mt; i++)
594 ril_mt->index.push_back(ril_mt->index[i] + count[i]);
597 ril_mt->il.resize(ril_mt->index[nat_mt]);
599 /* Store the interactions */
600 nint_mt = low_make_reverse_ilist(ilist, atoms->atom, count, bConstr, bSettle, bBCheck,
601 ril_mt->index, ril_mt->il, bLinkToAllAtoms, TRUE);
605 ril_mt->numAtomsInMolecule = atoms->nr;
610 /*! \brief Generate the reverse topology */
611 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t* mtop,
618 gmx_reverse_top_t rt;
620 /* Should we include constraints (for SHAKE) in rt? */
621 rt.bConstr = bConstr;
622 rt.bSettle = bSettle;
623 rt.bBCheck = bBCheck;
625 rt.bInterAtomicInteractions = mtop->bIntermolecularInteractions;
626 rt.ril_mt.resize(mtop->moltype.size());
627 rt.ril_mt_tot_size = 0;
628 std::vector<int> nint_mt;
629 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
631 const gmx_moltype_t& molt = mtop->moltype[mt];
632 if (molt.atoms.nr > 1)
634 rt.bInterAtomicInteractions = true;
637 /* Make the atom to interaction list for this molecule type */
638 int numberOfInteractions = make_reverse_ilist(
639 molt.ilist, &molt.atoms, rt.bConstr, rt.bSettle, rt.bBCheck, FALSE, &rt.ril_mt[mt]);
640 nint_mt.push_back(numberOfInteractions);
642 rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
646 fprintf(debug, "The total size of the atom to interaction index is %d integers\n",
651 for (const gmx_molblock_t& molblock : mtop->molblock)
653 *nint += molblock.nmol * nint_mt[molblock.type];
656 /* Make an intermolecular reverse top, if necessary */
657 rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
658 if (rt.bIntermolecularInteractions)
660 t_atoms atoms_global;
662 atoms_global.nr = mtop->natoms;
663 atoms_global.atom = nullptr; /* Only used with virtual sites */
665 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist,
666 "We should have an ilist when intermolecular interactions are on");
668 *nint += make_reverse_ilist(*mtop->intermolecular_ilist, &atoms_global, rt.bConstr,
669 rt.bSettle, rt.bBCheck, FALSE, &rt.ril_intermol);
672 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
674 rt.ilsort = ilsortFE_UNSORTED;
678 rt.ilsort = ilsortNO_FE;
681 /* Make a molblock index for fast searching */
683 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
685 const gmx_molblock_t& molb = mtop->molblock[mb];
686 const int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
689 i += molb.nmol * numAtomsPerMol;
691 mbi.natoms_mol = numAtomsPerMol;
692 mbi.type = molb.type;
693 rt.mbi.push_back(mbi);
696 rt.th_work.resize(gmx_omp_nthreads_get(emntDomdec));
701 void dd_make_reverse_top(FILE* fplog,
703 const gmx_mtop_t* mtop,
704 const gmx_vsite_t* vsite,
705 const t_inputrec* ir,
710 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
713 /* If normal and/or settle constraints act only within charge groups,
714 * we can store them in the reverse top and simply assign them to domains.
715 * Otherwise we need to assign them to multiple domains and set up
716 * the parallel version constraint algorithm(s).
719 dd->reverse_top = new gmx_reverse_top_t;
721 make_reverse_top(mtop, ir->efep != efepNO, !dd->comm->systemInfo.haveSplitConstraints,
722 !dd->comm->systemInfo.haveSplitSettles, bBCheck, &dd->nbonded_global);
724 dd->haveExclusions = false;
725 for (const gmx_molblock_t& molb : mtop->molblock)
727 const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop->moltype[molb.type].excls);
728 // We checked above that max 1 exclusion means only self exclusions
729 if (maxNumExclusionsPerAtom > 1)
731 dd->haveExclusions = true;
735 if (vsite && vsite->numInterUpdategroupVsites > 0)
740 "There are %d inter update-group virtual sites,\n"
741 "will an extra communication step for selected coordinates and forces\n",
742 vsite->numInterUpdategroupVsites);
744 init_domdec_vsites(dd, vsite->numInterUpdategroupVsites);
747 if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
749 init_domdec_constraints(dd, mtop);
753 fprintf(fplog, "\n");
757 /*! \brief Store a vsite interaction at the end of \p il
759 * This routine is very similar to add_ifunc, but vsites interactions
760 * have more work to do than other kinds of interactions, and the
761 * complex way nral (and thus vector contents) depends on ftype
762 * confuses static analysis tools unless we fuse the vsite
763 * atom-indexing organization code with the ifunc-adding code, so that
764 * they can see that nral is the same value. */
765 static inline void add_ifunc_for_vsites(t_iatom* tiatoms,
766 const gmx_ga2la_t& ga2la,
772 const t_iatom* iatoms,
777 if (il->nr + 1 + nral > il->nalloc)
779 il->nalloc = over_alloc_large(il->nr + 1 + nral);
780 srenew(il->iatoms, il->nalloc);
782 liatoms = il->iatoms + il->nr;
786 tiatoms[0] = iatoms[0];
790 /* We know the local index of the first atom */
795 /* Convert later in make_local_vsites */
796 tiatoms[1] = -a_gl - 1;
799 for (int k = 2; k < 1 + nral; k++)
801 int ak_gl = a_gl + iatoms[k] - a_mol;
802 if (const int* homeIndex = ga2la.findHome(ak_gl))
804 tiatoms[k] = *homeIndex;
808 /* Copy the global index, convert later in make_local_vsites */
809 tiatoms[k] = -(ak_gl + 1);
811 // Note that ga2la_get_home always sets the third parameter if
814 for (int k = 0; k < 1 + nral; k++)
816 liatoms[k] = tiatoms[k];
820 /*! \brief Store a bonded interaction at the end of \p il */
821 static inline void add_ifunc(int nral, const t_iatom* tiatoms, t_ilist* il)
826 if (il->nr + 1 + nral > il->nalloc)
828 il->nalloc = over_alloc_large(il->nr + 1 + nral);
829 srenew(il->iatoms, il->nalloc);
831 liatoms = il->iatoms + il->nr;
832 for (k = 0; k <= nral; k++)
834 liatoms[k] = tiatoms[k];
839 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
840 static void add_posres(int mol,
842 int numAtomsInMolecule,
843 const gmx_molblock_t* molb,
845 const t_iparams* ip_in,
851 /* This position restraint has not been added yet,
852 * so it's index is the current number of position restraints.
854 n = idef->il[F_POSRES].nr / 2;
855 if (n + 1 > idef->iparams_posres_nalloc)
857 idef->iparams_posres_nalloc = over_alloc_dd(n + 1);
858 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
860 ip = &idef->iparams_posres[n];
861 /* Copy the force constants */
862 *ip = ip_in[iatoms[0]];
864 /* Get the position restraint coordinates from the molblock */
865 a_molb = mol * numAtomsInMolecule + a_mol;
866 GMX_ASSERT(a_molb < ssize(molb->posres_xA),
867 "We need a sufficient number of position restraint coordinates");
868 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
869 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
870 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
871 if (!molb->posres_xB.empty())
873 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
874 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
875 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
879 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
880 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
881 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
883 /* Set the parameter index for idef->iparams_posre */
887 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
888 static void add_fbposres(int mol,
890 int numAtomsInMolecule,
891 const gmx_molblock_t* molb,
893 const t_iparams* ip_in,
899 /* This flat-bottom position restraint has not been added yet,
900 * so it's index is the current number of position restraints.
902 n = idef->il[F_FBPOSRES].nr / 2;
903 if (n + 1 > idef->iparams_fbposres_nalloc)
905 idef->iparams_fbposres_nalloc = over_alloc_dd(n + 1);
906 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
908 ip = &idef->iparams_fbposres[n];
909 /* Copy the force constants */
910 *ip = ip_in[iatoms[0]];
912 /* Get the position restraint coordinats from the molblock */
913 a_molb = mol * numAtomsInMolecule + a_mol;
914 GMX_ASSERT(a_molb < ssize(molb->posres_xA),
915 "We need a sufficient number of position restraint coordinates");
916 /* Take reference positions from A position of normal posres */
917 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
918 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
919 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
921 /* Note: no B-type for flat-bottom posres */
923 /* Set the parameter index for idef->iparams_posre */
927 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
928 static void add_vsite(const gmx_ga2la_t& ga2la,
929 gmx::ArrayRef<const int> index,
930 gmx::ArrayRef<const int> rtil,
937 const t_iatom* iatoms,
941 t_iatom tiatoms[1 + MAXATOMLIST];
942 int j, ftype_r, nral_r;
944 /* Add this interaction to the local topology */
945 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
947 if (iatoms[1 + nral])
949 /* Check for recursion */
950 for (k = 2; k < 1 + nral; k++)
952 if ((iatoms[1 + nral] & (2 << k)) && (tiatoms[k] < 0))
954 /* This construction atoms is a vsite and not a home atom */
957 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n",
958 iatoms[k] + 1, a_mol + 1);
960 /* Find the vsite construction */
962 /* Check all interactions assigned to this atom */
963 j = index[iatoms[k]];
964 while (j < index[iatoms[k] + 1])
967 nral_r = NRAL(ftype_r);
968 if (interaction_function[ftype_r].flags & IF_VSITE)
970 /* Add this vsite (recursion) */
971 add_vsite(ga2la, index, rtil, ftype_r, nral_r, FALSE, -1,
972 a_gl + iatoms[k] - iatoms[1], iatoms[k], rtil.data() + j, idef);
974 j += 1 + nral_rt(ftype_r);
981 /*! \brief Returns the squared distance between atoms \p i and \p j */
982 static real dd_dist2(t_pbc* pbc_null, const rvec* x, const int i, int j)
988 pbc_dx_aiuc(pbc_null, x[i], x[j], dx);
992 rvec_sub(x[i], x[j], dx);
998 /*! \brief Append t_idef structures 1 to nsrc in src to *dest */
999 static void combine_idef(t_idef* dest, gmx::ArrayRef<const thread_work_t> src)
1003 for (ftype = 0; ftype < F_NRE; ftype++)
1006 for (gmx::index s = 1; s < src.ssize(); s++)
1008 n += src[s].idef.il[ftype].nr;
1014 ild = &dest->il[ftype];
1016 if (ild->nr + n > ild->nalloc)
1018 ild->nalloc = over_alloc_large(ild->nr + n);
1019 srenew(ild->iatoms, ild->nalloc);
1022 for (gmx::index s = 1; s < src.ssize(); s++)
1024 const t_ilist& ils = src[s].idef.il[ftype];
1026 for (int i = 0; i < ils.nr; i++)
1028 ild->iatoms[ild->nr + i] = ils.iatoms[i];
1034 /* Position restraints need an additional treatment */
1035 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1037 int nposres = dest->il[ftype].nr / 2;
1038 // TODO: Simplify this code using std::vector
1039 t_iparams*& iparams_dest =
1040 (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1041 int& posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc
1042 : dest->iparams_fbposres_nalloc);
1043 if (nposres > posres_nalloc)
1045 posres_nalloc = over_alloc_large(nposres);
1046 srenew(iparams_dest, posres_nalloc);
1049 /* Set nposres to the number of original position restraints in dest */
1050 for (gmx::index s = 1; s < src.ssize(); s++)
1052 nposres -= src[s].idef.il[ftype].nr / 2;
1055 for (gmx::index s = 1; s < src.ssize(); s++)
1057 const t_iparams* iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres
1058 : src[s].idef.iparams_fbposres);
1060 for (int i = 0; i < src[s].idef.il[ftype].nr / 2; i++)
1062 /* Correct the index into iparams_posres */
1063 dest->il[ftype].iatoms[nposres * 2] = nposres;
1064 /* Copy the position restraint force parameters */
1065 iparams_dest[nposres] = iparams_src[i];
1074 /*! \brief Check and when available assign bonded interactions for local atom i
1076 static inline void check_assign_interactions_atom(int i,
1080 int numAtomsInMolecule,
1081 gmx::ArrayRef<const int> index,
1082 gmx::ArrayRef<const int> rtil,
1083 gmx_bool bInterMolInteractions,
1086 const gmx_domdec_t* dd,
1087 const gmx_domdec_zones_t* zones,
1088 const gmx_molblock_t* molb,
1095 const t_iparams* ip_in,
1101 gmx::ArrayRef<const DDPairInteractionRanges> iZones = zones->iZones;
1106 t_iatom tiatoms[1 + MAXATOMLIST];
1108 const int ftype = rtil[j++];
1109 auto iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1110 const int nral = NRAL(ftype);
1111 if (interaction_function[ftype].flags & IF_VSITE)
1113 assert(!bInterMolInteractions);
1114 /* The vsite construction goes where the vsite itself is */
1117 add_vsite(*dd->ga2la, index, rtil, ftype, nral, TRUE, i, i_gl, i_mol, iatoms.data(), idef);
1126 tiatoms[0] = iatoms[0];
1130 assert(!bInterMolInteractions);
1131 /* Assign single-body interactions to the home zone */
1136 if (ftype == F_POSRES)
1138 add_posres(mol, i_mol, numAtomsInMolecule, molb, tiatoms, ip_in, idef);
1140 else if (ftype == F_FBPOSRES)
1142 add_fbposres(mol, i_mol, numAtomsInMolecule, molb, tiatoms, ip_in, idef);
1152 /* This is a two-body interaction, we can assign
1153 * analogous to the non-bonded assignments.
1157 if (!bInterMolInteractions)
1159 /* Get the global index using the offset in the molecule */
1160 k_gl = i_gl + iatoms[2] - i_mol;
1166 if (const auto* entry = dd->ga2la->find(k_gl))
1168 int kz = entry->cell;
1173 /* Check zone interaction assignments */
1174 bUse = ((iz < iZones.ssize() && iz <= kz && iZones[iz].jZoneRange.isInRange(kz))
1175 || (kz < iZones.ssize() && iz > kz && iZones[kz].jZoneRange.isInRange(iz)));
1178 GMX_ASSERT(ftype != F_CONSTR || (iz == 0 && kz == 0),
1179 "Constraint assigned here should only involve home atoms");
1182 tiatoms[2] = entry->la;
1183 /* If necessary check the cgcm distance */
1184 if (bRCheck2B && dd_dist2(pbc_null, cg_cm, tiatoms[1], tiatoms[2]) >= rc2)
1197 /* Assign this multi-body bonded interaction to
1198 * the local node if we have all the atoms involved
1199 * (local or communicated) and the minimum zone shift
1200 * in each dimension is zero, for dimensions
1201 * with 2 DD cells an extra check may be necessary.
1203 ivec k_zero, k_plus;
1209 for (k = 1; k <= nral && bUse; k++)
1212 if (!bInterMolInteractions)
1214 /* Get the global index using the offset in the molecule */
1215 k_gl = i_gl + iatoms[k] - i_mol;
1221 const auto* entry = dd->ga2la->find(k_gl);
1222 if (entry == nullptr || entry->cell >= zones->n)
1224 /* We do not have this atom of this interaction
1225 * locally, or it comes from more than one cell
1234 tiatoms[k] = entry->la;
1235 for (d = 0; d < DIM; d++)
1237 if (zones->shift[entry->cell][d] == 0)
1248 bUse = (bUse && (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1253 for (d = 0; (d < DIM && bUse); d++)
1255 /* Check if the cg_cm distance falls within
1256 * the cut-off to avoid possible multiple
1257 * assignments of bonded interactions.
1259 if (rcheck[d] && k_plus[d]
1260 && dd_dist2(pbc_null, cg_cm, tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1269 /* Add this interaction to the local topology */
1270 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1271 /* Sum so we can check in global_stat
1272 * if we have everything.
1274 if (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))
1284 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1286 * With thread parallelizing each thread acts on a different atom range:
1287 * at_start to at_end.
1289 static int make_bondeds_zone(gmx_domdec_t* dd,
1290 const gmx_domdec_zones_t* zones,
1291 const std::vector<gmx_molblock_t>& molb,
1298 const t_iparams* ip_in,
1301 const gmx::Range<int>& atomRange)
1303 int mb, mt, mol, i_mol;
1305 gmx_reverse_top_t* rt;
1308 rt = dd->reverse_top;
1310 bBCheck = rt->bBCheck;
1314 for (int i : atomRange)
1316 /* Get the global atom number */
1317 const int i_gl = dd->globalAtomIndices[i];
1318 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1319 /* Check all intramolecular interactions assigned to this atom */
1320 gmx::ArrayRef<const int> index = rt->ril_mt[mt].index;
1321 gmx::ArrayRef<const t_iatom> rtil = rt->ril_mt[mt].il;
1323 check_assign_interactions_atom(i, i_gl, mol, i_mol, rt->ril_mt[mt].numAtomsInMolecule,
1324 index, rtil, FALSE, index[i_mol], index[i_mol + 1], dd,
1325 zones, &molb[mb], bRCheckMB, rcheck, bRCheck2B, rc2,
1326 pbc_null, cg_cm, ip_in, idef, izone, bBCheck, &nbonded_local);
1329 if (rt->bIntermolecularInteractions)
1331 /* Check all intermolecular interactions assigned to this atom */
1332 index = rt->ril_intermol.index;
1333 rtil = rt->ril_intermol.il;
1335 check_assign_interactions_atom(i, i_gl, mol, i_mol, rt->ril_mt[mt].numAtomsInMolecule,
1336 index, rtil, TRUE, index[i_gl], index[i_gl + 1], dd, zones,
1337 &molb[mb], bRCheckMB, rcheck, bRCheck2B, rc2, pbc_null,
1338 cg_cm, ip_in, idef, izone, bBCheck, &nbonded_local);
1342 return nbonded_local;
1345 /*! \brief Set the exclusion data for i-zone \p iz */
1346 static void make_exclusions_zone(gmx_domdec_t* dd,
1347 gmx_domdec_zones_t* zones,
1348 const std::vector<gmx_moltype_t>& moltype,
1350 ListOfLists<int>* lexcls,
1354 const gmx::ArrayRef<const int> intermolecularExclusionGroup)
1356 const gmx_ga2la_t& ga2la = *dd->ga2la;
1358 const auto& jAtomRange = zones->iZones[iz].jAtomRange;
1360 const gmx::index oldNumLists = lexcls->ssize();
1362 std::vector<int> exclusionsForAtom;
1363 for (int at = at_start; at < at_end; at++)
1365 exclusionsForAtom.clear();
1367 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1369 int a_gl, mb, mt, mol, a_mol;
1371 /* Copy the exclusions from the global top */
1372 a_gl = dd->globalAtomIndices[at];
1373 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1374 const auto excls = moltype[mt].excls[a_mol];
1375 for (const int aj_mol : excls)
1377 if (const auto* jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1379 /* This check is not necessary, but it can reduce
1380 * the number of exclusions in the list, which in turn
1381 * can speed up the pair list construction a bit.
1383 if (jAtomRange.isInRange(jEntry->la))
1385 exclusionsForAtom.push_back(jEntry->la);
1391 bool isExcludedAtom = !intermolecularExclusionGroup.empty()
1392 && std::find(intermolecularExclusionGroup.begin(),
1393 intermolecularExclusionGroup.end(), dd->globalAtomIndices[at])
1394 != intermolecularExclusionGroup.end();
1398 for (int qmAtomGlobalIndex : intermolecularExclusionGroup)
1400 if (const auto* entry = dd->ga2la->find(qmAtomGlobalIndex))
1402 exclusionsForAtom.push_back(entry->la);
1407 /* Append the exclusions for this atom to the topology */
1408 lexcls->pushBack(exclusionsForAtom);
1412 lexcls->ssize() - oldNumLists == at_end - at_start,
1413 "The number of exclusion list should match the number of atoms in the range");
1416 /*! \brief Clear a t_idef struct */
1417 static void clear_idef(t_idef* idef)
1421 /* Clear the counts */
1422 for (ftype = 0; ftype < F_NRE; ftype++)
1424 idef->il[ftype].nr = 0;
1428 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1429 static int make_local_bondeds_excls(gmx_domdec_t* dd,
1430 gmx_domdec_zones_t* zones,
1431 const gmx_mtop_t* mtop,
1440 ListOfLists<int>* lexcls,
1447 gmx_reverse_top_t* rt;
1449 if (dd->reverse_top->bInterAtomicInteractions)
1451 nzone_bondeds = zones->n;
1455 /* Only single charge group (or atom) molecules, so interactions don't
1456 * cross zone boundaries and we only need to assign in the home zone.
1461 /* We only use exclusions from i-zones to i- and j-zones */
1462 const int numIZonesForExclusions = (dd->haveExclusions ? zones->iZones.size() : 0);
1464 rt = dd->reverse_top;
1468 /* Clear the counts */
1475 for (int izone = 0; izone < nzone_bondeds; izone++)
1477 cg0 = zones->cg_range[izone];
1478 cg1 = zones->cg_range[izone + 1];
1480 const int numThreads = rt->th_work.size();
1481 #pragma omp parallel for num_threads(numThreads) schedule(static)
1482 for (int thread = 0; thread < numThreads; thread++)
1489 cg0t = cg0 + ((cg1 - cg0) * thread) / numThreads;
1490 cg1t = cg0 + ((cg1 - cg0) * (thread + 1)) / numThreads;
1498 idef_t = &rt->th_work[thread].idef;
1502 rt->th_work[thread].nbonded = make_bondeds_zone(
1503 dd, zones, mtop->molblock, bRCheckMB, rcheck, bRCheck2B, rc2, pbc_null,
1504 cg_cm, idef->iparams, idef_t, izone, gmx::Range<int>(cg0t, cg1t));
1506 if (izone < numIZonesForExclusions)
1508 ListOfLists<int>* excl_t;
1511 // Thread 0 stores exclusions directly in the final storage
1516 // Threads > 0 store in temporary storage, starting at list index 0
1517 excl_t = &rt->th_work[thread].excl;
1521 /* No charge groups and no distance check required */
1522 make_exclusions_zone(dd, zones, mtop->moltype, cginfo, excl_t, izone, cg0t,
1523 cg1t, mtop->intermolecularExclusionGroup);
1526 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1529 if (rt->th_work.size() > 1)
1531 combine_idef(idef, rt->th_work);
1534 for (const thread_work_t& th_work : rt->th_work)
1536 nbonded_local += th_work.nbonded;
1539 if (izone < numIZonesForExclusions)
1541 for (std::size_t th = 1; th < rt->th_work.size(); th++)
1543 lexcls->appendListOfLists(rt->th_work[th].excl);
1545 for (const thread_work_t& th_work : rt->th_work)
1547 *excl_count += th_work.excl_count;
1554 fprintf(debug, "We have %d exclusions, check count %d\n", lexcls->numElements(), *excl_count);
1557 return nbonded_local;
1560 void dd_make_local_top(gmx_domdec_t* dd,
1561 gmx_domdec_zones_t* zones,
1568 const gmx_mtop_t& mtop,
1569 gmx_localtop_t* ltop)
1571 gmx_bool bRCheckMB, bRCheck2B;
1575 t_pbc pbc, *pbc_null = nullptr;
1579 fprintf(debug, "Making local topology\n");
1585 if (dd->reverse_top->bInterAtomicInteractions)
1587 /* We need to check to which cell bondeds should be assigned */
1588 rc = dd_cutoff_twobody(dd);
1591 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1594 /* Should we check cg_cm distances when assigning bonded interactions? */
1595 for (d = 0; d < DIM; d++)
1598 /* Only need to check for dimensions where the part of the box
1599 * that is not communicated is smaller than the cut-off.
1601 if (d < npbcdim && dd->numCells[d] > 1
1602 && (dd->numCells[d] - npulse[d]) * cellsize_min[d] < 2 * rc)
1604 if (dd->numCells[d] == 2)
1609 /* Check for interactions between two atoms,
1610 * where we can allow interactions up to the cut-off,
1611 * instead of up to the smallest cell dimension.
1617 fprintf(debug, "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n", d,
1618 cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
1621 if (bRCheckMB || bRCheck2B)
1625 pbc_null = set_pbc_dd(&pbc, fr->pbcType, dd->numCells, TRUE, box);
1634 dd->nbonded_local = make_local_bondeds_excls(dd, zones, &mtop, fr->cginfo.data(), bRCheckMB,
1635 rcheck, bRCheck2B, rc, pbc_null, cgcm_or_x,
1636 <op->idef, <op->excls, &nexcl);
1638 /* The ilist is not sorted yet,
1639 * we can only do this when we have the charge arrays.
1641 ltop->idef.ilsort = ilsortUNKNOWN;
1643 ltop->atomtypes = mtop.atomtypes;
1646 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop)
1648 if (dd->reverse_top->ilsort == ilsortNO_FE)
1650 ltop->idef.ilsort = ilsortNO_FE;
1654 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
1658 void dd_init_local_top(const gmx_mtop_t& top_global, gmx_localtop_t* top)
1660 /* TODO: Get rid of the const casts below, e.g. by using a reference */
1661 top->idef.ntypes = top_global.ffparams.numTypes();
1662 top->idef.atnr = top_global.ffparams.atnr;
1663 top->idef.functype = const_cast<t_functype*>(top_global.ffparams.functype.data());
1664 top->idef.iparams = const_cast<t_iparams*>(top_global.ffparams.iparams.data());
1665 top->idef.fudgeQQ = top_global.ffparams.fudgeQQ;
1666 top->idef.cmap_grid = new gmx_cmap_t;
1667 *top->idef.cmap_grid = top_global.ffparams.cmap_grid;
1669 top->idef.ilsort = ilsortUNKNOWN;
1670 top->useInDomainDecomp_ = true;
1673 void dd_init_local_state(gmx_domdec_t* dd, const t_state* state_global, t_state* state_local)
1675 int buf[NITEM_DD_INIT_LOCAL_STATE];
1679 buf[0] = state_global->flags;
1680 buf[1] = state_global->ngtc;
1681 buf[2] = state_global->nnhpres;
1682 buf[3] = state_global->nhchainlength;
1683 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
1685 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE * sizeof(int), buf);
1687 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
1688 init_dfhist_state(state_local, buf[4]);
1689 state_local->flags = buf[0];
1692 /*! \brief Check if a link is stored in \p link between charge groups \p cg_gl and \p cg_gl_j and if not so, store a link */
1693 static void check_link(t_blocka* link, int cg_gl, int cg_gl_j)
1699 for (k = link->index[cg_gl]; k < link->index[cg_gl + 1]; k++)
1701 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
1702 if (link->a[k] == cg_gl_j)
1709 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl + 1] + 1 > link->nalloc_a,
1710 "Inconsistent allocation of link");
1711 /* Add this charge group link */
1712 if (link->index[cg_gl + 1] + 1 > link->nalloc_a)
1714 link->nalloc_a = over_alloc_large(link->index[cg_gl + 1] + 1);
1715 srenew(link->a, link->nalloc_a);
1717 link->a[link->index[cg_gl + 1]] = cg_gl_j;
1718 link->index[cg_gl + 1]++;
1722 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb)
1725 cginfo_mb_t* cgi_mb;
1727 /* For each atom make a list of other atoms in the system
1728 * that a linked to it via bonded interactions
1729 * which are also stored in reverse_top.
1732 reverse_ilist_t ril_intermol;
1733 if (mtop.bIntermolecularInteractions)
1737 atoms.nr = mtop.natoms;
1738 atoms.atom = nullptr;
1740 GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
1741 "We should have an ilist when intermolecular interactions are on");
1743 make_reverse_ilist(*mtop.intermolecular_ilist, &atoms, FALSE, FALSE, FALSE, TRUE, &ril_intermol);
1747 snew(link->index, mtop.natoms + 1);
1754 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
1756 const gmx_molblock_t& molb = mtop.molblock[mb];
1761 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1762 /* Make a reverse ilist in which the interactions are linked
1763 * to all atoms, not only the first atom as in gmx_reverse_top.
1764 * The constraints are discarded here.
1766 reverse_ilist_t ril;
1767 make_reverse_ilist(molt.ilist, &molt.atoms, FALSE, FALSE, FALSE, TRUE, &ril);
1769 cgi_mb = &cginfo_mb[mb];
1772 for (mol = 0; mol < (mtop.bIntermolecularInteractions ? molb.nmol : 1); mol++)
1774 for (int a = 0; a < molt.atoms.nr; a++)
1776 int cg_gl = cg_offset + a;
1777 link->index[cg_gl + 1] = link->index[cg_gl];
1778 int i = ril.index[a];
1779 while (i < ril.index[a + 1])
1781 int ftype = ril.il[i++];
1782 int nral = NRAL(ftype);
1783 /* Skip the ifunc index */
1785 for (int j = 0; j < nral; j++)
1787 int aj = ril.il[i + j];
1790 check_link(link, cg_gl, cg_offset + aj);
1793 i += nral_rt(ftype);
1796 if (mtop.bIntermolecularInteractions)
1798 int i = ril_intermol.index[cg_gl];
1799 while (i < ril_intermol.index[cg_gl + 1])
1801 int ftype = ril_intermol.il[i++];
1802 int nral = NRAL(ftype);
1803 /* Skip the ifunc index */
1805 for (int j = 0; j < nral; j++)
1807 /* Here we assume we have no charge groups;
1808 * this has been checked above.
1810 int aj = ril_intermol.il[i + j];
1811 check_link(link, cg_gl, aj);
1813 i += nral_rt(ftype);
1816 if (link->index[cg_gl + 1] - link->index[cg_gl] > 0)
1818 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[a]);
1823 cg_offset += molt.atoms.nr;
1825 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - molt.atoms.nr];
1829 fprintf(debug, "molecule type '%s' %d atoms has %d atom links through bonded interac.\n",
1830 *molt.name, molt.atoms.nr, nlink_mol);
1833 if (molb.nmol > mol)
1835 /* Copy the data for the rest of the molecules in this block */
1836 link->nalloc_a += (molb.nmol - mol) * nlink_mol;
1837 srenew(link->a, link->nalloc_a);
1838 for (; mol < molb.nmol; mol++)
1840 for (int a = 0; a < molt.atoms.nr; a++)
1842 int cg_gl = cg_offset + a;
1843 link->index[cg_gl + 1] = link->index[cg_gl + 1 - molt.atoms.nr] + nlink_mol;
1844 for (int j = link->index[cg_gl]; j < link->index[cg_gl + 1]; j++)
1846 link->a[j] = link->a[j - nlink_mol] + molt.atoms.nr;
1848 if (link->index[cg_gl + 1] - link->index[cg_gl] > 0
1849 && cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
1851 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
1855 cg_offset += molt.atoms.nr;
1862 fprintf(debug, "Of the %d atoms %d are linked via bonded interactions\n", mtop.natoms, ncgi);
1874 } bonded_distance_t;
1876 /*! \brief Compare distance^2 \p r2 against the distance in \p bd and if larger store it along with \p ftype and atom indices \p a1 and \p a2 */
1877 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2, bonded_distance_t* bd)
1888 /*! \brief Set the distance, function type and atom indices for the longest distance between charge-groups of molecule type \p molt for two-body and multi-body bonded interactions */
1889 static void bonded_cg_distance_mol(const gmx_moltype_t* molt,
1893 bonded_distance_t* bd_2b,
1894 bonded_distance_t* bd_mb)
1896 for (int ftype = 0; ftype < F_NRE; ftype++)
1898 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
1900 const auto& il = molt->ilist[ftype];
1901 int nral = NRAL(ftype);
1904 for (int i = 0; i < il.size(); i += 1 + nral)
1906 for (int ai = 0; ai < nral; ai++)
1908 int atomI = il.iatoms[i + 1 + ai];
1909 for (int aj = ai + 1; aj < nral; aj++)
1911 int atomJ = il.iatoms[i + 1 + aj];
1914 real rij2 = distance2(cg_cm[atomI], cg_cm[atomJ]);
1916 update_max_bonded_distance(rij2, ftype, atomI, atomJ,
1917 (nral == 2) ? bd_2b : bd_mb);
1927 const auto& excls = molt->excls;
1928 for (gmx::index ai = 0; ai < excls.ssize(); ai++)
1930 for (const int aj : excls[ai])
1934 real rij2 = distance2(cg_cm[ai], cg_cm[aj]);
1936 /* There is no function type for exclusions, use -1 */
1937 update_max_bonded_distance(rij2, -1, ai, aj, bd_2b);
1944 /*! \brief Set the distance, function type and atom indices for the longest atom distance involved in intermolecular interactions for two-body and multi-body bonded interactions */
1945 static void bonded_distance_intermol(const InteractionLists& ilists_intermol,
1950 bonded_distance_t* bd_2b,
1951 bonded_distance_t* bd_mb)
1955 set_pbc(&pbc, pbcType, box);
1957 for (int ftype = 0; ftype < F_NRE; ftype++)
1959 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
1961 const auto& il = ilists_intermol[ftype];
1962 int nral = NRAL(ftype);
1964 /* No nral>1 check here, since intermol interactions always
1965 * have nral>=2 (and the code is also correct for nral=1).
1967 for (int i = 0; i < il.size(); i += 1 + nral)
1969 for (int ai = 0; ai < nral; ai++)
1971 int atom_i = il.iatoms[i + 1 + ai];
1973 for (int aj = ai + 1; aj < nral; aj++)
1978 int atom_j = il.iatoms[i + 1 + aj];
1980 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
1984 update_max_bonded_distance(rij2, ftype, atom_i, atom_j, (nral == 2) ? bd_2b : bd_mb);
1992 //! Returns whether \p molt has at least one virtual site
1993 static bool moltypeHasVsite(const gmx_moltype_t& molt)
1995 bool hasVsite = false;
1996 for (int i = 0; i < F_NRE; i++)
1998 if ((interaction_function[i].flags & IF_VSITE) && molt.ilist[i].size() > 0)
2007 //! Returns coordinates not broken over PBC for a molecule
2008 static void getWholeMoleculeCoordinates(const gmx_moltype_t* molt,
2009 const gmx_ffparams_t* ffparams,
2018 if (pbcType != PbcType::No)
2020 mk_mshift(nullptr, graph, pbcType, box, x);
2022 shift_x(graph, box, x, xs);
2023 /* By doing an extra mk_mshift the molecules that are broken
2024 * because they were e.g. imported from another software
2025 * will be made whole again. Such are the healing powers
2028 mk_mshift(nullptr, graph, pbcType, box, xs);
2032 /* We copy the coordinates so the original coordinates remain
2033 * unchanged, just to be 100% sure that we do not affect
2034 * binary reproducibility of simulations.
2037 for (i = 0; i < n; i++)
2039 copy_rvec(x[i], xs[i]);
2043 if (moltypeHasVsite(*molt))
2045 /* Convert to old, deprecated format */
2046 t_ilist ilist[F_NRE];
2047 for (int ftype = 0; ftype < F_NRE; ftype++)
2049 if (interaction_function[ftype].flags & IF_VSITE)
2051 ilist[ftype].nr = molt->ilist[ftype].size();
2052 ilist[ftype].iatoms = const_cast<int*>(molt->ilist[ftype].iatoms.data());
2056 construct_vsites(nullptr, xs, 0.0, nullptr, ffparams->iparams.data(), ilist, PbcType::No,
2057 TRUE, nullptr, nullptr);
2061 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
2062 const gmx_mtop_t* mtop,
2063 const t_inputrec* ir,
2070 gmx_bool bExclRequired;
2074 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2075 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2077 bExclRequired = inputrecExclForces(ir);
2082 for (const gmx_molblock_t& molb : mtop->molblock)
2084 const gmx_moltype_t& molt = mtop->moltype[molb.type];
2085 if (molt.atoms.nr == 1 || molb.nmol == 0)
2087 at_offset += molb.nmol * molt.atoms.nr;
2091 if (ir->pbcType != PbcType::No)
2093 mk_graph_moltype(molt, &graph);
2096 snew(xs, molt.atoms.nr);
2097 for (int mol = 0; mol < molb.nmol; mol++)
2099 getWholeMoleculeCoordinates(&molt, &mtop->ffparams, ir->pbcType, &graph, box,
2102 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2103 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2105 bonded_cg_distance_mol(&molt, bBCheck, bExclRequired, xs, &bd_mol_2b, &bd_mol_mb);
2107 /* Process the mol data adding the atom index offset */
2108 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype, at_offset + bd_mol_2b.a1,
2109 at_offset + bd_mol_2b.a2, &bd_2b);
2110 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype, at_offset + bd_mol_mb.a1,
2111 at_offset + bd_mol_mb.a2, &bd_mb);
2113 at_offset += molt.atoms.nr;
2116 if (ir->pbcType != PbcType::No)
2123 if (mtop->bIntermolecularInteractions)
2125 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist,
2126 "We should have an ilist when intermolecular interactions are on");
2128 bonded_distance_intermol(*mtop->intermolecular_ilist, bBCheck, x, ir->pbcType, box, &bd_2b, &bd_mb);
2131 *r_2b = sqrt(bd_2b.r2);
2132 *r_mb = sqrt(bd_mb.r2);
2134 if (*r_2b > 0 || *r_mb > 0)
2136 GMX_LOG(mdlog.info).appendText("Initial maximum distances in bonded interactions:");
2140 .appendTextFormatted(
2141 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d", *r_2b,
2142 (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2143 bd_2b.a1 + 1, bd_2b.a2 + 1);
2148 .appendTextFormatted(
2149 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d", *r_mb,
2150 interaction_function[bd_mb.ftype].longname, bd_mb.a1 + 1, bd_mb.a2 + 1);