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