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