2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2011,2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file defines functions used by the domdec module
39 * while managing the construction, use and error checking for
40 * topologies local to a DD rank.
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_domdec
55 #include "gromacs/compat/make_unique.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/ga2la.h"
59 #include "gromacs/gmxlib/chargegroup.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/forcerec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/vsite.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/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 Do we require all exclusions to be assigned?
124 bool bExclRequired = false;
125 //! \brief The maximum number of exclusions one atom can have
126 int n_excl_at_max = 0;
127 //! \brief Are there constraints in this revserse top?
128 bool bConstr = false;
129 //! \brief Are there settles in this revserse top?
130 bool bSettle = false;
131 //! \brief All bonded interactions have to be assigned?
132 bool bBCheck = false;
133 //! \brief Are there bondeds/exclusions between charge-groups?
134 bool bInterCGInteractions = false;
135 //! \brief Reverse ilist for all moltypes
136 std::vector<reverse_ilist_t> ril_mt;
137 //! \brief The size of ril_mt[?].index summed over all entries
138 int ril_mt_tot_size = 0;
139 //! \brief The sorting state of bondeds for free energy
140 int ilsort = ilsortUNKNOWN;
141 //! \brief molblock to global atom index for quick lookup of molblocks on atom index
142 std::vector<MolblockIndices> mbi;
144 //! \brief Do we have intermolecular interactions?
145 bool bIntermolecularInteractions = false;
146 //! \brief Intermolecular reverse ilist
147 reverse_ilist_t ril_intermol;
149 /* Work data structures for multi-threading */
150 //! \brief Thread work array for local topology generation
151 std::vector<thread_work_t> th_work;
155 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
156 static int nral_rt(int ftype)
161 if (interaction_function[ftype].flags & IF_VSITE)
163 /* With vsites the reverse topology contains
164 * two extra entries for PBC.
172 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
173 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
174 gmx_bool bConstr, gmx_bool bSettle)
176 return ((((interaction_function[ftype].flags & IF_BOND) != 0u) &&
177 ((interaction_function[ftype].flags & IF_VSITE) == 0u) &&
178 (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0u))) ||
179 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
180 (bSettle && ftype == F_SETTLE));
183 /*! \brief Help print error output when interactions are missing */
185 print_missing_interactions_mb(t_commrec *cr,
186 const gmx_reverse_top_t *rt,
187 const char *moltypename,
188 const reverse_ilist_t *ril,
189 int a_start, int a_end,
190 int nat_mol, int nmol,
194 int nril_mol = ril->index[nat_mol];
195 snew(assigned, nmol*nril_mol);
196 gmx::StringOutputStream stream;
197 gmx::TextWriter log(&stream);
199 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
200 for (int ftype = 0; ftype < F_NRE; ftype++)
202 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
204 int nral = NRAL(ftype);
205 const t_ilist *il = &idef->il[ftype];
206 const t_iatom *ia = il->iatoms;
207 for (int i = 0; i < il->nr; i += 1+nral)
209 int a0 = gatindex[ia[1]];
210 /* Check if this interaction is in
211 * the currently checked molblock.
213 if (a0 >= a_start && a0 < a_end)
215 int mol = (a0 - a_start)/nat_mol;
216 int a0_mol = (a0 - a_start) - mol*nat_mol;
217 int j_mol = ril->index[a0_mol];
219 while (j_mol < ril->index[a0_mol+1] && !found)
221 int j = mol*nril_mol + j_mol;
222 int ftype_j = ril->il[j_mol];
223 /* Here we need to check if this interaction has
224 * not already been assigned, since we could have
225 * multiply defined interactions.
227 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
230 /* Check the atoms */
232 for (int a = 0; a < nral; a++)
234 if (gatindex[ia[1+a]] !=
235 a_start + mol*nat_mol + ril->il[j_mol+2+a])
245 j_mol += 2 + nral_rt(ftype_j);
249 gmx_incons("Some interactions seem to be assigned multiple times");
257 gmx_sumi(nmol*nril_mol, assigned, cr);
261 for (int mol = 0; mol < nmol; mol++)
264 while (j_mol < nril_mol)
266 int ftype = ril->il[j_mol];
267 int nral = NRAL(ftype);
268 int j = mol*nril_mol + j_mol;
269 if (assigned[j] == 0 &&
270 !(interaction_function[ftype].flags & IF_VSITE))
272 if (DDMASTER(cr->dd))
276 log.writeLineFormatted("Molecule type '%s'", moltypename);
277 log.writeLineFormatted(
278 "the first %d missing interactions, except for exclusions:", nprint);
280 log.writeStringFormatted("%20s atoms",
281 interaction_function[ftype].longname);
283 for (a = 0; a < nral; a++)
285 log.writeStringFormatted("%5d", ril->il[j_mol+2+a]+1);
289 log.writeString(" ");
292 log.writeString(" global");
293 for (a = 0; a < nral; a++)
295 log.writeStringFormatted("%6d",
296 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
298 log.ensureLineBreak();
306 j_mol += 2 + nral_rt(ftype);
311 return stream.toString();
314 /*! \brief Help print error output when interactions are missing */
315 static void print_missing_interactions_atoms(const gmx::MDLogger &mdlog,
317 const gmx_mtop_t *mtop,
320 const gmx_reverse_top_t *rt = cr->dd->reverse_top;
322 /* Print the atoms in the missing interactions per molblock */
324 for (const gmx_molblock_t &molb : mtop->molblock)
326 const gmx_moltype_t &moltype = mtop->moltype[molb.type];
328 a_end = a_start + molb.nmol*moltype.atoms.nr;
330 GMX_LOG(mdlog.warning).appendText(
331 print_missing_interactions_mb(cr, rt,
333 &rt->ril_mt[molb.type],
334 a_start, a_end, moltype.atoms.nr,
340 void dd_print_missing_interactions(const gmx::MDLogger &mdlog,
343 const gmx_mtop_t *top_global,
344 const gmx_localtop_t *top_local,
345 const t_state *state_local)
347 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
353 GMX_LOG(mdlog.warning).appendText(
354 "Not all bonded interactions have been properly assigned to the domain decomposition cells");
356 ndiff_tot = local_count - dd->nbonded_global;
358 for (ftype = 0; ftype < F_NRE; ftype++)
361 cl[ftype] = top_local->idef.il[ftype].nr/(1+nral);
364 gmx_sumi(F_NRE, cl, cr);
368 GMX_LOG(mdlog.warning).appendText("A list of missing interactions:");
369 rest_global = dd->nbonded_global;
370 rest_local = local_count;
371 for (ftype = 0; ftype < F_NRE; ftype++)
373 /* In the reverse and local top all constraints are merged
374 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
375 * and add these constraints when doing F_CONSTR.
377 if (((interaction_function[ftype].flags & IF_BOND) &&
378 (dd->reverse_top->bBCheck
379 || !(interaction_function[ftype].flags & IF_LIMZERO)))
380 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
381 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
383 n = gmx_mtop_ftype_count(top_global, ftype);
384 if (ftype == F_CONSTR)
386 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
388 ndiff = cl[ftype] - n;
391 GMX_LOG(mdlog.warning).appendTextFormatted(
392 "%20s of %6d missing %6d",
393 interaction_function[ftype].longname, n, -ndiff);
396 rest_local -= cl[ftype];
400 ndiff = rest_local - rest_global;
403 GMX_LOG(mdlog.warning).appendTextFormatted(
404 "%20s of %6d missing %6d", "exclusions",
405 rest_global, -ndiff);
409 print_missing_interactions_atoms(mdlog, cr, top_global, &top_local->idef);
410 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
411 -1, state_local->x.rvec_array(), state_local->box);
413 std::string errorMessage;
417 errorMessage = "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
421 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));
423 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
426 /*! \brief Return global topology molecule information for global atom index \p i_gl */
427 static void global_atomnr_to_moltype_ind(const gmx_reverse_top_t *rt,
429 int *mb, int *mt, int *mol, int *i_mol)
431 const MolblockIndices *mbi = rt->mbi.data();
433 int end = rt->mbi.size(); /* exclusive */
436 /* binary search for molblock_ind */
440 if (i_gl >= mbi[mid].a_end)
444 else if (i_gl < mbi[mid].a_start)
458 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
459 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
462 /*! \brief Count the exclusions for all atoms in \p cgs */
463 static void count_excls(const t_block *cgs, const t_blocka *excls,
464 int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
466 int cg, at0, at1, at, excl, atj;
471 for (cg = 0; cg < cgs->nr; cg++)
473 at0 = cgs->index[cg];
474 at1 = cgs->index[cg+1];
475 for (at = at0; at < at1; at++)
477 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
479 atj = excls->a[excl];
483 if (atj < at0 || atj >= at1)
490 *n_excl_at_max = std::max(*n_excl_at_max,
491 excls->index[at+1] - excls->index[at]);
496 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
497 static int low_make_reverse_ilist(const InteractionLists &il_mt,
499 gmx::ArrayRef < const std::vector < int>> vsitePbc,
501 gmx_bool bConstr, gmx_bool bSettle,
503 gmx::ArrayRef<const int> r_index,
504 gmx::ArrayRef<int> r_il,
505 gmx_bool bLinkToAllAtoms,
508 int ftype, j, nlink, link;
513 for (ftype = 0; ftype < F_NRE; ftype++)
515 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
516 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
517 (bSettle && ftype == F_SETTLE))
519 const bool bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0u);
520 const int nral = NRAL(ftype);
521 const auto &il = il_mt[ftype];
522 for (int i = 0; i < il.size(); i += 1+nral)
524 const int* ia = il.iatoms.data() + i;
529 /* We don't need the virtual sites for the cg-links */
539 /* Couple to the first atom in the interaction */
542 for (link = 0; link < nlink; link++)
547 GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
548 GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
549 r_il[r_index[a]+count[a]] =
550 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
551 r_il[r_index[a]+count[a]+1] = ia[0];
552 for (j = 1; j < 1+nral; j++)
554 /* Store the molecular atom number */
555 r_il[r_index[a]+count[a]+1+j] = ia[j];
558 if (interaction_function[ftype].flags & IF_VSITE)
562 /* Add an entry to iatoms for storing
563 * which of the constructing atoms are
566 r_il[r_index[a]+count[a]+2+nral] = 0;
567 for (j = 2; j < 1+nral; j++)
569 if (atom[ia[j]].ptype == eptVSite)
571 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
574 /* Store vsite pbc atom in a second extra entry */
575 r_il[r_index[a]+count[a]+2+nral+1] =
576 (!vsitePbc.empty() ? vsitePbc[ftype-F_VSITE2][i/(1+nral)] : -2);
581 /* We do not count vsites since they are always
582 * uniquely assigned and can be assigned
583 * to multiple nodes with recursive vsites.
586 !(interaction_function[ftype].flags & IF_LIMZERO))
591 count[a] += 2 + nral_rt(ftype);
600 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
601 static int make_reverse_ilist(const InteractionLists &ilist,
602 const t_atoms *atoms,
603 gmx::ArrayRef < const std::vector < int>> vsitePbc,
604 gmx_bool bConstr, gmx_bool bSettle,
606 gmx_bool bLinkToAllAtoms,
607 reverse_ilist_t *ril_mt)
609 int nat_mt, *count, i, nint_mt;
611 /* Count the interactions */
614 low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
616 bConstr, bSettle, bBCheck,
617 gmx::EmptyArrayRef(), gmx::EmptyArrayRef(),
618 bLinkToAllAtoms, FALSE);
620 ril_mt->index.push_back(0);
621 for (i = 0; i < nat_mt; i++)
623 ril_mt->index.push_back(ril_mt->index[i] + count[i]);
626 ril_mt->il.resize(ril_mt->index[nat_mt]);
628 /* Store the interactions */
630 low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
632 bConstr, bSettle, bBCheck,
633 ril_mt->index, ril_mt->il,
634 bLinkToAllAtoms, TRUE);
638 ril_mt->numAtomsInMolecule = atoms->nr;
643 /*! \brief Generate the reverse topology */
644 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
645 gmx::ArrayRef<const VsitePbc> vsitePbcPerMoltype,
646 gmx_bool bConstr, gmx_bool bSettle,
647 gmx_bool bBCheck, int *nint)
649 gmx_reverse_top_t rt;
651 /* Should we include constraints (for SHAKE) in rt? */
652 rt.bConstr = bConstr;
653 rt.bSettle = bSettle;
654 rt.bBCheck = bBCheck;
656 rt.bInterCGInteractions = mtop->bIntermolecularInteractions;
657 rt.ril_mt.resize(mtop->moltype.size());
658 rt.ril_mt_tot_size = 0;
659 std::vector<int> nint_mt;
660 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
662 const gmx_moltype_t &molt = mtop->moltype[mt];
665 rt.bInterCGInteractions = true;
668 /* Make the atom to interaction list for this molecule type */
669 gmx::ArrayRef < const std::vector < int>> vsitePbc;
670 if (!vsitePbcPerMoltype.empty())
672 vsitePbc = gmx::makeConstArrayRef(vsitePbcPerMoltype[mt]);
674 int numberOfInteractions =
675 make_reverse_ilist(molt.ilist, &molt.atoms, vsitePbc,
676 rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
678 nint_mt.push_back(numberOfInteractions);
680 rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
684 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt.ril_mt_tot_size);
688 for (const gmx_molblock_t &molblock : mtop->molblock)
690 *nint += molblock.nmol*nint_mt[molblock.type];
693 /* Make an intermolecular reverse top, if necessary */
694 rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
695 if (rt.bIntermolecularInteractions)
697 t_atoms atoms_global;
699 atoms_global.nr = mtop->natoms;
700 atoms_global.atom = nullptr; /* Only used with virtual sites */
702 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
705 make_reverse_ilist(*mtop->intermolecular_ilist,
707 gmx::EmptyArrayRef(),
708 rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
712 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
714 rt.ilsort = ilsortFE_UNSORTED;
718 rt.ilsort = ilsortNO_FE;
721 /* Make a molblock index for fast searching */
723 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
725 const gmx_molblock_t &molb = mtop->molblock[mb];
726 const int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
729 i += molb.nmol*numAtomsPerMol;
731 mbi.natoms_mol = numAtomsPerMol;
732 mbi.type = molb.type;
733 rt.mbi.push_back(mbi);
736 rt.th_work.resize(gmx_omp_nthreads_get(emntDomdec));
737 if (!vsitePbcPerMoltype.empty())
739 for (thread_work_t &th_work : rt.th_work)
741 th_work.vsitePbc = gmx::compat::make_unique<VsitePbc>();
748 void dd_make_reverse_top(FILE *fplog,
749 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
750 const gmx_vsite_t *vsite,
751 const t_inputrec *ir, gmx_bool bBCheck)
755 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
758 /* If normal and/or settle constraints act only within charge groups,
759 * we can store them in the reverse top and simply assign them to domains.
760 * Otherwise we need to assign them to multiple domains and set up
761 * the parallel version constraint algorithm(s).
764 gmx::ArrayRef<const VsitePbc> vsitePbcPerMoltype;
767 vsitePbcPerMoltype = gmx::makeConstArrayRef(vsite->vsite_pbc_molt);
770 dd->reverse_top = new gmx_reverse_top_t;
772 make_reverse_top(mtop, ir->efep != efepNO, vsitePbcPerMoltype,
773 !dd->splitConstraints, !dd->splitSettles,
774 bBCheck, &dd->nbonded_global);
776 gmx_reverse_top_t *rt = dd->reverse_top;
778 /* With the Verlet scheme, exclusions are handled in the non-bonded
779 * kernels and only exclusions inside the cut-off lead to exclusion
780 * forces. Since each atom pair is treated at most once in the non-bonded
781 * kernels, it doesn't matter if the exclusions for the same atom pair
782 * appear multiple times in the exclusion list. In contrast, the, old,
783 * group cut-off scheme loops over a list of exclusions, so there each
784 * excluded pair should appear exactly once.
786 rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
787 inputrecExclForces(ir));
790 dd->n_intercg_excl = 0;
791 rt->n_excl_at_max = 0;
792 for (const gmx_molblock_t &molb : mtop->molblock)
794 int n_excl_mol, n_excl_icg, n_excl_at_max;
796 const gmx_moltype_t &molt = mtop->moltype[molb.type];
797 count_excls(&molt.cgs, &molt.excls,
798 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
799 nexcl += molb.nmol*n_excl_mol;
800 dd->n_intercg_excl += molb.nmol*n_excl_icg;
801 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
803 if (rt->bExclRequired)
805 dd->nbonded_global += nexcl;
806 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
808 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
809 "will use an extra communication step for exclusion forces for %s\n",
810 dd->n_intercg_excl, eel_names[ir->coulombtype]);
814 if (vsite && vsite->n_intercg_vsite > 0)
818 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
819 "will an extra communication step for selected coordinates and forces\n",
820 vsite->n_intercg_vsite);
822 init_domdec_vsites(dd, vsite->n_intercg_vsite);
825 if (dd->splitConstraints || dd->splitSettles)
827 init_domdec_constraints(dd, mtop);
831 fprintf(fplog, "\n");
835 /*! \brief Store a vsite interaction at the end of \p il
837 * This routine is very similar to add_ifunc, but vsites interactions
838 * have more work to do than other kinds of interactions, and the
839 * complex way nral (and thus vector contents) depends on ftype
840 * confuses static analysis tools unless we fuse the vsite
841 * atom-indexing organization code with the ifunc-adding code, so that
842 * they can see that nral is the same value. */
844 add_ifunc_for_vsites(t_iatom *tiatoms, const gmx_ga2la_t &ga2la,
845 int nral, gmx_bool bHomeA,
846 int a, int a_gl, int a_mol,
847 const t_iatom *iatoms,
852 if (il->nr+1+nral > il->nalloc)
854 il->nalloc = over_alloc_large(il->nr+1+nral);
855 srenew(il->iatoms, il->nalloc);
857 liatoms = il->iatoms + il->nr;
861 tiatoms[0] = iatoms[0];
865 /* We know the local index of the first atom */
870 /* Convert later in make_local_vsites */
871 tiatoms[1] = -a_gl - 1;
874 for (int k = 2; k < 1+nral; k++)
876 int ak_gl = a_gl + iatoms[k] - a_mol;
877 if (const int *homeIndex = ga2la.findHome(ak_gl))
879 tiatoms[k] = *homeIndex;
883 /* Copy the global index, convert later in make_local_vsites */
884 tiatoms[k] = -(ak_gl + 1);
886 // Note that ga2la_get_home always sets the third parameter if
889 for (int k = 0; k < 1+nral; k++)
891 liatoms[k] = tiatoms[k];
895 /*! \brief Store a bonded interaction at the end of \p il */
896 static inline void add_ifunc(int nral, const t_iatom *tiatoms, t_ilist *il)
901 if (il->nr+1+nral > il->nalloc)
903 il->nalloc = over_alloc_large(il->nr+1+nral);
904 srenew(il->iatoms, il->nalloc);
906 liatoms = il->iatoms + il->nr;
907 for (k = 0; k <= nral; k++)
909 liatoms[k] = tiatoms[k];
914 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
915 static void add_posres(int mol, int a_mol, int numAtomsInMolecule,
916 const gmx_molblock_t *molb,
917 t_iatom *iatoms, const t_iparams *ip_in,
923 /* This position restraint has not been added yet,
924 * so it's index is the current number of position restraints.
926 n = idef->il[F_POSRES].nr/2;
927 if (n+1 > idef->iparams_posres_nalloc)
929 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
930 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
932 ip = &idef->iparams_posres[n];
933 /* Copy the force constants */
934 *ip = ip_in[iatoms[0]];
936 /* Get the position restraint coordinates from the molblock */
937 a_molb = mol*numAtomsInMolecule + a_mol;
938 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
939 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
940 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
941 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
942 if (!molb->posres_xB.empty())
944 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
945 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
946 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
950 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
951 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
952 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
954 /* Set the parameter index for idef->iparams_posre */
958 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
959 static void add_fbposres(int mol, int a_mol, int numAtomsInMolecule,
960 const gmx_molblock_t *molb,
961 t_iatom *iatoms, const t_iparams *ip_in,
967 /* This flat-bottom position restraint has not been added yet,
968 * so it's index is the current number of position restraints.
970 n = idef->il[F_FBPOSRES].nr/2;
971 if (n+1 > idef->iparams_fbposres_nalloc)
973 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
974 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
976 ip = &idef->iparams_fbposres[n];
977 /* Copy the force constants */
978 *ip = ip_in[iatoms[0]];
980 /* Get the position restraint coordinats from the molblock */
981 a_molb = mol*numAtomsInMolecule + a_mol;
982 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
983 /* Take reference positions from A position of normal posres */
984 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
985 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
986 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
988 /* Note: no B-type for flat-bottom posres */
990 /* Set the parameter index for idef->iparams_posre */
994 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
995 static void add_vsite(const gmx_ga2la_t &ga2la,
996 gmx::ArrayRef<const int> index,
997 gmx::ArrayRef<const int> rtil,
999 gmx_bool bHomeA, int a, int a_gl, int a_mol,
1000 const t_iatom *iatoms,
1005 t_iatom tiatoms[1+MAXATOMLIST];
1006 int j, ftype_r, nral_r;
1008 /* Add this interaction to the local topology */
1009 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1013 std::vector<int> &vsitePbcFtype = (*vsitePbc)[ftype - c_ftypeVsiteStart];
1014 const int vsi = idef->il[ftype].nr/(1+nral) - 1;
1015 if (static_cast<size_t>(vsi) >= vsitePbcFtype.size())
1017 vsitePbcFtype.resize(vsi + 1);
1021 pbc_a_mol = iatoms[1+nral+1];
1024 /* The pbc flag is one of the following two options:
1025 * -2: vsite and all constructing atoms are within the same cg, no pbc
1026 * -1: vsite and its first constructing atom are in the same cg, do pbc
1028 vsitePbcFtype[vsi] = pbc_a_mol;
1032 /* Set the pbc atom for this vsite so we can make its pbc
1033 * identical to the rest of the atoms in its charge group.
1034 * Since the order of the atoms does not change within a charge
1035 * group, we do not need the global to local atom index.
1037 vsitePbcFtype[vsi] = a + pbc_a_mol - iatoms[1];
1042 /* This vsite is non-home (required for recursion),
1043 * and therefore there is no charge group to match pbc with.
1044 * But we always turn on full_pbc to assure that higher order
1045 * recursion works correctly.
1047 vsitePbcFtype[vsi] = -1;
1053 /* Check for recursion */
1054 for (k = 2; k < 1+nral; k++)
1056 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1058 /* This construction atoms is a vsite and not a home atom */
1061 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1063 /* Find the vsite construction */
1065 /* Check all interactions assigned to this atom */
1066 j = index[iatoms[k]];
1067 while (j < index[iatoms[k]+1])
1069 ftype_r = rtil[j++];
1070 nral_r = NRAL(ftype_r);
1071 if (interaction_function[ftype_r].flags & IF_VSITE)
1073 /* Add this vsite (recursion) */
1074 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1075 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1078 j += 1 + nral_r + 2;
1090 /*! \brief Build the index that maps each local atom to its local atom group */
1091 static void makeLocalAtomGroupsFromAtoms(gmx_domdec_t *dd)
1093 const gmx::RangePartitioning &atomGrouping = dd->atomGrouping();
1095 dd->localAtomGroupFromAtom.clear();
1097 for (size_t g = 0; g < dd->globalAtomGroupIndices.size(); g++)
1099 for (int gmx_unused a : atomGrouping.block(g))
1101 dd->localAtomGroupFromAtom.push_back(g);
1106 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1107 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1113 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1117 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1123 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1124 static void combine_blocka(t_blocka *dest,
1125 gmx::ArrayRef<const thread_work_t> src)
1127 int ni = src.back().excl.nr;
1129 for (const thread_work_t &th_work : src)
1131 na += th_work.excl.nra;
1133 if (ni + 1 > dest->nalloc_index)
1135 dest->nalloc_index = over_alloc_large(ni+1);
1136 srenew(dest->index, dest->nalloc_index);
1138 if (dest->nra + na > dest->nalloc_a)
1140 dest->nalloc_a = over_alloc_large(dest->nra+na);
1141 srenew(dest->a, dest->nalloc_a);
1143 for (gmx::index s = 1; s < src.size(); s++)
1145 for (int i = dest->nr + 1; i < src[s].excl.nr + 1; i++)
1147 dest->index[i] = dest->nra + src[s].excl.index[i];
1149 for (int i = 0; i < src[s].excl.nra; i++)
1151 dest->a[dest->nra+i] = src[s].excl.a[i];
1153 dest->nr = src[s].excl.nr;
1154 dest->nra += src[s].excl.nra;
1158 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1159 * virtual sites need special attention, as pbc info differs per vsite.
1161 static void combine_idef(t_idef *dest,
1162 gmx::ArrayRef<const thread_work_t> src,
1167 for (ftype = 0; ftype < F_NRE; ftype++)
1170 for (gmx::index s = 1; s < src.size(); s++)
1172 n += src[s].idef.il[ftype].nr;
1178 ild = &dest->il[ftype];
1180 if (ild->nr + n > ild->nalloc)
1182 ild->nalloc = over_alloc_large(ild->nr+n);
1183 srenew(ild->iatoms, ild->nalloc);
1187 (((interaction_function[ftype].flags & IF_VSITE) != 0u) &&
1188 vsite->vsite_pbc_loc);
1189 const int nral1 = 1 + NRAL(ftype);
1190 const int ftv = ftype - c_ftypeVsiteStart;
1192 for (gmx::index s = 1; s < src.size(); s++)
1194 const t_ilist &ils = src[s].idef.il[ftype];
1196 for (int i = 0; i < ils.nr; i++)
1198 ild->iatoms[ild->nr + i] = ils.iatoms[i];
1202 const std::vector<int> &pbcSrc = (*src[s].vsitePbc)[ftv];
1203 std::vector<int> &pbcDest = (*vsite->vsite_pbc_loc)[ftv];
1204 pbcDest.resize((ild->nr + ils.nr)/nral1);
1205 for (int i = 0; i < ils.nr; i += nral1)
1207 pbcDest[(ild->nr + i)/nral1] = pbcSrc[i/nral1];
1214 /* Position restraints need an additional treatment */
1215 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1217 int nposres = dest->il[ftype].nr/2;
1218 // TODO: Simplify this code using std::vector
1219 t_iparams * &iparams_dest = (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1220 int &posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc : dest->iparams_fbposres_nalloc);
1221 if (nposres > posres_nalloc)
1223 posres_nalloc = over_alloc_large(nposres);
1224 srenew(iparams_dest, posres_nalloc);
1227 /* Set nposres to the number of original position restraints in dest */
1228 for (gmx::index s = 1; s < src.size(); s++)
1230 nposres -= src[s].idef.il[ftype].nr/2;
1233 for (gmx::index s = 1; s < src.size(); s++)
1235 const t_iparams *iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1237 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1239 /* Correct the index into iparams_posres */
1240 dest->il[ftype].iatoms[nposres*2] = nposres;
1241 /* Copy the position restraint force parameters */
1242 iparams_dest[nposres] = iparams_src[i];
1251 /*! \brief Check and when available assign bonded interactions for local atom i
1254 check_assign_interactions_atom(int i, int i_gl,
1256 int numAtomsInMolecule,
1257 gmx::ArrayRef<const int> index,
1258 gmx::ArrayRef<const int> rtil,
1259 gmx_bool bInterMolInteractions,
1260 int ind_start, int ind_end,
1261 const gmx_domdec_t *dd,
1262 const gmx_domdec_zones_t *zones,
1263 const gmx_molblock_t *molb,
1264 gmx_bool bRCheckMB, const ivec rcheck, gmx_bool bRCheck2B,
1269 const t_iparams *ip_in,
1281 t_iatom tiatoms[1 + MAXATOMLIST];
1283 const int ftype = rtil[j++];
1284 auto iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1285 const int nral = NRAL(ftype);
1286 if (interaction_function[ftype].flags & IF_VSITE)
1288 assert(!bInterMolInteractions);
1289 /* The vsite construction goes where the vsite itself is */
1292 add_vsite(*dd->ga2la, index, rtil, ftype, nral,
1293 TRUE, i, i_gl, i_mol,
1294 iatoms.data(), idef, vsitePbc);
1303 tiatoms[0] = iatoms[0];
1307 assert(!bInterMolInteractions);
1308 /* Assign single-body interactions to the home zone */
1313 if (ftype == F_POSRES)
1315 add_posres(mol, i_mol, numAtomsInMolecule,
1316 molb, tiatoms, ip_in, idef);
1318 else if (ftype == F_FBPOSRES)
1320 add_fbposres(mol, i_mol, numAtomsInMolecule,
1321 molb, tiatoms, ip_in, idef);
1331 /* This is a two-body interaction, we can assign
1332 * analogous to the non-bonded assignments.
1336 if (!bInterMolInteractions)
1338 /* Get the global index using the offset in the molecule */
1339 k_gl = i_gl + iatoms[2] - i_mol;
1345 if (const auto *entry = dd->ga2la->find(k_gl))
1347 int kz = entry->cell;
1352 /* Check zone interaction assignments */
1353 bUse = ((iz < zones->nizone &&
1355 kz >= zones->izone[iz].j0 &&
1356 kz < zones->izone[iz].j1) ||
1357 (kz < zones->nizone &&
1359 iz >= zones->izone[kz].j0 &&
1360 iz < zones->izone[kz].j1));
1364 tiatoms[2] = entry->la;
1365 /* If necessary check the cgcm distance */
1367 dd_dist2(pbc_null, cg_cm, la2lc,
1368 tiatoms[1], tiatoms[2]) >= rc2)
1381 /* Assign this multi-body bonded interaction to
1382 * the local node if we have all the atoms involved
1383 * (local or communicated) and the minimum zone shift
1384 * in each dimension is zero, for dimensions
1385 * with 2 DD cells an extra check may be necessary.
1387 ivec k_zero, k_plus;
1393 for (k = 1; k <= nral && bUse; k++)
1396 if (!bInterMolInteractions)
1398 /* Get the global index using the offset in the molecule */
1399 k_gl = i_gl + iatoms[k] - i_mol;
1405 const auto *entry = dd->ga2la->find(k_gl);
1406 if (entry == nullptr || entry->cell >= zones->n)
1408 /* We do not have this atom of this interaction
1409 * locally, or it comes from more than one cell
1418 tiatoms[k] = entry->la;
1419 for (d = 0; d < DIM; d++)
1421 if (zones->shift[entry->cell][d] == 0)
1433 (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1438 for (d = 0; (d < DIM && bUse); d++)
1440 /* Check if the cg_cm distance falls within
1441 * the cut-off to avoid possible multiple
1442 * assignments of bonded interactions.
1446 dd_dist2(pbc_null, cg_cm, la2lc,
1447 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1456 /* Add this interaction to the local topology */
1457 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1458 /* Sum so we can check in global_stat
1459 * if we have everything.
1462 !(interaction_function[ftype].flags & IF_LIMZERO))
1472 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1474 * With thread parallelizing each thread acts on a different atom range:
1475 * at_start to at_end.
1477 static int make_bondeds_zone(gmx_domdec_t *dd,
1478 const gmx_domdec_zones_t *zones,
1479 const std::vector<gmx_molblock_t> &molb,
1480 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1482 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1483 const t_iparams *ip_in,
1487 gmx::RangePartitioning::Block atomRange)
1489 int mb, mt, mol, i_mol;
1491 gmx_reverse_top_t *rt;
1494 rt = dd->reverse_top;
1496 bBCheck = rt->bBCheck;
1500 for (int i : atomRange)
1502 /* Get the global atom number */
1503 const int i_gl = dd->globalAtomIndices[i];
1504 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1505 /* Check all intramolecular interactions assigned to this atom */
1506 gmx::ArrayRef<const int> index = rt->ril_mt[mt].index;
1507 gmx::ArrayRef<const t_iatom> rtil = rt->ril_mt[mt].il;
1509 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1510 rt->ril_mt[mt].numAtomsInMolecule,
1512 index[i_mol], index[i_mol+1],
1515 bRCheckMB, rcheck, bRCheck2B, rc2,
1526 if (rt->bIntermolecularInteractions)
1528 /* Check all intermolecular interactions assigned to this atom */
1529 index = rt->ril_intermol.index;
1530 rtil = rt->ril_intermol.il;
1532 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1533 rt->ril_mt[mt].numAtomsInMolecule,
1535 index[i_gl], index[i_gl + 1],
1538 bRCheckMB, rcheck, bRCheck2B, rc2,
1550 return nbonded_local;
1553 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1554 static void set_no_exclusions_zone(const gmx_domdec_t *dd,
1555 const gmx_domdec_zones_t *zones,
1559 const auto zone = dd->atomGrouping().subRange(zones->cg_range[iz],
1560 zones->cg_range[iz + 1]);
1564 lexcls->index[a + 1] = lexcls->nra;
1568 /*! \brief Set the exclusion data for i-zone \p iz
1570 * This is a legacy version for the group scheme of the same routine below.
1571 * Here charge groups and distance checks to ensure unique exclusions
1574 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1575 const std::vector<gmx_moltype_t> &moltype,
1576 gmx_bool bRCheck, real rc2,
1577 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1581 int cg_start, int cg_end)
1585 const t_blocka *excls;
1587 const gmx_ga2la_t &ga2la = *dd->ga2la;
1590 dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1591 zones->izone[iz].jcg1);
1593 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1595 /* We set the end index, but note that we might not start at zero here */
1596 lexcls->nr = dd->atomGrouping().subRange(0, cg_end).size();
1598 int n = lexcls->nra;
1600 for (int cg = cg_start; cg < cg_end; cg++)
1602 if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1604 lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1605 srenew(lexcls->a, lexcls->nalloc_a);
1607 const auto atomGroup = dd->atomGrouping().block(cg);
1608 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1609 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1611 /* Copy the exclusions from the global top */
1612 for (int la : atomGroup)
1614 lexcls->index[la] = n;
1615 int a_gl = dd->globalAtomIndices[la];
1617 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1618 excls = &moltype[mt].excls;
1619 for (int j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1621 int aj_mol = excls->a[j];
1622 /* This computation of jla is only correct intra-cg */
1623 int jla = la + aj_mol - a_mol;
1624 if (atomGroup.inRange(jla))
1626 /* This is an intra-cg exclusion. We can skip
1627 * the global indexing and distance checking.
1629 /* Intra-cg exclusions are only required
1630 * for the home zone.
1634 lexcls->a[n++] = jla;
1635 /* Check to avoid double counts */
1644 /* This is a inter-cg exclusion */
1645 /* Since exclusions are pair interactions,
1646 * just like non-bonded interactions,
1647 * they can be assigned properly up
1648 * to the DD cutoff (not cutoff_min as
1649 * for the other bonded interactions).
1651 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1653 if (iz == 0 && jEntry->cell == 0)
1655 lexcls->a[n++] = jEntry->la;
1656 /* Check to avoid double counts */
1657 if (jEntry->la > la)
1662 else if (jRange.inRange(jEntry->la) &&
1664 dd_dist2(pbc_null, cg_cm, la2lc, la, jEntry->la) < rc2))
1666 /* jla > la, since jRange.begin() > la */
1667 lexcls->a[n++] = jEntry->la;
1677 /* There are no inter-cg excls and this cg is self-excluded.
1678 * These exclusions are only required for zone 0,
1679 * since other zones do not see themselves.
1683 for (int la : atomGroup)
1685 lexcls->index[la] = n;
1686 for (int j : atomGroup)
1691 count += (atomGroup.size()*(atomGroup.size() - 1))/2;
1695 /* We don't need exclusions for this cg */
1696 for (int la : atomGroup)
1698 lexcls->index[la] = n;
1704 lexcls->index[lexcls->nr] = n;
1710 /*! \brief Set the exclusion data for i-zone \p iz */
1711 static void make_exclusions_zone(gmx_domdec_t *dd,
1712 gmx_domdec_zones_t *zones,
1713 const std::vector<gmx_moltype_t> &moltype,
1717 int at_start, int at_end)
1719 int n_excl_at_max, n, at;
1721 const gmx_ga2la_t &ga2la = *dd->ga2la;
1724 dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1725 zones->izone[iz].jcg1);
1727 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1729 /* We set the end index, but note that we might not start at zero here */
1730 lexcls->nr = at_end;
1733 for (at = at_start; at < at_end; at++)
1735 if (n + 1000 > lexcls->nalloc_a)
1737 lexcls->nalloc_a = over_alloc_large(n + 1000);
1738 srenew(lexcls->a, lexcls->nalloc_a);
1740 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1742 int a_gl, mb, mt, mol, a_mol, j;
1743 const t_blocka *excls;
1745 if (n + n_excl_at_max > lexcls->nalloc_a)
1747 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1748 srenew(lexcls->a, lexcls->nalloc_a);
1751 /* Copy the exclusions from the global top */
1752 lexcls->index[at] = n;
1753 a_gl = dd->globalAtomIndices[at];
1754 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1755 &mb, &mt, &mol, &a_mol);
1756 excls = &moltype[mt].excls;
1757 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1759 const int aj_mol = excls->a[j];
1761 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1763 /* This check is not necessary, but it can reduce
1764 * the number of exclusions in the list, which in turn
1765 * can speed up the pair list construction a bit.
1767 if (jRange.inRange(jEntry->la))
1769 lexcls->a[n++] = jEntry->la;
1776 /* We don't need exclusions for this atom */
1777 lexcls->index[at] = n;
1781 lexcls->index[lexcls->nr] = n;
1786 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1787 static void check_alloc_index(t_blocka *ba, int nindex_max)
1789 if (nindex_max+1 > ba->nalloc_index)
1791 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1792 srenew(ba->index, ba->nalloc_index);
1796 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1797 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1800 int nr = dd->atomGrouping().subRange(0, zones->izone[zones->nizone - 1].cg1).size();
1802 check_alloc_index(lexcls, nr);
1804 for (size_t thread = 1; thread < dd->reverse_top->th_work.size(); thread++)
1806 check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1810 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1811 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1814 const auto nonhomeIzonesAtomRange =
1815 dd->atomGrouping().subRange(zones->izone[0].cg1,
1816 zones->izone[zones->nizone - 1].cg1);
1818 if (dd->n_intercg_excl == 0)
1820 /* There are no exclusions involving non-home charge groups,
1821 * but we need to set the indices for neighborsearching.
1823 for (int la : nonhomeIzonesAtomRange)
1825 lexcls->index[la] = lexcls->nra;
1828 /* nr is only used to loop over the exclusions for Ewald and RF,
1829 * so we can set it to the number of home atoms for efficiency.
1831 lexcls->nr = nonhomeIzonesAtomRange.begin();
1835 lexcls->nr = nonhomeIzonesAtomRange.end();
1839 /*! \brief Clear a t_idef struct */
1840 static void clear_idef(t_idef *idef)
1844 /* Clear the counts */
1845 for (ftype = 0; ftype < F_NRE; ftype++)
1847 idef->il[ftype].nr = 0;
1851 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1852 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1853 gmx_domdec_zones_t *zones,
1854 const gmx_mtop_t *mtop,
1856 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1858 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1859 t_idef *idef, gmx_vsite_t *vsite,
1860 t_blocka *lexcls, int *excl_count)
1862 int nzone_bondeds, nzone_excl;
1863 int izone, cg0, cg1;
1867 gmx_reverse_top_t *rt;
1869 if (dd->reverse_top->bInterCGInteractions)
1871 nzone_bondeds = zones->n;
1875 /* Only single charge group (or atom) molecules, so interactions don't
1876 * cross zone boundaries and we only need to assign in the home zone.
1881 if (dd->n_intercg_excl > 0)
1883 /* We only use exclusions from i-zones to i- and j-zones */
1884 nzone_excl = zones->nizone;
1888 /* There are no inter-cg exclusions and only zone 0 sees itself */
1892 check_exclusions_alloc(dd, zones, lexcls);
1894 rt = dd->reverse_top;
1898 /* Clear the counts */
1906 for (izone = 0; izone < nzone_bondeds; izone++)
1908 cg0 = zones->cg_range[izone];
1909 cg1 = zones->cg_range[izone + 1];
1911 const int numThreads = rt->th_work.size();
1912 #pragma omp parallel for num_threads(numThreads) schedule(static)
1913 for (thread = 0; thread < numThreads; thread++)
1921 cg0t = cg0 + ((cg1 - cg0)* thread )/numThreads;
1922 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/numThreads;
1930 idef_t = &rt->th_work[thread].idef;
1934 VsitePbc *vsitePbc = nullptr;
1935 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1939 vsitePbc = vsite->vsite_pbc_loc.get();
1943 vsitePbc = rt->th_work[thread].vsitePbc.get();
1947 rt->th_work[thread].nbonded =
1948 make_bondeds_zone(dd, zones,
1950 bRCheckMB, rcheck, bRCheck2B, rc2,
1951 la2lc, pbc_null, cg_cm, idef->iparams,
1954 dd->atomGrouping().subRange(cg0t, cg1t));
1956 if (izone < nzone_excl)
1964 excl_t = &rt->th_work[thread].excl;
1969 if (dd->atomGrouping().allBlocksHaveSizeOne() &&
1972 /* No charge groups and no distance check required */
1973 make_exclusions_zone(dd, zones,
1974 mtop->moltype, cginfo,
1981 rt->th_work[thread].excl_count =
1982 make_exclusions_zone_cg(dd, zones,
1983 mtop->moltype, bRCheck2B, rc2,
1984 la2lc, pbc_null, cg_cm, cginfo,
1991 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1994 if (rt->th_work.size() > 1)
1996 combine_idef(idef, rt->th_work, vsite);
1999 for (const thread_work_t &th_work : rt->th_work)
2001 nbonded_local += th_work.nbonded;
2004 if (izone < nzone_excl)
2006 if (rt->th_work.size() > 1)
2008 combine_blocka(lexcls, rt->th_work);
2011 for (const thread_work_t &th_work : rt->th_work)
2013 *excl_count += th_work.excl_count;
2018 /* Some zones might not have exclusions, but some code still needs to
2019 * loop over the index, so we set the indices here.
2021 for (izone = nzone_excl; izone < zones->nizone; izone++)
2023 set_no_exclusions_zone(dd, zones, izone, lexcls);
2026 finish_local_exclusions(dd, zones, lexcls);
2029 fprintf(debug, "We have %d exclusions, check count %d\n",
2030 lexcls->nra, *excl_count);
2033 return nbonded_local;
2036 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
2038 lcgs->nr = dd->globalAtomGroupIndices.size();
2039 lcgs->index = dd->atomGrouping_.rawIndex().data();
2042 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
2043 int npbcdim, matrix box,
2044 rvec cellsize_min, const ivec npulse,
2048 const gmx_mtop_t *mtop, gmx_localtop_t *ltop)
2050 gmx_bool bRCheckMB, bRCheck2B;
2054 t_pbc pbc, *pbc_null = nullptr;
2058 fprintf(debug, "Making local topology\n");
2061 dd_make_local_cgs(dd, <op->cgs);
2066 if (dd->reverse_top->bInterCGInteractions)
2068 /* We need to check to which cell bondeds should be assigned */
2069 rc = dd_cutoff_twobody(dd);
2072 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2075 /* Should we check cg_cm distances when assigning bonded interactions? */
2076 for (d = 0; d < DIM; d++)
2079 /* Only need to check for dimensions where the part of the box
2080 * that is not communicated is smaller than the cut-off.
2082 if (d < npbcdim && dd->nc[d] > 1 &&
2083 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2090 /* Check for interactions between two atoms,
2091 * where we can allow interactions up to the cut-off,
2092 * instead of up to the smallest cell dimension.
2099 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n",
2100 d, cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
2103 if (bRCheckMB || bRCheck2B)
2105 makeLocalAtomGroupsFromAtoms(dd);
2108 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
2118 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
2119 bRCheckMB, rcheck, bRCheck2B, rc,
2120 dd->localAtomGroupFromAtom.data(),
2121 pbc_null, cgcm_or_x,
2123 <op->excls, &nexcl);
2125 /* The ilist is not sorted yet,
2126 * we can only do this when we have the charge arrays.
2128 ltop->idef.ilsort = ilsortUNKNOWN;
2130 if (dd->reverse_top->bExclRequired)
2132 dd->nbonded_local += nexcl;
2135 ltop->atomtypes = mtop->atomtypes;
2138 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
2139 gmx_localtop_t *ltop)
2141 if (dd->reverse_top->ilsort == ilsortNO_FE)
2143 ltop->idef.ilsort = ilsortNO_FE;
2147 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
2151 gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global)
2153 gmx_localtop_t *top;
2157 /* TODO: Get rid of the const casts below, e.g. by using a reference */
2158 top->idef.ntypes = top_global->ffparams.numTypes();
2159 top->idef.atnr = top_global->ffparams.atnr;
2160 top->idef.functype = const_cast<t_functype *>(top_global->ffparams.functype.data());
2161 top->idef.iparams = const_cast<t_iparams *>(top_global->ffparams.iparams.data());
2162 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
2163 top->idef.cmap_grid = new gmx_cmap_t;
2164 *top->idef.cmap_grid = top_global->ffparams.cmap_grid;
2166 top->idef.ilsort = ilsortUNKNOWN;
2171 void dd_init_local_state(gmx_domdec_t *dd,
2172 const t_state *state_global, t_state *state_local)
2174 int buf[NITEM_DD_INIT_LOCAL_STATE];
2178 buf[0] = state_global->flags;
2179 buf[1] = state_global->ngtc;
2180 buf[2] = state_global->nnhpres;
2181 buf[3] = state_global->nhchainlength;
2182 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
2184 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2186 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
2187 init_dfhist_state(state_local, buf[4]);
2188 state_local->flags = buf[0];
2191 /*! \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 */
2192 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2198 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2200 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2201 if (link->a[k] == cg_gl_j)
2208 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2209 "Inconsistent allocation of link");
2210 /* Add this charge group link */
2211 if (link->index[cg_gl+1]+1 > link->nalloc_a)
2213 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2214 srenew(link->a, link->nalloc_a);
2216 link->a[link->index[cg_gl+1]] = cg_gl_j;
2217 link->index[cg_gl+1]++;
2221 /*! \brief Return a vector of the charge group index for all atoms */
2222 static std::vector<int> make_at2cg(const t_block &cgs)
2224 std::vector<int> at2cg(cgs.index[cgs.nr]);
2225 for (int cg = 0; cg < cgs.nr; cg++)
2227 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2236 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
2237 cginfo_mb_t *cginfo_mb)
2239 gmx_bool bExclRequired;
2241 cginfo_mb_t *cgi_mb;
2243 /* For each charge group make a list of other charge groups
2244 * in the system that a linked to it via bonded interactions
2245 * which are also stored in reverse_top.
2248 bExclRequired = dd->reverse_top->bExclRequired;
2250 reverse_ilist_t ril_intermol;
2251 if (mtop->bIntermolecularInteractions)
2253 if (ncg_mtop(mtop) < mtop->natoms)
2255 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.");
2260 atoms.nr = mtop->natoms;
2261 atoms.atom = nullptr;
2263 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
2265 make_reverse_ilist(*mtop->intermolecular_ilist,
2267 gmx::EmptyArrayRef(),
2268 FALSE, FALSE, FALSE, TRUE, &ril_intermol);
2272 snew(link->index, ncg_mtop(mtop)+1);
2279 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
2281 const gmx_molblock_t &molb = mtop->molblock[mb];
2286 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2287 const t_block &cgs = molt.cgs;
2288 const t_blocka &excls = molt.excls;
2289 std::vector<int> a2c = make_at2cg(cgs);
2290 /* Make a reverse ilist in which the interactions are linked
2291 * to all atoms, not only the first atom as in gmx_reverse_top.
2292 * The constraints are discarded here.
2294 reverse_ilist_t ril;
2295 make_reverse_ilist(molt.ilist, &molt.atoms, gmx::EmptyArrayRef(),
2296 FALSE, FALSE, FALSE, TRUE, &ril);
2298 cgi_mb = &cginfo_mb[mb];
2301 for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb.nmol : 1); mol++)
2303 for (int cg = 0; cg < cgs.nr; cg++)
2305 int cg_gl = cg_offset + cg;
2306 link->index[cg_gl+1] = link->index[cg_gl];
2307 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2309 int i = ril.index[a];
2310 while (i < ril.index[a+1])
2312 int ftype = ril.il[i++];
2313 int nral = NRAL(ftype);
2314 /* Skip the ifunc index */
2316 for (int j = 0; j < nral; j++)
2318 int aj = ril.il[i + j];
2321 check_link(link, cg_gl, cg_offset+a2c[aj]);
2324 i += nral_rt(ftype);
2328 /* Exclusions always go both ways */
2329 for (int j = excls.index[a]; j < excls.index[a + 1]; j++)
2331 int aj = excls.a[j];
2334 check_link(link, cg_gl, cg_offset+a2c[aj]);
2339 if (mtop->bIntermolecularInteractions)
2341 int i = ril_intermol.index[a];
2342 while (i < ril_intermol.index[a+1])
2344 int ftype = ril_intermol.il[i++];
2345 int nral = NRAL(ftype);
2346 /* Skip the ifunc index */
2348 for (int j = 0; j < nral; j++)
2350 /* Here we assume we have no charge groups;
2351 * this has been checked above.
2353 int aj = ril_intermol.il[i + j];
2354 check_link(link, cg_gl, aj);
2356 i += nral_rt(ftype);
2360 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2362 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2367 cg_offset += cgs.nr;
2369 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - cgs.nr];
2373 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt.name, cgs.nr, nlink_mol);
2376 if (molb.nmol > mol)
2378 /* Copy the data for the rest of the molecules in this block */
2379 link->nalloc_a += (molb.nmol - mol)*nlink_mol;
2380 srenew(link->a, link->nalloc_a);
2381 for (; mol < molb.nmol; mol++)
2383 for (int cg = 0; cg < cgs.nr; cg++)
2385 int cg_gl = cg_offset + cg;
2386 link->index[cg_gl + 1] =
2387 link->index[cg_gl + 1 - cgs.nr] + nlink_mol;
2388 for (int j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2390 link->a[j] = link->a[j - nlink_mol] + cgs.nr;
2392 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2393 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2395 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2399 cg_offset += cgs.nr;
2406 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2417 } bonded_distance_t;
2419 /*! \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 */
2420 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2421 bonded_distance_t *bd)
2432 /*! \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 */
2433 static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
2434 const std::vector<int> &at2cg,
2435 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2436 bonded_distance_t *bd_2b,
2437 bonded_distance_t *bd_mb)
2439 for (int ftype = 0; ftype < F_NRE; ftype++)
2441 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2443 const auto &il = molt->ilist[ftype];
2444 int nral = NRAL(ftype);
2447 for (int i = 0; i < il.size(); i += 1+nral)
2449 for (int ai = 0; ai < nral; ai++)
2451 int cgi = at2cg[il.iatoms[i+1+ai]];
2452 for (int aj = ai + 1; aj < nral; aj++)
2454 int cgj = at2cg[il.iatoms[i+1+aj]];
2457 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2459 update_max_bonded_distance(rij2, ftype,
2462 (nral == 2) ? bd_2b : bd_mb);
2472 const t_blocka *excls = &molt->excls;
2473 for (int ai = 0; ai < excls->nr; ai++)
2475 int cgi = at2cg[ai];
2476 for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2478 int cgj = at2cg[excls->a[j]];
2481 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2483 /* There is no function type for exclusions, use -1 */
2484 update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2491 /*! \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 */
2492 static void bonded_distance_intermol(const InteractionLists &ilists_intermol,
2494 const rvec *x, int ePBC, const matrix box,
2495 bonded_distance_t *bd_2b,
2496 bonded_distance_t *bd_mb)
2500 set_pbc(&pbc, ePBC, box);
2502 for (int ftype = 0; ftype < F_NRE; ftype++)
2504 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2506 const auto &il = ilists_intermol[ftype];
2507 int nral = NRAL(ftype);
2509 /* No nral>1 check here, since intermol interactions always
2510 * have nral>=2 (and the code is also correct for nral=1).
2512 for (int i = 0; i < il.size(); i += 1+nral)
2514 for (int ai = 0; ai < nral; ai++)
2516 int atom_i = il.iatoms[i + 1 + ai];
2518 for (int aj = ai + 1; aj < nral; aj++)
2523 int atom_j = il.iatoms[i + 1 + aj];
2525 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2529 update_max_bonded_distance(rij2, ftype,
2531 (nral == 2) ? bd_2b : bd_mb);
2539 //! Returns whether \p molt has at least one virtual site
2540 static bool moltypeHasVsite(const gmx_moltype_t &molt)
2542 bool hasVsite = false;
2543 for (int i = 0; i < F_NRE; i++)
2545 if ((interaction_function[i].flags & IF_VSITE) &&
2546 molt.ilist[i].size() > 0)
2555 //! Compute charge group centers of mass for molecule \p molt
2556 static void get_cgcm_mol(const gmx_moltype_t *molt,
2557 const gmx_ffparams_t *ffparams,
2558 int ePBC, t_graph *graph, const matrix box,
2559 const rvec *x, rvec *xs, rvec *cg_cm)
2563 if (ePBC != epbcNONE)
2565 mk_mshift(nullptr, graph, ePBC, box, x);
2567 shift_x(graph, box, x, xs);
2568 /* By doing an extra mk_mshift the molecules that are broken
2569 * because they were e.g. imported from another software
2570 * will be made whole again. Such are the healing powers
2573 mk_mshift(nullptr, graph, ePBC, box, xs);
2577 /* We copy the coordinates so the original coordinates remain
2578 * unchanged, just to be 100% sure that we do not affect
2579 * binary reproducibility of simulations.
2581 n = molt->cgs.index[molt->cgs.nr];
2582 for (i = 0; i < n; i++)
2584 copy_rvec(x[i], xs[i]);
2588 if (moltypeHasVsite(*molt))
2590 /* Convert to old, deprecated format */
2591 t_ilist ilist[F_NRE];
2592 for (int ftype = 0; ftype < F_NRE; ftype++)
2594 if (interaction_function[ftype].flags & IF_VSITE)
2596 ilist[ftype].nr = molt->ilist[ftype].size();
2597 ilist[ftype].iatoms = const_cast<int *>(molt->ilist[ftype].iatoms.data());
2601 construct_vsites(nullptr, xs, 0.0, nullptr,
2602 ffparams->iparams.data(), ilist,
2603 epbcNONE, TRUE, nullptr, nullptr);
2606 calc_cgcm(nullptr, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2609 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
2610 const gmx_mtop_t *mtop,
2611 const t_inputrec *ir,
2612 const rvec *x, const matrix box,
2614 real *r_2b, real *r_mb)
2616 gmx_bool bExclRequired;
2620 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2621 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2623 bExclRequired = inputrecExclForces(ir);
2628 for (const gmx_molblock_t &molb : mtop->molblock)
2630 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2631 if (molt.cgs.nr == 1 || molb.nmol == 0)
2633 at_offset += molb.nmol*molt.atoms.nr;
2637 if (ir->ePBC != epbcNONE)
2639 mk_graph_moltype(molt, &graph);
2642 std::vector<int> at2cg = make_at2cg(molt.cgs);
2643 snew(xs, molt.atoms.nr);
2644 snew(cg_cm, molt.cgs.nr);
2645 for (int mol = 0; mol < molb.nmol; mol++)
2647 get_cgcm_mol(&molt, &mtop->ffparams, ir->ePBC, &graph, box,
2648 x+at_offset, xs, cg_cm);
2650 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2651 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2653 bonded_cg_distance_mol(&molt, at2cg, bBCheck, bExclRequired, cg_cm,
2654 &bd_mol_2b, &bd_mol_mb);
2656 /* Process the mol data adding the atom index offset */
2657 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2658 at_offset + bd_mol_2b.a1,
2659 at_offset + bd_mol_2b.a2,
2661 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2662 at_offset + bd_mol_mb.a1,
2663 at_offset + bd_mol_mb.a2,
2666 at_offset += molt.atoms.nr;
2670 if (ir->ePBC != epbcNONE)
2677 if (mtop->bIntermolecularInteractions)
2679 if (ncg_mtop(mtop) < mtop->natoms)
2681 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.");
2684 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
2686 bonded_distance_intermol(*mtop->intermolecular_ilist,
2692 *r_2b = sqrt(bd_2b.r2);
2693 *r_mb = sqrt(bd_mb.r2);
2695 if (*r_2b > 0 || *r_mb > 0)
2697 GMX_LOG(mdlog.info).appendText("Initial maximum inter charge-group distances:");
2700 GMX_LOG(mdlog.info).appendTextFormatted(
2701 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2702 *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2703 bd_2b.a1 + 1, bd_2b.a2 + 1);
2707 GMX_LOG(mdlog.info).appendTextFormatted(
2708 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2709 *r_mb, interaction_function[bd_mb.ftype].longname,
2710 bd_mb.a1 + 1, bd_mb.a2 + 1);