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.
39 #include "gromacs/utility/smalloc.h"
40 #include "types/commrec.h"
46 #include "nbnxn_cuda_data_mgmt.h"
49 #include "md_logging.h"
50 #include "pme_loadbal.h"
52 /* Parameters and setting for one PP-PME setup */
54 real rcut_coulomb; /* Coulomb cut-off */
55 real rlist; /* pair-list cut-off */
56 real rlistlong; /* LR pair-list cut-off */
57 int nstcalclr; /* frequency of evaluating long-range forces for group scheme */
58 real spacing; /* (largest) PME grid spacing */
59 ivec grid; /* the PME grid dimensions */
60 real grid_efficiency; /* ineffiency factor for non-uniform grids <= 1 */
61 real ewaldcoeff_q; /* Electrostatic Ewald coefficient */
62 real ewaldcoeff_lj; /* LJ Ewald coefficient, only for the call to send_switchgrid */
63 gmx_pme_t pmedata; /* the data structure used in the PME code */
64 int count; /* number of times this setup has been timed */
65 double cycles; /* the fastest time for this setup in cycles */
68 /* In the initial scan, step by grids that are at least a factor 0.8 coarser */
69 #define PME_LB_GRID_SCALE_FAC 0.8
70 /* In the initial scan, try to skip grids with uneven x/y/z spacing,
71 * checking if the "efficiency" is more than 5% worse than the previous grid.
73 #define PME_LB_GRID_EFFICIENCY_REL_FAC 1.05
74 /* Rerun up till 12% slower setups than the fastest up till now */
75 #define PME_LB_SLOW_FAC 1.12
76 /* If setups get more than 2% faster, do another round to avoid
77 * choosing a slower setup due to acceleration or fluctuations.
79 #define PME_LB_ACCEL_TOL 1.02
82 epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
85 const char *pmelblim_str[epmelblimNR] =
86 { "no", "box size", "domain decompostion", "PME grid restriction" };
88 struct pme_load_balancing {
89 int nstage; /* the current maximum number of stages */
91 real cut_spacing; /* the minimum cutoff / PME grid spacing ratio */
92 real rcut_vdw; /* Vdw cutoff (does not change) */
93 real rcut_coulomb_start; /* Initial electrostatics cutoff */
94 int nstcalclr_start; /* Initial electrostatics cutoff */
95 real rbuf_coulomb; /* the pairlist buffer size */
96 real rbuf_vdw; /* the pairlist buffer size */
97 matrix box_start; /* the initial simulation box */
98 int n; /* the count of setup as well as the allocation size */
99 pme_setup_t *setup; /* the PME+cutoff setups */
100 int cur; /* the current setup */
101 int fastest; /* fastest setup up till now */
102 int start; /* start of setup range to consider in stage>0 */
103 int end; /* end of setup range to consider in stage>0 */
104 int elimited; /* was the balancing limited, uses enum above */
105 int cutoff_scheme; /* Verlet or group cut-offs */
107 int stage; /* the current stage */
110 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
111 const t_inputrec *ir, matrix box,
112 const interaction_const_t *ic,
115 pme_load_balancing_t pme_lb;
121 /* Any number of stages >= 2 is supported */
124 pme_lb->cutoff_scheme = ir->cutoff_scheme;
126 if (pme_lb->cutoff_scheme == ecutsVERLET)
128 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
129 pme_lb->rbuf_vdw = pme_lb->rbuf_coulomb;
133 if (ic->rcoulomb > ic->rlist)
135 pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
139 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
141 if (ic->rvdw > ic->rlist)
143 pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
147 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
151 copy_mat(box, pme_lb->box_start);
152 if (ir->ePBC == epbcXY && ir->nwall == 2)
154 svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
158 snew(pme_lb->setup, pme_lb->n);
160 pme_lb->rcut_vdw = ic->rvdw;
161 pme_lb->rcut_coulomb_start = ir->rcoulomb;
162 pme_lb->nstcalclr_start = ir->nstcalclr;
165 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
166 pme_lb->setup[0].rlist = ic->rlist;
167 pme_lb->setup[0].rlistlong = ic->rlistlong;
168 pme_lb->setup[0].nstcalclr = ir->nstcalclr;
169 pme_lb->setup[0].grid[XX] = ir->nkx;
170 pme_lb->setup[0].grid[YY] = ir->nky;
171 pme_lb->setup[0].grid[ZZ] = ir->nkz;
172 pme_lb->setup[0].ewaldcoeff_q = ic->ewaldcoeff_q;
173 pme_lb->setup[0].ewaldcoeff_lj = ic->ewaldcoeff_lj;
175 pme_lb->setup[0].pmedata = pmedata;
178 for (d = 0; d < DIM; d++)
180 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
186 pme_lb->setup[0].spacing = spm;
188 if (ir->fourier_spacing > 0)
190 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
194 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
202 pme_lb->elimited = epmelblimNO;
207 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t pme_lb,
209 const gmx_domdec_t *dd)
212 int npmenodes_x, npmenodes_y;
214 real tmpr_coulomb, tmpr_vdw;
218 /* Try to add a new setup with next larger cut-off to the list */
220 srenew(pme_lb->setup, pme_lb->n);
221 set = &pme_lb->setup[pme_lb->n-1];
224 get_pme_nnodes(dd, &npmenodes_x, &npmenodes_y);
229 /* Avoid infinite while loop, which can occur at the minimum grid size.
230 * Note that in practice load balancing will stop before this point.
231 * The factor 2.1 allows for the extreme case in which only grids
232 * of powers of 2 are allowed (the current code supports more grids).
242 clear_ivec(set->grid);
243 sp = calc_grid(NULL, pme_lb->box_start,
244 fac*pme_lb->setup[pme_lb->cur].spacing,
249 /* As here we can't easily check if one of the PME nodes
250 * uses threading, we do a conservative grid check.
251 * This means we can't use pme_order or less grid lines
252 * per PME node along x, which is not a strong restriction.
254 gmx_pme_check_restrictions(pme_order,
255 set->grid[XX], set->grid[YY], set->grid[ZZ],
256 npmenodes_x, npmenodes_y,
261 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
263 set->rcut_coulomb = pme_lb->cut_spacing*sp;
264 if (set->rcut_coulomb < pme_lb->rcut_coulomb_start)
266 /* This is unlikely, but can happen when e.g. continuing from
267 * a checkpoint after equilibration where the box shrank a lot.
268 * We want to avoid rcoulomb getting smaller than rvdw
269 * and there might be more issues with decreasing rcoulomb.
271 set->rcut_coulomb = pme_lb->rcut_coulomb_start;
274 if (pme_lb->cutoff_scheme == ecutsVERLET)
276 set->rlist = set->rcut_coulomb + pme_lb->rbuf_coulomb;
277 /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
278 set->rlistlong = set->rlist;
282 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
283 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
284 set->rlist = min(tmpr_coulomb, tmpr_vdw);
285 set->rlistlong = max(tmpr_coulomb, tmpr_vdw);
287 /* Set the long-range update frequency */
288 if (set->rlist == set->rlistlong)
290 /* No long-range interactions if the short-/long-range cutoffs are identical */
293 else if (pme_lb->nstcalclr_start == 0 || pme_lb->nstcalclr_start == 1)
295 /* We were not doing long-range before, but now we are since rlist!=rlistlong */
300 /* We were already doing long-range interactions from the start */
301 if (pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
303 /* We were originally doing long-range VdW-only interactions.
304 * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
305 * but if the coulomb cutoff has become longer we should update the long-range
308 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
312 /* We were not doing any long-range interaction from the start,
313 * since it is not possible to do twin-range coulomb for the PME interaction.
321 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
322 set->grid_efficiency = 1;
323 for (d = 0; d < DIM; d++)
325 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
327 /* The Ewald coefficient is inversly proportional to the cut-off */
329 pme_lb->setup[0].ewaldcoeff_q*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
330 /* We set ewaldcoeff_lj in set, even when LJ-PME is not used */
332 pme_lb->setup[0].ewaldcoeff_lj*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
339 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
340 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
345 static void print_grid(FILE *fp_err, FILE *fp_log,
348 const pme_setup_t *set,
351 char buf[STRLEN], buft[STRLEN];
355 sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
361 sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
363 desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
367 fprintf(fp_err, "\r%s\n", buf);
371 fprintf(fp_log, "%s\n", buf);
375 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
377 /* In the initial stage only n is set; end is not set yet */
388 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
390 pme_load_balancing_t pme_lb)
392 char buf[STRLEN], sbuf[22];
394 sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
395 gmx_step_str(step, sbuf),
396 pmelblim_str[pme_lb->elimited],
397 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
400 fprintf(fp_err, "\r%s\n", buf);
404 fprintf(fp_log, "%s\n", buf);
408 static void switch_to_stage1(pme_load_balancing_t pme_lb)
411 while (pme_lb->start+1 < pme_lb->n &&
412 (pme_lb->setup[pme_lb->start].count == 0 ||
413 pme_lb->setup[pme_lb->start].cycles >
414 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
418 while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
423 pme_lb->end = pme_lb->n;
424 if (pme_lb->setup[pme_lb->end-1].count > 0 &&
425 pme_lb->setup[pme_lb->end-1].cycles >
426 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
433 /* Next we want to choose setup pme_lb->start, but as we will increase
434 * pme_ln->cur by one right after returning, we subtract 1 here.
436 pme_lb->cur = pme_lb->start - 1;
439 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
446 interaction_const_t *ic,
447 nonbonded_verlet_t *nbv,
454 char buf[STRLEN], sbuf[22];
456 gmx_bool bUsesSimpleTables = TRUE;
458 if (pme_lb->stage == pme_lb->nstage)
465 gmx_sumd(1, &cycles, cr);
466 cycles /= cr->nnodes;
469 set = &pme_lb->setup[pme_lb->cur];
472 rtab = ir->rlistlong + ir->tabext;
474 if (set->count % 2 == 1)
476 /* Skip the first cycle, because the first step after a switch
477 * is much slower due to allocation and/or caching effects.
482 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
483 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
487 set->cycles = cycles;
491 if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
492 pme_lb->stage == pme_lb->nstage - 1)
494 /* The performance went up a lot (due to e.g. DD load balancing).
495 * Add a stage, keep the minima, but rescan all setups.
501 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
502 "Increased the number stages to %d"
503 " and ignoring the previous performance\n",
504 set->grid[XX], set->grid[YY], set->grid[ZZ],
505 cycles*1e-6, set->cycles*1e-6, PME_LB_ACCEL_TOL,
509 set->cycles = min(set->cycles, cycles);
512 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
514 pme_lb->fastest = pme_lb->cur;
516 if (DOMAINDECOMP(cr))
518 /* We found a new fastest setting, ensure that with subsequent
519 * shorter cut-off's the dynamic load balancing does not make
520 * the use of the current cut-off impossible. This solution is
521 * a trade-off, as the PME load balancing and DD domain size
522 * load balancing can interact in complex ways.
523 * With the Verlet kernels, DD load imbalance will usually be
524 * mainly due to bonded interaction imbalance, which will often
525 * quickly push the domain boundaries beyond the limit for the
526 * optimal, PME load balanced, cut-off. But it could be that
527 * better overal performance can be obtained with a slightly
528 * shorter cut-off and better DD load balancing.
530 change_dd_dlb_cutoff_limit(cr);
533 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
535 /* Check in stage 0 if we should stop scanning grids.
536 * Stop when the time is more than SLOW_FAC longer than the fastest.
538 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
539 cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
541 pme_lb->n = pme_lb->cur + 1;
542 /* Done with scanning, go to stage 1 */
543 switch_to_stage1(pme_lb);
546 if (pme_lb->stage == 0)
550 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
554 if (pme_lb->cur+1 < pme_lb->n)
556 /* We had already generated the next setup */
561 /* Find the next setup */
562 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
566 pme_lb->elimited = epmelblimPMEGRID;
570 if (OK && ir->ePBC != epbcNONE)
572 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
573 <= max_cutoff2(ir->ePBC, state->box));
576 pme_lb->elimited = epmelblimBOX;
584 if (DOMAINDECOMP(cr))
586 OK = change_dd_cutoff(cr, state, ir,
587 pme_lb->setup[pme_lb->cur].rlistlong);
590 /* Failed: do not use this setup */
592 pme_lb->elimited = epmelblimDD;
598 /* We hit the upper limit for the cut-off,
599 * the setup should not go further than cur.
601 pme_lb->n = pme_lb->cur + 1;
602 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
603 /* Switch to the next stage */
604 switch_to_stage1(pme_lb);
608 !(pme_lb->setup[pme_lb->cur].grid[XX]*
609 pme_lb->setup[pme_lb->cur].grid[YY]*
610 pme_lb->setup[pme_lb->cur].grid[ZZ] <
611 gridsize_start*PME_LB_GRID_SCALE_FAC
613 pme_lb->setup[pme_lb->cur].grid_efficiency <
614 pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
617 if (pme_lb->stage > 0 && pme_lb->end == 1)
620 pme_lb->stage = pme_lb->nstage;
622 else if (pme_lb->stage > 0 && pme_lb->end > 1)
624 /* If stage = nstage-1:
625 * scan over all setups, rerunning only those setups
626 * which are not much slower than the fastest
633 if (pme_lb->cur == pme_lb->end)
636 pme_lb->cur = pme_lb->start;
639 while (pme_lb->stage == pme_lb->nstage - 1 &&
640 pme_lb->setup[pme_lb->cur].count > 0 &&
641 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
643 if (pme_lb->stage == pme_lb->nstage)
645 /* We are done optimizing, use the fastest setup we found */
646 pme_lb->cur = pme_lb->fastest;
650 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
652 OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistlong);
655 /* Failsafe solution */
656 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
662 pme_lb->end = pme_lb->cur;
663 pme_lb->cur = pme_lb->start;
664 pme_lb->elimited = epmelblimDD;
665 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
669 /* Change the Coulomb cut-off and the PME grid */
671 set = &pme_lb->setup[pme_lb->cur];
673 ic->rcoulomb = set->rcut_coulomb;
674 ic->rlist = set->rlist;
675 ic->rlistlong = set->rlistlong;
676 ir->nstcalclr = set->nstcalclr;
677 ic->ewaldcoeff_q = set->ewaldcoeff_q;
678 /* TODO: centralize the code that sets the potentials shifts */
679 if (ic->coulomb_modifier == eintmodPOTSHIFT)
681 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
683 if (EVDW_PME(ic->vdwtype))
685 /* We have PME for both Coulomb and VdW, set rvdw equal to rcoulomb */
686 ic->rvdw = set->rcut_coulomb;
687 ic->ewaldcoeff_lj = set->ewaldcoeff_lj;
688 if (ic->vdw_modifier == eintmodPOTSHIFT)
692 ic->dispersion_shift.cpot = -pow(ic->rvdw, -6.0);
693 ic->repulsion_shift.cpot = -pow(ic->rvdw, -12.0);
694 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
695 crc2 = sqr(ic->ewaldcoeff_lj*ic->rvdw);
696 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, -6.0);
700 bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
701 if (pme_lb->cutoff_scheme == ecutsVERLET &&
702 nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
704 nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv, ic);
706 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
707 * also sharing texture references. To keep the code simple, we don't
708 * treat texture references as shared resources, but this means that
709 * the coulomb_tab texture ref will get updated by multiple threads.
710 * Hence, to ensure that the non-bonded kernels don't start before all
711 * texture binding operations are finished, we need to wait for all ranks
712 * to arrive here before continuing.
714 * Note that we could omit this barrier if GPUs are not shared (or
715 * texture objects are used), but as this is initialization code, there
716 * is not point in complicating things.
718 #ifdef GMX_THREAD_MPI
723 #endif /* GMX_THREAD_MPI */
726 /* Usually we won't need the simple tables with GPUs.
727 * But we do with hybrid acceleration and with free energy.
728 * To avoid bugs, we always re-initialize the simple tables here.
730 init_interaction_const_tables(NULL, ic, bUsesSimpleTables, rtab);
732 if (cr->duty & DUTY_PME)
734 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
736 /* Generate a new PME data structure,
737 * copying part of the old pointers.
739 gmx_pme_reinit(&set->pmedata,
740 cr, pme_lb->setup[0].pmedata, ir,
743 *pmedata = set->pmedata;
747 /* Tell our PME-only node to switch grid */
748 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff_q, set->ewaldcoeff_lj);
753 print_grid(NULL, debug, "", "switched to", set, -1);
756 if (pme_lb->stage == pme_lb->nstage)
758 print_grid(fp_err, fp_log, "", "optimal", set, -1);
764 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
769 static int pme_grid_points(const pme_setup_t *setup)
771 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
774 static real pme_loadbal_rlist(const pme_setup_t *setup)
776 /* With the group cut-off scheme we can have twin-range either
777 * for Coulomb or for VdW, so we need a check here.
778 * With the Verlet cut-off scheme rlist=rlistlong.
780 if (setup->rcut_coulomb > setup->rlist)
782 return setup->rlistlong;
790 static void print_pme_loadbal_setting(FILE *fplog,
792 const pme_setup_t *setup)
795 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
797 setup->rcut_coulomb, pme_loadbal_rlist(setup),
798 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
799 setup->spacing, 1/setup->ewaldcoeff_q);
802 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
805 gmx_bool bNonBondedOnGPU)
807 double pp_ratio, grid_ratio;
809 pp_ratio = pow(pme_loadbal_rlist(&pme_lb->setup[pme_lb->cur])/pme_loadbal_rlist(&pme_lb->setup[0]), 3.0);
810 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
811 (double)pme_grid_points(&pme_lb->setup[0]);
813 fprintf(fplog, "\n");
814 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
815 fprintf(fplog, "\n");
816 /* Here we only warn when the optimal setting is the last one */
817 if (pme_lb->elimited != epmelblimNO &&
818 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
820 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
821 pmelblim_str[pme_lb->elimited]);
822 fprintf(fplog, " you might not have reached a good load balance.\n");
823 if (pme_lb->elimited == epmelblimDD)
825 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
827 fprintf(fplog, "\n");
829 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
830 fprintf(fplog, " particle-particle PME\n");
831 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
832 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
833 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
834 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
835 pp_ratio, grid_ratio);
836 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
838 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
840 md_print_warn(cr, fplog,
841 "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
842 " For better performance, use (more) PME ranks (mdrun -npme),\n"
843 " or if you are beyond the scaling limit, use fewer total ranks (or nodes).\n");
847 fprintf(fplog, "\n");
851 void pme_loadbal_done(pme_load_balancing_t pme_lb,
852 t_commrec *cr, FILE *fplog,
853 gmx_bool bNonBondedOnGPU)
855 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
857 print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
860 /* TODO: Here we should free all pointers in pme_lb,
861 * but as it contains pme data structures,
862 * we need to first make pme.c free all data.