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.
38 * \brief This file implements functions for domdec to use
39 * while managing inter-atomic constraints.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "domdec_constraints.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/legacyheaders/constr.h"
55 #include "gromacs/legacyheaders/gmx_ga2la.h"
56 #include "gromacs/legacyheaders/gmx_hash.h"
57 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
58 #include "gromacs/legacyheaders/macros.h"
59 #include "gromacs/legacyheaders/types/commrec.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/smalloc.h"
67 #include "domdec_specatomcomm.h"
69 /*! \brief Struct used during constraint setup with domain decomposition */
70 typedef struct gmx_domdec_constraints {
71 //! @cond Doxygen_Suppress
72 int *molb_con_offset; /**< Offset in the constraint array for each molblock */
73 int *molb_ncon_mol; /**< The number of constraints per molecule for each molblock */
75 int ncon; /**< The fully local and conneced constraints */
76 /* The global constraint number, only required for clearing gc_req */
77 int *con_gl; /**< Global constraint indices for local constraints */
78 int *con_nlocat; /**< Number of local atoms (2/1/0) for each constraint */
79 int con_nalloc; /**< Allocation size for \p con_gl and \p con_nlocat */
81 char *gc_req; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
82 gmx_hash_t ga2la; /**< Global to local communicated constraint atom only index */
84 /* Multi-threading stuff */
85 int nthread; /**< Number of threads used for DD constraint setup */
86 t_ilist *ils; /**< Constraint ilist working arrays, size \p nthread */
88 } gmx_domdec_constraints_t;
90 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
91 rvec *x0, rvec *x1, gmx_bool bX1IsCoord)
93 if (dd->constraint_comm)
95 dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1, bX1IsCoord);
99 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd)
103 return dd->constraints->con_nlocat;
111 void dd_clear_local_constraint_indices(gmx_domdec_t *dd)
113 gmx_domdec_constraints_t *dc;
116 dc = dd->constraints;
118 for (i = 0; i < dc->ncon; i++)
120 dc->gc_req[dc->con_gl[i]] = 0;
123 if (dd->constraint_comm)
125 gmx_hash_clear_and_optimize(dc->ga2la);
129 void dd_clear_local_vsite_indices(gmx_domdec_t *dd)
133 gmx_hash_clear_and_optimize(dd->ga2la_vsite);
137 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
138 static void walk_out(int con, int con_offset, int a, int offset, int nrec,
139 int ncon1, const t_iatom *ia1, const t_iatom *ia2,
140 const t_blocka *at2con,
141 const gmx_ga2la_t ga2la, gmx_bool bHomeConnect,
142 gmx_domdec_constraints_t *dc,
143 gmx_domdec_specat_comm_t *dcc,
147 int a1_gl, a2_gl, a_loc, i, coni, b;
150 if (dc->gc_req[con_offset+con] == 0)
152 /* Add this non-home constraint to the list */
153 if (dc->ncon+1 > dc->con_nalloc)
155 dc->con_nalloc = over_alloc_large(dc->ncon+1);
156 srenew(dc->con_gl, dc->con_nalloc);
157 srenew(dc->con_nlocat, dc->con_nalloc);
159 dc->con_gl[dc->ncon] = con_offset + con;
160 dc->con_nlocat[dc->ncon] = (bHomeConnect ? 1 : 0);
161 dc->gc_req[con_offset+con] = 1;
162 if (il_local->nr + 3 > il_local->nalloc)
164 il_local->nalloc = over_alloc_dd(il_local->nr+3);
165 srenew(il_local->iatoms, il_local->nalloc);
167 iap = constr_iatomptr(ncon1, ia1, ia2, con);
168 il_local->iatoms[il_local->nr++] = iap[0];
169 a1_gl = offset + iap[1];
170 a2_gl = offset + iap[2];
171 /* The following indexing code can probably be optizimed */
172 if (ga2la_get_home(ga2la, a1_gl, &a_loc))
174 il_local->iatoms[il_local->nr++] = a_loc;
178 /* We set this index later */
179 il_local->iatoms[il_local->nr++] = -a1_gl - 1;
181 if (ga2la_get_home(ga2la, a2_gl, &a_loc))
183 il_local->iatoms[il_local->nr++] = a_loc;
187 /* We set this index later */
188 il_local->iatoms[il_local->nr++] = -a2_gl - 1;
192 /* Check to not ask for the same atom more than once */
193 if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
196 /* Add this non-home atom to the list */
197 if (ireq->n+1 > ireq->nalloc)
199 ireq->nalloc = over_alloc_large(ireq->n+1);
200 srenew(ireq->ind, ireq->nalloc);
202 ireq->ind[ireq->n++] = offset + a;
203 /* Temporarily mark with -2, we get the index later */
204 gmx_hash_set(dc->ga2la, offset+a, -2);
209 for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
215 iap = constr_iatomptr(ncon1, ia1, ia2, coni);
224 if (!ga2la_get_home(ga2la, offset+b, &a_loc))
226 walk_out(coni, con_offset, b, offset, nrec-1,
227 ncon1, ia1, ia2, at2con,
228 ga2la, FALSE, dc, dcc, il_local, ireq);
235 /*! \brief Looks up SETTLE constraints for a range of charge-groups */
236 static void atoms_to_settles(gmx_domdec_t *dd,
237 const gmx_mtop_t *mtop,
239 const int **at2settle_mt,
240 int cg_start, int cg_end,
245 gmx_mtop_atomlookup_t alook;
248 int cg, a, a_gl, a_glsa, a_gls[3], a_locs[3];
249 int mb, molnr, a_mol, offset;
250 const gmx_molblock_t *molb;
258 alook = gmx_mtop_atomlookup_settle_init(mtop);
260 nral = NRAL(F_SETTLE);
262 for (cg = cg_start; cg < cg_end; cg++)
264 if (GET_CGINFO_SETTLE(cginfo[cg]))
266 for (a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
268 a_gl = dd->gatindex[a];
270 gmx_mtop_atomnr_to_molblock_ind(alook, a_gl, &mb, &molnr, &a_mol);
271 molb = &mtop->molblock[mb];
273 settle = at2settle_mt[molb->type][a_mol];
277 offset = a_gl - a_mol;
279 ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms;
283 for (sa = 0; sa < nral; sa++)
285 a_glsa = offset + ia1[settle*(1+nral)+1+sa];
287 a_home[sa] = ga2la_get_home(ga2la, a_glsa, &a_locs[sa]);
290 if (nlocal == 0 && a_gl == a_glsa)
300 if (ils_local->nr+1+nral > ils_local->nalloc)
302 ils_local->nalloc = over_alloc_dd(ils_local->nr+1+nral);
303 srenew(ils_local->iatoms, ils_local->nalloc);
306 ils_local->iatoms[ils_local->nr++] = ia1[settle*4];
308 for (sa = 0; sa < nral; sa++)
310 if (ga2la_get_home(ga2la, a_gls[sa], &a_locs[sa]))
312 ils_local->iatoms[ils_local->nr++] = a_locs[sa];
316 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
317 /* Add this non-home atom to the list */
318 if (ireq->n+1 > ireq->nalloc)
320 ireq->nalloc = over_alloc_large(ireq->n+1);
321 srenew(ireq->ind, ireq->nalloc);
323 ireq->ind[ireq->n++] = a_gls[sa];
324 /* A check on double atom requests is
325 * not required for settle.
335 gmx_mtop_atomlookup_destroy(alook);
338 /*! \brief Looks up constraint for the local atoms */
339 static void atoms_to_constraints(gmx_domdec_t *dd,
340 const gmx_mtop_t *mtop,
342 const t_blocka *at2con_mt, int nrec,
346 const t_blocka *at2con;
348 gmx_mtop_atomlookup_t alook;
350 gmx_molblock_t *molb;
351 t_iatom *ia1, *ia2, *iap;
352 int nhome, cg, a, a_gl, a_mol, a_loc, b_lo, offset, mb, molnr, b_mol, i, con, con_offset;
353 gmx_domdec_constraints_t *dc;
354 gmx_domdec_specat_comm_t *dcc;
356 dc = dd->constraints;
357 dcc = dd->constraint_comm;
361 alook = gmx_mtop_atomlookup_init(mtop);
364 for (cg = 0; cg < dd->ncg_home; cg++)
366 if (GET_CGINFO_CONSTR(cginfo[cg]))
368 for (a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
370 a_gl = dd->gatindex[a];
372 gmx_mtop_atomnr_to_molblock_ind(alook, a_gl, &mb, &molnr, &a_mol);
373 molb = &mtop->molblock[mb];
375 ncon1 = mtop->moltype[molb->type].ilist[F_CONSTR].nr/NRAL(F_SETTLE);
377 ia1 = mtop->moltype[molb->type].ilist[F_CONSTR].iatoms;
378 ia2 = mtop->moltype[molb->type].ilist[F_CONSTRNC].iatoms;
380 /* Calculate the global constraint number offset for the molecule.
381 * This is only required for the global index to make sure
382 * that we use each constraint only once.
385 dc->molb_con_offset[mb] + molnr*dc->molb_ncon_mol[mb];
387 /* The global atom number offset for this molecule */
388 offset = a_gl - a_mol;
389 at2con = &at2con_mt[molb->type];
390 for (i = at2con->index[a_mol]; i < at2con->index[a_mol+1]; i++)
393 iap = constr_iatomptr(ncon1, ia1, ia2, con);
402 if (ga2la_get_home(ga2la, offset+b_mol, &a_loc))
404 /* Add this fully home constraint at the first atom */
407 if (dc->ncon+1 > dc->con_nalloc)
409 dc->con_nalloc = over_alloc_large(dc->ncon+1);
410 srenew(dc->con_gl, dc->con_nalloc);
411 srenew(dc->con_nlocat, dc->con_nalloc);
413 dc->con_gl[dc->ncon] = con_offset + con;
414 dc->con_nlocat[dc->ncon] = 2;
415 if (ilc_local->nr + 3 > ilc_local->nalloc)
417 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
418 srenew(ilc_local->iatoms, ilc_local->nalloc);
421 ilc_local->iatoms[ilc_local->nr++] = iap[0];
422 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
423 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a );
430 /* We need the nrec constraints coupled to this constraint,
431 * so we need to walk out of the home cell by nrec+1 atoms,
432 * since already atom bg is not locally present.
433 * Therefore we call walk_out with nrec recursions to go
434 * after this first call.
436 walk_out(con, con_offset, b_mol, offset, nrec,
437 ncon1, ia1, ia2, at2con,
438 dd->ga2la, TRUE, dc, dcc, ilc_local, ireq);
445 gmx_mtop_atomlookup_destroy(alook);
450 "Constraints: home %3d border %3d atoms: %3d\n",
451 nhome, dc->ncon-nhome,
452 dd->constraint_comm ? ireq->n : 0);
456 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
457 const gmx_mtop_t *mtop,
459 gmx_constr_t constr, int nrec,
462 gmx_domdec_constraints_t *dc;
463 t_ilist *ilc_local, *ils_local;
465 const t_blocka *at2con_mt;
466 const int **at2settle_mt;
467 gmx_hash_t ga2la_specat;
471 // This code should not be called unless this condition is true,
472 // because that's the only time init_domdec_constraints is
474 GMX_RELEASE_ASSERT(dd->bInterCGcons || dd->bInterCGsettles, "dd_make_local_constraints called when there are no local constraints");
475 // ... and init_domdec_constraints always sets
476 // dd->constraint_comm...
477 GMX_RELEASE_ASSERT(dd->constraint_comm, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
478 // ... which static analysis needs to be reassured about, because
479 // otherwise, when dd->bInterCGsettles is
480 // true. dd->constraint_comm is unilaterally dereferenced before
481 // the call to atoms_to_settles.
483 dc = dd->constraints;
485 ilc_local = &il_local[F_CONSTR];
486 ils_local = &il_local[F_SETTLE];
490 if (dd->constraint_comm)
492 at2con_mt = atom2constraints_moltype(constr);
493 ireq = &dd->constraint_comm->ireq[0];
498 // Currently unreachable
503 if (dd->bInterCGsettles)
505 at2settle_mt = atom2settle_moltype(constr);
510 /* Settle works inside charge groups, we assigned them already */
514 if (at2settle_mt == NULL)
516 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
524 /* Do the constraints, if present, on the first thread.
525 * Do the settles on all other threads.
527 t0_set = ((at2con_mt != NULL && dc->nthread > 1) ? 1 : 0);
529 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
530 for (thread = 0; thread < dc->nthread; thread++)
532 if (at2con_mt && thread == 0)
534 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
538 if (thread >= t0_set)
544 /* Distribute the settle check+assignments over
545 * dc->nthread or dc->nthread-1 threads.
547 cg0 = (dd->ncg_home*(thread-t0_set ))/(dc->nthread-t0_set);
548 cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);
550 if (thread == t0_set)
556 ilst = &dc->ils[thread];
560 ireqt = &dd->constraint_comm->ireq[thread];
566 atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
572 /* Combine the generate settles and requested indices */
573 for (thread = 1; thread < dc->nthread; thread++)
581 ilst = &dc->ils[thread];
582 if (ils_local->nr + ilst->nr > ils_local->nalloc)
584 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
585 srenew(ils_local->iatoms, ils_local->nalloc);
587 for (ia = 0; ia < ilst->nr; ia++)
589 ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
591 ils_local->nr += ilst->nr;
594 ireqt = &dd->constraint_comm->ireq[thread];
595 if (ireq->n+ireqt->n > ireq->nalloc)
597 ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
598 srenew(ireq->ind, ireq->nalloc);
600 for (ia = 0; ia < ireqt->n; ia++)
602 ireq->ind[ireq->n+ia] = ireqt->ind[ia];
609 fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
613 if (dd->constraint_comm)
618 setup_specat_communication(dd, ireq, dd->constraint_comm,
619 dd->constraints->ga2la,
621 "constraint", " or lincs-order");
623 /* Fill in the missing indices */
624 ga2la_specat = dd->constraints->ga2la;
626 nral1 = 1 + NRAL(F_CONSTR);
627 for (i = 0; i < ilc_local->nr; i += nral1)
629 iap = ilc_local->iatoms + i;
630 for (j = 1; j < nral1; j++)
634 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
639 nral1 = 1 + NRAL(F_SETTLE);
640 for (i = 0; i < ils_local->nr; i += nral1)
642 iap = ils_local->iatoms + i;
643 for (j = 1; j < nral1; j++)
647 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
654 // Currently unreachable
661 void init_domdec_constraints(gmx_domdec_t *dd,
664 gmx_domdec_constraints_t *dc;
665 gmx_molblock_t *molb;
670 fprintf(debug, "Begin init_domdec_constraints\n");
673 snew(dd->constraints, 1);
674 dc = dd->constraints;
676 snew(dc->molb_con_offset, mtop->nmolblock);
677 snew(dc->molb_ncon_mol, mtop->nmolblock);
680 for (mb = 0; mb < mtop->nmolblock; mb++)
682 molb = &mtop->molblock[mb];
683 dc->molb_con_offset[mb] = ncon;
684 dc->molb_ncon_mol[mb] =
685 mtop->moltype[molb->type].ilist[F_CONSTR].nr/3 +
686 mtop->moltype[molb->type].ilist[F_CONSTRNC].nr/3;
687 ncon += molb->nmol*dc->molb_ncon_mol[mb];
692 snew(dc->gc_req, ncon);
693 for (c = 0; c < ncon; c++)
699 /* Use a hash table for the global to local index.
700 * The number of keys is a rough estimate, it will be optimized later.
702 dc->ga2la = gmx_hash_init(std::min(mtop->natoms/20,
703 mtop->natoms/(2*dd->nnodes)));
705 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
706 snew(dc->ils, dc->nthread);
708 dd->constraint_comm = specat_comm_init(dc->nthread);