2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "gromacs/utility/smalloc.h"
44 #include "types/commrec.h"
46 #include "domdec_network.h"
47 #include "mtop_util.h"
48 #include "gmx_ga2la.h"
50 #include "gmx_omp_nthreads.h"
53 #include "gmx_fatal.h"
68 typedef struct gmx_domdec_specat_comm {
69 /* The number of indices to receive during the setup */
71 /* The atoms to send */
72 gmx_specatsend_t spas[DIM][2];
82 /* The range in the local buffer(s) for received atoms */
86 /* The atom indices we need from the surrounding cells.
87 * We can gather the indices over nthread threads.
91 } gmx_domdec_specat_comm_t;
93 typedef struct gmx_domdec_constraints {
96 /* The fully local and connected constraints */
98 /* The global constraint number, only required for clearing gc_req */
102 /* Boolean that tells if a global constraint index has been requested */
104 /* Global to local communicated constraint atom only index */
107 /* Multi-threading stuff */
110 } gmx_domdec_constraints_t;
113 static void dd_move_f_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
114 rvec *f, rvec *fshift)
116 gmx_specatsend_t *spas;
118 int n, n0, n1, d, dim, dir, i;
121 gmx_bool bPBC, bScrew;
124 for (d = dd->ndim-1; d >= 0; d--)
129 /* Pulse the grid forward and backward */
130 spas = spac->spas[d];
135 /* Send and receive the coordinates */
136 dd_sendrecv2_rvec(dd, d,
137 f+n+n1, n0, vbuf, spas[0].nsend,
138 f+n, n1, vbuf+spas[0].nsend, spas[1].nsend);
139 for (dir = 0; dir < 2; dir++)
141 bPBC = ((dir == 0 && dd->ci[dim] == 0) ||
142 (dir == 1 && dd->ci[dim] == dd->nc[dim]-1));
143 bScrew = (bPBC && dd->bScrewPBC && dim == XX);
145 spas = &spac->spas[d][dir];
146 /* Sum the buffer into the required forces */
147 if (!bPBC || (!bScrew && fshift == NULL))
149 for (i = 0; i < spas->nsend; i++)
151 rvec_inc(f[spas->a[i]], *vbuf);
158 vis[dim] = (dir == 0 ? 1 : -1);
162 /* Sum and add to shift forces */
163 for (i = 0; i < spas->nsend; i++)
165 rvec_inc(f[spas->a[i]], *vbuf);
166 rvec_inc(fshift[is], *vbuf);
172 /* Rotate the forces */
173 for (i = 0; i < spas->nsend; i++)
175 f[spas->a[i]][XX] += (*vbuf)[XX];
176 f[spas->a[i]][YY] -= (*vbuf)[YY];
177 f[spas->a[i]][ZZ] -= (*vbuf)[ZZ];
180 rvec_inc(fshift[is], *vbuf);
190 /* Two cells, so we only need to communicate one way */
191 spas = &spac->spas[d][0];
193 /* Send and receive the coordinates */
194 dd_sendrecv_rvec(dd, d, dddirForward,
195 f+n, spas->nrecv, spac->vbuf, spas->nsend);
196 /* Sum the buffer into the required forces */
197 if (dd->bScrewPBC && dim == XX &&
199 dd->ci[dim] == dd->nc[dim]-1))
201 for (i = 0; i < spas->nsend; i++)
203 /* Rotate the force */
204 f[spas->a[i]][XX] += spac->vbuf[i][XX];
205 f[spas->a[i]][YY] -= spac->vbuf[i][YY];
206 f[spas->a[i]][ZZ] -= spac->vbuf[i][ZZ];
211 for (i = 0; i < spas->nsend; i++)
213 rvec_inc(f[spas->a[i]], spac->vbuf[i]);
220 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift)
224 dd_move_f_specat(dd, dd->vsite_comm, f, fshift);
228 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f)
234 for (i = dd->vsite_comm->at_start; i < dd->vsite_comm->at_end; i++)
241 static void dd_move_x_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
244 rvec *x1, gmx_bool bX1IsCoord)
246 gmx_specatsend_t *spas;
247 rvec *x, *vbuf, *rbuf;
248 int nvec, v, n, nn, ns0, ns1, nr0, nr1, nr, d, dim, dir, i;
249 gmx_bool bPBC, bScrew = FALSE;
250 rvec shift = {0, 0, 0};
259 for (d = 0; d < dd->ndim; d++)
264 /* Pulse the grid forward and backward */
266 for (dir = 0; dir < 2; dir++)
268 if (dir == 0 && dd->ci[dim] == 0)
271 bScrew = (dd->bScrewPBC && dim == XX);
272 copy_rvec(box[dim], shift);
274 else if (dir == 1 && dd->ci[dim] == dd->nc[dim]-1)
277 bScrew = (dd->bScrewPBC && dim == XX);
278 for (i = 0; i < DIM; i++)
280 shift[i] = -box[dim][i];
288 spas = &spac->spas[d][dir];
289 for (v = 0; v < nvec; v++)
291 x = (v == 0 ? x0 : x1);
292 /* Copy the required coordinates to the send buffer */
293 if (!bPBC || (v == 1 && !bX1IsCoord))
296 for (i = 0; i < spas->nsend; i++)
298 copy_rvec(x[spas->a[i]], *vbuf);
304 /* Shift coordinates */
305 for (i = 0; i < spas->nsend; i++)
307 rvec_add(x[spas->a[i]], shift, *vbuf);
313 /* Shift and rotate coordinates */
314 for (i = 0; i < spas->nsend; i++)
316 (*vbuf)[XX] = x[spas->a[i]][XX] + shift[XX];
317 (*vbuf)[YY] = box[YY][YY] - x[spas->a[i]][YY] + shift[YY];
318 (*vbuf)[ZZ] = box[ZZ][ZZ] - x[spas->a[i]][ZZ] + shift[ZZ];
324 /* Send and receive the coordinates */
325 spas = spac->spas[d];
332 dd_sendrecv2_rvec(dd, d,
333 spac->vbuf+ns0, ns1, x0+n, nr1,
334 spac->vbuf, ns0, x0+n+nr1, nr0);
338 /* Communicate both vectors in one buffer */
340 dd_sendrecv2_rvec(dd, d,
341 spac->vbuf+2*ns0, 2*ns1, rbuf, 2*nr1,
342 spac->vbuf, 2*ns0, rbuf+2*nr1, 2*nr0);
343 /* Split the buffer into the two vectors */
345 for (dir = 1; dir >= 0; dir--)
347 nr = spas[dir].nrecv;
348 for (v = 0; v < 2; v++)
350 x = (v == 0 ? x0 : x1);
351 for (i = 0; i < nr; i++)
353 copy_rvec(*rbuf, x[nn+i]);
364 spas = &spac->spas[d][0];
365 /* Copy the required coordinates to the send buffer */
367 for (v = 0; v < nvec; v++)
369 x = (v == 0 ? x0 : x1);
370 if (dd->bScrewPBC && dim == XX &&
371 (dd->ci[XX] == 0 || dd->ci[XX] == dd->nc[XX]-1))
373 /* Here we only perform the rotation, the rest of the pbc
374 * is handled in the constraint or viste routines.
376 for (i = 0; i < spas->nsend; i++)
378 (*vbuf)[XX] = x[spas->a[i]][XX];
379 (*vbuf)[YY] = box[YY][YY] - x[spas->a[i]][YY];
380 (*vbuf)[ZZ] = box[ZZ][ZZ] - x[spas->a[i]][ZZ];
386 for (i = 0; i < spas->nsend; i++)
388 copy_rvec(x[spas->a[i]], *vbuf);
393 /* Send and receive the coordinates */
396 dd_sendrecv_rvec(dd, d, dddirBackward,
397 spac->vbuf, spas->nsend, x0+n, spas->nrecv);
401 /* Communicate both vectors in one buffer */
403 dd_sendrecv_rvec(dd, d, dddirBackward,
404 spac->vbuf, 2*spas->nsend, rbuf, 2*spas->nrecv);
405 /* Split the buffer into the two vectors */
407 for (v = 0; v < 2; v++)
409 x = (v == 0 ? x0 : x1);
410 for (i = 0; i < nr; i++)
412 copy_rvec(*rbuf, x[n+i]);
422 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
423 rvec *x0, rvec *x1, gmx_bool bX1IsCoord)
425 if (dd->constraint_comm)
427 dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1, bX1IsCoord);
431 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x)
435 dd_move_x_specat(dd, dd->vsite_comm, box, x, NULL, FALSE);
439 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd)
443 return dd->constraints->con_nlocat;
451 void dd_clear_local_constraint_indices(gmx_domdec_t *dd)
453 gmx_domdec_constraints_t *dc;
456 dc = dd->constraints;
458 for (i = 0; i < dc->ncon; i++)
460 dc->gc_req[dc->con_gl[i]] = 0;
463 if (dd->constraint_comm)
465 gmx_hash_clear_and_optimize(dc->ga2la);
469 void dd_clear_local_vsite_indices(gmx_domdec_t *dd)
475 gmx_hash_clear_and_optimize(dd->ga2la_vsite);
479 static int setup_specat_communication(gmx_domdec_t *dd,
481 gmx_domdec_specat_comm_t *spac,
482 gmx_hash_t ga2la_specat,
485 const char *specat_type,
488 int nsend[2], nlast, nsend_zero[2] = {0, 0}, *nsend_ptr;
489 int d, dim, ndir, dir, nr, ns, i, nrecv_local, n0, start, indr, ind, buf[2];
490 int nat_tot_specat, nat_tot_prev, nalloc_old;
491 gmx_bool bPBC, bFirst;
492 gmx_specatsend_t *spas;
496 fprintf(debug, "Begin setup_specat_communication for %s\n", specat_type);
499 /* nsend[0]: the number of atoms requested by this node only,
500 * we communicate this for more efficients checks
501 * nsend[1]: the total number of requested atoms
506 for (d = dd->ndim-1; d >= 0; d--)
508 /* Pulse the grid forward and backward */
510 bPBC = (dim < dd->npbcdim);
511 if (dd->nc[dim] == 2)
513 /* Only 2 cells, so we only need to communicate once */
520 for (dir = 0; dir < ndir; dir++)
524 ((dir == 0 && dd->ci[dim] == dd->nc[dim] - 1) ||
525 (dir == 1 && dd->ci[dim] == 0)))
527 /* No pbc: the fist/last cell should not request atoms */
528 nsend_ptr = nsend_zero;
534 /* Communicate the number of indices */
535 dd_sendrecv_int(dd, d, dir == 0 ? dddirForward : dddirBackward,
536 nsend_ptr, 2, spac->nreq[d][dir], 2);
537 nr = spac->nreq[d][dir][1];
538 if (nlast+nr > ireq->nalloc)
540 ireq->nalloc = over_alloc_dd(nlast+nr);
541 srenew(ireq->ind, ireq->nalloc);
543 /* Communicate the indices */
544 dd_sendrecv_int(dd, d, dir == 0 ? dddirForward : dddirBackward,
545 ireq->ind, nsend_ptr[1], ireq->ind+nlast, nr);
552 fprintf(debug, "Communicated the counts\n");
555 /* Search for the requested atoms and communicate the indices we have */
556 nat_tot_specat = at_start;
558 for (d = 0; d < dd->ndim; d++)
561 /* Pulse the grid forward and backward */
562 if (dd->dim[d] >= dd->npbcdim || dd->nc[dd->dim[d]] > 2)
570 nat_tot_prev = nat_tot_specat;
571 for (dir = ndir-1; dir >= 0; dir--)
573 if (nat_tot_specat > spac->bSendAtom_nalloc)
575 nalloc_old = spac->bSendAtom_nalloc;
576 spac->bSendAtom_nalloc = over_alloc_dd(nat_tot_specat);
577 srenew(spac->bSendAtom, spac->bSendAtom_nalloc);
578 for (i = nalloc_old; i < spac->bSendAtom_nalloc; i++)
580 spac->bSendAtom[i] = FALSE;
583 spas = &spac->spas[d][dir];
584 n0 = spac->nreq[d][dir][0];
585 nr = spac->nreq[d][dir][1];
588 fprintf(debug, "dim=%d, dir=%d, searching for %d atoms\n",
594 for (i = 0; i < nr; i++)
596 indr = ireq->ind[start+i];
598 /* Check if this is a home atom and if so ind will be set */
599 if (!ga2la_get_home(dd->ga2la, indr, &ind))
601 /* Search in the communicated atoms */
602 ind = gmx_hash_get_minone(ga2la_specat, indr);
606 if (i < n0 || !spac->bSendAtom[ind])
608 if (spas->nsend+1 > spas->a_nalloc)
610 spas->a_nalloc = over_alloc_large(spas->nsend+1);
611 srenew(spas->a, spas->a_nalloc);
613 /* Store the local index so we know which coordinates
616 spas->a[spas->nsend] = ind;
617 spac->bSendAtom[ind] = TRUE;
618 if (spas->nsend+1 > spac->ibuf_nalloc)
620 spac->ibuf_nalloc = over_alloc_large(spas->nsend+1);
621 srenew(spac->ibuf, spac->ibuf_nalloc);
623 /* Store the global index so we can send it now */
624 spac->ibuf[spas->nsend] = indr;
634 /* Clear the local flags */
635 for (i = 0; i < spas->nsend; i++)
637 spac->bSendAtom[spas->a[i]] = FALSE;
639 /* Send and receive the number of indices to communicate */
640 nsend[1] = spas->nsend;
641 dd_sendrecv_int(dd, d, dir == 0 ? dddirBackward : dddirForward,
645 fprintf(debug, "Send to rank %d, %d (%d) indices, "
646 "receive from rank %d, %d (%d) indices\n",
647 dd->neighbor[d][1-dir], nsend[1], nsend[0],
648 dd->neighbor[d][dir], buf[1], buf[0]);
651 for (i = 0; i < spas->nsend; i++)
653 fprintf(debug, " %d", spac->ibuf[i]+1);
655 fprintf(debug, "\n");
658 nrecv_local += buf[0];
659 spas->nrecv = buf[1];
660 if (nat_tot_specat + spas->nrecv > dd->gatindex_nalloc)
662 dd->gatindex_nalloc =
663 over_alloc_dd(nat_tot_specat + spas->nrecv);
664 srenew(dd->gatindex, dd->gatindex_nalloc);
666 /* Send and receive the indices */
667 dd_sendrecv_int(dd, d, dir == 0 ? dddirBackward : dddirForward,
668 spac->ibuf, spas->nsend,
669 dd->gatindex+nat_tot_specat, spas->nrecv);
670 nat_tot_specat += spas->nrecv;
673 /* Allocate the x/f communication buffers */
674 ns = spac->spas[d][0].nsend;
675 nr = spac->spas[d][0].nrecv;
678 ns += spac->spas[d][1].nsend;
679 nr += spac->spas[d][1].nrecv;
681 if (vbuf_fac*ns > spac->vbuf_nalloc)
683 spac->vbuf_nalloc = over_alloc_dd(vbuf_fac*ns);
684 srenew(spac->vbuf, spac->vbuf_nalloc);
686 if (vbuf_fac == 2 && vbuf_fac*nr > spac->vbuf2_nalloc)
688 spac->vbuf2_nalloc = over_alloc_dd(vbuf_fac*nr);
689 srenew(spac->vbuf2, spac->vbuf2_nalloc);
692 /* Make a global to local index for the communication atoms */
693 for (i = nat_tot_prev; i < nat_tot_specat; i++)
695 gmx_hash_change_or_set(ga2la_specat, dd->gatindex[i], i);
699 /* Check that in the end we got the number of atoms we asked for */
700 if (nrecv_local != ireq->n)
704 fprintf(debug, "Requested %d, received %d (tot recv %d)\n",
705 ireq->n, nrecv_local, nat_tot_specat-at_start);
708 for (i = 0; i < ireq->n; i++)
710 ind = gmx_hash_get_minone(ga2la_specat, ireq->ind[i]);
711 fprintf(debug, " %s%d",
712 (ind >= 0) ? "" : "!",
715 fprintf(debug, "\n");
718 fprintf(stderr, "\nDD cell %d %d %d: Neighboring cells do not have atoms:",
719 dd->ci[XX], dd->ci[YY], dd->ci[ZZ]);
720 for (i = 0; i < ireq->n; i++)
722 if (gmx_hash_get_minone(ga2la_specat, ireq->ind[i]) < 0)
724 fprintf(stderr, " %d", ireq->ind[i]+1);
727 fprintf(stderr, "\n");
728 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.",
729 dd->ci[XX], dd->ci[YY], dd->ci[ZZ],
730 nrecv_local, ireq->n, specat_type,
731 specat_type, add_err,
732 dd->bGridJump ? " or use the -rcon option of mdrun" : "");
735 spac->at_start = at_start;
736 spac->at_end = nat_tot_specat;
740 fprintf(debug, "Done setup_specat_communication\n");
743 return nat_tot_specat;
746 static void walk_out(int con, int con_offset, int a, int offset, int nrec,
747 int ncon1, const t_iatom *ia1, const t_iatom *ia2,
748 const t_blocka *at2con,
749 const gmx_ga2la_t ga2la, gmx_bool bHomeConnect,
750 gmx_domdec_constraints_t *dc,
751 gmx_domdec_specat_comm_t *dcc,
755 int a1_gl, a2_gl, a_loc, i, coni, b;
758 if (dc->gc_req[con_offset+con] == 0)
760 /* Add this non-home constraint to the list */
761 if (dc->ncon+1 > dc->con_nalloc)
763 dc->con_nalloc = over_alloc_large(dc->ncon+1);
764 srenew(dc->con_gl, dc->con_nalloc);
765 srenew(dc->con_nlocat, dc->con_nalloc);
767 dc->con_gl[dc->ncon] = con_offset + con;
768 dc->con_nlocat[dc->ncon] = (bHomeConnect ? 1 : 0);
769 dc->gc_req[con_offset+con] = 1;
770 if (il_local->nr + 3 > il_local->nalloc)
772 il_local->nalloc = over_alloc_dd(il_local->nr+3);
773 srenew(il_local->iatoms, il_local->nalloc);
775 iap = constr_iatomptr(ncon1, ia1, ia2, con);
776 il_local->iatoms[il_local->nr++] = iap[0];
777 a1_gl = offset + iap[1];
778 a2_gl = offset + iap[2];
779 /* The following indexing code can probably be optizimed */
780 if (ga2la_get_home(ga2la, a1_gl, &a_loc))
782 il_local->iatoms[il_local->nr++] = a_loc;
786 /* We set this index later */
787 il_local->iatoms[il_local->nr++] = -a1_gl - 1;
789 if (ga2la_get_home(ga2la, a2_gl, &a_loc))
791 il_local->iatoms[il_local->nr++] = a_loc;
795 /* We set this index later */
796 il_local->iatoms[il_local->nr++] = -a2_gl - 1;
800 /* Check to not ask for the same atom more than once */
801 if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
804 /* Add this non-home atom to the list */
805 if (ireq->n+1 > ireq->nalloc)
807 ireq->nalloc = over_alloc_large(ireq->n+1);
808 srenew(ireq->ind, ireq->nalloc);
810 ireq->ind[ireq->n++] = offset + a;
811 /* Temporarily mark with -2, we get the index later */
812 gmx_hash_set(dc->ga2la, offset+a, -2);
817 for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
823 iap = constr_iatomptr(ncon1, ia1, ia2, coni);
832 if (!ga2la_get_home(ga2la, offset+b, &a_loc))
834 walk_out(coni, con_offset, b, offset, nrec-1,
835 ncon1, ia1, ia2, at2con,
836 ga2la, FALSE, dc, dcc, il_local, ireq);
843 static void atoms_to_settles(gmx_domdec_t *dd,
844 const gmx_mtop_t *mtop,
846 const int **at2settle_mt,
847 int cg_start, int cg_end,
852 gmx_mtop_atomlookup_t alook;
855 int cg, a, a_gl, a_glsa, a_gls[3], a_locs[3];
856 int mb, molnr, a_mol, offset;
857 const gmx_molblock_t *molb;
865 alook = gmx_mtop_atomlookup_settle_init(mtop);
867 nral = NRAL(F_SETTLE);
869 for (cg = cg_start; cg < cg_end; cg++)
871 if (GET_CGINFO_SETTLE(cginfo[cg]))
873 for (a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
875 a_gl = dd->gatindex[a];
877 gmx_mtop_atomnr_to_molblock_ind(alook, a_gl, &mb, &molnr, &a_mol);
878 molb = &mtop->molblock[mb];
880 settle = at2settle_mt[molb->type][a_mol];
884 offset = a_gl - a_mol;
886 ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms;
890 for (sa = 0; sa < nral; sa++)
892 a_glsa = offset + ia1[settle*(1+nral)+1+sa];
894 a_home[sa] = ga2la_get_home(ga2la, a_glsa, &a_locs[sa]);
897 if (nlocal == 0 && a_gl == a_glsa)
907 if (ils_local->nr+1+nral > ils_local->nalloc)
909 ils_local->nalloc = over_alloc_dd(ils_local->nr+1+nral);
910 srenew(ils_local->iatoms, ils_local->nalloc);
913 ils_local->iatoms[ils_local->nr++] = ia1[settle*4];
915 for (sa = 0; sa < nral; sa++)
917 if (ga2la_get_home(ga2la, a_gls[sa], &a_locs[sa]))
919 ils_local->iatoms[ils_local->nr++] = a_locs[sa];
923 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
924 /* Add this non-home atom to the list */
925 if (ireq->n+1 > ireq->nalloc)
927 ireq->nalloc = over_alloc_large(ireq->n+1);
928 srenew(ireq->ind, ireq->nalloc);
930 ireq->ind[ireq->n++] = a_gls[sa];
931 /* A check on double atom requests is
932 * not required for settle.
942 gmx_mtop_atomlookup_destroy(alook);
945 static void atoms_to_constraints(gmx_domdec_t *dd,
946 const gmx_mtop_t *mtop,
948 const t_blocka *at2con_mt, int nrec,
952 const t_blocka *at2con;
954 gmx_mtop_atomlookup_t alook;
956 gmx_molblock_t *molb;
957 t_iatom *ia1, *ia2, *iap;
958 int nhome, cg, a, a_gl, a_mol, a_loc, b_lo, offset, mb, molnr, b_mol, i, con, con_offset;
959 gmx_domdec_constraints_t *dc;
960 gmx_domdec_specat_comm_t *dcc;
962 dc = dd->constraints;
963 dcc = dd->constraint_comm;
967 alook = gmx_mtop_atomlookup_init(mtop);
970 for (cg = 0; cg < dd->ncg_home; cg++)
972 if (GET_CGINFO_CONSTR(cginfo[cg]))
974 for (a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
976 a_gl = dd->gatindex[a];
978 gmx_mtop_atomnr_to_molblock_ind(alook, a_gl, &mb, &molnr, &a_mol);
979 molb = &mtop->molblock[mb];
981 ncon1 = mtop->moltype[molb->type].ilist[F_CONSTR].nr/NRAL(F_SETTLE);
983 ia1 = mtop->moltype[molb->type].ilist[F_CONSTR].iatoms;
984 ia2 = mtop->moltype[molb->type].ilist[F_CONSTRNC].iatoms;
986 /* Calculate the global constraint number offset for the molecule.
987 * This is only required for the global index to make sure
988 * that we use each constraint only once.
991 dc->molb_con_offset[mb] + molnr*dc->molb_ncon_mol[mb];
993 /* The global atom number offset for this molecule */
994 offset = a_gl - a_mol;
995 at2con = &at2con_mt[molb->type];
996 for (i = at2con->index[a_mol]; i < at2con->index[a_mol+1]; i++)
999 iap = constr_iatomptr(ncon1, ia1, ia2, con);
1000 if (a_mol == iap[1])
1008 if (ga2la_get_home(ga2la, offset+b_mol, &a_loc))
1010 /* Add this fully home constraint at the first atom */
1013 if (dc->ncon+1 > dc->con_nalloc)
1015 dc->con_nalloc = over_alloc_large(dc->ncon+1);
1016 srenew(dc->con_gl, dc->con_nalloc);
1017 srenew(dc->con_nlocat, dc->con_nalloc);
1019 dc->con_gl[dc->ncon] = con_offset + con;
1020 dc->con_nlocat[dc->ncon] = 2;
1021 if (ilc_local->nr + 3 > ilc_local->nalloc)
1023 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
1024 srenew(ilc_local->iatoms, ilc_local->nalloc);
1027 ilc_local->iatoms[ilc_local->nr++] = iap[0];
1028 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
1029 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a );
1036 /* We need the nrec constraints coupled to this constraint,
1037 * so we need to walk out of the home cell by nrec+1 atoms,
1038 * since already atom bg is not locally present.
1039 * Therefore we call walk_out with nrec recursions to go
1040 * after this first call.
1042 walk_out(con, con_offset, b_mol, offset, nrec,
1043 ncon1, ia1, ia2, at2con,
1044 dd->ga2la, TRUE, dc, dcc, ilc_local, ireq);
1051 gmx_mtop_atomlookup_destroy(alook);
1056 "Constraints: home %3d border %3d atoms: %3d\n",
1057 nhome, dc->ncon-nhome,
1058 dd->constraint_comm ? ireq->n : 0);
1062 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
1063 const gmx_mtop_t *mtop,
1065 gmx_constr_t constr, int nrec,
1068 gmx_domdec_constraints_t *dc;
1069 t_ilist *ilc_local, *ils_local;
1071 const t_blocka *at2con_mt;
1072 const int **at2settle_mt;
1073 gmx_hash_t ga2la_specat;
1077 dc = dd->constraints;
1079 ilc_local = &il_local[F_CONSTR];
1080 ils_local = &il_local[F_SETTLE];
1084 if (dd->constraint_comm)
1086 at2con_mt = atom2constraints_moltype(constr);
1087 ireq = &dd->constraint_comm->ireq[0];
1096 if (dd->bInterCGsettles)
1098 at2settle_mt = atom2settle_moltype(constr);
1103 /* Settle works inside charge groups, we assigned them already */
1104 at2settle_mt = NULL;
1107 if (at2settle_mt == NULL)
1109 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
1117 /* Do the constraints, if present, on the first thread.
1118 * Do the settles on all other threads.
1120 t0_set = ((at2con_mt != NULL && dc->nthread > 1) ? 1 : 0);
1122 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
1123 for (thread = 0; thread < dc->nthread; thread++)
1125 if (at2con_mt && thread == 0)
1127 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
1131 if (thread >= t0_set)
1137 /* Distribute the settle check+assignments over
1138 * dc->nthread or dc->nthread-1 threads.
1140 cg0 = (dd->ncg_home*(thread-t0_set ))/(dc->nthread-t0_set);
1141 cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);
1143 if (thread == t0_set)
1149 ilst = &dc->ils[thread];
1153 ireqt = &dd->constraint_comm->ireq[thread];
1159 atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
1165 /* Combine the generate settles and requested indices */
1166 for (thread = 1; thread < dc->nthread; thread++)
1172 if (thread > t0_set)
1174 ilst = &dc->ils[thread];
1175 if (ils_local->nr + ilst->nr > ils_local->nalloc)
1177 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
1178 srenew(ils_local->iatoms, ils_local->nalloc);
1180 for (ia = 0; ia < ilst->nr; ia++)
1182 ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
1184 ils_local->nr += ilst->nr;
1187 ireqt = &dd->constraint_comm->ireq[thread];
1188 if (ireq->n+ireqt->n > ireq->nalloc)
1190 ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
1191 srenew(ireq->ind, ireq->nalloc);
1193 for (ia = 0; ia < ireqt->n; ia++)
1195 ireq->ind[ireq->n+ia] = ireqt->ind[ia];
1197 ireq->n += ireqt->n;
1202 fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
1206 if (dd->constraint_comm)
1211 setup_specat_communication(dd, ireq, dd->constraint_comm,
1212 dd->constraints->ga2la,
1214 "constraint", " or lincs-order");
1216 /* Fill in the missing indices */
1217 ga2la_specat = dd->constraints->ga2la;
1219 nral1 = 1 + NRAL(F_CONSTR);
1220 for (i = 0; i < ilc_local->nr; i += nral1)
1222 iap = ilc_local->iatoms + i;
1223 for (j = 1; j < nral1; j++)
1227 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
1232 nral1 = 1 + NRAL(F_SETTLE);
1233 for (i = 0; i < ils_local->nr; i += nral1)
1235 iap = ils_local->iatoms + i;
1236 for (j = 1; j < nral1; j++)
1240 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
1253 int dd_make_local_vsites(gmx_domdec_t *dd, int at_start, t_ilist *lil)
1255 gmx_domdec_specat_comm_t *spac;
1257 gmx_hash_t ga2la_specat;
1258 int ftype, nral, i, j, gat, a;
1263 spac = dd->vsite_comm;
1264 ireq = &spac->ireq[0];
1265 ga2la_specat = dd->ga2la_vsite;
1268 /* Loop over all the home vsites */
1269 for (ftype = 0; ftype < F_NRE; ftype++)
1271 if (interaction_function[ftype].flags & IF_VSITE)
1275 for (i = 0; i < lilf->nr; i += 1+nral)
1277 iatoms = lilf->iatoms + i;
1278 /* Check if we have the other atoms */
1279 for (j = 1; j < 1+nral; j++)
1283 /* This is not a home atom,
1284 * we need to ask our neighbors.
1287 /* Check to not ask for the same atom more than once */
1288 if (gmx_hash_get_minone(dd->ga2la_vsite, a) == -1)
1290 /* Add this non-home atom to the list */
1291 if (ireq->n+1 > ireq->nalloc)
1293 ireq->nalloc = over_alloc_large(ireq->n+1);
1294 srenew(ireq->ind, ireq->nalloc);
1296 ireq->ind[ireq->n++] = a;
1297 /* Temporarily mark with -2,
1298 * we get the index later.
1300 gmx_hash_set(ga2la_specat, a, -2);
1308 at_end = setup_specat_communication(dd, ireq, dd->vsite_comm, ga2la_specat,
1309 at_start, 1, "vsite", "");
1311 /* Fill in the missing indices */
1312 for (ftype = 0; ftype < F_NRE; ftype++)
1314 if (interaction_function[ftype].flags & IF_VSITE)
1318 for (i = 0; i < lilf->nr; i += 1+nral)
1320 iatoms = lilf->iatoms + i;
1321 for (j = 1; j < 1+nral; j++)
1325 iatoms[j] = gmx_hash_get_minone(ga2la_specat, -iatoms[j]-1);
1335 static gmx_domdec_specat_comm_t *specat_comm_init(int nthread)
1337 gmx_domdec_specat_comm_t *spac;
1340 spac->nthread = nthread;
1341 snew(spac->ireq, spac->nthread);
1346 void init_domdec_constraints(gmx_domdec_t *dd,
1349 gmx_domdec_constraints_t *dc;
1350 gmx_molblock_t *molb;
1355 fprintf(debug, "Begin init_domdec_constraints\n");
1358 snew(dd->constraints, 1);
1359 dc = dd->constraints;
1361 snew(dc->molb_con_offset, mtop->nmolblock);
1362 snew(dc->molb_ncon_mol, mtop->nmolblock);
1365 for (mb = 0; mb < mtop->nmolblock; mb++)
1367 molb = &mtop->molblock[mb];
1368 dc->molb_con_offset[mb] = ncon;
1369 dc->molb_ncon_mol[mb] =
1370 mtop->moltype[molb->type].ilist[F_CONSTR].nr/3 +
1371 mtop->moltype[molb->type].ilist[F_CONSTRNC].nr/3;
1372 ncon += molb->nmol*dc->molb_ncon_mol[mb];
1377 snew(dc->gc_req, ncon);
1378 for (c = 0; c < ncon; c++)
1384 /* Use a hash table for the global to local index.
1385 * The number of keys is a rough estimate, it will be optimized later.
1387 dc->ga2la = gmx_hash_init(min(mtop->natoms/20,
1388 mtop->natoms/(2*dd->nnodes)));
1390 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
1391 snew(dc->ils, dc->nthread);
1393 dd->constraint_comm = specat_comm_init(dc->nthread);
1396 void init_domdec_vsites(gmx_domdec_t *dd, int n_intercg_vsite)
1399 gmx_domdec_constraints_t *dc;
1403 fprintf(debug, "Begin init_domdec_vsites\n");
1406 /* Use a hash table for the global to local index.
1407 * The number of keys is a rough estimate, it will be optimized later.
1409 dd->ga2la_vsite = gmx_hash_init(min(n_intercg_vsite/20,
1410 n_intercg_vsite/(2*dd->nnodes)));
1412 dd->vsite_comm = specat_comm_init(1);