Replace rounding using int(x+0.5) with roundToInt
[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, 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->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->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->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->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->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     comm = dd->comm;
349
350     const int           ncd       = dd->nc[dim];
351
352     const bool          dimHasPbc = (dim < ddbox->npbcdim);
353
354     gmx::ArrayRef<real> cell_size = rowMaster->buf_ncd;
355
356     if (debug)
357     {
358         fprintf(debug, "enforce_limits: %d %d\n", range[0], range[1]);
359     }
360
361     /* First we need to check if the scaling does not make cells
362      * smaller than the smallest allowed size.
363      * We need to do this iteratively, since if a cell is too small,
364      * it needs to be enlarged, which makes all the other cells smaller,
365      * which could in turn make another cell smaller than allowed.
366      */
367     for (int i = range[0]; i < range[1]; i++)
368     {
369         rowMaster->isCellMin[i] = false;
370     }
371     int nmin = 0;
372     int nmin_old;
373     do
374     {
375         nmin_old = nmin;
376         /* We need the total for normalization */
377         real fac = 0;
378         for (int i = range[0]; i < range[1]; i++)
379         {
380             if (!rowMaster->isCellMin[i])
381             {
382                 fac += cell_size[i];
383             }
384         }
385         fac = ( region_size - nmin*cellsize_limit_f)/fac; /* substracting cells already set to cellsize_limit_f */
386         /* Determine the cell boundaries */
387         for (int i = range[0]; i < range[1]; i++)
388         {
389             if (!rowMaster->isCellMin[i])
390             {
391                 cell_size[i] *= fac;
392                 if (!dimHasPbc && (i == 0 || i == dd->nc[dim] - 1))
393                 {
394                     cellsize_limit_f_i = 0;
395                 }
396                 else
397                 {
398                     cellsize_limit_f_i = cellsize_limit_f;
399                 }
400                 if (cell_size[i] < cellsize_limit_f_i)
401                 {
402                     rowMaster->isCellMin[i] = true;
403                     cell_size[i]            = cellsize_limit_f_i;
404                     nmin++;
405                 }
406             }
407             rowMaster->cellFrac[i + 1] = rowMaster->cellFrac[i] + cell_size[i];
408         }
409     }
410     while (nmin > nmin_old);
411
412     const int i  = range[1] - 1;
413     cell_size[i] = rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i];
414     /* For this check we should not use DD_CELL_MARGIN,
415      * but a slightly smaller factor,
416      * since rounding could get use below the limit.
417      */
418     if (dimHasPbc &&
419         cell_size[i] < cellsize_limit_f*DD_CELL_MARGIN2/DD_CELL_MARGIN)
420     {
421         char buf[22];
422         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",
423                   gmx_step_str(step, buf),
424                   dim2char(dim), ddbox->box_size[dim], ddbox->skew_fac[dim],
425                   ncd, comm->cellsize_min[dim]);
426     }
427
428     rowMaster->dlbIsLimited = (nmin > 0) || (range[0] > 0) || (range[1] < ncd);
429
430     if (!bUniform)
431     {
432         /* Check if the boundary did not displace more than halfway
433          * each of the cells it bounds, as this could cause problems,
434          * especially when the differences between cell sizes are large.
435          * If changes are applied, they will not make cells smaller
436          * than the cut-off, as we check all the boundaries which
437          * might be affected by a change and if the old state was ok,
438          * the cells will at most be shrunk back to their old size.
439          */
440         for (int i = range[0]+1; i < range[1]; i++)
441         {
442             halfway = 0.5*(rowMaster->oldCellFrac[i] + rowMaster->oldCellFrac[i - 1]);
443             if (rowMaster->cellFrac[i] < halfway)
444             {
445                 rowMaster->cellFrac[i] = halfway;
446                 /* Check if the change also causes shifts of the next boundaries */
447                 for (int j = i + 1; j < range[1]; j++)
448                 {
449                     if (rowMaster->cellFrac[j] < rowMaster->cellFrac[j - 1] + cellsize_limit_f)
450                     {
451                         rowMaster->cellFrac[j] =  rowMaster->cellFrac[j - 1] + cellsize_limit_f;
452                     }
453                 }
454             }
455             halfway = 0.5*(rowMaster->oldCellFrac[i] + rowMaster->oldCellFrac[i + 1]);
456             if (rowMaster->cellFrac[i] > halfway)
457             {
458                 rowMaster->cellFrac[i] = halfway;
459                 /* Check if the change also causes shifts of the next boundaries */
460                 for (int j = i - 1; j >= range[0] + 1; j--)
461                 {
462                     if (rowMaster->cellFrac[j] > rowMaster->cellFrac[j + 1] - cellsize_limit_f)
463                     {
464                         rowMaster->cellFrac[j] = rowMaster->cellFrac[j + 1] - cellsize_limit_f;
465                     }
466                 }
467             }
468         }
469     }
470
471     /* nrange is defined as [lower, upper) range for new call to enforce_limits */
472     /* find highest violation of LimLo (a) and the following violation of LimHi (thus the lowest following) (b)
473      * then call enforce_limits for (oldb,a), (a,b). In the next step: (b,nexta). oldb and nexta can be the boundaries.
474      * for a and b nrange is used */
475     if (d > 0)
476     {
477         /* Take care of the staggering of the cell boundaries */
478         if (bUniform)
479         {
480             for (int i = range[0]; i < range[1]; i++)
481             {
482                 rowMaster->bounds[i].cellFracLowerMax = rowMaster->cellFrac[i];
483                 rowMaster->bounds[i].cellFracUpperMin = rowMaster->cellFrac[i + 1];
484             }
485         }
486         else
487         {
488             for (int i = range[0] + 1; i < range[1]; i++)
489             {
490                 const RowMaster::Bounds &bounds = rowMaster->bounds[i];
491
492                 bool                     bLimLo = (rowMaster->cellFrac[i] < bounds.boundMin);
493                 bool                     bLimHi = (rowMaster->cellFrac[i] > bounds.boundMax);
494                 if (bLimLo && bLimHi)
495                 {
496                     /* Both limits violated, try the best we can */
497                     /* For this case we split the original range (range) in two parts and care about the other limitiations in the next iteration. */
498                     rowMaster->cellFrac[i] = 0.5*(bounds.boundMin + bounds.boundMax);
499                     nrange[0]              = range[0];
500                     nrange[1]              = i;
501                     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
502
503                     nrange[0] = i;
504                     nrange[1] = range[1];
505                     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
506
507                     return;
508                 }
509                 else if (bLimLo)
510                 {
511                     /* rowMaster->cellFrac[i] = rowMaster->boundMin[i]; */
512                     nrange[1] = i;  /* only store violation location. There could be a LimLo violation following with an higher index */
513                     bLastHi   = FALSE;
514                 }
515                 else if (bLimHi && !bLastHi)
516                 {
517                     bLastHi = TRUE;
518                     if (nrange[1] < range[1])   /* found a LimLo before */
519                     {
520                         rowMaster->cellFrac[nrange[1]] = rowMaster->bounds[1].boundMin;
521                         dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
522                         nrange[0] = nrange[1];
523                     }
524                     rowMaster->cellFrac[i] = rowMaster->bounds[i].boundMax;
525                     nrange[1]              = i;
526                     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
527                     nrange[0] = i;
528                     nrange[1] = range[1];
529                 }
530             }
531             if (nrange[1] < range[1])   /* found last a LimLo */
532             {
533                 rowMaster->cellFrac[nrange[1]] = rowMaster->bounds[nrange[1]].boundMin;
534                 dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
535                 nrange[0] = nrange[1];
536                 nrange[1] = range[1];
537                 dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
538             }
539             else if (nrange[0] > range[0]) /* found at least one LimHi */
540             {
541                 dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, nrange);
542             }
543         }
544     }
545 }
546
547
548 static void set_dd_cell_sizes_dlb_root(gmx_domdec_t *dd,
549                                        int d, int dim, RowMaster *rowMaster,
550                                        const gmx_ddbox_t *ddbox,
551                                        gmx_bool bDynamicBox,
552                                        gmx_bool bUniform, int64_t step)
553 {
554     gmx_domdec_comm_t *comm    = dd->comm;
555     constexpr real     c_relax = 0.5;
556     int                range[] = { 0, 0 };
557
558     /* Convert the maximum change from the input percentage to a fraction */
559     const real          change_limit = comm->dlb_scale_lim*0.01;
560
561     const int           ncd          = dd->nc[dim];
562
563     const bool          bPBC = (dim < ddbox->npbcdim);
564
565     gmx::ArrayRef<real> cell_size = rowMaster->buf_ncd;
566
567     /* Store the original boundaries */
568     for (int i = 0; i < ncd + 1; i++)
569     {
570         rowMaster->oldCellFrac[i] = rowMaster->cellFrac[i];
571     }
572     if (bUniform)
573     {
574         for (int i = 0; i < ncd; i++)
575         {
576             cell_size[i] = 1.0/ncd;
577         }
578     }
579     else if (dd_load_count(comm) > 0)
580     {
581         real load_aver  = comm->load[d].sum_m/ncd;
582         real change_max = 0;
583         real load_i;
584         real change;
585         for (int i = 0; i < ncd; i++)
586         {
587             /* Determine the relative imbalance of cell i */
588             load_i         = comm->load[d].load[i*comm->load[d].nload+2];
589             real imbalance = (load_i - load_aver)/(load_aver > 0 ? load_aver : 1);
590             /* Determine the change of the cell size using underrelaxation */
591             change     = -c_relax*imbalance;
592             change_max = std::max(change_max, std::max(change, -change));
593         }
594         /* Limit the amount of scaling.
595          * We need to use the same rescaling for all cells in one row,
596          * otherwise the load balancing might not converge.
597          */
598         real sc = c_relax;
599         if (change_max > change_limit)
600         {
601             sc *= change_limit/change_max;
602         }
603         for (int i = 0; i < ncd; i++)
604         {
605             /* Determine the relative imbalance of cell i */
606             load_i         = comm->load[d].load[i*comm->load[d].nload+2];
607             real imbalance = (load_i - load_aver)/(load_aver > 0 ? load_aver : 1);
608             /* Determine the change of the cell size using underrelaxation */
609             change       = -sc*imbalance;
610             cell_size[i] = (rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i])*(1 + change);
611         }
612     }
613
614     real cellsize_limit_f  = cellsize_min_dlb(comm, d, dim)/ddbox->box_size[dim];
615     cellsize_limit_f      *= DD_CELL_MARGIN;
616     real dist_min_f_hard   = grid_jump_limit(comm, comm->cutoff, d)/ddbox->box_size[dim];
617     real dist_min_f        = dist_min_f_hard * DD_CELL_MARGIN;
618     if (ddbox->tric_dir[dim])
619     {
620         cellsize_limit_f /= ddbox->skew_fac[dim];
621         dist_min_f       /= ddbox->skew_fac[dim];
622     }
623     if (bDynamicBox && d > 0)
624     {
625         dist_min_f *= DD_PRES_SCALE_MARGIN;
626     }
627     if (d > 0 && !bUniform)
628     {
629         /* Make sure that the grid is not shifted too much */
630         for (int i = 1; i < ncd; i++)
631         {
632             const RowMaster::Bounds &boundsNeighbor = rowMaster->bounds[i - 1];
633             RowMaster::Bounds       &bounds         = rowMaster->bounds[i];
634
635             if (bounds.cellFracUpperMin - boundsNeighbor.cellFracLowerMax < 2 * dist_min_f_hard)
636             {
637                 gmx_incons("Inconsistent DD boundary staggering limits!");
638             }
639             bounds.boundMin = boundsNeighbor.cellFracLowerMax + dist_min_f;
640             real space       = rowMaster->cellFrac[i] - (boundsNeighbor.cellFracLowerMax + dist_min_f);
641             if (space > 0)
642             {
643                 bounds.boundMin += 0.5*space;
644             }
645             bounds.boundMax  = bounds.cellFracUpperMin - dist_min_f;
646             space            = rowMaster->cellFrac[i] - (bounds.cellFracUpperMin - dist_min_f);
647             if (space < 0)
648             {
649                 rowMaster->bounds[i].boundMax += 0.5*space;
650             }
651             if (debug)
652             {
653                 fprintf(debug,
654                         "dim %d boundary %d %.3f < %.3f < %.3f < %.3f < %.3f\n",
655                         d, i,
656                         boundsNeighbor.cellFracLowerMax + dist_min_f,
657                         bounds.boundMin, rowMaster->cellFrac[i], bounds.boundMax,
658                         bounds.cellFracUpperMin - dist_min_f);
659             }
660         }
661     }
662     range[1]                 = ncd;
663     rowMaster->cellFrac[0]   = 0;
664     rowMaster->cellFrac[ncd] = 1;
665     dd_cell_sizes_dlb_root_enforce_limits(dd, d, dim, rowMaster, ddbox, bUniform, step, cellsize_limit_f, range);
666
667
668     /* After the checks above, the cells should obey the cut-off
669      * restrictions, but it does not hurt to check.
670      */
671     for (int i = 0; i < ncd; i++)
672     {
673         if (debug)
674         {
675             fprintf(debug, "Relative bounds dim %d  cell %d: %f %f\n",
676                     dim, i, rowMaster->cellFrac[i], rowMaster->cellFrac[i + 1]);
677         }
678
679         if ((bPBC || (i != 0 && i != dd->nc[dim]-1)) &&
680             rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i] <
681             cellsize_limit_f/DD_CELL_MARGIN)
682         {
683             char buf[22];
684             fprintf(stderr,
685                     "\nWARNING step %s: direction %c, cell %d too small: %f\n",
686                     gmx_step_str(step, buf), dim2char(dim), i,
687                     (rowMaster->cellFrac[i + 1] - rowMaster->cellFrac[i])
688                     *ddbox->box_size[dim]*ddbox->skew_fac[dim]);
689         }
690     }
691
692     int pos = ncd + 1;
693     /* Store the cell boundaries of the lower dimensions at the end */
694     for (int d1 = 0; d1 < d; d1++)
695     {
696         rowMaster->cellFrac[pos++] = comm->cellsizesWithDlb[d1].fracLower;
697         rowMaster->cellFrac[pos++] = comm->cellsizesWithDlb[d1].fracUpper;
698     }
699
700     if (d < comm->npmedecompdim)
701     {
702         /* The master determines the maximum shift for
703          * the coordinate communication between separate PME nodes.
704          */
705         set_pme_maxshift(dd, &comm->ddpme[d], bUniform, ddbox, rowMaster->cellFrac.data());
706     }
707     rowMaster->cellFrac[pos++] = comm->ddpme[0].maxshift;
708     if (d >= 1)
709     {
710         rowMaster->cellFrac[pos++] = comm->ddpme[1].maxshift;
711     }
712 }
713
714 static void relative_to_absolute_cell_bounds(gmx_domdec_t      *dd,
715                                              const gmx_ddbox_t *ddbox,
716                                              int                dimind)
717 {
718     gmx_domdec_comm_t        *comm      = dd->comm;
719     const DDCellsizesWithDlb &cellsizes = comm->cellsizesWithDlb[dimind];
720
721     /* Set the cell dimensions */
722     int dim            = dd->dim[dimind];
723     comm->cell_x0[dim] = cellsizes.fracLower*ddbox->box_size[dim];
724     comm->cell_x1[dim] = cellsizes.fracUpper*ddbox->box_size[dim];
725     if (dim >= ddbox->nboundeddim)
726     {
727         comm->cell_x0[dim] += ddbox->box0[dim];
728         comm->cell_x1[dim] += ddbox->box0[dim];
729     }
730 }
731
732 static void distribute_dd_cell_sizes_dlb(gmx_domdec_t *dd,
733                                          int d, int dim,
734                                          gmx::ArrayRef<real> cellFracRow,
735                                          const gmx_ddbox_t *ddbox)
736 {
737     gmx_domdec_comm_t &comm = *dd->comm;
738
739 #if GMX_MPI
740     /* Each node would only need to know two fractions,
741      * but it is probably cheaper to broadcast the whole array.
742      */
743     MPI_Bcast(cellFracRow.data(), ddCellFractionBufferSize(dd, d)*sizeof(real), MPI_BYTE,
744               0, comm.mpi_comm_load[d]);
745 #endif
746     /* Copy the fractions for this dimension from the buffer */
747     comm.cellsizesWithDlb[d].fracLower = cellFracRow[dd->ci[dim]    ];
748     comm.cellsizesWithDlb[d].fracUpper = cellFracRow[dd->ci[dim] + 1];
749     /* The whole array was communicated, so set the buffer position */
750     int pos = dd->nc[dim] + 1;
751     for (int d1 = 0; d1 <= d; d1++)
752     {
753         if (d1 < d)
754         {
755             /* Copy the cell fractions of the lower dimensions */
756             comm.cellsizesWithDlb[d1].fracLower = cellFracRow[pos++];
757             comm.cellsizesWithDlb[d1].fracUpper = cellFracRow[pos++];
758         }
759         relative_to_absolute_cell_bounds(dd, ddbox, d1);
760     }
761     /* Convert the communicated shift from float to int */
762     comm.ddpme[0].maxshift = gmx::roundToInt(cellFracRow[pos++]);
763     if (d >= 1)
764     {
765         comm.ddpme[1].maxshift = gmx::roundToInt(cellFracRow[pos++]);
766     }
767 }
768
769 static void set_dd_cell_sizes_dlb_change(gmx_domdec_t *dd,
770                                          const gmx_ddbox_t *ddbox,
771                                          gmx_bool bDynamicBox,
772                                          gmx_bool bUniform, int64_t step)
773 {
774     for (int d = 0; d < dd->ndim; d++)
775     {
776         const int dim         = dd->dim[d];
777         bool      isRowMember = true;
778         bool      isRowRoot   = true;
779         for (int d1 = d; d1 < dd->ndim; d1++)
780         {
781             if (dd->ci[dd->dim[d1]] > 0)
782             {
783                 if (d1 != d)
784                 {
785                     isRowMember = false;
786                 }
787                 isRowRoot = false;
788             }
789         }
790         if (isRowMember)
791         {
792             DDCellsizesWithDlb  &cellsizes = dd->comm->cellsizesWithDlb[d];
793             gmx::ArrayRef<real>  cellFracRow;
794
795             if (isRowRoot)
796             {
797                 RowMaster *rowMaster = cellsizes.rowMaster.get();
798                 set_dd_cell_sizes_dlb_root(dd, d, dim, rowMaster,
799                                            ddbox, bDynamicBox, bUniform, step);
800                 cellFracRow = rowMaster->cellFrac;
801             }
802             else
803             {
804                 cellFracRow = cellsizes.fracRow;
805             }
806             distribute_dd_cell_sizes_dlb(dd, d, dim, cellFracRow, ddbox);
807         }
808     }
809 }
810
811 static void set_dd_cell_sizes_dlb_nochange(gmx_domdec_t      *dd,
812                                            const gmx_ddbox_t *ddbox)
813 {
814     int d;
815
816     /* This function assumes the box is static and should therefore
817      * not be called when the box has changed since the last
818      * call to dd_partition_system.
819      */
820     for (d = 0; d < dd->ndim; d++)
821     {
822         relative_to_absolute_cell_bounds(dd, ddbox, d);
823     }
824 }
825
826
827
828 static void set_dd_cell_sizes_dlb(gmx_domdec_t *dd,
829                                   const gmx_ddbox_t *ddbox, gmx_bool bDynamicBox,
830                                   gmx_bool bUniform, gmx_bool bDoDLB, int64_t step,
831                                   gmx_wallcycle_t wcycle)
832 {
833     gmx_domdec_comm_t *comm;
834     int                dim;
835
836     comm = dd->comm;
837
838     if (bDoDLB)
839     {
840         wallcycle_start(wcycle, ewcDDCOMMBOUND);
841         set_dd_cell_sizes_dlb_change(dd, ddbox, bDynamicBox, bUniform, step);
842         wallcycle_stop(wcycle, ewcDDCOMMBOUND);
843     }
844     else if (bDynamicBox)
845     {
846         set_dd_cell_sizes_dlb_nochange(dd, ddbox);
847     }
848
849     /* Set the dimensions for which no DD is used */
850     for (dim = 0; dim < DIM; dim++)
851     {
852         if (dd->nc[dim] == 1)
853         {
854             comm->cell_x0[dim] = 0;
855             comm->cell_x1[dim] = ddbox->box_size[dim];
856             if (dim >= ddbox->nboundeddim)
857             {
858                 comm->cell_x0[dim] += ddbox->box0[dim];
859                 comm->cell_x1[dim] += ddbox->box0[dim];
860             }
861         }
862     }
863 }
864
865 void set_dd_cell_sizes(gmx_domdec_t *dd,
866                        const gmx_ddbox_t *ddbox, gmx_bool bDynamicBox,
867                        gmx_bool bUniform, gmx_bool bDoDLB, int64_t step,
868                        gmx_wallcycle_t wcycle)
869 {
870     gmx_domdec_comm_t *comm = dd->comm;
871
872     /* Copy the old cell boundaries for the cg displacement check */
873     copy_rvec(comm->cell_x0, comm->old_cell_x0);
874     copy_rvec(comm->cell_x1, comm->old_cell_x1);
875
876     if (isDlbOn(comm))
877     {
878         if (DDMASTER(dd))
879         {
880             check_box_size(dd, ddbox);
881         }
882         set_dd_cell_sizes_dlb(dd, ddbox, bDynamicBox, bUniform, bDoDLB, step, wcycle);
883     }
884     else
885     {
886         ivec numPulses;
887         set_dd_cell_sizes_slb(dd, ddbox, setcellsizeslbLOCAL, numPulses);
888
889         /* Check if the change in cell size requires a different number
890          * of communication pulses and if so change the number.
891          */
892         for (int d = 0; d < dd->ndim; d++)
893         {
894             gmx_domdec_comm_dim_t &cd           = comm->cd[d];
895             int                    numPulsesDim = numPulses[dd->dim[d]];
896             if (cd.numPulses() != numPulsesDim)
897             {
898                 if (debug)
899                 {
900                     fprintf(debug, "Changing the number of halo communication pulses along dim %c from %d to %d\n",
901                             dim2char(dd->dim[d]), cd.numPulses(), numPulsesDim);
902                 }
903                 cd.ind.resize(numPulsesDim);
904             }
905         }
906     }
907
908     if (debug)
909     {
910         for (int d = 0; d < DIM; d++)
911         {
912             fprintf(debug, "cell_x[%d] %f - %f skew_fac %f\n",
913                     d, comm->cell_x0[d], comm->cell_x1[d], ddbox->skew_fac[d]);
914         }
915     }
916 }