2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief This file defines functions for mdrun to call to make a new
38 * domain decomposition, and check it.
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_domdec
46 #include "partition.h"
55 #include "gromacs/domdec/collect.h"
56 #include "gromacs/domdec/dlb.h"
57 #include "gromacs/domdec/dlbtiming.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_network.h"
60 #include "gromacs/domdec/ga2la.h"
61 #include "gromacs/domdec/localatomsetmanager.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/ewald/pme_pp.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/math/functions.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/mdlib/forcerec.h"
70 #include "gromacs/mdlib/gmx_omp_nthreads.h"
71 #include "gromacs/mdlib/mdatoms.h"
72 #include "gromacs/mdlib/nsgrid.h"
73 #include "gromacs/mdlib/vsite.h"
74 #include "gromacs/mdtypes/commrec.h"
75 #include "gromacs/mdtypes/forcerec.h"
76 #include "gromacs/mdtypes/inputrec.h"
77 #include "gromacs/mdtypes/md_enums.h"
78 #include "gromacs/mdtypes/mdatom.h"
79 #include "gromacs/mdtypes/nblist.h"
80 #include "gromacs/mdtypes/state.h"
81 #include "gromacs/nbnxm/nbnxm.h"
82 #include "gromacs/pulling/pull.h"
83 #include "gromacs/timing/wallcycle.h"
84 #include "gromacs/topology/mtop_util.h"
85 #include "gromacs/topology/topology.h"
86 #include "gromacs/utility/cstringutil.h"
87 #include "gromacs/utility/fatalerror.h"
88 #include "gromacs/utility/logger.h"
89 #include "gromacs/utility/real.h"
90 #include "gromacs/utility/smalloc.h"
91 #include "gromacs/utility/strconvert.h"
92 #include "gromacs/utility/stringstream.h"
93 #include "gromacs/utility/stringutil.h"
94 #include "gromacs/utility/textwriter.h"
97 #include "cellsizes.h"
98 #include "distribute.h"
99 #include "domdec_constraints.h"
100 #include "domdec_internal.h"
101 #include "domdec_vsite.h"
103 #include "redistribute.h"
106 /*! \brief Turn on DLB when the load imbalance causes this amount of total loss.
108 * There is a bit of overhead with DLB and it's difficult to achieve
109 * a load imbalance of less than 2% with DLB.
111 #define DD_PERF_LOSS_DLB_ON 0.02
113 //! Warn about imbalance due to PP or PP/PME load imbalance at this loss.
114 #define DD_PERF_LOSS_WARN 0.05
117 //! Debug helper printing a DD zone
118 static void print_ddzone(FILE* fp, int d, int i, int j, gmx_ddzone_t* zone)
121 "zone d0 %d d1 %d d2 %d min0 %6.3f max1 %6.3f mch0 %6.3f mch1 %6.3f p1_0 %6.3f p1_1 "
134 /*! \brief Using the home grid size as input in cell_ns_x0 and cell_ns_x1
135 * takes the extremes over all home and remote zones in the halo
136 * and returns the results in cell_ns_x0 and cell_ns_x1.
137 * Note: only used with the group cut-off scheme.
139 static void dd_move_cellx(gmx_domdec_t* dd, const gmx_ddbox_t* ddbox, rvec cell_ns_x0, rvec cell_ns_x1)
141 constexpr int c_ddZoneCommMaxNumZones = 5;
142 gmx_ddzone_t buf_s[c_ddZoneCommMaxNumZones];
143 gmx_ddzone_t buf_r[c_ddZoneCommMaxNumZones];
144 gmx_ddzone_t buf_e[c_ddZoneCommMaxNumZones];
145 gmx_domdec_comm_t* comm = dd->comm;
149 for (int d = 1; d < dd->ndim; d++)
151 int dim = dd->dim[d];
152 gmx_ddzone_t& zp = (d == 1) ? comm->zone_d1[0] : comm->zone_d2[0][0];
154 /* Copy the base sizes of the home zone */
155 zp.min0 = cell_ns_x0[dim];
156 zp.max1 = cell_ns_x1[dim];
157 zp.min1 = cell_ns_x1[dim];
158 zp.mch0 = cell_ns_x0[dim];
159 zp.mch1 = cell_ns_x1[dim];
160 zp.p1_0 = cell_ns_x0[dim];
161 zp.p1_1 = cell_ns_x1[dim];
165 gmx::ArrayRef<DDCellsizesWithDlb> cellsizes = comm->cellsizesWithDlb;
167 /* Loop backward over the dimensions and aggregate the extremes
170 for (int d = dd->ndim - 2; d >= 0; d--)
172 const int dim = dd->dim[d];
173 const bool applyPbc = (dim < ddbox->npbcdim);
175 /* Use an rvec to store two reals */
176 extr_s[d][0] = cellsizes[d + 1].fracLower;
177 extr_s[d][1] = cellsizes[d + 1].fracUpper;
178 extr_s[d][2] = cellsizes[d + 1].fracUpper;
181 GMX_ASSERT(pos < c_ddZoneCommMaxNumZones, "The buffers should be sufficiently large");
182 /* Store the extremes in the backward sending buffer,
183 * so they get updated separately from the forward communication.
185 for (int d1 = d; d1 < dd->ndim - 1; d1++)
187 gmx_ddzone_t& buf = buf_s[pos];
189 /* We invert the order to be able to use the same loop for buf_e */
190 buf.min0 = extr_s[d1][1];
191 buf.max1 = extr_s[d1][0];
192 buf.min1 = extr_s[d1][2];
195 /* Store the cell corner of the dimension we communicate along */
196 buf.p1_0 = comm->cell_x0[dim];
202 buf_s[pos] = (dd->ndim == 2) ? comm->zone_d1[0] : comm->zone_d2[0][0];
205 if (dd->ndim == 3 && d == 0)
207 buf_s[pos] = comm->zone_d2[0][1];
209 buf_s[pos] = comm->zone_d1[0];
213 /* We only need to communicate the extremes
214 * in the forward direction
216 int numPulses = comm->cd[d].numPulses();
220 /* Take the minimum to avoid double communication */
221 numPulsesMin = std::min(numPulses, dd->numCells[dim] - 1 - numPulses);
225 /* Without PBC we should really not communicate over
226 * the boundaries, but implementing that complicates
227 * the communication setup and therefore we simply
228 * do all communication, but ignore some data.
230 numPulsesMin = numPulses;
232 for (int pulse = 0; pulse < numPulsesMin; pulse++)
234 /* Communicate the extremes forward */
235 bool receiveValidData = (applyPbc || dd->ci[dim] > 0);
237 int numElements = dd->ndim - d - 1;
238 ddSendrecv(dd, d, dddirForward, extr_s + d, numElements, extr_r + d, numElements);
240 if (receiveValidData)
242 for (int d1 = d; d1 < dd->ndim - 1; d1++)
244 extr_s[d1][0] = std::max(extr_s[d1][0], extr_r[d1][0]);
245 extr_s[d1][1] = std::min(extr_s[d1][1], extr_r[d1][1]);
246 extr_s[d1][2] = std::min(extr_s[d1][2], extr_r[d1][2]);
251 const int numElementsInBuffer = pos;
252 for (int pulse = 0; pulse < numPulses; pulse++)
254 /* Communicate all the zone information backward */
255 bool receiveValidData = (applyPbc || dd->ci[dim] < dd->numCells[dim] - 1);
258 sizeof(gmx_ddzone_t) == c_ddzoneNumReals * sizeof(real),
259 "Here we expect gmx_ddzone_t to consist of c_ddzoneNumReals reals (only)");
261 int numReals = numElementsInBuffer * c_ddzoneNumReals;
265 gmx::arrayRefFromArray(&buf_s[0].min0, numReals),
266 gmx::arrayRefFromArray(&buf_r[0].min0, numReals));
271 for (int d1 = d + 1; d1 < dd->ndim; d1++)
273 /* Determine the decrease of maximum required
274 * communication height along d1 due to the distance along d,
275 * this avoids a lot of useless atom communication.
277 real dist_d = comm->cell_x1[dim] - buf_r[0].p1_0;
280 if (ddbox->tric_dir[dim])
282 /* c is the off-diagonal coupling between the cell planes
283 * along directions d and d1.
285 c = ddbox->v[dim][dd->dim[d1]][dim];
291 real det = (1 + c * c) * gmx::square(comm->systemInfo.cutoff) - dist_d * dist_d;
294 dh[d1] = comm->systemInfo.cutoff - (c * dist_d + std::sqrt(det)) / (1 + c * c);
298 /* A negative value signals out of range */
304 /* Accumulate the extremes over all pulses */
305 for (int i = 0; i < numElementsInBuffer; i++)
313 if (receiveValidData)
315 buf_e[i].min0 = std::min(buf_e[i].min0, buf_r[i].min0);
316 buf_e[i].max1 = std::max(buf_e[i].max1, buf_r[i].max1);
317 buf_e[i].min1 = std::min(buf_e[i].min1, buf_r[i].min1);
321 if (dd->ndim == 3 && d == 0 && i == numElementsInBuffer - 1)
329 if (receiveValidData && dh[d1] >= 0)
331 buf_e[i].mch0 = std::max(buf_e[i].mch0, buf_r[i].mch0 - dh[d1]);
332 buf_e[i].mch1 = std::max(buf_e[i].mch1, buf_r[i].mch1 - dh[d1]);
335 /* Copy the received buffer to the send buffer,
336 * to pass the data through with the next pulse.
340 if (((applyPbc || dd->ci[dim] + numPulses < dd->numCells[dim]) && pulse == numPulses - 1)
341 || (!applyPbc && dd->ci[dim] + 1 + pulse == dd->numCells[dim] - 1))
343 /* Store the extremes */
346 for (int d1 = d; d1 < dd->ndim - 1; d1++)
348 extr_s[d1][1] = std::min(extr_s[d1][1], buf_e[pos].min0);
349 extr_s[d1][0] = std::max(extr_s[d1][0], buf_e[pos].max1);
350 extr_s[d1][2] = std::min(extr_s[d1][2], buf_e[pos].min1);
354 if (d == 1 || (d == 0 && dd->ndim == 3))
356 for (int i = d; i < 2; i++)
358 comm->zone_d2[1 - d][i] = buf_e[pos];
364 comm->zone_d1[1] = buf_e[pos];
370 if (d == 1 || (d == 0 && dd->ndim == 3))
372 for (int i = d; i < 2; i++)
374 comm->zone_d2[1 - d][i].dataSet = 0;
379 comm->zone_d1[1].dataSet = 0;
387 int dim = dd->dim[1];
388 for (int i = 0; i < 2; i++)
390 if (comm->zone_d1[i].dataSet != 0)
394 print_ddzone(debug, 1, i, 0, &comm->zone_d1[i]);
396 cell_ns_x0[dim] = std::min(cell_ns_x0[dim], comm->zone_d1[i].min0);
397 cell_ns_x1[dim] = std::max(cell_ns_x1[dim], comm->zone_d1[i].max1);
403 int dim = dd->dim[2];
404 for (int i = 0; i < 2; i++)
406 for (int j = 0; j < 2; j++)
408 if (comm->zone_d2[i][j].dataSet != 0)
412 print_ddzone(debug, 2, i, j, &comm->zone_d2[i][j]);
414 cell_ns_x0[dim] = std::min(cell_ns_x0[dim], comm->zone_d2[i][j].min0);
415 cell_ns_x1[dim] = std::max(cell_ns_x1[dim], comm->zone_d2[i][j].max1);
420 for (int d = 1; d < dd->ndim; d++)
422 cellsizes[d].fracLowerMax = extr_s[d - 1][0];
423 cellsizes[d].fracUpperMin = extr_s[d - 1][1];
427 "Cell fraction d %d, max0 %f, min1 %f\n",
429 cellsizes[d].fracLowerMax,
430 cellsizes[d].fracUpperMin);
435 //! Sets the charge-group zones to be equal to the home zone.
436 static void set_zones_ncg_home(gmx_domdec_t* dd)
438 gmx_domdec_zones_t* zones;
441 zones = &dd->comm->zones;
443 zones->cg_range[0] = 0;
444 for (i = 1; i < zones->n + 1; i++)
446 zones->cg_range[i] = dd->ncg_home;
448 /* zone_ncg1[0] should always be equal to ncg_home */
449 dd->comm->zone_ncg1[0] = dd->ncg_home;
452 //! Restore atom groups for the charge groups.
453 static void restoreAtomGroups(gmx_domdec_t* dd, const t_state* state)
455 gmx::ArrayRef<const int> atomGroupsState = state->cg_gl;
457 std::vector<int>& globalAtomGroupIndices = dd->globalAtomGroupIndices;
459 globalAtomGroupIndices.resize(atomGroupsState.size());
461 /* Copy back the global charge group indices from state
462 * and rebuild the local charge group to atom index.
464 for (gmx::index i = 0; i < atomGroupsState.ssize(); i++)
466 globalAtomGroupIndices[i] = atomGroupsState[i];
469 dd->ncg_home = atomGroupsState.size();
470 dd->comm->atomRanges.setEnd(DDAtomRanges::Type::Home, atomGroupsState.ssize());
472 set_zones_ncg_home(dd);
475 //! Sets the cginfo structures.
476 static void dd_set_cginfo(gmx::ArrayRef<const int> index_gl, int cg0, int cg1, t_forcerec* fr)
480 gmx::ArrayRef<cginfo_mb_t> cginfo_mb = fr->cginfo_mb;
481 gmx::ArrayRef<int> cginfo = fr->cginfo;
483 for (int cg = cg0; cg < cg1; cg++)
485 cginfo[cg] = ddcginfo(cginfo_mb, index_gl[cg]);
490 //! Makes the mappings between global and local atom indices during DD repartioning.
491 static void make_dd_indices(gmx_domdec_t* dd, const int atomStart)
493 const int numZones = dd->comm->zones.n;
494 gmx::ArrayRef<const int> zone2cg = dd->comm->zones.cg_range;
495 gmx::ArrayRef<const int> zone_ncg1 = dd->comm->zone_ncg1;
496 gmx::ArrayRef<const int> globalAtomGroupIndices = dd->globalAtomGroupIndices;
498 std::vector<int>& globalAtomIndices = dd->globalAtomIndices;
499 gmx_ga2la_t& ga2la = *dd->ga2la;
501 if (zone2cg[1] != dd->ncg_home)
503 gmx_incons("dd->ncg_zone is not up to date");
506 /* Make the local to global and global to local atom index */
508 globalAtomIndices.resize(a);
509 for (int zone = 0; zone < numZones; zone++)
520 int cg1 = zone2cg[zone + 1];
521 int cg1_p1 = cg0 + zone_ncg1[zone];
523 for (int cg = cg0; cg < cg1; cg++)
528 /* Signal that this cg is from more than one pulse away */
531 int cg_gl = globalAtomGroupIndices[cg];
532 globalAtomIndices.push_back(cg_gl);
533 ga2la.insert(cg_gl, { a, zone1 });
539 //! Checks whether global and local atom indices are consistent.
540 static void check_index_consistency(const gmx_domdec_t* dd, int natoms_sys, const char* where)
544 const int numAtomsInZones = dd->comm->atomRanges.end(DDAtomRanges::Type::Zones);
546 if (dd->comm->ddSettings.DD_debug > 1)
548 std::vector<int> have(natoms_sys);
549 for (int a = 0; a < numAtomsInZones; a++)
551 int globalAtomIndex = dd->globalAtomIndices[a];
552 if (have[globalAtomIndex] > 0)
555 "DD rank %d: global atom %d occurs twice: index %d and %d\n",
558 have[globalAtomIndex],
563 have[globalAtomIndex] = a + 1;
568 std::vector<int> have(numAtomsInZones);
571 for (int i = 0; i < natoms_sys; i++)
573 if (const auto entry = dd->ga2la->find(i))
575 const int a = entry->la;
576 if (a >= numAtomsInZones)
579 "DD rank %d: global atom %d marked as local atom %d, which is larger than "
590 if (dd->globalAtomIndices[a] != i)
593 "DD rank %d: global atom %d marked as local atom %d, which has global "
598 dd->globalAtomIndices[a] + 1);
605 if (ngl != numAtomsInZones)
607 fprintf(stderr, "DD rank %d, %s: %d global atom indices, %d local atoms\n", dd->rank, where, ngl, numAtomsInZones);
609 for (int a = 0; a < numAtomsInZones; a++)
614 "DD rank %d, %s: local atom %d, global %d has no global index\n",
618 dd->globalAtomIndices[a] + 1);
624 gmx_fatal(FARGS, "DD rank %d, %s: %d atom(group) index inconsistencies", dd->rank, where, nerr);
628 //! Clear all DD global state indices
629 static void clearDDStateIndices(gmx_domdec_t* dd, const bool keepLocalAtomIndices)
631 gmx_ga2la_t& ga2la = *dd->ga2la;
633 if (!keepLocalAtomIndices)
635 /* Clear the whole list without the overhead of searching */
640 const int numAtomsInZones = dd->comm->atomRanges.end(DDAtomRanges::Type::Zones);
641 for (int i = 0; i < numAtomsInZones; i++)
643 ga2la.erase(dd->globalAtomIndices[i]);
647 dd_clear_local_vsite_indices(dd);
651 dd_clear_local_constraint_indices(dd);
655 bool check_grid_jump(int64_t step, const gmx_domdec_t* dd, real cutoff, const gmx_ddbox_t* ddbox, gmx_bool bFatal)
657 gmx_domdec_comm_t* comm = dd->comm;
658 bool invalid = false;
660 for (int d = 1; d < dd->ndim; d++)
662 const DDCellsizesWithDlb& cellsizes = comm->cellsizesWithDlb[d];
663 const int dim = dd->dim[d];
664 const real limit = grid_jump_limit(comm, cutoff, d);
665 real bfac = ddbox->box_size[dim];
666 if (ddbox->tric_dir[dim])
668 bfac *= ddbox->skew_fac[dim];
670 if ((cellsizes.fracUpper - cellsizes.fracLowerMax) * bfac < limit
671 || (cellsizes.fracLower - cellsizes.fracUpperMin) * bfac > -limit)
679 /* This error should never be triggered under normal
680 * circumstances, but you never know ...
683 "step %s: The domain decomposition grid has shifted too much in the "
684 "%c-direction around cell %d %d %d. This should not have happened. "
685 "Running with fewer ranks might avoid this issue.",
686 gmx_step_str(step, buf),
697 //! Return the duration of force calculations on this rank.
698 static float dd_force_load(gmx_domdec_comm_t* comm)
702 if (comm->ddSettings.eFlop)
705 if (comm->ddSettings.eFlop > 1)
707 load *= 1.0 + (comm->ddSettings.eFlop - 1) * (0.1 * rand() / RAND_MAX - 0.05);
712 load = comm->cycl[ddCyclF];
713 if (comm->cycl_n[ddCyclF] > 1)
715 /* Subtract the maximum of the last n cycle counts
716 * to get rid of possible high counts due to other sources,
717 * for instance system activity, that would otherwise
718 * affect the dynamic load balancing.
720 load -= comm->cycl_max[ddCyclF];
724 if (comm->cycl_n[ddCyclWaitGPU] && comm->nrank_gpu_shared > 1)
726 float gpu_wait, gpu_wait_sum;
728 gpu_wait = comm->cycl[ddCyclWaitGPU];
729 if (comm->cycl_n[ddCyclF] > 1)
731 /* We should remove the WaitGPU time of the same MD step
732 * as the one with the maximum F time, since the F time
733 * and the wait time are not independent.
734 * Furthermore, the step for the max F time should be chosen
735 * the same on all ranks that share the same GPU.
736 * But to keep the code simple, we remove the average instead.
737 * The main reason for artificially long times at some steps
738 * is spurious CPU activity or MPI time, so we don't expect
739 * that changes in the GPU wait time matter a lot here.
741 gpu_wait *= (comm->cycl_n[ddCyclF] - 1) / static_cast<float>(comm->cycl_n[ddCyclF]);
743 /* Sum the wait times over the ranks that share the same GPU */
744 MPI_Allreduce(&gpu_wait, &gpu_wait_sum, 1, MPI_FLOAT, MPI_SUM, comm->mpi_comm_gpu_shared);
745 /* Replace the wait time by the average over the ranks */
746 load += -gpu_wait + gpu_wait_sum / comm->nrank_gpu_shared;
754 //! Runs cell size checks and communicates the boundaries.
755 static void comm_dd_ns_cell_sizes(gmx_domdec_t* dd, gmx_ddbox_t* ddbox, rvec cell_ns_x0, rvec cell_ns_x1, int64_t step)
757 gmx_domdec_comm_t* comm;
762 for (dim_ind = 0; dim_ind < dd->ndim; dim_ind++)
764 dim = dd->dim[dim_ind];
766 /* Without PBC we don't have restrictions on the outer cells */
767 if (!(dim >= ddbox->npbcdim && (dd->ci[dim] == 0 || dd->ci[dim] == dd->numCells[dim] - 1))
769 && (comm->cell_x1[dim] - comm->cell_x0[dim]) * ddbox->skew_fac[dim] < comm->cellsize_min[dim])
773 "step %s: The %c-size (%f) times the triclinic skew factor (%f) is smaller "
774 "than the smallest allowed cell size (%f) for domain decomposition grid cell "
776 gmx_step_str(step, buf),
778 comm->cell_x1[dim] - comm->cell_x0[dim],
779 ddbox->skew_fac[dim],
780 dd->comm->cellsize_min[dim],
787 if ((isDlbOn(dd->comm) && dd->ndim > 1) || ddbox->nboundeddim < DIM)
789 /* Communicate the boundaries and update cell_ns_x0/1 */
790 dd_move_cellx(dd, ddbox, cell_ns_x0, cell_ns_x1);
791 if (isDlbOn(dd->comm) && dd->ndim > 1)
793 check_grid_jump(step, dd, dd->comm->systemInfo.cutoff, ddbox, TRUE);
798 //! Compute and communicate to determine the load distribution across PP ranks.
799 static void get_load_distribution(gmx_domdec_t* dd, gmx_wallcycle_t wcycle)
801 gmx_domdec_comm_t* comm;
803 float cell_frac = 0, sbuf[DD_NLOAD_MAX];
808 fprintf(debug, "get_load_distribution start\n");
811 wallcycle_start(wcycle, ewcDDCOMMLOAD);
815 bSepPME = (dd->pme_nodeid >= 0);
817 if (dd->ndim == 0 && bSepPME)
819 /* Without decomposition, but with PME nodes, we need the load */
820 comm->load[0].mdf = comm->cycl[ddCyclPPduringPME];
821 comm->load[0].pme = comm->cycl[ddCyclPME];
824 // Either we have DLB off, or we have it on and the array is large enough
825 GMX_ASSERT(!isDlbOn(dd->comm) || static_cast<int>(dd->comm->cellsizesWithDlb.size()) == dd->ndim,
826 "DLB cell sizes data not set up properly ");
827 for (int d = dd->ndim - 1; d >= 0; d--)
829 const int dim = dd->dim[d];
830 /* Check if we participate in the communication in this dimension */
831 if (d == dd->ndim - 1 || (dd->ci[dd->dim[d + 1]] == 0 && dd->ci[dd->dim[dd->ndim - 1]] == 0))
833 load = &comm->load[d];
834 if (isDlbOn(dd->comm))
836 cell_frac = comm->cellsizesWithDlb[d].fracUpper - comm->cellsizesWithDlb[d].fracLower;
839 if (d == dd->ndim - 1)
841 sbuf[pos++] = dd_force_load(comm);
842 sbuf[pos++] = sbuf[0];
843 if (isDlbOn(dd->comm))
845 sbuf[pos++] = sbuf[0];
846 sbuf[pos++] = cell_frac;
849 sbuf[pos++] = comm->cellsizesWithDlb[d].fracLowerMax;
850 sbuf[pos++] = comm->cellsizesWithDlb[d].fracUpperMin;
855 sbuf[pos++] = comm->cycl[ddCyclPPduringPME];
856 sbuf[pos++] = comm->cycl[ddCyclPME];
861 sbuf[pos++] = comm->load[d + 1].sum;
862 sbuf[pos++] = comm->load[d + 1].max;
863 if (isDlbOn(dd->comm))
865 sbuf[pos++] = comm->load[d + 1].sum_m;
866 sbuf[pos++] = comm->load[d + 1].cvol_min * cell_frac;
867 sbuf[pos++] = comm->load[d + 1].flags;
870 sbuf[pos++] = comm->cellsizesWithDlb[d].fracLowerMax;
871 sbuf[pos++] = comm->cellsizesWithDlb[d].fracUpperMin;
876 sbuf[pos++] = comm->load[d + 1].mdf;
877 sbuf[pos++] = comm->load[d + 1].pme;
881 /* Communicate a row in DD direction d.
882 * The communicators are setup such that the root always has rank 0.
886 load->nload * sizeof(float),
889 load->nload * sizeof(float),
892 comm->mpi_comm_load[d]);
894 if (dd->ci[dim] == dd->master_ci[dim])
896 /* We are the master along this row, process this row */
897 RowMaster* rowMaster = nullptr;
901 rowMaster = comm->cellsizesWithDlb[d].rowMaster.get();
911 for (int i = 0; i < dd->numCells[dim]; i++)
913 load->sum += load->load[pos++];
914 load->max = std::max(load->max, load->load[pos]);
916 if (isDlbOn(dd->comm))
918 if (rowMaster->dlbIsLimited)
920 /* This direction could not be load balanced properly,
921 * therefore we need to use the maximum iso the average load.
923 load->sum_m = std::max(load->sum_m, load->load[pos]);
927 load->sum_m += load->load[pos];
930 load->cvol_min = std::min(load->cvol_min, load->load[pos]);
932 if (d < dd->ndim - 1)
934 load->flags = gmx::roundToInt(load->load[pos++]);
938 rowMaster->bounds[i].cellFracLowerMax = load->load[pos++];
939 rowMaster->bounds[i].cellFracUpperMin = load->load[pos++];
944 load->mdf = std::max(load->mdf, load->load[pos]);
946 load->pme = std::max(load->pme, load->load[pos]);
950 if (isDlbOn(comm) && rowMaster->dlbIsLimited)
952 load->sum_m *= dd->numCells[dim];
953 load->flags |= (1 << d);
961 comm->nload += dd_load_count(comm);
962 comm->load_step += comm->cycl[ddCyclStep];
963 comm->load_sum += comm->load[0].sum;
964 comm->load_max += comm->load[0].max;
967 for (int d = 0; d < dd->ndim; d++)
969 if (comm->load[0].flags & (1 << d))
977 comm->load_mdf += comm->load[0].mdf;
978 comm->load_pme += comm->load[0].pme;
982 wallcycle_stop(wcycle, ewcDDCOMMLOAD);
986 fprintf(debug, "get_load_distribution finished\n");
990 /*! \brief Return the relative performance loss on the total run time
991 * due to the force calculation load imbalance. */
992 static float dd_force_load_fraction(gmx_domdec_t* dd)
994 if (dd->comm->nload > 0 && dd->comm->load_step > 0)
996 return dd->comm->load_sum / (dd->comm->load_step * dd->nnodes);
1004 /*! \brief Return the relative performance loss on the total run time
1005 * due to the force calculation load imbalance. */
1006 static float dd_force_imb_perf_loss(gmx_domdec_t* dd)
1008 if (dd->comm->nload > 0 && dd->comm->load_step > 0)
1010 return (dd->comm->load_max * dd->nnodes - dd->comm->load_sum) / (dd->comm->load_step * dd->nnodes);
1018 //! Print load-balance report e.g. at the end of a run.
1019 static void print_dd_load_av(FILE* fplog, gmx_domdec_t* dd)
1021 gmx_domdec_comm_t* comm = dd->comm;
1023 /* Only the master rank prints loads and only if we measured loads */
1024 if (!DDMASTER(dd) || comm->nload == 0)
1030 int numPpRanks = dd->nnodes;
1031 int numPmeRanks = (comm->ddRankSetup.usePmeOnlyRanks ? comm->ddRankSetup.numRanksDoingPme : 0);
1032 int numRanks = numPpRanks + numPmeRanks;
1033 float lossFraction = 0;
1035 /* Print the average load imbalance and performance loss */
1036 if (dd->nnodes > 1 && comm->load_sum > 0)
1038 float imbalance = comm->load_max * numPpRanks / comm->load_sum - 1;
1039 lossFraction = dd_force_imb_perf_loss(dd);
1041 std::string msg = "\nDynamic load balancing report:\n";
1042 std::string dlbStateStr;
1044 switch (dd->comm->dlbState)
1046 case DlbState::offUser:
1047 dlbStateStr = "DLB was off during the run per user request.";
1049 case DlbState::offForever:
1050 /* Currectly this can happen due to performance loss observed, cell size
1051 * limitations or incompatibility with other settings observed during
1052 * determineInitialDlbState(). */
1053 dlbStateStr = "DLB got disabled because it was unsuitable to use.";
1055 case DlbState::offCanTurnOn:
1056 dlbStateStr = "DLB was off during the run due to low measured imbalance.";
1058 case DlbState::offTemporarilyLocked:
1060 "DLB was locked at the end of the run due to unfinished PP-PME "
1063 case DlbState::onCanTurnOff:
1064 dlbStateStr = "DLB was turned on during the run due to measured imbalance.";
1066 case DlbState::onUser:
1067 dlbStateStr = "DLB was permanently on during the run per user request.";
1069 default: GMX_ASSERT(false, "Undocumented DLB state");
1072 msg += " " + dlbStateStr + "\n";
1073 msg += gmx::formatString(" Average load imbalance: %.1f%%.\n", imbalance * 100);
1074 msg += gmx::formatString(
1075 " The balanceable part of the MD step is %d%%, load imbalance is computed from "
1077 gmx::roundToInt(dd_force_load_fraction(dd) * 100));
1078 msg += gmx::formatString(
1079 " Part of the total run time spent waiting due to load imbalance: %.1f%%.\n",
1080 lossFraction * 100);
1081 fprintf(fplog, "%s", msg.c_str());
1082 fprintf(stderr, "\n%s", msg.c_str());
1085 /* Print during what percentage of steps the load balancing was limited */
1086 bool dlbWasLimited = false;
1089 sprintf(buf, " Steps where the load balancing was limited by -rdd, -rcon and/or -dds:");
1090 for (int d = 0; d < dd->ndim; d++)
1092 int limitPercentage = (200 * comm->load_lim[d] + 1) / (2 * comm->nload);
1093 sprintf(buf + strlen(buf), " %c %d %%", dim2char(dd->dim[d]), limitPercentage);
1094 if (limitPercentage >= 50)
1096 dlbWasLimited = true;
1099 sprintf(buf + strlen(buf), "\n");
1100 fprintf(fplog, "%s", buf);
1101 fprintf(stderr, "%s", buf);
1104 /* Print the performance loss due to separate PME - PP rank imbalance */
1105 float lossFractionPme = 0;
1106 if (numPmeRanks > 0 && comm->load_mdf > 0 && comm->load_step > 0)
1108 float pmeForceRatio = comm->load_pme / comm->load_mdf;
1109 lossFractionPme = (comm->load_pme - comm->load_mdf) / comm->load_step;
1110 if (lossFractionPme <= 0)
1112 lossFractionPme *= numPmeRanks / static_cast<float>(numRanks);
1116 lossFractionPme *= numPpRanks / static_cast<float>(numRanks);
1118 sprintf(buf, " Average PME mesh/force load: %5.3f\n", pmeForceRatio);
1119 fprintf(fplog, "%s", buf);
1120 fprintf(stderr, "%s", buf);
1122 " Part of the total run time spent waiting due to PP/PME imbalance: %.1f %%\n",
1123 std::fabs(lossFractionPme) * 100);
1124 fprintf(fplog, "%s", buf);
1125 fprintf(stderr, "%s", buf);
1127 fprintf(fplog, "\n");
1128 fprintf(stderr, "\n");
1130 if ((lossFraction >= DD_PERF_LOSS_WARN) && (dd->comm->dlbState != DlbState::offTemporarilyLocked))
1132 std::string message = gmx::formatString(
1133 "NOTE: %.1f %% of the available CPU time was lost due to load imbalance\n"
1134 " in the domain decomposition.\n",
1135 lossFraction * 100);
1137 bool hadSuggestion = false;
1138 if (dd->comm->dlbState == DlbState::offUser)
1140 message += " You might want to allow dynamic load balancing (option -dlb auto.)\n";
1141 hadSuggestion = true;
1143 else if (dd->comm->dlbState == DlbState::offCanTurnOn)
1146 " Dynamic load balancing was automatically disabled, but it might be "
1147 "beneficial to manually tuning it on (option -dlb on.)\n";
1148 hadSuggestion = true;
1150 else if (dlbWasLimited)
1153 " You might want to decrease the cell size limit (options -rdd, -rcon "
1155 hadSuggestion = true;
1157 message += gmx::formatString(
1158 " You can %sconsider manually changing the decomposition (option -dd);\n"
1159 " e.g. by using fewer domains along the box dimension in which there is\n"
1160 " considerable inhomogeneity in the simulated system.",
1161 hadSuggestion ? "also " : "");
1163 fprintf(fplog, "%s\n", message.c_str());
1164 fprintf(stderr, "%s\n", message.c_str());
1166 if (numPmeRanks > 0 && std::fabs(lossFractionPme) >= DD_PERF_LOSS_WARN)
1169 "NOTE: %.1f %% performance was lost because the PME ranks\n"
1170 " had %s work to do than the PP ranks.\n"
1171 " You might want to %s the number of PME ranks\n"
1172 " or %s the cut-off and the grid spacing.\n",
1173 std::fabs(lossFractionPme * 100),
1174 (lossFractionPme < 0) ? "less" : "more",
1175 (lossFractionPme < 0) ? "decrease" : "increase",
1176 (lossFractionPme < 0) ? "decrease" : "increase");
1177 fprintf(fplog, "%s\n", buf);
1178 fprintf(stderr, "%s\n", buf);
1182 //! Return the minimum communication volume.
1183 static float dd_vol_min(gmx_domdec_t* dd)
1185 return dd->comm->load[0].cvol_min * dd->nnodes;
1188 //! Return the DD load flags.
1189 static int dd_load_flags(gmx_domdec_t* dd)
1191 return dd->comm->load[0].flags;
1194 //! Return the reported load imbalance in force calculations.
1195 static float dd_f_imbal(gmx_domdec_t* dd)
1197 if (dd->comm->load[0].sum > 0)
1199 return dd->comm->load[0].max * dd->nnodes / dd->comm->load[0].sum - 1.0F;
1203 /* Something is wrong in the cycle counting, report no load imbalance */
1208 //! Returns DD load balance report.
1209 static std::string dd_print_load(gmx_domdec_t* dd, int64_t step)
1211 gmx::StringOutputStream stream;
1212 gmx::TextWriter log(&stream);
1214 int flags = dd_load_flags(dd);
1217 log.writeString("DD load balancing is limited by minimum cell size in dimension");
1218 for (int d = 0; d < dd->ndim; d++)
1220 if (flags & (1 << d))
1222 log.writeStringFormatted(" %c", dim2char(dd->dim[d]));
1225 log.ensureLineBreak();
1227 log.writeString("DD step " + gmx::toString(step));
1228 if (isDlbOn(dd->comm))
1230 log.writeStringFormatted(" vol min/aver %5.3f%c", dd_vol_min(dd), flags ? '!' : ' ');
1234 log.writeStringFormatted(" load imb.: force %4.1f%%", dd_f_imbal(dd) * 100);
1236 if (dd->comm->cycl_n[ddCyclPME])
1238 log.writeStringFormatted(" pme mesh/force %5.3f", dd_pme_f_ratio(dd));
1240 log.ensureLineBreak();
1241 return stream.toString();
1244 //! Prints DD load balance report in mdrun verbose mode.
1245 static void dd_print_load_verbose(gmx_domdec_t* dd)
1247 if (isDlbOn(dd->comm))
1249 fprintf(stderr, "vol %4.2f%c ", dd_vol_min(dd), dd_load_flags(dd) ? '!' : ' ');
1253 fprintf(stderr, "imb F %2d%% ", gmx::roundToInt(dd_f_imbal(dd) * 100));
1255 if (dd->comm->cycl_n[ddCyclPME])
1257 fprintf(stderr, "pme/F %4.2f ", dd_pme_f_ratio(dd));
1261 //! Turns on dynamic load balancing if possible and needed.
1262 static void turn_on_dlb(const gmx::MDLogger& mdlog, gmx_domdec_t* dd, int64_t step)
1264 gmx_domdec_comm_t* comm = dd->comm;
1266 real cellsize_min = comm->cellsize_min[dd->dim[0]];
1267 for (int d = 1; d < dd->ndim; d++)
1269 cellsize_min = std::min(cellsize_min, comm->cellsize_min[dd->dim[d]]);
1272 /* Turn off DLB if we're too close to the cell size limit. */
1273 if (cellsize_min < comm->cellsize_limit * 1.05)
1276 .appendTextFormatted(
1277 "step %s Measured %.1f %% performance loss due to load imbalance, "
1278 "but the minimum cell size is smaller than 1.05 times the cell size limit. "
1279 "Will no longer try dynamic load balancing.",
1280 gmx::toString(step).c_str(),
1281 dd_force_imb_perf_loss(dd) * 100);
1283 comm->dlbState = DlbState::offForever;
1288 .appendTextFormatted(
1289 "step %s Turning on dynamic load balancing, because the performance loss due "
1290 "to load imbalance is %.1f %%.",
1291 gmx::toString(step).c_str(),
1292 dd_force_imb_perf_loss(dd) * 100);
1293 comm->dlbState = DlbState::onCanTurnOff;
1295 /* Store the non-DLB performance, so we can check if DLB actually
1296 * improves performance.
1298 GMX_RELEASE_ASSERT(comm->cycl_n[ddCyclStep] > 0,
1299 "When we turned on DLB, we should have measured cycles");
1300 comm->cyclesPerStepBeforeDLB = comm->cycl[ddCyclStep] / comm->cycl_n[ddCyclStep];
1304 /* We can set the required cell size info here,
1305 * so we do not need to communicate this.
1306 * The grid is completely uniform.
1308 for (int d = 0; d < dd->ndim; d++)
1310 RowMaster* rowMaster = comm->cellsizesWithDlb[d].rowMaster.get();
1314 comm->load[d].sum_m = comm->load[d].sum;
1316 int nc = dd->numCells[dd->dim[d]];
1317 for (int i = 0; i < nc; i++)
1319 rowMaster->cellFrac[i] = i / static_cast<real>(nc);
1322 rowMaster->bounds[i].cellFracLowerMax = i / static_cast<real>(nc);
1323 rowMaster->bounds[i].cellFracUpperMin = (i + 1) / static_cast<real>(nc);
1326 rowMaster->cellFrac[nc] = 1.0;
1331 //! Turns off dynamic load balancing (but leave it able to turn back on).
1332 static void turn_off_dlb(const gmx::MDLogger& mdlog, gmx_domdec_t* dd, int64_t step)
1336 "step " + gmx::toString(step)
1337 + " Turning off dynamic load balancing, because it is degrading performance.");
1338 dd->comm->dlbState = DlbState::offCanTurnOn;
1339 dd->comm->haveTurnedOffDlb = true;
1340 dd->comm->ddPartioningCountFirstDlbOff = dd->ddp_count;
1343 //! Turns off dynamic load balancing permanently.
1344 static void turn_off_dlb_forever(const gmx::MDLogger& mdlog, gmx_domdec_t* dd, int64_t step)
1346 GMX_RELEASE_ASSERT(dd->comm->dlbState == DlbState::offCanTurnOn,
1347 "Can only turn off DLB forever when it was in the can-turn-on state");
1350 "step " + gmx::toString(step)
1351 + " Will no longer try dynamic load balancing, as it degraded performance.");
1352 dd->comm->dlbState = DlbState::offForever;
1355 void set_dd_dlb_max_cutoff(t_commrec* cr, real cutoff)
1357 gmx_domdec_comm_t* comm;
1359 comm = cr->dd->comm;
1361 /* Turn on the DLB limiting (might have been on already) */
1362 comm->bPMELoadBalDLBLimits = TRUE;
1364 /* Change the cut-off limit */
1365 comm->PMELoadBal_max_cutoff = cutoff;
1370 "PME load balancing set a limit to the DLB staggering such that a %f cut-off will "
1371 "continue to fit\n",
1372 comm->PMELoadBal_max_cutoff);
1376 //! Merge atom buffers.
1377 static void merge_cg_buffers(int ncell,
1378 gmx_domdec_comm_dim_t* cd,
1381 gmx::ArrayRef<int> index_gl,
1383 gmx::ArrayRef<gmx::RVec> x,
1384 gmx::ArrayRef<const gmx::RVec> recv_vr,
1385 gmx::ArrayRef<cginfo_mb_t> cginfo_mb,
1386 gmx::ArrayRef<int> cginfo)
1388 gmx_domdec_ind_t *ind, *ind_p;
1389 int p, cell, c, cg, cg0, cg1, cg_gl;
1392 ind = &cd->ind[pulse];
1394 /* First correct the already stored data */
1395 shift = ind->nrecv[ncell];
1396 for (cell = ncell - 1; cell >= 0; cell--)
1398 shift -= ind->nrecv[cell];
1401 /* Move the cg's present from previous grid pulses */
1402 cg0 = ncg_cell[ncell + cell];
1403 cg1 = ncg_cell[ncell + cell + 1];
1404 for (cg = cg1 - 1; cg >= cg0; cg--)
1406 index_gl[cg + shift] = index_gl[cg];
1407 x[cg + shift] = x[cg];
1408 cginfo[cg + shift] = cginfo[cg];
1410 /* Correct the already stored send indices for the shift */
1411 for (p = 1; p <= pulse; p++)
1413 ind_p = &cd->ind[p];
1415 for (c = 0; c < cell; c++)
1417 cg0 += ind_p->nsend[c];
1419 cg1 = cg0 + ind_p->nsend[cell];
1420 for (cg = cg0; cg < cg1; cg++)
1422 ind_p->index[cg] += shift;
1428 /* Merge in the communicated buffers */
1431 for (cell = 0; cell < ncell; cell++)
1433 cg1 = ncg_cell[ncell + cell + 1] + shift;
1434 for (cg = 0; cg < ind->nrecv[cell]; cg++)
1436 /* Copy this atom from the buffer */
1437 index_gl[cg1] = recv_i[cg0];
1438 x[cg1] = recv_vr[cg0];
1439 /* Copy information */
1440 cg_gl = index_gl[cg1];
1441 cginfo[cg1] = ddcginfo(cginfo_mb, cg_gl);
1445 shift += ind->nrecv[cell];
1446 ncg_cell[ncell + cell + 1] = cg1;
1450 //! Makes a range partitioning for the atom groups wthin a cell
1451 static void make_cell2at_index(gmx_domdec_comm_dim_t* cd, int nzone, int atomGroupStart)
1453 /* Store the atom block boundaries for easy copying of communication buffers
1455 int g = atomGroupStart;
1456 for (int zone = 0; zone < nzone; zone++)
1458 for (gmx_domdec_ind_t& ind : cd->ind)
1460 ind.cell2at0[zone] = g;
1461 g += ind.nrecv[zone];
1462 ind.cell2at1[zone] = g;
1467 //! Returns whether a link is missing.
1468 static gmx_bool missing_link(const t_blocka& link, const int globalAtomIndex, const gmx_ga2la_t& ga2la)
1470 for (int i = link.index[globalAtomIndex]; i < link.index[globalAtomIndex + 1]; i++)
1472 if (!ga2la.findHome(link.a[i]))
1481 //! Domain corners for communication, a maximum of 4 i-zones see a j domain
1484 //! The corners for the non-bonded communication.
1486 //! Corner for rounding.
1488 //! Corners for rounding.
1490 //! Corners for bounded communication.
1492 //! Corner for rounding for bonded communication.
1496 //! Determine the corners of the domain(s) we are communicating with.
1497 static void set_dd_corners(const gmx_domdec_t* dd, int dim0, int dim1, int dim2, gmx_bool bDistMB, dd_corners_t* c)
1499 const gmx_domdec_comm_t* comm;
1500 const gmx_domdec_zones_t* zones;
1504 zones = &comm->zones;
1506 /* Keep the compiler happy */
1510 /* The first dimension is equal for all cells */
1511 c->c[0][0] = comm->cell_x0[dim0];
1514 c->bc[0] = c->c[0][0];
1519 /* This cell row is only seen from the first row */
1520 c->c[1][0] = comm->cell_x0[dim1];
1521 /* All rows can see this row */
1522 c->c[1][1] = comm->cell_x0[dim1];
1523 if (isDlbOn(dd->comm))
1525 c->c[1][1] = std::max(comm->cell_x0[dim1], comm->zone_d1[1].mch0);
1528 /* For the multi-body distance we need the maximum */
1529 c->bc[1] = std::max(comm->cell_x0[dim1], comm->zone_d1[1].p1_0);
1532 /* Set the upper-right corner for rounding */
1533 c->cr0 = comm->cell_x1[dim0];
1538 for (int j = 0; j < 4; j++)
1540 c->c[2][j] = comm->cell_x0[dim2];
1542 if (isDlbOn(dd->comm))
1544 /* Use the maximum of the i-cells that see a j-cell */
1545 for (const auto& iZone : zones->iZones)
1547 const int iZoneIndex = iZone.iZoneIndex;
1548 for (int jZone : iZone.jZoneRange)
1552 c->c[2][jZone - 4] = std::max(
1554 comm->zone_d2[zones->shift[iZoneIndex][dim0]][zones->shift[iZoneIndex][dim1]]
1561 /* For the multi-body distance we need the maximum */
1562 c->bc[2] = comm->cell_x0[dim2];
1563 for (int i = 0; i < 2; i++)
1565 for (int j = 0; j < 2; j++)
1567 c->bc[2] = std::max(c->bc[2], comm->zone_d2[i][j].p1_0);
1573 /* Set the upper-right corner for rounding */
1574 /* Cell (0,0,0) and cell (1,0,0) can see cell 4 (0,1,1)
1575 * Only cell (0,0,0) can see cell 7 (1,1,1)
1577 c->cr1[0] = comm->cell_x1[dim1];
1578 c->cr1[3] = comm->cell_x1[dim1];
1579 if (isDlbOn(dd->comm))
1581 c->cr1[0] = std::max(comm->cell_x1[dim1], comm->zone_d1[1].mch1);
1584 /* For the multi-body distance we need the maximum */
1585 c->bcr1 = std::max(comm->cell_x1[dim1], comm->zone_d1[1].p1_1);
1592 /*! \brief Add the atom groups we need to send in this pulse from this
1593 * zone to \p localAtomGroups and \p work. */
1594 static void get_zone_pulse_cgs(gmx_domdec_t* dd,
1599 gmx::ArrayRef<const int> globalAtomGroupIndices,
1608 bool distanceIsTriclinic,
1615 const dd_corners_t* c,
1616 const rvec sf2_round,
1617 gmx_bool bDistBonded,
1622 gmx::ArrayRef<const int> cginfo,
1623 std::vector<int>* localAtomGroups,
1624 dd_comm_setup_work_t* work)
1626 gmx_domdec_comm_t* comm;
1628 gmx_bool bDistMB_pulse;
1630 real r2, rb2, r, tric_sh;
1637 bScrew = (dd->unitCellInfo.haveScrewPBC && dim == XX);
1639 bDistMB_pulse = (bDistMB && bDistBonded);
1641 /* Unpack the work data */
1642 std::vector<int>& ibuf = work->atomGroupBuffer;
1643 std::vector<gmx::RVec>& vbuf = work->positionBuffer;
1647 for (cg = cg0; cg < cg1; cg++)
1651 if (!distanceIsTriclinic)
1653 /* Rectangular direction, easy */
1654 r = cg_cm[cg][dim] - c->c[dim_ind][zone];
1661 r = cg_cm[cg][dim] - c->bc[dim_ind];
1667 /* Rounding gives at most a 16% reduction
1668 * in communicated atoms
1670 if (dim_ind >= 1 && (zonei == 1 || zonei == 2))
1672 r = cg_cm[cg][dim0] - c->cr0;
1673 /* This is the first dimension, so always r >= 0 */
1680 if (dim_ind == 2 && (zonei == 2 || zonei == 3))
1682 r = cg_cm[cg][dim1] - c->cr1[zone];
1689 r = cg_cm[cg][dim1] - c->bcr1;
1699 /* Triclinic direction, more complicated */
1702 /* Rounding, conservative as the skew_fac multiplication
1703 * will slightly underestimate the distance.
1705 if (dim_ind >= 1 && (zonei == 1 || zonei == 2))
1707 rn[dim0] = cg_cm[cg][dim0] - c->cr0;
1708 for (i = dim0 + 1; i < DIM; i++)
1710 rn[dim0] -= cg_cm[cg][i] * v_0[i][dim0];
1712 r2 = rn[dim0] * rn[dim0] * sf2_round[dim0];
1715 rb[dim0] = rn[dim0];
1718 /* Take care that the cell planes along dim0 might not
1719 * be orthogonal to those along dim1 and dim2.
1721 for (i = 1; i <= dim_ind; i++)
1724 if (normal[dim0][dimd] > 0)
1726 rn[dimd] -= rn[dim0] * normal[dim0][dimd];
1729 rb[dimd] -= rb[dim0] * normal[dim0][dimd];
1734 if (dim_ind == 2 && (zonei == 2 || zonei == 3))
1736 GMX_ASSERT(dim1 >= 0 && dim1 < DIM, "Must have a valid dimension index");
1737 rn[dim1] += cg_cm[cg][dim1] - c->cr1[zone];
1739 for (i = dim1 + 1; i < DIM; i++)
1741 tric_sh -= cg_cm[cg][i] * v_1[i][dim1];
1743 rn[dim1] += tric_sh;
1746 r2 += rn[dim1] * rn[dim1] * sf2_round[dim1];
1747 /* Take care of coupling of the distances
1748 * to the planes along dim0 and dim1 through dim2.
1750 r2 -= rn[dim0] * rn[dim1] * skew_fac_01;
1751 /* Take care that the cell planes along dim1
1752 * might not be orthogonal to that along dim2.
1754 if (normal[dim1][dim2] > 0)
1756 rn[dim2] -= rn[dim1] * normal[dim1][dim2];
1761 rb[dim1] += cg_cm[cg][dim1] - c->bcr1 + tric_sh;
1764 rb2 += rb[dim1] * rb[dim1] * sf2_round[dim1];
1765 /* Take care of coupling of the distances
1766 * to the planes along dim0 and dim1 through dim2.
1768 rb2 -= rb[dim0] * rb[dim1] * skew_fac_01;
1769 /* Take care that the cell planes along dim1
1770 * might not be orthogonal to that along dim2.
1772 if (normal[dim1][dim2] > 0)
1774 rb[dim2] -= rb[dim1] * normal[dim1][dim2];
1779 /* The distance along the communication direction */
1780 rn[dim] += cg_cm[cg][dim] - c->c[dim_ind][zone];
1782 for (i = dim + 1; i < DIM; i++)
1784 tric_sh -= cg_cm[cg][i] * v_d[i][dim];
1789 r2 += rn[dim] * rn[dim] * skew_fac2_d;
1790 /* Take care of coupling of the distances
1791 * to the planes along dim0 and dim1 through dim2.
1793 if (dim_ind == 1 && zonei == 1)
1795 r2 -= rn[dim0] * rn[dim] * skew_fac_01;
1801 GMX_ASSERT(dim >= 0 && dim < DIM, "Must have a valid dimension index");
1802 rb[dim] += cg_cm[cg][dim] - c->bc[dim_ind] + tric_sh;
1805 rb2 += rb[dim] * rb[dim] * skew_fac2_d;
1806 /* Take care of coupling of the distances
1807 * to the planes along dim0 and dim1 through dim2.
1809 if (dim_ind == 1 && zonei == 1)
1811 rb2 -= rb[dim0] * rb[dim] * skew_fac_01;
1818 || (bDistBonded && ((bDistMB && rb2 < r_bcomm2) || (bDist2B && r2 < r_bcomm2))
1820 || (GET_CGINFO_BOND_INTER(cginfo[cg])
1821 && missing_link(*comm->bondedLinks, globalAtomGroupIndices[cg], *dd->ga2la)))))
1823 /* Store the local and global atom group indices and position */
1824 localAtomGroups->push_back(cg);
1825 ibuf.push_back(globalAtomGroupIndices[cg]);
1829 if (dd->ci[dim] == 0)
1831 /* Correct cg_cm for pbc */
1832 rvec_add(cg_cm[cg], box[dim], posPbc);
1835 posPbc[YY] = box[YY][YY] - posPbc[YY];
1836 posPbc[ZZ] = box[ZZ][ZZ] - posPbc[ZZ];
1841 copy_rvec(cg_cm[cg], posPbc);
1843 vbuf.emplace_back(posPbc[XX], posPbc[YY], posPbc[ZZ]);
1850 work->nsend_zone = nsend_z;
1854 static void clearCommSetupData(dd_comm_setup_work_t* work)
1856 work->localAtomGroupBuffer.clear();
1857 work->atomGroupBuffer.clear();
1858 work->positionBuffer.clear();
1860 work->nsend_zone = 0;
1863 //! Prepare DD communication.
1864 static void setup_dd_communication(gmx_domdec_t* dd, matrix box, gmx_ddbox_t* ddbox, t_forcerec* fr, t_state* state)
1866 int dim_ind, dim, dim0, dim1, dim2, dimd, nat_tot;
1867 int nzone, nzone_send, zone, zonei, cg0, cg1;
1869 int * zone_cg_range, pos_cg;
1870 gmx_domdec_comm_t* comm;
1871 gmx_domdec_zones_t* zones;
1872 gmx_domdec_comm_dim_t* cd;
1873 gmx_bool bBondComm, bDist2B, bDistMB, bDistBonded;
1874 dd_corners_t corners;
1875 rvec * normal, *v_d, *v_0 = nullptr, *v_1 = nullptr;
1876 real skew_fac2_d, skew_fac_01;
1881 fprintf(debug, "Setting up DD communication\n");
1886 if (comm->dth.empty())
1888 /* Initialize the thread data.
1889 * This can not be done in init_domain_decomposition,
1890 * as the numbers of threads is determined later.
1892 int numThreads = gmx_omp_nthreads_get(emntDomdec);
1893 comm->dth.resize(numThreads);
1896 bBondComm = comm->systemInfo.filterBondedCommunication;
1898 /* Do we need to determine extra distances for multi-body bondeds? */
1899 bDistMB = (comm->systemInfo.haveInterDomainMultiBodyBondeds && isDlbOn(dd->comm) && dd->ndim > 1);
1901 /* Do we need to determine extra distances for only two-body bondeds? */
1902 bDist2B = (bBondComm && !bDistMB);
1904 const real r_comm2 =
1905 gmx::square(domainToDomainIntoAtomToDomainCutoff(comm->systemInfo, comm->systemInfo.cutoff));
1906 const real r_bcomm2 =
1907 gmx::square(domainToDomainIntoAtomToDomainCutoff(comm->systemInfo, comm->cutoff_mbody));
1911 fprintf(debug, "bBondComm %s, r_bc %f\n", gmx::boolToString(bBondComm), std::sqrt(r_bcomm2));
1914 zones = &comm->zones;
1917 dim1 = (dd->ndim >= 2 ? dd->dim[1] : -1);
1918 dim2 = (dd->ndim >= 3 ? dd->dim[2] : -1);
1920 set_dd_corners(dd, dim0, dim1, dim2, bDistMB, &corners);
1922 /* Triclinic stuff */
1923 normal = ddbox->normal;
1927 v_0 = ddbox->v[dim0];
1928 if (ddbox->tric_dir[dim0] && ddbox->tric_dir[dim1])
1930 /* Determine the coupling coefficient for the distances
1931 * to the cell planes along dim0 and dim1 through dim2.
1932 * This is required for correct rounding.
1934 skew_fac_01 = ddbox->v[dim0][dim1 + 1][dim0] * ddbox->v[dim1][dim1 + 1][dim1];
1937 fprintf(debug, "\nskew_fac_01 %f\n", skew_fac_01);
1943 v_1 = ddbox->v[dim1];
1946 zone_cg_range = zones->cg_range.data();
1947 gmx::ArrayRef<cginfo_mb_t> cginfo_mb = fr->cginfo_mb;
1949 zone_cg_range[0] = 0;
1950 zone_cg_range[1] = dd->ncg_home;
1951 comm->zone_ncg1[0] = dd->ncg_home;
1952 pos_cg = dd->ncg_home;
1954 nat_tot = comm->atomRanges.numHomeAtoms();
1956 for (dim_ind = 0; dim_ind < dd->ndim; dim_ind++)
1958 dim = dd->dim[dim_ind];
1959 cd = &comm->cd[dim_ind];
1961 /* Check if we need to compute triclinic distances along this dim */
1962 bool distanceIsTriclinic = false;
1963 for (int i = 0; i <= dim_ind; i++)
1965 if (ddbox->tric_dir[dd->dim[i]])
1967 distanceIsTriclinic = true;
1971 if (dim >= ddbox->npbcdim && dd->ci[dim] == 0)
1973 /* No pbc in this dimension, the first node should not comm. */
1981 v_d = ddbox->v[dim];
1982 skew_fac2_d = gmx::square(ddbox->skew_fac[dim]);
1984 cd->receiveInPlace = true;
1985 for (int p = 0; p < cd->numPulses(); p++)
1987 /* Only atoms communicated in the first pulse are used
1988 * for multi-body bonded interactions or for bBondComm.
1990 bDistBonded = ((bDistMB || bDist2B) && p == 0);
1992 gmx_domdec_ind_t* ind = &cd->ind[p];
1994 /* Thread 0 writes in the global index array */
1996 clearCommSetupData(&comm->dth[0]);
1998 for (zone = 0; zone < nzone_send; zone++)
2000 if (dim_ind > 0 && distanceIsTriclinic)
2002 /* Determine slightly more optimized skew_fac's
2004 * This reduces the number of communicated atoms
2005 * by about 10% for 3D DD of rhombic dodecahedra.
2007 for (dimd = 0; dimd < dim; dimd++)
2009 sf2_round[dimd] = 1;
2010 if (ddbox->tric_dir[dimd])
2012 for (int i = dd->dim[dimd] + 1; i < DIM; i++)
2014 /* If we are shifted in dimension i
2015 * and the cell plane is tilted forward
2016 * in dimension i, skip this coupling.
2018 if (!(zones->shift[nzone + zone][i] && ddbox->v[dimd][i][dimd] >= 0))
2020 sf2_round[dimd] += gmx::square(ddbox->v[dimd][i][dimd]);
2023 sf2_round[dimd] = 1 / sf2_round[dimd];
2028 zonei = zone_perm[dim_ind][zone];
2031 /* Here we permutate the zones to obtain a convenient order
2032 * for neighbor searching
2034 cg0 = zone_cg_range[zonei];
2035 cg1 = zone_cg_range[zonei + 1];
2039 /* Look only at the cg's received in the previous grid pulse
2041 cg1 = zone_cg_range[nzone + zone + 1];
2042 cg0 = cg1 - cd->ind[p - 1].nrecv[zone];
2045 const int numThreads = gmx::ssize(comm->dth);
2046 #pragma omp parallel for num_threads(numThreads) schedule(static)
2047 for (int th = 0; th < numThreads; th++)
2051 dd_comm_setup_work_t& work = comm->dth[th];
2053 /* Retain data accumulated into buffers of thread 0 */
2056 clearCommSetupData(&work);
2059 int cg0_th = cg0 + ((cg1 - cg0) * th) / numThreads;
2060 int cg1_th = cg0 + ((cg1 - cg0) * (th + 1)) / numThreads;
2062 /* Get the cg's for this pulse in this zone */
2063 get_zone_pulse_cgs(dd,
2068 dd->globalAtomGroupIndices,
2077 distanceIsTriclinic,
2090 state->x.rvec_array(),
2092 th == 0 ? &ind->index : &work.localAtomGroupBuffer,
2095 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2098 std::vector<int>& atomGroups = comm->dth[0].atomGroupBuffer;
2099 std::vector<gmx::RVec>& positions = comm->dth[0].positionBuffer;
2100 ind->nsend[zone] = comm->dth[0].nsend_zone;
2101 /* Append data of threads>=1 to the communication buffers */
2102 for (int th = 1; th < numThreads; th++)
2104 const dd_comm_setup_work_t& dth = comm->dth[th];
2106 ind->index.insert(ind->index.end(),
2107 dth.localAtomGroupBuffer.begin(),
2108 dth.localAtomGroupBuffer.end());
2110 atomGroups.end(), dth.atomGroupBuffer.begin(), dth.atomGroupBuffer.end());
2112 positions.end(), dth.positionBuffer.begin(), dth.positionBuffer.end());
2113 comm->dth[0].nat += dth.nat;
2114 ind->nsend[zone] += dth.nsend_zone;
2117 /* Clear the counts in case we do not have pbc */
2118 for (zone = nzone_send; zone < nzone; zone++)
2120 ind->nsend[zone] = 0;
2122 ind->nsend[nzone] = ind->index.size();
2123 ind->nsend[nzone + 1] = comm->dth[0].nat;
2124 /* Communicate the number of cg's and atoms to receive */
2125 ddSendrecv(dd, dim_ind, dddirBackward, ind->nsend, nzone + 2, ind->nrecv, nzone + 2);
2129 /* We can receive in place if only the last zone is not empty */
2130 for (zone = 0; zone < nzone - 1; zone++)
2132 if (ind->nrecv[zone] > 0)
2134 cd->receiveInPlace = false;
2139 int receiveBufferSize = 0;
2140 if (!cd->receiveInPlace)
2142 receiveBufferSize = ind->nrecv[nzone];
2144 /* These buffer are actually only needed with in-place */
2145 DDBufferAccess<int> globalAtomGroupBuffer(comm->intBuffer, receiveBufferSize);
2146 DDBufferAccess<gmx::RVec> rvecBuffer(comm->rvecBuffer, receiveBufferSize);
2148 dd_comm_setup_work_t& work = comm->dth[0];
2150 /* Make space for the global cg indices */
2151 int numAtomGroupsNew = pos_cg + ind->nrecv[nzone];
2152 dd->globalAtomGroupIndices.resize(numAtomGroupsNew);
2153 /* Communicate the global cg indices */
2154 gmx::ArrayRef<int> integerBufferRef;
2155 if (cd->receiveInPlace)
2157 integerBufferRef = gmx::arrayRefFromArray(
2158 dd->globalAtomGroupIndices.data() + pos_cg, ind->nrecv[nzone]);
2162 integerBufferRef = globalAtomGroupBuffer.buffer;
2164 ddSendrecv<int>(dd, dim_ind, dddirBackward, work.atomGroupBuffer, integerBufferRef);
2166 /* Make space for cg_cm */
2167 dd_resize_atominfo_and_state(fr, state, pos_cg + ind->nrecv[nzone]);
2169 /* Communicate the coordinates */
2170 gmx::ArrayRef<gmx::RVec> rvecBufferRef;
2171 if (cd->receiveInPlace)
2173 rvecBufferRef = gmx::makeArrayRef(state->x).subArray(pos_cg, ind->nrecv[nzone]);
2177 rvecBufferRef = rvecBuffer.buffer;
2179 ddSendrecv<gmx::RVec>(dd, dim_ind, dddirBackward, work.positionBuffer, rvecBufferRef);
2181 /* Make the charge group index */
2182 if (cd->receiveInPlace)
2184 zone = (p == 0 ? 0 : nzone - 1);
2185 while (zone < nzone)
2187 for (int i = 0; i < ind->nrecv[zone]; i++)
2189 int globalAtomIndex = dd->globalAtomGroupIndices[pos_cg];
2190 fr->cginfo[pos_cg] = ddcginfo(cginfo_mb, globalAtomIndex);
2195 comm->zone_ncg1[nzone + zone] = ind->nrecv[zone];
2198 zone_cg_range[nzone + zone] = pos_cg;
2203 /* This part of the code is never executed with bBondComm. */
2204 merge_cg_buffers(nzone,
2208 dd->globalAtomGroupIndices,
2209 integerBufferRef.data(),
2214 pos_cg += ind->nrecv[nzone];
2216 nat_tot += ind->nrecv[nzone + 1];
2218 if (!cd->receiveInPlace)
2220 /* Store the atom block for easy copying of communication buffers */
2221 make_cell2at_index(cd, nzone, zone_cg_range[nzone]);
2226 comm->atomRanges.setEnd(DDAtomRanges::Type::Zones, nat_tot);
2230 /* We don't need to update cginfo, since that was alrady done above.
2231 * So we pass NULL for the forcerec.
2233 dd_set_cginfo(dd->globalAtomGroupIndices, dd->ncg_home, dd->globalAtomGroupIndices.size(), nullptr);
2238 fprintf(debug, "Finished setting up DD communication, zones:");
2239 for (c = 0; c < zones->n; c++)
2241 fprintf(debug, " %d", zones->cg_range[c + 1] - zones->cg_range[c]);
2243 fprintf(debug, "\n");
2247 //! Set boundaries for the charge group range.
2248 static void set_cg_boundaries(gmx_domdec_zones_t* zones)
2250 for (auto& iZone : zones->iZones)
2252 iZone.iAtomRange = gmx::Range<int>(0, zones->cg_range[iZone.iZoneIndex + 1]);
2253 iZone.jAtomRange = gmx::Range<int>(zones->cg_range[iZone.jZoneRange.begin()],
2254 zones->cg_range[iZone.jZoneRange.end()]);
2258 /*! \brief Set zone dimensions for zones \p zone_start to \p zone_end-1
2260 * Also sets the atom density for the home zone when \p zone_start=0.
2261 * For this \p numMovedChargeGroupsInHomeZone needs to be passed to tell
2262 * how many charge groups will move but are still part of the current range.
2263 * \todo When converting domdec to use proper classes, all these variables
2264 * should be private and a method should return the correct count
2265 * depending on an internal state.
2267 * \param[in,out] dd The domain decomposition struct
2268 * \param[in] box The box
2269 * \param[in] ddbox The domain decomposition box struct
2270 * \param[in] zone_start The start of the zone range to set sizes for
2271 * \param[in] zone_end The end of the zone range to set sizes for
2272 * \param[in] numMovedChargeGroupsInHomeZone The number of charge groups in the home zone that should moved but are still present in dd->comm->zones.cg_range
2274 static void set_zones_size(gmx_domdec_t* dd,
2276 const gmx_ddbox_t* ddbox,
2279 int numMovedChargeGroupsInHomeZone)
2281 gmx_domdec_comm_t* comm;
2282 gmx_domdec_zones_t* zones;
2291 zones = &comm->zones;
2293 /* Do we need to determine extra distances for multi-body bondeds? */
2294 bDistMB = (comm->systemInfo.haveInterDomainMultiBodyBondeds && isDlbOn(dd->comm) && dd->ndim > 1);
2296 for (z = zone_start; z < zone_end; z++)
2298 /* Copy cell limits to zone limits.
2299 * Valid for non-DD dims and non-shifted dims.
2301 copy_rvec(comm->cell_x0, zones->size[z].x0);
2302 copy_rvec(comm->cell_x1, zones->size[z].x1);
2305 for (d = 0; d < dd->ndim; d++)
2309 for (z = 0; z < zones->n; z++)
2311 /* With a staggered grid we have different sizes
2312 * for non-shifted dimensions.
2314 if (isDlbOn(dd->comm) && zones->shift[z][dim] == 0)
2318 zones->size[z].x0[dim] = comm->zone_d1[zones->shift[z][dd->dim[d - 1]]].min0;
2319 zones->size[z].x1[dim] = comm->zone_d1[zones->shift[z][dd->dim[d - 1]]].max1;
2323 zones->size[z].x0[dim] =
2324 comm->zone_d2[zones->shift[z][dd->dim[d - 2]]][zones->shift[z][dd->dim[d - 1]]]
2326 zones->size[z].x1[dim] =
2327 comm->zone_d2[zones->shift[z][dd->dim[d - 2]]][zones->shift[z][dd->dim[d - 1]]]
2333 rcs = comm->systemInfo.cutoff;
2334 rcmbs = comm->cutoff_mbody;
2335 if (ddbox->tric_dir[dim])
2337 rcs /= ddbox->skew_fac[dim];
2338 rcmbs /= ddbox->skew_fac[dim];
2341 /* Set the lower limit for the shifted zone dimensions */
2342 for (z = zone_start; z < zone_end; z++)
2344 if (zones->shift[z][dim] > 0)
2347 if (!isDlbOn(dd->comm) || d == 0)
2349 zones->size[z].x0[dim] = comm->cell_x1[dim];
2350 zones->size[z].x1[dim] = comm->cell_x1[dim] + rcs;
2354 /* Here we take the lower limit of the zone from
2355 * the lowest domain of the zone below.
2359 zones->size[z].x0[dim] = comm->zone_d1[zones->shift[z][dd->dim[d - 1]]].min1;
2365 zones->size[z].x0[dim] = zones->size[zone_perm[2][z - 4]].x0[dim];
2369 zones->size[z].x0[dim] =
2370 comm->zone_d2[zones->shift[z][dd->dim[d - 2]]][zones->shift[z][dd->dim[d - 1]]]
2374 /* A temporary limit, is updated below */
2375 zones->size[z].x1[dim] = zones->size[z].x0[dim];
2379 for (size_t zi = 0; zi < zones->iZones.size(); zi++)
2381 if (zones->shift[zi][dim] == 0)
2383 /* This takes the whole zone into account.
2384 * With multiple pulses this will lead
2385 * to a larger zone then strictly necessary.
2387 zones->size[z].x1[dim] = std::max(zones->size[z].x1[dim],
2388 zones->size[zi].x1[dim] + rcmbs);
2396 /* Loop over the i-zones to set the upper limit of each
2399 for (const auto& iZone : zones->iZones)
2401 const int zi = iZone.iZoneIndex;
2402 if (zones->shift[zi][dim] == 0)
2404 /* We should only use zones up to zone_end */
2405 const auto& jZoneRangeFull = iZone.jZoneRange;
2406 if (zone_end <= *jZoneRangeFull.begin())
2410 const gmx::Range<int> jZoneRange(*jZoneRangeFull.begin(),
2411 std::min(*jZoneRangeFull.end(), zone_end));
2412 for (int jZone : jZoneRange)
2414 if (zones->shift[jZone][dim] > 0)
2416 zones->size[jZone].x1[dim] =
2417 std::max(zones->size[jZone].x1[dim], zones->size[zi].x1[dim] + rcs);
2424 for (z = zone_start; z < zone_end; z++)
2426 /* Initialization only required to keep the compiler happy */
2427 rvec corner_min = { 0, 0, 0 }, corner_max = { 0, 0, 0 }, corner;
2430 /* To determine the bounding box for a zone we need to find
2431 * the extreme corners of 4, 2 or 1 corners.
2433 nc = 1 << (ddbox->nboundeddim - 1);
2435 for (c = 0; c < nc; c++)
2437 /* Set up a zone corner at x=0, ignoring trilinic couplings */
2441 corner[YY] = zones->size[z].x0[YY];
2445 corner[YY] = zones->size[z].x1[YY];
2449 corner[ZZ] = zones->size[z].x0[ZZ];
2453 corner[ZZ] = zones->size[z].x1[ZZ];
2455 if (dd->ndim == 1 && dd->dim[0] < ZZ && ZZ < dd->unitCellInfo.npbcdim
2456 && box[ZZ][1 - dd->dim[0]] != 0)
2458 /* With 1D domain decomposition the cg's are not in
2459 * the triclinic box, but triclinic x-y and rectangular y/x-z.
2460 * Shift the corner of the z-vector back to along the box
2461 * vector of dimension d, so it will later end up at 0 along d.
2462 * This can affect the location of this corner along dd->dim[0]
2463 * through the matrix operation below if box[d][dd->dim[0]]!=0.
2465 int d = 1 - dd->dim[0];
2467 corner[d] -= corner[ZZ] * box[ZZ][d] / box[ZZ][ZZ];
2469 /* Apply the triclinic couplings */
2470 for (i = YY; i < ddbox->npbcdim && i < DIM; i++)
2472 for (j = XX; j < i; j++)
2474 corner[j] += corner[i] * box[i][j] / box[i][i];
2479 copy_rvec(corner, corner_min);
2480 copy_rvec(corner, corner_max);
2484 for (i = 0; i < DIM; i++)
2486 corner_min[i] = std::min(corner_min[i], corner[i]);
2487 corner_max[i] = std::max(corner_max[i], corner[i]);
2491 /* Copy the extreme cornes without offset along x */
2492 for (i = 0; i < DIM; i++)
2494 zones->size[z].bb_x0[i] = corner_min[i];
2495 zones->size[z].bb_x1[i] = corner_max[i];
2497 /* Add the offset along x */
2498 zones->size[z].bb_x0[XX] += zones->size[z].x0[XX];
2499 zones->size[z].bb_x1[XX] += zones->size[z].x1[XX];
2502 if (zone_start == 0)
2505 for (dim = 0; dim < DIM; dim++)
2507 vol *= zones->size[0].x1[dim] - zones->size[0].x0[dim];
2510 (zones->cg_range[1] - zones->cg_range[0] - numMovedChargeGroupsInHomeZone) / vol;
2515 for (z = zone_start; z < zone_end; z++)
2518 "zone %d %6.3f - %6.3f %6.3f - %6.3f %6.3f - %6.3f\n",
2520 zones->size[z].x0[XX],
2521 zones->size[z].x1[XX],
2522 zones->size[z].x0[YY],
2523 zones->size[z].x1[YY],
2524 zones->size[z].x0[ZZ],
2525 zones->size[z].x1[ZZ]);
2527 "zone %d bb %6.3f - %6.3f %6.3f - %6.3f %6.3f - %6.3f\n",
2529 zones->size[z].bb_x0[XX],
2530 zones->size[z].bb_x1[XX],
2531 zones->size[z].bb_x0[YY],
2532 zones->size[z].bb_x1[YY],
2533 zones->size[z].bb_x0[ZZ],
2534 zones->size[z].bb_x1[ZZ]);
2539 /*! \brief Order data in \p dataToSort according to \p sort
2541 * Note: both buffers should have at least \p sort.size() elements.
2543 template<typename T>
2544 static void orderVector(gmx::ArrayRef<const gmx_cgsort_t> sort,
2545 gmx::ArrayRef<T> dataToSort,
2546 gmx::ArrayRef<T> sortBuffer)
2548 GMX_ASSERT(dataToSort.size() >= sort.size(), "The vector needs to be sufficiently large");
2549 GMX_ASSERT(sortBuffer.size() >= sort.size(),
2550 "The sorting buffer needs to be sufficiently large");
2552 /* Order the data into the temporary buffer */
2554 for (const gmx_cgsort_t& entry : sort)
2556 sortBuffer[i++] = dataToSort[entry.ind];
2559 /* Copy back to the original array */
2560 std::copy(sortBuffer.begin(), sortBuffer.begin() + sort.size(), dataToSort.begin());
2563 /*! \brief Order data in \p dataToSort according to \p sort
2565 * Note: \p vectorToSort should have at least \p sort.size() elements,
2566 * \p workVector is resized when it is too small.
2568 template<typename T>
2569 static void orderVector(gmx::ArrayRef<const gmx_cgsort_t> sort,
2570 gmx::ArrayRef<T> vectorToSort,
2571 std::vector<T>* workVector)
2573 if (gmx::index(workVector->size()) < sort.ssize())
2575 workVector->resize(sort.size());
2577 orderVector<T>(sort, vectorToSort, *workVector);
2580 //! Returns the sorting order for atoms based on the nbnxn grid order in sort
2581 static void dd_sort_order_nbnxn(const t_forcerec* fr, std::vector<gmx_cgsort_t>* sort)
2583 gmx::ArrayRef<const int> atomOrder = fr->nbv->getLocalAtomOrder();
2585 /* Using push_back() instead of this resize results in much slower code */
2586 sort->resize(atomOrder.size());
2587 gmx::ArrayRef<gmx_cgsort_t> buffer = *sort;
2588 size_t numSorted = 0;
2589 for (int i : atomOrder)
2593 /* The values of nsc and ind_gl are not used in this case */
2594 buffer[numSorted++].ind = i;
2597 sort->resize(numSorted);
2600 //! Returns the sorting state for DD.
2601 static void dd_sort_state(gmx_domdec_t* dd, t_forcerec* fr, t_state* state)
2603 gmx_domdec_sort_t* sort = dd->comm->sort.get();
2605 dd_sort_order_nbnxn(fr, &sort->sorted);
2607 /* We alloc with the old size, since cgindex is still old */
2608 DDBufferAccess<gmx::RVec> rvecBuffer(dd->comm->rvecBuffer, dd->ncg_home);
2610 /* Set the new home atom/charge group count */
2611 dd->ncg_home = sort->sorted.size();
2614 fprintf(debug, "Set the new home atom count to %d\n", dd->ncg_home);
2617 /* Reorder the state */
2618 gmx::ArrayRef<const gmx_cgsort_t> cgsort = sort->sorted;
2619 GMX_RELEASE_ASSERT(cgsort.ssize() == dd->ncg_home, "We should sort all the home atom groups");
2621 if (state->flags & enumValueToBitMask(StateEntry::X))
2623 orderVector(cgsort, makeArrayRef(state->x), rvecBuffer.buffer);
2625 if (state->flags & enumValueToBitMask(StateEntry::V))
2627 orderVector(cgsort, makeArrayRef(state->v), rvecBuffer.buffer);
2629 if (state->flags & enumValueToBitMask(StateEntry::Cgp))
2631 orderVector(cgsort, makeArrayRef(state->cg_p), rvecBuffer.buffer);
2634 /* Reorder the global cg index */
2635 orderVector<int>(cgsort, dd->globalAtomGroupIndices, &sort->intBuffer);
2636 /* Reorder the cginfo */
2637 orderVector<int>(cgsort, fr->cginfo, &sort->intBuffer);
2638 /* Set the home atom number */
2639 dd->comm->atomRanges.setEnd(DDAtomRanges::Type::Home, dd->ncg_home);
2641 /* The atoms are now exactly in grid order, update the grid order */
2642 fr->nbv->setLocalAtomOrder();
2645 //! Accumulates load statistics.
2646 static void add_dd_statistics(gmx_domdec_t* dd)
2648 gmx_domdec_comm_t* comm = dd->comm;
2650 for (int i = 0; i < static_cast<int>(DDAtomRanges::Type::Number); i++)
2652 auto range = static_cast<DDAtomRanges::Type>(i);
2653 comm->sum_nat[i] += comm->atomRanges.end(range) - comm->atomRanges.start(range);
2658 void reset_dd_statistics_counters(gmx_domdec_t* dd)
2660 gmx_domdec_comm_t* comm = dd->comm;
2662 /* Reset all the statistics and counters for total run counting */
2663 for (int i = 0; i < static_cast<int>(DDAtomRanges::Type::Number); i++)
2665 comm->sum_nat[i] = 0;
2669 comm->load_step = 0;
2672 clear_ivec(comm->load_lim);
2677 void print_dd_statistics(const t_commrec* cr, const t_inputrec& inputrec, FILE* fplog)
2679 gmx_domdec_comm_t* comm = cr->dd->comm;
2681 const int numRanges = static_cast<int>(DDAtomRanges::Type::Number);
2682 gmx_sumd(numRanges, comm->sum_nat, cr);
2684 if (fplog == nullptr)
2689 fprintf(fplog, "\n D O M A I N D E C O M P O S I T I O N S T A T I S T I C S\n\n");
2691 for (int i = static_cast<int>(DDAtomRanges::Type::Zones); i < numRanges; i++)
2693 auto range = static_cast<DDAtomRanges::Type>(i);
2694 double av = comm->sum_nat[i] / comm->ndecomp;
2697 case DDAtomRanges::Type::Zones:
2698 fprintf(fplog, " av. #atoms communicated per step for force: %d x %.1f\n", 2, av);
2700 case DDAtomRanges::Type::Vsites:
2701 if (cr->dd->vsite_comm)
2704 " av. #atoms communicated per step for vsites: %d x %.1f\n",
2705 (EEL_PME(inputrec.coulombtype)
2706 || inputrec.coulombtype == CoulombInteractionType::Ewald)
2712 case DDAtomRanges::Type::Constraints:
2713 if (cr->dd->constraint_comm)
2716 " av. #atoms communicated per step for LINCS: %d x %.1f\n",
2717 1 + inputrec.nLincsIter,
2721 default: gmx_incons(" Unknown type for DD statistics");
2724 fprintf(fplog, "\n");
2726 if (comm->ddSettings.recordLoad && EI_DYNAMICS(inputrec.eI))
2728 print_dd_load_av(fplog, cr->dd);
2732 //!\brief TODO Remove fplog when group scheme and charge groups are gone
2733 void dd_partition_system(FILE* fplog,
2734 const gmx::MDLogger& mdlog,
2736 const t_commrec* cr,
2737 gmx_bool bMasterState,
2739 t_state* state_global,
2740 const gmx_mtop_t& top_global,
2741 const t_inputrec& inputrec,
2742 gmx::ImdSession* imdSession,
2744 t_state* state_local,
2745 gmx::ForceBuffers* f,
2746 gmx::MDAtoms* mdAtoms,
2747 gmx_localtop_t* top_local,
2749 gmx::VirtualSitesHandler* vsite,
2750 gmx::Constraints* constr,
2752 gmx_wallcycle* wcycle,
2756 gmx_domdec_comm_t* comm;
2757 gmx_ddbox_t ddbox = { 0 };
2758 int64_t step_pcoupl;
2759 rvec cell_ns_x0, cell_ns_x1;
2760 int ncgindex_set, ncg_moved, nat_f_novirsum;
2761 gmx_bool bBoxChanged, bNStGlobalComm, bDoDLB, bCheckWhetherToTurnDlbOn, bLogLoad;
2767 wallcycle_start(wcycle, ewcDOMDEC);
2772 // TODO if the update code becomes accessible here, use
2773 // upd->deform for this logic.
2774 bBoxChanged = (bMasterState || inputrecDeform(&inputrec));
2775 if (inputrec.epc != PressureCoupling::No)
2777 /* With nstpcouple > 1 pressure coupling happens.
2778 * one step after calculating the pressure.
2779 * Box scaling happens at the end of the MD step,
2780 * after the DD partitioning.
2781 * We therefore have to do DLB in the first partitioning
2782 * after an MD step where P-coupling occurred.
2783 * We need to determine the last step in which p-coupling occurred.
2784 * MRS -- need to validate this for vv?
2786 int n = inputrec.nstpcouple;
2789 step_pcoupl = step - 1;
2793 step_pcoupl = ((step - 1) / n) * n + 1;
2795 if (step_pcoupl >= comm->partition_step)
2801 bNStGlobalComm = (step % nstglobalcomm == 0);
2809 /* Should we do dynamic load balacing this step?
2810 * Since it requires (possibly expensive) global communication,
2811 * we might want to do DLB less frequently.
2813 if (bBoxChanged || inputrec.epc != PressureCoupling::No)
2815 bDoDLB = bBoxChanged;
2819 bDoDLB = bNStGlobalComm;
2823 /* Check if we have recorded loads on the nodes */
2824 if (comm->ddSettings.recordLoad && dd_load_count(comm) > 0)
2826 bCheckWhetherToTurnDlbOn = dd_dlb_get_should_check_whether_to_turn_dlb_on(dd);
2828 /* Print load every nstlog, first and last step to the log file */
2829 bLogLoad = ((inputrec.nstlog > 0 && step % inputrec.nstlog == 0) || comm->n_load_collect == 0
2830 || (inputrec.nsteps >= 0
2831 && (step + inputrec.nstlist > inputrec.init_step + inputrec.nsteps)));
2833 /* Avoid extra communication due to verbose screen output
2834 * when nstglobalcomm is set.
2836 if (bDoDLB || bLogLoad || bCheckWhetherToTurnDlbOn
2837 || (bVerbose && (inputrec.nstlist == 0 || nstglobalcomm <= inputrec.nstlist)))
2839 get_load_distribution(dd, wcycle);
2844 GMX_LOG(mdlog.info).asParagraph().appendText(dd_print_load(dd, step - 1));
2848 dd_print_load_verbose(dd);
2851 comm->n_load_collect++;
2857 /* Add the measured cycles to the running average */
2858 const float averageFactor = 0.1F;
2859 comm->cyclesPerStepDlbExpAverage =
2860 (1 - averageFactor) * comm->cyclesPerStepDlbExpAverage
2861 + averageFactor * comm->cycl[ddCyclStep] / comm->cycl_n[ddCyclStep];
2863 if (comm->dlbState == DlbState::onCanTurnOff
2864 && dd->comm->n_load_have % c_checkTurnDlbOffInterval == c_checkTurnDlbOffInterval - 1)
2866 gmx_bool turnOffDlb;
2869 /* If the running averaged cycles with DLB are more
2870 * than before we turned on DLB, turn off DLB.
2871 * We will again run and check the cycles without DLB
2872 * and we can then decide if to turn off DLB forever.
2874 turnOffDlb = (comm->cyclesPerStepDlbExpAverage > comm->cyclesPerStepBeforeDLB);
2876 dd_bcast(dd, sizeof(turnOffDlb), &turnOffDlb);
2879 /* To turn off DLB, we need to redistribute the atoms */
2880 dd_collect_state(dd, state_local, state_global);
2881 bMasterState = TRUE;
2882 turn_off_dlb(mdlog, dd, step);
2886 else if (bCheckWhetherToTurnDlbOn)
2888 gmx_bool turnOffDlbForever = FALSE;
2889 gmx_bool turnOnDlb = FALSE;
2891 /* Since the timings are node dependent, the master decides */
2894 /* If we recently turned off DLB, we want to check if
2895 * performance is better without DLB. We want to do this
2896 * ASAP to minimize the chance that external factors
2897 * slowed down the DLB step are gone here and we
2898 * incorrectly conclude that DLB was causing the slowdown.
2899 * So we measure one nstlist block, no running average.
2901 if (comm->haveTurnedOffDlb
2902 && comm->cycl[ddCyclStep] / comm->cycl_n[ddCyclStep] < comm->cyclesPerStepDlbExpAverage)
2904 /* After turning off DLB we ran nstlist steps in fewer
2905 * cycles than with DLB. This likely means that DLB
2906 * in not benefical, but this could be due to a one
2907 * time unlucky fluctuation, so we require two such
2908 * observations in close succession to turn off DLB
2911 if (comm->dlbSlowerPartitioningCount > 0
2912 && dd->ddp_count < comm->dlbSlowerPartitioningCount + 10 * c_checkTurnDlbOnInterval)
2914 turnOffDlbForever = TRUE;
2916 comm->haveTurnedOffDlb = false;
2917 /* Register when we last measured DLB slowdown */
2918 comm->dlbSlowerPartitioningCount = dd->ddp_count;
2922 /* Here we check if the max PME rank load is more than 0.98
2923 * the max PP force load. If so, PP DLB will not help,
2924 * since we are (almost) limited by PME. Furthermore,
2925 * DLB will cause a significant extra x/f redistribution
2926 * cost on the PME ranks, which will then surely result
2927 * in lower total performance.
2929 if (comm->ddRankSetup.usePmeOnlyRanks && dd_pme_f_ratio(dd) > 1 - DD_PERF_LOSS_DLB_ON)
2935 turnOnDlb = (dd_force_imb_perf_loss(dd) >= DD_PERF_LOSS_DLB_ON);
2941 gmx_bool turnOffDlbForever;
2943 } bools{ turnOffDlbForever, turnOnDlb };
2944 dd_bcast(dd, sizeof(bools), &bools);
2945 if (bools.turnOffDlbForever)
2947 turn_off_dlb_forever(mdlog, dd, step);
2949 else if (bools.turnOnDlb)
2951 turn_on_dlb(mdlog, dd, step);
2956 comm->n_load_have++;
2962 /* Clear the old state */
2963 clearDDStateIndices(dd, false);
2966 auto xGlobal = positionsFromStatePointer(state_global);
2968 set_ddbox(*dd, true, DDMASTER(dd) ? state_global->box : nullptr, true, xGlobal, &ddbox);
2970 distributeState(mdlog, dd, top_global, state_global, ddbox, state_local);
2972 /* Ensure that we have space for the new distribution */
2973 dd_resize_atominfo_and_state(fr, state_local, dd->ncg_home);
2975 inc_nrnb(nrnb, eNR_CGCM, comm->atomRanges.numHomeAtoms());
2977 dd_set_cginfo(dd->globalAtomGroupIndices, 0, dd->ncg_home, fr);
2979 else if (state_local->ddp_count != dd->ddp_count)
2981 if (state_local->ddp_count > dd->ddp_count)
2984 "Internal inconsistency state_local->ddp_count (%d) > dd->ddp_count (%" PRId64
2986 state_local->ddp_count,
2990 if (state_local->ddp_count_cg_gl != state_local->ddp_count)
2993 "Internal inconsistency state_local->ddp_count_cg_gl (%d) != "
2994 "state_local->ddp_count (%d)",
2995 state_local->ddp_count_cg_gl,
2996 state_local->ddp_count);
2999 /* Clear the old state */
3000 clearDDStateIndices(dd, false);
3002 /* Restore the atom group indices from state_local */
3003 restoreAtomGroups(dd, state_local);
3004 make_dd_indices(dd, 0);
3005 ncgindex_set = dd->ncg_home;
3007 inc_nrnb(nrnb, eNR_CGCM, comm->atomRanges.numHomeAtoms());
3009 dd_set_cginfo(dd->globalAtomGroupIndices, 0, dd->ncg_home, fr);
3011 set_ddbox(*dd, bMasterState, state_local->box, true, state_local->x, &ddbox);
3013 bRedist = isDlbOn(comm);
3017 /* We have the full state, only redistribute the cgs */
3019 /* Clear the non-home indices */
3020 clearDDStateIndices(dd, true);
3023 /* To avoid global communication, we do not recompute the extent
3024 * of the system for dims without pbc. Therefore we need to copy
3025 * the previously computed values when we do not communicate.
3027 if (!bNStGlobalComm)
3029 copy_rvec(comm->box0, ddbox.box0);
3030 copy_rvec(comm->box_size, ddbox.box_size);
3032 set_ddbox(*dd, bMasterState, state_local->box, bNStGlobalComm, state_local->x, &ddbox);
3037 /* Copy needed for dim's without pbc when avoiding communication */
3038 copy_rvec(ddbox.box0, comm->box0);
3039 copy_rvec(ddbox.box_size, comm->box_size);
3041 set_dd_cell_sizes(dd, &ddbox, dd->unitCellInfo.ddBoxIsDynamic, bMasterState, bDoDLB, step, wcycle);
3043 if (comm->ddSettings.nstDDDumpGrid > 0 && step % comm->ddSettings.nstDDDumpGrid == 0)
3045 write_dd_grid_pdb("dd_grid", step, dd, state_local->box, &ddbox);
3048 if (comm->systemInfo.useUpdateGroups)
3050 comm->updateGroupsCog->addCogs(
3051 gmx::arrayRefFromArray(dd->globalAtomGroupIndices.data(), dd->ncg_home), state_local->x);
3054 /* Check if we should sort the charge groups */
3055 const bool bSortCG = (bMasterState || bRedist);
3057 /* When repartitioning we mark atom groups that will move to neighboring
3058 * DD cells, but we do not move them right away for performance reasons.
3059 * Thus we need to keep track of how many charge groups will move for
3060 * obtaining correct local charge group / atom counts.
3065 wallcycle_sub_start(wcycle, ewcsDD_REDIST);
3067 ncgindex_set = dd->ncg_home;
3068 dd_redistribute_cg(fplog, step, dd, ddbox.tric_dir, state_local, fr, nrnb, &ncg_moved);
3070 GMX_RELEASE_ASSERT(bSortCG, "Sorting is required after redistribution");
3072 if (comm->systemInfo.useUpdateGroups)
3074 comm->updateGroupsCog->addCogs(
3075 gmx::arrayRefFromArray(dd->globalAtomGroupIndices.data(), dd->ncg_home),
3079 wallcycle_sub_stop(wcycle, ewcsDD_REDIST);
3082 // TODO: Integrate this code in the nbnxm module
3083 get_nsgrid_boundaries(ddbox.nboundeddim,
3090 as_rvec_array(state_local->x.data()),
3097 comm_dd_ns_cell_sizes(dd, &ddbox, cell_ns_x0, cell_ns_x1, step);
3102 wallcycle_sub_start(wcycle, ewcsDD_GRID);
3104 /* Sort the state on charge group position.
3105 * This enables exact restarts from this step.
3106 * It also improves performance by about 15% with larger numbers
3107 * of atoms per node.
3110 /* Fill the ns grid with the home cell,
3111 * so we can sort with the indices.
3113 set_zones_ncg_home(dd);
3115 set_zones_size(dd, state_local->box, &ddbox, 0, 1, ncg_moved);
3117 nbnxn_put_on_grid(fr->nbv.get(),
3120 comm->zones.size[0].bb_x0,
3121 comm->zones.size[0].bb_x1,
3122 comm->updateGroupsCog.get(),
3123 { 0, dd->ncg_home },
3124 comm->zones.dens_zone0,
3128 bRedist ? comm->movedBuffer.data() : nullptr);
3132 fprintf(debug, "Step %s, sorting the %d home charge groups\n", gmx_step_str(step, sbuf), dd->ncg_home);
3134 dd_sort_state(dd, fr, state_local);
3136 /* After sorting and compacting we set the correct size */
3137 state_change_natoms(state_local, comm->atomRanges.numHomeAtoms());
3139 /* Rebuild all the indices */
3143 wallcycle_sub_stop(wcycle, ewcsDD_GRID);
3147 /* With the group scheme the sorting array is part of the DD state,
3148 * but it just got out of sync, so mark as invalid by emptying it.
3150 if (inputrec.cutoff_scheme == CutoffScheme::Group)
3152 comm->sort->sorted.clear();
3156 if (comm->systemInfo.useUpdateGroups)
3158 /* The update groups cog's are invalid after sorting
3159 * and need to be cleared before the next partitioning anyhow.
3161 comm->updateGroupsCog->clear();
3164 wallcycle_sub_start(wcycle, ewcsDD_SETUPCOMM);
3166 /* Set the induces for the home atoms */
3167 set_zones_ncg_home(dd);
3168 make_dd_indices(dd, ncgindex_set);
3170 /* Setup up the communication and communicate the coordinates */
3171 setup_dd_communication(dd, state_local->box, &ddbox, fr, state_local);
3173 /* Set the indices for the halo atoms */
3174 make_dd_indices(dd, dd->ncg_home);
3176 /* Set the charge group boundaries for neighbor searching */
3177 set_cg_boundaries(&comm->zones);
3179 /* When bSortCG=true, we have already set the size for zone 0 */
3180 set_zones_size(dd, state_local->box, &ddbox, bSortCG ? 1 : 0, comm->zones.n, 0);
3182 wallcycle_sub_stop(wcycle, ewcsDD_SETUPCOMM);
3185 write_dd_pdb("dd_home",step,"dump",top_global,cr,
3186 -1,state_local->x.rvec_array(),state_local->box);
3189 wallcycle_sub_start(wcycle, ewcsDD_MAKETOP);
3191 /* Extract a local topology from the global topology */
3192 for (int i = 0; i < dd->ndim; i++)
3194 np[dd->dim[i]] = comm->cd[i].numPulses();
3196 dd_make_local_top(dd,
3198 dd->unitCellInfo.npbcdim,
3203 state_local->x.rvec_array(),
3207 wallcycle_sub_stop(wcycle, ewcsDD_MAKETOP);
3209 wallcycle_sub_start(wcycle, ewcsDD_MAKECONSTR);
3211 /* Set up the special atom communication */
3212 int n = comm->atomRanges.end(DDAtomRanges::Type::Zones);
3213 for (int i = static_cast<int>(DDAtomRanges::Type::Zones) + 1;
3214 i < static_cast<int>(DDAtomRanges::Type::Number);
3217 auto range = static_cast<DDAtomRanges::Type>(i);
3220 case DDAtomRanges::Type::Vsites:
3221 if (vsite && vsite->numInterUpdategroupVirtualSites())
3223 n = dd_make_local_vsites(dd, n, top_local->idef.il);
3226 case DDAtomRanges::Type::Constraints:
3227 if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
3229 /* Only for inter-cg constraints we need special code */
3230 n = dd_make_local_constraints(dd,
3235 inputrec.nProjOrder,
3236 top_local->idef.il);
3239 default: gmx_incons("Unknown special atom type setup");
3241 comm->atomRanges.setEnd(range, n);
3244 wallcycle_sub_stop(wcycle, ewcsDD_MAKECONSTR);
3246 wallcycle_sub_start(wcycle, ewcsDD_TOPOTHER);
3248 /* Make space for the extra coordinates for virtual site
3249 * or constraint communication.
3251 state_local->natoms = comm->atomRanges.numAtomsTotal();
3253 state_change_natoms(state_local, state_local->natoms);
3255 if (vsite && vsite->numInterUpdategroupVirtualSites())
3257 nat_f_novirsum = comm->atomRanges.end(DDAtomRanges::Type::Vsites);
3261 if (EEL_FULL(inputrec.coulombtype) && dd->haveExclusions)
3263 nat_f_novirsum = comm->atomRanges.end(DDAtomRanges::Type::Zones);
3267 nat_f_novirsum = comm->atomRanges.numHomeAtoms();
3271 /* Set the number of atoms required for the force calculation.
3272 * Forces need to be constrained when doing energy
3273 * minimization. For simple simulations we could avoid some
3274 * allocation, zeroing and copying, but this is probably not worth
3275 * the complications and checking.
3277 forcerec_set_ranges(fr,
3278 comm->atomRanges.end(DDAtomRanges::Type::Zones),
3279 comm->atomRanges.end(DDAtomRanges::Type::Constraints),
3282 /* Update atom data for mdatoms and several algorithms */
3283 mdAlgorithmsSetupAtomData(cr, inputrec, top_global, top_local, fr, f, mdAtoms, constr, vsite, nullptr);
3285 auto mdatoms = mdAtoms->mdatoms();
3286 if (!thisRankHasDuty(cr, DUTY_PME))
3288 /* Send the charges and/or c6/sigmas to our PME only node */
3289 gmx_pme_send_parameters(cr,
3291 mdatoms->nChargePerturbed != 0,
3292 mdatoms->nTypePerturbed != 0,
3293 gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
3294 gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
3295 gmx::arrayRefFromArray(mdatoms->sqrt_c6A, mdatoms->nr),
3296 gmx::arrayRefFromArray(mdatoms->sqrt_c6B, mdatoms->nr),
3297 gmx::arrayRefFromArray(mdatoms->sigmaA, mdatoms->nr),
3298 gmx::arrayRefFromArray(mdatoms->sigmaB, mdatoms->nr),
3299 dd_pme_maxshift_x(*dd),
3300 dd_pme_maxshift_y(*dd));
3303 if (dd->atomSets != nullptr)
3305 /* Update the local atom sets */
3306 dd->atomSets->setIndicesInDomainDecomposition(*(dd->ga2la));
3309 // The pull group construction can need the atom sets updated above
3312 /* Update the local pull groups */
3313 dd_make_local_pull_groups(cr, pull_work);
3316 /* Update the local atoms to be communicated via the IMD protocol if bIMD is TRUE. */
3317 imdSession->dd_make_local_IMD_atoms(dd);
3319 add_dd_statistics(dd);
3321 /* Make sure we only count the cycles for this DD partitioning */
3322 clear_dd_cycle_counts(dd);
3324 /* Because the order of the atoms might have changed since
3325 * the last vsite construction, we need to communicate the constructing
3326 * atom coordinates again (for spreading the forces this MD step).
3328 dd_move_x_vsites(*dd, state_local->box, state_local->x.rvec_array());
3330 wallcycle_sub_stop(wcycle, ewcsDD_TOPOTHER);
3332 if (comm->ddSettings.nstDDDump > 0 && step % comm->ddSettings.nstDDDump == 0)
3334 dd_move_x(dd, state_local->box, state_local->x, nullWallcycle);
3335 write_dd_pdb("dd_dump",
3341 state_local->x.rvec_array(),
3345 /* Store the partitioning step */
3346 comm->partition_step = step;
3348 /* Increase the DD partitioning counter */
3350 /* The state currently matches this DD partitioning count, store it */
3351 state_local->ddp_count = dd->ddp_count;
3354 /* The DD master node knows the complete cg distribution,
3355 * store the count so we can possibly skip the cg info communication.
3357 comm->master_cg_ddp_count = (bSortCG ? 0 : dd->ddp_count);
3360 if (comm->ddSettings.DD_debug > 0)
3362 /* Set the env var GMX_DD_DEBUG if you suspect corrupted indices */
3363 check_index_consistency(dd, top_global.natoms, "after partitioning");
3366 wallcycle_stop(wcycle, ewcDOMDEC);
3369 /*! \brief Check whether bonded interactions are missing, if appropriate */
3370 void checkNumberOfBondedInteractions(const gmx::MDLogger& mdlog,
3372 int totalNumberOfBondedInteractions,
3373 const gmx_mtop_t& top_global,
3374 const gmx_localtop_t* top_local,
3375 gmx::ArrayRef<const gmx::RVec> x,
3377 bool* shouldCheckNumberOfBondedInteractions)
3379 if (*shouldCheckNumberOfBondedInteractions)
3381 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
3383 dd_print_missing_interactions(
3384 mdlog, cr, totalNumberOfBondedInteractions, top_global, top_local, x, box); // Does not return
3386 *shouldCheckNumberOfBondedInteractions = false;