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