2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2011,2012,2013,2014,2015, 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.
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.
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.
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.
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.
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.
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.
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_domdec
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_network.h"
55 #include "gromacs/legacyheaders/chargegroup.h"
56 #include "gromacs/legacyheaders/force.h"
57 #include "gromacs/legacyheaders/gmx_ga2la.h"
58 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
59 #include "gromacs/legacyheaders/names.h"
60 #include "gromacs/legacyheaders/network.h"
61 #include "gromacs/legacyheaders/typedefs.h"
62 #include "gromacs/legacyheaders/vsite.h"
63 #include "gromacs/legacyheaders/types/commrec.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/pbcutil/mshift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/topology/mtop_util.h"
68 #include "gromacs/topology/topsort.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/smalloc.h"
74 #include "domdec_constraints.h"
75 #include "domdec_internal.h"
76 #include "domdec_vsite.h"
78 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
79 #define NITEM_DD_INIT_LOCAL_STATE 5
82 int *index; /* Index for each atom into il */
83 int *il; /* ftype|type|a0|...|an|ftype|... */
93 /*! \brief Struct for thread local work data for local topology generation */
95 t_idef idef; /**< Partial local topology */
96 int **vsite_pbc; /**< vsite PBC structure */
97 int *vsite_pbc_nalloc; /**< Allocation sizes for vsite_pbc */
98 int nbonded; /**< The number of bondeds in this struct */
99 t_blocka excl; /**< List of exclusions */
100 int excl_count; /**< The total exclusion count for \p excl */
103 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
104 typedef struct gmx_reverse_top {
105 //! @cond Doxygen_Suppress
106 gmx_bool bExclRequired; /**< Do we require all exclusions to be assigned? */
107 int n_excl_at_max; /**< The maximum number of exclusions one atom can have */
108 gmx_bool bConstr; /**< Are there constraints in this revserse top? */
109 gmx_bool bSettle; /**< Are there settles in this revserse top? */
110 gmx_bool bBCheck; /**< All bonded interactions have to be assigned? */
111 gmx_bool bMultiCGmols; /**< Are the multi charge-group molecules? */
112 reverse_ilist_t *ril_mt; /**< Reverse ilist for all moltypes */
113 int ril_mt_tot_size; /**< The size of ril_mt[?].index summed over all entries */
114 int ilsort; /**< The sorting state of bondeds for free energy */
115 molblock_ind_t *mbi; /**< molblock to global atom index for quick lookup of molblocks on atom index */
116 int nmolblock; /**< The number of molblocks, size of \p mbi */
118 /* Work data structures for multi-threading */
119 int nthread; /**< The number of threads to be used */
120 thread_work_t *th_work; /**< Thread work array for local topology generation */
122 /* Pointers only used for an error message */
123 gmx_mtop_t *err_top_global; /**< Pointer to the global top, only used for error reporting */
124 gmx_localtop_t *err_top_local; /**< Pointer to the local top, only used for error reporting */
128 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
129 static int nral_rt(int ftype)
134 if (interaction_function[ftype].flags & IF_VSITE)
136 /* With vsites the reverse topology contains
137 * two extra entries for PBC.
145 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
146 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
147 gmx_bool bConstr, gmx_bool bSettle)
149 return (((interaction_function[ftype].flags & IF_BOND) &&
150 !(interaction_function[ftype].flags & IF_VSITE) &&
151 (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))) ||
152 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
153 (bSettle && ftype == F_SETTLE));
156 /*! \brief Print a header on error messages */
157 static void print_error_header(FILE *fplog, char *moltypename, int nprint)
159 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
160 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
162 "the first %d missing interactions, except for exclusions:\n",
165 "the first %d missing interactions, except for exclusions:\n",
169 /*! \brief Help print error output when interactions are missing */
170 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
171 gmx_reverse_top_t *rt,
173 reverse_ilist_t *ril,
174 int a_start, int a_end,
175 int nat_mol, int nmol,
178 int nril_mol, *assigned, *gatindex;
179 int ftype, ftype_j, nral, i, j_mol, j, a0, a0_mol, mol, a;
185 nril_mol = ril->index[nat_mol];
186 snew(assigned, nmol*nril_mol);
188 gatindex = cr->dd->gatindex;
189 for (ftype = 0; ftype < F_NRE; ftype++)
191 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
194 il = &idef->il[ftype];
196 for (i = 0; i < il->nr; i += 1+nral)
198 a0 = gatindex[ia[1]];
199 /* Check if this interaction is in
200 * the currently checked molblock.
202 if (a0 >= a_start && a0 < a_end)
204 mol = (a0 - a_start)/nat_mol;
205 a0_mol = (a0 - a_start) - mol*nat_mol;
206 j_mol = ril->index[a0_mol];
208 while (j_mol < ril->index[a0_mol+1] && !bFound)
210 j = mol*nril_mol + j_mol;
211 ftype_j = ril->il[j_mol];
212 /* Here we need to check if this interaction has
213 * not already been assigned, since we could have
214 * multiply defined interactions.
216 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
219 /* Check the atoms */
221 for (a = 0; a < nral; a++)
223 if (gatindex[ia[1+a]] !=
224 a_start + mol*nat_mol + ril->il[j_mol+2+a])
234 j_mol += 2 + nral_rt(ftype_j);
238 gmx_incons("Some interactions seem to be assigned multiple times");
246 gmx_sumi(nmol*nril_mol, assigned, cr);
250 for (mol = 0; mol < nmol; mol++)
253 while (j_mol < nril_mol)
255 ftype = ril->il[j_mol];
257 j = mol*nril_mol + j_mol;
258 if (assigned[j] == 0 &&
259 !(interaction_function[ftype].flags & IF_VSITE))
261 if (DDMASTER(cr->dd))
265 print_error_header(fplog, moltypename, nprint);
267 fprintf(fplog, "%20s atoms",
268 interaction_function[ftype].longname);
269 fprintf(stderr, "%20s atoms",
270 interaction_function[ftype].longname);
271 for (a = 0; a < nral; a++)
273 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
274 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
279 fprintf(stderr, " ");
282 fprintf(fplog, " global");
283 fprintf(stderr, " global");
284 for (a = 0; a < nral; a++)
286 fprintf(fplog, "%6d",
287 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
288 fprintf(stderr, "%6d",
289 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
291 fprintf(fplog, "\n");
292 fprintf(stderr, "\n");
300 j_mol += 2 + nral_rt(ftype);
307 /*! \brief Help print error output when interactions are missing */
308 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
309 gmx_mtop_t *mtop, t_idef *idef)
311 int mb, a_start, a_end;
312 gmx_molblock_t *molb;
313 gmx_reverse_top_t *rt;
315 rt = cr->dd->reverse_top;
317 /* Print the atoms in the missing interactions per molblock */
319 for (mb = 0; mb < mtop->nmolblock; mb++)
321 molb = &mtop->molblock[mb];
323 a_end = a_start + molb->nmol*molb->natoms_mol;
325 print_missing_interactions_mb(fplog, cr, rt,
326 *(mtop->moltype[molb->type].name),
327 &rt->ril_mt[molb->type],
328 a_start, a_end, molb->natoms_mol,
334 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr, int local_count, gmx_mtop_t *top_global, t_state *state_local)
336 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
340 gmx_mtop_t *err_top_global;
341 gmx_localtop_t *err_top_local;
345 err_top_global = dd->reverse_top->err_top_global;
346 err_top_local = dd->reverse_top->err_top_local;
350 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
354 ndiff_tot = local_count - dd->nbonded_global;
356 for (ftype = 0; ftype < F_NRE; ftype++)
359 cl[ftype] = err_top_local->idef.il[ftype].nr/(1+nral);
362 gmx_sumi(F_NRE, cl, cr);
368 fprintf(fplog, "\nA list of missing interactions:\n");
370 fprintf(stderr, "\nA list of missing interactions:\n");
371 rest_global = dd->nbonded_global;
372 rest_local = local_count;
373 for (ftype = 0; ftype < F_NRE; ftype++)
375 /* In the reverse and local top all constraints are merged
376 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
377 * and add these constraints when doing F_CONSTR.
379 if (((interaction_function[ftype].flags & IF_BOND) &&
380 (dd->reverse_top->bBCheck
381 || !(interaction_function[ftype].flags & IF_LIMZERO)))
382 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
383 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
385 n = gmx_mtop_ftype_count(err_top_global, ftype);
386 if (ftype == F_CONSTR)
388 n += gmx_mtop_ftype_count(err_top_global, F_CONSTRNC);
390 ndiff = cl[ftype] - n;
393 sprintf(buf, "%20s of %6d missing %6d",
394 interaction_function[ftype].longname, n, -ndiff);
397 fprintf(fplog, "%s\n", buf);
399 fprintf(stderr, "%s\n", buf);
402 rest_local -= cl[ftype];
406 ndiff = rest_local - rest_global;
409 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
410 rest_global, -ndiff);
413 fprintf(fplog, "%s\n", buf);
415 fprintf(stderr, "%s\n", buf);
419 print_missing_interactions_atoms(fplog, cr, err_top_global,
420 &err_top_local->idef);
421 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
422 -1, state_local->x, state_local->box);
427 gmx_incons("One or more interactions were multiple assigned in the domain decompostion");
431 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));
436 /*! \brief Return global topology molecule information for global atom index \p i_gl */
437 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
438 int *mb, int *mt, int *mol, int *i_mol)
440 molblock_ind_t *mbi = rt->mbi;
442 int end = rt->nmolblock; /* exclusive */
445 /* binary search for molblock_ind */
449 if (i_gl >= mbi[mid].a_end)
453 else if (i_gl < mbi[mid].a_start)
467 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
468 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
471 /*! \brief Count the exclusions for all atoms in \p cgs */
472 static void count_excls(const t_block *cgs, const t_blocka *excls,
473 int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
475 int cg, at0, at1, at, excl, atj;
480 for (cg = 0; cg < cgs->nr; cg++)
482 at0 = cgs->index[cg];
483 at1 = cgs->index[cg+1];
484 for (at = at0; at < at1; at++)
486 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
488 atj = excls->a[excl];
492 if (atj < at0 || atj >= at1)
499 *n_excl_at_max = std::max(*n_excl_at_max,
500 excls->index[at+1] - excls->index[at]);
505 /*! \brief Run the reverse ilist generation and store it when \p bAssign = TRUE */
506 static int low_make_reverse_ilist(t_ilist *il_mt, t_atom *atom,
509 gmx_bool bConstr, gmx_bool bSettle,
511 int *r_index, int *r_il,
512 gmx_bool bLinkToAllAtoms,
515 int ftype, nral, i, j, nlink, link;
523 for (ftype = 0; ftype < F_NRE; ftype++)
525 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
526 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
527 (bSettle && ftype == F_SETTLE))
529 bVSite = (interaction_function[ftype].flags & IF_VSITE);
532 for (i = 0; i < il->nr; i += 1+nral)
539 /* We don't need the virtual sites for the cg-links */
549 /* Couple to the first atom in the interaction */
552 for (link = 0; link < nlink; link++)
557 r_il[r_index[a]+count[a]] =
558 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
559 r_il[r_index[a]+count[a]+1] = ia[0];
560 for (j = 1; j < 1+nral; j++)
562 /* Store the molecular atom number */
563 r_il[r_index[a]+count[a]+1+j] = ia[j];
566 if (interaction_function[ftype].flags & IF_VSITE)
570 /* Add an entry to iatoms for storing
571 * which of the constructing atoms are
574 r_il[r_index[a]+count[a]+2+nral] = 0;
575 for (j = 2; j < 1+nral; j++)
577 if (atom[ia[j]].ptype == eptVSite)
579 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
582 /* Store vsite pbc atom in a second extra entry */
583 r_il[r_index[a]+count[a]+2+nral+1] =
584 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
589 /* We do not count vsites since they are always
590 * uniquely assigned and can be assigned
591 * to multiple nodes with recursive vsites.
594 !(interaction_function[ftype].flags & IF_LIMZERO))
599 count[a] += 2 + nral_rt(ftype);
608 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
609 static int make_reverse_ilist(gmx_moltype_t *molt,
611 gmx_bool bConstr, gmx_bool bSettle,
613 gmx_bool bLinkToAllAtoms,
614 reverse_ilist_t *ril_mt)
616 int nat_mt, *count, i, nint_mt;
618 /* Count the interactions */
619 nat_mt = molt->atoms.nr;
621 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
623 bConstr, bSettle, bBCheck, NULL, NULL,
624 bLinkToAllAtoms, FALSE);
626 snew(ril_mt->index, nat_mt+1);
627 ril_mt->index[0] = 0;
628 for (i = 0; i < nat_mt; i++)
630 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
633 snew(ril_mt->il, ril_mt->index[nat_mt]);
635 /* Store the interactions */
637 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
639 bConstr, bSettle, bBCheck,
640 ril_mt->index, ril_mt->il,
641 bLinkToAllAtoms, TRUE);
648 /*! \brief Destroys a reverse_ilist_t struct */
649 static void destroy_reverse_ilist(reverse_ilist_t *ril)
655 /*! \brief Generate the reverse topology */
656 static gmx_reverse_top_t *make_reverse_top(gmx_mtop_t *mtop, gmx_bool bFE,
657 int ***vsite_pbc_molt,
658 gmx_bool bConstr, gmx_bool bSettle,
659 gmx_bool bBCheck, int *nint)
662 gmx_reverse_top_t *rt;
669 /* Should we include constraints (for SHAKE) in rt? */
670 rt->bConstr = bConstr;
671 rt->bSettle = bSettle;
672 rt->bBCheck = bBCheck;
674 rt->bMultiCGmols = FALSE;
675 snew(nint_mt, mtop->nmoltype);
676 snew(rt->ril_mt, mtop->nmoltype);
677 rt->ril_mt_tot_size = 0;
678 for (mt = 0; mt < mtop->nmoltype; mt++)
680 molt = &mtop->moltype[mt];
681 if (molt->cgs.nr > 1)
683 rt->bMultiCGmols = TRUE;
686 /* Make the atom to interaction list for this molecule type */
688 make_reverse_ilist(molt, vsite_pbc_molt ? vsite_pbc_molt[mt] : NULL,
689 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
692 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt->atoms.nr];
696 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
700 for (mb = 0; mb < mtop->nmolblock; mb++)
702 *nint += mtop->molblock[mb].nmol*nint_mt[mtop->molblock[mb].type];
706 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
708 rt->ilsort = ilsortFE_UNSORTED;
712 rt->ilsort = ilsortNO_FE;
715 /* Make a molblock index for fast searching */
716 snew(rt->mbi, mtop->nmolblock);
717 rt->nmolblock = mtop->nmolblock;
719 for (mb = 0; mb < mtop->nmolblock; mb++)
721 rt->mbi[mb].a_start = i;
722 i += mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
723 rt->mbi[mb].a_end = i;
724 rt->mbi[mb].natoms_mol = mtop->molblock[mb].natoms_mol;
725 rt->mbi[mb].type = mtop->molblock[mb].type;
728 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
729 snew(rt->th_work, rt->nthread);
730 if (vsite_pbc_molt != NULL)
732 for (thread = 0; thread < rt->nthread; thread++)
734 snew(rt->th_work[thread].vsite_pbc, F_VSITEN-F_VSITE2+1);
735 snew(rt->th_work[thread].vsite_pbc_nalloc, F_VSITEN-F_VSITE2+1);
742 void dd_make_reverse_top(FILE *fplog,
743 gmx_domdec_t *dd, gmx_mtop_t *mtop,
745 t_inputrec *ir, gmx_bool bBCheck)
749 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
752 /* If normal and/or settle constraints act only within charge groups,
753 * we can store them in the reverse top and simply assign them to domains.
754 * Otherwise we need to assign them to multiple domains and set up
755 * the parallel version constraint algoirthm(s).
758 dd->reverse_top = make_reverse_top(mtop, ir->efep != efepNO,
759 vsite ? vsite->vsite_pbc_molt : NULL,
760 !dd->bInterCGcons, !dd->bInterCGsettles,
761 bBCheck, &dd->nbonded_global);
763 gmx_reverse_top_t *rt = dd->reverse_top;
765 if (rt->ril_mt_tot_size >= 200000 &&
767 mtop->nmolblock == 1 && mtop->molblock[0].nmol == 1)
769 /* mtop comes from a pre Gromacs 4 tpr file */
770 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";
773 fprintf(fplog, "\n%s\n\n", note);
777 fprintf(stderr, "\n%s\n\n", note);
781 /* With the Verlet scheme, exclusions are handled in the non-bonded
782 * kernels and only exclusions inside the cut-off lead to exclusion
783 * forces. Since each atom pair is treated at most once in the non-bonded
784 * kernels, it doesn't matter if the exclusions for the same atom pair
785 * appear multiple times in the exclusion list. In contrast, the, old,
786 * group cut-off scheme loops over a list of exclusions, so there each
787 * excluded pair should appear exactly once.
789 rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
790 IR_EXCL_FORCES(*ir));
795 dd->n_intercg_excl = 0;
796 rt->n_excl_at_max = 0;
797 for (mb = 0; mb < mtop->nmolblock; mb++)
799 gmx_molblock_t *molb;
801 int n_excl_mol, n_excl_icg, n_excl_at_max;
803 molb = &mtop->molblock[mb];
804 molt = &mtop->moltype[molb->type];
805 count_excls(&molt->cgs, &molt->excls,
806 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
807 nexcl += molb->nmol*n_excl_mol;
808 dd->n_intercg_excl += molb->nmol*n_excl_icg;
809 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
811 if (rt->bExclRequired)
813 dd->nbonded_global += nexcl;
814 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
816 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
817 "will use an extra communication step for exclusion forces for %s\n",
818 dd->n_intercg_excl, eel_names[ir->coulombtype]);
822 if (vsite && vsite->n_intercg_vsite > 0)
826 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
827 "will an extra communication step for selected coordinates and forces\n",
828 vsite->n_intercg_vsite);
830 init_domdec_vsites(dd, vsite->n_intercg_vsite);
833 if (dd->bInterCGcons || dd->bInterCGsettles)
835 init_domdec_constraints(dd, mtop);
839 fprintf(fplog, "\n");
843 /*! \brief Store a vsite interaction at the end of \p il
845 * This routine is very similar to add_ifunc, but vsites interactions
846 * have more work to do than other kinds of interactions, and the
847 * complex way nral (and thus vector contents) depends on ftype
848 * confuses static analysis tools unless we fuse the vsite
849 * atom-indexing organization code with the ifunc-adding code, so that
850 * they can see that nral is the same value. */
851 static gmx_inline void
852 add_ifunc_for_vsites(t_iatom *tiatoms, gmx_ga2la_t ga2la,
853 int nral, gmx_bool bHomeA,
854 int a, int a_gl, int a_mol,
855 const t_iatom *iatoms,
860 if (il->nr+1+nral > il->nalloc)
862 il->nalloc = over_alloc_large(il->nr+1+nral);
863 srenew(il->iatoms, il->nalloc);
865 liatoms = il->iatoms + il->nr;
869 tiatoms[0] = iatoms[0];
873 /* We know the local index of the first atom */
878 /* Convert later in make_local_vsites */
879 tiatoms[1] = -a_gl - 1;
882 for (int k = 2; k < 1+nral; k++)
884 int ak_gl = a_gl + iatoms[k] - a_mol;
885 if (!ga2la_get_home(ga2la, ak_gl, &tiatoms[k]))
887 /* Copy the global index, convert later in make_local_vsites */
888 tiatoms[k] = -(ak_gl + 1);
890 // Note that ga2la_get_home always sets the third parameter if
893 for (int k = 0; k < 1+nral; k++)
895 liatoms[k] = tiatoms[k];
899 /*! \brief Store a bonded interaction at the end of \p il */
900 static gmx_inline void add_ifunc(int nral, t_iatom *tiatoms, t_ilist *il)
905 if (il->nr+1+nral > il->nalloc)
907 il->nalloc = over_alloc_large(il->nr+1+nral);
908 srenew(il->iatoms, il->nalloc);
910 liatoms = il->iatoms + il->nr;
911 for (k = 0; k <= nral; k++)
913 liatoms[k] = tiatoms[k];
918 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
919 static void add_posres(int mol, int a_mol, const gmx_molblock_t *molb,
920 t_iatom *iatoms, const t_iparams *ip_in,
926 /* This position restraint has not been added yet,
927 * so it's index is the current number of position restraints.
929 n = idef->il[F_POSRES].nr/2;
930 if (n+1 > idef->iparams_posres_nalloc)
932 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
933 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
935 ip = &idef->iparams_posres[n];
936 /* Copy the force constants */
937 *ip = ip_in[iatoms[0]];
939 /* Get the position restraint coordinates from the molblock */
940 a_molb = mol*molb->natoms_mol + a_mol;
941 if (a_molb >= molb->nposres_xA)
943 gmx_incons("Not enough position restraint coordinates");
945 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
946 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
947 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
948 if (molb->nposres_xB > 0)
950 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
951 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
952 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
956 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
957 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
958 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
960 /* Set the parameter index for idef->iparams_posre */
964 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
965 static void add_fbposres(int mol, int a_mol, const gmx_molblock_t *molb,
966 t_iatom *iatoms, const t_iparams *ip_in,
972 /* This flat-bottom position restraint has not been added yet,
973 * so it's index is the current number of position restraints.
975 n = idef->il[F_FBPOSRES].nr/2;
976 if (n+1 > idef->iparams_fbposres_nalloc)
978 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
979 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
981 ip = &idef->iparams_fbposres[n];
982 /* Copy the force constants */
983 *ip = ip_in[iatoms[0]];
985 /* Get the position restriant coordinats from the molblock */
986 a_molb = mol*molb->natoms_mol + a_mol;
987 if (a_molb >= molb->nposres_xA)
989 gmx_incons("Not enough position restraint coordinates");
991 /* Take reference positions from A position of normal posres */
992 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
993 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
994 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
996 /* Note: no B-type for flat-bottom posres */
998 /* Set the parameter index for idef->iparams_posre */
1002 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
1003 static void add_vsite(gmx_ga2la_t ga2la, int *index, int *rtil,
1004 int ftype, int nral,
1005 gmx_bool bHomeA, int a, int a_gl, int a_mol,
1007 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
1009 int k, vsi, pbc_a_mol;
1010 t_iatom tiatoms[1+MAXATOMLIST];
1011 int j, ftype_r, nral_r;
1013 /* Add this interaction to the local topology */
1014 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1018 vsi = idef->il[ftype].nr/(1+nral) - 1;
1019 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
1021 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
1022 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
1026 pbc_a_mol = iatoms[1+nral+1];
1029 /* The pbc flag is one of the following two options:
1030 * -2: vsite and all constructing atoms are within the same cg, no pbc
1031 * -1: vsite and its first constructing atom are in the same cg, do pbc
1033 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
1037 /* Set the pbc atom for this vsite so we can make its pbc
1038 * identical to the rest of the atoms in its charge group.
1039 * Since the order of the atoms does not change within a charge
1040 * group, we do not need the global to local atom index.
1042 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
1047 /* This vsite is non-home (required for recursion),
1048 * and therefore there is no charge group to match pbc with.
1049 * But we always turn on full_pbc to assure that higher order
1050 * recursion works correctly.
1052 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
1058 /* Check for recursion */
1059 for (k = 2; k < 1+nral; k++)
1061 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1063 /* This construction atoms is a vsite and not a home atom */
1066 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1068 /* Find the vsite construction */
1070 /* Check all interactions assigned to this atom */
1071 j = index[iatoms[k]];
1072 while (j < index[iatoms[k]+1])
1074 ftype_r = rtil[j++];
1075 nral_r = NRAL(ftype_r);
1076 if (interaction_function[ftype_r].flags & IF_VSITE)
1078 /* Add this vsite (recursion) */
1079 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1080 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1081 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
1082 j += 1 + nral_r + 2;
1094 /*! \brief Update the local atom to local charge group index */
1095 static void make_la2lc(gmx_domdec_t *dd)
1097 int *cgindex, *la2lc, cg, a;
1099 cgindex = dd->cgindex;
1101 if (dd->nat_tot > dd->la2lc_nalloc)
1103 dd->la2lc_nalloc = over_alloc_dd(dd->nat_tot);
1104 snew(dd->la2lc, dd->la2lc_nalloc);
1108 /* Make the local atom to local cg index */
1109 for (cg = 0; cg < dd->ncg_tot; cg++)
1111 for (a = cgindex[cg]; a < cgindex[cg+1]; a++)
1118 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1119 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1125 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1129 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1135 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1136 static void combine_blocka(t_blocka *dest, const thread_work_t *src, int nsrc)
1140 ni = src[nsrc-1].excl.nr;
1142 for (s = 0; s < nsrc; s++)
1144 na += src[s].excl.nra;
1146 if (ni + 1 > dest->nalloc_index)
1148 dest->nalloc_index = over_alloc_large(ni+1);
1149 srenew(dest->index, dest->nalloc_index);
1151 if (dest->nra + na > dest->nalloc_a)
1153 dest->nalloc_a = over_alloc_large(dest->nra+na);
1154 srenew(dest->a, dest->nalloc_a);
1156 for (s = 1; s < nsrc; s++)
1158 for (i = dest->nr+1; i < src[s].excl.nr+1; i++)
1160 dest->index[i] = dest->nra + src[s].excl.index[i];
1162 for (i = 0; i < src[s].excl.nra; i++)
1164 dest->a[dest->nra+i] = src[s].excl.a[i];
1166 dest->nr = src[s].excl.nr;
1167 dest->nra += src[s].excl.nra;
1171 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1172 * virtual sites need special attention, as pbc info differs per vsite.
1174 static void combine_idef(t_idef *dest, const thread_work_t *src, int nsrc,
1179 for (ftype = 0; ftype < F_NRE; ftype++)
1184 for (s = 1; s < nsrc; s++)
1186 n += src[s].idef.il[ftype].nr;
1192 ild = &dest->il[ftype];
1194 if (ild->nr + n > ild->nalloc)
1196 ild->nalloc = over_alloc_large(ild->nr+n);
1197 srenew(ild->iatoms, ild->nalloc);
1201 int nral1 = 0, ftv = 0;
1203 vpbc = ((interaction_function[ftype].flags & IF_VSITE) &&
1204 vsite->vsite_pbc_loc != NULL);
1207 nral1 = 1 + NRAL(ftype);
1208 ftv = ftype - F_VSITE2;
1209 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1211 vsite->vsite_pbc_loc_nalloc[ftv] =
1212 over_alloc_large((ild->nr + n)/nral1);
1213 srenew(vsite->vsite_pbc_loc[ftv],
1214 vsite->vsite_pbc_loc_nalloc[ftv]);
1218 for (s = 1; s < nsrc; s++)
1223 ils = &src[s].idef.il[ftype];
1224 for (i = 0; i < ils->nr; i++)
1226 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1230 for (i = 0; i < ils->nr; i += nral1)
1232 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1233 src[s].vsite_pbc[ftv][i/nral1];
1242 /* Position restraints need an additional treatment */
1243 if (dest->il[F_POSRES].nr > 0)
1247 n = dest->il[F_POSRES].nr/2;
1248 if (n > dest->iparams_posres_nalloc)
1250 dest->iparams_posres_nalloc = over_alloc_large(n);
1251 srenew(dest->iparams_posres, dest->iparams_posres_nalloc);
1253 /* Set n to the number of original position restraints in dest */
1254 for (s = 1; s < nsrc; s++)
1256 n -= src[s].idef.il[F_POSRES].nr/2;
1258 for (s = 1; s < nsrc; s++)
1260 for (i = 0; i < src[s].idef.il[F_POSRES].nr/2; i++)
1262 /* Correct the index into iparams_posres */
1263 dest->il[F_POSRES].iatoms[n*2] = n;
1264 /* Copy the position restraint force parameters */
1265 dest->iparams_posres[n] = src[s].idef.iparams_posres[i];
1272 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1274 * With thread parallelizing each thread acts on a different atom range:
1275 * at_start to at_end.
1277 static int make_bondeds_zone(gmx_domdec_t *dd,
1278 const gmx_domdec_zones_t *zones,
1279 const gmx_molblock_t *molb,
1280 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1282 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1283 const t_iparams *ip_in,
1286 int *vsite_pbc_nalloc,
1288 int at_start, int at_end)
1290 int i, i_gl, mb, mt, mol, i_mol, j, ftype, nral, d, k;
1292 t_iatom *iatoms, tiatoms[1+MAXATOMLIST];
1293 gmx_bool bBCheck, bUse, bLocal;
1294 ivec k_zero, k_plus;
1299 const gmx_domdec_ns_ranges_t *izone;
1300 gmx_reverse_top_t *rt;
1303 nizone = zones->nizone;
1304 izone = zones->izone;
1306 rt = dd->reverse_top;
1308 bBCheck = rt->bBCheck;
1314 for (i = at_start; i < at_end; i++)
1316 /* Get the global atom number */
1317 i_gl = dd->gatindex[i];
1318 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1319 /* Check all interactions assigned to this atom */
1320 index = rt->ril_mt[mt].index;
1321 rtil = rt->ril_mt[mt].il;
1323 while (j < index[i_mol+1])
1328 if (ftype == F_SETTLE)
1330 /* Settles are only in the reverse top when they
1331 * operate within a charge group. So we can assign
1332 * them without checks. We do this only for performance
1333 * reasons; it could be handled by the code below.
1337 /* Home zone: add this settle to the local topology */
1338 tiatoms[0] = iatoms[0];
1340 tiatoms[2] = i + iatoms[2] - iatoms[1];
1341 tiatoms[3] = i + iatoms[3] - iatoms[1];
1342 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1347 else if (interaction_function[ftype].flags & IF_VSITE)
1349 /* The vsite construction goes where the vsite itself is */
1352 add_vsite(dd->ga2la, index, rtil, ftype, nral,
1353 TRUE, i, i_gl, i_mol,
1354 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1361 tiatoms[0] = iatoms[0];
1365 /* Assign single-body interactions to the home zone */
1370 if (ftype == F_POSRES)
1372 add_posres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1375 else if (ftype == F_FBPOSRES)
1377 add_fbposres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1388 /* This is a two-body interaction, we can assign
1389 * analogous to the non-bonded assignments.
1391 if (!ga2la_get(ga2la, i_gl+iatoms[2]-i_mol, &a_loc, &kz))
1401 /* Check zone interaction assignments */
1402 bUse = ((iz < nizone && iz <= kz &&
1403 izone[iz].j0 <= kz && kz < izone[iz].j1) ||
1404 (kz < nizone &&iz > kz &&
1405 izone[kz].j0 <= iz && iz < izone[kz].j1));
1410 /* If necessary check the cgcm distance */
1412 dd_dist2(pbc_null, cg_cm, la2lc,
1413 tiatoms[1], tiatoms[2]) >= rc2)
1422 /* Assign this multi-body bonded interaction to
1423 * the local node if we have all the atoms involved
1424 * (local or communicated) and the minimum zone shift
1425 * in each dimension is zero, for dimensions
1426 * with 2 DD cells an extra check may be necessary.
1431 for (k = 1; k <= nral && bUse; k++)
1433 bLocal = ga2la_get(ga2la, i_gl+iatoms[k]-i_mol,
1435 if (!bLocal || kz >= zones->n)
1437 /* We do not have this atom of this interaction
1438 * locally, or it comes from more than one cell
1446 for (d = 0; d < DIM; d++)
1448 if (zones->shift[kz][d] == 0)
1460 k_zero[XX] && k_zero[YY] && k_zero[ZZ]);
1463 for (d = 0; (d < DIM && bUse); d++)
1465 /* Check if the cg_cm distance falls within
1466 * the cut-off to avoid possible multiple
1467 * assignments of bonded interactions.
1471 dd_dist2(pbc_null, cg_cm, la2lc,
1472 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1481 /* Add this interaction to the local topology */
1482 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1483 /* Sum so we can check in global_stat
1484 * if we have everything.
1487 !(interaction_function[ftype].flags & IF_LIMZERO))
1497 return nbonded_local;
1500 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1501 static void set_no_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1502 int iz, t_blocka *lexcls)
1506 a0 = dd->cgindex[zones->cg_range[iz]];
1507 a1 = dd->cgindex[zones->cg_range[iz+1]];
1509 for (a = a0+1; a < a1+1; a++)
1511 lexcls->index[a] = lexcls->nra;
1515 /*! \brief Set the exclusion data for i-zone \p iz
1517 * This is a legacy version for the group scheme of the same routine below.
1518 * Here charge groups and distance checks to ensure unique exclusions
1521 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1522 const gmx_moltype_t *moltype,
1523 gmx_bool bRCheck, real rc2,
1524 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1528 int cg_start, int cg_end)
1530 int n_excl_at_max, n, count, jla0, jla1, jla;
1531 int cg, la0, la1, la, a_gl, mb, mt, mol, a_mol, j, aj_mol;
1532 const t_blocka *excls;
1538 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1539 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1541 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1543 /* We set the end index, but note that we might not start at zero here */
1544 lexcls->nr = dd->cgindex[cg_end];
1548 for (cg = cg_start; cg < cg_end; cg++)
1550 if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1552 lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1553 srenew(lexcls->a, lexcls->nalloc_a);
1555 la0 = dd->cgindex[cg];
1556 la1 = dd->cgindex[cg+1];
1557 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1558 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1560 /* Copy the exclusions from the global top */
1561 for (la = la0; la < la1; la++)
1563 lexcls->index[la] = n;
1564 a_gl = dd->gatindex[la];
1565 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1566 excls = &moltype[mt].excls;
1567 for (j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1569 aj_mol = excls->a[j];
1570 /* This computation of jla is only correct intra-cg */
1571 jla = la + aj_mol - a_mol;
1572 if (jla >= la0 && jla < la1)
1574 /* This is an intra-cg exclusion. We can skip
1575 * the global indexing and distance checking.
1577 /* Intra-cg exclusions are only required
1578 * for the home zone.
1582 lexcls->a[n++] = jla;
1583 /* Check to avoid double counts */
1592 /* This is a inter-cg exclusion */
1593 /* Since exclusions are pair interactions,
1594 * just like non-bonded interactions,
1595 * they can be assigned properly up
1596 * to the DD cutoff (not cutoff_min as
1597 * for the other bonded interactions).
1599 if (ga2la_get(ga2la, a_gl+aj_mol-a_mol, &jla, &cell))
1601 if (iz == 0 && cell == 0)
1603 lexcls->a[n++] = jla;
1604 /* Check to avoid double counts */
1610 else if (jla >= jla0 && jla < jla1 &&
1612 dd_dist2(pbc_null, cg_cm, la2lc, la, jla) < rc2))
1614 /* jla > la, since jla0 > la */
1615 lexcls->a[n++] = jla;
1625 /* There are no inter-cg excls and this cg is self-excluded.
1626 * These exclusions are only required for zone 0,
1627 * since other zones do not see themselves.
1631 for (la = la0; la < la1; la++)
1633 lexcls->index[la] = n;
1634 for (j = la0; j < la1; j++)
1639 count += ((la1 - la0)*(la1 - la0 - 1))/2;
1643 /* We don't need exclusions for this cg */
1644 for (la = la0; la < la1; la++)
1646 lexcls->index[la] = n;
1652 lexcls->index[lexcls->nr] = n;
1658 /*! \brief Set the exclusion data for i-zone \p iz */
1659 static void make_exclusions_zone(gmx_domdec_t *dd,
1660 gmx_domdec_zones_t *zones,
1661 const gmx_moltype_t *moltype,
1665 int at_start, int at_end)
1669 int n_excl_at_max, n, at;
1673 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1674 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1676 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1678 /* We set the end index, but note that we might not start at zero here */
1679 lexcls->nr = at_end;
1682 for (at = at_start; at < at_end; at++)
1684 if (n + 1000 > lexcls->nalloc_a)
1686 lexcls->nalloc_a = over_alloc_large(n + 1000);
1687 srenew(lexcls->a, lexcls->nalloc_a);
1689 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1691 int a_gl, mb, mt, mol, a_mol, j;
1692 const t_blocka *excls;
1694 if (n + n_excl_at_max > lexcls->nalloc_a)
1696 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1697 srenew(lexcls->a, lexcls->nalloc_a);
1700 /* Copy the exclusions from the global top */
1701 lexcls->index[at] = n;
1702 a_gl = dd->gatindex[at];
1703 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1704 &mb, &mt, &mol, &a_mol);
1705 excls = &moltype[mt].excls;
1706 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1708 int aj_mol, at_j, cell;
1710 aj_mol = excls->a[j];
1712 if (ga2la_get(ga2la, a_gl + aj_mol - a_mol, &at_j, &cell))
1714 /* This check is not necessary, but it can reduce
1715 * the number of exclusions in the list, which in turn
1716 * can speed up the pair list construction a bit.
1718 if (at_j >= jla0 && at_j < jla1)
1720 lexcls->a[n++] = at_j;
1727 /* We don't need exclusions for this atom */
1728 lexcls->index[at] = n;
1732 lexcls->index[lexcls->nr] = n;
1737 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1738 static void check_alloc_index(t_blocka *ba, int nindex_max)
1740 if (nindex_max+1 > ba->nalloc_index)
1742 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1743 srenew(ba->index, ba->nalloc_index);
1747 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1748 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1754 nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1756 check_alloc_index(lexcls, nr);
1758 for (thread = 1; thread < dd->reverse_top->nthread; thread++)
1760 check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1764 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1765 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1770 lexcls->nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1772 if (dd->n_intercg_excl == 0)
1774 /* There are no exclusions involving non-home charge groups,
1775 * but we need to set the indices for neighborsearching.
1777 la0 = dd->cgindex[zones->izone[0].cg1];
1778 for (la = la0; la < lexcls->nr; la++)
1780 lexcls->index[la] = lexcls->nra;
1783 /* nr is only used to loop over the exclusions for Ewald and RF,
1784 * so we can set it to the number of home atoms for efficiency.
1786 lexcls->nr = dd->cgindex[zones->izone[0].cg1];
1790 /*! \brief Clear a t_idef struct */
1791 static void clear_idef(t_idef *idef)
1795 /* Clear the counts */
1796 for (ftype = 0; ftype < F_NRE; ftype++)
1798 idef->il[ftype].nr = 0;
1802 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1803 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1804 gmx_domdec_zones_t *zones,
1805 const gmx_mtop_t *mtop,
1807 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1809 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1810 t_idef *idef, gmx_vsite_t *vsite,
1811 t_blocka *lexcls, int *excl_count)
1813 int nzone_bondeds, nzone_excl;
1818 gmx_reverse_top_t *rt;
1820 if (dd->reverse_top->bMultiCGmols)
1822 nzone_bondeds = zones->n;
1826 /* Only single charge group molecules, so interactions don't
1827 * cross zone boundaries and we only need to assign in the home zone.
1832 if (dd->n_intercg_excl > 0)
1834 /* We only use exclusions from i-zones to i- and j-zones */
1835 nzone_excl = zones->nizone;
1839 /* There are no inter-cg exclusions and only zone 0 sees itself */
1843 check_exclusions_alloc(dd, zones, lexcls);
1845 rt = dd->reverse_top;
1849 /* Clear the counts */
1857 for (iz = 0; iz < nzone_bondeds; iz++)
1859 cg0 = zones->cg_range[iz];
1860 cg1 = zones->cg_range[iz+1];
1862 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1863 for (thread = 0; thread < rt->nthread; thread++)
1868 int *vsite_pbc_nalloc;
1871 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
1872 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
1880 idef_t = &rt->th_work[thread].idef;
1884 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1888 vsite_pbc = vsite->vsite_pbc_loc;
1889 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
1893 vsite_pbc = rt->th_work[thread].vsite_pbc;
1894 vsite_pbc_nalloc = rt->th_work[thread].vsite_pbc_nalloc;
1900 vsite_pbc_nalloc = NULL;
1903 rt->th_work[thread].nbonded =
1904 make_bondeds_zone(dd, zones,
1906 bRCheckMB, rcheck, bRCheck2B, rc2,
1907 la2lc, pbc_null, cg_cm, idef->iparams,
1909 vsite_pbc, vsite_pbc_nalloc,
1911 dd->cgindex[cg0t], dd->cgindex[cg1t]);
1913 if (iz < nzone_excl)
1921 excl_t = &rt->th_work[thread].excl;
1926 if (dd->cgindex[dd->ncg_tot] == dd->ncg_tot &&
1929 /* No charge groups and no distance check required */
1930 make_exclusions_zone(dd, zones,
1931 mtop->moltype, cginfo,
1938 rt->th_work[thread].excl_count =
1939 make_exclusions_zone_cg(dd, zones,
1940 mtop->moltype, bRCheck2B, rc2,
1941 la2lc, pbc_null, cg_cm, cginfo,
1949 if (rt->nthread > 1)
1951 combine_idef(idef, rt->th_work, rt->nthread, vsite);
1954 for (thread = 0; thread < rt->nthread; thread++)
1956 nbonded_local += rt->th_work[thread].nbonded;
1959 if (iz < nzone_excl)
1961 if (rt->nthread > 1)
1963 combine_blocka(lexcls, rt->th_work, rt->nthread);
1966 for (thread = 0; thread < rt->nthread; thread++)
1968 *excl_count += rt->th_work[thread].excl_count;
1973 /* Some zones might not have exclusions, but some code still needs to
1974 * loop over the index, so we set the indices here.
1976 for (iz = nzone_excl; iz < zones->nizone; iz++)
1978 set_no_exclusions_zone(dd, zones, iz, lexcls);
1981 finish_local_exclusions(dd, zones, lexcls);
1984 fprintf(debug, "We have %d exclusions, check count %d\n",
1985 lexcls->nra, *excl_count);
1988 return nbonded_local;
1991 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
1993 lcgs->nr = dd->ncg_tot;
1994 lcgs->index = dd->cgindex;
1997 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1998 int npbcdim, matrix box,
1999 rvec cellsize_min, ivec npulse,
2003 gmx_mtop_t *mtop, gmx_localtop_t *ltop)
2005 gmx_bool bRCheckMB, bRCheck2B;
2009 t_pbc pbc, *pbc_null = NULL;
2013 fprintf(debug, "Making local topology\n");
2016 dd_make_local_cgs(dd, <op->cgs);
2021 if (dd->reverse_top->bMultiCGmols)
2023 /* We need to check to which cell bondeds should be assigned */
2024 rc = dd_cutoff_twobody(dd);
2027 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2030 /* Should we check cg_cm distances when assigning bonded interactions? */
2031 for (d = 0; d < DIM; d++)
2034 /* Only need to check for dimensions where the part of the box
2035 * that is not communicated is smaller than the cut-off.
2037 if (d < npbcdim && dd->nc[d] > 1 &&
2038 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2045 /* Check for interactions between two atoms,
2046 * where we can allow interactions up to the cut-off,
2047 * instead of up to the smallest cell dimension.
2054 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
2055 d, cellsize_min[d], d, rcheck[d], bRCheck2B);
2058 if (bRCheckMB || bRCheck2B)
2063 set_pbc_dd(&pbc, fr->ePBC, dd, TRUE, box);
2074 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
2075 bRCheckMB, rcheck, bRCheck2B, rc,
2077 pbc_null, cgcm_or_x,
2079 <op->excls, &nexcl);
2081 /* The ilist is not sorted yet,
2082 * we can only do this when we have the charge arrays.
2084 ltop->idef.ilsort = ilsortUNKNOWN;
2086 if (dd->reverse_top->bExclRequired)
2088 dd->nbonded_local += nexcl;
2090 forcerec_set_excl_load(fr, ltop);
2093 ltop->atomtypes = mtop->atomtypes;
2095 /* For an error message only */
2096 dd->reverse_top->err_top_global = mtop;
2097 dd->reverse_top->err_top_local = ltop;
2100 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
2101 gmx_localtop_t *ltop)
2103 if (dd->reverse_top->ilsort == ilsortNO_FE)
2105 ltop->idef.ilsort = ilsortNO_FE;
2109 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
2113 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global)
2115 gmx_localtop_t *top;
2120 top->idef.ntypes = top_global->ffparams.ntypes;
2121 top->idef.atnr = top_global->ffparams.atnr;
2122 top->idef.functype = top_global->ffparams.functype;
2123 top->idef.iparams = top_global->ffparams.iparams;
2124 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
2125 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
2127 for (i = 0; i < F_NRE; i++)
2129 top->idef.il[i].iatoms = NULL;
2130 top->idef.il[i].nalloc = 0;
2132 top->idef.ilsort = ilsortUNKNOWN;
2137 void dd_init_local_state(gmx_domdec_t *dd,
2138 t_state *state_global, t_state *state_local)
2140 int buf[NITEM_DD_INIT_LOCAL_STATE];
2144 buf[0] = state_global->flags;
2145 buf[1] = state_global->ngtc;
2146 buf[2] = state_global->nnhpres;
2147 buf[3] = state_global->nhchainlength;
2148 buf[4] = state_global->dfhist.nlambda;
2150 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2152 init_state(state_local, 0, buf[1], buf[2], buf[3], buf[4]);
2153 state_local->flags = buf[0];
2156 /*! \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 */
2157 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2163 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2165 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2166 if (link->a[k] == cg_gl_j)
2173 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2174 "Inconsistent allocation of link");
2175 /* Add this charge group link */
2176 if (link->index[cg_gl+1]+1 > link->nalloc_a)
2178 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2179 srenew(link->a, link->nalloc_a);
2181 link->a[link->index[cg_gl+1]] = cg_gl_j;
2182 link->index[cg_gl+1]++;
2186 /*! \brief Allocate and return an array of the charge group index for all atoms */
2187 static int *make_at2cg(t_block *cgs)
2191 snew(at2cg, cgs->index[cgs->nr]);
2192 for (cg = 0; cg < cgs->nr; cg++)
2194 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2203 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
2204 cginfo_mb_t *cginfo_mb)
2206 gmx_reverse_top_t *rt;
2207 int mb, cg_offset, cg, cg_gl, a, aj, i, j, ftype, nral, nlink_mol, mol, ncgi;
2208 gmx_molblock_t *molb;
2209 gmx_moltype_t *molt;
2213 reverse_ilist_t ril;
2215 cginfo_mb_t *cgi_mb;
2217 /* For each charge group make a list of other charge groups
2218 * in the system that a linked to it via bonded interactions
2219 * which are also stored in reverse_top.
2222 rt = dd->reverse_top;
2225 snew(link->index, ncg_mtop(mtop)+1);
2232 for (mb = 0; mb < mtop->nmolblock; mb++)
2234 molb = &mtop->molblock[mb];
2235 if (molb->nmol == 0)
2239 molt = &mtop->moltype[molb->type];
2241 excls = &molt->excls;
2242 a2c = make_at2cg(cgs);
2243 /* Make a reverse ilist in which the interactions are linked
2244 * to all atoms, not only the first atom as in gmx_reverse_top.
2245 * The constraints are discarded here.
2247 make_reverse_ilist(molt, NULL, FALSE, FALSE, FALSE, TRUE, &ril);
2249 cgi_mb = &cginfo_mb[mb];
2251 for (cg = 0; cg < cgs->nr; cg++)
2253 cg_gl = cg_offset + cg;
2254 link->index[cg_gl+1] = link->index[cg_gl];
2255 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2258 while (i < ril.index[a+1])
2260 ftype = ril.il[i++];
2262 /* Skip the ifunc index */
2264 for (j = 0; j < nral; j++)
2269 check_link(link, cg_gl, cg_offset+a2c[aj]);
2272 i += nral_rt(ftype);
2274 if (rt->bExclRequired)
2276 /* Exclusions always go both ways */
2277 for (j = excls->index[a]; j < excls->index[a+1]; j++)
2282 check_link(link, cg_gl, cg_offset+a2c[aj]);
2287 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2289 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2293 nlink_mol = link->index[cg_offset+cgs->nr] - link->index[cg_offset];
2295 cg_offset += cgs->nr;
2297 destroy_reverse_ilist(&ril);
2302 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt->name, cgs->nr, nlink_mol);
2307 /* Copy the data for the rest of the molecules in this block */
2308 link->nalloc_a += (molb->nmol - 1)*nlink_mol;
2309 srenew(link->a, link->nalloc_a);
2310 for (mol = 1; mol < molb->nmol; mol++)
2312 for (cg = 0; cg < cgs->nr; cg++)
2314 cg_gl = cg_offset + cg;
2315 link->index[cg_gl+1] =
2316 link->index[cg_gl+1-cgs->nr] + nlink_mol;
2317 for (j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2319 link->a[j] = link->a[j-nlink_mol] + cgs->nr;
2321 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2322 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2324 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2328 cg_offset += cgs->nr;
2335 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2341 /*! \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 */
2342 static void bonded_cg_distance_mol(gmx_moltype_t *molt, int *at2cg,
2343 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2344 real *r_2b, int *ft2b, int *a2_1, int *a2_2,
2345 real *r_mb, int *ftmb, int *am_1, int *am_2)
2347 int ftype, nral, i, j, ai, aj, cgi, cgj;
2350 real r2_2b, r2_mb, rij2;
2354 for (ftype = 0; ftype < F_NRE; ftype++)
2356 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2358 il = &molt->ilist[ftype];
2362 for (i = 0; i < il->nr; i += 1+nral)
2364 for (ai = 0; ai < nral; ai++)
2366 cgi = at2cg[il->iatoms[i+1+ai]];
2367 for (aj = 0; aj < nral; aj++)
2369 cgj = at2cg[il->iatoms[i+1+aj]];
2372 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2373 if (nral == 2 && rij2 > r2_2b)
2377 *a2_1 = il->iatoms[i+1+ai];
2378 *a2_2 = il->iatoms[i+1+aj];
2380 if (nral > 2 && rij2 > r2_mb)
2384 *am_1 = il->iatoms[i+1+ai];
2385 *am_2 = il->iatoms[i+1+aj];
2396 excls = &molt->excls;
2397 for (ai = 0; ai < excls->nr; ai++)
2400 for (j = excls->index[ai]; j < excls->index[ai+1]; j++)
2402 cgj = at2cg[excls->a[j]];
2405 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2415 *r_2b = sqrt(r2_2b);
2416 *r_mb = sqrt(r2_mb);
2419 //! Compute charge group centers of mass for molecule \p molt
2420 static void get_cgcm_mol(gmx_moltype_t *molt, gmx_ffparams_t *ffparams,
2421 int ePBC, t_graph *graph, matrix box,
2423 rvec *x, rvec *xs, rvec *cg_cm)
2427 if (ePBC != epbcNONE)
2429 mk_mshift(NULL, graph, ePBC, box, x);
2431 shift_x(graph, box, x, xs);
2432 /* By doing an extra mk_mshift the molecules that are broken
2433 * because they were e.g. imported from another software
2434 * will be made whole again. Such are the healing powers
2437 mk_mshift(NULL, graph, ePBC, box, xs);
2441 /* We copy the coordinates so the original coordinates remain
2442 * unchanged, just to be 100% sure that we do not affect
2443 * binary reproducibility of simulations.
2445 n = molt->cgs.index[molt->cgs.nr];
2446 for (i = 0; i < n; i++)
2448 copy_rvec(x[i], xs[i]);
2454 construct_vsites(vsite, xs, 0.0, NULL,
2455 ffparams->iparams, molt->ilist,
2456 epbcNONE, TRUE, NULL, NULL);
2459 calc_cgcm(NULL, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2462 //! Returns whether \p molt has a virtual site
2463 static int have_vsite_molt(gmx_moltype_t *molt)
2469 for (i = 0; i < F_NRE; i++)
2471 if ((interaction_function[i].flags & IF_VSITE) &&
2472 molt->ilist[i].nr > 0)
2481 void dd_bonded_cg_distance(FILE *fplog,
2483 t_inputrec *ir, rvec *x, matrix box,
2485 real *r_2b, real *r_mb)
2487 gmx_bool bExclRequired;
2488 int mb, at_offset, *at2cg, mol;
2491 gmx_molblock_t *molb;
2492 gmx_moltype_t *molt;
2494 real rmol_2b, rmol_mb;
2495 int ft2b = -1, a_2b_1 = -1, a_2b_2 = -1, ftmb = -1, a_mb_1 = -1, a_mb_2 = -1;
2496 int ftm2b = -1, amol_2b_1 = -1, amol_2b_2 = -1, ftmmb = -1, amol_mb_1 = -1, amol_mb_2 = -1;
2498 bExclRequired = IR_EXCL_FORCES(*ir);
2500 vsite = init_vsite(mtop, NULL, TRUE);
2505 for (mb = 0; mb < mtop->nmolblock; mb++)
2507 molb = &mtop->molblock[mb];
2508 molt = &mtop->moltype[molb->type];
2509 if (molt->cgs.nr == 1 || molb->nmol == 0)
2511 at_offset += molb->nmol*molt->atoms.nr;
2515 if (ir->ePBC != epbcNONE)
2517 mk_graph_ilist(NULL, molt->ilist, 0, molt->atoms.nr, FALSE, FALSE,
2521 at2cg = make_at2cg(&molt->cgs);
2522 snew(xs, molt->atoms.nr);
2523 snew(cg_cm, molt->cgs.nr);
2524 for (mol = 0; mol < molb->nmol; mol++)
2526 get_cgcm_mol(molt, &mtop->ffparams, ir->ePBC, &graph, box,
2527 have_vsite_molt(molt) ? vsite : NULL,
2528 x+at_offset, xs, cg_cm);
2530 bonded_cg_distance_mol(molt, at2cg, bBCheck, bExclRequired, cg_cm,
2531 &rmol_2b, &ftm2b, &amol_2b_1, &amol_2b_2,
2532 &rmol_mb, &ftmmb, &amol_mb_1, &amol_mb_2);
2533 if (rmol_2b > *r_2b)
2537 a_2b_1 = at_offset + amol_2b_1;
2538 a_2b_2 = at_offset + amol_2b_2;
2540 if (rmol_mb > *r_mb)
2544 a_mb_1 = at_offset + amol_mb_1;
2545 a_mb_2 = at_offset + amol_mb_2;
2548 at_offset += molt->atoms.nr;
2553 if (ir->ePBC != epbcNONE)
2560 /* We should have a vsite free routine, but here we can simply free */
2563 if (fplog && (ft2b >= 0 || ftmb >= 0))
2566 "Initial maximum inter charge-group distances:\n");
2570 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2571 *r_2b, interaction_function[ft2b].longname,
2572 a_2b_1+1, a_2b_2+1);
2577 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2578 *r_mb, interaction_function[ftmb].longname,
2579 a_mb_1+1, a_mb_2+1);