2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, 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.
42 #ifdef HAVE_SYS_TIME_H
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/tpxio.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/legacyheaders/calcgrid.h"
50 #include "gromacs/legacyheaders/checkpoint.h"
51 #include "gromacs/legacyheaders/inputrec.h"
52 #include "gromacs/legacyheaders/macros.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/perf_est.h"
55 #include "gromacs/legacyheaders/readinp.h"
56 #include "gromacs/legacyheaders/typedefs.h"
57 #include "gromacs/legacyheaders/types/commrec.h"
58 #include "gromacs/math/utilities.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/timing/walltime_accounting.h"
61 #include "gromacs/utility/baseversion.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/smalloc.h"
66 /* Enum for situations that can occur during log file parsing, the
67 * corresponding string entries can be found in do_the_tests() in
68 * const char* ParseLog[] */
69 /* TODO clean up CamelCasing of these enum names */
75 eParselogResetProblem,
79 eParselogLargePrimeFactor,
80 eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs,
89 int nPMEnodes; /* number of PME-only nodes used in this test */
90 int nx, ny, nz; /* DD grid */
91 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
92 double *Gcycles; /* This can contain more than one value if doing multiple tests */
96 float *PME_f_load; /* PME mesh/force load average*/
97 float PME_f_load_Av; /* Average average ;) ... */
98 char *mdrun_cmd_line; /* Mdrun command line used for this test */
104 int nr_inputfiles; /* The number of tpr and mdp input files */
105 gmx_int64_t orig_sim_steps; /* Number of steps to be done in the real simulation */
106 gmx_int64_t orig_init_step; /* Init step for the real simulation */
107 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
108 real *rvdw; /* The vdW radii */
109 real *rlist; /* Neighbourlist cutoff radius */
111 int *nkx, *nky, *nkz;
112 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
116 static void sep_line(FILE *fp)
118 fprintf(fp, "\n------------------------------------------------------------\n");
122 /* Wrapper for system calls */
123 static int gmx_system_call(char *command)
125 return ( system(command) );
129 /* Check if string starts with substring */
130 static gmx_bool str_starts(const char *string, const char *substring)
132 return ( strncmp(string, substring, strlen(substring)) == 0);
136 static void cleandata(t_perf *perfdata, int test_nr)
138 perfdata->Gcycles[test_nr] = 0.0;
139 perfdata->ns_per_day[test_nr] = 0.0;
140 perfdata->PME_f_load[test_nr] = 0.0;
146 static void remove_if_exists(const char *fn)
150 fprintf(stdout, "Deleting %s\n", fn);
156 static void finalize(const char *fn_out)
162 fp = fopen(fn_out, "r");
163 fprintf(stdout, "\n\n");
165 while (fgets(buf, STRLEN-1, fp) != NULL)
167 fprintf(stdout, "%s", buf);
170 fprintf(stdout, "\n\n");
175 eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr
178 static int parse_logfile(const char *logfile, const char *errfile,
179 t_perf *perfdata, int test_nr, int presteps, gmx_int64_t cpt_steps,
183 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
184 const char matchstrdd[] = "Domain decomposition grid";
185 const char matchstrcr[] = "resetting all time and cycle counters";
186 const char matchstrbal[] = "Average PME mesh/force load:";
187 const char matchstring[] = "R E A L C Y C L E A N D T I M E A C C O U N T I N G";
188 const char errSIG[] = "signal, stopping at the next";
190 float dum1, dum2, dum3, dum4;
193 gmx_int64_t resetsteps = -1;
194 gmx_bool bFoundResetStr = FALSE;
195 gmx_bool bResetChecked = FALSE;
198 if (!gmx_fexist(logfile))
200 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
201 cleandata(perfdata, test_nr);
202 return eParselogNotFound;
205 fp = fopen(logfile, "r");
206 perfdata->PME_f_load[test_nr] = -1.0;
207 perfdata->guessPME = -1;
209 iFound = eFoundNothing;
212 iFound = eFoundDDStr; /* Skip some case statements */
215 while (fgets(line, STRLEN, fp) != NULL)
217 /* Remove leading spaces */
220 /* Check for TERM and INT signals from user: */
221 if (strstr(line, errSIG) != NULL)
224 cleandata(perfdata, test_nr);
225 return eParselogTerm;
228 /* Check whether cycle resetting worked */
229 if (presteps > 0 && !bFoundResetStr)
231 if (strstr(line, matchstrcr) != NULL)
233 sprintf(dumstring, "step %s", "%"GMX_SCNd64);
234 sscanf(line, dumstring, &resetsteps);
235 bFoundResetStr = TRUE;
236 if (resetsteps == presteps+cpt_steps)
238 bResetChecked = TRUE;
242 sprintf(dumstring, "%"GMX_PRId64, resetsteps);
243 sprintf(dumstring2, "%"GMX_PRId64, presteps+cpt_steps);
244 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
245 " though they were supposed to be reset at step %s!\n",
246 dumstring, dumstring2);
251 /* Look for strings that appear in a certain order in the log file: */
255 /* Look for domain decomp grid and separate PME nodes: */
256 if (str_starts(line, matchstrdd))
258 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME ranks %d",
259 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
260 if (perfdata->nPMEnodes == -1)
262 perfdata->guessPME = npme;
264 else if (perfdata->nPMEnodes != npme)
266 gmx_fatal(FARGS, "PME ranks from command line and output file are not identical");
268 iFound = eFoundDDStr;
270 /* Catch a few errors that might have occured: */
271 else if (str_starts(line, "There is no domain decomposition for"))
274 return eParselogNoDDGrid;
276 else if (str_starts(line, "The number of ranks you selected"))
279 return eParselogLargePrimeFactor;
281 else if (str_starts(line, "reading tpx file"))
284 return eParselogTPXVersion;
286 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
289 return eParselogNotParallel;
293 /* Even after the "Domain decomposition grid" string was found,
294 * it could be that mdrun had to quit due to some error. */
295 if (str_starts(line, "Incorrect launch configuration: mismatching number of"))
298 return eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs;
300 else if (str_starts(line, "Some of the requested GPUs do not exist"))
303 return eParselogGpuProblem;
305 /* Look for PME mesh/force balance (not necessarily present, though) */
306 else if (str_starts(line, matchstrbal))
308 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
310 /* Look for matchstring */
311 else if (str_starts(line, matchstring))
313 iFound = eFoundAccountingStr;
316 case eFoundAccountingStr:
317 /* Already found matchstring - look for cycle data */
318 if (str_starts(line, "Total "))
320 sscanf(line, "Total %*f %lf", &(perfdata->Gcycles[test_nr]));
321 iFound = eFoundCycleStr;
325 /* Already found cycle data - look for remaining performance info and return */
326 if (str_starts(line, "Performance:"))
328 ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
329 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
330 perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
332 if (bResetChecked || presteps == 0)
338 return eParselogResetProblem;
345 /* Close the log file */
348 /* Check why there is no performance data in the log file.
349 * Did a fatal errors occur? */
350 if (gmx_fexist(errfile))
352 fp = fopen(errfile, "r");
353 while (fgets(line, STRLEN, fp) != NULL)
355 if (str_starts(line, "Fatal error:") )
357 if (fgets(line, STRLEN, fp) != NULL)
359 fprintf(stderr, "\nWARNING: An error occured during this benchmark:\n"
363 cleandata(perfdata, test_nr);
364 return eParselogFatal;
371 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
374 /* Giving up ... we could not find out why there is no performance data in
376 fprintf(stdout, "No performance data in log file.\n");
377 cleandata(perfdata, test_nr);
379 return eParselogNoPerfData;
383 static gmx_bool analyze_data(
392 int *index_tpr, /* OUT: Nr of mdp file with best settings */
393 int *npme_optimal) /* OUT: Optimal number of PME nodes */
396 int line = 0, line_win = -1;
397 int k_win = -1, i_win = -1, winPME;
398 double s = 0.0; /* standard deviation */
401 char str_PME_f_load[13];
402 gmx_bool bCanUseOrigTPR;
403 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
409 fprintf(fp, "Summary of successful runs:\n");
410 fprintf(fp, "Line tpr PME ranks Gcycles Av. Std.dev. ns/day PME/f");
413 fprintf(fp, " DD grid");
419 for (k = 0; k < ntprs; k++)
421 for (i = 0; i < ntests; i++)
423 /* Select the right dataset: */
424 pd = &(perfdata[k][i]);
426 pd->Gcycles_Av = 0.0;
427 pd->PME_f_load_Av = 0.0;
428 pd->ns_per_day_Av = 0.0;
430 if (pd->nPMEnodes == -1)
432 sprintf(strbuf, "(%3d)", pd->guessPME);
436 sprintf(strbuf, " ");
439 /* Get the average run time of a setting */
440 for (j = 0; j < nrepeats; j++)
442 pd->Gcycles_Av += pd->Gcycles[j];
443 pd->PME_f_load_Av += pd->PME_f_load[j];
445 pd->Gcycles_Av /= nrepeats;
446 pd->PME_f_load_Av /= nrepeats;
448 for (j = 0; j < nrepeats; j++)
450 if (pd->ns_per_day[j] > 0.0)
452 pd->ns_per_day_Av += pd->ns_per_day[j];
456 /* Somehow the performance number was not aquired for this run,
457 * therefor set the average to some negative value: */
458 pd->ns_per_day_Av = -1.0f*nrepeats;
462 pd->ns_per_day_Av /= nrepeats;
465 if (pd->PME_f_load_Av > 0.0)
467 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
471 sprintf(str_PME_f_load, "%s", " - ");
475 /* We assume we had a successful run if both averages are positive */
476 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
478 /* Output statistics if repeats were done */
481 /* Calculate the standard deviation */
483 for (j = 0; j < nrepeats; j++)
485 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
490 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
491 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
492 pd->ns_per_day_Av, str_PME_f_load);
495 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
499 /* Store the index of the best run found so far in 'winner': */
500 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
513 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
518 winPME = perfdata[k_win][i_win].nPMEnodes;
522 /* We stuck to a fixed number of PME-only nodes */
523 sprintf(strbuf, "settings No. %d", k_win);
527 /* We have optimized the number of PME-only nodes */
530 sprintf(strbuf, "%s", "the automatic number of PME ranks");
534 sprintf(strbuf, "%d PME ranks", winPME);
537 fprintf(fp, "Best performance was achieved with %s", strbuf);
538 if ((nrepeats > 1) && (ntests > 1))
540 fprintf(fp, " (see line %d)", line_win);
544 /* Only mention settings if they were modified: */
545 bRefinedCoul = !gmx_within_tol(info->rcoulomb[k_win], info->rcoulomb[0], GMX_REAL_EPS);
546 bRefinedVdW = !gmx_within_tol(info->rvdw[k_win], info->rvdw[0], GMX_REAL_EPS);
547 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
548 info->nky[k_win] == info->nky[0] &&
549 info->nkz[k_win] == info->nkz[0]);
551 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
553 fprintf(fp, "Optimized PME settings:\n");
554 bCanUseOrigTPR = FALSE;
558 bCanUseOrigTPR = TRUE;
563 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
568 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
573 fprintf(fp, " New Fourier grid xyz: %d %d %d (was %d %d %d)\n", info->nkx[k_win], info->nky[k_win], info->nkz[k_win],
574 info->nkx[0], info->nky[0], info->nkz[0]);
577 if (bCanUseOrigTPR && ntprs > 1)
579 fprintf(fp, "and original PME settings.\n");
584 /* Return the index of the mdp file that showed the highest performance
585 * and the optimal number of PME nodes */
587 *npme_optimal = winPME;
589 return bCanUseOrigTPR;
593 /* Get the commands we need to set up the runs from environment variables */
594 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
598 const char def_mpirun[] = "mpirun";
599 const char def_mdrun[] = "mdrun";
601 const char empty_mpirun[] = "";
603 /* Get the commands we need to set up the runs from environment variables */
606 if ( (cp = getenv("MPIRUN")) != NULL)
608 *cmd_mpirun = gmx_strdup(cp);
612 *cmd_mpirun = gmx_strdup(def_mpirun);
617 *cmd_mpirun = gmx_strdup(empty_mpirun);
620 if ( (cp = getenv("MDRUN" )) != NULL)
622 *cmd_mdrun = gmx_strdup(cp);
626 *cmd_mdrun = gmx_strdup(def_mdrun);
630 /* Check that the commands will run mdrun (perhaps via mpirun) by
631 * running a very quick test simulation. Requires MPI environment or
632 * GPU support to be available if applicable. */
633 /* TODO implement feature to parse the log file to get the list of
634 compatible GPUs from mdrun, if the user of gmx tune-pme has not
636 static void check_mdrun_works(gmx_bool bThreads,
637 const char *cmd_mpirun,
639 const char *cmd_mdrun,
640 gmx_bool bNeedGpuSupport)
642 char *command = NULL;
646 const char filename[] = "benchtest.log";
648 /* This string should always be identical to the one in copyrite.c,
649 * gmx_print_version_info() in the defined(GMX_MPI) section */
650 const char match_mpi[] = "MPI library: MPI";
651 const char match_mdrun[] = "Executable: ";
652 const char match_gpu[] = "GPU support: enabled";
653 gmx_bool bMdrun = FALSE;
654 gmx_bool bMPI = FALSE;
655 gmx_bool bHaveGpuSupport = FALSE;
657 /* Run a small test to see whether mpirun + mdrun work */
658 fprintf(stdout, "Making sure that mdrun can be executed. ");
661 snew(command, strlen(cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
662 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
666 snew(command, strlen(cmd_mpirun) + strlen(cmd_np) + strlen(cmd_mdrun) + strlen(filename) + 50);
667 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
669 fprintf(stdout, "Trying '%s' ... ", command);
670 make_backup(filename);
671 gmx_system_call(command);
673 /* Check if we find the characteristic string in the output: */
674 if (!gmx_fexist(filename))
676 gmx_fatal(FARGS, "Output from test run could not be found.");
679 fp = fopen(filename, "r");
680 /* We need to scan the whole output file, since sometimes the queuing system
681 * also writes stuff to stdout/err */
684 cp = fgets(line, STRLEN, fp);
687 if (str_starts(line, match_mdrun) )
691 if (str_starts(line, match_mpi) )
695 if (str_starts(line, match_gpu) )
697 bHaveGpuSupport = TRUE;
707 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
709 "seems to have been compiled with MPI instead.",
717 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
719 "seems to have been compiled without MPI support.",
726 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
730 if (bNeedGpuSupport && !bHaveGpuSupport)
732 gmx_fatal(FARGS, "The mdrun executable did not have the expected GPU support.");
735 fprintf(stdout, "passed.\n");
742 /*! \brief Helper struct so we can parse the string with eligible GPU
743 IDs outside do_the_tests. */
744 typedef struct eligible_gpu_ids
746 int n; /**< Length of ids */
747 int *ids; /**< Array of length n. NULL if no GPUs in use */
748 } t_eligible_gpu_ids;
750 /* Handles the no-GPU case by emitting an empty string. */
751 static char *make_gpu_id_command_line(int numRanks, int numPmeRanks, const t_eligible_gpu_ids *gpu_ids)
753 char *command_line, *flag = "-gpu_id ", *ptr;
756 /* Reserve enough room for the option name, enough single-digit
757 GPU ids (since that is currently all that is possible to use
758 with mdrun), and a terminating NULL. */
759 flag_length = strlen(flag);
760 snew(command_line, flag_length + numRanks + 1);
763 /* If the user has given no eligible GPU IDs, or we're trying the
764 * default behaviour, then there is nothing for g_tune_pme to give
765 * to mdrun -gpu_id */
766 if (gpu_ids->n > 0 && numPmeRanks > -1)
768 int numPpRanks, max_num_ranks_for_each_GPU;
771 /* Write the option flag */
775 numPpRanks = numRanks - numPmeRanks;
776 max_num_ranks_for_each_GPU = numPpRanks / gpu_ids->n;
777 if (max_num_ranks_for_each_GPU * gpu_ids->n != numPpRanks)
779 /* Some GPUs will receive more work than others, which
780 * we choose to be those with the lowest indices */
781 max_num_ranks_for_each_GPU++;
784 /* Loop over all eligible GPU ids */
785 for (gpu_id = 0, rank = 0; gpu_id < gpu_ids->n; gpu_id++)
787 int rank_for_this_GPU;
788 /* Loop over all PP ranks for GPU with ID gpu_id, building the
789 assignment string. */
790 for (rank_for_this_GPU = 0;
791 rank_for_this_GPU < max_num_ranks_for_each_GPU && rank < numPpRanks;
792 rank++, rank_for_this_GPU++)
794 *ptr = '0' + gpu_ids->ids[gpu_id];
804 static void launch_simulation(
805 gmx_bool bLaunch, /* Should the simulation be launched? */
806 FILE *fp, /* General log file */
807 gmx_bool bThreads, /* whether to use threads */
808 char *cmd_mpirun, /* Command for mpirun */
809 char *cmd_np, /* Switch for -np or -ntmpi or empty */
810 char *cmd_mdrun, /* Command for mdrun */
811 char *args_for_mdrun, /* Arguments for mdrun */
812 const char *simulation_tpr, /* This tpr will be simulated */
813 int nnodes, /* Number of ranks to use */
814 int nPMEnodes, /* Number of PME ranks to use */
815 const t_eligible_gpu_ids *gpu_ids) /* Struct containing GPU IDs for
816 * constructing mdrun command lines */
818 char *command, *cmd_gpu_ids;
821 /* Make enough space for the system call command,
822 * (200 extra chars for -npme ... etc. options should suffice): */
823 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+200);
825 cmd_gpu_ids = make_gpu_id_command_line(nnodes, nPMEnodes, gpu_ids);
827 /* Note that the -passall options requires args_for_mdrun to be at the end
828 * of the command line string */
831 sprintf(command, "%s%s-npme %d -s %s %s %s",
832 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids);
836 sprintf(command, "%s%s%s -npme %d -s %s %s %s",
837 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids);
840 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
844 /* Now the real thing! */
847 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
850 gmx_system_call(command);
855 static void modify_PMEsettings(
856 gmx_int64_t simsteps, /* Set this value as number of time steps */
857 gmx_int64_t init_step, /* Set this value as init_step */
858 const char *fn_best_tpr, /* tpr file with the best performance */
859 const char *fn_sim_tpr) /* name of tpr file to be launched */
867 read_tpx_state(fn_best_tpr, ir, &state, NULL, &mtop);
869 /* Reset nsteps and init_step to the value of the input .tpr file */
870 ir->nsteps = simsteps;
871 ir->init_step = init_step;
873 /* Write the tpr file which will be launched */
874 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%"GMX_PRId64);
875 fprintf(stdout, buf, ir->nsteps);
877 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
882 static gmx_bool can_scale_rvdw(int vdwtype)
884 return (evdwCUT == vdwtype ||
888 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
890 /* Make additional TPR files with more computational load for the
891 * direct space processors: */
892 static void make_benchmark_tprs(
893 const char *fn_sim_tpr, /* READ : User-provided tpr file */
894 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
895 gmx_int64_t benchsteps, /* Number of time steps for benchmark runs */
896 gmx_int64_t statesteps, /* Step counter in checkpoint file */
897 real rmin, /* Minimal Coulomb radius */
898 real rmax, /* Maximal Coulomb radius */
899 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
900 int *ntprs, /* No. of TPRs to write, each with a different
901 rcoulomb and fourierspacing */
902 t_inputinfo *info, /* Contains information about mdp file options */
903 FILE *fp) /* Write the output here */
909 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
912 gmx_bool bNote = FALSE;
913 real add; /* Add this to rcoul for the next test */
914 real fac = 1.0; /* Scaling factor for Coulomb radius */
915 real fourierspacing; /* Basic fourierspacing from tpr */
918 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
919 *ntprs > 1 ? "s" : "", "%"GMX_PRId64, benchsteps > 1 ? "s" : "");
920 fprintf(stdout, buf, benchsteps);
923 sprintf(buf, " (adding %s steps from checkpoint file)", "%"GMX_PRId64);
924 fprintf(stdout, buf, statesteps);
925 benchsteps += statesteps;
927 fprintf(stdout, ".\n");
931 read_tpx_state(fn_sim_tpr, ir, &state, NULL, &mtop);
933 /* Check if some kind of PME was chosen */
934 if (EEL_PME(ir->coulombtype) == FALSE)
936 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
940 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
941 if ( (ir->cutoff_scheme != ecutsVERLET) &&
942 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
944 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
945 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
947 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
948 else if (ir->rcoulomb > ir->rlist)
950 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
951 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
954 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
956 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
960 /* Reduce the number of steps for the benchmarks */
961 info->orig_sim_steps = ir->nsteps;
962 ir->nsteps = benchsteps;
963 /* We must not use init_step from the input tpr file for the benchmarks */
964 info->orig_init_step = ir->init_step;
967 /* For PME-switch potentials, keep the radial distance of the buffer region */
968 nlist_buffer = ir->rlist - ir->rcoulomb;
970 /* Determine length of triclinic box vectors */
971 for (d = 0; d < DIM; d++)
974 for (i = 0; i < DIM; i++)
976 box_size[d] += state.box[d][i]*state.box[d][i];
978 box_size[d] = sqrt(box_size[d]);
981 if (ir->fourier_spacing > 0)
983 info->fsx[0] = ir->fourier_spacing;
984 info->fsy[0] = ir->fourier_spacing;
985 info->fsz[0] = ir->fourier_spacing;
989 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
990 info->fsx[0] = box_size[XX]/ir->nkx;
991 info->fsy[0] = box_size[YY]/ir->nky;
992 info->fsz[0] = box_size[ZZ]/ir->nkz;
995 /* If no value for the fourierspacing was provided on the command line, we
996 * use the reconstruction from the tpr file */
997 if (ir->fourier_spacing > 0)
999 /* Use the spacing from the tpr */
1000 fourierspacing = ir->fourier_spacing;
1004 /* Use the maximum observed spacing */
1005 fourierspacing = max(max(info->fsx[0], info->fsy[0]), info->fsz[0]);
1008 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
1010 /* For performance comparisons the number of particles is useful to have */
1011 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
1013 /* Print information about settings of which some are potentially modified: */
1014 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
1015 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
1016 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
1017 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
1018 if (ir_vdw_switched(ir))
1020 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
1022 if (EPME_SWITCHED(ir->coulombtype))
1024 fprintf(fp, " rlist : %f nm\n", ir->rlist);
1026 if (ir->rlistlong != max_cutoff(ir->rvdw, ir->rcoulomb))
1028 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
1031 /* Print a descriptive line about the tpr settings tested */
1032 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
1033 fprintf(fp, " No. scaling rcoulomb");
1034 fprintf(fp, " nkx nky nkz");
1035 fprintf(fp, " spacing");
1036 if (can_scale_rvdw(ir->vdwtype))
1038 fprintf(fp, " rvdw");
1040 if (EPME_SWITCHED(ir->coulombtype))
1042 fprintf(fp, " rlist");
1044 if (ir->rlistlong != max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb)) )
1046 fprintf(fp, " rlistlong");
1048 fprintf(fp, " tpr file\n");
1050 /* Loop to create the requested number of tpr input files */
1051 for (j = 0; j < *ntprs; j++)
1053 /* The first .tpr is the provided one, just need to modify nsteps,
1054 * so skip the following block */
1057 /* Determine which Coulomb radii rc to use in the benchmarks */
1058 add = (rmax-rmin)/(*ntprs-1);
1059 if (gmx_within_tol(rmin, info->rcoulomb[0], GMX_REAL_EPS))
1061 ir->rcoulomb = rmin + j*add;
1063 else if (gmx_within_tol(rmax, info->rcoulomb[0], GMX_REAL_EPS))
1065 ir->rcoulomb = rmin + (j-1)*add;
1069 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
1070 add = (rmax-rmin)/(*ntprs-2);
1071 ir->rcoulomb = rmin + (j-1)*add;
1074 /* Determine the scaling factor fac */
1075 fac = ir->rcoulomb/info->rcoulomb[0];
1077 /* Scale the Fourier grid spacing */
1078 ir->nkx = ir->nky = ir->nkz = 0;
1079 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
1081 /* Adjust other radii since various conditions need to be fulfilled */
1082 if (eelPME == ir->coulombtype)
1084 /* plain PME, rcoulomb must be equal to rlist */
1085 ir->rlist = ir->rcoulomb;
1089 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1090 ir->rlist = ir->rcoulomb + nlist_buffer;
1093 if (bScaleRvdw && can_scale_rvdw(ir->vdwtype))
1095 if (ecutsVERLET == ir->cutoff_scheme ||
1096 evdwPME == ir->vdwtype)
1098 /* With either the Verlet cutoff-scheme or LJ-PME,
1099 the van der Waals radius must always equal the
1101 ir->rvdw = ir->rcoulomb;
1105 /* For vdw cutoff, rvdw >= rlist */
1106 ir->rvdw = max(info->rvdw[0], ir->rlist);
1110 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1112 } /* end of "if (j != 0)" */
1114 /* for j==0: Save the original settings
1115 * for j >0: Save modified radii and Fourier grids */
1116 info->rcoulomb[j] = ir->rcoulomb;
1117 info->rvdw[j] = ir->rvdw;
1118 info->nkx[j] = ir->nkx;
1119 info->nky[j] = ir->nky;
1120 info->nkz[j] = ir->nkz;
1121 info->rlist[j] = ir->rlist;
1122 info->rlistlong[j] = ir->rlistlong;
1123 info->fsx[j] = fac*fourierspacing;
1124 info->fsy[j] = fac*fourierspacing;
1125 info->fsz[j] = fac*fourierspacing;
1127 /* Write the benchmark tpr file */
1128 strncpy(fn_bench_tprs[j], fn_sim_tpr, strlen(fn_sim_tpr)-strlen(".tpr"));
1129 sprintf(buf, "_bench%.2d.tpr", j);
1130 strcat(fn_bench_tprs[j], buf);
1131 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1132 fprintf(stdout, "%"GMX_PRId64, ir->nsteps);
1135 fprintf(stdout, ", scaling factor %f\n", fac);
1139 fprintf(stdout, ", unmodified settings\n");
1142 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1144 /* Write information about modified tpr settings to log file */
1145 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1146 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1147 fprintf(fp, " %9f ", info->fsx[j]);
1148 if (can_scale_rvdw(ir->vdwtype))
1150 fprintf(fp, "%10f", ir->rvdw);
1152 if (EPME_SWITCHED(ir->coulombtype))
1154 fprintf(fp, "%10f", ir->rlist);
1156 if (info->rlistlong[0] != max_cutoff(info->rlist[0], max_cutoff(info->rvdw[0], info->rcoulomb[0])) )
1158 fprintf(fp, "%10f", ir->rlistlong);
1160 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1162 /* Make it clear to the user that some additional settings were modified */
1163 if (!gmx_within_tol(ir->rvdw, info->rvdw[0], GMX_REAL_EPS)
1164 || !gmx_within_tol(ir->rlistlong, info->rlistlong[0], GMX_REAL_EPS) )
1171 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1172 "other input settings were also changed (see table above).\n"
1173 "Please check if the modified settings are appropriate.\n");
1181 /* Rename the files we want to keep to some meaningful filename and
1182 * delete the rest */
1183 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1184 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1186 char numstring[STRLEN];
1187 char newfilename[STRLEN];
1188 const char *fn = NULL;
1193 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1195 for (i = 0; i < nfile; i++)
1197 opt = (char *)fnm[i].opt;
1198 if (strcmp(opt, "-p") == 0)
1200 /* do nothing; keep this file */
1203 else if (strcmp(opt, "-bg") == 0)
1205 /* Give the log file a nice name so one can later see which parameters were used */
1206 numstring[0] = '\0';
1209 sprintf(numstring, "_%d", nr);
1211 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1212 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1214 fprintf(stdout, "renaming log file to %s\n", newfilename);
1215 make_backup(newfilename);
1216 rename(opt2fn("-bg", nfile, fnm), newfilename);
1219 else if (strcmp(opt, "-err") == 0)
1221 /* This file contains the output of stderr. We want to keep it in
1222 * cases where there have been problems. */
1223 fn = opt2fn(opt, nfile, fnm);
1224 numstring[0] = '\0';
1227 sprintf(numstring, "_%d", nr);
1229 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1234 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1235 make_backup(newfilename);
1236 rename(fn, newfilename);
1240 fprintf(stdout, "Deleting %s\n", fn);
1245 /* Delete the files which are created for each benchmark run: (options -b*) */
1246 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1248 remove_if_exists(opt2fn(opt, nfile, fnm));
1255 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1258 /* Create a list of numbers of PME nodes to test */
1259 static void make_npme_list(
1260 const char *npmevalues_opt, /* Make a complete list with all
1261 * possibilities or a short list that keeps only
1262 * reasonable numbers of PME nodes */
1263 int *nentries, /* Number of entries we put in the nPMEnodes list */
1264 int *nPMEnodes[], /* Each entry contains the value for -npme */
1265 int nnodes, /* Total number of nodes to do the tests on */
1266 int minPMEnodes, /* Minimum number of PME nodes */
1267 int maxPMEnodes) /* Maximum number of PME nodes */
1270 int min_factor = 1; /* We request that npp and npme have this minimal
1271 * largest common factor (depends on npp) */
1272 int nlistmax; /* Max. list size */
1273 int nlist; /* Actual number of entries in list */
1277 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1278 if (0 == strcmp(npmevalues_opt, "all") )
1282 else if (0 == strcmp(npmevalues_opt, "subset") )
1284 eNPME = eNpmeSubset;
1286 else /* "auto" or "range" */
1292 else if (nnodes < 128)
1294 eNPME = eNpmeReduced;
1298 eNPME = eNpmeSubset;
1302 /* Calculate how many entries we could possibly have (in case of -npme all) */
1305 nlistmax = maxPMEnodes - minPMEnodes + 3;
1306 if (0 == minPMEnodes)
1316 /* Now make the actual list which is at most of size nlist */
1317 snew(*nPMEnodes, nlistmax);
1318 nlist = 0; /* start counting again, now the real entries in the list */
1319 for (i = 0; i < nlistmax - 2; i++)
1321 npme = maxPMEnodes - i;
1332 /* For 2d PME we want a common largest factor of at least the cube
1333 * root of the number of PP nodes */
1334 min_factor = (int) pow(npp, 1.0/3.0);
1337 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1340 if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1342 (*nPMEnodes)[nlist] = npme;
1346 /* We always test 0 PME nodes and the automatic number */
1347 *nentries = nlist + 2;
1348 (*nPMEnodes)[nlist ] = 0;
1349 (*nPMEnodes)[nlist+1] = -1;
1351 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1352 for (i = 0; i < *nentries-1; i++)
1354 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1356 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1360 /* Allocate memory to store the performance data */
1361 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1366 for (k = 0; k < ntprs; k++)
1368 snew(perfdata[k], datasets);
1369 for (i = 0; i < datasets; i++)
1371 for (j = 0; j < repeats; j++)
1373 snew(perfdata[k][i].Gcycles, repeats);
1374 snew(perfdata[k][i].ns_per_day, repeats);
1375 snew(perfdata[k][i].PME_f_load, repeats);
1382 /* Check for errors on mdrun -h */
1383 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp,
1384 const t_filenm *fnm, int nfile)
1386 const char *fn = NULL;
1387 char *command, *msg;
1391 snew(command, length + 15);
1392 snew(msg, length + 500);
1394 fprintf(stdout, "Making sure the benchmarks can be executed by running just 1 step...\n");
1395 sprintf(command, "%s -nsteps 1 -quiet", mdrun_cmd_line);
1396 fprintf(stdout, "Executing '%s' ...\n", command);
1397 ret = gmx_system_call(command);
1401 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1402 * get the error message from mdrun itself */
1404 "Cannot run the first benchmark simulation! Please check the error message of\n"
1405 "mdrun for the source of the problem. Did you provide a command line\n"
1406 "argument that neither gmx tune_pme nor mdrun understands? If you're\n"
1407 "sure your command line should work, you can bypass this check with \n"
1408 "gmx tune_pme -nocheck. The failing command was:\n"
1409 "\n%s\n\n", command);
1411 fprintf(stderr, "%s", msg);
1413 fprintf(fp, "%s", msg);
1417 fprintf(stdout, "Benchmarks can be executed!\n");
1419 /* Clean up the benchmark output files we just created */
1420 fprintf(stdout, "Cleaning up ...\n");
1421 remove_if_exists(opt2fn("-bc", nfile, fnm));
1422 remove_if_exists(opt2fn("-be", nfile, fnm));
1423 remove_if_exists(opt2fn("-bcpo", nfile, fnm));
1424 remove_if_exists(opt2fn("-bg", nfile, fnm));
1430 static void do_the_tests(
1431 FILE *fp, /* General g_tune_pme output file */
1432 char **tpr_names, /* Filenames of the input files to test */
1433 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1434 int minPMEnodes, /* Min fraction of nodes to use for PME */
1435 int npme_fixed, /* If >= -1, test fixed number of PME
1437 const char *npmevalues_opt, /* Which -npme values should be tested */
1438 t_perf **perfdata, /* Here the performace data is stored */
1439 int *pmeentries, /* Entries in the nPMEnodes list */
1440 int repeats, /* Repeat each test this often */
1441 int nnodes, /* Total number of nodes = nPP + nPME */
1442 int nr_tprs, /* Total number of tpr files to test */
1443 gmx_bool bThreads, /* Threads or MPI? */
1444 char *cmd_mpirun, /* mpirun command string */
1445 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1446 char *cmd_mdrun, /* mdrun command string */
1447 char *cmd_args_bench, /* arguments for mdrun in a string */
1448 const t_filenm *fnm, /* List of filenames from command line */
1449 int nfile, /* Number of files specified on the cmdl. */
1450 int presteps, /* DLB equilibration steps, is checked */
1451 gmx_int64_t cpt_steps, /* Time step counter in the checkpoint */
1452 gmx_bool bCheck, /* Check whether benchmark mdrun works */
1453 const t_eligible_gpu_ids *gpu_ids) /* Struct containing GPU IDs for
1454 * constructing mdrun command lines */
1456 int i, nr, k, ret, count = 0, totaltests;
1457 int *nPMEnodes = NULL;
1460 char *command, *cmd_stub;
1462 gmx_bool bResetProblem = FALSE;
1463 gmx_bool bFirst = TRUE;
1464 gmx_bool bUsingGpus = 0 < gpu_ids->n;
1466 /* This string array corresponds to the eParselog enum type at the start
1468 const char* ParseLog[] = {
1470 "Logfile not found!",
1471 "No timings, logfile truncated?",
1472 "Run was terminated.",
1473 "Counters were not reset properly.",
1474 "No DD grid found for these settings.",
1475 "TPX version conflict!",
1476 "mdrun was not started in parallel!",
1477 "Number of PP ranks has a prime factor that is too large.",
1478 "The number of PP ranks did not suit the number of GPUs.",
1479 "Some GPUs were not detected or are incompatible.",
1482 char str_PME_f_load[13];
1485 /* Allocate space for the mdrun command line. 100 extra characters should
1486 be more than enough for the -npme etcetera arguments */
1487 cmdline_length = strlen(cmd_mpirun)
1490 + strlen(cmd_args_bench)
1491 + strlen(tpr_names[0]) + 100;
1492 snew(command, cmdline_length);
1493 snew(cmd_stub, cmdline_length);
1495 /* Construct the part of the command line that stays the same for all tests: */
1498 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1502 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1505 /* Create a list of numbers of PME nodes to test */
1506 if (npme_fixed < -1)
1508 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1509 nnodes, minPMEnodes, maxPMEnodes);
1515 nPMEnodes[0] = npme_fixed;
1516 fprintf(stderr, "Will use a fixed number of %d PME-only ranks.\n", nPMEnodes[0]);
1521 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1523 finalize(opt2fn("-p", nfile, fnm));
1527 /* Allocate one dataset for each tpr input file: */
1528 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1530 /*****************************************/
1531 /* Main loop over all tpr files to test: */
1532 /*****************************************/
1533 totaltests = nr_tprs*(*pmeentries)*repeats;
1534 for (k = 0; k < nr_tprs; k++)
1536 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1537 fprintf(fp, "PME ranks Gcycles ns/day PME/f Remark\n");
1538 /* Loop over various numbers of PME nodes: */
1539 for (i = 0; i < *pmeentries; i++)
1541 char *cmd_gpu_ids = NULL;
1543 pd = &perfdata[k][i];
1545 cmd_gpu_ids = make_gpu_id_command_line(nnodes, nPMEnodes[i], gpu_ids);
1547 /* Loop over the repeats for each scenario: */
1548 for (nr = 0; nr < repeats; nr++)
1550 pd->nPMEnodes = nPMEnodes[i];
1552 /* Add -npme and -s to the command line and save it. Note that
1553 * the -passall (if set) options requires cmd_args_bench to be
1554 * at the end of the command line string */
1555 snew(pd->mdrun_cmd_line, cmdline_length);
1556 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s %s",
1557 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench, cmd_gpu_ids);
1559 /* To prevent that all benchmarks fail due to a show-stopper argument
1560 * on the mdrun command line, we make a quick check first.
1561 * This check can be turned off in cases where the automatically chosen
1562 * number of PME-only ranks leads to a number of PP ranks for which no
1563 * decomposition can be found (e.g. for large prime numbers) */
1564 if (bFirst && bCheck)
1566 /* TODO This check is really for a functional
1567 * .tpr, and if we need it, it should take place
1568 * for every .tpr, and the logic for it should be
1569 * immediately inside the loop over k, not in
1570 * this inner loop. */
1571 char *temporary_cmd_line;
1573 snew(temporary_cmd_line, cmdline_length);
1574 /* TODO -npme 0 is more likely to succeed at low
1575 parallelism than the default of -npme -1, but
1576 is more likely to fail at the scaling limit
1577 when the PP domains may be too small. "mpirun
1578 -np 1 mdrun" is probably a reasonable thing to
1579 do for this check, but it'll be easier to
1580 implement that after some refactoring of how
1581 the number of MPI ranks is managed. */
1582 sprintf(temporary_cmd_line, "%s-npme 0 -nb cpu -s %s %s",
1583 cmd_stub, tpr_names[k], cmd_args_bench);
1584 make_sure_it_runs(temporary_cmd_line, cmdline_length, fp, fnm, nfile);
1588 /* Do a benchmark simulation: */
1591 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1597 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1598 (100.0*count)/totaltests,
1599 k+1, nr_tprs, i+1, *pmeentries, buf);
1600 make_backup(opt2fn("-err", nfile, fnm));
1601 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1602 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1603 gmx_system_call(command);
1605 /* Collect the performance data from the log file; also check stderr
1606 * for fatal errors */
1607 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1608 pd, nr, presteps, cpt_steps, nnodes);
1609 if ((presteps > 0) && (ret == eParselogResetProblem))
1611 bResetProblem = TRUE;
1614 if (-1 == pd->nPMEnodes)
1616 sprintf(buf, "(%3d)", pd->guessPME);
1624 if (pd->PME_f_load[nr] > 0.0)
1626 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1630 sprintf(str_PME_f_load, "%s", " - ");
1633 /* Write the data we got to disk */
1634 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1635 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1636 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1638 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1644 /* Do some cleaning up and delete the files we do not need any more */
1645 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1647 /* If the first run with this number of processors already failed, do not try again: */
1648 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1650 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1651 count += repeats-(nr+1);
1654 } /* end of repeats loop */
1656 } /* end of -npme loop */
1657 } /* end of tpr file loop */
1662 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1670 static void check_input(
1677 real maxPMEfraction,
1678 real minPMEfraction,
1680 gmx_int64_t bench_nsteps,
1681 const t_filenm *fnm,
1691 /* Make sure the input file exists */
1692 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1694 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1697 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1698 if ( (0 == strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1700 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1701 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1704 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1707 gmx_fatal(FARGS, "Number of repeats < 0!");
1710 /* Check number of nodes */
1713 gmx_fatal(FARGS, "Number of ranks/threads must be a positive integer.");
1716 /* Automatically choose -ntpr if not set */
1726 /* Set a reasonable scaling factor for rcoulomb */
1729 *rmax = rcoulomb * 1.2;
1732 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1738 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1742 /* Make shure that rmin <= rcoulomb <= rmax */
1751 if (!(*rmin <= *rmax) )
1753 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1754 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1756 /* Add test scenarios if rmin or rmax were set */
1759 if (!gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1762 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1765 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1768 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1773 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1774 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) || !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1776 *ntprs = max(*ntprs, 2);
1779 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1780 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1782 *ntprs = max(*ntprs, 3);
1787 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1792 if (gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && gmx_within_tol(rcoulomb, *rmax, GMX_REAL_EPS)) /* We have just a single rc */
1794 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1795 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1796 "with correspondingly adjusted PME grid settings\n");
1801 /* Check whether max and min fraction are within required values */
1802 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1804 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1806 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1808 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1810 if (maxPMEfraction < minPMEfraction)
1812 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1815 /* Check whether the number of steps - if it was set - has a reasonable value */
1816 if (bench_nsteps < 0)
1818 gmx_fatal(FARGS, "Number of steps must be positive.");
1821 if (bench_nsteps > 10000 || bench_nsteps < 100)
1823 fprintf(stderr, "WARNING: steps=");
1824 fprintf(stderr, "%"GMX_PRId64, bench_nsteps);
1825 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1830 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1833 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1836 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1838 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1842 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1843 * only. We need to check whether the requested number of PME-only nodes
1845 if (npme_fixed > -1)
1847 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1848 if (2*npme_fixed > nnodes)
1850 gmx_fatal(FARGS, "Cannot have more than %d PME-only ranks for a total of %d ranks (you chose %d).\n",
1851 nnodes/2, nnodes, npme_fixed);
1853 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1855 fprintf(stderr, "WARNING: Only %g percent of the ranks are assigned as PME-only ranks.\n",
1856 100.0*((real)npme_fixed / (real)nnodes));
1858 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1860 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1861 " fixed number of PME-only ranks is requested with -fix.\n");
1867 /* Returns TRUE when "opt" is needed at launch time */
1868 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1870 if (0 == strncmp(opt, "-swap", 5))
1875 /* Apart from the input .tpr and the output log files we need all options that
1876 * were set on the command line and that do not start with -b */
1877 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2)
1878 || 0 == strncmp(opt, "-err", 4) || 0 == strncmp(opt, "-p", 2) )
1887 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1888 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1890 /* Apart from the input .tpr, all files starting with "-b" are for
1891 * _b_enchmark files exclusively */
1892 if (0 == strncmp(opt, "-s", 2))
1897 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2))
1899 if (!bOptional || bSet)
1916 if (bSet) /* These are additional input files like -cpi -ei */
1929 /* Adds 'buf' to 'str' */
1930 static void add_to_string(char **str, char *buf)
1935 len = strlen(*str) + strlen(buf) + 1;
1941 /* Create the command line for the benchmark as well as for the real run */
1942 static void create_command_line_snippets(
1943 gmx_bool bAppendFiles,
1944 gmx_bool bKeepAndNumCPT,
1945 gmx_bool bResetHWay,
1949 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1950 char *cmd_args_launch[], /* command line arguments for simulation run */
1951 char extra_args[]) /* Add this to the end of the command line */
1956 char strbuf[STRLEN];
1959 /* strlen needs at least '\0' as a string: */
1960 snew(*cmd_args_bench, 1);
1961 snew(*cmd_args_launch, 1);
1962 *cmd_args_launch[0] = '\0';
1963 *cmd_args_bench[0] = '\0';
1966 /*******************************************/
1967 /* 1. Process other command line arguments */
1968 /*******************************************/
1971 /* Add equilibration steps to benchmark options */
1972 sprintf(strbuf, "-resetstep %d ", presteps);
1973 add_to_string(cmd_args_bench, strbuf);
1975 /* These switches take effect only at launch time */
1976 if (FALSE == bAppendFiles)
1978 add_to_string(cmd_args_launch, "-noappend ");
1982 add_to_string(cmd_args_launch, "-cpnum ");
1986 add_to_string(cmd_args_launch, "-resethway ");
1989 /********************/
1990 /* 2. Process files */
1991 /********************/
1992 for (i = 0; i < nfile; i++)
1994 opt = (char *)fnm[i].opt;
1995 name = opt2fn(opt, nfile, fnm);
1997 /* Strbuf contains the options, now let's sort out where we need that */
1998 sprintf(strbuf, "%s %s ", opt, name);
2000 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
2002 /* All options starting with -b* need the 'b' removed,
2003 * therefore overwrite strbuf */
2004 if (0 == strncmp(opt, "-b", 2))
2006 sprintf(strbuf, "-%s %s ", &opt[2], name);
2009 add_to_string(cmd_args_bench, strbuf);
2012 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
2014 add_to_string(cmd_args_launch, strbuf);
2018 add_to_string(cmd_args_bench, extra_args);
2019 add_to_string(cmd_args_launch, extra_args);
2023 /* Set option opt */
2024 static void setopt(const char *opt, int nfile, t_filenm fnm[])
2028 for (i = 0; (i < nfile); i++)
2030 if (strcmp(opt, fnm[i].opt) == 0)
2032 fnm[i].flag |= ffSET;
2038 /* This routine inspects the tpr file and ...
2039 * 1. checks for output files that get triggered by a tpr option. These output
2040 * files are marked as 'set' to allow for a proper cleanup after each
2042 * 2. returns the PME:PP load ratio
2043 * 3. returns rcoulomb from the tpr */
2044 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
2046 gmx_bool bTpi; /* Is test particle insertion requested? */
2047 gmx_bool bFree; /* Is a free energy simulation requested? */
2048 gmx_bool bNM; /* Is a normal mode analysis requested? */
2049 gmx_bool bSwap; /* Is water/ion position swapping requested? */
2055 /* Check tpr file for options that trigger extra output files */
2056 read_tpx_state(opt2fn("-s", nfile, fnm), &ir, &state, NULL, &mtop);
2057 bFree = (efepNO != ir.efep );
2058 bNM = (eiNM == ir.eI );
2059 bSwap = (eswapNO != ir.eSwapCoords);
2060 bTpi = EI_TPI(ir.eI);
2062 /* Set these output files on the tuning command-line */
2065 setopt("-pf", nfile, fnm);
2066 setopt("-px", nfile, fnm);
2070 setopt("-dhdl", nfile, fnm);
2074 setopt("-tpi", nfile, fnm);
2075 setopt("-tpid", nfile, fnm);
2079 setopt("-mtx", nfile, fnm);
2083 setopt("-swap", nfile, fnm);
2086 *rcoulomb = ir.rcoulomb;
2088 /* Return the estimate for the number of PME nodes */
2089 return pme_load_estimate(&mtop, &ir, state.box);
2093 static void couple_files_options(int nfile, t_filenm fnm[])
2096 gmx_bool bSet, bBench;
2101 for (i = 0; i < nfile; i++)
2103 opt = (char *)fnm[i].opt;
2104 bSet = ((fnm[i].flag & ffSET) != 0);
2105 bBench = (0 == strncmp(opt, "-b", 2));
2107 /* Check optional files */
2108 /* If e.g. -eo is set, then -beo also needs to be set */
2109 if (is_optional(&fnm[i]) && bSet && !bBench)
2111 sprintf(buf, "-b%s", &opt[1]);
2112 setopt(buf, nfile, fnm);
2114 /* If -beo is set, then -eo also needs to be! */
2115 if (is_optional(&fnm[i]) && bSet && bBench)
2117 sprintf(buf, "-%s", &opt[2]);
2118 setopt(buf, nfile, fnm);
2124 #define BENCHSTEPS (1000)
2126 int gmx_tune_pme(int argc, char *argv[])
2128 const char *desc[] = {
2129 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of ranks, [THISMODULE] systematically",
2130 "times [gmx-mdrun] with various numbers of PME-only ranks and determines",
2131 "which setting is fastest. It will also test whether performance can",
2132 "be enhanced by shifting load from the reciprocal to the real space",
2133 "part of the Ewald sum. ",
2134 "Simply pass your [REF].tpr[ref] file to [THISMODULE] together with other options",
2135 "for [gmx-mdrun] as needed.[PAR]",
2136 "Which executables are used can be set in the environment variables",
2137 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
2138 "will be used as defaults. Note that for certain MPI frameworks you",
2139 "need to provide a machine- or hostfile. This can also be passed",
2140 "via the MPIRUN variable, e.g.[PAR]",
2141 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
2142 "Before doing the actual benchmark runs, [THISMODULE] will do a quick",
2143 "check whether mdrun works as expected with the provided parallel settings",
2144 "if the [TT]-check[tt] option is activated (the default).",
2145 "Please call [THISMODULE] with the normal options you would pass to",
2146 "[gmx-mdrun] and add [TT]-np[tt] for the number of ranks to perform the",
2147 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2148 "to repeat each test several times to get better statistics. [PAR]",
2149 "[THISMODULE] can test various real space / reciprocal space workloads",
2150 "for you. With [TT]-ntpr[tt] you control how many extra [REF].tpr[ref] files will be",
2151 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2152 "Typically, the first test (number 0) will be with the settings from the input",
2153 "[REF].tpr[ref] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2154 "specified by [TT]-rmax[tt] with a somewhat smaller PME grid at the same time. ",
2155 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2156 "The remaining [REF].tpr[ref] files will have equally-spaced Coulomb radii (and Fourier "
2157 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2158 "if you just seek the optimal number of PME-only ranks; in that case",
2159 "your input [REF].tpr[ref] file will remain unchanged.[PAR]",
2160 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2161 "MD systems. The dynamic load balancing needs about 100 time steps",
2162 "to adapt to local load imbalances, therefore the time step counters",
2163 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2164 "for a higher accuracy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2165 "From the 'DD' load imbalance entries in the md.log output file you",
2166 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2167 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2168 "After calling [gmx-mdrun] several times, detailed performance information",
2169 "is available in the output file [TT]perf.out[tt].",
2170 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2171 "(options [TT]-b*[tt]), these will be automatically deleted after each test.[PAR]",
2172 "If you want the simulation to be started automatically with the",
2173 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2174 "Basic support for GPU-enabled [TT]mdrun[tt] exists. Give a string containing the IDs",
2175 "of the GPUs that you wish to use in the optimization in the [TT]-gpu_id[tt]",
2176 "command-line argument. Unlike [TT]mdrun -gpu_id[tt], this does not imply a mapping",
2177 "but merely the eligible set. [TT]g_tune_pme[tt] will construct calls to",
2178 "mdrun that use this set appropriately, assuming that PP ranks with low indices",
2179 "should map to GPUs with low indices, and increasing both monotonically",
2180 "over the respective sets.[PAR]",
2185 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2186 real maxPMEfraction = 0.50;
2187 real minPMEfraction = 0.25;
2188 int maxPMEnodes, minPMEnodes;
2189 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2190 float guessPMEnodes;
2191 int npme_fixed = -2; /* If >= -1, use only this number
2192 * of PME-only nodes */
2194 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2195 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2196 gmx_bool bScaleRvdw = TRUE;
2197 gmx_int64_t bench_nsteps = BENCHSTEPS;
2198 gmx_int64_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2199 gmx_int64_t cpt_steps = 0; /* Step counter in .cpt input file */
2200 int presteps = 100; /* Do a full cycle reset after presteps steps */
2201 gmx_bool bOverwrite = FALSE, bKeepTPR;
2202 gmx_bool bLaunch = FALSE;
2203 char *ExtraArgs = NULL;
2204 char **tpr_names = NULL;
2205 const char *simulation_tpr = NULL;
2206 int best_npme, best_tpr;
2207 int sim_part = 1; /* For benchmarks with checkpoint files */
2211 /* Default program names if nothing else is found */
2212 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2213 char *cmd_args_bench, *cmd_args_launch;
2214 char *cmd_np = NULL;
2216 /* IDs of GPUs that are eligible for computation */
2217 char *eligible_gpu_ids = NULL;
2218 t_eligible_gpu_ids *gpu_ids = NULL;
2220 t_perf **perfdata = NULL;
2226 /* Print out how long the tuning took */
2229 static t_filenm fnm[] = {
2231 { efOUT, "-p", "perf", ffWRITE },
2232 { efLOG, "-err", "bencherr", ffWRITE },
2233 { efTPR, "-so", "tuned", ffWRITE },
2235 { efTPR, NULL, NULL, ffREAD },
2236 { efTRN, "-o", NULL, ffWRITE },
2237 { efCOMPRESSED, "-x", NULL, ffOPTWR },
2238 { efCPT, "-cpi", NULL, ffOPTRD },
2239 { efCPT, "-cpo", NULL, ffOPTWR },
2240 { efSTO, "-c", "confout", ffWRITE },
2241 { efEDR, "-e", "ener", ffWRITE },
2242 { efLOG, "-g", "md", ffWRITE },
2243 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2244 { efXVG, "-field", "field", ffOPTWR },
2245 { efXVG, "-table", "table", ffOPTRD },
2246 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2247 { efXVG, "-tablep", "tablep", ffOPTRD },
2248 { efXVG, "-tableb", "table", ffOPTRD },
2249 { efTRX, "-rerun", "rerun", ffOPTRD },
2250 { efXVG, "-tpi", "tpi", ffOPTWR },
2251 { efXVG, "-tpid", "tpidist", ffOPTWR },
2252 { efEDI, "-ei", "sam", ffOPTRD },
2253 { efXVG, "-eo", "edsam", ffOPTWR },
2254 { efXVG, "-devout", "deviatie", ffOPTWR },
2255 { efXVG, "-runav", "runaver", ffOPTWR },
2256 { efXVG, "-px", "pullx", ffOPTWR },
2257 { efXVG, "-pf", "pullf", ffOPTWR },
2258 { efXVG, "-ro", "rotation", ffOPTWR },
2259 { efLOG, "-ra", "rotangles", ffOPTWR },
2260 { efLOG, "-rs", "rotslabs", ffOPTWR },
2261 { efLOG, "-rt", "rottorque", ffOPTWR },
2262 { efMTX, "-mtx", "nm", ffOPTWR },
2263 { efNDX, "-dn", "dipole", ffOPTWR },
2264 { efXVG, "-swap", "swapions", ffOPTWR },
2265 /* Output files that are deleted after each benchmark run */
2266 { efTRN, "-bo", "bench", ffWRITE },
2267 { efXTC, "-bx", "bench", ffWRITE },
2268 { efCPT, "-bcpo", "bench", ffWRITE },
2269 { efSTO, "-bc", "bench", ffWRITE },
2270 { efEDR, "-be", "bench", ffWRITE },
2271 { efLOG, "-bg", "bench", ffWRITE },
2272 { efXVG, "-beo", "benchedo", ffOPTWR },
2273 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2274 { efXVG, "-bfield", "benchfld", ffOPTWR },
2275 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2276 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2277 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2278 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2279 { efXVG, "-bpx", "benchpx", ffOPTWR },
2280 { efXVG, "-bpf", "benchpf", ffOPTWR },
2281 { efXVG, "-bro", "benchrot", ffOPTWR },
2282 { efLOG, "-bra", "benchrota", ffOPTWR },
2283 { efLOG, "-brs", "benchrots", ffOPTWR },
2284 { efLOG, "-brt", "benchrott", ffOPTWR },
2285 { efMTX, "-bmtx", "benchn", ffOPTWR },
2286 { efNDX, "-bdn", "bench", ffOPTWR },
2287 { efXVG, "-bswap", "benchswp", ffOPTWR }
2290 gmx_bool bThreads = FALSE;
2294 const char *procstring[] =
2295 { NULL, "-np", "-n", "none", NULL };
2296 const char *npmevalues_opt[] =
2297 { NULL, "auto", "all", "subset", NULL };
2299 gmx_bool bAppendFiles = TRUE;
2300 gmx_bool bKeepAndNumCPT = FALSE;
2301 gmx_bool bResetCountersHalfWay = FALSE;
2302 gmx_bool bBenchmark = TRUE;
2303 gmx_bool bCheck = TRUE;
2305 output_env_t oenv = NULL;
2308 /***********************/
2309 /* g_tune_pme options: */
2310 /***********************/
2311 { "-np", FALSE, etINT, {&nnodes},
2312 "Number of ranks to run the tests on (must be > 2 for separate PME ranks)" },
2313 { "-npstring", FALSE, etENUM, {procstring},
2314 "Specify the number of ranks to [TT]$MPIRUN[tt] using this string"},
2315 { "-ntmpi", FALSE, etINT, {&nthreads},
2316 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2317 { "-r", FALSE, etINT, {&repeats},
2318 "Repeat each test this often" },
2319 { "-max", FALSE, etREAL, {&maxPMEfraction},
2320 "Max fraction of PME ranks to test with" },
2321 { "-min", FALSE, etREAL, {&minPMEfraction},
2322 "Min fraction of PME ranks to test with" },
2323 { "-npme", FALSE, etENUM, {npmevalues_opt},
2324 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2325 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2326 { "-fix", FALSE, etINT, {&npme_fixed},
2327 "If >= -1, do not vary the number of PME-only ranks, instead use this fixed value and only vary rcoulomb and the PME grid spacing."},
2328 { "-rmax", FALSE, etREAL, {&rmax},
2329 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2330 { "-rmin", FALSE, etREAL, {&rmin},
2331 "If >0, minimal rcoulomb for -ntpr>1" },
2332 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2333 "Scale rvdw along with rcoulomb"},
2334 { "-ntpr", FALSE, etINT, {&ntprs},
2335 "Number of [REF].tpr[ref] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2336 "If < 1, automatically choose the number of [REF].tpr[ref] files to test" },
2337 { "-steps", FALSE, etINT64, {&bench_nsteps},
2338 "Take timings for this many steps in the benchmark runs" },
2339 { "-resetstep", FALSE, etINT, {&presteps},
2340 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2341 { "-nsteps", FALSE, etINT64, {&new_sim_nsteps},
2342 "If non-negative, perform this many steps in the real run (overwrites nsteps from [REF].tpr[ref], add [REF].cpt[ref] steps)" },
2343 { "-launch", FALSE, etBOOL, {&bLaunch},
2344 "Launch the real simulation after optimization" },
2345 { "-bench", FALSE, etBOOL, {&bBenchmark},
2346 "Run the benchmarks or just create the input [REF].tpr[ref] files?" },
2347 { "-check", FALSE, etBOOL, {&bCheck},
2348 "Before the benchmark runs, check whether mdrun works in parallel" },
2349 { "-gpu_id", FALSE, etSTR, {&eligible_gpu_ids},
2350 "List of GPU device id-s that are eligible for use (unlike mdrun, does not imply any mapping)" },
2351 /******************/
2352 /* mdrun options: */
2353 /******************/
2354 /* We let g_tune_pme parse and understand these options, because we need to
2355 * prevent that they appear on the mdrun command line for the benchmarks */
2356 { "-append", FALSE, etBOOL, {&bAppendFiles},
2357 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2358 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2359 "Keep and number checkpoint files (launch only)" },
2360 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2361 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2364 #define NFILE asize(fnm)
2366 seconds = gmx_gettime();
2368 if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2369 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2375 /* Store the remaining unparsed command line entries in a string which
2376 * is then attached to the mdrun command line */
2378 ExtraArgs[0] = '\0';
2379 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2381 add_to_string(&ExtraArgs, argv[i]);
2382 add_to_string(&ExtraArgs, " ");
2385 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2388 if (opt2parg_bSet("-npstring", asize(pa), pa))
2390 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2395 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2397 /* and now we just set this; a bit of an ugly hack*/
2400 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2401 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2403 /* Automatically set -beo options if -eo is set etc. */
2404 couple_files_options(NFILE, fnm);
2406 /* Construct the command line arguments for benchmark runs
2407 * as well as for the simulation run */
2410 sprintf(bbuf, " -ntmpi %d ", nthreads);
2414 /* This string will be used for MPI runs and will appear after the
2415 * mpirun command. */
2416 if (strcmp(procstring[0], "none") != 0)
2418 sprintf(bbuf, " %s %d ", procstring[0], nnodes);
2428 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2429 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs);
2431 /* Prepare to use checkpoint file if requested */
2433 if (opt2bSet("-cpi", NFILE, fnm))
2435 const char *filename = opt2fn("-cpi", NFILE, fnm);
2437 read_checkpoint_part_and_step(filename,
2438 &cpt_sim_part, &cpt_steps);
2439 if (cpt_sim_part == 0)
2441 gmx_fatal(FARGS, "Checkpoint file %s could not be read!", filename);
2443 /* tune_pme will run the next part of the simulation */
2444 sim_part = cpt_sim_part + 1;
2447 /* Open performance output file and write header info */
2448 fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2450 /* Make a quick consistency check of command line parameters */
2451 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2452 maxPMEfraction, minPMEfraction, npme_fixed,
2453 bench_nsteps, fnm, NFILE, sim_part, presteps,
2455 /* Check any GPU IDs passed make sense, and fill the data structure for them */
2457 parse_digits_from_plain_string(eligible_gpu_ids, &gpu_ids->n, &gpu_ids->ids);
2459 /* Determine the maximum and minimum number of PME nodes to test,
2460 * the actual list of settings is build in do_the_tests(). */
2461 if ((nnodes > 2) && (npme_fixed < -1))
2463 if (0 == strcmp(npmevalues_opt[0], "auto"))
2465 /* Determine the npme range automatically based on the PME:PP load guess */
2466 if (guessPMEratio > 1.0)
2468 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2469 maxPMEnodes = nnodes/2;
2470 minPMEnodes = nnodes/2;
2474 /* PME : PP load is in the range 0..1, let's test around the guess */
2475 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2476 minPMEnodes = floor(0.7*guessPMEnodes);
2477 maxPMEnodes = ceil(1.6*guessPMEnodes);
2478 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2483 /* Determine the npme range based on user input */
2484 maxPMEnodes = floor(maxPMEfraction*nnodes);
2485 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2486 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2487 if (maxPMEnodes != minPMEnodes)
2489 fprintf(stdout, "- %d ", maxPMEnodes);
2491 fprintf(stdout, "PME-only ranks.\n Note that the automatic number of PME-only ranks and no separate PME ranks are always tested.\n");
2500 /* Get the commands we need to set up the runs from environment variables */
2501 get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2502 if (bBenchmark && repeats > 0)
2504 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun, NULL != eligible_gpu_ids);
2507 /* Print some header info to file */
2509 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2511 fprintf(fp, "%s for GROMACS %s\n", output_env_get_program_display_name(oenv),
2515 fprintf(fp, "Number of ranks : %d\n", nnodes);
2516 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2517 if (strcmp(procstring[0], "none") != 0)
2519 fprintf(fp, "Passing # of ranks via : %s\n", procstring[0]);
2523 fprintf(fp, "Not setting number of ranks in system call\n");
2528 fprintf(fp, "Number of threads : %d\n", nnodes);
2531 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2532 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2533 fprintf(fp, "Benchmark steps : ");
2534 fprintf(fp, "%"GMX_PRId64, bench_nsteps);
2536 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2539 fprintf(fp, "Checkpoint time step : ");
2540 fprintf(fp, "%"GMX_PRId64, cpt_steps);
2543 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2545 if (new_sim_nsteps >= 0)
2548 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2549 fprintf(stderr, "%"GMX_PRId64, new_sim_nsteps+cpt_steps);
2550 fprintf(stderr, " steps.\n");
2551 fprintf(fp, "Simulation steps : ");
2552 fprintf(fp, "%"GMX_PRId64, new_sim_nsteps);
2557 fprintf(fp, "Repeats for each test : %d\n", repeats);
2560 if (npme_fixed >= -1)
2562 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2565 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2566 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2568 /* Allocate memory for the inputinfo struct: */
2570 info->nr_inputfiles = ntprs;
2571 for (i = 0; i < ntprs; i++)
2573 snew(info->rcoulomb, ntprs);
2574 snew(info->rvdw, ntprs);
2575 snew(info->rlist, ntprs);
2576 snew(info->rlistlong, ntprs);
2577 snew(info->nkx, ntprs);
2578 snew(info->nky, ntprs);
2579 snew(info->nkz, ntprs);
2580 snew(info->fsx, ntprs);
2581 snew(info->fsy, ntprs);
2582 snew(info->fsz, ntprs);
2584 /* Make alternative tpr files to test: */
2585 snew(tpr_names, ntprs);
2586 for (i = 0; i < ntprs; i++)
2588 snew(tpr_names[i], STRLEN);
2591 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2592 * different grids could be found. */
2593 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2594 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2596 /********************************************************************************/
2597 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2598 /********************************************************************************/
2599 snew(perfdata, ntprs);
2602 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2603 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2604 cmd_args_bench, fnm, NFILE, presteps, cpt_steps, bCheck, gpu_ids);
2606 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2608 /* Analyse the results and give a suggestion for optimal settings: */
2609 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2610 repeats, info, &best_tpr, &best_npme);
2612 /* Take the best-performing tpr file and enlarge nsteps to original value */
2613 if (bKeepTPR && !bOverwrite)
2615 simulation_tpr = opt2fn("-s", NFILE, fnm);
2619 simulation_tpr = opt2fn("-so", NFILE, fnm);
2620 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2621 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2624 /* Let's get rid of the temporary benchmark input files */
2625 for (i = 0; i < ntprs; i++)
2627 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2628 remove(tpr_names[i]);
2631 /* Now start the real simulation if the user requested it ... */
2632 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2633 cmd_args_launch, simulation_tpr, nnodes, best_npme, gpu_ids);
2637 /* ... or simply print the performance results to screen: */
2640 finalize(opt2fn("-p", NFILE, fnm));