Make temperature and pressure coupling enums enum classes
[alexxy/gromacs.git] / src / gromacs / domdec / partition.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief This file defines functions for mdrun to call to make a new
38  * domain decomposition, and check it.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_domdec
42  */
43
44 #include "gmxpre.h"
45
46 #include "partition.h"
47
48 #include "config.h"
49
50 #include <cassert>
51 #include <cstdio>
52
53 #include <algorithm>
54
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"
95
96 #include "box.h"
97 #include "cellsizes.h"
98 #include "distribute.h"
99 #include "domdec_constraints.h"
100 #include "domdec_internal.h"
101 #include "domdec_vsite.h"
102 #include "dump.h"
103 #include "redistribute.h"
104 #include "utility.h"
105
106 /*! \brief Turn on DLB when the load imbalance causes this amount of total loss.
107  *
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.
110  */
111 #define DD_PERF_LOSS_DLB_ON 0.02
112
113 //! Warn about imbalance due to PP or PP/PME load imbalance at this loss.
114 #define DD_PERF_LOSS_WARN 0.05
115
116
117 //! Debug helper printing a DD zone
118 static void print_ddzone(FILE* fp, int d, int i, int j, gmx_ddzone_t* zone)
119 {
120     fprintf(fp,
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 "
122             "%6.3f\n",
123             d,
124             i,
125             j,
126             zone->min0,
127             zone->max1,
128             zone->mch0,
129             zone->mch0,
130             zone->p1_0,
131             zone->p1_1);
132 }
133
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.
138  */
139 static void dd_move_cellx(gmx_domdec_t* dd, const gmx_ddbox_t* ddbox, rvec cell_ns_x0, rvec cell_ns_x1)
140 {
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;
146
147     rvec extr_s[2];
148     rvec extr_r[2];
149     for (int d = 1; d < dd->ndim; d++)
150     {
151         int           dim = dd->dim[d];
152         gmx_ddzone_t& zp  = (d == 1) ? comm->zone_d1[0] : comm->zone_d2[0][0];
153
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];
162         zp.dataSet = 1;
163     }
164
165     gmx::ArrayRef<DDCellsizesWithDlb> cellsizes = comm->cellsizesWithDlb;
166
167     /* Loop backward over the dimensions and aggregate the extremes
168      * of the cell sizes.
169      */
170     for (int d = dd->ndim - 2; d >= 0; d--)
171     {
172         const int  dim      = dd->dim[d];
173         const bool applyPbc = (dim < ddbox->npbcdim);
174
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;
179
180         int pos = 0;
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.
184          */
185         for (int d1 = d; d1 < dd->ndim - 1; d1++)
186         {
187             gmx_ddzone_t& buf = buf_s[pos];
188
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];
193             buf.mch0 = 0;
194             buf.mch1 = 0;
195             /* Store the cell corner of the dimension we communicate along */
196             buf.p1_0    = comm->cell_x0[dim];
197             buf.p1_1    = 0;
198             buf.dataSet = 1;
199             pos++;
200         }
201
202         buf_s[pos] = (dd->ndim == 2) ? comm->zone_d1[0] : comm->zone_d2[0][0];
203         pos++;
204
205         if (dd->ndim == 3 && d == 0)
206         {
207             buf_s[pos] = comm->zone_d2[0][1];
208             pos++;
209             buf_s[pos] = comm->zone_d1[0];
210             pos++;
211         }
212
213         /* We only need to communicate the extremes
214          * in the forward direction
215          */
216         int numPulses = comm->cd[d].numPulses();
217         int numPulsesMin;
218         if (applyPbc)
219         {
220             /* Take the minimum to avoid double communication */
221             numPulsesMin = std::min(numPulses, dd->numCells[dim] - 1 - numPulses);
222         }
223         else
224         {
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.
229              */
230             numPulsesMin = numPulses;
231         }
232         for (int pulse = 0; pulse < numPulsesMin; pulse++)
233         {
234             /* Communicate the extremes forward */
235             bool receiveValidData = (applyPbc || dd->ci[dim] > 0);
236
237             int numElements = dd->ndim - d - 1;
238             ddSendrecv(dd, d, dddirForward, extr_s + d, numElements, extr_r + d, numElements);
239
240             if (receiveValidData)
241             {
242                 for (int d1 = d; d1 < dd->ndim - 1; d1++)
243                 {
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]);
247                 }
248             }
249         }
250
251         const int numElementsInBuffer = pos;
252         for (int pulse = 0; pulse < numPulses; pulse++)
253         {
254             /* Communicate all the zone information backward */
255             bool receiveValidData = (applyPbc || dd->ci[dim] < dd->numCells[dim] - 1);
256
257             static_assert(
258                     sizeof(gmx_ddzone_t) == c_ddzoneNumReals * sizeof(real),
259                     "Here we expect gmx_ddzone_t to consist of c_ddzoneNumReals reals (only)");
260
261             int numReals = numElementsInBuffer * c_ddzoneNumReals;
262             ddSendrecv(dd,
263                        d,
264                        dddirBackward,
265                        gmx::arrayRefFromArray(&buf_s[0].min0, numReals),
266                        gmx::arrayRefFromArray(&buf_r[0].min0, numReals));
267
268             rvec dh = { 0 };
269             if (pulse > 0)
270             {
271                 for (int d1 = d + 1; d1 < dd->ndim; d1++)
272                 {
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.
276                      */
277                     real dist_d = comm->cell_x1[dim] - buf_r[0].p1_0;
278
279                     real c;
280                     if (ddbox->tric_dir[dim])
281                     {
282                         /* c is the off-diagonal coupling between the cell planes
283                          * along directions d and d1.
284                          */
285                         c = ddbox->v[dim][dd->dim[d1]][dim];
286                     }
287                     else
288                     {
289                         c = 0;
290                     }
291                     real det = (1 + c * c) * gmx::square(comm->systemInfo.cutoff) - dist_d * dist_d;
292                     if (det > 0)
293                     {
294                         dh[d1] = comm->systemInfo.cutoff - (c * dist_d + std::sqrt(det)) / (1 + c * c);
295                     }
296                     else
297                     {
298                         /* A negative value signals out of range */
299                         dh[d1] = -1;
300                     }
301                 }
302             }
303
304             /* Accumulate the extremes over all pulses */
305             for (int i = 0; i < numElementsInBuffer; i++)
306             {
307                 if (pulse == 0)
308                 {
309                     buf_e[i] = buf_r[i];
310                 }
311                 else
312                 {
313                     if (receiveValidData)
314                     {
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);
318                     }
319
320                     int d1;
321                     if (dd->ndim == 3 && d == 0 && i == numElementsInBuffer - 1)
322                     {
323                         d1 = 1;
324                     }
325                     else
326                     {
327                         d1 = d + 1;
328                     }
329                     if (receiveValidData && dh[d1] >= 0)
330                     {
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]);
333                     }
334                 }
335                 /* Copy the received buffer to the send buffer,
336                  * to pass the data through with the next pulse.
337                  */
338                 buf_s[i] = buf_r[i];
339             }
340             if (((applyPbc || dd->ci[dim] + numPulses < dd->numCells[dim]) && pulse == numPulses - 1)
341                 || (!applyPbc && dd->ci[dim] + 1 + pulse == dd->numCells[dim] - 1))
342             {
343                 /* Store the extremes */
344                 int pos = 0;
345
346                 for (int d1 = d; d1 < dd->ndim - 1; d1++)
347                 {
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);
351                     pos++;
352                 }
353
354                 if (d == 1 || (d == 0 && dd->ndim == 3))
355                 {
356                     for (int i = d; i < 2; i++)
357                     {
358                         comm->zone_d2[1 - d][i] = buf_e[pos];
359                         pos++;
360                     }
361                 }
362                 if (d == 0)
363                 {
364                     comm->zone_d1[1] = buf_e[pos];
365                     pos++;
366                 }
367             }
368             else
369             {
370                 if (d == 1 || (d == 0 && dd->ndim == 3))
371                 {
372                     for (int i = d; i < 2; i++)
373                     {
374                         comm->zone_d2[1 - d][i].dataSet = 0;
375                     }
376                 }
377                 if (d == 0)
378                 {
379                     comm->zone_d1[1].dataSet = 0;
380                 }
381             }
382         }
383     }
384
385     if (dd->ndim >= 2)
386     {
387         int dim = dd->dim[1];
388         for (int i = 0; i < 2; i++)
389         {
390             if (comm->zone_d1[i].dataSet != 0)
391             {
392                 if (debug)
393                 {
394                     print_ddzone(debug, 1, i, 0, &comm->zone_d1[i]);
395                 }
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);
398             }
399         }
400     }
401     if (dd->ndim >= 3)
402     {
403         int dim = dd->dim[2];
404         for (int i = 0; i < 2; i++)
405         {
406             for (int j = 0; j < 2; j++)
407             {
408                 if (comm->zone_d2[i][j].dataSet != 0)
409                 {
410                     if (debug)
411                     {
412                         print_ddzone(debug, 2, i, j, &comm->zone_d2[i][j]);
413                     }
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);
416                 }
417             }
418         }
419     }
420     for (int d = 1; d < dd->ndim; d++)
421     {
422         cellsizes[d].fracLowerMax = extr_s[d - 1][0];
423         cellsizes[d].fracUpperMin = extr_s[d - 1][1];
424         if (debug)
425         {
426             fprintf(debug,
427                     "Cell fraction d %d, max0 %f, min1 %f\n",
428                     d,
429                     cellsizes[d].fracLowerMax,
430                     cellsizes[d].fracUpperMin);
431         }
432     }
433 }
434
435 //! Sets the charge-group zones to be equal to the home zone.
436 static void set_zones_ncg_home(gmx_domdec_t* dd)
437 {
438     gmx_domdec_zones_t* zones;
439     int                 i;
440
441     zones = &dd->comm->zones;
442
443     zones->cg_range[0] = 0;
444     for (i = 1; i < zones->n + 1; i++)
445     {
446         zones->cg_range[i] = dd->ncg_home;
447     }
448     /* zone_ncg1[0] should always be equal to ncg_home */
449     dd->comm->zone_ncg1[0] = dd->ncg_home;
450 }
451
452 //! Restore atom groups for the charge groups.
453 static void restoreAtomGroups(gmx_domdec_t* dd, const t_state* state)
454 {
455     gmx::ArrayRef<const int> atomGroupsState = state->cg_gl;
456
457     std::vector<int>& globalAtomGroupIndices = dd->globalAtomGroupIndices;
458
459     globalAtomGroupIndices.resize(atomGroupsState.size());
460
461     /* Copy back the global charge group indices from state
462      * and rebuild the local charge group to atom index.
463      */
464     for (gmx::index i = 0; i < atomGroupsState.ssize(); i++)
465     {
466         globalAtomGroupIndices[i] = atomGroupsState[i];
467     }
468
469     dd->ncg_home = atomGroupsState.size();
470     dd->comm->atomRanges.setEnd(DDAtomRanges::Type::Home, atomGroupsState.ssize());
471
472     set_zones_ncg_home(dd);
473 }
474
475 //! Sets the cginfo structures.
476 static void dd_set_cginfo(gmx::ArrayRef<const int> index_gl, int cg0, int cg1, t_forcerec* fr)
477 {
478     if (fr != nullptr)
479     {
480         gmx::ArrayRef<cginfo_mb_t> cginfo_mb = fr->cginfo_mb;
481         gmx::ArrayRef<int>         cginfo    = fr->cginfo;
482
483         for (int cg = cg0; cg < cg1; cg++)
484         {
485             cginfo[cg] = ddcginfo(cginfo_mb, index_gl[cg]);
486         }
487     }
488 }
489
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)
492 {
493     const int                numZones               = dd->comm->zones.n;
494     const int*               zone2cg                = dd->comm->zones.cg_range;
495     const int*               zone_ncg1              = dd->comm->zone_ncg1;
496     gmx::ArrayRef<const int> globalAtomGroupIndices = dd->globalAtomGroupIndices;
497
498     std::vector<int>& globalAtomIndices = dd->globalAtomIndices;
499     gmx_ga2la_t&      ga2la             = *dd->ga2la;
500
501     if (zone2cg[1] != dd->ncg_home)
502     {
503         gmx_incons("dd->ncg_zone is not up to date");
504     }
505
506     /* Make the local to global and global to local atom index */
507     int a = atomStart;
508     globalAtomIndices.resize(a);
509     for (int zone = 0; zone < numZones; zone++)
510     {
511         int cg0;
512         if (zone == 0)
513         {
514             cg0 = atomStart;
515         }
516         else
517         {
518             cg0 = zone2cg[zone];
519         }
520         int cg1    = zone2cg[zone + 1];
521         int cg1_p1 = cg0 + zone_ncg1[zone];
522
523         for (int cg = cg0; cg < cg1; cg++)
524         {
525             int zone1 = zone;
526             if (cg >= cg1_p1)
527             {
528                 /* Signal that this cg is from more than one pulse away */
529                 zone1 += numZones;
530             }
531             int cg_gl = globalAtomGroupIndices[cg];
532             globalAtomIndices.push_back(cg_gl);
533             ga2la.insert(cg_gl, { a, zone1 });
534             a++;
535         }
536     }
537 }
538
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)
541 {
542     int nerr = 0;
543
544     const int numAtomsInZones = dd->comm->atomRanges.end(DDAtomRanges::Type::Zones);
545
546     if (dd->comm->ddSettings.DD_debug > 1)
547     {
548         std::vector<int> have(natoms_sys);
549         for (int a = 0; a < numAtomsInZones; a++)
550         {
551             int globalAtomIndex = dd->globalAtomIndices[a];
552             if (have[globalAtomIndex] > 0)
553             {
554                 fprintf(stderr,
555                         "DD rank %d: global atom %d occurs twice: index %d and %d\n",
556                         dd->rank,
557                         globalAtomIndex + 1,
558                         have[globalAtomIndex],
559                         a + 1);
560             }
561             else
562             {
563                 have[globalAtomIndex] = a + 1;
564             }
565         }
566     }
567
568     std::vector<int> have(numAtomsInZones);
569
570     int ngl = 0;
571     for (int i = 0; i < natoms_sys; i++)
572     {
573         if (const auto entry = dd->ga2la->find(i))
574         {
575             const int a = entry->la;
576             if (a >= numAtomsInZones)
577             {
578                 fprintf(stderr,
579                         "DD rank %d: global atom %d marked as local atom %d, which is larger than "
580                         "nat_tot (%d)\n",
581                         dd->rank,
582                         i + 1,
583                         a + 1,
584                         numAtomsInZones);
585                 nerr++;
586             }
587             else
588             {
589                 have[a] = 1;
590                 if (dd->globalAtomIndices[a] != i)
591                 {
592                     fprintf(stderr,
593                             "DD rank %d: global atom %d marked as local atom %d, which has global "
594                             "atom index %d\n",
595                             dd->rank,
596                             i + 1,
597                             a + 1,
598                             dd->globalAtomIndices[a] + 1);
599                     nerr++;
600                 }
601             }
602             ngl++;
603         }
604     }
605     if (ngl != numAtomsInZones)
606     {
607         fprintf(stderr, "DD rank %d, %s: %d global atom indices, %d local atoms\n", dd->rank, where, ngl, numAtomsInZones);
608     }
609     for (int a = 0; a < numAtomsInZones; a++)
610     {
611         if (have[a] == 0)
612         {
613             fprintf(stderr,
614                     "DD rank %d, %s: local atom %d, global %d has no global index\n",
615                     dd->rank,
616                     where,
617                     a + 1,
618                     dd->globalAtomIndices[a] + 1);
619         }
620     }
621
622     if (nerr > 0)
623     {
624         gmx_fatal(FARGS, "DD rank %d, %s: %d atom(group) index inconsistencies", dd->rank, where, nerr);
625     }
626 }
627
628 //! Clear all DD global state indices
629 static void clearDDStateIndices(gmx_domdec_t* dd, const bool keepLocalAtomIndices)
630 {
631     gmx_ga2la_t& ga2la = *dd->ga2la;
632
633     if (!keepLocalAtomIndices)
634     {
635         /* Clear the whole list without the overhead of searching */
636         ga2la.clear();
637     }
638     else
639     {
640         const int numAtomsInZones = dd->comm->atomRanges.end(DDAtomRanges::Type::Zones);
641         for (int i = 0; i < numAtomsInZones; i++)
642         {
643             ga2la.erase(dd->globalAtomIndices[i]);
644         }
645     }
646
647     dd_clear_local_vsite_indices(dd);
648
649     if (dd->constraints)
650     {
651         dd_clear_local_constraint_indices(dd);
652     }
653 }
654
655 bool check_grid_jump(int64_t step, const gmx_domdec_t* dd, real cutoff, const gmx_ddbox_t* ddbox, gmx_bool bFatal)
656 {
657     gmx_domdec_comm_t* comm    = dd->comm;
658     bool               invalid = false;
659
660     for (int d = 1; d < dd->ndim; d++)
661     {
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])
667         {
668             bfac *= ddbox->skew_fac[dim];
669         }
670         if ((cellsizes.fracUpper - cellsizes.fracLowerMax) * bfac < limit
671             || (cellsizes.fracLower - cellsizes.fracUpperMin) * bfac > -limit)
672         {
673             invalid = true;
674
675             if (bFatal)
676             {
677                 char buf[22];
678
679                 /* This error should never be triggered under normal
680                  * circumstances, but you never know ...
681                  */
682                 gmx_fatal(FARGS,
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),
687                           dim2char(dim),
688                           dd->ci[XX],
689                           dd->ci[YY],
690                           dd->ci[ZZ]);
691             }
692         }
693     }
694
695     return invalid;
696 }
697 //! Return the duration of force calculations on this rank.
698 static float dd_force_load(gmx_domdec_comm_t* comm)
699 {
700     float load;
701
702     if (comm->ddSettings.eFlop)
703     {
704         load = comm->flop;
705         if (comm->ddSettings.eFlop > 1)
706         {
707             load *= 1.0 + (comm->ddSettings.eFlop - 1) * (0.1 * rand() / RAND_MAX - 0.05);
708         }
709     }
710     else
711     {
712         load = comm->cycl[ddCyclF];
713         if (comm->cycl_n[ddCyclF] > 1)
714         {
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.
719              */
720             load -= comm->cycl_max[ddCyclF];
721         }
722
723 #if GMX_MPI
724         if (comm->cycl_n[ddCyclWaitGPU] && comm->nrank_gpu_shared > 1)
725         {
726             float gpu_wait, gpu_wait_sum;
727
728             gpu_wait = comm->cycl[ddCyclWaitGPU];
729             if (comm->cycl_n[ddCyclF] > 1)
730             {
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.
740                  */
741                 gpu_wait *= (comm->cycl_n[ddCyclF] - 1) / static_cast<float>(comm->cycl_n[ddCyclF]);
742             }
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;
747         }
748 #endif
749     }
750
751     return load;
752 }
753
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)
756 {
757     gmx_domdec_comm_t* comm;
758     int                dim_ind, dim;
759
760     comm = dd->comm;
761
762     for (dim_ind = 0; dim_ind < dd->ndim; dim_ind++)
763     {
764         dim = dd->dim[dim_ind];
765
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))
768             && isDlbOn(comm)
769             && (comm->cell_x1[dim] - comm->cell_x0[dim]) * ddbox->skew_fac[dim] < comm->cellsize_min[dim])
770         {
771             char buf[22];
772             gmx_fatal(FARGS,
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 "
775                       "%d %d %d",
776                       gmx_step_str(step, buf),
777                       dim2char(dim),
778                       comm->cell_x1[dim] - comm->cell_x0[dim],
779                       ddbox->skew_fac[dim],
780                       dd->comm->cellsize_min[dim],
781                       dd->ci[XX],
782                       dd->ci[YY],
783                       dd->ci[ZZ]);
784         }
785     }
786
787     if ((isDlbOn(dd->comm) && dd->ndim > 1) || ddbox->nboundeddim < DIM)
788     {
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)
792         {
793             check_grid_jump(step, dd, dd->comm->systemInfo.cutoff, ddbox, TRUE);
794         }
795     }
796 }
797
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)
800 {
801     gmx_domdec_comm_t* comm;
802     domdec_load_t*     load;
803     float              cell_frac = 0, sbuf[DD_NLOAD_MAX];
804     gmx_bool           bSepPME;
805
806     if (debug)
807     {
808         fprintf(debug, "get_load_distribution start\n");
809     }
810
811     wallcycle_start(wcycle, ewcDDCOMMLOAD);
812
813     comm = dd->comm;
814
815     bSepPME = (dd->pme_nodeid >= 0);
816
817     if (dd->ndim == 0 && bSepPME)
818     {
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];
822     }
823
824     for (int d = dd->ndim - 1; d >= 0; d--)
825     {
826         const DDCellsizesWithDlb* cellsizes = (isDlbOn(dd->comm) ? &comm->cellsizesWithDlb[d] : nullptr);
827         const int                 dim       = dd->dim[d];
828         /* Check if we participate in the communication in this dimension */
829         if (d == dd->ndim - 1 || (dd->ci[dd->dim[d + 1]] == 0 && dd->ci[dd->dim[dd->ndim - 1]] == 0))
830         {
831             load = &comm->load[d];
832             if (isDlbOn(dd->comm))
833             {
834                 cell_frac = cellsizes->fracUpper - cellsizes->fracLower;
835             }
836             int pos = 0;
837             if (d == dd->ndim - 1)
838             {
839                 sbuf[pos++] = dd_force_load(comm);
840                 sbuf[pos++] = sbuf[0];
841                 if (isDlbOn(dd->comm))
842                 {
843                     sbuf[pos++] = sbuf[0];
844                     sbuf[pos++] = cell_frac;
845                     if (d > 0)
846                     {
847                         sbuf[pos++] = cellsizes->fracLowerMax;
848                         sbuf[pos++] = cellsizes->fracUpperMin;
849                     }
850                 }
851                 if (bSepPME)
852                 {
853                     sbuf[pos++] = comm->cycl[ddCyclPPduringPME];
854                     sbuf[pos++] = comm->cycl[ddCyclPME];
855                 }
856             }
857             else
858             {
859                 sbuf[pos++] = comm->load[d + 1].sum;
860                 sbuf[pos++] = comm->load[d + 1].max;
861                 if (isDlbOn(dd->comm))
862                 {
863                     sbuf[pos++] = comm->load[d + 1].sum_m;
864                     sbuf[pos++] = comm->load[d + 1].cvol_min * cell_frac;
865                     sbuf[pos++] = comm->load[d + 1].flags;
866                     if (d > 0)
867                     {
868                         sbuf[pos++] = cellsizes->fracLowerMax;
869                         sbuf[pos++] = cellsizes->fracUpperMin;
870                     }
871                 }
872                 if (bSepPME)
873                 {
874                     sbuf[pos++] = comm->load[d + 1].mdf;
875                     sbuf[pos++] = comm->load[d + 1].pme;
876                 }
877             }
878             load->nload = pos;
879             /* Communicate a row in DD direction d.
880              * The communicators are setup such that the root always has rank 0.
881              */
882 #if GMX_MPI
883             MPI_Gather(sbuf,
884                        load->nload * sizeof(float),
885                        MPI_BYTE,
886                        load->load,
887                        load->nload * sizeof(float),
888                        MPI_BYTE,
889                        0,
890                        comm->mpi_comm_load[d]);
891 #endif
892             if (dd->ci[dim] == dd->master_ci[dim])
893             {
894                 /* We are the master along this row, process this row */
895                 RowMaster* rowMaster = nullptr;
896
897                 if (isDlbOn(comm))
898                 {
899                     rowMaster = cellsizes->rowMaster.get();
900                 }
901                 load->sum      = 0;
902                 load->max      = 0;
903                 load->sum_m    = 0;
904                 load->cvol_min = 1;
905                 load->flags    = 0;
906                 load->mdf      = 0;
907                 load->pme      = 0;
908                 int pos        = 0;
909                 for (int i = 0; i < dd->numCells[dim]; i++)
910                 {
911                     load->sum += load->load[pos++];
912                     load->max = std::max(load->max, load->load[pos]);
913                     pos++;
914                     if (isDlbOn(dd->comm))
915                     {
916                         if (rowMaster->dlbIsLimited)
917                         {
918                             /* This direction could not be load balanced properly,
919                              * therefore we need to use the maximum iso the average load.
920                              */
921                             load->sum_m = std::max(load->sum_m, load->load[pos]);
922                         }
923                         else
924                         {
925                             load->sum_m += load->load[pos];
926                         }
927                         pos++;
928                         load->cvol_min = std::min(load->cvol_min, load->load[pos]);
929                         pos++;
930                         if (d < dd->ndim - 1)
931                         {
932                             load->flags = gmx::roundToInt(load->load[pos++]);
933                         }
934                         if (d > 0)
935                         {
936                             rowMaster->bounds[i].cellFracLowerMax = load->load[pos++];
937                             rowMaster->bounds[i].cellFracUpperMin = load->load[pos++];
938                         }
939                     }
940                     if (bSepPME)
941                     {
942                         load->mdf = std::max(load->mdf, load->load[pos]);
943                         pos++;
944                         load->pme = std::max(load->pme, load->load[pos]);
945                         pos++;
946                     }
947                 }
948                 if (isDlbOn(comm) && rowMaster->dlbIsLimited)
949                 {
950                     load->sum_m *= dd->numCells[dim];
951                     load->flags |= (1 << d);
952                 }
953             }
954         }
955     }
956
957     if (DDMASTER(dd))
958     {
959         comm->nload += dd_load_count(comm);
960         comm->load_step += comm->cycl[ddCyclStep];
961         comm->load_sum += comm->load[0].sum;
962         comm->load_max += comm->load[0].max;
963         if (isDlbOn(comm))
964         {
965             for (int d = 0; d < dd->ndim; d++)
966             {
967                 if (comm->load[0].flags & (1 << d))
968                 {
969                     comm->load_lim[d]++;
970                 }
971             }
972         }
973         if (bSepPME)
974         {
975             comm->load_mdf += comm->load[0].mdf;
976             comm->load_pme += comm->load[0].pme;
977         }
978     }
979
980     wallcycle_stop(wcycle, ewcDDCOMMLOAD);
981
982     if (debug)
983     {
984         fprintf(debug, "get_load_distribution finished\n");
985     }
986 }
987
988 /*! \brief Return the relative performance loss on the total run time
989  * due to the force calculation load imbalance. */
990 static float dd_force_load_fraction(gmx_domdec_t* dd)
991 {
992     if (dd->comm->nload > 0 && dd->comm->load_step > 0)
993     {
994         return dd->comm->load_sum / (dd->comm->load_step * dd->nnodes);
995     }
996     else
997     {
998         return 0;
999     }
1000 }
1001
1002 /*! \brief Return the relative performance loss on the total run time
1003  * due to the force calculation load imbalance. */
1004 static float dd_force_imb_perf_loss(gmx_domdec_t* dd)
1005 {
1006     if (dd->comm->nload > 0 && dd->comm->load_step > 0)
1007     {
1008         return (dd->comm->load_max * dd->nnodes - dd->comm->load_sum) / (dd->comm->load_step * dd->nnodes);
1009     }
1010     else
1011     {
1012         return 0;
1013     }
1014 }
1015
1016 //! Print load-balance report e.g. at the end of a run.
1017 static void print_dd_load_av(FILE* fplog, gmx_domdec_t* dd)
1018 {
1019     gmx_domdec_comm_t* comm = dd->comm;
1020
1021     /* Only the master rank prints loads and only if we measured loads */
1022     if (!DDMASTER(dd) || comm->nload == 0)
1023     {
1024         return;
1025     }
1026
1027     char buf[STRLEN];
1028     int  numPpRanks  = dd->nnodes;
1029     int  numPmeRanks = (comm->ddRankSetup.usePmeOnlyRanks ? comm->ddRankSetup.numRanksDoingPme : 0);
1030     int  numRanks    = numPpRanks + numPmeRanks;
1031     float lossFraction = 0;
1032
1033     /* Print the average load imbalance and performance loss */
1034     if (dd->nnodes > 1 && comm->load_sum > 0)
1035     {
1036         float imbalance = comm->load_max * numPpRanks / comm->load_sum - 1;
1037         lossFraction    = dd_force_imb_perf_loss(dd);
1038
1039         std::string msg = "\nDynamic load balancing report:\n";
1040         std::string dlbStateStr;
1041
1042         switch (dd->comm->dlbState)
1043         {
1044             case DlbState::offUser:
1045                 dlbStateStr = "DLB was off during the run per user request.";
1046                 break;
1047             case DlbState::offForever:
1048                 /* Currectly this can happen due to performance loss observed, cell size
1049                  * limitations or incompatibility with other settings observed during
1050                  * determineInitialDlbState(). */
1051                 dlbStateStr = "DLB got disabled because it was unsuitable to use.";
1052                 break;
1053             case DlbState::offCanTurnOn:
1054                 dlbStateStr = "DLB was off during the run due to low measured imbalance.";
1055                 break;
1056             case DlbState::offTemporarilyLocked:
1057                 dlbStateStr =
1058                         "DLB was locked at the end of the run due to unfinished PP-PME "
1059                         "balancing.";
1060                 break;
1061             case DlbState::onCanTurnOff:
1062                 dlbStateStr = "DLB was turned on during the run due to measured imbalance.";
1063                 break;
1064             case DlbState::onUser:
1065                 dlbStateStr = "DLB was permanently on during the run per user request.";
1066                 break;
1067             default: GMX_ASSERT(false, "Undocumented DLB state");
1068         }
1069
1070         msg += " " + dlbStateStr + "\n";
1071         msg += gmx::formatString(" Average load imbalance: %.1f%%.\n", imbalance * 100);
1072         msg += gmx::formatString(
1073                 " The balanceable part of the MD step is %d%%, load imbalance is computed from "
1074                 "this.\n",
1075                 gmx::roundToInt(dd_force_load_fraction(dd) * 100));
1076         msg += gmx::formatString(
1077                 " Part of the total run time spent waiting due to load imbalance: %.1f%%.\n",
1078                 lossFraction * 100);
1079         fprintf(fplog, "%s", msg.c_str());
1080         fprintf(stderr, "\n%s", msg.c_str());
1081     }
1082
1083     /* Print during what percentage of steps the  load balancing was limited */
1084     bool dlbWasLimited = false;
1085     if (isDlbOn(comm))
1086     {
1087         sprintf(buf, " Steps where the load balancing was limited by -rdd, -rcon and/or -dds:");
1088         for (int d = 0; d < dd->ndim; d++)
1089         {
1090             int limitPercentage = (200 * comm->load_lim[d] + 1) / (2 * comm->nload);
1091             sprintf(buf + strlen(buf), " %c %d %%", dim2char(dd->dim[d]), limitPercentage);
1092             if (limitPercentage >= 50)
1093             {
1094                 dlbWasLimited = true;
1095             }
1096         }
1097         sprintf(buf + strlen(buf), "\n");
1098         fprintf(fplog, "%s", buf);
1099         fprintf(stderr, "%s", buf);
1100     }
1101
1102     /* Print the performance loss due to separate PME - PP rank imbalance */
1103     float lossFractionPme = 0;
1104     if (numPmeRanks > 0 && comm->load_mdf > 0 && comm->load_step > 0)
1105     {
1106         float pmeForceRatio = comm->load_pme / comm->load_mdf;
1107         lossFractionPme     = (comm->load_pme - comm->load_mdf) / comm->load_step;
1108         if (lossFractionPme <= 0)
1109         {
1110             lossFractionPme *= numPmeRanks / static_cast<float>(numRanks);
1111         }
1112         else
1113         {
1114             lossFractionPme *= numPpRanks / static_cast<float>(numRanks);
1115         }
1116         sprintf(buf, " Average PME mesh/force load: %5.3f\n", pmeForceRatio);
1117         fprintf(fplog, "%s", buf);
1118         fprintf(stderr, "%s", buf);
1119         sprintf(buf,
1120                 " Part of the total run time spent waiting due to PP/PME imbalance: %.1f %%\n",
1121                 std::fabs(lossFractionPme) * 100);
1122         fprintf(fplog, "%s", buf);
1123         fprintf(stderr, "%s", buf);
1124     }
1125     fprintf(fplog, "\n");
1126     fprintf(stderr, "\n");
1127
1128     if (lossFraction >= DD_PERF_LOSS_WARN)
1129     {
1130         std::string message = gmx::formatString(
1131                 "NOTE: %.1f %% of the available CPU time was lost due to load imbalance\n"
1132                 "      in the domain decomposition.\n",
1133                 lossFraction * 100);
1134
1135         bool hadSuggestion = false;
1136         if (!isDlbOn(comm))
1137         {
1138             message += "      You might want to use dynamic load balancing (option -dlb.)\n";
1139             hadSuggestion = true;
1140         }
1141         else if (dlbWasLimited)
1142         {
1143             message +=
1144                     "      You might want to decrease the cell size limit (options -rdd, -rcon "
1145                     "and/or -dds).\n";
1146             hadSuggestion = true;
1147         }
1148         message += gmx::formatString(
1149                 "      You can %sconsider manually changing the decomposition (option -dd);\n"
1150                 "      e.g. by using fewer domains along the box dimension in which there is\n"
1151                 "      considerable inhomogeneity in the simulated system.",
1152                 hadSuggestion ? "also " : "");
1153
1154
1155         fprintf(fplog, "%s\n", message.c_str());
1156         fprintf(stderr, "%s\n", message.c_str());
1157     }
1158     if (numPmeRanks > 0 && std::fabs(lossFractionPme) >= DD_PERF_LOSS_WARN)
1159     {
1160         sprintf(buf,
1161                 "NOTE: %.1f %% performance was lost because the PME ranks\n"
1162                 "      had %s work to do than the PP ranks.\n"
1163                 "      You might want to %s the number of PME ranks\n"
1164                 "      or %s the cut-off and the grid spacing.\n",
1165                 std::fabs(lossFractionPme * 100),
1166                 (lossFractionPme < 0) ? "less" : "more",
1167                 (lossFractionPme < 0) ? "decrease" : "increase",
1168                 (lossFractionPme < 0) ? "decrease" : "increase");
1169         fprintf(fplog, "%s\n", buf);
1170         fprintf(stderr, "%s\n", buf);
1171     }
1172 }
1173
1174 //! Return the minimum communication volume.
1175 static float dd_vol_min(gmx_domdec_t* dd)
1176 {
1177     return dd->comm->load[0].cvol_min * dd->nnodes;
1178 }
1179
1180 //! Return the DD load flags.
1181 static int dd_load_flags(gmx_domdec_t* dd)
1182 {
1183     return dd->comm->load[0].flags;
1184 }
1185
1186 //! Return the reported load imbalance in force calculations.
1187 static float dd_f_imbal(gmx_domdec_t* dd)
1188 {
1189     if (dd->comm->load[0].sum > 0)
1190     {
1191         return dd->comm->load[0].max * dd->nnodes / dd->comm->load[0].sum - 1.0F;
1192     }
1193     else
1194     {
1195         /* Something is wrong in the cycle counting, report no load imbalance */
1196         return 0.0F;
1197     }
1198 }
1199
1200 //! Returns DD load balance report.
1201 static std::string dd_print_load(gmx_domdec_t* dd, int64_t step)
1202 {
1203     gmx::StringOutputStream stream;
1204     gmx::TextWriter         log(&stream);
1205
1206     int flags = dd_load_flags(dd);
1207     if (flags)
1208     {
1209         log.writeString("DD  load balancing is limited by minimum cell size in dimension");
1210         for (int d = 0; d < dd->ndim; d++)
1211         {
1212             if (flags & (1 << d))
1213             {
1214                 log.writeStringFormatted(" %c", dim2char(dd->dim[d]));
1215             }
1216         }
1217         log.ensureLineBreak();
1218     }
1219     log.writeString("DD  step " + gmx::toString(step));
1220     if (isDlbOn(dd->comm))
1221     {
1222         log.writeStringFormatted("  vol min/aver %5.3f%c", dd_vol_min(dd), flags ? '!' : ' ');
1223     }
1224     if (dd->nnodes > 1)
1225     {
1226         log.writeStringFormatted(" load imb.: force %4.1f%%", dd_f_imbal(dd) * 100);
1227     }
1228     if (dd->comm->cycl_n[ddCyclPME])
1229     {
1230         log.writeStringFormatted("  pme mesh/force %5.3f", dd_pme_f_ratio(dd));
1231     }
1232     log.ensureLineBreak();
1233     return stream.toString();
1234 }
1235
1236 //! Prints DD load balance report in mdrun verbose mode.
1237 static void dd_print_load_verbose(gmx_domdec_t* dd)
1238 {
1239     if (isDlbOn(dd->comm))
1240     {
1241         fprintf(stderr, "vol %4.2f%c ", dd_vol_min(dd), dd_load_flags(dd) ? '!' : ' ');
1242     }
1243     if (dd->nnodes > 1)
1244     {
1245         fprintf(stderr, "imb F %2d%% ", gmx::roundToInt(dd_f_imbal(dd) * 100));
1246     }
1247     if (dd->comm->cycl_n[ddCyclPME])
1248     {
1249         fprintf(stderr, "pme/F %4.2f ", dd_pme_f_ratio(dd));
1250     }
1251 }
1252
1253 //! Turns on dynamic load balancing if possible and needed.
1254 static void turn_on_dlb(const gmx::MDLogger& mdlog, gmx_domdec_t* dd, int64_t step)
1255 {
1256     gmx_domdec_comm_t* comm = dd->comm;
1257
1258     real cellsize_min = comm->cellsize_min[dd->dim[0]];
1259     for (int d = 1; d < dd->ndim; d++)
1260     {
1261         cellsize_min = std::min(cellsize_min, comm->cellsize_min[dd->dim[d]]);
1262     }
1263
1264     /* Turn off DLB if we're too close to the cell size limit. */
1265     if (cellsize_min < comm->cellsize_limit * 1.05)
1266     {
1267         GMX_LOG(mdlog.info)
1268                 .appendTextFormatted(
1269                         "step %s Measured %.1f %% performance loss due to load imbalance, "
1270                         "but the minimum cell size is smaller than 1.05 times the cell size limit. "
1271                         "Will no longer try dynamic load balancing.",
1272                         gmx::toString(step).c_str(),
1273                         dd_force_imb_perf_loss(dd) * 100);
1274
1275         comm->dlbState = DlbState::offForever;
1276         return;
1277     }
1278
1279     GMX_LOG(mdlog.info)
1280             .appendTextFormatted(
1281                     "step %s Turning on dynamic load balancing, because the performance loss due "
1282                     "to load imbalance is %.1f %%.",
1283                     gmx::toString(step).c_str(),
1284                     dd_force_imb_perf_loss(dd) * 100);
1285     comm->dlbState = DlbState::onCanTurnOff;
1286
1287     /* Store the non-DLB performance, so we can check if DLB actually
1288      * improves performance.
1289      */
1290     GMX_RELEASE_ASSERT(comm->cycl_n[ddCyclStep] > 0,
1291                        "When we turned on DLB, we should have measured cycles");
1292     comm->cyclesPerStepBeforeDLB = comm->cycl[ddCyclStep] / comm->cycl_n[ddCyclStep];
1293
1294     set_dlb_limits(dd);
1295
1296     /* We can set the required cell size info here,
1297      * so we do not need to communicate this.
1298      * The grid is completely uniform.
1299      */
1300     for (int d = 0; d < dd->ndim; d++)
1301     {
1302         RowMaster* rowMaster = comm->cellsizesWithDlb[d].rowMaster.get();
1303
1304         if (rowMaster)
1305         {
1306             comm->load[d].sum_m = comm->load[d].sum;
1307
1308             int nc = dd->numCells[dd->dim[d]];
1309             for (int i = 0; i < nc; i++)
1310             {
1311                 rowMaster->cellFrac[i] = i / static_cast<real>(nc);
1312                 if (d > 0)
1313                 {
1314                     rowMaster->bounds[i].cellFracLowerMax = i / static_cast<real>(nc);
1315                     rowMaster->bounds[i].cellFracUpperMin = (i + 1) / static_cast<real>(nc);
1316                 }
1317             }
1318             rowMaster->cellFrac[nc] = 1.0;
1319         }
1320     }
1321 }
1322
1323 //! Turns off dynamic load balancing (but leave it able to turn back on).
1324 static void turn_off_dlb(const gmx::MDLogger& mdlog, gmx_domdec_t* dd, int64_t step)
1325 {
1326     GMX_LOG(mdlog.info)
1327             .appendText(
1328                     "step " + gmx::toString(step)
1329                     + " Turning off dynamic load balancing, because it is degrading performance.");
1330     dd->comm->dlbState                     = DlbState::offCanTurnOn;
1331     dd->comm->haveTurnedOffDlb             = true;
1332     dd->comm->ddPartioningCountFirstDlbOff = dd->ddp_count;
1333 }
1334
1335 //! Turns off dynamic load balancing permanently.
1336 static void turn_off_dlb_forever(const gmx::MDLogger& mdlog, gmx_domdec_t* dd, int64_t step)
1337 {
1338     GMX_RELEASE_ASSERT(dd->comm->dlbState == DlbState::offCanTurnOn,
1339                        "Can only turn off DLB forever when it was in the can-turn-on state");
1340     GMX_LOG(mdlog.info)
1341             .appendText(
1342                     "step " + gmx::toString(step)
1343                     + " Will no longer try dynamic load balancing, as it degraded performance.");
1344     dd->comm->dlbState = DlbState::offForever;
1345 }
1346
1347 void set_dd_dlb_max_cutoff(t_commrec* cr, real cutoff)
1348 {
1349     gmx_domdec_comm_t* comm;
1350
1351     comm = cr->dd->comm;
1352
1353     /* Turn on the DLB limiting (might have been on already) */
1354     comm->bPMELoadBalDLBLimits = TRUE;
1355
1356     /* Change the cut-off limit */
1357     comm->PMELoadBal_max_cutoff = cutoff;
1358
1359     if (debug)
1360     {
1361         fprintf(debug,
1362                 "PME load balancing set a limit to the DLB staggering such that a %f cut-off will "
1363                 "continue to fit\n",
1364                 comm->PMELoadBal_max_cutoff);
1365     }
1366 }
1367
1368 //! Merge atom buffers.
1369 static void merge_cg_buffers(int                            ncell,
1370                              gmx_domdec_comm_dim_t*         cd,
1371                              int                            pulse,
1372                              int*                           ncg_cell,
1373                              gmx::ArrayRef<int>             index_gl,
1374                              const int*                     recv_i,
1375                              gmx::ArrayRef<gmx::RVec>       x,
1376                              gmx::ArrayRef<const gmx::RVec> recv_vr,
1377                              gmx::ArrayRef<cginfo_mb_t>     cginfo_mb,
1378                              gmx::ArrayRef<int>             cginfo)
1379 {
1380     gmx_domdec_ind_t *ind, *ind_p;
1381     int               p, cell, c, cg, cg0, cg1, cg_gl;
1382     int               shift;
1383
1384     ind = &cd->ind[pulse];
1385
1386     /* First correct the already stored data */
1387     shift = ind->nrecv[ncell];
1388     for (cell = ncell - 1; cell >= 0; cell--)
1389     {
1390         shift -= ind->nrecv[cell];
1391         if (shift > 0)
1392         {
1393             /* Move the cg's present from previous grid pulses */
1394             cg0 = ncg_cell[ncell + cell];
1395             cg1 = ncg_cell[ncell + cell + 1];
1396             for (cg = cg1 - 1; cg >= cg0; cg--)
1397             {
1398                 index_gl[cg + shift] = index_gl[cg];
1399                 x[cg + shift]        = x[cg];
1400                 cginfo[cg + shift]   = cginfo[cg];
1401             }
1402             /* Correct the already stored send indices for the shift */
1403             for (p = 1; p <= pulse; p++)
1404             {
1405                 ind_p = &cd->ind[p];
1406                 cg0   = 0;
1407                 for (c = 0; c < cell; c++)
1408                 {
1409                     cg0 += ind_p->nsend[c];
1410                 }
1411                 cg1 = cg0 + ind_p->nsend[cell];
1412                 for (cg = cg0; cg < cg1; cg++)
1413                 {
1414                     ind_p->index[cg] += shift;
1415                 }
1416             }
1417         }
1418     }
1419
1420     /* Merge in the communicated buffers */
1421     shift = 0;
1422     cg0   = 0;
1423     for (cell = 0; cell < ncell; cell++)
1424     {
1425         cg1 = ncg_cell[ncell + cell + 1] + shift;
1426         for (cg = 0; cg < ind->nrecv[cell]; cg++)
1427         {
1428             /* Copy this atom from the buffer */
1429             index_gl[cg1] = recv_i[cg0];
1430             x[cg1]        = recv_vr[cg0];
1431             /* Copy information */
1432             cg_gl       = index_gl[cg1];
1433             cginfo[cg1] = ddcginfo(cginfo_mb, cg_gl);
1434             cg0++;
1435             cg1++;
1436         }
1437         shift += ind->nrecv[cell];
1438         ncg_cell[ncell + cell + 1] = cg1;
1439     }
1440 }
1441
1442 //! Makes a range partitioning for the atom groups wthin a cell
1443 static void make_cell2at_index(gmx_domdec_comm_dim_t* cd, int nzone, int atomGroupStart)
1444 {
1445     /* Store the atom block boundaries for easy copying of communication buffers
1446      */
1447     int g = atomGroupStart;
1448     for (int zone = 0; zone < nzone; zone++)
1449     {
1450         for (gmx_domdec_ind_t& ind : cd->ind)
1451         {
1452             ind.cell2at0[zone] = g;
1453             g += ind.nrecv[zone];
1454             ind.cell2at1[zone] = g;
1455         }
1456     }
1457 }
1458
1459 //! Returns whether a link is missing.
1460 static gmx_bool missing_link(const t_blocka& link, const int globalAtomIndex, const gmx_ga2la_t& ga2la)
1461 {
1462     for (int i = link.index[globalAtomIndex]; i < link.index[globalAtomIndex + 1]; i++)
1463     {
1464         if (!ga2la.findHome(link.a[i]))
1465         {
1466             return true;
1467         }
1468     }
1469
1470     return false;
1471 }
1472
1473 //! Domain corners for communication, a maximum of 4 i-zones see a j domain
1474 typedef struct
1475 {
1476     //! The corners for the non-bonded communication.
1477     real c[DIM][4];
1478     //! Corner for rounding.
1479     real cr0;
1480     //! Corners for rounding.
1481     real cr1[4];
1482     //! Corners for bounded communication.
1483     real bc[DIM];
1484     //! Corner for rounding for bonded communication.
1485     real bcr1;
1486 } dd_corners_t;
1487
1488 //! Determine the corners of the domain(s) we are communicating with.
1489 static void set_dd_corners(const gmx_domdec_t* dd, int dim0, int dim1, int dim2, gmx_bool bDistMB, dd_corners_t* c)
1490 {
1491     const gmx_domdec_comm_t*  comm;
1492     const gmx_domdec_zones_t* zones;
1493
1494     comm = dd->comm;
1495
1496     zones = &comm->zones;
1497
1498     /* Keep the compiler happy */
1499     c->cr0  = 0;
1500     c->bcr1 = 0;
1501
1502     /* The first dimension is equal for all cells */
1503     c->c[0][0] = comm->cell_x0[dim0];
1504     if (bDistMB)
1505     {
1506         c->bc[0] = c->c[0][0];
1507     }
1508     if (dd->ndim >= 2)
1509     {
1510         dim1 = dd->dim[1];
1511         /* This cell row is only seen from the first row */
1512         c->c[1][0] = comm->cell_x0[dim1];
1513         /* All rows can see this row */
1514         c->c[1][1] = comm->cell_x0[dim1];
1515         if (isDlbOn(dd->comm))
1516         {
1517             c->c[1][1] = std::max(comm->cell_x0[dim1], comm->zone_d1[1].mch0);
1518             if (bDistMB)
1519             {
1520                 /* For the multi-body distance we need the maximum */
1521                 c->bc[1] = std::max(comm->cell_x0[dim1], comm->zone_d1[1].p1_0);
1522             }
1523         }
1524         /* Set the upper-right corner for rounding */
1525         c->cr0 = comm->cell_x1[dim0];
1526
1527         if (dd->ndim >= 3)
1528         {
1529             dim2 = dd->dim[2];
1530             for (int j = 0; j < 4; j++)
1531             {
1532                 c->c[2][j] = comm->cell_x0[dim2];
1533             }
1534             if (isDlbOn(dd->comm))
1535             {
1536                 /* Use the maximum of the i-cells that see a j-cell */
1537                 for (const auto& iZone : zones->iZones)
1538                 {
1539                     const int iZoneIndex = iZone.iZoneIndex;
1540                     for (int jZone : iZone.jZoneRange)
1541                     {
1542                         if (jZone >= 4)
1543                         {
1544                             c->c[2][jZone - 4] = std::max(
1545                                     c->c[2][jZone - 4],
1546                                     comm->zone_d2[zones->shift[iZoneIndex][dim0]][zones->shift[iZoneIndex][dim1]]
1547                                             .mch0);
1548                         }
1549                     }
1550                 }
1551                 if (bDistMB)
1552                 {
1553                     /* For the multi-body distance we need the maximum */
1554                     c->bc[2] = comm->cell_x0[dim2];
1555                     for (int i = 0; i < 2; i++)
1556                     {
1557                         for (int j = 0; j < 2; j++)
1558                         {
1559                             c->bc[2] = std::max(c->bc[2], comm->zone_d2[i][j].p1_0);
1560                         }
1561                     }
1562                 }
1563             }
1564
1565             /* Set the upper-right corner for rounding */
1566             /* Cell (0,0,0) and cell (1,0,0) can see cell 4 (0,1,1)
1567              * Only cell (0,0,0) can see cell 7 (1,1,1)
1568              */
1569             c->cr1[0] = comm->cell_x1[dim1];
1570             c->cr1[3] = comm->cell_x1[dim1];
1571             if (isDlbOn(dd->comm))
1572             {
1573                 c->cr1[0] = std::max(comm->cell_x1[dim1], comm->zone_d1[1].mch1);
1574                 if (bDistMB)
1575                 {
1576                     /* For the multi-body distance we need the maximum */
1577                     c->bcr1 = std::max(comm->cell_x1[dim1], comm->zone_d1[1].p1_1);
1578                 }
1579             }
1580         }
1581     }
1582 }
1583
1584 /*! \brief Add the atom groups we need to send in this pulse from this
1585  * zone to \p localAtomGroups and \p work. */
1586 static void get_zone_pulse_cgs(gmx_domdec_t*            dd,
1587                                int                      zonei,
1588                                int                      zone,
1589                                int                      cg0,
1590                                int                      cg1,
1591                                gmx::ArrayRef<const int> globalAtomGroupIndices,
1592                                int                      dim,
1593                                int                      dim_ind,
1594                                int                      dim0,
1595                                int                      dim1,
1596                                int                      dim2,
1597                                real                     r_comm2,
1598                                real                     r_bcomm2,
1599                                matrix                   box,
1600                                bool                     distanceIsTriclinic,
1601                                rvec*                    normal,
1602                                real                     skew_fac2_d,
1603                                real                     skew_fac_01,
1604                                rvec*                    v_d,
1605                                rvec*                    v_0,
1606                                rvec*                    v_1,
1607                                const dd_corners_t*      c,
1608                                const rvec               sf2_round,
1609                                gmx_bool                 bDistBonded,
1610                                gmx_bool                 bBondComm,
1611                                gmx_bool                 bDist2B,
1612                                gmx_bool                 bDistMB,
1613                                rvec*                    cg_cm,
1614                                gmx::ArrayRef<const int> cginfo,
1615                                std::vector<int>*        localAtomGroups,
1616                                dd_comm_setup_work_t*    work)
1617 {
1618     gmx_domdec_comm_t* comm;
1619     gmx_bool           bScrew;
1620     gmx_bool           bDistMB_pulse;
1621     int                cg, i;
1622     real               r2, rb2, r, tric_sh;
1623     rvec               rn, rb;
1624     int                dimd;
1625     int                nsend_z, nat;
1626
1627     comm = dd->comm;
1628
1629     bScrew = (dd->unitCellInfo.haveScrewPBC && dim == XX);
1630
1631     bDistMB_pulse = (bDistMB && bDistBonded);
1632
1633     /* Unpack the work data */
1634     std::vector<int>&       ibuf = work->atomGroupBuffer;
1635     std::vector<gmx::RVec>& vbuf = work->positionBuffer;
1636     nsend_z                      = 0;
1637     nat                          = work->nat;
1638
1639     for (cg = cg0; cg < cg1; cg++)
1640     {
1641         r2  = 0;
1642         rb2 = 0;
1643         if (!distanceIsTriclinic)
1644         {
1645             /* Rectangular direction, easy */
1646             r = cg_cm[cg][dim] - c->c[dim_ind][zone];
1647             if (r > 0)
1648             {
1649                 r2 += r * r;
1650             }
1651             if (bDistMB_pulse)
1652             {
1653                 r = cg_cm[cg][dim] - c->bc[dim_ind];
1654                 if (r > 0)
1655                 {
1656                     rb2 += r * r;
1657                 }
1658             }
1659             /* Rounding gives at most a 16% reduction
1660              * in communicated atoms
1661              */
1662             if (dim_ind >= 1 && (zonei == 1 || zonei == 2))
1663             {
1664                 r = cg_cm[cg][dim0] - c->cr0;
1665                 /* This is the first dimension, so always r >= 0 */
1666                 r2 += r * r;
1667                 if (bDistMB_pulse)
1668                 {
1669                     rb2 += r * r;
1670                 }
1671             }
1672             if (dim_ind == 2 && (zonei == 2 || zonei == 3))
1673             {
1674                 r = cg_cm[cg][dim1] - c->cr1[zone];
1675                 if (r > 0)
1676                 {
1677                     r2 += r * r;
1678                 }
1679                 if (bDistMB_pulse)
1680                 {
1681                     r = cg_cm[cg][dim1] - c->bcr1;
1682                     if (r > 0)
1683                     {
1684                         rb2 += r * r;
1685                     }
1686                 }
1687             }
1688         }
1689         else
1690         {
1691             /* Triclinic direction, more complicated */
1692             clear_rvec(rn);
1693             clear_rvec(rb);
1694             /* Rounding, conservative as the skew_fac multiplication
1695              * will slightly underestimate the distance.
1696              */
1697             if (dim_ind >= 1 && (zonei == 1 || zonei == 2))
1698             {
1699                 rn[dim0] = cg_cm[cg][dim0] - c->cr0;
1700                 for (i = dim0 + 1; i < DIM; i++)
1701                 {
1702                     rn[dim0] -= cg_cm[cg][i] * v_0[i][dim0];
1703                 }
1704                 r2 = rn[dim0] * rn[dim0] * sf2_round[dim0];
1705                 if (bDistMB_pulse)
1706                 {
1707                     rb[dim0] = rn[dim0];
1708                     rb2      = r2;
1709                 }
1710                 /* Take care that the cell planes along dim0 might not
1711                  * be orthogonal to those along dim1 and dim2.
1712                  */
1713                 for (i = 1; i <= dim_ind; i++)
1714                 {
1715                     dimd = dd->dim[i];
1716                     if (normal[dim0][dimd] > 0)
1717                     {
1718                         rn[dimd] -= rn[dim0] * normal[dim0][dimd];
1719                         if (bDistMB_pulse)
1720                         {
1721                             rb[dimd] -= rb[dim0] * normal[dim0][dimd];
1722                         }
1723                     }
1724                 }
1725             }
1726             if (dim_ind == 2 && (zonei == 2 || zonei == 3))
1727             {
1728                 GMX_ASSERT(dim1 >= 0 && dim1 < DIM, "Must have a valid dimension index");
1729                 rn[dim1] += cg_cm[cg][dim1] - c->cr1[zone];
1730                 tric_sh = 0;
1731                 for (i = dim1 + 1; i < DIM; i++)
1732                 {
1733                     tric_sh -= cg_cm[cg][i] * v_1[i][dim1];
1734                 }
1735                 rn[dim1] += tric_sh;
1736                 if (rn[dim1] > 0)
1737                 {
1738                     r2 += rn[dim1] * rn[dim1] * sf2_round[dim1];
1739                     /* Take care of coupling of the distances
1740                      * to the planes along dim0 and dim1 through dim2.
1741                      */
1742                     r2 -= rn[dim0] * rn[dim1] * skew_fac_01;
1743                     /* Take care that the cell planes along dim1
1744                      * might not be orthogonal to that along dim2.
1745                      */
1746                     if (normal[dim1][dim2] > 0)
1747                     {
1748                         rn[dim2] -= rn[dim1] * normal[dim1][dim2];
1749                     }
1750                 }
1751                 if (bDistMB_pulse)
1752                 {
1753                     rb[dim1] += cg_cm[cg][dim1] - c->bcr1 + tric_sh;
1754                     if (rb[dim1] > 0)
1755                     {
1756                         rb2 += rb[dim1] * rb[dim1] * sf2_round[dim1];
1757                         /* Take care of coupling of the distances
1758                          * to the planes along dim0 and dim1 through dim2.
1759                          */
1760                         rb2 -= rb[dim0] * rb[dim1] * skew_fac_01;
1761                         /* Take care that the cell planes along dim1
1762                          * might not be orthogonal to that along dim2.
1763                          */
1764                         if (normal[dim1][dim2] > 0)
1765                         {
1766                             rb[dim2] -= rb[dim1] * normal[dim1][dim2];
1767                         }
1768                     }
1769                 }
1770             }
1771             /* The distance along the communication direction */
1772             rn[dim] += cg_cm[cg][dim] - c->c[dim_ind][zone];
1773             tric_sh = 0;
1774             for (i = dim + 1; i < DIM; i++)
1775             {
1776                 tric_sh -= cg_cm[cg][i] * v_d[i][dim];
1777             }
1778             rn[dim] += tric_sh;
1779             if (rn[dim] > 0)
1780             {
1781                 r2 += rn[dim] * rn[dim] * skew_fac2_d;
1782                 /* Take care of coupling of the distances
1783                  * to the planes along dim0 and dim1 through dim2.
1784                  */
1785                 if (dim_ind == 1 && zonei == 1)
1786                 {
1787                     r2 -= rn[dim0] * rn[dim] * skew_fac_01;
1788                 }
1789             }
1790             if (bDistMB_pulse)
1791             {
1792                 clear_rvec(rb);
1793                 GMX_ASSERT(dim >= 0 && dim < DIM, "Must have a valid dimension index");
1794                 rb[dim] += cg_cm[cg][dim] - c->bc[dim_ind] + tric_sh;
1795                 if (rb[dim] > 0)
1796                 {
1797                     rb2 += rb[dim] * rb[dim] * skew_fac2_d;
1798                     /* Take care of coupling of the distances
1799                      * to the planes along dim0 and dim1 through dim2.
1800                      */
1801                     if (dim_ind == 1 && zonei == 1)
1802                     {
1803                         rb2 -= rb[dim0] * rb[dim] * skew_fac_01;
1804                     }
1805                 }
1806             }
1807         }
1808
1809         if (r2 < r_comm2
1810             || (bDistBonded && ((bDistMB && rb2 < r_bcomm2) || (bDist2B && r2 < r_bcomm2))
1811                 && (!bBondComm
1812                     || (GET_CGINFO_BOND_INTER(cginfo[cg])
1813                         && missing_link(*comm->bondedLinks, globalAtomGroupIndices[cg], *dd->ga2la)))))
1814         {
1815             /* Store the local and global atom group indices and position */
1816             localAtomGroups->push_back(cg);
1817             ibuf.push_back(globalAtomGroupIndices[cg]);
1818             nsend_z++;
1819
1820             rvec posPbc;
1821             if (dd->ci[dim] == 0)
1822             {
1823                 /* Correct cg_cm for pbc */
1824                 rvec_add(cg_cm[cg], box[dim], posPbc);
1825                 if (bScrew)
1826                 {
1827                     posPbc[YY] = box[YY][YY] - posPbc[YY];
1828                     posPbc[ZZ] = box[ZZ][ZZ] - posPbc[ZZ];
1829                 }
1830             }
1831             else
1832             {
1833                 copy_rvec(cg_cm[cg], posPbc);
1834             }
1835             vbuf.emplace_back(posPbc[XX], posPbc[YY], posPbc[ZZ]);
1836
1837             nat += 1;
1838         }
1839     }
1840
1841     work->nat        = nat;
1842     work->nsend_zone = nsend_z;
1843 }
1844
1845 //! Clear data.
1846 static void clearCommSetupData(dd_comm_setup_work_t* work)
1847 {
1848     work->localAtomGroupBuffer.clear();
1849     work->atomGroupBuffer.clear();
1850     work->positionBuffer.clear();
1851     work->nat        = 0;
1852     work->nsend_zone = 0;
1853 }
1854
1855 //! Prepare DD communication.
1856 static void setup_dd_communication(gmx_domdec_t* dd, matrix box, gmx_ddbox_t* ddbox, t_forcerec* fr, t_state* state)
1857 {
1858     int                    dim_ind, dim, dim0, dim1, dim2, dimd, nat_tot;
1859     int                    nzone, nzone_send, zone, zonei, cg0, cg1;
1860     int                    c;
1861     int *                  zone_cg_range, pos_cg;
1862     gmx_domdec_comm_t*     comm;
1863     gmx_domdec_zones_t*    zones;
1864     gmx_domdec_comm_dim_t* cd;
1865     gmx_bool               bBondComm, bDist2B, bDistMB, bDistBonded;
1866     dd_corners_t           corners;
1867     rvec *                 normal, *v_d, *v_0 = nullptr, *v_1 = nullptr;
1868     real                   skew_fac2_d, skew_fac_01;
1869     rvec                   sf2_round;
1870
1871     if (debug)
1872     {
1873         fprintf(debug, "Setting up DD communication\n");
1874     }
1875
1876     comm = dd->comm;
1877
1878     if (comm->dth.empty())
1879     {
1880         /* Initialize the thread data.
1881          * This can not be done in init_domain_decomposition,
1882          * as the numbers of threads is determined later.
1883          */
1884         int numThreads = gmx_omp_nthreads_get(emntDomdec);
1885         comm->dth.resize(numThreads);
1886     }
1887
1888     bBondComm = comm->systemInfo.filterBondedCommunication;
1889
1890     /* Do we need to determine extra distances for multi-body bondeds? */
1891     bDistMB = (comm->systemInfo.haveInterDomainMultiBodyBondeds && isDlbOn(dd->comm) && dd->ndim > 1);
1892
1893     /* Do we need to determine extra distances for only two-body bondeds? */
1894     bDist2B = (bBondComm && !bDistMB);
1895
1896     const real r_comm2 =
1897             gmx::square(domainToDomainIntoAtomToDomainCutoff(comm->systemInfo, comm->systemInfo.cutoff));
1898     const real r_bcomm2 =
1899             gmx::square(domainToDomainIntoAtomToDomainCutoff(comm->systemInfo, comm->cutoff_mbody));
1900
1901     if (debug)
1902     {
1903         fprintf(debug, "bBondComm %s, r_bc %f\n", gmx::boolToString(bBondComm), std::sqrt(r_bcomm2));
1904     }
1905
1906     zones = &comm->zones;
1907
1908     dim0 = dd->dim[0];
1909     dim1 = (dd->ndim >= 2 ? dd->dim[1] : -1);
1910     dim2 = (dd->ndim >= 3 ? dd->dim[2] : -1);
1911
1912     set_dd_corners(dd, dim0, dim1, dim2, bDistMB, &corners);
1913
1914     /* Triclinic stuff */
1915     normal      = ddbox->normal;
1916     skew_fac_01 = 0;
1917     if (dd->ndim >= 2)
1918     {
1919         v_0 = ddbox->v[dim0];
1920         if (ddbox->tric_dir[dim0] && ddbox->tric_dir[dim1])
1921         {
1922             /* Determine the coupling coefficient for the distances
1923              * to the cell planes along dim0 and dim1 through dim2.
1924              * This is required for correct rounding.
1925              */
1926             skew_fac_01 = ddbox->v[dim0][dim1 + 1][dim0] * ddbox->v[dim1][dim1 + 1][dim1];
1927             if (debug)
1928             {
1929                 fprintf(debug, "\nskew_fac_01 %f\n", skew_fac_01);
1930             }
1931         }
1932     }
1933     if (dd->ndim >= 3)
1934     {
1935         v_1 = ddbox->v[dim1];
1936     }
1937
1938     zone_cg_range                        = zones->cg_range;
1939     gmx::ArrayRef<cginfo_mb_t> cginfo_mb = fr->cginfo_mb;
1940
1941     zone_cg_range[0]   = 0;
1942     zone_cg_range[1]   = dd->ncg_home;
1943     comm->zone_ncg1[0] = dd->ncg_home;
1944     pos_cg             = dd->ncg_home;
1945
1946     nat_tot = comm->atomRanges.numHomeAtoms();
1947     nzone   = 1;
1948     for (dim_ind = 0; dim_ind < dd->ndim; dim_ind++)
1949     {
1950         dim = dd->dim[dim_ind];
1951         cd  = &comm->cd[dim_ind];
1952
1953         /* Check if we need to compute triclinic distances along this dim */
1954         bool distanceIsTriclinic = false;
1955         for (int i = 0; i <= dim_ind; i++)
1956         {
1957             if (ddbox->tric_dir[dd->dim[i]])
1958             {
1959                 distanceIsTriclinic = true;
1960             }
1961         }
1962
1963         if (dim >= ddbox->npbcdim && dd->ci[dim] == 0)
1964         {
1965             /* No pbc in this dimension, the first node should not comm. */
1966             nzone_send = 0;
1967         }
1968         else
1969         {
1970             nzone_send = nzone;
1971         }
1972
1973         v_d         = ddbox->v[dim];
1974         skew_fac2_d = gmx::square(ddbox->skew_fac[dim]);
1975
1976         cd->receiveInPlace = true;
1977         for (int p = 0; p < cd->numPulses(); p++)
1978         {
1979             /* Only atoms communicated in the first pulse are used
1980              * for multi-body bonded interactions or for bBondComm.
1981              */
1982             bDistBonded = ((bDistMB || bDist2B) && p == 0);
1983
1984             gmx_domdec_ind_t* ind = &cd->ind[p];
1985
1986             /* Thread 0 writes in the global index array */
1987             ind->index.clear();
1988             clearCommSetupData(&comm->dth[0]);
1989
1990             for (zone = 0; zone < nzone_send; zone++)
1991             {
1992                 if (dim_ind > 0 && distanceIsTriclinic)
1993                 {
1994                     /* Determine slightly more optimized skew_fac's
1995                      * for rounding.
1996                      * This reduces the number of communicated atoms
1997                      * by about 10% for 3D DD of rhombic dodecahedra.
1998                      */
1999                     for (dimd = 0; dimd < dim; dimd++)
2000                     {
2001                         sf2_round[dimd] = 1;
2002                         if (ddbox->tric_dir[dimd])
2003                         {
2004                             for (int i = dd->dim[dimd] + 1; i < DIM; i++)
2005                             {
2006                                 /* If we are shifted in dimension i
2007                                  * and the cell plane is tilted forward
2008                                  * in dimension i, skip this coupling.
2009                                  */
2010                                 if (!(zones->shift[nzone + zone][i] && ddbox->v[dimd][i][dimd] >= 0))
2011                                 {
2012                                     sf2_round[dimd] += gmx::square(ddbox->v[dimd][i][dimd]);
2013                                 }
2014                             }
2015                             sf2_round[dimd] = 1 / sf2_round[dimd];
2016                         }
2017                     }
2018                 }
2019
2020                 zonei = zone_perm[dim_ind][zone];
2021                 if (p == 0)
2022                 {
2023                     /* Here we permutate the zones to obtain a convenient order
2024                      * for neighbor searching
2025                      */
2026                     cg0 = zone_cg_range[zonei];
2027                     cg1 = zone_cg_range[zonei + 1];
2028                 }
2029                 else
2030                 {
2031                     /* Look only at the cg's received in the previous grid pulse
2032                      */
2033                     cg1 = zone_cg_range[nzone + zone + 1];
2034                     cg0 = cg1 - cd->ind[p - 1].nrecv[zone];
2035                 }
2036
2037                 const int numThreads = gmx::ssize(comm->dth);
2038 #pragma omp parallel for num_threads(numThreads) schedule(static)
2039                 for (int th = 0; th < numThreads; th++)
2040                 {
2041                     try
2042                     {
2043                         dd_comm_setup_work_t& work = comm->dth[th];
2044
2045                         /* Retain data accumulated into buffers of thread 0 */
2046                         if (th > 0)
2047                         {
2048                             clearCommSetupData(&work);
2049                         }
2050
2051                         int cg0_th = cg0 + ((cg1 - cg0) * th) / numThreads;
2052                         int cg1_th = cg0 + ((cg1 - cg0) * (th + 1)) / numThreads;
2053
2054                         /* Get the cg's for this pulse in this zone */
2055                         get_zone_pulse_cgs(dd,
2056                                            zonei,
2057                                            zone,
2058                                            cg0_th,
2059                                            cg1_th,
2060                                            dd->globalAtomGroupIndices,
2061                                            dim,
2062                                            dim_ind,
2063                                            dim0,
2064                                            dim1,
2065                                            dim2,
2066                                            r_comm2,
2067                                            r_bcomm2,
2068                                            box,
2069                                            distanceIsTriclinic,
2070                                            normal,
2071                                            skew_fac2_d,
2072                                            skew_fac_01,
2073                                            v_d,
2074                                            v_0,
2075                                            v_1,
2076                                            &corners,
2077                                            sf2_round,
2078                                            bDistBonded,
2079                                            bBondComm,
2080                                            bDist2B,
2081                                            bDistMB,
2082                                            state->x.rvec_array(),
2083                                            fr->cginfo,
2084                                            th == 0 ? &ind->index : &work.localAtomGroupBuffer,
2085                                            &work);
2086                     }
2087                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2088                 } // END
2089
2090                 std::vector<int>&       atomGroups = comm->dth[0].atomGroupBuffer;
2091                 std::vector<gmx::RVec>& positions  = comm->dth[0].positionBuffer;
2092                 ind->nsend[zone]                   = comm->dth[0].nsend_zone;
2093                 /* Append data of threads>=1 to the communication buffers */
2094                 for (int th = 1; th < numThreads; th++)
2095                 {
2096                     const dd_comm_setup_work_t& dth = comm->dth[th];
2097
2098                     ind->index.insert(ind->index.end(),
2099                                       dth.localAtomGroupBuffer.begin(),
2100                                       dth.localAtomGroupBuffer.end());
2101                     atomGroups.insert(
2102                             atomGroups.end(), dth.atomGroupBuffer.begin(), dth.atomGroupBuffer.end());
2103                     positions.insert(
2104                             positions.end(), dth.positionBuffer.begin(), dth.positionBuffer.end());
2105                     comm->dth[0].nat += dth.nat;
2106                     ind->nsend[zone] += dth.nsend_zone;
2107                 }
2108             }
2109             /* Clear the counts in case we do not have pbc */
2110             for (zone = nzone_send; zone < nzone; zone++)
2111             {
2112                 ind->nsend[zone] = 0;
2113             }
2114             ind->nsend[nzone]     = ind->index.size();
2115             ind->nsend[nzone + 1] = comm->dth[0].nat;
2116             /* Communicate the number of cg's and atoms to receive */
2117             ddSendrecv(dd, dim_ind, dddirBackward, ind->nsend, nzone + 2, ind->nrecv, nzone + 2);
2118
2119             if (p > 0)
2120             {
2121                 /* We can receive in place if only the last zone is not empty */
2122                 for (zone = 0; zone < nzone - 1; zone++)
2123                 {
2124                     if (ind->nrecv[zone] > 0)
2125                     {
2126                         cd->receiveInPlace = false;
2127                     }
2128                 }
2129             }
2130
2131             int receiveBufferSize = 0;
2132             if (!cd->receiveInPlace)
2133             {
2134                 receiveBufferSize = ind->nrecv[nzone];
2135             }
2136             /* These buffer are actually only needed with in-place */
2137             DDBufferAccess<int>       globalAtomGroupBuffer(comm->intBuffer, receiveBufferSize);
2138             DDBufferAccess<gmx::RVec> rvecBuffer(comm->rvecBuffer, receiveBufferSize);
2139
2140             dd_comm_setup_work_t& work = comm->dth[0];
2141
2142             /* Make space for the global cg indices */
2143             int numAtomGroupsNew = pos_cg + ind->nrecv[nzone];
2144             dd->globalAtomGroupIndices.resize(numAtomGroupsNew);
2145             /* Communicate the global cg indices */
2146             gmx::ArrayRef<int> integerBufferRef;
2147             if (cd->receiveInPlace)
2148             {
2149                 integerBufferRef = gmx::arrayRefFromArray(
2150                         dd->globalAtomGroupIndices.data() + pos_cg, ind->nrecv[nzone]);
2151             }
2152             else
2153             {
2154                 integerBufferRef = globalAtomGroupBuffer.buffer;
2155             }
2156             ddSendrecv<int>(dd, dim_ind, dddirBackward, work.atomGroupBuffer, integerBufferRef);
2157
2158             /* Make space for cg_cm */
2159             dd_resize_atominfo_and_state(fr, state, pos_cg + ind->nrecv[nzone]);
2160
2161             /* Communicate the coordinates */
2162             gmx::ArrayRef<gmx::RVec> rvecBufferRef;
2163             if (cd->receiveInPlace)
2164             {
2165                 rvecBufferRef = gmx::makeArrayRef(state->x).subArray(pos_cg, ind->nrecv[nzone]);
2166             }
2167             else
2168             {
2169                 rvecBufferRef = rvecBuffer.buffer;
2170             }
2171             ddSendrecv<gmx::RVec>(dd, dim_ind, dddirBackward, work.positionBuffer, rvecBufferRef);
2172
2173             /* Make the charge group index */
2174             if (cd->receiveInPlace)
2175             {
2176                 zone = (p == 0 ? 0 : nzone - 1);
2177                 while (zone < nzone)
2178                 {
2179                     for (int i = 0; i < ind->nrecv[zone]; i++)
2180                     {
2181                         int globalAtomIndex = dd->globalAtomGroupIndices[pos_cg];
2182                         fr->cginfo[pos_cg]  = ddcginfo(cginfo_mb, globalAtomIndex);
2183                         pos_cg++;
2184                     }
2185                     if (p == 0)
2186                     {
2187                         comm->zone_ncg1[nzone + zone] = ind->nrecv[zone];
2188                     }
2189                     zone++;
2190                     zone_cg_range[nzone + zone] = pos_cg;
2191                 }
2192             }
2193             else
2194             {
2195                 /* This part of the code is never executed with bBondComm. */
2196                 merge_cg_buffers(nzone,
2197                                  cd,
2198                                  p,
2199                                  zone_cg_range,
2200                                  dd->globalAtomGroupIndices,
2201                                  integerBufferRef.data(),
2202                                  state->x,
2203                                  rvecBufferRef,
2204                                  fr->cginfo_mb,
2205                                  fr->cginfo);
2206                 pos_cg += ind->nrecv[nzone];
2207             }
2208             nat_tot += ind->nrecv[nzone + 1];
2209         }
2210         if (!cd->receiveInPlace)
2211         {
2212             /* Store the atom block for easy copying of communication buffers */
2213             make_cell2at_index(cd, nzone, zone_cg_range[nzone]);
2214         }
2215         nzone += nzone;
2216     }
2217
2218     comm->atomRanges.setEnd(DDAtomRanges::Type::Zones, nat_tot);
2219
2220     if (!bBondComm)
2221     {
2222         /* We don't need to update cginfo, since that was alrady done above.
2223          * So we pass NULL for the forcerec.
2224          */
2225         dd_set_cginfo(dd->globalAtomGroupIndices, dd->ncg_home, dd->globalAtomGroupIndices.size(), nullptr);
2226     }
2227
2228     if (debug)
2229     {
2230         fprintf(debug, "Finished setting up DD communication, zones:");
2231         for (c = 0; c < zones->n; c++)
2232         {
2233             fprintf(debug, " %d", zones->cg_range[c + 1] - zones->cg_range[c]);
2234         }
2235         fprintf(debug, "\n");
2236     }
2237 }
2238
2239 //! Set boundaries for the charge group range.
2240 static void set_cg_boundaries(gmx_domdec_zones_t* zones)
2241 {
2242     for (auto& iZone : zones->iZones)
2243     {
2244         iZone.iAtomRange = gmx::Range<int>(0, zones->cg_range[iZone.iZoneIndex + 1]);
2245         iZone.jAtomRange = gmx::Range<int>(zones->cg_range[iZone.jZoneRange.begin()],
2246                                            zones->cg_range[iZone.jZoneRange.end()]);
2247     }
2248 }
2249
2250 /*! \brief Set zone dimensions for zones \p zone_start to \p zone_end-1
2251  *
2252  * Also sets the atom density for the home zone when \p zone_start=0.
2253  * For this \p numMovedChargeGroupsInHomeZone needs to be passed to tell
2254  * how many charge groups will move but are still part of the current range.
2255  * \todo When converting domdec to use proper classes, all these variables
2256  *       should be private and a method should return the correct count
2257  *       depending on an internal state.
2258  *
2259  * \param[in,out] dd          The domain decomposition struct
2260  * \param[in]     box         The box
2261  * \param[in]     ddbox       The domain decomposition box struct
2262  * \param[in]     zone_start  The start of the zone range to set sizes for
2263  * \param[in]     zone_end    The end of the zone range to set sizes for
2264  * \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
2265  */
2266 static void set_zones_size(gmx_domdec_t*      dd,
2267                            matrix             box,
2268                            const gmx_ddbox_t* ddbox,
2269                            int                zone_start,
2270                            int                zone_end,
2271                            int                numMovedChargeGroupsInHomeZone)
2272 {
2273     gmx_domdec_comm_t*  comm;
2274     gmx_domdec_zones_t* zones;
2275     gmx_bool            bDistMB;
2276     int                 z, d, dim;
2277     real                rcs, rcmbs;
2278     int                 i, j;
2279     real                vol;
2280
2281     comm = dd->comm;
2282
2283     zones = &comm->zones;
2284
2285     /* Do we need to determine extra distances for multi-body bondeds? */
2286     bDistMB = (comm->systemInfo.haveInterDomainMultiBodyBondeds && isDlbOn(dd->comm) && dd->ndim > 1);
2287
2288     for (z = zone_start; z < zone_end; z++)
2289     {
2290         /* Copy cell limits to zone limits.
2291          * Valid for non-DD dims and non-shifted dims.
2292          */
2293         copy_rvec(comm->cell_x0, zones->size[z].x0);
2294         copy_rvec(comm->cell_x1, zones->size[z].x1);
2295     }
2296
2297     for (d = 0; d < dd->ndim; d++)
2298     {
2299         dim = dd->dim[d];
2300
2301         for (z = 0; z < zones->n; z++)
2302         {
2303             /* With a staggered grid we have different sizes
2304              * for non-shifted dimensions.
2305              */
2306             if (isDlbOn(dd->comm) && zones->shift[z][dim] == 0)
2307             {
2308                 if (d == 1)
2309                 {
2310                     zones->size[z].x0[dim] = comm->zone_d1[zones->shift[z][dd->dim[d - 1]]].min0;
2311                     zones->size[z].x1[dim] = comm->zone_d1[zones->shift[z][dd->dim[d - 1]]].max1;
2312                 }
2313                 else if (d == 2)
2314                 {
2315                     zones->size[z].x0[dim] =
2316                             comm->zone_d2[zones->shift[z][dd->dim[d - 2]]][zones->shift[z][dd->dim[d - 1]]]
2317                                     .min0;
2318                     zones->size[z].x1[dim] =
2319                             comm->zone_d2[zones->shift[z][dd->dim[d - 2]]][zones->shift[z][dd->dim[d - 1]]]
2320                                     .max1;
2321                 }
2322             }
2323         }
2324
2325         rcs   = comm->systemInfo.cutoff;
2326         rcmbs = comm->cutoff_mbody;
2327         if (ddbox->tric_dir[dim])
2328         {
2329             rcs /= ddbox->skew_fac[dim];
2330             rcmbs /= ddbox->skew_fac[dim];
2331         }
2332
2333         /* Set the lower limit for the shifted zone dimensions */
2334         for (z = zone_start; z < zone_end; z++)
2335         {
2336             if (zones->shift[z][dim] > 0)
2337             {
2338                 dim = dd->dim[d];
2339                 if (!isDlbOn(dd->comm) || d == 0)
2340                 {
2341                     zones->size[z].x0[dim] = comm->cell_x1[dim];
2342                     zones->size[z].x1[dim] = comm->cell_x1[dim] + rcs;
2343                 }
2344                 else
2345                 {
2346                     /* Here we take the lower limit of the zone from
2347                      * the lowest domain of the zone below.
2348                      */
2349                     if (z < 4)
2350                     {
2351                         zones->size[z].x0[dim] = comm->zone_d1[zones->shift[z][dd->dim[d - 1]]].min1;
2352                     }
2353                     else
2354                     {
2355                         if (d == 1)
2356                         {
2357                             zones->size[z].x0[dim] = zones->size[zone_perm[2][z - 4]].x0[dim];
2358                         }
2359                         else
2360                         {
2361                             zones->size[z].x0[dim] =
2362                                     comm->zone_d2[zones->shift[z][dd->dim[d - 2]]][zones->shift[z][dd->dim[d - 1]]]
2363                                             .min1;
2364                         }
2365                     }
2366                     /* A temporary limit, is updated below */
2367                     zones->size[z].x1[dim] = zones->size[z].x0[dim];
2368
2369                     if (bDistMB)
2370                     {
2371                         for (size_t zi = 0; zi < zones->iZones.size(); zi++)
2372                         {
2373                             if (zones->shift[zi][dim] == 0)
2374                             {
2375                                 /* This takes the whole zone into account.
2376                                  * With multiple pulses this will lead
2377                                  * to a larger zone then strictly necessary.
2378                                  */
2379                                 zones->size[z].x1[dim] = std::max(zones->size[z].x1[dim],
2380                                                                   zones->size[zi].x1[dim] + rcmbs);
2381                             }
2382                         }
2383                     }
2384                 }
2385             }
2386         }
2387
2388         /* Loop over the i-zones to set the upper limit of each
2389          * j-zone they see.
2390          */
2391         for (const auto& iZone : zones->iZones)
2392         {
2393             const int zi = iZone.iZoneIndex;
2394             if (zones->shift[zi][dim] == 0)
2395             {
2396                 /* We should only use zones up to zone_end */
2397                 const auto& jZoneRangeFull = iZone.jZoneRange;
2398                 if (zone_end <= *jZoneRangeFull.begin())
2399                 {
2400                     continue;
2401                 }
2402                 const gmx::Range<int> jZoneRange(*jZoneRangeFull.begin(),
2403                                                  std::min(*jZoneRangeFull.end(), zone_end));
2404                 for (int jZone : jZoneRange)
2405                 {
2406                     if (zones->shift[jZone][dim] > 0)
2407                     {
2408                         zones->size[jZone].x1[dim] =
2409                                 std::max(zones->size[jZone].x1[dim], zones->size[zi].x1[dim] + rcs);
2410                     }
2411                 }
2412             }
2413         }
2414     }
2415
2416     for (z = zone_start; z < zone_end; z++)
2417     {
2418         /* Initialization only required to keep the compiler happy */
2419         rvec corner_min = { 0, 0, 0 }, corner_max = { 0, 0, 0 }, corner;
2420         int  nc, c;
2421
2422         /* To determine the bounding box for a zone we need to find
2423          * the extreme corners of 4, 2 or 1 corners.
2424          */
2425         nc = 1 << (ddbox->nboundeddim - 1);
2426
2427         for (c = 0; c < nc; c++)
2428         {
2429             /* Set up a zone corner at x=0, ignoring trilinic couplings */
2430             corner[XX] = 0;
2431             if ((c & 1) == 0)
2432             {
2433                 corner[YY] = zones->size[z].x0[YY];
2434             }
2435             else
2436             {
2437                 corner[YY] = zones->size[z].x1[YY];
2438             }
2439             if ((c & 2) == 0)
2440             {
2441                 corner[ZZ] = zones->size[z].x0[ZZ];
2442             }
2443             else
2444             {
2445                 corner[ZZ] = zones->size[z].x1[ZZ];
2446             }
2447             if (dd->ndim == 1 && dd->dim[0] < ZZ && ZZ < dd->unitCellInfo.npbcdim
2448                 && box[ZZ][1 - dd->dim[0]] != 0)
2449             {
2450                 /* With 1D domain decomposition the cg's are not in
2451                  * the triclinic box, but triclinic x-y and rectangular y/x-z.
2452                  * Shift the corner of the z-vector back to along the box
2453                  * vector of dimension d, so it will later end up at 0 along d.
2454                  * This can affect the location of this corner along dd->dim[0]
2455                  * through the matrix operation below if box[d][dd->dim[0]]!=0.
2456                  */
2457                 int d = 1 - dd->dim[0];
2458
2459                 corner[d] -= corner[ZZ] * box[ZZ][d] / box[ZZ][ZZ];
2460             }
2461             /* Apply the triclinic couplings */
2462             for (i = YY; i < ddbox->npbcdim && i < DIM; i++)
2463             {
2464                 for (j = XX; j < i; j++)
2465                 {
2466                     corner[j] += corner[i] * box[i][j] / box[i][i];
2467                 }
2468             }
2469             if (c == 0)
2470             {
2471                 copy_rvec(corner, corner_min);
2472                 copy_rvec(corner, corner_max);
2473             }
2474             else
2475             {
2476                 for (i = 0; i < DIM; i++)
2477                 {
2478                     corner_min[i] = std::min(corner_min[i], corner[i]);
2479                     corner_max[i] = std::max(corner_max[i], corner[i]);
2480                 }
2481             }
2482         }
2483         /* Copy the extreme cornes without offset along x */
2484         for (i = 0; i < DIM; i++)
2485         {
2486             zones->size[z].bb_x0[i] = corner_min[i];
2487             zones->size[z].bb_x1[i] = corner_max[i];
2488         }
2489         /* Add the offset along x */
2490         zones->size[z].bb_x0[XX] += zones->size[z].x0[XX];
2491         zones->size[z].bb_x1[XX] += zones->size[z].x1[XX];
2492     }
2493
2494     if (zone_start == 0)
2495     {
2496         vol = 1;
2497         for (dim = 0; dim < DIM; dim++)
2498         {
2499             vol *= zones->size[0].x1[dim] - zones->size[0].x0[dim];
2500         }
2501         zones->dens_zone0 =
2502                 (zones->cg_range[1] - zones->cg_range[0] - numMovedChargeGroupsInHomeZone) / vol;
2503     }
2504
2505     if (debug)
2506     {
2507         for (z = zone_start; z < zone_end; z++)
2508         {
2509             fprintf(debug,
2510                     "zone %d    %6.3f - %6.3f  %6.3f - %6.3f  %6.3f - %6.3f\n",
2511                     z,
2512                     zones->size[z].x0[XX],
2513                     zones->size[z].x1[XX],
2514                     zones->size[z].x0[YY],
2515                     zones->size[z].x1[YY],
2516                     zones->size[z].x0[ZZ],
2517                     zones->size[z].x1[ZZ]);
2518             fprintf(debug,
2519                     "zone %d bb %6.3f - %6.3f  %6.3f - %6.3f  %6.3f - %6.3f\n",
2520                     z,
2521                     zones->size[z].bb_x0[XX],
2522                     zones->size[z].bb_x1[XX],
2523                     zones->size[z].bb_x0[YY],
2524                     zones->size[z].bb_x1[YY],
2525                     zones->size[z].bb_x0[ZZ],
2526                     zones->size[z].bb_x1[ZZ]);
2527         }
2528     }
2529 }
2530
2531 /*! \brief Order data in \p dataToSort according to \p sort
2532  *
2533  * Note: both buffers should have at least \p sort.size() elements.
2534  */
2535 template<typename T>
2536 static void orderVector(gmx::ArrayRef<const gmx_cgsort_t> sort,
2537                         gmx::ArrayRef<T>                  dataToSort,
2538                         gmx::ArrayRef<T>                  sortBuffer)
2539 {
2540     GMX_ASSERT(dataToSort.size() >= sort.size(), "The vector needs to be sufficiently large");
2541     GMX_ASSERT(sortBuffer.size() >= sort.size(),
2542                "The sorting buffer needs to be sufficiently large");
2543
2544     /* Order the data into the temporary buffer */
2545     size_t i = 0;
2546     for (const gmx_cgsort_t& entry : sort)
2547     {
2548         sortBuffer[i++] = dataToSort[entry.ind];
2549     }
2550
2551     /* Copy back to the original array */
2552     std::copy(sortBuffer.begin(), sortBuffer.begin() + sort.size(), dataToSort.begin());
2553 }
2554
2555 /*! \brief Order data in \p dataToSort according to \p sort
2556  *
2557  * Note: \p vectorToSort should have at least \p sort.size() elements,
2558  *       \p workVector is resized when it is too small.
2559  */
2560 template<typename T>
2561 static void orderVector(gmx::ArrayRef<const gmx_cgsort_t> sort,
2562                         gmx::ArrayRef<T>                  vectorToSort,
2563                         std::vector<T>*                   workVector)
2564 {
2565     if (gmx::index(workVector->size()) < sort.ssize())
2566     {
2567         workVector->resize(sort.size());
2568     }
2569     orderVector<T>(sort, vectorToSort, *workVector);
2570 }
2571
2572 //! Returns the sorting order for atoms based on the nbnxn grid order in sort
2573 static void dd_sort_order_nbnxn(const t_forcerec* fr, std::vector<gmx_cgsort_t>* sort)
2574 {
2575     gmx::ArrayRef<const int> atomOrder = fr->nbv->getLocalAtomOrder();
2576
2577     /* Using push_back() instead of this resize results in much slower code */
2578     sort->resize(atomOrder.size());
2579     gmx::ArrayRef<gmx_cgsort_t> buffer    = *sort;
2580     size_t                      numSorted = 0;
2581     for (int i : atomOrder)
2582     {
2583         if (i >= 0)
2584         {
2585             /* The values of nsc and ind_gl are not used in this case */
2586             buffer[numSorted++].ind = i;
2587         }
2588     }
2589     sort->resize(numSorted);
2590 }
2591
2592 //! Returns the sorting state for DD.
2593 static void dd_sort_state(gmx_domdec_t* dd, t_forcerec* fr, t_state* state)
2594 {
2595     gmx_domdec_sort_t* sort = dd->comm->sort.get();
2596
2597     dd_sort_order_nbnxn(fr, &sort->sorted);
2598
2599     /* We alloc with the old size, since cgindex is still old */
2600     DDBufferAccess<gmx::RVec> rvecBuffer(dd->comm->rvecBuffer, dd->ncg_home);
2601
2602     /* Set the new home atom/charge group count */
2603     dd->ncg_home = sort->sorted.size();
2604     if (debug)
2605     {
2606         fprintf(debug, "Set the new home atom count to %d\n", dd->ncg_home);
2607     }
2608
2609     /* Reorder the state */
2610     gmx::ArrayRef<const gmx_cgsort_t> cgsort = sort->sorted;
2611     GMX_RELEASE_ASSERT(cgsort.ssize() == dd->ncg_home, "We should sort all the home atom groups");
2612
2613     if (state->flags & (1 << estX))
2614     {
2615         orderVector(cgsort, makeArrayRef(state->x), rvecBuffer.buffer);
2616     }
2617     if (state->flags & (1 << estV))
2618     {
2619         orderVector(cgsort, makeArrayRef(state->v), rvecBuffer.buffer);
2620     }
2621     if (state->flags & (1 << estCGP))
2622     {
2623         orderVector(cgsort, makeArrayRef(state->cg_p), rvecBuffer.buffer);
2624     }
2625
2626     /* Reorder the global cg index */
2627     orderVector<int>(cgsort, dd->globalAtomGroupIndices, &sort->intBuffer);
2628     /* Reorder the cginfo */
2629     orderVector<int>(cgsort, fr->cginfo, &sort->intBuffer);
2630     /* Set the home atom number */
2631     dd->comm->atomRanges.setEnd(DDAtomRanges::Type::Home, dd->ncg_home);
2632
2633     /* The atoms are now exactly in grid order, update the grid order */
2634     fr->nbv->setLocalAtomOrder();
2635 }
2636
2637 //! Accumulates load statistics.
2638 static void add_dd_statistics(gmx_domdec_t* dd)
2639 {
2640     gmx_domdec_comm_t* comm = dd->comm;
2641
2642     for (int i = 0; i < static_cast<int>(DDAtomRanges::Type::Number); i++)
2643     {
2644         auto range = static_cast<DDAtomRanges::Type>(i);
2645         comm->sum_nat[i] += comm->atomRanges.end(range) - comm->atomRanges.start(range);
2646     }
2647     comm->ndecomp++;
2648 }
2649
2650 void reset_dd_statistics_counters(gmx_domdec_t* dd)
2651 {
2652     gmx_domdec_comm_t* comm = dd->comm;
2653
2654     /* Reset all the statistics and counters for total run counting */
2655     for (int i = 0; i < static_cast<int>(DDAtomRanges::Type::Number); i++)
2656     {
2657         comm->sum_nat[i] = 0;
2658     }
2659     comm->ndecomp   = 0;
2660     comm->nload     = 0;
2661     comm->load_step = 0;
2662     comm->load_sum  = 0;
2663     comm->load_max  = 0;
2664     clear_ivec(comm->load_lim);
2665     comm->load_mdf = 0;
2666     comm->load_pme = 0;
2667 }
2668
2669 void print_dd_statistics(const t_commrec* cr, const t_inputrec* ir, FILE* fplog)
2670 {
2671     gmx_domdec_comm_t* comm = cr->dd->comm;
2672
2673     const int numRanges = static_cast<int>(DDAtomRanges::Type::Number);
2674     gmx_sumd(numRanges, comm->sum_nat, cr);
2675
2676     if (fplog == nullptr)
2677     {
2678         return;
2679     }
2680
2681     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");
2682
2683     for (int i = static_cast<int>(DDAtomRanges::Type::Zones); i < numRanges; i++)
2684     {
2685         auto   range = static_cast<DDAtomRanges::Type>(i);
2686         double av    = comm->sum_nat[i] / comm->ndecomp;
2687         switch (range)
2688         {
2689             case DDAtomRanges::Type::Zones:
2690                 fprintf(fplog, " av. #atoms communicated per step for force:  %d x %.1f\n", 2, av);
2691                 break;
2692             case DDAtomRanges::Type::Vsites:
2693                 if (cr->dd->vsite_comm)
2694                 {
2695                     fprintf(fplog,
2696                             " av. #atoms communicated per step for vsites: %d x %.1f\n",
2697                             (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD) ? 3 : 2,
2698                             av);
2699                 }
2700                 break;
2701             case DDAtomRanges::Type::Constraints:
2702                 if (cr->dd->constraint_comm)
2703                 {
2704                     fprintf(fplog,
2705                             " av. #atoms communicated per step for LINCS:  %d x %.1f\n",
2706                             1 + ir->nLincsIter,
2707                             av);
2708                 }
2709                 break;
2710             default: gmx_incons(" Unknown type for DD statistics");
2711         }
2712     }
2713     fprintf(fplog, "\n");
2714
2715     if (comm->ddSettings.recordLoad && EI_DYNAMICS(ir->eI))
2716     {
2717         print_dd_load_av(fplog, cr->dd);
2718     }
2719 }
2720
2721 //!\brief TODO Remove fplog when group scheme and charge groups are gone
2722 void dd_partition_system(FILE*                     fplog,
2723                          const gmx::MDLogger&      mdlog,
2724                          int64_t                   step,
2725                          const t_commrec*          cr,
2726                          gmx_bool                  bMasterState,
2727                          int                       nstglobalcomm,
2728                          t_state*                  state_global,
2729                          const gmx_mtop_t&         top_global,
2730                          const t_inputrec*         ir,
2731                          gmx::ImdSession*          imdSession,
2732                          pull_t*                   pull_work,
2733                          t_state*                  state_local,
2734                          gmx::ForceBuffers*        f,
2735                          gmx::MDAtoms*             mdAtoms,
2736                          gmx_localtop_t*           top_local,
2737                          t_forcerec*               fr,
2738                          gmx::VirtualSitesHandler* vsite,
2739                          gmx::Constraints*         constr,
2740                          t_nrnb*                   nrnb,
2741                          gmx_wallcycle*            wcycle,
2742                          gmx_bool                  bVerbose)
2743 {
2744     gmx_domdec_t*      dd;
2745     gmx_domdec_comm_t* comm;
2746     gmx_ddbox_t        ddbox = { 0 };
2747     int64_t            step_pcoupl;
2748     rvec               cell_ns_x0, cell_ns_x1;
2749     int                ncgindex_set, ncg_moved, nat_f_novirsum;
2750     gmx_bool           bBoxChanged, bNStGlobalComm, bDoDLB, bCheckWhetherToTurnDlbOn, bLogLoad;
2751     gmx_bool           bRedist;
2752     ivec               np;
2753     real               grid_density;
2754     char               sbuf[22];
2755
2756     wallcycle_start(wcycle, ewcDOMDEC);
2757
2758     dd   = cr->dd;
2759     comm = dd->comm;
2760
2761     // TODO if the update code becomes accessible here, use
2762     // upd->deform for this logic.
2763     bBoxChanged = (bMasterState || inputrecDeform(ir));
2764     if (ir->epc != PressureCoupling::No)
2765     {
2766         /* With nstpcouple > 1 pressure coupling happens.
2767          * one step after calculating the pressure.
2768          * Box scaling happens at the end of the MD step,
2769          * after the DD partitioning.
2770          * We therefore have to do DLB in the first partitioning
2771          * after an MD step where P-coupling occurred.
2772          * We need to determine the last step in which p-coupling occurred.
2773          * MRS -- need to validate this for vv?
2774          */
2775         int n = ir->nstpcouple;
2776         if (n == 1)
2777         {
2778             step_pcoupl = step - 1;
2779         }
2780         else
2781         {
2782             step_pcoupl = ((step - 1) / n) * n + 1;
2783         }
2784         if (step_pcoupl >= comm->partition_step)
2785         {
2786             bBoxChanged = TRUE;
2787         }
2788     }
2789
2790     bNStGlobalComm = (step % nstglobalcomm == 0);
2791
2792     if (!isDlbOn(comm))
2793     {
2794         bDoDLB = FALSE;
2795     }
2796     else
2797     {
2798         /* Should we do dynamic load balacing this step?
2799          * Since it requires (possibly expensive) global communication,
2800          * we might want to do DLB less frequently.
2801          */
2802         if (bBoxChanged || ir->epc != PressureCoupling::No)
2803         {
2804             bDoDLB = bBoxChanged;
2805         }
2806         else
2807         {
2808             bDoDLB = bNStGlobalComm;
2809         }
2810     }
2811
2812     /* Check if we have recorded loads on the nodes */
2813     if (comm->ddSettings.recordLoad && dd_load_count(comm) > 0)
2814     {
2815         bCheckWhetherToTurnDlbOn = dd_dlb_get_should_check_whether_to_turn_dlb_on(dd);
2816
2817         /* Print load every nstlog, first and last step to the log file */
2818         bLogLoad = ((ir->nstlog > 0 && step % ir->nstlog == 0) || comm->n_load_collect == 0
2819                     || (ir->nsteps >= 0 && (step + ir->nstlist > ir->init_step + ir->nsteps)));
2820
2821         /* Avoid extra communication due to verbose screen output
2822          * when nstglobalcomm is set.
2823          */
2824         if (bDoDLB || bLogLoad || bCheckWhetherToTurnDlbOn
2825             || (bVerbose && (ir->nstlist == 0 || nstglobalcomm <= ir->nstlist)))
2826         {
2827             get_load_distribution(dd, wcycle);
2828             if (DDMASTER(dd))
2829             {
2830                 if (bLogLoad)
2831                 {
2832                     GMX_LOG(mdlog.info).asParagraph().appendText(dd_print_load(dd, step - 1));
2833                 }
2834                 if (bVerbose)
2835                 {
2836                     dd_print_load_verbose(dd);
2837                 }
2838             }
2839             comm->n_load_collect++;
2840
2841             if (isDlbOn(comm))
2842             {
2843                 if (DDMASTER(dd))
2844                 {
2845                     /* Add the measured cycles to the running average */
2846                     const float averageFactor = 0.1F;
2847                     comm->cyclesPerStepDlbExpAverage =
2848                             (1 - averageFactor) * comm->cyclesPerStepDlbExpAverage
2849                             + averageFactor * comm->cycl[ddCyclStep] / comm->cycl_n[ddCyclStep];
2850                 }
2851                 if (comm->dlbState == DlbState::onCanTurnOff
2852                     && dd->comm->n_load_have % c_checkTurnDlbOffInterval == c_checkTurnDlbOffInterval - 1)
2853                 {
2854                     gmx_bool turnOffDlb;
2855                     if (DDMASTER(dd))
2856                     {
2857                         /* If the running averaged cycles with DLB are more
2858                          * than before we turned on DLB, turn off DLB.
2859                          * We will again run and check the cycles without DLB
2860                          * and we can then decide if to turn off DLB forever.
2861                          */
2862                         turnOffDlb = (comm->cyclesPerStepDlbExpAverage > comm->cyclesPerStepBeforeDLB);
2863                     }
2864                     dd_bcast(dd, sizeof(turnOffDlb), &turnOffDlb);
2865                     if (turnOffDlb)
2866                     {
2867                         /* To turn off DLB, we need to redistribute the atoms */
2868                         dd_collect_state(dd, state_local, state_global);
2869                         bMasterState = TRUE;
2870                         turn_off_dlb(mdlog, dd, step);
2871                     }
2872                 }
2873             }
2874             else if (bCheckWhetherToTurnDlbOn)
2875             {
2876                 gmx_bool turnOffDlbForever = FALSE;
2877                 gmx_bool turnOnDlb         = FALSE;
2878
2879                 /* Since the timings are node dependent, the master decides */
2880                 if (DDMASTER(dd))
2881                 {
2882                     /* If we recently turned off DLB, we want to check if
2883                      * performance is better without DLB. We want to do this
2884                      * ASAP to minimize the chance that external factors
2885                      * slowed down the DLB step are gone here and we
2886                      * incorrectly conclude that DLB was causing the slowdown.
2887                      * So we measure one nstlist block, no running average.
2888                      */
2889                     if (comm->haveTurnedOffDlb
2890                         && comm->cycl[ddCyclStep] / comm->cycl_n[ddCyclStep] < comm->cyclesPerStepDlbExpAverage)
2891                     {
2892                         /* After turning off DLB we ran nstlist steps in fewer
2893                          * cycles than with DLB. This likely means that DLB
2894                          * in not benefical, but this could be due to a one
2895                          * time unlucky fluctuation, so we require two such
2896                          * observations in close succession to turn off DLB
2897                          * forever.
2898                          */
2899                         if (comm->dlbSlowerPartitioningCount > 0
2900                             && dd->ddp_count < comm->dlbSlowerPartitioningCount + 10 * c_checkTurnDlbOnInterval)
2901                         {
2902                             turnOffDlbForever = TRUE;
2903                         }
2904                         comm->haveTurnedOffDlb = false;
2905                         /* Register when we last measured DLB slowdown */
2906                         comm->dlbSlowerPartitioningCount = dd->ddp_count;
2907                     }
2908                     else
2909                     {
2910                         /* Here we check if the max PME rank load is more than 0.98
2911                          * the max PP force load. If so, PP DLB will not help,
2912                          * since we are (almost) limited by PME. Furthermore,
2913                          * DLB will cause a significant extra x/f redistribution
2914                          * cost on the PME ranks, which will then surely result
2915                          * in lower total performance.
2916                          */
2917                         if (comm->ddRankSetup.usePmeOnlyRanks && dd_pme_f_ratio(dd) > 1 - DD_PERF_LOSS_DLB_ON)
2918                         {
2919                             turnOnDlb = FALSE;
2920                         }
2921                         else
2922                         {
2923                             turnOnDlb = (dd_force_imb_perf_loss(dd) >= DD_PERF_LOSS_DLB_ON);
2924                         }
2925                     }
2926                 }
2927                 struct
2928                 {
2929                     gmx_bool turnOffDlbForever;
2930                     gmx_bool turnOnDlb;
2931                 } bools{ turnOffDlbForever, turnOnDlb };
2932                 dd_bcast(dd, sizeof(bools), &bools);
2933                 if (bools.turnOffDlbForever)
2934                 {
2935                     turn_off_dlb_forever(mdlog, dd, step);
2936                 }
2937                 else if (bools.turnOnDlb)
2938                 {
2939                     turn_on_dlb(mdlog, dd, step);
2940                     bDoDLB = TRUE;
2941                 }
2942             }
2943         }
2944         comm->n_load_have++;
2945     }
2946
2947     bRedist = FALSE;
2948     if (bMasterState)
2949     {
2950         /* Clear the old state */
2951         clearDDStateIndices(dd, false);
2952         ncgindex_set = 0;
2953
2954         auto xGlobal = positionsFromStatePointer(state_global);
2955
2956         set_ddbox(*dd, true, DDMASTER(dd) ? state_global->box : nullptr, true, xGlobal, &ddbox);
2957
2958         distributeState(mdlog, dd, top_global, state_global, ddbox, state_local);
2959
2960         /* Ensure that we have space for the new distribution */
2961         dd_resize_atominfo_and_state(fr, state_local, dd->ncg_home);
2962
2963         inc_nrnb(nrnb, eNR_CGCM, comm->atomRanges.numHomeAtoms());
2964
2965         dd_set_cginfo(dd->globalAtomGroupIndices, 0, dd->ncg_home, fr);
2966     }
2967     else if (state_local->ddp_count != dd->ddp_count)
2968     {
2969         if (state_local->ddp_count > dd->ddp_count)
2970         {
2971             gmx_fatal(FARGS,
2972                       "Internal inconsistency state_local->ddp_count (%d) > dd->ddp_count (%" PRId64
2973                       ")",
2974                       state_local->ddp_count,
2975                       dd->ddp_count);
2976         }
2977
2978         if (state_local->ddp_count_cg_gl != state_local->ddp_count)
2979         {
2980             gmx_fatal(FARGS,
2981                       "Internal inconsistency state_local->ddp_count_cg_gl (%d) != "
2982                       "state_local->ddp_count (%d)",
2983                       state_local->ddp_count_cg_gl,
2984                       state_local->ddp_count);
2985         }
2986
2987         /* Clear the old state */
2988         clearDDStateIndices(dd, false);
2989
2990         /* Restore the atom group indices from state_local */
2991         restoreAtomGroups(dd, state_local);
2992         make_dd_indices(dd, 0);
2993         ncgindex_set = dd->ncg_home;
2994
2995         inc_nrnb(nrnb, eNR_CGCM, comm->atomRanges.numHomeAtoms());
2996
2997         dd_set_cginfo(dd->globalAtomGroupIndices, 0, dd->ncg_home, fr);
2998
2999         set_ddbox(*dd, bMasterState, state_local->box, true, state_local->x, &ddbox);
3000
3001         bRedist = isDlbOn(comm);
3002     }
3003     else
3004     {
3005         /* We have the full state, only redistribute the cgs */
3006
3007         /* Clear the non-home indices */
3008         clearDDStateIndices(dd, true);
3009         ncgindex_set = 0;
3010
3011         /* To avoid global communication, we do not recompute the extent
3012          * of the system for dims without pbc. Therefore we need to copy
3013          * the previously computed values when we do not communicate.
3014          */
3015         if (!bNStGlobalComm)
3016         {
3017             copy_rvec(comm->box0, ddbox.box0);
3018             copy_rvec(comm->box_size, ddbox.box_size);
3019         }
3020         set_ddbox(*dd, bMasterState, state_local->box, bNStGlobalComm, state_local->x, &ddbox);
3021
3022         bBoxChanged = TRUE;
3023         bRedist     = TRUE;
3024     }
3025     /* Copy needed for dim's without pbc when avoiding communication */
3026     copy_rvec(ddbox.box0, comm->box0);
3027     copy_rvec(ddbox.box_size, comm->box_size);
3028
3029     set_dd_cell_sizes(dd, &ddbox, dd->unitCellInfo.ddBoxIsDynamic, bMasterState, bDoDLB, step, wcycle);
3030
3031     if (comm->ddSettings.nstDDDumpGrid > 0 && step % comm->ddSettings.nstDDDumpGrid == 0)
3032     {
3033         write_dd_grid_pdb("dd_grid", step, dd, state_local->box, &ddbox);
3034     }
3035
3036     if (comm->systemInfo.useUpdateGroups)
3037     {
3038         comm->updateGroupsCog->addCogs(
3039                 gmx::arrayRefFromArray(dd->globalAtomGroupIndices.data(), dd->ncg_home), state_local->x);
3040     }
3041
3042     /* Check if we should sort the charge groups */
3043     const bool bSortCG = (bMasterState || bRedist);
3044
3045     /* When repartitioning we mark atom groups that will move to neighboring
3046      * DD cells, but we do not move them right away for performance reasons.
3047      * Thus we need to keep track of how many charge groups will move for
3048      * obtaining correct local charge group / atom counts.
3049      */
3050     ncg_moved = 0;
3051     if (bRedist)
3052     {
3053         wallcycle_sub_start(wcycle, ewcsDD_REDIST);
3054
3055         ncgindex_set = dd->ncg_home;
3056         dd_redistribute_cg(fplog, step, dd, ddbox.tric_dir, state_local, fr, nrnb, &ncg_moved);
3057
3058         GMX_RELEASE_ASSERT(bSortCG, "Sorting is required after redistribution");
3059
3060         if (comm->systemInfo.useUpdateGroups)
3061         {
3062             comm->updateGroupsCog->addCogs(
3063                     gmx::arrayRefFromArray(dd->globalAtomGroupIndices.data(), dd->ncg_home),
3064                     state_local->x);
3065         }
3066
3067         wallcycle_sub_stop(wcycle, ewcsDD_REDIST);
3068     }
3069
3070     // TODO: Integrate this code in the nbnxm module
3071     get_nsgrid_boundaries(ddbox.nboundeddim,
3072                           state_local->box,
3073                           dd,
3074                           &ddbox,
3075                           &comm->cell_x0,
3076                           &comm->cell_x1,
3077                           dd->ncg_home,
3078                           as_rvec_array(state_local->x.data()),
3079                           cell_ns_x0,
3080                           cell_ns_x1,
3081                           &grid_density);
3082
3083     if (bBoxChanged)
3084     {
3085         comm_dd_ns_cell_sizes(dd, &ddbox, cell_ns_x0, cell_ns_x1, step);
3086     }
3087
3088     if (bSortCG)
3089     {
3090         wallcycle_sub_start(wcycle, ewcsDD_GRID);
3091
3092         /* Sort the state on charge group position.
3093          * This enables exact restarts from this step.
3094          * It also improves performance by about 15% with larger numbers
3095          * of atoms per node.
3096          */
3097
3098         /* Fill the ns grid with the home cell,
3099          * so we can sort with the indices.
3100          */
3101         set_zones_ncg_home(dd);
3102
3103         set_zones_size(dd, state_local->box, &ddbox, 0, 1, ncg_moved);
3104
3105         nbnxn_put_on_grid(fr->nbv.get(),
3106                           state_local->box,
3107                           0,
3108                           comm->zones.size[0].bb_x0,
3109                           comm->zones.size[0].bb_x1,
3110                           comm->updateGroupsCog.get(),
3111                           { 0, dd->ncg_home },
3112                           comm->zones.dens_zone0,
3113                           fr->cginfo,
3114                           state_local->x,
3115                           ncg_moved,
3116                           bRedist ? comm->movedBuffer.data() : nullptr);
3117
3118         if (debug)
3119         {
3120             fprintf(debug, "Step %s, sorting the %d home charge groups\n", gmx_step_str(step, sbuf), dd->ncg_home);
3121         }
3122         dd_sort_state(dd, fr, state_local);
3123
3124         /* After sorting and compacting we set the correct size */
3125         state_change_natoms(state_local, comm->atomRanges.numHomeAtoms());
3126
3127         /* Rebuild all the indices */
3128         dd->ga2la->clear();
3129         ncgindex_set = 0;
3130
3131         wallcycle_sub_stop(wcycle, ewcsDD_GRID);
3132     }
3133     else
3134     {
3135         /* With the group scheme the sorting array is part of the DD state,
3136          * but it just got out of sync, so mark as invalid by emptying it.
3137          */
3138         if (ir->cutoff_scheme == ecutsGROUP)
3139         {
3140             comm->sort->sorted.clear();
3141         }
3142     }
3143
3144     if (comm->systemInfo.useUpdateGroups)
3145     {
3146         /* The update groups cog's are invalid after sorting
3147          * and need to be cleared before the next partitioning anyhow.
3148          */
3149         comm->updateGroupsCog->clear();
3150     }
3151
3152     wallcycle_sub_start(wcycle, ewcsDD_SETUPCOMM);
3153
3154     /* Set the induces for the home atoms */
3155     set_zones_ncg_home(dd);
3156     make_dd_indices(dd, ncgindex_set);
3157
3158     /* Setup up the communication and communicate the coordinates */
3159     setup_dd_communication(dd, state_local->box, &ddbox, fr, state_local);
3160
3161     /* Set the indices for the halo atoms */
3162     make_dd_indices(dd, dd->ncg_home);
3163
3164     /* Set the charge group boundaries for neighbor searching */
3165     set_cg_boundaries(&comm->zones);
3166
3167     /* When bSortCG=true, we have already set the size for zone 0 */
3168     set_zones_size(dd, state_local->box, &ddbox, bSortCG ? 1 : 0, comm->zones.n, 0);
3169
3170     wallcycle_sub_stop(wcycle, ewcsDD_SETUPCOMM);
3171
3172     /*
3173        write_dd_pdb("dd_home",step,"dump",top_global,cr,
3174                  -1,state_local->x.rvec_array(),state_local->box);
3175      */
3176
3177     wallcycle_sub_start(wcycle, ewcsDD_MAKETOP);
3178
3179     /* Extract a local topology from the global topology */
3180     for (int i = 0; i < dd->ndim; i++)
3181     {
3182         np[dd->dim[i]] = comm->cd[i].numPulses();
3183     }
3184     dd_make_local_top(dd,
3185                       &comm->zones,
3186                       dd->unitCellInfo.npbcdim,
3187                       state_local->box,
3188                       comm->cellsize_min,
3189                       np,
3190                       fr,
3191                       state_local->x.rvec_array(),
3192                       top_global,
3193                       top_local);
3194
3195     wallcycle_sub_stop(wcycle, ewcsDD_MAKETOP);
3196
3197     wallcycle_sub_start(wcycle, ewcsDD_MAKECONSTR);
3198
3199     /* Set up the special atom communication */
3200     int n = comm->atomRanges.end(DDAtomRanges::Type::Zones);
3201     for (int i = static_cast<int>(DDAtomRanges::Type::Zones) + 1;
3202          i < static_cast<int>(DDAtomRanges::Type::Number);
3203          i++)
3204     {
3205         auto range = static_cast<DDAtomRanges::Type>(i);
3206         switch (range)
3207         {
3208             case DDAtomRanges::Type::Vsites:
3209                 if (vsite && vsite->numInterUpdategroupVirtualSites())
3210                 {
3211                     n = dd_make_local_vsites(dd, n, top_local->idef.il);
3212                 }
3213                 break;
3214             case DDAtomRanges::Type::Constraints:
3215                 if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
3216                 {
3217                     /* Only for inter-cg constraints we need special code */
3218                     n = dd_make_local_constraints(dd,
3219                                                   n,
3220                                                   &top_global,
3221                                                   fr->cginfo.data(),
3222                                                   constr,
3223                                                   ir->nProjOrder,
3224                                                   top_local->idef.il);
3225                 }
3226                 break;
3227             default: gmx_incons("Unknown special atom type setup");
3228         }
3229         comm->atomRanges.setEnd(range, n);
3230     }
3231
3232     wallcycle_sub_stop(wcycle, ewcsDD_MAKECONSTR);
3233
3234     wallcycle_sub_start(wcycle, ewcsDD_TOPOTHER);
3235
3236     /* Make space for the extra coordinates for virtual site
3237      * or constraint communication.
3238      */
3239     state_local->natoms = comm->atomRanges.numAtomsTotal();
3240
3241     state_change_natoms(state_local, state_local->natoms);
3242
3243     if (vsite && vsite->numInterUpdategroupVirtualSites())
3244     {
3245         nat_f_novirsum = comm->atomRanges.end(DDAtomRanges::Type::Vsites);
3246     }
3247     else
3248     {
3249         if (EEL_FULL(ir->coulombtype) && dd->haveExclusions)
3250         {
3251             nat_f_novirsum = comm->atomRanges.end(DDAtomRanges::Type::Zones);
3252         }
3253         else
3254         {
3255             nat_f_novirsum = comm->atomRanges.numHomeAtoms();
3256         }
3257     }
3258
3259     /* Set the number of atoms required for the force calculation.
3260      * Forces need to be constrained when doing energy
3261      * minimization. For simple simulations we could avoid some
3262      * allocation, zeroing and copying, but this is probably not worth
3263      * the complications and checking.
3264      */
3265     forcerec_set_ranges(fr,
3266                         comm->atomRanges.end(DDAtomRanges::Type::Zones),
3267                         comm->atomRanges.end(DDAtomRanges::Type::Constraints),
3268                         nat_f_novirsum);
3269
3270     /* Update atom data for mdatoms and several algorithms */
3271     mdAlgorithmsSetupAtomData(cr, ir, top_global, top_local, fr, f, mdAtoms, constr, vsite, nullptr);
3272
3273     auto mdatoms = mdAtoms->mdatoms();
3274     if (!thisRankHasDuty(cr, DUTY_PME))
3275     {
3276         /* Send the charges and/or c6/sigmas to our PME only node */
3277         gmx_pme_send_parameters(cr,
3278                                 fr->ic,
3279                                 mdatoms->nChargePerturbed != 0,
3280                                 mdatoms->nTypePerturbed != 0,
3281                                 mdatoms->chargeA,
3282                                 mdatoms->chargeB,
3283                                 mdatoms->sqrt_c6A,
3284                                 mdatoms->sqrt_c6B,
3285                                 mdatoms->sigmaA,
3286                                 mdatoms->sigmaB,
3287                                 dd_pme_maxshift_x(dd),
3288                                 dd_pme_maxshift_y(dd));
3289     }
3290
3291     if (dd->atomSets != nullptr)
3292     {
3293         /* Update the local atom sets */
3294         dd->atomSets->setIndicesInDomainDecomposition(*(dd->ga2la));
3295     }
3296
3297     // The pull group construction can need the atom sets updated above
3298     if (ir->bPull)
3299     {
3300         /* Update the local pull groups */
3301         dd_make_local_pull_groups(cr, pull_work);
3302     }
3303
3304     /* Update the local atoms to be communicated via the IMD protocol if bIMD is TRUE. */
3305     imdSession->dd_make_local_IMD_atoms(dd);
3306
3307     add_dd_statistics(dd);
3308
3309     /* Make sure we only count the cycles for this DD partitioning */
3310     clear_dd_cycle_counts(dd);
3311
3312     /* Because the order of the atoms might have changed since
3313      * the last vsite construction, we need to communicate the constructing
3314      * atom coordinates again (for spreading the forces this MD step).
3315      */
3316     dd_move_x_vsites(*dd, state_local->box, state_local->x.rvec_array());
3317
3318     wallcycle_sub_stop(wcycle, ewcsDD_TOPOTHER);
3319
3320     if (comm->ddSettings.nstDDDump > 0 && step % comm->ddSettings.nstDDDump == 0)
3321     {
3322         dd_move_x(dd, state_local->box, state_local->x, nullWallcycle);
3323         write_dd_pdb("dd_dump",
3324                      step,
3325                      "dump",
3326                      &top_global,
3327                      cr,
3328                      -1,
3329                      state_local->x.rvec_array(),
3330                      state_local->box);
3331     }
3332
3333     /* Store the partitioning step */
3334     comm->partition_step = step;
3335
3336     /* Increase the DD partitioning counter */
3337     dd->ddp_count++;
3338     /* The state currently matches this DD partitioning count, store it */
3339     state_local->ddp_count = dd->ddp_count;
3340     if (bMasterState)
3341     {
3342         /* The DD master node knows the complete cg distribution,
3343          * store the count so we can possibly skip the cg info communication.
3344          */
3345         comm->master_cg_ddp_count = (bSortCG ? 0 : dd->ddp_count);
3346     }
3347
3348     if (comm->ddSettings.DD_debug > 0)
3349     {
3350         /* Set the env var GMX_DD_DEBUG if you suspect corrupted indices */
3351         check_index_consistency(dd, top_global.natoms, "after partitioning");
3352     }
3353
3354     wallcycle_stop(wcycle, ewcDOMDEC);
3355 }
3356
3357 /*! \brief Check whether bonded interactions are missing, if appropriate */
3358 void checkNumberOfBondedInteractions(const gmx::MDLogger&           mdlog,
3359                                      t_commrec*                     cr,
3360                                      int                            totalNumberOfBondedInteractions,
3361                                      const gmx_mtop_t*              top_global,
3362                                      const gmx_localtop_t*          top_local,
3363                                      gmx::ArrayRef<const gmx::RVec> x,
3364                                      const matrix                   box,
3365                                      bool* shouldCheckNumberOfBondedInteractions)
3366 {
3367     if (*shouldCheckNumberOfBondedInteractions)
3368     {
3369         if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
3370         {
3371             dd_print_missing_interactions(
3372                     mdlog, cr, totalNumberOfBondedInteractions, top_global, top_local, x, box); // Does not return
3373         }
3374         *shouldCheckNumberOfBondedInteractions = false;
3375     }
3376 }