Extract DDSystemInfo from gmx_domdec_comm_t
[alexxy/gromacs.git] / src / gromacs / domdec / cellsizes.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
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 Implements DD cell-size related functions.
38  *
39  * \author Berk Hess <hess@kth.se>
40  * \author Roland Schulz <roland.schulz@intel.com>
41  * \ingroup module_domdec
42  */
43
44 #include "gmxpre.h"
45
46 #include "cellsizes.h"
47
48 #include "config.h"
49
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdtypes/commrec.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55
56 #include "atomdistribution.h"
57 #include "domdec_internal.h"
58 #include "utility.h"
59
60 static void set_pme_maxshift(gmx_domdec_t *dd, gmx_ddpme_t *ddpme,
61                              gmx_bool bUniform, const gmx_ddbox_t *ddbox,
62                              const real *cellFrac)
63 {
64     gmx_domdec_comm_t *comm;
65     int                nc, ns, s;
66     int               *xmin, *xmax;
67     real               range, pme_boundary;
68     int                sh;
69
70     comm = dd->comm;
71     nc   = dd->nc[ddpme->dim];
72     ns   = ddpme->nslab;
73
74     if (!ddpme->dim_match)
75     {
76         /* PP decomposition is not along dim: the worst situation */
77         sh = ns/2;
78     }
79     else if (ns <= 3 || (bUniform && ns == nc))
80     {
81         /* The optimal situation */
82         sh = 1;
83     }
84     else
85     {
86         /* We need to check for all pme nodes which nodes they
87          * could possibly need to communicate with.
88          */
89         xmin = ddpme->pp_min;
90         xmax = ddpme->pp_max;
91         /* Allow for atoms to be maximally 2/3 times the cut-off
92          * out of their DD cell. This is a reasonable balance between
93          * between performance and support for most charge-group/cut-off
94          * combinations.
95          */
96         range  = 2.0/3.0*comm->systemInfo.cutoff/ddbox->box_size[ddpme->dim];
97         /* Avoid extra communication when we are exactly at a boundary */
98         range *= 0.999;
99
100         sh = 1;
101         for (s = 0; s < ns; s++)
102         {
103             /* PME slab s spreads atoms between box frac. s/ns and (s+1)/ns */
104             pme_boundary = static_cast<real>(s)/ns;
105             while (sh + 1 < ns &&
106                    ((s - (sh + 1) >= 0 &&
107                      cellFrac[xmax[s - (sh + 1)     ] + 1]     + range > pme_boundary) ||
108                     (s - (sh + 1) <  0 &&
109                      cellFrac[xmax[s - (sh + 1) + ns] + 1] - 1 + range > pme_boundary)))
110             {
111                 sh++;
112             }
113             pme_boundary = static_cast<real>(s+1)/ns;
114             while (sh + 1 < ns &&
115                    ((s + (sh + 1) <  ns &&
116                      cellFrac[xmin[s + (sh + 1)    ]  ]     - range < pme_boundary) ||
117                     (s + (sh + 1) >= ns &&
118                      cellFrac[xmin[s +(sh + 1) - ns]  ] + 1 - range < pme_boundary)))
119             {
120                 sh++;
121             }
122         }
123     }
124
125     ddpme->maxshift = sh;
126
127     if (debug)
128     {
129         fprintf(debug, "PME slab communication range for dim %d is %d\n",
130                 ddpme->dim, ddpme->maxshift);
131     }
132 }
133
134 static void check_box_size(const gmx_domdec_t *dd,
135                            const gmx_ddbox_t  *ddbox)
136 {
137     int d, dim;
138
139     for (d = 0; d < dd->ndim; d++)
140     {
141         dim = dd->dim[d];
142         if (dim < ddbox->nboundeddim &&
143             ddbox->box_size[dim]*ddbox->skew_fac[dim] <
144             dd->nc[dim]*dd->comm->cellsize_limit*DD_CELL_MARGIN)
145         {
146             gmx_fatal(FARGS, "The %c-size of the box (%f) times the triclinic skew factor (%f) is smaller than the number of DD cells (%d) times the smallest allowed cell size (%f)\n",
147                       dim2char(dim), ddbox->box_size[dim], ddbox->skew_fac[dim],
148                       dd->nc[dim], dd->comm->cellsize_limit);
149         }
150     }
151 }
152
153 real grid_jump_limit(const gmx_domdec_comm_t *comm,
154                      real                     cutoff,
155                      int                      dim_ind)
156 {
157     real grid_jump_limit;
158
159     /* The distance between the boundaries of cells at distance
160      * x+-1,y+-1 or y+-1,z+-1 is limited by the cut-off restrictions
161      * and by the fact that cells should not be shifted by more than
162      * half their size, such that cg's only shift by one cell
163      * at redecomposition.
164      */
165     grid_jump_limit = comm->cellsize_limit;
166     if (!comm->bVacDLBNoLimit)
167     {
168         if (comm->bPMELoadBalDLBLimits)
169         {
170             cutoff = std::max(cutoff, comm->PMELoadBal_max_cutoff);
171         }
172         grid_jump_limit = std::max(grid_jump_limit,
173                                    cutoff/comm->cd[dim_ind].numPulses());
174     }
175
176     return grid_jump_limit;
177 }
178
179 /* This function should be used for moving the domain boudaries during DLB,
180  * for obtaining the minimum cell size. It checks the initially set limit
181  * comm->cellsize_min, for bonded and initial non-bonded cut-offs,
182  * and, possibly, a longer cut-off limit set for PME load balancing.
183  */
184 static real cellsize_min_dlb(gmx_domdec_comm_t *comm, int dim_ind, int dim)
185 {
186     real cellsize_min;
187
188     cellsize_min = comm->cellsize_min[dim];
189
190     if (!comm->bVacDLBNoLimit)
191     {
192         /* The cut-off might have changed, e.g. by PME load balacning,
193          * from the value used to set comm->cellsize_min, so check it.
194          */
195         cellsize_min = std::max(cellsize_min, comm->systemInfo.cutoff/comm->cd[dim_ind].np_dlb);
196
197         if (comm->bPMELoadBalDLBLimits)
198         {
199             /* Check for the cut-off limit set by the PME load balancing */
200             cellsize_min = std::max(cellsize_min, comm->PMELoadBal_max_cutoff/comm->cd[dim_ind].np_dlb);
201         }
202     }
203
204     return cellsize_min;
205 }
206
207 /* Set the domain boundaries. Use for static (or no) load balancing,
208  * and also for the starting state for dynamic load balancing.
209  * setmode determine if and where the boundaries are stored, use enum above.
210  * Returns the number communication pulses in npulse.
211  */
212 gmx::ArrayRef < const std::vector < real>>
213 set_dd_cell_sizes_slb(gmx_domdec_t *dd, const gmx_ddbox_t *ddbox,
214                       int setmode, ivec npulse)
215 {
216     gmx_domdec_comm_t *comm = dd->comm;
217
218     gmx::ArrayRef < std::vector < real>> cell_x_master;
219     if (setmode == setcellsizeslbMASTER)
220     {
221         cell_x_master = dd->ma->cellSizesBuffer;
222     }
223
224     rvec cellsize_min;
225     for (int d = 0; d < DIM; d++)
226     {
227         cellsize_min[d] = ddbox->box_size[d]*ddbox->skew_fac[d];
228         npulse[d]       = 1;
229         if (dd->nc[d] == 1 || comm->slb_frac[d] == nullptr)
230         {
231             /* Uniform grid */
232             real cell_dx = ddbox->box_size[d]/dd->nc[d];
233             switch (setmode)
234             {
235                 case setcellsizeslbMASTER:
236                     for (int j = 0; j < dd->nc[d]+1; j++)
237                     {
238                         cell_x_master[d][j] = ddbox->box0[d] + j*cell_dx;
239                     }
240                     break;
241                 case setcellsizeslbLOCAL:
242                     comm->cell_x0[d] = ddbox->box0[d] + (dd->ci[d]  )*cell_dx;
243                     comm->cell_x1[d] = ddbox->box0[d] + (dd->ci[d]+1)*cell_dx;
244                     break;
245                 default:
246                     break;
247             }
248             real cellsize = cell_dx*ddbox->skew_fac[d];
249             while (cellsize*npulse[d] < comm->systemInfo.cutoff)
250             {
251                 npulse[d]++;
252             }
253             cellsize_min[d] = cellsize;
254         }
255         else
256         {
257             /* Statically load balanced grid */
258             /* Also when we are not doing a master distribution we determine
259              * all cell borders in a loop to obtain identical values
260              * to the master distribution case and to determine npulse.
261              */
262             gmx::ArrayRef<real> cell_x;
263             std::vector<real>   cell_x_buffer;
264             if (setmode == setcellsizeslbMASTER)
265             {
266                 cell_x = cell_x_master[d];
267             }
268             else
269             {
270                 cell_x_buffer.resize(dd->nc[d] + 1);
271                 cell_x = cell_x_buffer;
272             }
273             cell_x[0] = ddbox->box0[d];
274             for (int j = 0; j < dd->nc[d]; j++)
275             {
276                 real cell_dx  = ddbox->box_size[d]*comm->slb_frac[d][j];
277                 cell_x[j+1]   = cell_x[j] + cell_dx;
278                 real cellsize = cell_dx*ddbox->skew_fac[d];
279                 while (cellsize*npulse[d] < comm->systemInfo.cutoff &&
280                        npulse[d] < dd->nc[d]-1)
281                 {
282                     npulse[d]++;
283                 }
284                 cellsize_min[d] = std::min(cellsize_min[d], cellsize);
285             }
286             if (setmode == setcellsizeslbLOCAL)
287             {
288                 comm->cell_x0[d] = cell_x[dd->ci[d]];
289                 comm->cell_x1[d] = cell_x[dd->ci[d]+1];
290             }
291         }
292         /* The following limitation is to avoid that a cell would receive
293          * some of its own home charge groups back over the periodic boundary.
294          * Double charge groups cause trouble with the global indices.
295          */
296         if (d < ddbox->npbcdim &&
297             dd->nc[d] > 1 && npulse[d] >= dd->nc[d])
298         {
299             char error_string[STRLEN];
300
301             sprintf(error_string,
302                     "The box size in direction %c (%f) times the triclinic skew factor (%f) is too small for a cut-off of %f with %d domain decomposition cells, use 1 or more than %d %s or increase the box size in this direction",
303                     dim2char(d), ddbox->box_size[d], ddbox->skew_fac[d],
304                     comm->systemInfo.cutoff,
305                     dd->nc[d], dd->nc[d],
306                     dd->nnodes > dd->nc[d] ? "cells" : "ranks");
307
308             if (setmode == setcellsizeslbLOCAL)
309             {
310                 gmx_fatal_collective(FARGS, dd->mpi_comm_all, DDMASTER(dd),
311                                      "%s", error_string);
312             }
313             else
314             {
315                 gmx_fatal(FARGS, "%s", error_string);
316             }
317         }
318     }
319
320     if (!isDlbOn(comm))
321     {
322         copy_rvec(cellsize_min, comm->cellsize_min);
323     }
324
325     for (int d = 0; d < comm->npmedecompdim; d++)
326     {
327         set_pme_maxshift(dd, &comm->ddpme[d],
328                          comm->slb_frac[dd->dim[d]] == nullptr, ddbox,
329                          comm->ddpme[d].slb_dim_f);
330     }
331
332     return cell_x_master;
333 }
334
335
336 static void dd_cell_sizes_dlb_root_enforce_limits(gmx_domdec_t *dd,
337                                                   int d, int dim, RowMaster *rowMaster,
338                                                   const gmx_ddbox_t *ddbox,
339                                                   gmx_bool bUniform, int64_t step, real cellsize_limit_f, int range[])
340 {
341     gmx_domdec_comm_t *comm;
342     real               halfway, cellsize_limit_f_i, region_size;
343     gmx_bool           bLastHi  = FALSE;
344     int                nrange[] = {range[0], range[1]};
345
346     region_size = rowMaster->cellFrac[range[1]] - rowMaster->cellFrac[range[0]];
347
348     GMX_ASSERT(region_size >= (range[1] - range[0])*cellsize_limit_f,
349                "The region should fit all cells at minimum size");
350
351     comm = dd->comm;
352
353     const int           ncd       = dd->nc[dim];
354
355     const bool          dimHasPbc = (dim < ddbox->npbcdim);
356
357     gmx::ArrayRef<real> cell_size = rowMaster->buf_ncd;
358
359     if (debug)
360     {
361         fprintf(debug, "enforce_limits: %d %d\n", range[0], range[1]);
362     }
363
364     /* First we need to check if the scaling does not make cells
365      * smaller than the smallest allowed size.
366      * We need to do this iteratively, since if a cell is too small,
367      * it needs to be enlarged, which makes all the other cells smaller,
368      * which could in turn make another cell smaller than allowed.
369      */
370     for (int i = range[0]; i < range[1]; i++)
371     {
372         rowMaster->isCellMin[i] = false;
373     }
374     int nmin = 0;
375     int nmin_old;
376     do
377     {
378         nmin_old = nmin;
379         /* We need the total for normalization */
380         real fac = 0;
381         for (int i = range[0]; i < range[1]; i++)
382         {
383             if (!rowMaster->isCellMin[i])
384             {
385                 fac += cell_size[i];
386             }
387         }
388         /* This condition is ensured by the assertion at the end of the loop */
389         GMX_ASSERT(fac > 0, "We cannot have 0 size to work with");
390         fac = (region_size - nmin*cellsize_limit_f)/fac; /* substracting cells already set to cellsize_limit_f */
391         /* Determine the cell boundaries */
392         for (int i = range[0]; i < range[1]; i++)
393         {
394             if (!rowMaster->isCellMin[i])
395             {
396                 cell_size[i] *= fac;
397                 if (!dimHasPbc && (i == 0 || i == dd->nc[dim] - 1))
398                 {
399                     cellsize_limit_f_i = 0;
400                 }
401                 else
402                 {
403                     cellsize_limit_f_i = cellsize_limit_f;
404                 }
405                 if (cell_size[i] < cellsize_limit_f_i)
406                 {
407                     rowMaster->isCellMin[i] = true;
408                     cell_size[i]            = cellsize_limit_f_i;
409                     nmin++;
410                 }
411             }
412             rowMaster->cellFrac[i + 1] = rowMaster->cellFrac[i] + cell_size[i];
413         }
414
415         /* This is ensured by the assertion at the beginning of this function */
416         GMX_ASSERT(nmin < range[1] - range[0], "We can not have all cells limited");
417     }
418     while (nmin > nmin_old);
419
420     const int i  = range[1] - 1;
421     cell_size[i] = rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i];
422     /* For this check we should not use DD_CELL_MARGIN,
423      * but a slightly smaller factor,
424      * since rounding could get use below the limit.
425      */
426     if (dimHasPbc &&
427         cell_size[i] < cellsize_limit_f*DD_CELL_MARGIN2/DD_CELL_MARGIN)
428     {
429         char buf[22];
430         gmx_fatal(FARGS, "step %s: the dynamic load balancing could not balance dimension %c: box size %f, triclinic skew factor %f, #cells %d, minimum cell size %f\n",
431                   gmx_step_str(step, buf),
432                   dim2char(dim), ddbox->box_size[dim], ddbox->skew_fac[dim],
433                   ncd, comm->cellsize_min[dim]);
434     }
435
436     rowMaster->dlbIsLimited = (nmin > 0) || (range[0] > 0) || (range[1] < ncd);
437
438     if (!bUniform)
439     {
440         /* Check if the boundary did not displace more than halfway
441          * each of the cells it bounds, as this could cause problems,
442          * especially when the differences between cell sizes are large.
443          * If changes are applied, they will not make cells smaller
444          * than the cut-off, as we check all the boundaries which
445          * might be affected by a change and if the old state was ok,
446          * the cells will at most be shrunk back to their old size.
447          */
448         for (int i = range[0]+1; i < range[1]; i++)
449         {
450             halfway = 0.5*(rowMaster->oldCellFrac[i] + rowMaster->oldCellFrac[i - 1]);
451             if (rowMaster->cellFrac[i] < halfway)
452             {
453                 rowMaster->cellFrac[i] = halfway;
454                 /* Check if the change also causes shifts of the next boundaries */
455                 for (int j = i + 1; j < range[1]; j++)
456                 {
457                     if (rowMaster->cellFrac[j] < rowMaster->cellFrac[j - 1] + cellsize_limit_f)
458                     {
459                         rowMaster->cellFrac[j] =  rowMaster->cellFrac[j - 1] + cellsize_limit_f;
460                     }
461                 }
462             }
463             halfway = 0.5*(rowMaster->oldCellFrac[i] + rowMaster->oldCellFrac[i + 1]);
464             if (rowMaster->cellFrac[i] > halfway)
465             {
466                 rowMaster->cellFrac[i] = halfway;
467                 /* Check if the change also causes shifts of the next boundaries */
468                 for (int j = i - 1; j >= range[0] + 1; j--)
469                 {
470                     if (rowMaster->cellFrac[j] > rowMaster->cellFrac[j + 1] - cellsize_limit_f)
471                     {
472                         rowMaster->cellFrac[j] = rowMaster->cellFrac[j + 1] - cellsize_limit_f;
473                     }
474                 }
475             }
476         }
477     }
478
479     /* nrange is defined as [lower, upper) range for new call to enforce_limits */
480     /* find highest violation of LimLo (a) and the following violation of LimHi (thus the lowest following) (b)
481      * then call enforce_limits for (oldb,a), (a,b). In the next step: (b,nexta). oldb and nexta can be the boundaries.
482      * for a and b nrange is used */
483     if (d > 0)
484     {
485         /* Take care of the staggering of the cell boundaries */
486         if (bUniform)
487         {
488             for (int i = range[0]; i < range[1]; i++)
489             {
490                 rowMaster->bounds[i].cellFracLowerMax = rowMaster->cellFrac[i];
491                 rowMaster->bounds[i].cellFracUpperMin = rowMaster->cellFrac[i + 1];
492             }
493         }
494         else
495         {
496             for (int i = range[0] + 1; i < range[1]; i++)
497             {
498                 const RowMaster::Bounds &bounds = rowMaster->bounds[i];
499
500                 bool                     bLimLo = (rowMaster->cellFrac[i] < bounds.boundMin);
501                 bool                     bLimHi = (rowMaster->cellFrac[i] > bounds.boundMax);
502                 if (bLimLo && bLimHi)
503                 {
504                     /* Both limits violated, try the best we can */
505                     /* For this case we split the original range (range) in two parts and care about the other limitiations in the next iteration. */
506                     rowMaster->cellFrac[i] = 0.5*(bounds.boundMin + bounds.boundMax);
507                     nrange[0]              = range[0];
508                     nrange[1]              = i;
509                     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
510
511                     nrange[0] = i;
512                     nrange[1] = range[1];
513                     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
514
515                     return;
516                 }
517                 else if (bLimLo)
518                 {
519                     /* rowMaster->cellFrac[i] = rowMaster->boundMin[i]; */
520                     nrange[1] = i;  /* only store violation location. There could be a LimLo violation following with an higher index */
521                     bLastHi   = FALSE;
522                 }
523                 else if (bLimHi && !bLastHi)
524                 {
525                     bLastHi = TRUE;
526                     if (nrange[1] < range[1])   /* found a LimLo before */
527                     {
528                         rowMaster->cellFrac[nrange[1]] = rowMaster->bounds[nrange[1]].boundMin;
529                         dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
530                         nrange[0] = nrange[1];
531                     }
532                     rowMaster->cellFrac[i] = rowMaster->bounds[i].boundMax;
533                     nrange[1]              = i;
534                     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
535                     nrange[0] = i;
536                     nrange[1] = range[1];
537                 }
538             }
539             if (nrange[1] < range[1])   /* found last a LimLo */
540             {
541                 rowMaster->cellFrac[nrange[1]] = rowMaster->bounds[nrange[1]].boundMin;
542                 dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
543                 nrange[0] = nrange[1];
544                 nrange[1] = range[1];
545                 dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
546             }
547             else if (nrange[0] > range[0]) /* found at least one LimHi */
548             {
549                 dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
550             }
551         }
552     }
553 }
554
555
556 static void set_dd_cell_sizes_dlb_root(gmx_domdec_t *dd,
557                                        int d, int dim, RowMaster *rowMaster,
558                                        const gmx_ddbox_t *ddbox,
559                                        gmx_bool bDynamicBox,
560                                        gmx_bool bUniform, int64_t step)
561 {
562     gmx_domdec_comm_t *comm    = dd->comm;
563     constexpr real     c_relax = 0.5;
564     int                range[] = { 0, 0 };
565
566     /* Convert the maximum change from the input percentage to a fraction */
567     const real          change_limit = comm->dlb_scale_lim*0.01;
568
569     const int           ncd          = dd->nc[dim];
570
571     const bool          bPBC = (dim < ddbox->npbcdim);
572
573     gmx::ArrayRef<real> cell_size = rowMaster->buf_ncd;
574
575     /* Store the original boundaries */
576     for (int i = 0; i < ncd + 1; i++)
577     {
578         rowMaster->oldCellFrac[i] = rowMaster->cellFrac[i];
579     }
580     if (bUniform)
581     {
582         for (int i = 0; i < ncd; i++)
583         {
584             cell_size[i] = 1.0/ncd;
585         }
586     }
587     else if (dd_load_count(comm) > 0)
588     {
589         real load_aver  = comm->load[d].sum_m/ncd;
590         real change_max = 0;
591         real load_i;
592         real change;
593         for (int i = 0; i < ncd; i++)
594         {
595             /* Determine the relative imbalance of cell i */
596             load_i         = comm->load[d].load[i*comm->load[d].nload+2];
597             real imbalance = (load_i - load_aver)/(load_aver > 0 ? load_aver : 1);
598             /* Determine the change of the cell size using underrelaxation */
599             change     = -c_relax*imbalance;
600             change_max = std::max(change_max, std::max(change, -change));
601         }
602         /* Limit the amount of scaling.
603          * We need to use the same rescaling for all cells in one row,
604          * otherwise the load balancing might not converge.
605          */
606         real sc = c_relax;
607         if (change_max > change_limit)
608         {
609             sc *= change_limit/change_max;
610         }
611         for (int i = 0; i < ncd; i++)
612         {
613             /* Determine the relative imbalance of cell i */
614             load_i         = comm->load[d].load[i*comm->load[d].nload+2];
615             real imbalance = (load_i - load_aver)/(load_aver > 0 ? load_aver : 1);
616             /* Determine the change of the cell size using underrelaxation */
617             change       = -sc*imbalance;
618             cell_size[i] = (rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i])*(1 + change);
619         }
620     }
621
622     real cellsize_limit_f  = cellsize_min_dlb(comm, d, dim)/ddbox->box_size[dim];
623     cellsize_limit_f      *= DD_CELL_MARGIN;
624     real dist_min_f_hard   = grid_jump_limit(comm, comm->systemInfo.cutoff, d)/ddbox->box_size[dim];
625     real dist_min_f        = dist_min_f_hard * DD_CELL_MARGIN;
626     if (ddbox->tric_dir[dim])
627     {
628         cellsize_limit_f /= ddbox->skew_fac[dim];
629         dist_min_f       /= ddbox->skew_fac[dim];
630     }
631     if (bDynamicBox && d > 0)
632     {
633         dist_min_f *= DD_PRES_SCALE_MARGIN;
634     }
635     if (d > 0 && !bUniform)
636     {
637         /* Make sure that the grid is not shifted too much */
638         for (int i = 1; i < ncd; i++)
639         {
640             const RowMaster::Bounds &boundsNeighbor = rowMaster->bounds[i - 1];
641             RowMaster::Bounds       &bounds         = rowMaster->bounds[i];
642
643             if (bounds.cellFracUpperMin - boundsNeighbor.cellFracLowerMax < 2 * dist_min_f_hard)
644             {
645                 gmx_incons("Inconsistent DD boundary staggering limits!");
646             }
647             bounds.boundMin = boundsNeighbor.cellFracLowerMax + dist_min_f;
648             real space       = rowMaster->cellFrac[i] - (boundsNeighbor.cellFracLowerMax + dist_min_f);
649             if (space > 0)
650             {
651                 bounds.boundMin += 0.5*space;
652             }
653             bounds.boundMax  = bounds.cellFracUpperMin - dist_min_f;
654             space            = rowMaster->cellFrac[i] - (bounds.cellFracUpperMin - dist_min_f);
655             if (space < 0)
656             {
657                 rowMaster->bounds[i].boundMax += 0.5*space;
658             }
659             if (debug)
660             {
661                 fprintf(debug,
662                         "dim %d boundary %d %.3f < %.3f < %.3f < %.3f < %.3f\n",
663                         d, i,
664                         boundsNeighbor.cellFracLowerMax + dist_min_f,
665                         bounds.boundMin, rowMaster->cellFrac[i], bounds.boundMax,
666                         bounds.cellFracUpperMin - dist_min_f);
667             }
668         }
669     }
670     range[1]                 = ncd;
671     rowMaster->cellFrac[0]   = 0;
672     rowMaster->cellFrac[ncd] = 1;
673     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, range);
674
675
676     /* After the checks above, the cells should obey the cut-off
677      * restrictions, but it does not hurt to check.
678      */
679     for (int i = 0; i < ncd; i++)
680     {
681         if (debug)
682         {
683             fprintf(debug, "Relative bounds dim %d  cell %d: %f %f\n",
684                     dim, i, rowMaster->cellFrac[i], rowMaster->cellFrac[i + 1]);
685         }
686
687         if ((bPBC || (i != 0 && i != dd->nc[dim]-1)) &&
688             rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i] <
689             cellsize_limit_f/DD_CELL_MARGIN)
690         {
691             char buf[22];
692             fprintf(stderr,
693                     "\nWARNING step %s: direction %c, cell %d too small: %f\n",
694                     gmx_step_str(step, buf), dim2char(dim), i,
695                     (rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i])
696                     *ddbox->box_size[dim]*ddbox->skew_fac[dim]);
697         }
698     }
699
700     int pos = ncd + 1;
701     /* Store the cell boundaries of the lower dimensions at the end */
702     for (int d1 = 0; d1 < d; d1++)
703     {
704         rowMaster->cellFrac[pos++] = comm->cellsizesWithDlb[d1].fracLower;
705         rowMaster->cellFrac[pos++] = comm->cellsizesWithDlb[d1].fracUpper;
706     }
707
708     if (d < comm->npmedecompdim)
709     {
710         /* The master determines the maximum shift for
711          * the coordinate communication between separate PME nodes.
712          */
713         set_pme_maxshift(dd, &comm->ddpme[d], bUniform, ddbox, rowMaster->cellFrac.data());
714     }
715     rowMaster->cellFrac[pos++] = comm->ddpme[0].maxshift;
716     if (d >= 1)
717     {
718         rowMaster->cellFrac[pos++] = comm->ddpme[1].maxshift;
719     }
720 }
721
722 static void relative_to_absolute_cell_bounds(gmx_domdec_t      *dd,
723                                              const gmx_ddbox_t *ddbox,
724                                              int                dimind)
725 {
726     gmx_domdec_comm_t        *comm      = dd->comm;
727     const DDCellsizesWithDlb &cellsizes = comm->cellsizesWithDlb[dimind];
728
729     /* Set the cell dimensions */
730     int dim            = dd->dim[dimind];
731     comm->cell_x0[dim] = cellsizes.fracLower*ddbox->box_size[dim];
732     comm->cell_x1[dim] = cellsizes.fracUpper*ddbox->box_size[dim];
733     if (dim >= ddbox->nboundeddim)
734     {
735         comm->cell_x0[dim] += ddbox->box0[dim];
736         comm->cell_x1[dim] += ddbox->box0[dim];
737     }
738 }
739
740 static void distribute_dd_cell_sizes_dlb(gmx_domdec_t *dd,
741                                          int d, int dim,
742                                          gmx::ArrayRef<real> cellFracRow,
743                                          const gmx_ddbox_t *ddbox)
744 {
745     gmx_domdec_comm_t &comm = *dd->comm;
746
747 #if GMX_MPI
748     /* Each node would only need to know two fractions,
749      * but it is probably cheaper to broadcast the whole array.
750      */
751     MPI_Bcast(cellFracRow.data(), ddCellFractionBufferSize(dd, d)*sizeof(real), MPI_BYTE,
752               0, comm.mpi_comm_load[d]);
753 #endif
754     /* Copy the fractions for this dimension from the buffer */
755     comm.cellsizesWithDlb[d].fracLower = cellFracRow[dd->ci[dim]    ];
756     comm.cellsizesWithDlb[d].fracUpper = cellFracRow[dd->ci[dim] + 1];
757     /* The whole array was communicated, so set the buffer position */
758     int pos = dd->nc[dim] + 1;
759     for (int d1 = 0; d1 <= d; d1++)
760     {
761         if (d1 < d)
762         {
763             /* Copy the cell fractions of the lower dimensions */
764             comm.cellsizesWithDlb[d1].fracLower = cellFracRow[pos++];
765             comm.cellsizesWithDlb[d1].fracUpper = cellFracRow[pos++];
766         }
767         relative_to_absolute_cell_bounds(dd, ddbox, d1);
768     }
769     /* Convert the communicated shift from float to int */
770     comm.ddpme[0].maxshift = gmx::roundToInt(cellFracRow[pos++]);
771     if (d >= 1)
772     {
773         comm.ddpme[1].maxshift = gmx::roundToInt(cellFracRow[pos++]);
774     }
775 }
776
777 static void set_dd_cell_sizes_dlb_change(gmx_domdec_t *dd,
778                                          const gmx_ddbox_t *ddbox,
779                                          gmx_bool bDynamicBox,
780                                          gmx_bool bUniform, int64_t step)
781 {
782     for (int d = 0; d < dd->ndim; d++)
783     {
784         const int dim         = dd->dim[d];
785         bool      isRowMember = true;
786         bool      isRowRoot   = true;
787         for (int d1 = d; d1 < dd->ndim; d1++)
788         {
789             if (dd->ci[dd->dim[d1]] > 0)
790             {
791                 if (d1 != d)
792                 {
793                     isRowMember = false;
794                 }
795                 isRowRoot = false;
796             }
797         }
798         if (isRowMember)
799         {
800             DDCellsizesWithDlb  &cellsizes = dd->comm->cellsizesWithDlb[d];
801             gmx::ArrayRef<real>  cellFracRow;
802
803             if (isRowRoot)
804             {
805                 RowMaster *rowMaster = cellsizes.rowMaster.get();
806                 set_dd_cell_sizes_dlb_root(dd, d, dim, rowMaster,
807                                            ddbox, bDynamicBox, bUniform, step);
808                 cellFracRow = rowMaster->cellFrac;
809             }
810             else
811             {
812                 cellFracRow = cellsizes.fracRow;
813             }
814             distribute_dd_cell_sizes_dlb(dd, d, dim, cellFracRow, ddbox);
815         }
816     }
817 }
818
819 static void set_dd_cell_sizes_dlb_nochange(gmx_domdec_t      *dd,
820                                            const gmx_ddbox_t *ddbox)
821 {
822     int d;
823
824     /* This function assumes the box is static and should therefore
825      * not be called when the box has changed since the last
826      * call to dd_partition_system.
827      */
828     for (d = 0; d < dd->ndim; d++)
829     {
830         relative_to_absolute_cell_bounds(dd, ddbox, d);
831     }
832 }
833
834
835
836 static void set_dd_cell_sizes_dlb(gmx_domdec_t *dd,
837                                   const gmx_ddbox_t *ddbox, gmx_bool bDynamicBox,
838                                   gmx_bool bUniform, gmx_bool bDoDLB, int64_t step,
839                                   gmx_wallcycle_t wcycle)
840 {
841     gmx_domdec_comm_t *comm;
842     int                dim;
843
844     comm = dd->comm;
845
846     if (bDoDLB)
847     {
848         wallcycle_start(wcycle, ewcDDCOMMBOUND);
849         set_dd_cell_sizes_dlb_change(dd, ddbox, bDynamicBox, bUniform, step);
850         wallcycle_stop(wcycle, ewcDDCOMMBOUND);
851     }
852     else if (bDynamicBox)
853     {
854         set_dd_cell_sizes_dlb_nochange(dd, ddbox);
855     }
856
857     /* Set the dimensions for which no DD is used */
858     for (dim = 0; dim < DIM; dim++)
859     {
860         if (dd->nc[dim] == 1)
861         {
862             comm->cell_x0[dim] = 0;
863             comm->cell_x1[dim] = ddbox->box_size[dim];
864             if (dim >= ddbox->nboundeddim)
865             {
866                 comm->cell_x0[dim] += ddbox->box0[dim];
867                 comm->cell_x1[dim] += ddbox->box0[dim];
868             }
869         }
870     }
871 }
872
873 void set_dd_cell_sizes(gmx_domdec_t *dd,
874                        const gmx_ddbox_t *ddbox, gmx_bool bDynamicBox,
875                        gmx_bool bUniform, gmx_bool bDoDLB, int64_t step,
876                        gmx_wallcycle_t wcycle)
877 {
878     gmx_domdec_comm_t *comm = dd->comm;
879
880     /* Copy the old cell boundaries for the cg displacement check */
881     copy_rvec(comm->cell_x0, comm->old_cell_x0);
882     copy_rvec(comm->cell_x1, comm->old_cell_x1);
883
884     if (isDlbOn(comm))
885     {
886         if (DDMASTER(dd))
887         {
888             check_box_size(dd, ddbox);
889         }
890         set_dd_cell_sizes_dlb(dd, ddbox, bDynamicBox, bUniform, bDoDLB, step, wcycle);
891     }
892     else
893     {
894         ivec numPulses;
895         set_dd_cell_sizes_slb(dd, ddbox, setcellsizeslbLOCAL, numPulses);
896
897         /* Check if the change in cell size requires a different number
898          * of communication pulses and if so change the number.
899          */
900         for (int d = 0; d < dd->ndim; d++)
901         {
902             gmx_domdec_comm_dim_t &cd           = comm->cd[d];
903             int                    numPulsesDim = numPulses[dd->dim[d]];
904             if (cd.numPulses() != numPulsesDim)
905             {
906                 if (debug)
907                 {
908                     fprintf(debug, "Changing the number of halo communication pulses along dim %c from %d to %d\n",
909                             dim2char(dd->dim[d]), cd.numPulses(), numPulsesDim);
910                 }
911                 cd.ind.resize(numPulsesDim);
912             }
913         }
914     }
915
916     if (debug)
917     {
918         for (int d = 0; d < DIM; d++)
919         {
920             fprintf(debug, "cell_x[%d] %f - %f skew_fac %f\n",
921                     d, comm->cell_x0[d], comm->cell_x1[d], ddbox->skew_fac[d]);
922         }
923     }
924 }