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"
90 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
91 #define NITEM_DD_INIT_LOCAL_STATE 5
93 struct reverse_ilist_t
95 std::vector<int> index; /* Index for each atom into il */
96 std::vector<int> il; /* ftype|type|a0|...|an|ftype|... */
97 int numAtomsInMolecule; /* The number of atoms in this molecule */
100 struct MolblockIndices
108 /*! \brief Struct for thread local work data for local topology generation */
111 t_idef idef; /**< Partial local topology */
112 std::unique_ptr<VsitePbc> vsitePbc; /**< vsite PBC structure */
113 int nbonded; /**< The number of bondeds in this struct */
114 t_blocka excl; /**< List of exclusions */
115 int excl_count; /**< The total exclusion count for \p excl */
118 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
119 struct gmx_reverse_top_t
121 //! @cond Doxygen_Suppress
122 //! \brief Do we require all exclusions to be assigned?
123 bool bExclRequired = false;
124 //! \brief The maximum number of exclusions one atom can have
125 int n_excl_at_max = 0;
126 //! \brief Are there constraints in this revserse top?
127 bool bConstr = false;
128 //! \brief Are there settles in this revserse top?
129 bool bSettle = false;
130 //! \brief All bonded interactions have to be assigned?
131 bool bBCheck = false;
132 //! \brief Are there bondeds/exclusions between charge-groups?
133 bool bInterCGInteractions = false;
134 //! \brief Reverse ilist for all moltypes
135 std::vector<reverse_ilist_t> ril_mt;
136 //! \brief The size of ril_mt[?].index summed over all entries
137 int ril_mt_tot_size = 0;
138 //! \brief The sorting state of bondeds for free energy
139 int ilsort = ilsortUNKNOWN;
140 //! \brief molblock to global atom index for quick lookup of molblocks on atom index
141 std::vector<MolblockIndices> mbi;
143 //! \brief Do we have intermolecular interactions?
144 bool bIntermolecularInteractions = false;
145 //! \brief Intermolecular reverse ilist
146 reverse_ilist_t ril_intermol;
148 /* Work data structures for multi-threading */
149 //! \brief Thread work array for local topology generation
150 std::vector<thread_work_t> th_work;
154 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
155 static int nral_rt(int ftype)
160 if (interaction_function[ftype].flags & IF_VSITE)
162 /* With vsites the reverse topology contains
163 * two extra entries for PBC.
171 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
172 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
173 gmx_bool bConstr, gmx_bool bSettle)
175 return ((((interaction_function[ftype].flags & IF_BOND) != 0u) &&
176 ((interaction_function[ftype].flags & IF_VSITE) == 0u) &&
177 (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0u))) ||
178 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
179 (bSettle && ftype == F_SETTLE));
182 /*! \brief Help print error output when interactions are missing */
184 print_missing_interactions_mb(t_commrec *cr,
185 const gmx_reverse_top_t *rt,
186 const char *moltypename,
187 const reverse_ilist_t *ril,
188 int a_start, int a_end,
189 int nat_mol, int nmol,
193 int nril_mol = ril->index[nat_mol];
194 snew(assigned, nmol*nril_mol);
195 gmx::StringOutputStream stream;
196 gmx::TextWriter log(&stream);
198 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
199 for (int ftype = 0; ftype < F_NRE; ftype++)
201 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
203 int nral = NRAL(ftype);
204 const t_ilist *il = &idef->il[ftype];
205 const t_iatom *ia = il->iatoms;
206 for (int i = 0; i < il->nr; i += 1+nral)
208 int a0 = gatindex[ia[1]];
209 /* Check if this interaction is in
210 * the currently checked molblock.
212 if (a0 >= a_start && a0 < a_end)
214 int mol = (a0 - a_start)/nat_mol;
215 int a0_mol = (a0 - a_start) - mol*nat_mol;
216 int j_mol = ril->index[a0_mol];
218 while (j_mol < ril->index[a0_mol+1] && !found)
220 int j = mol*nril_mol + j_mol;
221 int ftype_j = ril->il[j_mol];
222 /* Here we need to check if this interaction has
223 * not already been assigned, since we could have
224 * multiply defined interactions.
226 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
229 /* Check the atoms */
231 for (int a = 0; a < nral; a++)
233 if (gatindex[ia[1+a]] !=
234 a_start + mol*nat_mol + ril->il[j_mol+2+a])
244 j_mol += 2 + nral_rt(ftype_j);
248 gmx_incons("Some interactions seem to be assigned multiple times");
256 gmx_sumi(nmol*nril_mol, assigned, cr);
260 for (int mol = 0; mol < nmol; mol++)
263 while (j_mol < nril_mol)
265 int ftype = ril->il[j_mol];
266 int nral = NRAL(ftype);
267 int j = mol*nril_mol + j_mol;
268 if (assigned[j] == 0 &&
269 !(interaction_function[ftype].flags & IF_VSITE))
271 if (DDMASTER(cr->dd))
275 log.writeLineFormatted("Molecule type '%s'", moltypename);
276 log.writeLineFormatted(
277 "the first %d missing interactions, except for exclusions:", nprint);
279 log.writeStringFormatted("%20s atoms",
280 interaction_function[ftype].longname);
282 for (a = 0; a < nral; a++)
284 log.writeStringFormatted("%5d", ril->il[j_mol+2+a]+1);
288 log.writeString(" ");
291 log.writeString(" global");
292 for (a = 0; a < nral; a++)
294 log.writeStringFormatted("%6d",
295 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
297 log.ensureLineBreak();
305 j_mol += 2 + nral_rt(ftype);
310 return stream.toString();
313 /*! \brief Help print error output when interactions are missing */
314 static void print_missing_interactions_atoms(const gmx::MDLogger &mdlog,
316 const gmx_mtop_t *mtop,
319 const gmx_reverse_top_t *rt = cr->dd->reverse_top;
321 /* Print the atoms in the missing interactions per molblock */
323 for (const gmx_molblock_t &molb : mtop->molblock)
325 const gmx_moltype_t &moltype = mtop->moltype[molb.type];
327 a_end = a_start + molb.nmol*moltype.atoms.nr;
329 GMX_LOG(mdlog.warning).appendText(
330 print_missing_interactions_mb(cr, rt,
332 &rt->ril_mt[molb.type],
333 a_start, a_end, moltype.atoms.nr,
339 void dd_print_missing_interactions(const gmx::MDLogger &mdlog,
342 const gmx_mtop_t *top_global,
343 const gmx_localtop_t *top_local,
344 const t_state *state_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,
410 -1, as_rvec_array(state_local->x.data()), state_local->box);
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,
498 gmx::ArrayRef < const std::vector < int>> vsitePbc,
500 gmx_bool bConstr, gmx_bool bSettle,
502 gmx::ArrayRef<const int> r_index,
503 gmx::ArrayRef<int> r_il,
504 gmx_bool bLinkToAllAtoms,
507 int ftype, j, nlink, link;
512 for (ftype = 0; ftype < F_NRE; ftype++)
514 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
515 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
516 (bSettle && ftype == F_SETTLE))
518 const bool bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0u);
519 const int nral = NRAL(ftype);
520 const auto &il = il_mt[ftype];
521 for (int i = 0; i < il.size(); i += 1+nral)
523 const int* ia = il.iatoms.data() + i;
528 /* We don't need the virtual sites for the cg-links */
538 /* Couple to the first atom in the interaction */
541 for (link = 0; link < nlink; link++)
546 GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
547 GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
548 r_il[r_index[a]+count[a]] =
549 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
550 r_il[r_index[a]+count[a]+1] = ia[0];
551 for (j = 1; j < 1+nral; j++)
553 /* Store the molecular atom number */
554 r_il[r_index[a]+count[a]+1+j] = ia[j];
557 if (interaction_function[ftype].flags & IF_VSITE)
561 /* Add an entry to iatoms for storing
562 * which of the constructing atoms are
565 r_il[r_index[a]+count[a]+2+nral] = 0;
566 for (j = 2; j < 1+nral; j++)
568 if (atom[ia[j]].ptype == eptVSite)
570 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
573 /* Store vsite pbc atom in a second extra entry */
574 r_il[r_index[a]+count[a]+2+nral+1] =
575 (!vsitePbc.empty() ? vsitePbc[ftype-F_VSITE2][i/(1+nral)] : -2);
580 /* We do not count vsites since they are always
581 * uniquely assigned and can be assigned
582 * to multiple nodes with recursive vsites.
585 !(interaction_function[ftype].flags & IF_LIMZERO))
590 count[a] += 2 + nral_rt(ftype);
599 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
600 static int make_reverse_ilist(const InteractionLists &ilist,
601 const t_atoms *atoms,
602 gmx::ArrayRef < const std::vector < int>> vsitePbc,
603 gmx_bool bConstr, gmx_bool bSettle,
605 gmx_bool bLinkToAllAtoms,
606 reverse_ilist_t *ril_mt)
608 int nat_mt, *count, i, nint_mt;
610 /* Count the interactions */
613 low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
615 bConstr, bSettle, bBCheck,
616 gmx::EmptyArrayRef(), gmx::EmptyArrayRef(),
617 bLinkToAllAtoms, FALSE);
619 ril_mt->index.push_back(0);
620 for (i = 0; i < nat_mt; i++)
622 ril_mt->index.push_back(ril_mt->index[i] + count[i]);
625 ril_mt->il.resize(ril_mt->index[nat_mt]);
627 /* Store the interactions */
629 low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
631 bConstr, bSettle, bBCheck,
632 ril_mt->index, ril_mt->il,
633 bLinkToAllAtoms, TRUE);
637 ril_mt->numAtomsInMolecule = atoms->nr;
642 /*! \brief Generate the reverse topology */
643 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
644 gmx::ArrayRef<const VsitePbc> vsitePbcPerMoltype,
645 gmx_bool bConstr, gmx_bool bSettle,
646 gmx_bool bBCheck, int *nint)
648 gmx_reverse_top_t rt;
650 /* Should we include constraints (for SHAKE) in rt? */
651 rt.bConstr = bConstr;
652 rt.bSettle = bSettle;
653 rt.bBCheck = bBCheck;
655 rt.bInterCGInteractions = mtop->bIntermolecularInteractions;
656 rt.ril_mt.resize(mtop->moltype.size());
657 rt.ril_mt_tot_size = 0;
658 std::vector<int> nint_mt;
659 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
661 const gmx_moltype_t &molt = mtop->moltype[mt];
664 rt.bInterCGInteractions = true;
667 /* Make the atom to interaction list for this molecule type */
668 gmx::ArrayRef < const std::vector < int>> vsitePbc;
669 if (!vsitePbcPerMoltype.empty())
671 vsitePbc = gmx::makeConstArrayRef(vsitePbcPerMoltype[mt]);
673 int numberOfInteractions =
674 make_reverse_ilist(molt.ilist, &molt.atoms, vsitePbc,
675 rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
677 nint_mt.push_back(numberOfInteractions);
679 rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
683 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt.ril_mt_tot_size);
687 for (const gmx_molblock_t &molblock : mtop->molblock)
689 *nint += molblock.nmol*nint_mt[molblock.type];
692 /* Make an intermolecular reverse top, if necessary */
693 rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
694 if (rt.bIntermolecularInteractions)
696 t_atoms atoms_global;
698 atoms_global.nr = mtop->natoms;
699 atoms_global.atom = nullptr; /* Only used with virtual sites */
701 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist.get(), "We should have an ilist when intermolecular interactions are on");
704 make_reverse_ilist(*mtop->intermolecular_ilist,
706 gmx::EmptyArrayRef(),
707 rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
711 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
713 rt.ilsort = ilsortFE_UNSORTED;
717 rt.ilsort = ilsortNO_FE;
720 /* Make a molblock index for fast searching */
722 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
724 const gmx_molblock_t &molb = mtop->molblock[mb];
725 const int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
728 i += molb.nmol*numAtomsPerMol;
730 mbi.natoms_mol = numAtomsPerMol;
731 mbi.type = molb.type;
732 rt.mbi.push_back(mbi);
735 rt.th_work.resize(gmx_omp_nthreads_get(emntDomdec));
736 if (!vsitePbcPerMoltype.empty())
738 for (thread_work_t &th_work : rt.th_work)
740 th_work.vsitePbc = gmx::compat::make_unique<VsitePbc>();
747 void dd_make_reverse_top(FILE *fplog,
748 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
749 const gmx_vsite_t *vsite,
750 const t_inputrec *ir, gmx_bool bBCheck)
754 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
757 /* If normal and/or settle constraints act only within charge groups,
758 * we can store them in the reverse top and simply assign them to domains.
759 * Otherwise we need to assign them to multiple domains and set up
760 * the parallel version constraint algorithm(s).
763 gmx::ArrayRef<const VsitePbc> vsitePbcPerMoltype;
766 vsitePbcPerMoltype = gmx::makeConstArrayRef(vsite->vsite_pbc_molt);
769 dd->reverse_top = new gmx_reverse_top_t;
771 make_reverse_top(mtop, ir->efep != efepNO, vsitePbcPerMoltype,
772 !dd->bInterCGcons, !dd->bInterCGsettles,
773 bBCheck, &dd->nbonded_global);
775 gmx_reverse_top_t *rt = dd->reverse_top;
777 /* With the Verlet scheme, exclusions are handled in the non-bonded
778 * kernels and only exclusions inside the cut-off lead to exclusion
779 * forces. Since each atom pair is treated at most once in the non-bonded
780 * kernels, it doesn't matter if the exclusions for the same atom pair
781 * appear multiple times in the exclusion list. In contrast, the, old,
782 * group cut-off scheme loops over a list of exclusions, so there each
783 * excluded pair should appear exactly once.
785 rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
786 inputrecExclForces(ir));
789 dd->n_intercg_excl = 0;
790 rt->n_excl_at_max = 0;
791 for (const gmx_molblock_t &molb : mtop->molblock)
793 int n_excl_mol, n_excl_icg, n_excl_at_max;
795 const gmx_moltype_t &molt = mtop->moltype[molb.type];
796 count_excls(&molt.cgs, &molt.excls,
797 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
798 nexcl += molb.nmol*n_excl_mol;
799 dd->n_intercg_excl += molb.nmol*n_excl_icg;
800 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
802 if (rt->bExclRequired)
804 dd->nbonded_global += nexcl;
805 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
807 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
808 "will use an extra communication step for exclusion forces for %s\n",
809 dd->n_intercg_excl, eel_names[ir->coulombtype]);
813 if (vsite && vsite->n_intercg_vsite > 0)
817 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
818 "will an extra communication step for selected coordinates and forces\n",
819 vsite->n_intercg_vsite);
821 init_domdec_vsites(dd, vsite->n_intercg_vsite);
824 if (dd->bInterCGcons || dd->bInterCGsettles)
826 init_domdec_constraints(dd, mtop);
830 fprintf(fplog, "\n");
834 /*! \brief Store a vsite interaction at the end of \p il
836 * This routine is very similar to add_ifunc, but vsites interactions
837 * have more work to do than other kinds of interactions, and the
838 * complex way nral (and thus vector contents) depends on ftype
839 * confuses static analysis tools unless we fuse the vsite
840 * atom-indexing organization code with the ifunc-adding code, so that
841 * they can see that nral is the same value. */
843 add_ifunc_for_vsites(t_iatom *tiatoms, const gmx_ga2la_t &ga2la,
844 int nral, gmx_bool bHomeA,
845 int a, int a_gl, int a_mol,
846 const t_iatom *iatoms,
851 if (il->nr+1+nral > il->nalloc)
853 il->nalloc = over_alloc_large(il->nr+1+nral);
854 srenew(il->iatoms, il->nalloc);
856 liatoms = il->iatoms + il->nr;
860 tiatoms[0] = iatoms[0];
864 /* We know the local index of the first atom */
869 /* Convert later in make_local_vsites */
870 tiatoms[1] = -a_gl - 1;
873 for (int k = 2; k < 1+nral; k++)
875 int ak_gl = a_gl + iatoms[k] - a_mol;
876 if (const int *homeIndex = ga2la.findHome(ak_gl))
878 tiatoms[k] = *homeIndex;
882 /* Copy the global index, convert later in make_local_vsites */
883 tiatoms[k] = -(ak_gl + 1);
885 // Note that ga2la_get_home always sets the third parameter if
888 for (int k = 0; k < 1+nral; k++)
890 liatoms[k] = tiatoms[k];
894 /*! \brief Store a bonded interaction at the end of \p il */
895 static inline void add_ifunc(int nral, const t_iatom *tiatoms, t_ilist *il)
900 if (il->nr+1+nral > il->nalloc)
902 il->nalloc = over_alloc_large(il->nr+1+nral);
903 srenew(il->iatoms, il->nalloc);
905 liatoms = il->iatoms + il->nr;
906 for (k = 0; k <= nral; k++)
908 liatoms[k] = tiatoms[k];
913 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
914 static void add_posres(int mol, int a_mol, int numAtomsInMolecule,
915 const gmx_molblock_t *molb,
916 t_iatom *iatoms, const t_iparams *ip_in,
922 /* This position restraint has not been added yet,
923 * so it's index is the current number of position restraints.
925 n = idef->il[F_POSRES].nr/2;
926 if (n+1 > idef->iparams_posres_nalloc)
928 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
929 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
931 ip = &idef->iparams_posres[n];
932 /* Copy the force constants */
933 *ip = ip_in[iatoms[0]];
935 /* Get the position restraint coordinates from the molblock */
936 a_molb = mol*numAtomsInMolecule + a_mol;
937 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
938 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
939 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
940 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
941 if (!molb->posres_xB.empty())
943 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
944 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
945 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
949 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
950 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
951 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
953 /* Set the parameter index for idef->iparams_posre */
957 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
958 static void add_fbposres(int mol, int a_mol, int numAtomsInMolecule,
959 const gmx_molblock_t *molb,
960 t_iatom *iatoms, const t_iparams *ip_in,
966 /* This flat-bottom position restraint has not been added yet,
967 * so it's index is the current number of position restraints.
969 n = idef->il[F_FBPOSRES].nr/2;
970 if (n+1 > idef->iparams_fbposres_nalloc)
972 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
973 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
975 ip = &idef->iparams_fbposres[n];
976 /* Copy the force constants */
977 *ip = ip_in[iatoms[0]];
979 /* Get the position restraint coordinats from the molblock */
980 a_molb = mol*numAtomsInMolecule + a_mol;
981 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
982 /* Take reference positions from A position of normal posres */
983 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
984 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
985 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
987 /* Note: no B-type for flat-bottom posres */
989 /* Set the parameter index for idef->iparams_posre */
993 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
994 static void add_vsite(const gmx_ga2la_t &ga2la,
995 gmx::ArrayRef<const int> index,
996 gmx::ArrayRef<const int> rtil,
998 gmx_bool bHomeA, int a, int a_gl, int a_mol,
999 const t_iatom *iatoms,
1004 t_iatom tiatoms[1+MAXATOMLIST];
1005 int j, ftype_r, nral_r;
1007 /* Add this interaction to the local topology */
1008 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1012 std::vector<int> &vsitePbcFtype = (*vsitePbc)[ftype - c_ftypeVsiteStart];
1013 const int vsi = idef->il[ftype].nr/(1+nral) - 1;
1014 if (static_cast<size_t>(vsi) >= vsitePbcFtype.size())
1016 vsitePbcFtype.resize(vsi + 1);
1020 pbc_a_mol = iatoms[1+nral+1];
1023 /* The pbc flag is one of the following two options:
1024 * -2: vsite and all constructing atoms are within the same cg, no pbc
1025 * -1: vsite and its first constructing atom are in the same cg, do pbc
1027 vsitePbcFtype[vsi] = pbc_a_mol;
1031 /* Set the pbc atom for this vsite so we can make its pbc
1032 * identical to the rest of the atoms in its charge group.
1033 * Since the order of the atoms does not change within a charge
1034 * group, we do not need the global to local atom index.
1036 vsitePbcFtype[vsi] = a + pbc_a_mol - iatoms[1];
1041 /* This vsite is non-home (required for recursion),
1042 * and therefore there is no charge group to match pbc with.
1043 * But we always turn on full_pbc to assure that higher order
1044 * recursion works correctly.
1046 vsitePbcFtype[vsi] = -1;
1052 /* Check for recursion */
1053 for (k = 2; k < 1+nral; k++)
1055 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1057 /* This construction atoms is a vsite and not a home atom */
1060 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1062 /* Find the vsite construction */
1064 /* Check all interactions assigned to this atom */
1065 j = index[iatoms[k]];
1066 while (j < index[iatoms[k]+1])
1068 ftype_r = rtil[j++];
1069 nral_r = NRAL(ftype_r);
1070 if (interaction_function[ftype_r].flags & IF_VSITE)
1072 /* Add this vsite (recursion) */
1073 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1074 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1077 j += 1 + nral_r + 2;
1089 /*! \brief Build the index that maps each local atom to its local atom group */
1090 static void makeLocalAtomGroupsFromAtoms(gmx_domdec_t *dd)
1092 const gmx::RangePartitioning &atomGrouping = dd->atomGrouping();
1094 dd->localAtomGroupFromAtom.clear();
1096 for (size_t g = 0; g < dd->globalAtomGroupIndices.size(); g++)
1098 for (int gmx_unused a : atomGrouping.block(g))
1100 dd->localAtomGroupFromAtom.push_back(g);
1105 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1106 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1112 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1116 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1122 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1123 static void combine_blocka(t_blocka *dest,
1124 gmx::ArrayRef<const thread_work_t> src)
1126 int ni = src.back().excl.nr;
1128 for (const thread_work_t &th_work : src)
1130 na += th_work.excl.nra;
1132 if (ni + 1 > dest->nalloc_index)
1134 dest->nalloc_index = over_alloc_large(ni+1);
1135 srenew(dest->index, dest->nalloc_index);
1137 if (dest->nra + na > dest->nalloc_a)
1139 dest->nalloc_a = over_alloc_large(dest->nra+na);
1140 srenew(dest->a, dest->nalloc_a);
1142 for (gmx::index s = 1; s < src.size(); s++)
1144 for (int i = dest->nr + 1; i < src[s].excl.nr + 1; i++)
1146 dest->index[i] = dest->nra + src[s].excl.index[i];
1148 for (int i = 0; i < src[s].excl.nra; i++)
1150 dest->a[dest->nra+i] = src[s].excl.a[i];
1152 dest->nr = src[s].excl.nr;
1153 dest->nra += src[s].excl.nra;
1157 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1158 * virtual sites need special attention, as pbc info differs per vsite.
1160 static void combine_idef(t_idef *dest,
1161 gmx::ArrayRef<const thread_work_t> src,
1166 for (ftype = 0; ftype < F_NRE; ftype++)
1169 for (gmx::index s = 1; s < src.size(); s++)
1171 n += src[s].idef.il[ftype].nr;
1177 ild = &dest->il[ftype];
1179 if (ild->nr + n > ild->nalloc)
1181 ild->nalloc = over_alloc_large(ild->nr+n);
1182 srenew(ild->iatoms, ild->nalloc);
1186 (((interaction_function[ftype].flags & IF_VSITE) != 0u) &&
1187 vsite->vsite_pbc_loc);
1188 const int nral1 = 1 + NRAL(ftype);
1189 const int ftv = ftype - c_ftypeVsiteStart;
1191 for (gmx::index s = 1; s < src.size(); s++)
1193 const t_ilist &ils = src[s].idef.il[ftype];
1195 for (int i = 0; i < ils.nr; i++)
1197 ild->iatoms[ild->nr + i] = ils.iatoms[i];
1201 const std::vector<int> &pbcSrc = (*src[s].vsitePbc)[ftv];
1202 std::vector<int> &pbcDest = (*vsite->vsite_pbc_loc)[ftv];
1203 pbcDest.resize((ild->nr + ils.nr)/nral1);
1204 for (int i = 0; i < ils.nr; i += nral1)
1206 pbcDest[(ild->nr + i)/nral1] = pbcSrc[i/nral1];
1213 /* Position restraints need an additional treatment */
1214 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1216 int nposres = dest->il[ftype].nr/2;
1217 // TODO: Simplify this code using std::vector
1218 t_iparams * &iparams_dest = (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1219 int &posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc : dest->iparams_fbposres_nalloc);
1220 if (nposres > posres_nalloc)
1222 posres_nalloc = over_alloc_large(nposres);
1223 srenew(iparams_dest, posres_nalloc);
1226 /* Set nposres to the number of original position restraints in dest */
1227 for (gmx::index s = 1; s < src.size(); s++)
1229 nposres -= src[s].idef.il[ftype].nr/2;
1232 for (gmx::index s = 1; s < src.size(); s++)
1234 const t_iparams *iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1236 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1238 /* Correct the index into iparams_posres */
1239 dest->il[ftype].iatoms[nposres*2] = nposres;
1240 /* Copy the position restraint force parameters */
1241 iparams_dest[nposres] = iparams_src[i];
1250 /*! \brief Check and when available assign bonded interactions for local atom i
1253 check_assign_interactions_atom(int i, int i_gl,
1255 int numAtomsInMolecule,
1256 gmx::ArrayRef<const int> index,
1257 gmx::ArrayRef<const int> rtil,
1258 gmx_bool bInterMolInteractions,
1259 int ind_start, int ind_end,
1260 const gmx_domdec_t *dd,
1261 const gmx_domdec_zones_t *zones,
1262 const gmx_molblock_t *molb,
1263 gmx_bool bRCheckMB, const ivec rcheck, gmx_bool bRCheck2B,
1268 const t_iparams *ip_in,
1282 t_iatom tiatoms[1 + MAXATOMLIST];
1285 gmx::ArrayRef<const t_iatom> iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1287 if (ftype == F_SETTLE)
1289 /* Settles are only in the reverse top when they
1290 * operate within a charge group. So we can assign
1291 * them without checks. We do this only for performance
1292 * reasons; it could be handled by the code below.
1296 /* Home zone: add this settle to the local topology */
1297 tiatoms[0] = iatoms[0];
1299 tiatoms[2] = i + iatoms[2] - iatoms[1];
1300 tiatoms[3] = i + iatoms[3] - iatoms[1];
1301 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1306 else if (interaction_function[ftype].flags & IF_VSITE)
1308 assert(!bInterMolInteractions);
1309 /* The vsite construction goes where the vsite itself is */
1312 add_vsite(*dd->ga2la, index, rtil, ftype, nral,
1313 TRUE, i, i_gl, i_mol,
1314 iatoms.data(), idef, vsitePbc);
1323 tiatoms[0] = iatoms[0];
1327 assert(!bInterMolInteractions);
1328 /* Assign single-body interactions to the home zone */
1333 if (ftype == F_POSRES)
1335 add_posres(mol, i_mol, numAtomsInMolecule,
1336 molb, tiatoms, ip_in, idef);
1338 else if (ftype == F_FBPOSRES)
1340 add_fbposres(mol, i_mol, numAtomsInMolecule,
1341 molb, tiatoms, ip_in, idef);
1351 /* This is a two-body interaction, we can assign
1352 * analogous to the non-bonded assignments.
1356 if (!bInterMolInteractions)
1358 /* Get the global index using the offset in the molecule */
1359 k_gl = i_gl + iatoms[2] - i_mol;
1365 if (const auto *entry = dd->ga2la->find(k_gl))
1367 int kz = entry->cell;
1372 /* Check zone interaction assignments */
1373 bUse = ((iz < zones->nizone &&
1375 kz >= zones->izone[iz].j0 &&
1376 kz < zones->izone[iz].j1) ||
1377 (kz < zones->nizone &&
1379 iz >= zones->izone[kz].j0 &&
1380 iz < zones->izone[kz].j1));
1384 tiatoms[2] = entry->la;
1385 /* If necessary check the cgcm distance */
1387 dd_dist2(pbc_null, cg_cm, la2lc,
1388 tiatoms[1], tiatoms[2]) >= rc2)
1401 /* Assign this multi-body bonded interaction to
1402 * the local node if we have all the atoms involved
1403 * (local or communicated) and the minimum zone shift
1404 * in each dimension is zero, for dimensions
1405 * with 2 DD cells an extra check may be necessary.
1407 ivec k_zero, k_plus;
1413 for (k = 1; k <= nral && bUse; k++)
1416 if (!bInterMolInteractions)
1418 /* Get the global index using the offset in the molecule */
1419 k_gl = i_gl + iatoms[k] - i_mol;
1425 const auto *entry = dd->ga2la->find(k_gl);
1426 if (entry == nullptr || entry->cell >= zones->n)
1428 /* We do not have this atom of this interaction
1429 * locally, or it comes from more than one cell
1438 tiatoms[k] = entry->la;
1439 for (d = 0; d < DIM; d++)
1441 if (zones->shift[entry->cell][d] == 0)
1453 (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1458 for (d = 0; (d < DIM && bUse); d++)
1460 /* Check if the cg_cm distance falls within
1461 * the cut-off to avoid possible multiple
1462 * assignments of bonded interactions.
1466 dd_dist2(pbc_null, cg_cm, la2lc,
1467 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1476 /* Add this interaction to the local topology */
1477 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1478 /* Sum so we can check in global_stat
1479 * if we have everything.
1482 !(interaction_function[ftype].flags & IF_LIMZERO))
1492 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1494 * With thread parallelizing each thread acts on a different atom range:
1495 * at_start to at_end.
1497 static int make_bondeds_zone(gmx_domdec_t *dd,
1498 const gmx_domdec_zones_t *zones,
1499 const std::vector<gmx_molblock_t> &molb,
1500 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1502 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1503 const t_iparams *ip_in,
1507 gmx::RangePartitioning::Block atomRange)
1509 int mb, mt, mol, i_mol;
1511 gmx_reverse_top_t *rt;
1514 rt = dd->reverse_top;
1516 bBCheck = rt->bBCheck;
1520 for (int i : atomRange)
1522 /* Get the global atom number */
1523 const int i_gl = dd->globalAtomIndices[i];
1524 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1525 /* Check all intramolecular interactions assigned to this atom */
1526 gmx::ArrayRef<const int> index = rt->ril_mt[mt].index;
1527 gmx::ArrayRef<const t_iatom> rtil = rt->ril_mt[mt].il;
1529 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1530 rt->ril_mt[mt].numAtomsInMolecule,
1532 index[i_mol], index[i_mol+1],
1535 bRCheckMB, rcheck, bRCheck2B, rc2,
1546 if (rt->bIntermolecularInteractions)
1548 /* Check all intermolecular interactions assigned to this atom */
1549 index = rt->ril_intermol.index;
1550 rtil = rt->ril_intermol.il;
1552 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1553 rt->ril_mt[mt].numAtomsInMolecule,
1555 index[i_gl], index[i_gl + 1],
1558 bRCheckMB, rcheck, bRCheck2B, rc2,
1570 return nbonded_local;
1573 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1574 static void set_no_exclusions_zone(const gmx_domdec_t *dd,
1575 const gmx_domdec_zones_t *zones,
1579 const auto zone = dd->atomGrouping().subRange(zones->cg_range[iz],
1580 zones->cg_range[iz + 1]);
1584 lexcls->index[a + 1] = lexcls->nra;
1588 /*! \brief Set the exclusion data for i-zone \p iz
1590 * This is a legacy version for the group scheme of the same routine below.
1591 * Here charge groups and distance checks to ensure unique exclusions
1594 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1595 const std::vector<gmx_moltype_t> &moltype,
1596 gmx_bool bRCheck, real rc2,
1597 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1601 int cg_start, int cg_end)
1605 const t_blocka *excls;
1607 const gmx_ga2la_t &ga2la = *dd->ga2la;
1610 dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1611 zones->izone[iz].jcg1);
1613 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1615 /* We set the end index, but note that we might not start at zero here */
1616 lexcls->nr = dd->atomGrouping().subRange(0, cg_end).size();
1618 int n = lexcls->nra;
1620 for (int cg = cg_start; cg < cg_end; cg++)
1622 if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1624 lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1625 srenew(lexcls->a, lexcls->nalloc_a);
1627 const auto atomGroup = dd->atomGrouping().block(cg);
1628 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1629 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1631 /* Copy the exclusions from the global top */
1632 for (int la : atomGroup)
1634 lexcls->index[la] = n;
1635 int a_gl = dd->globalAtomIndices[la];
1637 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1638 excls = &moltype[mt].excls;
1639 for (int j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1641 int aj_mol = excls->a[j];
1642 /* This computation of jla is only correct intra-cg */
1643 int jla = la + aj_mol - a_mol;
1644 if (atomGroup.inRange(jla))
1646 /* This is an intra-cg exclusion. We can skip
1647 * the global indexing and distance checking.
1649 /* Intra-cg exclusions are only required
1650 * for the home zone.
1654 lexcls->a[n++] = jla;
1655 /* Check to avoid double counts */
1664 /* This is a inter-cg exclusion */
1665 /* Since exclusions are pair interactions,
1666 * just like non-bonded interactions,
1667 * they can be assigned properly up
1668 * to the DD cutoff (not cutoff_min as
1669 * for the other bonded interactions).
1671 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1673 if (iz == 0 && jEntry->cell == 0)
1675 lexcls->a[n++] = jEntry->la;
1676 /* Check to avoid double counts */
1677 if (jEntry->la > la)
1682 else if (jRange.inRange(jEntry->la) &&
1684 dd_dist2(pbc_null, cg_cm, la2lc, la, jEntry->la) < rc2))
1686 /* jla > la, since jRange.begin() > la */
1687 lexcls->a[n++] = jEntry->la;
1697 /* There are no inter-cg excls and this cg is self-excluded.
1698 * These exclusions are only required for zone 0,
1699 * since other zones do not see themselves.
1703 for (int la : atomGroup)
1705 lexcls->index[la] = n;
1706 for (int j : atomGroup)
1711 count += (atomGroup.size()*(atomGroup.size() - 1))/2;
1715 /* We don't need exclusions for this cg */
1716 for (int la : atomGroup)
1718 lexcls->index[la] = n;
1724 lexcls->index[lexcls->nr] = n;
1730 /*! \brief Set the exclusion data for i-zone \p iz */
1731 static void make_exclusions_zone(gmx_domdec_t *dd,
1732 gmx_domdec_zones_t *zones,
1733 const std::vector<gmx_moltype_t> &moltype,
1737 int at_start, int at_end)
1739 int n_excl_at_max, n, at;
1741 const gmx_ga2la_t &ga2la = *dd->ga2la;
1744 dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1745 zones->izone[iz].jcg1);
1747 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1749 /* We set the end index, but note that we might not start at zero here */
1750 lexcls->nr = at_end;
1753 for (at = at_start; at < at_end; at++)
1755 if (n + 1000 > lexcls->nalloc_a)
1757 lexcls->nalloc_a = over_alloc_large(n + 1000);
1758 srenew(lexcls->a, lexcls->nalloc_a);
1760 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1762 int a_gl, mb, mt, mol, a_mol, j;
1763 const t_blocka *excls;
1765 if (n + n_excl_at_max > lexcls->nalloc_a)
1767 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1768 srenew(lexcls->a, lexcls->nalloc_a);
1771 /* Copy the exclusions from the global top */
1772 lexcls->index[at] = n;
1773 a_gl = dd->globalAtomIndices[at];
1774 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1775 &mb, &mt, &mol, &a_mol);
1776 excls = &moltype[mt].excls;
1777 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1779 const int aj_mol = excls->a[j];
1781 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1783 /* This check is not necessary, but it can reduce
1784 * the number of exclusions in the list, which in turn
1785 * can speed up the pair list construction a bit.
1787 if (jRange.inRange(jEntry->la))
1789 lexcls->a[n++] = jEntry->la;
1796 /* We don't need exclusions for this atom */
1797 lexcls->index[at] = n;
1801 lexcls->index[lexcls->nr] = n;
1806 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1807 static void check_alloc_index(t_blocka *ba, int nindex_max)
1809 if (nindex_max+1 > ba->nalloc_index)
1811 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1812 srenew(ba->index, ba->nalloc_index);
1816 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1817 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1820 int nr = dd->atomGrouping().subRange(0, zones->izone[zones->nizone - 1].cg1).size();
1822 check_alloc_index(lexcls, nr);
1824 for (size_t thread = 1; thread < dd->reverse_top->th_work.size(); thread++)
1826 check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1830 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1831 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1834 const auto nonhomeIzonesAtomRange =
1835 dd->atomGrouping().subRange(zones->izone[0].cg1,
1836 zones->izone[zones->nizone - 1].cg1);
1838 if (dd->n_intercg_excl == 0)
1840 /* There are no exclusions involving non-home charge groups,
1841 * but we need to set the indices for neighborsearching.
1843 for (int la : nonhomeIzonesAtomRange)
1845 lexcls->index[la] = lexcls->nra;
1848 /* nr is only used to loop over the exclusions for Ewald and RF,
1849 * so we can set it to the number of home atoms for efficiency.
1851 lexcls->nr = nonhomeIzonesAtomRange.begin();
1855 lexcls->nr = nonhomeIzonesAtomRange.end();
1859 /*! \brief Clear a t_idef struct */
1860 static void clear_idef(t_idef *idef)
1864 /* Clear the counts */
1865 for (ftype = 0; ftype < F_NRE; ftype++)
1867 idef->il[ftype].nr = 0;
1871 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1872 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1873 gmx_domdec_zones_t *zones,
1874 const gmx_mtop_t *mtop,
1876 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1878 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1879 t_idef *idef, gmx_vsite_t *vsite,
1880 t_blocka *lexcls, int *excl_count)
1882 int nzone_bondeds, nzone_excl;
1883 int izone, cg0, cg1;
1887 gmx_reverse_top_t *rt;
1889 if (dd->reverse_top->bInterCGInteractions)
1891 nzone_bondeds = zones->n;
1895 /* Only single charge group (or atom) molecules, so interactions don't
1896 * cross zone boundaries and we only need to assign in the home zone.
1901 if (dd->n_intercg_excl > 0)
1903 /* We only use exclusions from i-zones to i- and j-zones */
1904 nzone_excl = zones->nizone;
1908 /* There are no inter-cg exclusions and only zone 0 sees itself */
1912 check_exclusions_alloc(dd, zones, lexcls);
1914 rt = dd->reverse_top;
1918 /* Clear the counts */
1926 for (izone = 0; izone < nzone_bondeds; izone++)
1928 cg0 = zones->cg_range[izone];
1929 cg1 = zones->cg_range[izone + 1];
1931 const int numThreads = rt->th_work.size();
1932 #pragma omp parallel for num_threads(numThreads) schedule(static)
1933 for (thread = 0; thread < numThreads; thread++)
1941 cg0t = cg0 + ((cg1 - cg0)* thread )/numThreads;
1942 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/numThreads;
1950 idef_t = &rt->th_work[thread].idef;
1954 VsitePbc *vsitePbc = nullptr;
1955 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1959 vsitePbc = vsite->vsite_pbc_loc.get();
1963 vsitePbc = rt->th_work[thread].vsitePbc.get();
1967 rt->th_work[thread].nbonded =
1968 make_bondeds_zone(dd, zones,
1970 bRCheckMB, rcheck, bRCheck2B, rc2,
1971 la2lc, pbc_null, cg_cm, idef->iparams,
1974 dd->atomGrouping().subRange(cg0t, cg1t));
1976 if (izone < nzone_excl)
1984 excl_t = &rt->th_work[thread].excl;
1989 if (dd->atomGrouping().allBlocksHaveSizeOne() &&
1992 /* No charge groups and no distance check required */
1993 make_exclusions_zone(dd, zones,
1994 mtop->moltype, cginfo,
2001 rt->th_work[thread].excl_count =
2002 make_exclusions_zone_cg(dd, zones,
2003 mtop->moltype, bRCheck2B, rc2,
2004 la2lc, pbc_null, cg_cm, cginfo,
2011 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2014 if (rt->th_work.size() > 1)
2016 combine_idef(idef, rt->th_work, vsite);
2019 for (const thread_work_t &th_work : rt->th_work)
2021 nbonded_local += th_work.nbonded;
2024 if (izone < nzone_excl)
2026 if (rt->th_work.size() > 1)
2028 combine_blocka(lexcls, rt->th_work);
2031 for (const thread_work_t &th_work : rt->th_work)
2033 *excl_count += th_work.excl_count;
2038 /* Some zones might not have exclusions, but some code still needs to
2039 * loop over the index, so we set the indices here.
2041 for (izone = nzone_excl; izone < zones->nizone; izone++)
2043 set_no_exclusions_zone(dd, zones, izone, lexcls);
2046 finish_local_exclusions(dd, zones, lexcls);
2049 fprintf(debug, "We have %d exclusions, check count %d\n",
2050 lexcls->nra, *excl_count);
2053 return nbonded_local;
2056 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
2058 lcgs->nr = dd->globalAtomGroupIndices.size();
2059 lcgs->index = dd->atomGrouping_.rawIndex().data();
2062 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
2063 int npbcdim, matrix box,
2064 rvec cellsize_min, const ivec npulse,
2068 const gmx_mtop_t *mtop, gmx_localtop_t *ltop)
2070 gmx_bool bRCheckMB, bRCheck2B;
2074 t_pbc pbc, *pbc_null = nullptr;
2078 fprintf(debug, "Making local topology\n");
2081 dd_make_local_cgs(dd, <op->cgs);
2086 if (dd->reverse_top->bInterCGInteractions)
2088 /* We need to check to which cell bondeds should be assigned */
2089 rc = dd_cutoff_twobody(dd);
2092 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2095 /* Should we check cg_cm distances when assigning bonded interactions? */
2096 for (d = 0; d < DIM; d++)
2099 /* Only need to check for dimensions where the part of the box
2100 * that is not communicated is smaller than the cut-off.
2102 if (d < npbcdim && dd->nc[d] > 1 &&
2103 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2110 /* Check for interactions between two atoms,
2111 * where we can allow interactions up to the cut-off,
2112 * instead of up to the smallest cell dimension.
2119 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n",
2120 d, cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
2123 if (bRCheckMB || bRCheck2B)
2125 makeLocalAtomGroupsFromAtoms(dd);
2128 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
2138 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
2139 bRCheckMB, rcheck, bRCheck2B, rc,
2140 dd->localAtomGroupFromAtom.data(),
2141 pbc_null, cgcm_or_x,
2143 <op->excls, &nexcl);
2145 /* The ilist is not sorted yet,
2146 * we can only do this when we have the charge arrays.
2148 ltop->idef.ilsort = ilsortUNKNOWN;
2150 if (dd->reverse_top->bExclRequired)
2152 dd->nbonded_local += nexcl;
2155 ltop->atomtypes = mtop->atomtypes;
2158 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
2159 gmx_localtop_t *ltop)
2161 if (dd->reverse_top->ilsort == ilsortNO_FE)
2163 ltop->idef.ilsort = ilsortNO_FE;
2167 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
2171 gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global)
2173 gmx_localtop_t *top;
2177 top->idef.ntypes = top_global->ffparams.ntypes;
2178 top->idef.atnr = top_global->ffparams.atnr;
2179 top->idef.functype = top_global->ffparams.functype;
2180 top->idef.iparams = top_global->ffparams.iparams;
2181 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
2182 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
2184 top->idef.ilsort = ilsortUNKNOWN;
2189 void dd_init_local_state(gmx_domdec_t *dd,
2190 const t_state *state_global, t_state *state_local)
2192 int buf[NITEM_DD_INIT_LOCAL_STATE];
2196 buf[0] = state_global->flags;
2197 buf[1] = state_global->ngtc;
2198 buf[2] = state_global->nnhpres;
2199 buf[3] = state_global->nhchainlength;
2200 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
2202 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2204 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
2205 init_dfhist_state(state_local, buf[4]);
2206 state_local->flags = buf[0];
2209 /*! \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 */
2210 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2216 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2218 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2219 if (link->a[k] == cg_gl_j)
2226 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2227 "Inconsistent allocation of link");
2228 /* Add this charge group link */
2229 if (link->index[cg_gl+1]+1 > link->nalloc_a)
2231 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2232 srenew(link->a, link->nalloc_a);
2234 link->a[link->index[cg_gl+1]] = cg_gl_j;
2235 link->index[cg_gl+1]++;
2239 /*! \brief Return a vector of the charge group index for all atoms */
2240 static std::vector<int> make_at2cg(const t_block &cgs)
2242 std::vector<int> at2cg(cgs.index[cgs.nr]);
2243 for (int cg = 0; cg < cgs.nr; cg++)
2245 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2254 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
2255 cginfo_mb_t *cginfo_mb)
2257 gmx_bool bExclRequired;
2259 cginfo_mb_t *cgi_mb;
2261 /* For each charge group make a list of other charge groups
2262 * in the system that a linked to it via bonded interactions
2263 * which are also stored in reverse_top.
2266 bExclRequired = dd->reverse_top->bExclRequired;
2268 reverse_ilist_t ril_intermol;
2269 if (mtop->bIntermolecularInteractions)
2271 if (ncg_mtop(mtop) < mtop->natoms)
2273 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.");
2278 atoms.nr = mtop->natoms;
2279 atoms.atom = nullptr;
2281 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist.get(), "We should have an ilist when intermolecular interactions are on");
2283 make_reverse_ilist(*mtop->intermolecular_ilist,
2285 gmx::EmptyArrayRef(),
2286 FALSE, FALSE, FALSE, TRUE, &ril_intermol);
2290 snew(link->index, ncg_mtop(mtop)+1);
2297 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
2299 const gmx_molblock_t &molb = mtop->molblock[mb];
2304 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2305 const t_block &cgs = molt.cgs;
2306 const t_blocka &excls = molt.excls;
2307 std::vector<int> a2c = make_at2cg(cgs);
2308 /* Make a reverse ilist in which the interactions are linked
2309 * to all atoms, not only the first atom as in gmx_reverse_top.
2310 * The constraints are discarded here.
2312 reverse_ilist_t ril;
2313 make_reverse_ilist(molt.ilist, &molt.atoms, gmx::EmptyArrayRef(),
2314 FALSE, FALSE, FALSE, TRUE, &ril);
2316 cgi_mb = &cginfo_mb[mb];
2319 for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb.nmol : 1); mol++)
2321 for (int cg = 0; cg < cgs.nr; cg++)
2323 int cg_gl = cg_offset + cg;
2324 link->index[cg_gl+1] = link->index[cg_gl];
2325 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2327 int i = ril.index[a];
2328 while (i < ril.index[a+1])
2330 int ftype = ril.il[i++];
2331 int nral = NRAL(ftype);
2332 /* Skip the ifunc index */
2334 for (int j = 0; j < nral; j++)
2336 int aj = ril.il[i + j];
2339 check_link(link, cg_gl, cg_offset+a2c[aj]);
2342 i += nral_rt(ftype);
2346 /* Exclusions always go both ways */
2347 for (int j = excls.index[a]; j < excls.index[a + 1]; j++)
2349 int aj = excls.a[j];
2352 check_link(link, cg_gl, cg_offset+a2c[aj]);
2357 if (mtop->bIntermolecularInteractions)
2359 int i = ril_intermol.index[a];
2360 while (i < ril_intermol.index[a+1])
2362 int ftype = ril_intermol.il[i++];
2363 int nral = NRAL(ftype);
2364 /* Skip the ifunc index */
2366 for (int j = 0; j < nral; j++)
2368 /* Here we assume we have no charge groups;
2369 * this has been checked above.
2371 int aj = ril_intermol.il[i + j];
2372 check_link(link, cg_gl, aj);
2374 i += nral_rt(ftype);
2378 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2380 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2385 cg_offset += cgs.nr;
2387 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - cgs.nr];
2391 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt.name, cgs.nr, nlink_mol);
2394 if (molb.nmol > mol)
2396 /* Copy the data for the rest of the molecules in this block */
2397 link->nalloc_a += (molb.nmol - mol)*nlink_mol;
2398 srenew(link->a, link->nalloc_a);
2399 for (; mol < molb.nmol; mol++)
2401 for (int cg = 0; cg < cgs.nr; cg++)
2403 int cg_gl = cg_offset + cg;
2404 link->index[cg_gl + 1] =
2405 link->index[cg_gl + 1 - cgs.nr] + nlink_mol;
2406 for (int j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2408 link->a[j] = link->a[j - nlink_mol] + cgs.nr;
2410 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2411 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2413 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2417 cg_offset += cgs.nr;
2424 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2435 } bonded_distance_t;
2437 /*! \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 */
2438 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2439 bonded_distance_t *bd)
2450 /*! \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 */
2451 static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
2452 const std::vector<int> &at2cg,
2453 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2454 bonded_distance_t *bd_2b,
2455 bonded_distance_t *bd_mb)
2457 for (int ftype = 0; ftype < F_NRE; ftype++)
2459 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2461 const auto &il = molt->ilist[ftype];
2462 int nral = NRAL(ftype);
2465 for (int i = 0; i < il.size(); i += 1+nral)
2467 for (int ai = 0; ai < nral; ai++)
2469 int cgi = at2cg[il.iatoms[i+1+ai]];
2470 for (int aj = ai + 1; aj < nral; aj++)
2472 int cgj = at2cg[il.iatoms[i+1+aj]];
2475 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2477 update_max_bonded_distance(rij2, ftype,
2480 (nral == 2) ? bd_2b : bd_mb);
2490 const t_blocka *excls = &molt->excls;
2491 for (int ai = 0; ai < excls->nr; ai++)
2493 int cgi = at2cg[ai];
2494 for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2496 int cgj = at2cg[excls->a[j]];
2499 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2501 /* There is no function type for exclusions, use -1 */
2502 update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2509 /*! \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 */
2510 static void bonded_distance_intermol(const InteractionLists &ilists_intermol,
2512 const rvec *x, int ePBC, const matrix box,
2513 bonded_distance_t *bd_2b,
2514 bonded_distance_t *bd_mb)
2518 set_pbc(&pbc, ePBC, box);
2520 for (int ftype = 0; ftype < F_NRE; ftype++)
2522 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2524 const auto &il = ilists_intermol[ftype];
2525 int nral = NRAL(ftype);
2527 /* No nral>1 check here, since intermol interactions always
2528 * have nral>=2 (and the code is also correct for nral=1).
2530 for (int i = 0; i < il.size(); i += 1+nral)
2532 for (int ai = 0; ai < nral; ai++)
2534 int atom_i = il.iatoms[i + 1 + ai];
2536 for (int aj = ai + 1; aj < nral; aj++)
2541 int atom_j = il.iatoms[i + 1 + aj];
2543 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2547 update_max_bonded_distance(rij2, ftype,
2549 (nral == 2) ? bd_2b : bd_mb);
2557 //! Returns whether \p molt has at least one virtual site
2558 static bool moltypeHasVsite(const gmx_moltype_t &molt)
2560 bool hasVsite = false;
2561 for (int i = 0; i < F_NRE; i++)
2563 if ((interaction_function[i].flags & IF_VSITE) &&
2564 molt.ilist[i].size() > 0)
2573 //! Compute charge group centers of mass for molecule \p molt
2574 static void get_cgcm_mol(const gmx_moltype_t *molt,
2575 const gmx_ffparams_t *ffparams,
2576 int ePBC, t_graph *graph, const matrix box,
2577 const rvec *x, rvec *xs, rvec *cg_cm)
2581 if (ePBC != epbcNONE)
2583 mk_mshift(nullptr, graph, ePBC, box, x);
2585 shift_x(graph, box, x, xs);
2586 /* By doing an extra mk_mshift the molecules that are broken
2587 * because they were e.g. imported from another software
2588 * will be made whole again. Such are the healing powers
2591 mk_mshift(nullptr, graph, ePBC, box, xs);
2595 /* We copy the coordinates so the original coordinates remain
2596 * unchanged, just to be 100% sure that we do not affect
2597 * binary reproducibility of simulations.
2599 n = molt->cgs.index[molt->cgs.nr];
2600 for (i = 0; i < n; i++)
2602 copy_rvec(x[i], xs[i]);
2606 if (moltypeHasVsite(*molt))
2608 /* Convert to old, deprecated format */
2609 t_ilist ilist[F_NRE];
2610 for (int ftype = 0; ftype < F_NRE; ftype++)
2612 if (interaction_function[ftype].flags & IF_VSITE)
2614 ilist[ftype].nr = molt->ilist[ftype].size();
2615 ilist[ftype].iatoms = const_cast<int *>(molt->ilist[ftype].iatoms.data());
2619 construct_vsites(nullptr, xs, 0.0, nullptr,
2620 ffparams->iparams, ilist,
2621 epbcNONE, TRUE, nullptr, nullptr);
2624 calc_cgcm(nullptr, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2627 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
2628 const gmx_mtop_t *mtop,
2629 const t_inputrec *ir,
2630 const rvec *x, const matrix box,
2632 real *r_2b, real *r_mb)
2634 gmx_bool bExclRequired;
2638 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2639 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2641 bExclRequired = inputrecExclForces(ir);
2646 for (const gmx_molblock_t &molb : mtop->molblock)
2648 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2649 if (molt.cgs.nr == 1 || molb.nmol == 0)
2651 at_offset += molb.nmol*molt.atoms.nr;
2655 if (ir->ePBC != epbcNONE)
2657 mk_graph_moltype(molt, &graph);
2660 std::vector<int> at2cg = make_at2cg(molt.cgs);
2661 snew(xs, molt.atoms.nr);
2662 snew(cg_cm, molt.cgs.nr);
2663 for (int mol = 0; mol < molb.nmol; mol++)
2665 get_cgcm_mol(&molt, &mtop->ffparams, ir->ePBC, &graph, box,
2666 x+at_offset, xs, cg_cm);
2668 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2669 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2671 bonded_cg_distance_mol(&molt, at2cg, bBCheck, bExclRequired, cg_cm,
2672 &bd_mol_2b, &bd_mol_mb);
2674 /* Process the mol data adding the atom index offset */
2675 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2676 at_offset + bd_mol_2b.a1,
2677 at_offset + bd_mol_2b.a2,
2679 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2680 at_offset + bd_mol_mb.a1,
2681 at_offset + bd_mol_mb.a2,
2684 at_offset += molt.atoms.nr;
2688 if (ir->ePBC != epbcNONE)
2695 if (mtop->bIntermolecularInteractions)
2697 if (ncg_mtop(mtop) < mtop->natoms)
2699 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.");
2702 GMX_RELEASE_ASSERT(mtop->intermolecular_ilist.get(), "We should have an ilist when intermolecular interactions are on");
2704 bonded_distance_intermol(*mtop->intermolecular_ilist,
2710 *r_2b = sqrt(bd_2b.r2);
2711 *r_mb = sqrt(bd_mb.r2);
2713 if (*r_2b > 0 || *r_mb > 0)
2715 GMX_LOG(mdlog.info).appendText("Initial maximum inter charge-group distances:");
2718 GMX_LOG(mdlog.info).appendTextFormatted(
2719 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2720 *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2721 bd_2b.a1 + 1, bd_2b.a2 + 1);
2725 GMX_LOG(mdlog.info).appendTextFormatted(
2726 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2727 *r_mb, interaction_function[bd_mb.ftype].longname,
2728 bd_mb.a1 + 1, bd_mb.a2 + 1);