2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010 by the GROMACS development team.
5 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
6 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file implements functions for domdec to use
41 * while managing inter-atomic constraints.
43 * \author Berk Hess <hess@kth.se>
44 * \ingroup module_domdec
49 #include "domdec_constraints.h"
56 #include "gromacs/domdec/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/domdec/hashedmap.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/constr.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/topology/ifunc.h"
68 #include "gromacs/topology/mtop_lookup.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/listoflists.h"
74 #include "domdec_internal.h"
75 #include "domdec_specatomcomm.h"
77 using gmx::ListOfLists;
79 /*! \brief Struct used during constraint setup with domain decomposition */
80 struct gmx_domdec_constraints_t
82 //! @cond Doxygen_Suppress
83 std::vector<int> molb_con_offset; /**< Offset in the constraint array for each molblock */
84 std::vector<int> molb_ncon_mol; /**< The number of constraints per molecule for each molblock */
86 int ncon; /**< The fully local and conneced constraints */
87 /* The global constraint number, only required for clearing gc_req */
88 std::vector<int> con_gl; /**< Global constraint indices for local constraints */
89 std::vector<int> con_nlocat; /**< Number of local atoms (2/1/0) for each constraint */
91 std::vector<bool> gc_req; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
93 /* Hash table for keeping track of requests */
94 std::unique_ptr<gmx::HashedMap<int>> ga2la; /**< Global to local communicated constraint atom only index */
96 /* Multi-threading stuff */
97 int nthread; /**< Number of threads used for DD constraint setup */
98 std::vector<t_ilist> ils; /**< Constraint ilist working arrays, size \p nthread */
100 /* Buffers for requesting atoms */
101 std::vector<std::vector<int>> requestedGlobalAtomIndices; /**< Buffers for requesting global atom indices, one per thread */
106 void dd_move_x_constraints(gmx_domdec_t* dd,
108 gmx::ArrayRef<gmx::RVec> x0,
109 gmx::ArrayRef<gmx::RVec> x1,
112 if (dd->constraint_comm)
114 dd_move_x_specat(dd, dd->constraint_comm, box, as_rvec_array(x0.data()),
115 as_rvec_array(x1.data()), bX1IsCoord);
117 ddReopenBalanceRegionCpu(dd);
121 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd)
123 if (dd && dd->constraints)
125 return dd->constraints->con_nlocat;
133 void dd_clear_local_constraint_indices(gmx_domdec_t* dd)
135 gmx_domdec_constraints_t* dc = dd->constraints;
137 std::fill(dc->gc_req.begin(), dc->gc_req.end(), false);
139 if (dd->constraint_comm)
145 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
146 static void walk_out(int con,
151 gmx::ArrayRef<const int> ia1,
152 gmx::ArrayRef<const int> ia2,
153 const ListOfLists<int>& at2con,
154 const gmx_ga2la_t& ga2la,
155 gmx_bool bHomeConnect,
156 gmx_domdec_constraints_t* dc,
157 gmx_domdec_specat_comm_t* dcc,
159 std::vector<int>* ireq)
161 if (!dc->gc_req[con_offset + con])
163 /* Add this non-home constraint to the list */
164 dc->con_gl.push_back(con_offset + con);
165 dc->con_nlocat.push_back(bHomeConnect ? 1 : 0);
166 dc->gc_req[con_offset + con] = true;
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 const int* iap = constr_iatomptr(ia1, ia2, con);
173 il_local->iatoms[il_local->nr++] = iap[0];
174 const int a1_gl = offset + iap[1];
175 const int a2_gl = offset + iap[2];
176 /* The following indexing code can probably be optizimed */
177 if (const int* a_loc = ga2la.findHome(a1_gl))
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 (const int* a_loc = ga2la.findHome(a2_gl))
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 (!dc->ga2la->find(offset + a))
201 /* Add this non-home atom to the list */
202 ireq->push_back(offset + a);
203 /* Temporarily mark with -2, we get the index later */
204 dc->ga2la->insert(offset + a, -2);
209 /* Loop over the constraint connected to atom a */
210 for (const int coni : at2con[a])
215 const int* iap = constr_iatomptr(ia1, ia2, coni);
225 if (!ga2la.findHome(offset + b))
227 walk_out(coni, con_offset, b, offset, nrec - 1, ia1, ia2, at2con, ga2la, FALSE,
228 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 gmx::ArrayRef<const std::vector<int>> at2settle_mt,
243 std::vector<int>* ireq)
245 const gmx_ga2la_t& ga2la = *dd->ga2la;
246 int nral = NRAL(F_SETTLE);
249 for (int a = cg_start; a < cg_end; a++)
251 if (GET_CGINFO_SETTLE(cginfo[a]))
253 int a_gl = dd->globalAtomIndices[a];
255 mtopGetMolblockIndex(mtop, a_gl, &mb, nullptr, &a_mol);
257 const gmx_molblock_t* molb = &mtop->molblock[mb];
258 int settle = at2settle_mt[molb->type][a_mol];
262 int offset = a_gl - a_mol;
264 const int* ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms.data();
267 gmx_bool bAssign = FALSE;
269 for (int sa = 0; sa < nral; sa++)
271 int a_glsa = offset + ia1[settle * (1 + nral) + 1 + sa];
273 if (ga2la.findHome(a_glsa))
275 if (nlocal == 0 && a_gl == a_glsa)
285 if (ils_local->nr + 1 + nral > ils_local->nalloc)
287 ils_local->nalloc = over_alloc_dd(ils_local->nr + 1 + nral);
288 srenew(ils_local->iatoms, ils_local->nalloc);
291 ils_local->iatoms[ils_local->nr++] = ia1[settle * 4];
293 for (int sa = 0; sa < nral; sa++)
295 if (const int* a_loc = ga2la.findHome(a_gls[sa]))
297 ils_local->iatoms[ils_local->nr++] = *a_loc;
301 ils_local->iatoms[ils_local->nr++] = -a_gls[sa] - 1;
302 /* Add this non-home atom to the list */
303 ireq->push_back(a_gls[sa]);
304 /* A check on double atom requests is
305 * not required for settle.
315 /*! \brief Looks up constraint for the local atoms */
316 static void atoms_to_constraints(gmx_domdec_t* dd,
317 const gmx_mtop_t* mtop,
319 gmx::ArrayRef<const ListOfLists<int>> at2con_mt,
322 std::vector<int>* ireq)
324 gmx_domdec_constraints_t* dc = dd->constraints;
325 gmx_domdec_specat_comm_t* dcc = dd->constraint_comm;
327 const gmx_ga2la_t& ga2la = *dd->ga2la;
330 dc->con_nlocat.clear();
334 for (int a = 0; a < dd->ncg_home; a++)
336 if (GET_CGINFO_CONSTR(cginfo[a]))
338 int a_gl = dd->globalAtomIndices[a];
340 mtopGetMolblockIndex(mtop, a_gl, &mb, &molnr, &a_mol);
342 const gmx_molblock_t& molb = mtop->molblock[mb];
344 gmx::ArrayRef<const int> ia1 = mtop->moltype[molb.type].ilist[F_CONSTR].iatoms;
345 gmx::ArrayRef<const int> ia2 = mtop->moltype[molb.type].ilist[F_CONSTRNC].iatoms;
347 /* Calculate the global constraint number offset for the molecule.
348 * This is only required for the global index to make sure
349 * that we use each constraint only once.
351 const int con_offset = dc->molb_con_offset[mb] + molnr * dc->molb_ncon_mol[mb];
353 /* The global atom number offset for this molecule */
354 const int offset = a_gl - a_mol;
355 /* Loop over the constraints connected to atom a_mol in the molecule */
356 const auto& at2con = at2con_mt[molb.type];
357 for (const int con : at2con[a_mol])
359 const int* iap = constr_iatomptr(ia1, ia2, con);
369 if (const int* a_loc = ga2la.findHome(offset + b_mol))
371 /* Add this fully home constraint at the first atom */
374 dc->con_gl.push_back(con_offset + con);
375 dc->con_nlocat.push_back(2);
376 if (ilc_local->nr + 3 > ilc_local->nalloc)
378 ilc_local->nalloc = over_alloc_dd(ilc_local->nr + 3);
379 srenew(ilc_local->iatoms, ilc_local->nalloc);
381 const int b_lo = *a_loc;
382 ilc_local->iatoms[ilc_local->nr++] = iap[0];
383 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? a : b_lo);
384 ilc_local->iatoms[ilc_local->nr++] = (a_gl == iap[1] ? b_lo : a);
391 /* We need the nrec constraints coupled to this constraint,
392 * so we need to walk out of the home cell by nrec+1 atoms,
393 * since already atom bg is not locally present.
394 * Therefore we call walk_out with nrec recursions to go
395 * after this first call.
397 walk_out(con, con_offset, b_mol, offset, nrec, ia1, ia2, at2con, ga2la, TRUE,
398 dc, dcc, ilc_local, ireq);
404 GMX_ASSERT(dc->con_gl.size() == static_cast<size_t>(dc->ncon),
405 "con_gl size should match the number of constraints");
406 GMX_ASSERT(dc->con_nlocat.size() == static_cast<size_t>(dc->ncon),
407 "con_nlocat size should match the number of constraints");
411 fprintf(debug, "Constraints: home %3d border %3d atoms: %3zu\n", nhome, dc->ncon - nhome,
412 dd->constraint_comm ? ireq->size() : 0);
416 int dd_make_local_constraints(gmx_domdec_t* dd,
418 const struct gmx_mtop_t* mtop,
420 gmx::Constraints* constr,
424 gmx_domdec_constraints_t* dc;
425 t_ilist * ilc_local, *ils_local;
426 gmx::HashedMap<int>* ga2la_specat;
430 // This code should not be called unless this condition is true,
431 // because that's the only time init_domdec_constraints is
433 GMX_RELEASE_ASSERT(dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles,
434 "dd_make_local_constraints called when there are no local constraints");
435 // ... and init_domdec_constraints always sets
436 // dd->constraint_comm...
439 "Invalid use of dd_make_local_constraints before construction of constraint_comm");
440 // ... which static analysis needs to be reassured about, because
441 // otherwise, when dd->splitSettles is
442 // true. dd->constraint_comm is unilaterally dereferenced before
443 // the call to atoms_to_settles.
445 dc = dd->constraints;
447 ilc_local = &il_local[F_CONSTR];
448 ils_local = &il_local[F_SETTLE];
452 gmx::ArrayRef<const ListOfLists<int>> at2con_mt;
453 std::vector<int>* ireq = nullptr;
454 if (dd->constraint_comm)
456 // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
457 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
458 at2con_mt = constr->atom2constraints_moltype();
459 ireq = &dc->requestedGlobalAtomIndices[0];
463 gmx::ArrayRef<const std::vector<int>> at2settle_mt;
464 /* When settle works inside charge groups, we assigned them already */
465 if (dd->comm->systemInfo.haveSplitSettles)
467 // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
468 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
469 at2settle_mt = constr->atom2settle_moltype();
473 if (at2settle_mt.empty())
475 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec, ilc_local, ireq);
481 /* Do the constraints, if present, on the first thread.
482 * Do the settles on all other threads.
484 t0_set = ((!at2con_mt.empty() && dc->nthread > 1) ? 1 : 0);
486 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
487 for (int thread = 0; thread < dc->nthread; thread++)
491 if (!at2con_mt.empty() && thread == 0)
493 atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec, ilc_local, ireq);
496 if (thread >= t0_set)
501 /* Distribute the settle check+assignments over
502 * dc->nthread or dc->nthread-1 threads.
504 cg0 = (dd->ncg_home * (thread - t0_set)) / (dc->nthread - t0_set);
505 cg1 = (dd->ncg_home * (thread - t0_set + 1)) / (dc->nthread - t0_set);
507 if (thread == t0_set)
513 ilst = &dc->ils[thread];
517 std::vector<int>& ireqt = dc->requestedGlobalAtomIndices[thread];
523 atoms_to_settles(dd, mtop, cginfo, at2settle_mt, cg0, cg1, ilst, &ireqt);
526 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
529 /* Combine the generate settles and requested indices */
530 for (int thread = 1; thread < dc->nthread; thread++)
537 ilst = &dc->ils[thread];
538 if (ils_local->nr + ilst->nr > ils_local->nalloc)
540 ils_local->nalloc = over_alloc_large(ils_local->nr + ilst->nr);
541 srenew(ils_local->iatoms, ils_local->nalloc);
543 for (ia = 0; ia < ilst->nr; ia++)
545 ils_local->iatoms[ils_local->nr + ia] = ilst->iatoms[ia];
547 ils_local->nr += ilst->nr;
550 const std::vector<int>& ireqt = dc->requestedGlobalAtomIndices[thread];
551 ireq->insert(ireq->end(), ireqt.begin(), ireqt.end());
556 fprintf(debug, "Settles: total %3d\n", ils_local->nr / 4);
560 if (dd->constraint_comm)
564 at_end = setup_specat_communication(dd, ireq, dd->constraint_comm, dd->constraints->ga2la.get(),
565 at_start, 2, "constraint", " or lincs-order");
567 /* Fill in the missing indices */
568 ga2la_specat = dd->constraints->ga2la.get();
570 nral1 = 1 + NRAL(F_CONSTR);
571 for (i = 0; i < ilc_local->nr; i += nral1)
573 iap = ilc_local->iatoms + i;
574 for (j = 1; j < nral1; j++)
578 const int* a = ga2la_specat->find(-iap[j] - 1);
579 GMX_ASSERT(a, "We have checked before that this atom index has been set");
585 nral1 = 1 + NRAL(F_SETTLE);
586 for (i = 0; i < ils_local->nr; i += nral1)
588 iap = ils_local->iatoms + i;
589 for (j = 1; j < nral1; j++)
593 const int* a = ga2la_specat->find(-iap[j] - 1);
594 GMX_ASSERT(a, "We have checked before that this atom index has been set");
602 // Currently unreachable
609 void init_domdec_constraints(gmx_domdec_t* dd, const gmx_mtop_t* mtop)
611 gmx_domdec_constraints_t* dc;
612 const gmx_molblock_t* molb;
616 fprintf(debug, "Begin init_domdec_constraints\n");
619 dd->constraints = new gmx_domdec_constraints_t;
620 dc = dd->constraints;
622 dc->molb_con_offset.resize(mtop->molblock.size());
623 dc->molb_ncon_mol.resize(mtop->molblock.size());
626 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
628 molb = &mtop->molblock[mb];
629 dc->molb_con_offset[mb] = ncon;
630 dc->molb_ncon_mol[mb] = mtop->moltype[molb->type].ilist[F_CONSTR].size() / 3
631 + mtop->moltype[molb->type].ilist[F_CONSTRNC].size() / 3;
632 ncon += molb->nmol * dc->molb_ncon_mol[mb];
637 dc->gc_req.resize(ncon);
640 /* Use a hash table for the global to local index.
641 * The number of keys is a rough estimate, it will be optimized later.
643 int numKeysEstimate = std::min(mtop->natoms / 20, mtop->natoms / (2 * dd->nnodes));
644 dc->ga2la = std::make_unique<gmx::HashedMap<int>>(numKeysEstimate);
646 dc->nthread = gmx_omp_nthreads_get(emntDomdec);
647 dc->ils.resize(dc->nthread);
649 dd->constraint_comm = new gmx_domdec_specat_comm_t;
651 dc->requestedGlobalAtomIndices.resize(dc->nthread);