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,2019, 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"
54 #include "gromacs/domdec/dlbtiming.h"
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/domdec/ga2la.h"
58 #include "gromacs/domdec/hashedmap.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/constr.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/topology/ifunc.h"
66 #include "gromacs/topology/mtop_lookup.h"
67 #include "gromacs/utility/exceptions.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/listoflists.h"
72 #include "domdec_internal.h"
73 #include "domdec_specatomcomm.h"
75 using gmx::ListOfLists;
77 /*! \brief Struct used during constraint setup with domain decomposition */
78 struct gmx_domdec_constraints_t
80 //! @cond Doxygen_Suppress
81 std::vector<int> molb_con_offset; /**< Offset in the constraint array for each molblock */
82 std::vector<int> molb_ncon_mol; /**< The number of constraints per molecule for each molblock */
84 int ncon; /**< The fully local and conneced constraints */
85 /* The global constraint number, only required for clearing gc_req */
86 std::vector<int> con_gl; /**< Global constraint indices for local constraints */
87 std::vector<int> con_nlocat; /**< Number of local atoms (2/1/0) for each constraint */
89 std::vector<bool> gc_req; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
91 /* Hash table for keeping track of requests */
92 std::unique_ptr<gmx::HashedMap<int>> ga2la; /**< Global to local communicated constraint atom only index */
94 /* Multi-threading stuff */
95 int nthread; /**< Number of threads used for DD constraint setup */
96 std::vector<t_ilist> ils; /**< Constraint ilist working arrays, size \p nthread */
98 /* Buffers for requesting atoms */
99 std::vector<std::vector<int>> requestedGlobalAtomIndices; /**< Buffers for requesting global atom indices, one per thread */
104 void dd_move_x_constraints(gmx_domdec_t* dd, const matrix box, rvec* x0, rvec* x1, gmx_bool bX1IsCoord)
106 if (dd->constraint_comm)
108 dd_move_x_specat(dd, dd->constraint_comm, box, x0, x1, bX1IsCoord);
110 ddReopenBalanceRegionCpu(dd);
114 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd)
116 if (dd && dd->constraints)
118 return dd->constraints->con_nlocat;
126 void dd_clear_local_constraint_indices(gmx_domdec_t* dd)
128 gmx_domdec_constraints_t* dc = dd->constraints;
130 std::fill(dc->gc_req.begin(), dc->gc_req.end(), false);
132 if (dd->constraint_comm)
138 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
139 static void walk_out(int con,
144 gmx::ArrayRef<const int> ia1,
145 gmx::ArrayRef<const int> ia2,
146 const ListOfLists<int>& at2con,
147 const gmx_ga2la_t& ga2la,
148 gmx_bool bHomeConnect,
149 gmx_domdec_constraints_t* dc,
150 gmx_domdec_specat_comm_t* dcc,
152 std::vector<int>* ireq)
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 const int* iap = constr_iatomptr(ia1, ia2, con);
166 il_local->iatoms[il_local->nr++] = iap[0];
167 const int a1_gl = offset + iap[1];
168 const int a2_gl = offset + iap[2];
169 /* The following indexing code can probably be optizimed */
170 if (const int* a_loc = ga2la.findHome(a1_gl))
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 (const int* a_loc = ga2la.findHome(a2_gl))
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 (!dc->ga2la->find(offset + a))
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 dc->ga2la->insert(offset + a, -2);
202 /* Loop over the constraint connected to atom a */
203 for (const int coni : at2con[a])
208 const int* iap = constr_iatomptr(ia1, ia2, coni);
218 if (!ga2la.findHome(offset + b))
220 walk_out(coni, con_offset, b, offset, nrec - 1, ia1, ia2, at2con, ga2la, FALSE,
221 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 gmx::ArrayRef<const std::vector<int>> at2settle_mt,
236 std::vector<int>* ireq)
238 const gmx_ga2la_t& ga2la = *dd->ga2la;
239 int nral = NRAL(F_SETTLE);
242 for (int a = cg_start; a < cg_end; a++)
244 if (GET_CGINFO_SETTLE(cginfo[a]))
246 int a_gl = dd->globalAtomIndices[a];
248 mtopGetMolblockIndex(mtop, a_gl, &mb, nullptr, &a_mol);
250 const gmx_molblock_t* molb = &mtop->molblock[mb];
251 int settle = at2settle_mt[molb->type][a_mol];
255 int offset = a_gl - a_mol;
257 const int* ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms.data();
260 gmx_bool bAssign = FALSE;
262 for (int sa = 0; sa < nral; sa++)
264 int a_glsa = offset + ia1[settle * (1 + nral) + 1 + sa];
266 if (ga2la.findHome(a_glsa))
268 if (nlocal == 0 && a_gl == a_glsa)
278 if (ils_local->nr + 1 + nral > ils_local->nalloc)
280 ils_local->nalloc = over_alloc_dd(ils_local->nr + 1 + nral);
281 srenew(ils_local->iatoms, ils_local->nalloc);
284 ils_local->iatoms[ils_local->nr++] = ia1[settle * 4];
286 for (int sa = 0; sa < nral; sa++)
288 if (const int* a_loc = ga2la.findHome(a_gls[sa]))
290 ils_local->iatoms[ils_local->nr++] = *a_loc;
294 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
295 /* Add this non-home atom to the list */
296 ireq->push_back(a_gls[sa]);
297 /* A check on double atom requests is
298 * not required for settle.
308 /*! \brief Looks up constraint for the local atoms */
309 static void atoms_to_constraints(gmx_domdec_t* dd,
310 const gmx_mtop_t* mtop,
312 gmx::ArrayRef<const ListOfLists<int>> at2con_mt,
315 std::vector<int>* ireq)
317 gmx_domdec_constraints_t* dc = dd->constraints;
318 gmx_domdec_specat_comm_t* dcc = dd->constraint_comm;
320 const gmx_ga2la_t& ga2la = *dd->ga2la;
323 dc->con_nlocat.clear();
327 for (int a = 0; a < dd->ncg_home; a++)
329 if (GET_CGINFO_CONSTR(cginfo[a]))
331 int a_gl = dd->globalAtomIndices[a];
333 mtopGetMolblockIndex(mtop, a_gl, &mb, &molnr, &a_mol);
335 const gmx_molblock_t& molb = mtop->molblock[mb];
337 gmx::ArrayRef<const int> ia1 = mtop->moltype[molb.type].ilist[F_CONSTR].iatoms;
338 gmx::ArrayRef<const int> ia2 = mtop->moltype[molb.type].ilist[F_CONSTRNC].iatoms;
340 /* Calculate the global constraint number offset for the molecule.
341 * This is only required for the global index to make sure
342 * that we use each constraint only once.
344 const int con_offset = dc->molb_con_offset[mb] + molnr * dc->molb_ncon_mol[mb];
346 /* The global atom number offset for this molecule */
347 const int offset = a_gl - a_mol;
348 /* Loop over the constraints connected to atom a_mol in the molecule */
349 const auto& at2con = at2con_mt[molb.type];
350 for (const int con : at2con[a_mol])
352 const int* iap = constr_iatomptr(ia1, ia2, con);
362 if (const int* a_loc = ga2la.findHome(offset + b_mol))
364 /* Add this fully home constraint at the first atom */
367 dc->con_gl.push_back(con_offset + con);
368 dc->con_nlocat.push_back(2);
369 if (ilc_local->nr + 3 > ilc_local->nalloc)
371 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
372 srenew(ilc_local->iatoms, ilc_local->nalloc);
374 const int b_lo = *a_loc;
375 ilc_local->iatoms[ilc_local->nr++] = iap[0];
376 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
377 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a);
384 /* We need the nrec constraints coupled to this constraint,
385 * so we need to walk out of the home cell by nrec+1 atoms,
386 * since already atom bg is not locally present.
387 * Therefore we call walk_out with nrec recursions to go
388 * after this first call.
390 walk_out(con, con_offset, b_mol, offset, nrec, ia1, ia2, at2con, ga2la, TRUE,
391 dc, dcc, ilc_local, ireq);
397 GMX_ASSERT(dc->con_gl.size() == static_cast<size_t>(dc->ncon),
398 "con_gl size should match the number of constraints");
399 GMX_ASSERT(dc->con_nlocat.size() == static_cast<size_t>(dc->ncon),
400 "con_nlocat size should match the number of constraints");
404 fprintf(debug, "Constraints: home %3d border %3d atoms: %3zu\n", nhome, dc->ncon - nhome,
405 dd->constraint_comm ? ireq->size() : 0);
409 int dd_make_local_constraints(gmx_domdec_t* dd,
411 const struct gmx_mtop_t* mtop,
413 gmx::Constraints* constr,
417 gmx_domdec_constraints_t* dc;
418 t_ilist * ilc_local, *ils_local;
419 gmx::HashedMap<int>* ga2la_specat;
423 // This code should not be called unless this condition is true,
424 // because that's the only time init_domdec_constraints is
426 GMX_RELEASE_ASSERT(dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles,
427 "dd_make_local_constraints called when there are no local constraints");
428 // ... and init_domdec_constraints always sets
429 // dd->constraint_comm...
432 "Invalid use of dd_make_local_constraints before construction of constraint_comm");
433 // ... which static analysis needs to be reassured about, because
434 // otherwise, when dd->splitSettles is
435 // true. dd->constraint_comm is unilaterally dereferenced before
436 // the call to atoms_to_settles.
438 dc = dd->constraints;
440 ilc_local = &il_local[F_CONSTR];
441 ils_local = &il_local[F_SETTLE];
445 gmx::ArrayRef<const ListOfLists<int>> at2con_mt;
446 std::vector<int>* ireq = nullptr;
447 if (dd->constraint_comm)
449 // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
450 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
451 at2con_mt = constr->atom2constraints_moltype();
452 ireq = &dc->requestedGlobalAtomIndices[0];
456 gmx::ArrayRef<const std::vector<int>> at2settle_mt;
457 /* When settle works inside charge groups, we assigned them already */
458 if (dd->comm->systemInfo.haveSplitSettles)
460 // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
461 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
462 at2settle_mt = constr->atom2settle_moltype();
466 if (at2settle_mt.empty())
468 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec, ilc_local, ireq);
474 /* Do the constraints, if present, on the first thread.
475 * Do the settles on all other threads.
477 t0_set = ((!at2con_mt.empty() && dc->nthread > 1) ? 1 : 0);
479 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
480 for (int thread = 0; thread < dc->nthread; thread++)
484 if (!at2con_mt.empty() && thread == 0)
486 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec, ilc_local, ireq);
489 if (thread >= t0_set)
494 /* Distribute the settle check+assignments over
495 * dc->nthread or dc->nthread-1 threads.
497 cg0 = (dd->ncg_home * (thread - t0_set)) / (dc->nthread - t0_set);
498 cg1 = (dd->ncg_home * (thread - t0_set + 1)) / (dc->nthread - t0_set);
500 if (thread == t0_set)
506 ilst = &dc->ils[thread];
510 std::vector<int>& ireqt = dc->requestedGlobalAtomIndices[thread];
516 atoms_to_settles(dd, mtop, cginfo, at2settle_mt, cg0, cg1, ilst, &ireqt);
519 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
522 /* Combine the generate settles and requested indices */
523 for (int thread = 1; thread < dc->nthread; thread++)
530 ilst = &dc->ils[thread];
531 if (ils_local->nr + ilst->nr > ils_local->nalloc)
533 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
534 srenew(ils_local->iatoms, ils_local->nalloc);
536 for (ia = 0; ia < ilst->nr; ia++)
538 ils_local->iatoms[ils_local->nr + ia] = ilst->iatoms[ia];
540 ils_local->nr += ilst->nr;
543 const std::vector<int>& ireqt = dc->requestedGlobalAtomIndices[thread];
544 ireq->insert(ireq->end(), ireqt.begin(), ireqt.end());
549 fprintf(debug, "Settles: total %3d\n", ils_local->nr / 4);
553 if (dd->constraint_comm)
557 at_end = setup_specat_communication(dd, ireq, dd->constraint_comm, dd->constraints->ga2la.get(),
558 at_start, 2, "constraint", " or lincs-order");
560 /* Fill in the missing indices */
561 ga2la_specat = dd->constraints->ga2la.get();
563 nral1 = 1 + NRAL(F_CONSTR);
564 for (i = 0; i < ilc_local->nr; i += nral1)
566 iap = ilc_local->iatoms + i;
567 for (j = 1; j < nral1; j++)
571 const int* a = ga2la_specat->find(-iap[j] - 1);
572 GMX_ASSERT(a, "We have checked before that this atom index has been set");
578 nral1 = 1 + NRAL(F_SETTLE);
579 for (i = 0; i < ils_local->nr; i += nral1)
581 iap = ils_local->iatoms + i;
582 for (j = 1; j < nral1; j++)
586 const int* a = ga2la_specat->find(-iap[j] - 1);
587 GMX_ASSERT(a, "We have checked before that this atom index has been set");
595 // Currently unreachable
602 void init_domdec_constraints(gmx_domdec_t* dd, const gmx_mtop_t* mtop)
604 gmx_domdec_constraints_t* dc;
605 const gmx_molblock_t* molb;
609 fprintf(debug, "Begin init_domdec_constraints\n");
612 dd->constraints = new gmx_domdec_constraints_t;
613 dc = dd->constraints;
615 dc->molb_con_offset.resize(mtop->molblock.size());
616 dc->molb_ncon_mol.resize(mtop->molblock.size());
619 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
621 molb = &mtop->molblock[mb];
622 dc->molb_con_offset[mb] = ncon;
623 dc->molb_ncon_mol[mb] = mtop->moltype[molb->type].ilist[F_CONSTR].size() / 3
624 + mtop->moltype[molb->type].ilist[F_CONSTRNC].size() / 3;
625 ncon += molb->nmol * dc->molb_ncon_mol[mb];
630 dc->gc_req.resize(ncon);
633 /* Use a hash table for the global to local index.
634 * The number of keys is a rough estimate, it will be optimized later.
636 int numKeysEstimate = std::min(mtop->natoms / 20, mtop->natoms / (2 * dd->nnodes));
637 dc->ga2la = std::make_unique<gmx::HashedMap<int>>(numKeysEstimate);
639 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
640 dc->ils.resize(dc->nthread);
642 dd->constraint_comm = new gmx_domdec_specat_comm_t;
644 dc->requestedGlobalAtomIndices.resize(dc->nthread);