Added trivial const qualifiers
[gromacs.git] / src / gromacs / domdec / domdec_topology.cpp
blob76fa2786dbcd89cd5b03bcda1c453bc402bcbc81
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 <cassert>
49 #include <cstdlib>
50 #include <cstring>
52 #include <algorithm>
53 #include <string>
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_network.h"
57 #include "gromacs/domdec/ga2la.h"
58 #include "gromacs/gmxlib/chargegroup.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdlib/forcerec.h"
62 #include "gromacs/mdlib/gmx_omp_nthreads.h"
63 #include "gromacs/mdlib/vsite.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/inputrec.h"
66 #include "gromacs/mdtypes/md_enums.h"
67 #include "gromacs/mdtypes/mdatom.h"
68 #include "gromacs/mdtypes/state.h"
69 #include "gromacs/pbcutil/mshift.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/topology/ifunc.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "gromacs/topology/topsort.h"
74 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/utility/exceptions.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/smalloc.h"
79 #include "gromacs/utility/strconvert.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
89 struct reverse_ilist_t
91 int *index; /* Index for each atom into il */
92 int *il; /* ftype|type|a0|...|an|ftype|... */
93 int numAtomsInMolecule; /* The number of atoms in this molecule */
96 typedef struct {
97 int a_start;
98 int a_end;
99 int natoms_mol;
100 int type;
101 } molblock_ind_t;
103 /*! \brief Struct for thread local work data for local topology generation */
104 typedef struct {
105 t_idef idef; /**< Partial local topology */
106 int **vsite_pbc; /**< vsite PBC structure */
107 int *vsite_pbc_nalloc; /**< Allocation sizes for vsite_pbc */
108 int nbonded; /**< The number of bondeds in this struct */
109 t_blocka excl; /**< List of exclusions */
110 int excl_count; /**< The total exclusion count for \p excl */
111 } thread_work_t;
113 /*! \brief Struct for the reverse topology: links bonded interactions to atomsx */
114 struct gmx_reverse_top_t
116 //! @cond Doxygen_Suppress
117 gmx_bool bExclRequired; /**< Do we require all exclusions to be assigned? */
118 int n_excl_at_max; /**< The maximum number of exclusions one atom can have */
119 gmx_bool bConstr; /**< Are there constraints in this revserse top? */
120 gmx_bool bSettle; /**< Are there settles in this revserse top? */
121 gmx_bool bBCheck; /**< All bonded interactions have to be assigned? */
122 gmx_bool bInterCGInteractions; /**< Are there bondeds/exclusions between charge-groups? */
123 reverse_ilist_t *ril_mt; /**< Reverse ilist for all moltypes */
124 int ril_mt_tot_size; /**< The size of ril_mt[?].index summed over all entries */
125 int ilsort; /**< The sorting state of bondeds for free energy */
126 molblock_ind_t *mbi; /**< molblock to global atom index for quick lookup of molblocks on atom index */
127 int nmolblock; /**< The number of molblocks, size of \p mbi */
129 gmx_bool bIntermolecularInteractions; /**< Do we have intermolecular interactions? */
130 reverse_ilist_t ril_intermol; /**< Intermolecular reverse ilist */
132 /* Work data structures for multi-threading */
133 int nthread; /**< The number of threads to be used */
134 thread_work_t *th_work; /**< Thread work array for local topology generation */
135 //! @endcond
138 /*! \brief Returns the number of atom entries for il in gmx_reverse_top_t */
139 static int nral_rt(int ftype)
141 int nral;
143 nral = NRAL(ftype);
144 if (interaction_function[ftype].flags & IF_VSITE)
146 /* With vsites the reverse topology contains
147 * two extra entries for PBC.
149 nral += 2;
152 return nral;
155 /*! \brief Return whether interactions of type \p ftype need to be assigned exactly once */
156 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
157 gmx_bool bConstr, gmx_bool bSettle)
159 return ((((interaction_function[ftype].flags & IF_BOND) != 0u) &&
160 ((interaction_function[ftype].flags & IF_VSITE) == 0u) &&
161 (bBCheck || ((interaction_function[ftype].flags & IF_LIMZERO) == 0u))) ||
162 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
163 (bSettle && ftype == F_SETTLE));
166 /*! \brief Print a header on error messages */
167 static void print_error_header(FILE *fplog, const char *moltypename, int nprint)
169 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
170 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
171 fprintf(fplog,
172 "the first %d missing interactions, except for exclusions:\n",
173 nprint);
174 fprintf(stderr,
175 "the first %d missing interactions, except for exclusions:\n",
176 nprint);
179 /*! \brief Help print error output when interactions are missing */
180 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
181 const gmx_reverse_top_t *rt,
182 const char *moltypename,
183 const reverse_ilist_t *ril,
184 int a_start, int a_end,
185 int nat_mol, int nmol,
186 const t_idef *idef)
188 int *assigned;
189 int nril_mol = ril->index[nat_mol];
190 snew(assigned, nmol*nril_mol);
192 gmx::ArrayRef<const int> gatindex = cr->dd->globalAtomIndices;
193 for (int ftype = 0; ftype < F_NRE; ftype++)
195 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
197 int nral = NRAL(ftype);
198 const t_ilist *il = &idef->il[ftype];
199 const t_iatom *ia = il->iatoms;
200 for (int i = 0; i < il->nr; i += 1+nral)
202 int a0 = gatindex[ia[1]];
203 /* Check if this interaction is in
204 * the currently checked molblock.
206 if (a0 >= a_start && a0 < a_end)
208 int mol = (a0 - a_start)/nat_mol;
209 int a0_mol = (a0 - a_start) - mol*nat_mol;
210 int j_mol = ril->index[a0_mol];
211 bool found = false;
212 while (j_mol < ril->index[a0_mol+1] && !found)
214 int j = mol*nril_mol + j_mol;
215 int ftype_j = ril->il[j_mol];
216 /* Here we need to check if this interaction has
217 * not already been assigned, since we could have
218 * multiply defined interactions.
220 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
221 assigned[j] == 0)
223 /* Check the atoms */
224 found = true;
225 for (int a = 0; a < nral; a++)
227 if (gatindex[ia[1+a]] !=
228 a_start + mol*nat_mol + ril->il[j_mol+2+a])
230 found = false;
233 if (found)
235 assigned[j] = 1;
238 j_mol += 2 + nral_rt(ftype_j);
240 if (!found)
242 gmx_incons("Some interactions seem to be assigned multiple times");
245 ia += 1 + nral;
250 gmx_sumi(nmol*nril_mol, assigned, cr);
252 int nprint = 10;
253 int i = 0;
254 for (int mol = 0; mol < nmol; mol++)
256 int j_mol = 0;
257 while (j_mol < nril_mol)
259 int ftype = ril->il[j_mol];
260 int nral = NRAL(ftype);
261 int j = mol*nril_mol + j_mol;
262 if (assigned[j] == 0 &&
263 !(interaction_function[ftype].flags & IF_VSITE))
265 if (DDMASTER(cr->dd))
267 if (i == 0)
269 print_error_header(fplog, moltypename, nprint);
271 fprintf(fplog, "%20s atoms",
272 interaction_function[ftype].longname);
273 fprintf(stderr, "%20s atoms",
274 interaction_function[ftype].longname);
275 int a;
276 for (a = 0; a < nral; a++)
278 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
279 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
281 while (a < 4)
283 fprintf(fplog, " ");
284 fprintf(stderr, " ");
285 a++;
287 fprintf(fplog, " global");
288 fprintf(stderr, " global");
289 for (a = 0; a < nral; a++)
291 fprintf(fplog, "%6d",
292 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
293 fprintf(stderr, "%6d",
294 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
296 fprintf(fplog, "\n");
297 fprintf(stderr, "\n");
299 i++;
300 if (i >= nprint)
302 break;
305 j_mol += 2 + nral_rt(ftype);
309 sfree(assigned);
312 /*! \brief Help print error output when interactions are missing */
313 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
314 const gmx_mtop_t *mtop,
315 const t_idef *idef)
317 const gmx_reverse_top_t *rt = cr->dd->reverse_top;
319 /* Print the atoms in the missing interactions per molblock */
320 int a_end = 0;
321 for (const gmx_molblock_t &molb : mtop->molblock)
323 const gmx_moltype_t &moltype = mtop->moltype[molb.type];
324 int a_start = a_end;
325 a_end = a_start + molb.nmol*moltype.atoms.nr;
327 print_missing_interactions_mb(fplog, cr, rt,
328 *(moltype.name),
329 &rt->ril_mt[molb.type],
330 a_start, a_end, moltype.atoms.nr,
331 molb.nmol,
332 idef);
336 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
337 int local_count,
338 const gmx_mtop_t *top_global,
339 const gmx_localtop_t *top_local,
340 const t_state *state_local)
342 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
343 int ftype, nral;
344 char buf[STRLEN];
345 gmx_domdec_t *dd;
347 dd = cr->dd;
349 if (fplog)
351 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
352 fflush(fplog);
355 ndiff_tot = local_count - dd->nbonded_global;
357 for (ftype = 0; ftype < F_NRE; ftype++)
359 nral = NRAL(ftype);
360 cl[ftype] = top_local->idef.il[ftype].nr/(1+nral);
363 gmx_sumi(F_NRE, cl, cr);
365 if (DDMASTER(dd))
367 if (fplog)
369 fprintf(fplog, "\nA list of missing interactions:\n");
371 fprintf(stderr, "\nA list of missing interactions:\n");
372 rest_global = dd->nbonded_global;
373 rest_local = local_count;
374 for (ftype = 0; ftype < F_NRE; ftype++)
376 /* In the reverse and local top all constraints are merged
377 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
378 * and add these constraints when doing F_CONSTR.
380 if (((interaction_function[ftype].flags & IF_BOND) &&
381 (dd->reverse_top->bBCheck
382 || !(interaction_function[ftype].flags & IF_LIMZERO)))
383 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
384 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
386 n = gmx_mtop_ftype_count(top_global, ftype);
387 if (ftype == F_CONSTR)
389 n += gmx_mtop_ftype_count(top_global, F_CONSTRNC);
391 ndiff = cl[ftype] - n;
392 if (ndiff != 0)
394 sprintf(buf, "%20s of %6d missing %6d",
395 interaction_function[ftype].longname, n, -ndiff);
396 if (fplog)
398 fprintf(fplog, "%s\n", buf);
400 fprintf(stderr, "%s\n", buf);
402 rest_global -= n;
403 rest_local -= cl[ftype];
407 ndiff = rest_local - rest_global;
408 if (ndiff != 0)
410 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
411 rest_global, -ndiff);
412 if (fplog)
414 fprintf(fplog, "%s\n", buf);
416 fprintf(stderr, "%s\n", buf);
420 print_missing_interactions_atoms(fplog, cr, top_global, &top_local->idef);
421 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
422 -1, as_rvec_array(state_local->x.data()), state_local->box);
424 std::string errorMessage;
426 if (ndiff_tot > 0)
428 errorMessage = "One or more interactions were assigned to multiple domains of the domain decompostion. Please report this bug.";
430 else
432 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));
434 gmx_fatal_collective(FARGS, cr->mpi_comm_mygroup, MASTER(cr), "%s", errorMessage.c_str());
437 /*! \brief Return global topology molecule information for global atom index \p i_gl */
438 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
439 int *mb, int *mt, int *mol, int *i_mol)
441 molblock_ind_t *mbi = rt->mbi;
442 int start = 0;
443 int end = rt->nmolblock; /* exclusive */
444 int mid;
446 /* binary search for molblock_ind */
447 while (TRUE)
449 mid = (start+end)/2;
450 if (i_gl >= mbi[mid].a_end)
452 start = mid+1;
454 else if (i_gl < mbi[mid].a_start)
456 end = mid;
458 else
460 break;
464 *mb = mid;
465 mbi += mid;
467 *mt = mbi->type;
468 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
469 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
472 /*! \brief Count the exclusions for all atoms in \p cgs */
473 static void count_excls(const t_block *cgs, const t_blocka *excls,
474 int *n_excl, int *n_intercg_excl, int *n_excl_at_max)
476 int cg, at0, at1, at, excl, atj;
478 *n_excl = 0;
479 *n_intercg_excl = 0;
480 *n_excl_at_max = 0;
481 for (cg = 0; cg < cgs->nr; cg++)
483 at0 = cgs->index[cg];
484 at1 = cgs->index[cg+1];
485 for (at = at0; at < at1; at++)
487 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
489 atj = excls->a[excl];
490 if (atj > at)
492 (*n_excl)++;
493 if (atj < at0 || atj >= at1)
495 (*n_intercg_excl)++;
500 *n_excl_at_max = std::max(*n_excl_at_max,
501 excls->index[at+1] - excls->index[at]);
506 /*! \brief Run the reverse ilist generation and store it in r_il when \p bAssign = TRUE */
507 static int low_make_reverse_ilist(const t_ilist *il_mt, const t_atom *atom,
508 const int * const * vsite_pbc,
509 int *count,
510 gmx_bool bConstr, gmx_bool bSettle,
511 gmx_bool bBCheck,
512 const int *r_index, int *r_il,
513 gmx_bool bLinkToAllAtoms,
514 gmx_bool bAssign)
516 int ftype, nral, i, j, nlink, link;
517 const t_ilist *il;
518 const t_iatom *ia;
519 int a;
520 int nint;
521 gmx_bool bVSite;
523 nint = 0;
524 for (ftype = 0; ftype < F_NRE; ftype++)
526 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
527 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
528 (bSettle && ftype == F_SETTLE))
530 bVSite = ((interaction_function[ftype].flags & IF_VSITE) != 0u);
531 nral = NRAL(ftype);
532 il = &il_mt[ftype];
533 for (i = 0; i < il->nr; i += 1+nral)
535 ia = il->iatoms + i;
536 if (bLinkToAllAtoms)
538 if (bVSite)
540 /* We don't need the virtual sites for the cg-links */
541 nlink = 0;
543 else
545 nlink = nral;
548 else
550 /* Couple to the first atom in the interaction */
551 nlink = 1;
553 for (link = 0; link < nlink; link++)
555 a = ia[1+link];
556 if (bAssign)
558 assert(r_il); //with bAssign not allowed to be null
559 assert(r_index);
560 r_il[r_index[a]+count[a]] =
561 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
562 r_il[r_index[a]+count[a]+1] = ia[0];
563 for (j = 1; j < 1+nral; j++)
565 /* Store the molecular atom number */
566 r_il[r_index[a]+count[a]+1+j] = ia[j];
569 if (interaction_function[ftype].flags & IF_VSITE)
571 if (bAssign)
573 /* Add an entry to iatoms for storing
574 * which of the constructing atoms are
575 * vsites again.
577 r_il[r_index[a]+count[a]+2+nral] = 0;
578 for (j = 2; j < 1+nral; j++)
580 if (atom[ia[j]].ptype == eptVSite)
582 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
585 /* Store vsite pbc atom in a second extra entry */
586 r_il[r_index[a]+count[a]+2+nral+1] =
587 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
590 else
592 /* We do not count vsites since they are always
593 * uniquely assigned and can be assigned
594 * to multiple nodes with recursive vsites.
596 if (bBCheck ||
597 !(interaction_function[ftype].flags & IF_LIMZERO))
599 nint++;
602 count[a] += 2 + nral_rt(ftype);
608 return nint;
611 /*! \brief Make the reverse ilist: a list of bonded interactions linked to atoms */
612 static int make_reverse_ilist(const t_ilist *ilist,
613 const t_atoms *atoms,
614 const int * const * vsite_pbc,
615 gmx_bool bConstr, gmx_bool bSettle,
616 gmx_bool bBCheck,
617 gmx_bool bLinkToAllAtoms,
618 reverse_ilist_t *ril_mt)
620 int nat_mt, *count, i, nint_mt;
622 /* Count the interactions */
623 nat_mt = atoms->nr;
624 snew(count, nat_mt);
625 low_make_reverse_ilist(ilist, atoms->atom, vsite_pbc,
626 count,
627 bConstr, bSettle, bBCheck, nullptr, nullptr,
628 bLinkToAllAtoms, FALSE);
630 snew(ril_mt->index, nat_mt+1);
631 ril_mt->index[0] = 0;
632 for (i = 0; i < nat_mt; i++)
634 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
635 count[i] = 0;
637 snew(ril_mt->il, ril_mt->index[nat_mt]);
639 /* Store the interactions */
640 nint_mt =
641 low_make_reverse_ilist(ilist, atoms->atom, vsite_pbc,
642 count,
643 bConstr, bSettle, bBCheck,
644 ril_mt->index, ril_mt->il,
645 bLinkToAllAtoms, TRUE);
647 sfree(count);
649 ril_mt->numAtomsInMolecule = atoms->nr;
651 return nint_mt;
654 /*! \brief Destroys a reverse_ilist_t struct */
655 static void destroy_reverse_ilist(reverse_ilist_t *ril)
657 sfree(ril->index);
658 sfree(ril->il);
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)
667 gmx_reverse_top_t *rt;
668 int *nint_mt;
669 int thread;
671 snew(rt, 1);
673 /* Should we include constraints (for SHAKE) in rt? */
674 rt->bConstr = bConstr;
675 rt->bSettle = bSettle;
676 rt->bBCheck = bBCheck;
678 rt->bInterCGInteractions = mtop->bIntermolecularInteractions;
679 snew(nint_mt, mtop->moltype.size());
680 snew(rt->ril_mt, mtop->moltype.size());
681 rt->ril_mt_tot_size = 0;
682 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
684 const gmx_moltype_t &molt = mtop->moltype[mt];
685 if (molt.cgs.nr > 1)
687 rt->bInterCGInteractions = TRUE;
690 /* Make the atom to interaction list for this molecule type */
691 nint_mt[mt] =
692 make_reverse_ilist(molt.ilist, &molt.atoms,
693 vsite_pbc_molt ? vsite_pbc_molt[mt] : nullptr,
694 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
695 &rt->ril_mt[mt]);
697 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt.atoms.nr];
699 if (debug)
701 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
704 *nint = 0;
705 for (const gmx_molblock_t &molblock : mtop->molblock)
707 *nint += molblock.nmol*nint_mt[molblock.type];
709 sfree(nint_mt);
711 /* Make an intermolecular reverse top, if necessary */
712 rt->bIntermolecularInteractions = mtop->bIntermolecularInteractions;
713 if (rt->bIntermolecularInteractions)
715 t_atoms atoms_global;
717 rt->ril_intermol.index = nullptr;
718 rt->ril_intermol.il = nullptr;
720 atoms_global.nr = mtop->natoms;
721 atoms_global.atom = nullptr; /* Only used with virtual sites */
723 *nint +=
724 make_reverse_ilist(mtop->intermolecular_ilist, &atoms_global,
725 nullptr,
726 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
727 &rt->ril_intermol);
730 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
732 rt->ilsort = ilsortFE_UNSORTED;
734 else
736 rt->ilsort = ilsortNO_FE;
739 /* Make a molblock index for fast searching */
740 snew(rt->mbi, mtop->molblock.size());
741 rt->nmolblock = mtop->molblock.size();
742 int i = 0;
743 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
745 const gmx_molblock_t &molb = mtop->molblock[mb];
746 int numAtomsPerMol = mtop->moltype[molb.type].atoms.nr;
747 rt->mbi[mb].a_start = i;
748 i += molb.nmol*numAtomsPerMol;
749 rt->mbi[mb].a_end = i;
750 rt->mbi[mb].natoms_mol = numAtomsPerMol;
751 rt->mbi[mb].type = molb.type;
754 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
755 snew(rt->th_work, rt->nthread);
756 if (vsite_pbc_molt != nullptr)
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);
765 return rt;
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)
773 if (fplog)
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 : nullptr,
786 !dd->bInterCGcons, !dd->bInterCGsettles,
787 bBCheck, &dd->nbonded_global);
789 gmx_reverse_top_t *rt = dd->reverse_top;
791 /* With the Verlet scheme, exclusions are handled in the non-bonded
792 * kernels and only exclusions inside the cut-off lead to exclusion
793 * forces. Since each atom pair is treated at most once in the non-bonded
794 * kernels, it doesn't matter if the exclusions for the same atom pair
795 * appear multiple times in the exclusion list. In contrast, the, old,
796 * group cut-off scheme loops over a list of exclusions, so there each
797 * excluded pair should appear exactly once.
799 rt->bExclRequired = (ir->cutoff_scheme == ecutsGROUP &&
800 inputrecExclForces(ir));
802 int nexcl = 0;
803 dd->n_intercg_excl = 0;
804 rt->n_excl_at_max = 0;
805 for (const gmx_molblock_t &molb : mtop->molblock)
807 int n_excl_mol, n_excl_icg, n_excl_at_max;
809 const gmx_moltype_t &molt = mtop->moltype[molb.type];
810 count_excls(&molt.cgs, &molt.excls,
811 &n_excl_mol, &n_excl_icg, &n_excl_at_max);
812 nexcl += molb.nmol*n_excl_mol;
813 dd->n_intercg_excl += molb.nmol*n_excl_icg;
814 rt->n_excl_at_max = std::max(rt->n_excl_at_max, n_excl_at_max);
816 if (rt->bExclRequired)
818 dd->nbonded_global += nexcl;
819 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
821 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
822 "will use an extra communication step for exclusion forces for %s\n",
823 dd->n_intercg_excl, eel_names[ir->coulombtype]);
827 if (vsite && vsite->n_intercg_vsite > 0)
829 if (fplog)
831 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
832 "will an extra communication step for selected coordinates and forces\n",
833 vsite->n_intercg_vsite);
835 init_domdec_vsites(dd, vsite->n_intercg_vsite);
838 if (dd->bInterCGcons || dd->bInterCGsettles)
840 init_domdec_constraints(dd, mtop);
842 if (fplog)
844 fprintf(fplog, "\n");
848 /*! \brief Store a vsite interaction at the end of \p il
850 * This routine is very similar to add_ifunc, but vsites interactions
851 * have more work to do than other kinds of interactions, and the
852 * complex way nral (and thus vector contents) depends on ftype
853 * confuses static analysis tools unless we fuse the vsite
854 * atom-indexing organization code with the ifunc-adding code, so that
855 * they can see that nral is the same value. */
856 static inline void
857 add_ifunc_for_vsites(t_iatom *tiatoms, const gmx_ga2la_t &ga2la,
858 int nral, gmx_bool bHomeA,
859 int a, int a_gl, int a_mol,
860 const t_iatom *iatoms,
861 t_ilist *il)
863 t_iatom *liatoms;
865 if (il->nr+1+nral > il->nalloc)
867 il->nalloc = over_alloc_large(il->nr+1+nral);
868 srenew(il->iatoms, il->nalloc);
870 liatoms = il->iatoms + il->nr;
871 il->nr += 1 + nral;
873 /* Copy the type */
874 tiatoms[0] = iatoms[0];
876 if (bHomeA)
878 /* We know the local index of the first atom */
879 tiatoms[1] = a;
881 else
883 /* Convert later in make_local_vsites */
884 tiatoms[1] = -a_gl - 1;
887 for (int k = 2; k < 1+nral; k++)
889 int ak_gl = a_gl + iatoms[k] - a_mol;
890 if (const int *homeIndex = ga2la.findHome(ak_gl))
892 tiatoms[k] = *homeIndex;
894 else
896 /* Copy the global index, convert later in make_local_vsites */
897 tiatoms[k] = -(ak_gl + 1);
899 // Note that ga2la_get_home always sets the third parameter if
900 // it returns TRUE
902 for (int k = 0; k < 1+nral; k++)
904 liatoms[k] = tiatoms[k];
908 /*! \brief Store a bonded interaction at the end of \p il */
909 static inline void add_ifunc(int nral, const t_iatom *tiatoms, t_ilist *il)
911 t_iatom *liatoms;
912 int k;
914 if (il->nr+1+nral > il->nalloc)
916 il->nalloc = over_alloc_large(il->nr+1+nral);
917 srenew(il->iatoms, il->nalloc);
919 liatoms = il->iatoms + il->nr;
920 for (k = 0; k <= nral; k++)
922 liatoms[k] = tiatoms[k];
924 il->nr += 1 + nral;
927 /*! \brief Store a position restraint in idef and iatoms, complex because the parameters are different for each entry */
928 static void add_posres(int mol, int a_mol, int numAtomsInMolecule,
929 const gmx_molblock_t *molb,
930 t_iatom *iatoms, const t_iparams *ip_in,
931 t_idef *idef)
933 int n, a_molb;
934 t_iparams *ip;
936 /* This position restraint has not been added yet,
937 * so it's index is the current number of position restraints.
939 n = idef->il[F_POSRES].nr/2;
940 if (n+1 > idef->iparams_posres_nalloc)
942 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
943 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
945 ip = &idef->iparams_posres[n];
946 /* Copy the force constants */
947 *ip = ip_in[iatoms[0]];
949 /* Get the position restraint coordinates from the molblock */
950 a_molb = mol*numAtomsInMolecule + a_mol;
951 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
952 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
953 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
954 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
955 if (!molb->posres_xB.empty())
957 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
958 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
959 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
961 else
963 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
964 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
965 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
967 /* Set the parameter index for idef->iparams_posre */
968 iatoms[0] = n;
971 /*! \brief Store a flat-bottomed position restraint in idef and iatoms, complex because the parameters are different for each entry */
972 static void add_fbposres(int mol, int a_mol, int numAtomsInMolecule,
973 const gmx_molblock_t *molb,
974 t_iatom *iatoms, const t_iparams *ip_in,
975 t_idef *idef)
977 int n, a_molb;
978 t_iparams *ip;
980 /* This flat-bottom position restraint has not been added yet,
981 * so it's index is the current number of position restraints.
983 n = idef->il[F_FBPOSRES].nr/2;
984 if (n+1 > idef->iparams_fbposres_nalloc)
986 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
987 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
989 ip = &idef->iparams_fbposres[n];
990 /* Copy the force constants */
991 *ip = ip_in[iatoms[0]];
993 /* Get the position restraint coordinats from the molblock */
994 a_molb = mol*numAtomsInMolecule + a_mol;
995 GMX_ASSERT(a_molb < static_cast<int>(molb->posres_xA.size()), "We need a sufficient number of position restraint coordinates");
996 /* Take reference positions from A position of normal posres */
997 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
998 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
999 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
1001 /* Note: no B-type for flat-bottom posres */
1003 /* Set the parameter index for idef->iparams_posre */
1004 iatoms[0] = n;
1007 /*! \brief Store a virtual site interaction, complex because of PBC and recursion */
1008 static void add_vsite(const gmx_ga2la_t &ga2la,
1009 const int *index, const int *rtil,
1010 int ftype, int nral,
1011 gmx_bool bHomeA, int a, int a_gl, int a_mol,
1012 const t_iatom *iatoms,
1013 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
1015 int k, vsi, pbc_a_mol;
1016 t_iatom tiatoms[1+MAXATOMLIST];
1017 int j, ftype_r, nral_r;
1019 /* Add this interaction to the local topology */
1020 add_ifunc_for_vsites(tiatoms, ga2la, nral, bHomeA, a, a_gl, a_mol, iatoms, &idef->il[ftype]);
1022 if (vsite_pbc)
1024 vsi = idef->il[ftype].nr/(1+nral) - 1;
1025 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
1027 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
1028 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
1030 if (bHomeA)
1032 pbc_a_mol = iatoms[1+nral+1];
1033 if (pbc_a_mol < 0)
1035 /* The pbc flag is one of the following two options:
1036 * -2: vsite and all constructing atoms are within the same cg, no pbc
1037 * -1: vsite and its first constructing atom are in the same cg, do pbc
1039 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
1041 else
1043 /* Set the pbc atom for this vsite so we can make its pbc
1044 * identical to the rest of the atoms in its charge group.
1045 * Since the order of the atoms does not change within a charge
1046 * group, we do not need the global to local atom index.
1048 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
1051 else
1053 /* This vsite is non-home (required for recursion),
1054 * and therefore there is no charge group to match pbc with.
1055 * But we always turn on full_pbc to assure that higher order
1056 * recursion works correctly.
1058 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
1062 if (iatoms[1+nral])
1064 /* Check for recursion */
1065 for (k = 2; k < 1+nral; k++)
1067 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
1069 /* This construction atoms is a vsite and not a home atom */
1070 if (gmx_debug_at)
1072 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
1074 /* Find the vsite construction */
1076 /* Check all interactions assigned to this atom */
1077 j = index[iatoms[k]];
1078 while (j < index[iatoms[k]+1])
1080 ftype_r = rtil[j++];
1081 nral_r = NRAL(ftype_r);
1082 if (interaction_function[ftype_r].flags & IF_VSITE)
1084 /* Add this vsite (recursion) */
1085 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
1086 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
1087 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
1088 j += 1 + nral_r + 2;
1090 else
1092 j += 1 + nral_r;
1100 /*! \brief Build the index that maps each local atom to its local atom group */
1101 static void makeLocalAtomGroupsFromAtoms(gmx_domdec_t *dd)
1103 const gmx::RangePartitioning &atomGrouping = dd->atomGrouping();
1105 dd->localAtomGroupFromAtom.clear();
1107 for (size_t g = 0; g < dd->globalAtomGroupIndices.size(); g++)
1109 for (int gmx_unused a : atomGrouping.block(g))
1111 dd->localAtomGroupFromAtom.push_back(g);
1116 /*! \brief Returns the squared distance between charge groups \p i and \p j */
1117 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1119 rvec dx;
1121 if (pbc_null)
1123 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1125 else
1127 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1130 return norm2(dx);
1133 /*! \brief Append t_blocka block structures 1 to nsrc in src to *dest */
1134 static void combine_blocka(t_blocka *dest, const thread_work_t *src, int nsrc)
1136 int ni, na, s, i;
1138 ni = src[nsrc-1].excl.nr;
1139 na = 0;
1140 for (s = 0; s < nsrc; s++)
1142 na += src[s].excl.nra;
1144 if (ni + 1 > dest->nalloc_index)
1146 dest->nalloc_index = over_alloc_large(ni+1);
1147 srenew(dest->index, dest->nalloc_index);
1149 if (dest->nra + na > dest->nalloc_a)
1151 dest->nalloc_a = over_alloc_large(dest->nra+na);
1152 srenew(dest->a, dest->nalloc_a);
1154 for (s = 1; s < nsrc; s++)
1156 for (i = dest->nr+1; i < src[s].excl.nr+1; i++)
1158 dest->index[i] = dest->nra + src[s].excl.index[i];
1160 for (i = 0; i < src[s].excl.nra; i++)
1162 dest->a[dest->nra+i] = src[s].excl.a[i];
1164 dest->nr = src[s].excl.nr;
1165 dest->nra += src[s].excl.nra;
1169 /*! \brief Append t_idef structures 1 to nsrc in src to *dest,
1170 * virtual sites need special attention, as pbc info differs per vsite.
1172 static void combine_idef(t_idef *dest, const thread_work_t *src, int nsrc,
1173 gmx_vsite_t *vsite)
1175 int ftype;
1177 for (ftype = 0; ftype < F_NRE; ftype++)
1179 int n, s;
1181 n = 0;
1182 for (s = 1; s < nsrc; s++)
1184 n += src[s].idef.il[ftype].nr;
1186 if (n > 0)
1188 t_ilist *ild;
1190 ild = &dest->il[ftype];
1192 if (ild->nr + n > ild->nalloc)
1194 ild->nalloc = over_alloc_large(ild->nr+n);
1195 srenew(ild->iatoms, ild->nalloc);
1198 gmx_bool vpbc;
1199 int nral1 = 0, ftv = 0;
1201 vpbc = (((interaction_function[ftype].flags & IF_VSITE) != 0u) &&
1202 vsite->vsite_pbc_loc != nullptr);
1203 if (vpbc)
1205 nral1 = 1 + NRAL(ftype);
1206 ftv = ftype - F_VSITE2;
1207 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1209 vsite->vsite_pbc_loc_nalloc[ftv] =
1210 over_alloc_large((ild->nr + n)/nral1);
1211 srenew(vsite->vsite_pbc_loc[ftv],
1212 vsite->vsite_pbc_loc_nalloc[ftv]);
1216 for (s = 1; s < nsrc; s++)
1218 const t_ilist *ils;
1219 int i;
1221 ils = &src[s].idef.il[ftype];
1222 for (i = 0; i < ils->nr; i++)
1224 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1226 if (vpbc)
1228 for (i = 0; i < ils->nr; i += nral1)
1230 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1231 src[s].vsite_pbc[ftv][i/nral1];
1235 ild->nr += ils->nr;
1238 /* Position restraints need an additional treatment */
1239 if (ftype == F_POSRES || ftype == F_FBPOSRES)
1241 int nposres = dest->il[ftype].nr/2;
1242 // TODO: Simplify this code using std::vector
1243 t_iparams * &iparams_dest = (ftype == F_POSRES ? dest->iparams_posres : dest->iparams_fbposres);
1244 int &posres_nalloc = (ftype == F_POSRES ? dest->iparams_posres_nalloc : dest->iparams_fbposres_nalloc);
1245 if (nposres > posres_nalloc)
1247 posres_nalloc = over_alloc_large(nposres);
1248 srenew(iparams_dest, posres_nalloc);
1251 /* Set nposres to the number of original position restraints in dest */
1252 for (int s = 1; s < nsrc; s++)
1254 nposres -= src[s].idef.il[ftype].nr/2;
1257 for (int s = 1; s < nsrc; s++)
1259 const t_iparams *iparams_src = (ftype == F_POSRES ? src[s].idef.iparams_posres : src[s].idef.iparams_fbposres);
1261 for (int i = 0; i < src[s].idef.il[ftype].nr/2; i++)
1263 /* Correct the index into iparams_posres */
1264 dest->il[ftype].iatoms[nposres*2] = nposres;
1265 /* Copy the position restraint force parameters */
1266 iparams_dest[nposres] = iparams_src[i];
1267 nposres++;
1275 /*! \brief Check and when available assign bonded interactions for local atom i
1277 static inline void
1278 check_assign_interactions_atom(int i, int i_gl,
1279 int mol, int i_mol,
1280 int numAtomsInMolecule,
1281 const int *index, const int *rtil,
1282 gmx_bool bInterMolInteractions,
1283 int ind_start, int ind_end,
1284 const gmx_domdec_t *dd,
1285 const gmx_domdec_zones_t *zones,
1286 const gmx_molblock_t *molb,
1287 gmx_bool bRCheckMB, const ivec rcheck, gmx_bool bRCheck2B,
1288 real rc2,
1289 int *la2lc,
1290 t_pbc *pbc_null,
1291 rvec *cg_cm,
1292 const t_iparams *ip_in,
1293 t_idef *idef,
1294 int **vsite_pbc, int *vsite_pbc_nalloc,
1295 int iz,
1296 gmx_bool bBCheck,
1297 int *nbonded_local)
1299 int j;
1301 j = ind_start;
1302 while (j < ind_end)
1304 int ftype;
1305 const t_iatom *iatoms;
1306 int nral;
1307 t_iatom tiatoms[1 + MAXATOMLIST];
1309 ftype = rtil[j++];
1310 iatoms = rtil + j;
1311 nral = NRAL(ftype);
1312 if (ftype == F_SETTLE)
1314 /* Settles are only in the reverse top when they
1315 * operate within a charge group. So we can assign
1316 * them without checks. We do this only for performance
1317 * reasons; it could be handled by the code below.
1319 if (iz == 0)
1321 /* Home zone: add this settle to the local topology */
1322 tiatoms[0] = iatoms[0];
1323 tiatoms[1] = i;
1324 tiatoms[2] = i + iatoms[2] - iatoms[1];
1325 tiatoms[3] = i + iatoms[3] - iatoms[1];
1326 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1327 (*nbonded_local)++;
1329 j += 1 + nral;
1331 else if (interaction_function[ftype].flags & IF_VSITE)
1333 assert(!bInterMolInteractions);
1334 /* The vsite construction goes where the vsite itself is */
1335 if (iz == 0)
1337 add_vsite(*dd->ga2la, index, rtil, ftype, nral,
1338 TRUE, i, i_gl, i_mol,
1339 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1341 j += 1 + nral + 2;
1343 else
1345 gmx_bool bUse;
1347 /* Copy the type */
1348 tiatoms[0] = iatoms[0];
1350 if (nral == 1)
1352 assert(!bInterMolInteractions);
1353 /* Assign single-body interactions to the home zone */
1354 if (iz == 0)
1356 bUse = TRUE;
1357 tiatoms[1] = i;
1358 if (ftype == F_POSRES)
1360 add_posres(mol, i_mol, numAtomsInMolecule,
1361 molb, tiatoms, ip_in, idef);
1363 else if (ftype == F_FBPOSRES)
1365 add_fbposres(mol, i_mol, numAtomsInMolecule,
1366 molb, tiatoms, ip_in, idef);
1369 else
1371 bUse = FALSE;
1374 else if (nral == 2)
1376 /* This is a two-body interaction, we can assign
1377 * analogous to the non-bonded assignments.
1379 int k_gl;
1381 if (!bInterMolInteractions)
1383 /* Get the global index using the offset in the molecule */
1384 k_gl = i_gl + iatoms[2] - i_mol;
1386 else
1388 k_gl = iatoms[2];
1390 if (const auto *entry = dd->ga2la->find(k_gl))
1392 int kz = entry->cell;
1393 if (kz >= zones->n)
1395 kz -= zones->n;
1397 /* Check zone interaction assignments */
1398 bUse = ((iz < zones->nizone &&
1399 iz <= kz &&
1400 kz >= zones->izone[iz].j0 &&
1401 kz < zones->izone[iz].j1) ||
1402 (kz < zones->nizone &&
1403 iz > kz &&
1404 iz >= zones->izone[kz].j0 &&
1405 iz < zones->izone[kz].j1));
1406 if (bUse)
1408 tiatoms[1] = i;
1409 tiatoms[2] = entry->la;
1410 /* If necessary check the cgcm distance */
1411 if (bRCheck2B &&
1412 dd_dist2(pbc_null, cg_cm, la2lc,
1413 tiatoms[1], tiatoms[2]) >= rc2)
1415 bUse = FALSE;
1419 else
1421 bUse = false;
1424 else
1426 /* Assign this multi-body bonded interaction to
1427 * the local node if we have all the atoms involved
1428 * (local or communicated) and the minimum zone shift
1429 * in each dimension is zero, for dimensions
1430 * with 2 DD cells an extra check may be necessary.
1432 ivec k_zero, k_plus;
1433 int k;
1435 bUse = TRUE;
1436 clear_ivec(k_zero);
1437 clear_ivec(k_plus);
1438 for (k = 1; k <= nral && bUse; k++)
1440 int k_gl;
1441 if (!bInterMolInteractions)
1443 /* Get the global index using the offset in the molecule */
1444 k_gl = i_gl + iatoms[k] - i_mol;
1446 else
1448 k_gl = iatoms[k];
1450 const auto *entry = dd->ga2la->find(k_gl);
1451 if (entry == nullptr || entry->cell >= zones->n)
1453 /* We do not have this atom of this interaction
1454 * locally, or it comes from more than one cell
1455 * away.
1457 bUse = FALSE;
1459 else
1461 int d;
1463 tiatoms[k] = entry->la;
1464 for (d = 0; d < DIM; d++)
1466 if (zones->shift[entry->cell][d] == 0)
1468 k_zero[d] = k;
1470 else
1472 k_plus[d] = k;
1477 bUse = (bUse &&
1478 (k_zero[XX] != 0) && (k_zero[YY] != 0) && (k_zero[ZZ] != 0));
1479 if (bRCheckMB)
1481 int d;
1483 for (d = 0; (d < DIM && bUse); d++)
1485 /* Check if the cg_cm distance falls within
1486 * the cut-off to avoid possible multiple
1487 * assignments of bonded interactions.
1489 if (rcheck[d] &&
1490 k_plus[d] &&
1491 dd_dist2(pbc_null, cg_cm, la2lc,
1492 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1494 bUse = FALSE;
1499 if (bUse)
1501 /* Add this interaction to the local topology */
1502 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1503 /* Sum so we can check in global_stat
1504 * if we have everything.
1506 if (bBCheck ||
1507 !(interaction_function[ftype].flags & IF_LIMZERO))
1509 (*nbonded_local)++;
1512 j += 1 + nral;
1517 /*! \brief This function looks up and assigns bonded interactions for zone iz.
1519 * With thread parallelizing each thread acts on a different atom range:
1520 * at_start to at_end.
1522 static int make_bondeds_zone(gmx_domdec_t *dd,
1523 const gmx_domdec_zones_t *zones,
1524 const std::vector<gmx_molblock_t> &molb,
1525 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1526 real rc2,
1527 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1528 const t_iparams *ip_in,
1529 t_idef *idef,
1530 int **vsite_pbc,
1531 int *vsite_pbc_nalloc,
1532 int izone,
1533 gmx::RangePartitioning::Block atomRange)
1535 int mb, mt, mol, i_mol;
1536 int *index, *rtil;
1537 gmx_bool bBCheck;
1538 gmx_reverse_top_t *rt;
1539 int nbonded_local;
1541 rt = dd->reverse_top;
1543 bBCheck = rt->bBCheck;
1545 nbonded_local = 0;
1547 for (int i : atomRange)
1549 /* Get the global atom number */
1550 const int i_gl = dd->globalAtomIndices[i];
1551 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1552 /* Check all intramolecular interactions assigned to this atom */
1553 index = rt->ril_mt[mt].index;
1554 rtil = rt->ril_mt[mt].il;
1556 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1557 rt->ril_mt[mt].numAtomsInMolecule,
1558 index, rtil, FALSE,
1559 index[i_mol], index[i_mol+1],
1560 dd, zones,
1561 &molb[mb],
1562 bRCheckMB, rcheck, bRCheck2B, rc2,
1563 la2lc,
1564 pbc_null,
1565 cg_cm,
1566 ip_in,
1567 idef, vsite_pbc, vsite_pbc_nalloc,
1568 izone,
1569 bBCheck,
1570 &nbonded_local);
1573 if (rt->bIntermolecularInteractions)
1575 /* Check all intermolecular interactions assigned to this atom */
1576 index = rt->ril_intermol.index;
1577 rtil = rt->ril_intermol.il;
1579 check_assign_interactions_atom(i, i_gl, mol, i_mol,
1580 rt->ril_mt[mt].numAtomsInMolecule,
1581 index, rtil, TRUE,
1582 index[i_gl], index[i_gl + 1],
1583 dd, zones,
1584 &molb[mb],
1585 bRCheckMB, rcheck, bRCheck2B, rc2,
1586 la2lc,
1587 pbc_null,
1588 cg_cm,
1589 ip_in,
1590 idef, vsite_pbc, vsite_pbc_nalloc,
1591 izone,
1592 bBCheck,
1593 &nbonded_local);
1597 return nbonded_local;
1600 /*! \brief Set the exclusion data for i-zone \p iz for the case of no exclusions */
1601 static void set_no_exclusions_zone(const gmx_domdec_t *dd,
1602 const gmx_domdec_zones_t *zones,
1603 int iz,
1604 t_blocka *lexcls)
1606 const auto zone = dd->atomGrouping().subRange(zones->cg_range[iz],
1607 zones->cg_range[iz + 1]);
1609 for (int a : zone)
1611 lexcls->index[a + 1] = lexcls->nra;
1615 /*! \brief Set the exclusion data for i-zone \p iz
1617 * This is a legacy version for the group scheme of the same routine below.
1618 * Here charge groups and distance checks to ensure unique exclusions
1619 * are supported.
1621 static int make_exclusions_zone_cg(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1622 const std::vector<gmx_moltype_t> &moltype,
1623 gmx_bool bRCheck, real rc2,
1624 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1625 const int *cginfo,
1626 t_blocka *lexcls,
1627 int iz,
1628 int cg_start, int cg_end)
1630 int n_excl_at_max;
1631 int mb, mt, mol;
1632 const t_blocka *excls;
1634 const gmx_ga2la_t &ga2la = *dd->ga2la;
1636 const auto jRange =
1637 dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1638 zones->izone[iz].jcg1);
1640 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1642 /* We set the end index, but note that we might not start at zero here */
1643 lexcls->nr = dd->atomGrouping().subRange(0, cg_end).size();
1645 int n = lexcls->nra;
1646 int count = 0;
1647 for (int cg = cg_start; cg < cg_end; cg++)
1649 if (n + (cg_end - cg_start)*n_excl_at_max > lexcls->nalloc_a)
1651 lexcls->nalloc_a = over_alloc_large(n + (cg_end - cg_start)*n_excl_at_max);
1652 srenew(lexcls->a, lexcls->nalloc_a);
1654 const auto atomGroup = dd->atomGrouping().block(cg);
1655 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1656 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1658 /* Copy the exclusions from the global top */
1659 for (int la : atomGroup)
1661 lexcls->index[la] = n;
1662 int a_gl = dd->globalAtomIndices[la];
1663 int a_mol;
1664 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1665 excls = &moltype[mt].excls;
1666 for (int j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1668 int aj_mol = excls->a[j];
1669 /* This computation of jla is only correct intra-cg */
1670 int jla = la + aj_mol - a_mol;
1671 if (atomGroup.inRange(jla))
1673 /* This is an intra-cg exclusion. We can skip
1674 * the global indexing and distance checking.
1676 /* Intra-cg exclusions are only required
1677 * for the home zone.
1679 if (iz == 0)
1681 lexcls->a[n++] = jla;
1682 /* Check to avoid double counts */
1683 if (jla > la)
1685 count++;
1689 else
1691 /* This is a inter-cg exclusion */
1692 /* Since exclusions are pair interactions,
1693 * just like non-bonded interactions,
1694 * they can be assigned properly up
1695 * to the DD cutoff (not cutoff_min as
1696 * for the other bonded interactions).
1698 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
1700 if (iz == 0 && jEntry->cell == 0)
1702 lexcls->a[n++] = jEntry->la;
1703 /* Check to avoid double counts */
1704 if (jEntry->la > la)
1706 count++;
1709 else if (jRange.inRange(jEntry->la) &&
1710 (!bRCheck ||
1711 dd_dist2(pbc_null, cg_cm, la2lc, la, jEntry->la) < rc2))
1713 /* jla > la, since jRange.begin() > la */
1714 lexcls->a[n++] = jEntry->la;
1715 count++;
1722 else
1724 /* There are no inter-cg excls and this cg is self-excluded.
1725 * These exclusions are only required for zone 0,
1726 * since other zones do not see themselves.
1728 if (iz == 0)
1730 for (int la : atomGroup)
1732 lexcls->index[la] = n;
1733 for (int j : atomGroup)
1735 lexcls->a[n++] = j;
1738 count += (atomGroup.size()*(atomGroup.size() - 1))/2;
1740 else
1742 /* We don't need exclusions for this cg */
1743 for (int la : atomGroup)
1745 lexcls->index[la] = n;
1751 lexcls->index[lexcls->nr] = n;
1752 lexcls->nra = n;
1754 return count;
1757 /*! \brief Set the exclusion data for i-zone \p iz */
1758 static void make_exclusions_zone(gmx_domdec_t *dd,
1759 gmx_domdec_zones_t *zones,
1760 const std::vector<gmx_moltype_t> &moltype,
1761 const int *cginfo,
1762 t_blocka *lexcls,
1763 int iz,
1764 int at_start, int at_end)
1766 int n_excl_at_max, n, at;
1768 const gmx_ga2la_t &ga2la = *dd->ga2la;
1770 const auto jRange =
1771 dd->atomGrouping().subRange(zones->izone[iz].jcg0,
1772 zones->izone[iz].jcg1);
1774 n_excl_at_max = dd->reverse_top->n_excl_at_max;
1776 /* We set the end index, but note that we might not start at zero here */
1777 lexcls->nr = at_end;
1779 n = lexcls->nra;
1780 for (at = at_start; at < at_end; at++)
1782 if (n + 1000 > lexcls->nalloc_a)
1784 lexcls->nalloc_a = over_alloc_large(n + 1000);
1785 srenew(lexcls->a, lexcls->nalloc_a);
1787 if (GET_CGINFO_EXCL_INTER(cginfo[at]))
1789 int a_gl, mb, mt, mol, a_mol, j;
1790 const t_blocka *excls;
1792 if (n + n_excl_at_max > lexcls->nalloc_a)
1794 lexcls->nalloc_a = over_alloc_large(n + n_excl_at_max);
1795 srenew(lexcls->a, lexcls->nalloc_a);
1798 /* Copy the exclusions from the global top */
1799 lexcls->index[at] = n;
1800 a_gl = dd->globalAtomIndices[at];
1801 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl,
1802 &mb, &mt, &mol, &a_mol);
1803 excls = &moltype[mt].excls;
1804 for (j = excls->index[a_mol]; j < excls->index[a_mol + 1]; j++)
1806 const int aj_mol = excls->a[j];
1808 if (const auto *jEntry = ga2la.find(a_gl + aj_mol - a_mol))
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 (jRange.inRange(jEntry->la))
1816 lexcls->a[n++] = jEntry->la;
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->atomGrouping().subRange(0, zones->izone[zones->nizone - 1].cg1).size();
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 const auto nonhomeIzonesAtomRange =
1862 dd->atomGrouping().subRange(zones->izone[0].cg1,
1863 zones->izone[zones->nizone - 1].cg1);
1865 if (dd->n_intercg_excl == 0)
1867 /* There are no exclusions involving non-home charge groups,
1868 * but we need to set the indices for neighborsearching.
1870 for (int la : nonhomeIzonesAtomRange)
1872 lexcls->index[la] = lexcls->nra;
1875 /* nr is only used to loop over the exclusions for Ewald and RF,
1876 * so we can set it to the number of home atoms for efficiency.
1878 lexcls->nr = nonhomeIzonesAtomRange.begin();
1880 else
1882 lexcls->nr = nonhomeIzonesAtomRange.end();
1886 /*! \brief Clear a t_idef struct */
1887 static void clear_idef(t_idef *idef)
1889 int ftype;
1891 /* Clear the counts */
1892 for (ftype = 0; ftype < F_NRE; ftype++)
1894 idef->il[ftype].nr = 0;
1898 /*! \brief Generate and store all required local bonded interactions in \p idef and local exclusions in \p lexcls */
1899 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1900 gmx_domdec_zones_t *zones,
1901 const gmx_mtop_t *mtop,
1902 const int *cginfo,
1903 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1904 real rc,
1905 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1906 t_idef *idef, gmx_vsite_t *vsite,
1907 t_blocka *lexcls, int *excl_count)
1909 int nzone_bondeds, nzone_excl;
1910 int izone, cg0, cg1;
1911 real rc2;
1912 int nbonded_local;
1913 int thread;
1914 gmx_reverse_top_t *rt;
1916 if (dd->reverse_top->bInterCGInteractions)
1918 nzone_bondeds = zones->n;
1920 else
1922 /* Only single charge group (or atom) molecules, so interactions don't
1923 * cross zone boundaries and we only need to assign in the home zone.
1925 nzone_bondeds = 1;
1928 if (dd->n_intercg_excl > 0)
1930 /* We only use exclusions from i-zones to i- and j-zones */
1931 nzone_excl = zones->nizone;
1933 else
1935 /* There are no inter-cg exclusions and only zone 0 sees itself */
1936 nzone_excl = 1;
1939 check_exclusions_alloc(dd, zones, lexcls);
1941 rt = dd->reverse_top;
1943 rc2 = rc*rc;
1945 /* Clear the counts */
1946 clear_idef(idef);
1947 nbonded_local = 0;
1949 lexcls->nr = 0;
1950 lexcls->nra = 0;
1951 *excl_count = 0;
1953 for (izone = 0; izone < nzone_bondeds; izone++)
1955 cg0 = zones->cg_range[izone];
1956 cg1 = zones->cg_range[izone + 1];
1958 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1959 for (thread = 0; thread < rt->nthread; thread++)
1963 int cg0t, cg1t;
1964 t_idef *idef_t;
1965 int **vsite_pbc;
1966 int *vsite_pbc_nalloc;
1967 t_blocka *excl_t;
1969 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
1970 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
1972 if (thread == 0)
1974 idef_t = idef;
1976 else
1978 idef_t = &rt->th_work[thread].idef;
1979 clear_idef(idef_t);
1982 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1984 if (thread == 0)
1986 vsite_pbc = vsite->vsite_pbc_loc;
1987 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
1989 else
1991 vsite_pbc = rt->th_work[thread].vsite_pbc;
1992 vsite_pbc_nalloc = rt->th_work[thread].vsite_pbc_nalloc;
1995 else
1997 vsite_pbc = nullptr;
1998 vsite_pbc_nalloc = nullptr;
2001 rt->th_work[thread].nbonded =
2002 make_bondeds_zone(dd, zones,
2003 mtop->molblock,
2004 bRCheckMB, rcheck, bRCheck2B, rc2,
2005 la2lc, pbc_null, cg_cm, idef->iparams,
2006 idef_t,
2007 vsite_pbc, vsite_pbc_nalloc,
2008 izone,
2009 dd->atomGrouping().subRange(cg0t, cg1t));
2011 if (izone < nzone_excl)
2013 if (thread == 0)
2015 excl_t = lexcls;
2017 else
2019 excl_t = &rt->th_work[thread].excl;
2020 excl_t->nr = 0;
2021 excl_t->nra = 0;
2024 if (dd->atomGrouping().allBlocksHaveSizeOne() &&
2025 !rt->bExclRequired)
2027 /* No charge groups and no distance check required */
2028 make_exclusions_zone(dd, zones,
2029 mtop->moltype, cginfo,
2030 excl_t,
2031 izone,
2032 cg0t, cg1t);
2034 else
2036 rt->th_work[thread].excl_count =
2037 make_exclusions_zone_cg(dd, zones,
2038 mtop->moltype, bRCheck2B, rc2,
2039 la2lc, pbc_null, cg_cm, cginfo,
2040 excl_t,
2041 izone,
2042 cg0t, cg1t);
2046 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2049 if (rt->nthread > 1)
2051 combine_idef(idef, rt->th_work, rt->nthread, vsite);
2054 for (thread = 0; thread < rt->nthread; thread++)
2056 nbonded_local += rt->th_work[thread].nbonded;
2059 if (izone < nzone_excl)
2061 if (rt->nthread > 1)
2063 combine_blocka(lexcls, rt->th_work, rt->nthread);
2066 for (thread = 0; thread < rt->nthread; thread++)
2068 *excl_count += rt->th_work[thread].excl_count;
2073 /* Some zones might not have exclusions, but some code still needs to
2074 * loop over the index, so we set the indices here.
2076 for (izone = nzone_excl; izone < zones->nizone; izone++)
2078 set_no_exclusions_zone(dd, zones, izone, lexcls);
2081 finish_local_exclusions(dd, zones, lexcls);
2082 if (debug)
2084 fprintf(debug, "We have %d exclusions, check count %d\n",
2085 lexcls->nra, *excl_count);
2088 return nbonded_local;
2091 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
2093 lcgs->nr = dd->globalAtomGroupIndices.size();
2094 lcgs->index = dd->atomGrouping_.rawIndex().data();
2097 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
2098 int npbcdim, matrix box,
2099 rvec cellsize_min, const ivec npulse,
2100 t_forcerec *fr,
2101 rvec *cgcm_or_x,
2102 gmx_vsite_t *vsite,
2103 const gmx_mtop_t *mtop, gmx_localtop_t *ltop)
2105 gmx_bool bRCheckMB, bRCheck2B;
2106 real rc = -1;
2107 ivec rcheck;
2108 int d, nexcl;
2109 t_pbc pbc, *pbc_null = nullptr;
2111 if (debug)
2113 fprintf(debug, "Making local topology\n");
2116 dd_make_local_cgs(dd, &ltop->cgs);
2118 bRCheckMB = FALSE;
2119 bRCheck2B = FALSE;
2121 if (dd->reverse_top->bInterCGInteractions)
2123 /* We need to check to which cell bondeds should be assigned */
2124 rc = dd_cutoff_twobody(dd);
2125 if (debug)
2127 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
2130 /* Should we check cg_cm distances when assigning bonded interactions? */
2131 for (d = 0; d < DIM; d++)
2133 rcheck[d] = FALSE;
2134 /* Only need to check for dimensions where the part of the box
2135 * that is not communicated is smaller than the cut-off.
2137 if (d < npbcdim && dd->nc[d] > 1 &&
2138 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
2140 if (dd->nc[d] == 2)
2142 rcheck[d] = TRUE;
2143 bRCheckMB = TRUE;
2145 /* Check for interactions between two atoms,
2146 * where we can allow interactions up to the cut-off,
2147 * instead of up to the smallest cell dimension.
2149 bRCheck2B = TRUE;
2151 if (debug)
2153 fprintf(debug,
2154 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %s\n",
2155 d, cellsize_min[d], d, rcheck[d], gmx::boolToString(bRCheck2B));
2158 if (bRCheckMB || bRCheck2B)
2160 makeLocalAtomGroupsFromAtoms(dd);
2161 if (fr->bMolPBC)
2163 pbc_null = set_pbc_dd(&pbc, fr->ePBC, dd->nc, TRUE, box);
2165 else
2167 pbc_null = nullptr;
2172 dd->nbonded_local =
2173 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
2174 bRCheckMB, rcheck, bRCheck2B, rc,
2175 dd->localAtomGroupFromAtom.data(),
2176 pbc_null, cgcm_or_x,
2177 &ltop->idef, vsite,
2178 &ltop->excls, &nexcl);
2180 /* The ilist is not sorted yet,
2181 * we can only do this when we have the charge arrays.
2183 ltop->idef.ilsort = ilsortUNKNOWN;
2185 if (dd->reverse_top->bExclRequired)
2187 dd->nbonded_local += nexcl;
2190 ltop->atomtypes = mtop->atomtypes;
2193 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
2194 gmx_localtop_t *ltop)
2196 if (dd->reverse_top->ilsort == ilsortNO_FE)
2198 ltop->idef.ilsort = ilsortNO_FE;
2200 else
2202 gmx_sort_ilist_fe(&ltop->idef, mdatoms->chargeA, mdatoms->chargeB);
2206 gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global)
2208 gmx_localtop_t *top;
2209 int i;
2211 snew(top, 1);
2213 top->idef.ntypes = top_global->ffparams.ntypes;
2214 top->idef.atnr = top_global->ffparams.atnr;
2215 top->idef.functype = top_global->ffparams.functype;
2216 top->idef.iparams = top_global->ffparams.iparams;
2217 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
2218 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
2220 for (i = 0; i < F_NRE; i++)
2222 top->idef.il[i].iatoms = nullptr;
2223 top->idef.il[i].nalloc = 0;
2225 top->idef.ilsort = ilsortUNKNOWN;
2227 return top;
2230 void dd_init_local_state(gmx_domdec_t *dd,
2231 const t_state *state_global, t_state *state_local)
2233 int buf[NITEM_DD_INIT_LOCAL_STATE];
2235 if (DDMASTER(dd))
2237 buf[0] = state_global->flags;
2238 buf[1] = state_global->ngtc;
2239 buf[2] = state_global->nnhpres;
2240 buf[3] = state_global->nhchainlength;
2241 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
2243 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
2245 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
2246 init_dfhist_state(state_local, buf[4]);
2247 state_local->flags = buf[0];
2250 /*! \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 */
2251 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
2253 int k;
2254 gmx_bool bFound;
2256 bFound = FALSE;
2257 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
2259 GMX_RELEASE_ASSERT(link->a, "Inconsistent NULL pointer while making charge-group links");
2260 if (link->a[k] == cg_gl_j)
2262 bFound = TRUE;
2265 if (!bFound)
2267 GMX_RELEASE_ASSERT(link->a || link->index[cg_gl+1]+1 > link->nalloc_a,
2268 "Inconsistent allocation of link");
2269 /* Add this charge group link */
2270 if (link->index[cg_gl+1]+1 > link->nalloc_a)
2272 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
2273 srenew(link->a, link->nalloc_a);
2275 link->a[link->index[cg_gl+1]] = cg_gl_j;
2276 link->index[cg_gl+1]++;
2280 /*! \brief Return a vector of the charge group index for all atoms */
2281 static std::vector<int> make_at2cg(const t_block &cgs)
2283 std::vector<int> at2cg(cgs.index[cgs.nr]);
2284 for (int cg = 0; cg < cgs.nr; cg++)
2286 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2288 at2cg[a] = cg;
2292 return at2cg;
2295 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
2296 cginfo_mb_t *cginfo_mb)
2298 gmx_bool bExclRequired;
2299 reverse_ilist_t ril, ril_intermol;
2300 t_blocka *link;
2301 cginfo_mb_t *cgi_mb;
2303 /* For each charge group make a list of other charge groups
2304 * in the system that a linked to it via bonded interactions
2305 * which are also stored in reverse_top.
2308 bExclRequired = dd->reverse_top->bExclRequired;
2310 if (mtop->bIntermolecularInteractions)
2312 if (ncg_mtop(mtop) < mtop->natoms)
2314 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.");
2317 t_atoms atoms;
2319 atoms.nr = mtop->natoms;
2320 atoms.atom = nullptr;
2322 make_reverse_ilist(mtop->intermolecular_ilist, &atoms,
2323 nullptr, FALSE, FALSE, FALSE, TRUE, &ril_intermol);
2326 snew(link, 1);
2327 snew(link->index, ncg_mtop(mtop)+1);
2328 link->nalloc_a = 0;
2329 link->a = nullptr;
2331 link->index[0] = 0;
2332 int cg_offset = 0;
2333 int ncgi = 0;
2334 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
2336 const gmx_molblock_t &molb = mtop->molblock[mb];
2337 if (molb.nmol == 0)
2339 continue;
2341 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2342 const t_block &cgs = molt.cgs;
2343 const t_blocka &excls = molt.excls;
2344 std::vector<int> a2c = make_at2cg(cgs);
2345 /* Make a reverse ilist in which the interactions are linked
2346 * to all atoms, not only the first atom as in gmx_reverse_top.
2347 * The constraints are discarded here.
2349 make_reverse_ilist(molt.ilist, &molt.atoms,
2350 nullptr, FALSE, FALSE, FALSE, TRUE, &ril);
2352 cgi_mb = &cginfo_mb[mb];
2354 int mol;
2355 for (mol = 0; mol < (mtop->bIntermolecularInteractions ? molb.nmol : 1); mol++)
2357 for (int cg = 0; cg < cgs.nr; cg++)
2359 int cg_gl = cg_offset + cg;
2360 link->index[cg_gl+1] = link->index[cg_gl];
2361 for (int a = cgs.index[cg]; a < cgs.index[cg + 1]; a++)
2363 int i = ril.index[a];
2364 while (i < ril.index[a+1])
2366 int ftype = ril.il[i++];
2367 int nral = NRAL(ftype);
2368 /* Skip the ifunc index */
2369 i++;
2370 for (int j = 0; j < nral; j++)
2372 int aj = ril.il[i + j];
2373 if (a2c[aj] != cg)
2375 check_link(link, cg_gl, cg_offset+a2c[aj]);
2378 i += nral_rt(ftype);
2380 if (bExclRequired)
2382 /* Exclusions always go both ways */
2383 for (int j = excls.index[a]; j < excls.index[a + 1]; j++)
2385 int aj = excls.a[j];
2386 if (a2c[aj] != cg)
2388 check_link(link, cg_gl, cg_offset+a2c[aj]);
2393 if (mtop->bIntermolecularInteractions)
2395 int i = ril_intermol.index[a];
2396 while (i < ril_intermol.index[a+1])
2398 int ftype = ril_intermol.il[i++];
2399 int nral = NRAL(ftype);
2400 /* Skip the ifunc index */
2401 i++;
2402 for (int j = 0; j < nral; j++)
2404 /* Here we assume we have no charge groups;
2405 * this has been checked above.
2407 int aj = ril_intermol.il[i + j];
2408 check_link(link, cg_gl, aj);
2410 i += nral_rt(ftype);
2414 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2416 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2417 ncgi++;
2421 cg_offset += cgs.nr;
2423 int nlink_mol = link->index[cg_offset] - link->index[cg_offset - cgs.nr];
2425 destroy_reverse_ilist(&ril);
2427 if (debug)
2429 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt.name, cgs.nr, nlink_mol);
2432 if (molb.nmol > mol)
2434 /* Copy the data for the rest of the molecules in this block */
2435 link->nalloc_a += (molb.nmol - mol)*nlink_mol;
2436 srenew(link->a, link->nalloc_a);
2437 for (; mol < molb.nmol; mol++)
2439 for (int cg = 0; cg < cgs.nr; cg++)
2441 int cg_gl = cg_offset + cg;
2442 link->index[cg_gl + 1] =
2443 link->index[cg_gl + 1 - cgs.nr] + nlink_mol;
2444 for (int j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2446 link->a[j] = link->a[j - nlink_mol] + cgs.nr;
2448 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2449 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2451 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2452 ncgi++;
2455 cg_offset += cgs.nr;
2460 if (mtop->bIntermolecularInteractions)
2462 destroy_reverse_ilist(&ril_intermol);
2465 if (debug)
2467 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2470 return link;
2473 typedef struct {
2474 real r2;
2475 int ftype;
2476 int a1;
2477 int a2;
2478 } bonded_distance_t;
2480 /*! \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 */
2481 static void update_max_bonded_distance(real r2, int ftype, int a1, int a2,
2482 bonded_distance_t *bd)
2484 if (r2 > bd->r2)
2486 bd->r2 = r2;
2487 bd->ftype = ftype;
2488 bd->a1 = a1;
2489 bd->a2 = a2;
2493 /*! \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 */
2494 static void bonded_cg_distance_mol(const gmx_moltype_t *molt,
2495 const std::vector<int> &at2cg,
2496 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2497 bonded_distance_t *bd_2b,
2498 bonded_distance_t *bd_mb)
2500 for (int ftype = 0; ftype < F_NRE; ftype++)
2502 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2504 const t_ilist *il = &molt->ilist[ftype];
2505 int nral = NRAL(ftype);
2506 if (nral > 1)
2508 for (int i = 0; i < il->nr; i += 1+nral)
2510 for (int ai = 0; ai < nral; ai++)
2512 int cgi = at2cg[il->iatoms[i+1+ai]];
2513 for (int aj = ai + 1; aj < nral; aj++)
2515 int cgj = at2cg[il->iatoms[i+1+aj]];
2516 if (cgi != cgj)
2518 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2520 update_max_bonded_distance(rij2, ftype,
2521 il->iatoms[i+1+ai],
2522 il->iatoms[i+1+aj],
2523 (nral == 2) ? bd_2b : bd_mb);
2531 if (bExcl)
2533 const t_blocka *excls = &molt->excls;
2534 for (int ai = 0; ai < excls->nr; ai++)
2536 int cgi = at2cg[ai];
2537 for (int j = excls->index[ai]; j < excls->index[ai+1]; j++)
2539 int cgj = at2cg[excls->a[j]];
2540 if (cgi != cgj)
2542 real rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2544 /* There is no function type for exclusions, use -1 */
2545 update_max_bonded_distance(rij2, -1, ai, excls->a[j], bd_2b);
2552 /*! \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 */
2553 static void bonded_distance_intermol(const t_ilist *ilists_intermol,
2554 gmx_bool bBCheck,
2555 const rvec *x, int ePBC, const matrix box,
2556 bonded_distance_t *bd_2b,
2557 bonded_distance_t *bd_mb)
2559 t_pbc pbc;
2561 set_pbc(&pbc, ePBC, box);
2563 for (int ftype = 0; ftype < F_NRE; ftype++)
2565 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2567 const t_ilist *il = &ilists_intermol[ftype];
2568 int nral = NRAL(ftype);
2570 /* No nral>1 check here, since intermol interactions always
2571 * have nral>=2 (and the code is also correct for nral=1).
2573 for (int i = 0; i < il->nr; i += 1+nral)
2575 for (int ai = 0; ai < nral; ai++)
2577 int atom_i = il->iatoms[i + 1 + ai];
2579 for (int aj = ai + 1; aj < nral; aj++)
2581 rvec dx;
2582 real rij2;
2584 int atom_j = il->iatoms[i + 1 + aj];
2586 pbc_dx(&pbc, x[atom_i], x[atom_j], dx);
2588 rij2 = norm2(dx);
2590 update_max_bonded_distance(rij2, ftype,
2591 atom_i, atom_j,
2592 (nral == 2) ? bd_2b : bd_mb);
2600 //! Returns whether \p molt has at least one virtual site
2601 static bool moltypeHasVsite(const gmx_moltype_t &molt)
2603 bool hasVsite = false;
2604 for (int i = 0; i < F_NRE; i++)
2606 if ((interaction_function[i].flags & IF_VSITE) &&
2607 molt.ilist[i].nr > 0)
2609 hasVsite = true;
2613 return hasVsite;
2616 //! Compute charge group centers of mass for molecule \p molt
2617 static void get_cgcm_mol(const gmx_moltype_t *molt,
2618 const gmx_ffparams_t *ffparams,
2619 int ePBC, t_graph *graph, const matrix box,
2620 const rvec *x, rvec *xs, rvec *cg_cm)
2622 int n, i;
2624 if (ePBC != epbcNONE)
2626 mk_mshift(nullptr, graph, ePBC, box, x);
2628 shift_x(graph, box, x, xs);
2629 /* By doing an extra mk_mshift the molecules that are broken
2630 * because they were e.g. imported from another software
2631 * will be made whole again. Such are the healing powers
2632 * of GROMACS.
2634 mk_mshift(nullptr, graph, ePBC, box, xs);
2636 else
2638 /* We copy the coordinates so the original coordinates remain
2639 * unchanged, just to be 100% sure that we do not affect
2640 * binary reproducibility of simulations.
2642 n = molt->cgs.index[molt->cgs.nr];
2643 for (i = 0; i < n; i++)
2645 copy_rvec(x[i], xs[i]);
2649 if (moltypeHasVsite(*molt))
2651 construct_vsites(nullptr, xs, 0.0, nullptr,
2652 ffparams->iparams, molt->ilist,
2653 epbcNONE, TRUE, nullptr, nullptr);
2656 calc_cgcm(nullptr, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2659 void dd_bonded_cg_distance(FILE *fplog,
2660 const gmx_mtop_t *mtop,
2661 const t_inputrec *ir,
2662 const rvec *x, const matrix box,
2663 gmx_bool bBCheck,
2664 real *r_2b, real *r_mb)
2666 gmx_bool bExclRequired;
2667 int at_offset;
2668 t_graph graph;
2669 rvec *xs, *cg_cm;
2670 bonded_distance_t bd_2b = { 0, -1, -1, -1 };
2671 bonded_distance_t bd_mb = { 0, -1, -1, -1 };
2673 bExclRequired = inputrecExclForces(ir);
2675 *r_2b = 0;
2676 *r_mb = 0;
2677 at_offset = 0;
2678 for (const gmx_molblock_t &molb : mtop->molblock)
2680 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2681 if (molt.cgs.nr == 1 || molb.nmol == 0)
2683 at_offset += molb.nmol*molt.atoms.nr;
2685 else
2687 if (ir->ePBC != epbcNONE)
2689 mk_graph_ilist(nullptr, molt.ilist, 0, molt.atoms.nr, FALSE, FALSE,
2690 &graph);
2693 std::vector<int> at2cg = make_at2cg(molt.cgs);
2694 snew(xs, molt.atoms.nr);
2695 snew(cg_cm, molt.cgs.nr);
2696 for (int mol = 0; mol < molb.nmol; mol++)
2698 get_cgcm_mol(&molt, &mtop->ffparams, ir->ePBC, &graph, box,
2699 x+at_offset, xs, cg_cm);
2701 bonded_distance_t bd_mol_2b = { 0, -1, -1, -1 };
2702 bonded_distance_t bd_mol_mb = { 0, -1, -1, -1 };
2704 bonded_cg_distance_mol(&molt, at2cg, bBCheck, bExclRequired, cg_cm,
2705 &bd_mol_2b, &bd_mol_mb);
2707 /* Process the mol data adding the atom index offset */
2708 update_max_bonded_distance(bd_mol_2b.r2, bd_mol_2b.ftype,
2709 at_offset + bd_mol_2b.a1,
2710 at_offset + bd_mol_2b.a2,
2711 &bd_2b);
2712 update_max_bonded_distance(bd_mol_mb.r2, bd_mol_mb.ftype,
2713 at_offset + bd_mol_mb.a1,
2714 at_offset + bd_mol_mb.a2,
2715 &bd_mb);
2717 at_offset += molt.atoms.nr;
2719 sfree(cg_cm);
2720 sfree(xs);
2721 if (ir->ePBC != epbcNONE)
2723 done_graph(&graph);
2728 if (mtop->bIntermolecularInteractions)
2730 if (ncg_mtop(mtop) < mtop->natoms)
2732 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.");
2735 bonded_distance_intermol(mtop->intermolecular_ilist,
2736 bBCheck,
2737 x, ir->ePBC, box,
2738 &bd_2b, &bd_mb);
2741 *r_2b = sqrt(bd_2b.r2);
2742 *r_mb = sqrt(bd_mb.r2);
2744 if (fplog && (*r_2b > 0 || *r_mb > 0))
2746 fprintf(fplog,
2747 "Initial maximum inter charge-group distances:\n");
2748 if (*r_2b > 0)
2750 fprintf(fplog,
2751 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2752 *r_2b, (bd_2b.ftype >= 0) ? interaction_function[bd_2b.ftype].longname : "Exclusion",
2753 bd_2b.a1 + 1, bd_2b.a2 + 1);
2755 if (*r_mb > 0)
2757 fprintf(fplog,
2758 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2759 *r_mb, interaction_function[bd_mb.ftype].longname,
2760 bd_mb.a1 + 1, bd_mb.a2 + 1);