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
47 #include "gromacs/commandline/pargs.h"
49 #include "types/commrec.h"
50 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/utility/cstringutil.h"
57 #include "checkpoint.h"
63 #include "gromacs/timing/walltime_accounting.h"
64 #include "gromacs/math/utilities.h"
67 /* Enum for situations that can occur during log file parsing, the
68 * corresponding string entries can be found in do_the_tests() in
69 * const char* ParseLog[] */
75 eParselogResetProblem,
79 eParselogLargePrimeFactor,
87 int nPMEnodes; /* number of PME-only nodes used in this test */
88 int nx, ny, nz; /* DD grid */
89 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
90 double *Gcycles; /* This can contain more than one value if doing multiple tests */
94 float *PME_f_load; /* PME mesh/force load average*/
95 float PME_f_load_Av; /* Average average ;) ... */
96 char *mdrun_cmd_line; /* Mdrun command line used for this test */
102 int nr_inputfiles; /* The number of tpr and mdp input files */
103 gmx_int64_t orig_sim_steps; /* Number of steps to be done in the real simulation */
104 gmx_int64_t orig_init_step; /* Init step for the real simulation */
105 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
106 real *rvdw; /* The vdW radii */
107 real *rlist; /* Neighbourlist cutoff radius */
109 int *nkx, *nky, *nkz;
110 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
114 static void sep_line(FILE *fp)
116 fprintf(fp, "\n------------------------------------------------------------\n");
120 /* Wrapper for system calls */
121 static int gmx_system_call(char *command)
124 gmx_fatal(FARGS, "No calls to system(3) supported on this platform. Attempted to call:\n'%s'\n", command);
126 return ( system(command) );
131 /* Check if string starts with substring */
132 static gmx_bool str_starts(const char *string, const char *substring)
134 return ( strncmp(string, substring, strlen(substring)) == 0);
138 static void cleandata(t_perf *perfdata, int test_nr)
140 perfdata->Gcycles[test_nr] = 0.0;
141 perfdata->ns_per_day[test_nr] = 0.0;
142 perfdata->PME_f_load[test_nr] = 0.0;
148 static gmx_bool is_equal(real a, real b)
150 real diff, eps = 1.0e-7;
171 static void remove_if_exists(const char *fn)
175 fprintf(stdout, "Deleting %s\n", fn);
181 static void finalize(const char *fn_out)
187 fp = fopen(fn_out, "r");
188 fprintf(stdout, "\n\n");
190 while (fgets(buf, STRLEN-1, fp) != NULL)
192 fprintf(stdout, "%s", buf);
195 fprintf(stdout, "\n\n");
200 eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr
203 static int parse_logfile(const char *logfile, const char *errfile,
204 t_perf *perfdata, int test_nr, int presteps, gmx_int64_t cpt_steps,
208 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
209 const char matchstrdd[] = "Domain decomposition grid";
210 const char matchstrcr[] = "resetting all time and cycle counters";
211 const char matchstrbal[] = "Average PME mesh/force load:";
212 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";
213 const char errSIG[] = "signal, stopping at the next";
215 float dum1, dum2, dum3, dum4;
218 gmx_int64_t resetsteps = -1;
219 gmx_bool bFoundResetStr = FALSE;
220 gmx_bool bResetChecked = FALSE;
223 if (!gmx_fexist(logfile))
225 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
226 cleandata(perfdata, test_nr);
227 return eParselogNotFound;
230 fp = fopen(logfile, "r");
231 perfdata->PME_f_load[test_nr] = -1.0;
232 perfdata->guessPME = -1;
234 iFound = eFoundNothing;
237 iFound = eFoundDDStr; /* Skip some case statements */
240 while (fgets(line, STRLEN, fp) != NULL)
242 /* Remove leading spaces */
245 /* Check for TERM and INT signals from user: */
246 if (strstr(line, errSIG) != NULL)
249 cleandata(perfdata, test_nr);
250 return eParselogTerm;
253 /* Check whether cycle resetting worked */
254 if (presteps > 0 && !bFoundResetStr)
256 if (strstr(line, matchstrcr) != NULL)
258 sprintf(dumstring, "step %s", "%"GMX_SCNd64);
259 sscanf(line, dumstring, &resetsteps);
260 bFoundResetStr = TRUE;
261 if (resetsteps == presteps+cpt_steps)
263 bResetChecked = TRUE;
267 sprintf(dumstring, "%"GMX_PRId64, resetsteps);
268 sprintf(dumstring2, "%"GMX_PRId64, presteps+cpt_steps);
269 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
270 " though they were supposed to be reset at step %s!\n",
271 dumstring, dumstring2);
276 /* Look for strings that appear in a certain order in the log file: */
280 /* Look for domain decomp grid and separate PME nodes: */
281 if (str_starts(line, matchstrdd))
283 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME nodes %d",
284 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
285 if (perfdata->nPMEnodes == -1)
287 perfdata->guessPME = npme;
289 else if (perfdata->nPMEnodes != npme)
291 gmx_fatal(FARGS, "PME nodes from command line and output file are not identical");
293 iFound = eFoundDDStr;
295 /* Catch a few errors that might have occured: */
296 else if (str_starts(line, "There is no domain decomposition for"))
299 return eParselogNoDDGrid;
301 else if (str_starts(line, "The number of nodes you selected"))
304 return eParselogLargePrimeFactor;
306 else if (str_starts(line, "reading tpx file"))
309 return eParselogTPXVersion;
311 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
314 return eParselogNotParallel;
318 /* Look for PME mesh/force balance (not necessarily present, though) */
319 if (str_starts(line, matchstrbal))
321 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
323 /* Look for matchstring */
324 if (str_starts(line, matchstring))
326 iFound = eFoundAccountingStr;
329 case eFoundAccountingStr:
330 /* Already found matchstring - look for cycle data */
331 if (str_starts(line, "Total "))
333 sscanf(line, "Total %lf", &(perfdata->Gcycles[test_nr]));
334 iFound = eFoundCycleStr;
338 /* Already found cycle data - look for remaining performance info and return */
339 if (str_starts(line, "Performance:"))
341 ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
342 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
343 perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
345 if (bResetChecked || presteps == 0)
351 return eParselogResetProblem;
358 /* Close the log file */
361 /* Check why there is no performance data in the log file.
362 * Did a fatal errors occur? */
363 if (gmx_fexist(errfile))
365 fp = fopen(errfile, "r");
366 while (fgets(line, STRLEN, fp) != NULL)
368 if (str_starts(line, "Fatal error:") )
370 if (fgets(line, STRLEN, fp) != NULL)
372 fprintf(stderr, "\nWARNING: An error occured during this benchmark:\n"
376 cleandata(perfdata, test_nr);
377 return eParselogFatal;
384 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
387 /* Giving up ... we could not find out why there is no performance data in
389 fprintf(stdout, "No performance data in log file.\n");
390 cleandata(perfdata, test_nr);
392 return eParselogNoPerfData;
396 static gmx_bool analyze_data(
405 int *index_tpr, /* OUT: Nr of mdp file with best settings */
406 int *npme_optimal) /* OUT: Optimal number of PME nodes */
409 int line = 0, line_win = -1;
410 int k_win = -1, i_win = -1, winPME;
411 double s = 0.0; /* standard deviation */
414 char str_PME_f_load[13];
415 gmx_bool bCanUseOrigTPR;
416 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
422 fprintf(fp, "Summary of successful runs:\n");
423 fprintf(fp, "Line tpr PME nodes Gcycles Av. Std.dev. ns/day PME/f");
426 fprintf(fp, " DD grid");
432 for (k = 0; k < ntprs; k++)
434 for (i = 0; i < ntests; i++)
436 /* Select the right dataset: */
437 pd = &(perfdata[k][i]);
439 pd->Gcycles_Av = 0.0;
440 pd->PME_f_load_Av = 0.0;
441 pd->ns_per_day_Av = 0.0;
443 if (pd->nPMEnodes == -1)
445 sprintf(strbuf, "(%3d)", pd->guessPME);
449 sprintf(strbuf, " ");
452 /* Get the average run time of a setting */
453 for (j = 0; j < nrepeats; j++)
455 pd->Gcycles_Av += pd->Gcycles[j];
456 pd->PME_f_load_Av += pd->PME_f_load[j];
458 pd->Gcycles_Av /= nrepeats;
459 pd->PME_f_load_Av /= nrepeats;
461 for (j = 0; j < nrepeats; j++)
463 if (pd->ns_per_day[j] > 0.0)
465 pd->ns_per_day_Av += pd->ns_per_day[j];
469 /* Somehow the performance number was not aquired for this run,
470 * therefor set the average to some negative value: */
471 pd->ns_per_day_Av = -1.0f*nrepeats;
475 pd->ns_per_day_Av /= nrepeats;
478 if (pd->PME_f_load_Av > 0.0)
480 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
484 sprintf(str_PME_f_load, "%s", " - ");
488 /* We assume we had a successful run if both averages are positive */
489 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
491 /* Output statistics if repeats were done */
494 /* Calculate the standard deviation */
496 for (j = 0; j < nrepeats; j++)
498 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
503 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
504 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
505 pd->ns_per_day_Av, str_PME_f_load);
508 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
512 /* Store the index of the best run found so far in 'winner': */
513 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
526 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
531 winPME = perfdata[k_win][i_win].nPMEnodes;
535 /* We stuck to a fixed number of PME-only nodes */
536 sprintf(strbuf, "settings No. %d", k_win);
540 /* We have optimized the number of PME-only nodes */
543 sprintf(strbuf, "%s", "the automatic number of PME nodes");
547 sprintf(strbuf, "%d PME nodes", winPME);
550 fprintf(fp, "Best performance was achieved with %s", strbuf);
551 if ((nrepeats > 1) && (ntests > 1))
553 fprintf(fp, " (see line %d)", line_win);
557 /* Only mention settings if they were modified: */
558 bRefinedCoul = !is_equal(info->rcoulomb[k_win], info->rcoulomb[0]);
559 bRefinedVdW = !is_equal(info->rvdw[k_win], info->rvdw[0] );
560 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
561 info->nky[k_win] == info->nky[0] &&
562 info->nkz[k_win] == info->nkz[0]);
564 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
566 fprintf(fp, "Optimized PME settings:\n");
567 bCanUseOrigTPR = FALSE;
571 bCanUseOrigTPR = TRUE;
576 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
581 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
586 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],
587 info->nkx[0], info->nky[0], info->nkz[0]);
590 if (bCanUseOrigTPR && ntprs > 1)
592 fprintf(fp, "and original PME settings.\n");
597 /* Return the index of the mdp file that showed the highest performance
598 * and the optimal number of PME nodes */
600 *npme_optimal = winPME;
602 return bCanUseOrigTPR;
606 /* Get the commands we need to set up the runs from environment variables */
607 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
611 const char def_mpirun[] = "mpirun";
612 const char def_mdrun[] = "mdrun";
614 const char empty_mpirun[] = "";
616 /* Get the commands we need to set up the runs from environment variables */
619 if ( (cp = getenv("MPIRUN")) != NULL)
621 *cmd_mpirun = strdup(cp);
625 *cmd_mpirun = strdup(def_mpirun);
630 *cmd_mpirun = strdup(empty_mpirun);
633 if ( (cp = getenv("MDRUN" )) != NULL)
635 *cmd_mdrun = strdup(cp);
639 *cmd_mdrun = strdup(def_mdrun);
643 /* Check that the commands will run mdrun (perhaps via mpirun) by
644 * running a very quick test simulation. Requires MPI environment to
645 * be available if applicable. */
646 static void check_mdrun_works(gmx_bool bThreads,
647 const char *cmd_mpirun,
649 const char *cmd_mdrun)
651 char *command = NULL;
655 const char filename[] = "benchtest.log";
657 /* This string should always be identical to the one in copyrite.c,
658 * gmx_print_version_info() in the defined(GMX_MPI) section */
659 const char match_mpi[] = "MPI library: MPI";
660 const char match_mdrun[] = "Executable: ";
661 gmx_bool bMdrun = FALSE;
662 gmx_bool bMPI = FALSE;
664 /* Run a small test to see whether mpirun + mdrun work */
665 fprintf(stdout, "Making sure that mdrun can be executed. ");
668 snew(command, strlen(cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
669 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
673 snew(command, strlen(cmd_mpirun) + strlen(cmd_np) + strlen(cmd_mdrun) + strlen(filename) + 50);
674 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
676 fprintf(stdout, "Trying '%s' ... ", command);
677 make_backup(filename);
678 gmx_system_call(command);
680 /* Check if we find the characteristic string in the output: */
681 if (!gmx_fexist(filename))
683 gmx_fatal(FARGS, "Output from test run could not be found.");
686 fp = fopen(filename, "r");
687 /* We need to scan the whole output file, since sometimes the queuing system
688 * also writes stuff to stdout/err */
691 cp = fgets(line, STRLEN, fp);
694 if (str_starts(line, match_mdrun) )
698 if (str_starts(line, match_mpi) )
710 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
712 "seems to have been compiled with MPI instead.",
720 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
722 "seems to have been compiled without MPI support.",
729 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
733 fprintf(stdout, "passed.\n");
741 static void launch_simulation(
742 gmx_bool bLaunch, /* Should the simulation be launched? */
743 FILE *fp, /* General log file */
744 gmx_bool bThreads, /* whether to use threads */
745 char *cmd_mpirun, /* Command for mpirun */
746 char *cmd_np, /* Switch for -np or -ntmpi or empty */
747 char *cmd_mdrun, /* Command for mdrun */
748 char *args_for_mdrun, /* Arguments for mdrun */
749 const char *simulation_tpr, /* This tpr will be simulated */
750 int nPMEnodes) /* Number of PME nodes to use */
755 /* Make enough space for the system call command,
756 * (100 extra chars for -npme ... etc. options should suffice): */
757 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+100);
759 /* Note that the -passall options requires args_for_mdrun to be at the end
760 * of the command line string */
763 sprintf(command, "%s%s-npme %d -s %s %s",
764 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun);
768 sprintf(command, "%s%s%s -npme %d -s %s %s",
769 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun);
772 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
776 /* Now the real thing! */
779 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
782 gmx_system_call(command);
787 static void modify_PMEsettings(
788 gmx_int64_t simsteps, /* Set this value as number of time steps */
789 gmx_int64_t init_step, /* Set this value as init_step */
790 const char *fn_best_tpr, /* tpr file with the best performance */
791 const char *fn_sim_tpr) /* name of tpr file to be launched */
799 read_tpx_state(fn_best_tpr, ir, &state, NULL, &mtop);
801 /* Reset nsteps and init_step to the value of the input .tpr file */
802 ir->nsteps = simsteps;
803 ir->init_step = init_step;
805 /* Write the tpr file which will be launched */
806 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%"GMX_PRId64);
807 fprintf(stdout, buf, ir->nsteps);
809 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
815 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
817 /* Make additional TPR files with more computational load for the
818 * direct space processors: */
819 static void make_benchmark_tprs(
820 const char *fn_sim_tpr, /* READ : User-provided tpr file */
821 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
822 gmx_int64_t benchsteps, /* Number of time steps for benchmark runs */
823 gmx_int64_t statesteps, /* Step counter in checkpoint file */
824 real rmin, /* Minimal Coulomb radius */
825 real rmax, /* Maximal Coulomb radius */
826 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
827 int *ntprs, /* No. of TPRs to write, each with a different
828 rcoulomb and fourierspacing */
829 t_inputinfo *info, /* Contains information about mdp file options */
830 FILE *fp) /* Write the output here */
836 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
839 gmx_bool bNote = FALSE;
840 real add; /* Add this to rcoul for the next test */
841 real fac = 1.0; /* Scaling factor for Coulomb radius */
842 real fourierspacing; /* Basic fourierspacing from tpr */
845 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
846 *ntprs > 1 ? "s" : "", "%"GMX_PRId64, benchsteps > 1 ? "s" : "");
847 fprintf(stdout, buf, benchsteps);
850 sprintf(buf, " (adding %s steps from checkpoint file)", "%"GMX_PRId64);
851 fprintf(stdout, buf, statesteps);
852 benchsteps += statesteps;
854 fprintf(stdout, ".\n");
858 read_tpx_state(fn_sim_tpr, ir, &state, NULL, &mtop);
860 /* Check if some kind of PME was chosen */
861 if (EEL_PME(ir->coulombtype) == FALSE)
863 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
867 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
868 if ( (ir->cutoff_scheme != ecutsVERLET) &&
869 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
871 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
872 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
874 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
875 else if (ir->rcoulomb > ir->rlist)
877 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
878 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
881 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
883 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
887 /* Reduce the number of steps for the benchmarks */
888 info->orig_sim_steps = ir->nsteps;
889 ir->nsteps = benchsteps;
890 /* We must not use init_step from the input tpr file for the benchmarks */
891 info->orig_init_step = ir->init_step;
894 /* For PME-switch potentials, keep the radial distance of the buffer region */
895 nlist_buffer = ir->rlist - ir->rcoulomb;
897 /* Determine length of triclinic box vectors */
898 for (d = 0; d < DIM; d++)
901 for (i = 0; i < DIM; i++)
903 box_size[d] += state.box[d][i]*state.box[d][i];
905 box_size[d] = sqrt(box_size[d]);
908 if (ir->fourier_spacing > 0)
910 info->fsx[0] = ir->fourier_spacing;
911 info->fsy[0] = ir->fourier_spacing;
912 info->fsz[0] = ir->fourier_spacing;
916 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
917 info->fsx[0] = box_size[XX]/ir->nkx;
918 info->fsy[0] = box_size[YY]/ir->nky;
919 info->fsz[0] = box_size[ZZ]/ir->nkz;
922 /* If no value for the fourierspacing was provided on the command line, we
923 * use the reconstruction from the tpr file */
924 if (ir->fourier_spacing > 0)
926 /* Use the spacing from the tpr */
927 fourierspacing = ir->fourier_spacing;
931 /* Use the maximum observed spacing */
932 fourierspacing = max(max(info->fsx[0], info->fsy[0]), info->fsz[0]);
935 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
937 /* For performance comparisons the number of particles is useful to have */
938 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
940 /* Print information about settings of which some are potentially modified: */
941 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
942 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
943 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
944 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
945 if (ir_vdw_switched(ir))
947 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
949 if (EPME_SWITCHED(ir->coulombtype))
951 fprintf(fp, " rlist : %f nm\n", ir->rlist);
953 if (ir->rlistlong != max_cutoff(ir->rvdw, ir->rcoulomb))
955 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
958 /* Print a descriptive line about the tpr settings tested */
959 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
960 fprintf(fp, " No. scaling rcoulomb");
961 fprintf(fp, " nkx nky nkz");
962 fprintf(fp, " spacing");
963 if (evdwCUT == ir->vdwtype)
965 fprintf(fp, " rvdw");
967 if (EPME_SWITCHED(ir->coulombtype))
969 fprintf(fp, " rlist");
971 if (ir->rlistlong != max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb)) )
973 fprintf(fp, " rlistlong");
975 fprintf(fp, " tpr file\n");
977 /* Loop to create the requested number of tpr input files */
978 for (j = 0; j < *ntprs; j++)
980 /* The first .tpr is the provided one, just need to modify nsteps,
981 * so skip the following block */
984 /* Determine which Coulomb radii rc to use in the benchmarks */
985 add = (rmax-rmin)/(*ntprs-1);
986 if (is_equal(rmin, info->rcoulomb[0]))
988 ir->rcoulomb = rmin + j*add;
990 else if (is_equal(rmax, info->rcoulomb[0]))
992 ir->rcoulomb = rmin + (j-1)*add;
996 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
997 add = (rmax-rmin)/(*ntprs-2);
998 ir->rcoulomb = rmin + (j-1)*add;
1001 /* Determine the scaling factor fac */
1002 fac = ir->rcoulomb/info->rcoulomb[0];
1004 /* Scale the Fourier grid spacing */
1005 ir->nkx = ir->nky = ir->nkz = 0;
1006 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
1008 /* Adjust other radii since various conditions need to be fulfilled */
1009 if (eelPME == ir->coulombtype)
1011 /* plain PME, rcoulomb must be equal to rlist */
1012 ir->rlist = ir->rcoulomb;
1016 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1017 ir->rlist = ir->rcoulomb + nlist_buffer;
1020 if (bScaleRvdw && evdwCUT == ir->vdwtype)
1022 if (ecutsVERLET == ir->cutoff_scheme)
1024 /* With Verlet, the van der Waals radius must always equal the Coulomb radius */
1025 ir->rvdw = ir->rcoulomb;
1029 /* For vdw cutoff, rvdw >= rlist */
1030 ir->rvdw = max(info->rvdw[0], ir->rlist);
1034 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1036 } /* end of "if (j != 0)" */
1038 /* for j==0: Save the original settings
1039 * for j >0: Save modified radii and Fourier grids */
1040 info->rcoulomb[j] = ir->rcoulomb;
1041 info->rvdw[j] = ir->rvdw;
1042 info->nkx[j] = ir->nkx;
1043 info->nky[j] = ir->nky;
1044 info->nkz[j] = ir->nkz;
1045 info->rlist[j] = ir->rlist;
1046 info->rlistlong[j] = ir->rlistlong;
1047 info->fsx[j] = fac*fourierspacing;
1048 info->fsy[j] = fac*fourierspacing;
1049 info->fsz[j] = fac*fourierspacing;
1051 /* Write the benchmark tpr file */
1052 strncpy(fn_bench_tprs[j], fn_sim_tpr, strlen(fn_sim_tpr)-strlen(".tpr"));
1053 sprintf(buf, "_bench%.2d.tpr", j);
1054 strcat(fn_bench_tprs[j], buf);
1055 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1056 fprintf(stdout, "%"GMX_PRId64, ir->nsteps);
1059 fprintf(stdout, ", scaling factor %f\n", fac);
1063 fprintf(stdout, ", unmodified settings\n");
1066 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1068 /* Write information about modified tpr settings to log file */
1069 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1070 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1071 fprintf(fp, " %9f ", info->fsx[j]);
1072 if (evdwCUT == ir->vdwtype)
1074 fprintf(fp, "%10f", ir->rvdw);
1076 if (EPME_SWITCHED(ir->coulombtype))
1078 fprintf(fp, "%10f", ir->rlist);
1080 if (info->rlistlong[0] != max_cutoff(info->rlist[0], max_cutoff(info->rvdw[0], info->rcoulomb[0])) )
1082 fprintf(fp, "%10f", ir->rlistlong);
1084 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1086 /* Make it clear to the user that some additional settings were modified */
1087 if (!is_equal(ir->rvdw, info->rvdw[0])
1088 || !is_equal(ir->rlistlong, info->rlistlong[0]) )
1095 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1096 "other input settings were also changed (see table above).\n"
1097 "Please check if the modified settings are appropriate.\n");
1105 /* Rename the files we want to keep to some meaningful filename and
1106 * delete the rest */
1107 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1108 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1110 char numstring[STRLEN];
1111 char newfilename[STRLEN];
1112 const char *fn = NULL;
1117 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1119 for (i = 0; i < nfile; i++)
1121 opt = (char *)fnm[i].opt;
1122 if (strcmp(opt, "-p") == 0)
1124 /* do nothing; keep this file */
1127 else if (strcmp(opt, "-bg") == 0)
1129 /* Give the log file a nice name so one can later see which parameters were used */
1130 numstring[0] = '\0';
1133 sprintf(numstring, "_%d", nr);
1135 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1136 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1138 fprintf(stdout, "renaming log file to %s\n", newfilename);
1139 make_backup(newfilename);
1140 rename(opt2fn("-bg", nfile, fnm), newfilename);
1143 else if (strcmp(opt, "-err") == 0)
1145 /* This file contains the output of stderr. We want to keep it in
1146 * cases where there have been problems. */
1147 fn = opt2fn(opt, nfile, fnm);
1148 numstring[0] = '\0';
1151 sprintf(numstring, "_%d", nr);
1153 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1158 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1159 make_backup(newfilename);
1160 rename(fn, newfilename);
1164 fprintf(stdout, "Deleting %s\n", fn);
1169 /* Delete the files which are created for each benchmark run: (options -b*) */
1170 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1172 remove_if_exists(opt2fn(opt, nfile, fnm));
1179 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1182 /* Create a list of numbers of PME nodes to test */
1183 static void make_npme_list(
1184 const char *npmevalues_opt, /* Make a complete list with all
1185 * possibilities or a short list that keeps only
1186 * reasonable numbers of PME nodes */
1187 int *nentries, /* Number of entries we put in the nPMEnodes list */
1188 int *nPMEnodes[], /* Each entry contains the value for -npme */
1189 int nnodes, /* Total number of nodes to do the tests on */
1190 int minPMEnodes, /* Minimum number of PME nodes */
1191 int maxPMEnodes) /* Maximum number of PME nodes */
1194 int min_factor = 1; /* We request that npp and npme have this minimal
1195 * largest common factor (depends on npp) */
1196 int nlistmax; /* Max. list size */
1197 int nlist; /* Actual number of entries in list */
1201 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1202 if (0 == strcmp(npmevalues_opt, "all") )
1206 else if (0 == strcmp(npmevalues_opt, "subset") )
1208 eNPME = eNpmeSubset;
1210 else /* "auto" or "range" */
1216 else if (nnodes < 128)
1218 eNPME = eNpmeReduced;
1222 eNPME = eNpmeSubset;
1226 /* Calculate how many entries we could possibly have (in case of -npme all) */
1229 nlistmax = maxPMEnodes - minPMEnodes + 3;
1230 if (0 == minPMEnodes)
1240 /* Now make the actual list which is at most of size nlist */
1241 snew(*nPMEnodes, nlistmax);
1242 nlist = 0; /* start counting again, now the real entries in the list */
1243 for (i = 0; i < nlistmax - 2; i++)
1245 npme = maxPMEnodes - i;
1256 /* For 2d PME we want a common largest factor of at least the cube
1257 * root of the number of PP nodes */
1258 min_factor = (int) pow(npp, 1.0/3.0);
1261 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1264 if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1266 (*nPMEnodes)[nlist] = npme;
1270 /* We always test 0 PME nodes and the automatic number */
1271 *nentries = nlist + 2;
1272 (*nPMEnodes)[nlist ] = 0;
1273 (*nPMEnodes)[nlist+1] = -1;
1275 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1276 for (i = 0; i < *nentries-1; i++)
1278 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1280 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1284 /* Allocate memory to store the performance data */
1285 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1290 for (k = 0; k < ntprs; k++)
1292 snew(perfdata[k], datasets);
1293 for (i = 0; i < datasets; i++)
1295 for (j = 0; j < repeats; j++)
1297 snew(perfdata[k][i].Gcycles, repeats);
1298 snew(perfdata[k][i].ns_per_day, repeats);
1299 snew(perfdata[k][i].PME_f_load, repeats);
1306 /* Check for errors on mdrun -h */
1307 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp,
1308 const t_filenm *fnm, int nfile)
1310 const char *fn = NULL;
1311 char *command, *msg;
1315 snew(command, length + 15);
1316 snew(msg, length + 500);
1318 fprintf(stdout, "Making sure the benchmarks can be executed by running just 1 step...\n");
1319 sprintf(command, "%s -nsteps 1 -quiet", mdrun_cmd_line);
1320 fprintf(stdout, "Executing '%s' ...\n", command);
1321 ret = gmx_system_call(command);
1325 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1326 * get the error message from mdrun itself */
1327 sprintf(msg, "Cannot run the benchmark simulations! Please check the error message of\n"
1328 "mdrun for the source of the problem. Did you provide a command line\n"
1329 "argument that neither g_tune_pme nor mdrun understands? Offending command:\n"
1330 "\n%s\n\n", command);
1332 fprintf(stderr, "%s", msg);
1334 fprintf(fp, "%s", msg);
1338 fprintf(stdout, "Benchmarks can be executed!\n");
1340 /* Clean up the benchmark output files we just created */
1341 fprintf(stdout, "Cleaning up ...\n");
1342 remove_if_exists(opt2fn("-bc", nfile, fnm));
1343 remove_if_exists(opt2fn("-be", nfile, fnm));
1344 remove_if_exists(opt2fn("-bcpo", nfile, fnm));
1345 remove_if_exists(opt2fn("-bg", nfile, fnm));
1352 static void do_the_tests(
1353 FILE *fp, /* General g_tune_pme output file */
1354 char **tpr_names, /* Filenames of the input files to test */
1355 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1356 int minPMEnodes, /* Min fraction of nodes to use for PME */
1357 int npme_fixed, /* If >= -1, test fixed number of PME
1359 const char *npmevalues_opt, /* Which -npme values should be tested */
1360 t_perf **perfdata, /* Here the performace data is stored */
1361 int *pmeentries, /* Entries in the nPMEnodes list */
1362 int repeats, /* Repeat each test this often */
1363 int nnodes, /* Total number of nodes = nPP + nPME */
1364 int nr_tprs, /* Total number of tpr files to test */
1365 gmx_bool bThreads, /* Threads or MPI? */
1366 char *cmd_mpirun, /* mpirun command string */
1367 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1368 char *cmd_mdrun, /* mdrun command string */
1369 char *cmd_args_bench, /* arguments for mdrun in a string */
1370 const t_filenm *fnm, /* List of filenames from command line */
1371 int nfile, /* Number of files specified on the cmdl. */
1372 int presteps, /* DLB equilibration steps, is checked */
1373 gmx_int64_t cpt_steps) /* Time step counter in the checkpoint */
1375 int i, nr, k, ret, count = 0, totaltests;
1376 int *nPMEnodes = NULL;
1379 char *command, *cmd_stub;
1381 gmx_bool bResetProblem = FALSE;
1382 gmx_bool bFirst = TRUE;
1385 /* This string array corresponds to the eParselog enum type at the start
1387 const char* ParseLog[] = {
1389 "Logfile not found!",
1390 "No timings, logfile truncated?",
1391 "Run was terminated.",
1392 "Counters were not reset properly.",
1393 "No DD grid found for these settings.",
1394 "TPX version conflict!",
1395 "mdrun was not started in parallel!",
1396 "Number of PP nodes has a prime factor that is too large.",
1399 char str_PME_f_load[13];
1402 /* Allocate space for the mdrun command line. 100 extra characters should
1403 be more than enough for the -npme etcetera arguments */
1404 cmdline_length = strlen(cmd_mpirun)
1407 + strlen(cmd_args_bench)
1408 + strlen(tpr_names[0]) + 100;
1409 snew(command, cmdline_length);
1410 snew(cmd_stub, cmdline_length);
1412 /* Construct the part of the command line that stays the same for all tests: */
1415 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1419 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1422 /* Create a list of numbers of PME nodes to test */
1423 if (npme_fixed < -1)
1425 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1426 nnodes, minPMEnodes, maxPMEnodes);
1432 nPMEnodes[0] = npme_fixed;
1433 fprintf(stderr, "Will use a fixed number of %d PME-only nodes.\n", nPMEnodes[0]);
1438 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1440 finalize(opt2fn("-p", nfile, fnm));
1444 /* Allocate one dataset for each tpr input file: */
1445 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1447 /*****************************************/
1448 /* Main loop over all tpr files to test: */
1449 /*****************************************/
1450 totaltests = nr_tprs*(*pmeentries)*repeats;
1451 for (k = 0; k < nr_tprs; k++)
1453 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1454 fprintf(fp, "PME nodes Gcycles ns/day PME/f Remark\n");
1455 /* Loop over various numbers of PME nodes: */
1456 for (i = 0; i < *pmeentries; i++)
1458 pd = &perfdata[k][i];
1460 /* Loop over the repeats for each scenario: */
1461 for (nr = 0; nr < repeats; nr++)
1463 pd->nPMEnodes = nPMEnodes[i];
1465 /* Add -npme and -s to the command line and save it. Note that
1466 * the -passall (if set) options requires cmd_args_bench to be
1467 * at the end of the command line string */
1468 snew(pd->mdrun_cmd_line, cmdline_length);
1469 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s",
1470 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench);
1472 /* To prevent that all benchmarks fail due to a show-stopper argument
1473 * on the mdrun command line, we make a quick check first */
1476 make_sure_it_runs(pd->mdrun_cmd_line, cmdline_length, fp, fnm, nfile);
1480 /* Do a benchmark simulation: */
1483 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1489 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1490 (100.0*count)/totaltests,
1491 k+1, nr_tprs, i+1, *pmeentries, buf);
1492 make_backup(opt2fn("-err", nfile, fnm));
1493 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1494 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1495 gmx_system_call(command);
1497 /* Collect the performance data from the log file; also check stderr
1498 * for fatal errors */
1499 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1500 pd, nr, presteps, cpt_steps, nnodes);
1501 if ((presteps > 0) && (ret == eParselogResetProblem))
1503 bResetProblem = TRUE;
1506 if (-1 == pd->nPMEnodes)
1508 sprintf(buf, "(%3d)", pd->guessPME);
1516 if (pd->PME_f_load[nr] > 0.0)
1518 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1522 sprintf(str_PME_f_load, "%s", " - ");
1525 /* Write the data we got to disk */
1526 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1527 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1528 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1530 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1536 /* Do some cleaning up and delete the files we do not need any more */
1537 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1539 /* If the first run with this number of processors already failed, do not try again: */
1540 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1542 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1543 count += repeats-(nr+1);
1546 } /* end of repeats loop */
1547 } /* end of -npme loop */
1548 } /* end of tpr file loop */
1553 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1561 static void check_input(
1568 real maxPMEfraction,
1569 real minPMEfraction,
1571 gmx_int64_t bench_nsteps,
1572 const t_filenm *fnm,
1582 /* Make sure the input file exists */
1583 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1585 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1588 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1589 if ( (0 == strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1591 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1592 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1595 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1598 gmx_fatal(FARGS, "Number of repeats < 0!");
1601 /* Check number of nodes */
1604 gmx_fatal(FARGS, "Number of nodes/threads must be a positive integer.");
1607 /* Automatically choose -ntpr if not set */
1617 /* Set a reasonable scaling factor for rcoulomb */
1620 *rmax = rcoulomb * 1.2;
1623 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1629 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1633 /* Make shure that rmin <= rcoulomb <= rmax */
1642 if (!(*rmin <= *rmax) )
1644 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1645 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1647 /* Add test scenarios if rmin or rmax were set */
1650 if (!is_equal(*rmin, rcoulomb) && (*ntprs == 1) )
1653 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1656 if (!is_equal(*rmax, rcoulomb) && (*ntprs == 1) )
1659 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1664 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1665 if (!is_equal(*rmax, rcoulomb) || !is_equal(*rmin, rcoulomb) )
1667 *ntprs = max(*ntprs, 2);
1670 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1671 if (!is_equal(*rmax, rcoulomb) && !is_equal(*rmin, rcoulomb) )
1673 *ntprs = max(*ntprs, 3);
1678 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1683 if (is_equal(*rmin, rcoulomb) && is_equal(rcoulomb, *rmax)) /* We have just a single rc */
1685 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1686 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1687 "with correspondingly adjusted PME grid settings\n");
1692 /* Check whether max and min fraction are within required values */
1693 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1695 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1697 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1699 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1701 if (maxPMEfraction < minPMEfraction)
1703 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1706 /* Check whether the number of steps - if it was set - has a reasonable value */
1707 if (bench_nsteps < 0)
1709 gmx_fatal(FARGS, "Number of steps must be positive.");
1712 if (bench_nsteps > 10000 || bench_nsteps < 100)
1714 fprintf(stderr, "WARNING: steps=");
1715 fprintf(stderr, "%"GMX_PRId64, bench_nsteps);
1716 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1721 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1724 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1727 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1729 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1733 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1734 * only. We need to check whether the requested number of PME-only nodes
1736 if (npme_fixed > -1)
1738 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1739 if (2*npme_fixed > nnodes)
1741 gmx_fatal(FARGS, "Cannot have more than %d PME-only nodes for a total of %d nodes (you chose %d).\n",
1742 nnodes/2, nnodes, npme_fixed);
1744 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1746 fprintf(stderr, "WARNING: Only %g percent of the nodes are assigned as PME-only nodes.\n",
1747 100.0*((real)npme_fixed / (real)nnodes));
1749 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1751 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1752 " fixed number of PME-only nodes is requested with -fix.\n");
1758 /* Returns TRUE when "opt" is needed at launch time */
1759 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1761 if (0 == strncmp(opt, "-swap", 5))
1766 /* Apart from the input .tpr and the output log files we need all options that
1767 * were set on the command line and that do not start with -b */
1768 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2)
1769 || 0 == strncmp(opt, "-err", 4) || 0 == strncmp(opt, "-p", 2) )
1778 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1779 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1781 /* Apart from the input .tpr, all files starting with "-b" are for
1782 * _b_enchmark files exclusively */
1783 if (0 == strncmp(opt, "-s", 2))
1788 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2))
1790 if (!bOptional || bSet)
1807 if (bSet) /* These are additional input files like -cpi -ei */
1820 /* Adds 'buf' to 'str' */
1821 static void add_to_string(char **str, char *buf)
1826 len = strlen(*str) + strlen(buf) + 1;
1832 /* Create the command line for the benchmark as well as for the real run */
1833 static void create_command_line_snippets(
1834 gmx_bool bAppendFiles,
1835 gmx_bool bKeepAndNumCPT,
1836 gmx_bool bResetHWay,
1840 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1841 char *cmd_args_launch[], /* command line arguments for simulation run */
1842 char extra_args[]) /* Add this to the end of the command line */
1847 char strbuf[STRLEN];
1850 /* strlen needs at least '\0' as a string: */
1851 snew(*cmd_args_bench, 1);
1852 snew(*cmd_args_launch, 1);
1853 *cmd_args_launch[0] = '\0';
1854 *cmd_args_bench[0] = '\0';
1857 /*******************************************/
1858 /* 1. Process other command line arguments */
1859 /*******************************************/
1862 /* Add equilibration steps to benchmark options */
1863 sprintf(strbuf, "-resetstep %d ", presteps);
1864 add_to_string(cmd_args_bench, strbuf);
1866 /* These switches take effect only at launch time */
1867 if (FALSE == bAppendFiles)
1869 add_to_string(cmd_args_launch, "-noappend ");
1873 add_to_string(cmd_args_launch, "-cpnum ");
1877 add_to_string(cmd_args_launch, "-resethway ");
1880 /********************/
1881 /* 2. Process files */
1882 /********************/
1883 for (i = 0; i < nfile; i++)
1885 opt = (char *)fnm[i].opt;
1886 name = opt2fn(opt, nfile, fnm);
1888 /* Strbuf contains the options, now let's sort out where we need that */
1889 sprintf(strbuf, "%s %s ", opt, name);
1891 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1893 /* All options starting with -b* need the 'b' removed,
1894 * therefore overwrite strbuf */
1895 if (0 == strncmp(opt, "-b", 2))
1897 sprintf(strbuf, "-%s %s ", &opt[2], name);
1900 add_to_string(cmd_args_bench, strbuf);
1903 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
1905 add_to_string(cmd_args_launch, strbuf);
1909 add_to_string(cmd_args_bench, extra_args);
1910 add_to_string(cmd_args_launch, extra_args);
1914 /* Set option opt */
1915 static void setopt(const char *opt, int nfile, t_filenm fnm[])
1919 for (i = 0; (i < nfile); i++)
1921 if (strcmp(opt, fnm[i].opt) == 0)
1923 fnm[i].flag |= ffSET;
1929 /* This routine inspects the tpr file and ...
1930 * 1. checks for output files that get triggered by a tpr option. These output
1931 * files are marked as 'set' to allow for a proper cleanup after each
1933 * 2. returns the PME:PP load ratio
1934 * 3. returns rcoulomb from the tpr */
1935 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1937 gmx_bool bPull; /* Is pulling requested in .tpr file? */
1938 gmx_bool bTpi; /* Is test particle insertion requested? */
1939 gmx_bool bFree; /* Is a free energy simulation requested? */
1940 gmx_bool bNM; /* Is a normal mode analysis requested? */
1941 gmx_bool bSwap; /* Is water/ion position swapping requested? */
1947 /* Check tpr file for options that trigger extra output files */
1948 read_tpx_state(opt2fn("-s", nfile, fnm), &ir, &state, NULL, &mtop);
1949 bPull = (epullNO != ir.ePull);
1950 bFree = (efepNO != ir.efep );
1951 bNM = (eiNM == ir.eI );
1952 bSwap = (eswapNO != ir.eSwapCoords);
1953 bTpi = EI_TPI(ir.eI);
1955 /* Set these output files on the tuning command-line */
1958 setopt("-pf", nfile, fnm);
1959 setopt("-px", nfile, fnm);
1963 setopt("-dhdl", nfile, fnm);
1967 setopt("-tpi", nfile, fnm);
1968 setopt("-tpid", nfile, fnm);
1972 setopt("-mtx", nfile, fnm);
1976 setopt("-swap", nfile, fnm);
1979 *rcoulomb = ir.rcoulomb;
1981 /* Return the estimate for the number of PME nodes */
1982 return pme_load_estimate(&mtop, &ir, state.box);
1986 static void couple_files_options(int nfile, t_filenm fnm[])
1989 gmx_bool bSet, bBench;
1994 for (i = 0; i < nfile; i++)
1996 opt = (char *)fnm[i].opt;
1997 bSet = ((fnm[i].flag & ffSET) != 0);
1998 bBench = (0 == strncmp(opt, "-b", 2));
2000 /* Check optional files */
2001 /* If e.g. -eo is set, then -beo also needs to be set */
2002 if (is_optional(&fnm[i]) && bSet && !bBench)
2004 sprintf(buf, "-b%s", &opt[1]);
2005 setopt(buf, nfile, fnm);
2007 /* If -beo is set, then -eo also needs to be! */
2008 if (is_optional(&fnm[i]) && bSet && bBench)
2010 sprintf(buf, "-%s", &opt[2]);
2011 setopt(buf, nfile, fnm);
2017 #define BENCHSTEPS (1000)
2019 int gmx_tune_pme(int argc, char *argv[])
2021 const char *desc[] = {
2022 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of processors/threads, [THISMODULE] systematically",
2023 "times [gmx-mdrun] with various numbers of PME-only nodes and determines",
2024 "which setting is fastest. It will also test whether performance can",
2025 "be enhanced by shifting load from the reciprocal to the real space",
2026 "part of the Ewald sum. ",
2027 "Simply pass your [TT].tpr[tt] file to [THISMODULE] together with other options",
2028 "for [gmx-mdrun] as needed.[PAR]",
2029 "Which executables are used can be set in the environment variables",
2030 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
2031 "will be used as defaults. Note that for certain MPI frameworks you",
2032 "need to provide a machine- or hostfile. This can also be passed",
2033 "via the MPIRUN variable, e.g.[PAR]",
2034 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
2035 "Please call [THISMODULE] with the normal options you would pass to",
2036 "[gmx-mdrun] and add [TT]-np[tt] for the number of processors to perform the",
2037 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2038 "to repeat each test several times to get better statistics. [PAR]",
2039 "[THISMODULE] can test various real space / reciprocal space workloads",
2040 "for you. With [TT]-ntpr[tt] you control how many extra [TT].tpr[tt] files will be",
2041 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2042 "Typically, the first test (number 0) will be with the settings from the input",
2043 "[TT].tpr[tt] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2044 "specified by [TT]-rmax[tt] with a somwhat smaller PME grid at the same time. ",
2045 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2046 "The remaining [TT].tpr[tt] files will have equally-spaced Coulomb radii (and Fourier "
2047 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2048 "if you just seek the optimal number of PME-only nodes; in that case",
2049 "your input [TT].tpr[tt] file will remain unchanged.[PAR]",
2050 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2051 "MD systems. The dynamic load balancing needs about 100 time steps",
2052 "to adapt to local load imbalances, therefore the time step counters",
2053 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2054 "for a higher accuarcy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2055 "From the 'DD' load imbalance entries in the md.log output file you",
2056 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2057 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2058 "After calling [gmx-mdrun] several times, detailed performance information",
2059 "is available in the output file [TT]perf.out[tt].",
2060 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2061 "(options [TT]-b[tt]*), these will be automatically deleted after each test.[PAR]",
2062 "If you want the simulation to be started automatically with the",
2063 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2068 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2069 real maxPMEfraction = 0.50;
2070 real minPMEfraction = 0.25;
2071 int maxPMEnodes, minPMEnodes;
2072 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2073 float guessPMEnodes;
2074 int npme_fixed = -2; /* If >= -1, use only this number
2075 * of PME-only nodes */
2077 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2078 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2079 gmx_bool bScaleRvdw = TRUE;
2080 gmx_int64_t bench_nsteps = BENCHSTEPS;
2081 gmx_int64_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2082 gmx_int64_t cpt_steps = 0; /* Step counter in .cpt input file */
2083 int presteps = 100; /* Do a full cycle reset after presteps steps */
2084 gmx_bool bOverwrite = FALSE, bKeepTPR;
2085 gmx_bool bLaunch = FALSE;
2086 char *ExtraArgs = NULL;
2087 char **tpr_names = NULL;
2088 const char *simulation_tpr = NULL;
2089 int best_npme, best_tpr;
2090 int sim_part = 1; /* For benchmarks with checkpoint files */
2093 /* Default program names if nothing else is found */
2094 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2095 char *cmd_args_bench, *cmd_args_launch;
2096 char *cmd_np = NULL;
2098 t_perf **perfdata = NULL;
2104 /* Print out how long the tuning took */
2107 static t_filenm fnm[] = {
2109 { efOUT, "-p", "perf", ffWRITE },
2110 { efLOG, "-err", "bencherr", ffWRITE },
2111 { efTPX, "-so", "tuned", ffWRITE },
2113 { efTPX, NULL, NULL, ffREAD },
2114 { efTRN, "-o", NULL, ffWRITE },
2115 { efCOMPRESSED, "-x", NULL, ffOPTWR },
2116 { efCPT, "-cpi", NULL, ffOPTRD },
2117 { efCPT, "-cpo", NULL, ffOPTWR },
2118 { efSTO, "-c", "confout", ffWRITE },
2119 { efEDR, "-e", "ener", ffWRITE },
2120 { efLOG, "-g", "md", ffWRITE },
2121 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2122 { efXVG, "-field", "field", ffOPTWR },
2123 { efXVG, "-table", "table", ffOPTRD },
2124 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2125 { efXVG, "-tablep", "tablep", ffOPTRD },
2126 { efXVG, "-tableb", "table", ffOPTRD },
2127 { efTRX, "-rerun", "rerun", ffOPTRD },
2128 { efXVG, "-tpi", "tpi", ffOPTWR },
2129 { efXVG, "-tpid", "tpidist", ffOPTWR },
2130 { efEDI, "-ei", "sam", ffOPTRD },
2131 { efXVG, "-eo", "edsam", ffOPTWR },
2132 { efXVG, "-devout", "deviatie", ffOPTWR },
2133 { efXVG, "-runav", "runaver", ffOPTWR },
2134 { efXVG, "-px", "pullx", ffOPTWR },
2135 { efXVG, "-pf", "pullf", ffOPTWR },
2136 { efXVG, "-ro", "rotation", ffOPTWR },
2137 { efLOG, "-ra", "rotangles", ffOPTWR },
2138 { efLOG, "-rs", "rotslabs", ffOPTWR },
2139 { efLOG, "-rt", "rottorque", ffOPTWR },
2140 { efMTX, "-mtx", "nm", ffOPTWR },
2141 { efNDX, "-dn", "dipole", ffOPTWR },
2142 { efXVG, "-swap", "swapions", ffOPTWR },
2143 /* Output files that are deleted after each benchmark run */
2144 { efTRN, "-bo", "bench", ffWRITE },
2145 { efXTC, "-bx", "bench", ffWRITE },
2146 { efCPT, "-bcpo", "bench", ffWRITE },
2147 { efSTO, "-bc", "bench", ffWRITE },
2148 { efEDR, "-be", "bench", ffWRITE },
2149 { efLOG, "-bg", "bench", ffWRITE },
2150 { efXVG, "-beo", "benchedo", ffOPTWR },
2151 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2152 { efXVG, "-bfield", "benchfld", ffOPTWR },
2153 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2154 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2155 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2156 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2157 { efXVG, "-bpx", "benchpx", ffOPTWR },
2158 { efXVG, "-bpf", "benchpf", ffOPTWR },
2159 { efXVG, "-bro", "benchrot", ffOPTWR },
2160 { efLOG, "-bra", "benchrota", ffOPTWR },
2161 { efLOG, "-brs", "benchrots", ffOPTWR },
2162 { efLOG, "-brt", "benchrott", ffOPTWR },
2163 { efMTX, "-bmtx", "benchn", ffOPTWR },
2164 { efNDX, "-bdn", "bench", ffOPTWR },
2165 { efXVG, "-bswap", "benchswp", ffOPTWR }
2168 gmx_bool bThreads = FALSE;
2172 const char *procstring[] =
2173 { NULL, "-np", "-n", "none", NULL };
2174 const char *npmevalues_opt[] =
2175 { NULL, "auto", "all", "subset", NULL };
2177 gmx_bool bAppendFiles = TRUE;
2178 gmx_bool bKeepAndNumCPT = FALSE;
2179 gmx_bool bResetCountersHalfWay = FALSE;
2180 gmx_bool bBenchmark = TRUE;
2182 output_env_t oenv = NULL;
2185 /***********************/
2186 /* g_tune_pme options: */
2187 /***********************/
2188 { "-np", FALSE, etINT, {&nnodes},
2189 "Number of nodes to run the tests on (must be > 2 for separate PME nodes)" },
2190 { "-npstring", FALSE, etENUM, {procstring},
2191 "Specify the number of processors to [TT]$MPIRUN[tt] using this string"},
2192 { "-ntmpi", FALSE, etINT, {&nthreads},
2193 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2194 { "-r", FALSE, etINT, {&repeats},
2195 "Repeat each test this often" },
2196 { "-max", FALSE, etREAL, {&maxPMEfraction},
2197 "Max fraction of PME nodes to test with" },
2198 { "-min", FALSE, etREAL, {&minPMEfraction},
2199 "Min fraction of PME nodes to test with" },
2200 { "-npme", FALSE, etENUM, {npmevalues_opt},
2201 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2202 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2203 { "-fix", FALSE, etINT, {&npme_fixed},
2204 "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."},
2205 { "-rmax", FALSE, etREAL, {&rmax},
2206 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2207 { "-rmin", FALSE, etREAL, {&rmin},
2208 "If >0, minimal rcoulomb for -ntpr>1" },
2209 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2210 "Scale rvdw along with rcoulomb"},
2211 { "-ntpr", FALSE, etINT, {&ntprs},
2212 "Number of [TT].tpr[tt] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2213 "If < 1, automatically choose the number of [TT].tpr[tt] files to test" },
2214 { "-steps", FALSE, etINT64, {&bench_nsteps},
2215 "Take timings for this many steps in the benchmark runs" },
2216 { "-resetstep", FALSE, etINT, {&presteps},
2217 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2218 { "-simsteps", FALSE, etINT64, {&new_sim_nsteps},
2219 "If non-negative, perform this many steps in the real run (overwrites nsteps from [TT].tpr[tt], add [TT].cpt[tt] steps)" },
2220 { "-launch", FALSE, etBOOL, {&bLaunch},
2221 "Launch the real simulation after optimization" },
2222 { "-bench", FALSE, etBOOL, {&bBenchmark},
2223 "Run the benchmarks or just create the input [TT].tpr[tt] files?" },
2224 /******************/
2225 /* mdrun options: */
2226 /******************/
2227 /* We let g_tune_pme parse and understand these options, because we need to
2228 * prevent that they appear on the mdrun command line for the benchmarks */
2229 { "-append", FALSE, etBOOL, {&bAppendFiles},
2230 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2231 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2232 "Keep and number checkpoint files (launch only)" },
2233 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2234 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2237 #define NFILE asize(fnm)
2239 seconds = gmx_gettime();
2241 if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2242 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2248 /* Store the remaining unparsed command line entries in a string which
2249 * is then attached to the mdrun command line */
2251 ExtraArgs[0] = '\0';
2252 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2254 add_to_string(&ExtraArgs, argv[i]);
2255 add_to_string(&ExtraArgs, " ");
2258 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2261 if (opt2parg_bSet("-npstring", asize(pa), pa))
2263 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2268 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2270 /* and now we just set this; a bit of an ugly hack*/
2273 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2274 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2276 /* Automatically set -beo options if -eo is set etc. */
2277 couple_files_options(NFILE, fnm);
2279 /* Construct the command line arguments for benchmark runs
2280 * as well as for the simulation run */
2283 sprintf(bbuf, " -ntmpi %d ", nthreads);
2287 /* This string will be used for MPI runs and will appear after the
2288 * mpirun command. */
2289 if (strcmp(procstring[0], "none") != 0)
2291 sprintf(bbuf, " %s %d ", procstring[0], nnodes);
2301 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2302 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs);
2304 /* Read in checkpoint file if requested */
2306 if (opt2bSet("-cpi", NFILE, fnm))
2309 cr->duty = DUTY_PP; /* makes the following routine happy */
2310 read_checkpoint_simulation_part(opt2fn("-cpi", NFILE, fnm),
2311 &sim_part, &cpt_steps, cr,
2312 FALSE, NFILE, fnm, NULL, NULL);
2315 /* sim_part will now be 1 if no checkpoint file was found */
2318 gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi", NFILE, fnm));
2322 /* Open performance output file and write header info */
2323 fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2325 /* Make a quick consistency check of command line parameters */
2326 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2327 maxPMEfraction, minPMEfraction, npme_fixed,
2328 bench_nsteps, fnm, NFILE, sim_part, presteps,
2331 /* Determine the maximum and minimum number of PME nodes to test,
2332 * the actual list of settings is build in do_the_tests(). */
2333 if ((nnodes > 2) && (npme_fixed < -1))
2335 if (0 == strcmp(npmevalues_opt[0], "auto"))
2337 /* Determine the npme range automatically based on the PME:PP load guess */
2338 if (guessPMEratio > 1.0)
2340 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2341 maxPMEnodes = nnodes/2;
2342 minPMEnodes = nnodes/2;
2346 /* PME : PP load is in the range 0..1, let's test around the guess */
2347 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2348 minPMEnodes = floor(0.7*guessPMEnodes);
2349 maxPMEnodes = ceil(1.6*guessPMEnodes);
2350 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2355 /* Determine the npme range based on user input */
2356 maxPMEnodes = floor(maxPMEfraction*nnodes);
2357 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2358 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2359 if (maxPMEnodes != minPMEnodes)
2361 fprintf(stdout, "- %d ", maxPMEnodes);
2363 fprintf(stdout, "PME-only nodes.\n Note that the automatic number of PME-only nodes and no separate PME nodes are always tested.\n");
2372 /* Get the commands we need to set up the runs from environment variables */
2373 get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2374 if (bBenchmark && repeats > 0)
2376 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun);
2379 /* Print some header info to file */
2381 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2383 fprintf(fp, "%s for Gromacs %s\n", ShortProgram(), GromacsVersion());
2386 fprintf(fp, "Number of nodes : %d\n", nnodes);
2387 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2388 if (strcmp(procstring[0], "none") != 0)
2390 fprintf(fp, "Passing # of nodes via : %s\n", procstring[0]);
2394 fprintf(fp, "Not setting number of nodes in system call\n");
2399 fprintf(fp, "Number of threads : %d\n", nnodes);
2402 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2403 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2404 fprintf(fp, "Benchmark steps : ");
2405 fprintf(fp, "%"GMX_PRId64, bench_nsteps);
2407 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2410 fprintf(fp, "Checkpoint time step : ");
2411 fprintf(fp, "%"GMX_PRId64, cpt_steps);
2414 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2416 if (new_sim_nsteps >= 0)
2419 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2420 fprintf(stderr, "%"GMX_PRId64, new_sim_nsteps+cpt_steps);
2421 fprintf(stderr, " steps.\n");
2422 fprintf(fp, "Simulation steps : ");
2423 fprintf(fp, "%"GMX_PRId64, new_sim_nsteps);
2428 fprintf(fp, "Repeats for each test : %d\n", repeats);
2431 if (npme_fixed >= -1)
2433 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2436 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2437 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2439 /* Allocate memory for the inputinfo struct: */
2441 info->nr_inputfiles = ntprs;
2442 for (i = 0; i < ntprs; i++)
2444 snew(info->rcoulomb, ntprs);
2445 snew(info->rvdw, ntprs);
2446 snew(info->rlist, ntprs);
2447 snew(info->rlistlong, ntprs);
2448 snew(info->nkx, ntprs);
2449 snew(info->nky, ntprs);
2450 snew(info->nkz, ntprs);
2451 snew(info->fsx, ntprs);
2452 snew(info->fsy, ntprs);
2453 snew(info->fsz, ntprs);
2455 /* Make alternative tpr files to test: */
2456 snew(tpr_names, ntprs);
2457 for (i = 0; i < ntprs; i++)
2459 snew(tpr_names[i], STRLEN);
2462 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2463 * different grids could be found. */
2464 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2465 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2467 /********************************************************************************/
2468 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2469 /********************************************************************************/
2470 snew(perfdata, ntprs);
2473 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2474 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2475 cmd_args_bench, fnm, NFILE, presteps, cpt_steps);
2477 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2479 /* Analyse the results and give a suggestion for optimal settings: */
2480 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2481 repeats, info, &best_tpr, &best_npme);
2483 /* Take the best-performing tpr file and enlarge nsteps to original value */
2484 if (bKeepTPR && !bOverwrite)
2486 simulation_tpr = opt2fn("-s", NFILE, fnm);
2490 simulation_tpr = opt2fn("-so", NFILE, fnm);
2491 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2492 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2495 /* Let's get rid of the temporary benchmark input files */
2496 for (i = 0; i < ntprs; i++)
2498 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2499 remove(tpr_names[i]);
2502 /* Now start the real simulation if the user requested it ... */
2503 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2504 cmd_args_launch, simulation_tpr, best_npme);
2508 /* ... or simply print the performance results to screen: */
2511 finalize(opt2fn("-p", NFILE, fnm));