97587f933a6fc22b8e4419d07b3b37c121ad09f1
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_topology.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2006 - 2014, The GROMACS development team.
5  * Copyright (c) 2015,2016,2017,2018,2019, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *
39  * \brief This file defines functions used by the domdec module
40  * while managing the construction, use and error checking for
41  * topologies local to a DD rank.
42  *
43  * \author Berk Hess <hess@kth.se>
44  * \ingroup module_domdec
45  */
46
47 #include "gmxpre.h"
48
49 #include <cassert>
50 #include <cstdlib>
51 #include <cstring>
52
53 #include <algorithm>
54 #include <string>
55
56 #include "gromacs/compat/make_unique.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/gmxlib/chargegroup.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/forcerec.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdlib/vsite.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/ifunc.h"
74 #include "gromacs/topology/mtop_util.h"
75 #include "gromacs/topology/topsort.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/exceptions.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/gmxassert.h"
80 #include "gromacs/utility/logger.h"
81 #include "gromacs/utility/smalloc.h"
82 #include "gromacs/utility/strconvert.h"
83 #include "gromacs/utility/stringstream.h"
84 #include "gromacs/utility/stringutil.h"
85 #include "gromacs/utility/textwriter.h"
86
87 #include "domdec_constraints.h"
88 #include "domdec_internal.h"
89 #include "domdec_vsite.h"
90 #include "dump.h"
91
92 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
93 #define NITEM_DD_INIT_LOCAL_STATE 5
94
95 struct reverse_ilist_t
96 {
97     std::vector<int> index;              /* Index for each atom into il          */
98     std::vector<int> il;                 /* ftype|type|a0|...|an|ftype|...       */
99     int              numAtomsInMolecule; /* The number of atoms in this molecule */
100 };
101
102 struct MolblockIndices
103 {
104     int  a_start;
105     int  a_end;
106     int  natoms_mol;
107     int  type;
108 };
109
110 /*! \brief Struct for thread local work data for local topology generation */
111 struct thread_work_t
112 {
113     t_idef                    idef;       /**< Partial local topology */
114     std::unique_ptr<VsitePbc> vsitePbc;   /**< vsite PBC structure */
115     int                       nbonded;    /**< The number of bondeds in this struct */
116     t_blocka                  excl;       /**< List of exclusions */
117     int                       excl_count; /**< The total exclusion count for \p excl */
118 };
119
120 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
121 struct gmx_reverse_top_t
122 {
123     //! @cond Doxygen_Suppress
124     //! \brief Do we require all exclusions to be assigned?
125     bool                         bExclRequired = false;
126     //! \brief The maximum number of exclusions one atom can have
127     int                          n_excl_at_max = 0;
128     //! \brief Are there constraints in this revserse top?
129     bool                         bConstr = false;
130     //! \brief Are there settles in this revserse top?
131     bool                         bSettle = false;
132     //! \brief All bonded interactions have to be assigned?
133     bool                         bBCheck = false;
134     //! \brief Are there bondeds/exclusions between charge-groups?
135     bool                         bInterCGInteractions = false;
136     //! \brief Reverse ilist for all moltypes
137     std::vector<reverse_ilist_t> ril_mt;
138     //! \brief The size of ril_mt[?].index summed over all entries
139     int                          ril_mt_tot_size = 0;
140     //! \brief The sorting state of bondeds for free energy
141     int                          ilsort = ilsortUNKNOWN;
142     //! \brief molblock to global atom index for quick lookup of molblocks on atom index
143     std::vector<MolblockIndices> mbi;
144
145     //! \brief Do we have intermolecular interactions?
146     bool             bIntermolecularInteractions = false;
147     //! \brief Intermolecular reverse ilist
148     reverse_ilist_t  ril_intermol;
149
150     /* Work data structures for multi-threading */
151     //! \brief Thread work array for local topology generation
152     std::vector<thread_work_t> th_work;
153     //! @endcond
154 };
155
156 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
157 static int nral_rt(int ftype)
158 {
159     int nral;
160
161     nral = NRAL(ftype);
162     if (interaction_function[ftype].flags & IF_VSITE)
163     {
164         /* With vsites the reverse topology contains
165          * two extra entries for PBC.
166          */
167         nral += 2;
168     }
169
170     return nral;
171 }
172
173 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
174 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
175                                gmx_bool bConstr, gmx_bool bSettle)
176 {
177     return ((((interaction_function[ftype].flags & IF_BOND) != 0u) &&
178              ((interaction_function[ftype].flags & IF_VSITE) == 0u) &&
179              (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0u))) ||
180             (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
181             (bSettle && ftype == F_SETTLE));
182 }
183
184 /*! \brief Help print error output when interactions are missing */
185 static std::string
186 print_missing_interactions_mb(t_commrec *cr,
187                               const gmx_reverse_top_t *rt,
188                               const char *moltypename,
189                               const reverse_ilist_t *ril,
190                               int a_start, int a_end,
191                               int nat_mol, int nmol,
192                               const t_idef *idef)
193 {
194     int                     *assigned;
195     int                      nril_mol = ril->index[nat_mol];
196     snew(assigned, nmol*nril_mol);
197     gmx::StringOutputStream  stream;
198     gmx::TextWriter          log(&stream);
199
200     gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
201     for (int ftype = 0; ftype < F_NRE; ftype++)
202     {
203         if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
204         {
205             int            nral = NRAL(ftype);
206             const t_ilist *il   = &idef->il[ftype];
207             const t_iatom *ia   = il->iatoms;
208             for (int i = 0; i < il->nr; i += 1+nral)
209             {
210                 int a0 = gatindex[ia[1]];
211                 /* Check if this interaction is in
212                  * the currently checked molblock.
213                  */
214                 if (a0 >= a_start && a0 < a_end)
215                 {
216                     int  mol    = (a0 - a_start)/nat_mol;
217                     int  a0_mol = (a0 - a_start) - mol*nat_mol;
218                     int  j_mol  = ril->index[a0_mol];
219                     bool found  = false;
220                     while (j_mol < ril->index[a0_mol+1] && !found)
221                     {
222                         int j       = mol*nril_mol + j_mol;
223                         int ftype_j = ril->il[j_mol];
224                         /* Here we need to check if this interaction has
225                          * not already been assigned, since we could have
226                          * multiply defined interactions.
227                          */
228                         if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
229                             assigned[j] == 0)
230                         {
231                             /* Check the atoms */
232                             found = true;
233                             for (int a = 0; a < nral; a++)
234                             {
235                                 if (gatindex[ia[1+a]] !=
236                                     a_start + mol*nat_mol + ril->il[j_mol+2+a])
237                                 {
238                                     found = false;
239                                 }
240                             }
241                             if (found)
242                             {
243                                 assigned[j] = 1;
244                             }
245                         }
246                         j_mol += 2 + nral_rt(ftype_j);
247                     }
248                     if (!found)
249                     {
250                         gmx_incons("Some interactions seem to be assigned multiple times");
251                     }
252                 }
253                 ia += 1 + nral;
254             }
255         }
256     }
257
258     gmx_sumi(nmol*nril_mol, assigned, cr);
259
260     int nprint = 10;
261     int i      = 0;
262     for (int mol = 0; mol < nmol; mol++)
263     {
264         int j_mol = 0;
265         while (j_mol < nril_mol)
266         {
267             int ftype = ril->il[j_mol];
268             int nral  = NRAL(ftype);
269             int j     = mol*nril_mol + j_mol;
270             if (assigned[j] == 0 &&
271                 !(interaction_function[ftype].flags & IF_VSITE))
272             {
273                 if (DDMASTER(cr->dd))
274                 {
275                     if (i == 0)
276                     {
277                         log.writeLineFormatted("Molecule type '%s'", moltypename);
278                         log.writeLineFormatted(
279                                 "the first %d missing interactions, except for exclusions:", nprint);
280                     }
281                     log.writeStringFormatted("%20s atoms",
282                                              interaction_function[ftype].longname);
283                     int a;
284                     for (a = 0; a < nral; a++)
285                     {
286                         log.writeStringFormatted("%5d", ril->il[j_mol+2+a]+1);
287                     }
288                     while (a < 4)
289                     {
290                         log.writeString("     ");
291                         a++;
292                     }
293                     log.writeString(" global");
294                     for (a = 0; a < nral; a++)
295                     {
296                         log.writeStringFormatted("%6d",
297                                                  a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
298                     }
299                     log.ensureLineBreak();
300                 }
301                 i++;
302                 if (i >= nprint)
303                 {
304                     break;
305                 }
306             }
307             j_mol += 2 + nral_rt(ftype);
308         }
309     }
310
311     sfree(assigned);
312     return stream.toString();
313 }
314
315 /*! \brief Help print error output when interactions are missing */
316 static void print_missing_interactions_atoms(const gmx::MDLogger &mdlog,
317                                              t_commrec           *cr,
318                                              const gmx_mtop_t    *mtop,
319                                              const t_idef        *idef)
320 {
321     const gmx_reverse_top_t *rt = cr->dd->reverse_top;
322
323     /* Print the atoms in the missing interactions per molblock */
324     int a_end = 0;
325     for (const gmx_molblock_t &molb :  mtop->molblock)
326     {
327         const gmx_moltype_t &moltype  = mtop->moltype[molb.type];
328         int                  a_start  = a_end;
329         a_end                        = a_start + molb.nmol*moltype.atoms.nr;
330
331         GMX_LOG(mdlog.warning).appendText(
332                 print_missing_interactions_mb(cr, rt,
333                                               *(moltype.name),
334                                               &rt->ril_mt[molb.type],
335                                               a_start, a_end, moltype.atoms.nr,
336                                               molb.nmol,
337                                               idef));
338     }
339 }
340
341 void dd_print_missing_interactions(const gmx::MDLogger  &mdlog,
342                                    t_commrec            *cr,
343                                    int                   local_count,
344                                    const gmx_mtop_t     *top_global,
345                                    const gmx_localtop_t *top_local,
346                                    const t_state        *state_local)
347 {
348     int             ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
349     int             ftype, nral;
350     gmx_domdec_t   *dd;
351
352     dd = cr->dd;
353
354     GMX_LOG(mdlog.warning).appendText(
355             "Not all bonded interactions have been properly assigned to the domain decomposition cells");
356
357     ndiff_tot = local_count - dd->nbonded_global;
358
359     for (ftype = 0; ftype < F_NRE; ftype++)
360     {
361         nral      = NRAL(ftype);
362         cl[ftype] = top_local->idef.il[ftype].nr/(1+nral);
363     }
364
365     gmx_sumi(F_NRE, cl, cr);
366
367     if (DDMASTER(dd))
368     {
369         GMX_LOG(mdlog.warning).appendText("A list of missing interactions:");
370         rest_global = dd->nbonded_global;
371         rest_local  = local_count;
372         for (ftype = 0; ftype < F_NRE; ftype++)
373         {
374             /* In the reverse and local top all constraints are merged
375              * into F_CONSTR. So in the if statement we skip F_CONSTRNC
376              * and add these constraints when doing F_CONSTR.
377              */
378             if (((interaction_function[ftype].flags & IF_BOND) &&
379                  (dd->reverse_top->bBCheck
380                   || !(interaction_function[ftype].flags & IF_LIMZERO)))
381                 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
382                 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
383             {
384                 n    = gmx_mtop_ftype_count(top_global, ftype);
385                 if (ftype == F_CONSTR)
386                 {
387                     n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
388                 }
389                 ndiff = cl[ftype] - n;
390                 if (ndiff != 0)
391                 {
392                     GMX_LOG(mdlog.warning).appendTextFormatted(
393                             "%20s of %6d missing %6d",
394                             interaction_function[ftype].longname, n, -ndiff);
395                 }
396                 rest_global -= n;
397                 rest_local  -= cl[ftype];
398             }
399         }
400
401         ndiff = rest_local - rest_global;
402         if (ndiff != 0)
403         {
404             GMX_LOG(mdlog.warning).appendTextFormatted(
405                     "%20s of %6d missing %6d", "exclusions",
406                     rest_global, -ndiff);
407         }
408     }
409
410     print_missing_interactions_atoms(mdlog, cr, top_global, &top_local->idef);
411     write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
412                  -1, state_local->x.rvec_array(), state_local->box);
413
414     std::string errorMessage;
415
416     if (ndiff_tot > 0)
417     {
418         errorMessage = "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
419     }
420     else
421     {
422         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     }
424     gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
425 }
426
427 /*! \brief Return global topology molecule information for global atom index \p i_gl */
428 static void global_atomnr_to_moltype_ind(const gmx_reverse_top_t *rt,
429                                          int i_gl,
430                                          int *mb, int *mt, int *mol, int *i_mol)
431 {
432     const MolblockIndices *mbi   = rt->mbi.data();
433     int                    start = 0;
434     int                    end   = rt->mbi.size(); /* exclusive */
435     int                    mid;
436
437     /* binary search for molblock_ind */
438     while (TRUE)
439     {
440         mid = (start+end)/2;
441         if (i_gl >= mbi[mid].a_end)
442         {
443             start = mid+1;
444         }
445         else if (i_gl < mbi[mid].a_start)
446         {
447             end = mid;
448         }
449         else
450         {
451             break;
452         }
453     }
454
455     *mb  = mid;
456     mbi += mid;
457
458     *mt    = mbi->type;
459     *mol   = (i_gl - mbi->a_start) / mbi->natoms_mol;
460     *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
461 }
462
463 /*! \brief Count the exclusions for all atoms in \p cgs */
464 static void count_excls(const t_block *cgs, const t_blocka *excls,
465                         int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
466 {
467     int cg, at0, at1, at, excl, atj;
468
469     *n_excl         = 0;
470     *n_intercg_excl = 0;
471     *n_excl_at_max  = 0;
472     for (cg = 0; cg < cgs->nr; cg++)
473     {
474         at0 = cgs->index[cg];
475         at1 = cgs->index[cg+1];
476         for (at = at0; at < at1; at++)
477         {
478             for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
479             {
480                 atj = excls->a[excl];
481                 if (atj > at)
482                 {
483                     (*n_excl)++;
484                     if (atj < at0 || atj >= at1)
485                     {
486                         (*n_intercg_excl)++;
487                     }
488                 }
489             }
490
491             *n_excl_at_max = std::max(*n_excl_at_max,
492                                       excls->index[at+1] - excls->index[at]);
493         }
494     }
495 }
496
497 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
498 static int low_make_reverse_ilist(const InteractionLists &il_mt,
499                                   const t_atom *atom,
500                                   gmx::ArrayRef < const std::vector < int>> vsitePbc,
501                                   int *count,
502                                   gmx_bool bConstr, gmx_bool bSettle,
503                                   gmx_bool bBCheck,
504                                   gmx::ArrayRef<const int> r_index,
505                                   gmx::ArrayRef<int>       r_il,
506                                   gmx_bool bLinkToAllAtoms,
507                                   gmx_bool bAssign)
508 {
509     int            ftype, j, nlink, link;
510     int            a;
511     int            nint;
512
513     nint = 0;
514     for (ftype = 0; ftype < F_NRE; ftype++)
515     {
516         if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
517             (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
518             (bSettle && ftype == F_SETTLE))
519         {
520             const bool  bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0u);
521             const int   nral   = NRAL(ftype);
522             const auto &il     = il_mt[ftype];
523             for (int i = 0; i < il.size(); i += 1+nral)
524             {
525                 const int* ia = il.iatoms.data() + i;
526                 if (bLinkToAllAtoms)
527                 {
528                     if (bVSite)
529                     {
530                         /* We don't need the virtual sites for the cg-links */
531                         nlink = 0;
532                     }
533                     else
534                     {
535                         nlink = nral;
536                     }
537                 }
538                 else
539                 {
540                     /* Couple to the first atom in the interaction */
541                     nlink = 1;
542                 }
543                 for (link = 0; link < nlink; link++)
544                 {
545                     a = ia[1+link];
546                     if (bAssign)
547                     {
548                         GMX_ASSERT(!r_il.empty(), "with bAssign not allowed to be empty");
549                         GMX_ASSERT(!r_index.empty(), "with bAssign not allowed to be empty");
550                         r_il[r_index[a]+count[a]] =
551                             (ftype == F_CONSTRNC ? F_CONSTR : ftype);
552                         r_il[r_index[a]+count[a]+1] = ia[0];
553                         for (j = 1; j < 1+nral; j++)
554                         {
555                             /* Store the molecular atom number */
556                             r_il[r_index[a]+count[a]+1+j] = ia[j];
557                         }
558                     }
559                     if (interaction_function[ftype].flags & IF_VSITE)
560                     {
561                         if (bAssign)
562                         {
563                             /* Add an entry to iatoms for storing
564                              * which of the constructing atoms are
565                              * vsites again.
566                              */
567                             r_il[r_index[a]+count[a]+2+nral] = 0;
568                             for (j = 2; j < 1+nral; j++)
569                             {
570                                 if (atom[ia[j]].ptype == eptVSite)
571                                 {
572                                     r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
573                                 }
574                             }
575                             /* Store vsite pbc atom in a second extra entry */
576                             r_il[r_index[a]+count[a]+2+nral+1] =
577                                 (!vsitePbc.empty() ? vsitePbc[ftype-F_VSITE2][i/(1+nral)] : -2);
578                         }
579                     }
580                     else
581                     {
582                         /* We do not count vsites since they are always
583                          * uniquely assigned and can be assigned
584                          * to multiple nodes with recursive vsites.
585                          */
586                         if (bBCheck ||
587                             !(interaction_function[ftype].flags & IF_LIMZERO))
588                         {
589                             nint++;
590                         }
591                     }
592                     count[a] += 2 + nral_rt(ftype);
593                 }
594             }
595         }
596     }
597
598     return nint;
599 }
600
601 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
602 static int make_reverse_ilist(const InteractionLists &ilist,
603                               const t_atoms *atoms,
604                               gmx::ArrayRef < const std::vector < int>> vsitePbc,
605                               gmx_bool bConstr, gmx_bool bSettle,
606                               gmx_bool bBCheck,
607                               gmx_bool bLinkToAllAtoms,
608                               reverse_ilist_t *ril_mt)
609 {
610     int nat_mt, *count, i, nint_mt;
611
612     /* Count the interactions */
613     nat_mt = atoms->nr;
614     snew(count, nat_mt);
615     low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
616                            count,
617                            bConstr, bSettle, bBCheck,
618                            gmx::EmptyArrayRef(), gmx::EmptyArrayRef(),
619                            bLinkToAllAtoms, FALSE);
620
621     ril_mt->index.push_back(0);
622     for (i = 0; i < nat_mt; i++)
623     {
624         ril_mt->index.push_back(ril_mt->index[i] + count[i]);
625         count[i] = 0;
626     }
627     ril_mt->il.resize(ril_mt->index[nat_mt]);
628
629     /* Store the interactions */
630     nint_mt =
631         low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
632                                count,
633                                bConstr, bSettle, bBCheck,
634                                ril_mt->index, ril_mt->il,
635                                bLinkToAllAtoms, TRUE);
636
637     sfree(count);
638
639     ril_mt->numAtomsInMolecule = atoms->nr;
640
641     return nint_mt;
642 }
643
644 /*! \brief Generate the reverse topology */
645 static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
646                                           gmx::ArrayRef<const VsitePbc> vsitePbcPerMoltype,
647                                           gmx_bool bConstr, gmx_bool bSettle,
648                                           gmx_bool bBCheck, int *nint)
649 {
650     gmx_reverse_top_t  rt;
651
652     /* Should we include constraints (for SHAKE) in rt? */
653     rt.bConstr = bConstr;
654     rt.bSettle = bSettle;
655     rt.bBCheck = bBCheck;
656
657     rt.bInterCGInteractions = mtop->bIntermolecularInteractions;
658     rt.ril_mt.resize(mtop->moltype.size());
659     rt.ril_mt_tot_size = 0;
660     std::vector<int> nint_mt;
661     for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
662     {
663         const gmx_moltype_t &molt = mtop->moltype[mt];
664         if (molt.cgs.nr > 1)
665         {
666             rt.bInterCGInteractions = true;
667         }
668
669         /* Make the atom to interaction list for this molecule type */
670         gmx::ArrayRef < const std::vector < int>> vsitePbc;
671         if (!vsitePbcPerMoltype.empty())
672         {
673             vsitePbc = gmx::makeConstArrayRef(vsitePbcPerMoltype[mt]);
674         }
675         int numberOfInteractions =
676             make_reverse_ilist(molt.ilist, &molt.atoms, vsitePbc,
677                                rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
678                                &rt.ril_mt[mt]);
679         nint_mt.push_back(numberOfInteractions);
680
681         rt.ril_mt_tot_size += rt.ril_mt[mt].index[molt.atoms.nr];
682     }
683     if (debug)
684     {
685         fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt.ril_mt_tot_size);
686     }
687
688     *nint = 0;
689     for (const gmx_molblock_t &molblock : mtop->molblock)
690     {
691         *nint += molblock.nmol*nint_mt[molblock.type];
692     }
693
694     /* Make an intermolecular reverse top, if necessary */
695     rt.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
696     if (rt.bIntermolecularInteractions)
697     {
698         t_atoms atoms_global;
699
700         atoms_global.nr   = mtop->natoms;
701         atoms_global.atom = nullptr; /* Only used with virtual sites */
702
703         GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
704
705         *nint +=
706             make_reverse_ilist(*mtop->intermolecular_ilist,
707                                &atoms_global,
708                                gmx::EmptyArrayRef(),
709                                rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
710                                &rt.ril_intermol);
711     }
712
713     if (bFE && gmx_mtop_bondeds_free_energy(mtop))
714     {
715         rt.ilsort = ilsortFE_UNSORTED;
716     }
717     else
718     {
719         rt.ilsort = ilsortNO_FE;
720     }
721
722     /* Make a molblock index for fast searching */
723     int i         = 0;
724     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
725     {
726         const gmx_molblock_t &molb           = mtop->molblock[mb];
727         const int             numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
728         MolblockIndices       mbi;
729         mbi.a_start                          = i;
730         i                                   += molb.nmol*numAtomsPerMol;
731         mbi.a_end                            = i;
732         mbi.natoms_mol                       = numAtomsPerMol;
733         mbi.type                             = molb.type;
734         rt.mbi.push_back(mbi);
735     }
736
737     rt.th_work.resize(gmx_omp_nthreads_get(emntDomdec));
738     if (!vsitePbcPerMoltype.empty())
739     {
740         for (thread_work_t &th_work : rt.th_work)
741         {
742             th_work.vsitePbc = gmx::compat::make_unique<VsitePbc>();
743         }
744     }
745
746     return rt;
747 }
748
749 void dd_make_reverse_top(FILE *fplog,
750                          gmx_domdec_t *dd, const gmx_mtop_t *mtop,
751                          const gmx_vsite_t *vsite,
752                          const t_inputrec *ir, gmx_bool bBCheck)
753 {
754     if (fplog)
755     {
756         fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
757     }
758
759     /* If normal and/or settle constraints act only within charge groups,
760      * we can store them in the reverse top and simply assign them to domains.
761      * Otherwise we need to assign them to multiple domains and set up
762      * the parallel version constraint algorithm(s).
763      */
764
765     gmx::ArrayRef<const VsitePbc> vsitePbcPerMoltype;
766     if (vsite)
767     {
768         vsitePbcPerMoltype = gmx::makeConstArrayRef(vsite->vsite_pbc_molt);
769     }
770
771     dd->reverse_top  = new gmx_reverse_top_t;
772     *dd->reverse_top =
773         make_reverse_top(mtop, ir->efep != efepNO, vsitePbcPerMoltype,
774                          !dd->splitConstraints, !dd->splitSettles,
775                          bBCheck, &dd->nbonded_global);
776
777     gmx_reverse_top_t *rt = dd->reverse_top;
778
779     /* With the Verlet scheme, exclusions are handled in the non-bonded
780      * kernels and only exclusions inside the cut-off lead to exclusion
781      * forces. Since each atom pair is treated at most once in the non-bonded
782      * kernels, it doesn't matter if the exclusions for the same atom pair
783      * appear multiple times in the exclusion list. In contrast, the, old,
784      * group cut-off scheme loops over a list of exclusions, so there each
785      * excluded pair should appear exactly once.
786      */
787     rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
788                          inputrecExclForces(ir));
789
790     int nexcl          = 0;
791     dd->n_intercg_excl = 0;
792     rt->n_excl_at_max  = 0;
793     for (const gmx_molblock_t &molb : mtop->molblock)
794     {
795         int                  n_excl_mol, n_excl_icg, n_excl_at_max;
796
797         const gmx_moltype_t &molt = mtop->moltype[molb.type];
798         count_excls(&molt.cgs, &molt.excls,
799                     &n_excl_mol, &n_excl_icg, &n_excl_at_max);
800         nexcl              += molb.nmol*n_excl_mol;
801         dd->n_intercg_excl += molb.nmol*n_excl_icg;
802         rt->n_excl_at_max   = std::max(rt->n_excl_at_max, n_excl_at_max);
803     }
804     if (rt->bExclRequired)
805     {
806         dd->nbonded_global += nexcl;
807         if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
808         {
809             fprintf(fplog, "There are %d inter charge-group exclusions,\n"
810                     "will use an extra communication step for exclusion forces for %s\n",
811                     dd->n_intercg_excl, eel_names[ir->coulombtype]);
812         }
813     }
814
815     if (vsite && vsite->n_intercg_vsite > 0)
816     {
817         if (fplog)
818         {
819             fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
820                     "will an extra communication step for selected coordinates and forces\n",
821                     vsite->n_intercg_vsite);
822         }
823         init_domdec_vsites(dd, vsite->n_intercg_vsite);
824     }
825
826     if (dd->splitConstraints || dd->splitSettles)
827     {
828         init_domdec_constraints(dd, mtop);
829     }
830     if (fplog)
831     {
832         fprintf(fplog, "\n");
833     }
834 }
835
836 /*! \brief Store a vsite interaction at the end of \p il
837  *
838  * This routine is very similar to add_ifunc, but vsites interactions
839  * have more work to do than other kinds of interactions, and the
840  * complex way nral (and thus vector contents) depends on ftype
841  * confuses static analysis tools unless we fuse the vsite
842  * atom-indexing organization code with the ifunc-adding code, so that
843  * they can see that nral is the same value. */
844 static inline void
845 add_ifunc_for_vsites(t_iatom *tiatoms, const gmx_ga2la_t &ga2la,
846                      int nral, gmx_bool bHomeA,
847                      int a, int a_gl, int a_mol,
848                      const t_iatom *iatoms,
849                      t_ilist *il)
850 {
851     t_iatom *liatoms;
852
853     if (il->nr+1+nral > il->nalloc)
854     {
855         il->nalloc = over_alloc_large(il->nr+1+nral);
856         srenew(il->iatoms, il->nalloc);
857     }
858     liatoms = il->iatoms + il->nr;
859     il->nr += 1 + nral;
860
861     /* Copy the type */
862     tiatoms[0] = iatoms[0];
863
864     if (bHomeA)
865     {
866         /* We know the local index of the first atom */
867         tiatoms[1] = a;
868     }
869     else
870     {
871         /* Convert later in make_local_vsites */
872         tiatoms[1] = -a_gl - 1;
873     }
874
875     for (int k = 2; k < 1+nral; k++)
876     {
877         int ak_gl = a_gl + iatoms[k] - a_mol;
878         if (const int *homeIndex = ga2la.findHome(ak_gl))
879         {
880             tiatoms[k] = *homeIndex;
881         }
882         else
883         {
884             /* Copy the global index, convert later in make_local_vsites */
885             tiatoms[k] = -(ak_gl + 1);
886         }
887         // Note that ga2la_get_home always sets the third parameter if
888         // it returns TRUE
889     }
890     for (int k = 0; k < 1+nral; k++)
891     {
892         liatoms[k] = tiatoms[k];
893     }
894 }
895
896 /*! \brief Store a bonded interaction at the end of \p il */
897 static inline void add_ifunc(int nral, const t_iatom *tiatoms, t_ilist *il)
898 {
899     t_iatom *liatoms;
900     int      k;
901
902     if (il->nr+1+nral > il->nalloc)
903     {
904         il->nalloc = over_alloc_large(il->nr+1+nral);
905         srenew(il->iatoms, il->nalloc);
906     }
907     liatoms = il->iatoms + il->nr;
908     for (k = 0; k <= nral; k++)
909     {
910         liatoms[k] = tiatoms[k];
911     }
912     il->nr += 1 + nral;
913 }
914
915 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
916 static void add_posres(int mol, int a_mol, int numAtomsInMolecule,
917                        const gmx_molblock_t *molb,
918                        t_iatom *iatoms, const t_iparams *ip_in,
919                        t_idef *idef)
920 {
921     int        n, a_molb;
922     t_iparams *ip;
923
924     /* This position restraint has not been added yet,
925      * so it's index is the current number of position restraints.
926      */
927     n = idef->il[F_POSRES].nr/2;
928     if (n+1 > idef->iparams_posres_nalloc)
929     {
930         idef->iparams_posres_nalloc = over_alloc_dd(n+1);
931         srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
932     }
933     ip = &idef->iparams_posres[n];
934     /* Copy the force constants */
935     *ip = ip_in[iatoms[0]];
936
937     /* Get the position restraint coordinates from the molblock */
938     a_molb = mol*numAtomsInMolecule + a_mol;
939     GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
940     ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
941     ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
942     ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
943     if (!molb->posres_xB.empty())
944     {
945         ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
946         ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
947         ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
948     }
949     else
950     {
951         ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
952         ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
953         ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
954     }
955     /* Set the parameter index for idef->iparams_posre */
956     iatoms[0] = n;
957 }
958
959 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
960 static void add_fbposres(int mol, int a_mol, int numAtomsInMolecule,
961                          const gmx_molblock_t *molb,
962                          t_iatom *iatoms, const t_iparams *ip_in,
963                          t_idef *idef)
964 {
965     int        n, a_molb;
966     t_iparams *ip;
967
968     /* This flat-bottom position restraint has not been added yet,
969      * so it's index is the current number of position restraints.
970      */
971     n = idef->il[F_FBPOSRES].nr/2;
972     if (n+1 > idef->iparams_fbposres_nalloc)
973     {
974         idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
975         srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
976     }
977     ip = &idef->iparams_fbposres[n];
978     /* Copy the force constants */
979     *ip = ip_in[iatoms[0]];
980
981     /* Get the position restraint coordinats from the molblock */
982     a_molb = mol*numAtomsInMolecule + a_mol;
983     GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
984     /* Take reference positions from A position of normal posres */
985     ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
986     ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
987     ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
988
989     /* Note: no B-type for flat-bottom posres */
990
991     /* Set the parameter index for idef->iparams_posre */
992     iatoms[0] = n;
993 }
994
995 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
996 static void add_vsite(const gmx_ga2la_t &ga2la,
997                       gmx::ArrayRef<const int> index,
998                       gmx::ArrayRef<const int> rtil,
999                       int ftype, int nral,
1000                       gmx_bool bHomeA, int a, int a_gl, int a_mol,
1001                       const t_iatom *iatoms,
1002                       t_idef *idef,
1003                       VsitePbc *vsitePbc)
1004 {
1005     int     k, pbc_a_mol;
1006     t_iatom tiatoms[1+MAXATOMLIST];
1007     int     j, ftype_r, nral_r;
1008
1009     /* Add this interaction to the local topology */
1010     add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1011
1012     if (vsitePbc)
1013     {
1014         std::vector<int> &vsitePbcFtype = (*vsitePbc)[ftype - c_ftypeVsiteStart];
1015         const int         vsi           = idef->il[ftype].nr/(1+nral) - 1;
1016         if (static_cast<size_t>(vsi) >= vsitePbcFtype.size())
1017         {
1018             vsitePbcFtype.resize(vsi + 1);
1019         }
1020         if (bHomeA)
1021         {
1022             pbc_a_mol = iatoms[1+nral+1];
1023             if (pbc_a_mol < 0)
1024             {
1025                 /* The pbc flag is one of the following two options:
1026                  * -2: vsite and all constructing atoms are within the same cg, no pbc
1027                  * -1: vsite and its first constructing atom are in the same cg, do pbc
1028                  */
1029                 vsitePbcFtype[vsi] = pbc_a_mol;
1030             }
1031             else
1032             {
1033                 /* Set the pbc atom for this vsite so we can make its pbc
1034                  * identical to the rest of the atoms in its charge group.
1035                  * Since the order of the atoms does not change within a charge
1036                  * group, we do not need the global to local atom index.
1037                  */
1038                 vsitePbcFtype[vsi] = a + pbc_a_mol - iatoms[1];
1039             }
1040         }
1041         else
1042         {
1043             /* This vsite is non-home (required for recursion),
1044              * and therefore there is no charge group to match pbc with.
1045              * But we always turn on full_pbc to assure that higher order
1046              * recursion works correctly.
1047              */
1048             vsitePbcFtype[vsi] = -1;
1049         }
1050     }
1051
1052     if (iatoms[1+nral])
1053     {
1054         /* Check for recursion */
1055         for (k = 2; k < 1+nral; k++)
1056         {
1057             if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1058             {
1059                 /* This construction atoms is a vsite and not a home atom */
1060                 if (gmx_debug_at)
1061                 {
1062                     fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1063                 }
1064                 /* Find the vsite construction */
1065
1066                 /* Check all interactions assigned to this atom */
1067                 j = index[iatoms[k]];
1068                 while (j < index[iatoms[k]+1])
1069                 {
1070                     ftype_r = rtil[j++];
1071                     nral_r  = NRAL(ftype_r);
1072                     if (interaction_function[ftype_r].flags & IF_VSITE)
1073                     {
1074                         /* Add this vsite (recursion) */
1075                         add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1076                                   FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1077                                   rtil.data() + j,
1078                                   idef, vsitePbc);
1079                         j += 1 + nral_r + 2;
1080                     }
1081                     else
1082                     {
1083                         j += 1 + nral_r;
1084                     }
1085                 }
1086             }
1087         }
1088     }
1089 }
1090
1091 /*! \brief Build the index that maps each local atom to its local atom group */
1092 static void makeLocalAtomGroupsFromAtoms(gmx_domdec_t *dd)
1093 {
1094     const gmx::RangePartitioning &atomGrouping = dd->atomGrouping();
1095
1096     dd->localAtomGroupFromAtom.clear();
1097
1098     for (size_t g = 0; g < dd->globalAtomGroupIndices.size(); g++)
1099     {
1100         for (int gmx_unused a : atomGrouping.block(g))
1101         {
1102             dd->localAtomGroupFromAtom.push_back(g);
1103         }
1104     }
1105 }
1106
1107 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1108 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1109 {
1110     rvec dx;
1111
1112     if (pbc_null)
1113     {
1114         pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1115     }
1116     else
1117     {
1118         rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1119     }
1120
1121     return norm2(dx);
1122 }
1123
1124 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1125 static void combine_blocka(t_blocka                           *dest,
1126                            gmx::ArrayRef<const thread_work_t>  src)
1127 {
1128     int ni = src.back().excl.nr;
1129     int na = 0;
1130     for (const thread_work_t &th_work : src)
1131     {
1132         na += th_work.excl.nra;
1133     }
1134     if (ni + 1 > dest->nalloc_index)
1135     {
1136         dest->nalloc_index = over_alloc_large(ni+1);
1137         srenew(dest->index, dest->nalloc_index);
1138     }
1139     if (dest->nra + na > dest->nalloc_a)
1140     {
1141         dest->nalloc_a = over_alloc_large(dest->nra+na);
1142         srenew(dest->a, dest->nalloc_a);
1143     }
1144     for (gmx::index s = 1; s < src.size(); s++)
1145     {
1146         for (int i = dest->nr + 1; i < src[s].excl.nr + 1; i++)
1147         {
1148             dest->index[i] = dest->nra + src[s].excl.index[i];
1149         }
1150         for (int i = 0; i < src[s].excl.nra; i++)
1151         {
1152             dest->a[dest->nra+i] = src[s].excl.a[i];
1153         }
1154         dest->nr   = src[s].excl.nr;
1155         dest->nra += src[s].excl.nra;
1156     }
1157 }
1158
1159 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1160  * virtual sites need special attention, as pbc info differs per vsite.
1161  */
1162 static void combine_idef(t_idef                             *dest,
1163                          gmx::ArrayRef<const thread_work_t>  src,
1164                          gmx_vsite_t                        *vsite)
1165 {
1166     int ftype;
1167
1168     for (ftype = 0; ftype < F_NRE; ftype++)
1169     {
1170         int n = 0;
1171         for (gmx::index s = 1; s < src.size(); s++)
1172         {
1173             n += src[s].idef.il[ftype].nr;
1174         }
1175         if (n > 0)
1176         {
1177             t_ilist *ild;
1178
1179             ild = &dest->il[ftype];
1180
1181             if (ild->nr + n > ild->nalloc)
1182             {
1183                 ild->nalloc = over_alloc_large(ild->nr+n);
1184                 srenew(ild->iatoms, ild->nalloc);
1185             }
1186
1187             const bool vpbc  =
1188                 (((interaction_function[ftype].flags & IF_VSITE) != 0u) &&
1189                  vsite->vsite_pbc_loc);
1190             const int  nral1 = 1 + NRAL(ftype);
1191             const int  ftv   = ftype - c_ftypeVsiteStart;
1192
1193             for (gmx::index s = 1; s < src.size(); s++)
1194             {
1195                 const t_ilist &ils = src[s].idef.il[ftype];
1196
1197                 for (int i = 0; i < ils.nr; i++)
1198                 {
1199                     ild->iatoms[ild->nr + i] = ils.iatoms[i];
1200                 }
1201                 if (vpbc)
1202                 {
1203                     const std::vector<int> &pbcSrc  = (*src[s].vsitePbc)[ftv];
1204                     std::vector<int>       &pbcDest = (*vsite->vsite_pbc_loc)[ftv];
1205                     pbcDest.resize((ild->nr + ils.nr)/nral1);
1206                     for (int i = 0; i < ils.nr; i += nral1)
1207                     {
1208                         pbcDest[(ild->nr + i)/nral1] = pbcSrc[i/nral1];
1209                     }
1210                 }
1211
1212                 ild->nr += ils.nr;
1213             }
1214
1215             /* Position restraints need an additional treatment */
1216             if (ftype == F_POSRES || ftype == F_FBPOSRES)
1217             {
1218                 int          nposres       = dest->il[ftype].nr/2;
1219                 // TODO: Simplify this code using std::vector
1220                 t_iparams * &iparams_dest  = (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1221                 int         &posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc : dest->iparams_fbposres_nalloc);
1222                 if (nposres > posres_nalloc)
1223                 {
1224                     posres_nalloc = over_alloc_large(nposres);
1225                     srenew(iparams_dest, posres_nalloc);
1226                 }
1227
1228                 /* Set nposres to the number of original position restraints in dest */
1229                 for (gmx::index s = 1; s < src.size(); s++)
1230                 {
1231                     nposres -= src[s].idef.il[ftype].nr/2;
1232                 }
1233
1234                 for (gmx::index s = 1; s < src.size(); s++)
1235                 {
1236                     const t_iparams *iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1237
1238                     for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1239                     {
1240                         /* Correct the index into iparams_posres */
1241                         dest->il[ftype].iatoms[nposres*2] = nposres;
1242                         /* Copy the position restraint force parameters */
1243                         iparams_dest[nposres]             = iparams_src[i];
1244                         nposres++;
1245                     }
1246                 }
1247             }
1248         }
1249     }
1250 }
1251
1252 /*! \brief Check and when available assign bonded interactions for local atom i
1253  */
1254 static inline void
1255 check_assign_interactions_atom(int i, int i_gl,
1256                                int mol, int i_mol,
1257                                int numAtomsInMolecule,
1258                                gmx::ArrayRef<const int> index,
1259                                gmx::ArrayRef<const int> rtil,
1260                                gmx_bool bInterMolInteractions,
1261                                int ind_start, int ind_end,
1262                                const gmx_domdec_t *dd,
1263                                const gmx_domdec_zones_t *zones,
1264                                const gmx_molblock_t *molb,
1265                                gmx_bool bRCheckMB, const ivec rcheck, gmx_bool bRCheck2B,
1266                                real rc2,
1267                                int *la2lc,
1268                                t_pbc *pbc_null,
1269                                rvec *cg_cm,
1270                                const t_iparams *ip_in,
1271                                t_idef *idef,
1272                                VsitePbc *vsitePbc,
1273                                int iz,
1274                                gmx_bool bBCheck,
1275                                int *nbonded_local)
1276 {
1277     int j;
1278
1279     j = ind_start;
1280     while (j < ind_end)
1281     {
1282         t_iatom   tiatoms[1 + MAXATOMLIST];
1283
1284         const int ftype  = rtil[j++];
1285         auto      iatoms = gmx::constArrayRefFromArray(rtil.data() + j, rtil.size() - j);
1286         const int nral   = NRAL(ftype);
1287         if (interaction_function[ftype].flags & IF_VSITE)
1288         {
1289             assert(!bInterMolInteractions);
1290             /* The vsite construction goes where the vsite itself is */
1291             if (iz == 0)
1292             {
1293                 add_vsite(*dd->ga2la, index, rtil, ftype, nral,
1294                           TRUE, i, i_gl, i_mol,
1295                           iatoms.data(), idef, vsitePbc);
1296             }
1297             j += 1 + nral + 2;
1298         }
1299         else
1300         {
1301             gmx_bool bUse;
1302
1303             /* Copy the type */
1304             tiatoms[0] = iatoms[0];
1305
1306             if (nral == 1)
1307             {
1308                 assert(!bInterMolInteractions);
1309                 /* Assign single-body interactions to the home zone */
1310                 if (iz == 0)
1311                 {
1312                     bUse       = TRUE;
1313                     tiatoms[1] = i;
1314                     if (ftype == F_POSRES)
1315                     {
1316                         add_posres(mol, i_mol, numAtomsInMolecule,
1317                                    molb, tiatoms, ip_in, idef);
1318                     }
1319                     else if (ftype == F_FBPOSRES)
1320                     {
1321                         add_fbposres(mol, i_mol, numAtomsInMolecule,
1322                                      molb, tiatoms, ip_in, idef);
1323                     }
1324                 }
1325                 else
1326                 {
1327                     bUse = FALSE;
1328                 }
1329             }
1330             else if (nral == 2)
1331             {
1332                 /* This is a two-body interaction, we can assign
1333                  * analogous to the non-bonded assignments.
1334                  */
1335                 int k_gl;
1336
1337                 if (!bInterMolInteractions)
1338                 {
1339                     /* Get the global index using the offset in the molecule */
1340                     k_gl = i_gl + iatoms[2] - i_mol;
1341                 }
1342                 else
1343                 {
1344                     k_gl = iatoms[2];
1345                 }
1346                 if (const auto *entry = dd->ga2la->find(k_gl))
1347                 {
1348                     int kz = entry->cell;
1349                     if (kz >= zones->n)
1350                     {
1351                         kz -= zones->n;
1352                     }
1353                     /* Check zone interaction assignments */
1354                     bUse = ((iz < zones->nizone &&
1355                              iz <= kz &&
1356                              kz >= zones->izone[iz].j0 &&
1357                              kz <  zones->izone[iz].j1) ||
1358                             (kz < zones->nizone &&
1359                                   iz > kz &&
1360                              iz >= zones->izone[kz].j0 &&
1361                              iz <  zones->izone[kz].j1));
1362                     if (bUse)
1363                     {
1364                         GMX_ASSERT(ftype != F_CONSTR || (iz == 0 && kz == 0),
1365                                    "Constraint assigned here should only involve home atoms");
1366
1367                         tiatoms[1] = i;
1368                         tiatoms[2] = entry->la;
1369                         /* If necessary check the cgcm distance */
1370                         if (bRCheck2B &&
1371                             dd_dist2(pbc_null, cg_cm, la2lc,
1372                                      tiatoms[1], tiatoms[2]) >= rc2)
1373                         {
1374                             bUse = FALSE;
1375                         }
1376                     }
1377                 }
1378                 else
1379                 {
1380                     bUse = false;
1381                 }
1382             }
1383             else
1384             {
1385                 /* Assign this multi-body bonded interaction to
1386                  * the local node if we have all the atoms involved
1387                  * (local or communicated) and the minimum zone shift
1388                  * in each dimension is zero, for dimensions
1389                  * with 2 DD cells an extra check may be necessary.
1390                  */
1391                 ivec k_zero, k_plus;
1392                 int  k;
1393
1394                 bUse = TRUE;
1395                 clear_ivec(k_zero);
1396                 clear_ivec(k_plus);
1397                 for (k = 1; k <= nral && bUse; k++)
1398                 {
1399                     int k_gl;
1400                     if (!bInterMolInteractions)
1401                     {
1402                         /* Get the global index using the offset in the molecule */
1403                         k_gl = i_gl + iatoms[k] - i_mol;
1404                     }
1405                     else
1406                     {
1407                         k_gl = iatoms[k];
1408                     }
1409                     const auto *entry = dd->ga2la->find(k_gl);
1410                     if (entry == nullptr || entry->cell >= zones->n)
1411                     {
1412                         /* We do not have this atom of this interaction
1413                          * locally, or it comes from more than one cell
1414                          * away.
1415                          */
1416                         bUse = FALSE;
1417                     }
1418                     else
1419                     {
1420                         int d;
1421
1422                         tiatoms[k] = entry->la;
1423                         for (d = 0; d < DIM; d++)
1424                         {
1425                             if (zones->shift[entry->cell][d] == 0)
1426                             {
1427                                 k_zero[d] = k;
1428                             }
1429                             else
1430                             {
1431                                 k_plus[d] = k;
1432                             }
1433                         }
1434                     }
1435                 }
1436                 bUse = (bUse &&
1437                         (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1438                 if (bRCheckMB)
1439                 {
1440                     int d;
1441
1442                     for (d = 0; (d < DIM && bUse); d++)
1443                     {
1444                         /* Check if the cg_cm distance falls within
1445                          * the cut-off to avoid possible multiple
1446                          * assignments of bonded interactions.
1447                          */
1448                         if (rcheck[d] &&
1449                             k_plus[d] &&
1450                             dd_dist2(pbc_null, cg_cm, la2lc,
1451                                      tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1452                         {
1453                             bUse = FALSE;
1454                         }
1455                     }
1456                 }
1457             }
1458             if (bUse)
1459             {
1460                 /* Add this interaction to the local topology */
1461                 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1462                 /* Sum so we can check in global_stat
1463                  * if we have everything.
1464                  */
1465                 if (bBCheck ||
1466                     !(interaction_function[ftype].flags & IF_LIMZERO))
1467                 {
1468                     (*nbonded_local)++;
1469                 }
1470             }
1471             j += 1 + nral;
1472         }
1473     }
1474 }
1475
1476 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1477  *
1478  * With thread parallelizing each thread acts on a different atom range:
1479  * at_start to at_end.
1480  */
1481 static int make_bondeds_zone(gmx_domdec_t *dd,
1482                              const gmx_domdec_zones_t *zones,
1483                              const std::vector<gmx_molblock_t> &molb,
1484                              gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1485                              real rc2,
1486                              int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1487                              const t_iparams *ip_in,
1488                              t_idef *idef,
1489                              VsitePbc *vsitePbc,
1490                              int izone,
1491                              gmx::RangePartitioning::Block atomRange)
1492 {
1493     int                mb, mt, mol, i_mol;
1494     gmx_bool           bBCheck;
1495     gmx_reverse_top_t *rt;
1496     int                nbonded_local;
1497
1498     rt = dd->reverse_top;
1499
1500     bBCheck = rt->bBCheck;
1501
1502     nbonded_local = 0;
1503
1504     for (int i : atomRange)
1505     {
1506         /* Get the global atom number */
1507         const int i_gl = dd->globalAtomIndices[i];
1508         global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1509         /* Check all intramolecular interactions assigned to this atom */
1510         gmx::ArrayRef<const int>       index = rt->ril_mt[mt].index;
1511         gmx::ArrayRef<const t_iatom>   rtil  = rt->ril_mt[mt].il;
1512
1513         check_assign_interactions_atom(i, i_gl, mol, i_mol,
1514                                        rt->ril_mt[mt].numAtomsInMolecule,
1515                                        index, rtil, FALSE,
1516                                        index[i_mol], index[i_mol+1],
1517                                        dd, zones,
1518                                        &molb[mb],
1519                                        bRCheckMB, rcheck, bRCheck2B, rc2,
1520                                        la2lc,
1521                                        pbc_null,
1522                                        cg_cm,
1523                                        ip_in,
1524                                        idef, vsitePbc,
1525                                        izone,
1526                                        bBCheck,
1527                                        &nbonded_local);
1528
1529
1530         if (rt->bIntermolecularInteractions)
1531         {
1532             /* Check all intermolecular interactions assigned to this atom */
1533             index = rt->ril_intermol.index;
1534             rtil  = rt->ril_intermol.il;
1535
1536             check_assign_interactions_atom(i, i_gl, mol, i_mol,
1537                                            rt->ril_mt[mt].numAtomsInMolecule,
1538                                            index, rtil, TRUE,
1539                                            index[i_gl], index[i_gl + 1],
1540                                            dd, zones,
1541                                            &molb[mb],
1542                                            bRCheckMB, rcheck, bRCheck2B, rc2,
1543                                            la2lc,
1544                                            pbc_null,
1545                                            cg_cm,
1546                                            ip_in,
1547                                            idef, vsitePbc,
1548                                            izone,
1549                                            bBCheck,
1550                                            &nbonded_local);
1551         }
1552     }
1553
1554     return nbonded_local;
1555 }
1556
1557 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1558 static void set_no_exclusions_zone(const gmx_domdec_t       *dd,
1559                                    const gmx_domdec_zones_t *zones,
1560                                    int                       iz,
1561                                    t_blocka                 *lexcls)
1562 {
1563     const auto zone = dd->atomGrouping().subRange(zones->cg_range[iz],
1564                                                   zones->cg_range[iz + 1]);
1565
1566     for (int a : zone)
1567     {
1568         lexcls->index[a + 1] = lexcls->nra;
1569     }
1570 }
1571
1572 /*! \brief Set the exclusion data for i-zone \p iz
1573  *
1574  * This is a legacy version for the group scheme of the same routine below.
1575  * Here charge groups and distance checks to ensure unique exclusions
1576  * are supported.
1577  */
1578 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1579                                    const std::vector<gmx_moltype_t> &moltype,
1580                                    gmx_bool bRCheck, real rc2,
1581                                    int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1582                                    const int *cginfo,
1583                                    t_blocka *lexcls,
1584                                    int iz,
1585                                    int cg_start, int cg_end)
1586 {
1587     int                n_excl_at_max;
1588     int                mb, mt, mol;
1589     const t_blocka    *excls;
1590
1591     const gmx_ga2la_t &ga2la  = *dd->ga2la;
1592
1593     const auto         jRange =
1594         dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1595                                     zones->izone[iz].jcg1);
1596
1597     n_excl_at_max = dd->reverse_top->n_excl_at_max;
1598
1599     /* We set the end index, but note that we might not start at zero here */
1600     lexcls->nr = dd->atomGrouping().subRange(0, cg_end).size();
1601
1602     int n     = lexcls->nra;
1603     int count = 0;
1604     for (int cg = cg_start; cg < cg_end; cg++)
1605     {
1606         if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1607         {
1608             lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1609             srenew(lexcls->a, lexcls->nalloc_a);
1610         }
1611         const auto atomGroup = dd->atomGrouping().block(cg);
1612         if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1613             !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1614         {
1615             /* Copy the exclusions from the global top */
1616             for (int la : atomGroup)
1617             {
1618                 lexcls->index[la] = n;
1619                 int a_gl          = dd->globalAtomIndices[la];
1620                 int a_mol;
1621                 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1622                 excls = &moltype[mt].excls;
1623                 for (int j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1624                 {
1625                     int aj_mol = excls->a[j];
1626                     /* This computation of jla is only correct intra-cg */
1627                     int jla = la + aj_mol - a_mol;
1628                     if (atomGroup.inRange(jla))
1629                     {
1630                         /* This is an intra-cg exclusion. We can skip
1631                          *  the global indexing and distance checking.
1632                          */
1633                         /* Intra-cg exclusions are only required
1634                          * for the home zone.
1635                          */
1636                         if (iz == 0)
1637                         {
1638                             lexcls->a[n++] = jla;
1639                             /* Check to avoid double counts */
1640                             if (jla > la)
1641                             {
1642                                 count++;
1643                             }
1644                         }
1645                     }
1646                     else
1647                     {
1648                         /* This is a inter-cg exclusion */
1649                         /* Since exclusions are pair interactions,
1650                          * just like non-bonded interactions,
1651                          * they can be assigned properly up
1652                          * to the DD cutoff (not cutoff_min as
1653                          * for the other bonded interactions).
1654                          */
1655                         if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1656                         {
1657                             if (iz == 0 && jEntry->cell == 0)
1658                             {
1659                                 lexcls->a[n++] = jEntry->la;
1660                                 /* Check to avoid double counts */
1661                                 if (jEntry->la > la)
1662                                 {
1663                                     count++;
1664                                 }
1665                             }
1666                             else if (jRange.inRange(jEntry->la) &&
1667                                      (!bRCheck ||
1668                                       dd_dist2(pbc_null, cg_cm, la2lc, la, jEntry->la) < rc2))
1669                             {
1670                                 /* jla > la, since jRange.begin() > la */
1671                                 lexcls->a[n++] = jEntry->la;
1672                                 count++;
1673                             }
1674                         }
1675                     }
1676                 }
1677             }
1678         }
1679         else
1680         {
1681             /* There are no inter-cg excls and this cg is self-excluded.
1682              * These exclusions are only required for zone 0,
1683              * since other zones do not see themselves.
1684              */
1685             if (iz == 0)
1686             {
1687                 for (int la : atomGroup)
1688                 {
1689                     lexcls->index[la] = n;
1690                     for (int j : atomGroup)
1691                     {
1692                         lexcls->a[n++] = j;
1693                     }
1694                 }
1695                 count += (atomGroup.size()*(atomGroup.size() - 1))/2;
1696             }
1697             else
1698             {
1699                 /* We don't need exclusions for this cg */
1700                 for (int la : atomGroup)
1701                 {
1702                     lexcls->index[la] = n;
1703                 }
1704             }
1705         }
1706     }
1707
1708     lexcls->index[lexcls->nr] = n;
1709     lexcls->nra               = n;
1710
1711     return count;
1712 }
1713
1714 /*! \brief Set the exclusion data for i-zone \p iz */
1715 static void make_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1716                                  const std::vector<gmx_moltype_t> &moltype,
1717                                  const int *cginfo, t_blocka *lexcls, int iz,
1718                                  int at_start, int at_end,
1719                                  const gmx::ArrayRef<const int> intermolecularExclusionGroup)
1720 {
1721     int                n_excl_at_max, n, at;
1722
1723     const gmx_ga2la_t &ga2la  = *dd->ga2la;
1724
1725     const auto         jRange =
1726         dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1727                                     zones->izone[iz].jcg1);
1728
1729     n_excl_at_max = dd->reverse_top->n_excl_at_max;
1730
1731     /* We set the end index, but note that we might not start at zero here */
1732     lexcls->nr = at_end;
1733
1734     n = lexcls->nra;
1735     for (at = at_start; at < at_end; at++)
1736     {
1737         if (n + 1000 > lexcls->nalloc_a)
1738         {
1739             lexcls->nalloc_a = over_alloc_large(n + 1000);
1740             srenew(lexcls->a, lexcls->nalloc_a);
1741         }
1742
1743         if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1744         {
1745             int             a_gl, mb, mt, mol, a_mol, j;
1746             const t_blocka *excls;
1747
1748             if (n + n_excl_at_max > lexcls->nalloc_a)
1749             {
1750                 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1751                 srenew(lexcls->a, lexcls->nalloc_a);
1752             }
1753
1754             /* Copy the exclusions from the global top */
1755             lexcls->index[at] = n;
1756             a_gl              = dd->globalAtomIndices[at];
1757             global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1758                                          &mb, &mt, &mol, &a_mol);
1759             excls = &moltype[mt].excls;
1760             for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1761             {
1762                 const int aj_mol = excls->a[j];
1763
1764                 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1765                 {
1766                     /* This check is not necessary, but it can reduce
1767                      * the number of exclusions in the list, which in turn
1768                      * can speed up the pair list construction a bit.
1769                      */
1770                     if (jRange.inRange(jEntry->la))
1771                     {
1772                         lexcls->a[n++] = jEntry->la;
1773                     }
1774                 }
1775             }
1776         }
1777         else
1778         {
1779             /* We don't need exclusions for this atom */
1780             lexcls->index[at] = n;
1781         }
1782
1783         bool isExcludedAtom = !intermolecularExclusionGroup.empty() &&
1784             std::find(intermolecularExclusionGroup.begin(),
1785                       intermolecularExclusionGroup.end(),
1786                       dd->globalAtomIndices[at]) !=
1787             intermolecularExclusionGroup.end();
1788
1789         if (isExcludedAtom)
1790         {
1791             if (n + intermolecularExclusionGroup.size() > lexcls->nalloc_a)
1792             {
1793                 lexcls->nalloc_a =
1794                     over_alloc_large(n + intermolecularExclusionGroup.size());
1795                 srenew(lexcls->a, lexcls->nalloc_a);
1796             }
1797             for (int qmAtomGlobalIndex : intermolecularExclusionGroup)
1798             {
1799                 if (const auto *entry = dd->ga2la->find(qmAtomGlobalIndex))
1800                 {
1801                     lexcls->a[n++] = entry->la;
1802                 }
1803             }
1804         }
1805     }
1806
1807     lexcls->index[lexcls->nr] = n;
1808     lexcls->nra               = n;
1809 }
1810
1811
1812 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1813 static void check_alloc_index(t_blocka *ba, int nindex_max)
1814 {
1815     if (nindex_max+1 > ba->nalloc_index)
1816     {
1817         ba->nalloc_index = over_alloc_dd(nindex_max+1);
1818         srenew(ba->index, ba->nalloc_index);
1819     }
1820 }
1821
1822 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1823 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1824                                    t_blocka *lexcls)
1825 {
1826     int nr = dd->atomGrouping().subRange(0, zones->izone[zones->nizone - 1].cg1).size();
1827
1828     check_alloc_index(lexcls, nr);
1829
1830     for (size_t thread = 1; thread < dd->reverse_top->th_work.size(); thread++)
1831     {
1832         check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1833     }
1834 }
1835
1836 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1837 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1838                                     t_blocka *lexcls)
1839 {
1840     const auto nonhomeIzonesAtomRange =
1841         dd->atomGrouping().subRange(zones->izone[0].cg1,
1842                                     zones->izone[zones->nizone - 1].cg1);
1843
1844     if (dd->n_intercg_excl == 0)
1845     {
1846         /* There are no exclusions involving non-home charge groups,
1847          * but we need to set the indices for neighborsearching.
1848          */
1849         for (int la : nonhomeIzonesAtomRange)
1850         {
1851             lexcls->index[la] = lexcls->nra;
1852         }
1853
1854         /* nr is only used to loop over the exclusions for Ewald and RF,
1855          * so we can set it to the number of home atoms for efficiency.
1856          */
1857         lexcls->nr = nonhomeIzonesAtomRange.begin();
1858     }
1859     else
1860     {
1861         lexcls->nr = nonhomeIzonesAtomRange.end();
1862     }
1863 }
1864
1865 /*! \brief Clear a t_idef struct */
1866 static void clear_idef(t_idef *idef)
1867 {
1868     int  ftype;
1869
1870     /* Clear the counts */
1871     for (ftype = 0; ftype < F_NRE; ftype++)
1872     {
1873         idef->il[ftype].nr = 0;
1874     }
1875 }
1876
1877 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1878 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1879                                     gmx_domdec_zones_t *zones,
1880                                     const gmx_mtop_t *mtop,
1881                                     const int *cginfo,
1882                                     gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1883                                     real rc,
1884                                     int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1885                                     t_idef *idef, gmx_vsite_t *vsite,
1886                                     t_blocka *lexcls, int *excl_count)
1887 {
1888     int                nzone_bondeds, nzone_excl;
1889     int                izone, cg0, cg1;
1890     real               rc2;
1891     int                nbonded_local;
1892     int                thread;
1893     gmx_reverse_top_t *rt;
1894
1895     if (dd->reverse_top->bInterCGInteractions)
1896     {
1897         nzone_bondeds = zones->n;
1898     }
1899     else
1900     {
1901         /* Only single charge group (or atom) molecules, so interactions don't
1902          * cross zone boundaries and we only need to assign in the home zone.
1903          */
1904         nzone_bondeds = 1;
1905     }
1906
1907     if (dd->n_intercg_excl > 0)
1908     {
1909         /* We only use exclusions from i-zones to i- and j-zones */
1910         nzone_excl = zones->nizone;
1911     }
1912     else
1913     {
1914         /* There are no inter-cg exclusions and only zone 0 sees itself */
1915         nzone_excl = 1;
1916     }
1917
1918     check_exclusions_alloc(dd, zones, lexcls);
1919
1920     rt = dd->reverse_top;
1921
1922     rc2 = rc*rc;
1923
1924     /* Clear the counts */
1925     clear_idef(idef);
1926     nbonded_local = 0;
1927
1928     lexcls->nr    = 0;
1929     lexcls->nra   = 0;
1930     *excl_count   = 0;
1931
1932     for (izone = 0; izone < nzone_bondeds; izone++)
1933     {
1934         cg0 = zones->cg_range[izone];
1935         cg1 = zones->cg_range[izone + 1];
1936
1937         const int numThreads = rt->th_work.size();
1938 #pragma omp parallel for num_threads(numThreads) schedule(static)
1939         for (thread = 0; thread < numThreads; thread++)
1940         {
1941             try
1942             {
1943                 int       cg0t, cg1t;
1944                 t_idef   *idef_t;
1945                 t_blocka *excl_t;
1946
1947                 cg0t = cg0 + ((cg1 - cg0)* thread   )/numThreads;
1948                 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/numThreads;
1949
1950                 if (thread == 0)
1951                 {
1952                     idef_t = idef;
1953                 }
1954                 else
1955                 {
1956                     idef_t = &rt->th_work[thread].idef;
1957                     clear_idef(idef_t);
1958                 }
1959
1960                 VsitePbc *vsitePbc = nullptr;
1961                 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1962                 {
1963                     if (thread == 0)
1964                     {
1965                         vsitePbc = vsite->vsite_pbc_loc.get();
1966                     }
1967                     else
1968                     {
1969                         vsitePbc = rt->th_work[thread].vsitePbc.get();
1970                     }
1971                 }
1972
1973                 rt->th_work[thread].nbonded =
1974                     make_bondeds_zone(dd, zones,
1975                                       mtop->molblock,
1976                                       bRCheckMB, rcheck, bRCheck2B, rc2,
1977                                       la2lc, pbc_null, cg_cm, idef->iparams,
1978                                       idef_t, vsitePbc,
1979                                       izone,
1980                                       dd->atomGrouping().subRange(cg0t, cg1t));
1981
1982                 if (izone < nzone_excl)
1983                 {
1984                     if (thread == 0)
1985                     {
1986                         excl_t = lexcls;
1987                     }
1988                     else
1989                     {
1990                         excl_t      = &rt->th_work[thread].excl;
1991                         excl_t->nr  = 0;
1992                         excl_t->nra = 0;
1993                     }
1994
1995                     if (dd->atomGrouping().allBlocksHaveSizeOne() &&
1996                         !rt->bExclRequired)
1997                     {
1998                         /* No charge groups and no distance check required */
1999                         make_exclusions_zone(dd, zones, mtop->moltype, cginfo,
2000                                              excl_t, izone, cg0t,
2001                                              cg1t,
2002                                              mtop->intermolecularExclusionGroup);
2003                     }
2004                     else
2005                     {
2006                         rt->th_work[thread].excl_count =
2007                             make_exclusions_zone_cg(dd, zones,
2008                                                     mtop->moltype, bRCheck2B, rc2,
2009                                                     la2lc, pbc_null, cg_cm, cginfo,
2010                                                     excl_t,
2011                                                     izone,
2012                                                     cg0t, cg1t);
2013                     }
2014                 }
2015             }
2016             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2017         }
2018
2019         if (rt->th_work.size() > 1)
2020         {
2021             combine_idef(idef, rt->th_work, vsite);
2022         }
2023
2024         for (const thread_work_t &th_work : rt->th_work)
2025         {
2026             nbonded_local += th_work.nbonded;
2027         }
2028
2029         if (izone < nzone_excl)
2030         {
2031             if (rt->th_work.size() > 1)
2032             {
2033                 combine_blocka(lexcls, rt->th_work);
2034             }
2035
2036             for (const thread_work_t &th_work : rt->th_work)
2037             {
2038                 *excl_count += th_work.excl_count;
2039             }
2040         }
2041     }
2042
2043     /* Some zones might not have exclusions, but some code still needs to
2044      * loop over the index, so we set the indices here.
2045      */
2046     for (izone = nzone_excl; izone < zones->nizone; izone++)
2047     {
2048         set_no_exclusions_zone(dd, zones, izone, lexcls);
2049     }
2050
2051     finish_local_exclusions(dd, zones, lexcls);
2052     if (debug)
2053     {
2054         fprintf(debug, "We have %d exclusions, check count %d\n",
2055                 lexcls->nra, *excl_count);
2056     }
2057
2058     return nbonded_local;
2059 }
2060
2061 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
2062 {
2063     lcgs->nr    = dd->globalAtomGroupIndices.size();
2064     lcgs->index = dd->atomGrouping_.rawIndex().data();
2065 }
2066
2067 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
2068                        int npbcdim, matrix box,
2069                        rvec cellsize_min, const ivec npulse,
2070                        t_forcerec *fr,
2071                        rvec *cgcm_or_x,
2072                        gmx_vsite_t *vsite,
2073                        const gmx_mtop_t &mtop, gmx_localtop_t *ltop)
2074 {
2075     gmx_bool bRCheckMB, bRCheck2B;
2076     real     rc = -1;
2077     ivec     rcheck;
2078     int      d, nexcl;
2079     t_pbc    pbc, *pbc_null = nullptr;
2080
2081     if (debug)
2082     {
2083         fprintf(debug, "Making local topology\n");
2084     }
2085
2086     dd_make_local_cgs(dd, &ltop->cgs);
2087
2088     bRCheckMB   = FALSE;
2089     bRCheck2B   = FALSE;
2090
2091     if (dd->reverse_top->bInterCGInteractions)
2092     {
2093         /* We need to check to which cell bondeds should be assigned */
2094         rc = dd_cutoff_twobody(dd);
2095         if (debug)
2096         {
2097             fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2098         }
2099
2100         /* Should we check cg_cm distances when assigning bonded interactions? */
2101         for (d = 0; d < DIM; d++)
2102         {
2103             rcheck[d] = FALSE;
2104             /* Only need to check for dimensions where the part of the box
2105              * that is not communicated is smaller than the cut-off.
2106              */
2107             if (d < npbcdim && dd->nc[d] > 1 &&
2108                 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2109             {
2110                 if (dd->nc[d] == 2)
2111                 {
2112                     rcheck[d] = TRUE;
2113                     bRCheckMB = TRUE;
2114                 }
2115                 /* Check for interactions between two atoms,
2116                  * where we can allow interactions up to the cut-off,
2117                  * instead of up to the smallest cell dimension.
2118                  */
2119                 bRCheck2B = TRUE;
2120             }
2121             if (debug)
2122             {
2123                 fprintf(debug,
2124                         "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n",
2125                         d, cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
2126             }
2127         }
2128         if (bRCheckMB || bRCheck2B)
2129         {
2130             makeLocalAtomGroupsFromAtoms(dd);
2131             if (fr->bMolPBC)
2132             {
2133                 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
2134             }
2135             else
2136             {
2137                 pbc_null = nullptr;
2138             }
2139         }
2140     }
2141
2142     dd->nbonded_local =
2143         make_local_bondeds_excls(dd, zones, &mtop, fr->cginfo,
2144                                  bRCheckMB, rcheck, bRCheck2B, rc,
2145                                  dd->localAtomGroupFromAtom.data(),
2146                                  pbc_null, cgcm_or_x,
2147                                  &ltop->idef, vsite,
2148                                  &ltop->excls, &nexcl);
2149
2150     /* The ilist is not sorted yet,
2151      * we can only do this when we have the charge arrays.
2152      */
2153     ltop->idef.ilsort = ilsortUNKNOWN;
2154
2155     if (dd->reverse_top->bExclRequired)
2156     {
2157         dd->nbonded_local += nexcl;
2158     }
2159
2160     ltop->atomtypes  = mtop.atomtypes;
2161 }
2162
2163 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
2164                        gmx_localtop_t *ltop)
2165 {
2166     if (dd->reverse_top->ilsort == ilsortNO_FE)
2167     {
2168         ltop->idef.ilsort = ilsortNO_FE;
2169     }
2170     else
2171     {
2172         gmx_sort_ilist_fe(&ltop->idef, mdatoms->chargeA, mdatoms->chargeB);
2173     }
2174 }
2175
2176 void dd_init_local_top(const gmx_mtop_t &top_global,
2177                        gmx_localtop_t   *top)
2178 {
2179     /* TODO: Get rid of the const casts below, e.g. by using a reference */
2180     top->idef.ntypes     = top_global.ffparams.numTypes();
2181     top->idef.atnr       = top_global.ffparams.atnr;
2182     top->idef.functype   = const_cast<t_functype *>(top_global.ffparams.functype.data());
2183     top->idef.iparams    = const_cast<t_iparams *>(top_global.ffparams.iparams.data());
2184     top->idef.fudgeQQ    = top_global.ffparams.fudgeQQ;
2185     top->idef.cmap_grid  = new gmx_cmap_t;
2186     *top->idef.cmap_grid = top_global.ffparams.cmap_grid;
2187
2188     top->idef.ilsort        = ilsortUNKNOWN;
2189     top->useInDomainDecomp_ = true;
2190 }
2191
2192 void dd_init_local_state(gmx_domdec_t *dd,
2193                          const t_state *state_global, t_state *state_local)
2194 {
2195     int buf[NITEM_DD_INIT_LOCAL_STATE];
2196
2197     if (DDMASTER(dd))
2198     {
2199         buf[0] = state_global->flags;
2200         buf[1] = state_global->ngtc;
2201         buf[2] = state_global->nnhpres;
2202         buf[3] = state_global->nhchainlength;
2203         buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
2204     }
2205     dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2206
2207     init_gtc_state(state_local, buf[1], buf[2], buf[3]);
2208     init_dfhist_state(state_local, buf[4]);
2209     state_local->flags = buf[0];
2210 }
2211
2212 /*! \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 */
2213 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2214 {
2215     int      k;
2216     gmx_bool bFound;
2217
2218     bFound = FALSE;
2219     for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2220     {
2221         GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2222         if (link->a[k] == cg_gl_j)
2223         {
2224             bFound = TRUE;
2225         }
2226     }
2227     if (!bFound)
2228     {
2229         GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2230                            "Inconsistent allocation of link");
2231         /* Add this charge group link */
2232         if (link->index[cg_gl+1]+1 > link->nalloc_a)
2233         {
2234             link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2235             srenew(link->a, link->nalloc_a);
2236         }
2237         link->a[link->index[cg_gl+1]] = cg_gl_j;
2238         link->index[cg_gl+1]++;
2239     }
2240 }
2241
2242 /*! \brief Return a vector of the charge group index for all atoms */
2243 static std::vector<int> make_at2cg(const t_block &cgs)
2244 {
2245     std::vector<int> at2cg(cgs.index[cgs.nr]);
2246     for (int cg = 0; cg < cgs.nr; cg++)
2247     {
2248         for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2249         {
2250             at2cg[a] = cg;
2251         }
2252     }
2253
2254     return at2cg;
2255 }
2256
2257 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
2258                                   cginfo_mb_t *cginfo_mb)
2259 {
2260     gmx_bool            bExclRequired;
2261     t_blocka           *link;
2262     cginfo_mb_t        *cgi_mb;
2263
2264     /* For each charge group make a list of other charge groups
2265      * in the system that a linked to it via bonded interactions
2266      * which are also stored in reverse_top.
2267      */
2268
2269     bExclRequired = dd->reverse_top->bExclRequired;
2270
2271     reverse_ilist_t ril_intermol;
2272     if (mtop->bIntermolecularInteractions)
2273     {
2274         if (ncg_mtop(mtop) < mtop->natoms)
2275         {
2276             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.");
2277         }
2278
2279         t_atoms atoms;
2280
2281         atoms.nr   = mtop->natoms;
2282         atoms.atom = nullptr;
2283
2284         GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
2285
2286         make_reverse_ilist(*mtop->intermolecular_ilist,
2287                            &atoms,
2288                            gmx::EmptyArrayRef(),
2289                            FALSE, FALSE, FALSE, TRUE, &ril_intermol);
2290     }
2291
2292     snew(link, 1);
2293     snew(link->index, ncg_mtop(mtop)+1);
2294     link->nalloc_a = 0;
2295     link->a        = nullptr;
2296
2297     link->index[0] = 0;
2298     int cg_offset  = 0;
2299     int ncgi       = 0;
2300     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
2301     {
2302         const gmx_molblock_t &molb = mtop->molblock[mb];
2303         if (molb.nmol == 0)
2304         {
2305             continue;
2306         }
2307         const gmx_moltype_t &molt  = mtop->moltype[molb.type];
2308         const t_block       &cgs   = molt.cgs;
2309         const t_blocka      &excls = molt.excls;
2310         std::vector<int>     a2c   = make_at2cg(cgs);
2311         /* Make a reverse ilist in which the interactions are linked
2312          * to all atoms, not only the first atom as in gmx_reverse_top.
2313          * The constraints are discarded here.
2314          */
2315         reverse_ilist_t ril;
2316         make_reverse_ilist(molt.ilist, &molt.atoms, gmx::EmptyArrayRef(),
2317                            FALSE, FALSE, FALSE, TRUE, &ril);
2318
2319         cgi_mb = &cginfo_mb[mb];
2320
2321         int mol;
2322         for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb.nmol : 1); mol++)
2323         {
2324             for (int cg = 0; cg < cgs.nr; cg++)
2325             {
2326                 int cg_gl            = cg_offset + cg;
2327                 link->index[cg_gl+1] = link->index[cg_gl];
2328                 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2329                 {
2330                     int i = ril.index[a];
2331                     while (i < ril.index[a+1])
2332                     {
2333                         int ftype = ril.il[i++];
2334                         int nral  = NRAL(ftype);
2335                         /* Skip the ifunc index */
2336                         i++;
2337                         for (int j = 0; j < nral; j++)
2338                         {
2339                             int aj = ril.il[i + j];
2340                             if (a2c[aj] != cg)
2341                             {
2342                                 check_link(link, cg_gl, cg_offset+a2c[aj]);
2343                             }
2344                         }
2345                         i += nral_rt(ftype);
2346                     }
2347                     if (bExclRequired)
2348                     {
2349                         /* Exclusions always go both ways */
2350                         for (int j = excls.index[a]; j < excls.index[a + 1]; j++)
2351                         {
2352                             int aj = excls.a[j];
2353                             if (a2c[aj] != cg)
2354                             {
2355                                 check_link(link, cg_gl, cg_offset+a2c[aj]);
2356                             }
2357                         }
2358                     }
2359
2360                     if (mtop->bIntermolecularInteractions)
2361                     {
2362                         int i = ril_intermol.index[a];
2363                         while (i < ril_intermol.index[a+1])
2364                         {
2365                             int ftype = ril_intermol.il[i++];
2366                             int nral  = NRAL(ftype);
2367                             /* Skip the ifunc index */
2368                             i++;
2369                             for (int j = 0; j < nral; j++)
2370                             {
2371                                 /* Here we assume we have no charge groups;
2372                                  * this has been checked above.
2373                                  */
2374                                 int aj = ril_intermol.il[i + j];
2375                                 check_link(link, cg_gl, aj);
2376                             }
2377                             i += nral_rt(ftype);
2378                         }
2379                     }
2380                 }
2381                 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2382                 {
2383                     SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2384                     ncgi++;
2385                 }
2386             }
2387
2388             cg_offset += cgs.nr;
2389         }
2390         int nlink_mol = link->index[cg_offset] - link->index[cg_offset - cgs.nr];
2391
2392         if (debug)
2393         {
2394             fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt.name, cgs.nr, nlink_mol);
2395         }
2396
2397         if (molb.nmol > mol)
2398         {
2399             /* Copy the data for the rest of the molecules in this block */
2400             link->nalloc_a += (molb.nmol - mol)*nlink_mol;
2401             srenew(link->a, link->nalloc_a);
2402             for (; mol < molb.nmol; mol++)
2403             {
2404                 for (int cg = 0; cg < cgs.nr; cg++)
2405                 {
2406                     int cg_gl              = cg_offset + cg;
2407                     link->index[cg_gl + 1] =
2408                         link->index[cg_gl + 1 - cgs.nr] + nlink_mol;
2409                     for (int j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2410                     {
2411                         link->a[j] = link->a[j - nlink_mol] + cgs.nr;
2412                     }
2413                     if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2414                         cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2415                     {
2416                         SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2417                         ncgi++;
2418                     }
2419                 }
2420                 cg_offset += cgs.nr;
2421             }
2422         }
2423     }
2424
2425     if (debug)
2426     {
2427         fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2428     }
2429
2430     return link;
2431 }
2432
2433 typedef struct {
2434     real r2;
2435     int  ftype;
2436     int  a1;
2437     int  a2;
2438 } bonded_distance_t;
2439
2440 /*! \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 */
2441 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2442                                        bonded_distance_t *bd)
2443 {
2444     if (r2 > bd->r2)
2445     {
2446         bd->r2    = r2;
2447         bd->ftype = ftype;
2448         bd->a1    = a1;
2449         bd->a2    = a2;
2450     }
2451 }
2452
2453 /*! \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 */
2454 static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
2455                                    const std::vector<int> &at2cg,
2456                                    gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2457                                    bonded_distance_t *bd_2b,
2458                                    bonded_distance_t *bd_mb)
2459 {
2460     for (int ftype = 0; ftype < F_NRE; ftype++)
2461     {
2462         if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2463         {
2464             const auto    &il   = molt->ilist[ftype];
2465             int            nral = NRAL(ftype);
2466             if (nral > 1)
2467             {
2468                 for (int i = 0; i < il.size(); i += 1+nral)
2469                 {
2470                     for (int ai = 0; ai < nral; ai++)
2471                     {
2472                         int cgi = at2cg[il.iatoms[i+1+ai]];
2473                         for (int aj = ai + 1; aj < nral; aj++)
2474                         {
2475                             int cgj = at2cg[il.iatoms[i+1+aj]];
2476                             if (cgi != cgj)
2477                             {
2478                                 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2479
2480                                 update_max_bonded_distance(rij2, ftype,
2481                                                            il.iatoms[i+1+ai],
2482                                                            il.iatoms[i+1+aj],
2483                                                            (nral == 2) ? bd_2b : bd_mb);
2484                             }
2485                         }
2486                     }
2487                 }
2488             }
2489         }
2490     }
2491     if (bExcl)
2492     {
2493         const t_blocka *excls = &molt->excls;
2494         for (int ai = 0; ai < excls->nr; ai++)
2495         {
2496             int cgi = at2cg[ai];
2497             for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2498             {
2499                 int cgj = at2cg[excls->a[j]];
2500                 if (cgi != cgj)
2501                 {
2502                     real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2503
2504                     /* There is no function type for exclusions, use -1 */
2505                     update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2506                 }
2507             }
2508         }
2509     }
2510 }
2511
2512 /*! \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 */
2513 static void bonded_distance_intermol(const InteractionLists &ilists_intermol,
2514                                      gmx_bool bBCheck,
2515                                      const rvec *x, int ePBC, const matrix box,
2516                                      bonded_distance_t *bd_2b,
2517                                      bonded_distance_t *bd_mb)
2518 {
2519     t_pbc pbc;
2520
2521     set_pbc(&pbc, ePBC, box);
2522
2523     for (int ftype = 0; ftype < F_NRE; ftype++)
2524     {
2525         if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2526         {
2527             const auto    &il   = ilists_intermol[ftype];
2528             int            nral = NRAL(ftype);
2529
2530             /* No nral>1 check here, since intermol interactions always
2531              * have nral>=2 (and the code is also correct for nral=1).
2532              */
2533             for (int i = 0; i < il.size(); i += 1+nral)
2534             {
2535                 for (int ai = 0; ai < nral; ai++)
2536                 {
2537                     int atom_i = il.iatoms[i + 1 + ai];
2538
2539                     for (int aj = ai + 1; aj < nral; aj++)
2540                     {
2541                         rvec dx;
2542                         real rij2;
2543
2544                         int  atom_j = il.iatoms[i + 1 + aj];
2545
2546                         pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2547
2548                         rij2 = norm2(dx);
2549
2550                         update_max_bonded_distance(rij2, ftype,
2551                                                    atom_i, atom_j,
2552                                                    (nral == 2) ? bd_2b : bd_mb);
2553                     }
2554                 }
2555             }
2556         }
2557     }
2558 }
2559
2560 //! Returns whether \p molt has at least one virtual site
2561 static bool moltypeHasVsite(const gmx_moltype_t &molt)
2562 {
2563     bool hasVsite = false;
2564     for (int i = 0; i < F_NRE; i++)
2565     {
2566         if ((interaction_function[i].flags & IF_VSITE) &&
2567             molt.ilist[i].size() > 0)
2568         {
2569             hasVsite = true;
2570         }
2571     }
2572
2573     return hasVsite;
2574 }
2575
2576 //! Compute charge group centers of mass for molecule \p molt
2577 static void get_cgcm_mol(const gmx_moltype_t *molt,
2578                          const gmx_ffparams_t *ffparams,
2579                          int ePBC, t_graph *graph, const matrix box,
2580                          const rvec *x, rvec *xs, rvec *cg_cm)
2581 {
2582     int n, i;
2583
2584     if (ePBC != epbcNONE)
2585     {
2586         mk_mshift(nullptr, graph, ePBC, box, x);
2587
2588         shift_x(graph, box, x, xs);
2589         /* By doing an extra mk_mshift the molecules that are broken
2590          * because they were e.g. imported from another software
2591          * will be made whole again. Such are the healing powers
2592          * of GROMACS.
2593          */
2594         mk_mshift(nullptr, graph, ePBC, box, xs);
2595     }
2596     else
2597     {
2598         /* We copy the coordinates so the original coordinates remain
2599          * unchanged, just to be 100% sure that we do not affect
2600          * binary reproducibility of simulations.
2601          */
2602         n = molt->cgs.index[molt->cgs.nr];
2603         for (i = 0; i < n; i++)
2604         {
2605             copy_rvec(x[i], xs[i]);
2606         }
2607     }
2608
2609     if (moltypeHasVsite(*molt))
2610     {
2611         /* Convert to old, deprecated format */
2612         t_ilist ilist[F_NRE];
2613         for (int ftype = 0; ftype < F_NRE; ftype++)
2614         {
2615             if (interaction_function[ftype].flags & IF_VSITE)
2616             {
2617                 ilist[ftype].nr     = molt->ilist[ftype].size();
2618                 ilist[ftype].iatoms = const_cast<int *>(molt->ilist[ftype].iatoms.data());
2619             }
2620         }
2621
2622         construct_vsites(nullptr, xs, 0.0, nullptr,
2623                          ffparams->iparams.data(), ilist,
2624                          epbcNONE, TRUE, nullptr, nullptr);
2625     }
2626
2627     calc_cgcm(nullptr, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2628 }
2629
2630 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
2631                            const gmx_mtop_t *mtop,
2632                            const t_inputrec *ir,
2633                            const rvec *x, const matrix box,
2634                            gmx_bool bBCheck,
2635                            real *r_2b, real *r_mb)
2636 {
2637     gmx_bool           bExclRequired;
2638     int                at_offset;
2639     t_graph            graph;
2640     rvec              *xs, *cg_cm;
2641     bonded_distance_t  bd_2b = { 0, -1, -1, -1 };
2642     bonded_distance_t  bd_mb = { 0, -1, -1, -1 };
2643
2644     bExclRequired = inputrecExclForces(ir);
2645
2646     *r_2b     = 0;
2647     *r_mb     = 0;
2648     at_offset = 0;
2649     for (const gmx_molblock_t &molb : mtop->molblock)
2650     {
2651         const gmx_moltype_t &molt = mtop->moltype[molb.type];
2652         if (molt.cgs.nr == 1 || molb.nmol == 0)
2653         {
2654             at_offset += molb.nmol*molt.atoms.nr;
2655         }
2656         else
2657         {
2658             if (ir->ePBC != epbcNONE)
2659             {
2660                 mk_graph_moltype(molt, &graph);
2661             }
2662
2663             std::vector<int> at2cg = make_at2cg(molt.cgs);
2664             snew(xs, molt.atoms.nr);
2665             snew(cg_cm, molt.cgs.nr);
2666             for (int mol = 0; mol < molb.nmol; mol++)
2667             {
2668                 get_cgcm_mol(&molt, &mtop->ffparams, ir->ePBC, &graph, box,
2669                              x+at_offset, xs, cg_cm);
2670
2671                 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2672                 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2673
2674                 bonded_cg_distance_mol(&molt, at2cg, bBCheck, bExclRequired, cg_cm,
2675                                        &bd_mol_2b, &bd_mol_mb);
2676
2677                 /* Process the mol data adding the atom index offset */
2678                 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2679                                            at_offset + bd_mol_2b.a1,
2680                                            at_offset + bd_mol_2b.a2,
2681                                            &bd_2b);
2682                 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2683                                            at_offset + bd_mol_mb.a1,
2684                                            at_offset + bd_mol_mb.a2,
2685                                            &bd_mb);
2686
2687                 at_offset += molt.atoms.nr;
2688             }
2689             sfree(cg_cm);
2690             sfree(xs);
2691             if (ir->ePBC != epbcNONE)
2692             {
2693                 done_graph(&graph);
2694             }
2695         }
2696     }
2697
2698     if (mtop->bIntermolecularInteractions)
2699     {
2700         if (ncg_mtop(mtop) < mtop->natoms)
2701         {
2702             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.");
2703         }
2704
2705         GMX_RELEASE_ASSERT(mtop->intermolecular_ilist, "We should have an ilist when intermolecular interactions are on");
2706
2707         bonded_distance_intermol(*mtop->intermolecular_ilist,
2708                                  bBCheck,
2709                                  x, ir->ePBC, box,
2710                                  &bd_2b, &bd_mb);
2711     }
2712
2713     *r_2b = sqrt(bd_2b.r2);
2714     *r_mb = sqrt(bd_mb.r2);
2715
2716     if (*r_2b > 0 || *r_mb > 0)
2717     {
2718         GMX_LOG(mdlog.info).appendText("Initial maximum distances in bonded interactions:");
2719         if (*r_2b > 0)
2720         {
2721             GMX_LOG(mdlog.info).appendTextFormatted(
2722                     "    two-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2723                     *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2724                     bd_2b.a1 + 1, bd_2b.a2 + 1);
2725         }
2726         if (*r_mb > 0)
2727         {
2728             GMX_LOG(mdlog.info).appendTextFormatted(
2729                     "  multi-body bonded interactions: %5.3f nm, %s, atoms %d %d",
2730                     *r_mb, interaction_function[bd_mb.ftype].longname,
2731                     bd_mb.a1 + 1, bd_mb.a2 + 1);
2732         }
2733     }
2734 }