2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005,2006,2007,2008,2009 by the GROMACS development team.
5 * Copyright (c) 2010,2011,2012,2013,2014 by the GROMACS development team.
6 * Copyright (c) 2015,2016,2017,2018,2019 by the GROMACS development team.
7 * Copyright (c) 2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
56 #include "gromacs/domdec/builder.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/computemultibodycutoffs.h"
59 #include "gromacs/domdec/dlb.h"
60 #include "gromacs/domdec/dlbtiming.h"
61 #include "gromacs/domdec/domdec_network.h"
62 #include "gromacs/domdec/ga2la.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localtopologychecker.h"
65 #include "gromacs/domdec/options.h"
66 #include "gromacs/domdec/partition.h"
67 #include "gromacs/ewald/pme.h"
68 #include "gromacs/domdec/reversetopology.h"
69 #include "gromacs/gmxlib/network.h"
70 #include "gromacs/gmxlib/nrnb.h"
71 #include "gromacs/gpu_utils/device_stream_manager.h"
72 #include "gromacs/gpu_utils/gpu_utils.h"
73 #include "gromacs/hardware/hw_info.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/math/vectypes.h"
76 #include "gromacs/mdlib/calc_verletbuf.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/constraintrange.h"
79 #include "gromacs/mdlib/updategroups.h"
80 #include "gromacs/mdlib/vcm.h"
81 #include "gromacs/mdlib/vsite.h"
82 #include "gromacs/mdrun/mdmodules.h"
83 #include "gromacs/mdtypes/commrec.h"
84 #include "gromacs/mdtypes/forceoutput.h"
85 #include "gromacs/mdtypes/forcerec.h"
86 #include "gromacs/mdtypes/inputrec.h"
87 #include "gromacs/mdtypes/mdrunoptions.h"
88 #include "gromacs/mdtypes/state.h"
89 #include "gromacs/pbcutil/ishift.h"
90 #include "gromacs/pbcutil/pbc.h"
91 #include "gromacs/pulling/pull.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/topology/block.h"
94 #include "gromacs/topology/idef.h"
95 #include "gromacs/topology/ifunc.h"
96 #include "gromacs/topology/mtop_lookup.h"
97 #include "gromacs/topology/mtop_util.h"
98 #include "gromacs/topology/topology.h"
99 #include "gromacs/utility/basedefinitions.h"
100 #include "gromacs/utility/basenetwork.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/enumerationhelpers.h"
103 #include "gromacs/utility/exceptions.h"
104 #include "gromacs/utility/fatalerror.h"
105 #include "gromacs/utility/gmxmpi.h"
106 #include "gromacs/utility/logger.h"
107 #include "gromacs/utility/mdmodulesnotifiers.h"
108 #include "gromacs/utility/real.h"
109 #include "gromacs/utility/smalloc.h"
110 #include "gromacs/utility/strconvert.h"
111 #include "gromacs/utility/stringstream.h"
112 #include "gromacs/utility/stringutil.h"
113 #include "gromacs/utility/textwriter.h"
115 #include "atomdistribution.h"
117 #include "cellsizes.h"
118 #include "distribute.h"
119 #include "domdec_constraints.h"
120 #include "domdec_internal.h"
121 #include "domdec_setup.h"
122 #include "domdec_vsite.h"
123 #include "redistribute.h"
126 // TODO remove this when moving domdec into gmx namespace
128 using gmx::DdRankOrder;
129 using gmx::DlbOption;
130 using gmx::DomdecOptions;
131 using gmx::RangePartitioning;
133 static const char* enumValueToString(DlbState enumValue)
135 static constexpr gmx::EnumerationArray<DlbState, const char*> dlbStateNames = {
136 "off", "off", "auto", "locked", "on", "on"
138 return dlbStateNames[enumValue];
141 /* The DD zone order */
142 static const ivec dd_zo[DD_MAXZONE] = { { 0, 0, 0 }, { 1, 0, 0 }, { 1, 1, 0 }, { 0, 1, 0 },
143 { 0, 1, 1 }, { 0, 0, 1 }, { 1, 0, 1 }, { 1, 1, 1 } };
145 /* The non-bonded zone-pair setup for domain decomposition
146 * The first number is the i-zone, the second number the first j-zone seen by
147 * this i-zone, the third number the last+1 j-zone seen by this i-zone.
148 * As is, this is for 3D decomposition, where there are 4 i-zones.
149 * With 2D decomposition use only the first 2 i-zones and a last+1 j-zone of 4.
150 * With 1D decomposition use only the first i-zone and a last+1 j-zone of 2.
152 static const int ddNonbondedZonePairRanges[DD_MAXIZONE][3] = { { 0, 0, 8 },
157 static void ddindex2xyz(const ivec nc, int ind, ivec xyz)
159 xyz[XX] = ind / (nc[YY] * nc[ZZ]);
160 xyz[YY] = (ind / nc[ZZ]) % nc[YY];
161 xyz[ZZ] = ind % nc[ZZ];
164 static int ddcoord2ddnodeid(gmx_domdec_t* dd, ivec c)
168 const CartesianRankSetup& cartSetup = dd->comm->cartesianRankSetup;
169 const int ddindex = dd_index(dd->numCells, c);
170 if (cartSetup.bCartesianPP_PME)
172 ddnodeid = cartSetup.ddindex2ddnodeid[ddindex];
174 else if (cartSetup.bCartesianPP)
177 MPI_Cart_rank(dd->mpi_comm_all, c, &ddnodeid);
188 int ddglatnr(const gmx_domdec_t* dd, int i)
198 if (i >= dd->comm->atomRanges.numAtomsTotal())
201 "glatnr called with %d, which is larger than the local number of atoms (%d)",
203 dd->comm->atomRanges.numAtomsTotal());
205 atnr = dd->globalAtomIndices[i] + 1;
211 void dd_store_state(const gmx_domdec_t& dd, t_state* state)
213 if (state->ddp_count != dd.ddp_count)
215 gmx_incons("The MD state does not match the domain decomposition state");
218 state->cg_gl.resize(dd.numHomeAtoms);
219 for (int i = 0; i < dd.numHomeAtoms; i++)
221 state->cg_gl[i] = dd.globalAtomGroupIndices[i];
224 state->ddp_count_cg_gl = dd.ddp_count;
227 gmx_domdec_zones_t* domdec_zones(gmx_domdec_t* dd)
229 return &dd->comm->zones;
232 int dd_numAtomsZones(const gmx_domdec_t& dd)
234 return dd.comm->atomRanges.end(DDAtomRanges::Type::Zones);
237 int dd_numHomeAtoms(const gmx_domdec_t& dd)
239 return dd.comm->atomRanges.numHomeAtoms();
242 int dd_natoms_mdatoms(const gmx_domdec_t& dd)
244 /* We currently set mdatoms entries for all atoms:
245 * local + non-local + communicated for vsite + constraints
248 return dd.comm->atomRanges.numAtomsTotal();
251 int dd_natoms_vsite(const gmx_domdec_t& dd)
253 return dd.comm->atomRanges.end(DDAtomRanges::Type::Vsites);
256 void dd_get_constraint_range(const gmx_domdec_t& dd, int* at_start, int* at_end)
258 *at_start = dd.comm->atomRanges.start(DDAtomRanges::Type::Constraints);
259 *at_end = dd.comm->atomRanges.end(DDAtomRanges::Type::Constraints);
262 void dd_move_x(gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle)
264 wallcycle_start(wcycle, WallCycleCounter::MoveX);
266 rvec shift = { 0, 0, 0 };
268 gmx_domdec_comm_t* comm = dd->comm;
271 int nat_tot = comm->atomRanges.numHomeAtoms();
272 for (int d = 0; d < dd->ndim; d++)
274 const bool bPBC = (dd->ci[dd->dim[d]] == 0);
275 const bool bScrew = (bPBC && dd->unitCellInfo.haveScrewPBC && dd->dim[d] == XX);
278 copy_rvec(box[dd->dim[d]], shift);
280 gmx_domdec_comm_dim_t* cd = &comm->cd[d];
281 for (const gmx_domdec_ind_t& ind : cd->ind)
283 DDBufferAccess<gmx::RVec> sendBufferAccess(comm->rvecBuffer, ind.nsend[nzone + 1]);
284 gmx::ArrayRef<gmx::RVec>& sendBuffer = sendBufferAccess.buffer;
288 for (int j : ind.index)
290 sendBuffer[n] = x[j];
296 for (int j : ind.index)
298 /* We need to shift the coordinates */
299 for (int d = 0; d < DIM; d++)
301 sendBuffer[n][d] = x[j][d] + shift[d];
308 for (int j : ind.index)
311 sendBuffer[n][XX] = x[j][XX] + shift[XX];
313 * This operation requires a special shift force
314 * treatment, which is performed in calc_vir.
316 sendBuffer[n][YY] = box[YY][YY] - x[j][YY];
317 sendBuffer[n][ZZ] = box[ZZ][ZZ] - x[j][ZZ];
322 DDBufferAccess<gmx::RVec> receiveBufferAccess(
323 comm->rvecBuffer2, cd->receiveInPlace ? 0 : ind.nrecv[nzone + 1]);
325 gmx::ArrayRef<gmx::RVec> receiveBuffer;
326 if (cd->receiveInPlace)
328 receiveBuffer = gmx::arrayRefFromArray(x.data() + nat_tot, ind.nrecv[nzone + 1]);
332 receiveBuffer = receiveBufferAccess.buffer;
334 /* Send and receive the coordinates */
335 ddSendrecv(dd, d, dddirBackward, sendBuffer, receiveBuffer);
337 if (!cd->receiveInPlace)
340 for (int zone = 0; zone < nzone; zone++)
342 for (int i = ind.cell2at0[zone]; i < ind.cell2at1[zone]; i++)
344 x[i] = receiveBuffer[j++];
348 nat_tot += ind.nrecv[nzone + 1];
353 wallcycle_stop(wcycle, WallCycleCounter::MoveX);
356 void dd_move_f(gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle)
358 wallcycle_start(wcycle, WallCycleCounter::MoveF);
360 gmx::ArrayRef<gmx::RVec> f = forceWithShiftForces->force();
361 gmx::ArrayRef<gmx::RVec> fshift = forceWithShiftForces->shiftForces();
363 gmx_domdec_comm_t& comm = *dd->comm;
364 int nzone = comm.zones.n / 2;
365 int nat_tot = comm.atomRanges.end(DDAtomRanges::Type::Zones);
366 for (int d = dd->ndim - 1; d >= 0; d--)
368 /* Only forces in domains near the PBC boundaries need to
369 consider PBC in the treatment of fshift */
370 const bool shiftForcesNeedPbc =
371 (forceWithShiftForces->computeVirial() && dd->ci[dd->dim[d]] == 0);
372 const bool applyScrewPbc =
373 (shiftForcesNeedPbc && dd->unitCellInfo.haveScrewPBC && dd->dim[d] == XX);
374 /* Determine which shift vector we need */
375 ivec vis = { 0, 0, 0 };
377 const int is = gmx::ivecToShiftIndex(vis);
379 /* Loop over the pulses */
380 const gmx_domdec_comm_dim_t& cd = comm.cd[d];
381 for (int p = cd.numPulses() - 1; p >= 0; p--)
383 const gmx_domdec_ind_t& ind = cd.ind[p];
384 DDBufferAccess<gmx::RVec> receiveBufferAccess(comm.rvecBuffer, ind.nsend[nzone + 1]);
385 gmx::ArrayRef<gmx::RVec>& receiveBuffer = receiveBufferAccess.buffer;
387 nat_tot -= ind.nrecv[nzone + 1];
389 DDBufferAccess<gmx::RVec> sendBufferAccess(
390 comm.rvecBuffer2, cd.receiveInPlace ? 0 : ind.nrecv[nzone + 1]);
392 gmx::ArrayRef<gmx::RVec> sendBuffer;
393 if (cd.receiveInPlace)
395 sendBuffer = gmx::arrayRefFromArray(f.data() + nat_tot, ind.nrecv[nzone + 1]);
399 sendBuffer = sendBufferAccess.buffer;
401 for (int zone = 0; zone < nzone; zone++)
403 for (int i = ind.cell2at0[zone]; i < ind.cell2at1[zone]; i++)
405 sendBuffer[j++] = f[i];
409 /* Communicate the forces */
410 ddSendrecv(dd, d, dddirForward, sendBuffer, receiveBuffer);
411 /* Add the received forces */
413 if (!shiftForcesNeedPbc)
415 for (int j : ind.index)
417 for (int d = 0; d < DIM; d++)
419 f[j][d] += receiveBuffer[n][d];
424 else if (!applyScrewPbc)
426 for (int j : ind.index)
428 for (int d = 0; d < DIM; d++)
430 f[j][d] += receiveBuffer[n][d];
432 /* Add this force to the shift force */
433 for (int d = 0; d < DIM; d++)
435 fshift[is][d] += receiveBuffer[n][d];
442 for (int j : ind.index)
444 /* Rotate the force */
445 f[j][XX] += receiveBuffer[n][XX];
446 f[j][YY] -= receiveBuffer[n][YY];
447 f[j][ZZ] -= receiveBuffer[n][ZZ];
448 if (shiftForcesNeedPbc)
450 /* Add this force to the shift force */
451 for (int d = 0; d < DIM; d++)
453 fshift[is][d] += receiveBuffer[n][d];
462 wallcycle_stop(wcycle, WallCycleCounter::MoveF);
465 /* Convenience function for extracting a real buffer from an rvec buffer
467 * To reduce the number of temporary communication buffers and avoid
468 * cache polution, we reuse gmx::RVec buffers for storing reals.
469 * This functions return a real buffer reference with the same number
470 * of elements as the gmx::RVec buffer (so 1/3 of the size in bytes).
472 static gmx::ArrayRef<real> realArrayRefFromRvecArrayRef(gmx::ArrayRef<gmx::RVec> arrayRef)
474 return gmx::arrayRefFromArray(as_rvec_array(arrayRef.data())[0], arrayRef.size());
477 void dd_atom_spread_real(gmx_domdec_t* dd, real v[])
479 gmx_domdec_comm_t* comm = dd->comm;
482 int nat_tot = comm->atomRanges.numHomeAtoms();
483 for (int d = 0; d < dd->ndim; d++)
485 gmx_domdec_comm_dim_t* cd = &comm->cd[d];
486 for (const gmx_domdec_ind_t& ind : cd->ind)
488 /* Note: We provision for RVec instead of real, so a factor of 3
489 * more than needed. The buffer actually already has this size
490 * and we pass a plain pointer below, so this does not matter.
492 DDBufferAccess<gmx::RVec> sendBufferAccess(comm->rvecBuffer, ind.nsend[nzone + 1]);
493 gmx::ArrayRef<real> sendBuffer = realArrayRefFromRvecArrayRef(sendBufferAccess.buffer);
495 for (int j : ind.index)
497 sendBuffer[n++] = v[j];
500 DDBufferAccess<gmx::RVec> receiveBufferAccess(
501 comm->rvecBuffer2, cd->receiveInPlace ? 0 : ind.nrecv[nzone + 1]);
503 gmx::ArrayRef<real> receiveBuffer;
504 if (cd->receiveInPlace)
506 receiveBuffer = gmx::arrayRefFromArray(v + nat_tot, ind.nrecv[nzone + 1]);
510 receiveBuffer = realArrayRefFromRvecArrayRef(receiveBufferAccess.buffer);
512 /* Send and receive the data */
513 ddSendrecv(dd, d, dddirBackward, sendBuffer, receiveBuffer);
514 if (!cd->receiveInPlace)
517 for (int zone = 0; zone < nzone; zone++)
519 for (int i = ind.cell2at0[zone]; i < ind.cell2at1[zone]; i++)
521 v[i] = receiveBuffer[j++];
525 nat_tot += ind.nrecv[nzone + 1];
531 void dd_atom_sum_real(gmx_domdec_t* dd, real v[])
533 gmx_domdec_comm_t* comm = dd->comm;
535 int nzone = comm->zones.n / 2;
536 int nat_tot = comm->atomRanges.end(DDAtomRanges::Type::Zones);
537 for (int d = dd->ndim - 1; d >= 0; d--)
539 gmx_domdec_comm_dim_t* cd = &comm->cd[d];
540 for (int p = cd->numPulses() - 1; p >= 0; p--)
542 const gmx_domdec_ind_t& ind = cd->ind[p];
544 /* Note: We provision for RVec instead of real, so a factor of 3
545 * more than needed. The buffer actually already has this size
546 * and we typecast, so this works as intended.
548 DDBufferAccess<gmx::RVec> receiveBufferAccess(comm->rvecBuffer, ind.nsend[nzone + 1]);
549 gmx::ArrayRef<real> receiveBuffer = realArrayRefFromRvecArrayRef(receiveBufferAccess.buffer);
550 nat_tot -= ind.nrecv[nzone + 1];
552 DDBufferAccess<gmx::RVec> sendBufferAccess(
553 comm->rvecBuffer2, cd->receiveInPlace ? 0 : ind.nrecv[nzone + 1]);
555 gmx::ArrayRef<real> sendBuffer;
556 if (cd->receiveInPlace)
558 sendBuffer = gmx::arrayRefFromArray(v + nat_tot, ind.nrecv[nzone + 1]);
562 sendBuffer = realArrayRefFromRvecArrayRef(sendBufferAccess.buffer);
564 for (int zone = 0; zone < nzone; zone++)
566 for (int i = ind.cell2at0[zone]; i < ind.cell2at1[zone]; i++)
568 sendBuffer[j++] = v[i];
572 /* Communicate the forces */
573 ddSendrecv(dd, d, dddirForward, sendBuffer, receiveBuffer);
574 /* Add the received forces */
576 for (int j : ind.index)
578 v[j] += receiveBuffer[n];
586 real dd_cutoff_multibody(const gmx_domdec_t* dd)
588 const gmx_domdec_comm_t& comm = *dd->comm;
589 const DDSystemInfo& systemInfo = comm.systemInfo;
592 if (systemInfo.haveInterDomainMultiBodyBondeds)
594 if (comm.cutoff_mbody > 0)
596 r = comm.cutoff_mbody;
600 /* cutoff_mbody=0 means we do not have DLB */
601 r = comm.cellsize_min[dd->dim[0]];
602 for (int di = 1; di < dd->ndim; di++)
604 r = std::min(r, comm.cellsize_min[dd->dim[di]]);
606 if (comm.systemInfo.filterBondedCommunication)
608 r = std::max(r, comm.cutoff_mbody);
612 r = std::min(r, systemInfo.cutoff);
620 real dd_cutoff_twobody(const gmx_domdec_t* dd)
622 const real r_mb = dd_cutoff_multibody(dd);
624 return std::max(dd->comm->systemInfo.cutoff, r_mb);
628 static void dd_cart_coord2pmecoord(const DDRankSetup& ddRankSetup,
629 const CartesianRankSetup& cartSetup,
633 const int nc = ddRankSetup.numPPCells[cartSetup.cartpmedim];
634 const int ntot = cartSetup.ntot[cartSetup.cartpmedim];
635 copy_ivec(coord, coord_pme);
636 coord_pme[cartSetup.cartpmedim] =
637 nc + (coord[cartSetup.cartpmedim] * (ntot - nc) + (ntot - nc) / 2) / nc;
640 /* Returns the PME rank index in 0...npmenodes-1 for the PP cell with index ddCellIndex */
641 static int ddindex2pmeindex(const DDRankSetup& ddRankSetup, const int ddCellIndex)
643 const int npp = ddRankSetup.numPPRanks;
644 const int npme = ddRankSetup.numRanksDoingPme;
646 /* Here we assign a PME node to communicate with this DD node
647 * by assuming that the major index of both is x.
648 * We add npme/2 to obtain an even distribution.
650 return (ddCellIndex * npme + npme / 2) / npp;
653 static std::vector<int> dd_interleaved_pme_ranks(const DDRankSetup& ddRankSetup)
655 std::vector<int> pmeRanks(ddRankSetup.numRanksDoingPme);
658 for (int i = 0; i < ddRankSetup.numPPRanks; i++)
660 const int p0 = ddindex2pmeindex(ddRankSetup, i);
661 const int p1 = ddindex2pmeindex(ddRankSetup, i + 1);
662 if (i + 1 == ddRankSetup.numPPRanks || p1 > p0)
666 fprintf(debug, "pme_rank[%d] = %d\n", n, i + 1 + n);
668 pmeRanks[n] = i + 1 + n;
676 static int gmx_ddcoord2pmeindex(const gmx_domdec_t& dd, int x, int y, int z)
683 const int slab = ddindex2pmeindex(dd.comm->ddRankSetup, dd_index(dd.numCells, coords));
688 static int ddcoord2simnodeid(const t_commrec* cr, int x, int y, int z)
690 const CartesianRankSetup& cartSetup = cr->dd->comm->cartesianRankSetup;
691 ivec coords = { x, y, z };
694 if (cartSetup.bCartesianPP_PME)
697 MPI_Cart_rank(cr->mpi_comm_mysim, coords, &nodeid);
702 int ddindex = dd_index(cr->dd->numCells, coords);
703 if (cartSetup.bCartesianPP)
705 nodeid = cartSetup.ddindex2simnodeid[ddindex];
709 const DDRankSetup& rankSetup = cr->dd->comm->ddRankSetup;
710 if (rankSetup.rankOrder != DdRankOrder::pp_pme && rankSetup.usePmeOnlyRanks)
712 nodeid = ddindex + gmx_ddcoord2pmeindex(*cr->dd, x, y, z);
724 static int dd_simnode2pmenode(const DDRankSetup& ddRankSetup,
725 const CartesianRankSetup& cartSetup,
726 gmx::ArrayRef<const int> pmeRanks,
727 const t_commrec gmx_unused* cr,
728 const int sim_nodeid)
732 /* This assumes a uniform x domain decomposition grid cell size */
733 if (cartSetup.bCartesianPP_PME)
736 ivec coord, coord_pme;
737 MPI_Cart_coords(cr->mpi_comm_mysim, sim_nodeid, DIM, coord);
738 if (coord[cartSetup.cartpmedim] < ddRankSetup.numPPCells[cartSetup.cartpmedim])
740 /* This is a PP rank */
741 dd_cart_coord2pmecoord(ddRankSetup, cartSetup, coord, coord_pme);
742 MPI_Cart_rank(cr->mpi_comm_mysim, coord_pme, &pmenode);
746 else if (cartSetup.bCartesianPP)
748 if (sim_nodeid < ddRankSetup.numPPRanks)
750 pmenode = ddRankSetup.numPPRanks + ddindex2pmeindex(ddRankSetup, sim_nodeid);
755 /* This assumes DD cells with identical x coordinates
756 * are numbered sequentially.
758 if (pmeRanks.empty())
760 if (sim_nodeid < ddRankSetup.numPPRanks)
762 /* The DD index equals the nodeid */
763 pmenode = ddRankSetup.numPPRanks + ddindex2pmeindex(ddRankSetup, sim_nodeid);
769 while (sim_nodeid > pmeRanks[i])
773 if (sim_nodeid < pmeRanks[i])
775 pmenode = pmeRanks[i];
783 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd)
787 return { dd->comm->ddRankSetup.npmenodes_x, dd->comm->ddRankSetup.npmenodes_y };
795 std::vector<int> get_pme_ddranks(const t_commrec* cr, const int pmenodeid)
797 const DDRankSetup& ddRankSetup = cr->dd->comm->ddRankSetup;
798 const CartesianRankSetup& cartSetup = cr->dd->comm->cartesianRankSetup;
799 GMX_RELEASE_ASSERT(ddRankSetup.usePmeOnlyRanks,
800 "This function should only be called when PME-only ranks are in use");
801 std::vector<int> ddranks;
802 ddranks.reserve((ddRankSetup.numPPRanks + ddRankSetup.numRanksDoingPme - 1) / ddRankSetup.numRanksDoingPme);
804 for (int x = 0; x < ddRankSetup.numPPCells[XX]; x++)
806 for (int y = 0; y < ddRankSetup.numPPCells[YY]; y++)
808 for (int z = 0; z < ddRankSetup.numPPCells[ZZ]; z++)
810 if (cartSetup.bCartesianPP_PME)
812 ivec coord = { x, y, z };
814 dd_cart_coord2pmecoord(ddRankSetup, cartSetup, coord, coord_pme);
815 if (cr->dd->ci[XX] == coord_pme[XX] && cr->dd->ci[YY] == coord_pme[YY]
816 && cr->dd->ci[ZZ] == coord_pme[ZZ])
818 ddranks.push_back(ddcoord2simnodeid(cr, x, y, z));
823 /* The slab corresponds to the nodeid in the PME group */
824 if (gmx_ddcoord2pmeindex(*cr->dd, x, y, z) == pmenodeid)
826 ddranks.push_back(ddcoord2simnodeid(cr, x, y, z));
835 static gmx_bool receive_vir_ener(const gmx_domdec_t* dd, gmx::ArrayRef<const int> pmeRanks, const t_commrec* cr)
837 bool bReceive = true;
839 const DDRankSetup& ddRankSetup = dd->comm->ddRankSetup;
840 if (ddRankSetup.usePmeOnlyRanks)
842 const CartesianRankSetup& cartSetup = dd->comm->cartesianRankSetup;
843 if (cartSetup.bCartesianPP_PME)
846 int pmenode = dd_simnode2pmenode(ddRankSetup, cartSetup, pmeRanks, cr, cr->sim_nodeid);
848 MPI_Cart_coords(cr->mpi_comm_mysim, cr->sim_nodeid, DIM, coords);
849 coords[cartSetup.cartpmedim]++;
850 if (coords[cartSetup.cartpmedim] < dd->numCells[cartSetup.cartpmedim])
853 MPI_Cart_rank(cr->mpi_comm_mysim, coords, &rank);
854 if (dd_simnode2pmenode(ddRankSetup, cartSetup, pmeRanks, cr, rank) == pmenode)
856 /* This is not the last PP node for pmenode */
863 "Without MPI we should not have Cartesian PP-PME with #PMEnodes < #DDnodes");
868 int pmenode = dd_simnode2pmenode(ddRankSetup, cartSetup, pmeRanks, cr, cr->sim_nodeid);
869 if (cr->sim_nodeid + 1 < cr->nnodes
870 && dd_simnode2pmenode(ddRankSetup, cartSetup, pmeRanks, cr, cr->sim_nodeid + 1) == pmenode)
872 /* This is not the last PP node for pmenode */
881 static void set_slb_pme_dim_f(gmx_domdec_t* dd, int dim, real** dim_f)
883 gmx_domdec_comm_t* comm = dd->comm;
885 snew(*dim_f, dd->numCells[dim] + 1);
887 for (int i = 1; i < dd->numCells[dim]; i++)
889 if (comm->slb_frac[dim])
891 (*dim_f)[i] = (*dim_f)[i - 1] + comm->slb_frac[dim][i - 1];
895 (*dim_f)[i] = static_cast<real>(i) / static_cast<real>(dd->numCells[dim]);
898 (*dim_f)[dd->numCells[dim]] = 1;
901 static void init_ddpme(gmx_domdec_t* dd, gmx_ddpme_t* ddpme, int dimind)
903 const DDRankSetup& ddRankSetup = dd->comm->ddRankSetup;
905 if (dimind == 0 && dd->dim[0] == YY && ddRankSetup.npmenodes_x == 1)
913 ddpme->dim_match = (ddpme->dim == dd->dim[dimind]);
915 ddpme->nslab = (ddpme->dim == 0 ? ddRankSetup.npmenodes_x : ddRankSetup.npmenodes_y);
917 if (ddpme->nslab <= 1)
922 const int nso = ddRankSetup.numRanksDoingPme / ddpme->nslab;
923 /* Determine for each PME slab the PP location range for dimension dim */
924 snew(ddpme->pp_min, ddpme->nslab);
925 snew(ddpme->pp_max, ddpme->nslab);
926 for (int slab = 0; slab < ddpme->nslab; slab++)
928 ddpme->pp_min[slab] = dd->numCells[dd->dim[dimind]] - 1;
929 ddpme->pp_max[slab] = 0;
931 for (int i = 0; i < dd->nnodes; i++)
934 ddindex2xyz(dd->numCells, i, xyz);
935 /* For y only use our y/z slab.
936 * This assumes that the PME x grid size matches the DD grid size.
938 if (dimind == 0 || xyz[XX] == dd->ci[XX])
940 const int pmeindex = ddindex2pmeindex(ddRankSetup, i);
941 const int slab = (dimind == 0) ? (pmeindex / nso) : (pmeindex % ddpme->nslab);
942 ddpme->pp_min[slab] = std::min(ddpme->pp_min[slab], xyz[dimind]);
943 ddpme->pp_max[slab] = std::max(ddpme->pp_max[slab], xyz[dimind]);
947 set_slb_pme_dim_f(dd, ddpme->dim, &ddpme->slb_dim_f);
950 int dd_pme_maxshift_x(const gmx_domdec_t& dd)
952 const DDRankSetup& ddRankSetup = dd.comm->ddRankSetup;
954 if (ddRankSetup.ddpme[0].dim == XX)
956 return ddRankSetup.ddpme[0].maxshift;
964 int dd_pme_maxshift_y(const gmx_domdec_t& dd)
966 const DDRankSetup& ddRankSetup = dd.comm->ddRankSetup;
968 if (ddRankSetup.ddpme[0].dim == YY)
970 return ddRankSetup.ddpme[0].maxshift;
972 else if (ddRankSetup.npmedecompdim >= 2 && ddRankSetup.ddpme[1].dim == YY)
974 return ddRankSetup.ddpme[1].maxshift;
982 bool ddUsesUpdateGroups(const gmx_domdec_t& dd)
984 return dd.comm->systemInfo.useUpdateGroups;
987 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl)
989 /* Note that the cycles value can be incorrect, either 0 or some
990 * extremely large value, when our thread migrated to another core
991 * with an unsynchronized cycle counter. If this happens less often
992 * that once per nstlist steps, this will not cause issues, since
993 * we later subtract the maximum value from the sum over nstlist steps.
994 * A zero count will slightly lower the total, but that's a small effect.
995 * Note that the main purpose of the subtraction of the maximum value
996 * is to avoid throwing off the load balancing when stalls occur due
997 * e.g. system activity or network congestion.
999 dd->comm->cycl[ddCycl] += cycles;
1000 dd->comm->cycl_n[ddCycl]++;
1001 if (cycles > dd->comm->cycl_max[ddCycl])
1003 dd->comm->cycl_max[ddCycl] = cycles;
1008 static void make_load_communicator(gmx_domdec_t* dd, int dim_ind, ivec loc)
1010 MPI_Comm c_row = MPI_COMM_NULL;
1012 bool bPartOfGroup = false;
1014 const int dim = dd->dim[dim_ind];
1015 copy_ivec(loc, loc_c);
1016 for (int i = 0; i < dd->numCells[dim]; i++)
1019 const int rank = dd_index(dd->numCells, loc_c);
1020 if (rank == dd->rank)
1022 /* This process is part of the group */
1023 bPartOfGroup = TRUE;
1026 MPI_Comm_split(dd->mpi_comm_all, bPartOfGroup ? 0 : MPI_UNDEFINED, dd->rank, &c_row);
1029 dd->comm->mpi_comm_load[dim_ind] = c_row;
1030 if (!isDlbDisabled(dd->comm))
1032 DDCellsizesWithDlb& cellsizes = dd->comm->cellsizesWithDlb[dim_ind];
1034 if (dd->ci[dim] == dd->master_ci[dim])
1036 /* This is the root process of this row */
1037 cellsizes.rowMaster = std::make_unique<RowMaster>();
1039 RowMaster& rowMaster = *cellsizes.rowMaster;
1040 rowMaster.cellFrac.resize(ddCellFractionBufferSize(dd, dim_ind));
1041 rowMaster.oldCellFrac.resize(dd->numCells[dim] + 1);
1042 rowMaster.isCellMin.resize(dd->numCells[dim]);
1045 rowMaster.bounds.resize(dd->numCells[dim]);
1047 rowMaster.buf_ncd.resize(dd->numCells[dim]);
1051 /* This is not a root process, we only need to receive cell_f */
1052 cellsizes.fracRow.resize(ddCellFractionBufferSize(dd, dim_ind));
1055 if (dd->ci[dim] == dd->master_ci[dim])
1057 snew(dd->comm->load[dim_ind].load, dd->numCells[dim] * DD_NLOAD_MAX);
1063 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id)
1066 gmx_domdec_t* dd = cr->dd;
1068 if (!thisRankHasDuty(cr, DUTY_PP) || gpu_id < 0)
1070 /* Only ranks with short-ranged tasks (currently) use GPUs.
1071 * If we don't have GPUs assigned, there are no resources to share.
1076 if (cr->nnodes == 1)
1078 dd->comm->nrank_gpu_shared = 1;
1083 const int physicalnode_id_hash = gmx_physicalnode_id_hash();
1087 fprintf(debug, "dd_setup_dd_dlb_gpu_sharing:\n");
1088 fprintf(debug, "DD PP rank %d physical node hash %d gpu_id %d\n", dd->rank, physicalnode_id_hash, gpu_id);
1090 /* Split the PP communicator over the physical nodes */
1091 /* TODO: See if we should store this (before), as it's also used for
1092 * for the nodecomm summation.
1094 // TODO PhysicalNodeCommunicator could be extended/used to handle
1095 // the need for per-node per-group communicators.
1096 MPI_Comm mpi_comm_pp_physicalnode;
1097 MPI_Comm_split(dd->mpi_comm_all, physicalnode_id_hash, dd->rank, &mpi_comm_pp_physicalnode);
1098 MPI_Comm_split(mpi_comm_pp_physicalnode, gpu_id, dd->rank, &dd->comm->mpi_comm_gpu_shared);
1099 MPI_Comm_free(&mpi_comm_pp_physicalnode);
1100 MPI_Comm_size(dd->comm->mpi_comm_gpu_shared, &dd->comm->nrank_gpu_shared);
1104 fprintf(debug, "nrank_gpu_shared %d\n", dd->comm->nrank_gpu_shared);
1107 /* Note that some ranks could share a GPU, while others don't */
1109 if (dd->comm->nrank_gpu_shared == 1)
1111 MPI_Comm_free(&dd->comm->mpi_comm_gpu_shared);
1114 GMX_UNUSED_VALUE(cr);
1115 GMX_UNUSED_VALUE(gpu_id);
1119 static void make_load_communicators(gmx_domdec_t gmx_unused* dd)
1126 fprintf(debug, "Making load communicators\n");
1129 dd->comm->load = new domdec_load_t[std::max(dd->ndim, 1)];
1130 snew(dd->comm->mpi_comm_load, std::max(dd->ndim, 1));
1138 make_load_communicator(dd, 0, loc);
1141 const int dim0 = dd->dim[0];
1142 for (int i = 0; i < dd->numCells[dim0]; i++)
1145 make_load_communicator(dd, 1, loc);
1150 const int dim0 = dd->dim[0];
1151 for (int i = 0; i < dd->numCells[dim0]; i++)
1154 const int dim1 = dd->dim[1];
1155 for (int j = 0; j < dd->numCells[dim1]; j++)
1158 make_load_communicator(dd, 2, loc);
1165 fprintf(debug, "Finished making load communicators\n");
1170 /*! \brief Sets up the relation between neighboring domains and zones */
1171 static void setup_neighbor_relations(gmx_domdec_t* dd)
1174 GMX_ASSERT((dd->ndim >= 0) && (dd->ndim <= DIM), "Must have valid number of dimensions for DD");
1176 for (int d = 0; d < dd->ndim; d++)
1178 const int dim = dd->dim[d];
1179 copy_ivec(dd->ci, tmp);
1180 tmp[dim] = (tmp[dim] + 1) % dd->numCells[dim];
1181 dd->neighbor[d][0] = ddcoord2ddnodeid(dd, tmp);
1182 copy_ivec(dd->ci, tmp);
1183 tmp[dim] = (tmp[dim] - 1 + dd->numCells[dim]) % dd->numCells[dim];
1184 dd->neighbor[d][1] = ddcoord2ddnodeid(dd, tmp);
1188 "DD rank %d neighbor ranks in dir %d are + %d - %d\n",
1192 dd->neighbor[d][1]);
1196 int nzone = (1 << dd->ndim);
1197 int nizone = (1 << std::max(dd->ndim - 1, 0));
1198 assert(nizone >= 1 && nizone <= DD_MAXIZONE);
1200 gmx_domdec_zones_t* zones = &dd->comm->zones;
1202 for (int i = 0; i < nzone; i++)
1205 clear_ivec(zones->shift[i]);
1206 for (int d = 0; d < dd->ndim; d++)
1208 zones->shift[i][dd->dim[d]] = dd_zo[i][m++];
1213 for (int i = 0; i < nzone; i++)
1215 for (int d = 0; d < DIM; d++)
1217 s[d] = dd->ci[d] - zones->shift[i][d];
1220 s[d] += dd->numCells[d];
1222 else if (s[d] >= dd->numCells[d])
1224 s[d] -= dd->numCells[d];
1228 for (int iZoneIndex = 0; iZoneIndex < nizone; iZoneIndex++)
1231 ddNonbondedZonePairRanges[iZoneIndex][0] == iZoneIndex,
1232 "The first element for each ddNonbondedZonePairRanges should match its index");
1234 DDPairInteractionRanges iZone;
1235 iZone.iZoneIndex = iZoneIndex;
1236 /* dd_zp3 is for 3D decomposition, for fewer dimensions use only
1237 * j-zones up to nzone.
1239 iZone.jZoneRange = gmx::Range<int>(std::min(ddNonbondedZonePairRanges[iZoneIndex][1], nzone),
1240 std::min(ddNonbondedZonePairRanges[iZoneIndex][2], nzone));
1241 for (int dim = 0; dim < DIM; dim++)
1243 if (dd->numCells[dim] == 1)
1245 /* All shifts should be allowed */
1246 iZone.shift0[dim] = -1;
1247 iZone.shift1[dim] = 1;
1251 /* Determine the min/max j-zone shift wrt the i-zone */
1252 iZone.shift0[dim] = 1;
1253 iZone.shift1[dim] = -1;
1254 for (int jZone : iZone.jZoneRange)
1256 int shift_diff = zones->shift[jZone][dim] - zones->shift[iZoneIndex][dim];
1257 if (shift_diff < iZone.shift0[dim])
1259 iZone.shift0[dim] = shift_diff;
1261 if (shift_diff > iZone.shift1[dim])
1263 iZone.shift1[dim] = shift_diff;
1269 zones->iZones.push_back(iZone);
1272 if (!isDlbDisabled(dd->comm))
1274 dd->comm->cellsizesWithDlb.resize(dd->ndim);
1277 if (dd->comm->ddSettings.recordLoad)
1279 make_load_communicators(dd);
1283 static void make_pp_communicator(const gmx::MDLogger& mdlog,
1285 t_commrec gmx_unused* cr,
1286 bool gmx_unused reorder)
1289 gmx_domdec_comm_t* comm = dd->comm;
1290 CartesianRankSetup& cartSetup = comm->cartesianRankSetup;
1292 if (cartSetup.bCartesianPP)
1294 /* Set up cartesian communication for the particle-particle part */
1296 .appendTextFormatted("Will use a Cartesian communicator: %d x %d x %d",
1302 for (int i = 0; i < DIM; i++)
1306 MPI_Comm comm_cart = MPI_COMM_NULL;
1307 MPI_Cart_create(cr->mpi_comm_mygroup, DIM, dd->numCells, periods, static_cast<int>(reorder), &comm_cart);
1308 /* We overwrite the old communicator with the new cartesian one */
1309 cr->mpi_comm_mygroup = comm_cart;
1312 dd->mpi_comm_all = cr->mpi_comm_mygroup;
1313 MPI_Comm_rank(dd->mpi_comm_all, &dd->rank);
1315 if (cartSetup.bCartesianPP_PME)
1317 /* Since we want to use the original cartesian setup for sim,
1318 * and not the one after split, we need to make an index.
1320 cartSetup.ddindex2ddnodeid.resize(dd->nnodes);
1321 cartSetup.ddindex2ddnodeid[dd_index(dd->numCells, dd->ci)] = dd->rank;
1322 gmx_sumi(dd->nnodes, cartSetup.ddindex2ddnodeid.data(), cr);
1323 /* Get the rank of the DD master,
1324 * above we made sure that the master node is a PP node.
1326 int rank = MASTER(cr) ? dd->rank : 0;
1327 MPI_Allreduce(&rank, &dd->masterrank, 1, MPI_INT, MPI_SUM, dd->mpi_comm_all);
1329 else if (cartSetup.bCartesianPP)
1331 if (!comm->ddRankSetup.usePmeOnlyRanks)
1333 /* The PP communicator is also
1334 * the communicator for this simulation
1336 cr->mpi_comm_mysim = cr->mpi_comm_mygroup;
1338 cr->nodeid = dd->rank;
1340 MPI_Cart_coords(dd->mpi_comm_all, dd->rank, DIM, dd->ci);
1342 /* We need to make an index to go from the coordinates
1343 * to the nodeid of this simulation.
1345 cartSetup.ddindex2simnodeid.resize(dd->nnodes);
1346 std::vector<int> buf(dd->nnodes);
1347 if (thisRankHasDuty(cr, DUTY_PP))
1349 buf[dd_index(dd->numCells, dd->ci)] = cr->sim_nodeid;
1351 /* Communicate the ddindex to simulation nodeid index */
1352 MPI_Allreduce(buf.data(), cartSetup.ddindex2simnodeid.data(), dd->nnodes, MPI_INT, MPI_SUM, cr->mpi_comm_mysim);
1354 /* Determine the master coordinates and rank.
1355 * The DD master should be the same node as the master of this sim.
1357 for (int i = 0; i < dd->nnodes; i++)
1359 if (cartSetup.ddindex2simnodeid[i] == 0)
1361 ddindex2xyz(dd->numCells, i, dd->master_ci);
1362 MPI_Cart_rank(dd->mpi_comm_all, dd->master_ci, &dd->masterrank);
1367 fprintf(debug, "The master rank is %d\n", dd->masterrank);
1372 /* No Cartesian communicators */
1373 /* We use the rank in dd->comm->all as DD index */
1374 ddindex2xyz(dd->numCells, dd->rank, dd->ci);
1375 /* The simulation master nodeid is 0, so the DD master rank is also 0 */
1377 clear_ivec(dd->master_ci);
1382 .appendTextFormatted("Domain decomposition rank %d, coordinates %d %d %d\n",
1390 "Domain decomposition rank %d, coordinates %d %d %d\n\n",
1398 static void receive_ddindex2simnodeid(gmx_domdec_t* dd, t_commrec* cr)
1401 CartesianRankSetup& cartSetup = dd->comm->cartesianRankSetup;
1403 if (!cartSetup.bCartesianPP_PME && cartSetup.bCartesianPP)
1405 cartSetup.ddindex2simnodeid.resize(dd->nnodes);
1406 std::vector<int> buf(dd->nnodes);
1407 if (thisRankHasDuty(cr, DUTY_PP))
1409 buf[dd_index(dd->numCells, dd->ci)] = cr->sim_nodeid;
1411 /* Communicate the ddindex to simulation nodeid index */
1412 MPI_Allreduce(buf.data(), cartSetup.ddindex2simnodeid.data(), dd->nnodes, MPI_INT, MPI_SUM, cr->mpi_comm_mysim);
1415 GMX_UNUSED_VALUE(dd);
1416 GMX_UNUSED_VALUE(cr);
1420 static CartesianRankSetup split_communicator(const gmx::MDLogger& mdlog,
1422 const DdRankOrder ddRankOrder,
1423 bool gmx_unused reorder,
1424 const DDRankSetup& ddRankSetup,
1426 std::vector<int>* pmeRanks)
1428 CartesianRankSetup cartSetup;
1430 cartSetup.bCartesianPP = (ddRankOrder == DdRankOrder::cartesian);
1431 cartSetup.bCartesianPP_PME = false;
1433 const ivec& numDDCells = ddRankSetup.numPPCells;
1434 /* Initially we set ntot to the number of PP cells */
1435 copy_ivec(numDDCells, cartSetup.ntot);
1437 if (cartSetup.bCartesianPP)
1439 const int numDDCellsTot = ddRankSetup.numPPRanks;
1441 for (int i = 1; i < DIM; i++)
1443 bDiv[i] = ((ddRankSetup.numRanksDoingPme * numDDCells[i]) % numDDCellsTot == 0);
1445 if (bDiv[YY] || bDiv[ZZ])
1447 cartSetup.bCartesianPP_PME = TRUE;
1448 /* If we have 2D PME decomposition, which is always in x+y,
1449 * we stack the PME only nodes in z.
1450 * Otherwise we choose the direction that provides the thinnest slab
1451 * of PME only nodes as this will have the least effect
1452 * on the PP communication.
1453 * But for the PME communication the opposite might be better.
1455 if (bDiv[ZZ] && (ddRankSetup.npmenodes_y > 1 || !bDiv[YY] || numDDCells[YY] > numDDCells[ZZ]))
1457 cartSetup.cartpmedim = ZZ;
1461 cartSetup.cartpmedim = YY;
1463 cartSetup.ntot[cartSetup.cartpmedim] +=
1464 (ddRankSetup.numRanksDoingPme * numDDCells[cartSetup.cartpmedim]) / numDDCellsTot;
1469 .appendTextFormatted(
1470 "Number of PME-only ranks (%d) is not a multiple of nx*ny (%d*%d) or "
1472 ddRankSetup.numRanksDoingPme,
1478 .appendText("Will not use a Cartesian communicator for PP <-> PME\n");
1482 if (cartSetup.bCartesianPP_PME)
1489 .appendTextFormatted(
1490 "Will use a Cartesian communicator for PP <-> PME: %d x %d x %d",
1493 cartSetup.ntot[ZZ]);
1495 for (int i = 0; i < DIM; i++)
1499 MPI_Comm comm_cart = MPI_COMM_NULL;
1500 MPI_Cart_create(cr->mpi_comm_mysim, DIM, cartSetup.ntot, periods, static_cast<int>(reorder), &comm_cart);
1501 MPI_Comm_rank(comm_cart, &rank);
1502 if (MASTER(cr) && rank != 0)
1504 gmx_fatal(FARGS, "MPI rank 0 was renumbered by MPI_Cart_create, we do not allow this");
1507 /* With this assigment we loose the link to the original communicator
1508 * which will usually be MPI_COMM_WORLD, unless have multisim.
1510 cr->mpi_comm_mysim = comm_cart;
1511 cr->sim_nodeid = rank;
1513 MPI_Cart_coords(cr->mpi_comm_mysim, cr->sim_nodeid, DIM, ddCellIndex);
1516 .appendTextFormatted("Cartesian rank %d, coordinates %d %d %d\n",
1522 if (ddCellIndex[cartSetup.cartpmedim] < numDDCells[cartSetup.cartpmedim])
1526 if (!ddRankSetup.usePmeOnlyRanks
1527 || ddCellIndex[cartSetup.cartpmedim] >= numDDCells[cartSetup.cartpmedim])
1529 cr->duty = DUTY_PME;
1532 /* Split the sim communicator into PP and PME only nodes */
1533 MPI_Comm_split(cr->mpi_comm_mysim,
1534 getThisRankDuties(cr),
1535 dd_index(cartSetup.ntot, ddCellIndex),
1536 &cr->mpi_comm_mygroup);
1537 MPI_Comm_size(cr->mpi_comm_mygroup, &cr->sizeOfMyGroupCommunicator);
1539 GMX_UNUSED_VALUE(ddCellIndex);
1544 switch (ddRankOrder)
1546 case DdRankOrder::pp_pme:
1547 GMX_LOG(mdlog.info).appendText("Order of the ranks: PP first, PME last");
1549 case DdRankOrder::interleave:
1550 /* Interleave the PP-only and PME-only ranks */
1551 GMX_LOG(mdlog.info).appendText("Interleaving PP and PME ranks");
1552 *pmeRanks = dd_interleaved_pme_ranks(ddRankSetup);
1554 case DdRankOrder::cartesian: break;
1555 default: gmx_fatal(FARGS, "Invalid ddRankOrder=%d", static_cast<int>(ddRankOrder));
1558 if (dd_simnode2pmenode(ddRankSetup, cartSetup, *pmeRanks, cr, cr->sim_nodeid) == -1)
1560 cr->duty = DUTY_PME;
1567 /* Split the sim communicator into PP and PME only nodes */
1568 MPI_Comm_split(cr->mpi_comm_mysim, getThisRankDuties(cr), cr->nodeid, &cr->mpi_comm_mygroup);
1569 MPI_Comm_size(cr->mpi_comm_mygroup, &cr->sizeOfMyGroupCommunicator);
1570 MPI_Comm_rank(cr->mpi_comm_mygroup, &cr->nodeid);
1575 .appendTextFormatted("This rank does only %s work.\n",
1576 thisRankHasDuty(cr, DUTY_PP) ? "particle-particle" : "PME-mesh");
1581 /*! \brief Makes the PP communicator and the PME communicator, when needed
1583 * Returns the Cartesian rank setup.
1584 * Sets \p cr->mpi_comm_mygroup
1585 * For PP ranks, sets the DD PP cell index in \p ddCellIndex.
1586 * With separate PME ranks in interleaved order, set the PME ranks in \p pmeRanks.
1588 static CartesianRankSetup makeGroupCommunicators(const gmx::MDLogger& mdlog,
1589 const DDSettings& ddSettings,
1590 const DdRankOrder ddRankOrder,
1591 const DDRankSetup& ddRankSetup,
1594 std::vector<int>* pmeRanks)
1596 CartesianRankSetup cartSetup;
1598 // As a default, both group and sim communicators are equal to the default communicator
1599 cr->mpi_comm_mygroup = cr->mpiDefaultCommunicator;
1600 cr->mpi_comm_mysim = cr->mpiDefaultCommunicator;
1601 cr->nnodes = cr->sizeOfDefaultCommunicator;
1602 cr->nodeid = cr->rankInDefaultCommunicator;
1603 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1605 if (ddRankSetup.usePmeOnlyRanks)
1607 /* Split the communicator into a PP and PME part */
1608 cartSetup = split_communicator(
1609 mdlog, cr, ddRankOrder, ddSettings.useCartesianReorder, ddRankSetup, ddCellIndex, pmeRanks);
1613 /* All nodes do PP and PME */
1614 /* We do not require separate communicators */
1615 cartSetup.bCartesianPP = false;
1616 cartSetup.bCartesianPP_PME = false;
1622 /*! \brief For PP ranks, sets or makes the communicator
1624 * For PME ranks get the rank id.
1625 * For PP only ranks, sets the PME-only rank.
1627 static void setupGroupCommunication(const gmx::MDLogger& mdlog,
1628 const DDSettings& ddSettings,
1629 gmx::ArrayRef<const int> pmeRanks,
1631 const int numAtomsInSystem,
1634 const DDRankSetup& ddRankSetup = dd->comm->ddRankSetup;
1635 const CartesianRankSetup& cartSetup = dd->comm->cartesianRankSetup;
1637 if (thisRankHasDuty(cr, DUTY_PP))
1641 /* Copy or make a new PP communicator */
1643 /* We (possibly) reordered the nodes in split_communicator,
1644 * so it is no longer required in make_pp_communicator.
1646 const bool useCartesianReorder =
1647 (ddSettings.useCartesianReorder && !cartSetup.bCartesianPP_PME);
1649 make_pp_communicator(mdlog, dd, cr, useCartesianReorder);
1654 receive_ddindex2simnodeid(dd, cr);
1657 if (!thisRankHasDuty(cr, DUTY_PME))
1659 /* Set up the commnuication to our PME node */
1660 dd->pme_nodeid = dd_simnode2pmenode(ddRankSetup, cartSetup, pmeRanks, cr, cr->sim_nodeid);
1661 dd->pme_receive_vir_ener = receive_vir_ener(dd, pmeRanks, cr);
1665 "My pme_nodeid %d receive ener %s\n",
1667 gmx::boolToString(dd->pme_receive_vir_ener));
1672 dd->pme_nodeid = -1;
1675 /* We can not use DDMASTER(dd), because dd->masterrank is set later */
1678 dd->ma = std::make_unique<AtomDistribution>(dd->numCells, numAtomsInSystem, numAtomsInSystem);
1682 static real* get_slb_frac(const gmx::MDLogger& mdlog, const char* dir, int nc, const char* size_string)
1684 real* slb_frac = nullptr;
1685 if (nc > 1 && size_string != nullptr)
1687 GMX_LOG(mdlog.info).appendTextFormatted("Using static load balancing for the %s direction", dir);
1690 for (int i = 0; i < nc; i++)
1694 sscanf(size_string, "%20lf%n", &dbl, &n);
1698 "Incorrect or not enough DD cell size entries for direction %s: '%s'",
1707 std::string relativeCellSizes = "Relative cell sizes:";
1708 for (int i = 0; i < nc; i++)
1711 relativeCellSizes += gmx::formatString(" %5.3f", slb_frac[i]);
1713 GMX_LOG(mdlog.info).appendText(relativeCellSizes);
1719 static int multi_body_bondeds_count(const gmx_mtop_t& mtop)
1722 for (const auto ilists : IListRange(mtop))
1724 for (auto& ilist : extractILists(ilists.list(), IF_BOND))
1726 if (NRAL(ilist.functionType) > 2)
1728 n += ilists.nmol() * (ilist.iatoms.size() / ilistStride(ilist));
1736 static int dd_getenv(const gmx::MDLogger& mdlog, const char* env_var, int def)
1739 char* val = getenv(env_var);
1742 if (sscanf(val, "%20d", &nst) <= 0)
1746 GMX_LOG(mdlog.info).appendTextFormatted("Found env.var. %s = %s, using value %d", env_var, val, nst);
1752 static void check_dd_restrictions(const gmx_domdec_t* dd, const t_inputrec& inputrec, const gmx::MDLogger& mdlog)
1754 if (inputrec.pbcType == PbcType::Screw
1755 && (dd->numCells[XX] == 1 || dd->numCells[YY] > 1 || dd->numCells[ZZ] > 1))
1758 "With pbc=%s can only do domain decomposition in the x-direction",
1759 c_pbcTypeNames[inputrec.pbcType].c_str());
1762 if (inputrec.nstlist == 0)
1764 gmx_fatal(FARGS, "Domain decomposition does not work with nstlist=0");
1767 if (inputrec.comm_mode == ComRemovalAlgorithm::Angular && inputrec.pbcType != PbcType::No)
1769 GMX_LOG(mdlog.warning)
1771 "comm-mode angular will give incorrect results when the comm group "
1772 "partially crosses a periodic boundary");
1776 static real average_cellsize_min(const gmx_ddbox_t& ddbox, const ivec numDomains)
1778 real r = ddbox.box_size[XX];
1779 for (int d = 0; d < DIM; d++)
1781 if (numDomains[d] > 1)
1783 /* Check using the initial average cell size */
1784 r = std::min(r, ddbox.box_size[d] * ddbox.skew_fac[d] / numDomains[d]);
1791 /*! \brief Depending on the DLB initial value return the DLB switched off state or issue an error.
1793 static DlbState forceDlbOffOrBail(DlbState cmdlineDlbState,
1794 const std::string& reasonStr,
1795 const gmx::MDLogger& mdlog)
1797 std::string dlbNotSupportedErr = "Dynamic load balancing requested, but ";
1798 std::string dlbDisableNote = "NOTE: disabling dynamic load balancing as ";
1800 if (cmdlineDlbState == DlbState::onUser)
1802 gmx_fatal(FARGS, "%s", (dlbNotSupportedErr + reasonStr).c_str());
1804 else if (cmdlineDlbState == DlbState::offCanTurnOn)
1806 GMX_LOG(mdlog.info).appendText(dlbDisableNote + reasonStr);
1808 return DlbState::offForever;
1811 /*! \brief Return the dynamic load balancer's initial state based on initial conditions and user inputs.
1813 * This function parses the parameters of "-dlb" command line option setting
1814 * corresponding state values. Then it checks the consistency of the determined
1815 * state with other run parameters and settings. As a result, the initial state
1816 * may be altered or an error may be thrown if incompatibility of options is detected.
1818 * \param [in] mdlog Logger.
1819 * \param [in] dlbOption Enum value for the DLB option.
1820 * \param [in] bRecordLoad True if the load balancer is recording load information.
1821 * \param [in] mdrunOptions Options for mdrun.
1822 * \param [in] inputrec Pointer mdrun to input parameters.
1823 * \param [in] directGpuCommUsedWithGpuUpdate Direct GPU halo exchange and GPU update enabled
1824 * \returns DLB initial/startup state.
1826 static DlbState determineInitialDlbState(const gmx::MDLogger& mdlog,
1827 DlbOption dlbOption,
1828 gmx_bool bRecordLoad,
1829 const gmx::MdrunOptions& mdrunOptions,
1830 const t_inputrec& inputrec,
1831 const bool directGpuCommUsedWithGpuUpdate)
1833 DlbState dlbState = DlbState::offCanTurnOn;
1837 case DlbOption::turnOnWhenUseful: dlbState = DlbState::offCanTurnOn; break;
1838 case DlbOption::no: dlbState = DlbState::offUser; break;
1839 case DlbOption::yes: dlbState = DlbState::onUser; break;
1840 default: gmx_incons("Invalid dlbOption enum value");
1843 // P2P GPU comm + GPU update leads to case in which we enqueue async work for multiple timesteps
1844 // DLB needs to be disabled in that case
1845 if (directGpuCommUsedWithGpuUpdate)
1847 std::string reasonStr =
1848 "it is not supported with GPU direct communication + GPU update enabled.";
1849 return forceDlbOffOrBail(dlbState, reasonStr, mdlog);
1852 /* Reruns don't support DLB: bail or override auto mode */
1853 if (mdrunOptions.rerun)
1855 std::string reasonStr = "it is not supported in reruns.";
1856 return forceDlbOffOrBail(dlbState, reasonStr, mdlog);
1859 /* Unsupported integrators */
1860 if (!EI_DYNAMICS(inputrec.eI))
1863 gmx::formatString("it is only supported with dynamics, not with integrator '%s'.",
1864 enumValueToString(inputrec.eI));
1865 return forceDlbOffOrBail(dlbState, reasonStr, mdlog);
1868 /* Without cycle counters we can't time work to balance on */
1871 std::string reasonStr =
1872 "cycle counters unsupported or not enabled in the operating system kernel.";
1873 return forceDlbOffOrBail(dlbState, reasonStr, mdlog);
1876 if (mdrunOptions.reproducible)
1878 std::string reasonStr = "you started a reproducible run.";
1881 case DlbState::offUser: break;
1882 case DlbState::offForever:
1883 GMX_RELEASE_ASSERT(false, "DlbState::offForever is not a valid initial state");
1885 case DlbState::offCanTurnOn: return forceDlbOffOrBail(dlbState, reasonStr, mdlog);
1886 case DlbState::onCanTurnOff:
1887 GMX_RELEASE_ASSERT(false, "DlbState::offCanTurnOff is not a valid initial state");
1889 case DlbState::onUser:
1890 return forceDlbOffOrBail(
1893 + " In load balanced runs binary reproducibility cannot be "
1898 "Death horror: undefined case (%d) for load balancing choice",
1899 static_cast<int>(dlbState));
1906 static gmx_domdec_comm_t* init_dd_comm()
1908 gmx_domdec_comm_t* comm = new gmx_domdec_comm_t;
1910 comm->n_load_have = 0;
1911 comm->n_load_collect = 0;
1913 comm->haveTurnedOffDlb = false;
1915 for (int i = 0; i < static_cast<int>(DDAtomRanges::Type::Number); i++)
1917 comm->sum_nat[i] = 0;
1921 comm->load_step = 0;
1924 clear_ivec(comm->load_lim);
1928 /* This should be replaced by a unique pointer */
1929 comm->balanceRegion = ddBalanceRegionAllocate();
1934 static void setupUpdateGroups(const gmx::MDLogger& mdlog,
1935 const gmx_mtop_t& mtop,
1936 ArrayRef<const RangePartitioning> updateGroupingsPerMoleculeType,
1937 const bool useUpdateGroups,
1938 const real maxUpdateGroupRadius,
1939 DDSystemInfo* systemInfo)
1941 systemInfo->updateGroupingsPerMoleculeType = updateGroupingsPerMoleculeType;
1942 systemInfo->useUpdateGroups = useUpdateGroups;
1943 systemInfo->maxUpdateGroupRadius = maxUpdateGroupRadius;
1945 if (systemInfo->useUpdateGroups)
1947 int numUpdateGroups = 0;
1948 for (const auto& molblock : mtop.molblock)
1950 numUpdateGroups += molblock.nmol
1951 * systemInfo->updateGroupingsPerMoleculeType[molblock.type].numBlocks();
1955 .appendTextFormatted(
1956 "Using update groups, nr %d, average size %.1f atoms, max. radius %.3f "
1959 mtop.natoms / static_cast<double>(numUpdateGroups),
1960 systemInfo->maxUpdateGroupRadius);
1964 UnitCellInfo::UnitCellInfo(const t_inputrec& ir) :
1965 npbcdim(numPbcDimensions(ir.pbcType)),
1966 numBoundedDimensions(inputrec2nboundeddim(&ir)),
1967 ddBoxIsDynamic(numBoundedDimensions < DIM || inputrecDynamicBox(&ir)),
1968 haveScrewPBC(ir.pbcType == PbcType::Screw)
1972 /* Returns whether molecules are always whole, i.e. not broken by PBC */
1973 static bool moleculesAreAlwaysWhole(const gmx_mtop_t& mtop,
1974 const bool useUpdateGroups,
1975 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingsPerMoleculeType)
1977 if (useUpdateGroups)
1979 GMX_RELEASE_ASSERT(updateGroupingsPerMoleculeType.size() == mtop.moltype.size(),
1980 "Need one grouping per moltype");
1981 for (size_t mol = 0; mol < mtop.moltype.size(); mol++)
1983 if (updateGroupingsPerMoleculeType[mol].numBlocks() > 1)
1991 for (const auto& moltype : mtop.moltype)
1993 if (moltype.atoms.nr > 1)
2003 /*! \brief Generate the simulation system information */
2004 static DDSystemInfo getSystemInfo(const gmx::MDLogger& mdlog,
2006 MPI_Comm communicator,
2007 const DomdecOptions& options,
2008 const gmx_mtop_t& mtop,
2009 const t_inputrec& ir,
2011 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType,
2012 const bool useUpdateGroups,
2013 const real maxUpdateGroupRadius,
2014 gmx::ArrayRef<const gmx::RVec> xGlobal)
2016 const real tenPercentMargin = 1.1;
2018 DDSystemInfo systemInfo;
2021 mdlog, mtop, updateGroupingPerMoleculeType, useUpdateGroups, maxUpdateGroupRadius, &systemInfo);
2023 systemInfo.moleculesAreAlwaysWhole = moleculesAreAlwaysWhole(
2024 mtop, systemInfo.useUpdateGroups, systemInfo.updateGroupingsPerMoleculeType);
2025 systemInfo.haveInterDomainBondeds =
2026 (!systemInfo.moleculesAreAlwaysWhole || mtop.bIntermolecularInteractions);
2027 systemInfo.haveInterDomainMultiBodyBondeds =
2028 (systemInfo.haveInterDomainBondeds && multi_body_bondeds_count(mtop) > 0);
2030 if (systemInfo.useUpdateGroups)
2032 systemInfo.mayHaveSplitConstraints = false;
2033 systemInfo.mayHaveSplitSettles = false;
2037 systemInfo.mayHaveSplitConstraints = (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0
2038 || gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0);
2039 systemInfo.mayHaveSplitSettles = (gmx_mtop_ftype_count(mtop, F_SETTLE) > 0);
2044 /* Set the cut-off to some very large value,
2045 * so we don't need if statements everywhere in the code.
2046 * We use sqrt, since the cut-off is squared in some places.
2048 systemInfo.cutoff = GMX_CUTOFF_INF;
2052 systemInfo.cutoff = atomToAtomIntoDomainToDomainCutoff(systemInfo, ir.rlist);
2054 systemInfo.minCutoffForMultiBody = 0;
2056 /* Determine the minimum cell size limit, affected by many factors */
2057 systemInfo.cellsizeLimit = 0;
2058 systemInfo.filterBondedCommunication = false;
2060 /* We do not allow home atoms to move beyond the neighboring domain
2061 * between domain decomposition steps, which limits the cell size.
2062 * Get an estimate of cell size limit due to atom displacement.
2063 * In most cases this is a large overestimate, because it assumes
2064 * non-interaction atoms.
2065 * We set the chance to 1 in a trillion steps.
2067 constexpr real c_chanceThatAtomMovesBeyondDomain = 1e-12;
2068 const real limitForAtomDisplacement = minCellSizeForAtomDisplacement(
2069 mtop, ir, systemInfo.updateGroupingsPerMoleculeType, c_chanceThatAtomMovesBeyondDomain);
2070 GMX_LOG(mdlog.info).appendTextFormatted("Minimum cell size due to atom displacement: %.3f nm", limitForAtomDisplacement);
2072 systemInfo.cellsizeLimit = std::max(systemInfo.cellsizeLimit, limitForAtomDisplacement);
2074 /* TODO: PME decomposition currently requires atoms not to be more than
2075 * 2/3 of comm->cutoff, which is >=rlist, outside of their domain.
2076 * In nearly all cases, limitForAtomDisplacement will be smaller
2077 * than 2/3*rlist, so the PME requirement is satisfied.
2078 * But it would be better for both correctness and performance
2079 * to use limitForAtomDisplacement instead of 2/3*comm->cutoff.
2080 * Note that we would need to improve the pairlist buffer case.
2083 if (systemInfo.haveInterDomainBondeds)
2085 if (options.minimumCommunicationRange > 0)
2087 systemInfo.minCutoffForMultiBody =
2088 atomToAtomIntoDomainToDomainCutoff(systemInfo, options.minimumCommunicationRange);
2089 if (options.useBondedCommunication)
2091 systemInfo.filterBondedCommunication =
2092 (systemInfo.minCutoffForMultiBody > systemInfo.cutoff);
2096 systemInfo.cutoff = std::max(systemInfo.cutoff, systemInfo.minCutoffForMultiBody);
2099 else if (ir.bPeriodicMols)
2101 /* Can not easily determine the required cut-off */
2102 GMX_LOG(mdlog.warning)
2104 "NOTE: Periodic molecules are present in this system. Because of this, "
2105 "the domain decomposition algorithm cannot easily determine the "
2106 "minimum cell size that it requires for treating bonded interactions. "
2107 "Instead, domain decomposition will assume that half the non-bonded "
2108 "cut-off will be a suitable lower bound.");
2109 systemInfo.minCutoffForMultiBody = systemInfo.cutoff / 2;
2116 if (ddRole == DDRole::Master)
2118 dd_bonded_cg_distance(mdlog, mtop, ir, xGlobal, box, options.ddBondedChecking, &r_2b, &r_mb);
2120 gmx_bcast(sizeof(r_2b), &r_2b, communicator);
2121 gmx_bcast(sizeof(r_mb), &r_mb, communicator);
2123 /* We use an initial margin of 10% for the minimum cell size,
2124 * except when we are just below the non-bonded cut-off.
2126 if (options.useBondedCommunication)
2128 if (std::max(r_2b, r_mb) > systemInfo.cutoff)
2130 const real r_bonded = std::max(r_2b, r_mb);
2131 systemInfo.minCutoffForMultiBody = tenPercentMargin * r_bonded;
2132 /* This is the (only) place where we turn on the filtering */
2133 systemInfo.filterBondedCommunication = true;
2137 const real r_bonded = r_mb;
2138 systemInfo.minCutoffForMultiBody =
2139 std::min(tenPercentMargin * r_bonded, systemInfo.cutoff);
2141 /* We determine cutoff_mbody later */
2142 systemInfo.increaseMultiBodyCutoff = true;
2146 /* No special bonded communication,
2147 * simply increase the DD cut-off.
2149 systemInfo.minCutoffForMultiBody = tenPercentMargin * std::max(r_2b, r_mb);
2150 systemInfo.cutoff = std::max(systemInfo.cutoff, systemInfo.minCutoffForMultiBody);
2154 .appendTextFormatted("Minimum cell size due to bonded interactions: %.3f nm",
2155 systemInfo.minCutoffForMultiBody);
2157 systemInfo.cellsizeLimit = std::max(systemInfo.cellsizeLimit, systemInfo.minCutoffForMultiBody);
2160 systemInfo.constraintCommunicationRange = 0;
2161 if (systemInfo.mayHaveSplitConstraints && options.constraintCommunicationRange <= 0)
2163 /* There is a cell size limit due to the constraints (P-LINCS) */
2164 systemInfo.constraintCommunicationRange = gmx::constr_r_max(mdlog, &mtop, &ir);
2166 .appendTextFormatted("Estimated maximum distance required for P-LINCS: %.3f nm",
2167 systemInfo.constraintCommunicationRange);
2168 if (systemInfo.constraintCommunicationRange > systemInfo.cellsizeLimit)
2172 "This distance will limit the DD cell size, you can override this with "
2176 else if (options.constraintCommunicationRange > 0)
2178 /* Here we do not check for dd->splitConstraints.
2179 * because one can also set a cell size limit for virtual sites only
2180 * and at this point we don't know yet if there are intercg v-sites.
2183 .appendTextFormatted("User supplied maximum distance required for P-LINCS: %.3f nm",
2184 options.constraintCommunicationRange);
2185 systemInfo.constraintCommunicationRange = options.constraintCommunicationRange;
2187 systemInfo.cellsizeLimit = std::max(systemInfo.cellsizeLimit, systemInfo.constraintCommunicationRange);
2192 /*! \brief Exit with a fatal error if the DDGridSetup cannot be
2194 static void checkDDGridSetup(const DDGridSetup& ddGridSetup,
2196 MPI_Comm communicator,
2198 const DomdecOptions& options,
2199 const DDSettings& ddSettings,
2200 const DDSystemInfo& systemInfo,
2201 const real cellsizeLimit,
2202 const gmx_ddbox_t& ddbox)
2204 if (options.numCells[XX] <= 0 && (ddGridSetup.numDomains[XX] == 0))
2206 const bool bC = (systemInfo.mayHaveSplitConstraints
2207 && systemInfo.constraintCommunicationRange > systemInfo.minCutoffForMultiBody);
2208 std::string message =
2209 gmx::formatString("Change the number of ranks or mdrun option %s%s%s",
2210 !bC ? "-rdd" : "-rcon",
2211 ddSettings.initialDlbState != DlbState::offUser ? " or -dds" : "",
2212 bC ? " or your LINCS settings" : "");
2214 gmx_fatal_collective(FARGS,
2216 ddRole == DDRole::Master,
2217 "There is no domain decomposition for %d ranks that is compatible "
2218 "with the given box and a minimum cell size of %g nm\n"
2220 "Look in the log file for details on the domain decomposition",
2221 numNodes - ddGridSetup.numPmeOnlyRanks,
2226 const real acs = average_cellsize_min(ddbox, ddGridSetup.numDomains);
2227 if (acs < cellsizeLimit)
2229 if (options.numCells[XX] <= 0)
2233 "dd_choose_grid() should return a grid that satisfies the cell size limits");
2237 gmx_fatal_collective(
2240 ddRole == DDRole::Master,
2241 "The initial cell size (%f) is smaller than the cell size limit (%f), change "
2242 "options -dd, -rdd or -rcon, see the log file for details",
2248 const int numPPRanks =
2249 ddGridSetup.numDomains[XX] * ddGridSetup.numDomains[YY] * ddGridSetup.numDomains[ZZ];
2250 if (numNodes - numPPRanks != ddGridSetup.numPmeOnlyRanks)
2252 gmx_fatal_collective(FARGS,
2254 ddRole == DDRole::Master,
2255 "The size of the domain decomposition grid (%d) does not match the "
2256 "number of PP ranks (%d). The total number of ranks is %d",
2258 numNodes - ddGridSetup.numPmeOnlyRanks,
2261 if (ddGridSetup.numPmeOnlyRanks > numPPRanks)
2263 gmx_fatal_collective(FARGS,
2265 ddRole == DDRole::Master,
2266 "The number of separate PME ranks (%d) is larger than the number of "
2267 "PP ranks (%d), this is not supported.",
2268 ddGridSetup.numPmeOnlyRanks,
2273 /*! \brief Set the cell size and interaction limits, as well as the DD grid */
2274 static DDRankSetup getDDRankSetup(const gmx::MDLogger& mdlog,
2276 const DdRankOrder rankOrder,
2277 const DDGridSetup& ddGridSetup,
2278 const t_inputrec& ir)
2281 .appendTextFormatted("Domain decomposition grid %d x %d x %d, separate PME ranks %d",
2282 ddGridSetup.numDomains[XX],
2283 ddGridSetup.numDomains[YY],
2284 ddGridSetup.numDomains[ZZ],
2285 ddGridSetup.numPmeOnlyRanks);
2287 DDRankSetup ddRankSetup;
2289 ddRankSetup.rankOrder = rankOrder;
2291 ddRankSetup.numPPRanks = numNodes - ddGridSetup.numPmeOnlyRanks;
2292 copy_ivec(ddGridSetup.numDomains, ddRankSetup.numPPCells);
2294 ddRankSetup.usePmeOnlyRanks = (ddGridSetup.numPmeOnlyRanks > 0);
2295 if (ddRankSetup.usePmeOnlyRanks)
2297 ddRankSetup.numRanksDoingPme = ddGridSetup.numPmeOnlyRanks;
2301 ddRankSetup.numRanksDoingPme =
2302 ddGridSetup.numDomains[XX] * ddGridSetup.numDomains[YY] * ddGridSetup.numDomains[ZZ];
2305 if (EEL_PME(ir.coulombtype) || EVDW_PME(ir.vdwtype))
2307 /* The following choices should match those
2308 * in comm_cost_est in domdec_setup.c.
2309 * Note that here the checks have to take into account
2310 * that the decomposition might occur in a different order than xyz
2311 * (for instance through the env.var. GMX_DD_ORDER_ZYX),
2312 * in which case they will not match those in comm_cost_est,
2313 * but since that is mainly for testing purposes that's fine.
2315 if (ddGridSetup.numDDDimensions >= 2 && ddGridSetup.ddDimensions[0] == XX
2316 && ddGridSetup.ddDimensions[1] == YY
2317 && ddRankSetup.numRanksDoingPme > ddGridSetup.numDomains[XX]
2318 && ddRankSetup.numRanksDoingPme % ddGridSetup.numDomains[XX] == 0
2319 && getenv("GMX_PMEONEDD") == nullptr)
2321 ddRankSetup.npmedecompdim = 2;
2322 ddRankSetup.npmenodes_x = ddGridSetup.numDomains[XX];
2323 ddRankSetup.npmenodes_y = ddRankSetup.numRanksDoingPme / ddRankSetup.npmenodes_x;
2327 /* In case nc is 1 in both x and y we could still choose to
2328 * decompose pme in y instead of x, but we use x for simplicity.
2330 ddRankSetup.npmedecompdim = 1;
2331 if (ddGridSetup.ddDimensions[0] == YY)
2333 ddRankSetup.npmenodes_x = 1;
2334 ddRankSetup.npmenodes_y = ddRankSetup.numRanksDoingPme;
2338 ddRankSetup.npmenodes_x = ddRankSetup.numRanksDoingPme;
2339 ddRankSetup.npmenodes_y = 1;
2343 .appendTextFormatted("PME domain decomposition: %d x %d x %d",
2344 ddRankSetup.npmenodes_x,
2345 ddRankSetup.npmenodes_y,
2350 ddRankSetup.npmedecompdim = 0;
2351 ddRankSetup.npmenodes_x = 0;
2352 ddRankSetup.npmenodes_y = 0;
2358 /*! \brief Set the cell size and interaction limits */
2359 static void set_dd_limits(const gmx::MDLogger& mdlog,
2362 const DomdecOptions& options,
2363 const DDSettings& ddSettings,
2364 const DDSystemInfo& systemInfo,
2365 const DDGridSetup& ddGridSetup,
2366 const int numPPRanks,
2367 const gmx_mtop_t& mtop,
2368 const t_inputrec& ir,
2369 const gmx_ddbox_t& ddbox)
2371 gmx_domdec_comm_t* comm = dd->comm;
2372 comm->ddSettings = ddSettings;
2374 /* Initialize to GPU share count to 0, might change later */
2375 comm->nrank_gpu_shared = 0;
2377 comm->dlbState = comm->ddSettings.initialDlbState;
2378 dd_dlb_set_should_check_whether_to_turn_dlb_on(dd, TRUE);
2379 /* To consider turning DLB on after 2*nstlist steps we need to check
2380 * at partitioning count 3. Thus we need to increase the first count by 2.
2382 comm->ddPartioningCountFirstDlbOff += 2;
2384 comm->bPMELoadBalDLBLimits = FALSE;
2386 /* Allocate the charge group/atom sorting struct */
2387 comm->sort = std::make_unique<gmx_domdec_sort_t>();
2389 comm->systemInfo = systemInfo;
2391 if (systemInfo.useUpdateGroups)
2393 /* Note: We would like to use dd->nnodes for the atom count estimate,
2394 * but that is not yet available here. But this anyhow only
2395 * affect performance up to the second dd_partition_system call.
2397 const int homeAtomCountEstimate = mtop.natoms / numPPRanks;
2398 comm->updateGroupsCog = std::make_unique<gmx::UpdateGroupsCog>(
2399 mtop, systemInfo.updateGroupingsPerMoleculeType, maxReferenceTemperature(ir), homeAtomCountEstimate);
2402 /* Set the DD setup given by ddGridSetup */
2403 copy_ivec(ddGridSetup.numDomains, dd->numCells);
2404 dd->ndim = ddGridSetup.numDDDimensions;
2405 copy_ivec(ddGridSetup.ddDimensions, dd->dim);
2407 dd->nnodes = dd->numCells[XX] * dd->numCells[YY] * dd->numCells[ZZ];
2409 snew(comm->slb_frac, DIM);
2410 if (isDlbDisabled(comm))
2412 comm->slb_frac[XX] = get_slb_frac(mdlog, "x", dd->numCells[XX], options.cellSizeX);
2413 comm->slb_frac[YY] = get_slb_frac(mdlog, "y", dd->numCells[YY], options.cellSizeY);
2414 comm->slb_frac[ZZ] = get_slb_frac(mdlog, "z", dd->numCells[ZZ], options.cellSizeZ);
2417 /* Set the multi-body cut-off and cellsize limit for DLB */
2418 comm->cutoff_mbody = systemInfo.minCutoffForMultiBody;
2419 comm->cellsize_limit = systemInfo.cellsizeLimit;
2420 if (systemInfo.haveInterDomainBondeds && systemInfo.increaseMultiBodyCutoff)
2422 if (systemInfo.filterBondedCommunication || !isDlbDisabled(comm))
2424 /* Set the bonded communication distance to halfway
2425 * the minimum and the maximum,
2426 * since the extra communication cost is nearly zero.
2428 real acs = average_cellsize_min(ddbox, dd->numCells);
2429 comm->cutoff_mbody = 0.5 * (systemInfo.minCutoffForMultiBody + acs);
2430 if (!isDlbDisabled(comm))
2432 /* Check if this does not limit the scaling */
2433 comm->cutoff_mbody = std::min(comm->cutoff_mbody, options.dlbScaling * acs);
2435 if (!systemInfo.filterBondedCommunication)
2437 /* Without bBondComm do not go beyond the n.b. cut-off */
2438 comm->cutoff_mbody = std::min(comm->cutoff_mbody, systemInfo.cutoff);
2439 if (comm->cellsize_limit >= systemInfo.cutoff)
2441 /* We don't loose a lot of efficieny
2442 * when increasing it to the n.b. cut-off.
2443 * It can even be slightly faster, because we need
2444 * less checks for the communication setup.
2446 comm->cutoff_mbody = systemInfo.cutoff;
2449 /* Check if we did not end up below our original limit */
2450 comm->cutoff_mbody = std::max(comm->cutoff_mbody, systemInfo.minCutoffForMultiBody);
2452 if (comm->cutoff_mbody > comm->cellsize_limit)
2454 comm->cellsize_limit = comm->cutoff_mbody;
2457 /* Without DLB and cutoff_mbody<cutoff, cutoff_mbody is dynamic */
2463 "Bonded atom communication beyond the cut-off: %s\n"
2464 "cellsize limit %f\n",
2465 gmx::boolToString(systemInfo.filterBondedCommunication),
2466 comm->cellsize_limit);
2469 if (ddRole == DDRole::Master)
2471 check_dd_restrictions(dd, ir, mdlog);
2475 static void writeSettings(gmx::TextWriter* log,
2477 const gmx_mtop_t& mtop,
2478 const t_inputrec& ir,
2479 gmx_bool bDynLoadBal,
2481 const gmx_ddbox_t* ddbox)
2483 gmx_domdec_comm_t* comm = dd->comm;
2487 log->writeString("The maximum number of communication pulses is:");
2488 for (int d = 0; d < dd->ndim; d++)
2490 log->writeStringFormatted(" %c %d", dim2char(dd->dim[d]), comm->cd[d].np_dlb);
2492 log->ensureLineBreak();
2493 log->writeLineFormatted("The minimum size for domain decomposition cells is %.3f nm",
2494 comm->cellsize_limit);
2495 log->writeLineFormatted("The requested allowed shrink of DD cells (option -dds) is: %.2f", dlb_scale);
2496 log->writeString("The allowed shrink of domain decomposition cells is:");
2497 for (int d = 0; d < DIM; d++)
2499 if (dd->numCells[d] > 1)
2502 (d >= ddbox->npbcdim && dd->numCells[d] == 2)
2504 : comm->cellsize_min_dlb[d]
2505 / (ddbox->box_size[d] * ddbox->skew_fac[d] / dd->numCells[d]);
2506 log->writeStringFormatted(" %c %.2f", dim2char(d), shrink);
2509 log->ensureLineBreak();
2514 set_dd_cell_sizes_slb(dd, ddbox, setcellsizeslbPULSE_ONLY, np);
2515 log->writeString("The initial number of communication pulses is:");
2516 for (int d = 0; d < dd->ndim; d++)
2518 log->writeStringFormatted(" %c %d", dim2char(dd->dim[d]), np[dd->dim[d]]);
2520 log->ensureLineBreak();
2521 log->writeString("The initial domain decomposition cell size is:");
2522 for (int d = 0; d < DIM; d++)
2524 if (dd->numCells[d] > 1)
2526 log->writeStringFormatted(" %c %.2f nm", dim2char(d), dd->comm->cellsize_min[d]);
2529 log->ensureLineBreak();
2533 const bool haveInterDomainVsites =
2534 (countInterUpdategroupVsites(mtop, comm->systemInfo.updateGroupingsPerMoleculeType) != 0);
2536 if (comm->systemInfo.haveInterDomainBondeds || haveInterDomainVsites
2537 || comm->systemInfo.mayHaveSplitConstraints || comm->systemInfo.mayHaveSplitSettles)
2539 std::string decompUnits;
2540 if (comm->systemInfo.useUpdateGroups)
2542 decompUnits = "atom groups";
2546 decompUnits = "atoms";
2549 log->writeLineFormatted("The maximum allowed distance for %s involved in interactions is:",
2550 decompUnits.c_str());
2551 log->writeLineFormatted(
2552 "%40s %-7s %6.3f nm", "non-bonded interactions", "", comm->systemInfo.cutoff);
2557 limit = dd->comm->cellsize_limit;
2561 if (dd->unitCellInfo.ddBoxIsDynamic)
2564 "(the following are initial values, they could change due to box "
2567 limit = dd->comm->cellsize_min[XX];
2568 for (int d = 1; d < DIM; d++)
2570 limit = std::min(limit, dd->comm->cellsize_min[d]);
2574 if (comm->systemInfo.haveInterDomainBondeds)
2576 log->writeLineFormatted("%40s %-7s %6.3f nm",
2577 "two-body bonded interactions",
2579 std::max(comm->systemInfo.cutoff, comm->cutoff_mbody));
2580 log->writeLineFormatted("%40s %-7s %6.3f nm",
2581 "multi-body bonded interactions",
2583 (comm->systemInfo.filterBondedCommunication || isDlbOn(dd->comm))
2584 ? comm->cutoff_mbody
2585 : std::min(comm->systemInfo.cutoff, limit));
2587 if (haveInterDomainVsites)
2589 log->writeLineFormatted("%40s %-7s %6.3f nm", "virtual site constructions", "(-rcon)", limit);
2591 if (comm->systemInfo.mayHaveSplitConstraints || comm->systemInfo.mayHaveSplitSettles)
2593 std::string separation =
2594 gmx::formatString("atoms separated by up to %d constraints", 1 + ir.nProjOrder);
2595 log->writeLineFormatted("%40s %-7s %6.3f nm\n", separation.c_str(), "(-rcon)", limit);
2597 log->ensureLineBreak();
2601 static void logSettings(const gmx::MDLogger& mdlog,
2603 const gmx_mtop_t& mtop,
2604 const t_inputrec& ir,
2606 const gmx_ddbox_t* ddbox)
2608 gmx::StringOutputStream stream;
2609 gmx::TextWriter log(&stream);
2610 writeSettings(&log, dd, mtop, ir, isDlbOn(dd->comm), dlb_scale, ddbox);
2611 if (dd->comm->dlbState == DlbState::offCanTurnOn)
2614 log.ensureEmptyLine();
2616 "When dynamic load balancing gets turned on, these settings will change to:");
2618 writeSettings(&log, dd, mtop, ir, true, dlb_scale, ddbox);
2620 GMX_LOG(mdlog.info).asParagraph().appendText(stream.toString());
2623 static void set_cell_limits_dlb(const gmx::MDLogger& mdlog,
2626 const t_inputrec& inputrec,
2627 const gmx_ddbox_t* ddbox)
2630 int npulse_d_max = 0;
2633 gmx_domdec_comm_t* comm = dd->comm;
2635 bool bNoCutOff = (inputrec.rvdw == 0 || inputrec.rcoulomb == 0);
2637 /* Determine the maximum number of comm. pulses in one dimension */
2639 comm->cellsize_limit = std::max(comm->cellsize_limit, comm->cutoff_mbody);
2641 /* Determine the maximum required number of grid pulses */
2642 if (comm->cellsize_limit >= comm->systemInfo.cutoff)
2644 /* Only a single pulse is required */
2647 else if (!bNoCutOff && comm->cellsize_limit > 0)
2649 /* We round down slightly here to avoid overhead due to the latency
2650 * of extra communication calls when the cut-off
2651 * would be only slightly longer than the cell size.
2652 * Later cellsize_limit is redetermined,
2653 * so we can not miss interactions due to this rounding.
2655 npulse = static_cast<int>(0.96 + comm->systemInfo.cutoff / comm->cellsize_limit);
2659 /* There is no cell size limit */
2660 npulse = std::max(dd->numCells[XX] - 1, std::max(dd->numCells[YY] - 1, dd->numCells[ZZ] - 1));
2663 if (!bNoCutOff && npulse > 1)
2665 /* See if we can do with less pulses, based on dlb_scale */
2667 for (int d = 0; d < dd->ndim; d++)
2669 int dim = dd->dim[d];
2670 npulse_d = static_cast<int>(
2672 + dd->numCells[dim] * comm->systemInfo.cutoff
2673 / (ddbox->box_size[dim] * ddbox->skew_fac[dim] * dlb_scale));
2674 npulse_d_max = std::max(npulse_d_max, npulse_d);
2676 npulse = std::min(npulse, npulse_d_max);
2679 /* This env var can override npulse */
2680 const int ddPulseEnv = dd_getenv(mdlog, "GMX_DD_NPULSE", 0);
2683 npulse = ddPulseEnv;
2687 comm->bVacDLBNoLimit = (inputrec.pbcType == PbcType::No);
2688 for (int d = 0; d < dd->ndim; d++)
2690 comm->cd[d].np_dlb = std::min(npulse, dd->numCells[dd->dim[d]] - 1);
2691 comm->maxpulse = std::max(comm->maxpulse, comm->cd[d].np_dlb);
2692 if (comm->cd[d].np_dlb < dd->numCells[dd->dim[d]] - 1)
2694 comm->bVacDLBNoLimit = FALSE;
2698 /* cellsize_limit is set for LINCS in init_domain_decomposition */
2699 if (!comm->bVacDLBNoLimit)
2701 comm->cellsize_limit = std::max(comm->cellsize_limit, comm->systemInfo.cutoff / comm->maxpulse);
2703 comm->cellsize_limit = std::max(comm->cellsize_limit, comm->cutoff_mbody);
2704 /* Set the minimum cell size for each DD dimension */
2705 for (int d = 0; d < dd->ndim; d++)
2707 if (comm->bVacDLBNoLimit || comm->cd[d].np_dlb * comm->cellsize_limit >= comm->systemInfo.cutoff)
2709 comm->cellsize_min_dlb[dd->dim[d]] = comm->cellsize_limit;
2713 comm->cellsize_min_dlb[dd->dim[d]] = comm->systemInfo.cutoff / comm->cd[d].np_dlb;
2716 if (comm->cutoff_mbody <= 0)
2718 comm->cutoff_mbody = std::min(comm->systemInfo.cutoff, comm->cellsize_limit);
2726 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd)
2728 return dd.comm->systemInfo.moleculesAreAlwaysWhole;
2731 bool dd_bonded_molpbc(const gmx_domdec_t& dd, PbcType pbcType)
2733 /* If each molecule is a single charge group
2734 * or we use domain decomposition for each periodic dimension,
2735 * we do not need to take pbc into account for the bonded interactions.
2737 return (pbcType != PbcType::No && dd.comm->systemInfo.haveInterDomainBondeds
2738 && !(dd.numCells[XX] > 1 && dd.numCells[YY] > 1
2739 && (dd.numCells[ZZ] > 1 || pbcType == PbcType::XY)));
2742 /*! \brief Sets grid size limits and PP-PME setup, prints settings to log */
2743 static void set_ddgrid_parameters(const gmx::MDLogger& mdlog,
2746 const gmx_mtop_t& mtop,
2747 const t_inputrec& inputrec,
2748 const gmx_ddbox_t* ddbox)
2750 gmx_domdec_comm_t* comm = dd->comm;
2751 DDRankSetup& ddRankSetup = comm->ddRankSetup;
2753 if (EEL_PME(inputrec.coulombtype) || EVDW_PME(inputrec.vdwtype))
2755 init_ddpme(dd, &ddRankSetup.ddpme[0], 0);
2756 if (ddRankSetup.npmedecompdim >= 2)
2758 init_ddpme(dd, &ddRankSetup.ddpme[1], 1);
2763 ddRankSetup.numRanksDoingPme = 0;
2764 if (dd->pme_nodeid >= 0)
2766 gmx_fatal_collective(FARGS,
2769 "Can not have separate PME ranks without PME electrostatics");
2775 fprintf(debug, "The DD cut-off is %f\n", comm->systemInfo.cutoff);
2777 if (!isDlbDisabled(comm))
2779 set_cell_limits_dlb(mdlog, dd, dlb_scale, inputrec, ddbox);
2782 logSettings(mdlog, dd, mtop, inputrec, dlb_scale, ddbox);
2784 const real vol_frac = (inputrec.pbcType == PbcType::No)
2785 ? (1 - 1 / static_cast<double>(dd->nnodes))
2786 : ((1 + comm_box_frac(dd->numCells, comm->systemInfo.cutoff, *ddbox))
2787 / static_cast<double>(dd->nnodes));
2790 fprintf(debug, "Volume fraction for all DD zones: %f\n", vol_frac);
2792 int natoms_tot = mtop.natoms;
2794 dd->ga2la = new gmx_ga2la_t(natoms_tot, static_cast<int>(vol_frac * natoms_tot));
2797 /*! \brief Get some important DD parameters which can be modified by env.vars */
2798 static DDSettings getDDSettings(const gmx::MDLogger& mdlog,
2799 const DomdecOptions& options,
2800 const gmx::MdrunOptions& mdrunOptions,
2801 const t_inputrec& ir,
2802 const bool directGpuCommUsedWithGpuUpdate)
2804 DDSettings ddSettings;
2806 ddSettings.useSendRecv2 = (dd_getenv(mdlog, "GMX_DD_USE_SENDRECV2", 0) != 0);
2807 ddSettings.dlb_scale_lim = dd_getenv(mdlog, "GMX_DLB_MAX_BOX_SCALING", 10);
2808 ddSettings.useDDOrderZYX = bool(dd_getenv(mdlog, "GMX_DD_ORDER_ZYX", 0));
2809 ddSettings.useCartesianReorder = bool(dd_getenv(mdlog, "GMX_NO_CART_REORDER", 1));
2810 ddSettings.eFlop = dd_getenv(mdlog, "GMX_DLB_BASED_ON_FLOPS", 0);
2811 const int recload = dd_getenv(mdlog, "GMX_DD_RECORD_LOAD", 1);
2812 ddSettings.nstDDDump = dd_getenv(mdlog, "GMX_DD_NST_DUMP", 0);
2813 ddSettings.nstDDDumpGrid = dd_getenv(mdlog, "GMX_DD_NST_DUMP_GRID", 0);
2814 ddSettings.DD_debug = dd_getenv(mdlog, "GMX_DD_DEBUG", 0);
2816 if (ddSettings.useSendRecv2)
2820 "Will use two sequential MPI_Sendrecv calls instead of two simultaneous "
2821 "non-blocking MPI_Irecv and MPI_Isend pairs for constraint and vsite "
2825 if (ddSettings.eFlop)
2827 GMX_LOG(mdlog.info).appendText("Will load balance based on FLOP count");
2828 ddSettings.recordLoad = true;
2832 ddSettings.recordLoad = (wallcycle_have_counter() && recload > 0);
2835 ddSettings.initialDlbState = determineInitialDlbState(
2836 mdlog, options.dlbOption, ddSettings.recordLoad, mdrunOptions, ir, directGpuCommUsedWithGpuUpdate);
2838 .appendTextFormatted("Dynamic load balancing: %s",
2839 enumValueToString(ddSettings.initialDlbState));
2844 gmx_domdec_t::gmx_domdec_t(const t_inputrec& ir) : unitCellInfo(ir) {}
2846 gmx_domdec_t::~gmx_domdec_t() = default;
2851 // TODO once the functionality stablizes, move this class and
2852 // supporting functionality into builder.cpp
2853 /*! \brief Impl class for DD builder */
2854 class DomainDecompositionBuilder::Impl
2858 Impl(const MDLogger& mdlog,
2860 const DomdecOptions& options,
2861 const MdrunOptions& mdrunOptions,
2862 const gmx_mtop_t& mtop,
2863 const t_inputrec& ir,
2864 const MDModulesNotifiers& notifiers,
2866 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType,
2867 bool useUpdateGroups,
2868 real maxUpdateGroupRadius,
2869 ArrayRef<const RVec> xGlobal,
2870 bool useGpuForNonbonded,
2872 bool directGpuCommUsedWithGpuUpdate);
2874 //! Build the resulting DD manager
2875 gmx_domdec_t* build(LocalAtomSetManager* atomSets,
2876 const gmx_localtop_t& localTopology,
2877 const t_state& localState,
2878 ObservablesReducerBuilder* observablesReducerBuilder);
2880 //! Objects used in constructing and configuring DD
2883 const MDLogger& mdlog_;
2884 //! Communication object
2886 //! User-supplied options configuring DD behavior
2887 const DomdecOptions options_;
2888 //! Global system topology
2889 const gmx_mtop_t& mtop_;
2890 //! User input values from the tpr file
2891 const t_inputrec& ir_;
2892 //! MdModules object
2893 const MDModulesNotifiers& notifiers_;
2896 //! Internal objects used in constructing DD
2898 //! Settings combined from the user input
2899 DDSettings ddSettings_;
2900 //! Information derived from the simulation system
2901 DDSystemInfo systemInfo_;
2903 gmx_ddbox_t ddbox_ = { 0 };
2904 //! Organization of the DD grids
2905 DDGridSetup ddGridSetup_;
2906 //! Organzation of the DD ranks
2907 DDRankSetup ddRankSetup_;
2908 //! Number of DD cells in each dimension
2909 ivec ddCellIndex_ = { 0, 0, 0 };
2910 //! IDs of PME-only ranks
2911 std::vector<int> pmeRanks_;
2912 //! Contains a valid Cartesian-communicator-based setup, or defaults.
2913 CartesianRankSetup cartSetup_;
2917 DomainDecompositionBuilder::Impl::Impl(const MDLogger& mdlog,
2919 const DomdecOptions& options,
2920 const MdrunOptions& mdrunOptions,
2921 const gmx_mtop_t& mtop,
2922 const t_inputrec& ir,
2923 const MDModulesNotifiers& notifiers,
2925 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType,
2926 const bool useUpdateGroups,
2927 const real maxUpdateGroupRadius,
2928 ArrayRef<const RVec> xGlobal,
2929 bool useGpuForNonbonded,
2931 bool directGpuCommUsedWithGpuUpdate) :
2932 mdlog_(mdlog), cr_(cr), options_(options), mtop_(mtop), ir_(ir), notifiers_(notifiers)
2934 GMX_LOG(mdlog_.info).appendTextFormatted("\nInitializing Domain Decomposition on %d ranks", cr_->sizeOfDefaultCommunicator);
2936 ddSettings_ = getDDSettings(mdlog_, options_, mdrunOptions, ir_, directGpuCommUsedWithGpuUpdate);
2938 if (ddSettings_.eFlop > 1)
2940 /* Ensure that we have different random flop counts on different ranks */
2941 srand(1 + cr_->rankInDefaultCommunicator);
2944 systemInfo_ = getSystemInfo(mdlog_,
2945 MASTER(cr_) ? DDRole::Master : DDRole::Agent,
2946 cr->mpiDefaultCommunicator,
2951 updateGroupingPerMoleculeType,
2953 maxUpdateGroupRadius,
2956 const int numRanksRequested = cr_->sizeOfDefaultCommunicator;
2957 const bool checkForLargePrimeFactors = (options_.numCells[0] <= 0);
2960 /* Checks for ability to use PME-only ranks */
2961 auto separatePmeRanksPermitted = checkForSeparatePmeRanks(
2962 notifiers_, options_, numRanksRequested, useGpuForNonbonded, useGpuForPme);
2964 /* Checks for validity of requested Ranks setup */
2965 checkForValidRankCountRequests(numRanksRequested,
2966 EEL_PME(ir_.coulombtype) | EVDW_PME(ir_.vdwtype),
2967 options_.numPmeRanks,
2968 separatePmeRanksPermitted,
2969 checkForLargePrimeFactors);
2971 // DD grid setup uses a more different cell size limit for
2972 // automated setup than the one in systemInfo_. The latter is used
2973 // in set_dd_limits() to configure DLB, for example.
2974 const real gridSetupCellsizeLimit =
2975 getDDGridSetupCellSizeLimit(mdlog_,
2976 !isDlbDisabled(ddSettings_.initialDlbState),
2977 options_.dlbScaling,
2979 systemInfo_.cellsizeLimit);
2980 ddGridSetup_ = getDDGridSetup(mdlog_,
2981 MASTER(cr_) ? DDRole::Master : DDRole::Agent,
2982 cr->mpiDefaultCommunicator,
2987 gridSetupCellsizeLimit,
2990 separatePmeRanksPermitted,
2994 checkDDGridSetup(ddGridSetup_,
2995 MASTER(cr_) ? DDRole::Master : DDRole::Agent,
2996 cr->mpiDefaultCommunicator,
2997 cr->sizeOfDefaultCommunicator,
3001 gridSetupCellsizeLimit,
3004 cr_->npmenodes = ddGridSetup_.numPmeOnlyRanks;
3006 ddRankSetup_ = getDDRankSetup(
3007 mdlog_, cr_->sizeOfDefaultCommunicator, options_.rankOrder, ddGridSetup_, ir_);
3009 /* Generate the group communicator, also decides the duty of each rank */
3010 cartSetup_ = makeGroupCommunicators(
3011 mdlog_, ddSettings_, options_.rankOrder, ddRankSetup_, cr_, ddCellIndex_, &pmeRanks_);
3014 gmx_domdec_t* DomainDecompositionBuilder::Impl::build(LocalAtomSetManager* atomSets,
3015 const gmx_localtop_t& localTopology,
3016 const t_state& localState,
3017 ObservablesReducerBuilder* observablesReducerBuilder)
3019 gmx_domdec_t* dd = new gmx_domdec_t(ir_);
3021 copy_ivec(ddCellIndex_, dd->ci);
3023 dd->comm = init_dd_comm();
3025 dd->comm->ddRankSetup = ddRankSetup_;
3026 dd->comm->cartesianRankSetup = cartSetup_;
3028 set_dd_limits(mdlog_,
3029 MASTER(cr_) ? DDRole::Master : DDRole::Agent,
3035 ddRankSetup_.numPPRanks,
3040 setupGroupCommunication(mdlog_, ddSettings_, pmeRanks_, cr_, mtop_.natoms, dd);
3042 if (thisRankHasDuty(cr_, DUTY_PP))
3044 set_ddgrid_parameters(mdlog_, dd, options_.dlbScaling, mtop_, ir_, &ddbox_);
3046 setup_neighbor_relations(dd);
3049 /* Set overallocation to avoid frequent reallocation of arrays */
3050 set_over_alloc_dd(true);
3052 dd->atomSets = atomSets;
3054 dd->localTopologyChecker = std::make_unique<LocalTopologyChecker>(
3055 mdlog_, cr_, mtop_, localTopology, localState, dd->comm->systemInfo.useUpdateGroups, observablesReducerBuilder);
3060 DomainDecompositionBuilder::DomainDecompositionBuilder(const MDLogger& mdlog,
3062 const DomdecOptions& options,
3063 const MdrunOptions& mdrunOptions,
3064 const gmx_mtop_t& mtop,
3065 const t_inputrec& ir,
3066 const MDModulesNotifiers& notifiers,
3068 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType,
3069 const bool useUpdateGroups,
3070 const real maxUpdateGroupRadius,
3071 ArrayRef<const RVec> xGlobal,
3072 const bool useGpuForNonbonded,
3073 const bool useGpuForPme,
3074 const bool directGpuCommUsedWithGpuUpdate) :
3075 impl_(new Impl(mdlog,
3083 updateGroupingPerMoleculeType,
3085 maxUpdateGroupRadius,
3089 directGpuCommUsedWithGpuUpdate))
3093 gmx_domdec_t* DomainDecompositionBuilder::build(LocalAtomSetManager* atomSets,
3094 const gmx_localtop_t& localTopology,
3095 const t_state& localState,
3096 ObservablesReducerBuilder* observablesReducerBuilder)
3098 return impl_->build(atomSets, localTopology, localState, observablesReducerBuilder);
3101 DomainDecompositionBuilder::~DomainDecompositionBuilder() = default;
3105 static gmx_bool test_dd_cutoff(const t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested)
3108 int LocallyLimited = 0;
3110 const auto* dd = cr->dd;
3112 set_ddbox(*dd, false, box, true, x, &ddbox);
3116 for (int d = 0; d < dd->ndim; d++)
3118 const int dim = dd->dim[d];
3120 real inv_cell_size = DD_CELL_MARGIN * dd->numCells[dim] / ddbox.box_size[dim];
3121 if (dd->unitCellInfo.ddBoxIsDynamic)
3123 inv_cell_size *= DD_PRES_SCALE_MARGIN;
3126 const int np = 1 + static_cast<int>(cutoffRequested * inv_cell_size * ddbox.skew_fac[dim]);
3128 if (!isDlbDisabled(dd->comm) && (dim < ddbox.npbcdim) && (dd->comm->cd[d].np_dlb > 0))
3130 if (np > dd->comm->cd[d].np_dlb)
3135 /* If a current local cell size is smaller than the requested
3136 * cut-off, we could still fix it, but this gets very complicated.
3137 * Without fixing here, we might actually need more checks.
3139 real cellSizeAlongDim =
3140 (dd->comm->cell_x1[dim] - dd->comm->cell_x0[dim]) * ddbox.skew_fac[dim];
3141 if (cellSizeAlongDim * dd->comm->cd[d].np_dlb < cutoffRequested)
3148 if (!isDlbDisabled(dd->comm))
3150 /* If DLB is not active yet, we don't need to check the grid jumps.
3151 * Actually we shouldn't, because then the grid jump data is not set.
3153 if (isDlbOn(dd->comm) && gmx::check_grid_jump(0, dd, cutoffRequested, &ddbox, FALSE))
3158 gmx_sumi(1, &LocallyLimited, cr);
3160 if (LocallyLimited > 0)
3169 bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested)
3171 bool bCutoffAllowed = test_dd_cutoff(cr, box, x, cutoffRequested);
3175 cr->dd->comm->systemInfo.cutoff = cutoffRequested;
3178 return bCutoffAllowed;
3181 void constructGpuHaloExchange(const gmx::MDLogger& mdlog,
3182 const t_commrec& cr,
3183 const gmx::DeviceStreamManager& deviceStreamManager,
3184 gmx_wallcycle* wcycle)
3186 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
3187 "Local non-bonded stream should be valid when using"
3188 "GPU halo exchange.");
3189 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
3190 "Non-local non-bonded stream should be valid when using "
3191 "GPU halo exchange.");
3193 if (cr.dd->gpuHaloExchange[0].empty())
3195 GMX_LOG(mdlog.warning)
3197 .appendTextFormatted(
3198 "NOTE: Activating the 'GPU halo exchange' feature, enabled "
3200 "GMX_GPU_DD_COMMS environment variable.");
3203 for (int d = 0; d < cr.dd->ndim; d++)
3205 for (int pulse = cr.dd->gpuHaloExchange[d].size(); pulse < cr.dd->comm->cd[d].numPulses(); pulse++)
3207 cr.dd->gpuHaloExchange[d].push_back(std::make_unique<gmx::GpuHaloExchange>(
3210 cr.mpi_comm_mygroup,
3211 deviceStreamManager.context(),
3212 deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal),
3213 deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal),
3220 void reinitGpuHaloExchange(const t_commrec& cr,
3221 const DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
3222 const DeviceBuffer<gmx::RVec> d_forcesBuffer)
3224 for (int d = 0; d < cr.dd->ndim; d++)
3226 for (int pulse = 0; pulse < cr.dd->comm->cd[d].numPulses(); pulse++)
3228 cr.dd->gpuHaloExchange[d][pulse]->reinitHalo(d_coordinatesBuffer, d_forcesBuffer);
3233 void communicateGpuHaloCoordinates(const t_commrec& cr,
3235 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
3237 for (int d = 0; d < cr.dd->ndim; d++)
3239 for (int pulse = 0; pulse < cr.dd->comm->cd[d].numPulses(); pulse++)
3241 cr.dd->gpuHaloExchange[d][pulse]->communicateHaloCoordinates(box, coordinatesReadyOnDeviceEvent);
3246 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces)
3248 for (int d = cr.dd->ndim - 1; d >= 0; d--)
3250 for (int pulse = cr.dd->comm->cd[d].numPulses() - 1; pulse >= 0; pulse--)
3252 cr.dd->gpuHaloExchange[d][pulse]->communicateHaloForces(accumulateForces);
3257 void dd_init_local_state(const gmx_domdec_t& dd, const t_state* state_global, t_state* state_local)
3259 std::array<int, 5> buf;
3263 buf[0] = state_global->flags;
3264 buf[1] = state_global->ngtc;
3265 buf[2] = state_global->nnhpres;
3266 buf[3] = state_global->nhchainlength;
3267 buf[4] = state_global->dfhist ? state_global->dfhist->nlambda : 0;
3269 dd_bcast(&dd, buf.size() * sizeof(int), buf.data());
3271 init_gtc_state(state_local, buf[1], buf[2], buf[3]);
3272 init_dfhist_state(state_local, buf[4]);
3273 state_local->flags = buf[0];
3276 void putUpdateGroupAtomsInSamePeriodicImage(const gmx_domdec_t& dd,
3277 const gmx_mtop_t& mtop,
3279 gmx::ArrayRef<gmx::RVec> positions)
3282 for (const gmx_molblock_t& molblock : mtop.molblock)
3284 const auto& updateGrouping = dd.comm->systemInfo.updateGroupingsPerMoleculeType[molblock.type];
3286 for (int mol = 0; mol < molblock.nmol; mol++)
3288 for (int g = 0; g < updateGrouping.numBlocks(); g++)
3290 const auto& block = updateGrouping.block(g);
3291 const int atomBegin = atomOffset + block.begin();
3292 const int atomEnd = atomOffset + block.end();
3293 for (int a = atomBegin + 1; a < atomEnd; a++)
3295 // Make sure that atoms in the same update group
3296 // are in the same periodic image after restarts.
3297 for (int d = DIM - 1; d >= 0; d--)
3299 while (positions[a][d] - positions[atomBegin][d] > 0.5_real * box[d][d])
3301 positions[a] -= box[d];
3303 while (positions[a][d] - positions[atomBegin][d] < -0.5_real * box[d][d])
3305 positions[a] += box[d];
3310 atomOffset += updateGrouping.fullRange().end();