2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2011,2012,2013,2014,2015,2016, 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
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/gmxlib/chargegroup.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/force.h"
64 #include "gromacs/mdlib/forcerec.h"
65 #include "gromacs/mdlib/gmx_omp_nthreads.h"
66 #include "gromacs/mdlib/vsite.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/topology/topsort.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/exceptions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/gmxassert.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/stringutil.h"
82 #include "domdec_constraints.h"
83 #include "domdec_internal.h"
84 #include "domdec_vsite.h"
86 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
87 #define NITEM_DD_INIT_LOCAL_STATE 5
90 int *index; /* Index for each atom into il */
91 int *il; /* ftype|type|a0|...|an|ftype|... */
101 /*! \brief Struct for thread local work data for local topology generation */
103 t_idef idef; /**< Partial local topology */
104 int **vsite_pbc; /**< vsite PBC structure */
105 int *vsite_pbc_nalloc; /**< Allocation sizes for vsite_pbc */
106 int nbonded; /**< The number of bondeds in this struct */
107 t_blocka excl; /**< List of exclusions */
108 int excl_count; /**< The total exclusion count for \p excl */
111 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
112 struct gmx_reverse_top_t
114 //! @cond Doxygen_Suppress
115 gmx_bool bExclRequired; /**< Do we require all exclusions to be assigned? */
116 int n_excl_at_max; /**< The maximum number of exclusions one atom can have */
117 gmx_bool bConstr; /**< Are there constraints in this revserse top? */
118 gmx_bool bSettle; /**< Are there settles in this revserse top? */
119 gmx_bool bBCheck; /**< All bonded interactions have to be assigned? */
120 gmx_bool bInterCGInteractions; /**< Are there bondeds/exclusions between charge-groups? */
121 reverse_ilist_t *ril_mt; /**< Reverse ilist for all moltypes */
122 int ril_mt_tot_size; /**< The size of ril_mt[?].index summed over all entries */
123 int ilsort; /**< The sorting state of bondeds for free energy */
124 molblock_ind_t *mbi; /**< molblock to global atom index for quick lookup of molblocks on atom index */
125 int nmolblock; /**< The number of molblocks, size of \p mbi */
127 gmx_bool bIntermolecularInteractions; /**< Do we have intermolecular interactions? */
128 reverse_ilist_t ril_intermol; /**< Intermolecular reverse ilist */
130 /* Work data structures for multi-threading */
131 int nthread; /**< The number of threads to be used */
132 thread_work_t *th_work; /**< Thread work array for local topology generation */
136 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
137 static int nral_rt(int ftype)
142 if (interaction_function[ftype].flags & IF_VSITE)
144 /* With vsites the reverse topology contains
145 * two extra entries for PBC.
153 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
154 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
155 gmx_bool bConstr, gmx_bool bSettle)
157 return (((interaction_function[ftype].flags & IF_BOND) &&
158 !(interaction_function[ftype].flags & IF_VSITE) &&
159 (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))) ||
160 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
161 (bSettle && ftype == F_SETTLE));
164 /*! \brief Print a header on error messages */
165 static void print_error_header(FILE *fplog, const char *moltypename, int nprint)
167 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
168 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
170 "the first %d missing interactions, except for exclusions:\n",
173 "the first %d missing interactions, except for exclusions:\n",
177 /*! \brief Help print error output when interactions are missing */
178 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
179 const gmx_reverse_top_t *rt,
180 const char *moltypename,
181 const reverse_ilist_t *ril,
182 int a_start, int a_end,
183 int nat_mol, int nmol,
187 int nril_mol = ril->index[nat_mol];
188 snew(assigned, nmol*nril_mol);
190 int *gatindex = cr->dd->gatindex;
191 for (int ftype = 0; ftype < F_NRE; ftype++)
193 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
195 int nral = NRAL(ftype);
196 const t_ilist *il = &idef->il[ftype];
197 const t_iatom *ia = il->iatoms;
198 for (int i = 0; i < il->nr; i += 1+nral)
200 int a0 = gatindex[ia[1]];
201 /* Check if this interaction is in
202 * the currently checked molblock.
204 if (a0 >= a_start && a0 < a_end)
206 int mol = (a0 - a_start)/nat_mol;
207 int a0_mol = (a0 - a_start) - mol*nat_mol;
208 int j_mol = ril->index[a0_mol];
210 while (j_mol < ril->index[a0_mol+1] && !found)
212 int j = mol*nril_mol + j_mol;
213 int ftype_j = ril->il[j_mol];
214 /* Here we need to check if this interaction has
215 * not already been assigned, since we could have
216 * multiply defined interactions.
218 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
221 /* Check the atoms */
223 for (int a = 0; a < nral; a++)
225 if (gatindex[ia[1+a]] !=
226 a_start + mol*nat_mol + ril->il[j_mol+2+a])
236 j_mol += 2 + nral_rt(ftype_j);
240 gmx_incons("Some interactions seem to be assigned multiple times");
248 gmx_sumi(nmol*nril_mol, assigned, cr);
252 for (int mol = 0; mol < nmol; mol++)
255 while (j_mol < nril_mol)
257 int ftype = ril->il[j_mol];
258 int nral = NRAL(ftype);
259 int j = mol*nril_mol + j_mol;
260 if (assigned[j] == 0 &&
261 !(interaction_function[ftype].flags & IF_VSITE))
263 if (DDMASTER(cr->dd))
267 print_error_header(fplog, moltypename, nprint);
269 fprintf(fplog, "%20s atoms",
270 interaction_function[ftype].longname);
271 fprintf(stderr, "%20s atoms",
272 interaction_function[ftype].longname);
274 for (a = 0; a < nral; a++)
276 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
277 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
282 fprintf(stderr, " ");
285 fprintf(fplog, " global");
286 fprintf(stderr, " global");
287 for (a = 0; a < nral; a++)
289 fprintf(fplog, "%6d",
290 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
291 fprintf(stderr, "%6d",
292 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
294 fprintf(fplog, "\n");
295 fprintf(stderr, "\n");
303 j_mol += 2 + nral_rt(ftype);
310 /*! \brief Help print error output when interactions are missing */
311 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
312 const gmx_mtop_t *mtop,
315 int mb, a_start, a_end;
316 const gmx_molblock_t *molb;
317 const gmx_reverse_top_t *rt;
319 rt = cr->dd->reverse_top;
321 /* Print the atoms in the missing interactions per molblock */
323 for (mb = 0; mb < mtop->nmolblock; mb++)
325 molb = &mtop->molblock[mb];
327 a_end = a_start + molb->nmol*molb->natoms_mol;
329 print_missing_interactions_mb(fplog, cr, rt,
330 *(mtop->moltype[molb->type].name),
331 &rt->ril_mt[molb->type],
332 a_start, a_end, molb->natoms_mol,
338 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
340 const gmx_mtop_t *top_global,
341 const gmx_localtop_t *top_local,
342 t_state *state_local)
344 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
353 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
357 ndiff_tot = local_count - dd->nbonded_global;
359 for (ftype = 0; ftype < F_NRE; ftype++)
362 cl[ftype] = top_local->idef.il[ftype].nr/(1+nral);
365 gmx_sumi(F_NRE, cl, cr);
371 fprintf(fplog, "\nA list of missing interactions:\n");
373 fprintf(stderr, "\nA list of missing interactions:\n");
374 rest_global = dd->nbonded_global;
375 rest_local = local_count;
376 for (ftype = 0; ftype < F_NRE; ftype++)
378 /* In the reverse and local top all constraints are merged
379 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
380 * and add these constraints when doing F_CONSTR.
382 if (((interaction_function[ftype].flags & IF_BOND) &&
383 (dd->reverse_top->bBCheck
384 || !(interaction_function[ftype].flags & IF_LIMZERO)))
385 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
386 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
388 n = gmx_mtop_ftype_count(top_global, ftype);
389 if (ftype == F_CONSTR)
391 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
393 ndiff = cl[ftype] - n;
396 sprintf(buf, "%20s of %6d missing %6d",
397 interaction_function[ftype].longname, n, -ndiff);
400 fprintf(fplog, "%s\n", buf);
402 fprintf(stderr, "%s\n", buf);
405 rest_local -= cl[ftype];
409 ndiff = rest_local - rest_global;
412 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
413 rest_global, -ndiff);
416 fprintf(fplog, "%s\n", buf);
418 fprintf(stderr, "%s\n", buf);
422 print_missing_interactions_atoms(fplog, cr, top_global, &top_local->idef);
423 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
424 -1, as_rvec_array(state_local->x.data()), state_local->box);
426 std::string errorMessage;
430 errorMessage = "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
434 errorMessage = gmx::formatString("%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(dd), dd_cutoff_twobody(dd));
436 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), errorMessage.c_str());
439 /*! \brief Return global topology molecule information for global atom index \p i_gl */
440 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
441 int *mb, int *mt, int *mol, int *i_mol)
443 molblock_ind_t *mbi = rt->mbi;
445 int end = rt->nmolblock; /* exclusive */
448 /* binary search for molblock_ind */
452 if (i_gl >= mbi[mid].a_end)
456 else if (i_gl < mbi[mid].a_start)
470 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
471 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
474 /*! \brief Count the exclusions for all atoms in \p cgs */
475 static void count_excls(const t_block *cgs, const t_blocka *excls,
476 int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
478 int cg, at0, at1, at, excl, atj;
483 for (cg = 0; cg < cgs->nr; cg++)
485 at0 = cgs->index[cg];
486 at1 = cgs->index[cg+1];
487 for (at = at0; at < at1; at++)
489 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
491 atj = excls->a[excl];
495 if (atj < at0 || atj >= at1)
502 *n_excl_at_max = std::max(*n_excl_at_max,
503 excls->index[at+1] - excls->index[at]);
508 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
509 static int low_make_reverse_ilist(const t_ilist *il_mt, const t_atom *atom,
510 const int * const * vsite_pbc,
512 gmx_bool bConstr, gmx_bool bSettle,
514 int *r_index, int *r_il,
515 gmx_bool bLinkToAllAtoms,
518 int ftype, nral, i, j, nlink, link;
526 for (ftype = 0; ftype < F_NRE; ftype++)
528 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
529 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
530 (bSettle && ftype == F_SETTLE))
532 bVSite = (interaction_function[ftype].flags & IF_VSITE);
535 for (i = 0; i < il->nr; i += 1+nral)
542 /* We don't need the virtual sites for the cg-links */
552 /* Couple to the first atom in the interaction */
555 for (link = 0; link < nlink; link++)
560 assert(r_il); //with bAssign not allowed to be null
562 r_il[r_index[a]+count[a]] =
563 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
564 r_il[r_index[a]+count[a]+1] = ia[0];
565 for (j = 1; j < 1+nral; j++)
567 /* Store the molecular atom number */
568 r_il[r_index[a]+count[a]+1+j] = ia[j];
571 if (interaction_function[ftype].flags & IF_VSITE)
575 /* Add an entry to iatoms for storing
576 * which of the constructing atoms are
579 r_il[r_index[a]+count[a]+2+nral] = 0;
580 for (j = 2; j < 1+nral; j++)
582 if (atom[ia[j]].ptype == eptVSite)
584 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
587 /* Store vsite pbc atom in a second extra entry */
588 r_il[r_index[a]+count[a]+2+nral+1] =
589 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
594 /* We do not count vsites since they are always
595 * uniquely assigned and can be assigned
596 * to multiple nodes with recursive vsites.
599 !(interaction_function[ftype].flags & IF_LIMZERO))
604 count[a] += 2 + nral_rt(ftype);
613 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
614 static int make_reverse_ilist(const t_ilist *ilist,
615 const t_atoms *atoms,
616 const int * const * vsite_pbc,
617 gmx_bool bConstr, gmx_bool bSettle,
619 gmx_bool bLinkToAllAtoms,
620 reverse_ilist_t *ril_mt)
622 int nat_mt, *count, i, nint_mt;
624 /* Count the interactions */
627 low_make_reverse_ilist(ilist, atoms->atom, vsite_pbc,
629 bConstr, bSettle, bBCheck, NULL, NULL,
630 bLinkToAllAtoms, FALSE);
632 snew(ril_mt->index, nat_mt+1);
633 ril_mt->index[0] = 0;
634 for (i = 0; i < nat_mt; i++)
636 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
639 snew(ril_mt->il, ril_mt->index[nat_mt]);
641 /* Store the interactions */
643 low_make_reverse_ilist(ilist, atoms->atom, vsite_pbc,
645 bConstr, bSettle, bBCheck,
646 ril_mt->index, ril_mt->il,
647 bLinkToAllAtoms, TRUE);
654 /*! \brief Destroys a reverse_ilist_t struct */
655 static void destroy_reverse_ilist(reverse_ilist_t *ril)
661 /*! \brief Generate the reverse topology */
662 static gmx_reverse_top_t *make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
663 const int * const * const * vsite_pbc_molt,
664 gmx_bool bConstr, gmx_bool bSettle,
665 gmx_bool bBCheck, int *nint)
668 gmx_reverse_top_t *rt;
675 /* Should we include constraints (for SHAKE) in rt? */
676 rt->bConstr = bConstr;
677 rt->bSettle = bSettle;
678 rt->bBCheck = bBCheck;
680 rt->bInterCGInteractions = mtop->bIntermolecularInteractions;
681 snew(nint_mt, mtop->nmoltype);
682 snew(rt->ril_mt, mtop->nmoltype);
683 rt->ril_mt_tot_size = 0;
684 for (mt = 0; mt < mtop->nmoltype; mt++)
686 molt = &mtop->moltype[mt];
687 if (molt->cgs.nr > 1)
689 rt->bInterCGInteractions = TRUE;
692 /* Make the atom to interaction list for this molecule type */
694 make_reverse_ilist(molt->ilist, &molt->atoms,
695 vsite_pbc_molt ? vsite_pbc_molt[mt] : NULL,
696 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
699 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt->atoms.nr];
703 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
707 for (mb = 0; mb < mtop->nmolblock; mb++)
709 *nint += mtop->molblock[mb].nmol*nint_mt[mtop->molblock[mb].type];
713 /* Make an intermolecular reverse top, if necessary */
714 rt->bIntermolecularInteractions = mtop->bIntermolecularInteractions;
715 if (rt->bIntermolecularInteractions)
717 t_atoms atoms_global;
719 rt->ril_intermol.index = NULL;
720 rt->ril_intermol.il = NULL;
722 atoms_global.nr = mtop->natoms;
723 atoms_global.atom = NULL; /* Only used with virtual sites */
726 make_reverse_ilist(mtop->intermolecular_ilist, &atoms_global,
728 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
732 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
734 rt->ilsort = ilsortFE_UNSORTED;
738 rt->ilsort = ilsortNO_FE;
741 /* Make a molblock index for fast searching */
742 snew(rt->mbi, mtop->nmolblock);
743 rt->nmolblock = mtop->nmolblock;
745 for (mb = 0; mb < mtop->nmolblock; mb++)
747 rt->mbi[mb].a_start = i;
748 i += mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
749 rt->mbi[mb].a_end = i;
750 rt->mbi[mb].natoms_mol = mtop->molblock[mb].natoms_mol;
751 rt->mbi[mb].type = mtop->molblock[mb].type;
754 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
755 snew(rt->th_work, rt->nthread);
756 if (vsite_pbc_molt != NULL)
758 for (thread = 0; thread < rt->nthread; thread++)
760 snew(rt->th_work[thread].vsite_pbc, F_VSITEN-F_VSITE2+1);
761 snew(rt->th_work[thread].vsite_pbc_nalloc, F_VSITEN-F_VSITE2+1);
768 void dd_make_reverse_top(FILE *fplog,
769 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
770 const gmx_vsite_t *vsite,
771 const t_inputrec *ir, gmx_bool bBCheck)
775 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
778 /* If normal and/or settle constraints act only within charge groups,
779 * we can store them in the reverse top and simply assign them to domains.
780 * Otherwise we need to assign them to multiple domains and set up
781 * the parallel version constraint algorithm(s).
784 dd->reverse_top = make_reverse_top(mtop, ir->efep != efepNO,
785 vsite ? vsite->vsite_pbc_molt : NULL,
786 !dd->bInterCGcons, !dd->bInterCGsettles,
787 bBCheck, &dd->nbonded_global);
789 gmx_reverse_top_t *rt = dd->reverse_top;
791 if (rt->ril_mt_tot_size >= 200000 &&
793 mtop->nmolblock == 1 && mtop->molblock[0].nmol == 1)
795 /* mtop comes from a pre Gromacs 4 tpr file */
796 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";
799 fprintf(fplog, "\n%s\n\n", note);
803 fprintf(stderr, "\n%s\n\n", note);
807 /* With the Verlet scheme, exclusions are handled in the non-bonded
808 * kernels and only exclusions inside the cut-off lead to exclusion
809 * forces. Since each atom pair is treated at most once in the non-bonded
810 * kernels, it doesn't matter if the exclusions for the same atom pair
811 * appear multiple times in the exclusion list. In contrast, the, old,
812 * group cut-off scheme loops over a list of exclusions, so there each
813 * excluded pair should appear exactly once.
815 rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
816 inputrecExclForces(ir));
821 dd->n_intercg_excl = 0;
822 rt->n_excl_at_max = 0;
823 for (mb = 0; mb < mtop->nmolblock; mb++)
825 gmx_molblock_t *molb;
827 int n_excl_mol, n_excl_icg, n_excl_at_max;
829 molb = &mtop->molblock[mb];
830 molt = &mtop->moltype[molb->type];
831 count_excls(&molt->cgs, &molt->excls,
832 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
833 nexcl += molb->nmol*n_excl_mol;
834 dd->n_intercg_excl += molb->nmol*n_excl_icg;
835 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
837 if (rt->bExclRequired)
839 dd->nbonded_global += nexcl;
840 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
842 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
843 "will use an extra communication step for exclusion forces for %s\n",
844 dd->n_intercg_excl, eel_names[ir->coulombtype]);
848 if (vsite && vsite->n_intercg_vsite > 0)
852 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
853 "will an extra communication step for selected coordinates and forces\n",
854 vsite->n_intercg_vsite);
856 init_domdec_vsites(dd, vsite->n_intercg_vsite);
859 if (dd->bInterCGcons || dd->bInterCGsettles)
861 init_domdec_constraints(dd, mtop);
865 fprintf(fplog, "\n");
869 /*! \brief Store a vsite interaction at the end of \p il
871 * This routine is very similar to add_ifunc, but vsites interactions
872 * have more work to do than other kinds of interactions, and the
873 * complex way nral (and thus vector contents) depends on ftype
874 * confuses static analysis tools unless we fuse the vsite
875 * atom-indexing organization code with the ifunc-adding code, so that
876 * they can see that nral is the same value. */
877 static gmx_inline void
878 add_ifunc_for_vsites(t_iatom *tiatoms, gmx_ga2la_t *ga2la,
879 int nral, gmx_bool bHomeA,
880 int a, int a_gl, int a_mol,
881 const t_iatom *iatoms,
886 if (il->nr+1+nral > il->nalloc)
888 il->nalloc = over_alloc_large(il->nr+1+nral);
889 srenew(il->iatoms, il->nalloc);
891 liatoms = il->iatoms + il->nr;
895 tiatoms[0] = iatoms[0];
899 /* We know the local index of the first atom */
904 /* Convert later in make_local_vsites */
905 tiatoms[1] = -a_gl - 1;
908 for (int k = 2; k < 1+nral; k++)
910 int ak_gl = a_gl + iatoms[k] - a_mol;
911 if (!ga2la_get_home(ga2la, ak_gl, &tiatoms[k]))
913 /* Copy the global index, convert later in make_local_vsites */
914 tiatoms[k] = -(ak_gl + 1);
916 // Note that ga2la_get_home always sets the third parameter if
919 for (int k = 0; k < 1+nral; k++)
921 liatoms[k] = tiatoms[k];
925 /*! \brief Store a bonded interaction at the end of \p il */
926 static gmx_inline void add_ifunc(int nral, t_iatom *tiatoms, t_ilist *il)
931 if (il->nr+1+nral > il->nalloc)
933 il->nalloc = over_alloc_large(il->nr+1+nral);
934 srenew(il->iatoms, il->nalloc);
936 liatoms = il->iatoms + il->nr;
937 for (k = 0; k <= nral; k++)
939 liatoms[k] = tiatoms[k];
944 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
945 static void add_posres(int mol, int a_mol, const gmx_molblock_t *molb,
946 t_iatom *iatoms, const t_iparams *ip_in,
952 /* This position restraint has not been added yet,
953 * so it's index is the current number of position restraints.
955 n = idef->il[F_POSRES].nr/2;
956 if (n+1 > idef->iparams_posres_nalloc)
958 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
959 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
961 ip = &idef->iparams_posres[n];
962 /* Copy the force constants */
963 *ip = ip_in[iatoms[0]];
965 /* Get the position restraint coordinates from the molblock */
966 a_molb = mol*molb->natoms_mol + a_mol;
967 if (a_molb >= molb->nposres_xA)
969 gmx_incons("Not enough position restraint coordinates");
971 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
972 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
973 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
974 if (molb->nposres_xB > 0)
976 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
977 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
978 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
982 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
983 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
984 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
986 /* Set the parameter index for idef->iparams_posre */
990 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
991 static void add_fbposres(int mol, int a_mol, const gmx_molblock_t *molb,
992 t_iatom *iatoms, const t_iparams *ip_in,
998 /* This flat-bottom position restraint has not been added yet,
999 * so it's index is the current number of position restraints.
1001 n = idef->il[F_FBPOSRES].nr/2;
1002 if (n+1 > idef->iparams_fbposres_nalloc)
1004 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
1005 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
1007 ip = &idef->iparams_fbposres[n];
1008 /* Copy the force constants */
1009 *ip = ip_in[iatoms[0]];
1011 /* Get the position restriant coordinats from the molblock */
1012 a_molb = mol*molb->natoms_mol + a_mol;
1013 if (a_molb >= molb->nposres_xA)
1015 gmx_incons("Not enough position restraint coordinates");
1017 /* Take reference positions from A position of normal posres */
1018 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
1019 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
1020 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
1022 /* Note: no B-type for flat-bottom posres */
1024 /* Set the parameter index for idef->iparams_posre */
1028 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
1029 static void add_vsite(gmx_ga2la_t *ga2la, const int *index, const int *rtil,
1030 int ftype, int nral,
1031 gmx_bool bHomeA, int a, int a_gl, int a_mol,
1032 const t_iatom *iatoms,
1033 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
1035 int k, vsi, pbc_a_mol;
1036 t_iatom tiatoms[1+MAXATOMLIST];
1037 int j, ftype_r, nral_r;
1039 /* Add this interaction to the local topology */
1040 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1044 vsi = idef->il[ftype].nr/(1+nral) - 1;
1045 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
1047 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
1048 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
1052 pbc_a_mol = iatoms[1+nral+1];
1055 /* The pbc flag is one of the following two options:
1056 * -2: vsite and all constructing atoms are within the same cg, no pbc
1057 * -1: vsite and its first constructing atom are in the same cg, do pbc
1059 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
1063 /* Set the pbc atom for this vsite so we can make its pbc
1064 * identical to the rest of the atoms in its charge group.
1065 * Since the order of the atoms does not change within a charge
1066 * group, we do not need the global to local atom index.
1068 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
1073 /* This vsite is non-home (required for recursion),
1074 * and therefore there is no charge group to match pbc with.
1075 * But we always turn on full_pbc to assure that higher order
1076 * recursion works correctly.
1078 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
1084 /* Check for recursion */
1085 for (k = 2; k < 1+nral; k++)
1087 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1089 /* This construction atoms is a vsite and not a home atom */
1092 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1094 /* Find the vsite construction */
1096 /* Check all interactions assigned to this atom */
1097 j = index[iatoms[k]];
1098 while (j < index[iatoms[k]+1])
1100 ftype_r = rtil[j++];
1101 nral_r = NRAL(ftype_r);
1102 if (interaction_function[ftype_r].flags & IF_VSITE)
1104 /* Add this vsite (recursion) */
1105 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1106 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1107 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
1108 j += 1 + nral_r + 2;
1120 /*! \brief Update the local atom to local charge group index */
1121 static void make_la2lc(gmx_domdec_t *dd)
1123 int *cgindex, *la2lc, cg, a;
1125 cgindex = dd->cgindex;
1127 if (dd->nat_tot > dd->la2lc_nalloc)
1129 dd->la2lc_nalloc = over_alloc_dd(dd->nat_tot);
1130 snew(dd->la2lc, dd->la2lc_nalloc);
1134 /* Make the local atom to local cg index */
1135 for (cg = 0; cg < dd->ncg_tot; cg++)
1137 for (a = cgindex[cg]; a < cgindex[cg+1]; a++)
1144 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1145 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1151 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1155 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1161 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1162 static void combine_blocka(t_blocka *dest, const thread_work_t *src, int nsrc)
1166 ni = src[nsrc-1].excl.nr;
1168 for (s = 0; s < nsrc; s++)
1170 na += src[s].excl.nra;
1172 if (ni + 1 > dest->nalloc_index)
1174 dest->nalloc_index = over_alloc_large(ni+1);
1175 srenew(dest->index, dest->nalloc_index);
1177 if (dest->nra + na > dest->nalloc_a)
1179 dest->nalloc_a = over_alloc_large(dest->nra+na);
1180 srenew(dest->a, dest->nalloc_a);
1182 for (s = 1; s < nsrc; s++)
1184 for (i = dest->nr+1; i < src[s].excl.nr+1; i++)
1186 dest->index[i] = dest->nra + src[s].excl.index[i];
1188 for (i = 0; i < src[s].excl.nra; i++)
1190 dest->a[dest->nra+i] = src[s].excl.a[i];
1192 dest->nr = src[s].excl.nr;
1193 dest->nra += src[s].excl.nra;
1197 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1198 * virtual sites need special attention, as pbc info differs per vsite.
1200 static void combine_idef(t_idef *dest, const thread_work_t *src, int nsrc,
1205 for (ftype = 0; ftype < F_NRE; ftype++)
1210 for (s = 1; s < nsrc; s++)
1212 n += src[s].idef.il[ftype].nr;
1218 ild = &dest->il[ftype];
1220 if (ild->nr + n > ild->nalloc)
1222 ild->nalloc = over_alloc_large(ild->nr+n);
1223 srenew(ild->iatoms, ild->nalloc);
1227 int nral1 = 0, ftv = 0;
1229 vpbc = ((interaction_function[ftype].flags & IF_VSITE) &&
1230 vsite->vsite_pbc_loc != NULL);
1233 nral1 = 1 + NRAL(ftype);
1234 ftv = ftype - F_VSITE2;
1235 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1237 vsite->vsite_pbc_loc_nalloc[ftv] =
1238 over_alloc_large((ild->nr + n)/nral1);
1239 srenew(vsite->vsite_pbc_loc[ftv],
1240 vsite->vsite_pbc_loc_nalloc[ftv]);
1244 for (s = 1; s < nsrc; s++)
1249 ils = &src[s].idef.il[ftype];
1250 for (i = 0; i < ils->nr; i++)
1252 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1256 for (i = 0; i < ils->nr; i += nral1)
1258 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1259 src[s].vsite_pbc[ftv][i/nral1];
1266 /* Position restraints need an additional treatment */
1267 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1269 int nposres = dest->il[ftype].nr/2;
1270 if (nposres > dest->iparams_posres_nalloc)
1272 dest->iparams_posres_nalloc = over_alloc_large(nposres);
1273 srenew(dest->iparams_posres, dest->iparams_posres_nalloc);
1275 /* Set nposres to the number of original position restraints in dest */
1276 for (int s = 1; s < nsrc; s++)
1278 nposres -= src[s].idef.il[ftype].nr/2;
1280 for (int s = 1; s < nsrc; s++)
1282 if (ftype == F_POSRES)
1284 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1286 /* Correct the index into iparams_posres */
1287 dest->il[ftype].iatoms[nposres*2] = nposres;
1288 /* Copy the position restraint force parameters */
1289 dest->iparams_posres[nposres] = src[s].idef.iparams_posres[i];
1295 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1297 /* Correct the index into iparams_fbposres */
1298 dest->il[ftype].iatoms[nposres*2] = nposres;
1299 /* Copy the position restraint force parameters */
1300 dest->iparams_fbposres[nposres] = src[s].idef.iparams_fbposres[i];
1310 /*! \brief Check and when available assign bonded interactions for local atom i
1312 static gmx_inline void
1313 check_assign_interactions_atom(int i, int i_gl,
1315 const int *index, const int *rtil,
1316 gmx_bool bInterMolInteractions,
1317 int ind_start, int ind_end,
1318 const gmx_domdec_t *dd,
1319 const gmx_domdec_zones_t *zones,
1320 const gmx_molblock_t *molb,
1321 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1326 const t_iparams *ip_in,
1328 int **vsite_pbc, int *vsite_pbc_nalloc,
1339 const t_iatom *iatoms;
1341 t_iatom tiatoms[1 + MAXATOMLIST];
1346 if (ftype == F_SETTLE)
1348 /* Settles are only in the reverse top when they
1349 * operate within a charge group. So we can assign
1350 * them without checks. We do this only for performance
1351 * reasons; it could be handled by the code below.
1355 /* Home zone: add this settle to the local topology */
1356 tiatoms[0] = iatoms[0];
1358 tiatoms[2] = i + iatoms[2] - iatoms[1];
1359 tiatoms[3] = i + iatoms[3] - iatoms[1];
1360 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1365 else if (interaction_function[ftype].flags & IF_VSITE)
1367 assert(!bInterMolInteractions);
1368 /* The vsite construction goes where the vsite itself is */
1371 add_vsite(dd->ga2la, index, rtil, ftype, nral,
1372 TRUE, i, i_gl, i_mol,
1373 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1382 tiatoms[0] = iatoms[0];
1386 assert(!bInterMolInteractions);
1387 /* Assign single-body interactions to the home zone */
1392 if (ftype == F_POSRES)
1394 add_posres(mol, i_mol, molb, tiatoms, ip_in,
1397 else if (ftype == F_FBPOSRES)
1399 add_fbposres(mol, i_mol, molb, tiatoms, ip_in,
1410 /* This is a two-body interaction, we can assign
1411 * analogous to the non-bonded assignments.
1413 int k_gl, a_loc, kz;
1415 if (!bInterMolInteractions)
1417 /* Get the global index using the offset in the molecule */
1418 k_gl = i_gl + iatoms[2] - i_mol;
1424 if (!ga2la_get(dd->ga2la, k_gl, &a_loc, &kz))
1434 /* Check zone interaction assignments */
1435 bUse = ((iz < zones->nizone &&
1437 kz >= zones->izone[iz].j0 &&
1438 kz < zones->izone[iz].j1) ||
1439 (kz < zones->nizone &&
1441 iz >= zones->izone[kz].j0 &&
1442 iz < zones->izone[kz].j1));
1447 /* If necessary check the cgcm distance */
1449 dd_dist2(pbc_null, cg_cm, la2lc,
1450 tiatoms[1], tiatoms[2]) >= rc2)
1459 /* Assign this multi-body bonded interaction to
1460 * the local node if we have all the atoms involved
1461 * (local or communicated) and the minimum zone shift
1462 * in each dimension is zero, for dimensions
1463 * with 2 DD cells an extra check may be necessary.
1465 ivec k_zero, k_plus;
1471 for (k = 1; k <= nral && bUse; k++)
1477 if (!bInterMolInteractions)
1479 /* Get the global index using the offset in the molecule */
1480 k_gl = i_gl + iatoms[k] - i_mol;
1486 bLocal = ga2la_get(dd->ga2la, k_gl, &a_loc, &kz);
1487 if (!bLocal || kz >= zones->n)
1489 /* We do not have this atom of this interaction
1490 * locally, or it comes from more than one cell
1500 for (d = 0; d < DIM; d++)
1502 if (zones->shift[kz][d] == 0)
1514 k_zero[XX] && k_zero[YY] && k_zero[ZZ]);
1519 for (d = 0; (d < DIM && bUse); d++)
1521 /* Check if the cg_cm distance falls within
1522 * the cut-off to avoid possible multiple
1523 * assignments of bonded interactions.
1527 dd_dist2(pbc_null, cg_cm, la2lc,
1528 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1537 /* Add this interaction to the local topology */
1538 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1539 /* Sum so we can check in global_stat
1540 * if we have everything.
1543 !(interaction_function[ftype].flags & IF_LIMZERO))
1553 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1555 * With thread parallelizing each thread acts on a different atom range:
1556 * at_start to at_end.
1558 static int make_bondeds_zone(gmx_domdec_t *dd,
1559 const gmx_domdec_zones_t *zones,
1560 const gmx_molblock_t *molb,
1561 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1563 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1564 const t_iparams *ip_in,
1567 int *vsite_pbc_nalloc,
1569 int at_start, int at_end)
1571 int i, i_gl, mb, mt, mol, i_mol;
1574 gmx_reverse_top_t *rt;
1577 rt = dd->reverse_top;
1579 bBCheck = rt->bBCheck;
1583 for (i = at_start; i < at_end; i++)
1585 /* Get the global atom number */
1586 i_gl = dd->gatindex[i];
1587 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1588 /* Check all intramolecular interactions assigned to this atom */
1589 index = rt->ril_mt[mt].index;
1590 rtil = rt->ril_mt[mt].il;
1592 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1594 index[i_mol], index[i_mol+1],
1597 bRCheckMB, rcheck, bRCheck2B, rc2,
1602 idef, vsite_pbc, vsite_pbc_nalloc,
1608 if (rt->bIntermolecularInteractions)
1610 /* Check all intermolecular interactions assigned to this atom */
1611 index = rt->ril_intermol.index;
1612 rtil = rt->ril_intermol.il;
1614 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1616 index[i_gl], index[i_gl + 1],
1619 bRCheckMB, rcheck, bRCheck2B, rc2,
1624 idef, vsite_pbc, vsite_pbc_nalloc,
1631 return nbonded_local;
1634 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1635 static void set_no_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1636 int iz, t_blocka *lexcls)
1640 a0 = dd->cgindex[zones->cg_range[iz]];
1641 a1 = dd->cgindex[zones->cg_range[iz+1]];
1643 for (a = a0+1; a < a1+1; a++)
1645 lexcls->index[a] = lexcls->nra;
1649 /*! \brief Set the exclusion data for i-zone \p iz
1651 * This is a legacy version for the group scheme of the same routine below.
1652 * Here charge groups and distance checks to ensure unique exclusions
1655 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1656 const gmx_moltype_t *moltype,
1657 gmx_bool bRCheck, real rc2,
1658 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1662 int cg_start, int cg_end)
1664 int n_excl_at_max, n, count, jla0, jla1, jla;
1665 int cg, la0, la1, la, a_gl, mb, mt, mol, a_mol, j, aj_mol;
1666 const t_blocka *excls;
1672 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1673 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1675 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1677 /* We set the end index, but note that we might not start at zero here */
1678 lexcls->nr = dd->cgindex[cg_end];
1682 for (cg = cg_start; cg < cg_end; cg++)
1684 if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1686 lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1687 srenew(lexcls->a, lexcls->nalloc_a);
1689 la0 = dd->cgindex[cg];
1690 la1 = dd->cgindex[cg+1];
1691 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1692 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1694 /* Copy the exclusions from the global top */
1695 for (la = la0; la < la1; la++)
1697 lexcls->index[la] = n;
1698 a_gl = dd->gatindex[la];
1699 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1700 excls = &moltype[mt].excls;
1701 for (j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1703 aj_mol = excls->a[j];
1704 /* This computation of jla is only correct intra-cg */
1705 jla = la + aj_mol - a_mol;
1706 if (jla >= la0 && jla < la1)
1708 /* This is an intra-cg exclusion. We can skip
1709 * the global indexing and distance checking.
1711 /* Intra-cg exclusions are only required
1712 * for the home zone.
1716 lexcls->a[n++] = jla;
1717 /* Check to avoid double counts */
1726 /* This is a inter-cg exclusion */
1727 /* Since exclusions are pair interactions,
1728 * just like non-bonded interactions,
1729 * they can be assigned properly up
1730 * to the DD cutoff (not cutoff_min as
1731 * for the other bonded interactions).
1733 if (ga2la_get(ga2la, a_gl+aj_mol-a_mol, &jla, &cell))
1735 if (iz == 0 && cell == 0)
1737 lexcls->a[n++] = jla;
1738 /* Check to avoid double counts */
1744 else if (jla >= jla0 && jla < jla1 &&
1746 dd_dist2(pbc_null, cg_cm, la2lc, la, jla) < rc2))
1748 /* jla > la, since jla0 > la */
1749 lexcls->a[n++] = jla;
1759 /* There are no inter-cg excls and this cg is self-excluded.
1760 * These exclusions are only required for zone 0,
1761 * since other zones do not see themselves.
1765 for (la = la0; la < la1; la++)
1767 lexcls->index[la] = n;
1768 for (j = la0; j < la1; j++)
1773 count += ((la1 - la0)*(la1 - la0 - 1))/2;
1777 /* We don't need exclusions for this cg */
1778 for (la = la0; la < la1; la++)
1780 lexcls->index[la] = n;
1786 lexcls->index[lexcls->nr] = n;
1792 /*! \brief Set the exclusion data for i-zone \p iz */
1793 static void make_exclusions_zone(gmx_domdec_t *dd,
1794 gmx_domdec_zones_t *zones,
1795 const gmx_moltype_t *moltype,
1799 int at_start, int at_end)
1803 int n_excl_at_max, n, at;
1807 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1808 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1810 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1812 /* We set the end index, but note that we might not start at zero here */
1813 lexcls->nr = at_end;
1816 for (at = at_start; at < at_end; at++)
1818 if (n + 1000 > lexcls->nalloc_a)
1820 lexcls->nalloc_a = over_alloc_large(n + 1000);
1821 srenew(lexcls->a, lexcls->nalloc_a);
1823 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1825 int a_gl, mb, mt, mol, a_mol, j;
1826 const t_blocka *excls;
1828 if (n + n_excl_at_max > lexcls->nalloc_a)
1830 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1831 srenew(lexcls->a, lexcls->nalloc_a);
1834 /* Copy the exclusions from the global top */
1835 lexcls->index[at] = n;
1836 a_gl = dd->gatindex[at];
1837 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1838 &mb, &mt, &mol, &a_mol);
1839 excls = &moltype[mt].excls;
1840 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1842 int aj_mol, at_j, cell;
1844 aj_mol = excls->a[j];
1846 if (ga2la_get(ga2la, a_gl + aj_mol - a_mol, &at_j, &cell))
1848 /* This check is not necessary, but it can reduce
1849 * the number of exclusions in the list, which in turn
1850 * can speed up the pair list construction a bit.
1852 if (at_j >= jla0 && at_j < jla1)
1854 lexcls->a[n++] = at_j;
1861 /* We don't need exclusions for this atom */
1862 lexcls->index[at] = n;
1866 lexcls->index[lexcls->nr] = n;
1871 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1872 static void check_alloc_index(t_blocka *ba, int nindex_max)
1874 if (nindex_max+1 > ba->nalloc_index)
1876 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1877 srenew(ba->index, ba->nalloc_index);
1881 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1882 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1888 nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1890 check_alloc_index(lexcls, nr);
1892 for (thread = 1; thread < dd->reverse_top->nthread; thread++)
1894 check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1898 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1899 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1904 lexcls->nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1906 if (dd->n_intercg_excl == 0)
1908 /* There are no exclusions involving non-home charge groups,
1909 * but we need to set the indices for neighborsearching.
1911 la0 = dd->cgindex[zones->izone[0].cg1];
1912 for (la = la0; la < lexcls->nr; la++)
1914 lexcls->index[la] = lexcls->nra;
1917 /* nr is only used to loop over the exclusions for Ewald and RF,
1918 * so we can set it to the number of home atoms for efficiency.
1920 lexcls->nr = dd->cgindex[zones->izone[0].cg1];
1924 /*! \brief Clear a t_idef struct */
1925 static void clear_idef(t_idef *idef)
1929 /* Clear the counts */
1930 for (ftype = 0; ftype < F_NRE; ftype++)
1932 idef->il[ftype].nr = 0;
1936 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1937 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1938 gmx_domdec_zones_t *zones,
1939 const gmx_mtop_t *mtop,
1941 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1943 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1944 t_idef *idef, gmx_vsite_t *vsite,
1945 t_blocka *lexcls, int *excl_count)
1947 int nzone_bondeds, nzone_excl;
1948 int izone, cg0, cg1;
1952 gmx_reverse_top_t *rt;
1954 if (dd->reverse_top->bInterCGInteractions)
1956 nzone_bondeds = zones->n;
1960 /* Only single charge group (or atom) molecules, so interactions don't
1961 * cross zone boundaries and we only need to assign in the home zone.
1966 if (dd->n_intercg_excl > 0)
1968 /* We only use exclusions from i-zones to i- and j-zones */
1969 nzone_excl = zones->nizone;
1973 /* There are no inter-cg exclusions and only zone 0 sees itself */
1977 check_exclusions_alloc(dd, zones, lexcls);
1979 rt = dd->reverse_top;
1983 /* Clear the counts */
1991 for (izone = 0; izone < nzone_bondeds; izone++)
1993 cg0 = zones->cg_range[izone];
1994 cg1 = zones->cg_range[izone + 1];
1996 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1997 for (thread = 0; thread < rt->nthread; thread++)
2004 int *vsite_pbc_nalloc;
2007 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
2008 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
2016 idef_t = &rt->th_work[thread].idef;
2020 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
2024 vsite_pbc = vsite->vsite_pbc_loc;
2025 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
2029 vsite_pbc = rt->th_work[thread].vsite_pbc;
2030 vsite_pbc_nalloc = rt->th_work[thread].vsite_pbc_nalloc;
2036 vsite_pbc_nalloc = NULL;
2039 rt->th_work[thread].nbonded =
2040 make_bondeds_zone(dd, zones,
2042 bRCheckMB, rcheck, bRCheck2B, rc2,
2043 la2lc, pbc_null, cg_cm, idef->iparams,
2045 vsite_pbc, vsite_pbc_nalloc,
2047 dd->cgindex[cg0t], dd->cgindex[cg1t]);
2049 if (izone < nzone_excl)
2057 excl_t = &rt->th_work[thread].excl;
2062 if (dd->cgindex[dd->ncg_tot] == dd->ncg_tot &&
2065 /* No charge groups and no distance check required */
2066 make_exclusions_zone(dd, zones,
2067 mtop->moltype, cginfo,
2074 rt->th_work[thread].excl_count =
2075 make_exclusions_zone_cg(dd, zones,
2076 mtop->moltype, bRCheck2B, rc2,
2077 la2lc, pbc_null, cg_cm, cginfo,
2084 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2087 if (rt->nthread > 1)
2089 combine_idef(idef, rt->th_work, rt->nthread, vsite);
2092 for (thread = 0; thread < rt->nthread; thread++)
2094 nbonded_local += rt->th_work[thread].nbonded;
2097 if (izone < nzone_excl)
2099 if (rt->nthread > 1)
2101 combine_blocka(lexcls, rt->th_work, rt->nthread);
2104 for (thread = 0; thread < rt->nthread; thread++)
2106 *excl_count += rt->th_work[thread].excl_count;
2111 /* Some zones might not have exclusions, but some code still needs to
2112 * loop over the index, so we set the indices here.
2114 for (izone = nzone_excl; izone < zones->nizone; izone++)
2116 set_no_exclusions_zone(dd, zones, izone, lexcls);
2119 finish_local_exclusions(dd, zones, lexcls);
2122 fprintf(debug, "We have %d exclusions, check count %d\n",
2123 lexcls->nra, *excl_count);
2126 return nbonded_local;
2129 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
2131 lcgs->nr = dd->ncg_tot;
2132 lcgs->index = dd->cgindex;
2135 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
2136 int npbcdim, matrix box,
2137 rvec cellsize_min, ivec npulse,
2141 const gmx_mtop_t *mtop, gmx_localtop_t *ltop)
2143 gmx_bool bRCheckMB, bRCheck2B;
2147 t_pbc pbc, *pbc_null = NULL;
2151 fprintf(debug, "Making local topology\n");
2154 dd_make_local_cgs(dd, <op->cgs);
2159 if (dd->reverse_top->bInterCGInteractions)
2161 /* We need to check to which cell bondeds should be assigned */
2162 rc = dd_cutoff_twobody(dd);
2165 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2168 /* Should we check cg_cm distances when assigning bonded interactions? */
2169 for (d = 0; d < DIM; d++)
2172 /* Only need to check for dimensions where the part of the box
2173 * that is not communicated is smaller than the cut-off.
2175 if (d < npbcdim && dd->nc[d] > 1 &&
2176 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2183 /* Check for interactions between two atoms,
2184 * where we can allow interactions up to the cut-off,
2185 * instead of up to the smallest cell dimension.
2192 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
2193 d, cellsize_min[d], d, rcheck[d], bRCheck2B);
2196 if (bRCheckMB || bRCheck2B)
2201 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
2211 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
2212 bRCheckMB, rcheck, bRCheck2B, rc,
2214 pbc_null, cgcm_or_x,
2216 <op->excls, &nexcl);
2218 /* The ilist is not sorted yet,
2219 * we can only do this when we have the charge arrays.
2221 ltop->idef.ilsort = ilsortUNKNOWN;
2223 if (dd->reverse_top->bExclRequired)
2225 dd->nbonded_local += nexcl;
2228 ltop->atomtypes = mtop->atomtypes;
2231 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
2232 gmx_localtop_t *ltop)
2234 if (dd->reverse_top->ilsort == ilsortNO_FE)
2236 ltop->idef.ilsort = ilsortNO_FE;
2240 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
2244 gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global)
2246 gmx_localtop_t *top;
2251 top->idef.ntypes = top_global->ffparams.ntypes;
2252 top->idef.atnr = top_global->ffparams.atnr;
2253 top->idef.functype = top_global->ffparams.functype;
2254 top->idef.iparams = top_global->ffparams.iparams;
2255 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
2256 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
2258 for (i = 0; i < F_NRE; i++)
2260 top->idef.il[i].iatoms = NULL;
2261 top->idef.il[i].nalloc = 0;
2263 top->idef.ilsort = ilsortUNKNOWN;
2268 void dd_init_local_state(gmx_domdec_t *dd,
2269 t_state *state_global, t_state *state_local)
2271 int buf[NITEM_DD_INIT_LOCAL_STATE];
2275 buf[0] = state_global->flags;
2276 buf[1] = state_global->ngtc;
2277 buf[2] = state_global->nnhpres;
2278 buf[3] = state_global->nhchainlength;
2279 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
2281 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2283 init_state(state_local, 0, buf[1], buf[2], buf[3], buf[4]);
2284 state_local->flags = buf[0];
2287 /*! \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 */
2288 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2294 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2296 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2297 if (link->a[k] == cg_gl_j)
2304 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2305 "Inconsistent allocation of link");
2306 /* Add this charge group link */
2307 if (link->index[cg_gl+1]+1 > link->nalloc_a)
2309 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2310 srenew(link->a, link->nalloc_a);
2312 link->a[link->index[cg_gl+1]] = cg_gl_j;
2313 link->index[cg_gl+1]++;
2317 /*! \brief Allocate and return an array of the charge group index for all atoms */
2318 static int *make_at2cg(t_block *cgs)
2322 snew(at2cg, cgs->index[cgs->nr]);
2323 for (cg = 0; cg < cgs->nr; cg++)
2325 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2334 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
2335 cginfo_mb_t *cginfo_mb)
2337 gmx_bool bExclRequired;
2338 int mb, cg_offset, cg, cg_gl, a, aj, i, j, ftype, nral, nlink_mol, mol, ncgi;
2339 gmx_molblock_t *molb;
2340 gmx_moltype_t *molt;
2344 reverse_ilist_t ril, ril_intermol;
2346 cginfo_mb_t *cgi_mb;
2348 /* For each charge group make a list of other charge groups
2349 * in the system that a linked to it via bonded interactions
2350 * which are also stored in reverse_top.
2353 bExclRequired = dd->reverse_top->bExclRequired;
2355 if (mtop->bIntermolecularInteractions)
2357 if (ncg_mtop(mtop) < mtop->natoms)
2359 gmx_fatal(FARGS, "The combination of intermolecular interactions, charge groups and domain decomposition is not supported. Use cutoff-scheme=Verlet (which removes the charge groups) or run without domain decomposition.");
2364 atoms.nr = mtop->natoms;
2367 make_reverse_ilist(mtop->intermolecular_ilist, &atoms,
2368 NULL, FALSE, FALSE, FALSE, TRUE, &ril_intermol);
2372 snew(link->index, ncg_mtop(mtop)+1);
2379 for (mb = 0; mb < mtop->nmolblock; mb++)
2381 molb = &mtop->molblock[mb];
2382 if (molb->nmol == 0)
2386 molt = &mtop->moltype[molb->type];
2388 excls = &molt->excls;
2389 a2c = make_at2cg(cgs);
2390 /* Make a reverse ilist in which the interactions are linked
2391 * to all atoms, not only the first atom as in gmx_reverse_top.
2392 * The constraints are discarded here.
2394 make_reverse_ilist(molt->ilist, &molt->atoms,
2395 NULL, FALSE, FALSE, FALSE, TRUE, &ril);
2397 cgi_mb = &cginfo_mb[mb];
2399 for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb->nmol : 1); mol++)
2401 for (cg = 0; cg < cgs->nr; cg++)
2403 cg_gl = cg_offset + cg;
2404 link->index[cg_gl+1] = link->index[cg_gl];
2405 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2408 while (i < ril.index[a+1])
2410 ftype = ril.il[i++];
2412 /* Skip the ifunc index */
2414 for (j = 0; j < nral; j++)
2419 check_link(link, cg_gl, cg_offset+a2c[aj]);
2422 i += nral_rt(ftype);
2426 /* Exclusions always go both ways */
2427 for (j = excls->index[a]; j < excls->index[a+1]; j++)
2432 check_link(link, cg_gl, cg_offset+a2c[aj]);
2437 if (mtop->bIntermolecularInteractions)
2439 i = ril_intermol.index[a];
2440 while (i < ril_intermol.index[a+1])
2442 ftype = ril_intermol.il[i++];
2444 /* Skip the ifunc index */
2446 for (j = 0; j < nral; j++)
2448 /* Here we assume we have no charge groups;
2449 * this has been checked above.
2451 aj = ril_intermol.il[i+j];
2452 check_link(link, cg_gl, aj);
2454 i += nral_rt(ftype);
2458 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2460 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2465 cg_offset += cgs->nr;
2467 nlink_mol = link->index[cg_offset] - link->index[cg_offset-cgs->nr];
2469 destroy_reverse_ilist(&ril);
2474 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt->name, cgs->nr, nlink_mol);
2477 if (molb->nmol > mol)
2479 /* Copy the data for the rest of the molecules in this block */
2480 link->nalloc_a += (molb->nmol - mol)*nlink_mol;
2481 srenew(link->a, link->nalloc_a);
2482 for (; mol < molb->nmol; mol++)
2484 for (cg = 0; cg < cgs->nr; cg++)
2486 cg_gl = cg_offset + cg;
2487 link->index[cg_gl+1] =
2488 link->index[cg_gl+1-cgs->nr] + nlink_mol;
2489 for (j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2491 link->a[j] = link->a[j-nlink_mol] + cgs->nr;
2493 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2494 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2496 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2500 cg_offset += cgs->nr;
2505 if (mtop->bIntermolecularInteractions)
2507 destroy_reverse_ilist(&ril_intermol);
2512 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2523 } bonded_distance_t;
2525 /*! \brief Compare distance^2 \p r2 against the distance in \p bd and if larger store it along with \p ftype and atom indices \p a1 and \p a2 */
2526 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2527 bonded_distance_t *bd)
2538 /*! \brief Set the distance, function type and atom indices for the longest distance between charge-groups of molecule type \p molt for two-body and multi-body bonded interactions */
2539 static void bonded_cg_distance_mol(gmx_moltype_t *molt, int *at2cg,
2540 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2541 bonded_distance_t *bd_2b,
2542 bonded_distance_t *bd_mb)
2544 for (int ftype = 0; ftype < F_NRE; ftype++)
2546 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2548 const t_ilist *il = &molt->ilist[ftype];
2549 int nral = NRAL(ftype);
2552 for (int i = 0; i < il->nr; i += 1+nral)
2554 for (int ai = 0; ai < nral; ai++)
2556 int cgi = at2cg[il->iatoms[i+1+ai]];
2557 for (int aj = ai + 1; aj < nral; aj++)
2559 int cgj = at2cg[il->iatoms[i+1+aj]];
2562 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2564 update_max_bonded_distance(rij2, ftype,
2567 (nral == 2) ? bd_2b : bd_mb);
2577 const t_blocka *excls = &molt->excls;
2578 for (int ai = 0; ai < excls->nr; ai++)
2580 int cgi = at2cg[ai];
2581 for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2583 int cgj = at2cg[excls->a[j]];
2586 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2588 /* There is no function type for exclusions, use -1 */
2589 update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2596 /*! \brief Set the distance, function type and atom indices for the longest atom distance involved in intermolecular interactions for two-body and multi-body bonded interactions */
2597 static void bonded_distance_intermol(const t_ilist *ilists_intermol,
2599 const rvec *x, int ePBC, matrix box,
2600 bonded_distance_t *bd_2b,
2601 bonded_distance_t *bd_mb)
2605 set_pbc(&pbc, ePBC, box);
2607 for (int ftype = 0; ftype < F_NRE; ftype++)
2609 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2611 const t_ilist *il = &ilists_intermol[ftype];
2612 int nral = NRAL(ftype);
2614 /* No nral>1 check here, since intermol interactions always
2615 * have nral>=2 (and the code is also correct for nral=1).
2617 for (int i = 0; i < il->nr; i += 1+nral)
2619 for (int ai = 0; ai < nral; ai++)
2621 int atom_i = il->iatoms[i + 1 + ai];
2623 for (int aj = ai + 1; aj < nral; aj++)
2628 int atom_j = il->iatoms[i + 1 + aj];
2630 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2634 update_max_bonded_distance(rij2, ftype,
2636 (nral == 2) ? bd_2b : bd_mb);
2644 //! Compute charge group centers of mass for molecule \p molt
2645 static void get_cgcm_mol(const gmx_moltype_t *molt,
2646 const gmx_ffparams_t *ffparams,
2647 int ePBC, t_graph *graph, matrix box,
2648 const gmx_vsite_t *vsite,
2649 const rvec *x, rvec *xs, rvec *cg_cm)
2653 if (ePBC != epbcNONE)
2655 mk_mshift(NULL, graph, ePBC, box, x);
2657 shift_x(graph, box, x, xs);
2658 /* By doing an extra mk_mshift the molecules that are broken
2659 * because they were e.g. imported from another software
2660 * will be made whole again. Such are the healing powers
2663 mk_mshift(NULL, graph, ePBC, box, xs);
2667 /* We copy the coordinates so the original coordinates remain
2668 * unchanged, just to be 100% sure that we do not affect
2669 * binary reproducibility of simulations.
2671 n = molt->cgs.index[molt->cgs.nr];
2672 for (i = 0; i < n; i++)
2674 copy_rvec(x[i], xs[i]);
2680 construct_vsites(vsite, xs, 0.0, NULL,
2681 ffparams->iparams, molt->ilist,
2682 epbcNONE, TRUE, NULL, NULL);
2685 calc_cgcm(NULL, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2688 //! Returns whether \p molt has a virtual site
2689 static int have_vsite_molt(gmx_moltype_t *molt)
2695 for (i = 0; i < F_NRE; i++)
2697 if ((interaction_function[i].flags & IF_VSITE) &&
2698 molt->ilist[i].nr > 0)
2707 void dd_bonded_cg_distance(FILE *fplog,
2708 const gmx_mtop_t *mtop,
2709 const t_inputrec *ir,
2710 const rvec *x, matrix box,
2712 real *r_2b, real *r_mb)
2714 gmx_bool bExclRequired;
2715 int mb, at_offset, *at2cg, mol;
2718 gmx_molblock_t *molb;
2719 gmx_moltype_t *molt;
2721 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2722 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2724 bExclRequired = inputrecExclForces(ir);
2726 vsite = init_vsite(mtop, NULL, TRUE);
2731 for (mb = 0; mb < mtop->nmolblock; mb++)
2733 molb = &mtop->molblock[mb];
2734 molt = &mtop->moltype[molb->type];
2735 if (molt->cgs.nr == 1 || molb->nmol == 0)
2737 at_offset += molb->nmol*molt->atoms.nr;
2741 if (ir->ePBC != epbcNONE)
2743 mk_graph_ilist(NULL, molt->ilist, 0, molt->atoms.nr, FALSE, FALSE,
2747 at2cg = make_at2cg(&molt->cgs);
2748 snew(xs, molt->atoms.nr);
2749 snew(cg_cm, molt->cgs.nr);
2750 for (mol = 0; mol < molb->nmol; mol++)
2752 get_cgcm_mol(molt, &mtop->ffparams, ir->ePBC, &graph, box,
2753 have_vsite_molt(molt) ? vsite : NULL,
2754 x+at_offset, xs, cg_cm);
2756 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2757 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2759 bonded_cg_distance_mol(molt, at2cg, bBCheck, bExclRequired, cg_cm,
2760 &bd_mol_2b, &bd_mol_mb);
2762 /* Process the mol data adding the atom index offset */
2763 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2764 at_offset + bd_mol_2b.a1,
2765 at_offset + bd_mol_2b.a2,
2767 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2768 at_offset + bd_mol_mb.a1,
2769 at_offset + bd_mol_mb.a2,
2772 at_offset += molt->atoms.nr;
2777 if (ir->ePBC != epbcNONE)
2784 /* We should have a vsite free routine, but here we can simply free */
2787 if (mtop->bIntermolecularInteractions)
2789 if (ncg_mtop(mtop) < mtop->natoms)
2791 gmx_fatal(FARGS, "The combination of intermolecular interactions, charge groups and domain decomposition is not supported. Use cutoff-scheme=Verlet (which removes the charge groups) or run without domain decomposition.");
2794 bonded_distance_intermol(mtop->intermolecular_ilist,
2800 *r_2b = sqrt(bd_2b.r2);
2801 *r_mb = sqrt(bd_mb.r2);
2803 if (fplog && (*r_2b > 0 || *r_mb > 0))
2806 "Initial maximum inter charge-group distances:\n");
2810 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2811 *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2812 bd_2b.a1 + 1, bd_2b.a2 + 1);
2817 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2818 *r_mb, interaction_function[bd_mb.ftype].longname,
2819 bd_mb.a1 + 1, bd_mb.a2 + 1);