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