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/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/inputrec.h"
66 #include "gromacs/mdtypes/md_enums.h"
67 #include "gromacs/mdtypes/mdatom.h"
68 #include "gromacs/mdtypes/state.h"
69 #include "gromacs/pbcutil/mshift.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/topology/ifunc.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "gromacs/topology/topsort.h"
74 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/utility/exceptions.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/listoflists.h"
79 #include "gromacs/utility/logger.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/strconvert.h"
82 #include "gromacs/utility/stringstream.h"
83 #include "gromacs/utility/stringutil.h"
84 #include "gromacs/utility/textwriter.h"
86 #include "domdec_constraints.h"
87 #include "domdec_internal.h"
88 #include "domdec_vsite.h"
91 using gmx::ListOfLists;
93 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
94 #define NITEM_DD_INIT_LOCAL_STATE 5
96 struct reverse_ilist_t
98 std::vector<int> index; /* Index for each atom into il */
99 std::vector<int> il; /* ftype|type|a0|...|an|ftype|... */
100 int numAtomsInMolecule; /* The number of atoms in this molecule */
103 struct MolblockIndices
111 /*! \brief Struct for thread local work data for local topology generation */
114 t_idef idef; /**< Partial local topology */
115 std::unique_ptr<VsitePbc> vsitePbc; /**< vsite PBC structure */
116 int nbonded; /**< The number of bondeds in this struct */
117 ListOfLists<int> excl; /**< List of exclusions */
118 int excl_count; /**< The total exclusion count for \p excl */
121 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
122 struct gmx_reverse_top_t
124 //! @cond Doxygen_Suppress
125 //! \brief Are there constraints in this revserse top?
126 bool bConstr = false;
127 //! \brief Are there settles in this revserse top?
128 bool bSettle = false;
129 //! \brief All bonded interactions have to be assigned?
130 bool bBCheck = false;
131 //! \brief Are there bondeds/exclusions between atoms?
132 bool bInterAtomicInteractions = false;
133 //! \brief Reverse ilist for all moltypes
134 std::vector<reverse_ilist_t> ril_mt;
135 //! \brief The size of ril_mt[?].index summed over all entries
136 int ril_mt_tot_size = 0;
137 //! \brief The sorting state of bondeds for free energy
138 int ilsort = ilsortUNKNOWN;
139 //! \brief molblock to global atom index for quick lookup of molblocks on atom index
140 std::vector<MolblockIndices> mbi;
142 //! \brief Do we have intermolecular interactions?
143 bool bIntermolecularInteractions = false;
144 //! \brief Intermolecular reverse ilist
145 reverse_ilist_t ril_intermol;
147 /* Work data structures for multi-threading */
148 //! \brief Thread work array for local topology generation
149 std::vector<thread_work_t> th_work;
153 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
154 static int nral_rt(int ftype)
159 if (interaction_function[ftype].flags & IF_VSITE)
161 /* With vsites the reverse topology contains an extra entry
162 * for storing if constructing atoms are vsites.
170 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
171 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck, gmx_bool bConstr, gmx_bool bSettle)
173 return ((((interaction_function[ftype].flags & IF_BOND) != 0U)
174 && ((interaction_function[ftype].flags & IF_VSITE) == 0U)
175 && (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U)))
176 || (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) || (bSettle && ftype == F_SETTLE));
179 /*! \brief Help print error output when interactions are missing */
180 static std::string print_missing_interactions_mb(t_commrec* cr,
181 const gmx_reverse_top_t* rt,
182 const char* moltypename,
183 const reverse_ilist_t* ril,
191 int nril_mol = ril->index[nat_mol];
192 snew(assigned, nmol * nril_mol);
193 gmx::StringOutputStream stream;
194 gmx::TextWriter log(&stream);
196 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
197 for (int ftype = 0; ftype < F_NRE; ftype++)
199 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
201 int nral = NRAL(ftype);
202 const t_ilist* il = &idef->il[ftype];
203 const t_iatom* ia = il->iatoms;
204 for (int i = 0; i < il->nr; i += 1 + nral)
206 int a0 = gatindex[ia[1]];
207 /* Check if this interaction is in
208 * the currently checked molblock.
210 if (a0 >= a_start && a0 < a_end)
212 int mol = (a0 - a_start) / nat_mol;
213 int a0_mol = (a0 - a_start) - mol * nat_mol;
214 int j_mol = ril->index[a0_mol];
216 while (j_mol < ril->index[a0_mol + 1] && !found)
218 int j = mol * nril_mol + j_mol;
219 int ftype_j = ril->il[j_mol];
220 /* Here we need to check if this interaction has
221 * not already been assigned, since we could have
222 * multiply defined interactions.
224 if (ftype == ftype_j && ia[0] == ril->il[j_mol + 1] && assigned[j] == 0)
226 /* Check the atoms */
228 for (int a = 0; a < nral; a++)
230 if (gatindex[ia[1 + a]] != a_start + mol * nat_mol + ril->il[j_mol + 2 + a])
240 j_mol += 2 + nral_rt(ftype_j);
244 gmx_incons("Some interactions seem to be assigned multiple times");
252 gmx_sumi(nmol * nril_mol, assigned, cr);
256 for (int mol = 0; mol < nmol; mol++)
259 while (j_mol < nril_mol)
261 int ftype = ril->il[j_mol];
262 int nral = NRAL(ftype);
263 int j = mol * nril_mol + j_mol;
264 if (assigned[j] == 0 && !(interaction_function[ftype].flags & IF_VSITE))
266 if (DDMASTER(cr->dd))
270 log.writeLineFormatted("Molecule type '%s'", moltypename);
271 log.writeLineFormatted(
272 "the first %d missing interactions, except for exclusions:", nprint);
274 log.writeStringFormatted("%20s atoms", interaction_function[ftype].longname);
276 for (a = 0; a < nral; a++)
278 log.writeStringFormatted("%5d", ril->il[j_mol + 2 + a] + 1);
282 log.writeString(" ");
285 log.writeString(" global");
286 for (a = 0; a < nral; a++)
288 log.writeStringFormatted(
289 "%6d", a_start + mol * nat_mol + ril->il[j_mol + 2 + a] + 1);
291 log.ensureLineBreak();
299 j_mol += 2 + nral_rt(ftype);
304 return stream.toString();
307 /*! \brief Help print error output when interactions are missing */
308 static void print_missing_interactions_atoms(const gmx::MDLogger& mdlog,
310 const gmx_mtop_t* mtop,
313 const gmx_reverse_top_t* rt = cr->dd->reverse_top;
315 /* Print the atoms in the missing interactions per molblock */
317 for (const gmx_molblock_t& molb : mtop->molblock)
319 const gmx_moltype_t& moltype = mtop->moltype[molb.type];
321 a_end = a_start + molb.nmol * moltype.atoms.nr;
323 GMX_LOG(mdlog.warning)
324 .appendText(print_missing_interactions_mb(cr, rt, *(moltype.name),
325 &rt->ril_mt[molb.type], a_start, a_end,
326 moltype.atoms.nr, molb.nmol, idef));
330 void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
333 const gmx_mtop_t* top_global,
334 const gmx_localtop_t* top_local,
338 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
344 GMX_LOG(mdlog.warning)
346 "Not all bonded interactions have been properly assigned to the domain "
347 "decomposition cells");
349 ndiff_tot = local_count - dd->nbonded_global;
351 for (ftype = 0; ftype < F_NRE; ftype++)
354 cl[ftype] = top_local->idef.il[ftype].nr / (1 + nral);
357 gmx_sumi(F_NRE, cl, cr);
361 GMX_LOG(mdlog.warning).appendText("A list of missing interactions:");
362 rest_global = dd->nbonded_global;
363 rest_local = local_count;
364 for (ftype = 0; ftype < F_NRE; ftype++)
366 /* In the reverse and local top all constraints are merged
367 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
368 * and add these constraints when doing F_CONSTR.
370 if (((interaction_function[ftype].flags & IF_BOND)
371 && (dd->reverse_top->bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO)))
372 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
373 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
375 n = gmx_mtop_ftype_count(top_global, ftype);
376 if (ftype == F_CONSTR)
378 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
380 ndiff = cl[ftype] - n;
383 GMX_LOG(mdlog.warning)
384 .appendTextFormatted("%20s of %6d missing %6d",
385 interaction_function[ftype].longname, n, -ndiff);
388 rest_local -= cl[ftype];
392 ndiff = rest_local - rest_global;
395 GMX_LOG(mdlog.warning).appendTextFormatted("%20s of %6d missing %6d", "exclusions", rest_global, -ndiff);
399 print_missing_interactions_atoms(mdlog, cr, top_global, &top_local->idef);
400 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr, -1, x, box);
402 std::string errorMessage;
407 "One or more interactions were assigned to multiple domains of the domain "
408 "decompostion. Please report this bug.";
412 errorMessage = gmx::formatString(
413 "%d of the %d bonded interactions could not be calculated because some atoms "
414 "involved moved further apart than the multi-body cut-off distance (%g nm) or the "
415 "two-body cut-off distance (%g nm), see option -rdd, for pairs and tabulated bonds "
416 "also see option -ddcheck",
417 -ndiff_tot, cr->dd->nbonded_global, dd_cutoff_multibody(dd), dd_cutoff_twobody(dd));
419 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
422 /*! \brief Return global topology molecule information for global atom index \p i_gl */
423 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)
425 const MolblockIndices* mbi = rt->mbi.data();
427 int end = rt->mbi.size(); /* exclusive */
430 /* binary search for molblock_ind */
433 mid = (start + end) / 2;
434 if (i_gl >= mbi[mid].a_end)
438 else if (i_gl < mbi[mid].a_start)
452 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
453 *i_mol = (i_gl - mbi->a_start) - (*mol) * mbi->natoms_mol;
456 /*! \brief Returns the maximum number of exclusions per atom */
457 static int getMaxNumExclusionsPerAtom(const ListOfLists<int>& excls)
460 for (gmx::index at = 0; at < excls.ssize(); at++)
462 const auto list = excls[at];
463 const int numExcls = list.ssize();
465 GMX_RELEASE_ASSERT(numExcls != 1 || list[0] == at,
466 "With 1 exclusion we expect a self-exclusion");
468 maxNumExcls = std::max(maxNumExcls, numExcls);
474 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
475 static int low_make_reverse_ilist(const InteractionLists& il_mt,
481 gmx::ArrayRef<const int> r_index,
482 gmx::ArrayRef<int> r_il,
483 gmx_bool bLinkToAllAtoms,
486 int ftype, j, nlink, link;
491 for (ftype = 0; ftype < F_NRE; ftype++)
493 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE))
494 || (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) || (bSettle && ftype == F_SETTLE))
496 const bool bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
497 const int nral = NRAL(ftype);
498 const auto& il = il_mt[ftype];
499 for (int i = 0; i < il.size(); i += 1 + nral)
501 const int* ia = il.iatoms.data() + i;
506 /* We don't need the virtual sites for the cg-links */
516 /* Couple to the first atom in the interaction */
519 for (link = 0; link < nlink; link++)
524 GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
525 GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
526 r_il[r_index[a] + count[a]] = (ftype == F_CONSTRNC ? F_CONSTR : ftype);
527 r_il[r_index[a] + count[a] + 1] = ia[0];
528 for (j = 1; j < 1 + nral; j++)
530 /* Store the molecular atom number */
531 r_il[r_index[a] + count[a] + 1 + j] = ia[j];
534 if (interaction_function[ftype].flags & IF_VSITE)
538 /* Add an entry to iatoms for storing
539 * which of the constructing atoms are
542 r_il[r_index[a] + count[a] + 2 + nral] = 0;
543 for (j = 2; j < 1 + nral; j++)
545 if (atom[ia[j]].ptype == eptVSite)
547 r_il[r_index[a] + count[a] + 2 + nral] |= (2 << j);
554 /* We do not count vsites since they are always
555 * uniquely assigned and can be assigned
556 * to multiple nodes with recursive vsites.
558 if (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))
563 count[a] += 2 + nral_rt(ftype);
572 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
573 static int make_reverse_ilist(const InteractionLists& ilist,
574 const t_atoms* atoms,
578 gmx_bool bLinkToAllAtoms,
579 reverse_ilist_t* ril_mt)
581 int nat_mt, *count, i, nint_mt;
583 /* Count the interactions */
586 low_make_reverse_ilist(ilist, atoms->atom, count, bConstr, bSettle, bBCheck, {}, {},
587 bLinkToAllAtoms, FALSE);
589 ril_mt->index.push_back(0);
590 for (i = 0; i < nat_mt; i++)
592 ril_mt->index.push_back(ril_mt->index[i] + count[i]);
595 ril_mt->il.resize(ril_mt->index[nat_mt]);
597 /* Store the interactions */
598 nint_mt = low_make_reverse_ilist(ilist, atoms->atom, count, bConstr, bSettle, bBCheck,
599 ril_mt->index, ril_mt->il, bLinkToAllAtoms, TRUE);
603 ril_mt->numAtomsInMolecule = atoms->nr;
608 /*! \brief Generate the reverse topology */
609 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t* mtop,
616 gmx_reverse_top_t rt;
618 /* Should we include constraints (for SHAKE) in rt? */
619 rt.bConstr = bConstr;
620 rt.bSettle = bSettle;
621 rt.bBCheck = bBCheck;
623 rt.bInterAtomicInteractions = mtop->bIntermolecularInteractions;
624 rt.ril_mt.resize(mtop->moltype.size());
625 rt.ril_mt_tot_size = 0;
626 std::vector<int> nint_mt;
627 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
629 const gmx_moltype_t& molt = mtop->moltype[mt];
630 if (molt.atoms.nr > 1)
632 rt.bInterAtomicInteractions = true;
635 /* Make the atom to interaction list for this molecule type */
636 int numberOfInteractions = make_reverse_ilist(
637 molt.ilist, &molt.atoms, rt.bConstr, rt.bSettle, rt.bBCheck, FALSE, &rt.ril_mt[mt]);
638 nint_mt.push_back(numberOfInteractions);
640 rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
644 fprintf(debug, "The total size of the atom to interaction index is %d integers\n",
649 for (const gmx_molblock_t& molblock : mtop->molblock)
651 *nint += molblock.nmol * nint_mt[molblock.type];
654 /* Make an intermolecular reverse top, if necessary */
655 rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
656 if (rt.bIntermolecularInteractions)
658 t_atoms atoms_global;
660 atoms_global.nr = mtop->natoms;
661 atoms_global.atom = nullptr; /* Only used with virtual sites */
663 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist,
664 "We should have an ilist when intermolecular interactions are on");
666 *nint += make_reverse_ilist(*mtop->intermolecular_ilist, &atoms_global, rt.bConstr,
667 rt.bSettle, rt.bBCheck, FALSE, &rt.ril_intermol);
670 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
672 rt.ilsort = ilsortFE_UNSORTED;
676 rt.ilsort = ilsortNO_FE;
679 /* Make a molblock index for fast searching */
681 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
683 const gmx_molblock_t& molb = mtop->molblock[mb];
684 const int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
687 i += molb.nmol * numAtomsPerMol;
689 mbi.natoms_mol = numAtomsPerMol;
690 mbi.type = molb.type;
691 rt.mbi.push_back(mbi);
694 rt.th_work.resize(gmx_omp_nthreads_get(emntDomdec));
699 void dd_make_reverse_top(FILE* fplog,
701 const gmx_mtop_t* mtop,
702 const gmx_vsite_t* vsite,
703 const t_inputrec* ir,
708 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
711 /* If normal and/or settle constraints act only within charge groups,
712 * we can store them in the reverse top and simply assign them to domains.
713 * Otherwise we need to assign them to multiple domains and set up
714 * the parallel version constraint algorithm(s).
717 dd->reverse_top = new gmx_reverse_top_t;
719 make_reverse_top(mtop, ir->efep != efepNO, !dd->comm->systemInfo.haveSplitConstraints,
720 !dd->comm->systemInfo.haveSplitSettles, bBCheck, &dd->nbonded_global);
722 dd->haveExclusions = false;
723 for (const gmx_molblock_t& molb : mtop->molblock)
725 const int maxNumExclusionsPerAtom = getMaxNumExclusionsPerAtom(mtop->moltype[molb.type].excls);
726 // We checked above that max 1 exclusion means only self exclusions
727 if (maxNumExclusionsPerAtom > 1)
729 dd->haveExclusions = true;
733 if (vsite && vsite->numInterUpdategroupVsites > 0)
738 "There are %d inter update-group virtual sites,\n"
739 "will an extra communication step for selected coordinates and forces\n",
740 vsite->numInterUpdategroupVsites);
742 init_domdec_vsites(dd, vsite->numInterUpdategroupVsites);
745 if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
747 init_domdec_constraints(dd, mtop);
751 fprintf(fplog, "\n");
755 /*! \brief Store a vsite interaction at the end of \p il
757 * This routine is very similar to add_ifunc, but vsites interactions
758 * have more work to do than other kinds of interactions, and the
759 * complex way nral (and thus vector contents) depends on ftype
760 * confuses static analysis tools unless we fuse the vsite
761 * atom-indexing organization code with the ifunc-adding code, so that
762 * they can see that nral is the same value. */
763 static inline void add_ifunc_for_vsites(t_iatom* tiatoms,
764 const gmx_ga2la_t& ga2la,
770 const t_iatom* iatoms,
775 if (il->nr + 1 + nral > il->nalloc)
777 il->nalloc = over_alloc_large(il->nr + 1 + nral);
778 srenew(il->iatoms, il->nalloc);
780 liatoms = il->iatoms + il->nr;
784 tiatoms[0] = iatoms[0];
788 /* We know the local index of the first atom */
793 /* Convert later in make_local_vsites */
794 tiatoms[1] = -a_gl - 1;
797 for (int k = 2; k < 1 + nral; k++)
799 int ak_gl = a_gl + iatoms[k] - a_mol;
800 if (const int* homeIndex = ga2la.findHome(ak_gl))
802 tiatoms[k] = *homeIndex;
806 /* Copy the global index, convert later in make_local_vsites */
807 tiatoms[k] = -(ak_gl + 1);
809 // Note that ga2la_get_home always sets the third parameter if
812 for (int k = 0; k < 1 + nral; k++)
814 liatoms[k] = tiatoms[k];
818 /*! \brief Store a bonded interaction at the end of \p il */
819 static inline void add_ifunc(int nral, const t_iatom* tiatoms, t_ilist* il)
824 if (il->nr + 1 + nral > il->nalloc)
826 il->nalloc = over_alloc_large(il->nr + 1 + nral);
827 srenew(il->iatoms, il->nalloc);
829 liatoms = il->iatoms + il->nr;
830 for (k = 0; k <= nral; k++)
832 liatoms[k] = tiatoms[k];
837 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
838 static void add_posres(int mol,
840 int numAtomsInMolecule,
841 const gmx_molblock_t* molb,
843 const t_iparams* ip_in,
849 /* This position restraint has not been added yet,
850 * so it's index is the current number of position restraints.
852 n = idef->il[F_POSRES].nr / 2;
853 if (n + 1 > idef->iparams_posres_nalloc)
855 idef->iparams_posres_nalloc = over_alloc_dd(n + 1);
856 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
858 ip = &idef->iparams_posres[n];
859 /* Copy the force constants */
860 *ip = ip_in[iatoms[0]];
862 /* Get the position restraint coordinates from the molblock */
863 a_molb = mol * numAtomsInMolecule + a_mol;
864 GMX_ASSERT(a_molb < ssize(molb->posres_xA),
865 "We need a sufficient number of position restraint coordinates");
866 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
867 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
868 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
869 if (!molb->posres_xB.empty())
871 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
872 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
873 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
877 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
878 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
879 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
881 /* Set the parameter index for idef->iparams_posre */
885 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
886 static void add_fbposres(int mol,
888 int numAtomsInMolecule,
889 const gmx_molblock_t* molb,
891 const t_iparams* ip_in,
897 /* This flat-bottom position restraint has not been added yet,
898 * so it's index is the current number of position restraints.
900 n = idef->il[F_FBPOSRES].nr / 2;
901 if (n + 1 > idef->iparams_fbposres_nalloc)
903 idef->iparams_fbposres_nalloc = over_alloc_dd(n + 1);
904 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
906 ip = &idef->iparams_fbposres[n];
907 /* Copy the force constants */
908 *ip = ip_in[iatoms[0]];
910 /* Get the position restraint coordinats from the molblock */
911 a_molb = mol * numAtomsInMolecule + a_mol;
912 GMX_ASSERT(a_molb < ssize(molb->posres_xA),
913 "We need a sufficient number of position restraint coordinates");
914 /* Take reference positions from A position of normal posres */
915 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
916 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
917 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
919 /* Note: no B-type for flat-bottom posres */
921 /* Set the parameter index for idef->iparams_posre */
925 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
926 static void add_vsite(const gmx_ga2la_t& ga2la,
927 gmx::ArrayRef<const int> index,
928 gmx::ArrayRef<const int> rtil,
935 const t_iatom* iatoms,
939 t_iatom tiatoms[1 + MAXATOMLIST];
940 int j, ftype_r, nral_r;
942 /* Add this interaction to the local topology */
943 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
945 if (iatoms[1 + nral])
947 /* Check for recursion */
948 for (k = 2; k < 1 + nral; k++)
950 if ((iatoms[1 + nral] & (2 << k)) && (tiatoms[k] < 0))
952 /* This construction atoms is a vsite and not a home atom */
955 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n",
956 iatoms[k] + 1, a_mol + 1);
958 /* Find the vsite construction */
960 /* Check all interactions assigned to this atom */
961 j = index[iatoms[k]];
962 while (j < index[iatoms[k] + 1])
965 nral_r = NRAL(ftype_r);
966 if (interaction_function[ftype_r].flags & IF_VSITE)
968 /* Add this vsite (recursion) */
969 add_vsite(ga2la, index, rtil, ftype_r, nral_r, FALSE, -1,
970 a_gl + iatoms[k] - iatoms[1], iatoms[k], rtil.data() + j, idef);
972 j += 1 + nral_rt(ftype_r);
979 /*! \brief Returns the squared distance between atoms \p i and \p j */
980 static real dd_dist2(t_pbc* pbc_null, const rvec* x, const int i, int j)
986 pbc_dx_aiuc(pbc_null, x[i], x[j], dx);
990 rvec_sub(x[i], x[j], dx);
996 /*! \brief Append t_idef structures 1 to nsrc in src to *dest */
997 static void combine_idef(t_idef* dest, gmx::ArrayRef<const thread_work_t> src)
1001 for (ftype = 0; ftype < F_NRE; ftype++)
1004 for (gmx::index s = 1; s < src.ssize(); s++)
1006 n += src[s].idef.il[ftype].nr;
1012 ild = &dest->il[ftype];
1014 if (ild->nr + n > ild->nalloc)
1016 ild->nalloc = over_alloc_large(ild->nr + n);
1017 srenew(ild->iatoms, ild->nalloc);
1020 for (gmx::index s = 1; s < src.ssize(); s++)
1022 const t_ilist& ils = src[s].idef.il[ftype];
1024 for (int i = 0; i < ils.nr; i++)
1026 ild->iatoms[ild->nr + i] = ils.iatoms[i];
1032 /* Position restraints need an additional treatment */
1033 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1035 int nposres = dest->il[ftype].nr / 2;
1036 // TODO: Simplify this code using std::vector
1037 t_iparams*& iparams_dest =
1038 (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1039 int& posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc
1040 : dest->iparams_fbposres_nalloc);
1041 if (nposres > posres_nalloc)
1043 posres_nalloc = over_alloc_large(nposres);
1044 srenew(iparams_dest, posres_nalloc);
1047 /* Set nposres to the number of original position restraints in dest */
1048 for (gmx::index s = 1; s < src.ssize(); s++)
1050 nposres -= src[s].idef.il[ftype].nr / 2;
1053 for (gmx::index s = 1; s < src.ssize(); s++)
1055 const t_iparams* iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres
1056 : src[s].idef.iparams_fbposres);
1058 for (int i = 0; i < src[s].idef.il[ftype].nr / 2; i++)
1060 /* Correct the index into iparams_posres */
1061 dest->il[ftype].iatoms[nposres * 2] = nposres;
1062 /* Copy the position restraint force parameters */
1063 iparams_dest[nposres] = iparams_src[i];
1072 /*! \brief Check and when available assign bonded interactions for local atom i
1074 static inline void check_assign_interactions_atom(int i,
1078 int numAtomsInMolecule,
1079 gmx::ArrayRef<const int> index,
1080 gmx::ArrayRef<const int> rtil,
1081 gmx_bool bInterMolInteractions,
1084 const gmx_domdec_t* dd,
1085 const gmx_domdec_zones_t* zones,
1086 const gmx_molblock_t* molb,
1093 const t_iparams* ip_in,
1099 gmx::ArrayRef<const DDPairInteractionRanges> iZones = zones->iZones;
1104 t_iatom tiatoms[1 + MAXATOMLIST];
1106 const int ftype = rtil[j++];
1107 auto iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1108 const int nral = NRAL(ftype);
1109 if (interaction_function[ftype].flags & IF_VSITE)
1111 assert(!bInterMolInteractions);
1112 /* The vsite construction goes where the vsite itself is */
1115 add_vsite(*dd->ga2la, index, rtil, ftype, nral, TRUE, i, i_gl, i_mol, iatoms.data(), idef);
1124 tiatoms[0] = iatoms[0];
1128 assert(!bInterMolInteractions);
1129 /* Assign single-body interactions to the home zone */
1134 if (ftype == F_POSRES)
1136 add_posres(mol, i_mol, numAtomsInMolecule, molb, tiatoms, ip_in, idef);
1138 else if (ftype == F_FBPOSRES)
1140 add_fbposres(mol, i_mol, numAtomsInMolecule, molb, tiatoms, ip_in, idef);
1150 /* This is a two-body interaction, we can assign
1151 * analogous to the non-bonded assignments.
1155 if (!bInterMolInteractions)
1157 /* Get the global index using the offset in the molecule */
1158 k_gl = i_gl + iatoms[2] - i_mol;
1164 if (const auto* entry = dd->ga2la->find(k_gl))
1166 int kz = entry->cell;
1171 /* Check zone interaction assignments */
1172 bUse = ((iz < iZones.ssize() && iz <= kz && iZones[iz].jZoneRange.isInRange(kz))
1173 || (kz < iZones.ssize() && iz > kz && iZones[kz].jZoneRange.isInRange(iz)));
1176 GMX_ASSERT(ftype != F_CONSTR || (iz == 0 && kz == 0),
1177 "Constraint assigned here should only involve home atoms");
1180 tiatoms[2] = entry->la;
1181 /* If necessary check the cgcm distance */
1182 if (bRCheck2B && dd_dist2(pbc_null, cg_cm, tiatoms[1], tiatoms[2]) >= rc2)
1195 /* Assign this multi-body bonded interaction to
1196 * the local node if we have all the atoms involved
1197 * (local or communicated) and the minimum zone shift
1198 * in each dimension is zero, for dimensions
1199 * with 2 DD cells an extra check may be necessary.
1201 ivec k_zero, k_plus;
1207 for (k = 1; k <= nral && bUse; k++)
1210 if (!bInterMolInteractions)
1212 /* Get the global index using the offset in the molecule */
1213 k_gl = i_gl + iatoms[k] - i_mol;
1219 const auto* entry = dd->ga2la->find(k_gl);
1220 if (entry == nullptr || entry->cell >= zones->n)
1222 /* We do not have this atom of this interaction
1223 * locally, or it comes from more than one cell
1232 tiatoms[k] = entry->la;
1233 for (d = 0; d < DIM; d++)
1235 if (zones->shift[entry->cell][d] == 0)
1246 bUse = (bUse && (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1251 for (d = 0; (d < DIM && bUse); d++)
1253 /* Check if the cg_cm distance falls within
1254 * the cut-off to avoid possible multiple
1255 * assignments of bonded interactions.
1257 if (rcheck[d] && k_plus[d]
1258 && dd_dist2(pbc_null, cg_cm, tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1267 /* Add this interaction to the local topology */
1268 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1269 /* Sum so we can check in global_stat
1270 * if we have everything.
1272 if (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))
1282 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1284 * With thread parallelizing each thread acts on a different atom range:
1285 * at_start to at_end.
1287 static int make_bondeds_zone(gmx_domdec_t* dd,
1288 const gmx_domdec_zones_t* zones,
1289 const std::vector<gmx_molblock_t>& molb,
1296 const t_iparams* ip_in,
1299 const gmx::Range<int>& atomRange)
1301 int mb, mt, mol, i_mol;
1303 gmx_reverse_top_t* rt;
1306 rt = dd->reverse_top;
1308 bBCheck = rt->bBCheck;
1312 for (int i : atomRange)
1314 /* Get the global atom number */
1315 const int i_gl = dd->globalAtomIndices[i];
1316 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1317 /* Check all intramolecular interactions assigned to this atom */
1318 gmx::ArrayRef<const int> index = rt->ril_mt[mt].index;
1319 gmx::ArrayRef<const t_iatom> rtil = rt->ril_mt[mt].il;
1321 check_assign_interactions_atom(i, i_gl, mol, i_mol, rt->ril_mt[mt].numAtomsInMolecule,
1322 index, rtil, FALSE, index[i_mol], index[i_mol + 1], dd,
1323 zones, &molb[mb], bRCheckMB, rcheck, bRCheck2B, rc2,
1324 pbc_null, cg_cm, ip_in, idef, izone, bBCheck, &nbonded_local);
1327 if (rt->bIntermolecularInteractions)
1329 /* Check all intermolecular interactions assigned to this atom */
1330 index = rt->ril_intermol.index;
1331 rtil = rt->ril_intermol.il;
1333 check_assign_interactions_atom(i, i_gl, mol, i_mol, rt->ril_mt[mt].numAtomsInMolecule,
1334 index, rtil, TRUE, index[i_gl], index[i_gl + 1], dd, zones,
1335 &molb[mb], bRCheckMB, rcheck, bRCheck2B, rc2, pbc_null,
1336 cg_cm, ip_in, idef, izone, bBCheck, &nbonded_local);
1340 return nbonded_local;
1343 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1344 static void set_no_exclusions_zone(const gmx_domdec_zones_t* zones, int iz, ListOfLists<int>* lexcls)
1346 for (int a = zones->cg_range[iz]; a < zones->cg_range[iz + 1]; a++)
1348 lexcls->pushBack({});
1352 /*! \brief Set the exclusion data for i-zone \p iz */
1353 static void make_exclusions_zone(gmx_domdec_t* dd,
1354 gmx_domdec_zones_t* zones,
1355 const std::vector<gmx_moltype_t>& moltype,
1357 ListOfLists<int>* lexcls,
1361 const gmx::ArrayRef<const int> intermolecularExclusionGroup)
1363 const gmx_ga2la_t& ga2la = *dd->ga2la;
1365 const auto& jAtomRange = zones->iZones[iz].jAtomRange;
1367 const gmx::index oldNumLists = lexcls->ssize();
1369 std::vector<int> exclusionsForAtom;
1370 for (int at = at_start; at < at_end; at++)
1372 exclusionsForAtom.clear();
1374 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1376 int a_gl, mb, mt, mol, a_mol;
1378 /* Copy the exclusions from the global top */
1379 a_gl = dd->globalAtomIndices[at];
1380 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1381 const auto excls = moltype[mt].excls[a_mol];
1382 for (const int aj_mol : excls)
1384 if (const auto* jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1386 /* This check is not necessary, but it can reduce
1387 * the number of exclusions in the list, which in turn
1388 * can speed up the pair list construction a bit.
1390 if (jAtomRange.isInRange(jEntry->la))
1392 exclusionsForAtom.push_back(jEntry->la);
1398 bool isExcludedAtom = !intermolecularExclusionGroup.empty()
1399 && std::find(intermolecularExclusionGroup.begin(),
1400 intermolecularExclusionGroup.end(), dd->globalAtomIndices[at])
1401 != intermolecularExclusionGroup.end();
1405 for (int qmAtomGlobalIndex : intermolecularExclusionGroup)
1407 if (const auto* entry = dd->ga2la->find(qmAtomGlobalIndex))
1409 exclusionsForAtom.push_back(entry->la);
1414 /* Append the exclusions for this atom to the topology */
1415 lexcls->pushBack(exclusionsForAtom);
1419 lexcls->ssize() - oldNumLists == at_end - at_start,
1420 "The number of exclusion list should match the number of atoms in the range");
1423 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1424 static void finish_local_exclusions(gmx_domdec_t* dd, gmx_domdec_zones_t* zones, ListOfLists<int>* lexcls)
1426 const gmx::Range<int> nonhomeIzonesAtomRange(zones->iZones[0].iAtomRange.end(),
1427 zones->iZones.back().iAtomRange.end());
1429 if (!dd->haveExclusions)
1431 /* There are no exclusions involving non-home charge groups,
1432 * but we need to set the indices for neighborsearching.
1434 for (int gmx_unused la : nonhomeIzonesAtomRange)
1436 lexcls->pushBack({});
1441 /*! \brief Clear a t_idef struct */
1442 static void clear_idef(t_idef* idef)
1446 /* Clear the counts */
1447 for (ftype = 0; ftype < F_NRE; ftype++)
1449 idef->il[ftype].nr = 0;
1453 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1454 static int make_local_bondeds_excls(gmx_domdec_t* dd,
1455 gmx_domdec_zones_t* zones,
1456 const gmx_mtop_t* mtop,
1465 ListOfLists<int>* lexcls,
1468 int nzone_bondeds, nzone_excl;
1472 gmx_reverse_top_t* rt;
1474 if (dd->reverse_top->bInterAtomicInteractions)
1476 nzone_bondeds = zones->n;
1480 /* Only single charge group (or atom) molecules, so interactions don't
1481 * cross zone boundaries and we only need to assign in the home zone.
1486 if (dd->haveExclusions)
1488 /* We only use exclusions from i-zones to i- and j-zones */
1489 nzone_excl = zones->iZones.size();
1493 /* There are no exclusions and only zone 0 sees itself */
1497 rt = dd->reverse_top;
1501 /* Clear the counts */
1508 for (int izone = 0; izone < nzone_bondeds; izone++)
1510 cg0 = zones->cg_range[izone];
1511 cg1 = zones->cg_range[izone + 1];
1513 const int numThreads = rt->th_work.size();
1514 #pragma omp parallel for num_threads(numThreads) schedule(static)
1515 for (int thread = 0; thread < numThreads; thread++)
1522 cg0t = cg0 + ((cg1 - cg0) * thread) / numThreads;
1523 cg1t = cg0 + ((cg1 - cg0) * (thread + 1)) / numThreads;
1531 idef_t = &rt->th_work[thread].idef;
1535 rt->th_work[thread].nbonded = make_bondeds_zone(
1536 dd, zones, mtop->molblock, bRCheckMB, rcheck, bRCheck2B, rc2, pbc_null,
1537 cg_cm, idef->iparams, idef_t, izone, gmx::Range<int>(cg0t, cg1t));
1539 if (izone < nzone_excl)
1541 ListOfLists<int>* excl_t;
1544 // Thread 0 stores exclusions directly in the final storage
1549 // Threads > 0 store in temporary storage, starting at list index 0
1550 excl_t = &rt->th_work[thread].excl;
1554 /* No charge groups and no distance check required */
1555 make_exclusions_zone(dd, zones, mtop->moltype, cginfo, excl_t, izone, cg0t,
1556 cg1t, mtop->intermolecularExclusionGroup);
1559 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1562 if (rt->th_work.size() > 1)
1564 combine_idef(idef, rt->th_work);
1567 for (const thread_work_t& th_work : rt->th_work)
1569 nbonded_local += th_work.nbonded;
1572 if (izone < nzone_excl)
1574 for (std::size_t th = 1; th < rt->th_work.size(); th++)
1576 lexcls->appendListOfLists(rt->th_work[th].excl);
1578 for (const thread_work_t& th_work : rt->th_work)
1580 *excl_count += th_work.excl_count;
1585 /* Some zones might not have exclusions, but some code still needs to
1586 * loop over the index, so we set the indices here.
1588 for (size_t iZone = nzone_excl; iZone < zones->iZones.size(); iZone++)
1590 set_no_exclusions_zone(zones, iZone, lexcls);
1593 finish_local_exclusions(dd, zones, lexcls);
1596 fprintf(debug, "We have %d exclusions, check count %d\n", lexcls->numElements(), *excl_count);
1599 return nbonded_local;
1602 void dd_make_local_top(gmx_domdec_t* dd,
1603 gmx_domdec_zones_t* zones,
1610 const gmx_mtop_t& mtop,
1611 gmx_localtop_t* ltop)
1613 gmx_bool bRCheckMB, bRCheck2B;
1617 t_pbc pbc, *pbc_null = nullptr;
1621 fprintf(debug, "Making local topology\n");
1627 if (dd->reverse_top->bInterAtomicInteractions)
1629 /* We need to check to which cell bondeds should be assigned */
1630 rc = dd_cutoff_twobody(dd);
1633 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1636 /* Should we check cg_cm distances when assigning bonded interactions? */
1637 for (d = 0; d < DIM; d++)
1640 /* Only need to check for dimensions where the part of the box
1641 * that is not communicated is smaller than the cut-off.
1643 if (d < npbcdim && dd->numCells[d] > 1
1644 && (dd->numCells[d] - npulse[d]) * cellsize_min[d] < 2 * rc)
1646 if (dd->numCells[d] == 2)
1651 /* Check for interactions between two atoms,
1652 * where we can allow interactions up to the cut-off,
1653 * instead of up to the smallest cell dimension.
1659 fprintf(debug, "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n", d,
1660 cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
1663 if (bRCheckMB || bRCheck2B)
1667 pbc_null = set_pbc_dd(&pbc, fr->pbcType, dd->numCells, TRUE, box);
1676 dd->nbonded_local = make_local_bondeds_excls(dd, zones, &mtop, fr->cginfo.data(), bRCheckMB,
1677 rcheck, bRCheck2B, rc, pbc_null, cgcm_or_x,
1678 <op->idef, <op->excls, &nexcl);
1680 /* The ilist is not sorted yet,
1681 * we can only do this when we have the charge arrays.
1683 ltop->idef.ilsort = ilsortUNKNOWN;
1685 ltop->atomtypes = mtop.atomtypes;
1688 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop)
1690 if (dd->reverse_top->ilsort == ilsortNO_FE)
1692 ltop->idef.ilsort = ilsortNO_FE;
1696 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
1700 void dd_init_local_top(const gmx_mtop_t& top_global, gmx_localtop_t* top)
1702 /* TODO: Get rid of the const casts below, e.g. by using a reference */
1703 top->idef.ntypes = top_global.ffparams.numTypes();
1704 top->idef.atnr = top_global.ffparams.atnr;
1705 top->idef.functype = const_cast<t_functype*>(top_global.ffparams.functype.data());
1706 top->idef.iparams = const_cast<t_iparams*>(top_global.ffparams.iparams.data());
1707 top->idef.fudgeQQ = top_global.ffparams.fudgeQQ;
1708 top->idef.cmap_grid = new gmx_cmap_t;
1709 *top->idef.cmap_grid = top_global.ffparams.cmap_grid;
1711 top->idef.ilsort = ilsortUNKNOWN;
1712 top->useInDomainDecomp_ = true;
1715 void dd_init_local_state(gmx_domdec_t* dd, const t_state* state_global, t_state* state_local)
1717 int buf[NITEM_DD_INIT_LOCAL_STATE];
1721 buf[0] = state_global->flags;
1722 buf[1] = state_global->ngtc;
1723 buf[2] = state_global->nnhpres;
1724 buf[3] = state_global->nhchainlength;
1725 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
1727 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE * sizeof(int), buf);
1729 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
1730 init_dfhist_state(state_local, buf[4]);
1731 state_local->flags = buf[0];
1734 /*! \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 */
1735 static void check_link(t_blocka* link, int cg_gl, int cg_gl_j)
1741 for (k = link->index[cg_gl]; k < link->index[cg_gl + 1]; k++)
1743 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
1744 if (link->a[k] == cg_gl_j)
1751 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl + 1] + 1 > link->nalloc_a,
1752 "Inconsistent allocation of link");
1753 /* Add this charge group link */
1754 if (link->index[cg_gl + 1] + 1 > link->nalloc_a)
1756 link->nalloc_a = over_alloc_large(link->index[cg_gl + 1] + 1);
1757 srenew(link->a, link->nalloc_a);
1759 link->a[link->index[cg_gl + 1]] = cg_gl_j;
1760 link->index[cg_gl + 1]++;
1764 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb)
1767 cginfo_mb_t* cgi_mb;
1769 /* For each atom make a list of other atoms in the system
1770 * that a linked to it via bonded interactions
1771 * which are also stored in reverse_top.
1774 reverse_ilist_t ril_intermol;
1775 if (mtop.bIntermolecularInteractions)
1779 atoms.nr = mtop.natoms;
1780 atoms.atom = nullptr;
1782 GMX_RELEASE_ASSERT(mtop.intermolecular_ilist,
1783 "We should have an ilist when intermolecular interactions are on");
1785 make_reverse_ilist(*mtop.intermolecular_ilist, &atoms, FALSE, FALSE, FALSE, TRUE, &ril_intermol);
1789 snew(link->index, mtop.natoms + 1);
1796 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
1798 const gmx_molblock_t& molb = mtop.molblock[mb];
1803 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1804 /* Make a reverse ilist in which the interactions are linked
1805 * to all atoms, not only the first atom as in gmx_reverse_top.
1806 * The constraints are discarded here.
1808 reverse_ilist_t ril;
1809 make_reverse_ilist(molt.ilist, &molt.atoms, FALSE, FALSE, FALSE, TRUE, &ril);
1811 cgi_mb = &cginfo_mb[mb];
1814 for (mol = 0; mol < (mtop.bIntermolecularInteractions ? molb.nmol : 1); mol++)
1816 for (int a = 0; a < molt.atoms.nr; a++)
1818 int cg_gl = cg_offset + a;
1819 link->index[cg_gl + 1] = link->index[cg_gl];
1820 int i = ril.index[a];
1821 while (i < ril.index[a + 1])
1823 int ftype = ril.il[i++];
1824 int nral = NRAL(ftype);
1825 /* Skip the ifunc index */
1827 for (int j = 0; j < nral; j++)
1829 int aj = ril.il[i + j];
1832 check_link(link, cg_gl, cg_offset + aj);
1835 i += nral_rt(ftype);
1838 if (mtop.bIntermolecularInteractions)
1840 int i = ril_intermol.index[cg_gl];
1841 while (i < ril_intermol.index[cg_gl + 1])
1843 int ftype = ril_intermol.il[i++];
1844 int nral = NRAL(ftype);
1845 /* Skip the ifunc index */
1847 for (int j = 0; j < nral; j++)
1849 /* Here we assume we have no charge groups;
1850 * this has been checked above.
1852 int aj = ril_intermol.il[i + j];
1853 check_link(link, cg_gl, aj);
1855 i += nral_rt(ftype);
1858 if (link->index[cg_gl + 1] - link->index[cg_gl] > 0)
1860 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[a]);
1865 cg_offset += molt.atoms.nr;
1867 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - molt.atoms.nr];
1871 fprintf(debug, "molecule type '%s' %d atoms has %d atom links through bonded interac.\n",
1872 *molt.name, molt.atoms.nr, nlink_mol);
1875 if (molb.nmol > mol)
1877 /* Copy the data for the rest of the molecules in this block */
1878 link->nalloc_a += (molb.nmol - mol) * nlink_mol;
1879 srenew(link->a, link->nalloc_a);
1880 for (; mol < molb.nmol; mol++)
1882 for (int a = 0; a < molt.atoms.nr; a++)
1884 int cg_gl = cg_offset + a;
1885 link->index[cg_gl + 1] = link->index[cg_gl + 1 - molt.atoms.nr] + nlink_mol;
1886 for (int j = link->index[cg_gl]; j < link->index[cg_gl + 1]; j++)
1888 link->a[j] = link->a[j - nlink_mol] + molt.atoms.nr;
1890 if (link->index[cg_gl + 1] - link->index[cg_gl] > 0
1891 && cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
1893 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
1897 cg_offset += molt.atoms.nr;
1904 fprintf(debug, "Of the %d atoms %d are linked via bonded interactions\n", mtop.natoms, ncgi);
1916 } bonded_distance_t;
1918 /*! \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 */
1919 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2, bonded_distance_t* bd)
1930 /*! \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 */
1931 static void bonded_cg_distance_mol(const gmx_moltype_t* molt,
1935 bonded_distance_t* bd_2b,
1936 bonded_distance_t* bd_mb)
1938 for (int ftype = 0; ftype < F_NRE; ftype++)
1940 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
1942 const auto& il = molt->ilist[ftype];
1943 int nral = NRAL(ftype);
1946 for (int i = 0; i < il.size(); i += 1 + nral)
1948 for (int ai = 0; ai < nral; ai++)
1950 int atomI = il.iatoms[i + 1 + ai];
1951 for (int aj = ai + 1; aj < nral; aj++)
1953 int atomJ = il.iatoms[i + 1 + aj];
1956 real rij2 = distance2(cg_cm[atomI], cg_cm[atomJ]);
1958 update_max_bonded_distance(rij2, ftype, atomI, atomJ,
1959 (nral == 2) ? bd_2b : bd_mb);
1969 const auto& excls = molt->excls;
1970 for (gmx::index ai = 0; ai < excls.ssize(); ai++)
1972 for (const int aj : excls[ai])
1976 real rij2 = distance2(cg_cm[ai], cg_cm[aj]);
1978 /* There is no function type for exclusions, use -1 */
1979 update_max_bonded_distance(rij2, -1, ai, aj, bd_2b);
1986 /*! \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 */
1987 static void bonded_distance_intermol(const InteractionLists& ilists_intermol,
1992 bonded_distance_t* bd_2b,
1993 bonded_distance_t* bd_mb)
1997 set_pbc(&pbc, pbcType, box);
1999 for (int ftype = 0; ftype < F_NRE; ftype++)
2001 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2003 const auto& il = ilists_intermol[ftype];
2004 int nral = NRAL(ftype);
2006 /* No nral>1 check here, since intermol interactions always
2007 * have nral>=2 (and the code is also correct for nral=1).
2009 for (int i = 0; i < il.size(); i += 1 + nral)
2011 for (int ai = 0; ai < nral; ai++)
2013 int atom_i = il.iatoms[i + 1 + ai];
2015 for (int aj = ai + 1; aj < nral; aj++)
2020 int atom_j = il.iatoms[i + 1 + aj];
2022 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2026 update_max_bonded_distance(rij2, ftype, atom_i, atom_j, (nral == 2) ? bd_2b : bd_mb);
2034 //! Returns whether \p molt has at least one virtual site
2035 static bool moltypeHasVsite(const gmx_moltype_t& molt)
2037 bool hasVsite = false;
2038 for (int i = 0; i < F_NRE; i++)
2040 if ((interaction_function[i].flags & IF_VSITE) && molt.ilist[i].size() > 0)
2049 //! Returns coordinates not broken over PBC for a molecule
2050 static void getWholeMoleculeCoordinates(const gmx_moltype_t* molt,
2051 const gmx_ffparams_t* ffparams,
2060 if (pbcType != PbcType::No)
2062 mk_mshift(nullptr, graph, pbcType, box, x);
2064 shift_x(graph, box, x, xs);
2065 /* By doing an extra mk_mshift the molecules that are broken
2066 * because they were e.g. imported from another software
2067 * will be made whole again. Such are the healing powers
2070 mk_mshift(nullptr, graph, pbcType, box, xs);
2074 /* We copy the coordinates so the original coordinates remain
2075 * unchanged, just to be 100% sure that we do not affect
2076 * binary reproducibility of simulations.
2079 for (i = 0; i < n; i++)
2081 copy_rvec(x[i], xs[i]);
2085 if (moltypeHasVsite(*molt))
2087 /* Convert to old, deprecated format */
2088 t_ilist ilist[F_NRE];
2089 for (int ftype = 0; ftype < F_NRE; ftype++)
2091 if (interaction_function[ftype].flags & IF_VSITE)
2093 ilist[ftype].nr = molt->ilist[ftype].size();
2094 ilist[ftype].iatoms = const_cast<int*>(molt->ilist[ftype].iatoms.data());
2098 construct_vsites(nullptr, xs, 0.0, nullptr, ffparams->iparams.data(), ilist, PbcType::No,
2099 TRUE, nullptr, nullptr);
2103 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
2104 const gmx_mtop_t* mtop,
2105 const t_inputrec* ir,
2112 gmx_bool bExclRequired;
2116 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2117 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2119 bExclRequired = inputrecExclForces(ir);
2124 for (const gmx_molblock_t& molb : mtop->molblock)
2126 const gmx_moltype_t& molt = mtop->moltype[molb.type];
2127 if (molt.atoms.nr == 1 || molb.nmol == 0)
2129 at_offset += molb.nmol * molt.atoms.nr;
2133 if (ir->pbcType != PbcType::No)
2135 mk_graph_moltype(molt, &graph);
2138 snew(xs, molt.atoms.nr);
2139 for (int mol = 0; mol < molb.nmol; mol++)
2141 getWholeMoleculeCoordinates(&molt, &mtop->ffparams, ir->pbcType, &graph, box,
2144 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2145 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2147 bonded_cg_distance_mol(&molt, bBCheck, bExclRequired, xs, &bd_mol_2b, &bd_mol_mb);
2149 /* Process the mol data adding the atom index offset */
2150 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype, at_offset + bd_mol_2b.a1,
2151 at_offset + bd_mol_2b.a2, &bd_2b);
2152 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype, at_offset + bd_mol_mb.a1,
2153 at_offset + bd_mol_mb.a2, &bd_mb);
2155 at_offset += molt.atoms.nr;
2158 if (ir->pbcType != PbcType::No)
2165 if (mtop->bIntermolecularInteractions)
2167 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist,
2168 "We should have an ilist when intermolecular interactions are on");
2170 bonded_distance_intermol(*mtop->intermolecular_ilist, bBCheck, x, ir->pbcType, box, &bd_2b, &bd_mb);
2173 *r_2b = sqrt(bd_2b.r2);
2174 *r_mb = sqrt(bd_mb.r2);
2176 if (*r_2b > 0 || *r_mb > 0)
2178 GMX_LOG(mdlog.info).appendText("Initial maximum distances in bonded interactions:");
2182 .appendTextFormatted(
2183 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d", *r_2b,
2184 (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2185 bd_2b.a1 + 1, bd_2b.a2 + 1);
2190 .appendTextFormatted(
2191 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d", *r_mb,
2192 interaction_function[bd_mb.ftype].longname, bd_mb.a1 + 1, bd_mb.a2 + 1);