2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, 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.
45 #include "nbnxn_cuda_data_mgmt.h"
48 #include "md_logging.h"
49 #include "pme_loadbal.h"
51 /* Parameters and setting for one PP-PME setup */
53 real rcut_coulomb; /* Coulomb cut-off */
54 real rlist; /* pair-list cut-off */
55 real rlistlong; /* LR pair-list cut-off */
56 int nstcalclr; /* frequency of evaluating long-range forces for group scheme */
57 real spacing; /* (largest) PME grid spacing */
58 ivec grid; /* the PME grid dimensions */
59 real grid_efficiency; /* ineffiency factor for non-uniform grids <= 1 */
60 real ewaldcoeff; /* the Ewald coefficient */
61 gmx_pme_t pmedata; /* the data structure used in the PME code */
63 int count; /* number of times this setup has been timed */
64 double cycles; /* the fastest time for this setup in cycles */
67 /* In the initial scan, step by grids that are at least a factor 0.8 coarser */
68 #define PME_LB_GRID_SCALE_FAC 0.8
69 /* In the initial scan, try to skip grids with uneven x/y/z spacing,
70 * checking if the "efficiency" is more than 5% worse than the previous grid.
72 #define PME_LB_GRID_EFFICIENCY_REL_FAC 1.05
73 /* Rerun up till 12% slower setups than the fastest up till now */
74 #define PME_LB_SLOW_FAC 1.12
75 /* If setups get more than 2% faster, do another round to avoid
76 * choosing a slower setup due to acceleration or fluctuations.
78 #define PME_LB_ACCEL_TOL 1.02
81 epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimPMEGRID, epmelblimNR
84 const char *pmelblim_str[epmelblimNR] =
85 { "no", "box size", "domain decompostion", "PME grid restriction" };
87 struct pme_load_balancing {
88 int nstage; /* the current maximum number of stages */
90 real cut_spacing; /* the minimum cutoff / PME grid spacing ratio */
91 real rcut_vdw; /* Vdw cutoff (does not change) */
92 real rcut_coulomb_start; /* Initial electrostatics cutoff */
93 int nstcalclr_start; /* Initial electrostatics cutoff */
94 real rbuf_coulomb; /* the pairlist buffer size */
95 real rbuf_vdw; /* the pairlist buffer size */
96 matrix box_start; /* the initial simulation box */
97 int n; /* the count of setup as well as the allocation size */
98 pme_setup_t *setup; /* the PME+cutoff setups */
99 int cur; /* the current setup */
100 int fastest; /* fastest setup up till now */
101 int start; /* start of setup range to consider in stage>0 */
102 int end; /* end of setup range to consider in stage>0 */
103 int elimited; /* was the balancing limited, uses enum above */
104 int cutoff_scheme; /* Verlet or group cut-offs */
106 int stage; /* the current stage */
109 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
110 const t_inputrec *ir, matrix box,
111 const interaction_const_t *ic,
114 pme_load_balancing_t pme_lb;
120 /* Any number of stages >= 2 is supported */
123 pme_lb->cutoff_scheme = ir->cutoff_scheme;
125 if (pme_lb->cutoff_scheme == ecutsVERLET)
127 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
128 pme_lb->rbuf_vdw = pme_lb->rbuf_coulomb;
132 if (ic->rcoulomb > ic->rlist)
134 pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
138 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
140 if (ic->rvdw > ic->rlist)
142 pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
146 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
150 copy_mat(box, pme_lb->box_start);
151 if (ir->ePBC == epbcXY && ir->nwall == 2)
153 svmul(ir->wall_ewald_zfac, pme_lb->box_start[ZZ], pme_lb->box_start[ZZ]);
157 snew(pme_lb->setup, pme_lb->n);
159 pme_lb->rcut_vdw = ic->rvdw;
160 pme_lb->rcut_coulomb_start = ir->rcoulomb;
161 pme_lb->nstcalclr_start = ir->nstcalclr;
164 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
165 pme_lb->setup[0].rlist = ic->rlist;
166 pme_lb->setup[0].rlistlong = ic->rlistlong;
167 pme_lb->setup[0].nstcalclr = ir->nstcalclr;
168 pme_lb->setup[0].grid[XX] = ir->nkx;
169 pme_lb->setup[0].grid[YY] = ir->nky;
170 pme_lb->setup[0].grid[ZZ] = ir->nkz;
171 pme_lb->setup[0].ewaldcoeff = ic->ewaldcoeff;
173 pme_lb->setup[0].pmedata = pmedata;
176 for (d = 0; d < DIM; d++)
178 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
184 pme_lb->setup[0].spacing = spm;
186 if (ir->fourier_spacing > 0)
188 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
192 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
200 pme_lb->elimited = epmelblimNO;
205 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t pme_lb,
207 const gmx_domdec_t *dd)
210 int npmenodes_x, npmenodes_y;
212 real tmpr_coulomb, tmpr_vdw;
216 /* Try to add a new setup with next larger cut-off to the list */
218 srenew(pme_lb->setup, pme_lb->n);
219 set = &pme_lb->setup[pme_lb->n-1];
222 get_pme_nnodes(dd, &npmenodes_x, &npmenodes_y);
227 /* Avoid infinite while loop, which can occur at the minimum grid size.
228 * Note that in practice load balancing will stop before this point.
229 * The factor 2.1 allows for the extreme case in which only grids
230 * of powers of 2 are allowed (the current code supports more grids).
240 clear_ivec(set->grid);
241 sp = calc_grid(NULL, pme_lb->box_start,
242 fac*pme_lb->setup[pme_lb->cur].spacing,
247 /* As here we can't easily check if one of the PME nodes
248 * uses threading, we do a conservative grid check.
249 * This means we can't use pme_order or less grid lines
250 * per PME node along x, which is not a strong restriction.
252 gmx_pme_check_restrictions(pme_order,
253 set->grid[XX], set->grid[YY], set->grid[ZZ],
254 npmenodes_x, npmenodes_y,
259 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing || !grid_ok);
261 set->rcut_coulomb = pme_lb->cut_spacing*sp;
263 if (pme_lb->cutoff_scheme == ecutsVERLET)
265 set->rlist = set->rcut_coulomb + pme_lb->rbuf_coulomb;
266 /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
267 set->rlistlong = set->rlist;
271 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
272 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
273 set->rlist = min(tmpr_coulomb, tmpr_vdw);
274 set->rlistlong = max(tmpr_coulomb, tmpr_vdw);
276 /* Set the long-range update frequency */
277 if (set->rlist == set->rlistlong)
279 /* No long-range interactions if the short-/long-range cutoffs are identical */
282 else if (pme_lb->nstcalclr_start == 0 || pme_lb->nstcalclr_start == 1)
284 /* We were not doing long-range before, but now we are since rlist!=rlistlong */
289 /* We were already doing long-range interactions from the start */
290 if (pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
292 /* We were originally doing long-range VdW-only interactions.
293 * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
294 * but if the coulomb cutoff has become longer we should update the long-range
297 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
301 /* We were not doing any long-range interaction from the start,
302 * since it is not possible to do twin-range coulomb for the PME interaction.
310 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
311 set->grid_efficiency = 1;
312 for (d = 0; d < DIM; d++)
314 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
316 /* The Ewald coefficient is inversly proportional to the cut-off */
318 pme_lb->setup[0].ewaldcoeff*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
325 fprintf(debug, "PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
326 set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb);
331 static void print_grid(FILE *fp_err, FILE *fp_log,
334 const pme_setup_t *set,
337 char buf[STRLEN], buft[STRLEN];
341 sprintf(buft, ": %.1f M-cycles", cycles*1e-6);
347 sprintf(buf, "%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
349 desc, set->grid[XX], set->grid[YY], set->grid[ZZ], set->rcut_coulomb,
353 fprintf(fp_err, "\r%s\n", buf);
357 fprintf(fp_log, "%s\n", buf);
361 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
363 /* In the initial stage only n is set; end is not set yet */
374 static void print_loadbal_limited(FILE *fp_err, FILE *fp_log,
375 gmx_large_int_t step,
376 pme_load_balancing_t pme_lb)
378 char buf[STRLEN], sbuf[22];
380 sprintf(buf, "step %4s: the %s limits the PME load balancing to a coulomb cut-off of %.3f",
381 gmx_step_str(step, sbuf),
382 pmelblim_str[pme_lb->elimited],
383 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
386 fprintf(fp_err, "\r%s\n", buf);
390 fprintf(fp_log, "%s\n", buf);
394 static void switch_to_stage1(pme_load_balancing_t pme_lb)
397 while (pme_lb->start+1 < pme_lb->n &&
398 (pme_lb->setup[pme_lb->start].count == 0 ||
399 pme_lb->setup[pme_lb->start].cycles >
400 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
404 while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
409 pme_lb->end = pme_lb->n;
410 if (pme_lb->setup[pme_lb->end-1].count > 0 &&
411 pme_lb->setup[pme_lb->end-1].cycles >
412 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
419 /* Next we want to choose setup pme_lb->start, but as we will increase
420 * pme_ln->cur by one right after returning, we subtract 1 here.
422 pme_lb->cur = pme_lb->start - 1;
425 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
432 interaction_const_t *ic,
433 nonbonded_verlet_t *nbv,
435 gmx_large_int_t step)
440 char buf[STRLEN], sbuf[22];
442 gmx_bool bUsesSimpleTables = TRUE;
444 if (pme_lb->stage == pme_lb->nstage)
451 gmx_sumd(1, &cycles, cr);
452 cycles /= cr->nnodes;
455 set = &pme_lb->setup[pme_lb->cur];
458 rtab = ir->rlistlong + ir->tabext;
460 if (set->count % 2 == 1)
462 /* Skip the first cycle, because the first step after a switch
463 * is much slower due to allocation and/or caching effects.
468 sprintf(buf, "step %4s: ", gmx_step_str(step, sbuf));
469 print_grid(fp_err, fp_log, buf, "timed with", set, cycles);
473 set->cycles = cycles;
477 if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
478 pme_lb->stage == pme_lb->nstage - 1)
480 /* The performance went up a lot (due to e.g. DD load balancing).
481 * Add a stage, keep the minima, but rescan all setups.
487 fprintf(debug, "The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
488 "Increased the number stages to %d"
489 " and ignoring the previous performance\n",
490 set->grid[XX], set->grid[YY], set->grid[ZZ],
491 cycles*1e-6, set->cycles*1e-6, PME_LB_ACCEL_TOL,
495 set->cycles = min(set->cycles, cycles);
498 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
500 pme_lb->fastest = pme_lb->cur;
502 if (DOMAINDECOMP(cr))
504 /* We found a new fastest setting, ensure that with subsequent
505 * shorter cut-off's the dynamic load balancing does not make
506 * the use of the current cut-off impossible. This solution is
507 * a trade-off, as the PME load balancing and DD domain size
508 * load balancing can interact in complex ways.
509 * With the Verlet kernels, DD load imbalance will usually be
510 * mainly due to bonded interaction imbalance, which will often
511 * quickly push the domain boundaries beyond the limit for the
512 * optimal, PME load balanced, cut-off. But it could be that
513 * better overal performance can be obtained with a slightly
514 * shorter cut-off and better DD load balancing.
516 change_dd_dlb_cutoff_limit(cr);
519 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
521 /* Check in stage 0 if we should stop scanning grids.
522 * Stop when the time is more than SLOW_FAC longer than the fastest.
524 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
525 cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
527 pme_lb->n = pme_lb->cur + 1;
528 /* Done with scanning, go to stage 1 */
529 switch_to_stage1(pme_lb);
532 if (pme_lb->stage == 0)
536 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
540 if (pme_lb->cur+1 < pme_lb->n)
542 /* We had already generated the next setup */
547 /* Find the next setup */
548 OK = pme_loadbal_increase_cutoff(pme_lb, ir->pme_order, cr->dd);
552 pme_lb->elimited = epmelblimPMEGRID;
556 if (OK && ir->ePBC != epbcNONE)
558 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
559 <= max_cutoff2(ir->ePBC, state->box));
562 pme_lb->elimited = epmelblimBOX;
570 if (DOMAINDECOMP(cr))
572 OK = change_dd_cutoff(cr, state, ir,
573 pme_lb->setup[pme_lb->cur].rlistlong);
576 /* Failed: do not use this setup */
578 pme_lb->elimited = epmelblimDD;
584 /* We hit the upper limit for the cut-off,
585 * the setup should not go further than cur.
587 pme_lb->n = pme_lb->cur + 1;
588 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
589 /* Switch to the next stage */
590 switch_to_stage1(pme_lb);
594 !(pme_lb->setup[pme_lb->cur].grid[XX]*
595 pme_lb->setup[pme_lb->cur].grid[YY]*
596 pme_lb->setup[pme_lb->cur].grid[ZZ] <
597 gridsize_start*PME_LB_GRID_SCALE_FAC
599 pme_lb->setup[pme_lb->cur].grid_efficiency <
600 pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
603 if (pme_lb->stage > 0 && pme_lb->end == 1)
606 pme_lb->stage = pme_lb->nstage;
608 else if (pme_lb->stage > 0 && pme_lb->end > 1)
610 /* If stage = nstage-1:
611 * scan over all setups, rerunning only those setups
612 * which are not much slower than the fastest
619 if (pme_lb->cur == pme_lb->end)
622 pme_lb->cur = pme_lb->start;
625 while (pme_lb->stage == pme_lb->nstage - 1 &&
626 pme_lb->setup[pme_lb->cur].count > 0 &&
627 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
629 if (pme_lb->stage == pme_lb->nstage)
631 /* We are done optimizing, use the fastest setup we found */
632 pme_lb->cur = pme_lb->fastest;
636 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
638 OK = change_dd_cutoff(cr, state, ir, pme_lb->setup[pme_lb->cur].rlistlong);
641 /* Failsafe solution */
642 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
648 pme_lb->end = pme_lb->cur;
649 pme_lb->cur = pme_lb->start;
650 pme_lb->elimited = epmelblimDD;
651 print_loadbal_limited(fp_err, fp_log, step, pme_lb);
655 /* Change the Coulomb cut-off and the PME grid */
657 set = &pme_lb->setup[pme_lb->cur];
659 ic->rcoulomb = set->rcut_coulomb;
660 ic->rlist = set->rlist;
661 ic->rlistlong = set->rlistlong;
662 ir->nstcalclr = set->nstcalclr;
663 ic->ewaldcoeff = set->ewaldcoeff;
665 bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
666 if (pme_lb->cutoff_scheme == ecutsVERLET &&
667 nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
669 nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv, ic);
671 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
672 * also sharing texture references. To keep the code simple, we don't
673 * treat texture references as shared resources, but this means that
674 * the coulomb_tab texture ref will get updated by multiple threads.
675 * Hence, to ensure that the non-bonded kernels don't start before all
676 * texture binding operations are finished, we need to wait for all ranks
677 * to arrive here before continuing.
679 * Note that we could omit this barrier if GPUs are not shared (or
680 * texture objects are used), but as this is initialization code, there
681 * is not point in complicating things.
683 #ifdef GMX_THREAD_MPI
688 #endif /* GMX_THREAD_MPI */
692 init_interaction_const_tables(NULL, ic, bUsesSimpleTables,
696 if (pme_lb->cutoff_scheme == ecutsVERLET && nbv->ngrp > 1)
698 init_interaction_const_tables(NULL, ic, bUsesSimpleTables,
702 if (cr->duty & DUTY_PME)
704 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
706 /* Generate a new PME data structure,
707 * copying part of the old pointers.
709 gmx_pme_reinit(&set->pmedata,
710 cr, pme_lb->setup[0].pmedata, ir,
713 *pmedata = set->pmedata;
717 /* Tell our PME-only node to switch grid */
718 gmx_pme_send_switchgrid(cr, set->grid, set->ewaldcoeff);
723 print_grid(NULL, debug, "", "switched to", set, -1);
726 if (pme_lb->stage == pme_lb->nstage)
728 print_grid(fp_err, fp_log, "", "optimal", set, -1);
734 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
739 static int pme_grid_points(const pme_setup_t *setup)
741 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
744 static real pme_loadbal_rlist(const pme_setup_t *setup)
746 /* With the group cut-off scheme we can have twin-range either
747 * for Coulomb or for VdW, so we need a check here.
748 * With the Verlet cut-off scheme rlist=rlistlong.
750 if (setup->rcut_coulomb > setup->rlist)
752 return setup->rlistlong;
760 static void print_pme_loadbal_setting(FILE *fplog,
762 const pme_setup_t *setup)
765 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
767 setup->rcut_coulomb, pme_loadbal_rlist(setup),
768 setup->grid[XX], setup->grid[YY], setup->grid[ZZ],
769 setup->spacing, 1/setup->ewaldcoeff);
772 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
775 gmx_bool bNonBondedOnGPU)
777 double pp_ratio, grid_ratio;
779 pp_ratio = pow(pme_loadbal_rlist(&pme_lb->setup[pme_lb->cur])/pme_loadbal_rlist(&pme_lb->setup[0]), 3.0);
780 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
781 (double)pme_grid_points(&pme_lb->setup[0]);
783 fprintf(fplog, "\n");
784 fprintf(fplog, " P P - P M E L O A D B A L A N C I N G\n");
785 fprintf(fplog, "\n");
786 /* Here we only warn when the optimal setting is the last one */
787 if (pme_lb->elimited != epmelblimNO &&
788 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
790 fprintf(fplog, " NOTE: The PP/PME load balancing was limited by the %s,\n",
791 pmelblim_str[pme_lb->elimited]);
792 fprintf(fplog, " you might not have reached a good load balance.\n");
793 if (pme_lb->elimited == epmelblimDD)
795 fprintf(fplog, " Try different mdrun -dd settings or lower the -dds value.\n");
797 fprintf(fplog, "\n");
799 fprintf(fplog, " PP/PME load balancing changed the cut-off and PME settings:\n");
800 fprintf(fplog, " particle-particle PME\n");
801 fprintf(fplog, " rcoulomb rlist grid spacing 1/beta\n");
802 print_pme_loadbal_setting(fplog, "initial", &pme_lb->setup[0]);
803 print_pme_loadbal_setting(fplog, "final", &pme_lb->setup[pme_lb->cur]);
804 fprintf(fplog, " cost-ratio %4.2f %4.2f\n",
805 pp_ratio, grid_ratio);
806 fprintf(fplog, " (note that these numbers concern only part of the total PP and PME load)\n");
808 if (pp_ratio > 1.5 && !bNonBondedOnGPU)
810 md_print_warn(cr, fplog,
811 "NOTE: PME load balancing increased the non-bonded workload by more than 50%%.\n"
812 " For better performance use (more) PME nodes (mdrun -npme),\n"
813 " or in case you are beyond the scaling limit, use less nodes in total.\n");
817 fprintf(fplog, "\n");
821 void pme_loadbal_done(pme_load_balancing_t pme_lb,
822 t_commrec *cr, FILE *fplog,
823 gmx_bool bNonBondedOnGPU)
825 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
827 print_pme_loadbal_settings(pme_lb, cr, fplog, bNonBondedOnGPU);
830 /* TODO: Here we should free all pointers in pme_lb,
831 * but as it contains pme data structures,
832 * we need to first make pme.c free all data.