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