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