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.
44 #include "domdec_network.h"
49 #include "chargegroup.h"
50 #include "gromacs/gmxlib/topsort.h"
51 #include "mtop_util.h"
54 #include "gmx_ga2la.h"
56 #include "gmx_omp_nthreads.h"
58 /* for dd_init_local_state */
59 #define NITEM_DD_INIT_LOCAL_STATE 5
62 int *index; /* Index for each atom into il */
63 int *il; /* ftype|type|a0|...|an|ftype|... */
64 } gmx_reverse_ilist_t;
73 typedef struct gmx_reverse_top {
74 gmx_bool bExclRequired; /* Do we require all exclusions to be assigned? */
75 gmx_bool bConstr; /* Are there constraints in this revserse top? */
76 gmx_bool bSettle; /* Are there settles in this revserse top? */
77 gmx_bool bBCheck; /* All bonded interactions have to be assigned? */
78 gmx_bool bMultiCGmols; /* Are the multi charge-group molecules? */
79 gmx_reverse_ilist_t *ril_mt; /* Reverse ilist for all moltypes */
81 int ilsort; /* The sorting state of bondeds for free energy */
82 gmx_molblock_ind_t *mbi;
85 /* Work data structures for multi-threading */
89 int **vsite_pbc_nalloc;
91 t_blocka *excl_thread;
92 int *excl_count_thread;
94 /* Pointers only used for an error message */
95 gmx_mtop_t *err_top_global;
96 gmx_localtop_t *err_top_local;
99 static int nral_rt(int ftype)
101 /* Returns the number of atom entries for il in gmx_reverse_top_t */
105 if (interaction_function[ftype].flags & IF_VSITE)
107 /* With vsites the reverse topology contains
108 * two extra entries for PBC.
116 /* This function tells which interactions need to be assigned exactly once */
117 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
118 gmx_bool bConstr, gmx_bool bSettle)
120 return (((interaction_function[ftype].flags & IF_BOND) &&
121 !(interaction_function[ftype].flags & IF_VSITE) &&
122 (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))) ||
123 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
124 (bSettle && ftype == F_SETTLE));
127 static void print_error_header(FILE *fplog, char *moltypename, int nprint)
129 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
130 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
132 "the first %d missing interactions, except for exclusions:\n",
135 "the first %d missing interactions, except for exclusions:\n",
139 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
140 gmx_reverse_top_t *rt,
142 gmx_reverse_ilist_t *ril,
143 int a_start, int a_end,
144 int nat_mol, int nmol,
147 int nril_mol, *assigned, *gatindex;
148 int ftype, ftype_j, nral, i, j_mol, j, k, a0, a0_mol, mol, a, a_gl;
154 nril_mol = ril->index[nat_mol];
155 snew(assigned, nmol*nril_mol);
157 gatindex = cr->dd->gatindex;
158 for (ftype = 0; ftype < F_NRE; ftype++)
160 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
163 il = &idef->il[ftype];
165 for (i = 0; i < il->nr; i += 1+nral)
167 a0 = gatindex[ia[1]];
168 /* Check if this interaction is in
169 * the currently checked molblock.
171 if (a0 >= a_start && a0 < a_end)
173 mol = (a0 - a_start)/nat_mol;
174 a0_mol = (a0 - a_start) - mol*nat_mol;
175 j_mol = ril->index[a0_mol];
177 while (j_mol < ril->index[a0_mol+1] && !bFound)
179 j = mol*nril_mol + j_mol;
180 ftype_j = ril->il[j_mol];
181 /* Here we need to check if this interaction has
182 * not already been assigned, since we could have
183 * multiply defined interactions.
185 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
188 /* Check the atoms */
190 for (a = 0; a < nral; a++)
192 if (gatindex[ia[1+a]] !=
193 a_start + mol*nat_mol + ril->il[j_mol+2+a])
203 j_mol += 2 + nral_rt(ftype_j);
207 gmx_incons("Some interactions seem to be assigned multiple times");
215 gmx_sumi(nmol*nril_mol, assigned, cr);
219 for (mol = 0; mol < nmol; mol++)
222 while (j_mol < nril_mol)
224 ftype = ril->il[j_mol];
226 j = mol*nril_mol + j_mol;
227 if (assigned[j] == 0 &&
228 !(interaction_function[ftype].flags & IF_VSITE))
230 if (DDMASTER(cr->dd))
234 print_error_header(fplog, moltypename, nprint);
236 fprintf(fplog, "%20s atoms",
237 interaction_function[ftype].longname);
238 fprintf(stderr, "%20s atoms",
239 interaction_function[ftype].longname);
240 for (a = 0; a < nral; a++)
242 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
243 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
248 fprintf(stderr, " ");
251 fprintf(fplog, " global");
252 fprintf(stderr, " global");
253 for (a = 0; a < nral; a++)
255 fprintf(fplog, "%6d",
256 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
257 fprintf(stderr, "%6d",
258 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
260 fprintf(fplog, "\n");
261 fprintf(stderr, "\n");
269 j_mol += 2 + nral_rt(ftype);
276 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
277 gmx_mtop_t *mtop, t_idef *idef)
279 int mb, a_start, a_end;
280 gmx_molblock_t *molb;
281 gmx_reverse_top_t *rt;
283 rt = cr->dd->reverse_top;
285 /* Print the atoms in the missing interactions per molblock */
287 for (mb = 0; mb < mtop->nmolblock; mb++)
289 molb = &mtop->molblock[mb];
291 a_end = a_start + molb->nmol*molb->natoms_mol;
293 print_missing_interactions_mb(fplog, cr, rt,
294 *(mtop->moltype[molb->type].name),
295 &rt->ril_mt[molb->type],
296 a_start, a_end, molb->natoms_mol,
302 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr, int local_count, gmx_mtop_t *top_global, t_state *state_local)
304 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
308 gmx_mtop_t *err_top_global;
309 gmx_localtop_t *err_top_local;
313 err_top_global = dd->reverse_top->err_top_global;
314 err_top_local = dd->reverse_top->err_top_local;
318 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
322 ndiff_tot = local_count - dd->nbonded_global;
324 for (ftype = 0; ftype < F_NRE; ftype++)
327 cl[ftype] = err_top_local->idef.il[ftype].nr/(1+nral);
330 gmx_sumi(F_NRE, cl, cr);
334 fprintf(fplog, "\nA list of missing interactions:\n");
335 fprintf(stderr, "\nA list of missing interactions:\n");
336 rest_global = dd->nbonded_global;
337 rest_local = local_count;
338 for (ftype = 0; ftype < F_NRE; ftype++)
340 /* In the reverse and local top all constraints are merged
341 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
342 * and add these constraints when doing F_CONSTR.
344 if (((interaction_function[ftype].flags & IF_BOND) &&
345 (dd->reverse_top->bBCheck
346 || !(interaction_function[ftype].flags & IF_LIMZERO)))
347 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
348 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
351 n = gmx_mtop_ftype_count(err_top_global, ftype);
352 if (ftype == F_CONSTR)
354 n += gmx_mtop_ftype_count(err_top_global, F_CONSTRNC);
356 ndiff = cl[ftype] - n;
359 sprintf(buf, "%20s of %6d missing %6d",
360 interaction_function[ftype].longname, n, -ndiff);
361 fprintf(fplog, "%s\n", buf);
362 fprintf(stderr, "%s\n", buf);
365 rest_local -= cl[ftype];
369 ndiff = rest_local - rest_global;
372 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
373 rest_global, -ndiff);
374 fprintf(fplog, "%s\n", buf);
375 fprintf(stderr, "%s\n", buf);
379 print_missing_interactions_atoms(fplog, cr, err_top_global,
380 &err_top_local->idef);
381 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
382 -1, state_local->x, state_local->box);
387 gmx_incons("One or more interactions were multiple assigned in the domain decompostion");
391 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));
396 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
397 int *mb, int *mt, int *mol, int *i_mol)
402 gmx_molblock_ind_t *mbi = rt->mbi;
404 int end = rt->nmolblock; /* exclusive */
407 /* binary search for molblock_ind */
411 if (i_gl >= mbi[mid].a_end)
415 else if (i_gl < mbi[mid].a_start)
429 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
430 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
433 static int count_excls(t_block *cgs, t_blocka *excls, int *n_intercg_excl)
435 int n, n_inter, cg, at0, at1, at, excl, atj;
439 for (cg = 0; cg < cgs->nr; cg++)
441 at0 = cgs->index[cg];
442 at1 = cgs->index[cg+1];
443 for (at = at0; at < at1; at++)
445 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
447 atj = excls->a[excl];
451 if (atj < at0 || atj >= at1)
463 static int low_make_reverse_ilist(t_ilist *il_mt, t_atom *atom,
466 gmx_bool bConstr, gmx_bool bSettle,
468 int *r_index, int *r_il,
469 gmx_bool bLinkToAllAtoms,
472 int ftype, nral, i, j, nlink, link;
480 for (ftype = 0; ftype < F_NRE; ftype++)
482 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
483 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
484 (bSettle && ftype == F_SETTLE))
486 bVSite = (interaction_function[ftype].flags & IF_VSITE);
490 for (i = 0; i < il->nr; i += 1+nral)
497 /* We don't need the virtual sites for the cg-links */
507 /* Couple to the first atom in the interaction */
510 for (link = 0; link < nlink; link++)
515 r_il[r_index[a]+count[a]] =
516 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
517 r_il[r_index[a]+count[a]+1] = ia[0];
518 for (j = 1; j < 1+nral; j++)
520 /* Store the molecular atom number */
521 r_il[r_index[a]+count[a]+1+j] = ia[j];
524 if (interaction_function[ftype].flags & IF_VSITE)
528 /* Add an entry to iatoms for storing
529 * which of the constructing atoms are
532 r_il[r_index[a]+count[a]+2+nral] = 0;
533 for (j = 2; j < 1+nral; j++)
535 if (atom[ia[j]].ptype == eptVSite)
537 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
540 /* Store vsite pbc atom in a second extra entry */
541 r_il[r_index[a]+count[a]+2+nral+1] =
542 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
547 /* We do not count vsites since they are always
548 * uniquely assigned and can be assigned
549 * to multiple nodes with recursive vsites.
552 !(interaction_function[ftype].flags & IF_LIMZERO))
557 count[a] += 2 + nral_rt(ftype);
566 static int make_reverse_ilist(gmx_moltype_t *molt,
568 gmx_bool bConstr, gmx_bool bSettle,
570 gmx_bool bLinkToAllAtoms,
571 gmx_reverse_ilist_t *ril_mt)
573 int nat_mt, *count, i, nint_mt;
575 /* Count the interactions */
576 nat_mt = molt->atoms.nr;
578 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
580 bConstr, bSettle, bBCheck, NULL, NULL,
581 bLinkToAllAtoms, FALSE);
583 snew(ril_mt->index, nat_mt+1);
584 ril_mt->index[0] = 0;
585 for (i = 0; i < nat_mt; i++)
587 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
590 snew(ril_mt->il, ril_mt->index[nat_mt]);
592 /* Store the interactions */
594 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
596 bConstr, bSettle, bBCheck,
597 ril_mt->index, ril_mt->il,
598 bLinkToAllAtoms, TRUE);
605 static void destroy_reverse_ilist(gmx_reverse_ilist_t *ril)
611 static gmx_reverse_top_t *make_reverse_top(gmx_mtop_t *mtop, gmx_bool bFE,
612 int ***vsite_pbc_molt,
613 gmx_bool bConstr, gmx_bool bSettle,
614 gmx_bool bBCheck, int *nint)
617 gmx_reverse_top_t *rt;
624 /* Should we include constraints (for SHAKE) in rt? */
625 rt->bConstr = bConstr;
626 rt->bSettle = bSettle;
627 rt->bBCheck = bBCheck;
629 rt->bMultiCGmols = FALSE;
630 snew(nint_mt, mtop->nmoltype);
631 snew(rt->ril_mt, mtop->nmoltype);
632 rt->ril_mt_tot_size = 0;
633 for (mt = 0; mt < mtop->nmoltype; mt++)
635 molt = &mtop->moltype[mt];
636 if (molt->cgs.nr > 1)
638 rt->bMultiCGmols = TRUE;
641 /* Make the atom to interaction list for this molecule type */
643 make_reverse_ilist(molt, vsite_pbc_molt ? vsite_pbc_molt[mt] : NULL,
644 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
647 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt->atoms.nr];
651 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
655 for (mb = 0; mb < mtop->nmolblock; mb++)
657 *nint += mtop->molblock[mb].nmol*nint_mt[mtop->molblock[mb].type];
661 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
663 rt->ilsort = ilsortFE_UNSORTED;
667 rt->ilsort = ilsortNO_FE;
670 /* Make a molblock index for fast searching */
671 snew(rt->mbi, mtop->nmolblock);
672 rt->nmolblock = mtop->nmolblock;
674 for (mb = 0; mb < mtop->nmolblock; mb++)
676 rt->mbi[mb].a_start = i;
677 i += mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
678 rt->mbi[mb].a_end = i;
679 rt->mbi[mb].natoms_mol = mtop->molblock[mb].natoms_mol;
680 rt->mbi[mb].type = mtop->molblock[mb].type;
683 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
684 snew(rt->idef_thread, rt->nthread);
685 if (vsite_pbc_molt != NULL)
687 snew(rt->vsite_pbc, rt->nthread);
688 snew(rt->vsite_pbc_nalloc, rt->nthread);
689 for (thread = 0; thread < rt->nthread; thread++)
691 snew(rt->vsite_pbc[thread], F_VSITEN-F_VSITE2+1);
692 snew(rt->vsite_pbc_nalloc[thread], F_VSITEN-F_VSITE2+1);
695 snew(rt->nbonded_thread, rt->nthread);
696 snew(rt->excl_thread, rt->nthread);
697 snew(rt->excl_count_thread, rt->nthread);
702 void dd_make_reverse_top(FILE *fplog,
703 gmx_domdec_t *dd, gmx_mtop_t *mtop,
705 t_inputrec *ir, gmx_bool bBCheck)
707 int mb, n_recursive_vsite, nexcl, nexcl_icg, a;
708 gmx_molblock_t *molb;
713 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
716 /* If normal and/or settle constraints act only within charge groups,
717 * we can store them in the reverse top and simply assign them to domains.
718 * Otherwise we need to assign them to multiple domains and set up
719 * the parallel version constraint algoirthm(s).
722 dd->reverse_top = make_reverse_top(mtop, ir->efep != efepNO,
723 vsite ? vsite->vsite_pbc_molt : NULL,
724 !dd->bInterCGcons, !dd->bInterCGsettles,
725 bBCheck, &dd->nbonded_global);
727 if (dd->reverse_top->ril_mt_tot_size >= 200000 &&
729 mtop->nmolblock == 1 && mtop->molblock[0].nmol == 1)
731 /* mtop comes from a pre Gromacs 4 tpr file */
732 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";
735 fprintf(fplog, "\n%s\n\n", note);
739 fprintf(stderr, "\n%s\n\n", note);
743 dd->reverse_top->bExclRequired = IR_EXCL_FORCES(*ir);
746 dd->n_intercg_excl = 0;
747 for (mb = 0; mb < mtop->nmolblock; mb++)
749 molb = &mtop->molblock[mb];
750 molt = &mtop->moltype[molb->type];
751 nexcl += molb->nmol*count_excls(&molt->cgs, &molt->excls, &nexcl_icg);
752 dd->n_intercg_excl += molb->nmol*nexcl_icg;
754 if (dd->reverse_top->bExclRequired)
756 dd->nbonded_global += nexcl;
757 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
759 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
760 "will use an extra communication step for exclusion forces for %s\n",
761 dd->n_intercg_excl, eel_names[ir->coulombtype]);
765 if (vsite && vsite->n_intercg_vsite > 0)
769 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
770 "will an extra communication step for selected coordinates and forces\n",
771 vsite->n_intercg_vsite);
773 init_domdec_vsites(dd, vsite->n_intercg_vsite);
776 if (dd->bInterCGcons || dd->bInterCGsettles)
778 init_domdec_constraints(dd, mtop);
782 fprintf(fplog, "\n");
786 static gmx_inline void add_ifunc(int nral, t_iatom *tiatoms, t_ilist *il)
791 if (il->nr+1+nral > il->nalloc)
793 il->nalloc = over_alloc_large(il->nr+1+nral);
794 srenew(il->iatoms, il->nalloc);
796 liatoms = il->iatoms + il->nr;
797 for (k = 0; k <= nral; k++)
799 liatoms[k] = tiatoms[k];
804 static void add_posres(int mol, int a_mol, const gmx_molblock_t *molb,
805 t_iatom *iatoms, const t_iparams *ip_in,
811 /* This position restraint has not been added yet,
812 * so it's index is the current number of position restraints.
814 n = idef->il[F_POSRES].nr/2;
815 if (n+1 > idef->iparams_posres_nalloc)
817 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
818 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
820 ip = &idef->iparams_posres[n];
821 /* Copy the force constants */
822 *ip = ip_in[iatoms[0]];
824 /* Get the position restraint coordinates from the molblock */
825 a_molb = mol*molb->natoms_mol + a_mol;
826 if (a_molb >= molb->nposres_xA)
828 gmx_incons("Not enough position restraint coordinates");
830 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
831 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
832 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
833 if (molb->nposres_xB > 0)
835 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
836 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
837 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
841 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
842 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
843 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
845 /* Set the parameter index for idef->iparams_posre */
849 static void add_fbposres(int mol, int a_mol, const gmx_molblock_t *molb,
850 t_iatom *iatoms, const t_iparams *ip_in,
856 /* This flat-bottom position restraint has not been added yet,
857 * so it's index is the current number of position restraints.
859 n = idef->il[F_FBPOSRES].nr/2;
860 if (n+1 > idef->iparams_fbposres_nalloc)
862 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
863 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
865 ip = &idef->iparams_fbposres[n];
866 /* Copy the force constants */
867 *ip = ip_in[iatoms[0]];
869 /* Get the position restriant coordinats from the molblock */
870 a_molb = mol*molb->natoms_mol + a_mol;
871 if (a_molb >= molb->nposres_xA)
873 gmx_incons("Not enough position restraint coordinates");
875 /* Take reference positions from A position of normal posres */
876 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
877 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
878 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
880 /* Note: no B-type for flat-bottom posres */
882 /* Set the parameter index for idef->iparams_posre */
886 static void add_vsite(gmx_ga2la_t ga2la, int *index, int *rtil,
888 gmx_bool bHomeA, int a, int a_gl, int a_mol,
890 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
892 int k, ak_gl, vsi, pbc_a_mol;
893 t_iatom tiatoms[1+MAXATOMLIST], *iatoms_r;
894 int j, ftype_r, nral_r;
897 tiatoms[0] = iatoms[0];
901 /* We know the local index of the first atom */
906 /* Convert later in make_local_vsites */
907 tiatoms[1] = -a_gl - 1;
910 for (k = 2; k < 1+nral; k++)
912 ak_gl = a_gl + iatoms[k] - a_mol;
913 if (!ga2la_get_home(ga2la, ak_gl, &tiatoms[k]))
915 /* Copy the global index, convert later in make_local_vsites */
916 tiatoms[k] = -(ak_gl + 1);
920 /* Add this interaction to the local topology */
921 add_ifunc(nral, tiatoms, &idef->il[ftype]);
924 vsi = idef->il[ftype].nr/(1+nral) - 1;
925 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
927 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
928 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
932 pbc_a_mol = iatoms[1+nral+1];
935 /* The pbc flag is one of the following two options:
936 * -2: vsite and all constructing atoms are within the same cg, no pbc
937 * -1: vsite and its first constructing atom are in the same cg, do pbc
939 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
943 /* Set the pbc atom for this vsite so we can make its pbc
944 * identical to the rest of the atoms in its charge group.
945 * Since the order of the atoms does not change within a charge
946 * group, we do not need the global to local atom index.
948 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
953 /* This vsite is non-home (required for recursion),
954 * and therefore there is no charge group to match pbc with.
955 * But we always turn on full_pbc to assure that higher order
956 * recursion works correctly.
958 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
964 /* Check for recursion */
965 for (k = 2; k < 1+nral; k++)
967 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
969 /* This construction atoms is a vsite and not a home atom */
972 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
974 /* Find the vsite construction */
976 /* Check all interactions assigned to this atom */
977 j = index[iatoms[k]];
978 while (j < index[iatoms[k]+1])
981 nral_r = NRAL(ftype_r);
982 if (interaction_function[ftype_r].flags & IF_VSITE)
984 /* Add this vsite (recursion) */
985 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
986 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
987 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
1000 static void make_la2lc(gmx_domdec_t *dd)
1002 int *cgindex, *la2lc, cg, a;
1004 cgindex = dd->cgindex;
1006 if (dd->nat_tot > dd->la2lc_nalloc)
1008 dd->la2lc_nalloc = over_alloc_dd(dd->nat_tot);
1009 snew(dd->la2lc, dd->la2lc_nalloc);
1013 /* Make the local atom to local cg index */
1014 for (cg = 0; cg < dd->ncg_tot; cg++)
1016 for (a = cgindex[cg]; a < cgindex[cg+1]; a++)
1023 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1029 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1033 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1039 /* Append the nsrc t_blocka block structures in src to *dest */
1040 static void combine_blocka(t_blocka *dest, const t_blocka *src, int nsrc)
1044 ni = src[nsrc-1].nr;
1046 for (s = 0; s < nsrc; s++)
1050 if (ni + 1 > dest->nalloc_index)
1052 dest->nalloc_index = over_alloc_large(ni+1);
1053 srenew(dest->index, dest->nalloc_index);
1055 if (dest->nra + na > dest->nalloc_a)
1057 dest->nalloc_a = over_alloc_large(dest->nra+na);
1058 srenew(dest->a, dest->nalloc_a);
1060 for (s = 0; s < nsrc; s++)
1062 for (i = dest->nr+1; i < src[s].nr+1; i++)
1064 dest->index[i] = dest->nra + src[s].index[i];
1066 for (i = 0; i < src[s].nra; i++)
1068 dest->a[dest->nra+i] = src[s].a[i];
1070 dest->nr = src[s].nr;
1071 dest->nra += src[s].nra;
1075 /* Append the nsrc t_idef structures in src to *dest,
1076 * virtual sites need special attention, as pbc info differs per vsite.
1078 static void combine_idef(t_idef *dest, const t_idef *src, int nsrc,
1079 gmx_vsite_t *vsite, int ***vsite_pbc_t)
1085 int nral1 = 0, ftv = 0;
1087 for (ftype = 0; ftype < F_NRE; ftype++)
1090 for (s = 0; s < nsrc; s++)
1092 n += src[s].il[ftype].nr;
1096 ild = &dest->il[ftype];
1098 if (ild->nr + n > ild->nalloc)
1100 ild->nalloc = over_alloc_large(ild->nr+n);
1101 srenew(ild->iatoms, ild->nalloc);
1104 vpbc = ((interaction_function[ftype].flags & IF_VSITE) &&
1105 vsite->vsite_pbc_loc != NULL);
1108 nral1 = 1 + NRAL(ftype);
1109 ftv = ftype - F_VSITE2;
1110 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1112 vsite->vsite_pbc_loc_nalloc[ftv] =
1113 over_alloc_large((ild->nr + n)/nral1);
1114 srenew(vsite->vsite_pbc_loc[ftv],
1115 vsite->vsite_pbc_loc_nalloc[ftv]);
1119 for (s = 0; s < nsrc; s++)
1121 ils = &src[s].il[ftype];
1122 for (i = 0; i < ils->nr; i++)
1124 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1128 for (i = 0; i < ils->nr; i += nral1)
1130 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1131 vsite_pbc_t[s][ftv][i/nral1];
1140 /* Position restraints need an additional treatment */
1141 if (dest->il[F_POSRES].nr > 0)
1143 n = dest->il[F_POSRES].nr/2;
1144 if (n > dest->iparams_posres_nalloc)
1146 dest->iparams_posres_nalloc = over_alloc_large(n);
1147 srenew(dest->iparams_posres, dest->iparams_posres_nalloc);
1149 /* Set n to the number of original position restraints in dest */
1150 for (s = 0; s < nsrc; s++)
1152 n -= src[s].il[F_POSRES].nr/2;
1154 for (s = 0; s < nsrc; s++)
1156 for (i = 0; i < src[s].il[F_POSRES].nr/2; i++)
1158 /* Correct the index into iparams_posres */
1159 dest->il[F_POSRES].iatoms[n*2] = n;
1160 /* Copy the position restraint force parameters */
1161 dest->iparams_posres[n] = src[s].iparams_posres[i];
1168 /* This function looks up and assigns bonded interactions for zone iz.
1169 * With thread parallelizing each thread acts on a different atom range:
1170 * at_start to at_end.
1172 static int make_bondeds_zone(gmx_domdec_t *dd,
1173 const gmx_domdec_zones_t *zones,
1174 const gmx_molblock_t *molb,
1175 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1177 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1178 const t_iparams *ip_in,
1181 int *vsite_pbc_nalloc,
1183 int at_start, int at_end)
1185 int i, i_gl, mb, mt, mol, i_mol, j, ftype, nral, d, k;
1187 t_iatom *iatoms, tiatoms[1+MAXATOMLIST];
1188 gmx_bool bBCheck, bUse, bLocal;
1189 ivec k_zero, k_plus;
1194 const gmx_domdec_ns_ranges_t *izone;
1195 gmx_reverse_top_t *rt;
1198 nizone = zones->nizone;
1199 izone = zones->izone;
1201 rt = dd->reverse_top;
1203 bBCheck = rt->bBCheck;
1209 for (i = at_start; i < at_end; i++)
1211 /* Get the global atom number */
1212 i_gl = dd->gatindex[i];
1213 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1214 /* Check all interactions assigned to this atom */
1215 index = rt->ril_mt[mt].index;
1216 rtil = rt->ril_mt[mt].il;
1218 while (j < index[i_mol+1])
1223 if (ftype == F_SETTLE)
1225 /* Settles are only in the reverse top when they
1226 * operate within a charge group. So we can assign
1227 * them without checks. We do this only for performance
1228 * reasons; it could be handled by the code below.
1232 /* Home zone: add this settle to the local topology */
1233 tiatoms[0] = iatoms[0];
1235 tiatoms[2] = i + iatoms[2] - iatoms[1];
1236 tiatoms[3] = i + iatoms[3] - iatoms[1];
1237 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1242 else if (interaction_function[ftype].flags & IF_VSITE)
1244 /* The vsite construction goes where the vsite itself is */
1247 add_vsite(dd->ga2la, index, rtil, ftype, nral,
1248 TRUE, i, i_gl, i_mol,
1249 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1256 tiatoms[0] = iatoms[0];
1260 /* Assign single-body interactions to the home zone */
1265 if (ftype == F_POSRES)
1267 add_posres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1270 else if (ftype == F_FBPOSRES)
1272 add_fbposres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1283 /* This is a two-body interaction, we can assign
1284 * analogous to the non-bonded assignments.
1286 if (!ga2la_get(ga2la, i_gl+iatoms[2]-i_mol, &a_loc, &kz))
1296 /* Check zone interaction assignments */
1297 bUse = ((iz < nizone && iz <= kz &&
1298 izone[iz].j0 <= kz && kz < izone[iz].j1) ||
1299 (kz < nizone &&iz > kz &&
1300 izone[kz].j0 <= iz && iz < izone[kz].j1));
1305 /* If necessary check the cgcm distance */
1307 dd_dist2(pbc_null, cg_cm, la2lc,
1308 tiatoms[1], tiatoms[2]) >= rc2)
1317 /* Assign this multi-body bonded interaction to
1318 * the local node if we have all the atoms involved
1319 * (local or communicated) and the minimum zone shift
1320 * in each dimension is zero, for dimensions
1321 * with 2 DD cells an extra check may be necessary.
1326 for (k = 1; k <= nral && bUse; k++)
1328 bLocal = ga2la_get(ga2la, i_gl+iatoms[k]-i_mol,
1330 if (!bLocal || kz >= zones->n)
1332 /* We do not have this atom of this interaction
1333 * locally, or it comes from more than one cell
1341 for (d = 0; d < DIM; d++)
1343 if (zones->shift[kz][d] == 0)
1355 k_zero[XX] && k_zero[YY] && k_zero[ZZ]);
1358 for (d = 0; (d < DIM && bUse); d++)
1360 /* Check if the cg_cm distance falls within
1361 * the cut-off to avoid possible multiple
1362 * assignments of bonded interactions.
1366 dd_dist2(pbc_null, cg_cm, la2lc,
1367 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1376 /* Add this interaction to the local topology */
1377 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1378 /* Sum so we can check in global_stat
1379 * if we have everything.
1382 !(interaction_function[ftype].flags & IF_LIMZERO))
1392 return nbonded_local;
1395 static void set_no_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1396 int iz, t_blocka *lexcls)
1400 a0 = dd->cgindex[zones->cg_range[iz]];
1401 a1 = dd->cgindex[zones->cg_range[iz+1]];
1403 for (a = a0+1; a < a1+1; a++)
1405 lexcls->index[a] = lexcls->nra;
1409 static int make_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1410 const gmx_moltype_t *moltype,
1411 gmx_bool bRCheck, real rc2,
1412 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1416 int cg_start, int cg_end)
1418 int nizone, n, count, jla0, jla1, jla;
1419 int cg, la0, la1, la, a_gl, mb, mt, mol, a_mol, j, aj_mol;
1420 const t_blocka *excls;
1427 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1428 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1430 /* We set the end index, but note that we might not start at zero here */
1431 lexcls->nr = dd->cgindex[cg_end];
1435 for (cg = cg_start; cg < cg_end; cg++)
1437 /* Here we assume the number of exclusions in one charge group
1438 * is never larger than 1000.
1440 if (n+1000 > lexcls->nalloc_a)
1442 lexcls->nalloc_a = over_alloc_large(n+1000);
1443 srenew(lexcls->a, lexcls->nalloc_a);
1445 la0 = dd->cgindex[cg];
1446 la1 = dd->cgindex[cg+1];
1447 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1448 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1450 /* Copy the exclusions from the global top */
1451 for (la = la0; la < la1; la++)
1453 lexcls->index[la] = n;
1454 a_gl = dd->gatindex[la];
1455 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1456 excls = &moltype[mt].excls;
1457 for (j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1459 aj_mol = excls->a[j];
1460 /* This computation of jla is only correct intra-cg */
1461 jla = la + aj_mol - a_mol;
1462 if (jla >= la0 && jla < la1)
1464 /* This is an intra-cg exclusion. We can skip
1465 * the global indexing and distance checking.
1467 /* Intra-cg exclusions are only required
1468 * for the home zone.
1472 lexcls->a[n++] = jla;
1473 /* Check to avoid double counts */
1482 /* This is a inter-cg exclusion */
1483 /* Since exclusions are pair interactions,
1484 * just like non-bonded interactions,
1485 * they can be assigned properly up
1486 * to the DD cutoff (not cutoff_min as
1487 * for the other bonded interactions).
1489 if (ga2la_get(ga2la, a_gl+aj_mol-a_mol, &jla, &cell))
1491 if (iz == 0 && cell == 0)
1493 lexcls->a[n++] = jla;
1494 /* Check to avoid double counts */
1500 else if (jla >= jla0 && jla < jla1 &&
1502 dd_dist2(pbc_null, cg_cm, la2lc, la, jla) < rc2))
1504 /* jla > la, since jla0 > la */
1505 lexcls->a[n++] = jla;
1515 /* There are no inter-cg excls and this cg is self-excluded.
1516 * These exclusions are only required for zone 0,
1517 * since other zones do not see themselves.
1521 for (la = la0; la < la1; la++)
1523 lexcls->index[la] = n;
1524 for (j = la0; j < la1; j++)
1529 count += ((la1 - la0)*(la1 - la0 - 1))/2;
1533 /* We don't need exclusions for this cg */
1534 for (la = la0; la < la1; la++)
1536 lexcls->index[la] = n;
1542 lexcls->index[lexcls->nr] = n;
1548 static void check_alloc_index(t_blocka *ba, int nindex_max)
1550 if (nindex_max+1 > ba->nalloc_index)
1552 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1553 srenew(ba->index, ba->nalloc_index);
1557 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1563 nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1565 check_alloc_index(lexcls, nr);
1567 for (thread = 1; thread < dd->reverse_top->nthread; thread++)
1569 check_alloc_index(&dd->reverse_top->excl_thread[thread], nr);
1573 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1578 lexcls->nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1580 if (dd->n_intercg_excl == 0)
1582 /* There are no exclusions involving non-home charge groups,
1583 * but we need to set the indices for neighborsearching.
1585 la0 = dd->cgindex[zones->izone[0].cg1];
1586 for (la = la0; la < lexcls->nr; la++)
1588 lexcls->index[la] = lexcls->nra;
1591 /* nr is only used to loop over the exclusions for Ewald and RF,
1592 * so we can set it to the number of home atoms for efficiency.
1594 lexcls->nr = dd->cgindex[zones->izone[0].cg1];
1598 static void clear_idef(t_idef *idef)
1602 /* Clear the counts */
1603 for (ftype = 0; ftype < F_NRE; ftype++)
1605 idef->il[ftype].nr = 0;
1609 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1610 gmx_domdec_zones_t *zones,
1611 const gmx_mtop_t *mtop,
1613 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1615 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1616 t_idef *idef, gmx_vsite_t *vsite,
1617 t_blocka *lexcls, int *excl_count)
1619 int nzone_bondeds, nzone_excl;
1624 gmx_reverse_top_t *rt;
1626 if (dd->reverse_top->bMultiCGmols)
1628 nzone_bondeds = zones->n;
1632 /* Only single charge group molecules, so interactions don't
1633 * cross zone boundaries and we only need to assign in the home zone.
1638 if (dd->n_intercg_excl > 0)
1640 /* We only use exclusions from i-zones to i- and j-zones */
1641 nzone_excl = zones->nizone;
1645 /* There are no inter-cg exclusions and only zone 0 sees itself */
1649 check_exclusions_alloc(dd, zones, lexcls);
1651 rt = dd->reverse_top;
1655 /* Clear the counts */
1663 for (iz = 0; iz < nzone_bondeds; iz++)
1665 cg0 = zones->cg_range[iz];
1666 cg1 = zones->cg_range[iz+1];
1668 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1669 for (thread = 0; thread < rt->nthread; thread++)
1675 int *vsite_pbc_nalloc;
1678 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
1679 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
1687 idef_t = &rt->idef_thread[thread];
1691 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1695 vsite_pbc = vsite->vsite_pbc_loc;
1696 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
1700 vsite_pbc = rt->vsite_pbc[thread];
1701 vsite_pbc_nalloc = rt->vsite_pbc_nalloc[thread];
1707 vsite_pbc_nalloc = NULL;
1710 rt->nbonded_thread[thread] =
1711 make_bondeds_zone(dd, zones,
1713 bRCheckMB, rcheck, bRCheck2B, rc2,
1714 la2lc, pbc_null, cg_cm, idef->iparams,
1716 vsite_pbc, vsite_pbc_nalloc,
1718 dd->cgindex[cg0t], dd->cgindex[cg1t]);
1720 if (iz < nzone_excl)
1728 excl_t = &rt->excl_thread[thread];
1733 rt->excl_count_thread[thread] =
1734 make_exclusions_zone(dd, zones,
1735 mtop->moltype, bRCheck2B, rc2,
1736 la2lc, pbc_null, cg_cm, cginfo,
1743 if (rt->nthread > 1)
1745 combine_idef(idef, rt->idef_thread+1, rt->nthread-1,
1746 vsite, rt->vsite_pbc+1);
1749 for (thread = 0; thread < rt->nthread; thread++)
1751 nbonded_local += rt->nbonded_thread[thread];
1754 if (iz < nzone_excl)
1756 if (rt->nthread > 1)
1758 combine_blocka(lexcls, rt->excl_thread+1, rt->nthread-1);
1761 for (thread = 0; thread < rt->nthread; thread++)
1763 *excl_count += rt->excl_count_thread[thread];
1768 /* Some zones might not have exclusions, but some code still needs to
1769 * loop over the index, so we set the indices here.
1771 for (iz = nzone_excl; iz < zones->nizone; iz++)
1773 set_no_exclusions_zone(dd, zones, iz, lexcls);
1776 finish_local_exclusions(dd, zones, lexcls);
1779 fprintf(debug, "We have %d exclusions, check count %d\n",
1780 lexcls->nra, *excl_count);
1783 return nbonded_local;
1786 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
1788 lcgs->nr = dd->ncg_tot;
1789 lcgs->index = dd->cgindex;
1792 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1793 int npbcdim, matrix box,
1794 rvec cellsize_min, ivec npulse,
1798 gmx_mtop_t *mtop, gmx_localtop_t *ltop)
1800 gmx_bool bUniqueExcl, bRCheckMB, bRCheck2B, bRCheckExcl;
1804 t_pbc pbc, *pbc_null = NULL;
1808 fprintf(debug, "Making local topology\n");
1811 dd_make_local_cgs(dd, <op->cgs);
1815 bRCheckExcl = FALSE;
1817 if (dd->reverse_top->bMultiCGmols)
1819 /* We need to check to which cell bondeds should be assigned */
1820 rc = dd_cutoff_twobody(dd);
1823 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1826 /* Should we check cg_cm distances when assigning bonded interactions? */
1827 for (d = 0; d < DIM; d++)
1830 /* Only need to check for dimensions where the part of the box
1831 * that is not communicated is smaller than the cut-off.
1833 if (d < npbcdim && dd->nc[d] > 1 &&
1834 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
1841 /* Check for interactions between two atoms,
1842 * where we can allow interactions up to the cut-off,
1843 * instead of up to the smallest cell dimension.
1850 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
1851 d, cellsize_min[d], d, rcheck[d], bRCheck2B);
1854 if (dd->reverse_top->bExclRequired)
1856 bRCheckExcl = bRCheck2B;
1860 /* If we don't have forces on exclusions,
1861 * we don't care about exclusions being assigned mulitple times.
1863 bRCheckExcl = FALSE;
1865 if (bRCheckMB || bRCheck2B)
1870 set_pbc_dd(&pbc, fr->ePBC, dd, TRUE, box);
1881 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
1882 bRCheckMB, rcheck, bRCheck2B, rc,
1884 pbc_null, cgcm_or_x,
1886 <op->excls, &nexcl);
1888 /* The ilist is not sorted yet,
1889 * we can only do this when we have the charge arrays.
1891 ltop->idef.ilsort = ilsortUNKNOWN;
1893 if (dd->reverse_top->bExclRequired)
1895 dd->nbonded_local += nexcl;
1897 forcerec_set_excl_load(fr, ltop);
1900 ltop->atomtypes = mtop->atomtypes;
1902 /* For an error message only */
1903 dd->reverse_top->err_top_global = mtop;
1904 dd->reverse_top->err_top_local = ltop;
1907 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
1908 gmx_localtop_t *ltop)
1910 if (dd->reverse_top->ilsort == ilsortNO_FE)
1912 ltop->idef.ilsort = ilsortNO_FE;
1916 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
1920 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global)
1922 gmx_localtop_t *top;
1927 top->idef.ntypes = top_global->ffparams.ntypes;
1928 top->idef.atnr = top_global->ffparams.atnr;
1929 top->idef.functype = top_global->ffparams.functype;
1930 top->idef.iparams = top_global->ffparams.iparams;
1931 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
1932 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
1934 for (i = 0; i < F_NRE; i++)
1936 top->idef.il[i].iatoms = NULL;
1937 top->idef.il[i].nalloc = 0;
1939 top->idef.ilsort = ilsortUNKNOWN;
1944 void dd_init_local_state(gmx_domdec_t *dd,
1945 t_state *state_global, t_state *state_local)
1947 int buf[NITEM_DD_INIT_LOCAL_STATE];
1951 buf[0] = state_global->flags;
1952 buf[1] = state_global->ngtc;
1953 buf[2] = state_global->nnhpres;
1954 buf[3] = state_global->nhchainlength;
1955 buf[4] = state_global->dfhist.nlambda;
1957 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
1959 init_state(state_local, 0, buf[1], buf[2], buf[3], buf[4]);
1960 state_local->flags = buf[0];
1963 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
1969 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
1971 if (link->a[k] == cg_gl_j)
1978 /* Add this charge group link */
1979 if (link->index[cg_gl+1]+1 > link->nalloc_a)
1981 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
1982 srenew(link->a, link->nalloc_a);
1984 link->a[link->index[cg_gl+1]] = cg_gl_j;
1985 link->index[cg_gl+1]++;
1989 static int *make_at2cg(t_block *cgs)
1993 snew(at2cg, cgs->index[cgs->nr]);
1994 for (cg = 0; cg < cgs->nr; cg++)
1996 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2005 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
2006 cginfo_mb_t *cginfo_mb)
2008 gmx_reverse_top_t *rt;
2009 int mb, cg_offset, cg, cg_gl, a, aj, i, j, ftype, nral, nlink_mol, mol, ncgi;
2010 gmx_molblock_t *molb;
2011 gmx_moltype_t *molt;
2015 gmx_reverse_ilist_t ril;
2017 cginfo_mb_t *cgi_mb;
2019 /* For each charge group make a list of other charge groups
2020 * in the system that a linked to it via bonded interactions
2021 * which are also stored in reverse_top.
2024 rt = dd->reverse_top;
2027 snew(link->index, ncg_mtop(mtop)+1);
2034 for (mb = 0; mb < mtop->nmolblock; mb++)
2036 molb = &mtop->molblock[mb];
2037 if (molb->nmol == 0)
2041 molt = &mtop->moltype[molb->type];
2043 excls = &molt->excls;
2044 a2c = make_at2cg(cgs);
2045 /* Make a reverse ilist in which the interactions are linked
2046 * to all atoms, not only the first atom as in gmx_reverse_top.
2047 * The constraints are discarded here.
2049 make_reverse_ilist(molt, NULL, FALSE, FALSE, FALSE, TRUE, &ril);
2051 cgi_mb = &cginfo_mb[mb];
2053 for (cg = 0; cg < cgs->nr; cg++)
2055 cg_gl = cg_offset + cg;
2056 link->index[cg_gl+1] = link->index[cg_gl];
2057 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2060 while (i < ril.index[a+1])
2062 ftype = ril.il[i++];
2064 /* Skip the ifunc index */
2066 for (j = 0; j < nral; j++)
2071 check_link(link, cg_gl, cg_offset+a2c[aj]);
2074 i += nral_rt(ftype);
2076 if (rt->bExclRequired)
2078 /* Exclusions always go both ways */
2079 for (j = excls->index[a]; j < excls->index[a+1]; j++)
2084 check_link(link, cg_gl, cg_offset+a2c[aj]);
2089 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2091 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2095 nlink_mol = link->index[cg_offset+cgs->nr] - link->index[cg_offset];
2097 cg_offset += cgs->nr;
2099 destroy_reverse_ilist(&ril);
2104 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt->name, cgs->nr, nlink_mol);
2109 /* Copy the data for the rest of the molecules in this block */
2110 link->nalloc_a += (molb->nmol - 1)*nlink_mol;
2111 srenew(link->a, link->nalloc_a);
2112 for (mol = 1; mol < molb->nmol; mol++)
2114 for (cg = 0; cg < cgs->nr; cg++)
2116 cg_gl = cg_offset + cg;
2117 link->index[cg_gl+1] =
2118 link->index[cg_gl+1-cgs->nr] + nlink_mol;
2119 for (j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2121 link->a[j] = link->a[j-nlink_mol] + cgs->nr;
2123 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2124 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2126 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2130 cg_offset += cgs->nr;
2137 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2143 static void bonded_cg_distance_mol(gmx_moltype_t *molt, int *at2cg,
2144 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2145 real *r_2b, int *ft2b, int *a2_1, int *a2_2,
2146 real *r_mb, int *ftmb, int *am_1, int *am_2)
2148 int ftype, nral, i, j, ai, aj, cgi, cgj;
2151 real r2_2b, r2_mb, rij2;
2155 for (ftype = 0; ftype < F_NRE; ftype++)
2157 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2159 il = &molt->ilist[ftype];
2163 for (i = 0; i < il->nr; i += 1+nral)
2165 for (ai = 0; ai < nral; ai++)
2167 cgi = at2cg[il->iatoms[i+1+ai]];
2168 for (aj = 0; aj < nral; aj++)
2170 cgj = at2cg[il->iatoms[i+1+aj]];
2173 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2174 if (nral == 2 && rij2 > r2_2b)
2178 *a2_1 = il->iatoms[i+1+ai];
2179 *a2_2 = il->iatoms[i+1+aj];
2181 if (nral > 2 && rij2 > r2_mb)
2185 *am_1 = il->iatoms[i+1+ai];
2186 *am_2 = il->iatoms[i+1+aj];
2197 excls = &molt->excls;
2198 for (ai = 0; ai < excls->nr; ai++)
2201 for (j = excls->index[ai]; j < excls->index[ai+1]; j++)
2203 cgj = at2cg[excls->a[j]];
2206 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2216 *r_2b = sqrt(r2_2b);
2217 *r_mb = sqrt(r2_mb);
2220 static void get_cgcm_mol(gmx_moltype_t *molt, gmx_ffparams_t *ffparams,
2221 int ePBC, t_graph *graph, matrix box,
2223 rvec *x, rvec *xs, rvec *cg_cm)
2227 if (ePBC != epbcNONE)
2229 mk_mshift(NULL, graph, ePBC, box, x);
2231 shift_x(graph, box, x, xs);
2232 /* By doing an extra mk_mshift the molecules that are broken
2233 * because they were e.g. imported from another software
2234 * will be made whole again. Such are the healing powers
2237 mk_mshift(NULL, graph, ePBC, box, xs);
2241 /* We copy the coordinates so the original coordinates remain
2242 * unchanged, just to be 100% sure that we do not affect
2243 * binary reproducibility of simulations.
2245 n = molt->cgs.index[molt->cgs.nr];
2246 for (i = 0; i < n; i++)
2248 copy_rvec(x[i], xs[i]);
2254 construct_vsites(vsite, xs, 0.0, NULL,
2255 ffparams->iparams, molt->ilist,
2256 epbcNONE, TRUE, NULL, NULL);
2259 calc_cgcm(NULL, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2262 static int have_vsite_molt(gmx_moltype_t *molt)
2268 for (i = 0; i < F_NRE; i++)
2270 if ((interaction_function[i].flags & IF_VSITE) &&
2271 molt->ilist[i].nr > 0)
2280 void dd_bonded_cg_distance(FILE *fplog,
2282 t_inputrec *ir, rvec *x, matrix box,
2284 real *r_2b, real *r_mb)
2286 gmx_bool bExclRequired;
2287 int mb, cg_offset, at_offset, *at2cg, mol;
2290 gmx_molblock_t *molb;
2291 gmx_moltype_t *molt;
2293 real rmol_2b, rmol_mb;
2294 int ft2b = -1, a_2b_1 = -1, a_2b_2 = -1, ftmb = -1, a_mb_1 = -1, a_mb_2 = -1;
2295 int ftm2b = -1, amol_2b_1 = -1, amol_2b_2 = -1, ftmmb = -1, amol_mb_1 = -1, amol_mb_2 = -1;
2297 bExclRequired = IR_EXCL_FORCES(*ir);
2299 vsite = init_vsite(mtop, NULL, TRUE);
2305 for (mb = 0; mb < mtop->nmolblock; mb++)
2307 molb = &mtop->molblock[mb];
2308 molt = &mtop->moltype[molb->type];
2309 if (molt->cgs.nr == 1 || molb->nmol == 0)
2311 cg_offset += molb->nmol*molt->cgs.nr;
2312 at_offset += molb->nmol*molt->atoms.nr;
2316 if (ir->ePBC != epbcNONE)
2318 mk_graph_ilist(NULL, molt->ilist, 0, molt->atoms.nr, FALSE, FALSE,
2322 at2cg = make_at2cg(&molt->cgs);
2323 snew(xs, molt->atoms.nr);
2324 snew(cg_cm, molt->cgs.nr);
2325 for (mol = 0; mol < molb->nmol; mol++)
2327 get_cgcm_mol(molt, &mtop->ffparams, ir->ePBC, &graph, box,
2328 have_vsite_molt(molt) ? vsite : NULL,
2329 x+at_offset, xs, cg_cm);
2331 bonded_cg_distance_mol(molt, at2cg, bBCheck, bExclRequired, cg_cm,
2332 &rmol_2b, &ftm2b, &amol_2b_1, &amol_2b_2,
2333 &rmol_mb, &ftmmb, &amol_mb_1, &amol_mb_2);
2334 if (rmol_2b > *r_2b)
2338 a_2b_1 = at_offset + amol_2b_1;
2339 a_2b_2 = at_offset + amol_2b_2;
2341 if (rmol_mb > *r_mb)
2345 a_mb_1 = at_offset + amol_mb_1;
2346 a_mb_2 = at_offset + amol_mb_2;
2349 cg_offset += molt->cgs.nr;
2350 at_offset += molt->atoms.nr;
2355 if (ir->ePBC != epbcNONE)
2362 /* We should have a vsite free routine, but here we can simply free */
2365 if (fplog && (ft2b >= 0 || ftmb >= 0))
2368 "Initial maximum inter charge-group distances:\n");
2372 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2373 *r_2b, interaction_function[ft2b].longname,
2374 a_2b_1+1, a_2b_2+1);
2379 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2380 *r_mb, interaction_function[ftmb].longname,
2381 a_mb_1+1, a_mb_2+1);