2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2011, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "nbnxn_cuda_data_mgmt.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; /* the Ewald coefficient */
62 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
81 enum { epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimNR };
83 const char *pmelblim_str[epmelblimNR] =
84 { "no", "box size", "domain decompostion" };
86 struct pme_load_balancing {
87 int nstage; /* the current maximum number of stages */
89 real cut_spacing; /* the minimum cutoff / PME grid spacing ratio */
90 real rcut_vdw; /* Vdw cutoff (does not change) */
91 real rcut_coulomb_start; /* Initial electrostatics cutoff */
92 int nstcalclr_start; /* Initial electrostatics cutoff */
93 real rbuf_coulomb; /* the pairlist buffer size */
94 real rbuf_vdw; /* the pairlist buffer size */
95 matrix box_start; /* the initial simulation box */
96 int n; /* the count of setup as well as the allocation size */
97 pme_setup_t *setup; /* the PME+cutoff setups */
98 int cur; /* the current setup */
99 int fastest; /* fastest setup up till now */
100 int start; /* start of setup range to consider in stage>0 */
101 int end; /* end of setup range to consider in stage>0 */
102 int elimited; /* was the balancing limited, uses enum above */
103 int cutoff_scheme; /* Verlet or group cut-offs */
105 int stage; /* the current stage */
108 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
109 const t_inputrec *ir,matrix box,
110 const interaction_const_t *ic,
113 pme_load_balancing_t pme_lb;
119 /* Any number of stages >= 2 is supported */
122 pme_lb->cutoff_scheme = ir->cutoff_scheme;
124 if(pme_lb->cutoff_scheme == ecutsVERLET)
126 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
127 pme_lb->rbuf_vdw = pme_lb->rbuf_coulomb;
131 if(ic->rcoulomb > ic->rlist)
133 pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
137 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
139 if(ic->rvdw > ic->rlist)
141 pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
145 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
149 copy_mat(box,pme_lb->box_start);
150 if (ir->ePBC==epbcXY && ir->nwall==2)
152 svmul(ir->wall_ewald_zfac,pme_lb->box_start[ZZ],pme_lb->box_start[ZZ]);
156 snew(pme_lb->setup,pme_lb->n);
158 pme_lb->rcut_vdw = ic->rvdw;
159 pme_lb->rcut_coulomb_start = ir->rcoulomb;
160 pme_lb->nstcalclr_start = ir->nstcalclr;
163 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
164 pme_lb->setup[0].rlist = ic->rlist;
165 pme_lb->setup[0].rlistlong = ic->rlistlong;
166 pme_lb->setup[0].nstcalclr = ir->nstcalclr;
167 pme_lb->setup[0].grid[XX] = ir->nkx;
168 pme_lb->setup[0].grid[YY] = ir->nky;
169 pme_lb->setup[0].grid[ZZ] = ir->nkz;
170 pme_lb->setup[0].ewaldcoeff = ic->ewaldcoeff;
172 pme_lb->setup[0].pmedata = pmedata;
177 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
183 pme_lb->setup[0].spacing = spm;
185 if (ir->fourier_spacing > 0)
187 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
191 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
199 pme_lb->elimited = epmelblimNO;
204 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t pme_lb,
209 real tmpr_coulomb,tmpr_vdw;
212 /* Try to add a new setup with next larger cut-off to the list */
214 srenew(pme_lb->setup,pme_lb->n);
215 set = &pme_lb->setup[pme_lb->n-1];
222 clear_ivec(set->grid);
223 sp = calc_grid(NULL,pme_lb->box_start,
224 fac*pme_lb->setup[pme_lb->cur].spacing,
229 /* In parallel we can't have grids smaller than 2*pme_order,
230 * and we would anyhow not gain much speed at these grid sizes.
234 if (set->grid[d] <= 2*pme_order)
242 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing);
244 set->rcut_coulomb = pme_lb->cut_spacing*sp;
246 if(pme_lb->cutoff_scheme == ecutsVERLET)
248 set->rlist = set->rcut_coulomb + pme_lb->rbuf_coulomb;
249 /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
250 set->rlistlong = set->rlist;
254 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
255 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
256 set->rlist = min(tmpr_coulomb,tmpr_vdw);
257 set->rlistlong = max(tmpr_coulomb,tmpr_vdw);
259 /* Set the long-range update frequency */
260 if(set->rlist == set->rlistlong)
262 /* No long-range interactions if the short-/long-range cutoffs are identical */
265 else if(pme_lb->nstcalclr_start==0 || pme_lb->nstcalclr_start==1)
267 /* We were not doing long-range before, but now we are since rlist!=rlistlong */
272 /* We were already doing long-range interactions from the start */
273 if(pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
275 /* We were originally doing long-range VdW-only interactions.
276 * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
277 * but if the coulomb cutoff has become longer we should update the long-range
280 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
284 /* We were not doing any long-range interaction from the start,
285 * since it is not possible to do twin-range coulomb for the PME interaction.
293 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
294 set->grid_efficiency = 1;
297 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
299 /* The Ewald coefficient is inversly proportional to the cut-off */
301 pme_lb->setup[0].ewaldcoeff*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
308 fprintf(debug,"PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
309 set->grid[XX],set->grid[YY],set->grid[ZZ],set->rcut_coulomb);
314 static void print_grid(FILE *fp_err,FILE *fp_log,
317 const pme_setup_t *set,
320 char buf[STRLEN],buft[STRLEN];
324 sprintf(buft,": %.1f M-cycles",cycles*1e-6);
330 sprintf(buf,"%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
332 desc,set->grid[XX],set->grid[YY],set->grid[ZZ],set->rcut_coulomb,
336 fprintf(fp_err,"\r%s\n",buf);
340 fprintf(fp_log,"%s\n",buf);
344 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
346 /* In the initial stage only n is set; end is not set yet */
357 static void print_loadbal_limited(FILE *fp_err,FILE *fp_log,
358 gmx_large_int_t step,
359 pme_load_balancing_t pme_lb)
361 char buf[STRLEN],sbuf[22];
363 sprintf(buf,"step %4s: the %s limited the PME load balancing to a coulomb cut-off of %.3f",
364 gmx_step_str(step,sbuf),
365 pmelblim_str[pme_lb->elimited],
366 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
369 fprintf(fp_err,"\r%s\n",buf);
373 fprintf(fp_log,"%s\n",buf);
377 static void switch_to_stage1(pme_load_balancing_t pme_lb)
380 while (pme_lb->start+1 < pme_lb->n &&
381 (pme_lb->setup[pme_lb->start].count == 0 ||
382 pme_lb->setup[pme_lb->start].cycles >
383 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
387 while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
392 pme_lb->end = pme_lb->n;
393 if (pme_lb->setup[pme_lb->end-1].count > 0 &&
394 pme_lb->setup[pme_lb->end-1].cycles >
395 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
402 /* Next we want to choose setup pme_lb->start, but as we will increase
403 * pme_ln->cur by one right after returning, we subtract 1 here.
405 pme_lb->cur = pme_lb->start - 1;
408 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
415 interaction_const_t *ic,
416 nonbonded_verlet_t *nbv,
418 gmx_large_int_t step)
423 char buf[STRLEN],sbuf[22];
425 gmx_bool bUsesSimpleTables = TRUE;
427 if (pme_lb->stage == pme_lb->nstage)
434 gmx_sumd(1,&cycles,cr);
435 cycles /= cr->nnodes;
438 set = &pme_lb->setup[pme_lb->cur];
441 rtab = ir->rlistlong + ir->tabext;
443 if (set->count % 2 == 1)
445 /* Skip the first cycle, because the first step after a switch
446 * is much slower due to allocation and/or caching effects.
451 sprintf(buf, "step %4s: ", gmx_step_str(step,sbuf));
452 print_grid(fp_err,fp_log,buf,"timed with",set,cycles);
456 set->cycles = cycles;
460 if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
461 pme_lb->stage == pme_lb->nstage - 1)
463 /* The performance went up a lot (due to e.g. DD load balancing).
464 * Add a stage, keep the minima, but rescan all setups.
470 fprintf(debug,"The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
471 "Increased the number stages to %d"
472 " and ignoring the previous performance\n",
473 set->grid[XX],set->grid[YY],set->grid[ZZ],
474 cycles*1e-6,set->cycles*1e-6,PME_LB_ACCEL_TOL,
478 set->cycles = min(set->cycles,cycles);
481 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
483 pme_lb->fastest = pme_lb->cur;
485 if (DOMAINDECOMP(cr))
487 /* We found a new fastest setting, ensure that with subsequent
488 * shorter cut-off's the dynamic load balancing does not make
489 * the use of the current cut-off impossible. This solution is
490 * a trade-off, as the PME load balancing and DD domain size
491 * load balancing can interact in complex ways.
492 * With the Verlet kernels, DD load imbalance will usually be
493 * mainly due to bonded interaction imbalance, which will often
494 * quickly push the domain boundaries beyond the limit for the
495 * optimal, PME load balanced, cut-off. But it could be that
496 * better overal performance can be obtained with a slightly
497 * shorter cut-off and better DD load balancing.
499 change_dd_dlb_cutoff_limit(cr);
502 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
504 /* Check in stage 0 if we should stop scanning grids.
505 * Stop when the time is more than SLOW_FAC longer than the fastest.
507 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
508 cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
510 pme_lb->n = pme_lb->cur + 1;
511 /* Done with scanning, go to stage 1 */
512 switch_to_stage1(pme_lb);
515 if (pme_lb->stage == 0)
519 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
523 if (pme_lb->cur+1 < pme_lb->n)
525 /* We had already generated the next setup */
530 /* Find the next setup */
531 OK = pme_loadbal_increase_cutoff(pme_lb,ir->pme_order);
534 if (OK && ir->ePBC != epbcNONE)
536 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
537 <= max_cutoff2(ir->ePBC,state->box));
540 pme_lb->elimited = epmelblimBOX;
548 if (DOMAINDECOMP(cr))
550 OK = change_dd_cutoff(cr,state,ir,
551 pme_lb->setup[pme_lb->cur].rlistlong);
554 /* Failed: do not use this setup */
556 pme_lb->elimited = epmelblimDD;
562 /* We hit the upper limit for the cut-off,
563 * the setup should not go further than cur.
565 pme_lb->n = pme_lb->cur + 1;
566 print_loadbal_limited(fp_err,fp_log,step,pme_lb);
567 /* Switch to the next stage */
568 switch_to_stage1(pme_lb);
572 !(pme_lb->setup[pme_lb->cur].grid[XX]*
573 pme_lb->setup[pme_lb->cur].grid[YY]*
574 pme_lb->setup[pme_lb->cur].grid[ZZ] <
575 gridsize_start*PME_LB_GRID_SCALE_FAC
577 pme_lb->setup[pme_lb->cur].grid_efficiency <
578 pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
581 if (pme_lb->stage > 0 && pme_lb->end == 1)
584 pme_lb->stage = pme_lb->nstage;
586 else if (pme_lb->stage > 0 && pme_lb->end > 1)
588 /* If stage = nstage-1:
589 * scan over all setups, rerunning only those setups
590 * which are not much slower than the fastest
597 if (pme_lb->cur == pme_lb->end)
600 pme_lb->cur = pme_lb->start;
603 while (pme_lb->stage == pme_lb->nstage - 1 &&
604 pme_lb->setup[pme_lb->cur].count > 0 &&
605 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
607 if (pme_lb->stage == pme_lb->nstage)
609 /* We are done optimizing, use the fastest setup we found */
610 pme_lb->cur = pme_lb->fastest;
614 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
616 OK = change_dd_cutoff(cr,state,ir,pme_lb->setup[pme_lb->cur].rlistlong);
619 /* Failsafe solution */
620 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
626 pme_lb->end = pme_lb->cur;
627 pme_lb->cur = pme_lb->start;
628 pme_lb->elimited = epmelblimDD;
629 print_loadbal_limited(fp_err,fp_log,step,pme_lb);
633 /* Change the Coulomb cut-off and the PME grid */
635 set = &pme_lb->setup[pme_lb->cur];
637 ic->rcoulomb = set->rcut_coulomb;
638 ic->rlist = set->rlist;
639 ic->rlistlong = set->rlistlong;
640 ir->nstcalclr = set->nstcalclr;
641 ic->ewaldcoeff = set->ewaldcoeff;
643 bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
644 if (pme_lb->cutoff_scheme == ecutsVERLET &&
645 nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
647 nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv,ic);
651 init_interaction_const_tables(NULL,ic,bUsesSimpleTables,
655 if (pme_lb->cutoff_scheme == ecutsVERLET && nbv->ngrp > 1)
657 init_interaction_const_tables(NULL,ic,bUsesSimpleTables,
661 if (cr->duty & DUTY_PME)
663 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
665 /* Generate a new PME data structure,
666 * copying part of the old pointers.
668 gmx_pme_reinit(&set->pmedata,
669 cr,pme_lb->setup[0].pmedata,ir,
672 *pmedata = set->pmedata;
676 /* Tell our PME-only node to switch grid */
677 gmx_pme_send_switch(cr, set->grid, set->ewaldcoeff);
682 print_grid(NULL,debug,"","switched to",set,-1);
685 if (pme_lb->stage == pme_lb->nstage)
687 print_grid(fp_err,fp_log,"","optimal",set,-1);
693 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
698 static int pme_grid_points(const pme_setup_t *setup)
700 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
703 static void print_pme_loadbal_setting(FILE *fplog,
705 const pme_setup_t *setup)
708 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
710 setup->rcut_coulomb,setup->rlist,
711 setup->grid[XX],setup->grid[YY],setup->grid[ZZ],
712 setup->spacing,1/setup->ewaldcoeff);
715 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
718 double pp_ratio,grid_ratio;
720 pp_ratio = pow(pme_lb->setup[pme_lb->cur].rlist/pme_lb->setup[0].rlistlong,3.0);
721 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
722 (double)pme_grid_points(&pme_lb->setup[0]);
725 fprintf(fplog," P P - P M E L O A D B A L A N C I N G\n");
727 /* Here we only warn when the optimal setting is the last one */
728 if (pme_lb->elimited != epmelblimNO &&
729 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
731 fprintf(fplog," NOTE: The PP/PME load balancing was limited by the %s,\n",
732 pmelblim_str[pme_lb->elimited]);
733 fprintf(fplog," you might not have reached a good load balance.\n");
734 if (pme_lb->elimited == epmelblimDD)
736 fprintf(fplog," Try different mdrun -dd settings or lower the -dds value.\n");
740 fprintf(fplog," PP/PME load balancing changed the cut-off and PME settings:\n");
741 fprintf(fplog," particle-particle PME\n");
742 fprintf(fplog," rcoulomb rlist grid spacing 1/beta\n");
743 print_pme_loadbal_setting(fplog,"initial",&pme_lb->setup[0]);
744 print_pme_loadbal_setting(fplog,"final" ,&pme_lb->setup[pme_lb->cur]);
745 fprintf(fplog," cost-ratio %4.2f %4.2f\n",
746 pp_ratio,grid_ratio);
747 fprintf(fplog," (note that these numbers concern only part of the total PP and PME load)\n");
751 void pme_loadbal_done(pme_load_balancing_t pme_lb, FILE *fplog)
753 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
755 print_pme_loadbal_settings(pme_lb,fplog);
758 /* TODO: Here we should free all pointers in pme_lb,
759 * but as it contains pme data structures,
760 * we need to first make pme.c free all data.