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