2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2012,2013, 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.
45 #include "domdec_network.h"
46 #include "mtop_util.h"
47 #include "gmx_ga2la.h"
49 #include "gmx_omp_nthreads.h"
65 typedef struct gmx_domdec_specat_comm {
66 /* The number of indices to receive during the setup */
68 /* The atoms to send */
69 gmx_specatsend_t spas[DIM][2];
79 /* The range in the local buffer(s) for received atoms */
83 /* The atom indices we need from the surrounding cells.
84 * We can gather the indices over nthread threads.
88 } gmx_domdec_specat_comm_t;
90 typedef struct gmx_domdec_constraints {
93 /* The fully local and connected constraints */
95 /* The global constraint number, only required for clearing gc_req */
99 /* Boolean that tells if a global constraint index has been requested */
101 /* Global to local communicated constraint atom only index */
104 /* Multi-threading stuff */
107 } gmx_domdec_constraints_t;
110 static void dd_move_f_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
111 rvec *f, rvec *fshift)
113 gmx_specatsend_t *spas;
115 int n, n0, n1, d, dim, dir, i;
118 gmx_bool bPBC, bScrew;
121 for (d = dd->ndim-1; d >= 0; d--)
126 /* Pulse the grid forward and backward */
127 spas = spac->spas[d];
132 /* Send and receive the coordinates */
133 dd_sendrecv2_rvec(dd, d,
134 f+n+n1, n0, vbuf, spas[0].nsend,
135 f+n, n1, vbuf+spas[0].nsend, spas[1].nsend);
136 for (dir = 0; dir < 2; dir++)
138 bPBC = ((dir == 0 && dd->ci[dim] == 0) ||
139 (dir == 1 && dd->ci[dim] == dd->nc[dim]-1));
140 bScrew = (bPBC && dd->bScrewPBC && dim == XX);
142 spas = &spac->spas[d][dir];
143 /* Sum the buffer into the required forces */
144 if (!bPBC || (!bScrew && fshift == NULL))
146 for (i = 0; i < spas->nsend; i++)
148 rvec_inc(f[spas->a[i]], *vbuf);
155 vis[dim] = (dir == 0 ? 1 : -1);
159 /* Sum and add to shift forces */
160 for (i = 0; i < spas->nsend; i++)
162 rvec_inc(f[spas->a[i]], *vbuf);
163 rvec_inc(fshift[is], *vbuf);
169 /* Rotate the forces */
170 for (i = 0; i < spas->nsend; i++)
172 f[spas->a[i]][XX] += (*vbuf)[XX];
173 f[spas->a[i]][YY] -= (*vbuf)[YY];
174 f[spas->a[i]][ZZ] -= (*vbuf)[ZZ];
177 rvec_inc(fshift[is], *vbuf);
187 /* Two cells, so we only need to communicate one way */
188 spas = &spac->spas[d][0];
190 /* Send and receive the coordinates */
191 dd_sendrecv_rvec(dd, d, dddirForward,
192 f+n, spas->nrecv, spac->vbuf, spas->nsend);
193 /* Sum the buffer into the required forces */
194 if (dd->bScrewPBC && dim == XX &&
196 dd->ci[dim] == dd->nc[dim]-1))
198 for (i = 0; i < spas->nsend; i++)
200 /* Rotate the force */
201 f[spas->a[i]][XX] += spac->vbuf[i][XX];
202 f[spas->a[i]][YY] -= spac->vbuf[i][YY];
203 f[spas->a[i]][ZZ] -= spac->vbuf[i][ZZ];
208 for (i = 0; i < spas->nsend; i++)
210 rvec_inc(f[spas->a[i]], spac->vbuf[i]);
217 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift)
221 dd_move_f_specat(dd, dd->vsite_comm, f, fshift);
225 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f)
231 for (i = dd->vsite_comm->at_start; i < dd->vsite_comm->at_end; i++)
238 static void dd_move_x_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
239 matrix box, rvec *x0, rvec *x1)
241 gmx_specatsend_t *spas;
242 rvec *x, *vbuf, *rbuf;
243 int nvec, v, n, nn, ns0, ns1, nr0, nr1, nr, d, dim, dir, i;
244 gmx_bool bPBC, bScrew = FALSE;
245 rvec shift = {0, 0, 0};
254 for (d = 0; d < dd->ndim; d++)
259 /* Pulse the grid forward and backward */
261 for (dir = 0; dir < 2; dir++)
263 if (dir == 0 && dd->ci[dim] == 0)
266 bScrew = (dd->bScrewPBC && dim == XX);
267 copy_rvec(box[dim], shift);
269 else if (dir == 1 && dd->ci[dim] == dd->nc[dim]-1)
272 bScrew = (dd->bScrewPBC && dim == XX);
273 for (i = 0; i < DIM; i++)
275 shift[i] = -box[dim][i];
283 spas = &spac->spas[d][dir];
284 for (v = 0; v < nvec; v++)
286 x = (v == 0 ? x0 : x1);
287 /* Copy the required coordinates to the send buffer */
291 for (i = 0; i < spas->nsend; i++)
293 copy_rvec(x[spas->a[i]], *vbuf);
299 /* Shift coordinates */
300 for (i = 0; i < spas->nsend; i++)
302 rvec_add(x[spas->a[i]], shift, *vbuf);
308 /* Shift and rotate coordinates */
309 for (i = 0; i < spas->nsend; i++)
311 (*vbuf)[XX] = x[spas->a[i]][XX] + shift[XX];
312 (*vbuf)[YY] = box[YY][YY] - x[spas->a[i]][YY] + shift[YY];
313 (*vbuf)[ZZ] = box[ZZ][ZZ] - x[spas->a[i]][ZZ] + shift[ZZ];
319 /* Send and receive the coordinates */
320 spas = spac->spas[d];
327 dd_sendrecv2_rvec(dd, d,
328 spac->vbuf+ns0, ns1, x0+n, nr1,
329 spac->vbuf, ns0, x0+n+nr1, nr0);
333 /* Communicate both vectors in one buffer */
335 dd_sendrecv2_rvec(dd, d,
336 spac->vbuf+2*ns0, 2*ns1, rbuf, 2*nr1,
337 spac->vbuf, 2*ns0, rbuf+2*nr1, 2*nr0);
338 /* Split the buffer into the two vectors */
340 for (dir = 1; dir >= 0; dir--)
342 nr = spas[dir].nrecv;
343 for (v = 0; v < 2; v++)
345 x = (v == 0 ? x0 : x1);
346 for (i = 0; i < nr; i++)
348 copy_rvec(*rbuf, x[nn+i]);
359 spas = &spac->spas[d][0];
360 /* Copy the required coordinates to the send buffer */
362 for (v = 0; v < nvec; v++)
364 x = (v == 0 ? x0 : x1);
365 if (dd->bScrewPBC && dim == XX &&
366 (dd->ci[XX] == 0 || dd->ci[XX] == dd->nc[XX]-1))
368 /* Here we only perform the rotation, the rest of the pbc
369 * is handled in the constraint or viste routines.
371 for (i = 0; i < spas->nsend; i++)
373 (*vbuf)[XX] = x[spas->a[i]][XX];
374 (*vbuf)[YY] = box[YY][YY] - x[spas->a[i]][YY];
375 (*vbuf)[ZZ] = box[ZZ][ZZ] - x[spas->a[i]][ZZ];
381 for (i = 0; i < spas->nsend; i++)
383 copy_rvec(x[spas->a[i]], *vbuf);
388 /* Send and receive the coordinates */
391 dd_sendrecv_rvec(dd, d, dddirBackward,
392 spac->vbuf, spas->nsend, x0+n, spas->nrecv);
396 /* Communicate both vectors in one buffer */
398 dd_sendrecv_rvec(dd, d, dddirBackward,
399 spac->vbuf, 2*spas->nsend, rbuf, 2*spas->nrecv);
400 /* Split the buffer into the two vectors */
402 for (v = 0; v < 2; v++)
404 x = (v == 0 ? x0 : x1);
405 for (i = 0; i < nr; i++)
407 copy_rvec(*rbuf, x[n+i]);
417 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box, rvec *x0, rvec *x1)
419 if (dd->constraint_comm)
421 dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1);
425 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x)
429 dd_move_x_specat(dd, dd->vsite_comm, box, x, NULL);
433 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd)
437 return dd->constraints->con_nlocat;
445 void dd_clear_local_constraint_indices(gmx_domdec_t *dd)
447 gmx_domdec_constraints_t *dc;
450 dc = dd->constraints;
452 for (i = 0; i < dc->ncon; i++)
454 dc->gc_req[dc->con_gl[i]] = 0;
457 if (dd->constraint_comm)
459 gmx_hash_clear_and_optimize(dc->ga2la);
463 void dd_clear_local_vsite_indices(gmx_domdec_t *dd)
469 gmx_hash_clear_and_optimize(dd->ga2la_vsite);
473 static int setup_specat_communication(gmx_domdec_t *dd,
475 gmx_domdec_specat_comm_t *spac,
476 gmx_hash_t ga2la_specat,
479 const char *specat_type,
482 int nsend[2], nlast, nsend_zero[2] = {0, 0}, *nsend_ptr;
483 int d, dim, ndir, dir, nr, ns, i, nrecv_local, n0, start, indr, ind, buf[2];
484 int nat_tot_specat, nat_tot_prev, nalloc_old;
485 gmx_bool bPBC, bFirst;
486 gmx_specatsend_t *spas;
490 fprintf(debug, "Begin setup_specat_communication for %s\n", specat_type);
493 /* nsend[0]: the number of atoms requested by this node only,
494 * we communicate this for more efficients checks
495 * nsend[1]: the total number of requested atoms
500 for (d = dd->ndim-1; d >= 0; d--)
502 /* Pulse the grid forward and backward */
504 bPBC = (dim < dd->npbcdim);
505 if (dd->nc[dim] == 2)
507 /* Only 2 cells, so we only need to communicate once */
514 for (dir = 0; dir < ndir; dir++)
518 ((dir == 0 && dd->ci[dim] == dd->nc[dim] - 1) ||
519 (dir == 1 && dd->ci[dim] == 0)))
521 /* No pbc: the fist/last cell should not request atoms */
522 nsend_ptr = nsend_zero;
528 /* Communicate the number of indices */
529 dd_sendrecv_int(dd, d, dir == 0 ? dddirForward : dddirBackward,
530 nsend_ptr, 2, spac->nreq[d][dir], 2);
531 nr = spac->nreq[d][dir][1];
532 if (nlast+nr > ireq->nalloc)
534 ireq->nalloc = over_alloc_dd(nlast+nr);
535 srenew(ireq->ind, ireq->nalloc);
537 /* Communicate the indices */
538 dd_sendrecv_int(dd, d, dir == 0 ? dddirForward : dddirBackward,
539 ireq->ind, nsend_ptr[1], ireq->ind+nlast, nr);
546 fprintf(debug, "Communicated the counts\n");
549 /* Search for the requested atoms and communicate the indices we have */
550 nat_tot_specat = at_start;
552 for (d = 0; d < dd->ndim; d++)
555 /* Pulse the grid forward and backward */
556 if (dd->dim[d] >= dd->npbcdim || dd->nc[dd->dim[d]] > 2)
564 nat_tot_prev = nat_tot_specat;
565 for (dir = ndir-1; dir >= 0; dir--)
567 if (nat_tot_specat > spac->bSendAtom_nalloc)
569 nalloc_old = spac->bSendAtom_nalloc;
570 spac->bSendAtom_nalloc = over_alloc_dd(nat_tot_specat);
571 srenew(spac->bSendAtom, spac->bSendAtom_nalloc);
572 for (i = nalloc_old; i < spac->bSendAtom_nalloc; i++)
574 spac->bSendAtom[i] = FALSE;
577 spas = &spac->spas[d][dir];
578 n0 = spac->nreq[d][dir][0];
579 nr = spac->nreq[d][dir][1];
582 fprintf(debug, "dim=%d, dir=%d, searching for %d atoms\n",
588 for (i = 0; i < nr; i++)
590 indr = ireq->ind[start+i];
592 /* Check if this is a home atom and if so ind will be set */
593 if (!ga2la_get_home(dd->ga2la, indr, &ind))
595 /* Search in the communicated atoms */
596 ind = gmx_hash_get_minone(ga2la_specat, indr);
600 if (i < n0 || !spac->bSendAtom[ind])
602 if (spas->nsend+1 > spas->a_nalloc)
604 spas->a_nalloc = over_alloc_large(spas->nsend+1);
605 srenew(spas->a, spas->a_nalloc);
607 /* Store the local index so we know which coordinates
610 spas->a[spas->nsend] = ind;
611 spac->bSendAtom[ind] = TRUE;
612 if (spas->nsend+1 > spac->ibuf_nalloc)
614 spac->ibuf_nalloc = over_alloc_large(spas->nsend+1);
615 srenew(spac->ibuf, spac->ibuf_nalloc);
617 /* Store the global index so we can send it now */
618 spac->ibuf[spas->nsend] = indr;
628 /* Clear the local flags */
629 for (i = 0; i < spas->nsend; i++)
631 spac->bSendAtom[spas->a[i]] = FALSE;
633 /* Send and receive the number of indices to communicate */
634 nsend[1] = spas->nsend;
635 dd_sendrecv_int(dd, d, dir == 0 ? dddirBackward : dddirForward,
639 fprintf(debug, "Send to node %d, %d (%d) indices, "
640 "receive from node %d, %d (%d) indices\n",
641 dd->neighbor[d][1-dir], nsend[1], nsend[0],
642 dd->neighbor[d][dir], buf[1], buf[0]);
645 for (i = 0; i < spas->nsend; i++)
647 fprintf(debug, " %d", spac->ibuf[i]+1);
649 fprintf(debug, "\n");
652 nrecv_local += buf[0];
653 spas->nrecv = buf[1];
654 if (nat_tot_specat + spas->nrecv > dd->gatindex_nalloc)
656 dd->gatindex_nalloc =
657 over_alloc_dd(nat_tot_specat + spas->nrecv);
658 srenew(dd->gatindex, dd->gatindex_nalloc);
660 /* Send and receive the indices */
661 dd_sendrecv_int(dd, d, dir == 0 ? dddirBackward : dddirForward,
662 spac->ibuf, spas->nsend,
663 dd->gatindex+nat_tot_specat, spas->nrecv);
664 nat_tot_specat += spas->nrecv;
667 /* Allocate the x/f communication buffers */
668 ns = spac->spas[d][0].nsend;
669 nr = spac->spas[d][0].nrecv;
672 ns += spac->spas[d][1].nsend;
673 nr += spac->spas[d][1].nrecv;
675 if (vbuf_fac*ns > spac->vbuf_nalloc)
677 spac->vbuf_nalloc = over_alloc_dd(vbuf_fac*ns);
678 srenew(spac->vbuf, spac->vbuf_nalloc);
680 if (vbuf_fac == 2 && vbuf_fac*nr > spac->vbuf2_nalloc)
682 spac->vbuf2_nalloc = over_alloc_dd(vbuf_fac*nr);
683 srenew(spac->vbuf2, spac->vbuf2_nalloc);
686 /* Make a global to local index for the communication atoms */
687 for (i = nat_tot_prev; i < nat_tot_specat; i++)
689 gmx_hash_change_or_set(ga2la_specat, dd->gatindex[i], i);
693 /* Check that in the end we got the number of atoms we asked for */
694 if (nrecv_local != ireq->n)
698 fprintf(debug, "Requested %d, received %d (tot recv %d)\n",
699 ireq->n, nrecv_local, nat_tot_specat-at_start);
702 for (i = 0; i < ireq->n; i++)
704 ind = gmx_hash_get_minone(ga2la_specat, ireq->ind[i]);
705 fprintf(debug, " %s%d",
706 (ind >= 0) ? "" : "!",
709 fprintf(debug, "\n");
712 fprintf(stderr, "\nDD cell %d %d %d: Neighboring cells do not have atoms:",
713 dd->ci[XX], dd->ci[YY], dd->ci[ZZ]);
714 for (i = 0; i < ireq->n; i++)
716 if (gmx_hash_get_minone(ga2la_specat, ireq->ind[i]) < 0)
718 fprintf(stderr, " %d", ireq->ind[i]+1);
721 fprintf(stderr, "\n");
722 gmx_fatal(FARGS, "DD cell %d %d %d could only obtain %d of the %d atoms that are connected via %ss from the neighboring cells. This probably means your %s lengths are too long compared to the domain decomposition cell size. Decrease the number of domain decomposition grid cells%s%s.",
723 dd->ci[XX], dd->ci[YY], dd->ci[ZZ],
724 nrecv_local, ireq->n, specat_type,
725 specat_type, add_err,
726 dd->bGridJump ? " or use the -rcon option of mdrun" : "");
729 spac->at_start = at_start;
730 spac->at_end = nat_tot_specat;
734 fprintf(debug, "Done setup_specat_communication\n");
737 return nat_tot_specat;
740 static void walk_out(int con, int con_offset, int a, int offset, int nrec,
741 int ncon1, const t_iatom *ia1, const t_iatom *ia2,
742 const t_blocka *at2con,
743 const gmx_ga2la_t ga2la, gmx_bool bHomeConnect,
744 gmx_domdec_constraints_t *dc,
745 gmx_domdec_specat_comm_t *dcc,
749 int a1_gl, a2_gl, a_loc, i, coni, b;
752 if (dc->gc_req[con_offset+con] == 0)
754 /* Add this non-home constraint to the list */
755 if (dc->ncon+1 > dc->con_nalloc)
757 dc->con_nalloc = over_alloc_large(dc->ncon+1);
758 srenew(dc->con_gl, dc->con_nalloc);
759 srenew(dc->con_nlocat, dc->con_nalloc);
761 dc->con_gl[dc->ncon] = con_offset + con;
762 dc->con_nlocat[dc->ncon] = (bHomeConnect ? 1 : 0);
763 dc->gc_req[con_offset+con] = 1;
764 if (il_local->nr + 3 > il_local->nalloc)
766 il_local->nalloc = over_alloc_dd(il_local->nr+3);
767 srenew(il_local->iatoms, il_local->nalloc);
769 iap = constr_iatomptr(ncon1, ia1, ia2, con);
770 il_local->iatoms[il_local->nr++] = iap[0];
771 a1_gl = offset + iap[1];
772 a2_gl = offset + iap[2];
773 /* The following indexing code can probably be optizimed */
774 if (ga2la_get_home(ga2la, a1_gl, &a_loc))
776 il_local->iatoms[il_local->nr++] = a_loc;
780 /* We set this index later */
781 il_local->iatoms[il_local->nr++] = -a1_gl - 1;
783 if (ga2la_get_home(ga2la, a2_gl, &a_loc))
785 il_local->iatoms[il_local->nr++] = a_loc;
789 /* We set this index later */
790 il_local->iatoms[il_local->nr++] = -a2_gl - 1;
794 /* Check to not ask for the same atom more than once */
795 if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
798 /* Add this non-home atom to the list */
799 if (ireq->n+1 > ireq->nalloc)
801 ireq->nalloc = over_alloc_large(ireq->n+1);
802 srenew(ireq->ind, ireq->nalloc);
804 ireq->ind[ireq->n++] = offset + a;
805 /* Temporarily mark with -2, we get the index later */
806 gmx_hash_set(dc->ga2la, offset+a, -2);
811 for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
817 iap = constr_iatomptr(ncon1, ia1, ia2, coni);
826 if (!ga2la_get_home(ga2la, offset+b, &a_loc))
828 walk_out(coni, con_offset, b, offset, nrec-1,
829 ncon1, ia1, ia2, at2con,
830 ga2la, FALSE, dc, dcc, il_local, ireq);
837 static void atoms_to_settles(gmx_domdec_t *dd,
838 const gmx_mtop_t *mtop,
840 const int **at2settle_mt,
841 int cg_start, int cg_end,
846 gmx_mtop_atomlookup_t alook;
849 int cg, a, a_gl, a_glsa, a_gls[3], a_locs[3];
850 int mb, molnr, a_mol, offset;
851 const gmx_molblock_t *molb;
859 alook = gmx_mtop_atomlookup_settle_init(mtop);
861 nral = NRAL(F_SETTLE);
863 for (cg = cg_start; cg < cg_end; cg++)
865 if (GET_CGINFO_SETTLE(cginfo[cg]))
867 for (a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
869 a_gl = dd->gatindex[a];
871 gmx_mtop_atomnr_to_molblock_ind(alook, a_gl, &mb, &molnr, &a_mol);
872 molb = &mtop->molblock[mb];
874 settle = at2settle_mt[molb->type][a_mol];
878 offset = a_gl - a_mol;
880 ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms;
884 for (sa = 0; sa < nral; sa++)
886 a_glsa = offset + ia1[settle*(1+nral)+1+sa];
888 a_home[sa] = ga2la_get_home(ga2la, a_glsa, &a_locs[sa]);
891 if (nlocal == 0 && a_gl == a_glsa)
901 if (ils_local->nr+1+nral > ils_local->nalloc)
903 ils_local->nalloc = over_alloc_dd(ils_local->nr+1+nral);
904 srenew(ils_local->iatoms, ils_local->nalloc);
907 ils_local->iatoms[ils_local->nr++] = ia1[settle*4];
909 for (sa = 0; sa < nral; sa++)
911 if (ga2la_get_home(ga2la, a_gls[sa], &a_locs[sa]))
913 ils_local->iatoms[ils_local->nr++] = a_locs[sa];
917 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
918 /* Add this non-home atom to the list */
919 if (ireq->n+1 > ireq->nalloc)
921 ireq->nalloc = over_alloc_large(ireq->n+1);
922 srenew(ireq->ind, ireq->nalloc);
924 ireq->ind[ireq->n++] = a_gls[sa];
925 /* A check on double atom requests is
926 * not required for settle.
936 gmx_mtop_atomlookup_destroy(alook);
939 static void atoms_to_constraints(gmx_domdec_t *dd,
940 const gmx_mtop_t *mtop,
942 const t_blocka *at2con_mt, int nrec,
946 const t_blocka *at2con;
948 gmx_mtop_atomlookup_t alook;
950 gmx_molblock_t *molb;
951 t_iatom *ia1, *ia2, *iap;
952 int nhome, cg, a, a_gl, a_mol, a_loc, b_lo, offset, mb, molnr, b_mol, i, con, con_offset;
953 gmx_domdec_constraints_t *dc;
954 gmx_domdec_specat_comm_t *dcc;
956 dc = dd->constraints;
957 dcc = dd->constraint_comm;
961 alook = gmx_mtop_atomlookup_init(mtop);
964 for (cg = 0; cg < dd->ncg_home; cg++)
966 if (GET_CGINFO_CONSTR(cginfo[cg]))
968 for (a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
970 a_gl = dd->gatindex[a];
972 gmx_mtop_atomnr_to_molblock_ind(alook, a_gl, &mb, &molnr, &a_mol);
973 molb = &mtop->molblock[mb];
975 ncon1 = mtop->moltype[molb->type].ilist[F_CONSTR].nr/NRAL(F_SETTLE);
977 ia1 = mtop->moltype[molb->type].ilist[F_CONSTR].iatoms;
978 ia2 = mtop->moltype[molb->type].ilist[F_CONSTRNC].iatoms;
980 /* Calculate the global constraint number offset for the molecule.
981 * This is only required for the global index to make sure
982 * that we use each constraint only once.
985 dc->molb_con_offset[mb] + molnr*dc->molb_ncon_mol[mb];
987 /* The global atom number offset for this molecule */
988 offset = a_gl - a_mol;
989 at2con = &at2con_mt[molb->type];
990 for (i = at2con->index[a_mol]; i < at2con->index[a_mol+1]; i++)
993 iap = constr_iatomptr(ncon1, ia1, ia2, con);
1002 if (ga2la_get_home(ga2la, offset+b_mol, &a_loc))
1004 /* Add this fully home constraint at the first atom */
1007 if (dc->ncon+1 > dc->con_nalloc)
1009 dc->con_nalloc = over_alloc_large(dc->ncon+1);
1010 srenew(dc->con_gl, dc->con_nalloc);
1011 srenew(dc->con_nlocat, dc->con_nalloc);
1013 dc->con_gl[dc->ncon] = con_offset + con;
1014 dc->con_nlocat[dc->ncon] = 2;
1015 if (ilc_local->nr + 3 > ilc_local->nalloc)
1017 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
1018 srenew(ilc_local->iatoms, ilc_local->nalloc);
1021 ilc_local->iatoms[ilc_local->nr++] = iap[0];
1022 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
1023 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a );
1030 /* We need the nrec constraints coupled to this constraint,
1031 * so we need to walk out of the home cell by nrec+1 atoms,
1032 * since already atom bg is not locally present.
1033 * Therefore we call walk_out with nrec recursions to go
1034 * after this first call.
1036 walk_out(con, con_offset, b_mol, offset, nrec,
1037 ncon1, ia1, ia2, at2con,
1038 dd->ga2la, TRUE, dc, dcc, ilc_local, ireq);
1045 gmx_mtop_atomlookup_destroy(alook);
1050 "Constraints: home %3d border %3d atoms: %3d\n",
1051 nhome, dc->ncon-nhome,
1052 dd->constraint_comm ? ireq->n : 0);
1056 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
1057 const gmx_mtop_t *mtop,
1059 gmx_constr_t constr, int nrec,
1062 gmx_domdec_constraints_t *dc;
1063 t_ilist *ilc_local, *ils_local;
1065 const t_blocka *at2con_mt;
1066 const int **at2settle_mt;
1067 gmx_hash_t ga2la_specat;
1071 dc = dd->constraints;
1073 ilc_local = &il_local[F_CONSTR];
1074 ils_local = &il_local[F_SETTLE];
1078 if (dd->constraint_comm)
1080 at2con_mt = atom2constraints_moltype(constr);
1081 ireq = &dd->constraint_comm->ireq[0];
1090 if (dd->bInterCGsettles)
1092 at2settle_mt = atom2settle_moltype(constr);
1097 /* Settle works inside charge groups, we assigned them already */
1098 at2settle_mt = NULL;
1101 if (at2settle_mt == NULL)
1103 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
1111 /* Do the constraints, if present, on the first thread.
1112 * Do the settles on all other threads.
1114 t0_set = ((at2con_mt != NULL && dc->nthread > 1) ? 1 : 0);
1116 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
1117 for (thread = 0; thread < dc->nthread; thread++)
1119 if (at2con_mt && thread == 0)
1121 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
1125 if (thread >= t0_set)
1131 /* Distribute the settle check+assignments over
1132 * dc->nthread or dc->nthread-1 threads.
1134 cg0 = (dd->ncg_home*(thread-t0_set ))/(dc->nthread-t0_set);
1135 cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);
1137 if (thread == t0_set)
1143 ilst = &dc->ils[thread];
1147 ireqt = &dd->constraint_comm->ireq[thread];
1153 atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
1159 /* Combine the generate settles and requested indices */
1160 for (thread = 1; thread < dc->nthread; thread++)
1166 if (thread > t0_set)
1168 ilst = &dc->ils[thread];
1169 if (ils_local->nr + ilst->nr > ils_local->nalloc)
1171 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
1172 srenew(ils_local->iatoms, ils_local->nalloc);
1174 for (ia = 0; ia < ilst->nr; ia++)
1176 ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
1178 ils_local->nr += ilst->nr;
1181 ireqt = &dd->constraint_comm->ireq[thread];
1182 if (ireq->n+ireqt->n > ireq->nalloc)
1184 ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
1185 srenew(ireq->ind, ireq->nalloc);
1187 for (ia = 0; ia < ireqt->n; ia++)
1189 ireq->ind[ireq->n+ia] = ireqt->ind[ia];
1191 ireq->n += ireqt->n;
1196 fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
1200 if (dd->constraint_comm)
1205 setup_specat_communication(dd, ireq, dd->constraint_comm,
1206 dd->constraints->ga2la,
1208 "constraint", " or lincs-order");
1210 /* Fill in the missing indices */
1211 ga2la_specat = dd->constraints->ga2la;
1213 nral1 = 1 + NRAL(F_CONSTR);
1214 for (i = 0; i < ilc_local->nr; i += nral1)
1216 iap = ilc_local->iatoms + i;
1217 for (j = 1; j < nral1; j++)
1221 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
1226 nral1 = 1 + NRAL(F_SETTLE);
1227 for (i = 0; i < ils_local->nr; i += nral1)
1229 iap = ils_local->iatoms + i;
1230 for (j = 1; j < nral1; j++)
1234 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
1247 int dd_make_local_vsites(gmx_domdec_t *dd, int at_start, t_ilist *lil)
1249 gmx_domdec_specat_comm_t *spac;
1251 gmx_hash_t ga2la_specat;
1252 int ftype, nral, i, j, gat, a;
1257 spac = dd->vsite_comm;
1258 ireq = &spac->ireq[0];
1259 ga2la_specat = dd->ga2la_vsite;
1262 /* Loop over all the home vsites */
1263 for (ftype = 0; ftype < F_NRE; ftype++)
1265 if (interaction_function[ftype].flags & IF_VSITE)
1269 for (i = 0; i < lilf->nr; i += 1+nral)
1271 iatoms = lilf->iatoms + i;
1272 /* Check if we have the other atoms */
1273 for (j = 1; j < 1+nral; j++)
1277 /* This is not a home atom,
1278 * we need to ask our neighbors.
1281 /* Check to not ask for the same atom more than once */
1282 if (gmx_hash_get_minone(dd->ga2la_vsite, a) == -1)
1284 /* Add this non-home atom to the list */
1285 if (ireq->n+1 > ireq->nalloc)
1287 ireq->nalloc = over_alloc_large(ireq->n+1);
1288 srenew(ireq->ind, ireq->nalloc);
1290 ireq->ind[ireq->n++] = a;
1291 /* Temporarily mark with -2,
1292 * we get the index later.
1294 gmx_hash_set(ga2la_specat, a, -2);
1302 at_end = setup_specat_communication(dd, ireq, dd->vsite_comm, ga2la_specat,
1303 at_start, 1, "vsite", "");
1305 /* Fill in the missing indices */
1306 for (ftype = 0; ftype < F_NRE; ftype++)
1308 if (interaction_function[ftype].flags & IF_VSITE)
1312 for (i = 0; i < lilf->nr; i += 1+nral)
1314 iatoms = lilf->iatoms + i;
1315 for (j = 1; j < 1+nral; j++)
1319 iatoms[j] = gmx_hash_get_minone(ga2la_specat, -iatoms[j]-1);
1329 static gmx_domdec_specat_comm_t *specat_comm_init(int nthread)
1331 gmx_domdec_specat_comm_t *spac;
1334 spac->nthread = nthread;
1335 snew(spac->ireq, spac->nthread);
1340 void init_domdec_constraints(gmx_domdec_t *dd,
1343 gmx_domdec_constraints_t *dc;
1344 gmx_molblock_t *molb;
1349 fprintf(debug, "Begin init_domdec_constraints\n");
1352 snew(dd->constraints, 1);
1353 dc = dd->constraints;
1355 snew(dc->molb_con_offset, mtop->nmolblock);
1356 snew(dc->molb_ncon_mol, mtop->nmolblock);
1359 for (mb = 0; mb < mtop->nmolblock; mb++)
1361 molb = &mtop->molblock[mb];
1362 dc->molb_con_offset[mb] = ncon;
1363 dc->molb_ncon_mol[mb] =
1364 mtop->moltype[molb->type].ilist[F_CONSTR].nr/3 +
1365 mtop->moltype[molb->type].ilist[F_CONSTRNC].nr/3;
1366 ncon += molb->nmol*dc->molb_ncon_mol[mb];
1371 snew(dc->gc_req, ncon);
1372 for (c = 0; c < ncon; c++)
1378 /* Use a hash table for the global to local index.
1379 * The number of keys is a rough estimate, it will be optimized later.
1381 dc->ga2la = gmx_hash_init(min(mtop->natoms/20,
1382 mtop->natoms/(2*dd->nnodes)));
1384 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
1385 snew(dc->ils, dc->nthread);
1387 dd->constraint_comm = specat_comm_init(dc->nthread);
1390 void init_domdec_vsites(gmx_domdec_t *dd, int n_intercg_vsite)
1393 gmx_domdec_constraints_t *dc;
1397 fprintf(debug, "Begin init_domdec_vsites\n");
1400 /* Use a hash table for the global to local index.
1401 * The number of keys is a rough estimate, it will be optimized later.
1403 dd->ga2la_vsite = gmx_hash_init(min(n_intercg_vsite/20,
1404 n_intercg_vsite/(2*dd->nnodes)));
1406 dd->vsite_comm = specat_comm_init(1);