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