Partial conversion of gmx_domdec_t to C++
[gromacs.git] / src / gromacs / domdec / domdec_topology.cpp
blob129b87b784c9482ed0debba2c145d861ea245fbe
1 /*
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,2017,2018, 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.
36 /*! \internal \file
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
46 #include "gmxpre.h"
48 #include <assert.h>
49 #include <stdlib.h>
50 #include <string.h>
52 #include <cassert>
54 #include <algorithm>
55 #include <string>
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/forcerec.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdlib/vsite.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/ifunc.h"
74 #include "gromacs/topology/mtop_util.h"
75 #include "gromacs/topology/topsort.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/exceptions.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/gmxassert.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/stringutil.h"
83 #include "domdec_constraints.h"
84 #include "domdec_internal.h"
85 #include "domdec_vsite.h"
87 /*! \brief The number of integer item in the local state, used for broadcasting of the state */
88 #define NITEM_DD_INIT_LOCAL_STATE 5
90 struct reverse_ilist_t
92 int *index; /* Index for each atom into il */
93 int *il; /* ftype|type|a0|...|an|ftype|... */
94 int numAtomsInMolecule; /* The number of atoms in this molecule */
97 typedef struct {
98 int a_start;
99 int a_end;
100 int natoms_mol;
101 int type;
102 } molblock_ind_t;
104 /*! \brief Struct for thread local work data for local topology generation */
105 typedef struct {
106 t_idef idef; /**< Partial local topology */
107 int **vsite_pbc; /**< vsite PBC structure */
108 int *vsite_pbc_nalloc; /**< Allocation sizes for vsite_pbc */
109 int nbonded; /**< The number of bondeds in this struct */
110 t_blocka excl; /**< List of exclusions */
111 int excl_count; /**< The total exclusion count for \p excl */
112 } thread_work_t;
114 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
115 struct gmx_reverse_top_t
117 //! @cond Doxygen_Suppress
118 gmx_bool bExclRequired; /**< Do we require all exclusions to be assigned? */
119 int n_excl_at_max; /**< The maximum number of exclusions one atom can have */
120 gmx_bool bConstr; /**< Are there constraints in this revserse top? */
121 gmx_bool bSettle; /**< Are there settles in this revserse top? */
122 gmx_bool bBCheck; /**< All bonded interactions have to be assigned? */
123 gmx_bool bInterCGInteractions; /**< Are there bondeds/exclusions between charge-groups? */
124 reverse_ilist_t *ril_mt; /**< Reverse ilist for all moltypes */
125 int ril_mt_tot_size; /**< The size of ril_mt[?].index summed over all entries */
126 int ilsort; /**< The sorting state of bondeds for free energy */
127 molblock_ind_t *mbi; /**< molblock to global atom index for quick lookup of molblocks on atom index */
128 int nmolblock; /**< The number of molblocks, size of \p mbi */
130 gmx_bool bIntermolecularInteractions; /**< Do we have intermolecular interactions? */
131 reverse_ilist_t ril_intermol; /**< Intermolecular reverse ilist */
133 /* Work data structures for multi-threading */
134 int nthread; /**< The number of threads to be used */
135 thread_work_t *th_work; /**< Thread work array for local topology generation */
136 //! @endcond
139 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
140 static int nral_rt(int ftype)
142 int nral;
144 nral = NRAL(ftype);
145 if (interaction_function[ftype].flags & IF_VSITE)
147 /* With vsites the reverse topology contains
148 * two extra entries for PBC.
150 nral += 2;
153 return nral;
156 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
157 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
158 gmx_bool bConstr, gmx_bool bSettle)
160 return (((interaction_function[ftype].flags & IF_BOND) &&
161 !(interaction_function[ftype].flags & IF_VSITE) &&
162 (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))) ||
163 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
164 (bSettle && ftype == F_SETTLE));
167 /*! \brief Print a header on error messages */
168 static void print_error_header(FILE *fplog, const char *moltypename, int nprint)
170 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
171 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
172 fprintf(fplog,
173 "the first %d missing interactions, except for exclusions:\n",
174 nprint);
175 fprintf(stderr,
176 "the first %d missing interactions, except for exclusions:\n",
177 nprint);
180 /*! \brief Help print error output when interactions are missing */
181 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
182 const gmx_reverse_top_t *rt,
183 const char *moltypename,
184 const reverse_ilist_t *ril,
185 int a_start, int a_end,
186 int nat_mol, int nmol,
187 const t_idef *idef)
189 int *assigned;
190 int nril_mol = ril->index[nat_mol];
191 snew(assigned, nmol*nril_mol);
193 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
194 for (int ftype = 0; ftype < F_NRE; ftype++)
196 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
198 int nral = NRAL(ftype);
199 const t_ilist *il = &idef->il[ftype];
200 const t_iatom *ia = il->iatoms;
201 for (int i = 0; i < il->nr; i += 1+nral)
203 int a0 = gatindex[ia[1]];
204 /* Check if this interaction is in
205 * the currently checked molblock.
207 if (a0 >= a_start && a0 < a_end)
209 int mol = (a0 - a_start)/nat_mol;
210 int a0_mol = (a0 - a_start) - mol*nat_mol;
211 int j_mol = ril->index[a0_mol];
212 bool found = false;
213 while (j_mol < ril->index[a0_mol+1] && !found)
215 int j = mol*nril_mol + j_mol;
216 int ftype_j = ril->il[j_mol];
217 /* Here we need to check if this interaction has
218 * not already been assigned, since we could have
219 * multiply defined interactions.
221 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
222 assigned[j] == 0)
224 /* Check the atoms */
225 found = true;
226 for (int a = 0; a < nral; a++)
228 if (gatindex[ia[1+a]] !=
229 a_start + mol*nat_mol + ril->il[j_mol+2+a])
231 found = false;
234 if (found)
236 assigned[j] = 1;
239 j_mol += 2 + nral_rt(ftype_j);
241 if (!found)
243 gmx_incons("Some interactions seem to be assigned multiple times");
246 ia += 1 + nral;
251 gmx_sumi(nmol*nril_mol, assigned, cr);
253 int nprint = 10;
254 int i = 0;
255 for (int mol = 0; mol < nmol; mol++)
257 int j_mol = 0;
258 while (j_mol < nril_mol)
260 int ftype = ril->il[j_mol];
261 int nral = NRAL(ftype);
262 int j = mol*nril_mol + j_mol;
263 if (assigned[j] == 0 &&
264 !(interaction_function[ftype].flags & IF_VSITE))
266 if (DDMASTER(cr->dd))
268 if (i == 0)
270 print_error_header(fplog, moltypename, nprint);
272 fprintf(fplog, "%20s atoms",
273 interaction_function[ftype].longname);
274 fprintf(stderr, "%20s atoms",
275 interaction_function[ftype].longname);
276 int a;
277 for (a = 0; a < nral; a++)
279 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
280 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
282 while (a < 4)
284 fprintf(fplog, " ");
285 fprintf(stderr, " ");
286 a++;
288 fprintf(fplog, " global");
289 fprintf(stderr, " global");
290 for (a = 0; a < nral; a++)
292 fprintf(fplog, "%6d",
293 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
294 fprintf(stderr, "%6d",
295 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
297 fprintf(fplog, "\n");
298 fprintf(stderr, "\n");
300 i++;
301 if (i >= nprint)
303 break;
306 j_mol += 2 + nral_rt(ftype);
310 sfree(assigned);
313 /*! \brief Help print error output when interactions are missing */
314 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
315 const gmx_mtop_t *mtop,
316 const t_idef *idef)
318 const gmx_reverse_top_t *rt = cr->dd->reverse_top;
320 /* Print the atoms in the missing interactions per molblock */
321 int a_end = 0;
322 for (const gmx_molblock_t &molb : mtop->molblock)
324 const gmx_moltype_t &moltype = mtop->moltype[molb.type];
325 int a_start = a_end;
326 a_end = a_start + molb.nmol*moltype.atoms.nr;
328 print_missing_interactions_mb(fplog, cr, rt,
329 *(moltype.name),
330 &rt->ril_mt[molb.type],
331 a_start, a_end, moltype.atoms.nr,
332 molb.nmol,
333 idef);
337 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
338 int local_count,
339 const gmx_mtop_t *top_global,
340 const gmx_localtop_t *top_local,
341 const t_state *state_local)
343 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
344 int ftype, nral;
345 char buf[STRLEN];
346 gmx_domdec_t *dd;
348 dd = cr->dd;
350 if (fplog)
352 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
353 fflush(fplog);
356 ndiff_tot = local_count - dd->nbonded_global;
358 for (ftype = 0; ftype < F_NRE; ftype++)
360 nral = NRAL(ftype);
361 cl[ftype] = top_local->idef.il[ftype].nr/(1+nral);
364 gmx_sumi(F_NRE, cl, cr);
366 if (DDMASTER(dd))
368 if (fplog)
370 fprintf(fplog, "\nA list of missing interactions:\n");
372 fprintf(stderr, "\nA list of missing interactions:\n");
373 rest_global = dd->nbonded_global;
374 rest_local = local_count;
375 for (ftype = 0; ftype < F_NRE; ftype++)
377 /* In the reverse and local top all constraints are merged
378 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
379 * and add these constraints when doing F_CONSTR.
381 if (((interaction_function[ftype].flags & IF_BOND) &&
382 (dd->reverse_top->bBCheck
383 || !(interaction_function[ftype].flags & IF_LIMZERO)))
384 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
385 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
387 n = gmx_mtop_ftype_count(top_global, ftype);
388 if (ftype == F_CONSTR)
390 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
392 ndiff = cl[ftype] - n;
393 if (ndiff != 0)
395 sprintf(buf, "%20s of %6d missing %6d",
396 interaction_function[ftype].longname, n, -ndiff);
397 if (fplog)
399 fprintf(fplog, "%s\n", buf);
401 fprintf(stderr, "%s\n", buf);
403 rest_global -= n;
404 rest_local -= cl[ftype];
408 ndiff = rest_local - rest_global;
409 if (ndiff != 0)
411 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
412 rest_global, -ndiff);
413 if (fplog)
415 fprintf(fplog, "%s\n", buf);
417 fprintf(stderr, "%s\n", buf);
421 print_missing_interactions_atoms(fplog, cr, top_global, &top_local->idef);
422 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
423 -1, as_rvec_array(state_local->x.data()), state_local->box);
425 std::string errorMessage;
427 if (ndiff_tot > 0)
429 errorMessage = "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
431 else
433 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));
435 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), errorMessage.c_str());
438 /*! \brief Return global topology molecule information for global atom index \p i_gl */
439 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
440 int *mb, int *mt, int *mol, int *i_mol)
442 molblock_ind_t *mbi = rt->mbi;
443 int start = 0;
444 int end = rt->nmolblock; /* exclusive */
445 int mid;
447 /* binary search for molblock_ind */
448 while (TRUE)
450 mid = (start+end)/2;
451 if (i_gl >= mbi[mid].a_end)
453 start = mid+1;
455 else if (i_gl < mbi[mid].a_start)
457 end = mid;
459 else
461 break;
465 *mb = mid;
466 mbi += mid;
468 *mt = mbi->type;
469 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
470 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
473 /*! \brief Count the exclusions for all atoms in \p cgs */
474 static void count_excls(const t_block *cgs, const t_blocka *excls,
475 int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
477 int cg, at0, at1, at, excl, atj;
479 *n_excl = 0;
480 *n_intercg_excl = 0;
481 *n_excl_at_max = 0;
482 for (cg = 0; cg < cgs->nr; cg++)
484 at0 = cgs->index[cg];
485 at1 = cgs->index[cg+1];
486 for (at = at0; at < at1; at++)
488 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
490 atj = excls->a[excl];
491 if (atj > at)
493 (*n_excl)++;
494 if (atj < at0 || atj >= at1)
496 (*n_intercg_excl)++;
501 *n_excl_at_max = std::max(*n_excl_at_max,
502 excls->index[at+1] - excls->index[at]);
507 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
508 static int low_make_reverse_ilist(const t_ilist *il_mt, const t_atom *atom,
509 const int * const * vsite_pbc,
510 int *count,
511 gmx_bool bConstr, gmx_bool bSettle,
512 gmx_bool bBCheck,
513 int *r_index, int *r_il,
514 gmx_bool bLinkToAllAtoms,
515 gmx_bool bAssign)
517 int ftype, nral, i, j, nlink, link;
518 const t_ilist *il;
519 const t_iatom *ia;
520 int a;
521 int nint;
522 gmx_bool bVSite;
524 nint = 0;
525 for (ftype = 0; ftype < F_NRE; ftype++)
527 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
528 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
529 (bSettle && ftype == F_SETTLE))
531 bVSite = (interaction_function[ftype].flags & IF_VSITE);
532 nral = NRAL(ftype);
533 il = &il_mt[ftype];
534 for (i = 0; i < il->nr; i += 1+nral)
536 ia = il->iatoms + i;
537 if (bLinkToAllAtoms)
539 if (bVSite)
541 /* We don't need the virtual sites for the cg-links */
542 nlink = 0;
544 else
546 nlink = nral;
549 else
551 /* Couple to the first atom in the interaction */
552 nlink = 1;
554 for (link = 0; link < nlink; link++)
556 a = ia[1+link];
557 if (bAssign)
559 assert(r_il); //with bAssign not allowed to be null
560 assert(r_index);
561 r_il[r_index[a]+count[a]] =
562 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
563 r_il[r_index[a]+count[a]+1] = ia[0];
564 for (j = 1; j < 1+nral; j++)
566 /* Store the molecular atom number */
567 r_il[r_index[a]+count[a]+1+j] = ia[j];
570 if (interaction_function[ftype].flags & IF_VSITE)
572 if (bAssign)
574 /* Add an entry to iatoms for storing
575 * which of the constructing atoms are
576 * vsites again.
578 r_il[r_index[a]+count[a]+2+nral] = 0;
579 for (j = 2; j < 1+nral; j++)
581 if (atom[ia[j]].ptype == eptVSite)
583 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
586 /* Store vsite pbc atom in a second extra entry */
587 r_il[r_index[a]+count[a]+2+nral+1] =
588 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
591 else
593 /* We do not count vsites since they are always
594 * uniquely assigned and can be assigned
595 * to multiple nodes with recursive vsites.
597 if (bBCheck ||
598 !(interaction_function[ftype].flags & IF_LIMZERO))
600 nint++;
603 count[a] += 2 + nral_rt(ftype);
609 return nint;
612 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
613 static int make_reverse_ilist(const t_ilist *ilist,
614 const t_atoms *atoms,
615 const int * const * vsite_pbc,
616 gmx_bool bConstr, gmx_bool bSettle,
617 gmx_bool bBCheck,
618 gmx_bool bLinkToAllAtoms,
619 reverse_ilist_t *ril_mt)
621 int nat_mt, *count, i, nint_mt;
623 /* Count the interactions */
624 nat_mt = atoms->nr;
625 snew(count, nat_mt);
626 low_make_reverse_ilist(ilist, atoms->atom, vsite_pbc,
627 count,
628 bConstr, bSettle, bBCheck, nullptr, nullptr,
629 bLinkToAllAtoms, FALSE);
631 snew(ril_mt->index, nat_mt+1);
632 ril_mt->index[0] = 0;
633 for (i = 0; i < nat_mt; i++)
635 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
636 count[i] = 0;
638 snew(ril_mt->il, ril_mt->index[nat_mt]);
640 /* Store the interactions */
641 nint_mt =
642 low_make_reverse_ilist(ilist, atoms->atom, vsite_pbc,
643 count,
644 bConstr, bSettle, bBCheck,
645 ril_mt->index, ril_mt->il,
646 bLinkToAllAtoms, TRUE);
648 sfree(count);
650 ril_mt->numAtomsInMolecule = atoms->nr;
652 return nint_mt;
655 /*! \brief Destroys a reverse_ilist_t struct */
656 static void destroy_reverse_ilist(reverse_ilist_t *ril)
658 sfree(ril->index);
659 sfree(ril->il);
662 /*! \brief Generate the reverse topology */
663 static gmx_reverse_top_t *make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
664 const int * const * const * vsite_pbc_molt,
665 gmx_bool bConstr, gmx_bool bSettle,
666 gmx_bool bBCheck, int *nint)
668 gmx_reverse_top_t *rt;
669 int *nint_mt;
670 int thread;
672 snew(rt, 1);
674 /* Should we include constraints (for SHAKE) in rt? */
675 rt->bConstr = bConstr;
676 rt->bSettle = bSettle;
677 rt->bBCheck = bBCheck;
679 rt->bInterCGInteractions = mtop->bIntermolecularInteractions;
680 snew(nint_mt, mtop->moltype.size());
681 snew(rt->ril_mt, mtop->moltype.size());
682 rt->ril_mt_tot_size = 0;
683 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
685 const gmx_moltype_t &molt = mtop->moltype[mt];
686 if (molt.cgs.nr > 1)
688 rt->bInterCGInteractions = TRUE;
691 /* Make the atom to interaction list for this molecule type */
692 nint_mt[mt] =
693 make_reverse_ilist(molt.ilist, &molt.atoms,
694 vsite_pbc_molt ? vsite_pbc_molt[mt] : nullptr,
695 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
696 &rt->ril_mt[mt]);
698 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt.atoms.nr];
700 if (debug)
702 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
705 *nint = 0;
706 for (const gmx_molblock_t &molblock : mtop->molblock)
708 *nint += molblock.nmol*nint_mt[molblock.type];
710 sfree(nint_mt);
712 /* Make an intermolecular reverse top, if necessary */
713 rt->bIntermolecularInteractions = mtop->bIntermolecularInteractions;
714 if (rt->bIntermolecularInteractions)
716 t_atoms atoms_global;
718 rt->ril_intermol.index = nullptr;
719 rt->ril_intermol.il = nullptr;
721 atoms_global.nr = mtop->natoms;
722 atoms_global.atom = nullptr; /* Only used with virtual sites */
724 *nint +=
725 make_reverse_ilist(mtop->intermolecular_ilist, &atoms_global,
726 nullptr,
727 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
728 &rt->ril_intermol);
731 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
733 rt->ilsort = ilsortFE_UNSORTED;
735 else
737 rt->ilsort = ilsortNO_FE;
740 /* Make a molblock index for fast searching */
741 snew(rt->mbi, mtop->molblock.size());
742 rt->nmolblock = mtop->molblock.size();
743 int i = 0;
744 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
746 const gmx_molblock_t &molb = mtop->molblock[mb];
747 int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
748 rt->mbi[mb].a_start = i;
749 i += molb.nmol*numAtomsPerMol;
750 rt->mbi[mb].a_end = i;
751 rt->mbi[mb].natoms_mol = numAtomsPerMol;
752 rt->mbi[mb].type = molb.type;
755 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
756 snew(rt->th_work, rt->nthread);
757 if (vsite_pbc_molt != nullptr)
759 for (thread = 0; thread < rt->nthread; thread++)
761 snew(rt->th_work[thread].vsite_pbc, F_VSITEN-F_VSITE2+1);
762 snew(rt->th_work[thread].vsite_pbc_nalloc, F_VSITEN-F_VSITE2+1);
766 return rt;
769 void dd_make_reverse_top(FILE *fplog,
770 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
771 const gmx_vsite_t *vsite,
772 const t_inputrec *ir, gmx_bool bBCheck)
774 if (fplog)
776 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
779 /* If normal and/or settle constraints act only within charge groups,
780 * we can store them in the reverse top and simply assign them to domains.
781 * Otherwise we need to assign them to multiple domains and set up
782 * the parallel version constraint algorithm(s).
785 dd->reverse_top = make_reverse_top(mtop, ir->efep != efepNO,
786 vsite ? vsite->vsite_pbc_molt : nullptr,
787 !dd->bInterCGcons, !dd->bInterCGsettles,
788 bBCheck, &dd->nbonded_global);
790 gmx_reverse_top_t *rt = dd->reverse_top;
792 /* With the Verlet scheme, exclusions are handled in the non-bonded
793 * kernels and only exclusions inside the cut-off lead to exclusion
794 * forces. Since each atom pair is treated at most once in the non-bonded
795 * kernels, it doesn't matter if the exclusions for the same atom pair
796 * appear multiple times in the exclusion list. In contrast, the, old,
797 * group cut-off scheme loops over a list of exclusions, so there each
798 * excluded pair should appear exactly once.
800 rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
801 inputrecExclForces(ir));
803 int nexcl = 0;
804 dd->n_intercg_excl = 0;
805 rt->n_excl_at_max = 0;
806 for (const gmx_molblock_t &molb : mtop->molblock)
808 int n_excl_mol, n_excl_icg, n_excl_at_max;
810 const gmx_moltype_t &molt = mtop->moltype[molb.type];
811 count_excls(&molt.cgs, &molt.excls,
812 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
813 nexcl += molb.nmol*n_excl_mol;
814 dd->n_intercg_excl += molb.nmol*n_excl_icg;
815 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
817 if (rt->bExclRequired)
819 dd->nbonded_global += nexcl;
820 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
822 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
823 "will use an extra communication step for exclusion forces for %s\n",
824 dd->n_intercg_excl, eel_names[ir->coulombtype]);
828 if (vsite && vsite->n_intercg_vsite > 0)
830 if (fplog)
832 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
833 "will an extra communication step for selected coordinates and forces\n",
834 vsite->n_intercg_vsite);
836 init_domdec_vsites(dd, vsite->n_intercg_vsite);
839 if (dd->bInterCGcons || dd->bInterCGsettles)
841 init_domdec_constraints(dd, mtop);
843 if (fplog)
845 fprintf(fplog, "\n");
849 /*! \brief Store a vsite interaction at the end of \p il
851 * This routine is very similar to add_ifunc, but vsites interactions
852 * have more work to do than other kinds of interactions, and the
853 * complex way nral (and thus vector contents) depends on ftype
854 * confuses static analysis tools unless we fuse the vsite
855 * atom-indexing organization code with the ifunc-adding code, so that
856 * they can see that nral is the same value. */
857 static inline void
858 add_ifunc_for_vsites(t_iatom *tiatoms, gmx_ga2la_t *ga2la,
859 int nral, gmx_bool bHomeA,
860 int a, int a_gl, int a_mol,
861 const t_iatom *iatoms,
862 t_ilist *il)
864 t_iatom *liatoms;
866 if (il->nr+1+nral > il->nalloc)
868 il->nalloc = over_alloc_large(il->nr+1+nral);
869 srenew(il->iatoms, il->nalloc);
871 liatoms = il->iatoms + il->nr;
872 il->nr += 1 + nral;
874 /* Copy the type */
875 tiatoms[0] = iatoms[0];
877 if (bHomeA)
879 /* We know the local index of the first atom */
880 tiatoms[1] = a;
882 else
884 /* Convert later in make_local_vsites */
885 tiatoms[1] = -a_gl - 1;
888 for (int k = 2; k < 1+nral; k++)
890 int ak_gl = a_gl + iatoms[k] - a_mol;
891 if (!ga2la_get_home(ga2la, ak_gl, &tiatoms[k]))
893 /* Copy the global index, convert later in make_local_vsites */
894 tiatoms[k] = -(ak_gl + 1);
896 // Note that ga2la_get_home always sets the third parameter if
897 // it returns TRUE
899 for (int k = 0; k < 1+nral; k++)
901 liatoms[k] = tiatoms[k];
905 /*! \brief Store a bonded interaction at the end of \p il */
906 static inline void add_ifunc(int nral, t_iatom *tiatoms, t_ilist *il)
908 t_iatom *liatoms;
909 int k;
911 if (il->nr+1+nral > il->nalloc)
913 il->nalloc = over_alloc_large(il->nr+1+nral);
914 srenew(il->iatoms, il->nalloc);
916 liatoms = il->iatoms + il->nr;
917 for (k = 0; k <= nral; k++)
919 liatoms[k] = tiatoms[k];
921 il->nr += 1 + nral;
924 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
925 static void add_posres(int mol, int a_mol, int numAtomsInMolecule,
926 const gmx_molblock_t *molb,
927 t_iatom *iatoms, const t_iparams *ip_in,
928 t_idef *idef)
930 int n, a_molb;
931 t_iparams *ip;
933 /* This position restraint has not been added yet,
934 * so it's index is the current number of position restraints.
936 n = idef->il[F_POSRES].nr/2;
937 if (n+1 > idef->iparams_posres_nalloc)
939 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
940 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
942 ip = &idef->iparams_posres[n];
943 /* Copy the force constants */
944 *ip = ip_in[iatoms[0]];
946 /* Get the position restraint coordinates from the molblock */
947 a_molb = mol*numAtomsInMolecule + a_mol;
948 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
949 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
950 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
951 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
952 if (!molb->posres_xB.empty())
954 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
955 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
956 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
958 else
960 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
961 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
962 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
964 /* Set the parameter index for idef->iparams_posre */
965 iatoms[0] = n;
968 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
969 static void add_fbposres(int mol, int a_mol, int numAtomsInMolecule,
970 const gmx_molblock_t *molb,
971 t_iatom *iatoms, const t_iparams *ip_in,
972 t_idef *idef)
974 int n, a_molb;
975 t_iparams *ip;
977 /* This flat-bottom position restraint has not been added yet,
978 * so it's index is the current number of position restraints.
980 n = idef->il[F_FBPOSRES].nr/2;
981 if (n+1 > idef->iparams_fbposres_nalloc)
983 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
984 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
986 ip = &idef->iparams_fbposres[n];
987 /* Copy the force constants */
988 *ip = ip_in[iatoms[0]];
990 /* Get the position restraint coordinats from the molblock */
991 a_molb = mol*numAtomsInMolecule + a_mol;
992 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
993 /* Take reference positions from A position of normal posres */
994 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
995 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
996 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
998 /* Note: no B-type for flat-bottom posres */
1000 /* Set the parameter index for idef->iparams_posre */
1001 iatoms[0] = n;
1004 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
1005 static void add_vsite(gmx_ga2la_t *ga2la, const int *index, const int *rtil,
1006 int ftype, int nral,
1007 gmx_bool bHomeA, int a, int a_gl, int a_mol,
1008 const t_iatom *iatoms,
1009 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
1011 int k, vsi, pbc_a_mol;
1012 t_iatom tiatoms[1+MAXATOMLIST];
1013 int j, ftype_r, nral_r;
1015 /* Add this interaction to the local topology */
1016 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1018 if (vsite_pbc)
1020 vsi = idef->il[ftype].nr/(1+nral) - 1;
1021 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
1023 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
1024 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
1026 if (bHomeA)
1028 pbc_a_mol = iatoms[1+nral+1];
1029 if (pbc_a_mol < 0)
1031 /* The pbc flag is one of the following two options:
1032 * -2: vsite and all constructing atoms are within the same cg, no pbc
1033 * -1: vsite and its first constructing atom are in the same cg, do pbc
1035 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
1037 else
1039 /* Set the pbc atom for this vsite so we can make its pbc
1040 * identical to the rest of the atoms in its charge group.
1041 * Since the order of the atoms does not change within a charge
1042 * group, we do not need the global to local atom index.
1044 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
1047 else
1049 /* This vsite is non-home (required for recursion),
1050 * and therefore there is no charge group to match pbc with.
1051 * But we always turn on full_pbc to assure that higher order
1052 * recursion works correctly.
1054 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
1058 if (iatoms[1+nral])
1060 /* Check for recursion */
1061 for (k = 2; k < 1+nral; k++)
1063 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1065 /* This construction atoms is a vsite and not a home atom */
1066 if (gmx_debug_at)
1068 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1070 /* Find the vsite construction */
1072 /* Check all interactions assigned to this atom */
1073 j = index[iatoms[k]];
1074 while (j < index[iatoms[k]+1])
1076 ftype_r = rtil[j++];
1077 nral_r = NRAL(ftype_r);
1078 if (interaction_function[ftype_r].flags & IF_VSITE)
1080 /* Add this vsite (recursion) */
1081 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1082 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1083 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
1084 j += 1 + nral_r + 2;
1086 else
1088 j += 1 + nral_r;
1096 /*! \brief Build the index that maps each local atom to its local atom group */
1097 static void makeLocalAtomGroupFromAtom(gmx_domdec_t *dd)
1099 const gmx::BlockRanges &atomGroups = dd->atomGroups();
1101 dd->localAtomGroupFromAtom.clear();
1103 for (size_t g = 0; g < dd->globalAtomGroupIndices.size(); g++)
1105 for (int a = atomGroups.index[g]; a < atomGroups.index[g + 1]; a++)
1107 dd->localAtomGroupFromAtom.push_back(g);
1112 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1113 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1115 rvec dx;
1117 if (pbc_null)
1119 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1121 else
1123 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1126 return norm2(dx);
1129 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1130 static void combine_blocka(t_blocka *dest, const thread_work_t *src, int nsrc)
1132 int ni, na, s, i;
1134 ni = src[nsrc-1].excl.nr;
1135 na = 0;
1136 for (s = 0; s < nsrc; s++)
1138 na += src[s].excl.nra;
1140 if (ni + 1 > dest->nalloc_index)
1142 dest->nalloc_index = over_alloc_large(ni+1);
1143 srenew(dest->index, dest->nalloc_index);
1145 if (dest->nra + na > dest->nalloc_a)
1147 dest->nalloc_a = over_alloc_large(dest->nra+na);
1148 srenew(dest->a, dest->nalloc_a);
1150 for (s = 1; s < nsrc; s++)
1152 for (i = dest->nr+1; i < src[s].excl.nr+1; i++)
1154 dest->index[i] = dest->nra + src[s].excl.index[i];
1156 for (i = 0; i < src[s].excl.nra; i++)
1158 dest->a[dest->nra+i] = src[s].excl.a[i];
1160 dest->nr = src[s].excl.nr;
1161 dest->nra += src[s].excl.nra;
1165 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1166 * virtual sites need special attention, as pbc info differs per vsite.
1168 static void combine_idef(t_idef *dest, const thread_work_t *src, int nsrc,
1169 gmx_vsite_t *vsite)
1171 int ftype;
1173 for (ftype = 0; ftype < F_NRE; ftype++)
1175 int n, s;
1177 n = 0;
1178 for (s = 1; s < nsrc; s++)
1180 n += src[s].idef.il[ftype].nr;
1182 if (n > 0)
1184 t_ilist *ild;
1186 ild = &dest->il[ftype];
1188 if (ild->nr + n > ild->nalloc)
1190 ild->nalloc = over_alloc_large(ild->nr+n);
1191 srenew(ild->iatoms, ild->nalloc);
1194 gmx_bool vpbc;
1195 int nral1 = 0, ftv = 0;
1197 vpbc = ((interaction_function[ftype].flags & IF_VSITE) &&
1198 vsite->vsite_pbc_loc != nullptr);
1199 if (vpbc)
1201 nral1 = 1 + NRAL(ftype);
1202 ftv = ftype - F_VSITE2;
1203 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1205 vsite->vsite_pbc_loc_nalloc[ftv] =
1206 over_alloc_large((ild->nr + n)/nral1);
1207 srenew(vsite->vsite_pbc_loc[ftv],
1208 vsite->vsite_pbc_loc_nalloc[ftv]);
1212 for (s = 1; s < nsrc; s++)
1214 const t_ilist *ils;
1215 int i;
1217 ils = &src[s].idef.il[ftype];
1218 for (i = 0; i < ils->nr; i++)
1220 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1222 if (vpbc)
1224 for (i = 0; i < ils->nr; i += nral1)
1226 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1227 src[s].vsite_pbc[ftv][i/nral1];
1231 ild->nr += ils->nr;
1234 /* Position restraints need an additional treatment */
1235 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1237 int nposres = dest->il[ftype].nr/2;
1238 // TODO: Simplify this code using std::vector
1239 t_iparams * &iparams_dest = (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1240 int &posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc : dest->iparams_fbposres_nalloc);
1241 if (nposres > posres_nalloc)
1243 posres_nalloc = over_alloc_large(nposres);
1244 srenew(iparams_dest, posres_nalloc);
1247 /* Set nposres to the number of original position restraints in dest */
1248 for (int s = 1; s < nsrc; s++)
1250 nposres -= src[s].idef.il[ftype].nr/2;
1253 for (int s = 1; s < nsrc; s++)
1255 const t_iparams *iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1257 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1259 /* Correct the index into iparams_posres */
1260 dest->il[ftype].iatoms[nposres*2] = nposres;
1261 /* Copy the position restraint force parameters */
1262 iparams_dest[nposres] = iparams_src[i];
1263 nposres++;
1271 /*! \brief Check and when available assign bonded interactions for local atom i
1273 static inline void
1274 check_assign_interactions_atom(int i, int i_gl,
1275 int mol, int i_mol,
1276 int numAtomsInMolecule,
1277 const int *index, const int *rtil,
1278 gmx_bool bInterMolInteractions,
1279 int ind_start, int ind_end,
1280 const gmx_domdec_t *dd,
1281 const gmx_domdec_zones_t *zones,
1282 const gmx_molblock_t *molb,
1283 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1284 real rc2,
1285 int *la2lc,
1286 t_pbc *pbc_null,
1287 rvec *cg_cm,
1288 const t_iparams *ip_in,
1289 t_idef *idef,
1290 int **vsite_pbc, int *vsite_pbc_nalloc,
1291 int iz,
1292 gmx_bool bBCheck,
1293 int *nbonded_local)
1295 int j;
1297 j = ind_start;
1298 while (j < ind_end)
1300 int ftype;
1301 const t_iatom *iatoms;
1302 int nral;
1303 t_iatom tiatoms[1 + MAXATOMLIST];
1305 ftype = rtil[j++];
1306 iatoms = rtil + j;
1307 nral = NRAL(ftype);
1308 if (ftype == F_SETTLE)
1310 /* Settles are only in the reverse top when they
1311 * operate within a charge group. So we can assign
1312 * them without checks. We do this only for performance
1313 * reasons; it could be handled by the code below.
1315 if (iz == 0)
1317 /* Home zone: add this settle to the local topology */
1318 tiatoms[0] = iatoms[0];
1319 tiatoms[1] = i;
1320 tiatoms[2] = i + iatoms[2] - iatoms[1];
1321 tiatoms[3] = i + iatoms[3] - iatoms[1];
1322 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1323 (*nbonded_local)++;
1325 j += 1 + nral;
1327 else if (interaction_function[ftype].flags & IF_VSITE)
1329 assert(!bInterMolInteractions);
1330 /* The vsite construction goes where the vsite itself is */
1331 if (iz == 0)
1333 add_vsite(dd->ga2la, index, rtil, ftype, nral,
1334 TRUE, i, i_gl, i_mol,
1335 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1337 j += 1 + nral + 2;
1339 else
1341 gmx_bool bUse;
1343 /* Copy the type */
1344 tiatoms[0] = iatoms[0];
1346 if (nral == 1)
1348 assert(!bInterMolInteractions);
1349 /* Assign single-body interactions to the home zone */
1350 if (iz == 0)
1352 bUse = TRUE;
1353 tiatoms[1] = i;
1354 if (ftype == F_POSRES)
1356 add_posres(mol, i_mol, numAtomsInMolecule,
1357 molb, tiatoms, ip_in, idef);
1359 else if (ftype == F_FBPOSRES)
1361 add_fbposres(mol, i_mol, numAtomsInMolecule,
1362 molb, tiatoms, ip_in, idef);
1365 else
1367 bUse = FALSE;
1370 else if (nral == 2)
1372 /* This is a two-body interaction, we can assign
1373 * analogous to the non-bonded assignments.
1375 int k_gl, a_loc, kz;
1377 if (!bInterMolInteractions)
1379 /* Get the global index using the offset in the molecule */
1380 k_gl = i_gl + iatoms[2] - i_mol;
1382 else
1384 k_gl = iatoms[2];
1386 if (!ga2la_get(dd->ga2la, k_gl, &a_loc, &kz))
1388 bUse = FALSE;
1390 else
1392 if (kz >= zones->n)
1394 kz -= zones->n;
1396 /* Check zone interaction assignments */
1397 bUse = ((iz < zones->nizone &&
1398 iz <= kz &&
1399 kz >= zones->izone[iz].j0 &&
1400 kz < zones->izone[iz].j1) ||
1401 (kz < zones->nizone &&
1402 iz > kz &&
1403 iz >= zones->izone[kz].j0 &&
1404 iz < zones->izone[kz].j1));
1405 if (bUse)
1407 tiatoms[1] = i;
1408 tiatoms[2] = a_loc;
1409 /* If necessary check the cgcm distance */
1410 if (bRCheck2B &&
1411 dd_dist2(pbc_null, cg_cm, la2lc,
1412 tiatoms[1], tiatoms[2]) >= rc2)
1414 bUse = FALSE;
1419 else
1421 /* Assign this multi-body bonded interaction to
1422 * the local node if we have all the atoms involved
1423 * (local or communicated) and the minimum zone shift
1424 * in each dimension is zero, for dimensions
1425 * with 2 DD cells an extra check may be necessary.
1427 ivec k_zero, k_plus;
1428 int k;
1430 bUse = TRUE;
1431 clear_ivec(k_zero);
1432 clear_ivec(k_plus);
1433 for (k = 1; k <= nral && bUse; k++)
1435 gmx_bool bLocal;
1436 int k_gl, a_loc;
1437 int kz;
1439 if (!bInterMolInteractions)
1441 /* Get the global index using the offset in the molecule */
1442 k_gl = i_gl + iatoms[k] - i_mol;
1444 else
1446 k_gl = iatoms[k];
1448 bLocal = ga2la_get(dd->ga2la, k_gl, &a_loc, &kz);
1449 if (!bLocal || kz >= zones->n)
1451 /* We do not have this atom of this interaction
1452 * locally, or it comes from more than one cell
1453 * away.
1455 bUse = FALSE;
1457 else
1459 int d;
1461 tiatoms[k] = a_loc;
1462 for (d = 0; d < DIM; d++)
1464 if (zones->shift[kz][d] == 0)
1466 k_zero[d] = k;
1468 else
1470 k_plus[d] = k;
1475 bUse = (bUse &&
1476 k_zero[XX] && k_zero[YY] && k_zero[ZZ]);
1477 if (bRCheckMB)
1479 int d;
1481 for (d = 0; (d < DIM && bUse); d++)
1483 /* Check if the cg_cm distance falls within
1484 * the cut-off to avoid possible multiple
1485 * assignments of bonded interactions.
1487 if (rcheck[d] &&
1488 k_plus[d] &&
1489 dd_dist2(pbc_null, cg_cm, la2lc,
1490 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1492 bUse = FALSE;
1497 if (bUse)
1499 /* Add this interaction to the local topology */
1500 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1501 /* Sum so we can check in global_stat
1502 * if we have everything.
1504 if (bBCheck ||
1505 !(interaction_function[ftype].flags & IF_LIMZERO))
1507 (*nbonded_local)++;
1510 j += 1 + nral;
1515 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1517 * With thread parallelizing each thread acts on a different atom range:
1518 * at_start to at_end.
1520 static int make_bondeds_zone(gmx_domdec_t *dd,
1521 const gmx_domdec_zones_t *zones,
1522 const std::vector<gmx_molblock_t> &molb,
1523 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1524 real rc2,
1525 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1526 const t_iparams *ip_in,
1527 t_idef *idef,
1528 int **vsite_pbc,
1529 int *vsite_pbc_nalloc,
1530 int izone,
1531 int at_start, int at_end)
1533 int mb, mt, mol, i_mol;
1534 int *index, *rtil;
1535 gmx_bool bBCheck;
1536 gmx_reverse_top_t *rt;
1537 int nbonded_local;
1539 rt = dd->reverse_top;
1541 bBCheck = rt->bBCheck;
1543 nbonded_local = 0;
1545 for (int i = at_start; i < at_end; i++)
1547 /* Get the global atom number */
1548 const int i_gl = dd->globalAtomIndices[i];
1549 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1550 /* Check all intramolecular interactions assigned to this atom */
1551 index = rt->ril_mt[mt].index;
1552 rtil = rt->ril_mt[mt].il;
1554 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1555 rt->ril_mt[mt].numAtomsInMolecule,
1556 index, rtil, FALSE,
1557 index[i_mol], index[i_mol+1],
1558 dd, zones,
1559 &molb[mb],
1560 bRCheckMB, rcheck, bRCheck2B, rc2,
1561 la2lc,
1562 pbc_null,
1563 cg_cm,
1564 ip_in,
1565 idef, vsite_pbc, vsite_pbc_nalloc,
1566 izone,
1567 bBCheck,
1568 &nbonded_local);
1571 if (rt->bIntermolecularInteractions)
1573 /* Check all intermolecular interactions assigned to this atom */
1574 index = rt->ril_intermol.index;
1575 rtil = rt->ril_intermol.il;
1577 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1578 rt->ril_mt[mt].numAtomsInMolecule,
1579 index, rtil, TRUE,
1580 index[i_gl], index[i_gl + 1],
1581 dd, zones,
1582 &molb[mb],
1583 bRCheckMB, rcheck, bRCheck2B, rc2,
1584 la2lc,
1585 pbc_null,
1586 cg_cm,
1587 ip_in,
1588 idef, vsite_pbc, vsite_pbc_nalloc,
1589 izone,
1590 bBCheck,
1591 &nbonded_local);
1595 return nbonded_local;
1598 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1599 static void set_no_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1600 int iz, t_blocka *lexcls)
1602 int a0 = dd->atomGroups().index[zones->cg_range[iz]];
1603 int a1 = dd->atomGroups().index[zones->cg_range[iz + 1]];
1605 for (int a = a0 + 1; a < a1 + 1; a++)
1607 lexcls->index[a] = lexcls->nra;
1611 /*! \brief Set the exclusion data for i-zone \p iz
1613 * This is a legacy version for the group scheme of the same routine below.
1614 * Here charge groups and distance checks to ensure unique exclusions
1615 * are supported.
1617 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1618 const std::vector<gmx_moltype_t> &moltype,
1619 gmx_bool bRCheck, real rc2,
1620 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1621 const int *cginfo,
1622 t_blocka *lexcls,
1623 int iz,
1624 int cg_start, int cg_end)
1626 int n_excl_at_max;
1627 int mb, mt, mol;
1628 const t_blocka *excls;
1629 gmx_ga2la_t *ga2la;
1630 int cell;
1632 ga2la = dd->ga2la;
1634 int jla0 = dd->atomGroups().index[zones->izone[iz].jcg0];
1635 int jla1 = dd->atomGroups().index[zones->izone[iz].jcg1];
1637 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1639 /* We set the end index, but note that we might not start at zero here */
1640 lexcls->nr = dd->atomGroups().index[cg_end];
1642 int n = lexcls->nra;
1643 int count = 0;
1644 for (int cg = cg_start; cg < cg_end; cg++)
1646 if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1648 lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1649 srenew(lexcls->a, lexcls->nalloc_a);
1651 int la0 = dd->atomGroups().index[cg];
1652 int la1 = dd->atomGroups().index[cg + 1];
1653 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1654 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1656 /* Copy the exclusions from the global top */
1657 for (int la = la0; la < la1; la++)
1659 lexcls->index[la] = n;
1660 int a_gl = dd->globalAtomIndices[la];
1661 int a_mol;
1662 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1663 excls = &moltype[mt].excls;
1664 for (int j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1666 int aj_mol = excls->a[j];
1667 /* This computation of jla is only correct intra-cg */
1668 int jla = la + aj_mol - a_mol;
1669 if (jla >= la0 && jla < la1)
1671 /* This is an intra-cg exclusion. We can skip
1672 * the global indexing and distance checking.
1674 /* Intra-cg exclusions are only required
1675 * for the home zone.
1677 if (iz == 0)
1679 lexcls->a[n++] = jla;
1680 /* Check to avoid double counts */
1681 if (jla > la)
1683 count++;
1687 else
1689 /* This is a inter-cg exclusion */
1690 /* Since exclusions are pair interactions,
1691 * just like non-bonded interactions,
1692 * they can be assigned properly up
1693 * to the DD cutoff (not cutoff_min as
1694 * for the other bonded interactions).
1696 if (ga2la_get(ga2la, a_gl+aj_mol-a_mol, &jla, &cell))
1698 if (iz == 0 && cell == 0)
1700 lexcls->a[n++] = jla;
1701 /* Check to avoid double counts */
1702 if (jla > la)
1704 count++;
1707 else if (jla >= jla0 && jla < jla1 &&
1708 (!bRCheck ||
1709 dd_dist2(pbc_null, cg_cm, la2lc, la, jla) < rc2))
1711 /* jla > la, since jla0 > la */
1712 lexcls->a[n++] = jla;
1713 count++;
1720 else
1722 /* There are no inter-cg excls and this cg is self-excluded.
1723 * These exclusions are only required for zone 0,
1724 * since other zones do not see themselves.
1726 if (iz == 0)
1728 for (int la = la0; la < la1; la++)
1730 lexcls->index[la] = n;
1731 for (int j = la0; j < la1; j++)
1733 lexcls->a[n++] = j;
1736 count += ((la1 - la0)*(la1 - la0 - 1))/2;
1738 else
1740 /* We don't need exclusions for this cg */
1741 for (int la = la0; la < la1; la++)
1743 lexcls->index[la] = n;
1749 lexcls->index[lexcls->nr] = n;
1750 lexcls->nra = n;
1752 return count;
1755 /*! \brief Set the exclusion data for i-zone \p iz */
1756 static void make_exclusions_zone(gmx_domdec_t *dd,
1757 gmx_domdec_zones_t *zones,
1758 const std::vector<gmx_moltype_t> &moltype,
1759 const int *cginfo,
1760 t_blocka *lexcls,
1761 int iz,
1762 int at_start, int at_end)
1764 gmx_ga2la_t *ga2la;
1765 int n_excl_at_max, n, at;
1767 ga2la = dd->ga2la;
1769 int jla0 = dd->atomGroups().index[zones->izone[iz].jcg0];
1770 int jla1 = dd->atomGroups().index[zones->izone[iz].jcg1];
1772 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1774 /* We set the end index, but note that we might not start at zero here */
1775 lexcls->nr = at_end;
1777 n = lexcls->nra;
1778 for (at = at_start; at < at_end; at++)
1780 if (n + 1000 > lexcls->nalloc_a)
1782 lexcls->nalloc_a = over_alloc_large(n + 1000);
1783 srenew(lexcls->a, lexcls->nalloc_a);
1785 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1787 int a_gl, mb, mt, mol, a_mol, j;
1788 const t_blocka *excls;
1790 if (n + n_excl_at_max > lexcls->nalloc_a)
1792 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1793 srenew(lexcls->a, lexcls->nalloc_a);
1796 /* Copy the exclusions from the global top */
1797 lexcls->index[at] = n;
1798 a_gl = dd->globalAtomIndices[at];
1799 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1800 &mb, &mt, &mol, &a_mol);
1801 excls = &moltype[mt].excls;
1802 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1804 int aj_mol, at_j, cell;
1806 aj_mol = excls->a[j];
1808 if (ga2la_get(ga2la, a_gl + aj_mol - a_mol, &at_j, &cell))
1810 /* This check is not necessary, but it can reduce
1811 * the number of exclusions in the list, which in turn
1812 * can speed up the pair list construction a bit.
1814 if (at_j >= jla0 && at_j < jla1)
1816 lexcls->a[n++] = at_j;
1821 else
1823 /* We don't need exclusions for this atom */
1824 lexcls->index[at] = n;
1828 lexcls->index[lexcls->nr] = n;
1829 lexcls->nra = n;
1833 /*! \brief Ensure we have enough space in \p ba for \p nindex_max indices */
1834 static void check_alloc_index(t_blocka *ba, int nindex_max)
1836 if (nindex_max+1 > ba->nalloc_index)
1838 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1839 srenew(ba->index, ba->nalloc_index);
1843 /*! \brief Ensure that we have enough space for exclusion storate in \p lexcls */
1844 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1845 t_blocka *lexcls)
1847 int nr = dd->atomGroups().index[zones->izone[zones->nizone - 1].cg1];
1849 check_alloc_index(lexcls, nr);
1851 for (int thread = 1; thread < dd->reverse_top->nthread; thread++)
1853 check_alloc_index(&dd->reverse_top->th_work[thread].excl, nr);
1857 /*! \brief Set the total count indexes for the local exclusions, needed by several functions */
1858 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1859 t_blocka *lexcls)
1861 lexcls->nr = dd->atomGroups().index[zones->izone[zones->nizone - 1].cg1];
1863 if (dd->n_intercg_excl == 0)
1865 /* There are no exclusions involving non-home charge groups,
1866 * but we need to set the indices for neighborsearching.
1868 int la0 = dd->atomGroups().index[zones->izone[0].cg1];
1869 for (int la = la0; la < lexcls->nr; la++)
1871 lexcls->index[la] = lexcls->nra;
1874 /* nr is only used to loop over the exclusions for Ewald and RF,
1875 * so we can set it to the number of home atoms for efficiency.
1877 lexcls->nr = dd->atomGroups().index[zones->izone[0].cg1];
1881 /*! \brief Clear a t_idef struct */
1882 static void clear_idef(t_idef *idef)
1884 int ftype;
1886 /* Clear the counts */
1887 for (ftype = 0; ftype < F_NRE; ftype++)
1889 idef->il[ftype].nr = 0;
1893 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1894 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1895 gmx_domdec_zones_t *zones,
1896 const gmx_mtop_t *mtop,
1897 const int *cginfo,
1898 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1899 real rc,
1900 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1901 t_idef *idef, gmx_vsite_t *vsite,
1902 t_blocka *lexcls, int *excl_count)
1904 int nzone_bondeds, nzone_excl;
1905 int izone, cg0, cg1;
1906 real rc2;
1907 int nbonded_local;
1908 int thread;
1909 gmx_reverse_top_t *rt;
1911 if (dd->reverse_top->bInterCGInteractions)
1913 nzone_bondeds = zones->n;
1915 else
1917 /* Only single charge group (or atom) molecules, so interactions don't
1918 * cross zone boundaries and we only need to assign in the home zone.
1920 nzone_bondeds = 1;
1923 if (dd->n_intercg_excl > 0)
1925 /* We only use exclusions from i-zones to i- and j-zones */
1926 nzone_excl = zones->nizone;
1928 else
1930 /* There are no inter-cg exclusions and only zone 0 sees itself */
1931 nzone_excl = 1;
1934 check_exclusions_alloc(dd, zones, lexcls);
1936 rt = dd->reverse_top;
1938 rc2 = rc*rc;
1940 /* Clear the counts */
1941 clear_idef(idef);
1942 nbonded_local = 0;
1944 lexcls->nr = 0;
1945 lexcls->nra = 0;
1946 *excl_count = 0;
1948 for (izone = 0; izone < nzone_bondeds; izone++)
1950 cg0 = zones->cg_range[izone];
1951 cg1 = zones->cg_range[izone + 1];
1953 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1954 for (thread = 0; thread < rt->nthread; thread++)
1958 int cg0t, cg1t;
1959 t_idef *idef_t;
1960 int **vsite_pbc;
1961 int *vsite_pbc_nalloc;
1962 t_blocka *excl_t;
1964 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
1965 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
1967 if (thread == 0)
1969 idef_t = idef;
1971 else
1973 idef_t = &rt->th_work[thread].idef;
1974 clear_idef(idef_t);
1977 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1979 if (thread == 0)
1981 vsite_pbc = vsite->vsite_pbc_loc;
1982 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
1984 else
1986 vsite_pbc = rt->th_work[thread].vsite_pbc;
1987 vsite_pbc_nalloc = rt->th_work[thread].vsite_pbc_nalloc;
1990 else
1992 vsite_pbc = nullptr;
1993 vsite_pbc_nalloc = nullptr;
1996 rt->th_work[thread].nbonded =
1997 make_bondeds_zone(dd, zones,
1998 mtop->molblock,
1999 bRCheckMB, rcheck, bRCheck2B, rc2,
2000 la2lc, pbc_null, cg_cm, idef->iparams,
2001 idef_t,
2002 vsite_pbc, vsite_pbc_nalloc,
2003 izone,
2004 dd->atomGroups().index[cg0t],
2005 dd->atomGroups().index[cg1t]);
2007 if (izone < nzone_excl)
2009 if (thread == 0)
2011 excl_t = lexcls;
2013 else
2015 excl_t = &rt->th_work[thread].excl;
2016 excl_t->nr = 0;
2017 excl_t->nra = 0;
2020 int numAtomGroups = dd->globalAtomGroupIndices.size();
2021 if (dd->atomGroups().index[numAtomGroups] == numAtomGroups &&
2022 !rt->bExclRequired)
2024 /* No charge groups and no distance check required */
2025 make_exclusions_zone(dd, zones,
2026 mtop->moltype, cginfo,
2027 excl_t,
2028 izone,
2029 cg0t, cg1t);
2031 else
2033 rt->th_work[thread].excl_count =
2034 make_exclusions_zone_cg(dd, zones,
2035 mtop->moltype, bRCheck2B, rc2,
2036 la2lc, pbc_null, cg_cm, cginfo,
2037 excl_t,
2038 izone,
2039 cg0t, cg1t);
2043 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2046 if (rt->nthread > 1)
2048 combine_idef(idef, rt->th_work, rt->nthread, vsite);
2051 for (thread = 0; thread < rt->nthread; thread++)
2053 nbonded_local += rt->th_work[thread].nbonded;
2056 if (izone < nzone_excl)
2058 if (rt->nthread > 1)
2060 combine_blocka(lexcls, rt->th_work, rt->nthread);
2063 for (thread = 0; thread < rt->nthread; thread++)
2065 *excl_count += rt->th_work[thread].excl_count;
2070 /* Some zones might not have exclusions, but some code still needs to
2071 * loop over the index, so we set the indices here.
2073 for (izone = nzone_excl; izone < zones->nizone; izone++)
2075 set_no_exclusions_zone(dd, zones, izone, lexcls);
2078 finish_local_exclusions(dd, zones, lexcls);
2079 if (debug)
2081 fprintf(debug, "We have %d exclusions, check count %d\n",
2082 lexcls->nra, *excl_count);
2085 return nbonded_local;
2088 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
2090 lcgs->nr = dd->globalAtomGroupIndices.size();
2091 lcgs->index = dd->atomGroups_.index.data();
2094 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
2095 int npbcdim, matrix box,
2096 rvec cellsize_min, ivec npulse,
2097 t_forcerec *fr,
2098 rvec *cgcm_or_x,
2099 gmx_vsite_t *vsite,
2100 const gmx_mtop_t *mtop, gmx_localtop_t *ltop)
2102 gmx_bool bRCheckMB, bRCheck2B;
2103 real rc = -1;
2104 ivec rcheck;
2105 int d, nexcl;
2106 t_pbc pbc, *pbc_null = nullptr;
2108 if (debug)
2110 fprintf(debug, "Making local topology\n");
2113 dd_make_local_cgs(dd, &ltop->cgs);
2115 bRCheckMB = FALSE;
2116 bRCheck2B = FALSE;
2118 if (dd->reverse_top->bInterCGInteractions)
2120 /* We need to check to which cell bondeds should be assigned */
2121 rc = dd_cutoff_twobody(dd);
2122 if (debug)
2124 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2127 /* Should we check cg_cm distances when assigning bonded interactions? */
2128 for (d = 0; d < DIM; d++)
2130 rcheck[d] = FALSE;
2131 /* Only need to check for dimensions where the part of the box
2132 * that is not communicated is smaller than the cut-off.
2134 if (d < npbcdim && dd->nc[d] > 1 &&
2135 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2137 if (dd->nc[d] == 2)
2139 rcheck[d] = TRUE;
2140 bRCheckMB = TRUE;
2142 /* Check for interactions between two atoms,
2143 * where we can allow interactions up to the cut-off,
2144 * instead of up to the smallest cell dimension.
2146 bRCheck2B = TRUE;
2148 if (debug)
2150 fprintf(debug,
2151 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
2152 d, cellsize_min[d], d, rcheck[d], bRCheck2B);
2155 if (bRCheckMB || bRCheck2B)
2157 makeLocalAtomGroupFromAtom(dd);
2158 if (fr->bMolPBC)
2160 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
2162 else
2164 pbc_null = nullptr;
2169 dd->nbonded_local =
2170 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
2171 bRCheckMB, rcheck, bRCheck2B, rc,
2172 dd->localAtomGroupFromAtom.data(),
2173 pbc_null, cgcm_or_x,
2174 &ltop->idef, vsite,
2175 &ltop->excls, &nexcl);
2177 /* The ilist is not sorted yet,
2178 * we can only do this when we have the charge arrays.
2180 ltop->idef.ilsort = ilsortUNKNOWN;
2182 if (dd->reverse_top->bExclRequired)
2184 dd->nbonded_local += nexcl;
2187 ltop->atomtypes = mtop->atomtypes;
2190 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
2191 gmx_localtop_t *ltop)
2193 if (dd->reverse_top->ilsort == ilsortNO_FE)
2195 ltop->idef.ilsort = ilsortNO_FE;
2197 else
2199 gmx_sort_ilist_fe(&ltop->idef, mdatoms->chargeA, mdatoms->chargeB);
2203 gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global)
2205 gmx_localtop_t *top;
2206 int i;
2208 snew(top, 1);
2210 top->idef.ntypes = top_global->ffparams.ntypes;
2211 top->idef.atnr = top_global->ffparams.atnr;
2212 top->idef.functype = top_global->ffparams.functype;
2213 top->idef.iparams = top_global->ffparams.iparams;
2214 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
2215 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
2217 for (i = 0; i < F_NRE; i++)
2219 top->idef.il[i].iatoms = nullptr;
2220 top->idef.il[i].nalloc = 0;
2222 top->idef.ilsort = ilsortUNKNOWN;
2224 return top;
2227 void dd_init_local_state(gmx_domdec_t *dd,
2228 t_state *state_global, t_state *state_local)
2230 int buf[NITEM_DD_INIT_LOCAL_STATE];
2232 if (DDMASTER(dd))
2234 buf[0] = state_global->flags;
2235 buf[1] = state_global->ngtc;
2236 buf[2] = state_global->nnhpres;
2237 buf[3] = state_global->nhchainlength;
2238 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
2240 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2242 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
2243 init_dfhist_state(state_local, buf[4]);
2244 state_local->flags = buf[0];
2247 /*! \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 */
2248 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2250 int k;
2251 gmx_bool bFound;
2253 bFound = FALSE;
2254 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2256 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2257 if (link->a[k] == cg_gl_j)
2259 bFound = TRUE;
2262 if (!bFound)
2264 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2265 "Inconsistent allocation of link");
2266 /* Add this charge group link */
2267 if (link->index[cg_gl+1]+1 > link->nalloc_a)
2269 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2270 srenew(link->a, link->nalloc_a);
2272 link->a[link->index[cg_gl+1]] = cg_gl_j;
2273 link->index[cg_gl+1]++;
2277 /*! \brief Return a vector of the charge group index for all atoms */
2278 static std::vector<int> make_at2cg(const t_block &cgs)
2280 std::vector<int> at2cg(cgs.index[cgs.nr]);
2281 for (int cg = 0; cg < cgs.nr; cg++)
2283 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2285 at2cg[a] = cg;
2289 return at2cg;
2292 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
2293 cginfo_mb_t *cginfo_mb)
2295 gmx_bool bExclRequired;
2296 reverse_ilist_t ril, ril_intermol;
2297 t_blocka *link;
2298 cginfo_mb_t *cgi_mb;
2300 /* For each charge group make a list of other charge groups
2301 * in the system that a linked to it via bonded interactions
2302 * which are also stored in reverse_top.
2305 bExclRequired = dd->reverse_top->bExclRequired;
2307 if (mtop->bIntermolecularInteractions)
2309 if (ncg_mtop(mtop) < mtop->natoms)
2311 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.");
2314 t_atoms atoms;
2316 atoms.nr = mtop->natoms;
2317 atoms.atom = nullptr;
2319 make_reverse_ilist(mtop->intermolecular_ilist, &atoms,
2320 nullptr, FALSE, FALSE, FALSE, TRUE, &ril_intermol);
2323 snew(link, 1);
2324 snew(link->index, ncg_mtop(mtop)+1);
2325 link->nalloc_a = 0;
2326 link->a = nullptr;
2328 link->index[0] = 0;
2329 int cg_offset = 0;
2330 int ncgi = 0;
2331 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
2333 const gmx_molblock_t &molb = mtop->molblock[mb];
2334 if (molb.nmol == 0)
2336 continue;
2338 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2339 const t_block &cgs = molt.cgs;
2340 const t_blocka &excls = molt.excls;
2341 std::vector<int> a2c = make_at2cg(cgs);
2342 /* Make a reverse ilist in which the interactions are linked
2343 * to all atoms, not only the first atom as in gmx_reverse_top.
2344 * The constraints are discarded here.
2346 make_reverse_ilist(molt.ilist, &molt.atoms,
2347 nullptr, FALSE, FALSE, FALSE, TRUE, &ril);
2349 cgi_mb = &cginfo_mb[mb];
2351 int mol;
2352 for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb.nmol : 1); mol++)
2354 for (int cg = 0; cg < cgs.nr; cg++)
2356 int cg_gl = cg_offset + cg;
2357 link->index[cg_gl+1] = link->index[cg_gl];
2358 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2360 int i = ril.index[a];
2361 while (i < ril.index[a+1])
2363 int ftype = ril.il[i++];
2364 int nral = NRAL(ftype);
2365 /* Skip the ifunc index */
2366 i++;
2367 for (int j = 0; j < nral; j++)
2369 int aj = ril.il[i + j];
2370 if (a2c[aj] != cg)
2372 check_link(link, cg_gl, cg_offset+a2c[aj]);
2375 i += nral_rt(ftype);
2377 if (bExclRequired)
2379 /* Exclusions always go both ways */
2380 for (int j = excls.index[a]; j < excls.index[a + 1]; j++)
2382 int aj = excls.a[j];
2383 if (a2c[aj] != cg)
2385 check_link(link, cg_gl, cg_offset+a2c[aj]);
2390 if (mtop->bIntermolecularInteractions)
2392 int i = ril_intermol.index[a];
2393 while (i < ril_intermol.index[a+1])
2395 int ftype = ril_intermol.il[i++];
2396 int nral = NRAL(ftype);
2397 /* Skip the ifunc index */
2398 i++;
2399 for (int j = 0; j < nral; j++)
2401 /* Here we assume we have no charge groups;
2402 * this has been checked above.
2404 int aj = ril_intermol.il[i + j];
2405 check_link(link, cg_gl, aj);
2407 i += nral_rt(ftype);
2411 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2413 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2414 ncgi++;
2418 cg_offset += cgs.nr;
2420 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - cgs.nr];
2422 destroy_reverse_ilist(&ril);
2424 if (debug)
2426 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt.name, cgs.nr, nlink_mol);
2429 if (molb.nmol > mol)
2431 /* Copy the data for the rest of the molecules in this block */
2432 link->nalloc_a += (molb.nmol - mol)*nlink_mol;
2433 srenew(link->a, link->nalloc_a);
2434 for (; mol < molb.nmol; mol++)
2436 for (int cg = 0; cg < cgs.nr; cg++)
2438 int cg_gl = cg_offset + cg;
2439 link->index[cg_gl + 1] =
2440 link->index[cg_gl + 1 - cgs.nr] + nlink_mol;
2441 for (int j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2443 link->a[j] = link->a[j - nlink_mol] + cgs.nr;
2445 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2446 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2448 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2449 ncgi++;
2452 cg_offset += cgs.nr;
2457 if (mtop->bIntermolecularInteractions)
2459 destroy_reverse_ilist(&ril_intermol);
2462 if (debug)
2464 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2467 return link;
2470 typedef struct {
2471 real r2;
2472 int ftype;
2473 int a1;
2474 int a2;
2475 } bonded_distance_t;
2477 /*! \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 */
2478 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2479 bonded_distance_t *bd)
2481 if (r2 > bd->r2)
2483 bd->r2 = r2;
2484 bd->ftype = ftype;
2485 bd->a1 = a1;
2486 bd->a2 = a2;
2490 /*! \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 */
2491 static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
2492 const std::vector<int> &at2cg,
2493 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2494 bonded_distance_t *bd_2b,
2495 bonded_distance_t *bd_mb)
2497 for (int ftype = 0; ftype < F_NRE; ftype++)
2499 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2501 const t_ilist *il = &molt->ilist[ftype];
2502 int nral = NRAL(ftype);
2503 if (nral > 1)
2505 for (int i = 0; i < il->nr; i += 1+nral)
2507 for (int ai = 0; ai < nral; ai++)
2509 int cgi = at2cg[il->iatoms[i+1+ai]];
2510 for (int aj = ai + 1; aj < nral; aj++)
2512 int cgj = at2cg[il->iatoms[i+1+aj]];
2513 if (cgi != cgj)
2515 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2517 update_max_bonded_distance(rij2, ftype,
2518 il->iatoms[i+1+ai],
2519 il->iatoms[i+1+aj],
2520 (nral == 2) ? bd_2b : bd_mb);
2528 if (bExcl)
2530 const t_blocka *excls = &molt->excls;
2531 for (int ai = 0; ai < excls->nr; ai++)
2533 int cgi = at2cg[ai];
2534 for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2536 int cgj = at2cg[excls->a[j]];
2537 if (cgi != cgj)
2539 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2541 /* There is no function type for exclusions, use -1 */
2542 update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2549 /*! \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 */
2550 static void bonded_distance_intermol(const t_ilist *ilists_intermol,
2551 gmx_bool bBCheck,
2552 const rvec *x, int ePBC, const matrix box,
2553 bonded_distance_t *bd_2b,
2554 bonded_distance_t *bd_mb)
2556 t_pbc pbc;
2558 set_pbc(&pbc, ePBC, box);
2560 for (int ftype = 0; ftype < F_NRE; ftype++)
2562 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2564 const t_ilist *il = &ilists_intermol[ftype];
2565 int nral = NRAL(ftype);
2567 /* No nral>1 check here, since intermol interactions always
2568 * have nral>=2 (and the code is also correct for nral=1).
2570 for (int i = 0; i < il->nr; i += 1+nral)
2572 for (int ai = 0; ai < nral; ai++)
2574 int atom_i = il->iatoms[i + 1 + ai];
2576 for (int aj = ai + 1; aj < nral; aj++)
2578 rvec dx;
2579 real rij2;
2581 int atom_j = il->iatoms[i + 1 + aj];
2583 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2585 rij2 = norm2(dx);
2587 update_max_bonded_distance(rij2, ftype,
2588 atom_i, atom_j,
2589 (nral == 2) ? bd_2b : bd_mb);
2597 //! Returns whether \p molt has at least one virtual site
2598 static bool moltypeHasVsite(const gmx_moltype_t &molt)
2600 bool hasVsite = false;
2601 for (int i = 0; i < F_NRE; i++)
2603 if ((interaction_function[i].flags & IF_VSITE) &&
2604 molt.ilist[i].nr > 0)
2606 hasVsite = true;
2610 return hasVsite;
2613 //! Compute charge group centers of mass for molecule \p molt
2614 static void get_cgcm_mol(const gmx_moltype_t *molt,
2615 const gmx_ffparams_t *ffparams,
2616 int ePBC, t_graph *graph, const matrix box,
2617 const rvec *x, rvec *xs, rvec *cg_cm)
2619 int n, i;
2621 if (ePBC != epbcNONE)
2623 mk_mshift(nullptr, graph, ePBC, box, x);
2625 shift_x(graph, box, x, xs);
2626 /* By doing an extra mk_mshift the molecules that are broken
2627 * because they were e.g. imported from another software
2628 * will be made whole again. Such are the healing powers
2629 * of GROMACS.
2631 mk_mshift(nullptr, graph, ePBC, box, xs);
2633 else
2635 /* We copy the coordinates so the original coordinates remain
2636 * unchanged, just to be 100% sure that we do not affect
2637 * binary reproducibility of simulations.
2639 n = molt->cgs.index[molt->cgs.nr];
2640 for (i = 0; i < n; i++)
2642 copy_rvec(x[i], xs[i]);
2646 if (moltypeHasVsite(*molt))
2648 construct_vsites(nullptr, xs, 0.0, nullptr,
2649 ffparams->iparams, molt->ilist,
2650 epbcNONE, TRUE, nullptr, nullptr);
2653 calc_cgcm(nullptr, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2656 void dd_bonded_cg_distance(FILE *fplog,
2657 const gmx_mtop_t *mtop,
2658 const t_inputrec *ir,
2659 const rvec *x, const matrix box,
2660 gmx_bool bBCheck,
2661 real *r_2b, real *r_mb)
2663 gmx_bool bExclRequired;
2664 int at_offset;
2665 t_graph graph;
2666 rvec *xs, *cg_cm;
2667 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2668 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2670 bExclRequired = inputrecExclForces(ir);
2672 *r_2b = 0;
2673 *r_mb = 0;
2674 at_offset = 0;
2675 for (const gmx_molblock_t &molb : mtop->molblock)
2677 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2678 if (molt.cgs.nr == 1 || molb.nmol == 0)
2680 at_offset += molb.nmol*molt.atoms.nr;
2682 else
2684 if (ir->ePBC != epbcNONE)
2686 mk_graph_ilist(nullptr, molt.ilist, 0, molt.atoms.nr, FALSE, FALSE,
2687 &graph);
2690 std::vector<int> at2cg = make_at2cg(molt.cgs);
2691 snew(xs, molt.atoms.nr);
2692 snew(cg_cm, molt.cgs.nr);
2693 for (int mol = 0; mol < molb.nmol; mol++)
2695 get_cgcm_mol(&molt, &mtop->ffparams, ir->ePBC, &graph, box,
2696 x+at_offset, xs, cg_cm);
2698 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2699 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2701 bonded_cg_distance_mol(&molt, at2cg, bBCheck, bExclRequired, cg_cm,
2702 &bd_mol_2b, &bd_mol_mb);
2704 /* Process the mol data adding the atom index offset */
2705 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2706 at_offset + bd_mol_2b.a1,
2707 at_offset + bd_mol_2b.a2,
2708 &bd_2b);
2709 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2710 at_offset + bd_mol_mb.a1,
2711 at_offset + bd_mol_mb.a2,
2712 &bd_mb);
2714 at_offset += molt.atoms.nr;
2716 sfree(cg_cm);
2717 sfree(xs);
2718 if (ir->ePBC != epbcNONE)
2720 done_graph(&graph);
2725 if (mtop->bIntermolecularInteractions)
2727 if (ncg_mtop(mtop) < mtop->natoms)
2729 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.");
2732 bonded_distance_intermol(mtop->intermolecular_ilist,
2733 bBCheck,
2734 x, ir->ePBC, box,
2735 &bd_2b, &bd_mb);
2738 *r_2b = sqrt(bd_2b.r2);
2739 *r_mb = sqrt(bd_mb.r2);
2741 if (fplog && (*r_2b > 0 || *r_mb > 0))
2743 fprintf(fplog,
2744 "Initial maximum inter charge-group distances:\n");
2745 if (*r_2b > 0)
2747 fprintf(fplog,
2748 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2749 *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2750 bd_2b.a1 + 1, bd_2b.a2 + 1);
2752 if (*r_mb > 0)
2754 fprintf(fplog,
2755 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2756 *r_mb, interaction_function[bd_mb.ftype].longname,
2757 bd_mb.a1 + 1, bd_mb.a2 + 1);