2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
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.
41 #include "types/commrec.h"
43 #include "domdec_network.h"
46 #include "chargegroup.h"
48 #include "gmx_ga2la.h"
50 #include "gmx_omp_nthreads.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/mshift.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/topsort.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/smalloc.h"
61 /* for dd_init_local_state */
62 #define NITEM_DD_INIT_LOCAL_STATE 5
65 int *index; /* Index for each atom into il */
66 int *il; /* ftype|type|a0|...|an|ftype|... */
67 } gmx_reverse_ilist_t;
76 typedef struct gmx_reverse_top {
77 gmx_bool bExclRequired; /* Do we require all exclusions to be assigned? */
78 gmx_bool bConstr; /* Are there constraints in this revserse top? */
79 gmx_bool bSettle; /* Are there settles in this revserse top? */
80 gmx_bool bBCheck; /* All bonded interactions have to be assigned? */
81 gmx_bool bMultiCGmols; /* Are the multi charge-group molecules? */
82 gmx_reverse_ilist_t *ril_mt; /* Reverse ilist for all moltypes */
84 int ilsort; /* The sorting state of bondeds for free energy */
85 gmx_molblock_ind_t *mbi;
88 /* Work data structures for multi-threading */
92 int **vsite_pbc_nalloc;
94 t_blocka *excl_thread;
95 int *excl_count_thread;
97 /* Pointers only used for an error message */
98 gmx_mtop_t *err_top_global;
99 gmx_localtop_t *err_top_local;
102 static int nral_rt(int ftype)
104 /* Returns the number of atom entries for il in gmx_reverse_top_t */
108 if (interaction_function[ftype].flags & IF_VSITE)
110 /* With vsites the reverse topology contains
111 * two extra entries for PBC.
119 /* This function tells which interactions need to be assigned exactly once */
120 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
121 gmx_bool bConstr, gmx_bool bSettle)
123 return (((interaction_function[ftype].flags & IF_BOND) &&
124 !(interaction_function[ftype].flags & IF_VSITE) &&
125 (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))) ||
126 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
127 (bSettle && ftype == F_SETTLE));
130 static void print_error_header(FILE *fplog, char *moltypename, int nprint)
132 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
133 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
135 "the first %d missing interactions, except for exclusions:\n",
138 "the first %d missing interactions, except for exclusions:\n",
142 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
143 gmx_reverse_top_t *rt,
145 gmx_reverse_ilist_t *ril,
146 int a_start, int a_end,
147 int nat_mol, int nmol,
150 int nril_mol, *assigned, *gatindex;
151 int ftype, ftype_j, nral, i, j_mol, j, k, a0, a0_mol, mol, a, a_gl;
157 nril_mol = ril->index[nat_mol];
158 snew(assigned, nmol*nril_mol);
160 gatindex = cr->dd->gatindex;
161 for (ftype = 0; ftype < F_NRE; ftype++)
163 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
166 il = &idef->il[ftype];
168 for (i = 0; i < il->nr; i += 1+nral)
170 a0 = gatindex[ia[1]];
171 /* Check if this interaction is in
172 * the currently checked molblock.
174 if (a0 >= a_start && a0 < a_end)
176 mol = (a0 - a_start)/nat_mol;
177 a0_mol = (a0 - a_start) - mol*nat_mol;
178 j_mol = ril->index[a0_mol];
180 while (j_mol < ril->index[a0_mol+1] && !bFound)
182 j = mol*nril_mol + j_mol;
183 ftype_j = ril->il[j_mol];
184 /* Here we need to check if this interaction has
185 * not already been assigned, since we could have
186 * multiply defined interactions.
188 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
191 /* Check the atoms */
193 for (a = 0; a < nral; a++)
195 if (gatindex[ia[1+a]] !=
196 a_start + mol*nat_mol + ril->il[j_mol+2+a])
206 j_mol += 2 + nral_rt(ftype_j);
210 gmx_incons("Some interactions seem to be assigned multiple times");
218 gmx_sumi(nmol*nril_mol, assigned, cr);
222 for (mol = 0; mol < nmol; mol++)
225 while (j_mol < nril_mol)
227 ftype = ril->il[j_mol];
229 j = mol*nril_mol + j_mol;
230 if (assigned[j] == 0 &&
231 !(interaction_function[ftype].flags & IF_VSITE))
233 if (DDMASTER(cr->dd))
237 print_error_header(fplog, moltypename, nprint);
239 fprintf(fplog, "%20s atoms",
240 interaction_function[ftype].longname);
241 fprintf(stderr, "%20s atoms",
242 interaction_function[ftype].longname);
243 for (a = 0; a < nral; a++)
245 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
246 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
251 fprintf(stderr, " ");
254 fprintf(fplog, " global");
255 fprintf(stderr, " global");
256 for (a = 0; a < nral; a++)
258 fprintf(fplog, "%6d",
259 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
260 fprintf(stderr, "%6d",
261 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
263 fprintf(fplog, "\n");
264 fprintf(stderr, "\n");
272 j_mol += 2 + nral_rt(ftype);
279 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
280 gmx_mtop_t *mtop, t_idef *idef)
282 int mb, a_start, a_end;
283 gmx_molblock_t *molb;
284 gmx_reverse_top_t *rt;
286 rt = cr->dd->reverse_top;
288 /* Print the atoms in the missing interactions per molblock */
290 for (mb = 0; mb < mtop->nmolblock; mb++)
292 molb = &mtop->molblock[mb];
294 a_end = a_start + molb->nmol*molb->natoms_mol;
296 print_missing_interactions_mb(fplog, cr, rt,
297 *(mtop->moltype[molb->type].name),
298 &rt->ril_mt[molb->type],
299 a_start, a_end, molb->natoms_mol,
305 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr, int local_count, gmx_mtop_t *top_global, t_state *state_local)
307 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
311 gmx_mtop_t *err_top_global;
312 gmx_localtop_t *err_top_local;
316 err_top_global = dd->reverse_top->err_top_global;
317 err_top_local = dd->reverse_top->err_top_local;
321 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
325 ndiff_tot = local_count - dd->nbonded_global;
327 for (ftype = 0; ftype < F_NRE; ftype++)
330 cl[ftype] = err_top_local->idef.il[ftype].nr/(1+nral);
333 gmx_sumi(F_NRE, cl, cr);
337 fprintf(fplog, "\nA list of missing interactions:\n");
338 fprintf(stderr, "\nA list of missing interactions:\n");
339 rest_global = dd->nbonded_global;
340 rest_local = local_count;
341 for (ftype = 0; ftype < F_NRE; ftype++)
343 /* In the reverse and local top all constraints are merged
344 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
345 * and add these constraints when doing F_CONSTR.
347 if (((interaction_function[ftype].flags & IF_BOND) &&
348 (dd->reverse_top->bBCheck
349 || !(interaction_function[ftype].flags & IF_LIMZERO)))
350 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
351 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
354 n = gmx_mtop_ftype_count(err_top_global, ftype);
355 if (ftype == F_CONSTR)
357 n += gmx_mtop_ftype_count(err_top_global, F_CONSTRNC);
359 ndiff = cl[ftype] - n;
362 sprintf(buf, "%20s of %6d missing %6d",
363 interaction_function[ftype].longname, n, -ndiff);
364 fprintf(fplog, "%s\n", buf);
365 fprintf(stderr, "%s\n", buf);
368 rest_local -= cl[ftype];
372 ndiff = rest_local - rest_global;
375 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
376 rest_global, -ndiff);
377 fprintf(fplog, "%s\n", buf);
378 fprintf(stderr, "%s\n", buf);
382 print_missing_interactions_atoms(fplog, cr, err_top_global,
383 &err_top_local->idef);
384 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
385 -1, state_local->x, state_local->box);
390 gmx_incons("One or more interactions were multiple assigned in the domain decompostion");
394 gmx_fatal(FARGS, "%d of the %d bonded interactions could not be calculated because some atoms involved moved further apart than the multi-body cut-off distance (%g nm) or the two-body cut-off distance (%g nm), see option -rdd, for pairs and tabulated bonds also see option -ddcheck", -ndiff_tot, cr->dd->nbonded_global, dd_cutoff_mbody(cr->dd), dd_cutoff_twobody(cr->dd));
399 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
400 int *mb, int *mt, int *mol, int *i_mol)
405 gmx_molblock_ind_t *mbi = rt->mbi;
407 int end = rt->nmolblock; /* exclusive */
410 /* binary search for molblock_ind */
414 if (i_gl >= mbi[mid].a_end)
418 else if (i_gl < mbi[mid].a_start)
432 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
433 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
436 static int count_excls(t_block *cgs, t_blocka *excls, int *n_intercg_excl)
438 int n, n_inter, cg, at0, at1, at, excl, atj;
442 for (cg = 0; cg < cgs->nr; cg++)
444 at0 = cgs->index[cg];
445 at1 = cgs->index[cg+1];
446 for (at = at0; at < at1; at++)
448 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
450 atj = excls->a[excl];
454 if (atj < at0 || atj >= at1)
466 static int low_make_reverse_ilist(t_ilist *il_mt, t_atom *atom,
469 gmx_bool bConstr, gmx_bool bSettle,
471 int *r_index, int *r_il,
472 gmx_bool bLinkToAllAtoms,
475 int ftype, nral, i, j, nlink, link;
483 for (ftype = 0; ftype < F_NRE; ftype++)
485 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
486 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
487 (bSettle && ftype == F_SETTLE))
489 bVSite = (interaction_function[ftype].flags & IF_VSITE);
493 for (i = 0; i < il->nr; i += 1+nral)
500 /* We don't need the virtual sites for the cg-links */
510 /* Couple to the first atom in the interaction */
513 for (link = 0; link < nlink; link++)
518 r_il[r_index[a]+count[a]] =
519 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
520 r_il[r_index[a]+count[a]+1] = ia[0];
521 for (j = 1; j < 1+nral; j++)
523 /* Store the molecular atom number */
524 r_il[r_index[a]+count[a]+1+j] = ia[j];
527 if (interaction_function[ftype].flags & IF_VSITE)
531 /* Add an entry to iatoms for storing
532 * which of the constructing atoms are
535 r_il[r_index[a]+count[a]+2+nral] = 0;
536 for (j = 2; j < 1+nral; j++)
538 if (atom[ia[j]].ptype == eptVSite)
540 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
543 /* Store vsite pbc atom in a second extra entry */
544 r_il[r_index[a]+count[a]+2+nral+1] =
545 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
550 /* We do not count vsites since they are always
551 * uniquely assigned and can be assigned
552 * to multiple nodes with recursive vsites.
555 !(interaction_function[ftype].flags & IF_LIMZERO))
560 count[a] += 2 + nral_rt(ftype);
569 static int make_reverse_ilist(gmx_moltype_t *molt,
571 gmx_bool bConstr, gmx_bool bSettle,
573 gmx_bool bLinkToAllAtoms,
574 gmx_reverse_ilist_t *ril_mt)
576 int nat_mt, *count, i, nint_mt;
578 /* Count the interactions */
579 nat_mt = molt->atoms.nr;
581 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
583 bConstr, bSettle, bBCheck, NULL, NULL,
584 bLinkToAllAtoms, FALSE);
586 snew(ril_mt->index, nat_mt+1);
587 ril_mt->index[0] = 0;
588 for (i = 0; i < nat_mt; i++)
590 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
593 snew(ril_mt->il, ril_mt->index[nat_mt]);
595 /* Store the interactions */
597 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
599 bConstr, bSettle, bBCheck,
600 ril_mt->index, ril_mt->il,
601 bLinkToAllAtoms, TRUE);
608 static void destroy_reverse_ilist(gmx_reverse_ilist_t *ril)
614 static gmx_reverse_top_t *make_reverse_top(gmx_mtop_t *mtop, gmx_bool bFE,
615 int ***vsite_pbc_molt,
616 gmx_bool bConstr, gmx_bool bSettle,
617 gmx_bool bBCheck, int *nint)
620 gmx_reverse_top_t *rt;
627 /* Should we include constraints (for SHAKE) in rt? */
628 rt->bConstr = bConstr;
629 rt->bSettle = bSettle;
630 rt->bBCheck = bBCheck;
632 rt->bMultiCGmols = FALSE;
633 snew(nint_mt, mtop->nmoltype);
634 snew(rt->ril_mt, mtop->nmoltype);
635 rt->ril_mt_tot_size = 0;
636 for (mt = 0; mt < mtop->nmoltype; mt++)
638 molt = &mtop->moltype[mt];
639 if (molt->cgs.nr > 1)
641 rt->bMultiCGmols = TRUE;
644 /* Make the atom to interaction list for this molecule type */
646 make_reverse_ilist(molt, vsite_pbc_molt ? vsite_pbc_molt[mt] : NULL,
647 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
650 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt->atoms.nr];
654 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
658 for (mb = 0; mb < mtop->nmolblock; mb++)
660 *nint += mtop->molblock[mb].nmol*nint_mt[mtop->molblock[mb].type];
664 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
666 rt->ilsort = ilsortFE_UNSORTED;
670 rt->ilsort = ilsortNO_FE;
673 /* Make a molblock index for fast searching */
674 snew(rt->mbi, mtop->nmolblock);
675 rt->nmolblock = mtop->nmolblock;
677 for (mb = 0; mb < mtop->nmolblock; mb++)
679 rt->mbi[mb].a_start = i;
680 i += mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
681 rt->mbi[mb].a_end = i;
682 rt->mbi[mb].natoms_mol = mtop->molblock[mb].natoms_mol;
683 rt->mbi[mb].type = mtop->molblock[mb].type;
686 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
687 snew(rt->idef_thread, rt->nthread);
688 if (vsite_pbc_molt != NULL)
690 snew(rt->vsite_pbc, rt->nthread);
691 snew(rt->vsite_pbc_nalloc, rt->nthread);
692 for (thread = 0; thread < rt->nthread; thread++)
694 snew(rt->vsite_pbc[thread], F_VSITEN-F_VSITE2+1);
695 snew(rt->vsite_pbc_nalloc[thread], F_VSITEN-F_VSITE2+1);
698 snew(rt->nbonded_thread, rt->nthread);
699 snew(rt->excl_thread, rt->nthread);
700 snew(rt->excl_count_thread, rt->nthread);
705 void dd_make_reverse_top(FILE *fplog,
706 gmx_domdec_t *dd, gmx_mtop_t *mtop,
708 t_inputrec *ir, gmx_bool bBCheck)
710 int mb, n_recursive_vsite, nexcl, nexcl_icg, a;
711 gmx_molblock_t *molb;
716 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
719 /* If normal and/or settle constraints act only within charge groups,
720 * we can store them in the reverse top and simply assign them to domains.
721 * Otherwise we need to assign them to multiple domains and set up
722 * the parallel version constraint algoirthm(s).
725 dd->reverse_top = make_reverse_top(mtop, ir->efep != efepNO,
726 vsite ? vsite->vsite_pbc_molt : NULL,
727 !dd->bInterCGcons, !dd->bInterCGsettles,
728 bBCheck, &dd->nbonded_global);
730 if (dd->reverse_top->ril_mt_tot_size >= 200000 &&
732 mtop->nmolblock == 1 && mtop->molblock[0].nmol == 1)
734 /* mtop comes from a pre Gromacs 4 tpr file */
735 const char *note = "NOTE: The tpr file used for this simulation is in an old format, for less memory usage and possibly more performance create a new tpr file with an up to date version of grompp";
738 fprintf(fplog, "\n%s\n\n", note);
742 fprintf(stderr, "\n%s\n\n", note);
746 dd->reverse_top->bExclRequired = IR_EXCL_FORCES(*ir);
749 dd->n_intercg_excl = 0;
750 for (mb = 0; mb < mtop->nmolblock; mb++)
752 molb = &mtop->molblock[mb];
753 molt = &mtop->moltype[molb->type];
754 nexcl += molb->nmol*count_excls(&molt->cgs, &molt->excls, &nexcl_icg);
755 dd->n_intercg_excl += molb->nmol*nexcl_icg;
757 if (dd->reverse_top->bExclRequired)
759 dd->nbonded_global += nexcl;
760 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
762 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
763 "will use an extra communication step for exclusion forces for %s\n",
764 dd->n_intercg_excl, eel_names[ir->coulombtype]);
768 if (vsite && vsite->n_intercg_vsite > 0)
772 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
773 "will an extra communication step for selected coordinates and forces\n",
774 vsite->n_intercg_vsite);
776 init_domdec_vsites(dd, vsite->n_intercg_vsite);
779 if (dd->bInterCGcons || dd->bInterCGsettles)
781 init_domdec_constraints(dd, mtop);
785 fprintf(fplog, "\n");
789 static gmx_inline void add_ifunc(int nral, t_iatom *tiatoms, t_ilist *il)
794 if (il->nr+1+nral > il->nalloc)
796 il->nalloc = over_alloc_large(il->nr+1+nral);
797 srenew(il->iatoms, il->nalloc);
799 liatoms = il->iatoms + il->nr;
800 for (k = 0; k <= nral; k++)
802 liatoms[k] = tiatoms[k];
807 static void add_posres(int mol, int a_mol, const gmx_molblock_t *molb,
808 t_iatom *iatoms, const t_iparams *ip_in,
814 /* This position restraint has not been added yet,
815 * so it's index is the current number of position restraints.
817 n = idef->il[F_POSRES].nr/2;
818 if (n+1 > idef->iparams_posres_nalloc)
820 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
821 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
823 ip = &idef->iparams_posres[n];
824 /* Copy the force constants */
825 *ip = ip_in[iatoms[0]];
827 /* Get the position restraint coordinates from the molblock */
828 a_molb = mol*molb->natoms_mol + a_mol;
829 if (a_molb >= molb->nposres_xA)
831 gmx_incons("Not enough position restraint coordinates");
833 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
834 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
835 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
836 if (molb->nposres_xB > 0)
838 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
839 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
840 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
844 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
845 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
846 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
848 /* Set the parameter index for idef->iparams_posre */
852 static void add_fbposres(int mol, int a_mol, const gmx_molblock_t *molb,
853 t_iatom *iatoms, const t_iparams *ip_in,
859 /* This flat-bottom position restraint has not been added yet,
860 * so it's index is the current number of position restraints.
862 n = idef->il[F_FBPOSRES].nr/2;
863 if (n+1 > idef->iparams_fbposres_nalloc)
865 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
866 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
868 ip = &idef->iparams_fbposres[n];
869 /* Copy the force constants */
870 *ip = ip_in[iatoms[0]];
872 /* Get the position restriant coordinats from the molblock */
873 a_molb = mol*molb->natoms_mol + a_mol;
874 if (a_molb >= molb->nposres_xA)
876 gmx_incons("Not enough position restraint coordinates");
878 /* Take reference positions from A position of normal posres */
879 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
880 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
881 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
883 /* Note: no B-type for flat-bottom posres */
885 /* Set the parameter index for idef->iparams_posre */
889 static void add_vsite(gmx_ga2la_t ga2la, int *index, int *rtil,
891 gmx_bool bHomeA, int a, int a_gl, int a_mol,
893 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
895 int k, ak_gl, vsi, pbc_a_mol;
896 t_iatom tiatoms[1+MAXATOMLIST], *iatoms_r;
897 int j, ftype_r, nral_r;
900 tiatoms[0] = iatoms[0];
904 /* We know the local index of the first atom */
909 /* Convert later in make_local_vsites */
910 tiatoms[1] = -a_gl - 1;
913 for (k = 2; k < 1+nral; k++)
915 ak_gl = a_gl + iatoms[k] - a_mol;
916 if (!ga2la_get_home(ga2la, ak_gl, &tiatoms[k]))
918 /* Copy the global index, convert later in make_local_vsites */
919 tiatoms[k] = -(ak_gl + 1);
923 /* Add this interaction to the local topology */
924 add_ifunc(nral, tiatoms, &idef->il[ftype]);
927 vsi = idef->il[ftype].nr/(1+nral) - 1;
928 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
930 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
931 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
935 pbc_a_mol = iatoms[1+nral+1];
938 /* The pbc flag is one of the following two options:
939 * -2: vsite and all constructing atoms are within the same cg, no pbc
940 * -1: vsite and its first constructing atom are in the same cg, do pbc
942 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
946 /* Set the pbc atom for this vsite so we can make its pbc
947 * identical to the rest of the atoms in its charge group.
948 * Since the order of the atoms does not change within a charge
949 * group, we do not need the global to local atom index.
951 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
956 /* This vsite is non-home (required for recursion),
957 * and therefore there is no charge group to match pbc with.
958 * But we always turn on full_pbc to assure that higher order
959 * recursion works correctly.
961 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
967 /* Check for recursion */
968 for (k = 2; k < 1+nral; k++)
970 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
972 /* This construction atoms is a vsite and not a home atom */
975 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
977 /* Find the vsite construction */
979 /* Check all interactions assigned to this atom */
980 j = index[iatoms[k]];
981 while (j < index[iatoms[k]+1])
984 nral_r = NRAL(ftype_r);
985 if (interaction_function[ftype_r].flags & IF_VSITE)
987 /* Add this vsite (recursion) */
988 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
989 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
990 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
1003 static void make_la2lc(gmx_domdec_t *dd)
1005 int *cgindex, *la2lc, cg, a;
1007 cgindex = dd->cgindex;
1009 if (dd->nat_tot > dd->la2lc_nalloc)
1011 dd->la2lc_nalloc = over_alloc_dd(dd->nat_tot);
1012 snew(dd->la2lc, dd->la2lc_nalloc);
1016 /* Make the local atom to local cg index */
1017 for (cg = 0; cg < dd->ncg_tot; cg++)
1019 for (a = cgindex[cg]; a < cgindex[cg+1]; a++)
1026 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1032 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1036 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1042 /* Append the nsrc t_blocka block structures in src to *dest */
1043 static void combine_blocka(t_blocka *dest, const t_blocka *src, int nsrc)
1047 ni = src[nsrc-1].nr;
1049 for (s = 0; s < nsrc; s++)
1053 if (ni + 1 > dest->nalloc_index)
1055 dest->nalloc_index = over_alloc_large(ni+1);
1056 srenew(dest->index, dest->nalloc_index);
1058 if (dest->nra + na > dest->nalloc_a)
1060 dest->nalloc_a = over_alloc_large(dest->nra+na);
1061 srenew(dest->a, dest->nalloc_a);
1063 for (s = 0; s < nsrc; s++)
1065 for (i = dest->nr+1; i < src[s].nr+1; i++)
1067 dest->index[i] = dest->nra + src[s].index[i];
1069 for (i = 0; i < src[s].nra; i++)
1071 dest->a[dest->nra+i] = src[s].a[i];
1073 dest->nr = src[s].nr;
1074 dest->nra += src[s].nra;
1078 /* Append the nsrc t_idef structures in src to *dest,
1079 * virtual sites need special attention, as pbc info differs per vsite.
1081 static void combine_idef(t_idef *dest, const t_idef *src, int nsrc,
1082 gmx_vsite_t *vsite, int ***vsite_pbc_t)
1088 int nral1 = 0, ftv = 0;
1090 for (ftype = 0; ftype < F_NRE; ftype++)
1093 for (s = 0; s < nsrc; s++)
1095 n += src[s].il[ftype].nr;
1099 ild = &dest->il[ftype];
1101 if (ild->nr + n > ild->nalloc)
1103 ild->nalloc = over_alloc_large(ild->nr+n);
1104 srenew(ild->iatoms, ild->nalloc);
1107 vpbc = ((interaction_function[ftype].flags & IF_VSITE) &&
1108 vsite->vsite_pbc_loc != NULL);
1111 nral1 = 1 + NRAL(ftype);
1112 ftv = ftype - F_VSITE2;
1113 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1115 vsite->vsite_pbc_loc_nalloc[ftv] =
1116 over_alloc_large((ild->nr + n)/nral1);
1117 srenew(vsite->vsite_pbc_loc[ftv],
1118 vsite->vsite_pbc_loc_nalloc[ftv]);
1122 for (s = 0; s < nsrc; s++)
1124 ils = &src[s].il[ftype];
1125 for (i = 0; i < ils->nr; i++)
1127 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1131 for (i = 0; i < ils->nr; i += nral1)
1133 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1134 vsite_pbc_t[s][ftv][i/nral1];
1143 /* Position restraints need an additional treatment */
1144 if (dest->il[F_POSRES].nr > 0)
1146 n = dest->il[F_POSRES].nr/2;
1147 if (n > dest->iparams_posres_nalloc)
1149 dest->iparams_posres_nalloc = over_alloc_large(n);
1150 srenew(dest->iparams_posres, dest->iparams_posres_nalloc);
1152 /* Set n to the number of original position restraints in dest */
1153 for (s = 0; s < nsrc; s++)
1155 n -= src[s].il[F_POSRES].nr/2;
1157 for (s = 0; s < nsrc; s++)
1159 for (i = 0; i < src[s].il[F_POSRES].nr/2; i++)
1161 /* Correct the index into iparams_posres */
1162 dest->il[F_POSRES].iatoms[n*2] = n;
1163 /* Copy the position restraint force parameters */
1164 dest->iparams_posres[n] = src[s].iparams_posres[i];
1171 /* This function looks up and assigns bonded interactions for zone iz.
1172 * With thread parallelizing each thread acts on a different atom range:
1173 * at_start to at_end.
1175 static int make_bondeds_zone(gmx_domdec_t *dd,
1176 const gmx_domdec_zones_t *zones,
1177 const gmx_molblock_t *molb,
1178 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1180 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1181 const t_iparams *ip_in,
1184 int *vsite_pbc_nalloc,
1186 int at_start, int at_end)
1188 int i, i_gl, mb, mt, mol, i_mol, j, ftype, nral, d, k;
1190 t_iatom *iatoms, tiatoms[1+MAXATOMLIST];
1191 gmx_bool bBCheck, bUse, bLocal;
1192 ivec k_zero, k_plus;
1197 const gmx_domdec_ns_ranges_t *izone;
1198 gmx_reverse_top_t *rt;
1201 nizone = zones->nizone;
1202 izone = zones->izone;
1204 rt = dd->reverse_top;
1206 bBCheck = rt->bBCheck;
1212 for (i = at_start; i < at_end; i++)
1214 /* Get the global atom number */
1215 i_gl = dd->gatindex[i];
1216 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1217 /* Check all interactions assigned to this atom */
1218 index = rt->ril_mt[mt].index;
1219 rtil = rt->ril_mt[mt].il;
1221 while (j < index[i_mol+1])
1226 if (ftype == F_SETTLE)
1228 /* Settles are only in the reverse top when they
1229 * operate within a charge group. So we can assign
1230 * them without checks. We do this only for performance
1231 * reasons; it could be handled by the code below.
1235 /* Home zone: add this settle to the local topology */
1236 tiatoms[0] = iatoms[0];
1238 tiatoms[2] = i + iatoms[2] - iatoms[1];
1239 tiatoms[3] = i + iatoms[3] - iatoms[1];
1240 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1245 else if (interaction_function[ftype].flags & IF_VSITE)
1247 /* The vsite construction goes where the vsite itself is */
1250 add_vsite(dd->ga2la, index, rtil, ftype, nral,
1251 TRUE, i, i_gl, i_mol,
1252 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1259 tiatoms[0] = iatoms[0];
1263 /* Assign single-body interactions to the home zone */
1268 if (ftype == F_POSRES)
1270 add_posres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1273 else if (ftype == F_FBPOSRES)
1275 add_fbposres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1286 /* This is a two-body interaction, we can assign
1287 * analogous to the non-bonded assignments.
1289 if (!ga2la_get(ga2la, i_gl+iatoms[2]-i_mol, &a_loc, &kz))
1299 /* Check zone interaction assignments */
1300 bUse = ((iz < nizone && iz <= kz &&
1301 izone[iz].j0 <= kz && kz < izone[iz].j1) ||
1302 (kz < nizone &&iz > kz &&
1303 izone[kz].j0 <= iz && iz < izone[kz].j1));
1308 /* If necessary check the cgcm distance */
1310 dd_dist2(pbc_null, cg_cm, la2lc,
1311 tiatoms[1], tiatoms[2]) >= rc2)
1320 /* Assign this multi-body bonded interaction to
1321 * the local node if we have all the atoms involved
1322 * (local or communicated) and the minimum zone shift
1323 * in each dimension is zero, for dimensions
1324 * with 2 DD cells an extra check may be necessary.
1329 for (k = 1; k <= nral && bUse; k++)
1331 bLocal = ga2la_get(ga2la, i_gl+iatoms[k]-i_mol,
1333 if (!bLocal || kz >= zones->n)
1335 /* We do not have this atom of this interaction
1336 * locally, or it comes from more than one cell
1344 for (d = 0; d < DIM; d++)
1346 if (zones->shift[kz][d] == 0)
1358 k_zero[XX] && k_zero[YY] && k_zero[ZZ]);
1361 for (d = 0; (d < DIM && bUse); d++)
1363 /* Check if the cg_cm distance falls within
1364 * the cut-off to avoid possible multiple
1365 * assignments of bonded interactions.
1369 dd_dist2(pbc_null, cg_cm, la2lc,
1370 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1379 /* Add this interaction to the local topology */
1380 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1381 /* Sum so we can check in global_stat
1382 * if we have everything.
1385 !(interaction_function[ftype].flags & IF_LIMZERO))
1395 return nbonded_local;
1398 static void set_no_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1399 int iz, t_blocka *lexcls)
1403 a0 = dd->cgindex[zones->cg_range[iz]];
1404 a1 = dd->cgindex[zones->cg_range[iz+1]];
1406 for (a = a0+1; a < a1+1; a++)
1408 lexcls->index[a] = lexcls->nra;
1412 static int make_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1413 const gmx_moltype_t *moltype,
1414 gmx_bool bRCheck, real rc2,
1415 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1419 int cg_start, int cg_end)
1421 int nizone, n, count, jla0, jla1, jla;
1422 int cg, la0, la1, la, a_gl, mb, mt, mol, a_mol, j, aj_mol;
1423 const t_blocka *excls;
1430 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1431 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1433 /* We set the end index, but note that we might not start at zero here */
1434 lexcls->nr = dd->cgindex[cg_end];
1438 for (cg = cg_start; cg < cg_end; cg++)
1440 /* Here we assume the number of exclusions in one charge group
1441 * is never larger than 1000.
1443 if (n+1000 > lexcls->nalloc_a)
1445 lexcls->nalloc_a = over_alloc_large(n+1000);
1446 srenew(lexcls->a, lexcls->nalloc_a);
1448 la0 = dd->cgindex[cg];
1449 la1 = dd->cgindex[cg+1];
1450 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1451 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1453 /* Copy the exclusions from the global top */
1454 for (la = la0; la < la1; la++)
1456 lexcls->index[la] = n;
1457 a_gl = dd->gatindex[la];
1458 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1459 excls = &moltype[mt].excls;
1460 for (j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1462 aj_mol = excls->a[j];
1463 /* This computation of jla is only correct intra-cg */
1464 jla = la + aj_mol - a_mol;
1465 if (jla >= la0 && jla < la1)
1467 /* This is an intra-cg exclusion. We can skip
1468 * the global indexing and distance checking.
1470 /* Intra-cg exclusions are only required
1471 * for the home zone.
1475 lexcls->a[n++] = jla;
1476 /* Check to avoid double counts */
1485 /* This is a inter-cg exclusion */
1486 /* Since exclusions are pair interactions,
1487 * just like non-bonded interactions,
1488 * they can be assigned properly up
1489 * to the DD cutoff (not cutoff_min as
1490 * for the other bonded interactions).
1492 if (ga2la_get(ga2la, a_gl+aj_mol-a_mol, &jla, &cell))
1494 if (iz == 0 && cell == 0)
1496 lexcls->a[n++] = jla;
1497 /* Check to avoid double counts */
1503 else if (jla >= jla0 && jla < jla1 &&
1505 dd_dist2(pbc_null, cg_cm, la2lc, la, jla) < rc2))
1507 /* jla > la, since jla0 > la */
1508 lexcls->a[n++] = jla;
1518 /* There are no inter-cg excls and this cg is self-excluded.
1519 * These exclusions are only required for zone 0,
1520 * since other zones do not see themselves.
1524 for (la = la0; la < la1; la++)
1526 lexcls->index[la] = n;
1527 for (j = la0; j < la1; j++)
1532 count += ((la1 - la0)*(la1 - la0 - 1))/2;
1536 /* We don't need exclusions for this cg */
1537 for (la = la0; la < la1; la++)
1539 lexcls->index[la] = n;
1545 lexcls->index[lexcls->nr] = n;
1551 static void check_alloc_index(t_blocka *ba, int nindex_max)
1553 if (nindex_max+1 > ba->nalloc_index)
1555 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1556 srenew(ba->index, ba->nalloc_index);
1560 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1566 nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1568 check_alloc_index(lexcls, nr);
1570 for (thread = 1; thread < dd->reverse_top->nthread; thread++)
1572 check_alloc_index(&dd->reverse_top->excl_thread[thread], nr);
1576 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1581 lexcls->nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1583 if (dd->n_intercg_excl == 0)
1585 /* There are no exclusions involving non-home charge groups,
1586 * but we need to set the indices for neighborsearching.
1588 la0 = dd->cgindex[zones->izone[0].cg1];
1589 for (la = la0; la < lexcls->nr; la++)
1591 lexcls->index[la] = lexcls->nra;
1594 /* nr is only used to loop over the exclusions for Ewald and RF,
1595 * so we can set it to the number of home atoms for efficiency.
1597 lexcls->nr = dd->cgindex[zones->izone[0].cg1];
1601 static void clear_idef(t_idef *idef)
1605 /* Clear the counts */
1606 for (ftype = 0; ftype < F_NRE; ftype++)
1608 idef->il[ftype].nr = 0;
1612 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1613 gmx_domdec_zones_t *zones,
1614 const gmx_mtop_t *mtop,
1616 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1618 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1619 t_idef *idef, gmx_vsite_t *vsite,
1620 t_blocka *lexcls, int *excl_count)
1622 int nzone_bondeds, nzone_excl;
1627 gmx_reverse_top_t *rt;
1629 if (dd->reverse_top->bMultiCGmols)
1631 nzone_bondeds = zones->n;
1635 /* Only single charge group molecules, so interactions don't
1636 * cross zone boundaries and we only need to assign in the home zone.
1641 if (dd->n_intercg_excl > 0)
1643 /* We only use exclusions from i-zones to i- and j-zones */
1644 nzone_excl = zones->nizone;
1648 /* There are no inter-cg exclusions and only zone 0 sees itself */
1652 check_exclusions_alloc(dd, zones, lexcls);
1654 rt = dd->reverse_top;
1658 /* Clear the counts */
1666 for (iz = 0; iz < nzone_bondeds; iz++)
1668 cg0 = zones->cg_range[iz];
1669 cg1 = zones->cg_range[iz+1];
1671 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1672 for (thread = 0; thread < rt->nthread; thread++)
1678 int *vsite_pbc_nalloc;
1681 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
1682 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
1690 idef_t = &rt->idef_thread[thread];
1694 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1698 vsite_pbc = vsite->vsite_pbc_loc;
1699 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
1703 vsite_pbc = rt->vsite_pbc[thread];
1704 vsite_pbc_nalloc = rt->vsite_pbc_nalloc[thread];
1710 vsite_pbc_nalloc = NULL;
1713 rt->nbonded_thread[thread] =
1714 make_bondeds_zone(dd, zones,
1716 bRCheckMB, rcheck, bRCheck2B, rc2,
1717 la2lc, pbc_null, cg_cm, idef->iparams,
1719 vsite_pbc, vsite_pbc_nalloc,
1721 dd->cgindex[cg0t], dd->cgindex[cg1t]);
1723 if (iz < nzone_excl)
1731 excl_t = &rt->excl_thread[thread];
1736 rt->excl_count_thread[thread] =
1737 make_exclusions_zone(dd, zones,
1738 mtop->moltype, bRCheck2B, rc2,
1739 la2lc, pbc_null, cg_cm, cginfo,
1746 if (rt->nthread > 1)
1748 combine_idef(idef, rt->idef_thread+1, rt->nthread-1,
1749 vsite, rt->vsite_pbc+1);
1752 for (thread = 0; thread < rt->nthread; thread++)
1754 nbonded_local += rt->nbonded_thread[thread];
1757 if (iz < nzone_excl)
1759 if (rt->nthread > 1)
1761 combine_blocka(lexcls, rt->excl_thread+1, rt->nthread-1);
1764 for (thread = 0; thread < rt->nthread; thread++)
1766 *excl_count += rt->excl_count_thread[thread];
1771 /* Some zones might not have exclusions, but some code still needs to
1772 * loop over the index, so we set the indices here.
1774 for (iz = nzone_excl; iz < zones->nizone; iz++)
1776 set_no_exclusions_zone(dd, zones, iz, lexcls);
1779 finish_local_exclusions(dd, zones, lexcls);
1782 fprintf(debug, "We have %d exclusions, check count %d\n",
1783 lexcls->nra, *excl_count);
1786 return nbonded_local;
1789 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
1791 lcgs->nr = dd->ncg_tot;
1792 lcgs->index = dd->cgindex;
1795 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1796 int npbcdim, matrix box,
1797 rvec cellsize_min, ivec npulse,
1801 gmx_mtop_t *mtop, gmx_localtop_t *ltop)
1803 gmx_bool bUniqueExcl, bRCheckMB, bRCheck2B, bRCheckExcl;
1807 t_pbc pbc, *pbc_null = NULL;
1811 fprintf(debug, "Making local topology\n");
1814 dd_make_local_cgs(dd, <op->cgs);
1818 bRCheckExcl = FALSE;
1820 if (dd->reverse_top->bMultiCGmols)
1822 /* We need to check to which cell bondeds should be assigned */
1823 rc = dd_cutoff_twobody(dd);
1826 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1829 /* Should we check cg_cm distances when assigning bonded interactions? */
1830 for (d = 0; d < DIM; d++)
1833 /* Only need to check for dimensions where the part of the box
1834 * that is not communicated is smaller than the cut-off.
1836 if (d < npbcdim && dd->nc[d] > 1 &&
1837 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
1844 /* Check for interactions between two atoms,
1845 * where we can allow interactions up to the cut-off,
1846 * instead of up to the smallest cell dimension.
1853 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
1854 d, cellsize_min[d], d, rcheck[d], bRCheck2B);
1857 if (dd->reverse_top->bExclRequired)
1859 bRCheckExcl = bRCheck2B;
1863 /* If we don't have forces on exclusions,
1864 * we don't care about exclusions being assigned mulitple times.
1866 bRCheckExcl = FALSE;
1868 if (bRCheckMB || bRCheck2B)
1873 set_pbc_dd(&pbc, fr->ePBC, dd, TRUE, box);
1884 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
1885 bRCheckMB, rcheck, bRCheck2B, rc,
1887 pbc_null, cgcm_or_x,
1889 <op->excls, &nexcl);
1891 /* The ilist is not sorted yet,
1892 * we can only do this when we have the charge arrays.
1894 ltop->idef.ilsort = ilsortUNKNOWN;
1896 if (dd->reverse_top->bExclRequired)
1898 dd->nbonded_local += nexcl;
1900 forcerec_set_excl_load(fr, ltop);
1903 ltop->atomtypes = mtop->atomtypes;
1905 /* For an error message only */
1906 dd->reverse_top->err_top_global = mtop;
1907 dd->reverse_top->err_top_local = ltop;
1910 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
1911 gmx_localtop_t *ltop)
1913 if (dd->reverse_top->ilsort == ilsortNO_FE)
1915 ltop->idef.ilsort = ilsortNO_FE;
1919 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
1923 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global)
1925 gmx_localtop_t *top;
1930 top->idef.ntypes = top_global->ffparams.ntypes;
1931 top->idef.atnr = top_global->ffparams.atnr;
1932 top->idef.functype = top_global->ffparams.functype;
1933 top->idef.iparams = top_global->ffparams.iparams;
1934 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
1935 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
1937 for (i = 0; i < F_NRE; i++)
1939 top->idef.il[i].iatoms = NULL;
1940 top->idef.il[i].nalloc = 0;
1942 top->idef.ilsort = ilsortUNKNOWN;
1947 void dd_init_local_state(gmx_domdec_t *dd,
1948 t_state *state_global, t_state *state_local)
1950 int buf[NITEM_DD_INIT_LOCAL_STATE];
1954 buf[0] = state_global->flags;
1955 buf[1] = state_global->ngtc;
1956 buf[2] = state_global->nnhpres;
1957 buf[3] = state_global->nhchainlength;
1958 buf[4] = state_global->dfhist.nlambda;
1960 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
1962 init_state(state_local, 0, buf[1], buf[2], buf[3], buf[4]);
1963 state_local->flags = buf[0];
1966 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
1972 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
1974 if (link->a[k] == cg_gl_j)
1981 /* Add this charge group link */
1982 if (link->index[cg_gl+1]+1 > link->nalloc_a)
1984 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
1985 srenew(link->a, link->nalloc_a);
1987 link->a[link->index[cg_gl+1]] = cg_gl_j;
1988 link->index[cg_gl+1]++;
1992 static int *make_at2cg(t_block *cgs)
1996 snew(at2cg, cgs->index[cgs->nr]);
1997 for (cg = 0; cg < cgs->nr; cg++)
1999 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2008 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
2009 cginfo_mb_t *cginfo_mb)
2011 gmx_reverse_top_t *rt;
2012 int mb, cg_offset, cg, cg_gl, a, aj, i, j, ftype, nral, nlink_mol, mol, ncgi;
2013 gmx_molblock_t *molb;
2014 gmx_moltype_t *molt;
2018 gmx_reverse_ilist_t ril;
2020 cginfo_mb_t *cgi_mb;
2022 /* For each charge group make a list of other charge groups
2023 * in the system that a linked to it via bonded interactions
2024 * which are also stored in reverse_top.
2027 rt = dd->reverse_top;
2030 snew(link->index, ncg_mtop(mtop)+1);
2037 for (mb = 0; mb < mtop->nmolblock; mb++)
2039 molb = &mtop->molblock[mb];
2040 if (molb->nmol == 0)
2044 molt = &mtop->moltype[molb->type];
2046 excls = &molt->excls;
2047 a2c = make_at2cg(cgs);
2048 /* Make a reverse ilist in which the interactions are linked
2049 * to all atoms, not only the first atom as in gmx_reverse_top.
2050 * The constraints are discarded here.
2052 make_reverse_ilist(molt, NULL, FALSE, FALSE, FALSE, TRUE, &ril);
2054 cgi_mb = &cginfo_mb[mb];
2056 for (cg = 0; cg < cgs->nr; cg++)
2058 cg_gl = cg_offset + cg;
2059 link->index[cg_gl+1] = link->index[cg_gl];
2060 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2063 while (i < ril.index[a+1])
2065 ftype = ril.il[i++];
2067 /* Skip the ifunc index */
2069 for (j = 0; j < nral; j++)
2074 check_link(link, cg_gl, cg_offset+a2c[aj]);
2077 i += nral_rt(ftype);
2079 if (rt->bExclRequired)
2081 /* Exclusions always go both ways */
2082 for (j = excls->index[a]; j < excls->index[a+1]; j++)
2087 check_link(link, cg_gl, cg_offset+a2c[aj]);
2092 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2094 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2098 nlink_mol = link->index[cg_offset+cgs->nr] - link->index[cg_offset];
2100 cg_offset += cgs->nr;
2102 destroy_reverse_ilist(&ril);
2107 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt->name, cgs->nr, nlink_mol);
2112 /* Copy the data for the rest of the molecules in this block */
2113 link->nalloc_a += (molb->nmol - 1)*nlink_mol;
2114 srenew(link->a, link->nalloc_a);
2115 for (mol = 1; mol < molb->nmol; mol++)
2117 for (cg = 0; cg < cgs->nr; cg++)
2119 cg_gl = cg_offset + cg;
2120 link->index[cg_gl+1] =
2121 link->index[cg_gl+1-cgs->nr] + nlink_mol;
2122 for (j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2124 link->a[j] = link->a[j-nlink_mol] + cgs->nr;
2126 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2127 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2129 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2133 cg_offset += cgs->nr;
2140 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2146 static void bonded_cg_distance_mol(gmx_moltype_t *molt, int *at2cg,
2147 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2148 real *r_2b, int *ft2b, int *a2_1, int *a2_2,
2149 real *r_mb, int *ftmb, int *am_1, int *am_2)
2151 int ftype, nral, i, j, ai, aj, cgi, cgj;
2154 real r2_2b, r2_mb, rij2;
2158 for (ftype = 0; ftype < F_NRE; ftype++)
2160 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2162 il = &molt->ilist[ftype];
2166 for (i = 0; i < il->nr; i += 1+nral)
2168 for (ai = 0; ai < nral; ai++)
2170 cgi = at2cg[il->iatoms[i+1+ai]];
2171 for (aj = 0; aj < nral; aj++)
2173 cgj = at2cg[il->iatoms[i+1+aj]];
2176 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2177 if (nral == 2 && rij2 > r2_2b)
2181 *a2_1 = il->iatoms[i+1+ai];
2182 *a2_2 = il->iatoms[i+1+aj];
2184 if (nral > 2 && rij2 > r2_mb)
2188 *am_1 = il->iatoms[i+1+ai];
2189 *am_2 = il->iatoms[i+1+aj];
2200 excls = &molt->excls;
2201 for (ai = 0; ai < excls->nr; ai++)
2204 for (j = excls->index[ai]; j < excls->index[ai+1]; j++)
2206 cgj = at2cg[excls->a[j]];
2209 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2219 *r_2b = sqrt(r2_2b);
2220 *r_mb = sqrt(r2_mb);
2223 static void get_cgcm_mol(gmx_moltype_t *molt, gmx_ffparams_t *ffparams,
2224 int ePBC, t_graph *graph, matrix box,
2226 rvec *x, rvec *xs, rvec *cg_cm)
2230 if (ePBC != epbcNONE)
2232 mk_mshift(NULL, graph, ePBC, box, x);
2234 shift_x(graph, box, x, xs);
2235 /* By doing an extra mk_mshift the molecules that are broken
2236 * because they were e.g. imported from another software
2237 * will be made whole again. Such are the healing powers
2240 mk_mshift(NULL, graph, ePBC, box, xs);
2244 /* We copy the coordinates so the original coordinates remain
2245 * unchanged, just to be 100% sure that we do not affect
2246 * binary reproducibility of simulations.
2248 n = molt->cgs.index[molt->cgs.nr];
2249 for (i = 0; i < n; i++)
2251 copy_rvec(x[i], xs[i]);
2257 construct_vsites(vsite, xs, 0.0, NULL,
2258 ffparams->iparams, molt->ilist,
2259 epbcNONE, TRUE, NULL, NULL);
2262 calc_cgcm(NULL, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2265 static int have_vsite_molt(gmx_moltype_t *molt)
2271 for (i = 0; i < F_NRE; i++)
2273 if ((interaction_function[i].flags & IF_VSITE) &&
2274 molt->ilist[i].nr > 0)
2283 void dd_bonded_cg_distance(FILE *fplog,
2285 t_inputrec *ir, rvec *x, matrix box,
2287 real *r_2b, real *r_mb)
2289 gmx_bool bExclRequired;
2290 int mb, cg_offset, at_offset, *at2cg, mol;
2293 gmx_molblock_t *molb;
2294 gmx_moltype_t *molt;
2296 real rmol_2b, rmol_mb;
2297 int ft2b = -1, a_2b_1 = -1, a_2b_2 = -1, ftmb = -1, a_mb_1 = -1, a_mb_2 = -1;
2298 int ftm2b = -1, amol_2b_1 = -1, amol_2b_2 = -1, ftmmb = -1, amol_mb_1 = -1, amol_mb_2 = -1;
2300 bExclRequired = IR_EXCL_FORCES(*ir);
2302 vsite = init_vsite(mtop, NULL, TRUE);
2308 for (mb = 0; mb < mtop->nmolblock; mb++)
2310 molb = &mtop->molblock[mb];
2311 molt = &mtop->moltype[molb->type];
2312 if (molt->cgs.nr == 1 || molb->nmol == 0)
2314 cg_offset += molb->nmol*molt->cgs.nr;
2315 at_offset += molb->nmol*molt->atoms.nr;
2319 if (ir->ePBC != epbcNONE)
2321 mk_graph_ilist(NULL, molt->ilist, 0, molt->atoms.nr, FALSE, FALSE,
2325 at2cg = make_at2cg(&molt->cgs);
2326 snew(xs, molt->atoms.nr);
2327 snew(cg_cm, molt->cgs.nr);
2328 for (mol = 0; mol < molb->nmol; mol++)
2330 get_cgcm_mol(molt, &mtop->ffparams, ir->ePBC, &graph, box,
2331 have_vsite_molt(molt) ? vsite : NULL,
2332 x+at_offset, xs, cg_cm);
2334 bonded_cg_distance_mol(molt, at2cg, bBCheck, bExclRequired, cg_cm,
2335 &rmol_2b, &ftm2b, &amol_2b_1, &amol_2b_2,
2336 &rmol_mb, &ftmmb, &amol_mb_1, &amol_mb_2);
2337 if (rmol_2b > *r_2b)
2341 a_2b_1 = at_offset + amol_2b_1;
2342 a_2b_2 = at_offset + amol_2b_2;
2344 if (rmol_mb > *r_mb)
2348 a_mb_1 = at_offset + amol_mb_1;
2349 a_mb_2 = at_offset + amol_mb_2;
2352 cg_offset += molt->cgs.nr;
2353 at_offset += molt->atoms.nr;
2358 if (ir->ePBC != epbcNONE)
2365 /* We should have a vsite free routine, but here we can simply free */
2368 if (fplog && (ft2b >= 0 || ftmb >= 0))
2371 "Initial maximum inter charge-group distances:\n");
2375 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2376 *r_2b, interaction_function[ft2b].longname,
2377 a_2b_1+1, a_2b_2+1);
2382 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2383 *r_mb, interaction_function[ftmb].longname,
2384 a_mb_1+1, a_mb_2+1);