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, 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/chargegroup.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/forcerec.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/mshift.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/topology/ifunc.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/topology/topsort.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/exceptions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/gmxassert.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 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
92 #define NITEM_DD_INIT_LOCAL_STATE 5
94 struct reverse_ilist_t
96 std::vector<int> index; /* Index for each atom into il */
97 std::vector<int> il; /* ftype|type|a0|...|an|ftype|... */
98 int numAtomsInMolecule; /* The number of atoms in this molecule */
101 struct MolblockIndices
109 /*! \brief Struct for thread local work data for local topology generation */
112 t_idef idef; /**< Partial local topology */
113 std::unique_ptr<VsitePbc> vsitePbc; /**< vsite PBC structure */
114 int nbonded; /**< The number of bondeds in this struct */
115 t_blocka excl; /**< List of exclusions */
116 int excl_count; /**< The total exclusion count for \p excl */
119 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
120 struct gmx_reverse_top_t
122 //! @cond Doxygen_Suppress
123 //! \brief The maximum number of exclusions one atom can have
124 int n_excl_at_max = 0;
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 charge-groups?
132 bool bInterCGInteractions = 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,
172 gmx_bool bConstr, gmx_bool bSettle)
174 return ((((interaction_function[ftype].flags & IF_BOND) != 0U) &&
175 ((interaction_function[ftype].flags & IF_VSITE) == 0U) &&
176 (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0U))) ||
177 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
178 (bSettle && ftype == F_SETTLE));
181 /*! \brief Help print error output when interactions are missing */
183 print_missing_interactions_mb(t_commrec *cr,
184 const gmx_reverse_top_t *rt,
185 const char *moltypename,
186 const reverse_ilist_t *ril,
187 int a_start, int a_end,
188 int nat_mol, int nmol,
192 int nril_mol = ril->index[nat_mol];
193 snew(assigned, nmol*nril_mol);
194 gmx::StringOutputStream stream;
195 gmx::TextWriter log(&stream);
197 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
198 for (int ftype = 0; ftype < F_NRE; ftype++)
200 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
202 int nral = NRAL(ftype);
203 const t_ilist *il = &idef->il[ftype];
204 const t_iatom *ia = il->iatoms;
205 for (int i = 0; i < il->nr; i += 1+nral)
207 int a0 = gatindex[ia[1]];
208 /* Check if this interaction is in
209 * the currently checked molblock.
211 if (a0 >= a_start && a0 < a_end)
213 int mol = (a0 - a_start)/nat_mol;
214 int a0_mol = (a0 - a_start) - mol*nat_mol;
215 int j_mol = ril->index[a0_mol];
217 while (j_mol < ril->index[a0_mol+1] && !found)
219 int j = mol*nril_mol + j_mol;
220 int ftype_j = ril->il[j_mol];
221 /* Here we need to check if this interaction has
222 * not already been assigned, since we could have
223 * multiply defined interactions.
225 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
228 /* Check the atoms */
230 for (int a = 0; a < nral; a++)
232 if (gatindex[ia[1+a]] !=
233 a_start + mol*nat_mol + ril->il[j_mol+2+a])
243 j_mol += 2 + nral_rt(ftype_j);
247 gmx_incons("Some interactions seem to be assigned multiple times");
255 gmx_sumi(nmol*nril_mol, assigned, cr);
259 for (int mol = 0; mol < nmol; mol++)
262 while (j_mol < nril_mol)
264 int ftype = ril->il[j_mol];
265 int nral = NRAL(ftype);
266 int j = mol*nril_mol + j_mol;
267 if (assigned[j] == 0 &&
268 !(interaction_function[ftype].flags & IF_VSITE))
270 if (DDMASTER(cr->dd))
274 log.writeLineFormatted("Molecule type '%s'", moltypename);
275 log.writeLineFormatted(
276 "the first %d missing interactions, except for exclusions:", nprint);
278 log.writeStringFormatted("%20s atoms",
279 interaction_function[ftype].longname);
281 for (a = 0; a < nral; a++)
283 log.writeStringFormatted("%5d", ril->il[j_mol+2+a]+1);
287 log.writeString(" ");
290 log.writeString(" global");
291 for (a = 0; a < nral; a++)
293 log.writeStringFormatted("%6d",
294 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
296 log.ensureLineBreak();
304 j_mol += 2 + nral_rt(ftype);
309 return stream.toString();
312 /*! \brief Help print error output when interactions are missing */
313 static void print_missing_interactions_atoms(const gmx::MDLogger &mdlog,
315 const gmx_mtop_t *mtop,
318 const gmx_reverse_top_t *rt = cr->dd->reverse_top;
320 /* Print the atoms in the missing interactions per molblock */
322 for (const gmx_molblock_t &molb : mtop->molblock)
324 const gmx_moltype_t &moltype = mtop->moltype[molb.type];
326 a_end = a_start + molb.nmol*moltype.atoms.nr;
328 GMX_LOG(mdlog.warning).appendText(
329 print_missing_interactions_mb(cr, rt,
331 &rt->ril_mt[molb.type],
332 a_start, a_end, moltype.atoms.nr,
338 void dd_print_missing_interactions(const gmx::MDLogger &mdlog,
341 const gmx_mtop_t *top_global,
342 const gmx_localtop_t *top_local,
346 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
352 GMX_LOG(mdlog.warning).appendText(
353 "Not all bonded interactions have been properly assigned to the domain decomposition cells");
355 ndiff_tot = local_count - dd->nbonded_global;
357 for (ftype = 0; ftype < F_NRE; ftype++)
360 cl[ftype] = top_local->idef.il[ftype].nr/(1+nral);
363 gmx_sumi(F_NRE, cl, cr);
367 GMX_LOG(mdlog.warning).appendText("A list of missing interactions:");
368 rest_global = dd->nbonded_global;
369 rest_local = local_count;
370 for (ftype = 0; ftype < F_NRE; ftype++)
372 /* In the reverse and local top all constraints are merged
373 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
374 * and add these constraints when doing F_CONSTR.
376 if (((interaction_function[ftype].flags & IF_BOND) &&
377 (dd->reverse_top->bBCheck
378 || !(interaction_function[ftype].flags & IF_LIMZERO)))
379 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
380 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
382 n = gmx_mtop_ftype_count(top_global, ftype);
383 if (ftype == F_CONSTR)
385 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
387 ndiff = cl[ftype] - n;
390 GMX_LOG(mdlog.warning).appendTextFormatted(
391 "%20s of %6d missing %6d",
392 interaction_function[ftype].longname, n, -ndiff);
395 rest_local -= cl[ftype];
399 ndiff = rest_local - rest_global;
402 GMX_LOG(mdlog.warning).appendTextFormatted(
403 "%20s of %6d missing %6d", "exclusions",
404 rest_global, -ndiff);
408 print_missing_interactions_atoms(mdlog, cr, top_global, &top_local->idef);
409 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
412 std::string errorMessage;
416 errorMessage = "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
420 errorMessage = gmx::formatString("%d of the %d bonded interactions could not be calculated because some atoms involved moved further apart than the multi-body cut-off distance (%g nm) or the two-body cut-off distance (%g nm), see option -rdd, for pairs and tabulated bonds also see option -ddcheck", -ndiff_tot, cr->dd->nbonded_global, dd_cutoff_multibody(dd), dd_cutoff_twobody(dd));
422 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
425 /*! \brief Return global topology molecule information for global atom index \p i_gl */
426 static void global_atomnr_to_moltype_ind(const gmx_reverse_top_t *rt,
428 int *mb, int *mt, int *mol, int *i_mol)
430 const MolblockIndices *mbi = rt->mbi.data();
432 int end = rt->mbi.size(); /* exclusive */
435 /* binary search for molblock_ind */
439 if (i_gl >= mbi[mid].a_end)
443 else if (i_gl < mbi[mid].a_start)
457 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
458 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
461 /*! \brief Count the exclusions for all atoms in \p cgs */
462 static void count_excls(const t_block *cgs, const t_blocka *excls,
463 int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
465 int cg, at0, at1, at, excl, atj;
470 for (cg = 0; cg < cgs->nr; cg++)
472 at0 = cgs->index[cg];
473 at1 = cgs->index[cg+1];
474 for (at = at0; at < at1; at++)
476 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
478 atj = excls->a[excl];
482 if (atj < at0 || atj >= at1)
489 *n_excl_at_max = std::max(*n_excl_at_max,
490 excls->index[at+1] - excls->index[at]);
495 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
496 static int low_make_reverse_ilist(const InteractionLists &il_mt,
499 gmx_bool bConstr, gmx_bool bSettle,
501 gmx::ArrayRef<const int> r_index,
502 gmx::ArrayRef<int> r_il,
503 gmx_bool bLinkToAllAtoms,
506 int ftype, j, nlink, link;
511 for (ftype = 0; ftype < F_NRE; ftype++)
513 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
514 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
515 (bSettle && ftype == F_SETTLE))
517 const bool bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0U);
518 const int nral = NRAL(ftype);
519 const auto &il = il_mt[ftype];
520 for (int i = 0; i < il.size(); i += 1+nral)
522 const int* ia = il.iatoms.data() + i;
527 /* We don't need the virtual sites for the cg-links */
537 /* Couple to the first atom in the interaction */
540 for (link = 0; link < nlink; link++)
545 GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
546 GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
547 r_il[r_index[a]+count[a]] =
548 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
549 r_il[r_index[a]+count[a]+1] = ia[0];
550 for (j = 1; j < 1+nral; j++)
552 /* Store the molecular atom number */
553 r_il[r_index[a]+count[a]+1+j] = ia[j];
556 if (interaction_function[ftype].flags & IF_VSITE)
560 /* Add an entry to iatoms for storing
561 * which of the constructing atoms are
564 r_il[r_index[a]+count[a]+2+nral] = 0;
565 for (j = 2; j < 1+nral; j++)
567 if (atom[ia[j]].ptype == eptVSite)
569 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
576 /* We do not count vsites since they are always
577 * uniquely assigned and can be assigned
578 * to multiple nodes with recursive vsites.
581 !(interaction_function[ftype].flags & IF_LIMZERO))
586 count[a] += 2 + nral_rt(ftype);
595 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
596 static int make_reverse_ilist(const InteractionLists &ilist,
597 const t_atoms *atoms,
598 gmx_bool bConstr, gmx_bool bSettle,
600 gmx_bool bLinkToAllAtoms,
601 reverse_ilist_t *ril_mt)
603 int nat_mt, *count, i, nint_mt;
605 /* Count the interactions */
608 low_make_reverse_ilist(ilist, atoms->atom,
610 bConstr, bSettle, bBCheck,
612 bLinkToAllAtoms, FALSE);
614 ril_mt->index.push_back(0);
615 for (i = 0; i < nat_mt; i++)
617 ril_mt->index.push_back(ril_mt->index[i] + count[i]);
620 ril_mt->il.resize(ril_mt->index[nat_mt]);
622 /* Store the interactions */
624 low_make_reverse_ilist(ilist, atoms->atom,
626 bConstr, bSettle, bBCheck,
627 ril_mt->index, ril_mt->il,
628 bLinkToAllAtoms, TRUE);
632 ril_mt->numAtomsInMolecule = atoms->nr;
637 /*! \brief Generate the reverse topology */
638 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
639 gmx_bool bConstr, gmx_bool bSettle,
640 gmx_bool bBCheck, int *nint)
642 gmx_reverse_top_t rt;
644 /* Should we include constraints (for SHAKE) in rt? */
645 rt.bConstr = bConstr;
646 rt.bSettle = bSettle;
647 rt.bBCheck = bBCheck;
649 rt.bInterCGInteractions = mtop->bIntermolecularInteractions;
650 rt.ril_mt.resize(mtop->moltype.size());
651 rt.ril_mt_tot_size = 0;
652 std::vector<int> nint_mt;
653 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
655 const gmx_moltype_t &molt = mtop->moltype[mt];
658 rt.bInterCGInteractions = true;
661 /* Make the atom to interaction list for this molecule type */
662 int numberOfInteractions =
663 make_reverse_ilist(molt.ilist, &molt.atoms,
664 rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
666 nint_mt.push_back(numberOfInteractions);
668 rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
672 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt.ril_mt_tot_size);
676 for (const gmx_molblock_t &molblock : mtop->molblock)
678 *nint += molblock.nmol*nint_mt[molblock.type];
681 /* Make an intermolecular reverse top, if necessary */
682 rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
683 if (rt.bIntermolecularInteractions)
685 t_atoms atoms_global;
687 atoms_global.nr = mtop->natoms;
688 atoms_global.atom = nullptr; /* Only used with virtual sites */
690 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
693 make_reverse_ilist(*mtop->intermolecular_ilist,
695 rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
699 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
701 rt.ilsort = ilsortFE_UNSORTED;
705 rt.ilsort = ilsortNO_FE;
708 /* Make a molblock index for fast searching */
710 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
712 const gmx_molblock_t &molb = mtop->molblock[mb];
713 const int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
716 i += molb.nmol*numAtomsPerMol;
718 mbi.natoms_mol = numAtomsPerMol;
719 mbi.type = molb.type;
720 rt.mbi.push_back(mbi);
723 rt.th_work.resize(gmx_omp_nthreads_get(emntDomdec));
728 void dd_make_reverse_top(FILE *fplog,
729 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
730 const gmx_vsite_t *vsite,
731 const t_inputrec *ir, gmx_bool bBCheck)
735 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
738 /* If normal and/or settle constraints act only within charge groups,
739 * we can store them in the reverse top and simply assign them to domains.
740 * Otherwise we need to assign them to multiple domains and set up
741 * the parallel version constraint algorithm(s).
744 dd->reverse_top = new gmx_reverse_top_t;
746 make_reverse_top(mtop, ir->efep != efepNO,
747 !dd->comm->systemInfo.haveSplitConstraints,
748 !dd->comm->systemInfo.haveSplitSettles,
749 bBCheck, &dd->nbonded_global);
751 gmx_reverse_top_t *rt = dd->reverse_top;
753 dd->n_intercg_excl = 0;
754 rt->n_excl_at_max = 0;
755 for (const gmx_molblock_t &molb : mtop->molblock)
757 int n_excl_mol, n_excl_icg, n_excl_at_max;
759 const gmx_moltype_t &molt = mtop->moltype[molb.type];
760 count_excls(&molt.cgs, &molt.excls,
761 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
762 dd->n_intercg_excl += molb.nmol*n_excl_icg;
763 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
766 if (vsite && vsite->numInterUpdategroupVsites > 0)
770 fprintf(fplog, "There are %d inter update-group virtual sites,\n"
771 "will an extra communication step for selected coordinates and forces\n",
772 vsite->numInterUpdategroupVsites);
774 init_domdec_vsites(dd, vsite->numInterUpdategroupVsites);
777 if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
779 init_domdec_constraints(dd, mtop);
783 fprintf(fplog, "\n");
787 /*! \brief Store a vsite interaction at the end of \p il
789 * This routine is very similar to add_ifunc, but vsites interactions
790 * have more work to do than other kinds of interactions, and the
791 * complex way nral (and thus vector contents) depends on ftype
792 * confuses static analysis tools unless we fuse the vsite
793 * atom-indexing organization code with the ifunc-adding code, so that
794 * they can see that nral is the same value. */
796 add_ifunc_for_vsites(t_iatom *tiatoms, const gmx_ga2la_t &ga2la,
797 int nral, gmx_bool bHomeA,
798 int a, int a_gl, int a_mol,
799 const t_iatom *iatoms,
804 if (il->nr+1+nral > il->nalloc)
806 il->nalloc = over_alloc_large(il->nr+1+nral);
807 srenew(il->iatoms, il->nalloc);
809 liatoms = il->iatoms + il->nr;
813 tiatoms[0] = iatoms[0];
817 /* We know the local index of the first atom */
822 /* Convert later in make_local_vsites */
823 tiatoms[1] = -a_gl - 1;
826 for (int k = 2; k < 1+nral; k++)
828 int ak_gl = a_gl + iatoms[k] - a_mol;
829 if (const int *homeIndex = ga2la.findHome(ak_gl))
831 tiatoms[k] = *homeIndex;
835 /* Copy the global index, convert later in make_local_vsites */
836 tiatoms[k] = -(ak_gl + 1);
838 // Note that ga2la_get_home always sets the third parameter if
841 for (int k = 0; k < 1+nral; k++)
843 liatoms[k] = tiatoms[k];
847 /*! \brief Store a bonded interaction at the end of \p il */
848 static inline void add_ifunc(int nral, const t_iatom *tiatoms, t_ilist *il)
853 if (il->nr+1+nral > il->nalloc)
855 il->nalloc = over_alloc_large(il->nr+1+nral);
856 srenew(il->iatoms, il->nalloc);
858 liatoms = il->iatoms + il->nr;
859 for (k = 0; k <= nral; k++)
861 liatoms[k] = tiatoms[k];
866 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
867 static void add_posres(int mol, int a_mol, int numAtomsInMolecule,
868 const gmx_molblock_t *molb,
869 t_iatom *iatoms, const t_iparams *ip_in,
875 /* This position restraint has not been added yet,
876 * so it's index is the current number of position restraints.
878 n = idef->il[F_POSRES].nr/2;
879 if (n+1 > idef->iparams_posres_nalloc)
881 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
882 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
884 ip = &idef->iparams_posres[n];
885 /* Copy the force constants */
886 *ip = ip_in[iatoms[0]];
888 /* Get the position restraint coordinates from the molblock */
889 a_molb = mol*numAtomsInMolecule + a_mol;
890 GMX_ASSERT(a_molb < ssize(molb->posres_xA), "We need a sufficient number of position restraint coordinates");
891 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
892 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
893 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
894 if (!molb->posres_xB.empty())
896 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
897 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
898 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
902 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
903 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
904 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
906 /* Set the parameter index for idef->iparams_posre */
910 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
911 static void add_fbposres(int mol, int a_mol, int numAtomsInMolecule,
912 const gmx_molblock_t *molb,
913 t_iatom *iatoms, const t_iparams *ip_in,
919 /* This flat-bottom position restraint has not been added yet,
920 * so it's index is the current number of position restraints.
922 n = idef->il[F_FBPOSRES].nr/2;
923 if (n+1 > idef->iparams_fbposres_nalloc)
925 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
926 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
928 ip = &idef->iparams_fbposres[n];
929 /* Copy the force constants */
930 *ip = ip_in[iatoms[0]];
932 /* Get the position restraint coordinats from the molblock */
933 a_molb = mol*numAtomsInMolecule + a_mol;
934 GMX_ASSERT(a_molb < ssize(molb->posres_xA), "We need a sufficient number of position restraint coordinates");
935 /* Take reference positions from A position of normal posres */
936 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
937 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
938 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
940 /* Note: no B-type for flat-bottom posres */
942 /* Set the parameter index for idef->iparams_posre */
946 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
947 static void add_vsite(const gmx_ga2la_t &ga2la,
948 gmx::ArrayRef<const int> index,
949 gmx::ArrayRef<const int> rtil,
951 gmx_bool bHomeA, int a, int a_gl, int a_mol,
952 const t_iatom *iatoms,
956 t_iatom tiatoms[1+MAXATOMLIST];
957 int j, ftype_r, nral_r;
959 /* Add this interaction to the local topology */
960 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
964 /* Check for recursion */
965 for (k = 2; k < 1+nral; k++)
967 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
969 /* This construction atoms is a vsite and not a home atom */
972 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
974 /* Find the vsite construction */
976 /* Check all interactions assigned to this atom */
977 j = index[iatoms[k]];
978 while (j < index[iatoms[k]+1])
981 nral_r = NRAL(ftype_r);
982 if (interaction_function[ftype_r].flags & IF_VSITE)
984 /* Add this vsite (recursion) */
985 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
986 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
990 j += 1 + nral_rt(ftype_r);
997 /*! \brief Returns the squared distance between atoms \p i and \p j */
998 static real dd_dist2(t_pbc *pbc_null, const rvec *x, const int i, int j)
1004 pbc_dx_aiuc(pbc_null, x[i], x[j], dx);
1008 rvec_sub(x[i], x[j], dx);
1014 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1015 static void combine_blocka(t_blocka *dest,
1016 gmx::ArrayRef<const thread_work_t> src)
1018 int ni = src.back().excl.nr;
1020 for (const thread_work_t &th_work : src)
1022 na += th_work.excl.nra;
1024 if (ni + 1 > dest->nalloc_index)
1026 dest->nalloc_index = over_alloc_large(ni+1);
1027 srenew(dest->index, dest->nalloc_index);
1029 if (dest->nra + na > dest->nalloc_a)
1031 dest->nalloc_a = over_alloc_large(dest->nra+na);
1032 srenew(dest->a, dest->nalloc_a);
1034 for (gmx::index s = 1; s < src.ssize(); s++)
1036 for (int i = dest->nr + 1; i < src[s].excl.nr + 1; i++)
1038 dest->index[i] = dest->nra + src[s].excl.index[i];
1040 for (int i = 0; i < src[s].excl.nra; i++)
1042 dest->a[dest->nra+i] = src[s].excl.a[i];
1044 dest->nr = src[s].excl.nr;
1045 dest->nra += src[s].excl.nra;
1049 /*! \brief Append t_idef structures 1 to nsrc in src to *dest */
1050 static void combine_idef(t_idef *dest,
1051 gmx::ArrayRef<const thread_work_t> src)
1055 for (ftype = 0; ftype < F_NRE; ftype++)
1058 for (gmx::index s = 1; s < src.ssize(); s++)
1060 n += src[s].idef.il[ftype].nr;
1066 ild = &dest->il[ftype];
1068 if (ild->nr + n > ild->nalloc)
1070 ild->nalloc = over_alloc_large(ild->nr+n);
1071 srenew(ild->iatoms, ild->nalloc);
1074 for (gmx::index s = 1; s < src.ssize(); s++)
1076 const t_ilist &ils = src[s].idef.il[ftype];
1078 for (int i = 0; i < ils.nr; i++)
1080 ild->iatoms[ild->nr + i] = ils.iatoms[i];
1086 /* Position restraints need an additional treatment */
1087 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1089 int nposres = dest->il[ftype].nr/2;
1090 // TODO: Simplify this code using std::vector
1091 t_iparams * &iparams_dest = (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1092 int &posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc : dest->iparams_fbposres_nalloc);
1093 if (nposres > posres_nalloc)
1095 posres_nalloc = over_alloc_large(nposres);
1096 srenew(iparams_dest, posres_nalloc);
1099 /* Set nposres to the number of original position restraints in dest */
1100 for (gmx::index s = 1; s < src.ssize(); s++)
1102 nposres -= src[s].idef.il[ftype].nr/2;
1105 for (gmx::index s = 1; s < src.ssize(); s++)
1107 const t_iparams *iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1109 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1111 /* Correct the index into iparams_posres */
1112 dest->il[ftype].iatoms[nposres*2] = nposres;
1113 /* Copy the position restraint force parameters */
1114 iparams_dest[nposres] = iparams_src[i];
1123 /*! \brief Check and when available assign bonded interactions for local atom i
1126 check_assign_interactions_atom(int i, int i_gl,
1128 int numAtomsInMolecule,
1129 gmx::ArrayRef<const int> index,
1130 gmx::ArrayRef<const int> rtil,
1131 gmx_bool bInterMolInteractions,
1132 int ind_start, int ind_end,
1133 const gmx_domdec_t *dd,
1134 const gmx_domdec_zones_t *zones,
1135 const gmx_molblock_t *molb,
1136 gmx_bool bRCheckMB, const ivec rcheck, gmx_bool bRCheck2B,
1140 const t_iparams *ip_in,
1151 t_iatom tiatoms[1 + MAXATOMLIST];
1153 const int ftype = rtil[j++];
1154 auto iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1155 const int nral = NRAL(ftype);
1156 if (interaction_function[ftype].flags & IF_VSITE)
1158 assert(!bInterMolInteractions);
1159 /* The vsite construction goes where the vsite itself is */
1162 add_vsite(*dd->ga2la, index, rtil, ftype, nral,
1163 TRUE, i, i_gl, i_mol,
1164 iatoms.data(), idef);
1173 tiatoms[0] = iatoms[0];
1177 assert(!bInterMolInteractions);
1178 /* Assign single-body interactions to the home zone */
1183 if (ftype == F_POSRES)
1185 add_posres(mol, i_mol, numAtomsInMolecule,
1186 molb, tiatoms, ip_in, idef);
1188 else if (ftype == F_FBPOSRES)
1190 add_fbposres(mol, i_mol, numAtomsInMolecule,
1191 molb, tiatoms, ip_in, idef);
1201 /* This is a two-body interaction, we can assign
1202 * analogous to the non-bonded assignments.
1206 if (!bInterMolInteractions)
1208 /* Get the global index using the offset in the molecule */
1209 k_gl = i_gl + iatoms[2] - i_mol;
1215 if (const auto *entry = dd->ga2la->find(k_gl))
1217 int kz = entry->cell;
1222 /* Check zone interaction assignments */
1223 bUse = ((iz < zones->nizone &&
1225 kz >= zones->izone[iz].j0 &&
1226 kz < zones->izone[iz].j1) ||
1227 (kz < zones->nizone &&
1229 iz >= zones->izone[kz].j0 &&
1230 iz < zones->izone[kz].j1));
1233 GMX_ASSERT(ftype != F_CONSTR || (iz == 0 && kz == 0),
1234 "Constraint assigned here should only involve home atoms");
1237 tiatoms[2] = entry->la;
1238 /* If necessary check the cgcm distance */
1240 dd_dist2(pbc_null, cg_cm,
1241 tiatoms[1], tiatoms[2]) >= rc2)
1254 /* Assign this multi-body bonded interaction to
1255 * the local node if we have all the atoms involved
1256 * (local or communicated) and the minimum zone shift
1257 * in each dimension is zero, for dimensions
1258 * with 2 DD cells an extra check may be necessary.
1260 ivec k_zero, k_plus;
1266 for (k = 1; k <= nral && bUse; k++)
1269 if (!bInterMolInteractions)
1271 /* Get the global index using the offset in the molecule */
1272 k_gl = i_gl + iatoms[k] - i_mol;
1278 const auto *entry = dd->ga2la->find(k_gl);
1279 if (entry == nullptr || entry->cell >= zones->n)
1281 /* We do not have this atom of this interaction
1282 * locally, or it comes from more than one cell
1291 tiatoms[k] = entry->la;
1292 for (d = 0; d < DIM; d++)
1294 if (zones->shift[entry->cell][d] == 0)
1306 (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1311 for (d = 0; (d < DIM && bUse); d++)
1313 /* Check if the cg_cm distance falls within
1314 * the cut-off to avoid possible multiple
1315 * assignments of bonded interactions.
1319 dd_dist2(pbc_null, cg_cm,
1320 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1329 /* Add this interaction to the local topology */
1330 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1331 /* Sum so we can check in global_stat
1332 * if we have everything.
1335 !(interaction_function[ftype].flags & IF_LIMZERO))
1345 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1347 * With thread parallelizing each thread acts on a different atom range:
1348 * at_start to at_end.
1350 static int make_bondeds_zone(gmx_domdec_t *dd,
1351 const gmx_domdec_zones_t *zones,
1352 const std::vector<gmx_molblock_t> &molb,
1353 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1355 t_pbc *pbc_null, rvec *cg_cm,
1356 const t_iparams *ip_in,
1359 gmx::RangePartitioning::Block atomRange)
1361 int mb, mt, mol, i_mol;
1363 gmx_reverse_top_t *rt;
1366 rt = dd->reverse_top;
1368 bBCheck = rt->bBCheck;
1372 for (int i : atomRange)
1374 /* Get the global atom number */
1375 const int i_gl = dd->globalAtomIndices[i];
1376 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1377 /* Check all intramolecular interactions assigned to this atom */
1378 gmx::ArrayRef<const int> index = rt->ril_mt[mt].index;
1379 gmx::ArrayRef<const t_iatom> rtil = rt->ril_mt[mt].il;
1381 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1382 rt->ril_mt[mt].numAtomsInMolecule,
1384 index[i_mol], index[i_mol+1],
1387 bRCheckMB, rcheck, bRCheck2B, rc2,
1397 if (rt->bIntermolecularInteractions)
1399 /* Check all intermolecular interactions assigned to this atom */
1400 index = rt->ril_intermol.index;
1401 rtil = rt->ril_intermol.il;
1403 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1404 rt->ril_mt[mt].numAtomsInMolecule,
1406 index[i_gl], index[i_gl + 1],
1409 bRCheckMB, rcheck, bRCheck2B, rc2,
1420 return nbonded_local;
1423 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1424 static void set_no_exclusions_zone(const gmx_domdec_zones_t *zones,
1428 for (int a = zones->cg_range[iz]; a < zones->cg_range[iz + 1]; a++)
1430 lexcls->index[a + 1] = lexcls->nra;
1434 /*! \brief Set the exclusion data for i-zone \p iz */
1435 static void make_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1436 const std::vector<gmx_moltype_t> &moltype,
1437 const int *cginfo, t_blocka *lexcls, int iz,
1438 int at_start, int at_end,
1439 const gmx::ArrayRef<const int> intermolecularExclusionGroup)
1441 int n_excl_at_max, n, at;
1443 const gmx_ga2la_t &ga2la = *dd->ga2la;
1445 // TODO: Replace this by a more standard range
1446 const gmx::RangePartitioning::Block jRange(zones->izone[iz].jcg0,
1447 zones->izone[iz].jcg1);
1449 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1451 /* We set the end index, but note that we might not start at zero here */
1452 lexcls->nr = at_end;
1455 for (at = at_start; at < at_end; at++)
1457 if (n + 1000 > lexcls->nalloc_a)
1459 lexcls->nalloc_a = over_alloc_large(n + 1000);
1460 srenew(lexcls->a, lexcls->nalloc_a);
1463 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1465 int a_gl, mb, mt, mol, a_mol, j;
1466 const t_blocka *excls;
1468 if (n + n_excl_at_max > lexcls->nalloc_a)
1470 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1471 srenew(lexcls->a, lexcls->nalloc_a);
1474 /* Copy the exclusions from the global top */
1475 lexcls->index[at] = n;
1476 a_gl = dd->globalAtomIndices[at];
1477 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1478 &mb, &mt, &mol, &a_mol);
1479 excls = &moltype[mt].excls;
1480 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1482 const int aj_mol = excls->a[j];
1484 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1486 /* This check is not necessary, but it can reduce
1487 * the number of exclusions in the list, which in turn
1488 * can speed up the pair list construction a bit.
1490 if (jRange.inRange(jEntry->la))
1492 lexcls->a[n++] = jEntry->la;
1499 /* We don't need exclusions for this atom */
1500 lexcls->index[at] = n;
1503 bool isExcludedAtom = !intermolecularExclusionGroup.empty() &&
1504 std::find(intermolecularExclusionGroup.begin(),
1505 intermolecularExclusionGroup.end(),
1506 dd->globalAtomIndices[at]) !=
1507 intermolecularExclusionGroup.end();
1511 if (n + intermolecularExclusionGroup.ssize() > lexcls->nalloc_a)
1514 over_alloc_large(n + intermolecularExclusionGroup.size());
1515 srenew(lexcls->a, lexcls->nalloc_a);
1517 for (int qmAtomGlobalIndex : intermolecularExclusionGroup)
1519 if (const auto *entry = dd->ga2la->find(qmAtomGlobalIndex))
1521 lexcls->a[n++] = entry->la;
1527 lexcls->index[lexcls->nr] = n;
1532 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1533 static void check_alloc_index(t_blocka *ba, int nindex_max)
1535 if (nindex_max+1 > ba->nalloc_index)
1537 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1538 srenew(ba->index, ba->nalloc_index);
1542 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1543 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1546 const int nr = zones->izone[zones->nizone - 1].cg1;
1548 check_alloc_index(lexcls, nr);
1550 for (size_t thread = 1; thread < dd->reverse_top->th_work.size(); thread++)
1552 check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1556 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1557 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1560 // TODO: Replace this by a more standard range
1561 const gmx::RangePartitioning::Block nonhomeIzonesAtomRange(zones->izone[0].cg1,
1562 zones->izone[zones->nizone - 1].cg1);
1564 if (dd->n_intercg_excl == 0)
1566 /* There are no exclusions involving non-home charge groups,
1567 * but we need to set the indices for neighborsearching.
1569 for (int la : nonhomeIzonesAtomRange)
1571 lexcls->index[la] = lexcls->nra;
1574 /* nr is only used to loop over the exclusions for Ewald and RF,
1575 * so we can set it to the number of home atoms for efficiency.
1577 lexcls->nr = nonhomeIzonesAtomRange.begin();
1581 lexcls->nr = nonhomeIzonesAtomRange.end();
1585 /*! \brief Clear a t_idef struct */
1586 static void clear_idef(t_idef *idef)
1590 /* Clear the counts */
1591 for (ftype = 0; ftype < F_NRE; ftype++)
1593 idef->il[ftype].nr = 0;
1597 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1598 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1599 gmx_domdec_zones_t *zones,
1600 const gmx_mtop_t *mtop,
1602 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1604 t_pbc *pbc_null, rvec *cg_cm,
1606 t_blocka *lexcls, int *excl_count)
1608 int nzone_bondeds, nzone_excl;
1609 int izone, cg0, cg1;
1612 gmx_reverse_top_t *rt;
1614 if (dd->reverse_top->bInterCGInteractions)
1616 nzone_bondeds = zones->n;
1620 /* Only single charge group (or atom) molecules, so interactions don't
1621 * cross zone boundaries and we only need to assign in the home zone.
1626 if (dd->n_intercg_excl > 0)
1628 /* We only use exclusions from i-zones to i- and j-zones */
1629 nzone_excl = zones->nizone;
1633 /* There are no inter-cg exclusions and only zone 0 sees itself */
1637 check_exclusions_alloc(dd, zones, lexcls);
1639 rt = dd->reverse_top;
1643 /* Clear the counts */
1651 for (izone = 0; izone < nzone_bondeds; izone++)
1653 cg0 = zones->cg_range[izone];
1654 cg1 = zones->cg_range[izone + 1];
1656 const int numThreads = rt->th_work.size();
1657 #pragma omp parallel for num_threads(numThreads) schedule(static)
1658 for (int thread = 0; thread < numThreads; thread++)
1666 cg0t = cg0 + ((cg1 - cg0)* thread )/numThreads;
1667 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/numThreads;
1675 idef_t = &rt->th_work[thread].idef;
1679 rt->th_work[thread].nbonded =
1680 make_bondeds_zone(dd, zones,
1682 bRCheckMB, rcheck, bRCheck2B, rc2,
1683 pbc_null, cg_cm, idef->iparams,
1686 gmx::RangePartitioning::Block(cg0t, cg1t));
1688 if (izone < nzone_excl)
1696 excl_t = &rt->th_work[thread].excl;
1701 /* No charge groups and no distance check required */
1702 make_exclusions_zone(dd, zones, mtop->moltype, cginfo,
1703 excl_t, izone, cg0t,
1705 mtop->intermolecularExclusionGroup);
1708 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1711 if (rt->th_work.size() > 1)
1713 combine_idef(idef, rt->th_work);
1716 for (const thread_work_t &th_work : rt->th_work)
1718 nbonded_local += th_work.nbonded;
1721 if (izone < nzone_excl)
1723 if (rt->th_work.size() > 1)
1725 combine_blocka(lexcls, rt->th_work);
1728 for (const thread_work_t &th_work : rt->th_work)
1730 *excl_count += th_work.excl_count;
1735 /* Some zones might not have exclusions, but some code still needs to
1736 * loop over the index, so we set the indices here.
1738 for (izone = nzone_excl; izone < zones->nizone; izone++)
1740 set_no_exclusions_zone(zones, izone, lexcls);
1743 finish_local_exclusions(dd, zones, lexcls);
1746 fprintf(debug, "We have %d exclusions, check count %d\n",
1747 lexcls->nra, *excl_count);
1750 return nbonded_local;
1753 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1754 int npbcdim, matrix box,
1755 rvec cellsize_min, const ivec npulse,
1758 const gmx_mtop_t &mtop, gmx_localtop_t *ltop)
1760 gmx_bool bRCheckMB, bRCheck2B;
1764 t_pbc pbc, *pbc_null = nullptr;
1768 fprintf(debug, "Making local topology\n");
1774 if (dd->reverse_top->bInterCGInteractions)
1776 /* We need to check to which cell bondeds should be assigned */
1777 rc = dd_cutoff_twobody(dd);
1780 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1783 /* Should we check cg_cm distances when assigning bonded interactions? */
1784 for (d = 0; d < DIM; d++)
1787 /* Only need to check for dimensions where the part of the box
1788 * that is not communicated is smaller than the cut-off.
1790 if (d < npbcdim && dd->nc[d] > 1 &&
1791 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
1798 /* Check for interactions between two atoms,
1799 * where we can allow interactions up to the cut-off,
1800 * instead of up to the smallest cell dimension.
1807 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n",
1808 d, cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
1811 if (bRCheckMB || bRCheck2B)
1815 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
1825 make_local_bondeds_excls(dd, zones, &mtop, fr->cginfo.data(),
1826 bRCheckMB, rcheck, bRCheck2B, rc,
1827 pbc_null, cgcm_or_x,
1829 <op->excls, &nexcl);
1831 /* The ilist is not sorted yet,
1832 * we can only do this when we have the charge arrays.
1834 ltop->idef.ilsort = ilsortUNKNOWN;
1836 ltop->atomtypes = mtop.atomtypes;
1839 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
1840 gmx_localtop_t *ltop)
1842 if (dd->reverse_top->ilsort == ilsortNO_FE)
1844 ltop->idef.ilsort = ilsortNO_FE;
1848 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
1852 void dd_init_local_top(const gmx_mtop_t &top_global,
1853 gmx_localtop_t *top)
1855 /* TODO: Get rid of the const casts below, e.g. by using a reference */
1856 top->idef.ntypes = top_global.ffparams.numTypes();
1857 top->idef.atnr = top_global.ffparams.atnr;
1858 top->idef.functype = const_cast<t_functype *>(top_global.ffparams.functype.data());
1859 top->idef.iparams = const_cast<t_iparams *>(top_global.ffparams.iparams.data());
1860 top->idef.fudgeQQ = top_global.ffparams.fudgeQQ;
1861 top->idef.cmap_grid = new gmx_cmap_t;
1862 *top->idef.cmap_grid = top_global.ffparams.cmap_grid;
1864 top->idef.ilsort = ilsortUNKNOWN;
1865 top->useInDomainDecomp_ = true;
1868 void dd_init_local_state(gmx_domdec_t *dd,
1869 const t_state *state_global, t_state *state_local)
1871 int buf[NITEM_DD_INIT_LOCAL_STATE];
1875 buf[0] = state_global->flags;
1876 buf[1] = state_global->ngtc;
1877 buf[2] = state_global->nnhpres;
1878 buf[3] = state_global->nhchainlength;
1879 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
1881 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
1883 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
1884 init_dfhist_state(state_local, buf[4]);
1885 state_local->flags = buf[0];
1888 /*! \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 */
1889 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
1895 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
1897 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
1898 if (link->a[k] == cg_gl_j)
1905 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
1906 "Inconsistent allocation of link");
1907 /* Add this charge group link */
1908 if (link->index[cg_gl+1]+1 > link->nalloc_a)
1910 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
1911 srenew(link->a, link->nalloc_a);
1913 link->a[link->index[cg_gl+1]] = cg_gl_j;
1914 link->index[cg_gl+1]++;
1918 /*! \brief Return a vector of the charge group index for all atoms */
1919 static std::vector<int> make_at2cg(const t_block &cgs)
1921 std::vector<int> at2cg(cgs.index[cgs.nr]);
1922 for (int cg = 0; cg < cgs.nr; cg++)
1924 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
1933 t_blocka *makeBondedLinks(const gmx_mtop_t *mtop,
1934 cginfo_mb_t *cginfo_mb)
1937 cginfo_mb_t *cgi_mb;
1939 /* For each charge group make a list of other charge groups
1940 * in the system that a linked to it via bonded interactions
1941 * which are also stored in reverse_top.
1944 reverse_ilist_t ril_intermol;
1945 if (mtop->bIntermolecularInteractions)
1947 if (ncg_mtop(mtop) < mtop->natoms)
1949 gmx_fatal(FARGS, "The combination of intermolecular interactions, charge groups and domain decomposition is not supported. Use cutoff-scheme=Verlet (which removes the charge groups) or run without domain decomposition.");
1954 atoms.nr = mtop->natoms;
1955 atoms.atom = nullptr;
1957 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
1959 make_reverse_ilist(*mtop->intermolecular_ilist,
1961 FALSE, FALSE, FALSE, TRUE, &ril_intermol);
1965 snew(link->index, ncg_mtop(mtop)+1);
1972 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
1974 const gmx_molblock_t &molb = mtop->molblock[mb];
1979 const gmx_moltype_t &molt = mtop->moltype[molb.type];
1980 const t_block &cgs = molt.cgs;
1981 std::vector<int> a2c = make_at2cg(cgs);
1982 /* Make a reverse ilist in which the interactions are linked
1983 * to all atoms, not only the first atom as in gmx_reverse_top.
1984 * The constraints are discarded here.
1986 reverse_ilist_t ril;
1987 make_reverse_ilist(molt.ilist, &molt.atoms,
1988 FALSE, FALSE, FALSE, TRUE, &ril);
1990 cgi_mb = &cginfo_mb[mb];
1993 for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb.nmol : 1); mol++)
1995 for (int cg = 0; cg < cgs.nr; cg++)
1997 int cg_gl = cg_offset + cg;
1998 link->index[cg_gl+1] = link->index[cg_gl];
1999 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2001 int i = ril.index[a];
2002 while (i < ril.index[a+1])
2004 int ftype = ril.il[i++];
2005 int nral = NRAL(ftype);
2006 /* Skip the ifunc index */
2008 for (int j = 0; j < nral; j++)
2010 int aj = ril.il[i + j];
2013 check_link(link, cg_gl, cg_offset+a2c[aj]);
2016 i += nral_rt(ftype);
2019 if (mtop->bIntermolecularInteractions)
2021 int i = ril_intermol.index[a];
2022 while (i < ril_intermol.index[a+1])
2024 int ftype = ril_intermol.il[i++];
2025 int nral = NRAL(ftype);
2026 /* Skip the ifunc index */
2028 for (int j = 0; j < nral; j++)
2030 /* Here we assume we have no charge groups;
2031 * this has been checked above.
2033 int aj = ril_intermol.il[i + j];
2034 check_link(link, cg_gl, aj);
2036 i += nral_rt(ftype);
2040 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2042 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2047 cg_offset += cgs.nr;
2049 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - cgs.nr];
2053 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt.name, cgs.nr, nlink_mol);
2056 if (molb.nmol > mol)
2058 /* Copy the data for the rest of the molecules in this block */
2059 link->nalloc_a += (molb.nmol - mol)*nlink_mol;
2060 srenew(link->a, link->nalloc_a);
2061 for (; mol < molb.nmol; mol++)
2063 for (int cg = 0; cg < cgs.nr; cg++)
2065 int cg_gl = cg_offset + cg;
2066 link->index[cg_gl + 1] =
2067 link->index[cg_gl + 1 - cgs.nr] + nlink_mol;
2068 for (int j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2070 link->a[j] = link->a[j - nlink_mol] + cgs.nr;
2072 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2073 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2075 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2079 cg_offset += cgs.nr;
2086 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2097 } bonded_distance_t;
2099 /*! \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 */
2100 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2101 bonded_distance_t *bd)
2112 /*! \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 */
2113 static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
2114 const std::vector<int> &at2cg,
2115 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2116 bonded_distance_t *bd_2b,
2117 bonded_distance_t *bd_mb)
2119 for (int ftype = 0; ftype < F_NRE; ftype++)
2121 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2123 const auto &il = molt->ilist[ftype];
2124 int nral = NRAL(ftype);
2127 for (int i = 0; i < il.size(); i += 1+nral)
2129 for (int ai = 0; ai < nral; ai++)
2131 int cgi = at2cg[il.iatoms[i+1+ai]];
2132 for (int aj = ai + 1; aj < nral; aj++)
2134 int cgj = at2cg[il.iatoms[i+1+aj]];
2137 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2139 update_max_bonded_distance(rij2, ftype,
2142 (nral == 2) ? bd_2b : bd_mb);
2152 const t_blocka *excls = &molt->excls;
2153 for (int ai = 0; ai < excls->nr; ai++)
2155 int cgi = at2cg[ai];
2156 for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2158 int cgj = at2cg[excls->a[j]];
2161 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2163 /* There is no function type for exclusions, use -1 */
2164 update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2171 /*! \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 */
2172 static void bonded_distance_intermol(const InteractionLists &ilists_intermol,
2174 const rvec *x, int ePBC, const matrix box,
2175 bonded_distance_t *bd_2b,
2176 bonded_distance_t *bd_mb)
2180 set_pbc(&pbc, ePBC, box);
2182 for (int ftype = 0; ftype < F_NRE; ftype++)
2184 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2186 const auto &il = ilists_intermol[ftype];
2187 int nral = NRAL(ftype);
2189 /* No nral>1 check here, since intermol interactions always
2190 * have nral>=2 (and the code is also correct for nral=1).
2192 for (int i = 0; i < il.size(); i += 1+nral)
2194 for (int ai = 0; ai < nral; ai++)
2196 int atom_i = il.iatoms[i + 1 + ai];
2198 for (int aj = ai + 1; aj < nral; aj++)
2203 int atom_j = il.iatoms[i + 1 + aj];
2205 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2209 update_max_bonded_distance(rij2, ftype,
2211 (nral == 2) ? bd_2b : bd_mb);
2219 //! Returns whether \p molt has at least one virtual site
2220 static bool moltypeHasVsite(const gmx_moltype_t &molt)
2222 bool hasVsite = false;
2223 for (int i = 0; i < F_NRE; i++)
2225 if ((interaction_function[i].flags & IF_VSITE) &&
2226 molt.ilist[i].size() > 0)
2235 //! Compute charge group centers of mass for molecule \p molt
2236 static void get_cgcm_mol(const gmx_moltype_t *molt,
2237 const gmx_ffparams_t *ffparams,
2238 int ePBC, t_graph *graph, const matrix box,
2239 const rvec *x, rvec *xs, rvec *cg_cm)
2243 if (ePBC != epbcNONE)
2245 mk_mshift(nullptr, graph, ePBC, box, x);
2247 shift_x(graph, box, x, xs);
2248 /* By doing an extra mk_mshift the molecules that are broken
2249 * because they were e.g. imported from another software
2250 * will be made whole again. Such are the healing powers
2253 mk_mshift(nullptr, graph, ePBC, box, xs);
2257 /* We copy the coordinates so the original coordinates remain
2258 * unchanged, just to be 100% sure that we do not affect
2259 * binary reproducibility of simulations.
2261 n = molt->cgs.index[molt->cgs.nr];
2262 for (i = 0; i < n; i++)
2264 copy_rvec(x[i], xs[i]);
2268 if (moltypeHasVsite(*molt))
2270 /* Convert to old, deprecated format */
2271 t_ilist ilist[F_NRE];
2272 for (int ftype = 0; ftype < F_NRE; ftype++)
2274 if (interaction_function[ftype].flags & IF_VSITE)
2276 ilist[ftype].nr = molt->ilist[ftype].size();
2277 ilist[ftype].iatoms = const_cast<int *>(molt->ilist[ftype].iatoms.data());
2281 construct_vsites(nullptr, xs, 0.0, nullptr,
2282 ffparams->iparams.data(), ilist,
2283 epbcNONE, TRUE, nullptr, nullptr);
2286 calc_cgcm(nullptr, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2289 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
2290 const gmx_mtop_t *mtop,
2291 const t_inputrec *ir,
2292 const rvec *x, const matrix box,
2294 real *r_2b, real *r_mb)
2296 gmx_bool bExclRequired;
2300 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2301 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2303 bExclRequired = inputrecExclForces(ir);
2308 for (const gmx_molblock_t &molb : mtop->molblock)
2310 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2311 if (molt.cgs.nr == 1 || molb.nmol == 0)
2313 at_offset += molb.nmol*molt.atoms.nr;
2317 if (ir->ePBC != epbcNONE)
2319 mk_graph_moltype(molt, &graph);
2322 std::vector<int> at2cg = make_at2cg(molt.cgs);
2323 snew(xs, molt.atoms.nr);
2324 snew(cg_cm, molt.cgs.nr);
2325 for (int mol = 0; mol < molb.nmol; mol++)
2327 get_cgcm_mol(&molt, &mtop->ffparams, ir->ePBC, &graph, box,
2328 x+at_offset, xs, cg_cm);
2330 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2331 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2333 bonded_cg_distance_mol(&molt, at2cg, bBCheck, bExclRequired, cg_cm,
2334 &bd_mol_2b, &bd_mol_mb);
2336 /* Process the mol data adding the atom index offset */
2337 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2338 at_offset + bd_mol_2b.a1,
2339 at_offset + bd_mol_2b.a2,
2341 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2342 at_offset + bd_mol_mb.a1,
2343 at_offset + bd_mol_mb.a2,
2346 at_offset += molt.atoms.nr;
2350 if (ir->ePBC != epbcNONE)
2357 if (mtop->bIntermolecularInteractions)
2359 if (ncg_mtop(mtop) < mtop->natoms)
2361 gmx_fatal(FARGS, "The combination of intermolecular interactions, charge groups and domain decomposition is not supported. Use cutoff-scheme=Verlet (which removes the charge groups) or run without domain decomposition.");
2364 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
2366 bonded_distance_intermol(*mtop->intermolecular_ilist,
2372 *r_2b = sqrt(bd_2b.r2);
2373 *r_mb = sqrt(bd_mb.r2);
2375 if (*r_2b > 0 || *r_mb > 0)
2377 GMX_LOG(mdlog.info).appendText("Initial maximum distances in bonded interactions:");
2380 GMX_LOG(mdlog.info).appendTextFormatted(
2381 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2382 *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2383 bd_2b.a1 + 1, bd_2b.a2 + 1);
2387 GMX_LOG(mdlog.info).appendTextFormatted(
2388 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2389 *r_mb, interaction_function[bd_mb.ftype].longname,
2390 bd_mb.a1 + 1, bd_mb.a2 + 1);