2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
44 #ifdef HAVE_SYS_TIME_H
60 #include "checkpoint.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_large_int_t orig_sim_steps; /* Number of steps to be done in the real simulation */
104 gmx_large_int_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_large_int_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_large_int_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_large_int_pfmt);
250 sscanf(line, dumstring, &resetsteps);
251 bFoundResetStr = TRUE;
252 if (resetsteps == presteps+cpt_steps)
254 bResetChecked = TRUE;
258 sprintf(dumstring, gmx_large_int_pfmt, resetsteps);
259 sprintf(dumstring2, gmx_large_int_pfmt, 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_np[],
599 char *cmd_mdrun[], int repeats)
603 const char def_mpirun[] = "mpirun";
604 const char def_mdrun[] = "mdrun";
606 const char empty_mpirun[] = "";
608 /* Get the commands we need to set up the runs from environment variables */
611 if ( (cp = getenv("MPIRUN")) != NULL)
613 *cmd_mpirun = strdup(cp);
617 *cmd_mpirun = strdup(def_mpirun);
622 *cmd_mpirun = strdup(empty_mpirun);
625 if ( (cp = getenv("MDRUN" )) != NULL)
627 *cmd_mdrun = strdup(cp);
631 *cmd_mdrun = strdup(def_mdrun);
635 /* Check that the commands will run mdrun (perhaps via mpirun) by
636 * running a very quick test simulation. Requires MPI environment to
637 * be available if applicable. */
638 static void check_mdrun_works(gmx_bool bThreads,
639 const char *cmd_mpirun,
641 const char *cmd_mdrun)
643 char *command = NULL;
647 const char filename[] = "benchtest.log";
649 /* This string should always be identical to the one in copyrite.c,
650 * gmx_print_version_info() in the defined(GMX_MPI) section */
651 const char match_mpi[] = "MPI library: MPI";
652 const char match_mdrun[] = "Program: ";
653 gmx_bool bMdrun = FALSE;
654 gmx_bool bMPI = FALSE;
656 /* Run a small test to see whether mpirun + mdrun work */
657 fprintf(stdout, "Making sure that mdrun can be executed. ");
660 snew(command, strlen(cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
661 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
665 snew(command, strlen(cmd_mpirun) + strlen(cmd_np) + strlen(cmd_mdrun) + strlen(filename) + 50);
666 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
668 fprintf(stdout, "Trying '%s' ... ", command);
669 make_backup(filename);
670 gmx_system_call(command);
672 /* Check if we find the characteristic string in the output: */
673 if (!gmx_fexist(filename))
675 gmx_fatal(FARGS, "Output from test run could not be found.");
678 fp = fopen(filename, "r");
679 /* We need to scan the whole output file, since sometimes the queuing system
680 * also writes stuff to stdout/err */
683 cp = fgets(line, STRLEN, fp);
686 if (str_starts(line, match_mdrun) )
690 if (str_starts(line, match_mpi) )
702 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
704 "seems to have been compiled with MPI instead.",
712 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
714 "seems to have been compiled without MPI support.",
721 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
725 fprintf(stdout, "passed.\n");
733 static void launch_simulation(
734 gmx_bool bLaunch, /* Should the simulation be launched? */
735 FILE *fp, /* General log file */
736 gmx_bool bThreads, /* whether to use threads */
737 char *cmd_mpirun, /* Command for mpirun */
738 char *cmd_np, /* Switch for -np or -ntmpi or empty */
739 char *cmd_mdrun, /* Command for mdrun */
740 char *args_for_mdrun, /* Arguments for mdrun */
741 const char *simulation_tpr, /* This tpr will be simulated */
742 int nPMEnodes) /* Number of PME nodes to use */
747 /* Make enough space for the system call command,
748 * (100 extra chars for -npme ... etc. options should suffice): */
749 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+100);
751 /* Note that the -passall options requires args_for_mdrun to be at the end
752 * of the command line string */
755 sprintf(command, "%s%s-npme %d -s %s %s",
756 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun);
760 sprintf(command, "%s%s%s -npme %d -s %s %s",
761 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun);
764 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
768 /* Now the real thing! */
771 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
774 gmx_system_call(command);
779 static void modify_PMEsettings(
780 gmx_large_int_t simsteps, /* Set this value as number of time steps */
781 gmx_large_int_t init_step, /* Set this value as init_step */
782 const char *fn_best_tpr, /* tpr file with the best performance */
783 const char *fn_sim_tpr) /* name of tpr file to be launched */
791 read_tpx_state(fn_best_tpr, ir, &state, NULL, &mtop);
793 /* Reset nsteps and init_step to the value of the input .tpr file */
794 ir->nsteps = simsteps;
795 ir->init_step = init_step;
797 /* Write the tpr file which will be launched */
798 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, gmx_large_int_pfmt);
799 fprintf(stdout, buf, ir->nsteps);
801 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
807 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
809 /* Make additional TPR files with more computational load for the
810 * direct space processors: */
811 static void make_benchmark_tprs(
812 const char *fn_sim_tpr, /* READ : User-provided tpr file */
813 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
814 gmx_large_int_t benchsteps, /* Number of time steps for benchmark runs */
815 gmx_large_int_t statesteps, /* Step counter in checkpoint file */
816 real rmin, /* Minimal Coulomb radius */
817 real rmax, /* Maximal Coulomb radius */
818 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
819 int *ntprs, /* No. of TPRs to write, each with a different
820 rcoulomb and fourierspacing */
821 t_inputinfo *info, /* Contains information about mdp file options */
822 FILE *fp) /* Write the output here */
828 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
831 gmx_bool bNote = FALSE;
832 real add; /* Add this to rcoul for the next test */
833 real fac = 1.0; /* Scaling factor for Coulomb radius */
834 real fourierspacing; /* Basic fourierspacing from tpr */
837 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
838 *ntprs > 1 ? "s" : "", gmx_large_int_pfmt, benchsteps > 1 ? "s" : "");
839 fprintf(stdout, buf, benchsteps);
842 sprintf(buf, " (adding %s steps from checkpoint file)", gmx_large_int_pfmt);
843 fprintf(stdout, buf, statesteps);
844 benchsteps += statesteps;
846 fprintf(stdout, ".\n");
850 read_tpx_state(fn_sim_tpr, ir, &state, NULL, &mtop);
852 /* Check if some kind of PME was chosen */
853 if (EEL_PME(ir->coulombtype) == FALSE)
855 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
859 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
860 if ( (ir->cutoff_scheme != ecutsVERLET) &&
861 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
863 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
864 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
866 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
867 else if (ir->rcoulomb > ir->rlist)
869 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
870 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
873 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
875 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
879 /* Reduce the number of steps for the benchmarks */
880 info->orig_sim_steps = ir->nsteps;
881 ir->nsteps = benchsteps;
882 /* We must not use init_step from the input tpr file for the benchmarks */
883 info->orig_init_step = ir->init_step;
886 /* For PME-switch potentials, keep the radial distance of the buffer region */
887 nlist_buffer = ir->rlist - ir->rcoulomb;
889 /* Determine length of triclinic box vectors */
890 for (d = 0; d < DIM; d++)
893 for (i = 0; i < DIM; i++)
895 box_size[d] += state.box[d][i]*state.box[d][i];
897 box_size[d] = sqrt(box_size[d]);
900 if (ir->fourier_spacing > 0)
902 info->fsx[0] = ir->fourier_spacing;
903 info->fsy[0] = ir->fourier_spacing;
904 info->fsz[0] = ir->fourier_spacing;
908 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
909 info->fsx[0] = box_size[XX]/ir->nkx;
910 info->fsy[0] = box_size[YY]/ir->nky;
911 info->fsz[0] = box_size[ZZ]/ir->nkz;
914 /* If no value for the fourierspacing was provided on the command line, we
915 * use the reconstruction from the tpr file */
916 if (ir->fourier_spacing > 0)
918 /* Use the spacing from the tpr */
919 fourierspacing = ir->fourier_spacing;
923 /* Use the maximum observed spacing */
924 fourierspacing = max(max(info->fsx[0], info->fsy[0]), info->fsz[0]);
927 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
929 /* For performance comparisons the number of particles is useful to have */
930 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
932 /* Print information about settings of which some are potentially modified: */
933 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
934 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
935 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
936 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
937 if (EVDW_SWITCHED(ir->vdwtype))
939 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
941 if (EPME_SWITCHED(ir->coulombtype))
943 fprintf(fp, " rlist : %f nm\n", ir->rlist);
945 if (ir->rlistlong != max_cutoff(ir->rvdw, ir->rcoulomb))
947 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
950 /* Print a descriptive line about the tpr settings tested */
951 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
952 fprintf(fp, " No. scaling rcoulomb");
953 fprintf(fp, " nkx nky nkz");
954 fprintf(fp, " spacing");
955 if (evdwCUT == ir->vdwtype)
957 fprintf(fp, " rvdw");
959 if (EPME_SWITCHED(ir->coulombtype))
961 fprintf(fp, " rlist");
963 if (ir->rlistlong != max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb)) )
965 fprintf(fp, " rlistlong");
967 fprintf(fp, " tpr file\n");
969 /* Loop to create the requested number of tpr input files */
970 for (j = 0; j < *ntprs; j++)
972 /* The first .tpr is the provided one, just need to modify nsteps,
973 * so skip the following block */
976 /* Determine which Coulomb radii rc to use in the benchmarks */
977 add = (rmax-rmin)/(*ntprs-1);
978 if (is_equal(rmin, info->rcoulomb[0]))
980 ir->rcoulomb = rmin + j*add;
982 else if (is_equal(rmax, info->rcoulomb[0]))
984 ir->rcoulomb = rmin + (j-1)*add;
988 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
989 add = (rmax-rmin)/(*ntprs-2);
990 ir->rcoulomb = rmin + (j-1)*add;
993 /* Determine the scaling factor fac */
994 fac = ir->rcoulomb/info->rcoulomb[0];
996 /* Scale the Fourier grid spacing */
997 ir->nkx = ir->nky = ir->nkz = 0;
998 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
1000 /* Adjust other radii since various conditions need to be fulfilled */
1001 if (eelPME == ir->coulombtype)
1003 /* plain PME, rcoulomb must be equal to rlist */
1004 ir->rlist = ir->rcoulomb;
1008 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1009 ir->rlist = ir->rcoulomb + nlist_buffer;
1012 if (bScaleRvdw && evdwCUT == ir->vdwtype)
1014 if ( ecutsVERLET == ir->cutoff_scheme)
1016 /* With Verlet, the van der Waals radius must always equal the Coulomb radius */
1017 ir->rvdw = ir->rcoulomb;
1021 /* For vdw cutoff, rvdw >= rlist */
1022 ir->rvdw = max(info->rvdw[0], ir->rlist);
1026 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1028 } /* end of "if (j != 0)" */
1030 /* for j==0: Save the original settings
1031 * for j >0: Save modified radii and Fourier grids */
1032 info->rcoulomb[j] = ir->rcoulomb;
1033 info->rvdw[j] = ir->rvdw;
1034 info->nkx[j] = ir->nkx;
1035 info->nky[j] = ir->nky;
1036 info->nkz[j] = ir->nkz;
1037 info->rlist[j] = ir->rlist;
1038 info->rlistlong[j] = ir->rlistlong;
1039 info->fsx[j] = fac*fourierspacing;
1040 info->fsy[j] = fac*fourierspacing;
1041 info->fsz[j] = fac*fourierspacing;
1043 /* Write the benchmark tpr file */
1044 strncpy(fn_bench_tprs[j], fn_sim_tpr, strlen(fn_sim_tpr)-strlen(".tpr"));
1045 sprintf(buf, "_bench%.2d.tpr", j);
1046 strcat(fn_bench_tprs[j], buf);
1047 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1048 fprintf(stdout, gmx_large_int_pfmt, ir->nsteps);
1051 fprintf(stdout, ", scaling factor %f\n", fac);
1055 fprintf(stdout, ", unmodified settings\n");
1058 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1060 /* Write information about modified tpr settings to log file */
1061 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1062 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1063 fprintf(fp, " %9f ", info->fsx[j]);
1064 if (evdwCUT == ir->vdwtype)
1066 fprintf(fp, "%10f", ir->rvdw);
1068 if (EPME_SWITCHED(ir->coulombtype))
1070 fprintf(fp, "%10f", ir->rlist);
1072 if (info->rlistlong[0] != max_cutoff(info->rlist[0], max_cutoff(info->rvdw[0], info->rcoulomb[0])) )
1074 fprintf(fp, "%10f", ir->rlistlong);
1076 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1078 /* Make it clear to the user that some additional settings were modified */
1079 if (!is_equal(ir->rvdw, info->rvdw[0])
1080 || !is_equal(ir->rlistlong, info->rlistlong[0]) )
1087 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1088 "other input settings were also changed (see table above).\n"
1089 "Please check if the modified settings are appropriate.\n");
1097 /* Rename the files we want to keep to some meaningful filename and
1098 * delete the rest */
1099 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1100 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1102 char numstring[STRLEN];
1103 char newfilename[STRLEN];
1104 const char *fn = NULL;
1109 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1111 for (i = 0; i < nfile; i++)
1113 opt = (char *)fnm[i].opt;
1114 if (strcmp(opt, "-p") == 0)
1116 /* do nothing; keep this file */
1119 else if (strcmp(opt, "-bg") == 0)
1121 /* Give the log file a nice name so one can later see which parameters were used */
1122 numstring[0] = '\0';
1125 sprintf(numstring, "_%d", nr);
1127 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1128 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1130 fprintf(stdout, "renaming log file to %s\n", newfilename);
1131 make_backup(newfilename);
1132 rename(opt2fn("-bg", nfile, fnm), newfilename);
1135 else if (strcmp(opt, "-err") == 0)
1137 /* This file contains the output of stderr. We want to keep it in
1138 * cases where there have been problems. */
1139 fn = opt2fn(opt, nfile, fnm);
1140 numstring[0] = '\0';
1143 sprintf(numstring, "_%d", nr);
1145 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1150 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1151 make_backup(newfilename);
1152 rename(fn, newfilename);
1156 fprintf(stdout, "Deleting %s\n", fn);
1161 /* Delete the files which are created for each benchmark run: (options -b*) */
1162 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1164 fn = opt2fn(opt, nfile, fnm);
1167 fprintf(stdout, "Deleting %s\n", fn);
1175 /* Returns the largest common factor of n1 and n2 */
1176 static int largest_common_factor(int n1, int n2)
1181 for (factor = nmax; factor > 0; factor--)
1183 if (0 == (n1 % factor) && 0 == (n2 % factor) )
1188 return 0; /* one for the compiler */
1192 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1195 /* Create a list of numbers of PME nodes to test */
1196 static void make_npme_list(
1197 const char *npmevalues_opt, /* Make a complete list with all
1198 * possibilities or a short list that keeps only
1199 * reasonable numbers of PME nodes */
1200 int *nentries, /* Number of entries we put in the nPMEnodes list */
1201 int *nPMEnodes[], /* Each entry contains the value for -npme */
1202 int nnodes, /* Total number of nodes to do the tests on */
1203 int minPMEnodes, /* Minimum number of PME nodes */
1204 int maxPMEnodes) /* Maximum number of PME nodes */
1207 int min_factor = 1; /* We request that npp and npme have this minimal
1208 * largest common factor (depends on npp) */
1209 int nlistmax; /* Max. list size */
1210 int nlist; /* Actual number of entries in list */
1214 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1215 if (0 == strcmp(npmevalues_opt, "all") )
1219 else if (0 == strcmp(npmevalues_opt, "subset") )
1221 eNPME = eNpmeSubset;
1223 else /* "auto" or "range" */
1229 else if (nnodes < 128)
1231 eNPME = eNpmeReduced;
1235 eNPME = eNpmeSubset;
1239 /* Calculate how many entries we could possibly have (in case of -npme all) */
1242 nlistmax = maxPMEnodes - minPMEnodes + 3;
1243 if (0 == minPMEnodes)
1253 /* Now make the actual list which is at most of size nlist */
1254 snew(*nPMEnodes, nlistmax);
1255 nlist = 0; /* start counting again, now the real entries in the list */
1256 for (i = 0; i < nlistmax - 2; i++)
1258 npme = maxPMEnodes - i;
1269 /* For 2d PME we want a common largest factor of at least the cube
1270 * root of the number of PP nodes */
1271 min_factor = (int) pow(npp, 1.0/3.0);
1274 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1277 if (largest_common_factor(npp, npme) >= min_factor)
1279 (*nPMEnodes)[nlist] = npme;
1283 /* We always test 0 PME nodes and the automatic number */
1284 *nentries = nlist + 2;
1285 (*nPMEnodes)[nlist ] = 0;
1286 (*nPMEnodes)[nlist+1] = -1;
1288 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1289 for (i = 0; i < *nentries-1; i++)
1291 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1293 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1297 /* Allocate memory to store the performance data */
1298 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1303 for (k = 0; k < ntprs; k++)
1305 snew(perfdata[k], datasets);
1306 for (i = 0; i < datasets; i++)
1308 for (j = 0; j < repeats; j++)
1310 snew(perfdata[k][i].Gcycles, repeats);
1311 snew(perfdata[k][i].ns_per_day, repeats);
1312 snew(perfdata[k][i].PME_f_load, repeats);
1319 /* Check for errors on mdrun -h */
1320 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp)
1322 char *command, *msg;
1326 snew(command, length + 15);
1327 snew(msg, length + 500);
1329 fprintf(stdout, "Making shure the benchmarks can be executed ...\n");
1330 sprintf(command, "%s-h -quiet", mdrun_cmd_line);
1331 ret = gmx_system_call(command);
1335 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1336 * get the error message from mdrun itself */
1337 sprintf(msg, "Cannot run the benchmark simulations! Please check the error message of\n"
1338 "mdrun for the source of the problem. Did you provide a command line\n"
1339 "argument that neither g_tune_pme nor mdrun understands? Offending command:\n"
1340 "\n%s\n\n", command);
1342 fprintf(stderr, "%s", msg);
1344 fprintf(fp, "%s", msg);
1354 static void do_the_tests(
1355 FILE *fp, /* General g_tune_pme output file */
1356 char **tpr_names, /* Filenames of the input files to test */
1357 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1358 int minPMEnodes, /* Min fraction of nodes to use for PME */
1359 int npme_fixed, /* If >= -1, test fixed number of PME
1361 const char *npmevalues_opt, /* Which -npme values should be tested */
1362 t_perf **perfdata, /* Here the performace data is stored */
1363 int *pmeentries, /* Entries in the nPMEnodes list */
1364 int repeats, /* Repeat each test this often */
1365 int nnodes, /* Total number of nodes = nPP + nPME */
1366 int nr_tprs, /* Total number of tpr files to test */
1367 gmx_bool bThreads, /* Threads or MPI? */
1368 char *cmd_mpirun, /* mpirun command string */
1369 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1370 char *cmd_mdrun, /* mdrun command string */
1371 char *cmd_args_bench, /* arguments for mdrun in a string */
1372 const t_filenm *fnm, /* List of filenames from command line */
1373 int nfile, /* Number of files specified on the cmdl. */
1374 int presteps, /* DLB equilibration steps, is checked */
1375 gmx_large_int_t cpt_steps) /* Time step counter in the checkpoint */
1377 int i, nr, k, ret, count = 0, totaltests;
1378 int *nPMEnodes = NULL;
1381 char *command, *cmd_stub;
1383 gmx_bool bResetProblem = FALSE;
1384 gmx_bool bFirst = TRUE;
1387 /* This string array corresponds to the eParselog enum type at the start
1389 const char* ParseLog[] = {
1391 "Logfile not found!",
1392 "No timings, logfile truncated?",
1393 "Run was terminated.",
1394 "Counters were not reset properly.",
1395 "No DD grid found for these settings.",
1396 "TPX version conflict!",
1397 "mdrun was not started in parallel!",
1398 "Number of PP nodes has a prime factor that is too large.",
1401 char str_PME_f_load[13];
1404 /* Allocate space for the mdrun command line. 100 extra characters should
1405 be more than enough for the -npme etcetera arguments */
1406 cmdline_length = strlen(cmd_mpirun)
1409 + strlen(cmd_args_bench)
1410 + strlen(tpr_names[0]) + 100;
1411 snew(command, cmdline_length);
1412 snew(cmd_stub, cmdline_length);
1414 /* Construct the part of the command line that stays the same for all tests: */
1417 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1421 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1424 /* Create a list of numbers of PME nodes to test */
1425 if (npme_fixed < -1)
1427 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1428 nnodes, minPMEnodes, maxPMEnodes);
1434 nPMEnodes[0] = npme_fixed;
1435 fprintf(stderr, "Will use a fixed number of %d PME-only nodes.\n", nPMEnodes[0]);
1440 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1442 finalize(opt2fn("-p", nfile, fnm));
1446 /* Allocate one dataset for each tpr input file: */
1447 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1449 /*****************************************/
1450 /* Main loop over all tpr files to test: */
1451 /*****************************************/
1452 totaltests = nr_tprs*(*pmeentries)*repeats;
1453 for (k = 0; k < nr_tprs; k++)
1455 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1456 fprintf(fp, "PME nodes Gcycles ns/day PME/f Remark\n");
1457 /* Loop over various numbers of PME nodes: */
1458 for (i = 0; i < *pmeentries; i++)
1460 pd = &perfdata[k][i];
1462 /* Loop over the repeats for each scenario: */
1463 for (nr = 0; nr < repeats; nr++)
1465 pd->nPMEnodes = nPMEnodes[i];
1467 /* Add -npme and -s to the command line and save it. Note that
1468 * the -passall (if set) options requires cmd_args_bench to be
1469 * at the end of the command line string */
1470 snew(pd->mdrun_cmd_line, cmdline_length);
1471 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s",
1472 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench);
1474 /* To prevent that all benchmarks fail due to a show-stopper argument
1475 * on the mdrun command line, we make a quick check with mdrun -h first */
1478 make_sure_it_runs(pd->mdrun_cmd_line, cmdline_length, fp);
1482 /* Do a benchmark simulation: */
1485 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1491 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1492 (100.0*count)/totaltests,
1493 k+1, nr_tprs, i+1, *pmeentries, buf);
1494 make_backup(opt2fn("-err", nfile, fnm));
1495 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1496 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1497 gmx_system_call(command);
1499 /* Collect the performance data from the log file; also check stderr
1500 * for fatal errors */
1501 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1502 pd, nr, presteps, cpt_steps, nnodes);
1503 if ((presteps > 0) && (ret == eParselogResetProblem))
1505 bResetProblem = TRUE;
1508 if (-1 == pd->nPMEnodes)
1510 sprintf(buf, "(%3d)", pd->guessPME);
1518 if (pd->PME_f_load[nr] > 0.0)
1520 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1524 sprintf(str_PME_f_load, "%s", " - ");
1527 /* Write the data we got to disk */
1528 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1529 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1530 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1532 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1538 /* Do some cleaning up and delete the files we do not need any more */
1539 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1541 /* If the first run with this number of processors already failed, do not try again: */
1542 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1544 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1545 count += repeats-(nr+1);
1548 } /* end of repeats loop */
1549 } /* end of -npme loop */
1550 } /* end of tpr file loop */
1555 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1563 static void check_input(
1570 real maxPMEfraction,
1571 real minPMEfraction,
1573 gmx_large_int_t bench_nsteps,
1574 const t_filenm *fnm,
1584 /* Make sure the input file exists */
1585 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1587 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1590 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1591 if ( (0 == strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1593 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1594 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1597 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1600 gmx_fatal(FARGS, "Number of repeats < 0!");
1603 /* Check number of nodes */
1606 gmx_fatal(FARGS, "Number of nodes/threads must be a positive integer.");
1609 /* Automatically choose -ntpr if not set */
1619 /* Set a reasonable scaling factor for rcoulomb */
1622 *rmax = rcoulomb * 1.2;
1625 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1631 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1635 /* Make shure that rmin <= rcoulomb <= rmax */
1644 if (!(*rmin <= *rmax) )
1646 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1647 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1649 /* Add test scenarios if rmin or rmax were set */
1652 if (!is_equal(*rmin, rcoulomb) && (*ntprs == 1) )
1655 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1658 if (!is_equal(*rmax, rcoulomb) && (*ntprs == 1) )
1661 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1666 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1667 if (!is_equal(*rmax, rcoulomb) || !is_equal(*rmin, rcoulomb) )
1669 *ntprs = max(*ntprs, 2);
1672 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1673 if (!is_equal(*rmax, rcoulomb) && !is_equal(*rmin, rcoulomb) )
1675 *ntprs = max(*ntprs, 3);
1680 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1685 if (is_equal(*rmin, rcoulomb) && is_equal(rcoulomb, *rmax)) /* We have just a single rc */
1687 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1688 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1689 "with correspondingly adjusted PME grid settings\n");
1694 /* Check whether max and min fraction are within required values */
1695 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1697 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1699 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1701 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1703 if (maxPMEfraction < minPMEfraction)
1705 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1708 /* Check whether the number of steps - if it was set - has a reasonable value */
1709 if (bench_nsteps < 0)
1711 gmx_fatal(FARGS, "Number of steps must be positive.");
1714 if (bench_nsteps > 10000 || bench_nsteps < 100)
1716 fprintf(stderr, "WARNING: steps=");
1717 fprintf(stderr, gmx_large_int_pfmt, bench_nsteps);
1718 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1723 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1726 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1729 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1731 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1735 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1736 * only. We need to check whether the requested number of PME-only nodes
1738 if (npme_fixed > -1)
1740 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1741 if (2*npme_fixed > nnodes)
1743 gmx_fatal(FARGS, "Cannot have more than %d PME-only nodes for a total of %d nodes (you chose %d).\n",
1744 nnodes/2, nnodes, npme_fixed);
1746 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1748 fprintf(stderr, "WARNING: Only %g percent of the nodes are assigned as PME-only nodes.\n",
1749 100.0*((real)npme_fixed / (real)nnodes));
1751 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1753 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1754 " fixed number of PME-only nodes is requested with -fix.\n");
1760 /* Returns TRUE when "opt" is needed at launch time */
1761 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1763 /* Apart from the input .tpr and the output log files we need all options that
1764 * were set on the command line and that do not start with -b */
1765 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2)
1766 || 0 == strncmp(opt, "-err", 4) || 0 == strncmp(opt, "-p", 2) )
1775 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1776 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1778 /* Apart from the input .tpr, all files starting with "-b" are for
1779 * _b_enchmark files exclusively */
1780 if (0 == strncmp(opt, "-s", 2))
1785 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2))
1787 if (!bOptional || bSet)
1804 if (bSet) /* These are additional input files like -cpi -ei */
1817 /* Adds 'buf' to 'str' */
1818 static void add_to_string(char **str, char *buf)
1823 len = strlen(*str) + strlen(buf) + 1;
1829 /* Create the command line for the benchmark as well as for the real run */
1830 static void create_command_line_snippets(
1831 gmx_bool bAppendFiles,
1832 gmx_bool bKeepAndNumCPT,
1833 gmx_bool bResetHWay,
1837 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1838 char *cmd_args_launch[], /* command line arguments for simulation run */
1839 char extra_args[]) /* Add this to the end of the command line */
1844 char strbuf[STRLEN];
1847 /* strlen needs at least '\0' as a string: */
1848 snew(*cmd_args_bench, 1);
1849 snew(*cmd_args_launch, 1);
1850 *cmd_args_launch[0] = '\0';
1851 *cmd_args_bench[0] = '\0';
1854 /*******************************************/
1855 /* 1. Process other command line arguments */
1856 /*******************************************/
1859 /* Add equilibration steps to benchmark options */
1860 sprintf(strbuf, "-resetstep %d ", presteps);
1861 add_to_string(cmd_args_bench, strbuf);
1863 /* These switches take effect only at launch time */
1864 if (FALSE == bAppendFiles)
1866 add_to_string(cmd_args_launch, "-noappend ");
1870 add_to_string(cmd_args_launch, "-cpnum ");
1874 add_to_string(cmd_args_launch, "-resethway ");
1877 /********************/
1878 /* 2. Process files */
1879 /********************/
1880 for (i = 0; i < nfile; i++)
1882 opt = (char *)fnm[i].opt;
1883 name = opt2fn(opt, nfile, fnm);
1885 /* Strbuf contains the options, now let's sort out where we need that */
1886 sprintf(strbuf, "%s %s ", opt, name);
1888 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1890 /* All options starting with -b* need the 'b' removed,
1891 * therefore overwrite strbuf */
1892 if (0 == strncmp(opt, "-b", 2))
1894 sprintf(strbuf, "-%s %s ", &opt[2], name);
1897 add_to_string(cmd_args_bench, strbuf);
1900 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
1902 add_to_string(cmd_args_launch, strbuf);
1906 add_to_string(cmd_args_bench, extra_args);
1907 add_to_string(cmd_args_launch, extra_args);
1911 /* Set option opt */
1912 static void setopt(const char *opt, int nfile, t_filenm fnm[])
1916 for (i = 0; (i < nfile); i++)
1918 if (strcmp(opt, fnm[i].opt) == 0)
1920 fnm[i].flag |= ffSET;
1926 /* This routine inspects the tpr file and ...
1927 * 1. checks for output files that get triggered by a tpr option. These output
1928 * files are marked as 'set' to allow for a proper cleanup after each
1930 * 2. returns the PME:PP load ratio
1931 * 3. returns rcoulomb from the tpr */
1932 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1934 gmx_bool bPull; /* Is pulling requested in .tpr file? */
1935 gmx_bool bTpi; /* Is test particle insertion requested? */
1936 gmx_bool bFree; /* Is a free energy simulation requested? */
1937 gmx_bool bNM; /* Is a normal mode analysis requested? */
1943 /* Check tpr file for options that trigger extra output files */
1944 read_tpx_state(opt2fn("-s", nfile, fnm), &ir, &state, NULL, &mtop);
1945 bPull = (epullNO != ir.ePull);
1946 bFree = (efepNO != ir.efep );
1947 bNM = (eiNM == ir.eI );
1948 bTpi = EI_TPI(ir.eI);
1950 /* Set these output files on the tuning command-line */
1953 setopt("-pf", nfile, fnm);
1954 setopt("-px", nfile, fnm);
1958 setopt("-dhdl", nfile, fnm);
1962 setopt("-tpi", nfile, fnm);
1963 setopt("-tpid", nfile, fnm);
1967 setopt("-mtx", nfile, fnm);
1970 *rcoulomb = ir.rcoulomb;
1972 /* Return the estimate for the number of PME nodes */
1973 return pme_load_estimate(&mtop, &ir, state.box);
1977 static void couple_files_options(int nfile, t_filenm fnm[])
1980 gmx_bool bSet, bBench;
1985 for (i = 0; i < nfile; i++)
1987 opt = (char *)fnm[i].opt;
1988 bSet = ((fnm[i].flag & ffSET) != 0);
1989 bBench = (0 == strncmp(opt, "-b", 2));
1991 /* Check optional files */
1992 /* If e.g. -eo is set, then -beo also needs to be set */
1993 if (is_optional(&fnm[i]) && bSet && !bBench)
1995 sprintf(buf, "-b%s", &opt[1]);
1996 setopt(buf, nfile, fnm);
1998 /* If -beo is set, then -eo also needs to be! */
1999 if (is_optional(&fnm[i]) && bSet && bBench)
2001 sprintf(buf, "-%s", &opt[2]);
2002 setopt(buf, nfile, fnm);
2008 static double gettime()
2010 #ifdef HAVE_GETTIMEOFDAY
2014 gettimeofday(&t, NULL);
2016 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
2022 seconds = time(NULL);
2029 #define BENCHSTEPS (1000)
2031 int gmx_tune_pme(int argc, char *argv[])
2033 const char *desc[] = {
2034 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of processors/threads, this program systematically",
2035 "times [TT]mdrun[tt] with various numbers of PME-only nodes and determines",
2036 "which setting is fastest. It will also test whether performance can",
2037 "be enhanced by shifting load from the reciprocal to the real space",
2038 "part of the Ewald sum. ",
2039 "Simply pass your [TT].tpr[tt] file to [TT]g_tune_pme[tt] together with other options",
2040 "for [TT]mdrun[tt] as needed.[PAR]",
2041 "Which executables are used can be set in the environment variables",
2042 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
2043 "will be used as defaults. Note that for certain MPI frameworks you",
2044 "need to provide a machine- or hostfile. This can also be passed",
2045 "via the MPIRUN variable, e.g.[PAR]",
2046 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
2047 "Please call [TT]g_tune_pme[tt] with the normal options you would pass to",
2048 "[TT]mdrun[tt] and add [TT]-np[tt] for the number of processors to perform the",
2049 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2050 "to repeat each test several times to get better statistics. [PAR]",
2051 "[TT]g_tune_pme[tt] can test various real space / reciprocal space workloads",
2052 "for you. With [TT]-ntpr[tt] you control how many extra [TT].tpr[tt] files will be",
2053 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2054 "Typically, the first test (number 0) will be with the settings from the input",
2055 "[TT].tpr[tt] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2056 "specified by [TT]-rmax[tt] with a somwhat smaller PME grid at the same time. ",
2057 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2058 "The remaining [TT].tpr[tt] files will have equally-spaced Coulomb radii (and Fourier "
2059 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2060 "if you just seek the optimal number of PME-only nodes; in that case",
2061 "your input [TT].tpr[tt] file will remain unchanged.[PAR]",
2062 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2063 "MD systems. The dynamic load balancing needs about 100 time steps",
2064 "to adapt to local load imbalances, therefore the time step counters",
2065 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2066 "for a higher accuarcy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2067 "From the 'DD' load imbalance entries in the md.log output file you",
2068 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2069 "[TT]g_tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2070 "After calling [TT]mdrun[tt] several times, detailed performance information",
2071 "is available in the output file [TT]perf.out.[tt] ",
2072 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2073 "(options [TT]-b[tt]*), these will be automatically deleted after each test.[PAR]",
2074 "If you want the simulation to be started automatically with the",
2075 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2080 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2081 real maxPMEfraction = 0.50;
2082 real minPMEfraction = 0.25;
2083 int maxPMEnodes, minPMEnodes;
2084 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2085 float guessPMEnodes;
2086 int npme_fixed = -2; /* If >= -1, use only this number
2087 * of PME-only nodes */
2089 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2090 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2091 gmx_bool bScaleRvdw = TRUE;
2092 gmx_large_int_t bench_nsteps = BENCHSTEPS;
2093 gmx_large_int_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2094 gmx_large_int_t cpt_steps = 0; /* Step counter in .cpt input file */
2095 int presteps = 100; /* Do a full cycle reset after presteps steps */
2096 gmx_bool bOverwrite = FALSE, bKeepTPR;
2097 gmx_bool bLaunch = FALSE;
2098 char *ExtraArgs = NULL;
2099 char **tpr_names = NULL;
2100 const char *simulation_tpr = NULL;
2101 int best_npme, best_tpr;
2102 int sim_part = 1; /* For benchmarks with checkpoint files */
2105 /* Default program names if nothing else is found */
2106 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2107 char *cmd_args_bench, *cmd_args_launch;
2108 char *cmd_np = NULL;
2110 t_perf **perfdata = NULL;
2116 /* Print out how long the tuning took */
2119 static t_filenm fnm[] = {
2121 { efOUT, "-p", "perf", ffWRITE },
2122 { efLOG, "-err", "bencherr", ffWRITE },
2123 { efTPX, "-so", "tuned", ffWRITE },
2125 { efTPX, NULL, NULL, ffREAD },
2126 { efTRN, "-o", NULL, ffWRITE },
2127 { efXTC, "-x", NULL, ffOPTWR },
2128 { efCPT, "-cpi", NULL, ffOPTRD },
2129 { efCPT, "-cpo", NULL, ffOPTWR },
2130 { efSTO, "-c", "confout", ffWRITE },
2131 { efEDR, "-e", "ener", ffWRITE },
2132 { efLOG, "-g", "md", ffWRITE },
2133 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2134 { efXVG, "-field", "field", ffOPTWR },
2135 { efXVG, "-table", "table", ffOPTRD },
2136 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2137 { efXVG, "-tablep", "tablep", ffOPTRD },
2138 { efXVG, "-tableb", "table", ffOPTRD },
2139 { efTRX, "-rerun", "rerun", ffOPTRD },
2140 { efXVG, "-tpi", "tpi", ffOPTWR },
2141 { efXVG, "-tpid", "tpidist", ffOPTWR },
2142 { efEDI, "-ei", "sam", ffOPTRD },
2143 { efXVG, "-eo", "edsam", ffOPTWR },
2144 { efGCT, "-j", "wham", ffOPTRD },
2145 { efGCT, "-jo", "bam", ffOPTWR },
2146 { efXVG, "-ffout", "gct", ffOPTWR },
2147 { efXVG, "-devout", "deviatie", ffOPTWR },
2148 { efXVG, "-runav", "runaver", ffOPTWR },
2149 { efXVG, "-px", "pullx", ffOPTWR },
2150 { efXVG, "-pf", "pullf", ffOPTWR },
2151 { efXVG, "-ro", "rotation", ffOPTWR },
2152 { efLOG, "-ra", "rotangles", ffOPTWR },
2153 { efLOG, "-rs", "rotslabs", ffOPTWR },
2154 { efLOG, "-rt", "rottorque", ffOPTWR },
2155 { efMTX, "-mtx", "nm", ffOPTWR },
2156 { efNDX, "-dn", "dipole", ffOPTWR },
2157 /* Output files that are deleted after each benchmark run */
2158 { efTRN, "-bo", "bench", ffWRITE },
2159 { efXTC, "-bx", "bench", ffWRITE },
2160 { efCPT, "-bcpo", "bench", ffWRITE },
2161 { efSTO, "-bc", "bench", ffWRITE },
2162 { efEDR, "-be", "bench", ffWRITE },
2163 { efLOG, "-bg", "bench", ffWRITE },
2164 { efXVG, "-beo", "benchedo", ffOPTWR },
2165 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2166 { efXVG, "-bfield", "benchfld", ffOPTWR },
2167 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2168 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2169 { efGCT, "-bjo", "bench", ffOPTWR },
2170 { efXVG, "-bffout", "benchgct", ffOPTWR },
2171 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2172 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2173 { efXVG, "-bpx", "benchpx", ffOPTWR },
2174 { efXVG, "-bpf", "benchpf", ffOPTWR },
2175 { efXVG, "-bro", "benchrot", ffOPTWR },
2176 { efLOG, "-bra", "benchrota", ffOPTWR },
2177 { efLOG, "-brs", "benchrots", ffOPTWR },
2178 { efLOG, "-brt", "benchrott", ffOPTWR },
2179 { efMTX, "-bmtx", "benchn", ffOPTWR },
2180 { efNDX, "-bdn", "bench", ffOPTWR }
2183 gmx_bool bThreads = FALSE;
2187 const char *procstring[] =
2188 { NULL, "-np", "-n", "none", NULL };
2189 const char *npmevalues_opt[] =
2190 { NULL, "auto", "all", "subset", NULL };
2192 gmx_bool bAppendFiles = TRUE;
2193 gmx_bool bKeepAndNumCPT = FALSE;
2194 gmx_bool bResetCountersHalfWay = FALSE;
2195 gmx_bool bBenchmark = TRUE;
2197 output_env_t oenv = NULL;
2200 /***********************/
2201 /* g_tune_pme options: */
2202 /***********************/
2203 { "-np", FALSE, etINT, {&nnodes},
2204 "Number of nodes to run the tests on (must be > 2 for separate PME nodes)" },
2205 { "-npstring", FALSE, etENUM, {procstring},
2206 "Specify the number of processors to [TT]$MPIRUN[tt] using this string"},
2207 { "-ntmpi", FALSE, etINT, {&nthreads},
2208 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2209 { "-r", FALSE, etINT, {&repeats},
2210 "Repeat each test this often" },
2211 { "-max", FALSE, etREAL, {&maxPMEfraction},
2212 "Max fraction of PME nodes to test with" },
2213 { "-min", FALSE, etREAL, {&minPMEfraction},
2214 "Min fraction of PME nodes to test with" },
2215 { "-npme", FALSE, etENUM, {npmevalues_opt},
2216 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2217 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2218 { "-fix", FALSE, etINT, {&npme_fixed},
2219 "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."},
2220 { "-rmax", FALSE, etREAL, {&rmax},
2221 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2222 { "-rmin", FALSE, etREAL, {&rmin},
2223 "If >0, minimal rcoulomb for -ntpr>1" },
2224 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2225 "Scale rvdw along with rcoulomb"},
2226 { "-ntpr", FALSE, etINT, {&ntprs},
2227 "Number of [TT].tpr[tt] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2228 "If < 1, automatically choose the number of [TT].tpr[tt] files to test" },
2229 { "-steps", FALSE, etGMX_LARGE_INT, {&bench_nsteps},
2230 "Take timings for this many steps in the benchmark runs" },
2231 { "-resetstep", FALSE, etINT, {&presteps},
2232 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2233 { "-simsteps", FALSE, etGMX_LARGE_INT, {&new_sim_nsteps},
2234 "If non-negative, perform this many steps in the real run (overwrites nsteps from [TT].tpr[tt], add [TT].cpt[tt] steps)" },
2235 { "-launch", FALSE, etBOOL, {&bLaunch},
2236 "Launch the real simulation after optimization" },
2237 { "-bench", FALSE, etBOOL, {&bBenchmark},
2238 "Run the benchmarks or just create the input [TT].tpr[tt] files?" },
2239 /******************/
2240 /* mdrun options: */
2241 /******************/
2242 /* We let g_tune_pme parse and understand these options, because we need to
2243 * prevent that they appear on the mdrun command line for the benchmarks */
2244 { "-append", FALSE, etBOOL, {&bAppendFiles},
2245 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2246 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2247 "Keep and number checkpoint files (launch only)" },
2248 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2249 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2253 #define NFILE asize(fnm)
2255 CopyRight(stderr, argv[0]);
2257 seconds = gettime();
2259 parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2260 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2263 /* Store the remaining unparsed command line entries in a string which
2264 * is then attached to the mdrun command line */
2266 ExtraArgs[0] = '\0';
2267 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2269 add_to_string(&ExtraArgs, argv[i]);
2270 add_to_string(&ExtraArgs, " ");
2273 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2276 if (opt2parg_bSet("-npstring", asize(pa), pa))
2278 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2283 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2285 /* and now we just set this; a bit of an ugly hack*/
2288 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2289 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2291 /* Automatically set -beo options if -eo is set etc. */
2292 couple_files_options(NFILE, fnm);
2294 /* Construct the command line arguments for benchmark runs
2295 * as well as for the simulation run */
2298 sprintf(bbuf, " -ntmpi %d ", nthreads);
2302 /* This string will be used for MPI runs and will appear after the
2303 * mpirun command. */
2304 if (strcmp(procstring[0], "none") != 0)
2306 sprintf(bbuf, " %s %d ", procstring[0], nnodes);
2316 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2317 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs);
2319 /* Read in checkpoint file if requested */
2321 if (opt2bSet("-cpi", NFILE, fnm))
2324 cr->duty = DUTY_PP; /* makes the following routine happy */
2325 read_checkpoint_simulation_part(opt2fn("-cpi", NFILE, fnm),
2326 &sim_part, &cpt_steps, cr,
2327 FALSE, NFILE, fnm, NULL, NULL);
2330 /* sim_part will now be 1 if no checkpoint file was found */
2333 gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi", NFILE, fnm));
2337 /* Open performance output file and write header info */
2338 fp = ffopen(opt2fn("-p", NFILE, fnm), "w");
2340 /* Make a quick consistency check of command line parameters */
2341 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2342 maxPMEfraction, minPMEfraction, npme_fixed,
2343 bench_nsteps, fnm, NFILE, sim_part, presteps,
2346 /* Determine the maximum and minimum number of PME nodes to test,
2347 * the actual list of settings is build in do_the_tests(). */
2348 if ((nnodes > 2) && (npme_fixed < -1))
2350 if (0 == strcmp(npmevalues_opt[0], "auto"))
2352 /* Determine the npme range automatically based on the PME:PP load guess */
2353 if (guessPMEratio > 1.0)
2355 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2356 maxPMEnodes = nnodes/2;
2357 minPMEnodes = nnodes/2;
2361 /* PME : PP load is in the range 0..1, let's test around the guess */
2362 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2363 minPMEnodes = floor(0.7*guessPMEnodes);
2364 maxPMEnodes = ceil(1.6*guessPMEnodes);
2365 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2370 /* Determine the npme range based on user input */
2371 maxPMEnodes = floor(maxPMEfraction*nnodes);
2372 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2373 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2374 if (maxPMEnodes != minPMEnodes)
2376 fprintf(stdout, "- %d ", maxPMEnodes);
2378 fprintf(stdout, "PME-only nodes.\n Note that the automatic number of PME-only nodes and no separate PME nodes are always tested.\n");
2387 /* Get the commands we need to set up the runs from environment variables */
2388 get_program_paths(bThreads, &cmd_mpirun, cmd_np, &cmd_mdrun, repeats);
2389 if (bBenchmark && repeats > 0)
2391 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun);
2394 /* Print some header info to file */
2396 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2398 fprintf(fp, "%s for Gromacs %s\n", ShortProgram(), GromacsVersion());
2401 fprintf(fp, "Number of nodes : %d\n", nnodes);
2402 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2403 if (strcmp(procstring[0], "none") != 0)
2405 fprintf(fp, "Passing # of nodes via : %s\n", procstring[0]);
2409 fprintf(fp, "Not setting number of nodes in system call\n");
2414 fprintf(fp, "Number of threads : %d\n", nnodes);
2417 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2418 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2419 fprintf(fp, "Benchmark steps : ");
2420 fprintf(fp, gmx_large_int_pfmt, bench_nsteps);
2422 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2425 fprintf(fp, "Checkpoint time step : ");
2426 fprintf(fp, gmx_large_int_pfmt, cpt_steps);
2429 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2431 if (new_sim_nsteps >= 0)
2434 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2435 fprintf(stderr, gmx_large_int_pfmt, new_sim_nsteps+cpt_steps);
2436 fprintf(stderr, " steps.\n");
2437 fprintf(fp, "Simulation steps : ");
2438 fprintf(fp, gmx_large_int_pfmt, new_sim_nsteps);
2443 fprintf(fp, "Repeats for each test : %d\n", repeats);
2446 if (npme_fixed >= -1)
2448 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2451 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2452 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2454 /* Allocate memory for the inputinfo struct: */
2456 info->nr_inputfiles = ntprs;
2457 for (i = 0; i < ntprs; i++)
2459 snew(info->rcoulomb, ntprs);
2460 snew(info->rvdw, ntprs);
2461 snew(info->rlist, ntprs);
2462 snew(info->rlistlong, ntprs);
2463 snew(info->nkx, ntprs);
2464 snew(info->nky, ntprs);
2465 snew(info->nkz, ntprs);
2466 snew(info->fsx, ntprs);
2467 snew(info->fsy, ntprs);
2468 snew(info->fsz, ntprs);
2470 /* Make alternative tpr files to test: */
2471 snew(tpr_names, ntprs);
2472 for (i = 0; i < ntprs; i++)
2474 snew(tpr_names[i], STRLEN);
2477 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2478 * different grids could be found. */
2479 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2480 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2482 /********************************************************************************/
2483 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2484 /********************************************************************************/
2485 snew(perfdata, ntprs);
2488 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2489 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2490 cmd_args_bench, fnm, NFILE, presteps, cpt_steps);
2492 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gettime()-seconds)/60.0);
2494 /* Analyse the results and give a suggestion for optimal settings: */
2495 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2496 repeats, info, &best_tpr, &best_npme);
2498 /* Take the best-performing tpr file and enlarge nsteps to original value */
2499 if (bKeepTPR && !bOverwrite)
2501 simulation_tpr = opt2fn("-s", NFILE, fnm);
2505 simulation_tpr = opt2fn("-so", NFILE, fnm);
2506 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2507 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2510 /* Let's get rid of the temporary benchmark input files */
2511 for (i = 0; i < ntprs; i++)
2513 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2514 remove(tpr_names[i]);
2517 /* Now start the real simulation if the user requested it ... */
2518 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2519 cmd_args_launch, simulation_tpr, best_npme);
2523 /* ... or simply print the performance results to screen: */
2526 finalize(opt2fn("-p", NFILE, fnm));