1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2011, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
46 #include "nbnxn_cuda_data_mgmt.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
80 enum { epmelblimNO, epmelblimBOX, epmelblimDD, epmelblimNR };
82 const char *pmelblim_str[epmelblimNR] =
83 { "no", "box size", "domain decompostion" };
85 struct pme_load_balancing {
86 int nstage; /* the current maximum number of stages */
88 real cut_spacing; /* the minimum cutoff / PME grid spacing ratio */
89 real rcut_vdw; /* Vdw cutoff (does not change) */
90 real rcut_coulomb_start; /* Initial electrostatics cutoff */
91 int nstcalclr_start; /* Initial electrostatics cutoff */
92 real rbuf_coulomb; /* the pairlist buffer size */
93 real rbuf_vdw; /* the pairlist buffer size */
94 matrix box_start; /* the initial simulation box */
95 int n; /* the count of setup as well as the allocation size */
96 pme_setup_t *setup; /* the PME+cutoff setups */
97 int cur; /* the current setup */
98 int fastest; /* fastest setup up till now */
99 int start; /* start of setup range to consider in stage>0 */
100 int end; /* end of setup range to consider in stage>0 */
101 int elimited; /* was the balancing limited, uses enum above */
102 int cutoff_scheme; /* Verlet or group cut-offs */
104 int stage; /* the current stage */
107 void pme_loadbal_init(pme_load_balancing_t *pme_lb_p,
108 const t_inputrec *ir,matrix box,
109 const interaction_const_t *ic,
112 pme_load_balancing_t pme_lb;
118 /* Any number of stages >= 2 is supported */
121 pme_lb->cutoff_scheme = ir->cutoff_scheme;
123 if(pme_lb->cutoff_scheme == ecutsVERLET)
125 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
126 pme_lb->rbuf_vdw = pme_lb->rbuf_coulomb;
130 if(ic->rcoulomb > ic->rlist)
132 pme_lb->rbuf_coulomb = ic->rlistlong - ic->rcoulomb;
136 pme_lb->rbuf_coulomb = ic->rlist - ic->rcoulomb;
138 if(ic->rvdw > ic->rlist)
140 pme_lb->rbuf_vdw = ic->rlistlong - ic->rvdw;
144 pme_lb->rbuf_vdw = ic->rlist - ic->rvdw;
148 copy_mat(box,pme_lb->box_start);
149 if (ir->ePBC==epbcXY && ir->nwall==2)
151 svmul(ir->wall_ewald_zfac,pme_lb->box_start[ZZ],pme_lb->box_start[ZZ]);
155 snew(pme_lb->setup,pme_lb->n);
157 pme_lb->rcut_vdw = ic->rvdw;
158 pme_lb->rcut_coulomb_start = ir->rcoulomb;
159 pme_lb->nstcalclr_start = ir->nstcalclr;
162 pme_lb->setup[0].rcut_coulomb = ic->rcoulomb;
163 pme_lb->setup[0].rlist = ic->rlist;
164 pme_lb->setup[0].rlistlong = ic->rlistlong;
165 pme_lb->setup[0].nstcalclr = ir->nstcalclr;
166 pme_lb->setup[0].grid[XX] = ir->nkx;
167 pme_lb->setup[0].grid[YY] = ir->nky;
168 pme_lb->setup[0].grid[ZZ] = ir->nkz;
169 pme_lb->setup[0].ewaldcoeff = ic->ewaldcoeff;
171 pme_lb->setup[0].pmedata = pmedata;
176 sp = norm(pme_lb->box_start[d])/pme_lb->setup[0].grid[d];
182 pme_lb->setup[0].spacing = spm;
184 if (ir->fourier_spacing > 0)
186 pme_lb->cut_spacing = ir->rcoulomb/ir->fourier_spacing;
190 pme_lb->cut_spacing = ir->rcoulomb/pme_lb->setup[0].spacing;
198 pme_lb->elimited = epmelblimNO;
203 static gmx_bool pme_loadbal_increase_cutoff(pme_load_balancing_t pme_lb,
208 real tmpr_coulomb,tmpr_vdw;
211 /* Try to add a new setup with next larger cut-off to the list */
213 srenew(pme_lb->setup,pme_lb->n);
214 set = &pme_lb->setup[pme_lb->n-1];
221 clear_ivec(set->grid);
222 sp = calc_grid(NULL,pme_lb->box_start,
223 fac*pme_lb->setup[pme_lb->cur].spacing,
228 /* In parallel we can't have grids smaller than 2*pme_order,
229 * and we would anyhow not gain much speed at these grid sizes.
233 if (set->grid[d] <= 2*pme_order)
241 while (sp <= 1.001*pme_lb->setup[pme_lb->cur].spacing);
243 set->rcut_coulomb = pme_lb->cut_spacing*sp;
245 if(pme_lb->cutoff_scheme == ecutsVERLET)
247 set->rlist = set->rcut_coulomb + pme_lb->rbuf_coulomb;
248 /* We dont use LR lists with Verlet, but this avoids if-statements in further checks */
249 set->rlistlong = set->rlist;
253 tmpr_coulomb = set->rcut_coulomb + pme_lb->rbuf_coulomb;
254 tmpr_vdw = pme_lb->rcut_vdw + pme_lb->rbuf_vdw;
255 set->rlist = min(tmpr_coulomb,tmpr_vdw);
256 set->rlistlong = max(tmpr_coulomb,tmpr_vdw);
258 /* Set the long-range update frequency */
259 if(set->rlist == set->rlistlong)
261 /* No long-range interactions if the short-/long-range cutoffs are identical */
264 else if(pme_lb->nstcalclr_start==0 || pme_lb->nstcalclr_start==1)
266 /* We were not doing long-range before, but now we are since rlist!=rlistlong */
271 /* We were already doing long-range interactions from the start */
272 if(pme_lb->rcut_vdw > pme_lb->rcut_coulomb_start)
274 /* We were originally doing long-range VdW-only interactions.
275 * If rvdw is still longer than rcoulomb we keep the original nstcalclr,
276 * but if the coulomb cutoff has become longer we should update the long-range
279 set->nstcalclr = (tmpr_vdw > tmpr_coulomb) ? pme_lb->nstcalclr_start : 1;
283 /* We were not doing any long-range interaction from the start,
284 * since it is not possible to do twin-range coulomb for the PME interaction.
292 /* The grid efficiency is the size wrt a grid with uniform x/y/z spacing */
293 set->grid_efficiency = 1;
296 set->grid_efficiency *= (set->grid[d]*sp)/norm(pme_lb->box_start[d]);
298 /* The Ewald coefficient is inversly proportional to the cut-off */
300 pme_lb->setup[0].ewaldcoeff*pme_lb->setup[0].rcut_coulomb/set->rcut_coulomb;
307 fprintf(debug,"PME loadbal: grid %d %d %d, coulomb cutoff %f\n",
308 set->grid[XX],set->grid[YY],set->grid[ZZ],set->rcut_coulomb);
313 static void print_grid(FILE *fp_err,FILE *fp_log,
316 const pme_setup_t *set,
319 char buf[STRLEN],buft[STRLEN];
323 sprintf(buft,": %.1f M-cycles",cycles*1e-6);
329 sprintf(buf,"%-11s%10s pme grid %d %d %d, coulomb cutoff %.3f%s",
331 desc,set->grid[XX],set->grid[YY],set->grid[ZZ],set->rcut_coulomb,
335 fprintf(fp_err,"\r%s\n",buf);
339 fprintf(fp_log,"%s\n",buf);
343 static int pme_loadbal_end(pme_load_balancing_t pme_lb)
345 /* In the initial stage only n is set; end is not set yet */
356 static void print_loadbal_limited(FILE *fp_err,FILE *fp_log,
357 gmx_large_int_t step,
358 pme_load_balancing_t pme_lb)
360 char buf[STRLEN],sbuf[22];
362 sprintf(buf,"step %4s: the %s limited the PME load balancing to a coulomb cut-off of %.3f",
363 gmx_step_str(step,sbuf),
364 pmelblim_str[pme_lb->elimited],
365 pme_lb->setup[pme_loadbal_end(pme_lb)-1].rcut_coulomb);
368 fprintf(fp_err,"\r%s\n",buf);
372 fprintf(fp_log,"%s\n",buf);
376 static void switch_to_stage1(pme_load_balancing_t pme_lb)
379 while (pme_lb->start+1 < pme_lb->n &&
380 (pme_lb->setup[pme_lb->start].count == 0 ||
381 pme_lb->setup[pme_lb->start].cycles >
382 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC))
386 while (pme_lb->start > 0 && pme_lb->setup[pme_lb->start-1].cycles == 0)
391 pme_lb->end = pme_lb->n;
392 if (pme_lb->setup[pme_lb->end-1].count > 0 &&
393 pme_lb->setup[pme_lb->end-1].cycles >
394 pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
401 /* Next we want to choose setup pme_lb->start, but as we will increase
402 * pme_ln->cur by one right after returning, we subtract 1 here.
404 pme_lb->cur = pme_lb->start - 1;
407 gmx_bool pme_load_balance(pme_load_balancing_t pme_lb,
414 interaction_const_t *ic,
415 nonbonded_verlet_t *nbv,
417 gmx_large_int_t step)
422 char buf[STRLEN],sbuf[22];
424 gmx_bool bUsesSimpleTables = TRUE;
426 if (pme_lb->stage == pme_lb->nstage)
433 gmx_sumd(1,&cycles,cr);
434 cycles /= cr->nnodes;
437 set = &pme_lb->setup[pme_lb->cur];
440 rtab = ir->rlistlong + ir->tabext;
442 if (set->count % 2 == 1)
444 /* Skip the first cycle, because the first step after a switch
445 * is much slower due to allocation and/or caching effects.
450 sprintf(buf, "step %4s: ", gmx_step_str(step,sbuf));
451 print_grid(fp_err,fp_log,buf,"timed with",set,cycles);
455 set->cycles = cycles;
459 if (cycles*PME_LB_ACCEL_TOL < set->cycles &&
460 pme_lb->stage == pme_lb->nstage - 1)
462 /* The performance went up a lot (due to e.g. DD load balancing).
463 * Add a stage, keep the minima, but rescan all setups.
469 fprintf(debug,"The performance for grid %d %d %d went from %.3f to %.1f M-cycles, this is more than %f\n"
470 "Increased the number stages to %d"
471 " and ignoring the previous performance\n",
472 set->grid[XX],set->grid[YY],set->grid[ZZ],
473 cycles*1e-6,set->cycles*1e-6,PME_LB_ACCEL_TOL,
477 set->cycles = min(set->cycles,cycles);
480 if (set->cycles < pme_lb->setup[pme_lb->fastest].cycles)
482 pme_lb->fastest = pme_lb->cur;
484 cycles_fast = pme_lb->setup[pme_lb->fastest].cycles;
486 /* Check in stage 0 if we should stop scanning grids.
487 * Stop when the time is more than SLOW_FAC longer than the fastest.
489 if (pme_lb->stage == 0 && pme_lb->cur > 0 &&
490 cycles > pme_lb->setup[pme_lb->fastest].cycles*PME_LB_SLOW_FAC)
492 pme_lb->n = pme_lb->cur + 1;
493 /* Done with scanning, go to stage 1 */
494 switch_to_stage1(pme_lb);
497 if (pme_lb->stage == 0)
501 gridsize_start = set->grid[XX]*set->grid[YY]*set->grid[ZZ];
505 if (pme_lb->cur+1 < pme_lb->n)
507 /* We had already generated the next setup */
512 /* Find the next setup */
513 OK = pme_loadbal_increase_cutoff(pme_lb,ir->pme_order);
516 if (OK && ir->ePBC != epbcNONE)
518 OK = (sqr(pme_lb->setup[pme_lb->cur+1].rlistlong)
519 <= max_cutoff2(ir->ePBC,state->box));
522 pme_lb->elimited = epmelblimBOX;
530 if (DOMAINDECOMP(cr))
532 OK = change_dd_cutoff(cr,state,ir,
533 pme_lb->setup[pme_lb->cur].rlistlong);
536 /* Failed: do not use this setup */
538 pme_lb->elimited = epmelblimDD;
544 /* We hit the upper limit for the cut-off,
545 * the setup should not go further than cur.
547 pme_lb->n = pme_lb->cur + 1;
548 print_loadbal_limited(fp_err,fp_log,step,pme_lb);
549 /* Switch to the next stage */
550 switch_to_stage1(pme_lb);
554 !(pme_lb->setup[pme_lb->cur].grid[XX]*
555 pme_lb->setup[pme_lb->cur].grid[YY]*
556 pme_lb->setup[pme_lb->cur].grid[ZZ] <
557 gridsize_start*PME_LB_GRID_SCALE_FAC
559 pme_lb->setup[pme_lb->cur].grid_efficiency <
560 pme_lb->setup[pme_lb->cur-1].grid_efficiency*PME_LB_GRID_EFFICIENCY_REL_FAC));
563 if (pme_lb->stage > 0 && pme_lb->end == 1)
566 pme_lb->stage = pme_lb->nstage;
568 else if (pme_lb->stage > 0 && pme_lb->end > 1)
570 /* If stage = nstage-1:
571 * scan over all setups, rerunning only those setups
572 * which are not much slower than the fastest
579 if (pme_lb->cur == pme_lb->end)
582 pme_lb->cur = pme_lb->start;
585 while (pme_lb->stage == pme_lb->nstage - 1 &&
586 pme_lb->setup[pme_lb->cur].count > 0 &&
587 pme_lb->setup[pme_lb->cur].cycles > cycles_fast*PME_LB_SLOW_FAC);
589 if (pme_lb->stage == pme_lb->nstage)
591 /* We are done optimizing, use the fastest setup we found */
592 pme_lb->cur = pme_lb->fastest;
596 if (DOMAINDECOMP(cr) && pme_lb->stage > 0)
598 OK = change_dd_cutoff(cr,state,ir,pme_lb->setup[pme_lb->cur].rlistlong);
601 /* Failsafe solution */
602 if (pme_lb->cur > 1 && pme_lb->stage == pme_lb->nstage)
608 pme_lb->end = pme_lb->cur;
609 pme_lb->cur = pme_lb->start;
610 pme_lb->elimited = epmelblimDD;
611 print_loadbal_limited(fp_err,fp_log,step,pme_lb);
615 /* Change the Coulomb cut-off and the PME grid */
617 set = &pme_lb->setup[pme_lb->cur];
619 ic->rcoulomb = set->rcut_coulomb;
620 ic->rlist = set->rlist;
621 ic->rlistlong = set->rlistlong;
622 ir->nstcalclr = set->nstcalclr;
623 ic->ewaldcoeff = set->ewaldcoeff;
625 bUsesSimpleTables = uses_simple_tables(ir->cutoff_scheme, nbv, 0);
626 if (pme_lb->cutoff_scheme == ecutsVERLET && nbv->grp[0].kernel_type == nbk8x8x8_CUDA)
628 nbnxn_cuda_pme_loadbal_update_param(nbv->cu_nbv,ic);
632 init_interaction_const_tables(NULL,ic,bUsesSimpleTables,
636 if (pme_lb->cutoff_scheme == ecutsVERLET && nbv->ngrp > 1)
638 init_interaction_const_tables(NULL,ic,bUsesSimpleTables,
642 if (cr->duty & DUTY_PME)
644 if (pme_lb->setup[pme_lb->cur].pmedata == NULL)
646 /* Generate a new PME data structure,
647 * copying part of the old pointers.
649 gmx_pme_reinit(&set->pmedata,
650 cr,pme_lb->setup[0].pmedata,ir,
653 *pmedata = set->pmedata;
657 /* Tell our PME-only node to switch grid */
658 gmx_pme_send_switch(cr, set->grid, set->ewaldcoeff);
663 print_grid(NULL,debug,"","switched to",set,-1);
666 if (pme_lb->stage == pme_lb->nstage)
668 print_grid(fp_err,fp_log,"","optimal",set,-1);
674 void restart_pme_loadbal(pme_load_balancing_t pme_lb, int n)
679 static int pme_grid_points(const pme_setup_t *setup)
681 return setup->grid[XX]*setup->grid[YY]*setup->grid[ZZ];
684 static void print_pme_loadbal_setting(FILE *fplog,
686 const pme_setup_t *setup)
689 " %-7s %6.3f nm %6.3f nm %3d %3d %3d %5.3f nm %5.3f nm\n",
691 setup->rcut_coulomb,setup->rlist,
692 setup->grid[XX],setup->grid[YY],setup->grid[ZZ],
693 setup->spacing,1/setup->ewaldcoeff);
696 static void print_pme_loadbal_settings(pme_load_balancing_t pme_lb,
699 double pp_ratio,grid_ratio;
701 pp_ratio = pow(pme_lb->setup[pme_lb->cur].rlist/pme_lb->setup[0].rlistlong,3.0);
702 grid_ratio = pme_grid_points(&pme_lb->setup[pme_lb->cur])/
703 (double)pme_grid_points(&pme_lb->setup[0]);
706 fprintf(fplog," P P - P M E L O A D B A L A N C I N G\n");
708 /* Here we only warn when the optimal setting is the last one */
709 if (pme_lb->elimited != epmelblimNO &&
710 pme_lb->cur == pme_loadbal_end(pme_lb)-1)
712 fprintf(fplog," NOTE: The PP/PME load balancing was limited by the %s,\n",
713 pmelblim_str[pme_lb->elimited]);
714 fprintf(fplog," you might not have reached a good load balance.\n");
715 if (pme_lb->elimited == epmelblimDD)
717 fprintf(fplog," Try different mdrun -dd settings or lower the -dds value.\n");
721 fprintf(fplog," PP/PME load balancing changed the cut-off and PME settings:\n");
722 fprintf(fplog," particle-particle PME\n");
723 fprintf(fplog," rcoulomb rlist grid spacing 1/beta\n");
724 print_pme_loadbal_setting(fplog,"initial",&pme_lb->setup[0]);
725 print_pme_loadbal_setting(fplog,"final" ,&pme_lb->setup[pme_lb->cur]);
726 fprintf(fplog," cost-ratio %4.2f %4.2f\n",
727 pp_ratio,grid_ratio);
728 fprintf(fplog," (note that these numbers concern only part of the total PP and PME load)\n");
732 void pme_loadbal_done(pme_load_balancing_t pme_lb, FILE *fplog)
734 if (fplog != NULL && (pme_lb->cur > 0 || pme_lb->elimited != epmelblimNO))
736 print_pme_loadbal_settings(pme_lb,fplog);
739 /* TODO: Here we should free all pointers in pme_lb,
740 * but as it contains pme data structures,
741 * we need to first make pme.c free all data.