1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
16 * Gnomes, ROck Monsters And Chili Sauce
27 #include "domdec_network.h"
32 #include "chargegroup.h"
33 #include "gmx_random.h"
35 #include "mtop_util.h"
38 #include "gmx_ga2la.h"
40 #include "gmx_omp_nthreads.h"
42 /* for dd_init_local_state */
43 #define NITEM_DD_INIT_LOCAL_STATE 5
46 int *index; /* Index for each atom into il */
47 int *il; /* ftype|type|a0|...|an|ftype|... */
48 } gmx_reverse_ilist_t;
57 typedef struct gmx_reverse_top {
58 gmx_bool bExclRequired; /* Do we require all exclusions to be assigned? */
59 gmx_bool bConstr; /* Are there constraints in this revserse top? */
60 gmx_bool bSettle; /* Are there settles in this revserse top? */
61 gmx_bool bBCheck; /* All bonded interactions have to be assigned? */
62 gmx_bool bMultiCGmols; /* Are the multi charge-group molecules? */
63 gmx_reverse_ilist_t *ril_mt; /* Reverse ilist for all moltypes */
65 int ilsort; /* The sorting state of bondeds for free energy */
66 gmx_molblock_ind_t *mbi;
69 /* Work data structures for multi-threading */
73 int **vsite_pbc_nalloc;
75 t_blocka *excl_thread;
76 int *excl_count_thread;
78 /* Pointers only used for an error message */
79 gmx_mtop_t *err_top_global;
80 gmx_localtop_t *err_top_local;
83 static int nral_rt(int ftype)
85 /* Returns the number of atom entries for il in gmx_reverse_top_t */
89 if (interaction_function[ftype].flags & IF_VSITE)
91 /* With vsites the reverse topology contains
92 * two extra entries for PBC.
100 /* This function tells which interactions need to be assigned exactly once */
101 static gmx_bool dd_check_ftype(int ftype, gmx_bool bBCheck,
102 gmx_bool bConstr, gmx_bool bSettle)
104 return (((interaction_function[ftype].flags & IF_BOND) &&
105 !(interaction_function[ftype].flags & IF_VSITE) &&
106 (bBCheck || !(interaction_function[ftype].flags & IF_LIMZERO))) ||
107 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
108 (bSettle && ftype == F_SETTLE));
111 static void print_error_header(FILE *fplog, char *moltypename, int nprint)
113 fprintf(fplog, "\nMolecule type '%s'\n", moltypename);
114 fprintf(stderr, "\nMolecule type '%s'\n", moltypename);
116 "the first %d missing interactions, except for exclusions:\n",
119 "the first %d missing interactions, except for exclusions:\n",
123 static void print_missing_interactions_mb(FILE *fplog, t_commrec *cr,
124 gmx_reverse_top_t *rt,
126 gmx_reverse_ilist_t *ril,
127 int a_start, int a_end,
128 int nat_mol, int nmol,
131 int nril_mol, *assigned, *gatindex;
132 int ftype, ftype_j, nral, i, j_mol, j, k, a0, a0_mol, mol, a, a_gl;
138 nril_mol = ril->index[nat_mol];
139 snew(assigned, nmol*nril_mol);
141 gatindex = cr->dd->gatindex;
142 for (ftype = 0; ftype < F_NRE; ftype++)
144 if (dd_check_ftype(ftype, rt->bBCheck, rt->bConstr, rt->bSettle))
147 il = &idef->il[ftype];
149 for (i = 0; i < il->nr; i += 1+nral)
151 a0 = gatindex[ia[1]];
152 /* Check if this interaction is in
153 * the currently checked molblock.
155 if (a0 >= a_start && a0 < a_end)
157 mol = (a0 - a_start)/nat_mol;
158 a0_mol = (a0 - a_start) - mol*nat_mol;
159 j_mol = ril->index[a0_mol];
161 while (j_mol < ril->index[a0_mol+1] && !bFound)
163 j = mol*nril_mol + j_mol;
164 ftype_j = ril->il[j_mol];
165 /* Here we need to check if this interaction has
166 * not already been assigned, since we could have
167 * multiply defined interactions.
169 if (ftype == ftype_j && ia[0] == ril->il[j_mol+1] &&
172 /* Check the atoms */
174 for (a = 0; a < nral; a++)
176 if (gatindex[ia[1+a]] !=
177 a_start + mol*nat_mol + ril->il[j_mol+2+a])
187 j_mol += 2 + nral_rt(ftype_j);
191 gmx_incons("Some interactions seem to be assigned multiple times");
199 gmx_sumi(nmol*nril_mol, assigned, cr);
203 for (mol = 0; mol < nmol; mol++)
206 while (j_mol < nril_mol)
208 ftype = ril->il[j_mol];
210 j = mol*nril_mol + j_mol;
211 if (assigned[j] == 0 &&
212 !(interaction_function[ftype].flags & IF_VSITE))
214 if (DDMASTER(cr->dd))
218 print_error_header(fplog, moltypename, nprint);
220 fprintf(fplog, "%20s atoms",
221 interaction_function[ftype].longname);
222 fprintf(stderr, "%20s atoms",
223 interaction_function[ftype].longname);
224 for (a = 0; a < nral; a++)
226 fprintf(fplog, "%5d", ril->il[j_mol+2+a]+1);
227 fprintf(stderr, "%5d", ril->il[j_mol+2+a]+1);
232 fprintf(stderr, " ");
235 fprintf(fplog, " global");
236 fprintf(stderr, " global");
237 for (a = 0; a < nral; a++)
239 fprintf(fplog, "%6d",
240 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
241 fprintf(stderr, "%6d",
242 a_start+mol*nat_mol+ril->il[j_mol+2+a]+1);
244 fprintf(fplog, "\n");
245 fprintf(stderr, "\n");
253 j_mol += 2 + nral_rt(ftype);
260 static void print_missing_interactions_atoms(FILE *fplog, t_commrec *cr,
261 gmx_mtop_t *mtop, t_idef *idef)
263 int mb, a_start, a_end;
264 gmx_molblock_t *molb;
265 gmx_reverse_top_t *rt;
267 rt = cr->dd->reverse_top;
269 /* Print the atoms in the missing interactions per molblock */
271 for (mb = 0; mb < mtop->nmolblock; mb++)
273 molb = &mtop->molblock[mb];
275 a_end = a_start + molb->nmol*molb->natoms_mol;
277 print_missing_interactions_mb(fplog, cr, rt,
278 *(mtop->moltype[molb->type].name),
279 &rt->ril_mt[molb->type],
280 a_start, a_end, molb->natoms_mol,
286 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr, int local_count, gmx_mtop_t *top_global, t_state *state_local)
288 int ndiff_tot, cl[F_NRE], n, ndiff, rest_global, rest_local;
292 gmx_mtop_t *err_top_global;
293 gmx_localtop_t *err_top_local;
297 err_top_global = dd->reverse_top->err_top_global;
298 err_top_local = dd->reverse_top->err_top_local;
302 fprintf(fplog, "\nNot all bonded interactions have been properly assigned to the domain decomposition cells\n");
306 ndiff_tot = local_count - dd->nbonded_global;
308 for (ftype = 0; ftype < F_NRE; ftype++)
311 cl[ftype] = err_top_local->idef.il[ftype].nr/(1+nral);
314 gmx_sumi(F_NRE, cl, cr);
318 fprintf(fplog, "\nA list of missing interactions:\n");
319 fprintf(stderr, "\nA list of missing interactions:\n");
320 rest_global = dd->nbonded_global;
321 rest_local = local_count;
322 for (ftype = 0; ftype < F_NRE; ftype++)
324 /* In the reverse and local top all constraints are merged
325 * into F_CONSTR. So in the if statement we skip F_CONSTRNC
326 * and add these constraints when doing F_CONSTR.
328 if (((interaction_function[ftype].flags & IF_BOND) &&
329 (dd->reverse_top->bBCheck
330 || !(interaction_function[ftype].flags & IF_LIMZERO)))
331 || (dd->reverse_top->bConstr && ftype == F_CONSTR)
332 || (dd->reverse_top->bSettle && ftype == F_SETTLE))
335 n = gmx_mtop_ftype_count(err_top_global, ftype);
336 if (ftype == F_CONSTR)
338 n += gmx_mtop_ftype_count(err_top_global, F_CONSTRNC);
340 ndiff = cl[ftype] - n;
343 sprintf(buf, "%20s of %6d missing %6d",
344 interaction_function[ftype].longname, n, -ndiff);
345 fprintf(fplog, "%s\n", buf);
346 fprintf(stderr, "%s\n", buf);
349 rest_local -= cl[ftype];
353 ndiff = rest_local - rest_global;
356 sprintf(buf, "%20s of %6d missing %6d", "exclusions",
357 rest_global, -ndiff);
358 fprintf(fplog, "%s\n", buf);
359 fprintf(stderr, "%s\n", buf);
363 print_missing_interactions_atoms(fplog, cr, err_top_global,
364 &err_top_local->idef);
365 write_dd_pdb("dd_dump_err", 0, "dump", top_global, cr,
366 -1, state_local->x, state_local->box);
371 gmx_incons("One or more interactions were multiple assigned in the domain decompostion");
375 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));
380 static void global_atomnr_to_moltype_ind(gmx_reverse_top_t *rt, int i_gl,
381 int *mb, int *mt, int *mol, int *i_mol)
386 gmx_molblock_ind_t *mbi = rt->mbi;
388 int end = rt->nmolblock; /* exclusive */
391 /* binary search for molblock_ind */
395 if (i_gl >= mbi[mid].a_end)
399 else if (i_gl < mbi[mid].a_start)
413 *mol = (i_gl - mbi->a_start) / mbi->natoms_mol;
414 *i_mol = (i_gl - mbi->a_start) - (*mol)*mbi->natoms_mol;
417 static int count_excls(t_block *cgs, t_blocka *excls, int *n_intercg_excl)
419 int n, n_inter, cg, at0, at1, at, excl, atj;
423 for (cg = 0; cg < cgs->nr; cg++)
425 at0 = cgs->index[cg];
426 at1 = cgs->index[cg+1];
427 for (at = at0; at < at1; at++)
429 for (excl = excls->index[at]; excl < excls->index[at+1]; excl++)
431 atj = excls->a[excl];
435 if (atj < at0 || atj >= at1)
447 static int low_make_reverse_ilist(t_ilist *il_mt, t_atom *atom,
450 gmx_bool bConstr, gmx_bool bSettle,
452 int *r_index, int *r_il,
453 gmx_bool bLinkToAllAtoms,
456 int ftype, nral, i, j, nlink, link;
464 for (ftype = 0; ftype < F_NRE; ftype++)
466 if ((interaction_function[ftype].flags & (IF_BOND | IF_VSITE)) ||
467 (bConstr && (ftype == F_CONSTR || ftype == F_CONSTRNC)) ||
468 (bSettle && ftype == F_SETTLE))
470 bVSite = (interaction_function[ftype].flags & IF_VSITE);
474 for (i = 0; i < il->nr; i += 1+nral)
481 /* We don't need the virtual sites for the cg-links */
491 /* Couple to the first atom in the interaction */
494 for (link = 0; link < nlink; link++)
499 r_il[r_index[a]+count[a]] =
500 (ftype == F_CONSTRNC ? F_CONSTR : ftype);
501 r_il[r_index[a]+count[a]+1] = ia[0];
502 for (j = 1; j < 1+nral; j++)
504 /* Store the molecular atom number */
505 r_il[r_index[a]+count[a]+1+j] = ia[j];
508 if (interaction_function[ftype].flags & IF_VSITE)
512 /* Add an entry to iatoms for storing
513 * which of the constructing atoms are
516 r_il[r_index[a]+count[a]+2+nral] = 0;
517 for (j = 2; j < 1+nral; j++)
519 if (atom[ia[j]].ptype == eptVSite)
521 r_il[r_index[a]+count[a]+2+nral] |= (2<<j);
524 /* Store vsite pbc atom in a second extra entry */
525 r_il[r_index[a]+count[a]+2+nral+1] =
526 (vsite_pbc ? vsite_pbc[ftype-F_VSITE2][i/(1+nral)] : -2);
531 /* We do not count vsites since they are always
532 * uniquely assigned and can be assigned
533 * to multiple nodes with recursive vsites.
536 !(interaction_function[ftype].flags & IF_LIMZERO))
541 count[a] += 2 + nral_rt(ftype);
550 static int make_reverse_ilist(gmx_moltype_t *molt,
552 gmx_bool bConstr, gmx_bool bSettle,
554 gmx_bool bLinkToAllAtoms,
555 gmx_reverse_ilist_t *ril_mt)
557 int nat_mt, *count, i, nint_mt;
559 /* Count the interactions */
560 nat_mt = molt->atoms.nr;
562 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
564 bConstr, bSettle, bBCheck, NULL, NULL,
565 bLinkToAllAtoms, FALSE);
567 snew(ril_mt->index, nat_mt+1);
568 ril_mt->index[0] = 0;
569 for (i = 0; i < nat_mt; i++)
571 ril_mt->index[i+1] = ril_mt->index[i] + count[i];
574 snew(ril_mt->il, ril_mt->index[nat_mt]);
576 /* Store the interactions */
578 low_make_reverse_ilist(molt->ilist, molt->atoms.atom, vsite_pbc,
580 bConstr, bSettle, bBCheck,
581 ril_mt->index, ril_mt->il,
582 bLinkToAllAtoms, TRUE);
589 static void destroy_reverse_ilist(gmx_reverse_ilist_t *ril)
595 static gmx_reverse_top_t *make_reverse_top(gmx_mtop_t *mtop, gmx_bool bFE,
596 int ***vsite_pbc_molt,
597 gmx_bool bConstr, gmx_bool bSettle,
598 gmx_bool bBCheck, int *nint)
601 gmx_reverse_top_t *rt;
608 /* Should we include constraints (for SHAKE) in rt? */
609 rt->bConstr = bConstr;
610 rt->bSettle = bSettle;
611 rt->bBCheck = bBCheck;
613 rt->bMultiCGmols = FALSE;
614 snew(nint_mt, mtop->nmoltype);
615 snew(rt->ril_mt, mtop->nmoltype);
616 rt->ril_mt_tot_size = 0;
617 for (mt = 0; mt < mtop->nmoltype; mt++)
619 molt = &mtop->moltype[mt];
620 if (molt->cgs.nr > 1)
622 rt->bMultiCGmols = TRUE;
625 /* Make the atom to interaction list for this molecule type */
627 make_reverse_ilist(molt, vsite_pbc_molt ? vsite_pbc_molt[mt] : NULL,
628 rt->bConstr, rt->bSettle, rt->bBCheck, FALSE,
631 rt->ril_mt_tot_size += rt->ril_mt[mt].index[molt->atoms.nr];
635 fprintf(debug, "The total size of the atom to interaction index is %d integers\n", rt->ril_mt_tot_size);
639 for (mb = 0; mb < mtop->nmolblock; mb++)
641 *nint += mtop->molblock[mb].nmol*nint_mt[mtop->molblock[mb].type];
645 if (bFE && gmx_mtop_bondeds_free_energy(mtop))
647 rt->ilsort = ilsortFE_UNSORTED;
651 rt->ilsort = ilsortNO_FE;
654 /* Make a molblock index for fast searching */
655 snew(rt->mbi, mtop->nmolblock);
656 rt->nmolblock = mtop->nmolblock;
658 for (mb = 0; mb < mtop->nmolblock; mb++)
660 rt->mbi[mb].a_start = i;
661 i += mtop->molblock[mb].nmol*mtop->molblock[mb].natoms_mol;
662 rt->mbi[mb].a_end = i;
663 rt->mbi[mb].natoms_mol = mtop->molblock[mb].natoms_mol;
664 rt->mbi[mb].type = mtop->molblock[mb].type;
667 rt->nthread = gmx_omp_nthreads_get(emntDomdec);
668 snew(rt->idef_thread, rt->nthread);
669 if (vsite_pbc_molt != NULL)
671 snew(rt->vsite_pbc, rt->nthread);
672 snew(rt->vsite_pbc_nalloc, rt->nthread);
673 for (thread = 0; thread < rt->nthread; thread++)
675 snew(rt->vsite_pbc[thread], F_VSITEN-F_VSITE2+1);
676 snew(rt->vsite_pbc_nalloc[thread], F_VSITEN-F_VSITE2+1);
679 snew(rt->nbonded_thread, rt->nthread);
680 snew(rt->excl_thread, rt->nthread);
681 snew(rt->excl_count_thread, rt->nthread);
686 void dd_make_reverse_top(FILE *fplog,
687 gmx_domdec_t *dd, gmx_mtop_t *mtop,
688 gmx_vsite_t *vsite, gmx_constr_t constr,
689 t_inputrec *ir, gmx_bool bBCheck)
691 int mb, n_recursive_vsite, nexcl, nexcl_icg, a;
692 gmx_molblock_t *molb;
697 fprintf(fplog, "\nLinking all bonded interactions to atoms\n");
700 /* If normal and/or settle constraints act only within charge groups,
701 * we can store them in the reverse top and simply assign them to domains.
702 * Otherwise we need to assign them to multiple domains and set up
703 * the parallel version constraint algoirthm(s).
706 dd->reverse_top = make_reverse_top(mtop, ir->efep != efepNO,
707 vsite ? vsite->vsite_pbc_molt : NULL,
708 !dd->bInterCGcons, !dd->bInterCGsettles,
709 bBCheck, &dd->nbonded_global);
711 if (dd->reverse_top->ril_mt_tot_size >= 200000 &&
713 mtop->nmolblock == 1 && mtop->molblock[0].nmol == 1)
715 /* mtop comes from a pre Gromacs 4 tpr file */
716 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";
719 fprintf(fplog, "\n%s\n\n", note);
723 fprintf(stderr, "\n%s\n\n", note);
727 dd->reverse_top->bExclRequired = IR_EXCL_FORCES(*ir);
730 dd->n_intercg_excl = 0;
731 for (mb = 0; mb < mtop->nmolblock; mb++)
733 molb = &mtop->molblock[mb];
734 molt = &mtop->moltype[molb->type];
735 nexcl += molb->nmol*count_excls(&molt->cgs, &molt->excls, &nexcl_icg);
736 dd->n_intercg_excl += molb->nmol*nexcl_icg;
738 if (dd->reverse_top->bExclRequired)
740 dd->nbonded_global += nexcl;
741 if (EEL_FULL(ir->coulombtype) && dd->n_intercg_excl > 0 && fplog)
743 fprintf(fplog, "There are %d inter charge-group exclusions,\n"
744 "will use an extra communication step for exclusion forces for %s\n",
745 dd->n_intercg_excl, eel_names[ir->coulombtype]);
749 if (vsite && vsite->n_intercg_vsite > 0)
753 fprintf(fplog, "There are %d inter charge-group virtual sites,\n"
754 "will an extra communication step for selected coordinates and forces\n",
755 vsite->n_intercg_vsite);
757 init_domdec_vsites(dd, vsite->n_intercg_vsite);
760 if (dd->bInterCGcons || dd->bInterCGsettles)
762 init_domdec_constraints(dd, mtop, constr);
766 fprintf(fplog, "\n");
770 static inline void add_ifunc(int nral, t_iatom *tiatoms, t_ilist *il)
775 if (il->nr+1+nral > il->nalloc)
777 il->nalloc = over_alloc_large(il->nr+1+nral);
778 srenew(il->iatoms, il->nalloc);
780 liatoms = il->iatoms + il->nr;
781 for (k = 0; k <= nral; k++)
783 liatoms[k] = tiatoms[k];
788 static void add_posres(int mol, int a_mol, const gmx_molblock_t *molb,
789 t_iatom *iatoms, const t_iparams *ip_in,
795 /* This position restraint has not been added yet,
796 * so it's index is the current number of position restraints.
798 n = idef->il[F_POSRES].nr/2;
799 if (n+1 > idef->iparams_posres_nalloc)
801 idef->iparams_posres_nalloc = over_alloc_dd(n+1);
802 srenew(idef->iparams_posres, idef->iparams_posres_nalloc);
804 ip = &idef->iparams_posres[n];
805 /* Copy the force constants */
806 *ip = ip_in[iatoms[0]];
808 /* Get the position restraint coordinates from the molblock */
809 a_molb = mol*molb->natoms_mol + a_mol;
810 if (a_molb >= molb->nposres_xA)
812 gmx_incons("Not enough position restraint coordinates");
814 ip->posres.pos0A[XX] = molb->posres_xA[a_molb][XX];
815 ip->posres.pos0A[YY] = molb->posres_xA[a_molb][YY];
816 ip->posres.pos0A[ZZ] = molb->posres_xA[a_molb][ZZ];
817 if (molb->nposres_xB > 0)
819 ip->posres.pos0B[XX] = molb->posres_xB[a_molb][XX];
820 ip->posres.pos0B[YY] = molb->posres_xB[a_molb][YY];
821 ip->posres.pos0B[ZZ] = molb->posres_xB[a_molb][ZZ];
825 ip->posres.pos0B[XX] = ip->posres.pos0A[XX];
826 ip->posres.pos0B[YY] = ip->posres.pos0A[YY];
827 ip->posres.pos0B[ZZ] = ip->posres.pos0A[ZZ];
829 /* Set the parameter index for idef->iparams_posre */
833 static void add_fbposres(int mol, int a_mol, const gmx_molblock_t *molb,
834 t_iatom *iatoms, const t_iparams *ip_in,
840 /* This flat-bottom position restraint has not been added yet,
841 * so it's index is the current number of position restraints.
843 n = idef->il[F_FBPOSRES].nr/2;
844 if (n+1 > idef->iparams_fbposres_nalloc)
846 idef->iparams_fbposres_nalloc = over_alloc_dd(n+1);
847 srenew(idef->iparams_fbposres, idef->iparams_fbposres_nalloc);
849 ip = &idef->iparams_fbposres[n];
850 /* Copy the force constants */
851 *ip = ip_in[iatoms[0]];
853 /* Get the position restriant coordinats from the molblock */
854 a_molb = mol*molb->natoms_mol + a_mol;
855 if (a_molb >= molb->nposres_xA)
857 gmx_incons("Not enough position restraint coordinates");
859 /* Take reference positions from A position of normal posres */
860 ip->fbposres.pos0[XX] = molb->posres_xA[a_molb][XX];
861 ip->fbposres.pos0[YY] = molb->posres_xA[a_molb][YY];
862 ip->fbposres.pos0[ZZ] = molb->posres_xA[a_molb][ZZ];
864 /* Note: no B-type for flat-bottom posres */
866 /* Set the parameter index for idef->iparams_posre */
870 static void add_vsite(gmx_ga2la_t ga2la, int *index, int *rtil,
872 gmx_bool bHomeA, int a, int a_gl, int a_mol,
874 t_idef *idef, int **vsite_pbc, int *vsite_pbc_nalloc)
876 int k, ak_gl, vsi, pbc_a_mol;
877 t_iatom tiatoms[1+MAXATOMLIST], *iatoms_r;
878 int j, ftype_r, nral_r;
881 tiatoms[0] = iatoms[0];
885 /* We know the local index of the first atom */
890 /* Convert later in make_local_vsites */
891 tiatoms[1] = -a_gl - 1;
894 for (k = 2; k < 1+nral; k++)
896 ak_gl = a_gl + iatoms[k] - a_mol;
897 if (!ga2la_get_home(ga2la, ak_gl, &tiatoms[k]))
899 /* Copy the global index, convert later in make_local_vsites */
900 tiatoms[k] = -(ak_gl + 1);
904 /* Add this interaction to the local topology */
905 add_ifunc(nral, tiatoms, &idef->il[ftype]);
908 vsi = idef->il[ftype].nr/(1+nral) - 1;
909 if (vsi >= vsite_pbc_nalloc[ftype-F_VSITE2])
911 vsite_pbc_nalloc[ftype-F_VSITE2] = over_alloc_large(vsi+1);
912 srenew(vsite_pbc[ftype-F_VSITE2], vsite_pbc_nalloc[ftype-F_VSITE2]);
916 pbc_a_mol = iatoms[1+nral+1];
919 /* The pbc flag is one of the following two options:
920 * -2: vsite and all constructing atoms are within the same cg, no pbc
921 * -1: vsite and its first constructing atom are in the same cg, do pbc
923 vsite_pbc[ftype-F_VSITE2][vsi] = pbc_a_mol;
927 /* Set the pbc atom for this vsite so we can make its pbc
928 * identical to the rest of the atoms in its charge group.
929 * Since the order of the atoms does not change within a charge
930 * group, we do not need the global to local atom index.
932 vsite_pbc[ftype-F_VSITE2][vsi] = a + pbc_a_mol - iatoms[1];
937 /* This vsite is non-home (required for recursion),
938 * and therefore there is no charge group to match pbc with.
939 * But we always turn on full_pbc to assure that higher order
940 * recursion works correctly.
942 vsite_pbc[ftype-F_VSITE2][vsi] = -1;
948 /* Check for recursion */
949 for (k = 2; k < 1+nral; k++)
951 if ((iatoms[1+nral] & (2<<k)) && (tiatoms[k] < 0))
953 /* This construction atoms is a vsite and not a home atom */
956 fprintf(debug, "Constructing atom %d of vsite atom %d is a vsite and non-home\n", iatoms[k]+1, a_mol+1);
958 /* Find the vsite construction */
960 /* Check all interactions assigned to this atom */
961 j = index[iatoms[k]];
962 while (j < index[iatoms[k]+1])
965 nral_r = NRAL(ftype_r);
966 if (interaction_function[ftype_r].flags & IF_VSITE)
968 /* Add this vsite (recursion) */
969 add_vsite(ga2la, index, rtil, ftype_r, nral_r,
970 FALSE, -1, a_gl+iatoms[k]-iatoms[1], iatoms[k],
971 rtil+j, idef, vsite_pbc, vsite_pbc_nalloc);
984 static void make_la2lc(gmx_domdec_t *dd)
986 int *cgindex, *la2lc, cg, a;
988 cgindex = dd->cgindex;
990 if (dd->nat_tot > dd->la2lc_nalloc)
992 dd->la2lc_nalloc = over_alloc_dd(dd->nat_tot);
993 snew(dd->la2lc, dd->la2lc_nalloc);
997 /* Make the local atom to local cg index */
998 for (cg = 0; cg < dd->ncg_tot; cg++)
1000 for (a = cgindex[cg]; a < cgindex[cg+1]; a++)
1007 static real dd_dist2(t_pbc *pbc_null, rvec *cg_cm, const int *la2lc, int i, int j)
1013 pbc_dx_aiuc(pbc_null, cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1017 rvec_sub(cg_cm[la2lc[i]], cg_cm[la2lc[j]], dx);
1023 /* Append the nsrc t_blocka block structures in src to *dest */
1024 static void combine_blocka(t_blocka *dest, const t_blocka *src, int nsrc)
1028 ni = src[nsrc-1].nr;
1030 for (s = 0; s < nsrc; s++)
1034 if (ni + 1 > dest->nalloc_index)
1036 dest->nalloc_index = over_alloc_large(ni+1);
1037 srenew(dest->index, dest->nalloc_index);
1039 if (dest->nra + na > dest->nalloc_a)
1041 dest->nalloc_a = over_alloc_large(dest->nra+na);
1042 srenew(dest->a, dest->nalloc_a);
1044 for (s = 0; s < nsrc; s++)
1046 for (i = dest->nr+1; i < src[s].nr+1; i++)
1048 dest->index[i] = dest->nra + src[s].index[i];
1050 for (i = 0; i < src[s].nra; i++)
1052 dest->a[dest->nra+i] = src[s].a[i];
1054 dest->nr = src[s].nr;
1055 dest->nra += src[s].nra;
1059 /* Append the nsrc t_idef structures in src to *dest,
1060 * virtual sites need special attention, as pbc info differs per vsite.
1062 static void combine_idef(t_idef *dest, const t_idef *src, int nsrc,
1063 gmx_vsite_t *vsite, int ***vsite_pbc_t)
1069 int nral1 = 0, ftv = 0;
1071 for (ftype = 0; ftype < F_NRE; ftype++)
1074 for (s = 0; s < nsrc; s++)
1076 n += src[s].il[ftype].nr;
1080 ild = &dest->il[ftype];
1082 if (ild->nr + n > ild->nalloc)
1084 ild->nalloc = over_alloc_large(ild->nr+n);
1085 srenew(ild->iatoms, ild->nalloc);
1088 vpbc = ((interaction_function[ftype].flags & IF_VSITE) &&
1089 vsite->vsite_pbc_loc != NULL);
1092 nral1 = 1 + NRAL(ftype);
1093 ftv = ftype - F_VSITE2;
1094 if ((ild->nr + n)/nral1 > vsite->vsite_pbc_loc_nalloc[ftv])
1096 vsite->vsite_pbc_loc_nalloc[ftv] =
1097 over_alloc_large((ild->nr + n)/nral1);
1098 srenew(vsite->vsite_pbc_loc[ftv],
1099 vsite->vsite_pbc_loc_nalloc[ftv]);
1103 for (s = 0; s < nsrc; s++)
1105 ils = &src[s].il[ftype];
1106 for (i = 0; i < ils->nr; i++)
1108 ild->iatoms[ild->nr+i] = ils->iatoms[i];
1112 for (i = 0; i < ils->nr; i += nral1)
1114 vsite->vsite_pbc_loc[ftv][(ild->nr+i)/nral1] =
1115 vsite_pbc_t[s][ftv][i/nral1];
1124 /* Position restraints need an additional treatment */
1125 if (dest->il[F_POSRES].nr > 0)
1127 n = dest->il[F_POSRES].nr/2;
1128 if (n > dest->iparams_posres_nalloc)
1130 dest->iparams_posres_nalloc = over_alloc_large(n);
1131 srenew(dest->iparams_posres, dest->iparams_posres_nalloc);
1133 /* Set n to the number of original position restraints in dest */
1134 for (s = 0; s < nsrc; s++)
1136 n -= src[s].il[F_POSRES].nr/2;
1138 for (s = 0; s < nsrc; s++)
1140 for (i = 0; i < src[s].il[F_POSRES].nr/2; i++)
1142 /* Correct the index into iparams_posres */
1143 dest->il[F_POSRES].iatoms[n*2] = n;
1144 /* Copy the position restraint force parameters */
1145 dest->iparams_posres[n] = src[s].iparams_posres[i];
1152 /* This function looks up and assigns bonded interactions for zone iz.
1153 * With thread parallelizing each thread acts on a different atom range:
1154 * at_start to at_end.
1156 static int make_bondeds_zone(gmx_domdec_t *dd,
1157 const gmx_domdec_zones_t *zones,
1158 const gmx_molblock_t *molb,
1159 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1161 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1162 const t_iparams *ip_in,
1163 t_idef *idef, gmx_vsite_t *vsite,
1165 int *vsite_pbc_nalloc,
1167 int at_start, int at_end)
1169 int i, i_gl, mb, mt, mol, i_mol, j, ftype, nral, d, k;
1171 t_iatom *iatoms, tiatoms[1+MAXATOMLIST];
1172 gmx_bool bBCheck, bUse, bLocal;
1173 ivec k_zero, k_plus;
1178 const gmx_domdec_ns_ranges_t *izone;
1179 gmx_reverse_top_t *rt;
1182 nizone = zones->nizone;
1183 izone = zones->izone;
1185 rt = dd->reverse_top;
1187 bBCheck = rt->bBCheck;
1193 for (i = at_start; i < at_end; i++)
1195 /* Get the global atom number */
1196 i_gl = dd->gatindex[i];
1197 global_atomnr_to_moltype_ind(rt, i_gl, &mb, &mt, &mol, &i_mol);
1198 /* Check all interactions assigned to this atom */
1199 index = rt->ril_mt[mt].index;
1200 rtil = rt->ril_mt[mt].il;
1202 while (j < index[i_mol+1])
1207 if (ftype == F_SETTLE)
1209 /* Settles are only in the reverse top when they
1210 * operate within a charge group. So we can assign
1211 * them without checks. We do this only for performance
1212 * reasons; it could be handled by the code below.
1216 /* Home zone: add this settle to the local topology */
1217 tiatoms[0] = iatoms[0];
1219 tiatoms[2] = i + iatoms[2] - iatoms[1];
1220 tiatoms[3] = i + iatoms[3] - iatoms[1];
1221 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1226 else if (interaction_function[ftype].flags & IF_VSITE)
1228 /* The vsite construction goes where the vsite itself is */
1231 add_vsite(dd->ga2la, index, rtil, ftype, nral,
1232 TRUE, i, i_gl, i_mol,
1233 iatoms, idef, vsite_pbc, vsite_pbc_nalloc);
1240 tiatoms[0] = iatoms[0];
1244 /* Assign single-body interactions to the home zone */
1249 if (ftype == F_POSRES)
1251 add_posres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1254 else if (ftype == F_FBPOSRES)
1256 add_fbposres(mol, i_mol, &molb[mb], tiatoms, ip_in,
1267 /* This is a two-body interaction, we can assign
1268 * analogous to the non-bonded assignments.
1270 if (!ga2la_get(ga2la, i_gl+iatoms[2]-i_mol, &a_loc, &kz))
1280 /* Check zone interaction assignments */
1281 bUse = ((iz < nizone && iz <= kz &&
1282 izone[iz].j0 <= kz && kz < izone[iz].j1) ||
1283 (kz < nizone &&iz > kz &&
1284 izone[kz].j0 <= iz && iz < izone[kz].j1));
1289 /* If necessary check the cgcm distance */
1291 dd_dist2(pbc_null, cg_cm, la2lc,
1292 tiatoms[1], tiatoms[2]) >= rc2)
1301 /* Assign this multi-body bonded interaction to
1302 * the local node if we have all the atoms involved
1303 * (local or communicated) and the minimum zone shift
1304 * in each dimension is zero, for dimensions
1305 * with 2 DD cells an extra check may be necessary.
1310 for (k = 1; k <= nral && bUse; k++)
1312 bLocal = ga2la_get(ga2la, i_gl+iatoms[k]-i_mol,
1314 if (!bLocal || kz >= zones->n)
1316 /* We do not have this atom of this interaction
1317 * locally, or it comes from more than one cell
1325 for (d = 0; d < DIM; d++)
1327 if (zones->shift[kz][d] == 0)
1339 k_zero[XX] && k_zero[YY] && k_zero[ZZ]);
1342 for (d = 0; (d < DIM && bUse); d++)
1344 /* Check if the cg_cm distance falls within
1345 * the cut-off to avoid possible multiple
1346 * assignments of bonded interactions.
1350 dd_dist2(pbc_null, cg_cm, la2lc,
1351 tiatoms[k_zero[d]], tiatoms[k_plus[d]]) >= rc2)
1360 /* Add this interaction to the local topology */
1361 add_ifunc(nral, tiatoms, &idef->il[ftype]);
1362 /* Sum so we can check in global_stat
1363 * if we have everything.
1366 !(interaction_function[ftype].flags & IF_LIMZERO))
1376 return nbonded_local;
1379 static void set_no_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1380 int iz, t_blocka *lexcls)
1384 a0 = dd->cgindex[zones->cg_range[iz]];
1385 a1 = dd->cgindex[zones->cg_range[iz+1]];
1387 for (a = a0+1; a < a1+1; a++)
1389 lexcls->index[a] = lexcls->nra;
1393 static int make_exclusions_zone(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1394 const gmx_moltype_t *moltype,
1395 gmx_bool bRCheck, real rc2,
1396 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1400 int cg_start, int cg_end)
1402 int nizone, n, count, jla0, jla1, jla;
1403 int cg, la0, la1, la, a_gl, mb, mt, mol, a_mol, j, aj_mol;
1404 const t_blocka *excls;
1411 jla0 = dd->cgindex[zones->izone[iz].jcg0];
1412 jla1 = dd->cgindex[zones->izone[iz].jcg1];
1414 /* We set the end index, but note that we might not start at zero here */
1415 lexcls->nr = dd->cgindex[cg_end];
1419 for (cg = cg_start; cg < cg_end; cg++)
1421 /* Here we assume the number of exclusions in one charge group
1422 * is never larger than 1000.
1424 if (n+1000 > lexcls->nalloc_a)
1426 lexcls->nalloc_a = over_alloc_large(n+1000);
1427 srenew(lexcls->a, lexcls->nalloc_a);
1429 la0 = dd->cgindex[cg];
1430 la1 = dd->cgindex[cg+1];
1431 if (GET_CGINFO_EXCL_INTER(cginfo[cg]) ||
1432 !GET_CGINFO_EXCL_INTRA(cginfo[cg]))
1434 /* Copy the exclusions from the global top */
1435 for (la = la0; la < la1; la++)
1437 lexcls->index[la] = n;
1438 a_gl = dd->gatindex[la];
1439 global_atomnr_to_moltype_ind(dd->reverse_top, a_gl, &mb, &mt, &mol, &a_mol);
1440 excls = &moltype[mt].excls;
1441 for (j = excls->index[a_mol]; j < excls->index[a_mol+1]; j++)
1443 aj_mol = excls->a[j];
1444 /* This computation of jla is only correct intra-cg */
1445 jla = la + aj_mol - a_mol;
1446 if (jla >= la0 && jla < la1)
1448 /* This is an intra-cg exclusion. We can skip
1449 * the global indexing and distance checking.
1451 /* Intra-cg exclusions are only required
1452 * for the home zone.
1456 lexcls->a[n++] = jla;
1457 /* Check to avoid double counts */
1466 /* This is a inter-cg exclusion */
1467 /* Since exclusions are pair interactions,
1468 * just like non-bonded interactions,
1469 * they can be assigned properly up
1470 * to the DD cutoff (not cutoff_min as
1471 * for the other bonded interactions).
1473 if (ga2la_get(ga2la, a_gl+aj_mol-a_mol, &jla, &cell))
1475 if (iz == 0 && cell == 0)
1477 lexcls->a[n++] = jla;
1478 /* Check to avoid double counts */
1484 else if (jla >= jla0 && jla < jla1 &&
1486 dd_dist2(pbc_null, cg_cm, la2lc, la, jla) < rc2))
1488 /* jla > la, since jla0 > la */
1489 lexcls->a[n++] = jla;
1499 /* There are no inter-cg excls and this cg is self-excluded.
1500 * These exclusions are only required for zone 0,
1501 * since other zones do not see themselves.
1505 for (la = la0; la < la1; la++)
1507 lexcls->index[la] = n;
1508 for (j = la0; j < la1; j++)
1513 count += ((la1 - la0)*(la1 - la0 - 1))/2;
1517 /* We don't need exclusions for this cg */
1518 for (la = la0; la < la1; la++)
1520 lexcls->index[la] = n;
1526 lexcls->index[lexcls->nr] = n;
1532 static void check_alloc_index(t_blocka *ba, int nindex_max)
1534 if (nindex_max+1 > ba->nalloc_index)
1536 ba->nalloc_index = over_alloc_dd(nindex_max+1);
1537 srenew(ba->index, ba->nalloc_index);
1541 static void check_exclusions_alloc(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1547 nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1549 check_alloc_index(lexcls, nr);
1551 for (thread = 1; thread < dd->reverse_top->nthread; thread++)
1553 check_alloc_index(&dd->reverse_top->excl_thread[thread], nr);
1557 static void finish_local_exclusions(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1562 lexcls->nr = dd->cgindex[zones->izone[zones->nizone-1].cg1];
1564 if (dd->n_intercg_excl == 0)
1566 /* There are no exclusions involving non-home charge groups,
1567 * but we need to set the indices for neighborsearching.
1569 la0 = dd->cgindex[zones->izone[0].cg1];
1570 for (la = la0; la < lexcls->nr; la++)
1572 lexcls->index[la] = lexcls->nra;
1575 /* nr is only used to loop over the exclusions for Ewald and RF,
1576 * so we can set it to the number of home atoms for efficiency.
1578 lexcls->nr = dd->cgindex[zones->izone[0].cg1];
1582 static void clear_idef(t_idef *idef)
1586 /* Clear the counts */
1587 for (ftype = 0; ftype < F_NRE; ftype++)
1589 idef->il[ftype].nr = 0;
1593 static int make_local_bondeds_excls(gmx_domdec_t *dd,
1594 gmx_domdec_zones_t *zones,
1595 const gmx_mtop_t *mtop,
1597 gmx_bool bRCheckMB, ivec rcheck, gmx_bool bRCheck2B,
1599 int *la2lc, t_pbc *pbc_null, rvec *cg_cm,
1600 t_idef *idef, gmx_vsite_t *vsite,
1601 t_blocka *lexcls, int *excl_count)
1603 int nzone_bondeds, nzone_excl;
1608 gmx_reverse_top_t *rt;
1610 if (dd->reverse_top->bMultiCGmols)
1612 nzone_bondeds = zones->n;
1616 /* Only single charge group molecules, so interactions don't
1617 * cross zone boundaries and we only need to assign in the home zone.
1622 if (dd->n_intercg_excl > 0)
1624 /* We only use exclusions from i-zones to i- and j-zones */
1625 nzone_excl = zones->nizone;
1629 /* There are no inter-cg exclusions and only zone 0 sees itself */
1633 check_exclusions_alloc(dd, zones, lexcls);
1635 rt = dd->reverse_top;
1639 /* Clear the counts */
1647 for (iz = 0; iz < nzone_bondeds; iz++)
1649 cg0 = zones->cg_range[iz];
1650 cg1 = zones->cg_range[iz+1];
1652 #pragma omp parallel for num_threads(rt->nthread) schedule(static)
1653 for (thread = 0; thread < rt->nthread; thread++)
1659 int *vsite_pbc_nalloc;
1662 cg0t = cg0 + ((cg1 - cg0)* thread )/rt->nthread;
1663 cg1t = cg0 + ((cg1 - cg0)*(thread+1))/rt->nthread;
1671 idef_t = &rt->idef_thread[thread];
1675 if (vsite && vsite->bHaveChargeGroups && vsite->n_intercg_vsite > 0)
1679 vsite_pbc = vsite->vsite_pbc_loc;
1680 vsite_pbc_nalloc = vsite->vsite_pbc_loc_nalloc;
1684 vsite_pbc = rt->vsite_pbc[thread];
1685 vsite_pbc_nalloc = rt->vsite_pbc_nalloc[thread];
1691 vsite_pbc_nalloc = NULL;
1694 rt->nbonded_thread[thread] =
1695 make_bondeds_zone(dd, zones,
1697 bRCheckMB, rcheck, bRCheck2B, rc2,
1698 la2lc, pbc_null, cg_cm, idef->iparams,
1700 vsite, vsite_pbc, vsite_pbc_nalloc,
1702 dd->cgindex[cg0t], dd->cgindex[cg1t]);
1704 if (iz < nzone_excl)
1712 excl_t = &rt->excl_thread[thread];
1717 rt->excl_count_thread[thread] =
1718 make_exclusions_zone(dd, zones,
1719 mtop->moltype, bRCheck2B, rc2,
1720 la2lc, pbc_null, cg_cm, cginfo,
1727 if (rt->nthread > 1)
1729 combine_idef(idef, rt->idef_thread+1, rt->nthread-1,
1730 vsite, rt->vsite_pbc+1);
1733 for (thread = 0; thread < rt->nthread; thread++)
1735 nbonded_local += rt->nbonded_thread[thread];
1738 if (iz < nzone_excl)
1740 if (rt->nthread > 1)
1742 combine_blocka(lexcls, rt->excl_thread+1, rt->nthread-1);
1745 for (thread = 0; thread < rt->nthread; thread++)
1747 *excl_count += rt->excl_count_thread[thread];
1752 /* Some zones might not have exclusions, but some code still needs to
1753 * loop over the index, so we set the indices here.
1755 for (iz = nzone_excl; iz < zones->nizone; iz++)
1757 set_no_exclusions_zone(dd, zones, iz, lexcls);
1760 finish_local_exclusions(dd, zones, lexcls);
1763 fprintf(debug, "We have %d exclusions, check count %d\n",
1764 lexcls->nra, *excl_count);
1767 return nbonded_local;
1770 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs)
1772 lcgs->nr = dd->ncg_tot;
1773 lcgs->index = dd->cgindex;
1776 void dd_make_local_top(FILE *fplog,
1777 gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
1778 int npbcdim, matrix box,
1779 rvec cellsize_min, ivec npulse,
1783 gmx_mtop_t *mtop, gmx_localtop_t *ltop)
1785 gmx_bool bUniqueExcl, bRCheckMB, bRCheck2B, bRCheckExcl;
1789 t_pbc pbc, *pbc_null = NULL;
1793 fprintf(debug, "Making local topology\n");
1796 dd_make_local_cgs(dd, <op->cgs);
1800 bRCheckExcl = FALSE;
1802 if (dd->reverse_top->bMultiCGmols)
1804 /* We need to check to which cell bondeds should be assigned */
1805 rc = dd_cutoff_twobody(dd);
1808 fprintf(debug, "Two-body bonded cut-off distance is %g\n", rc);
1811 /* Should we check cg_cm distances when assigning bonded interactions? */
1812 for (d = 0; d < DIM; d++)
1815 /* Only need to check for dimensions where the part of the box
1816 * that is not communicated is smaller than the cut-off.
1818 if (d < npbcdim && dd->nc[d] > 1 &&
1819 (dd->nc[d] - npulse[d])*cellsize_min[d] < 2*rc)
1826 /* Check for interactions between two atoms,
1827 * where we can allow interactions up to the cut-off,
1828 * instead of up to the smallest cell dimension.
1835 "dim %d cellmin %f bonded rcheck[%d] = %d, bRCheck2B = %d\n",
1836 d, cellsize_min[d], d, rcheck[d], bRCheck2B);
1839 if (dd->reverse_top->bExclRequired)
1841 bRCheckExcl = bRCheck2B;
1845 /* If we don't have forces on exclusions,
1846 * we don't care about exclusions being assigned mulitple times.
1848 bRCheckExcl = FALSE;
1850 if (bRCheckMB || bRCheck2B)
1855 set_pbc_dd(&pbc, fr->ePBC, dd, TRUE, box);
1866 make_local_bondeds_excls(dd, zones, mtop, fr->cginfo,
1867 bRCheckMB, rcheck, bRCheck2B, rc,
1869 pbc_null, cgcm_or_x,
1871 <op->excls, &nexcl);
1873 /* The ilist is not sorted yet,
1874 * we can only do this when we have the charge arrays.
1876 ltop->idef.ilsort = ilsortUNKNOWN;
1878 if (dd->reverse_top->bExclRequired)
1880 dd->nbonded_local += nexcl;
1882 forcerec_set_excl_load(fr, ltop, NULL);
1885 ltop->atomtypes = mtop->atomtypes;
1887 /* For an error message only */
1888 dd->reverse_top->err_top_global = mtop;
1889 dd->reverse_top->err_top_local = ltop;
1892 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
1893 gmx_localtop_t *ltop)
1895 if (dd->reverse_top->ilsort == ilsortNO_FE)
1897 ltop->idef.ilsort = ilsortNO_FE;
1901 gmx_sort_ilist_fe(<op->idef, mdatoms->chargeA, mdatoms->chargeB);
1905 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global)
1907 gmx_localtop_t *top;
1912 top->idef.ntypes = top_global->ffparams.ntypes;
1913 top->idef.atnr = top_global->ffparams.atnr;
1914 top->idef.functype = top_global->ffparams.functype;
1915 top->idef.iparams = top_global->ffparams.iparams;
1916 top->idef.fudgeQQ = top_global->ffparams.fudgeQQ;
1917 top->idef.cmap_grid = top_global->ffparams.cmap_grid;
1919 for (i = 0; i < F_NRE; i++)
1921 top->idef.il[i].iatoms = NULL;
1922 top->idef.il[i].nalloc = 0;
1924 top->idef.ilsort = ilsortUNKNOWN;
1929 void dd_init_local_state(gmx_domdec_t *dd,
1930 t_state *state_global, t_state *state_local)
1932 int buf[NITEM_DD_INIT_LOCAL_STATE];
1936 buf[0] = state_global->flags;
1937 buf[1] = state_global->ngtc;
1938 buf[2] = state_global->nnhpres;
1939 buf[3] = state_global->nhchainlength;
1940 buf[4] = state_global->dfhist.nlambda;
1942 dd_bcast(dd, NITEM_DD_INIT_LOCAL_STATE*sizeof(int), buf);
1944 init_state(state_local, 0, buf[1], buf[2], buf[3], buf[4]);
1945 state_local->flags = buf[0];
1947 /* With Langevin Dynamics we need to make proper storage space
1948 * in the global and local state for the random numbers.
1950 if (state_local->flags & (1<<estLD_RNG))
1952 if (DDMASTER(dd) && state_global->nrngi > 1)
1954 state_global->nrng = dd->nnodes*gmx_rng_n();
1955 srenew(state_global->ld_rng, state_global->nrng);
1957 state_local->nrng = gmx_rng_n();
1958 snew(state_local->ld_rng, state_local->nrng);
1960 if (state_local->flags & (1<<estLD_RNGI))
1962 if (DDMASTER(dd) && state_global->nrngi > 1)
1964 state_global->nrngi = dd->nnodes;
1965 srenew(state_global->ld_rngi, state_global->nrngi);
1967 snew(state_local->ld_rngi, 1);
1971 static void check_link(t_blocka *link, int cg_gl, int cg_gl_j)
1977 for (k = link->index[cg_gl]; k < link->index[cg_gl+1]; k++)
1979 if (link->a[k] == cg_gl_j)
1986 /* Add this charge group link */
1987 if (link->index[cg_gl+1]+1 > link->nalloc_a)
1989 link->nalloc_a = over_alloc_large(link->index[cg_gl+1]+1);
1990 srenew(link->a, link->nalloc_a);
1992 link->a[link->index[cg_gl+1]] = cg_gl_j;
1993 link->index[cg_gl+1]++;
1997 static int *make_at2cg(t_block *cgs)
2001 snew(at2cg, cgs->index[cgs->nr]);
2002 for (cg = 0; cg < cgs->nr; cg++)
2004 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2013 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
2014 cginfo_mb_t *cginfo_mb)
2016 gmx_reverse_top_t *rt;
2017 int mb, cg_offset, cg, cg_gl, a, aj, i, j, ftype, nral, nlink_mol, mol, ncgi;
2018 gmx_molblock_t *molb;
2019 gmx_moltype_t *molt;
2023 gmx_reverse_ilist_t ril;
2025 cginfo_mb_t *cgi_mb;
2027 /* For each charge group make a list of other charge groups
2028 * in the system that a linked to it via bonded interactions
2029 * which are also stored in reverse_top.
2032 rt = dd->reverse_top;
2035 snew(link->index, ncg_mtop(mtop)+1);
2042 for (mb = 0; mb < mtop->nmolblock; mb++)
2044 molb = &mtop->molblock[mb];
2045 if (molb->nmol == 0)
2049 molt = &mtop->moltype[molb->type];
2051 excls = &molt->excls;
2052 a2c = make_at2cg(cgs);
2053 /* Make a reverse ilist in which the interactions are linked
2054 * to all atoms, not only the first atom as in gmx_reverse_top.
2055 * The constraints are discarded here.
2057 make_reverse_ilist(molt, NULL, FALSE, FALSE, FALSE, TRUE, &ril);
2059 cgi_mb = &cginfo_mb[mb];
2061 for (cg = 0; cg < cgs->nr; cg++)
2063 cg_gl = cg_offset + cg;
2064 link->index[cg_gl+1] = link->index[cg_gl];
2065 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
2068 while (i < ril.index[a+1])
2070 ftype = ril.il[i++];
2072 /* Skip the ifunc index */
2074 for (j = 0; j < nral; j++)
2079 check_link(link, cg_gl, cg_offset+a2c[aj]);
2082 i += nral_rt(ftype);
2084 if (rt->bExclRequired)
2086 /* Exclusions always go both ways */
2087 for (j = excls->index[a]; j < excls->index[a+1]; j++)
2092 check_link(link, cg_gl, cg_offset+a2c[aj]);
2097 if (link->index[cg_gl+1] - link->index[cg_gl] > 0)
2099 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg]);
2103 nlink_mol = link->index[cg_offset+cgs->nr] - link->index[cg_offset];
2105 cg_offset += cgs->nr;
2107 destroy_reverse_ilist(&ril);
2112 fprintf(debug, "molecule type '%s' %d cgs has %d cg links through bonded interac.\n", *molt->name, cgs->nr, nlink_mol);
2117 /* Copy the data for the rest of the molecules in this block */
2118 link->nalloc_a += (molb->nmol - 1)*nlink_mol;
2119 srenew(link->a, link->nalloc_a);
2120 for (mol = 1; mol < molb->nmol; mol++)
2122 for (cg = 0; cg < cgs->nr; cg++)
2124 cg_gl = cg_offset + cg;
2125 link->index[cg_gl+1] =
2126 link->index[cg_gl+1-cgs->nr] + nlink_mol;
2127 for (j = link->index[cg_gl]; j < link->index[cg_gl+1]; j++)
2129 link->a[j] = link->a[j-nlink_mol] + cgs->nr;
2131 if (link->index[cg_gl+1] - link->index[cg_gl] > 0 &&
2132 cg_gl - cgi_mb->cg_start < cgi_mb->cg_mod)
2134 SET_CGINFO_BOND_INTER(cgi_mb->cginfo[cg_gl - cgi_mb->cg_start]);
2138 cg_offset += cgs->nr;
2145 fprintf(debug, "Of the %d charge groups %d are linked via bonded interactions\n", ncg_mtop(mtop), ncgi);
2151 static void bonded_cg_distance_mol(gmx_moltype_t *molt, int *at2cg,
2152 gmx_bool bBCheck, gmx_bool bExcl, rvec *cg_cm,
2153 real *r_2b, int *ft2b, int *a2_1, int *a2_2,
2154 real *r_mb, int *ftmb, int *am_1, int *am_2)
2156 int ftype, nral, i, j, ai, aj, cgi, cgj;
2159 real r2_2b, r2_mb, rij2;
2163 for (ftype = 0; ftype < F_NRE; ftype++)
2165 if (dd_check_ftype(ftype, bBCheck, FALSE, FALSE))
2167 il = &molt->ilist[ftype];
2171 for (i = 0; i < il->nr; i += 1+nral)
2173 for (ai = 0; ai < nral; ai++)
2175 cgi = at2cg[il->iatoms[i+1+ai]];
2176 for (aj = 0; aj < nral; aj++)
2178 cgj = at2cg[il->iatoms[i+1+aj]];
2181 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2182 if (nral == 2 && rij2 > r2_2b)
2186 *a2_1 = il->iatoms[i+1+ai];
2187 *a2_2 = il->iatoms[i+1+aj];
2189 if (nral > 2 && rij2 > r2_mb)
2193 *am_1 = il->iatoms[i+1+ai];
2194 *am_2 = il->iatoms[i+1+aj];
2205 excls = &molt->excls;
2206 for (ai = 0; ai < excls->nr; ai++)
2209 for (j = excls->index[ai]; j < excls->index[ai+1]; j++)
2211 cgj = at2cg[excls->a[j]];
2214 rij2 = distance2(cg_cm[cgi], cg_cm[cgj]);
2224 *r_2b = sqrt(r2_2b);
2225 *r_mb = sqrt(r2_mb);
2228 static void get_cgcm_mol(gmx_moltype_t *molt, gmx_ffparams_t *ffparams,
2229 int ePBC, t_graph *graph, matrix box,
2231 rvec *x, rvec *xs, rvec *cg_cm)
2235 if (ePBC != epbcNONE)
2237 mk_mshift(NULL, graph, ePBC, box, x);
2239 shift_x(graph, box, x, xs);
2240 /* By doing an extra mk_mshift the molecules that are broken
2241 * because they were e.g. imported from another software
2242 * will be made whole again. Such are the healing powers
2245 mk_mshift(NULL, graph, ePBC, box, xs);
2249 /* We copy the coordinates so the original coordinates remain
2250 * unchanged, just to be 100% sure that we do not affect
2251 * binary reproducibility of simulations.
2253 n = molt->cgs.index[molt->cgs.nr];
2254 for (i = 0; i < n; i++)
2256 copy_rvec(x[i], xs[i]);
2262 construct_vsites(NULL, vsite, xs, NULL, 0.0, NULL,
2263 ffparams->iparams, molt->ilist,
2264 epbcNONE, TRUE, NULL, NULL, NULL);
2267 calc_cgcm(NULL, 0, molt->cgs.nr, &molt->cgs, xs, cg_cm);
2270 static int have_vsite_molt(gmx_moltype_t *molt)
2276 for (i = 0; i < F_NRE; i++)
2278 if ((interaction_function[i].flags & IF_VSITE) &&
2279 molt->ilist[i].nr > 0)
2288 void dd_bonded_cg_distance(FILE *fplog,
2289 gmx_domdec_t *dd, gmx_mtop_t *mtop,
2290 t_inputrec *ir, rvec *x, matrix box,
2292 real *r_2b, real *r_mb)
2294 gmx_bool bExclRequired;
2295 int mb, cg_offset, at_offset, *at2cg, mol;
2298 gmx_molblock_t *molb;
2299 gmx_moltype_t *molt;
2301 real rmol_2b, rmol_mb;
2302 int ft2b = -1, a_2b_1 = -1, a_2b_2 = -1, ftmb = -1, a_mb_1 = -1, a_mb_2 = -1;
2303 int ftm2b = -1, amol_2b_1 = -1, amol_2b_2 = -1, ftmmb = -1, amol_mb_1 = -1, amol_mb_2 = -1;
2305 bExclRequired = IR_EXCL_FORCES(*ir);
2307 vsite = init_vsite(mtop, NULL, TRUE);
2313 for (mb = 0; mb < mtop->nmolblock; mb++)
2315 molb = &mtop->molblock[mb];
2316 molt = &mtop->moltype[molb->type];
2317 if (molt->cgs.nr == 1 || molb->nmol == 0)
2319 cg_offset += molb->nmol*molt->cgs.nr;
2320 at_offset += molb->nmol*molt->atoms.nr;
2324 if (ir->ePBC != epbcNONE)
2326 mk_graph_ilist(NULL, molt->ilist, 0, molt->atoms.nr, FALSE, FALSE,
2330 at2cg = make_at2cg(&molt->cgs);
2331 snew(xs, molt->atoms.nr);
2332 snew(cg_cm, molt->cgs.nr);
2333 for (mol = 0; mol < molb->nmol; mol++)
2335 get_cgcm_mol(molt, &mtop->ffparams, ir->ePBC, &graph, box,
2336 have_vsite_molt(molt) ? vsite : NULL,
2337 x+at_offset, xs, cg_cm);
2339 bonded_cg_distance_mol(molt, at2cg, bBCheck, bExclRequired, cg_cm,
2340 &rmol_2b, &ftm2b, &amol_2b_1, &amol_2b_2,
2341 &rmol_mb, &ftmmb, &amol_mb_1, &amol_mb_2);
2342 if (rmol_2b > *r_2b)
2346 a_2b_1 = at_offset + amol_2b_1;
2347 a_2b_2 = at_offset + amol_2b_2;
2349 if (rmol_mb > *r_mb)
2353 a_mb_1 = at_offset + amol_mb_1;
2354 a_mb_2 = at_offset + amol_mb_2;
2357 cg_offset += molt->cgs.nr;
2358 at_offset += molt->atoms.nr;
2363 if (ir->ePBC != epbcNONE)
2370 /* We should have a vsite free routine, but here we can simply free */
2373 if (fplog && (ft2b >= 0 || ftmb >= 0))
2376 "Initial maximum inter charge-group distances:\n");
2380 " two-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2381 *r_2b, interaction_function[ft2b].longname,
2382 a_2b_1+1, a_2b_2+1);
2387 " multi-body bonded interactions: %5.3f nm, %s, atoms %d %d\n",
2388 *r_mb, interaction_function[ftmb].longname,
2389 a_mb_1+1, a_mb_2+1);