2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
41 #ifdef HAVE_SYS_TIME_H
45 #include "gromacs/commandline/pargs.h"
47 #include "types/commrec.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/math/vec.h"
51 #include "gromacs/fileio/tpxio.h"
52 #include "gromacs/utility/cstringutil.h"
55 #include "checkpoint.h"
61 #include "gromacs/timing/walltime_accounting.h"
62 #include "gromacs/math/utilities.h"
64 #include "gromacs/utility/fatalerror.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[] */
74 eParselogResetProblem,
78 eParselogLargePrimeFactor,
86 int nPMEnodes; /* number of PME-only nodes used in this test */
87 int nx, ny, nz; /* DD grid */
88 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
89 double *Gcycles; /* This can contain more than one value if doing multiple tests */
93 float *PME_f_load; /* PME mesh/force load average*/
94 float PME_f_load_Av; /* Average average ;) ... */
95 char *mdrun_cmd_line; /* Mdrun command line used for this test */
101 int nr_inputfiles; /* The number of tpr and mdp input files */
102 gmx_int64_t orig_sim_steps; /* Number of steps to be done in the real simulation */
103 gmx_int64_t orig_init_step; /* Init step for the real simulation */
104 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
105 real *rvdw; /* The vdW radii */
106 real *rlist; /* Neighbourlist cutoff radius */
108 int *nkx, *nky, *nkz;
109 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
113 static void sep_line(FILE *fp)
115 fprintf(fp, "\n------------------------------------------------------------\n");
119 /* Wrapper for system calls */
120 static int gmx_system_call(char *command)
123 gmx_fatal(FARGS, "No calls to system(3) supported on this platform. Attempted to call:\n'%s'\n", command);
125 return ( system(command) );
130 /* Check if string starts with substring */
131 static gmx_bool str_starts(const char *string, const char *substring)
133 return ( strncmp(string, substring, strlen(substring)) == 0);
137 static void cleandata(t_perf *perfdata, int test_nr)
139 perfdata->Gcycles[test_nr] = 0.0;
140 perfdata->ns_per_day[test_nr] = 0.0;
141 perfdata->PME_f_load[test_nr] = 0.0;
147 static void remove_if_exists(const char *fn)
151 fprintf(stdout, "Deleting %s\n", fn);
157 static void finalize(const char *fn_out)
163 fp = fopen(fn_out, "r");
164 fprintf(stdout, "\n\n");
166 while (fgets(buf, STRLEN-1, fp) != NULL)
168 fprintf(stdout, "%s", buf);
171 fprintf(stdout, "\n\n");
176 eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr
179 static int parse_logfile(const char *logfile, const char *errfile,
180 t_perf *perfdata, int test_nr, int presteps, gmx_int64_t cpt_steps,
184 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
185 const char matchstrdd[] = "Domain decomposition grid";
186 const char matchstrcr[] = "resetting all time and cycle counters";
187 const char matchstrbal[] = "Average PME mesh/force load:";
188 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";
189 const char errSIG[] = "signal, stopping at the next";
191 float dum1, dum2, dum3, dum4;
194 gmx_int64_t resetsteps = -1;
195 gmx_bool bFoundResetStr = FALSE;
196 gmx_bool bResetChecked = FALSE;
199 if (!gmx_fexist(logfile))
201 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
202 cleandata(perfdata, test_nr);
203 return eParselogNotFound;
206 fp = fopen(logfile, "r");
207 perfdata->PME_f_load[test_nr] = -1.0;
208 perfdata->guessPME = -1;
210 iFound = eFoundNothing;
213 iFound = eFoundDDStr; /* Skip some case statements */
216 while (fgets(line, STRLEN, fp) != NULL)
218 /* Remove leading spaces */
221 /* Check for TERM and INT signals from user: */
222 if (strstr(line, errSIG) != NULL)
225 cleandata(perfdata, test_nr);
226 return eParselogTerm;
229 /* Check whether cycle resetting worked */
230 if (presteps > 0 && !bFoundResetStr)
232 if (strstr(line, matchstrcr) != NULL)
234 sprintf(dumstring, "step %s", "%"GMX_SCNd64);
235 sscanf(line, dumstring, &resetsteps);
236 bFoundResetStr = TRUE;
237 if (resetsteps == presteps+cpt_steps)
239 bResetChecked = TRUE;
243 sprintf(dumstring, "%"GMX_PRId64, resetsteps);
244 sprintf(dumstring2, "%"GMX_PRId64, presteps+cpt_steps);
245 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
246 " though they were supposed to be reset at step %s!\n",
247 dumstring, dumstring2);
252 /* Look for strings that appear in a certain order in the log file: */
256 /* Look for domain decomp grid and separate PME nodes: */
257 if (str_starts(line, matchstrdd))
259 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME nodes %d",
260 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
261 if (perfdata->nPMEnodes == -1)
263 perfdata->guessPME = npme;
265 else if (perfdata->nPMEnodes != npme)
267 gmx_fatal(FARGS, "PME nodes from command line and output file are not identical");
269 iFound = eFoundDDStr;
271 /* Catch a few errors that might have occured: */
272 else if (str_starts(line, "There is no domain decomposition for"))
275 return eParselogNoDDGrid;
277 else if (str_starts(line, "The number of nodes you selected"))
280 return eParselogLargePrimeFactor;
282 else if (str_starts(line, "reading tpx file"))
285 return eParselogTPXVersion;
287 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
290 return eParselogNotParallel;
294 /* Look for PME mesh/force balance (not necessarily present, though) */
295 if (str_starts(line, matchstrbal))
297 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
299 /* Look for matchstring */
300 if (str_starts(line, matchstring))
302 iFound = eFoundAccountingStr;
305 case eFoundAccountingStr:
306 /* Already found matchstring - look for cycle data */
307 if (str_starts(line, "Total "))
309 sscanf(line, "Total %*f %lf", &(perfdata->Gcycles[test_nr]));
310 iFound = eFoundCycleStr;
314 /* Already found cycle data - look for remaining performance info and return */
315 if (str_starts(line, "Performance:"))
317 ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
318 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
319 perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
321 if (bResetChecked || presteps == 0)
327 return eParselogResetProblem;
334 /* Close the log file */
337 /* Check why there is no performance data in the log file.
338 * Did a fatal errors occur? */
339 if (gmx_fexist(errfile))
341 fp = fopen(errfile, "r");
342 while (fgets(line, STRLEN, fp) != NULL)
344 if (str_starts(line, "Fatal error:") )
346 if (fgets(line, STRLEN, fp) != NULL)
348 fprintf(stderr, "\nWARNING: An error occured during this benchmark:\n"
352 cleandata(perfdata, test_nr);
353 return eParselogFatal;
360 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
363 /* Giving up ... we could not find out why there is no performance data in
365 fprintf(stdout, "No performance data in log file.\n");
366 cleandata(perfdata, test_nr);
368 return eParselogNoPerfData;
372 static gmx_bool analyze_data(
381 int *index_tpr, /* OUT: Nr of mdp file with best settings */
382 int *npme_optimal) /* OUT: Optimal number of PME nodes */
385 int line = 0, line_win = -1;
386 int k_win = -1, i_win = -1, winPME;
387 double s = 0.0; /* standard deviation */
390 char str_PME_f_load[13];
391 gmx_bool bCanUseOrigTPR;
392 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
398 fprintf(fp, "Summary of successful runs:\n");
399 fprintf(fp, "Line tpr PME nodes Gcycles Av. Std.dev. ns/day PME/f");
402 fprintf(fp, " DD grid");
408 for (k = 0; k < ntprs; k++)
410 for (i = 0; i < ntests; i++)
412 /* Select the right dataset: */
413 pd = &(perfdata[k][i]);
415 pd->Gcycles_Av = 0.0;
416 pd->PME_f_load_Av = 0.0;
417 pd->ns_per_day_Av = 0.0;
419 if (pd->nPMEnodes == -1)
421 sprintf(strbuf, "(%3d)", pd->guessPME);
425 sprintf(strbuf, " ");
428 /* Get the average run time of a setting */
429 for (j = 0; j < nrepeats; j++)
431 pd->Gcycles_Av += pd->Gcycles[j];
432 pd->PME_f_load_Av += pd->PME_f_load[j];
434 pd->Gcycles_Av /= nrepeats;
435 pd->PME_f_load_Av /= nrepeats;
437 for (j = 0; j < nrepeats; j++)
439 if (pd->ns_per_day[j] > 0.0)
441 pd->ns_per_day_Av += pd->ns_per_day[j];
445 /* Somehow the performance number was not aquired for this run,
446 * therefor set the average to some negative value: */
447 pd->ns_per_day_Av = -1.0f*nrepeats;
451 pd->ns_per_day_Av /= nrepeats;
454 if (pd->PME_f_load_Av > 0.0)
456 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
460 sprintf(str_PME_f_load, "%s", " - ");
464 /* We assume we had a successful run if both averages are positive */
465 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
467 /* Output statistics if repeats were done */
470 /* Calculate the standard deviation */
472 for (j = 0; j < nrepeats; j++)
474 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
479 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
480 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
481 pd->ns_per_day_Av, str_PME_f_load);
484 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
488 /* Store the index of the best run found so far in 'winner': */
489 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
502 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
507 winPME = perfdata[k_win][i_win].nPMEnodes;
511 /* We stuck to a fixed number of PME-only nodes */
512 sprintf(strbuf, "settings No. %d", k_win);
516 /* We have optimized the number of PME-only nodes */
519 sprintf(strbuf, "%s", "the automatic number of PME nodes");
523 sprintf(strbuf, "%d PME nodes", winPME);
526 fprintf(fp, "Best performance was achieved with %s", strbuf);
527 if ((nrepeats > 1) && (ntests > 1))
529 fprintf(fp, " (see line %d)", line_win);
533 /* Only mention settings if they were modified: */
534 bRefinedCoul = !gmx_within_tol(info->rcoulomb[k_win], info->rcoulomb[0], GMX_REAL_EPS);
535 bRefinedVdW = !gmx_within_tol(info->rvdw[k_win], info->rvdw[0], GMX_REAL_EPS);
536 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
537 info->nky[k_win] == info->nky[0] &&
538 info->nkz[k_win] == info->nkz[0]);
540 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
542 fprintf(fp, "Optimized PME settings:\n");
543 bCanUseOrigTPR = FALSE;
547 bCanUseOrigTPR = TRUE;
552 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
557 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
562 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],
563 info->nkx[0], info->nky[0], info->nkz[0]);
566 if (bCanUseOrigTPR && ntprs > 1)
568 fprintf(fp, "and original PME settings.\n");
573 /* Return the index of the mdp file that showed the highest performance
574 * and the optimal number of PME nodes */
576 *npme_optimal = winPME;
578 return bCanUseOrigTPR;
582 /* Get the commands we need to set up the runs from environment variables */
583 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
587 const char def_mpirun[] = "mpirun";
588 const char def_mdrun[] = "mdrun";
590 const char empty_mpirun[] = "";
592 /* Get the commands we need to set up the runs from environment variables */
595 if ( (cp = getenv("MPIRUN")) != NULL)
597 *cmd_mpirun = strdup(cp);
601 *cmd_mpirun = strdup(def_mpirun);
606 *cmd_mpirun = strdup(empty_mpirun);
609 if ( (cp = getenv("MDRUN" )) != NULL)
611 *cmd_mdrun = strdup(cp);
615 *cmd_mdrun = strdup(def_mdrun);
619 /* Check that the commands will run mdrun (perhaps via mpirun) by
620 * running a very quick test simulation. Requires MPI environment to
621 * be available if applicable. */
622 static void check_mdrun_works(gmx_bool bThreads,
623 const char *cmd_mpirun,
625 const char *cmd_mdrun)
627 char *command = NULL;
631 const char filename[] = "benchtest.log";
633 /* This string should always be identical to the one in copyrite.c,
634 * gmx_print_version_info() in the defined(GMX_MPI) section */
635 const char match_mpi[] = "MPI library: MPI";
636 const char match_mdrun[] = "Executable: ";
637 gmx_bool bMdrun = FALSE;
638 gmx_bool bMPI = FALSE;
640 /* Run a small test to see whether mpirun + mdrun work */
641 fprintf(stdout, "Making sure that mdrun can be executed. ");
644 snew(command, strlen(cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
645 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
649 snew(command, strlen(cmd_mpirun) + strlen(cmd_np) + strlen(cmd_mdrun) + strlen(filename) + 50);
650 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
652 fprintf(stdout, "Trying '%s' ... ", command);
653 make_backup(filename);
654 gmx_system_call(command);
656 /* Check if we find the characteristic string in the output: */
657 if (!gmx_fexist(filename))
659 gmx_fatal(FARGS, "Output from test run could not be found.");
662 fp = fopen(filename, "r");
663 /* We need to scan the whole output file, since sometimes the queuing system
664 * also writes stuff to stdout/err */
667 cp = fgets(line, STRLEN, fp);
670 if (str_starts(line, match_mdrun) )
674 if (str_starts(line, match_mpi) )
686 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
688 "seems to have been compiled with MPI instead.",
696 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
698 "seems to have been compiled without MPI support.",
705 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
709 fprintf(stdout, "passed.\n");
717 static void launch_simulation(
718 gmx_bool bLaunch, /* Should the simulation be launched? */
719 FILE *fp, /* General log file */
720 gmx_bool bThreads, /* whether to use threads */
721 char *cmd_mpirun, /* Command for mpirun */
722 char *cmd_np, /* Switch for -np or -ntmpi or empty */
723 char *cmd_mdrun, /* Command for mdrun */
724 char *args_for_mdrun, /* Arguments for mdrun */
725 const char *simulation_tpr, /* This tpr will be simulated */
726 int nPMEnodes) /* Number of PME nodes to use */
731 /* Make enough space for the system call command,
732 * (100 extra chars for -npme ... etc. options should suffice): */
733 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+100);
735 /* Note that the -passall options requires args_for_mdrun to be at the end
736 * of the command line string */
739 sprintf(command, "%s%s-npme %d -s %s %s",
740 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun);
744 sprintf(command, "%s%s%s -npme %d -s %s %s",
745 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun);
748 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
752 /* Now the real thing! */
755 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
758 gmx_system_call(command);
763 static void modify_PMEsettings(
764 gmx_int64_t simsteps, /* Set this value as number of time steps */
765 gmx_int64_t init_step, /* Set this value as init_step */
766 const char *fn_best_tpr, /* tpr file with the best performance */
767 const char *fn_sim_tpr) /* name of tpr file to be launched */
775 read_tpx_state(fn_best_tpr, ir, &state, NULL, &mtop);
777 /* Reset nsteps and init_step to the value of the input .tpr file */
778 ir->nsteps = simsteps;
779 ir->init_step = init_step;
781 /* Write the tpr file which will be launched */
782 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%"GMX_PRId64);
783 fprintf(stdout, buf, ir->nsteps);
785 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
791 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
793 /* Make additional TPR files with more computational load for the
794 * direct space processors: */
795 static void make_benchmark_tprs(
796 const char *fn_sim_tpr, /* READ : User-provided tpr file */
797 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
798 gmx_int64_t benchsteps, /* Number of time steps for benchmark runs */
799 gmx_int64_t statesteps, /* Step counter in checkpoint file */
800 real rmin, /* Minimal Coulomb radius */
801 real rmax, /* Maximal Coulomb radius */
802 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
803 int *ntprs, /* No. of TPRs to write, each with a different
804 rcoulomb and fourierspacing */
805 t_inputinfo *info, /* Contains information about mdp file options */
806 FILE *fp) /* Write the output here */
812 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
815 gmx_bool bNote = FALSE;
816 real add; /* Add this to rcoul for the next test */
817 real fac = 1.0; /* Scaling factor for Coulomb radius */
818 real fourierspacing; /* Basic fourierspacing from tpr */
821 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
822 *ntprs > 1 ? "s" : "", "%"GMX_PRId64, benchsteps > 1 ? "s" : "");
823 fprintf(stdout, buf, benchsteps);
826 sprintf(buf, " (adding %s steps from checkpoint file)", "%"GMX_PRId64);
827 fprintf(stdout, buf, statesteps);
828 benchsteps += statesteps;
830 fprintf(stdout, ".\n");
834 read_tpx_state(fn_sim_tpr, ir, &state, NULL, &mtop);
836 /* Check if some kind of PME was chosen */
837 if (EEL_PME(ir->coulombtype) == FALSE)
839 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
843 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
844 if ( (ir->cutoff_scheme != ecutsVERLET) &&
845 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
847 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
848 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
850 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
851 else if (ir->rcoulomb > ir->rlist)
853 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
854 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
857 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
859 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
863 /* Reduce the number of steps for the benchmarks */
864 info->orig_sim_steps = ir->nsteps;
865 ir->nsteps = benchsteps;
866 /* We must not use init_step from the input tpr file for the benchmarks */
867 info->orig_init_step = ir->init_step;
870 /* For PME-switch potentials, keep the radial distance of the buffer region */
871 nlist_buffer = ir->rlist - ir->rcoulomb;
873 /* Determine length of triclinic box vectors */
874 for (d = 0; d < DIM; d++)
877 for (i = 0; i < DIM; i++)
879 box_size[d] += state.box[d][i]*state.box[d][i];
881 box_size[d] = sqrt(box_size[d]);
884 if (ir->fourier_spacing > 0)
886 info->fsx[0] = ir->fourier_spacing;
887 info->fsy[0] = ir->fourier_spacing;
888 info->fsz[0] = ir->fourier_spacing;
892 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
893 info->fsx[0] = box_size[XX]/ir->nkx;
894 info->fsy[0] = box_size[YY]/ir->nky;
895 info->fsz[0] = box_size[ZZ]/ir->nkz;
898 /* If no value for the fourierspacing was provided on the command line, we
899 * use the reconstruction from the tpr file */
900 if (ir->fourier_spacing > 0)
902 /* Use the spacing from the tpr */
903 fourierspacing = ir->fourier_spacing;
907 /* Use the maximum observed spacing */
908 fourierspacing = max(max(info->fsx[0], info->fsy[0]), info->fsz[0]);
911 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
913 /* For performance comparisons the number of particles is useful to have */
914 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
916 /* Print information about settings of which some are potentially modified: */
917 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
918 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
919 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
920 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
921 if (ir_vdw_switched(ir))
923 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
925 if (EPME_SWITCHED(ir->coulombtype))
927 fprintf(fp, " rlist : %f nm\n", ir->rlist);
929 if (ir->rlistlong != max_cutoff(ir->rvdw, ir->rcoulomb))
931 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
934 /* Print a descriptive line about the tpr settings tested */
935 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
936 fprintf(fp, " No. scaling rcoulomb");
937 fprintf(fp, " nkx nky nkz");
938 fprintf(fp, " spacing");
939 if (evdwCUT == ir->vdwtype)
941 fprintf(fp, " rvdw");
943 if (EPME_SWITCHED(ir->coulombtype))
945 fprintf(fp, " rlist");
947 if (ir->rlistlong != max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb)) )
949 fprintf(fp, " rlistlong");
951 fprintf(fp, " tpr file\n");
953 /* Loop to create the requested number of tpr input files */
954 for (j = 0; j < *ntprs; j++)
956 /* The first .tpr is the provided one, just need to modify nsteps,
957 * so skip the following block */
960 /* Determine which Coulomb radii rc to use in the benchmarks */
961 add = (rmax-rmin)/(*ntprs-1);
962 if (gmx_within_tol(rmin, info->rcoulomb[0], GMX_REAL_EPS))
964 ir->rcoulomb = rmin + j*add;
966 else if (gmx_within_tol(rmax, info->rcoulomb[0], GMX_REAL_EPS))
968 ir->rcoulomb = rmin + (j-1)*add;
972 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
973 add = (rmax-rmin)/(*ntprs-2);
974 ir->rcoulomb = rmin + (j-1)*add;
977 /* Determine the scaling factor fac */
978 fac = ir->rcoulomb/info->rcoulomb[0];
980 /* Scale the Fourier grid spacing */
981 ir->nkx = ir->nky = ir->nkz = 0;
982 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
984 /* Adjust other radii since various conditions need to be fulfilled */
985 if (eelPME == ir->coulombtype)
987 /* plain PME, rcoulomb must be equal to rlist */
988 ir->rlist = ir->rcoulomb;
992 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
993 ir->rlist = ir->rcoulomb + nlist_buffer;
996 if (bScaleRvdw && evdwCUT == ir->vdwtype)
998 if (ecutsVERLET == ir->cutoff_scheme)
1000 /* With Verlet, the van der Waals radius must always equal the Coulomb radius */
1001 ir->rvdw = ir->rcoulomb;
1005 /* For vdw cutoff, rvdw >= rlist */
1006 ir->rvdw = max(info->rvdw[0], ir->rlist);
1010 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1012 } /* end of "if (j != 0)" */
1014 /* for j==0: Save the original settings
1015 * for j >0: Save modified radii and Fourier grids */
1016 info->rcoulomb[j] = ir->rcoulomb;
1017 info->rvdw[j] = ir->rvdw;
1018 info->nkx[j] = ir->nkx;
1019 info->nky[j] = ir->nky;
1020 info->nkz[j] = ir->nkz;
1021 info->rlist[j] = ir->rlist;
1022 info->rlistlong[j] = ir->rlistlong;
1023 info->fsx[j] = fac*fourierspacing;
1024 info->fsy[j] = fac*fourierspacing;
1025 info->fsz[j] = fac*fourierspacing;
1027 /* Write the benchmark tpr file */
1028 strncpy(fn_bench_tprs[j], fn_sim_tpr, strlen(fn_sim_tpr)-strlen(".tpr"));
1029 sprintf(buf, "_bench%.2d.tpr", j);
1030 strcat(fn_bench_tprs[j], buf);
1031 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1032 fprintf(stdout, "%"GMX_PRId64, ir->nsteps);
1035 fprintf(stdout, ", scaling factor %f\n", fac);
1039 fprintf(stdout, ", unmodified settings\n");
1042 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1044 /* Write information about modified tpr settings to log file */
1045 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1046 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1047 fprintf(fp, " %9f ", info->fsx[j]);
1048 if (evdwCUT == ir->vdwtype)
1050 fprintf(fp, "%10f", ir->rvdw);
1052 if (EPME_SWITCHED(ir->coulombtype))
1054 fprintf(fp, "%10f", ir->rlist);
1056 if (info->rlistlong[0] != max_cutoff(info->rlist[0], max_cutoff(info->rvdw[0], info->rcoulomb[0])) )
1058 fprintf(fp, "%10f", ir->rlistlong);
1060 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1062 /* Make it clear to the user that some additional settings were modified */
1063 if (!gmx_within_tol(ir->rvdw, info->rvdw[0], GMX_REAL_EPS)
1064 || !gmx_within_tol(ir->rlistlong, info->rlistlong[0], GMX_REAL_EPS) )
1071 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1072 "other input settings were also changed (see table above).\n"
1073 "Please check if the modified settings are appropriate.\n");
1081 /* Rename the files we want to keep to some meaningful filename and
1082 * delete the rest */
1083 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1084 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1086 char numstring[STRLEN];
1087 char newfilename[STRLEN];
1088 const char *fn = NULL;
1093 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1095 for (i = 0; i < nfile; i++)
1097 opt = (char *)fnm[i].opt;
1098 if (strcmp(opt, "-p") == 0)
1100 /* do nothing; keep this file */
1103 else if (strcmp(opt, "-bg") == 0)
1105 /* Give the log file a nice name so one can later see which parameters were used */
1106 numstring[0] = '\0';
1109 sprintf(numstring, "_%d", nr);
1111 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1112 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1114 fprintf(stdout, "renaming log file to %s\n", newfilename);
1115 make_backup(newfilename);
1116 rename(opt2fn("-bg", nfile, fnm), newfilename);
1119 else if (strcmp(opt, "-err") == 0)
1121 /* This file contains the output of stderr. We want to keep it in
1122 * cases where there have been problems. */
1123 fn = opt2fn(opt, nfile, fnm);
1124 numstring[0] = '\0';
1127 sprintf(numstring, "_%d", nr);
1129 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1134 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1135 make_backup(newfilename);
1136 rename(fn, newfilename);
1140 fprintf(stdout, "Deleting %s\n", fn);
1145 /* Delete the files which are created for each benchmark run: (options -b*) */
1146 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1148 remove_if_exists(opt2fn(opt, nfile, fnm));
1155 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1158 /* Create a list of numbers of PME nodes to test */
1159 static void make_npme_list(
1160 const char *npmevalues_opt, /* Make a complete list with all
1161 * possibilities or a short list that keeps only
1162 * reasonable numbers of PME nodes */
1163 int *nentries, /* Number of entries we put in the nPMEnodes list */
1164 int *nPMEnodes[], /* Each entry contains the value for -npme */
1165 int nnodes, /* Total number of nodes to do the tests on */
1166 int minPMEnodes, /* Minimum number of PME nodes */
1167 int maxPMEnodes) /* Maximum number of PME nodes */
1170 int min_factor = 1; /* We request that npp and npme have this minimal
1171 * largest common factor (depends on npp) */
1172 int nlistmax; /* Max. list size */
1173 int nlist; /* Actual number of entries in list */
1177 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1178 if (0 == strcmp(npmevalues_opt, "all") )
1182 else if (0 == strcmp(npmevalues_opt, "subset") )
1184 eNPME = eNpmeSubset;
1186 else /* "auto" or "range" */
1192 else if (nnodes < 128)
1194 eNPME = eNpmeReduced;
1198 eNPME = eNpmeSubset;
1202 /* Calculate how many entries we could possibly have (in case of -npme all) */
1205 nlistmax = maxPMEnodes - minPMEnodes + 3;
1206 if (0 == minPMEnodes)
1216 /* Now make the actual list which is at most of size nlist */
1217 snew(*nPMEnodes, nlistmax);
1218 nlist = 0; /* start counting again, now the real entries in the list */
1219 for (i = 0; i < nlistmax - 2; i++)
1221 npme = maxPMEnodes - i;
1232 /* For 2d PME we want a common largest factor of at least the cube
1233 * root of the number of PP nodes */
1234 min_factor = (int) pow(npp, 1.0/3.0);
1237 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1240 if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1242 (*nPMEnodes)[nlist] = npme;
1246 /* We always test 0 PME nodes and the automatic number */
1247 *nentries = nlist + 2;
1248 (*nPMEnodes)[nlist ] = 0;
1249 (*nPMEnodes)[nlist+1] = -1;
1251 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1252 for (i = 0; i < *nentries-1; i++)
1254 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1256 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1260 /* Allocate memory to store the performance data */
1261 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1266 for (k = 0; k < ntprs; k++)
1268 snew(perfdata[k], datasets);
1269 for (i = 0; i < datasets; i++)
1271 for (j = 0; j < repeats; j++)
1273 snew(perfdata[k][i].Gcycles, repeats);
1274 snew(perfdata[k][i].ns_per_day, repeats);
1275 snew(perfdata[k][i].PME_f_load, repeats);
1282 /* Check for errors on mdrun -h */
1283 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp,
1284 const t_filenm *fnm, int nfile)
1286 const char *fn = NULL;
1287 char *command, *msg;
1291 snew(command, length + 15);
1292 snew(msg, length + 500);
1294 fprintf(stdout, "Making sure the benchmarks can be executed by running just 1 step...\n");
1295 sprintf(command, "%s -nsteps 1 -quiet", mdrun_cmd_line);
1296 fprintf(stdout, "Executing '%s' ...\n", command);
1297 ret = gmx_system_call(command);
1301 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1302 * get the error message from mdrun itself */
1303 sprintf(msg, "Cannot run the benchmark simulations! Please check the error message of\n"
1304 "mdrun for the source of the problem. Did you provide a command line\n"
1305 "argument that neither g_tune_pme nor mdrun understands? Offending command:\n"
1306 "\n%s\n\n", command);
1308 fprintf(stderr, "%s", msg);
1310 fprintf(fp, "%s", msg);
1314 fprintf(stdout, "Benchmarks can be executed!\n");
1316 /* Clean up the benchmark output files we just created */
1317 fprintf(stdout, "Cleaning up ...\n");
1318 remove_if_exists(opt2fn("-bc", nfile, fnm));
1319 remove_if_exists(opt2fn("-be", nfile, fnm));
1320 remove_if_exists(opt2fn("-bcpo", nfile, fnm));
1321 remove_if_exists(opt2fn("-bg", nfile, fnm));
1328 static void do_the_tests(
1329 FILE *fp, /* General g_tune_pme output file */
1330 char **tpr_names, /* Filenames of the input files to test */
1331 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1332 int minPMEnodes, /* Min fraction of nodes to use for PME */
1333 int npme_fixed, /* If >= -1, test fixed number of PME
1335 const char *npmevalues_opt, /* Which -npme values should be tested */
1336 t_perf **perfdata, /* Here the performace data is stored */
1337 int *pmeentries, /* Entries in the nPMEnodes list */
1338 int repeats, /* Repeat each test this often */
1339 int nnodes, /* Total number of nodes = nPP + nPME */
1340 int nr_tprs, /* Total number of tpr files to test */
1341 gmx_bool bThreads, /* Threads or MPI? */
1342 char *cmd_mpirun, /* mpirun command string */
1343 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1344 char *cmd_mdrun, /* mdrun command string */
1345 char *cmd_args_bench, /* arguments for mdrun in a string */
1346 const t_filenm *fnm, /* List of filenames from command line */
1347 int nfile, /* Number of files specified on the cmdl. */
1348 int presteps, /* DLB equilibration steps, is checked */
1349 gmx_int64_t cpt_steps) /* Time step counter in the checkpoint */
1351 int i, nr, k, ret, count = 0, totaltests;
1352 int *nPMEnodes = NULL;
1355 char *command, *cmd_stub;
1357 gmx_bool bResetProblem = FALSE;
1358 gmx_bool bFirst = TRUE;
1361 /* This string array corresponds to the eParselog enum type at the start
1363 const char* ParseLog[] = {
1365 "Logfile not found!",
1366 "No timings, logfile truncated?",
1367 "Run was terminated.",
1368 "Counters were not reset properly.",
1369 "No DD grid found for these settings.",
1370 "TPX version conflict!",
1371 "mdrun was not started in parallel!",
1372 "Number of PP nodes has a prime factor that is too large.",
1375 char str_PME_f_load[13];
1378 /* Allocate space for the mdrun command line. 100 extra characters should
1379 be more than enough for the -npme etcetera arguments */
1380 cmdline_length = strlen(cmd_mpirun)
1383 + strlen(cmd_args_bench)
1384 + strlen(tpr_names[0]) + 100;
1385 snew(command, cmdline_length);
1386 snew(cmd_stub, cmdline_length);
1388 /* Construct the part of the command line that stays the same for all tests: */
1391 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1395 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1398 /* Create a list of numbers of PME nodes to test */
1399 if (npme_fixed < -1)
1401 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1402 nnodes, minPMEnodes, maxPMEnodes);
1408 nPMEnodes[0] = npme_fixed;
1409 fprintf(stderr, "Will use a fixed number of %d PME-only nodes.\n", nPMEnodes[0]);
1414 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1416 finalize(opt2fn("-p", nfile, fnm));
1420 /* Allocate one dataset for each tpr input file: */
1421 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1423 /*****************************************/
1424 /* Main loop over all tpr files to test: */
1425 /*****************************************/
1426 totaltests = nr_tprs*(*pmeentries)*repeats;
1427 for (k = 0; k < nr_tprs; k++)
1429 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1430 fprintf(fp, "PME nodes Gcycles ns/day PME/f Remark\n");
1431 /* Loop over various numbers of PME nodes: */
1432 for (i = 0; i < *pmeentries; i++)
1434 pd = &perfdata[k][i];
1436 /* Loop over the repeats for each scenario: */
1437 for (nr = 0; nr < repeats; nr++)
1439 pd->nPMEnodes = nPMEnodes[i];
1441 /* Add -npme and -s to the command line and save it. Note that
1442 * the -passall (if set) options requires cmd_args_bench to be
1443 * at the end of the command line string */
1444 snew(pd->mdrun_cmd_line, cmdline_length);
1445 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s",
1446 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench);
1448 /* To prevent that all benchmarks fail due to a show-stopper argument
1449 * on the mdrun command line, we make a quick check first */
1452 make_sure_it_runs(pd->mdrun_cmd_line, cmdline_length, fp, fnm, nfile);
1456 /* Do a benchmark simulation: */
1459 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1465 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1466 (100.0*count)/totaltests,
1467 k+1, nr_tprs, i+1, *pmeentries, buf);
1468 make_backup(opt2fn("-err", nfile, fnm));
1469 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1470 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1471 gmx_system_call(command);
1473 /* Collect the performance data from the log file; also check stderr
1474 * for fatal errors */
1475 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1476 pd, nr, presteps, cpt_steps, nnodes);
1477 if ((presteps > 0) && (ret == eParselogResetProblem))
1479 bResetProblem = TRUE;
1482 if (-1 == pd->nPMEnodes)
1484 sprintf(buf, "(%3d)", pd->guessPME);
1492 if (pd->PME_f_load[nr] > 0.0)
1494 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1498 sprintf(str_PME_f_load, "%s", " - ");
1501 /* Write the data we got to disk */
1502 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1503 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1504 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1506 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1512 /* Do some cleaning up and delete the files we do not need any more */
1513 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1515 /* If the first run with this number of processors already failed, do not try again: */
1516 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1518 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1519 count += repeats-(nr+1);
1522 } /* end of repeats loop */
1523 } /* end of -npme loop */
1524 } /* end of tpr file loop */
1529 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1537 static void check_input(
1544 real maxPMEfraction,
1545 real minPMEfraction,
1547 gmx_int64_t bench_nsteps,
1548 const t_filenm *fnm,
1558 /* Make sure the input file exists */
1559 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1561 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1564 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1565 if ( (0 == strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1567 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1568 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1571 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1574 gmx_fatal(FARGS, "Number of repeats < 0!");
1577 /* Check number of nodes */
1580 gmx_fatal(FARGS, "Number of nodes/threads must be a positive integer.");
1583 /* Automatically choose -ntpr if not set */
1593 /* Set a reasonable scaling factor for rcoulomb */
1596 *rmax = rcoulomb * 1.2;
1599 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1605 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1609 /* Make shure that rmin <= rcoulomb <= rmax */
1618 if (!(*rmin <= *rmax) )
1620 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1621 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1623 /* Add test scenarios if rmin or rmax were set */
1626 if (!gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1629 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1632 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1635 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1640 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1641 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) || !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1643 *ntprs = max(*ntprs, 2);
1646 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1647 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1649 *ntprs = max(*ntprs, 3);
1654 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1659 if (gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && gmx_within_tol(rcoulomb, *rmax, GMX_REAL_EPS)) /* We have just a single rc */
1661 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1662 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1663 "with correspondingly adjusted PME grid settings\n");
1668 /* Check whether max and min fraction are within required values */
1669 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1671 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1673 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1675 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1677 if (maxPMEfraction < minPMEfraction)
1679 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1682 /* Check whether the number of steps - if it was set - has a reasonable value */
1683 if (bench_nsteps < 0)
1685 gmx_fatal(FARGS, "Number of steps must be positive.");
1688 if (bench_nsteps > 10000 || bench_nsteps < 100)
1690 fprintf(stderr, "WARNING: steps=");
1691 fprintf(stderr, "%"GMX_PRId64, bench_nsteps);
1692 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1697 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1700 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1703 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1705 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1709 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1710 * only. We need to check whether the requested number of PME-only nodes
1712 if (npme_fixed > -1)
1714 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1715 if (2*npme_fixed > nnodes)
1717 gmx_fatal(FARGS, "Cannot have more than %d PME-only nodes for a total of %d nodes (you chose %d).\n",
1718 nnodes/2, nnodes, npme_fixed);
1720 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1722 fprintf(stderr, "WARNING: Only %g percent of the nodes are assigned as PME-only nodes.\n",
1723 100.0*((real)npme_fixed / (real)nnodes));
1725 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1727 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1728 " fixed number of PME-only nodes is requested with -fix.\n");
1734 /* Returns TRUE when "opt" is needed at launch time */
1735 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1737 if (0 == strncmp(opt, "-swap", 5))
1742 /* Apart from the input .tpr and the output log files we need all options that
1743 * were set on the command line and that do not start with -b */
1744 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2)
1745 || 0 == strncmp(opt, "-err", 4) || 0 == strncmp(opt, "-p", 2) )
1754 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1755 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1757 /* Apart from the input .tpr, all files starting with "-b" are for
1758 * _b_enchmark files exclusively */
1759 if (0 == strncmp(opt, "-s", 2))
1764 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2))
1766 if (!bOptional || bSet)
1783 if (bSet) /* These are additional input files like -cpi -ei */
1796 /* Adds 'buf' to 'str' */
1797 static void add_to_string(char **str, char *buf)
1802 len = strlen(*str) + strlen(buf) + 1;
1808 /* Create the command line for the benchmark as well as for the real run */
1809 static void create_command_line_snippets(
1810 gmx_bool bAppendFiles,
1811 gmx_bool bKeepAndNumCPT,
1812 gmx_bool bResetHWay,
1816 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1817 char *cmd_args_launch[], /* command line arguments for simulation run */
1818 char extra_args[]) /* Add this to the end of the command line */
1823 char strbuf[STRLEN];
1826 /* strlen needs at least '\0' as a string: */
1827 snew(*cmd_args_bench, 1);
1828 snew(*cmd_args_launch, 1);
1829 *cmd_args_launch[0] = '\0';
1830 *cmd_args_bench[0] = '\0';
1833 /*******************************************/
1834 /* 1. Process other command line arguments */
1835 /*******************************************/
1838 /* Add equilibration steps to benchmark options */
1839 sprintf(strbuf, "-resetstep %d ", presteps);
1840 add_to_string(cmd_args_bench, strbuf);
1842 /* These switches take effect only at launch time */
1843 if (FALSE == bAppendFiles)
1845 add_to_string(cmd_args_launch, "-noappend ");
1849 add_to_string(cmd_args_launch, "-cpnum ");
1853 add_to_string(cmd_args_launch, "-resethway ");
1856 /********************/
1857 /* 2. Process files */
1858 /********************/
1859 for (i = 0; i < nfile; i++)
1861 opt = (char *)fnm[i].opt;
1862 name = opt2fn(opt, nfile, fnm);
1864 /* Strbuf contains the options, now let's sort out where we need that */
1865 sprintf(strbuf, "%s %s ", opt, name);
1867 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1869 /* All options starting with -b* need the 'b' removed,
1870 * therefore overwrite strbuf */
1871 if (0 == strncmp(opt, "-b", 2))
1873 sprintf(strbuf, "-%s %s ", &opt[2], name);
1876 add_to_string(cmd_args_bench, strbuf);
1879 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
1881 add_to_string(cmd_args_launch, strbuf);
1885 add_to_string(cmd_args_bench, extra_args);
1886 add_to_string(cmd_args_launch, extra_args);
1890 /* Set option opt */
1891 static void setopt(const char *opt, int nfile, t_filenm fnm[])
1895 for (i = 0; (i < nfile); i++)
1897 if (strcmp(opt, fnm[i].opt) == 0)
1899 fnm[i].flag |= ffSET;
1905 /* This routine inspects the tpr file and ...
1906 * 1. checks for output files that get triggered by a tpr option. These output
1907 * files are marked as 'set' to allow for a proper cleanup after each
1909 * 2. returns the PME:PP load ratio
1910 * 3. returns rcoulomb from the tpr */
1911 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1913 gmx_bool bPull; /* Is pulling requested in .tpr file? */
1914 gmx_bool bTpi; /* Is test particle insertion requested? */
1915 gmx_bool bFree; /* Is a free energy simulation requested? */
1916 gmx_bool bNM; /* Is a normal mode analysis requested? */
1917 gmx_bool bSwap; /* Is water/ion position swapping requested? */
1923 /* Check tpr file for options that trigger extra output files */
1924 read_tpx_state(opt2fn("-s", nfile, fnm), &ir, &state, NULL, &mtop);
1925 bPull = (epullNO != ir.ePull);
1926 bFree = (efepNO != ir.efep );
1927 bNM = (eiNM == ir.eI );
1928 bSwap = (eswapNO != ir.eSwapCoords);
1929 bTpi = EI_TPI(ir.eI);
1931 /* Set these output files on the tuning command-line */
1934 setopt("-pf", nfile, fnm);
1935 setopt("-px", nfile, fnm);
1939 setopt("-dhdl", nfile, fnm);
1943 setopt("-tpi", nfile, fnm);
1944 setopt("-tpid", nfile, fnm);
1948 setopt("-mtx", nfile, fnm);
1952 setopt("-swap", nfile, fnm);
1955 *rcoulomb = ir.rcoulomb;
1957 /* Return the estimate for the number of PME nodes */
1958 return pme_load_estimate(&mtop, &ir, state.box);
1962 static void couple_files_options(int nfile, t_filenm fnm[])
1965 gmx_bool bSet, bBench;
1970 for (i = 0; i < nfile; i++)
1972 opt = (char *)fnm[i].opt;
1973 bSet = ((fnm[i].flag & ffSET) != 0);
1974 bBench = (0 == strncmp(opt, "-b", 2));
1976 /* Check optional files */
1977 /* If e.g. -eo is set, then -beo also needs to be set */
1978 if (is_optional(&fnm[i]) && bSet && !bBench)
1980 sprintf(buf, "-b%s", &opt[1]);
1981 setopt(buf, nfile, fnm);
1983 /* If -beo is set, then -eo also needs to be! */
1984 if (is_optional(&fnm[i]) && bSet && bBench)
1986 sprintf(buf, "-%s", &opt[2]);
1987 setopt(buf, nfile, fnm);
1993 #define BENCHSTEPS (1000)
1995 int gmx_tune_pme(int argc, char *argv[])
1997 const char *desc[] = {
1998 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of processors/threads, [THISMODULE] systematically",
1999 "times [gmx-mdrun] with various numbers of PME-only nodes and determines",
2000 "which setting is fastest. It will also test whether performance can",
2001 "be enhanced by shifting load from the reciprocal to the real space",
2002 "part of the Ewald sum. ",
2003 "Simply pass your [TT].tpr[tt] file to [THISMODULE] together with other options",
2004 "for [gmx-mdrun] as needed.[PAR]",
2005 "Which executables are used can be set in the environment variables",
2006 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
2007 "will be used as defaults. Note that for certain MPI frameworks you",
2008 "need to provide a machine- or hostfile. This can also be passed",
2009 "via the MPIRUN variable, e.g.[PAR]",
2010 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
2011 "Please call [THISMODULE] with the normal options you would pass to",
2012 "[gmx-mdrun] and add [TT]-np[tt] for the number of processors to perform the",
2013 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2014 "to repeat each test several times to get better statistics. [PAR]",
2015 "[THISMODULE] can test various real space / reciprocal space workloads",
2016 "for you. With [TT]-ntpr[tt] you control how many extra [TT].tpr[tt] files will be",
2017 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2018 "Typically, the first test (number 0) will be with the settings from the input",
2019 "[TT].tpr[tt] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2020 "specified by [TT]-rmax[tt] with a somwhat smaller PME grid at the same time. ",
2021 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2022 "The remaining [TT].tpr[tt] files will have equally-spaced Coulomb radii (and Fourier "
2023 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2024 "if you just seek the optimal number of PME-only nodes; in that case",
2025 "your input [TT].tpr[tt] file will remain unchanged.[PAR]",
2026 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2027 "MD systems. The dynamic load balancing needs about 100 time steps",
2028 "to adapt to local load imbalances, therefore the time step counters",
2029 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2030 "for a higher accuarcy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2031 "From the 'DD' load imbalance entries in the md.log output file you",
2032 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2033 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2034 "After calling [gmx-mdrun] several times, detailed performance information",
2035 "is available in the output file [TT]perf.out[tt].",
2036 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2037 "(options [TT]-b[tt]*), these will be automatically deleted after each test.[PAR]",
2038 "If you want the simulation to be started automatically with the",
2039 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2044 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2045 real maxPMEfraction = 0.50;
2046 real minPMEfraction = 0.25;
2047 int maxPMEnodes, minPMEnodes;
2048 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2049 float guessPMEnodes;
2050 int npme_fixed = -2; /* If >= -1, use only this number
2051 * of PME-only nodes */
2053 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2054 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2055 gmx_bool bScaleRvdw = TRUE;
2056 gmx_int64_t bench_nsteps = BENCHSTEPS;
2057 gmx_int64_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2058 gmx_int64_t cpt_steps = 0; /* Step counter in .cpt input file */
2059 int presteps = 100; /* Do a full cycle reset after presteps steps */
2060 gmx_bool bOverwrite = FALSE, bKeepTPR;
2061 gmx_bool bLaunch = FALSE;
2062 char *ExtraArgs = NULL;
2063 char **tpr_names = NULL;
2064 const char *simulation_tpr = NULL;
2065 int best_npme, best_tpr;
2066 int sim_part = 1; /* For benchmarks with checkpoint files */
2069 /* Default program names if nothing else is found */
2070 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2071 char *cmd_args_bench, *cmd_args_launch;
2072 char *cmd_np = NULL;
2074 t_perf **perfdata = NULL;
2080 /* Print out how long the tuning took */
2083 static t_filenm fnm[] = {
2085 { efOUT, "-p", "perf", ffWRITE },
2086 { efLOG, "-err", "bencherr", ffWRITE },
2087 { efTPX, "-so", "tuned", ffWRITE },
2089 { efTPX, NULL, NULL, ffREAD },
2090 { efTRN, "-o", NULL, ffWRITE },
2091 { efCOMPRESSED, "-x", NULL, ffOPTWR },
2092 { efCPT, "-cpi", NULL, ffOPTRD },
2093 { efCPT, "-cpo", NULL, ffOPTWR },
2094 { efSTO, "-c", "confout", ffWRITE },
2095 { efEDR, "-e", "ener", ffWRITE },
2096 { efLOG, "-g", "md", ffWRITE },
2097 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2098 { efXVG, "-field", "field", ffOPTWR },
2099 { efXVG, "-table", "table", ffOPTRD },
2100 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2101 { efXVG, "-tablep", "tablep", ffOPTRD },
2102 { efXVG, "-tableb", "table", ffOPTRD },
2103 { efTRX, "-rerun", "rerun", ffOPTRD },
2104 { efXVG, "-tpi", "tpi", ffOPTWR },
2105 { efXVG, "-tpid", "tpidist", ffOPTWR },
2106 { efEDI, "-ei", "sam", ffOPTRD },
2107 { efXVG, "-eo", "edsam", ffOPTWR },
2108 { efXVG, "-devout", "deviatie", ffOPTWR },
2109 { efXVG, "-runav", "runaver", ffOPTWR },
2110 { efXVG, "-px", "pullx", ffOPTWR },
2111 { efXVG, "-pf", "pullf", ffOPTWR },
2112 { efXVG, "-ro", "rotation", ffOPTWR },
2113 { efLOG, "-ra", "rotangles", ffOPTWR },
2114 { efLOG, "-rs", "rotslabs", ffOPTWR },
2115 { efLOG, "-rt", "rottorque", ffOPTWR },
2116 { efMTX, "-mtx", "nm", ffOPTWR },
2117 { efNDX, "-dn", "dipole", ffOPTWR },
2118 { efXVG, "-swap", "swapions", ffOPTWR },
2119 /* Output files that are deleted after each benchmark run */
2120 { efTRN, "-bo", "bench", ffWRITE },
2121 { efXTC, "-bx", "bench", ffWRITE },
2122 { efCPT, "-bcpo", "bench", ffWRITE },
2123 { efSTO, "-bc", "bench", ffWRITE },
2124 { efEDR, "-be", "bench", ffWRITE },
2125 { efLOG, "-bg", "bench", ffWRITE },
2126 { efXVG, "-beo", "benchedo", ffOPTWR },
2127 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2128 { efXVG, "-bfield", "benchfld", ffOPTWR },
2129 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2130 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2131 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2132 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2133 { efXVG, "-bpx", "benchpx", ffOPTWR },
2134 { efXVG, "-bpf", "benchpf", ffOPTWR },
2135 { efXVG, "-bro", "benchrot", ffOPTWR },
2136 { efLOG, "-bra", "benchrota", ffOPTWR },
2137 { efLOG, "-brs", "benchrots", ffOPTWR },
2138 { efLOG, "-brt", "benchrott", ffOPTWR },
2139 { efMTX, "-bmtx", "benchn", ffOPTWR },
2140 { efNDX, "-bdn", "bench", ffOPTWR },
2141 { efXVG, "-bswap", "benchswp", ffOPTWR }
2144 gmx_bool bThreads = FALSE;
2148 const char *procstring[] =
2149 { NULL, "-np", "-n", "none", NULL };
2150 const char *npmevalues_opt[] =
2151 { NULL, "auto", "all", "subset", NULL };
2153 gmx_bool bAppendFiles = TRUE;
2154 gmx_bool bKeepAndNumCPT = FALSE;
2155 gmx_bool bResetCountersHalfWay = FALSE;
2156 gmx_bool bBenchmark = TRUE;
2158 output_env_t oenv = NULL;
2161 /***********************/
2162 /* g_tune_pme options: */
2163 /***********************/
2164 { "-np", FALSE, etINT, {&nnodes},
2165 "Number of nodes to run the tests on (must be > 2 for separate PME nodes)" },
2166 { "-npstring", FALSE, etENUM, {procstring},
2167 "Specify the number of processors to [TT]$MPIRUN[tt] using this string"},
2168 { "-ntmpi", FALSE, etINT, {&nthreads},
2169 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2170 { "-r", FALSE, etINT, {&repeats},
2171 "Repeat each test this often" },
2172 { "-max", FALSE, etREAL, {&maxPMEfraction},
2173 "Max fraction of PME nodes to test with" },
2174 { "-min", FALSE, etREAL, {&minPMEfraction},
2175 "Min fraction of PME nodes to test with" },
2176 { "-npme", FALSE, etENUM, {npmevalues_opt},
2177 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2178 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2179 { "-fix", FALSE, etINT, {&npme_fixed},
2180 "If >= -1, do not vary the number of PME-only nodes, instead use this fixed value and only vary rcoulomb and the PME grid spacing."},
2181 { "-rmax", FALSE, etREAL, {&rmax},
2182 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2183 { "-rmin", FALSE, etREAL, {&rmin},
2184 "If >0, minimal rcoulomb for -ntpr>1" },
2185 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2186 "Scale rvdw along with rcoulomb"},
2187 { "-ntpr", FALSE, etINT, {&ntprs},
2188 "Number of [TT].tpr[tt] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2189 "If < 1, automatically choose the number of [TT].tpr[tt] files to test" },
2190 { "-steps", FALSE, etINT64, {&bench_nsteps},
2191 "Take timings for this many steps in the benchmark runs" },
2192 { "-resetstep", FALSE, etINT, {&presteps},
2193 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2194 { "-simsteps", FALSE, etINT64, {&new_sim_nsteps},
2195 "If non-negative, perform this many steps in the real run (overwrites nsteps from [TT].tpr[tt], add [TT].cpt[tt] steps)" },
2196 { "-launch", FALSE, etBOOL, {&bLaunch},
2197 "Launch the real simulation after optimization" },
2198 { "-bench", FALSE, etBOOL, {&bBenchmark},
2199 "Run the benchmarks or just create the input [TT].tpr[tt] files?" },
2200 /******************/
2201 /* mdrun options: */
2202 /******************/
2203 /* We let g_tune_pme parse and understand these options, because we need to
2204 * prevent that they appear on the mdrun command line for the benchmarks */
2205 { "-append", FALSE, etBOOL, {&bAppendFiles},
2206 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2207 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2208 "Keep and number checkpoint files (launch only)" },
2209 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2210 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2213 #define NFILE asize(fnm)
2215 seconds = gmx_gettime();
2217 if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2218 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2224 /* Store the remaining unparsed command line entries in a string which
2225 * is then attached to the mdrun command line */
2227 ExtraArgs[0] = '\0';
2228 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2230 add_to_string(&ExtraArgs, argv[i]);
2231 add_to_string(&ExtraArgs, " ");
2234 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2237 if (opt2parg_bSet("-npstring", asize(pa), pa))
2239 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2244 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2246 /* and now we just set this; a bit of an ugly hack*/
2249 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2250 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2252 /* Automatically set -beo options if -eo is set etc. */
2253 couple_files_options(NFILE, fnm);
2255 /* Construct the command line arguments for benchmark runs
2256 * as well as for the simulation run */
2259 sprintf(bbuf, " -ntmpi %d ", nthreads);
2263 /* This string will be used for MPI runs and will appear after the
2264 * mpirun command. */
2265 if (strcmp(procstring[0], "none") != 0)
2267 sprintf(bbuf, " %s %d ", procstring[0], nnodes);
2277 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2278 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs);
2280 /* Read in checkpoint file if requested */
2282 if (opt2bSet("-cpi", NFILE, fnm))
2285 cr->duty = DUTY_PP; /* makes the following routine happy */
2286 read_checkpoint_simulation_part(opt2fn("-cpi", NFILE, fnm),
2287 &sim_part, &cpt_steps, cr,
2288 FALSE, NFILE, fnm, NULL, NULL);
2291 /* sim_part will now be 1 if no checkpoint file was found */
2294 gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi", NFILE, fnm));
2298 /* Open performance output file and write header info */
2299 fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2301 /* Make a quick consistency check of command line parameters */
2302 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2303 maxPMEfraction, minPMEfraction, npme_fixed,
2304 bench_nsteps, fnm, NFILE, sim_part, presteps,
2307 /* Determine the maximum and minimum number of PME nodes to test,
2308 * the actual list of settings is build in do_the_tests(). */
2309 if ((nnodes > 2) && (npme_fixed < -1))
2311 if (0 == strcmp(npmevalues_opt[0], "auto"))
2313 /* Determine the npme range automatically based on the PME:PP load guess */
2314 if (guessPMEratio > 1.0)
2316 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2317 maxPMEnodes = nnodes/2;
2318 minPMEnodes = nnodes/2;
2322 /* PME : PP load is in the range 0..1, let's test around the guess */
2323 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2324 minPMEnodes = floor(0.7*guessPMEnodes);
2325 maxPMEnodes = ceil(1.6*guessPMEnodes);
2326 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2331 /* Determine the npme range based on user input */
2332 maxPMEnodes = floor(maxPMEfraction*nnodes);
2333 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2334 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2335 if (maxPMEnodes != minPMEnodes)
2337 fprintf(stdout, "- %d ", maxPMEnodes);
2339 fprintf(stdout, "PME-only nodes.\n Note that the automatic number of PME-only nodes and no separate PME nodes are always tested.\n");
2348 /* Get the commands we need to set up the runs from environment variables */
2349 get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2350 if (bBenchmark && repeats > 0)
2352 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun);
2355 /* Print some header info to file */
2357 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2359 fprintf(fp, "%s for Gromacs %s\n", ShortProgram(), GromacsVersion());
2362 fprintf(fp, "Number of nodes : %d\n", nnodes);
2363 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2364 if (strcmp(procstring[0], "none") != 0)
2366 fprintf(fp, "Passing # of nodes via : %s\n", procstring[0]);
2370 fprintf(fp, "Not setting number of nodes in system call\n");
2375 fprintf(fp, "Number of threads : %d\n", nnodes);
2378 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2379 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2380 fprintf(fp, "Benchmark steps : ");
2381 fprintf(fp, "%"GMX_PRId64, bench_nsteps);
2383 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2386 fprintf(fp, "Checkpoint time step : ");
2387 fprintf(fp, "%"GMX_PRId64, cpt_steps);
2390 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2392 if (new_sim_nsteps >= 0)
2395 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2396 fprintf(stderr, "%"GMX_PRId64, new_sim_nsteps+cpt_steps);
2397 fprintf(stderr, " steps.\n");
2398 fprintf(fp, "Simulation steps : ");
2399 fprintf(fp, "%"GMX_PRId64, new_sim_nsteps);
2404 fprintf(fp, "Repeats for each test : %d\n", repeats);
2407 if (npme_fixed >= -1)
2409 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2412 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2413 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2415 /* Allocate memory for the inputinfo struct: */
2417 info->nr_inputfiles = ntprs;
2418 for (i = 0; i < ntprs; i++)
2420 snew(info->rcoulomb, ntprs);
2421 snew(info->rvdw, ntprs);
2422 snew(info->rlist, ntprs);
2423 snew(info->rlistlong, ntprs);
2424 snew(info->nkx, ntprs);
2425 snew(info->nky, ntprs);
2426 snew(info->nkz, ntprs);
2427 snew(info->fsx, ntprs);
2428 snew(info->fsy, ntprs);
2429 snew(info->fsz, ntprs);
2431 /* Make alternative tpr files to test: */
2432 snew(tpr_names, ntprs);
2433 for (i = 0; i < ntprs; i++)
2435 snew(tpr_names[i], STRLEN);
2438 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2439 * different grids could be found. */
2440 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2441 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2443 /********************************************************************************/
2444 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2445 /********************************************************************************/
2446 snew(perfdata, ntprs);
2449 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2450 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2451 cmd_args_bench, fnm, NFILE, presteps, cpt_steps);
2453 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2455 /* Analyse the results and give a suggestion for optimal settings: */
2456 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2457 repeats, info, &best_tpr, &best_npme);
2459 /* Take the best-performing tpr file and enlarge nsteps to original value */
2460 if (bKeepTPR && !bOverwrite)
2462 simulation_tpr = opt2fn("-s", NFILE, fnm);
2466 simulation_tpr = opt2fn("-so", NFILE, fnm);
2467 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2468 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2471 /* Let's get rid of the temporary benchmark input files */
2472 for (i = 0; i < ntprs; i++)
2474 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2475 remove(tpr_names[i]);
2478 /* Now start the real simulation if the user requested it ... */
2479 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2480 cmd_args_launch, simulation_tpr, best_npme);
2484 /* ... or simply print the performance results to screen: */
2487 finalize(opt2fn("-p", NFILE, fnm));