2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "types/commrec.h"
42 #include "nbnxn_cuda_data_mgmt.h"
45 #include "md_logging.h"
46 #include "pme_loadbal.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/smalloc.h"
53 /* Parameters and setting for one PP-PME setup */
55 real rcut_coulomb; /* Coulomb cut-off */
56 real rlist; /* pair-list cut-off */
57 real rlistlong; /* LR pair-list cut-off */
58 int nstcalclr; /* frequency of evaluating long-range forces for group scheme */
59 real spacing; /* (largest) PME grid spacing */
60 ivec grid; /* the PME grid dimensions */
61 real grid_efficiency; /* ineffiency factor for non-uniform grids <= 1 */
62 real ewaldcoeff_q; /* Electrostatic Ewald coefficient */
63 real ewaldcoeff_lj; /* LJ Ewald coefficient, only for the call to send_switchgrid */
64 gmx_pme_t pmedata; /* the data structure used in the PME code */
65 int count; /* number of times this setup has been timed */
66 double cycles; /* the fastest time for this setup in cycles */
69 /* In the initial scan, step by grids that are at least a factor 0.8 coarser */
70 #define PME_LB_GRID_SCALE_FAC 0.8
71 /* In the initial scan, try to skip grids with uneven x/y/z spacing,
72 * checking if the "efficiency" is more than 5% worse than the previous grid.
74 #define PME_LB_GRID_EFFICIENCY_REL_FAC 1.05
75 /* Rerun up till 12% slower setups than the fastest up till now */
76 #define PME_LB_SLOW_FAC 1.12
77 /* If setups get more than 2% faster, do another round to avoid
78 * choosing a slower setup due to acceleration or fluctuations.
80 #define PME_LB_ACCEL_TOL 1.02
83 epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
86 const char *pmelblim_str[epmelblimNR] =
87 { "no", "box size", "domain decompostion", "PME grid restriction" };
89 struct pme_load_balancing {
90 int nstage; /* the current maximum number of stages */
92 real cut_spacing; /* the minimum cutoff / PME grid spacing ratio */
93 real rcut_vdw; /* Vdw cutoff (does not change) */
94 real rcut_coulomb_start; /* Initial electrostatics cutoff */
95 int nstcalclr_start; /* Initial electrostatics cutoff */
96 real rbuf_coulomb; /* the pairlist buffer size */
97 real rbuf_vdw; /* the pairlist buffer size */
98 matrix box_start; /* the initial simulation box */
99 int n; /* the count of setup as well as the allocation size */
100 pme_setup_t *setup; /* the PME+cutoff setups */
101 int cur; /* the current setup */
102 int fastest; /* fastest setup up till now */
103 int start; /* start of setup range to consider in stage>0 */
104 int end; /* end of setup range to consider in stage>0 */
105 int elimited; /* was the balancing limited, uses enum above */
106 int cutoff_scheme; /* Verlet or group cut-offs */
108 int stage; /* the current stage */
111 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
112 const t_inputrec *ir, matrix box,
113 const interaction_const_t *ic,
116 pme_load_balancing_t pme_lb;
122 /* Any number of stages >= 2 is supported */
125 pme_lb->cutoff_scheme = ir->cutoff_scheme;
127 if (pme_lb->cutoff_scheme == ecutsVERLET)
129 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
130 pme_lb->rbuf_vdw = pme_lb->rbuf_coulomb;
134 if (ic->rcoulomb > ic->rlist)
136 pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
140 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
142 if (ic->rvdw > ic->rlist)
144 pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
148 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
152 copy_mat(box, pme_lb->box_start);
153 if (ir->ePBC == epbcXY && ir->nwall == 2)
155 svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
159 snew(pme_lb->setup, pme_lb->n);
161 pme_lb->rcut_vdw = ic->rvdw;
162 pme_lb->rcut_coulomb_start = ir->rcoulomb;
163 pme_lb->nstcalclr_start = ir->nstcalclr;
166 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
167 pme_lb->setup[0].rlist = ic->rlist;
168 pme_lb->setup[0].rlistlong = ic->rlistlong;
169 pme_lb->setup[0].nstcalclr = ir->nstcalclr;
170 pme_lb->setup[0].grid[XX] = ir->nkx;
171 pme_lb->setup[0].grid[YY] = ir->nky;
172 pme_lb->setup[0].grid[ZZ] = ir->nkz;
173 pme_lb->setup[0].ewaldcoeff_q = ic->ewaldcoeff_q;
174 pme_lb->setup[0].ewaldcoeff_lj = ic->ewaldcoeff_lj;
176 pme_lb->setup[0].pmedata = pmedata;
179 for (d = 0; d < DIM; d++)
181 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
187 pme_lb->setup[0].spacing = spm;
189 if (ir->fourier_spacing > 0)
191 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
195 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
203 pme_lb->elimited = epmelblimNO;
208 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t pme_lb,
210 const gmx_domdec_t *dd)
213 int npmenodes_x, npmenodes_y;
215 real tmpr_coulomb, tmpr_vdw;
219 /* Try to add a new setup with next larger cut-off to the list */
221 srenew(pme_lb->setup, pme_lb->n);
222 set = &pme_lb->setup[pme_lb->n-1];
225 get_pme_nnodes(dd, &npmenodes_x, &npmenodes_y);
230 /* Avoid infinite while loop, which can occur at the minimum grid size.
231 * Note that in practice load balancing will stop before this point.
232 * The factor 2.1 allows for the extreme case in which only grids
233 * of powers of 2 are allowed (the current code supports more grids).
243 clear_ivec(set->grid);
244 sp = calc_grid(NULL, pme_lb->box_start,
245 fac*pme_lb->setup[pme_lb->cur].spacing,
250 /* As here we can't easily check if one of the PME nodes
251 * uses threading, we do a conservative grid check.
252 * This means we can't use pme_order or less grid lines
253 * per PME node along x, which is not a strong restriction.
255 gmx_pme_check_restrictions(pme_order,
256 set->grid[XX], set->grid[YY], set->grid[ZZ],
257 npmenodes_x, npmenodes_y,
262 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
264 set->rcut_coulomb = pme_lb->cut_spacing*sp;
266 if (pme_lb->cutoff_scheme == ecutsVERLET)
268 set->rlist = set->rcut_coulomb + pme_lb->rbuf_coulomb;
269 /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
270 set->rlistlong = set->rlist;
274 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
275 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
276 set->rlist = min(tmpr_coulomb, tmpr_vdw);
277 set->rlistlong = max(tmpr_coulomb, tmpr_vdw);
279 /* Set the long-range update frequency */
280 if (set->rlist == set->rlistlong)
282 /* No long-range interactions if the short-/long-range cutoffs are identical */
285 else if (pme_lb->nstcalclr_start == 0 || pme_lb->nstcalclr_start == 1)
287 /* We were not doing long-range before, but now we are since rlist!=rlistlong */
292 /* We were already doing long-range interactions from the start */
293 if (pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
295 /* We were originally doing long-range VdW-only interactions.
296 * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
297 * but if the coulomb cutoff has become longer we should update the long-range
300 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
304 /* We were not doing any long-range interaction from the start,
305 * since it is not possible to do twin-range coulomb for the PME interaction.
313 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
314 set->grid_efficiency = 1;
315 for (d = 0; d < DIM; d++)
317 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
319 /* The Ewald coefficient is inversly proportional to the cut-off */
321 pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
322 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
324 pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
331 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
332 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
337 static void print_grid(FILE *fp_err, FILE *fp_log,
340 const pme_setup_t *set,
343 char buf[STRLEN], buft[STRLEN];
347 sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
353 sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
355 desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
359 fprintf(fp_err, "\r%s\n", buf);
363 fprintf(fp_log, "%s\n", buf);
367 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
369 /* In the initial stage only n is set; end is not set yet */
380 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
382 pme_load_balancing_t pme_lb)
384 char buf[STRLEN], sbuf[22];
386 sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
387 gmx_step_str(step, sbuf),
388 pmelblim_str[pme_lb->elimited],
389 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
392 fprintf(fp_err, "\r%s\n", buf);
396 fprintf(fp_log, "%s\n", buf);
400 static void switch_to_stage1(pme_load_balancing_t pme_lb)
403 while (pme_lb->start+1 < pme_lb->n &&
404 (pme_lb->setup[pme_lb->start].count == 0 ||
405 pme_lb->setup[pme_lb->start].cycles >
406 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
410 while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
415 pme_lb->end = pme_lb->n;
416 if (pme_lb->setup[pme_lb->end-1].count > 0 &&
417 pme_lb->setup[pme_lb->end-1].cycles >
418 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
425 /* Next we want to choose setup pme_lb->start, but as we will increase
426 * pme_ln->cur by one right after returning, we subtract 1 here.
428 pme_lb->cur = pme_lb->start - 1;
431 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
438 interaction_const_t *ic,
439 nonbonded_verlet_t *nbv,
446 char buf[STRLEN], sbuf[22];
448 gmx_bool bUsesSimpleTables = TRUE;
450 if (pme_lb->stage == pme_lb->nstage)
457 gmx_sumd(1, &cycles, cr);
458 cycles /= cr->nnodes;
461 set = &pme_lb->setup[pme_lb->cur];
464 rtab = ir->rlistlong + ir->tabext;
466 if (set->count % 2 == 1)
468 /* Skip the first cycle, because the first step after a switch
469 * is much slower due to allocation and/or caching effects.
474 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
475 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
479 set->cycles = cycles;
483 if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
484 pme_lb->stage == pme_lb->nstage - 1)
486 /* The performance went up a lot (due to e.g. DD load balancing).
487 * Add a stage, keep the minima, but rescan all setups.
493 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
494 "Increased the number stages to %d"
495 " and ignoring the previous performance\n",
496 set->grid[XX], set->grid[YY], set->grid[ZZ],
497 cycles*1e-6, set->cycles*1e-6, PME_LB_ACCEL_TOL,
501 set->cycles = min(set->cycles, cycles);
504 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
506 pme_lb->fastest = pme_lb->cur;
508 if (DOMAINDECOMP(cr))
510 /* We found a new fastest setting, ensure that with subsequent
511 * shorter cut-off's the dynamic load balancing does not make
512 * the use of the current cut-off impossible. This solution is
513 * a trade-off, as the PME load balancing and DD domain size
514 * load balancing can interact in complex ways.
515 * With the Verlet kernels, DD load imbalance will usually be
516 * mainly due to bonded interaction imbalance, which will often
517 * quickly push the domain boundaries beyond the limit for the
518 * optimal, PME load balanced, cut-off. But it could be that
519 * better overal performance can be obtained with a slightly
520 * shorter cut-off and better DD load balancing.
522 change_dd_dlb_cutoff_limit(cr);
525 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
527 /* Check in stage 0 if we should stop scanning grids.
528 * Stop when the time is more than SLOW_FAC longer than the fastest.
530 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
531 cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
533 pme_lb->n = pme_lb->cur + 1;
534 /* Done with scanning, go to stage 1 */
535 switch_to_stage1(pme_lb);
538 if (pme_lb->stage == 0)
542 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
546 if (pme_lb->cur+1 < pme_lb->n)
548 /* We had already generated the next setup */
553 /* Find the next setup */
554 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
558 pme_lb->elimited = epmelblimPMEGRID;
562 if (OK && ir->ePBC != epbcNONE)
564 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
565 <= max_cutoff2(ir->ePBC, state->box));
568 pme_lb->elimited = epmelblimBOX;
576 if (DOMAINDECOMP(cr))
578 OK = change_dd_cutoff(cr, state, ir,
579 pme_lb->setup[pme_lb->cur].rlistlong);
582 /* Failed: do not use this setup */
584 pme_lb->elimited = epmelblimDD;
590 /* We hit the upper limit for the cut-off,
591 * the setup should not go further than cur.
593 pme_lb->n = pme_lb->cur + 1;
594 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
595 /* Switch to the next stage */
596 switch_to_stage1(pme_lb);
600 !(pme_lb->setup[pme_lb->cur].grid[XX]*
601 pme_lb->setup[pme_lb->cur].grid[YY]*
602 pme_lb->setup[pme_lb->cur].grid[ZZ] <
603 gridsize_start*PME_LB_GRID_SCALE_FAC
605 pme_lb->setup[pme_lb->cur].grid_efficiency <
606 pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
609 if (pme_lb->stage > 0 && pme_lb->end == 1)
612 pme_lb->stage = pme_lb->nstage;
614 else if (pme_lb->stage > 0 && pme_lb->end > 1)
616 /* If stage = nstage-1:
617 * scan over all setups, rerunning only those setups
618 * which are not much slower than the fastest
625 if (pme_lb->cur == pme_lb->end)
628 pme_lb->cur = pme_lb->start;
631 while (pme_lb->stage == pme_lb->nstage - 1 &&
632 pme_lb->setup[pme_lb->cur].count > 0 &&
633 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
635 if (pme_lb->stage == pme_lb->nstage)
637 /* We are done optimizing, use the fastest setup we found */
638 pme_lb->cur = pme_lb->fastest;
642 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
644 OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistlong);
647 /* Failsafe solution */
648 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
654 pme_lb->end = pme_lb->cur;
655 pme_lb->cur = pme_lb->start;
656 pme_lb->elimited = epmelblimDD;
657 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
661 /* Change the Coulomb cut-off and the PME grid */
663 set = &pme_lb->setup[pme_lb->cur];
665 ic->rcoulomb = set->rcut_coulomb;
666 ic->rlist = set->rlist;
667 ic->rlistlong = set->rlistlong;
668 ir->nstcalclr = set->nstcalclr;
669 ic->ewaldcoeff_q = set->ewaldcoeff_q;
670 /* TODO: centralize the code that sets the potentials shifts */
671 if (ic->coulomb_modifier == eintmodPOTSHIFT)
673 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
675 if (EVDW_PME(ic->vdwtype))
677 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
678 ic->rvdw = set->rcut_coulomb;
679 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
680 if (ic->vdw_modifier == eintmodPOTSHIFT)
684 ic->dispersion_shift.cpot = -pow(ic->rvdw, -6.0);
685 ic->repulsion_shift.cpot = -pow(ic->rvdw, -12.0);
686 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
687 crc2 = sqr(ic->ewaldcoeff_lj*ic->rvdw);
688 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, -6.0);
692 bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
693 if (pme_lb->cutoff_scheme == ecutsVERLET &&
694 nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
696 nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv, ic);
698 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
699 * also sharing texture references. To keep the code simple, we don't
700 * treat texture references as shared resources, but this means that
701 * the coulomb_tab texture ref will get updated by multiple threads.
702 * Hence, to ensure that the non-bonded kernels don't start before all
703 * texture binding operations are finished, we need to wait for all ranks
704 * to arrive here before continuing.
706 * Note that we could omit this barrier if GPUs are not shared (or
707 * texture objects are used), but as this is initialization code, there
708 * is not point in complicating things.
710 #ifdef GMX_THREAD_MPI
715 #endif /* GMX_THREAD_MPI */
718 /* Usually we won't need the simple tables with GPUs.
719 * But we do with hybrid acceleration and with free energy.
720 * To avoid bugs, we always re-initialize the simple tables here.
722 init_interaction_const_tables(NULL, ic, bUsesSimpleTables, rtab);
724 if (cr->duty & DUTY_PME)
726 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
728 /* Generate a new PME data structure,
729 * copying part of the old pointers.
731 gmx_pme_reinit(&set->pmedata,
732 cr, pme_lb->setup[0].pmedata, ir,
735 *pmedata = set->pmedata;
739 /* Tell our PME-only node to switch grid */
740 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
745 print_grid(NULL, debug, "", "switched to", set, -1);
748 if (pme_lb->stage == pme_lb->nstage)
750 print_grid(fp_err, fp_log, "", "optimal", set, -1);
756 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
761 static int pme_grid_points(const pme_setup_t *setup)
763 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
766 static real pme_loadbal_rlist(const pme_setup_t *setup)
768 /* With the group cut-off scheme we can have twin-range either
769 * for Coulomb or for VdW, so we need a check here.
770 * With the Verlet cut-off scheme rlist=rlistlong.
772 if (setup->rcut_coulomb > setup->rlist)
774 return setup->rlistlong;
782 static void print_pme_loadbal_setting(FILE *fplog,
784 const pme_setup_t *setup)
787 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
789 setup->rcut_coulomb, pme_loadbal_rlist(setup),
790 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
791 setup->spacing, 1/setup->ewaldcoeff_q);
794 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
797 gmx_bool bNonBondedOnGPU)
799 double pp_ratio, grid_ratio;
801 pp_ratio = pow(pme_loadbal_rlist(&pme_lb->setup[pme_lb->cur])/pme_loadbal_rlist(&pme_lb->setup[0]), 3.0);
802 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
803 (double)pme_grid_points(&pme_lb->setup[0]);
805 fprintf(fplog, "\n");
806 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
807 fprintf(fplog, "\n");
808 /* Here we only warn when the optimal setting is the last one */
809 if (pme_lb->elimited != epmelblimNO &&
810 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
812 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
813 pmelblim_str[pme_lb->elimited]);
814 fprintf(fplog, " you might not have reached a good load balance.\n");
815 if (pme_lb->elimited == epmelblimDD)
817 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
819 fprintf(fplog, "\n");
821 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
822 fprintf(fplog, " particle-particle PME\n");
823 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
824 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
825 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
826 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
827 pp_ratio, grid_ratio);
828 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
830 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
832 md_print_warn(cr, fplog,
833 "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
834 " For better performance, use (more) PME ranks (mdrun -npme),\n"
835 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).\n");
839 fprintf(fplog, "\n");
843 void pme_loadbal_done(pme_load_balancing_t pme_lb,
844 t_commrec *cr, FILE *fplog,
845 gmx_bool bNonBondedOnGPU)
847 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
849 print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
852 /* TODO: Here we should free all pointers in pme_lb,
853 * but as it contains pme data structures,
854 * we need to first make pme.c free all data.