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,2018, 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/mdtypes/forcerec.h" // only for GET_CGINFO_*
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/topology/ifunc.h"
64 #include "gromacs/topology/mtop_lookup.h"
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxassert.h"
69 #include "domdec_specatomcomm.h"
70 #include "domdec_vsite.h"
73 /*! \brief Struct used during constraint setup with domain decomposition */
74 struct gmx_domdec_constraints_t
76 //! @cond Doxygen_Suppress
77 std::vector<int> molb_con_offset; /**< Offset in the constraint array for each molblock */
78 std::vector<int> molb_ncon_mol; /**< The number of constraints per molecule for each molblock */
80 int ncon; /**< The fully local and conneced constraints */
81 /* The global constraint number, only required for clearing gc_req */
82 std::vector<int> con_gl; /**< Global constraint indices for local constraints */
83 std::vector<int> con_nlocat; /**< Number of local atoms (2/1/0) for each constraint */
85 std::vector<bool> gc_req; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
86 gmx_hash_t *ga2la; /**< Global to local communicated constraint atom only index, TODO: get rid of the plain pointer */
88 /* Multi-threading stuff */
89 int nthread; /**< Number of threads used for DD constraint setup */
90 std::vector<t_ilist> ils; /**< Constraint ilist working arrays, size \p nthread */
92 /* Buffers for requesting atoms */
93 std::vector < std::vector < int>> requestedGlobalAtomIndices; /**< Buffers for requesting global atom indices, one per thread */
98 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
99 rvec *x0, rvec *x1, gmx_bool bX1IsCoord)
101 if (dd->constraint_comm)
103 dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1, bX1IsCoord);
105 ddReopenBalanceRegionCpu(dd);
109 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd)
111 if (dd && dd->constraints)
113 return dd->constraints->con_nlocat;
117 return gmx::EmptyArrayRef();
121 void dd_clear_local_constraint_indices(gmx_domdec_t *dd)
123 gmx_domdec_constraints_t *dc = dd->constraints;
125 std::fill(dc->gc_req.begin(), dc->gc_req.end(), false);
127 if (dd->constraint_comm)
129 gmx_hash_clear_and_optimize(dc->ga2la);
133 void dd_clear_local_vsite_indices(gmx_domdec_t *dd)
137 gmx_hash_clear_and_optimize(dd->ga2la_vsite);
141 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
142 static void walk_out(int con, int con_offset, int a, int offset, int nrec,
143 int ncon1, const t_iatom *ia1, const t_iatom *ia2,
144 const t_blocka *at2con,
145 const gmx_ga2la_t *ga2la, gmx_bool bHomeConnect,
146 gmx_domdec_constraints_t *dc,
147 gmx_domdec_specat_comm_t *dcc,
149 std::vector<int> *ireq)
151 int a1_gl, a2_gl, a_loc, i, coni, b;
154 if (!dc->gc_req[con_offset + con])
156 /* Add this non-home constraint to the list */
157 dc->con_gl.push_back(con_offset + con);
158 dc->con_nlocat.push_back(bHomeConnect ? 1 : 0);
159 dc->gc_req[con_offset + con] = true;
160 if (il_local->nr + 3 > il_local->nalloc)
162 il_local->nalloc = over_alloc_dd(il_local->nr+3);
163 srenew(il_local->iatoms, il_local->nalloc);
165 iap = constr_iatomptr(ncon1, ia1, ia2, con);
166 il_local->iatoms[il_local->nr++] = iap[0];
167 a1_gl = offset + iap[1];
168 a2_gl = offset + iap[2];
169 /* The following indexing code can probably be optizimed */
170 if (ga2la_get_home(ga2la, a1_gl, &a_loc))
172 il_local->iatoms[il_local->nr++] = a_loc;
176 /* We set this index later */
177 il_local->iatoms[il_local->nr++] = -a1_gl - 1;
179 if (ga2la_get_home(ga2la, a2_gl, &a_loc))
181 il_local->iatoms[il_local->nr++] = a_loc;
185 /* We set this index later */
186 il_local->iatoms[il_local->nr++] = -a2_gl - 1;
190 /* Check to not ask for the same atom more than once */
191 if (gmx_hash_get_minone(dc->ga2la, offset+a) == -1)
194 /* Add this non-home atom to the list */
195 ireq->push_back(offset + a);
196 /* Temporarily mark with -2, we get the index later */
197 gmx_hash_set(dc->ga2la, offset+a, -2);
202 for (i = at2con->index[a]; i < at2con->index[a+1]; i++)
208 iap = constr_iatomptr(ncon1, ia1, ia2, coni);
217 if (!ga2la_get_home(ga2la, offset+b, &a_loc))
219 walk_out(coni, con_offset, b, offset, nrec-1,
220 ncon1, ia1, ia2, at2con,
221 ga2la, FALSE, dc, dcc, il_local, ireq);
228 /*! \brief Looks up SETTLE constraints for a range of charge-groups */
229 static void atoms_to_settles(gmx_domdec_t *dd,
230 const gmx_mtop_t *mtop,
232 const int **at2settle_mt,
233 int cg_start, int cg_end,
235 std::vector<int> *ireq)
237 gmx_ga2la_t *ga2la = dd->ga2la;
238 int nral = NRAL(F_SETTLE);
241 for (int cg = cg_start; cg < cg_end; cg++)
243 if (GET_CGINFO_SETTLE(cginfo[cg]))
245 for (int a : dd->atomGrouping().block(cg))
247 int a_gl = dd->globalAtomIndices[a];
249 mtopGetMolblockIndex(mtop, a_gl, &mb, nullptr, &a_mol);
251 const gmx_molblock_t *molb = &mtop->molblock[mb];
252 int settle = at2settle_mt[molb->type][a_mol];
256 int offset = a_gl - a_mol;
258 t_iatom *ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms;
260 int a_gls[3], a_locs[3];
261 gmx_bool bAssign = FALSE;
263 for (int sa = 0; sa < nral; sa++)
265 int a_glsa = offset + ia1[settle*(1+nral)+1+sa];
267 if (ga2la_get_home(ga2la, a_glsa, &a_locs[sa]))
269 if (nlocal == 0 && a_gl == a_glsa)
279 if (ils_local->nr+1+nral > ils_local->nalloc)
281 ils_local->nalloc = over_alloc_dd(ils_local->nr+1+nral);
282 srenew(ils_local->iatoms, ils_local->nalloc);
285 ils_local->iatoms[ils_local->nr++] = ia1[settle*4];
287 for (int sa = 0; sa < nral; sa++)
289 if (ga2la_get_home(ga2la, a_gls[sa], &a_locs[sa]))
291 ils_local->iatoms[ils_local->nr++] = a_locs[sa];
295 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
296 /* Add this non-home atom to the list */
297 ireq->push_back(a_gls[sa]);
298 /* A check on double atom requests is
299 * not required for settle.
310 /*! \brief Looks up constraint for the local atoms */
311 static void atoms_to_constraints(gmx_domdec_t *dd,
312 const gmx_mtop_t *mtop,
314 gmx::ArrayRef<const t_blocka> at2con_mt, int nrec,
316 std::vector<int> *ireq)
318 const t_blocka *at2con;
320 t_iatom *ia1, *ia2, *iap;
321 int a_loc, b_lo, offset, b_mol, i, con, con_offset;
323 gmx_domdec_constraints_t *dc = dd->constraints;
324 gmx_domdec_specat_comm_t *dcc = dd->constraint_comm;
326 gmx_ga2la_t *ga2la = dd->ga2la;
329 dc->con_nlocat.clear();
333 for (int cg = 0; cg < dd->ncg_home; cg++)
335 if (GET_CGINFO_CONSTR(cginfo[cg]))
337 for (int a : dd->atomGrouping().block(cg))
339 int a_gl = dd->globalAtomIndices[a];
341 mtopGetMolblockIndex(mtop, a_gl, &mb, &molnr, &a_mol);
343 const gmx_molblock_t *molb = &mtop->molblock[mb];
345 ncon1 = mtop->moltype[molb->type].ilist[F_CONSTR].nr/NRAL(F_SETTLE);
347 ia1 = mtop->moltype[molb->type].ilist[F_CONSTR].iatoms;
348 ia2 = mtop->moltype[molb->type].ilist[F_CONSTRNC].iatoms;
350 /* Calculate the global constraint number offset for the molecule.
351 * This is only required for the global index to make sure
352 * that we use each constraint only once.
355 dc->molb_con_offset[mb] + molnr*dc->molb_ncon_mol[mb];
357 /* The global atom number offset for this molecule */
358 offset = a_gl - a_mol;
359 at2con = &at2con_mt[molb->type];
360 for (i = at2con->index[a_mol]; i < at2con->index[a_mol+1]; i++)
363 iap = constr_iatomptr(ncon1, ia1, ia2, con);
372 if (ga2la_get_home(ga2la, offset+b_mol, &a_loc))
374 /* Add this fully home constraint at the first atom */
377 dc->con_gl.push_back(con_offset + con);
378 dc->con_nlocat.push_back(2);
379 if (ilc_local->nr + 3 > ilc_local->nalloc)
381 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
382 srenew(ilc_local->iatoms, ilc_local->nalloc);
385 ilc_local->iatoms[ilc_local->nr++] = iap[0];
386 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
387 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a );
394 /* We need the nrec constraints coupled to this constraint,
395 * so we need to walk out of the home cell by nrec+1 atoms,
396 * since already atom bg is not locally present.
397 * Therefore we call walk_out with nrec recursions to go
398 * after this first call.
400 walk_out(con, con_offset, b_mol, offset, nrec,
401 ncon1, ia1, ia2, at2con,
402 dd->ga2la, TRUE, dc, dcc, ilc_local, ireq);
409 GMX_ASSERT(dc->con_gl.size() == static_cast<size_t>(dc->ncon), "con_gl size should match the number of constraints");
410 GMX_ASSERT(dc->con_nlocat.size() == static_cast<size_t>(dc->ncon), "con_nlocat size should match the number of constraints");
415 "Constraints: home %3d border %3d atoms: %3zu\n",
416 nhome, dc->ncon-nhome,
417 dd->constraint_comm ? ireq->size() : 0);
421 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
422 const struct gmx_mtop_t *mtop,
424 gmx::Constraints *constr, int nrec,
427 gmx_domdec_constraints_t *dc;
428 t_ilist *ilc_local, *ils_local;
429 std::vector<int> *ireq;
430 gmx::ArrayRef<const t_blocka> at2con_mt;
431 const int **at2settle_mt;
432 gmx_hash_t *ga2la_specat;
436 // This code should not be called unless this condition is true,
437 // because that's the only time init_domdec_constraints is
439 GMX_RELEASE_ASSERT(dd->bInterCGcons || dd->bInterCGsettles, "dd_make_local_constraints called when there are no local constraints");
440 // ... and init_domdec_constraints always sets
441 // dd->constraint_comm...
442 GMX_RELEASE_ASSERT(dd->constraint_comm, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
443 // ... which static analysis needs to be reassured about, because
444 // otherwise, when dd->bInterCGsettles is
445 // true. dd->constraint_comm is unilaterally dereferenced before
446 // the call to atoms_to_settles.
448 dc = dd->constraints;
450 ilc_local = &il_local[F_CONSTR];
451 ils_local = &il_local[F_SETTLE];
455 if (dd->constraint_comm)
457 // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
458 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
459 at2con_mt = constr->atom2constraints_moltype();
460 ireq = &dc->requestedGlobalAtomIndices[0];
465 // Currently unreachable
466 at2con_mt = gmx::EmptyArrayRef();
470 if (dd->bInterCGsettles)
472 // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
473 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
474 at2settle_mt = constr->atom2settle_moltype();
479 /* Settle works inside charge groups, we assigned them already */
480 at2settle_mt = nullptr;
483 if (at2settle_mt == nullptr)
485 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
493 /* Do the constraints, if present, on the first thread.
494 * Do the settles on all other threads.
496 t0_set = ((!at2con_mt.empty() && dc->nthread > 1) ? 1 : 0);
498 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
499 for (thread = 0; thread < dc->nthread; thread++)
503 if (!at2con_mt.empty() && thread == 0)
505 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec,
509 if (thread >= t0_set)
514 /* Distribute the settle check+assignments over
515 * dc->nthread or dc->nthread-1 threads.
517 cg0 = (dd->ncg_home*(thread-t0_set ))/(dc->nthread-t0_set);
518 cg1 = (dd->ncg_home*(thread-t0_set+1))/(dc->nthread-t0_set);
520 if (thread == t0_set)
526 ilst = &dc->ils[thread];
530 std::vector<int> &ireqt = dc->requestedGlobalAtomIndices[thread];
536 atoms_to_settles(dd, mtop, cginfo, at2settle_mt,
541 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
544 /* Combine the generate settles and requested indices */
545 for (thread = 1; thread < dc->nthread; thread++)
552 ilst = &dc->ils[thread];
553 if (ils_local->nr + ilst->nr > ils_local->nalloc)
555 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
556 srenew(ils_local->iatoms, ils_local->nalloc);
558 for (ia = 0; ia < ilst->nr; ia++)
560 ils_local->iatoms[ils_local->nr+ia] = ilst->iatoms[ia];
562 ils_local->nr += ilst->nr;
565 const std::vector<int> &ireqt = dc->requestedGlobalAtomIndices[thread];
566 ireq->insert(ireq->end(), ireqt.begin(), ireqt.end());
571 fprintf(debug, "Settles: total %3d\n", ils_local->nr/4);
575 if (dd->constraint_comm)
580 setup_specat_communication(dd, ireq, dd->constraint_comm,
581 dd->constraints->ga2la,
583 "constraint", " or lincs-order");
585 /* Fill in the missing indices */
586 ga2la_specat = dd->constraints->ga2la;
588 nral1 = 1 + NRAL(F_CONSTR);
589 for (i = 0; i < ilc_local->nr; i += nral1)
591 iap = ilc_local->iatoms + i;
592 for (j = 1; j < nral1; j++)
596 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
601 nral1 = 1 + NRAL(F_SETTLE);
602 for (i = 0; i < ils_local->nr; i += nral1)
604 iap = ils_local->iatoms + i;
605 for (j = 1; j < nral1; j++)
609 iap[j] = gmx_hash_get_minone(ga2la_specat, -iap[j]-1);
616 // Currently unreachable
623 void init_domdec_constraints(gmx_domdec_t *dd,
624 const gmx_mtop_t *mtop)
626 gmx_domdec_constraints_t *dc;
627 const gmx_molblock_t *molb;
631 fprintf(debug, "Begin init_domdec_constraints\n");
634 dd->constraints = new gmx_domdec_constraints_t;
635 dc = dd->constraints;
637 dc->molb_con_offset.resize(mtop->molblock.size());
638 dc->molb_ncon_mol.resize(mtop->molblock.size());
641 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
643 molb = &mtop->molblock[mb];
644 dc->molb_con_offset[mb] = ncon;
645 dc->molb_ncon_mol[mb] =
646 mtop->moltype[molb->type].ilist[F_CONSTR].nr/3 +
647 mtop->moltype[molb->type].ilist[F_CONSTRNC].nr/3;
648 ncon += molb->nmol*dc->molb_ncon_mol[mb];
653 dc->gc_req.resize(ncon);
656 /* Use a hash table for the global to local index.
657 * The number of keys is a rough estimate, it will be optimized later.
659 dc->ga2la = gmx_hash_init(std::min(mtop->natoms/20,
660 mtop->natoms/(2*dd->nnodes)));
662 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
663 dc->ils.resize(dc->nthread);
665 dd->constraint_comm = new gmx_domdec_specat_comm_t;
667 dc->requestedGlobalAtomIndices.resize(dc->nthread);