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