2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2012,2013,2014,2015,2016,2017, 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/dlbtiming.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/domdec/ga2la.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/constr.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/topology/mtop_lookup.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/smalloc.h"
68 #include "domdec_specatomcomm.h"
69 #include "domdec_vsite.h"
72 /*! \brief Struct used during constraint setup with domain decomposition */
73 struct gmx_domdec_constraints_t {
74 //! @cond Doxygen_Suppress
75 int *molb_con_offset; /**< Offset in the constraint array for each molblock */
76 int *molb_ncon_mol; /**< The number of constraints per molecule for each molblock */
78 int ncon; /**< The fully local and conneced constraints */
79 /* The global constraint number, only required for clearing gc_req */
80 int *con_gl; /**< Global constraint indices for local constraints */
81 int *con_nlocat; /**< Number of local atoms (2/1/0) for each constraint */
82 int con_nalloc; /**< Allocation size for \p con_gl and \p con_nlocat */
84 char *gc_req; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
85 gmx_hash_t *ga2la; /**< Global to local communicated constraint atom only index */
87 /* Multi-threading stuff */
88 int nthread; /**< Number of threads used for DD constraint setup */
89 t_ilist *ils; /**< Constraint ilist working arrays, size \p nthread */
93 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
94 rvec *x0, rvec *x1, gmx_bool bX1IsCoord)
96 if (dd->constraint_comm)
98 dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1, bX1IsCoord);
100 ddReopenBalanceRegionCpu(dd);
104 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd)
108 return dd->constraints->con_nlocat;
116 void dd_clear_local_constraint_indices(gmx_domdec_t *dd)
118 gmx_domdec_constraints_t *dc;
121 dc = dd->constraints;
123 for (i = 0; i < dc->ncon; i++)
125 dc->gc_req[dc->con_gl[i]] = 0;
128 if (dd->constraint_comm)
130 gmx_hash_clear_and_optimize(dc->ga2la);
134 void dd_clear_local_vsite_indices(gmx_domdec_t *dd)
138 gmx_hash_clear_and_optimize(dd->ga2la_vsite);
142 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
143 static void walk_out(int con, int con_offset, int a, int offset, int nrec,
144 int ncon1, const t_iatom *ia1, const t_iatom *ia2,
145 const t_blocka *at2con,
146 const gmx_ga2la_t *ga2la, gmx_bool bHomeConnect,
147 gmx_domdec_constraints_t *dc,
148 gmx_domdec_specat_comm_t *dcc,
152 int a1_gl, a2_gl, a_loc, i, coni, b;
155 if (dc->gc_req[con_offset+con] == 0)
157 /* Add this non-home constraint to the list */
158 if (dc->ncon+1 > dc->con_nalloc)
160 dc->con_nalloc = over_alloc_large(dc->ncon+1);
161 srenew(dc->con_gl, dc->con_nalloc);
162 srenew(dc->con_nlocat, dc->con_nalloc);
164 dc->con_gl[dc->ncon] = con_offset + con;
165 dc->con_nlocat[dc->ncon] = (bHomeConnect ? 1 : 0);
166 dc->gc_req[con_offset+con] = 1;
167 if (il_local->nr + 3 > il_local->nalloc)
169 il_local->nalloc = over_alloc_dd(il_local->nr+3);
170 srenew(il_local->iatoms, il_local->nalloc);
172 iap = constr_iatomptr(ncon1, ia1, ia2, con);
173 il_local->iatoms[il_local->nr++] = iap[0];
174 a1_gl = offset + iap[1];
175 a2_gl = offset + iap[2];
176 /* The following indexing code can probably be optizimed */
177 if (ga2la_get_home(ga2la, a1_gl, &a_loc))
179 il_local->iatoms[il_local->nr++] = a_loc;
183 /* We set this index later */
184 il_local->iatoms[il_local->nr++] = -a1_gl - 1;
186 if (ga2la_get_home(ga2la, a2_gl, &a_loc))
188 il_local->iatoms[il_local->nr++] = a_loc;
192 /* We set this index later */
193 il_local->iatoms[il_local->nr++] = -a2_gl - 1;
197 /* Check to not ask for the same atom more than once */
198 if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
201 /* Add this non-home atom to the list */
202 if (ireq->n+1 > ireq->nalloc)
204 ireq->nalloc = over_alloc_large(ireq->n+1);
205 srenew(ireq->ind, ireq->nalloc);
207 ireq->ind[ireq->n++] = offset + a;
208 /* Temporarily mark with -2, we get the index later */
209 gmx_hash_set(dc->ga2la, offset+a, -2);
214 for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
220 iap = constr_iatomptr(ncon1, ia1, ia2, coni);
229 if (!ga2la_get_home(ga2la, offset+b, &a_loc))
231 walk_out(coni, con_offset, b, offset, nrec-1,
232 ncon1, ia1, ia2, at2con,
233 ga2la, FALSE, dc, dcc, il_local, ireq);
240 /*! \brief Looks up SETTLE constraints for a range of charge-groups */
241 static void atoms_to_settles(gmx_domdec_t *dd,
242 const gmx_mtop_t *mtop,
244 const int **at2settle_mt,
245 int cg_start, int cg_end,
249 gmx_ga2la_t *ga2la = dd->ga2la;
250 int nral = NRAL(F_SETTLE);
253 for (int cg = cg_start; cg < cg_end; cg++)
255 if (GET_CGINFO_SETTLE(cginfo[cg]))
257 for (int a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
259 int a_gl = dd->gatindex[a];
261 mtopGetMolblockIndex(mtop, a_gl, &mb, nullptr, &a_mol);
263 const gmx_molblock_t *molb = &mtop->molblock[mb];
264 int settle = at2settle_mt[molb->type][a_mol];
268 int offset = a_gl - a_mol;
270 t_iatom *ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms;
272 int a_gls[3], a_locs[3];
273 gmx_bool bAssign = FALSE;
275 for (int sa = 0; sa < nral; sa++)
277 int a_glsa = offset + ia1[settle*(1+nral)+1+sa];
279 if (ga2la_get_home(ga2la, a_glsa, &a_locs[sa]))
281 if (nlocal == 0 && a_gl == a_glsa)
291 if (ils_local->nr+1+nral > ils_local->nalloc)
293 ils_local->nalloc = over_alloc_dd(ils_local->nr+1+nral);
294 srenew(ils_local->iatoms, ils_local->nalloc);
297 ils_local->iatoms[ils_local->nr++] = ia1[settle*4];
299 for (int sa = 0; sa < nral; sa++)
301 if (ga2la_get_home(ga2la, a_gls[sa], &a_locs[sa]))
303 ils_local->iatoms[ils_local->nr++] = a_locs[sa];
307 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
308 /* Add this non-home atom to the list */
309 if (ireq->n+1 > ireq->nalloc)
311 ireq->nalloc = over_alloc_large(ireq->n+1);
312 srenew(ireq->ind, ireq->nalloc);
314 ireq->ind[ireq->n++] = a_gls[sa];
315 /* A check on double atom requests is
316 * not required for settle.
327 /*! \brief Looks up constraint for the local atoms */
328 static void atoms_to_constraints(gmx_domdec_t *dd,
329 const gmx_mtop_t *mtop,
331 const t_blocka *at2con_mt, int nrec,
335 const t_blocka *at2con;
337 t_iatom *ia1, *ia2, *iap;
338 int a_loc, b_lo, offset, b_mol, i, con, con_offset;
340 gmx_domdec_constraints_t *dc = dd->constraints;
341 gmx_domdec_specat_comm_t *dcc = dd->constraint_comm;
343 gmx_ga2la_t *ga2la = dd->ga2la;
347 for (int cg = 0; cg < dd->ncg_home; cg++)
349 if (GET_CGINFO_CONSTR(cginfo[cg]))
351 for (int a = dd->cgindex[cg]; a < dd->cgindex[cg+1]; a++)
353 int a_gl = dd->gatindex[a];
355 mtopGetMolblockIndex(mtop, a_gl, &mb, &molnr, &a_mol);
357 const gmx_molblock_t *molb = &mtop->molblock[mb];
359 ncon1 = mtop->moltype[molb->type].ilist[F_CONSTR].nr/NRAL(F_SETTLE);
361 ia1 = mtop->moltype[molb->type].ilist[F_CONSTR].iatoms;
362 ia2 = mtop->moltype[molb->type].ilist[F_CONSTRNC].iatoms;
364 /* Calculate the global constraint number offset for the molecule.
365 * This is only required for the global index to make sure
366 * that we use each constraint only once.
369 dc->molb_con_offset[mb] + molnr*dc->molb_ncon_mol[mb];
371 /* The global atom number offset for this molecule */
372 offset = a_gl - a_mol;
373 at2con = &at2con_mt[molb->type];
374 for (i = at2con->index[a_mol]; i < at2con->index[a_mol+1]; i++)
377 iap = constr_iatomptr(ncon1, ia1, ia2, con);
386 if (ga2la_get_home(ga2la, offset+b_mol, &a_loc))
388 /* Add this fully home constraint at the first atom */
391 if (dc->ncon+1 > dc->con_nalloc)
393 dc->con_nalloc = over_alloc_large(dc->ncon+1);
394 srenew(dc->con_gl, dc->con_nalloc);
395 srenew(dc->con_nlocat, dc->con_nalloc);
397 dc->con_gl[dc->ncon] = con_offset + con;
398 dc->con_nlocat[dc->ncon] = 2;
399 if (ilc_local->nr + 3 > ilc_local->nalloc)
401 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
402 srenew(ilc_local->iatoms, ilc_local->nalloc);
405 ilc_local->iatoms[ilc_local->nr++] = iap[0];
406 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
407 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a );
414 /* We need the nrec constraints coupled to this constraint,
415 * so we need to walk out of the home cell by nrec+1 atoms,
416 * since already atom bg is not locally present.
417 * Therefore we call walk_out with nrec recursions to go
418 * after this first call.
420 walk_out(con, con_offset, b_mol, offset, nrec,
421 ncon1, ia1, ia2, at2con,
422 dd->ga2la, TRUE, dc, dcc, ilc_local, ireq);
432 "Constraints: home %3d border %3d atoms: %3d\n",
433 nhome, dc->ncon-nhome,
434 dd->constraint_comm ? ireq->n : 0);
438 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
439 const struct gmx_mtop_t *mtop,
441 gmx_constr_t constr, int nrec,
444 gmx_domdec_constraints_t *dc;
445 t_ilist *ilc_local, *ils_local;
447 const t_blocka *at2con_mt;
448 const int **at2settle_mt;
449 gmx_hash_t *ga2la_specat;
453 // This code should not be called unless this condition is true,
454 // because that's the only time init_domdec_constraints is
456 GMX_RELEASE_ASSERT(dd->bInterCGcons || dd->bInterCGsettles, "dd_make_local_constraints called when there are no local constraints");
457 // ... and init_domdec_constraints always sets
458 // dd->constraint_comm...
459 GMX_RELEASE_ASSERT(dd->constraint_comm, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
460 // ... which static analysis needs to be reassured about, because
461 // otherwise, when dd->bInterCGsettles is
462 // true. dd->constraint_comm is unilaterally dereferenced before
463 // the call to atoms_to_settles.
465 dc = dd->constraints;
467 ilc_local = &il_local[F_CONSTR];
468 ils_local = &il_local[F_SETTLE];
472 if (dd->constraint_comm)
474 at2con_mt = atom2constraints_moltype(constr);
475 ireq = &dd->constraint_comm->ireq[0];
480 // Currently unreachable
485 if (dd->bInterCGsettles)
487 at2settle_mt = atom2settle_moltype(constr);
492 /* Settle works inside charge groups, we assigned them already */
493 at2settle_mt = nullptr;
496 if (at2settle_mt == nullptr)
498 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
506 /* Do the constraints, if present, on the first thread.
507 * Do the settles on all other threads.
509 t0_set = ((at2con_mt != nullptr && dc->nthread > 1) ? 1 : 0);
511 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
512 for (thread = 0; thread < dc->nthread; thread++)
516 if (at2con_mt && thread == 0)
518 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
522 if (thread >= t0_set)
528 /* Distribute the settle check+assignments over
529 * dc->nthread or dc->nthread-1 threads.
531 cg0 = (dd->ncg_home*(thread-t0_set ))/(dc->nthread-t0_set);
532 cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);
534 if (thread == t0_set)
540 ilst = &dc->ils[thread];
544 ireqt = &dd->constraint_comm->ireq[thread];
550 atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
555 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
558 /* Combine the generate settles and requested indices */
559 for (thread = 1; thread < dc->nthread; thread++)
567 ilst = &dc->ils[thread];
568 if (ils_local->nr + ilst->nr > ils_local->nalloc)
570 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
571 srenew(ils_local->iatoms, ils_local->nalloc);
573 for (ia = 0; ia < ilst->nr; ia++)
575 ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
577 ils_local->nr += ilst->nr;
580 ireqt = &dd->constraint_comm->ireq[thread];
581 if (ireq->n+ireqt->n > ireq->nalloc)
583 ireq->nalloc = over_alloc_large(ireq->n+ireqt->n);
584 srenew(ireq->ind, ireq->nalloc);
586 for (ia = 0; ia < ireqt->n; ia++)
588 ireq->ind[ireq->n+ia] = ireqt->ind[ia];
595 fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
599 if (dd->constraint_comm)
604 setup_specat_communication(dd, ireq, dd->constraint_comm,
605 dd->constraints->ga2la,
607 "constraint", " or lincs-order");
609 /* Fill in the missing indices */
610 ga2la_specat = dd->constraints->ga2la;
612 nral1 = 1 + NRAL(F_CONSTR);
613 for (i = 0; i < ilc_local->nr; i += nral1)
615 iap = ilc_local->iatoms + i;
616 for (j = 1; j < nral1; j++)
620 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
625 nral1 = 1 + NRAL(F_SETTLE);
626 for (i = 0; i < ils_local->nr; i += nral1)
628 iap = ils_local->iatoms + i;
629 for (j = 1; j < nral1; j++)
633 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
640 // Currently unreachable
647 void init_domdec_constraints(gmx_domdec_t *dd,
648 const gmx_mtop_t *mtop)
650 gmx_domdec_constraints_t *dc;
651 const gmx_molblock_t *molb;
656 fprintf(debug, "Begin init_domdec_constraints\n");
659 snew(dd->constraints, 1);
660 dc = dd->constraints;
662 snew(dc->molb_con_offset, mtop->nmolblock);
663 snew(dc->molb_ncon_mol, mtop->nmolblock);
666 for (mb = 0; mb < mtop->nmolblock; mb++)
668 molb = &mtop->molblock[mb];
669 dc->molb_con_offset[mb] = ncon;
670 dc->molb_ncon_mol[mb] =
671 mtop->moltype[molb->type].ilist[F_CONSTR].nr/3 +
672 mtop->moltype[molb->type].ilist[F_CONSTRNC].nr/3;
673 ncon += molb->nmol*dc->molb_ncon_mol[mb];
678 snew(dc->gc_req, ncon);
679 for (c = 0; c < ncon; c++)
685 /* Use a hash table for the global to local index.
686 * The number of keys is a rough estimate, it will be optimized later.
688 dc->ga2la = gmx_hash_init(std::min(mtop->natoms/20,
689 mtop->natoms/(2*dd->nnodes)));
691 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
692 snew(dc->ils, dc->nthread);
694 dd->constraint_comm = specat_comm_init(dc->nthread);