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"
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 finalize(const char *fn_out)
177 fp = fopen(fn_out, "r");
178 fprintf(stdout, "\n\n");
180 while (fgets(buf, STRLEN-1, fp) != NULL)
182 fprintf(stdout, "%s", buf);
185 fprintf(stdout, "\n\n");
190 eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr
193 static int parse_logfile(const char *logfile, const char *errfile,
194 t_perf *perfdata, int test_nr, int presteps, gmx_int64_t cpt_steps,
198 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
199 const char matchstrdd[] = "Domain decomposition grid";
200 const char matchstrcr[] = "resetting all time and cycle counters";
201 const char matchstrbal[] = "Average PME mesh/force load:";
202 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";
203 const char errSIG[] = "signal, stopping at the next";
206 float dum1, dum2, dum3, dum4;
209 gmx_int64_t resetsteps = -1;
210 gmx_bool bFoundResetStr = FALSE;
211 gmx_bool bResetChecked = FALSE;
214 if (!gmx_fexist(logfile))
216 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
217 cleandata(perfdata, test_nr);
218 return eParselogNotFound;
221 fp = fopen(logfile, "r");
222 perfdata->PME_f_load[test_nr] = -1.0;
223 perfdata->guessPME = -1;
225 iFound = eFoundNothing;
228 iFound = eFoundDDStr; /* Skip some case statements */
231 while (fgets(line, STRLEN, fp) != NULL)
233 /* Remove leading spaces */
236 /* Check for TERM and INT signals from user: */
237 if (strstr(line, errSIG) != NULL)
240 cleandata(perfdata, test_nr);
241 return eParselogTerm;
244 /* Check whether cycle resetting worked */
245 if (presteps > 0 && !bFoundResetStr)
247 if (strstr(line, matchstrcr) != NULL)
249 sprintf(dumstring, "step %s", "%"GMX_SCNd64);
250 sscanf(line, dumstring, &resetsteps);
251 bFoundResetStr = TRUE;
252 if (resetsteps == presteps+cpt_steps)
254 bResetChecked = TRUE;
258 sprintf(dumstring, "%"GMX_PRId64, resetsteps);
259 sprintf(dumstring2, "%"GMX_PRId64, presteps+cpt_steps);
260 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
261 " though they were supposed to be reset at step %s!\n",
262 dumstring, dumstring2);
267 /* Look for strings that appear in a certain order in the log file: */
271 /* Look for domain decomp grid and separate PME nodes: */
272 if (str_starts(line, matchstrdd))
274 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME nodes %d",
275 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
276 if (perfdata->nPMEnodes == -1)
278 perfdata->guessPME = npme;
280 else if (perfdata->nPMEnodes != npme)
282 gmx_fatal(FARGS, "PME nodes from command line and output file are not identical");
284 iFound = eFoundDDStr;
286 /* Catch a few errors that might have occured: */
287 else if (str_starts(line, "There is no domain decomposition for"))
290 return eParselogNoDDGrid;
292 else if (str_starts(line, "The number of nodes you selected"))
295 return eParselogLargePrimeFactor;
297 else if (str_starts(line, "reading tpx file"))
300 return eParselogTPXVersion;
302 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
305 return eParselogNotParallel;
309 /* Look for PME mesh/force balance (not necessarily present, though) */
310 if (str_starts(line, matchstrbal))
312 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
314 /* Look for matchstring */
315 if (str_starts(line, matchstring))
317 iFound = eFoundAccountingStr;
320 case eFoundAccountingStr:
321 /* Already found matchstring - look for cycle data */
322 if (str_starts(line, "Total "))
324 sscanf(line, "Total %d %lf", &procs, &(perfdata->Gcycles[test_nr]));
325 iFound = eFoundCycleStr;
329 /* Already found cycle data - look for remaining performance info and return */
330 if (str_starts(line, "Performance:"))
332 ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
333 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
334 perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
336 if (bResetChecked || presteps == 0)
342 return eParselogResetProblem;
349 /* Close the log file */
352 /* Check why there is no performance data in the log file.
353 * Did a fatal errors occur? */
354 if (gmx_fexist(errfile))
356 fp = fopen(errfile, "r");
357 while (fgets(line, STRLEN, fp) != NULL)
359 if (str_starts(line, "Fatal error:") )
361 if (fgets(line, STRLEN, fp) != NULL)
363 fprintf(stderr, "\nWARNING: An error occured during this benchmark:\n"
367 cleandata(perfdata, test_nr);
368 return eParselogFatal;
375 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
378 /* Giving up ... we could not find out why there is no performance data in
380 fprintf(stdout, "No performance data in log file.\n");
381 cleandata(perfdata, test_nr);
383 return eParselogNoPerfData;
387 static gmx_bool analyze_data(
396 int *index_tpr, /* OUT: Nr of mdp file with best settings */
397 int *npme_optimal) /* OUT: Optimal number of PME nodes */
400 int line = 0, line_win = -1;
401 int k_win = -1, i_win = -1, winPME;
402 double s = 0.0; /* standard deviation */
405 char str_PME_f_load[13];
406 gmx_bool bCanUseOrigTPR;
407 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
413 fprintf(fp, "Summary of successful runs:\n");
414 fprintf(fp, "Line tpr PME nodes Gcycles Av. Std.dev. ns/day PME/f");
417 fprintf(fp, " DD grid");
423 for (k = 0; k < ntprs; k++)
425 for (i = 0; i < ntests; i++)
427 /* Select the right dataset: */
428 pd = &(perfdata[k][i]);
430 pd->Gcycles_Av = 0.0;
431 pd->PME_f_load_Av = 0.0;
432 pd->ns_per_day_Av = 0.0;
434 if (pd->nPMEnodes == -1)
436 sprintf(strbuf, "(%3d)", pd->guessPME);
440 sprintf(strbuf, " ");
443 /* Get the average run time of a setting */
444 for (j = 0; j < nrepeats; j++)
446 pd->Gcycles_Av += pd->Gcycles[j];
447 pd->PME_f_load_Av += pd->PME_f_load[j];
449 pd->Gcycles_Av /= nrepeats;
450 pd->PME_f_load_Av /= nrepeats;
452 for (j = 0; j < nrepeats; j++)
454 if (pd->ns_per_day[j] > 0.0)
456 pd->ns_per_day_Av += pd->ns_per_day[j];
460 /* Somehow the performance number was not aquired for this run,
461 * therefor set the average to some negative value: */
462 pd->ns_per_day_Av = -1.0f*nrepeats;
466 pd->ns_per_day_Av /= nrepeats;
469 if (pd->PME_f_load_Av > 0.0)
471 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
475 sprintf(str_PME_f_load, "%s", " - ");
479 /* We assume we had a successful run if both averages are positive */
480 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
482 /* Output statistics if repeats were done */
485 /* Calculate the standard deviation */
487 for (j = 0; j < nrepeats; j++)
489 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
494 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
495 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
496 pd->ns_per_day_Av, str_PME_f_load);
499 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
503 /* Store the index of the best run found so far in 'winner': */
504 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
517 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
522 winPME = perfdata[k_win][i_win].nPMEnodes;
526 /* We stuck to a fixed number of PME-only nodes */
527 sprintf(strbuf, "settings No. %d", k_win);
531 /* We have optimized the number of PME-only nodes */
534 sprintf(strbuf, "%s", "the automatic number of PME nodes");
538 sprintf(strbuf, "%d PME nodes", winPME);
541 fprintf(fp, "Best performance was achieved with %s", strbuf);
542 if ((nrepeats > 1) && (ntests > 1))
544 fprintf(fp, " (see line %d)", line_win);
548 /* Only mention settings if they were modified: */
549 bRefinedCoul = !is_equal(info->rcoulomb[k_win], info->rcoulomb[0]);
550 bRefinedVdW = !is_equal(info->rvdw[k_win], info->rvdw[0] );
551 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
552 info->nky[k_win] == info->nky[0] &&
553 info->nkz[k_win] == info->nkz[0]);
555 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
557 fprintf(fp, "Optimized PME settings:\n");
558 bCanUseOrigTPR = FALSE;
562 bCanUseOrigTPR = TRUE;
567 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
572 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
577 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],
578 info->nkx[0], info->nky[0], info->nkz[0]);
581 if (bCanUseOrigTPR && ntprs > 1)
583 fprintf(fp, "and original PME settings.\n");
588 /* Return the index of the mdp file that showed the highest performance
589 * and the optimal number of PME nodes */
591 *npme_optimal = winPME;
593 return bCanUseOrigTPR;
597 /* Get the commands we need to set up the runs from environment variables */
598 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
602 const char def_mpirun[] = "mpirun";
603 const char def_mdrun[] = "mdrun";
605 const char empty_mpirun[] = "";
607 /* Get the commands we need to set up the runs from environment variables */
610 if ( (cp = getenv("MPIRUN")) != NULL)
612 *cmd_mpirun = strdup(cp);
616 *cmd_mpirun = strdup(def_mpirun);
621 *cmd_mpirun = strdup(empty_mpirun);
624 if ( (cp = getenv("MDRUN" )) != NULL)
626 *cmd_mdrun = strdup(cp);
630 *cmd_mdrun = strdup(def_mdrun);
634 /* Check that the commands will run mdrun (perhaps via mpirun) by
635 * running a very quick test simulation. Requires MPI environment to
636 * be available if applicable. */
637 static void check_mdrun_works(gmx_bool bThreads,
638 const char *cmd_mpirun,
640 const char *cmd_mdrun)
642 char *command = NULL;
646 const char filename[] = "benchtest.log";
648 /* This string should always be identical to the one in copyrite.c,
649 * gmx_print_version_info() in the defined(GMX_MPI) section */
650 const char match_mpi[] = "MPI library: MPI";
651 const char match_mdrun[] = "Program: ";
652 gmx_bool bMdrun = FALSE;
653 gmx_bool bMPI = FALSE;
655 /* Run a small test to see whether mpirun + mdrun work */
656 fprintf(stdout, "Making sure that mdrun can be executed. ");
659 snew(command, strlen(cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
660 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
664 snew(command, strlen(cmd_mpirun) + strlen(cmd_np) + strlen(cmd_mdrun) + strlen(filename) + 50);
665 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
667 fprintf(stdout, "Trying '%s' ... ", command);
668 make_backup(filename);
669 gmx_system_call(command);
671 /* Check if we find the characteristic string in the output: */
672 if (!gmx_fexist(filename))
674 gmx_fatal(FARGS, "Output from test run could not be found.");
677 fp = fopen(filename, "r");
678 /* We need to scan the whole output file, since sometimes the queuing system
679 * also writes stuff to stdout/err */
682 cp = fgets(line, STRLEN, fp);
685 if (str_starts(line, match_mdrun) )
689 if (str_starts(line, match_mpi) )
701 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
703 "seems to have been compiled with MPI instead.",
711 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
713 "seems to have been compiled without MPI support.",
720 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
724 fprintf(stdout, "passed.\n");
732 static void launch_simulation(
733 gmx_bool bLaunch, /* Should the simulation be launched? */
734 FILE *fp, /* General log file */
735 gmx_bool bThreads, /* whether to use threads */
736 char *cmd_mpirun, /* Command for mpirun */
737 char *cmd_np, /* Switch for -np or -ntmpi or empty */
738 char *cmd_mdrun, /* Command for mdrun */
739 char *args_for_mdrun, /* Arguments for mdrun */
740 const char *simulation_tpr, /* This tpr will be simulated */
741 int nPMEnodes) /* Number of PME nodes to use */
746 /* Make enough space for the system call command,
747 * (100 extra chars for -npme ... etc. options should suffice): */
748 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+100);
750 /* Note that the -passall options requires args_for_mdrun to be at the end
751 * of the command line string */
754 sprintf(command, "%s%s-npme %d -s %s %s",
755 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun);
759 sprintf(command, "%s%s%s -npme %d -s %s %s",
760 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun);
763 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
767 /* Now the real thing! */
770 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
773 gmx_system_call(command);
778 static void modify_PMEsettings(
779 gmx_int64_t simsteps, /* Set this value as number of time steps */
780 gmx_int64_t init_step, /* Set this value as init_step */
781 const char *fn_best_tpr, /* tpr file with the best performance */
782 const char *fn_sim_tpr) /* name of tpr file to be launched */
790 read_tpx_state(fn_best_tpr, ir, &state, NULL, &mtop);
792 /* Reset nsteps and init_step to the value of the input .tpr file */
793 ir->nsteps = simsteps;
794 ir->init_step = init_step;
796 /* Write the tpr file which will be launched */
797 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%"GMX_PRId64);
798 fprintf(stdout, buf, ir->nsteps);
800 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
806 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
808 /* Make additional TPR files with more computational load for the
809 * direct space processors: */
810 static void make_benchmark_tprs(
811 const char *fn_sim_tpr, /* READ : User-provided tpr file */
812 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
813 gmx_int64_t benchsteps, /* Number of time steps for benchmark runs */
814 gmx_int64_t statesteps, /* Step counter in checkpoint file */
815 real rmin, /* Minimal Coulomb radius */
816 real rmax, /* Maximal Coulomb radius */
817 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
818 int *ntprs, /* No. of TPRs to write, each with a different
819 rcoulomb and fourierspacing */
820 t_inputinfo *info, /* Contains information about mdp file options */
821 FILE *fp) /* Write the output here */
827 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
830 gmx_bool bNote = FALSE;
831 real add; /* Add this to rcoul for the next test */
832 real fac = 1.0; /* Scaling factor for Coulomb radius */
833 real fourierspacing; /* Basic fourierspacing from tpr */
836 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
837 *ntprs > 1 ? "s" : "", "%"GMX_PRId64, benchsteps > 1 ? "s" : "");
838 fprintf(stdout, buf, benchsteps);
841 sprintf(buf, " (adding %s steps from checkpoint file)", "%"GMX_PRId64);
842 fprintf(stdout, buf, statesteps);
843 benchsteps += statesteps;
845 fprintf(stdout, ".\n");
849 read_tpx_state(fn_sim_tpr, ir, &state, NULL, &mtop);
851 /* Check if some kind of PME was chosen */
852 if (EEL_PME(ir->coulombtype) == FALSE)
854 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
858 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
859 if ( (ir->cutoff_scheme != ecutsVERLET) &&
860 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
862 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
863 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
865 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
866 else if (ir->rcoulomb > ir->rlist)
868 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
869 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
872 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
874 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
878 /* Reduce the number of steps for the benchmarks */
879 info->orig_sim_steps = ir->nsteps;
880 ir->nsteps = benchsteps;
881 /* We must not use init_step from the input tpr file for the benchmarks */
882 info->orig_init_step = ir->init_step;
885 /* For PME-switch potentials, keep the radial distance of the buffer region */
886 nlist_buffer = ir->rlist - ir->rcoulomb;
888 /* Determine length of triclinic box vectors */
889 for (d = 0; d < DIM; d++)
892 for (i = 0; i < DIM; i++)
894 box_size[d] += state.box[d][i]*state.box[d][i];
896 box_size[d] = sqrt(box_size[d]);
899 if (ir->fourier_spacing > 0)
901 info->fsx[0] = ir->fourier_spacing;
902 info->fsy[0] = ir->fourier_spacing;
903 info->fsz[0] = ir->fourier_spacing;
907 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
908 info->fsx[0] = box_size[XX]/ir->nkx;
909 info->fsy[0] = box_size[YY]/ir->nky;
910 info->fsz[0] = box_size[ZZ]/ir->nkz;
913 /* If no value for the fourierspacing was provided on the command line, we
914 * use the reconstruction from the tpr file */
915 if (ir->fourier_spacing > 0)
917 /* Use the spacing from the tpr */
918 fourierspacing = ir->fourier_spacing;
922 /* Use the maximum observed spacing */
923 fourierspacing = max(max(info->fsx[0], info->fsy[0]), info->fsz[0]);
926 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
928 /* For performance comparisons the number of particles is useful to have */
929 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
931 /* Print information about settings of which some are potentially modified: */
932 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
933 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
934 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
935 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
936 if (ir_vdw_switched(ir))
938 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
940 if (EPME_SWITCHED(ir->coulombtype))
942 fprintf(fp, " rlist : %f nm\n", ir->rlist);
944 if (ir->rlistlong != max_cutoff(ir->rvdw, ir->rcoulomb))
946 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
949 /* Print a descriptive line about the tpr settings tested */
950 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
951 fprintf(fp, " No. scaling rcoulomb");
952 fprintf(fp, " nkx nky nkz");
953 fprintf(fp, " spacing");
954 if (evdwCUT == ir->vdwtype)
956 fprintf(fp, " rvdw");
958 if (EPME_SWITCHED(ir->coulombtype))
960 fprintf(fp, " rlist");
962 if (ir->rlistlong != max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb)) )
964 fprintf(fp, " rlistlong");
966 fprintf(fp, " tpr file\n");
968 /* Loop to create the requested number of tpr input files */
969 for (j = 0; j < *ntprs; j++)
971 /* The first .tpr is the provided one, just need to modify nsteps,
972 * so skip the following block */
975 /* Determine which Coulomb radii rc to use in the benchmarks */
976 add = (rmax-rmin)/(*ntprs-1);
977 if (is_equal(rmin, info->rcoulomb[0]))
979 ir->rcoulomb = rmin + j*add;
981 else if (is_equal(rmax, info->rcoulomb[0]))
983 ir->rcoulomb = rmin + (j-1)*add;
987 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
988 add = (rmax-rmin)/(*ntprs-2);
989 ir->rcoulomb = rmin + (j-1)*add;
992 /* Determine the scaling factor fac */
993 fac = ir->rcoulomb/info->rcoulomb[0];
995 /* Scale the Fourier grid spacing */
996 ir->nkx = ir->nky = ir->nkz = 0;
997 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
999 /* Adjust other radii since various conditions need to be fulfilled */
1000 if (eelPME == ir->coulombtype)
1002 /* plain PME, rcoulomb must be equal to rlist */
1003 ir->rlist = ir->rcoulomb;
1007 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1008 ir->rlist = ir->rcoulomb + nlist_buffer;
1011 if (bScaleRvdw && evdwCUT == ir->vdwtype)
1013 if (ecutsVERLET == ir->cutoff_scheme)
1015 /* With Verlet, the van der Waals radius must always equal the Coulomb radius */
1016 ir->rvdw = ir->rcoulomb;
1020 /* For vdw cutoff, rvdw >= rlist */
1021 ir->rvdw = max(info->rvdw[0], ir->rlist);
1025 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1027 } /* end of "if (j != 0)" */
1029 /* for j==0: Save the original settings
1030 * for j >0: Save modified radii and Fourier grids */
1031 info->rcoulomb[j] = ir->rcoulomb;
1032 info->rvdw[j] = ir->rvdw;
1033 info->nkx[j] = ir->nkx;
1034 info->nky[j] = ir->nky;
1035 info->nkz[j] = ir->nkz;
1036 info->rlist[j] = ir->rlist;
1037 info->rlistlong[j] = ir->rlistlong;
1038 info->fsx[j] = fac*fourierspacing;
1039 info->fsy[j] = fac*fourierspacing;
1040 info->fsz[j] = fac*fourierspacing;
1042 /* Write the benchmark tpr file */
1043 strncpy(fn_bench_tprs[j], fn_sim_tpr, strlen(fn_sim_tpr)-strlen(".tpr"));
1044 sprintf(buf, "_bench%.2d.tpr", j);
1045 strcat(fn_bench_tprs[j], buf);
1046 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1047 fprintf(stdout, "%"GMX_PRId64, ir->nsteps);
1050 fprintf(stdout, ", scaling factor %f\n", fac);
1054 fprintf(stdout, ", unmodified settings\n");
1057 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1059 /* Write information about modified tpr settings to log file */
1060 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1061 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1062 fprintf(fp, " %9f ", info->fsx[j]);
1063 if (evdwCUT == ir->vdwtype)
1065 fprintf(fp, "%10f", ir->rvdw);
1067 if (EPME_SWITCHED(ir->coulombtype))
1069 fprintf(fp, "%10f", ir->rlist);
1071 if (info->rlistlong[0] != max_cutoff(info->rlist[0], max_cutoff(info->rvdw[0], info->rcoulomb[0])) )
1073 fprintf(fp, "%10f", ir->rlistlong);
1075 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1077 /* Make it clear to the user that some additional settings were modified */
1078 if (!is_equal(ir->rvdw, info->rvdw[0])
1079 || !is_equal(ir->rlistlong, info->rlistlong[0]) )
1086 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1087 "other input settings were also changed (see table above).\n"
1088 "Please check if the modified settings are appropriate.\n");
1096 /* Rename the files we want to keep to some meaningful filename and
1097 * delete the rest */
1098 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1099 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1101 char numstring[STRLEN];
1102 char newfilename[STRLEN];
1103 const char *fn = NULL;
1108 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1110 for (i = 0; i < nfile; i++)
1112 opt = (char *)fnm[i].opt;
1113 if (strcmp(opt, "-p") == 0)
1115 /* do nothing; keep this file */
1118 else if (strcmp(opt, "-bg") == 0)
1120 /* Give the log file a nice name so one can later see which parameters were used */
1121 numstring[0] = '\0';
1124 sprintf(numstring, "_%d", nr);
1126 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1127 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1129 fprintf(stdout, "renaming log file to %s\n", newfilename);
1130 make_backup(newfilename);
1131 rename(opt2fn("-bg", nfile, fnm), newfilename);
1134 else if (strcmp(opt, "-err") == 0)
1136 /* This file contains the output of stderr. We want to keep it in
1137 * cases where there have been problems. */
1138 fn = opt2fn(opt, nfile, fnm);
1139 numstring[0] = '\0';
1142 sprintf(numstring, "_%d", nr);
1144 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1149 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1150 make_backup(newfilename);
1151 rename(fn, newfilename);
1155 fprintf(stdout, "Deleting %s\n", fn);
1160 /* Delete the files which are created for each benchmark run: (options -b*) */
1161 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1163 fn = opt2fn(opt, nfile, fnm);
1166 fprintf(stdout, "Deleting %s\n", fn);
1175 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1178 /* Create a list of numbers of PME nodes to test */
1179 static void make_npme_list(
1180 const char *npmevalues_opt, /* Make a complete list with all
1181 * possibilities or a short list that keeps only
1182 * reasonable numbers of PME nodes */
1183 int *nentries, /* Number of entries we put in the nPMEnodes list */
1184 int *nPMEnodes[], /* Each entry contains the value for -npme */
1185 int nnodes, /* Total number of nodes to do the tests on */
1186 int minPMEnodes, /* Minimum number of PME nodes */
1187 int maxPMEnodes) /* Maximum number of PME nodes */
1190 int min_factor = 1; /* We request that npp and npme have this minimal
1191 * largest common factor (depends on npp) */
1192 int nlistmax; /* Max. list size */
1193 int nlist; /* Actual number of entries in list */
1197 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1198 if (0 == strcmp(npmevalues_opt, "all") )
1202 else if (0 == strcmp(npmevalues_opt, "subset") )
1204 eNPME = eNpmeSubset;
1206 else /* "auto" or "range" */
1212 else if (nnodes < 128)
1214 eNPME = eNpmeReduced;
1218 eNPME = eNpmeSubset;
1222 /* Calculate how many entries we could possibly have (in case of -npme all) */
1225 nlistmax = maxPMEnodes - minPMEnodes + 3;
1226 if (0 == minPMEnodes)
1236 /* Now make the actual list which is at most of size nlist */
1237 snew(*nPMEnodes, nlistmax);
1238 nlist = 0; /* start counting again, now the real entries in the list */
1239 for (i = 0; i < nlistmax - 2; i++)
1241 npme = maxPMEnodes - i;
1252 /* For 2d PME we want a common largest factor of at least the cube
1253 * root of the number of PP nodes */
1254 min_factor = (int) pow(npp, 1.0/3.0);
1257 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1260 if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1262 (*nPMEnodes)[nlist] = npme;
1266 /* We always test 0 PME nodes and the automatic number */
1267 *nentries = nlist + 2;
1268 (*nPMEnodes)[nlist ] = 0;
1269 (*nPMEnodes)[nlist+1] = -1;
1271 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1272 for (i = 0; i < *nentries-1; i++)
1274 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1276 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1280 /* Allocate memory to store the performance data */
1281 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1286 for (k = 0; k < ntprs; k++)
1288 snew(perfdata[k], datasets);
1289 for (i = 0; i < datasets; i++)
1291 for (j = 0; j < repeats; j++)
1293 snew(perfdata[k][i].Gcycles, repeats);
1294 snew(perfdata[k][i].ns_per_day, repeats);
1295 snew(perfdata[k][i].PME_f_load, repeats);
1302 /* Check for errors on mdrun -h */
1303 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp)
1305 char *command, *msg;
1309 snew(command, length + 15);
1310 snew(msg, length + 500);
1312 fprintf(stdout, "Making sure the benchmarks can be executed ...\n");
1313 /* FIXME: mdrun -h no longer actually does anything useful.
1314 * It unconditionally prints the help, ignoring all other options. */
1315 sprintf(command, "%s-h -quiet", mdrun_cmd_line);
1316 ret = gmx_system_call(command);
1320 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1321 * get the error message from mdrun itself */
1322 sprintf(msg, "Cannot run the benchmark simulations! Please check the error message of\n"
1323 "mdrun for the source of the problem. Did you provide a command line\n"
1324 "argument that neither g_tune_pme nor mdrun understands? Offending command:\n"
1325 "\n%s\n\n", command);
1327 fprintf(stderr, "%s", msg);
1329 fprintf(fp, "%s", msg);
1339 static void do_the_tests(
1340 FILE *fp, /* General g_tune_pme output file */
1341 char **tpr_names, /* Filenames of the input files to test */
1342 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1343 int minPMEnodes, /* Min fraction of nodes to use for PME */
1344 int npme_fixed, /* If >= -1, test fixed number of PME
1346 const char *npmevalues_opt, /* Which -npme values should be tested */
1347 t_perf **perfdata, /* Here the performace data is stored */
1348 int *pmeentries, /* Entries in the nPMEnodes list */
1349 int repeats, /* Repeat each test this often */
1350 int nnodes, /* Total number of nodes = nPP + nPME */
1351 int nr_tprs, /* Total number of tpr files to test */
1352 gmx_bool bThreads, /* Threads or MPI? */
1353 char *cmd_mpirun, /* mpirun command string */
1354 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1355 char *cmd_mdrun, /* mdrun command string */
1356 char *cmd_args_bench, /* arguments for mdrun in a string */
1357 const t_filenm *fnm, /* List of filenames from command line */
1358 int nfile, /* Number of files specified on the cmdl. */
1359 int presteps, /* DLB equilibration steps, is checked */
1360 gmx_int64_t cpt_steps) /* Time step counter in the checkpoint */
1362 int i, nr, k, ret, count = 0, totaltests;
1363 int *nPMEnodes = NULL;
1366 char *command, *cmd_stub;
1368 gmx_bool bResetProblem = FALSE;
1369 gmx_bool bFirst = TRUE;
1372 /* This string array corresponds to the eParselog enum type at the start
1374 const char* ParseLog[] = {
1376 "Logfile not found!",
1377 "No timings, logfile truncated?",
1378 "Run was terminated.",
1379 "Counters were not reset properly.",
1380 "No DD grid found for these settings.",
1381 "TPX version conflict!",
1382 "mdrun was not started in parallel!",
1383 "Number of PP nodes has a prime factor that is too large.",
1386 char str_PME_f_load[13];
1389 /* Allocate space for the mdrun command line. 100 extra characters should
1390 be more than enough for the -npme etcetera arguments */
1391 cmdline_length = strlen(cmd_mpirun)
1394 + strlen(cmd_args_bench)
1395 + strlen(tpr_names[0]) + 100;
1396 snew(command, cmdline_length);
1397 snew(cmd_stub, cmdline_length);
1399 /* Construct the part of the command line that stays the same for all tests: */
1402 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1406 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1409 /* Create a list of numbers of PME nodes to test */
1410 if (npme_fixed < -1)
1412 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1413 nnodes, minPMEnodes, maxPMEnodes);
1419 nPMEnodes[0] = npme_fixed;
1420 fprintf(stderr, "Will use a fixed number of %d PME-only nodes.\n", nPMEnodes[0]);
1425 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1427 finalize(opt2fn("-p", nfile, fnm));
1431 /* Allocate one dataset for each tpr input file: */
1432 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1434 /*****************************************/
1435 /* Main loop over all tpr files to test: */
1436 /*****************************************/
1437 totaltests = nr_tprs*(*pmeentries)*repeats;
1438 for (k = 0; k < nr_tprs; k++)
1440 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1441 fprintf(fp, "PME nodes Gcycles ns/day PME/f Remark\n");
1442 /* Loop over various numbers of PME nodes: */
1443 for (i = 0; i < *pmeentries; i++)
1445 pd = &perfdata[k][i];
1447 /* Loop over the repeats for each scenario: */
1448 for (nr = 0; nr < repeats; nr++)
1450 pd->nPMEnodes = nPMEnodes[i];
1452 /* Add -npme and -s to the command line and save it. Note that
1453 * the -passall (if set) options requires cmd_args_bench to be
1454 * at the end of the command line string */
1455 snew(pd->mdrun_cmd_line, cmdline_length);
1456 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s",
1457 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench);
1459 /* To prevent that all benchmarks fail due to a show-stopper argument
1460 * on the mdrun command line, we make a quick check with mdrun -h first */
1463 make_sure_it_runs(pd->mdrun_cmd_line, cmdline_length, fp);
1467 /* Do a benchmark simulation: */
1470 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1476 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1477 (100.0*count)/totaltests,
1478 k+1, nr_tprs, i+1, *pmeentries, buf);
1479 make_backup(opt2fn("-err", nfile, fnm));
1480 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1481 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1482 gmx_system_call(command);
1484 /* Collect the performance data from the log file; also check stderr
1485 * for fatal errors */
1486 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1487 pd, nr, presteps, cpt_steps, nnodes);
1488 if ((presteps > 0) && (ret == eParselogResetProblem))
1490 bResetProblem = TRUE;
1493 if (-1 == pd->nPMEnodes)
1495 sprintf(buf, "(%3d)", pd->guessPME);
1503 if (pd->PME_f_load[nr] > 0.0)
1505 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1509 sprintf(str_PME_f_load, "%s", " - ");
1512 /* Write the data we got to disk */
1513 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1514 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1515 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1517 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1523 /* Do some cleaning up and delete the files we do not need any more */
1524 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1526 /* If the first run with this number of processors already failed, do not try again: */
1527 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1529 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1530 count += repeats-(nr+1);
1533 } /* end of repeats loop */
1534 } /* end of -npme loop */
1535 } /* end of tpr file loop */
1540 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1548 static void check_input(
1555 real maxPMEfraction,
1556 real minPMEfraction,
1558 gmx_int64_t bench_nsteps,
1559 const t_filenm *fnm,
1569 /* Make sure the input file exists */
1570 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1572 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1575 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1576 if ( (0 == strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1578 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1579 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1582 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1585 gmx_fatal(FARGS, "Number of repeats < 0!");
1588 /* Check number of nodes */
1591 gmx_fatal(FARGS, "Number of nodes/threads must be a positive integer.");
1594 /* Automatically choose -ntpr if not set */
1604 /* Set a reasonable scaling factor for rcoulomb */
1607 *rmax = rcoulomb * 1.2;
1610 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1616 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1620 /* Make shure that rmin <= rcoulomb <= rmax */
1629 if (!(*rmin <= *rmax) )
1631 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1632 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1634 /* Add test scenarios if rmin or rmax were set */
1637 if (!is_equal(*rmin, rcoulomb) && (*ntprs == 1) )
1640 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1643 if (!is_equal(*rmax, rcoulomb) && (*ntprs == 1) )
1646 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1651 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1652 if (!is_equal(*rmax, rcoulomb) || !is_equal(*rmin, rcoulomb) )
1654 *ntprs = max(*ntprs, 2);
1657 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1658 if (!is_equal(*rmax, rcoulomb) && !is_equal(*rmin, rcoulomb) )
1660 *ntprs = max(*ntprs, 3);
1665 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1670 if (is_equal(*rmin, rcoulomb) && is_equal(rcoulomb, *rmax)) /* We have just a single rc */
1672 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1673 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1674 "with correspondingly adjusted PME grid settings\n");
1679 /* Check whether max and min fraction are within required values */
1680 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1682 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1684 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1686 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1688 if (maxPMEfraction < minPMEfraction)
1690 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1693 /* Check whether the number of steps - if it was set - has a reasonable value */
1694 if (bench_nsteps < 0)
1696 gmx_fatal(FARGS, "Number of steps must be positive.");
1699 if (bench_nsteps > 10000 || bench_nsteps < 100)
1701 fprintf(stderr, "WARNING: steps=");
1702 fprintf(stderr, "%"GMX_PRId64, bench_nsteps);
1703 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1708 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1711 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1714 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1716 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1720 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1721 * only. We need to check whether the requested number of PME-only nodes
1723 if (npme_fixed > -1)
1725 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1726 if (2*npme_fixed > nnodes)
1728 gmx_fatal(FARGS, "Cannot have more than %d PME-only nodes for a total of %d nodes (you chose %d).\n",
1729 nnodes/2, nnodes, npme_fixed);
1731 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1733 fprintf(stderr, "WARNING: Only %g percent of the nodes are assigned as PME-only nodes.\n",
1734 100.0*((real)npme_fixed / (real)nnodes));
1736 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1738 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1739 " fixed number of PME-only nodes is requested with -fix.\n");
1745 /* Returns TRUE when "opt" is needed at launch time */
1746 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1748 if (0 == strncmp(opt, "-swap", 5))
1753 /* Apart from the input .tpr and the output log files we need all options that
1754 * were set on the command line and that do not start with -b */
1755 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2)
1756 || 0 == strncmp(opt, "-err", 4) || 0 == strncmp(opt, "-p", 2) )
1765 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1766 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1768 /* Apart from the input .tpr, all files starting with "-b" are for
1769 * _b_enchmark files exclusively */
1770 if (0 == strncmp(opt, "-s", 2))
1775 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2))
1777 if (!bOptional || bSet)
1794 if (bSet) /* These are additional input files like -cpi -ei */
1807 /* Adds 'buf' to 'str' */
1808 static void add_to_string(char **str, char *buf)
1813 len = strlen(*str) + strlen(buf) + 1;
1819 /* Create the command line for the benchmark as well as for the real run */
1820 static void create_command_line_snippets(
1821 gmx_bool bAppendFiles,
1822 gmx_bool bKeepAndNumCPT,
1823 gmx_bool bResetHWay,
1827 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1828 char *cmd_args_launch[], /* command line arguments for simulation run */
1829 char extra_args[]) /* Add this to the end of the command line */
1834 char strbuf[STRLEN];
1837 /* strlen needs at least '\0' as a string: */
1838 snew(*cmd_args_bench, 1);
1839 snew(*cmd_args_launch, 1);
1840 *cmd_args_launch[0] = '\0';
1841 *cmd_args_bench[0] = '\0';
1844 /*******************************************/
1845 /* 1. Process other command line arguments */
1846 /*******************************************/
1849 /* Add equilibration steps to benchmark options */
1850 sprintf(strbuf, "-resetstep %d ", presteps);
1851 add_to_string(cmd_args_bench, strbuf);
1853 /* These switches take effect only at launch time */
1854 if (FALSE == bAppendFiles)
1856 add_to_string(cmd_args_launch, "-noappend ");
1860 add_to_string(cmd_args_launch, "-cpnum ");
1864 add_to_string(cmd_args_launch, "-resethway ");
1867 /********************/
1868 /* 2. Process files */
1869 /********************/
1870 for (i = 0; i < nfile; i++)
1872 opt = (char *)fnm[i].opt;
1873 name = opt2fn(opt, nfile, fnm);
1875 /* Strbuf contains the options, now let's sort out where we need that */
1876 sprintf(strbuf, "%s %s ", opt, name);
1878 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1880 /* All options starting with -b* need the 'b' removed,
1881 * therefore overwrite strbuf */
1882 if (0 == strncmp(opt, "-b", 2))
1884 sprintf(strbuf, "-%s %s ", &opt[2], name);
1887 add_to_string(cmd_args_bench, strbuf);
1890 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
1892 add_to_string(cmd_args_launch, strbuf);
1896 add_to_string(cmd_args_bench, extra_args);
1897 add_to_string(cmd_args_launch, extra_args);
1901 /* Set option opt */
1902 static void setopt(const char *opt, int nfile, t_filenm fnm[])
1906 for (i = 0; (i < nfile); i++)
1908 if (strcmp(opt, fnm[i].opt) == 0)
1910 fnm[i].flag |= ffSET;
1916 /* This routine inspects the tpr file and ...
1917 * 1. checks for output files that get triggered by a tpr option. These output
1918 * files are marked as 'set' to allow for a proper cleanup after each
1920 * 2. returns the PME:PP load ratio
1921 * 3. returns rcoulomb from the tpr */
1922 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1924 gmx_bool bPull; /* Is pulling requested in .tpr file? */
1925 gmx_bool bTpi; /* Is test particle insertion requested? */
1926 gmx_bool bFree; /* Is a free energy simulation requested? */
1927 gmx_bool bNM; /* Is a normal mode analysis requested? */
1928 gmx_bool bSwap; /* Is water/ion position swapping requested? */
1934 /* Check tpr file for options that trigger extra output files */
1935 read_tpx_state(opt2fn("-s", nfile, fnm), &ir, &state, NULL, &mtop);
1936 bPull = (epullNO != ir.ePull);
1937 bFree = (efepNO != ir.efep );
1938 bNM = (eiNM == ir.eI );
1939 bSwap = (eswapNO != ir.eSwapCoords);
1940 bTpi = EI_TPI(ir.eI);
1942 /* Set these output files on the tuning command-line */
1945 setopt("-pf", nfile, fnm);
1946 setopt("-px", nfile, fnm);
1950 setopt("-dhdl", nfile, fnm);
1954 setopt("-tpi", nfile, fnm);
1955 setopt("-tpid", nfile, fnm);
1959 setopt("-mtx", nfile, fnm);
1963 setopt("-swap", nfile, fnm);
1966 *rcoulomb = ir.rcoulomb;
1968 /* Return the estimate for the number of PME nodes */
1969 return pme_load_estimate(&mtop, &ir, state.box);
1973 static void couple_files_options(int nfile, t_filenm fnm[])
1976 gmx_bool bSet, bBench;
1981 for (i = 0; i < nfile; i++)
1983 opt = (char *)fnm[i].opt;
1984 bSet = ((fnm[i].flag & ffSET) != 0);
1985 bBench = (0 == strncmp(opt, "-b", 2));
1987 /* Check optional files */
1988 /* If e.g. -eo is set, then -beo also needs to be set */
1989 if (is_optional(&fnm[i]) && bSet && !bBench)
1991 sprintf(buf, "-b%s", &opt[1]);
1992 setopt(buf, nfile, fnm);
1994 /* If -beo is set, then -eo also needs to be! */
1995 if (is_optional(&fnm[i]) && bSet && bBench)
1997 sprintf(buf, "-%s", &opt[2]);
1998 setopt(buf, nfile, fnm);
2004 #define BENCHSTEPS (1000)
2006 int gmx_tune_pme(int argc, char *argv[])
2008 const char *desc[] = {
2009 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of processors/threads, [THISMODULE] systematically",
2010 "times [gmx-mdrun] with various numbers of PME-only nodes and determines",
2011 "which setting is fastest. It will also test whether performance can",
2012 "be enhanced by shifting load from the reciprocal to the real space",
2013 "part of the Ewald sum. ",
2014 "Simply pass your [TT].tpr[tt] file to [THISMODULE] together with other options",
2015 "for [gmx-mdrun] as needed.[PAR]",
2016 "Which executables are used can be set in the environment variables",
2017 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
2018 "will be used as defaults. Note that for certain MPI frameworks you",
2019 "need to provide a machine- or hostfile. This can also be passed",
2020 "via the MPIRUN variable, e.g.[PAR]",
2021 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
2022 "Please call [THISMODULE] with the normal options you would pass to",
2023 "[gmx-mdrun] and add [TT]-np[tt] for the number of processors to perform the",
2024 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2025 "to repeat each test several times to get better statistics. [PAR]",
2026 "[THISMODULE] can test various real space / reciprocal space workloads",
2027 "for you. With [TT]-ntpr[tt] you control how many extra [TT].tpr[tt] files will be",
2028 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2029 "Typically, the first test (number 0) will be with the settings from the input",
2030 "[TT].tpr[tt] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2031 "specified by [TT]-rmax[tt] with a somwhat smaller PME grid at the same time. ",
2032 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2033 "The remaining [TT].tpr[tt] files will have equally-spaced Coulomb radii (and Fourier "
2034 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2035 "if you just seek the optimal number of PME-only nodes; in that case",
2036 "your input [TT].tpr[tt] file will remain unchanged.[PAR]",
2037 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2038 "MD systems. The dynamic load balancing needs about 100 time steps",
2039 "to adapt to local load imbalances, therefore the time step counters",
2040 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2041 "for a higher accuarcy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2042 "From the 'DD' load imbalance entries in the md.log output file you",
2043 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2044 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2045 "After calling [gmx-mdrun] several times, detailed performance information",
2046 "is available in the output file [TT]perf.out[tt].",
2047 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2048 "(options [TT]-b[tt]*), these will be automatically deleted after each test.[PAR]",
2049 "If you want the simulation to be started automatically with the",
2050 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2055 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2056 real maxPMEfraction = 0.50;
2057 real minPMEfraction = 0.25;
2058 int maxPMEnodes, minPMEnodes;
2059 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2060 float guessPMEnodes;
2061 int npme_fixed = -2; /* If >= -1, use only this number
2062 * of PME-only nodes */
2064 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2065 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2066 gmx_bool bScaleRvdw = TRUE;
2067 gmx_int64_t bench_nsteps = BENCHSTEPS;
2068 gmx_int64_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2069 gmx_int64_t cpt_steps = 0; /* Step counter in .cpt input file */
2070 int presteps = 100; /* Do a full cycle reset after presteps steps */
2071 gmx_bool bOverwrite = FALSE, bKeepTPR;
2072 gmx_bool bLaunch = FALSE;
2073 char *ExtraArgs = NULL;
2074 char **tpr_names = NULL;
2075 const char *simulation_tpr = NULL;
2076 int best_npme, best_tpr;
2077 int sim_part = 1; /* For benchmarks with checkpoint files */
2080 /* Default program names if nothing else is found */
2081 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2082 char *cmd_args_bench, *cmd_args_launch;
2083 char *cmd_np = NULL;
2085 t_perf **perfdata = NULL;
2091 /* Print out how long the tuning took */
2094 static t_filenm fnm[] = {
2096 { efOUT, "-p", "perf", ffWRITE },
2097 { efLOG, "-err", "bencherr", ffWRITE },
2098 { efTPX, "-so", "tuned", ffWRITE },
2100 { efTPX, NULL, NULL, ffREAD },
2101 { efTRN, "-o", NULL, ffWRITE },
2102 { efCOMPRESSED, "-x", NULL, ffOPTWR },
2103 { efCPT, "-cpi", NULL, ffOPTRD },
2104 { efCPT, "-cpo", NULL, ffOPTWR },
2105 { efSTO, "-c", "confout", ffWRITE },
2106 { efEDR, "-e", "ener", ffWRITE },
2107 { efLOG, "-g", "md", ffWRITE },
2108 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2109 { efXVG, "-field", "field", ffOPTWR },
2110 { efXVG, "-table", "table", ffOPTRD },
2111 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2112 { efXVG, "-tablep", "tablep", ffOPTRD },
2113 { efXVG, "-tableb", "table", ffOPTRD },
2114 { efTRX, "-rerun", "rerun", ffOPTRD },
2115 { efXVG, "-tpi", "tpi", ffOPTWR },
2116 { efXVG, "-tpid", "tpidist", ffOPTWR },
2117 { efEDI, "-ei", "sam", ffOPTRD },
2118 { efXVG, "-eo", "edsam", ffOPTWR },
2119 { efXVG, "-devout", "deviatie", ffOPTWR },
2120 { efXVG, "-runav", "runaver", ffOPTWR },
2121 { efXVG, "-px", "pullx", ffOPTWR },
2122 { efXVG, "-pf", "pullf", ffOPTWR },
2123 { efXVG, "-ro", "rotation", ffOPTWR },
2124 { efLOG, "-ra", "rotangles", ffOPTWR },
2125 { efLOG, "-rs", "rotslabs", ffOPTWR },
2126 { efLOG, "-rt", "rottorque", ffOPTWR },
2127 { efMTX, "-mtx", "nm", ffOPTWR },
2128 { efNDX, "-dn", "dipole", ffOPTWR },
2129 { efXVG, "-swap", "swapions", ffOPTWR },
2130 /* Output files that are deleted after each benchmark run */
2131 { efTRN, "-bo", "bench", ffWRITE },
2132 { efXTC, "-bx", "bench", ffWRITE },
2133 { efCPT, "-bcpo", "bench", ffWRITE },
2134 { efSTO, "-bc", "bench", ffWRITE },
2135 { efEDR, "-be", "bench", ffWRITE },
2136 { efLOG, "-bg", "bench", ffWRITE },
2137 { efXVG, "-beo", "benchedo", ffOPTWR },
2138 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2139 { efXVG, "-bfield", "benchfld", ffOPTWR },
2140 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2141 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2142 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2143 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2144 { efXVG, "-bpx", "benchpx", ffOPTWR },
2145 { efXVG, "-bpf", "benchpf", ffOPTWR },
2146 { efXVG, "-bro", "benchrot", ffOPTWR },
2147 { efLOG, "-bra", "benchrota", ffOPTWR },
2148 { efLOG, "-brs", "benchrots", ffOPTWR },
2149 { efLOG, "-brt", "benchrott", ffOPTWR },
2150 { efMTX, "-bmtx", "benchn", ffOPTWR },
2151 { efNDX, "-bdn", "bench", ffOPTWR },
2152 { efXVG, "-bswap", "benchswp", ffOPTWR }
2155 gmx_bool bThreads = FALSE;
2159 const char *procstring[] =
2160 { NULL, "-np", "-n", "none", NULL };
2161 const char *npmevalues_opt[] =
2162 { NULL, "auto", "all", "subset", NULL };
2164 gmx_bool bAppendFiles = TRUE;
2165 gmx_bool bKeepAndNumCPT = FALSE;
2166 gmx_bool bResetCountersHalfWay = FALSE;
2167 gmx_bool bBenchmark = TRUE;
2169 output_env_t oenv = NULL;
2172 /***********************/
2173 /* g_tune_pme options: */
2174 /***********************/
2175 { "-np", FALSE, etINT, {&nnodes},
2176 "Number of nodes to run the tests on (must be > 2 for separate PME nodes)" },
2177 { "-npstring", FALSE, etENUM, {procstring},
2178 "Specify the number of processors to [TT]$MPIRUN[tt] using this string"},
2179 { "-ntmpi", FALSE, etINT, {&nthreads},
2180 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2181 { "-r", FALSE, etINT, {&repeats},
2182 "Repeat each test this often" },
2183 { "-max", FALSE, etREAL, {&maxPMEfraction},
2184 "Max fraction of PME nodes to test with" },
2185 { "-min", FALSE, etREAL, {&minPMEfraction},
2186 "Min fraction of PME nodes to test with" },
2187 { "-npme", FALSE, etENUM, {npmevalues_opt},
2188 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2189 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2190 { "-fix", FALSE, etINT, {&npme_fixed},
2191 "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."},
2192 { "-rmax", FALSE, etREAL, {&rmax},
2193 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2194 { "-rmin", FALSE, etREAL, {&rmin},
2195 "If >0, minimal rcoulomb for -ntpr>1" },
2196 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2197 "Scale rvdw along with rcoulomb"},
2198 { "-ntpr", FALSE, etINT, {&ntprs},
2199 "Number of [TT].tpr[tt] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2200 "If < 1, automatically choose the number of [TT].tpr[tt] files to test" },
2201 { "-steps", FALSE, etINT64, {&bench_nsteps},
2202 "Take timings for this many steps in the benchmark runs" },
2203 { "-resetstep", FALSE, etINT, {&presteps},
2204 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2205 { "-simsteps", FALSE, etINT64, {&new_sim_nsteps},
2206 "If non-negative, perform this many steps in the real run (overwrites nsteps from [TT].tpr[tt], add [TT].cpt[tt] steps)" },
2207 { "-launch", FALSE, etBOOL, {&bLaunch},
2208 "Launch the real simulation after optimization" },
2209 { "-bench", FALSE, etBOOL, {&bBenchmark},
2210 "Run the benchmarks or just create the input [TT].tpr[tt] files?" },
2211 /******************/
2212 /* mdrun options: */
2213 /******************/
2214 /* We let g_tune_pme parse and understand these options, because we need to
2215 * prevent that they appear on the mdrun command line for the benchmarks */
2216 { "-append", FALSE, etBOOL, {&bAppendFiles},
2217 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2218 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2219 "Keep and number checkpoint files (launch only)" },
2220 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2221 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2224 #define NFILE asize(fnm)
2226 seconds = gmx_gettime();
2228 if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2229 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2235 /* Store the remaining unparsed command line entries in a string which
2236 * is then attached to the mdrun command line */
2238 ExtraArgs[0] = '\0';
2239 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2241 add_to_string(&ExtraArgs, argv[i]);
2242 add_to_string(&ExtraArgs, " ");
2245 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2248 if (opt2parg_bSet("-npstring", asize(pa), pa))
2250 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2255 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2257 /* and now we just set this; a bit of an ugly hack*/
2260 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2261 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2263 /* Automatically set -beo options if -eo is set etc. */
2264 couple_files_options(NFILE, fnm);
2266 /* Construct the command line arguments for benchmark runs
2267 * as well as for the simulation run */
2270 sprintf(bbuf, " -ntmpi %d ", nthreads);
2274 /* This string will be used for MPI runs and will appear after the
2275 * mpirun command. */
2276 if (strcmp(procstring[0], "none") != 0)
2278 sprintf(bbuf, " %s %d ", procstring[0], nnodes);
2288 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2289 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs);
2291 /* Read in checkpoint file if requested */
2293 if (opt2bSet("-cpi", NFILE, fnm))
2296 cr->duty = DUTY_PP; /* makes the following routine happy */
2297 read_checkpoint_simulation_part(opt2fn("-cpi", NFILE, fnm),
2298 &sim_part, &cpt_steps, cr,
2299 FALSE, NFILE, fnm, NULL, NULL);
2302 /* sim_part will now be 1 if no checkpoint file was found */
2305 gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi", NFILE, fnm));
2309 /* Open performance output file and write header info */
2310 fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2312 /* Make a quick consistency check of command line parameters */
2313 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2314 maxPMEfraction, minPMEfraction, npme_fixed,
2315 bench_nsteps, fnm, NFILE, sim_part, presteps,
2318 /* Determine the maximum and minimum number of PME nodes to test,
2319 * the actual list of settings is build in do_the_tests(). */
2320 if ((nnodes > 2) && (npme_fixed < -1))
2322 if (0 == strcmp(npmevalues_opt[0], "auto"))
2324 /* Determine the npme range automatically based on the PME:PP load guess */
2325 if (guessPMEratio > 1.0)
2327 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2328 maxPMEnodes = nnodes/2;
2329 minPMEnodes = nnodes/2;
2333 /* PME : PP load is in the range 0..1, let's test around the guess */
2334 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2335 minPMEnodes = floor(0.7*guessPMEnodes);
2336 maxPMEnodes = ceil(1.6*guessPMEnodes);
2337 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2342 /* Determine the npme range based on user input */
2343 maxPMEnodes = floor(maxPMEfraction*nnodes);
2344 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2345 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2346 if (maxPMEnodes != minPMEnodes)
2348 fprintf(stdout, "- %d ", maxPMEnodes);
2350 fprintf(stdout, "PME-only nodes.\n Note that the automatic number of PME-only nodes and no separate PME nodes are always tested.\n");
2359 /* Get the commands we need to set up the runs from environment variables */
2360 get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2361 if (bBenchmark && repeats > 0)
2363 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun);
2366 /* Print some header info to file */
2368 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2370 fprintf(fp, "%s for Gromacs %s\n", ShortProgram(), GromacsVersion());
2373 fprintf(fp, "Number of nodes : %d\n", nnodes);
2374 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2375 if (strcmp(procstring[0], "none") != 0)
2377 fprintf(fp, "Passing # of nodes via : %s\n", procstring[0]);
2381 fprintf(fp, "Not setting number of nodes in system call\n");
2386 fprintf(fp, "Number of threads : %d\n", nnodes);
2389 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2390 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2391 fprintf(fp, "Benchmark steps : ");
2392 fprintf(fp, "%"GMX_PRId64, bench_nsteps);
2394 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2397 fprintf(fp, "Checkpoint time step : ");
2398 fprintf(fp, "%"GMX_PRId64, cpt_steps);
2401 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2403 if (new_sim_nsteps >= 0)
2406 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2407 fprintf(stderr, "%"GMX_PRId64, new_sim_nsteps+cpt_steps);
2408 fprintf(stderr, " steps.\n");
2409 fprintf(fp, "Simulation steps : ");
2410 fprintf(fp, "%"GMX_PRId64, new_sim_nsteps);
2415 fprintf(fp, "Repeats for each test : %d\n", repeats);
2418 if (npme_fixed >= -1)
2420 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2423 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2424 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2426 /* Allocate memory for the inputinfo struct: */
2428 info->nr_inputfiles = ntprs;
2429 for (i = 0; i < ntprs; i++)
2431 snew(info->rcoulomb, ntprs);
2432 snew(info->rvdw, ntprs);
2433 snew(info->rlist, ntprs);
2434 snew(info->rlistlong, ntprs);
2435 snew(info->nkx, ntprs);
2436 snew(info->nky, ntprs);
2437 snew(info->nkz, ntprs);
2438 snew(info->fsx, ntprs);
2439 snew(info->fsy, ntprs);
2440 snew(info->fsz, ntprs);
2442 /* Make alternative tpr files to test: */
2443 snew(tpr_names, ntprs);
2444 for (i = 0; i < ntprs; i++)
2446 snew(tpr_names[i], STRLEN);
2449 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2450 * different grids could be found. */
2451 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2452 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2454 /********************************************************************************/
2455 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2456 /********************************************************************************/
2457 snew(perfdata, ntprs);
2460 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2461 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2462 cmd_args_bench, fnm, NFILE, presteps, cpt_steps);
2464 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2466 /* Analyse the results and give a suggestion for optimal settings: */
2467 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2468 repeats, info, &best_tpr, &best_npme);
2470 /* Take the best-performing tpr file and enlarge nsteps to original value */
2471 if (bKeepTPR && !bOverwrite)
2473 simulation_tpr = opt2fn("-s", NFILE, fnm);
2477 simulation_tpr = opt2fn("-so", NFILE, fnm);
2478 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2479 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2482 /* Let's get rid of the temporary benchmark input files */
2483 for (i = 0; i < ntprs; i++)
2485 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2486 remove(tpr_names[i]);
2489 /* Now start the real simulation if the user requested it ... */
2490 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2491 cmd_args_launch, simulation_tpr, best_npme);
2495 /* ... or simply print the performance results to screen: */
2498 finalize(opt2fn("-p", NFILE, fnm));