int setup_specat_communication(gmx_domdec_t               *dd,
                               ind_req_t                  *ireq,
                               gmx_domdec_specat_comm_t   *spac,
                               gmx_hash_t                 *ga2la_specat,
                               int                         at_start,
                               int                         vbuf_fac,
                               const char                 *specat_type,
                               const char                 *add_err)
{
    int               nsend[2], nlast, nsend_zero[2] = {0, 0}, *nsend_ptr;
    int               d, dim, ndir, dir, nr, ns, i, nrecv_local, n0, start, indr, ind, buf[2];
    int               nat_tot_specat, nat_tot_prev, nalloc_old;
    gmx_bool          bPBC;
    gmx_specatsend_t *spas;

    if (debug)
    {
        fprintf(debug, "Begin setup_specat_communication for %s\n", specat_type);
    }

    /* nsend[0]: the number of atoms requested by this node only,
     *           we communicate this for more efficients checks
     * nsend[1]: the total number of requested atoms
     */
    nsend[0] = ireq->n;
    nsend[1] = nsend[0];
    nlast    = nsend[1];
    for (d = dd->ndim-1; d >= 0; d--)
    {
        /* Pulse the grid forward and backward */
        dim  = dd->dim[d];
        bPBC = (dim < dd->npbcdim);
        if (dd->nc[dim] == 2)
        {
            /* Only 2 cells, so we only need to communicate once */
            ndir = 1;
        }
        else
        {
            ndir = 2;
        }
        for (dir = 0; dir < ndir; dir++)
        {
            if (!bPBC &&
                dd->nc[dim] > 2 &&
                ((dir == 0 && dd->ci[dim] == dd->nc[dim] - 1) ||
                 (dir == 1 && dd->ci[dim] == 0)))
            {
                /* No pbc: the fist/last cell should not request atoms */
                nsend_ptr = nsend_zero;
            }
            else
            {
                nsend_ptr = nsend;
            }
            /* Communicate the number of indices */
            dd_sendrecv_int(dd, d, dir == 0 ? dddirForward : dddirBackward,
                            nsend_ptr, 2, spac->nreq[d][dir], 2);
            nr = spac->nreq[d][dir][1];
            if (nlast+nr > ireq->nalloc)
            {
                ireq->nalloc = over_alloc_dd(nlast+nr);
                srenew(ireq->ind, ireq->nalloc);
            }
            /* Communicate the indices */
            dd_sendrecv_int(dd, d, dir == 0 ? dddirForward : dddirBackward,
                            ireq->ind, nsend_ptr[1], ireq->ind+nlast, nr);
            nlast += nr;
        }
        nsend[1] = nlast;
    }
    if (debug)
    {
        fprintf(debug, "Communicated the counts\n");
    }

    /* Search for the requested atoms and communicate the indices we have */
    nat_tot_specat = at_start;
    nrecv_local    = 0;
    for (d = 0; d < dd->ndim; d++)
    {
        /* Pulse the grid forward and backward */
        if (dd->dim[d] >= dd->npbcdim || dd->nc[dd->dim[d]] > 2)
        {
            ndir = 2;
        }
        else
        {
            ndir = 1;
        }
        nat_tot_prev = nat_tot_specat;
        for (dir = ndir-1; dir >= 0; dir--)
        {
            if (nat_tot_specat > spac->bSendAtom_nalloc)
            {
                nalloc_old             = spac->bSendAtom_nalloc;
                spac->bSendAtom_nalloc = over_alloc_dd(nat_tot_specat);
                srenew(spac->bSendAtom, spac->bSendAtom_nalloc);
                for (i = nalloc_old; i < spac->bSendAtom_nalloc; i++)
                {
                    spac->bSendAtom[i] = FALSE;
                }
            }
            spas = &spac->spas[d][dir];
            n0   = spac->nreq[d][dir][0];
            nr   = spac->nreq[d][dir][1];
            if (debug)
            {
                fprintf(debug, "dim=%d, dir=%d, searching for %d atoms\n",
                        d, dir, nr);
            }
            start       = nlast - nr;
            spas->nsend = 0;
            nsend[0]    = 0;
            for (i = 0; i < nr; i++)
            {
                indr = ireq->ind[start+i];
                ind  = -1;
                /* Check if this is a home atom and if so ind will be set */
                if (!ga2la_get_home(dd->ga2la, indr, &ind))
                {
                    /* Search in the communicated atoms */
                    ind = gmx_hash_get_minone(ga2la_specat, indr);
                }
                if (ind >= 0)
                {
                    if (i < n0 || !spac->bSendAtom[ind])
                    {
                        if (spas->nsend+1 > spas->a_nalloc)
                        {
                            spas->a_nalloc = over_alloc_large(spas->nsend+1);
                            srenew(spas->a, spas->a_nalloc);
                        }
                        /* Store the local index so we know which coordinates
                         * to send out later.
                         */
                        spas->a[spas->nsend] = ind;
                        spac->bSendAtom[ind] = TRUE;
                        if (spas->nsend+1 > spac->ibuf_nalloc)
                        {
                            spac->ibuf_nalloc = over_alloc_large(spas->nsend+1);
                            srenew(spac->ibuf, spac->ibuf_nalloc);
                        }
                        /* Store the global index so we can send it now */
                        spac->ibuf[spas->nsend] = indr;
                        if (i < n0)
                        {
                            nsend[0]++;
                        }
                        spas->nsend++;
                    }
                }
            }
            nlast = start;
            /* Clear the local flags */
            for (i = 0; i < spas->nsend; i++)
            {
                spac->bSendAtom[spas->a[i]] = FALSE;
            }
            /* Send and receive the number of indices to communicate */
            nsend[1] = spas->nsend;
            dd_sendrecv_int(dd, d, dir == 0 ? dddirBackward : dddirForward,
                            nsend, 2, buf, 2);
            if (debug)
            {
                fprintf(debug, "Send to rank %d, %d (%d) indices, "
                        "receive from rank %d, %d (%d) indices\n",
                        dd->neighbor[d][1-dir], nsend[1], nsend[0],
                        dd->neighbor[d][dir], buf[1], buf[0]);
                if (gmx_debug_at)
                {
                    for (i = 0; i < spas->nsend; i++)
                    {
                        fprintf(debug, " %d", spac->ibuf[i]+1);
                    }
                    fprintf(debug, "\n");
                }
            }
            nrecv_local += buf[0];
            spas->nrecv  = buf[1];
            if (nat_tot_specat + spas->nrecv > dd->gatindex_nalloc)
            {
                dd->gatindex_nalloc =
                    over_alloc_dd(nat_tot_specat + spas->nrecv);
                srenew(dd->gatindex, dd->gatindex_nalloc);
            }
            /* Send and receive the indices */
            dd_sendrecv_int(dd, d, dir == 0 ? dddirBackward : dddirForward,
                            spac->ibuf, spas->nsend,
                            dd->gatindex+nat_tot_specat, spas->nrecv);
            nat_tot_specat += spas->nrecv;
        }

        /* Allocate the x/f communication buffers */
        ns = spac->spas[d][0].nsend;
        nr = spac->spas[d][0].nrecv;
        if (ndir == 2)
        {
            ns += spac->spas[d][1].nsend;
            nr += spac->spas[d][1].nrecv;
        }
        if (vbuf_fac*ns > spac->vbuf_nalloc)
        {
            spac->vbuf_nalloc = over_alloc_dd(vbuf_fac*ns);
            srenew(spac->vbuf, spac->vbuf_nalloc);
        }
        if (vbuf_fac == 2 && vbuf_fac*nr > spac->vbuf2_nalloc)
        {
            spac->vbuf2_nalloc = over_alloc_dd(vbuf_fac*nr);
            srenew(spac->vbuf2, spac->vbuf2_nalloc);
        }

        /* Make a global to local index for the communication atoms */
        for (i = nat_tot_prev; i < nat_tot_specat; i++)
        {
            gmx_hash_change_or_set(ga2la_specat, dd->gatindex[i], i);
        }
    }

    /* Check that in the end we got the number of atoms we asked for */
    if (nrecv_local != ireq->n)
    {
        if (debug)
        {
            fprintf(debug, "Requested %d, received %d (tot recv %d)\n",
                    ireq->n, nrecv_local, nat_tot_specat-at_start);
            if (gmx_debug_at)
            {
                for (i = 0; i < ireq->n; i++)
                {
                    ind = gmx_hash_get_minone(ga2la_specat, ireq->ind[i]);
                    fprintf(debug, " %s%d",
                            (ind >= 0) ? "" : "!",
                            ireq->ind[i]+1);
                }
                fprintf(debug, "\n");
            }
        }
        fprintf(stderr, "\nDD cell %d %d %d: Neighboring cells do not have atoms:",
                dd->ci[XX], dd->ci[YY], dd->ci[ZZ]);
        for (i = 0; i < ireq->n; i++)
        {
            if (gmx_hash_get_minone(ga2la_specat, ireq->ind[i]) < 0)
            {
                fprintf(stderr, " %d", ireq->ind[i]+1);
            }
        }
        fprintf(stderr, "\n");
        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.",
                  dd->ci[XX], dd->ci[YY], dd->ci[ZZ],
                  nrecv_local, ireq->n, specat_type,
                  specat_type, add_err,
                  dd_dlb_is_on(dd) ? " or use the -rcon option of mdrun" : "");
    }

    spac->at_start = at_start;
    spac->at_end   = nat_tot_specat;

    if (debug)
    {
        fprintf(debug, "Done setup_specat_communication\n");
    }

    return nat_tot_specat;
}
Esempio n. 2
0
int dd_make_local_vsites(gmx_domdec_t *dd, int at_start, t_ilist *lil)
{
    gmx_domdec_specat_comm_t *spac;
    ind_req_t                *ireq;
    gmx_hash_t                ga2la_specat;
    int  ftype, nral, i, j, a;
    t_ilist                  *lilf;
    t_iatom                  *iatoms;
    int  at_end;

    spac         = dd->vsite_comm;
    ireq         = &spac->ireq[0];
    ga2la_specat = dd->ga2la_vsite;

    ireq->n = 0;
    /* Loop over all the home vsites */
    for (ftype = 0; ftype < F_NRE; ftype++)
    {
        if (interaction_function[ftype].flags & IF_VSITE)
        {
            nral = NRAL(ftype);
            lilf = &lil[ftype];
            for (i = 0; i < lilf->nr; i += 1+nral)
            {
                iatoms = lilf->iatoms + i;
                /* Check if we have the other atoms */
                for (j = 1; j < 1+nral; j++)
                {
                    if (iatoms[j] < 0)
                    {
                        /* This is not a home atom,
                         * we need to ask our neighbors.
                         */
                        a = -iatoms[j] - 1;
                        /* Check to not ask for the same atom more than once */
                        if (gmx_hash_get_minone(dd->ga2la_vsite, a) == -1)
                        {
                            /* Add this non-home atom to the list */
                            if (ireq->n+1 > ireq->nalloc)
                            {
                                ireq->nalloc = over_alloc_large(ireq->n+1);
                                srenew(ireq->ind, ireq->nalloc);
                            }
                            ireq->ind[ireq->n++] = a;
                            /* Temporarily mark with -2,
                             * we get the index later.
                             */
                            gmx_hash_set(ga2la_specat, a, -2);
                        }
                    }
                }
            }
        }
    }

    at_end = setup_specat_communication(dd, ireq, dd->vsite_comm, ga2la_specat,
                                        at_start, 1, "vsite", "");

    /* Fill in the missing indices */
    for (ftype = 0; ftype < F_NRE; ftype++)
    {
        if (interaction_function[ftype].flags & IF_VSITE)
        {
            nral = NRAL(ftype);
            lilf = &lil[ftype];
            for (i = 0; i < lilf->nr; i += 1+nral)
            {
                iatoms = lilf->iatoms + i;
                for (j = 1; j < 1+nral; j++)
                {
                    if (iatoms[j] < 0)
                    {
                        iatoms[j] = gmx_hash_get_minone(ga2la_specat, -iatoms[j]-1);
                    }
                }
            }
        }
    }

    return at_end;
}
/*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
static void walk_out(int con, int con_offset, int a, int offset, int nrec,
                     int ncon1, const t_iatom *ia1, const t_iatom *ia2,
                     const t_blocka *at2con,
                     const gmx_ga2la_t *ga2la, gmx_bool bHomeConnect,
                     gmx_domdec_constraints_t *dc,
                     gmx_domdec_specat_comm_t *dcc,
                     t_ilist *il_local,
                     ind_req_t *ireq)
{
    int            a1_gl, a2_gl, a_loc, i, coni, b;
    const t_iatom *iap;

    if (dc->gc_req[con_offset+con] == 0)
    {
        /* Add this non-home constraint to the list */
        if (dc->ncon+1 > dc->con_nalloc)
        {
            dc->con_nalloc = over_alloc_large(dc->ncon+1);
            srenew(dc->con_gl, dc->con_nalloc);
            srenew(dc->con_nlocat, dc->con_nalloc);
        }
        dc->con_gl[dc->ncon]       = con_offset + con;
        dc->con_nlocat[dc->ncon]   = (bHomeConnect ? 1 : 0);
        dc->gc_req[con_offset+con] = 1;
        if (il_local->nr + 3 > il_local->nalloc)
        {
            il_local->nalloc = over_alloc_dd(il_local->nr+3);
            srenew(il_local->iatoms, il_local->nalloc);
        }
        iap = constr_iatomptr(ncon1, ia1, ia2, con);
        il_local->iatoms[il_local->nr++] = iap[0];
        a1_gl = offset + iap[1];
        a2_gl = offset + iap[2];
        /* The following indexing code can probably be optizimed */
        if (ga2la_get_home(ga2la, a1_gl, &a_loc))
        {
            il_local->iatoms[il_local->nr++] = a_loc;
        }
        else
        {
            /* We set this index later */
            il_local->iatoms[il_local->nr++] = -a1_gl - 1;
        }
        if (ga2la_get_home(ga2la, a2_gl, &a_loc))
        {
            il_local->iatoms[il_local->nr++] = a_loc;
        }
        else
        {
            /* We set this index later */
            il_local->iatoms[il_local->nr++] = -a2_gl - 1;
        }
        dc->ncon++;
    }
    /* Check to not ask for the same atom more than once */
    if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
    {
        assert(dcc);
        /* Add this non-home atom to the list */
        if (ireq->n+1 > ireq->nalloc)
        {
            ireq->nalloc = over_alloc_large(ireq->n+1);
            srenew(ireq->ind, ireq->nalloc);
        }
        ireq->ind[ireq->n++] = offset + a;
        /* Temporarily mark with -2, we get the index later */
        gmx_hash_set(dc->ga2la, offset+a, -2);
    }

    if (nrec > 0)
    {
        for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
        {
            coni = at2con->a[i];
            if (coni != con)
            {
                /* Walk further */
                iap = constr_iatomptr(ncon1, ia1, ia2, coni);
                if (a == iap[1])
                {
                    b = iap[2];
                }
                else
                {
                    b = iap[1];
                }
                if (!ga2la_get_home(ga2la, offset+b, &a_loc))
                {
                    walk_out(coni, con_offset, b, offset, nrec-1,
                             ncon1, ia1, ia2, at2con,
                             ga2la, FALSE, dc, dcc, il_local, ireq);
                }
            }
        }
    }
}
int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
                              const struct gmx_mtop_t *mtop,
                              const int *cginfo,
                              gmx_constr_t constr, int nrec,
                              t_ilist *il_local)
{
    gmx_domdec_constraints_t   *dc;
    t_ilist                    *ilc_local, *ils_local;
    ind_req_t                  *ireq;
    const t_blocka             *at2con_mt;
    const int                 **at2settle_mt;
    gmx_hash_t                 *ga2la_specat;
    int at_end, i, j;
    t_iatom                    *iap;

    // This code should not be called unless this condition is true,
    // because that's the only time init_domdec_constraints is
    // called...
    GMX_RELEASE_ASSERT(dd->bInterCGcons || dd->bInterCGsettles, "dd_make_local_constraints called when there are no local constraints");
    // ... and init_domdec_constraints always sets
    // dd->constraint_comm...
    GMX_RELEASE_ASSERT(dd->constraint_comm, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
    // ... which static analysis needs to be reassured about, because
    // otherwise, when dd->bInterCGsettles is
    // true. dd->constraint_comm is unilaterally dereferenced before
    // the call to atoms_to_settles.

    dc = dd->constraints;

    ilc_local = &il_local[F_CONSTR];
    ils_local = &il_local[F_SETTLE];

    dc->ncon      = 0;
    ilc_local->nr = 0;
    if (dd->constraint_comm)
    {
        at2con_mt = atom2constraints_moltype(constr);
        ireq      = &dd->constraint_comm->ireq[0];
        ireq->n   = 0;
    }
    else
    {
        // Currently unreachable
        at2con_mt = NULL;
        ireq      = NULL;
    }

    if (dd->bInterCGsettles)
    {
        at2settle_mt  = atom2settle_moltype(constr);
        ils_local->nr = 0;
    }
    else
    {
        /* Settle works inside charge groups, we assigned them already */
        at2settle_mt = NULL;
    }

    if (at2settle_mt == NULL)
    {
        atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
                             ilc_local, ireq);
    }
    else
    {
        int t0_set;
        int thread;

        /* Do the constraints, if present, on the first thread.
         * Do the settles on all other threads.
         */
        t0_set = ((at2con_mt != NULL && dc->nthread > 1) ? 1 : 0);

#pragma omp parallel for num_threads(dc->nthread) schedule(static)
        for (thread = 0; thread < dc->nthread; thread++)
        {
            try
            {
                if (at2con_mt && thread == 0)
                {
                    atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
                                         ilc_local, ireq);
                }

                if (thread >= t0_set)
                {
                    int        cg0, cg1;
                    t_ilist   *ilst;
                    ind_req_t *ireqt;

                    /* Distribute the settle check+assignments over
                     * dc->nthread or dc->nthread-1 threads.
                     */
                    cg0 = (dd->ncg_home*(thread-t0_set  ))/(dc->nthread-t0_set);
                    cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);

                    if (thread == t0_set)
                    {
                        ilst = ils_local;
                    }
                    else
                    {
                        ilst = &dc->ils[thread];
                    }
                    ilst->nr = 0;

                    ireqt = &dd->constraint_comm->ireq[thread];
                    if (thread > 0)
                    {
                        ireqt->n = 0;
                    }

                    atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
                                     cg0, cg1,
                                     ilst, ireqt);
                }
            }
            GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
        }

        /* Combine the generate settles and requested indices */
        for (thread = 1; thread < dc->nthread; thread++)
        {
            t_ilist   *ilst;
            ind_req_t *ireqt;
            int        ia;

            if (thread > t0_set)
            {
                ilst = &dc->ils[thread];
                if (ils_local->nr + ilst->nr > ils_local->nalloc)
                {
                    ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
                    srenew(ils_local->iatoms, ils_local->nalloc);
                }
                for (ia = 0; ia < ilst->nr; ia++)
                {
                    ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
                }
                ils_local->nr += ilst->nr;
            }

            ireqt = &dd->constraint_comm->ireq[thread];
            if (ireq->n+ireqt->n > ireq->nalloc)
            {
                ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
                srenew(ireq->ind, ireq->nalloc);
            }
            for (ia = 0; ia < ireqt->n; ia++)
            {
                ireq->ind[ireq->n+ia] = ireqt->ind[ia];
            }
            ireq->n += ireqt->n;
        }

        if (debug)
        {
            fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
        }
    }

    if (dd->constraint_comm)
    {
        int nral1;

        at_end =
            setup_specat_communication(dd, ireq, dd->constraint_comm,
                                       dd->constraints->ga2la,
                                       at_start, 2,
                                       "constraint", " or lincs-order");

        /* Fill in the missing indices */
        ga2la_specat = dd->constraints->ga2la;

        nral1 = 1 + NRAL(F_CONSTR);
        for (i = 0; i < ilc_local->nr; i += nral1)
        {
            iap = ilc_local->iatoms + i;
            for (j = 1; j < nral1; j++)
            {
                if (iap[j] < 0)
                {
                    iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
                }
            }
        }

        nral1 = 1 + NRAL(F_SETTLE);
        for (i = 0; i < ils_local->nr; i += nral1)
        {
            iap = ils_local->iatoms + i;
            for (j = 1; j < nral1; j++)
            {
                if (iap[j] < 0)
                {
                    iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
                }
            }
        }
    }
    else
    {
        // Currently unreachable
        at_end = at_start;
    }

    return at_end;
}
Esempio n. 5
0
int dd_make_local_constraints(gmx_domdec_t *dd,int at_start,
                              const gmx_mtop_t *mtop,
                              const int *cginfo,
                              gmx_constr_t constr,int nrec,
                              t_ilist *il_local)
{
    gmx_domdec_constraints_t *dc;
    t_ilist *ilc_local,*ils_local;
    ind_req_t *ireq;
    const t_blocka *at2con_mt;
    const int **at2settle_mt;
    gmx_hash_t ga2la_specat;
    int at_end,i,j;
    t_iatom *iap;
    
    dc = dd->constraints;
    
    ilc_local = &il_local[F_CONSTR];
    ils_local = &il_local[F_SETTLE];

    dc->ncon      = 0;
    ilc_local->nr = 0;
    if (dd->constraint_comm)
    {
        at2con_mt = atom2constraints_moltype(constr);
        ireq = &dd->constraint_comm->ireq[0];
        ireq->n = 0;
    }
    else
    {
        at2con_mt = NULL;
        ireq = NULL;
    }

    if (dd->bInterCGsettles)
    {
        at2settle_mt = atom2settle_moltype(constr);
        ils_local->nr = 0;
    }
    else
    {
        /* Settle works inside charge groups, we assigned them already */
        at2settle_mt = NULL;
    }

    if (at2settle_mt == NULL)
    {
        atoms_to_constraints(dd,mtop,cginfo,at2con_mt,nrec,
                             ilc_local,ireq);
    }
    else
    {
        int t0_set;
        int thread;

        /* Do the constraints, if present, on the first thread.
         * Do the settles on all other threads.
         */
        t0_set = ((at2con_mt != NULL && dc->nthread > 1) ? 1 : 0);

#pragma omp parallel for num_threads(dc->nthread) schedule(static)
        for(thread=0; thread<dc->nthread; thread++)
        {
            if (at2con_mt && thread == 0)
            {
                atoms_to_constraints(dd,mtop,cginfo,at2con_mt,nrec,
                                     ilc_local,ireq);
            }

            if (thread >= t0_set)
            {
                int cg0,cg1;
                t_ilist *ilst;
                ind_req_t *ireqt;

                /* Distribute the settle check+assignments over
                 * dc->nthread or dc->nthread-1 threads.
                 */
                cg0 = (dd->ncg_home*(thread-t0_set  ))/(dc->nthread-t0_set);
                cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);

                if (thread == t0_set)
                {
                    ilst = ils_local;
                }
                else
                {
                    ilst = &dc->ils[thread];
                }
                ilst->nr = 0;

                ireqt = &dd->constraint_comm->ireq[thread];
                if (thread > 0)
                {
                    ireqt->n = 0;
                }

                atoms_to_settles(dd,mtop,cginfo,at2settle_mt,
                                 cg0,cg1,
                                 ilst,ireqt);
            }
        }

        /* Combine the generate settles and requested indices */
        for(thread=1; thread<dc->nthread; thread++)
        {
            t_ilist *ilst;
            ind_req_t *ireqt;
            int ia;

            if (thread > t0_set)
            {
                ilst = &dc->ils[thread];
                if (ils_local->nr + ilst->nr > ils_local->nalloc)
                {
                    ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
                    srenew(ils_local->iatoms,ils_local->nalloc);
                }
                for(ia=0; ia<ilst->nr; ia++)
                {
                    ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
                }
                ils_local->nr += ilst->nr;
            }

            ireqt = &dd->constraint_comm->ireq[thread];
            if (ireq->n+ireqt->n > ireq->nalloc)
            {
                ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
                srenew(ireq->ind,ireq->nalloc);
            }
            for(ia=0; ia<ireqt->n; ia++)
            {
                ireq->ind[ireq->n+ia] = ireqt->ind[ia];
            }
            ireq->n += ireqt->n;
        }

        if (debug)
        {
            fprintf(debug,"Settles: total %3d\n",ils_local->nr/4);
        }
    }

    if (dd->constraint_comm) {
        int nral1;

        at_end =
            setup_specat_communication(dd,ireq,dd->constraint_comm,
                                       dd->constraints->ga2la,
                                       at_start,2,
                                       "constraint"," or lincs-order");
        
        /* Fill in the missing indices */
        ga2la_specat = dd->constraints->ga2la;

        nral1 = 1 + NRAL(F_CONSTR);
        for(i=0; i<ilc_local->nr; i+=nral1)
        {
            iap = ilc_local->iatoms + i;
            for(j=1; j<nral1; j++)
            {
                if (iap[j] < 0)
                {
                    iap[j] = gmx_hash_get_minone(ga2la_specat,-iap[j]-1);
                }
            }
        }

        nral1 = 1 + NRAL(F_SETTLE);
        for(i=0; i<ils_local->nr; i+=nral1)
        {
            iap = ils_local->iatoms + i;
            for(j=1; j<nral1; j++)
            {
                if (iap[j] < 0)
                {
                    iap[j] = gmx_hash_get_minone(ga2la_specat,-iap[j]-1);
                }
            }
        }
    }
    else
    {
        at_end = at_start;
    }
    
    return at_end;
}