3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
40 #ifdef HAVE_SYS_TIME_H
56 #include "checkpoint.h"
64 /* Enum for situations that can occur during log file parsing, the
65 * corresponding string entries can be found in do_the_tests() in
66 * const char* ParseLog[] */
72 eParselogResetProblem,
83 int nPMEnodes; /* number of PME-only nodes used in this test */
84 int nx, ny, nz; /* DD grid */
85 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
86 double *Gcycles; /* This can contain more than one value if doing multiple tests */
90 float *PME_f_load; /* PME mesh/force load average*/
91 float PME_f_load_Av; /* Average average ;) ... */
92 char *mdrun_cmd_line; /* Mdrun command line used for this test */
98 int nr_inputfiles; /* The number of tpr and mdp input files */
99 gmx_large_int_t orig_sim_steps; /* Number of steps to be done in the real simulation */
100 gmx_large_int_t orig_init_step; /* Init step for the real simulation */
101 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
102 real *rvdw; /* The vdW radii */
103 real *rlist; /* Neighbourlist cutoff radius */
105 int *nkx, *nky, *nkz;
106 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
110 static void sep_line(FILE *fp)
112 fprintf(fp, "\n------------------------------------------------------------\n");
116 /* Wrapper for system calls */
117 static int gmx_system_call(char *command)
120 gmx_fatal(FARGS,"No calls to system(3) supported on this platform. Attempted to call:\n'%s'\n",command);
122 return ( system(command) );
127 /* Check if string starts with substring */
128 static gmx_bool str_starts(const char *string, const char *substring)
130 return ( strncmp(string, substring, strlen(substring)) == 0);
134 static void cleandata(t_perf *perfdata, int test_nr)
136 perfdata->Gcycles[test_nr] = 0.0;
137 perfdata->ns_per_day[test_nr] = 0.0;
138 perfdata->PME_f_load[test_nr] = 0.0;
144 static gmx_bool is_equal(real a, real b)
146 real diff, eps=1.0e-7;
151 if (diff < 0.0) diff = -diff;
164 static void finalize(const char *fn_out)
170 fp = fopen(fn_out,"r");
171 fprintf(stdout,"\n\n");
173 while( fgets(buf,STRLEN-1,fp) != NULL )
175 fprintf(stdout,"%s",buf);
178 fprintf(stdout,"\n\n");
182 enum {eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr};
184 static int parse_logfile(const char *logfile, const char *errfile,
185 t_perf *perfdata, int test_nr, int presteps, gmx_large_int_t cpt_steps,
189 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
190 const char matchstrdd[]="Domain decomposition grid";
191 const char matchstrcr[]="resetting all time and cycle counters";
192 const char matchstrbal[]="Average PME mesh/force load:";
193 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";
194 const char errSIG[]="signal, stopping at the next";
197 float dum1,dum2,dum3,dum4;
200 gmx_large_int_t resetsteps=-1;
201 gmx_bool bFoundResetStr = FALSE;
202 gmx_bool bResetChecked = FALSE;
205 if (!gmx_fexist(logfile))
207 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
208 cleandata(perfdata, test_nr);
209 return eParselogNotFound;
212 fp = fopen(logfile, "r");
213 perfdata->PME_f_load[test_nr] = -1.0;
214 perfdata->guessPME = -1;
216 iFound = eFoundNothing;
219 iFound = eFoundDDStr; /* Skip some case statements */
222 while (fgets(line, STRLEN, fp) != NULL)
224 /* Remove leading spaces */
227 /* Check for TERM and INT signals from user: */
228 if ( strstr(line, errSIG) != NULL )
231 cleandata(perfdata, test_nr);
232 return eParselogTerm;
235 /* Check whether cycle resetting worked */
236 if (presteps > 0 && !bFoundResetStr)
238 if (strstr(line, matchstrcr) != NULL)
240 sprintf(dumstring, "step %s", gmx_large_int_pfmt);
241 sscanf(line, dumstring, &resetsteps);
242 bFoundResetStr = TRUE;
243 if (resetsteps == presteps+cpt_steps)
245 bResetChecked = TRUE;
249 sprintf(dumstring , gmx_large_int_pfmt, resetsteps);
250 sprintf(dumstring2, gmx_large_int_pfmt, presteps+cpt_steps);
251 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
252 " though they were supposed to be reset at step %s!\n",
253 dumstring, dumstring2);
258 /* Look for strings that appear in a certain order in the log file: */
262 /* Look for domain decomp grid and separate PME nodes: */
263 if (str_starts(line, matchstrdd))
265 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME nodes %d",
266 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
267 if (perfdata->nPMEnodes == -1)
269 perfdata->guessPME = npme;
271 else if (perfdata->nPMEnodes != npme)
273 gmx_fatal(FARGS, "PME nodes from command line and output file are not identical");
275 iFound = eFoundDDStr;
277 /* Catch a few errors that might have occured: */
278 else if (str_starts(line, "There is no domain decomposition for"))
281 return eParselogNoDDGrid;
283 else if (str_starts(line, "reading tpx file"))
286 return eParselogTPXVersion;
288 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
291 return eParselogNotParallel;
295 /* Look for PME mesh/force balance (not necessarily present, though) */
296 if (str_starts(line, matchstrbal))
298 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
300 /* Look for matchstring */
301 if (str_starts(line, matchstring))
303 iFound = eFoundAccountingStr;
306 case eFoundAccountingStr:
307 /* Already found matchstring - look for cycle data */
308 if (str_starts(line, "Total "))
310 sscanf(line,"Total %d %lf",&procs,&(perfdata->Gcycles[test_nr]));
311 iFound = eFoundCycleStr;
315 /* Already found cycle data - look for remaining performance info and return */
316 if (str_starts(line, "Performance:"))
318 ndum = sscanf(line,"%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
319 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
320 perfdata->ns_per_day[test_nr] = (ndum==5)? dum3 : dum1;
322 if (bResetChecked || presteps == 0)
328 return eParselogResetProblem;
335 /* Close the log file */
338 /* Check why there is no performance data in the log file.
339 * Did a fatal errors occur? */
340 if (gmx_fexist(errfile))
342 fp = fopen(errfile, "r");
343 while (fgets(line, STRLEN, fp) != NULL)
345 if ( str_starts(line, "Fatal error:") )
347 if (fgets(line, STRLEN, fp) != NULL)
349 fprintf(stderr, "\nWARNING: An error occured during this benchmark:\n"
353 cleandata(perfdata, test_nr);
354 return eParselogFatal;
361 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
364 /* Giving up ... we could not find out why there is no performance data in
366 fprintf(stdout, "No performance data in log file.\n");
367 cleandata(perfdata, test_nr);
369 return eParselogNoPerfData;
373 static gmx_bool analyze_data(
382 int *index_tpr, /* OUT: Nr of mdp file with best settings */
383 int *npme_optimal) /* OUT: Optimal number of PME nodes */
386 int line=0, line_win=-1;
387 int k_win=-1, i_win=-1, winPME;
388 double s=0.0; /* standard deviation */
391 char str_PME_f_load[13];
392 gmx_bool bCanUseOrigTPR;
393 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
399 fprintf(fp, "Summary of successful runs:\n");
400 fprintf(fp, "Line tpr PME nodes Gcycles Av. Std.dev. ns/day PME/f");
403 fprintf(fp, " DD grid");
409 for (k=0; k<ntprs; k++)
411 for (i=0; i<ntests; i++)
413 /* Select the right dataset: */
414 pd = &(perfdata[k][i]);
416 pd->Gcycles_Av = 0.0;
417 pd->PME_f_load_Av = 0.0;
418 pd->ns_per_day_Av = 0.0;
420 if (pd->nPMEnodes == -1)
422 sprintf(strbuf, "(%3d)", pd->guessPME);
426 sprintf(strbuf, " ");
429 /* Get the average run time of a setting */
430 for (j=0; j<nrepeats; j++)
432 pd->Gcycles_Av += pd->Gcycles[j];
433 pd->PME_f_load_Av += pd->PME_f_load[j];
435 pd->Gcycles_Av /= nrepeats;
436 pd->PME_f_load_Av /= nrepeats;
438 for (j=0; j<nrepeats; j++)
440 if (pd->ns_per_day[j] > 0.0)
442 pd->ns_per_day_Av += pd->ns_per_day[j];
446 /* Somehow the performance number was not aquired for this run,
447 * therefor set the average to some negative value: */
448 pd->ns_per_day_Av = -1.0f*nrepeats;
452 pd->ns_per_day_Av /= nrepeats;
455 if (pd->PME_f_load_Av > 0.0)
457 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
461 sprintf(str_PME_f_load, "%s", " - ");
465 /* We assume we had a successful run if both averages are positive */
466 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
468 /* Output statistics if repeats were done */
471 /* Calculate the standard deviation */
473 for (j=0; j<nrepeats; j++)
475 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
480 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
481 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
482 pd->ns_per_day_Av, str_PME_f_load);
485 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
489 /* Store the index of the best run found so far in 'winner': */
490 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
503 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
508 winPME = perfdata[k_win][i_win].nPMEnodes;
512 /* We stuck to a fixed number of PME-only nodes */
513 sprintf(strbuf, "settings No. %d", k_win);
517 /* We have optimized the number of PME-only nodes */
520 sprintf(strbuf, "%s", "the automatic number of PME nodes");
524 sprintf(strbuf, "%d PME nodes", winPME);
527 fprintf(fp, "Best performance was achieved with %s", strbuf);
528 if ((nrepeats > 1) && (ntests > 1))
530 fprintf(fp, " (see line %d)", line_win);
534 /* Only mention settings if they were modified: */
535 bRefinedCoul = !is_equal(info->rcoulomb[k_win], info->rcoulomb[0]);
536 bRefinedVdW = !is_equal(info->rvdw[k_win] , info->rvdw[0] );
537 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
538 info->nky[k_win] == info->nky[0] &&
539 info->nkz[k_win] == info->nkz[0]);
541 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
543 fprintf(fp, "Optimized PME settings:\n");
544 bCanUseOrigTPR = FALSE;
548 bCanUseOrigTPR = TRUE;
553 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
558 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
563 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],
564 info->nkx[0], info->nky[0], info->nkz[0]);
567 if (bCanUseOrigTPR && ntprs > 1)
569 fprintf(fp, "and original PME settings.\n");
574 /* Return the index of the mdp file that showed the highest performance
575 * and the optimal number of PME nodes */
577 *npme_optimal = winPME;
579 return bCanUseOrigTPR;
583 /* Get the commands we need to set up the runs from environment variables */
584 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char cmd_np[],
585 char *cmd_mdrun[], int repeats)
592 const char def_mpirun[] = "mpirun";
593 const char def_mdrun[] = "mdrun";
594 const char filename[] = "benchtest.log";
595 const char match_mpi[] = "NNODES=";
596 const char match_mdrun[]= "Program: ";
597 const char empty_mpirun[] = "";
598 gmx_bool bMdrun = FALSE;
599 gmx_bool bMPI = FALSE;
602 /* Get the commands we need to set up the runs from environment variables */
605 if ( (cp = getenv("MPIRUN")) != NULL)
607 *cmd_mpirun = strdup(cp);
611 *cmd_mpirun = strdup(def_mpirun);
616 *cmd_mpirun = strdup(empty_mpirun);
619 if ( (cp = getenv("MDRUN" )) != NULL )
621 *cmd_mdrun = strdup(cp);
625 *cmd_mdrun = strdup(def_mdrun);
629 /* If no simulations have to be performed, we are done here */
635 /* Run a small test to see whether mpirun + mdrun work */
636 fprintf(stdout, "Making sure that mdrun can be executed. ");
639 snew(command, strlen(*cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
640 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", *cmd_mdrun, cmd_np, filename);
644 snew(command, strlen(*cmd_mpirun) + strlen(cmd_np) + strlen(*cmd_mdrun) + strlen(filename) + 50);
645 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", *cmd_mpirun, cmd_np, *cmd_mdrun, filename);
647 fprintf(stdout, "Trying '%s' ... ", command);
648 make_backup(filename);
649 gmx_system_call(command);
651 /* Check if we find the characteristic string in the output: */
652 if (!gmx_fexist(filename))
653 gmx_fatal(FARGS, "Output from test run could not be found.");
655 fp = fopen(filename, "r");
656 /* We need to scan the whole output file, since sometimes the queuing system
657 * also writes stuff to stdout/err */
660 cp2=fgets(line, STRLEN, fp);
663 if ( str_starts(line, match_mdrun) )
667 if ( str_starts(line, match_mpi) )
679 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
681 "seems to have been compiled with MPI instead.",
689 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
691 "seems to have been compiled without MPI support.",
698 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
702 fprintf(stdout, "passed.\n");
710 static void launch_simulation(
711 gmx_bool bLaunch, /* Should the simulation be launched? */
712 FILE *fp, /* General log file */
713 gmx_bool bThreads, /* whether to use threads */
714 char *cmd_mpirun, /* Command for mpirun */
715 char *cmd_np, /* Switch for -np or -ntmpi or empty */
716 char *cmd_mdrun, /* Command for mdrun */
717 char *args_for_mdrun, /* Arguments for mdrun */
718 const char *simulation_tpr, /* This tpr will be simulated */
719 int nPMEnodes) /* Number of PME nodes to use */
724 /* Make enough space for the system call command,
725 * (100 extra chars for -npme ... etc. options should suffice): */
726 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+100);
728 /* Note that the -passall options requires args_for_mdrun to be at the end
729 * of the command line string */
732 sprintf(command, "%s%s-npme %d -s %s %s",
733 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun);
737 sprintf(command, "%s%s%s -npme %d -s %s %s",
738 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun);
741 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch? "Using":"Please use", command);
745 /* Now the real thing! */
748 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
751 gmx_system_call(command);
756 static void modify_PMEsettings(
757 gmx_large_int_t simsteps, /* Set this value as number of time steps */
758 gmx_large_int_t init_step, /* Set this value as init_step */
759 const char *fn_best_tpr, /* tpr file with the best performance */
760 const char *fn_sim_tpr) /* name of tpr file to be launched */
768 read_tpx_state(fn_best_tpr,ir,&state,NULL,&mtop);
770 /* Reset nsteps and init_step to the value of the input .tpr file */
771 ir->nsteps = simsteps;
772 ir->init_step = init_step;
774 /* Write the tpr file which will be launched */
775 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, gmx_large_int_pfmt);
776 fprintf(stdout,buf,ir->nsteps);
778 write_tpx_state(fn_sim_tpr,ir,&state,&mtop);
784 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
786 /* Make additional TPR files with more computational load for the
787 * direct space processors: */
788 static void make_benchmark_tprs(
789 const char *fn_sim_tpr, /* READ : User-provided tpr file */
790 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
791 gmx_large_int_t benchsteps, /* Number of time steps for benchmark runs */
792 gmx_large_int_t statesteps, /* Step counter in checkpoint file */
793 real rmin, /* Minimal Coulomb radius */
794 real rmax, /* Maximal Coulomb radius */
795 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
796 int *ntprs, /* No. of TPRs to write, each with a different
797 rcoulomb and fourierspacing */
798 t_inputinfo *info, /* Contains information about mdp file options */
799 FILE *fp) /* Write the output here */
805 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
808 gmx_bool bNote = FALSE;
809 real add; /* Add this to rcoul for the next test */
810 real fac = 1.0; /* Scaling factor for Coulomb radius */
811 real fourierspacing; /* Basic fourierspacing from tpr */
814 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
815 *ntprs > 1? "s":"", gmx_large_int_pfmt, benchsteps>1?"s":"");
816 fprintf(stdout, buf, benchsteps);
819 sprintf(buf, " (adding %s steps from checkpoint file)", gmx_large_int_pfmt);
820 fprintf(stdout, buf, statesteps);
821 benchsteps += statesteps;
823 fprintf(stdout, ".\n");
827 read_tpx_state(fn_sim_tpr,ir,&state,NULL,&mtop);
829 /* Check if some kind of PME was chosen */
830 if (EEL_PME(ir->coulombtype) == FALSE)
832 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
836 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
837 if ( (ir->cutoff_scheme != ecutsVERLET) &&
838 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
840 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
841 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
843 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
844 else if (ir->rcoulomb > ir->rlist)
846 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
847 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
850 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
852 fprintf(stdout,"NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
856 /* Reduce the number of steps for the benchmarks */
857 info->orig_sim_steps = ir->nsteps;
858 ir->nsteps = benchsteps;
859 /* We must not use init_step from the input tpr file for the benchmarks */
860 info->orig_init_step = ir->init_step;
863 /* For PME-switch potentials, keep the radial distance of the buffer region */
864 nlist_buffer = ir->rlist - ir->rcoulomb;
866 /* Determine length of triclinic box vectors */
872 box_size[d] += state.box[d][i]*state.box[d][i];
874 box_size[d] = sqrt(box_size[d]);
877 if (ir->fourier_spacing > 0)
879 info->fsx[0] = ir->fourier_spacing;
880 info->fsy[0] = ir->fourier_spacing;
881 info->fsz[0] = ir->fourier_spacing;
885 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
886 info->fsx[0] = box_size[XX]/ir->nkx;
887 info->fsy[0] = box_size[YY]/ir->nky;
888 info->fsz[0] = box_size[ZZ]/ir->nkz;
891 /* If no value for the fourierspacing was provided on the command line, we
892 * use the reconstruction from the tpr file */
893 if (ir->fourier_spacing > 0)
895 /* Use the spacing from the tpr */
896 fourierspacing = ir->fourier_spacing;
900 /* Use the maximum observed spacing */
901 fourierspacing = max(max(info->fsx[0],info->fsy[0]),info->fsz[0]);
904 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
906 /* For performance comparisons the number of particles is useful to have */
907 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
909 /* Print information about settings of which some are potentially modified: */
910 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
911 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
912 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
913 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
914 if (EVDW_SWITCHED(ir->vdwtype))
916 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
918 if (EPME_SWITCHED(ir->coulombtype))
920 fprintf(fp, " rlist : %f nm\n", ir->rlist);
922 if (ir->rlistlong != max_cutoff(ir->rvdw,ir->rcoulomb))
924 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
927 /* Print a descriptive line about the tpr settings tested */
928 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
929 fprintf(fp, " No. scaling rcoulomb");
930 fprintf(fp, " nkx nky nkz");
931 fprintf(fp, " spacing");
932 if (evdwCUT == ir->vdwtype)
934 fprintf(fp, " rvdw");
936 if (EPME_SWITCHED(ir->coulombtype))
938 fprintf(fp, " rlist");
940 if ( ir->rlistlong != max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb)) )
942 fprintf(fp, " rlistlong");
944 fprintf(fp, " tpr file\n");
946 /* Loop to create the requested number of tpr input files */
947 for (j = 0; j < *ntprs; j++)
949 /* The first .tpr is the provided one, just need to modify nsteps,
950 * so skip the following block */
953 /* Determine which Coulomb radii rc to use in the benchmarks */
954 add = (rmax-rmin)/(*ntprs-1);
955 if (is_equal(rmin,info->rcoulomb[0]))
957 ir->rcoulomb = rmin + j*add;
959 else if (is_equal(rmax,info->rcoulomb[0]))
961 ir->rcoulomb = rmin + (j-1)*add;
965 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
966 add = (rmax-rmin)/(*ntprs-2);
967 ir->rcoulomb = rmin + (j-1)*add;
970 /* Determine the scaling factor fac */
971 fac = ir->rcoulomb/info->rcoulomb[0];
973 /* Scale the Fourier grid spacing */
974 ir->nkx = ir->nky = ir->nkz = 0;
975 calc_grid(NULL,state.box,fourierspacing*fac,&ir->nkx,&ir->nky,&ir->nkz);
977 /* Adjust other radii since various conditions neet to be fulfilled */
978 if (eelPME == ir->coulombtype)
980 /* plain PME, rcoulomb must be equal to rlist */
981 ir->rlist = ir->rcoulomb;
985 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
986 ir->rlist = ir->rcoulomb + nlist_buffer;
989 if (bScaleRvdw && evdwCUT == ir->vdwtype)
991 /* For vdw cutoff, rvdw >= rlist */
992 ir->rvdw = max(info->rvdw[0], ir->rlist);
995 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
997 } /* end of "if (j != 0)" */
999 /* for j==0: Save the original settings
1000 * for j >0: Save modified radii and Fourier grids */
1001 info->rcoulomb[j] = ir->rcoulomb;
1002 info->rvdw[j] = ir->rvdw;
1003 info->nkx[j] = ir->nkx;
1004 info->nky[j] = ir->nky;
1005 info->nkz[j] = ir->nkz;
1006 info->rlist[j] = ir->rlist;
1007 info->rlistlong[j] = ir->rlistlong;
1008 info->fsx[j] = fac*fourierspacing;
1009 info->fsy[j] = fac*fourierspacing;
1010 info->fsz[j] = fac*fourierspacing;
1012 /* Write the benchmark tpr file */
1013 strncpy(fn_bench_tprs[j],fn_sim_tpr,strlen(fn_sim_tpr)-strlen(".tpr"));
1014 sprintf(buf, "_bench%.2d.tpr", j);
1015 strcat(fn_bench_tprs[j], buf);
1016 fprintf(stdout,"Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1017 fprintf(stdout, gmx_large_int_pfmt, ir->nsteps);
1020 fprintf(stdout,", scaling factor %f\n", fac);
1024 fprintf(stdout,", unmodified settings\n");
1027 write_tpx_state(fn_bench_tprs[j],ir,&state,&mtop);
1029 /* Write information about modified tpr settings to log file */
1030 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1031 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1032 fprintf(fp, " %9f ", info->fsx[j]);
1033 if (evdwCUT == ir->vdwtype)
1035 fprintf(fp, "%10f", ir->rvdw);
1037 if (EPME_SWITCHED(ir->coulombtype))
1039 fprintf(fp, "%10f", ir->rlist);
1041 if ( info->rlistlong[0] != max_cutoff(info->rlist[0],max_cutoff(info->rvdw[0],info->rcoulomb[0])) )
1043 fprintf(fp, "%10f", ir->rlistlong);
1045 fprintf(fp, " %-14s\n",fn_bench_tprs[j]);
1047 /* Make it clear to the user that some additional settings were modified */
1048 if ( !is_equal(ir->rvdw , info->rvdw[0])
1049 || !is_equal(ir->rlistlong, info->rlistlong[0]) )
1056 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1057 "other input settings were also changed (see table above).\n"
1058 "Please check if the modified settings are appropriate.\n");
1066 /* Rename the files we want to keep to some meaningful filename and
1067 * delete the rest */
1068 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1069 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1071 char numstring[STRLEN];
1072 char newfilename[STRLEN];
1073 const char *fn=NULL;
1078 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1080 for (i=0; i<nfile; i++)
1082 opt = (char *)fnm[i].opt;
1083 if ( strcmp(opt, "-p") == 0 )
1085 /* do nothing; keep this file */
1088 else if (strcmp(opt, "-bg") == 0)
1090 /* Give the log file a nice name so one can later see which parameters were used */
1091 numstring[0] = '\0';
1094 sprintf(numstring, "_%d", nr);
1096 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg",nfile,fnm), k, nnodes, nPMEnodes, numstring);
1097 if (gmx_fexist(opt2fn("-bg",nfile,fnm)))
1099 fprintf(stdout, "renaming log file to %s\n", newfilename);
1100 make_backup(newfilename);
1101 rename(opt2fn("-bg",nfile,fnm), newfilename);
1104 else if (strcmp(opt, "-err") == 0)
1106 /* This file contains the output of stderr. We want to keep it in
1107 * cases where there have been problems. */
1108 fn = opt2fn(opt, nfile, fnm);
1109 numstring[0] = '\0';
1112 sprintf(numstring, "_%d", nr);
1114 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1119 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1120 make_backup(newfilename);
1121 rename(fn, newfilename);
1125 fprintf(stdout, "Deleting %s\n", fn);
1130 /* Delete the files which are created for each benchmark run: (options -b*) */
1131 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt,nfile,fnm) || !is_optional(&fnm[i])) )
1133 fn = opt2fn(opt, nfile, fnm);
1136 fprintf(stdout, "Deleting %s\n", fn);
1144 /* Returns the largest common factor of n1 and n2 */
1145 static int largest_common_factor(int n1, int n2)
1150 for (factor=nmax; factor > 0; factor--)
1152 if ( 0==(n1 % factor) && 0==(n2 % factor) )
1157 return 0; /* one for the compiler */
1160 enum {eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr};
1162 /* Create a list of numbers of PME nodes to test */
1163 static void make_npme_list(
1164 const char *npmevalues_opt, /* Make a complete list with all
1165 * possibilities or a short list that keeps only
1166 * reasonable numbers of PME nodes */
1167 int *nentries, /* Number of entries we put in the nPMEnodes list */
1168 int *nPMEnodes[], /* Each entry contains the value for -npme */
1169 int nnodes, /* Total number of nodes to do the tests on */
1170 int minPMEnodes, /* Minimum number of PME nodes */
1171 int maxPMEnodes) /* Maximum number of PME nodes */
1174 int min_factor=1; /* We request that npp and npme have this minimal
1175 * largest common factor (depends on npp) */
1176 int nlistmax; /* Max. list size */
1177 int nlist; /* Actual number of entries in list */
1181 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1182 if ( 0 == strcmp(npmevalues_opt, "all") )
1186 else if ( 0 == strcmp(npmevalues_opt, "subset") )
1188 eNPME = eNpmeSubset;
1190 else /* "auto" or "range" */
1196 else if (nnodes < 128)
1198 eNPME = eNpmeReduced;
1202 eNPME = eNpmeSubset;
1206 /* Calculate how many entries we could possibly have (in case of -npme all) */
1209 nlistmax = maxPMEnodes - minPMEnodes + 3;
1210 if (0 == minPMEnodes)
1220 /* Now make the actual list which is at most of size nlist */
1221 snew(*nPMEnodes, nlistmax);
1222 nlist = 0; /* start counting again, now the real entries in the list */
1223 for (i = 0; i < nlistmax - 2; i++)
1225 npme = maxPMEnodes - i;
1236 /* For 2d PME we want a common largest factor of at least the cube
1237 * root of the number of PP nodes */
1238 min_factor = (int) pow(npp, 1.0/3.0);
1241 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1244 if (largest_common_factor(npp, npme) >= min_factor)
1246 (*nPMEnodes)[nlist] = npme;
1250 /* We always test 0 PME nodes and the automatic number */
1251 *nentries = nlist + 2;
1252 (*nPMEnodes)[nlist ] = 0;
1253 (*nPMEnodes)[nlist+1] = -1;
1255 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1256 for (i=0; i<*nentries-1; i++)
1258 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1260 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1264 /* Allocate memory to store the performance data */
1265 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1270 for (k=0; k<ntprs; k++)
1272 snew(perfdata[k], datasets);
1273 for (i=0; i<datasets; i++)
1275 for (j=0; j<repeats; j++)
1277 snew(perfdata[k][i].Gcycles , repeats);
1278 snew(perfdata[k][i].ns_per_day, repeats);
1279 snew(perfdata[k][i].PME_f_load, repeats);
1286 /* Check for errors on mdrun -h */
1287 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp)
1289 char *command, *msg;
1293 snew(command, length + 15);
1294 snew(msg , length + 500);
1296 fprintf(stdout, "Making shure the benchmarks can be executed ...\n");
1297 sprintf(command, "%s-h -quiet", mdrun_cmd_line);
1298 ret = gmx_system_call(command);
1302 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1303 * get the error message from mdrun itself */
1304 sprintf(msg, "Cannot run the benchmark simulations! Please check the error message of\n"
1305 "mdrun for the source of the problem. Did you provide a command line\n"
1306 "argument that neither g_tune_pme nor mdrun understands? Offending command:\n"
1307 "\n%s\n\n", command);
1309 fprintf(stderr, "%s", msg);
1311 fprintf(fp , "%s", msg);
1321 static void do_the_tests(
1322 FILE *fp, /* General g_tune_pme output file */
1323 char **tpr_names, /* Filenames of the input files to test */
1324 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1325 int minPMEnodes, /* Min fraction of nodes to use for PME */
1326 int npme_fixed, /* If >= -1, test fixed number of PME
1328 const char *npmevalues_opt, /* Which -npme values should be tested */
1329 t_perf **perfdata, /* Here the performace data is stored */
1330 int *pmeentries, /* Entries in the nPMEnodes list */
1331 int repeats, /* Repeat each test this often */
1332 int nnodes, /* Total number of nodes = nPP + nPME */
1333 int nr_tprs, /* Total number of tpr files to test */
1334 gmx_bool bThreads, /* Threads or MPI? */
1335 char *cmd_mpirun, /* mpirun command string */
1336 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1337 char *cmd_mdrun, /* mdrun command string */
1338 char *cmd_args_bench, /* arguments for mdrun in a string */
1339 const t_filenm *fnm, /* List of filenames from command line */
1340 int nfile, /* Number of files specified on the cmdl. */
1341 int presteps, /* DLB equilibration steps, is checked */
1342 gmx_large_int_t cpt_steps) /* Time step counter in the checkpoint */
1344 int i,nr,k,ret,count=0,totaltests;
1345 int *nPMEnodes=NULL;
1348 char *command, *cmd_stub;
1350 gmx_bool bResetProblem=FALSE;
1351 gmx_bool bFirst=TRUE;
1354 /* This string array corresponds to the eParselog enum type at the start
1356 const char* ParseLog[] = {"OK.",
1357 "Logfile not found!",
1358 "No timings, logfile truncated?",
1359 "Run was terminated.",
1360 "Counters were not reset properly.",
1361 "No DD grid found for these settings.",
1362 "TPX version conflict!",
1363 "mdrun was not started in parallel!",
1364 "An error occured." };
1365 char str_PME_f_load[13];
1368 /* Allocate space for the mdrun command line. 100 extra characters should
1369 be more than enough for the -npme etcetera arguments */
1370 cmdline_length = strlen(cmd_mpirun)
1373 + strlen(cmd_args_bench)
1374 + strlen(tpr_names[0]) + 100;
1375 snew(command , cmdline_length);
1376 snew(cmd_stub, cmdline_length);
1378 /* Construct the part of the command line that stays the same for all tests: */
1381 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1385 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1388 /* Create a list of numbers of PME nodes to test */
1389 if (npme_fixed < -1)
1391 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1392 nnodes, minPMEnodes, maxPMEnodes);
1398 nPMEnodes[0] = npme_fixed;
1399 fprintf(stderr, "Will use a fixed number of %d PME-only nodes.\n", nPMEnodes[0]);
1404 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1406 finalize(opt2fn("-p", nfile, fnm));
1410 /* Allocate one dataset for each tpr input file: */
1411 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1413 /*****************************************/
1414 /* Main loop over all tpr files to test: */
1415 /*****************************************/
1416 totaltests = nr_tprs*(*pmeentries)*repeats;
1417 for (k=0; k<nr_tprs;k++)
1419 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1420 fprintf(fp, "PME nodes Gcycles ns/day PME/f Remark\n");
1421 /* Loop over various numbers of PME nodes: */
1422 for (i = 0; i < *pmeentries; i++)
1424 pd = &perfdata[k][i];
1426 /* Loop over the repeats for each scenario: */
1427 for (nr = 0; nr < repeats; nr++)
1429 pd->nPMEnodes = nPMEnodes[i];
1431 /* Add -npme and -s to the command line and save it. Note that
1432 * the -passall (if set) options requires cmd_args_bench to be
1433 * at the end of the command line string */
1434 snew(pd->mdrun_cmd_line, cmdline_length);
1435 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s",
1436 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench);
1438 /* To prevent that all benchmarks fail due to a show-stopper argument
1439 * on the mdrun command line, we make a quick check with mdrun -h first */
1442 make_sure_it_runs(pd->mdrun_cmd_line, cmdline_length, fp);
1446 /* Do a benchmark simulation: */
1449 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1455 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1456 (100.0*count)/totaltests,
1457 k+1, nr_tprs, i+1, *pmeentries, buf);
1458 make_backup(opt2fn("-err",nfile,fnm));
1459 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err",nfile,fnm));
1460 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1461 gmx_system_call(command);
1463 /* Collect the performance data from the log file; also check stderr
1464 * for fatal errors */
1465 ret = parse_logfile(opt2fn("-bg",nfile,fnm), opt2fn("-err",nfile,fnm),
1466 pd, nr, presteps, cpt_steps, nnodes);
1467 if ((presteps > 0) && (ret == eParselogResetProblem))
1469 bResetProblem = TRUE;
1472 if (-1 == pd->nPMEnodes)
1474 sprintf(buf, "(%3d)", pd->guessPME);
1482 if (pd->PME_f_load[nr] > 0.0)
1484 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1488 sprintf(str_PME_f_load, "%s", " - ");
1491 /* Write the data we got to disk */
1492 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1493 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1494 if (! (ret==eParselogOK || ret==eParselogNoDDGrid || ret==eParselogNotFound) )
1496 fprintf(fp, " Check %s file for problems.", ret==eParselogFatal? "err":"log");
1502 /* Do some cleaning up and delete the files we do not need any more */
1503 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret==eParselogFatal);
1505 /* If the first run with this number of processors already failed, do not try again: */
1506 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1508 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1509 count += repeats-(nr+1);
1512 } /* end of repeats loop */
1513 } /* end of -npme loop */
1514 } /* end of tpr file loop */
1519 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1527 static void check_input(
1534 real maxPMEfraction,
1535 real minPMEfraction,
1537 gmx_large_int_t bench_nsteps,
1538 const t_filenm *fnm,
1548 /* Make sure the input file exists */
1549 if (!gmx_fexist(opt2fn("-s",nfile,fnm)))
1551 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s",nfile,fnm));
1554 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1555 if ( (0 == strcmp(opt2fn("-cpi",nfile,fnm), opt2fn("-bcpo",nfile,fnm)) ) && (sim_part > 1) )
1557 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1558 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1561 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1564 gmx_fatal(FARGS, "Number of repeats < 0!");
1567 /* Check number of nodes */
1570 gmx_fatal(FARGS, "Number of nodes/threads must be a positive integer.");
1573 /* Automatically choose -ntpr if not set */
1583 /* Set a reasonable scaling factor for rcoulomb */
1585 *rmax = rcoulomb * 1.2;
1587 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs==1?"":"s");
1593 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1597 /* Make shure that rmin <= rcoulomb <= rmax */
1606 if ( !(*rmin <= *rmax) )
1608 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1609 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1611 /* Add test scenarios if rmin or rmax were set */
1614 if ( !is_equal(*rmin, rcoulomb) && (*ntprs == 1) )
1617 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1620 if (!is_equal(*rmax, rcoulomb) && (*ntprs == 1) )
1623 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1628 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1629 if ( !is_equal(*rmax, rcoulomb) || !is_equal(*rmin, rcoulomb) )
1631 *ntprs = max(*ntprs, 2);
1634 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1635 if ( !is_equal(*rmax, rcoulomb) && !is_equal(*rmin, rcoulomb) )
1637 *ntprs = max(*ntprs, 3);
1642 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1647 if (is_equal(*rmin,rcoulomb) && is_equal(rcoulomb,*rmax)) /* We have just a single rc */
1649 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1650 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1651 "with correspondingly adjusted PME grid settings\n");
1656 /* Check whether max and min fraction are within required values */
1657 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1659 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1661 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1663 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1665 if (maxPMEfraction < minPMEfraction)
1667 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1670 /* Check whether the number of steps - if it was set - has a reasonable value */
1671 if (bench_nsteps < 0)
1673 gmx_fatal(FARGS, "Number of steps must be positive.");
1676 if (bench_nsteps > 10000 || bench_nsteps < 100)
1678 fprintf(stderr, "WARNING: steps=");
1679 fprintf(stderr, gmx_large_int_pfmt, bench_nsteps);
1680 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100)? "few" : "many");
1685 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1688 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1691 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1693 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1697 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1698 * only. We need to check whether the requested number of PME-only nodes
1700 if (npme_fixed > -1)
1702 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1703 if (2*npme_fixed > nnodes)
1705 gmx_fatal(FARGS, "Cannot have more than %d PME-only nodes for a total of %d nodes (you chose %d).\n",
1706 nnodes/2, nnodes, npme_fixed);
1708 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1710 fprintf(stderr, "WARNING: Only %g percent of the nodes are assigned as PME-only nodes.\n",
1711 100.0*((real)npme_fixed / (real)nnodes));
1713 if (opt2parg_bSet("-min",npargs,pa) || opt2parg_bSet("-max",npargs,pa))
1715 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1716 " fixed number of PME-only nodes is requested with -fix.\n");
1722 /* Returns TRUE when "opt" is needed at launch time */
1723 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1725 /* Apart from the input .tpr and the error log we need all options that were set
1726 * on the command line and that do not start with -b */
1727 if (0 == strncmp(opt,"-b", 2) || 0 == strncmp(opt,"-s", 2) || 0 == strncmp(opt,"-err", 4))
1736 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1737 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1739 /* Apart from the input .tpr, all files starting with "-b" are for
1740 * _b_enchmark files exclusively */
1741 if (0 == strncmp(opt,"-s", 2))
1746 if (0 == strncmp(opt,"-b", 2) || 0 == strncmp(opt,"-s", 2))
1748 if (!bOptional || bSet)
1765 if (bSet) /* These are additional input files like -cpi -ei */
1778 /* Adds 'buf' to 'str' */
1779 static void add_to_string(char **str, char *buf)
1784 len = strlen(*str) + strlen(buf) + 1;
1790 /* Create the command line for the benchmark as well as for the real run */
1791 static void create_command_line_snippets(
1792 gmx_bool bAppendFiles,
1793 gmx_bool bKeepAndNumCPT,
1794 gmx_bool bResetHWay,
1798 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1799 char *cmd_args_launch[], /* command line arguments for simulation run */
1800 char extra_args[]) /* Add this to the end of the command line */
1805 char strbuf[STRLEN];
1808 /* strlen needs at least '\0' as a string: */
1809 snew(*cmd_args_bench ,1);
1810 snew(*cmd_args_launch,1);
1811 *cmd_args_launch[0]='\0';
1812 *cmd_args_bench[0] ='\0';
1815 /*******************************************/
1816 /* 1. Process other command line arguments */
1817 /*******************************************/
1820 /* Add equilibration steps to benchmark options */
1821 sprintf(strbuf, "-resetstep %d ", presteps);
1822 add_to_string(cmd_args_bench, strbuf);
1824 /* These switches take effect only at launch time */
1825 if (FALSE == bAppendFiles)
1827 add_to_string(cmd_args_launch, "-noappend ");
1831 add_to_string(cmd_args_launch, "-cpnum ");
1835 add_to_string(cmd_args_launch, "-resethway ");
1838 /********************/
1839 /* 2. Process files */
1840 /********************/
1841 for (i=0; i<nfile; i++)
1843 opt = (char *)fnm[i].opt;
1844 name = opt2fn(opt,nfile,fnm);
1846 /* Strbuf contains the options, now let's sort out where we need that */
1847 sprintf(strbuf, "%s %s ", opt, name);
1849 if ( is_bench_file(opt, opt2bSet(opt,nfile,fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1851 /* All options starting with -b* need the 'b' removed,
1852 * therefore overwrite strbuf */
1853 if (0 == strncmp(opt, "-b", 2))
1855 sprintf(strbuf, "-%s %s ", &opt[2], name);
1858 add_to_string(cmd_args_bench, strbuf);
1861 if ( is_launch_file(opt,opt2bSet(opt,nfile,fnm)) )
1863 add_to_string(cmd_args_launch, strbuf);
1867 add_to_string(cmd_args_bench , extra_args);
1868 add_to_string(cmd_args_launch, extra_args);
1872 /* Set option opt */
1873 static void setopt(const char *opt,int nfile,t_filenm fnm[])
1877 for(i=0; (i<nfile); i++)
1879 if (strcmp(opt,fnm[i].opt)==0)
1881 fnm[i].flag |= ffSET;
1887 /* This routine inspects the tpr file and ...
1888 * 1. checks for output files that get triggered by a tpr option. These output
1889 * files are marked as 'set' to allow for a proper cleanup after each
1891 * 2. returns the PME:PP load ratio
1892 * 3. returns rcoulomb from the tpr */
1893 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1895 gmx_bool bPull; /* Is pulling requested in .tpr file? */
1896 gmx_bool bTpi; /* Is test particle insertion requested? */
1897 gmx_bool bFree; /* Is a free energy simulation requested? */
1898 gmx_bool bNM; /* Is a normal mode analysis requested? */
1904 /* Check tpr file for options that trigger extra output files */
1905 read_tpx_state(opt2fn("-s",nfile,fnm),&ir,&state,NULL,&mtop);
1906 bPull = (epullNO != ir.ePull);
1907 bFree = (efepNO != ir.efep );
1908 bNM = (eiNM == ir.eI );
1909 bTpi = EI_TPI(ir.eI);
1911 /* Set these output files on the tuning command-line */
1914 setopt("-pf" , nfile, fnm);
1915 setopt("-px" , nfile, fnm);
1919 setopt("-dhdl", nfile, fnm);
1923 setopt("-tpi" , nfile, fnm);
1924 setopt("-tpid", nfile, fnm);
1928 setopt("-mtx" , nfile, fnm);
1931 *rcoulomb = ir.rcoulomb;
1933 /* Return the estimate for the number of PME nodes */
1934 return pme_load_estimate(&mtop,&ir,state.box);
1938 static void couple_files_options(int nfile, t_filenm fnm[])
1941 gmx_bool bSet,bBench;
1946 for (i=0; i<nfile; i++)
1948 opt = (char *)fnm[i].opt;
1949 bSet = ((fnm[i].flag & ffSET) != 0);
1950 bBench = (0 == strncmp(opt,"-b", 2));
1952 /* Check optional files */
1953 /* If e.g. -eo is set, then -beo also needs to be set */
1954 if (is_optional(&fnm[i]) && bSet && !bBench)
1956 sprintf(buf, "-b%s", &opt[1]);
1957 setopt(buf,nfile,fnm);
1959 /* If -beo is set, then -eo also needs to be! */
1960 if (is_optional(&fnm[i]) && bSet && bBench)
1962 sprintf(buf, "-%s", &opt[2]);
1963 setopt(buf,nfile,fnm);
1969 static double gettime()
1971 #ifdef HAVE_GETTIMEOFDAY
1975 gettimeofday(&t,NULL);
1977 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
1983 seconds = time(NULL);
1990 #define BENCHSTEPS (1000)
1992 int gmx_tune_pme(int argc,char *argv[])
1994 const char *desc[] = {
1995 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of processors/threads, this program systematically",
1996 "times [TT]mdrun[tt] with various numbers of PME-only nodes and determines",
1997 "which setting is fastest. It will also test whether performance can",
1998 "be enhanced by shifting load from the reciprocal to the real space",
1999 "part of the Ewald sum. ",
2000 "Simply pass your [TT].tpr[tt] file to [TT]g_tune_pme[tt] together with other options",
2001 "for [TT]mdrun[tt] as needed.[PAR]",
2002 "Which executables are used can be set in the environment variables",
2003 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
2004 "will be used as defaults. Note that for certain MPI frameworks you",
2005 "need to provide a machine- or hostfile. This can also be passed",
2006 "via the MPIRUN variable, e.g.[PAR]",
2007 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
2008 "Please call [TT]g_tune_pme[tt] with the normal options you would pass to",
2009 "[TT]mdrun[tt] and add [TT]-np[tt] for the number of processors to perform the",
2010 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2011 "to repeat each test several times to get better statistics. [PAR]",
2012 "[TT]g_tune_pme[tt] can test various real space / reciprocal space workloads",
2013 "for you. With [TT]-ntpr[tt] you control how many extra [TT].tpr[tt] files will be",
2014 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2015 "Typically, the first test (number 0) will be with the settings from the input",
2016 "[TT].tpr[tt] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2017 "specified by [TT]-rmax[tt] with a somwhat smaller PME grid at the same time. ",
2018 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2019 "The remaining [TT].tpr[tt] files will have equally-spaced Coulomb radii (and Fourier "
2020 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2021 "if you just seek the optimal number of PME-only nodes; in that case",
2022 "your input [TT].tpr[tt] file will remain unchanged.[PAR]",
2023 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2024 "MD systems. The dynamic load balancing needs about 100 time steps",
2025 "to adapt to local load imbalances, therefore the time step counters",
2026 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2027 "for a higher accuarcy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2028 "From the 'DD' load imbalance entries in the md.log output file you",
2029 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2030 "[TT]g_tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2031 "After calling [TT]mdrun[tt] several times, detailed performance information",
2032 "is available in the output file [TT]perf.out.[tt] ",
2033 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2034 "(options [TT]-b[tt]*), these will be automatically deleted after each test.[PAR]",
2035 "If you want the simulation to be started automatically with the",
2036 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2041 int pmeentries=0; /* How many values for -npme do we actually test for each tpr file */
2042 real maxPMEfraction=0.50;
2043 real minPMEfraction=0.25;
2044 int maxPMEnodes, minPMEnodes;
2045 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2046 float guessPMEnodes;
2047 int npme_fixed=-2; /* If >= -1, use only this number
2048 * of PME-only nodes */
2050 real rmin=0.0,rmax=0.0; /* min and max value for rcoulomb if scaling is requested */
2051 real rcoulomb=-1.0; /* Coulomb radius as set in .tpr file */
2052 gmx_bool bScaleRvdw=TRUE;
2053 gmx_large_int_t bench_nsteps=BENCHSTEPS;
2054 gmx_large_int_t new_sim_nsteps=-1; /* -1 indicates: not set by the user */
2055 gmx_large_int_t cpt_steps=0; /* Step counter in .cpt input file */
2056 int presteps=100; /* Do a full cycle reset after presteps steps */
2057 gmx_bool bOverwrite=FALSE, bKeepTPR;
2058 gmx_bool bLaunch=FALSE;
2059 char *ExtraArgs=NULL;
2060 char **tpr_names=NULL;
2061 const char *simulation_tpr=NULL;
2062 int best_npme, best_tpr;
2063 int sim_part = 1; /* For benchmarks with checkpoint files */
2066 /* Default program names if nothing else is found */
2067 char *cmd_mpirun=NULL, *cmd_mdrun=NULL;
2068 char *cmd_args_bench, *cmd_args_launch;
2071 t_perf **perfdata=NULL;
2077 /* Print out how long the tuning took */
2080 static t_filenm fnm[] = {
2082 { efOUT, "-p", "perf", ffWRITE },
2083 { efLOG, "-err", "bencherr", ffWRITE },
2084 { efTPX, "-so", "tuned", ffWRITE },
2086 { efTPX, NULL, NULL, ffREAD },
2087 { efTRN, "-o", NULL, ffWRITE },
2088 { efXTC, "-x", NULL, ffOPTWR },
2089 { efCPT, "-cpi", NULL, ffOPTRD },
2090 { efCPT, "-cpo", NULL, ffOPTWR },
2091 { efSTO, "-c", "confout", ffWRITE },
2092 { efEDR, "-e", "ener", ffWRITE },
2093 { efLOG, "-g", "md", ffWRITE },
2094 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2095 { efXVG, "-field", "field", ffOPTWR },
2096 { efXVG, "-table", "table", ffOPTRD },
2097 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2098 { efXVG, "-tablep", "tablep", ffOPTRD },
2099 { efXVG, "-tableb", "table", ffOPTRD },
2100 { efTRX, "-rerun", "rerun", ffOPTRD },
2101 { efXVG, "-tpi", "tpi", ffOPTWR },
2102 { efXVG, "-tpid", "tpidist", ffOPTWR },
2103 { efEDI, "-ei", "sam", ffOPTRD },
2104 { efXVG, "-eo", "edsam", ffOPTWR },
2105 { efGCT, "-j", "wham", ffOPTRD },
2106 { efGCT, "-jo", "bam", ffOPTWR },
2107 { efXVG, "-ffout", "gct", ffOPTWR },
2108 { efXVG, "-devout", "deviatie", ffOPTWR },
2109 { efXVG, "-runav", "runaver", ffOPTWR },
2110 { efXVG, "-px", "pullx", ffOPTWR },
2111 { efXVG, "-pf", "pullf", ffOPTWR },
2112 { efXVG, "-ro", "rotation", ffOPTWR },
2113 { efLOG, "-ra", "rotangles",ffOPTWR },
2114 { efLOG, "-rs", "rotslabs", ffOPTWR },
2115 { efLOG, "-rt", "rottorque",ffOPTWR },
2116 { efMTX, "-mtx", "nm", ffOPTWR },
2117 { efNDX, "-dn", "dipole", ffOPTWR },
2118 /* Output files that are deleted after each benchmark run */
2119 { efTRN, "-bo", "bench", ffWRITE },
2120 { efXTC, "-bx", "bench", ffWRITE },
2121 { efCPT, "-bcpo", "bench", ffWRITE },
2122 { efSTO, "-bc", "bench", ffWRITE },
2123 { efEDR, "-be", "bench", ffWRITE },
2124 { efLOG, "-bg", "bench", ffWRITE },
2125 { efXVG, "-beo", "benchedo", ffOPTWR },
2126 { efXVG, "-bdhdl", "benchdhdl",ffOPTWR },
2127 { efXVG, "-bfield", "benchfld" ,ffOPTWR },
2128 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2129 { efXVG, "-btpid", "benchtpid",ffOPTWR },
2130 { efGCT, "-bjo", "bench", ffOPTWR },
2131 { efXVG, "-bffout", "benchgct", ffOPTWR },
2132 { efXVG, "-bdevout","benchdev", ffOPTWR },
2133 { efXVG, "-brunav", "benchrnav",ffOPTWR },
2134 { efXVG, "-bpx", "benchpx", ffOPTWR },
2135 { efXVG, "-bpf", "benchpf", ffOPTWR },
2136 { efXVG, "-bro", "benchrot", ffOPTWR },
2137 { efLOG, "-bra", "benchrota",ffOPTWR },
2138 { efLOG, "-brs", "benchrots",ffOPTWR },
2139 { efLOG, "-brt", "benchrott",ffOPTWR },
2140 { efMTX, "-bmtx", "benchn", ffOPTWR },
2141 { efNDX, "-bdn", "bench", ffOPTWR }
2144 gmx_bool bThreads = FALSE;
2148 const char *procstring[] =
2149 { NULL, "-np", "-n", "none", NULL };
2150 const char *npmevalues_opt[] =
2151 { NULL, "auto", "all", "subset", NULL };
2153 gmx_bool bAppendFiles=TRUE;
2154 gmx_bool bKeepAndNumCPT=FALSE;
2155 gmx_bool bResetCountersHalfWay=FALSE;
2156 gmx_bool bBenchmark=TRUE;
2158 output_env_t oenv=NULL;
2161 /***********************/
2162 /* g_tune_pme options: */
2163 /***********************/
2164 { "-np", FALSE, etINT, {&nnodes},
2165 "Number of nodes to run the tests on (must be > 2 for separate PME nodes)" },
2166 { "-npstring", FALSE, etENUM, {procstring},
2167 "Specify the number of processors to [TT]$MPIRUN[tt] using this string"},
2168 { "-ntmpi", FALSE, etINT, {&nthreads},
2169 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2170 { "-r", FALSE, etINT, {&repeats},
2171 "Repeat each test this often" },
2172 { "-max", FALSE, etREAL, {&maxPMEfraction},
2173 "Max fraction of PME nodes to test with" },
2174 { "-min", FALSE, etREAL, {&minPMEfraction},
2175 "Min fraction of PME nodes to test with" },
2176 { "-npme", FALSE, etENUM, {npmevalues_opt},
2177 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2178 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2179 { "-fix", FALSE, etINT, {&npme_fixed},
2180 "If >= -1, do not vary the number of PME-only nodes, instead use this fixed value and only vary rcoulomb and the PME grid spacing."},
2181 { "-rmax", FALSE, etREAL, {&rmax},
2182 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2183 { "-rmin", FALSE, etREAL, {&rmin},
2184 "If >0, minimal rcoulomb for -ntpr>1" },
2185 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2186 "Scale rvdw along with rcoulomb"},
2187 { "-ntpr", FALSE, etINT, {&ntprs},
2188 "Number of [TT].tpr[tt] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2189 "If < 1, automatically choose the number of [TT].tpr[tt] files to test" },
2190 { "-steps", FALSE, etGMX_LARGE_INT, {&bench_nsteps},
2191 "Take timings for this many steps in the benchmark runs" },
2192 { "-resetstep",FALSE, etINT, {&presteps},
2193 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2194 { "-simsteps", FALSE, etGMX_LARGE_INT, {&new_sim_nsteps},
2195 "If non-negative, perform this many steps in the real run (overwrites nsteps from [TT].tpr[tt], add [TT].cpt[tt] steps)" },
2196 { "-launch", FALSE, etBOOL, {&bLaunch},
2197 "Launch the real simulation after optimization" },
2198 { "-bench", FALSE, etBOOL, {&bBenchmark},
2199 "Run the benchmarks or just create the input [TT].tpr[tt] files?" },
2200 /******************/
2201 /* mdrun options: */
2202 /******************/
2203 /* We let g_tune_pme parse and understand these options, because we need to
2204 * prevent that they appear on the mdrun command line for the benchmarks */
2205 { "-append", FALSE, etBOOL, {&bAppendFiles},
2206 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2207 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2208 "Keep and number checkpoint files (launch only)" },
2209 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2210 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2214 #define NFILE asize(fnm)
2216 CopyRight(stderr,argv[0]);
2218 seconds = gettime();
2220 parse_common_args(&argc,argv,PCA_NOEXIT_ON_ARGS,
2221 NFILE,fnm,asize(pa),pa,asize(desc),desc,
2224 /* Store the remaining unparsed command line entries in a string which
2225 * is then attached to the mdrun command line */
2227 ExtraArgs[0] = '\0';
2228 for (i=1; i<argc; i++) /* argc will now be 1 if everything was understood */
2230 add_to_string(&ExtraArgs, argv[i]);
2231 add_to_string(&ExtraArgs, " ");
2234 if (opt2parg_bSet("-ntmpi",asize(pa),pa))
2237 if (opt2parg_bSet("-npstring",asize(pa),pa))
2239 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2244 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2246 /* and now we just set this; a bit of an ugly hack*/
2249 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2250 guessPMEratio = inspect_tpr(NFILE,fnm,&rcoulomb);
2252 /* Automatically set -beo options if -eo is set etc. */
2253 couple_files_options(NFILE,fnm);
2255 /* Construct the command line arguments for benchmark runs
2256 * as well as for the simulation run */
2259 sprintf(bbuf," -ntmpi %d ", nthreads);
2263 sprintf(bbuf," -np %d ", nnodes);
2268 create_command_line_snippets(bAppendFiles,bKeepAndNumCPT,bResetCountersHalfWay,presteps,
2269 NFILE,fnm,&cmd_args_bench, &cmd_args_launch, ExtraArgs);
2271 /* Read in checkpoint file if requested */
2273 if(opt2bSet("-cpi",NFILE,fnm))
2276 cr->duty=DUTY_PP; /* makes the following routine happy */
2277 read_checkpoint_simulation_part(opt2fn("-cpi",NFILE,fnm),
2278 &sim_part,&cpt_steps,cr,
2279 FALSE,NFILE,fnm,NULL,NULL);
2282 /* sim_part will now be 1 if no checkpoint file was found */
2285 gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi",NFILE,fnm));
2289 /* Open performance output file and write header info */
2290 fp = ffopen(opt2fn("-p",NFILE,fnm),"w");
2292 /* Make a quick consistency check of command line parameters */
2293 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2294 maxPMEfraction, minPMEfraction, npme_fixed,
2295 bench_nsteps, fnm, NFILE, sim_part, presteps,
2298 /* Determine the maximum and minimum number of PME nodes to test,
2299 * the actual list of settings is build in do_the_tests(). */
2300 if ((nnodes > 2) && (npme_fixed < -1))
2302 if (0 == strcmp(npmevalues_opt[0], "auto"))
2304 /* Determine the npme range automatically based on the PME:PP load guess */
2305 if (guessPMEratio > 1.0)
2307 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2308 maxPMEnodes=nnodes/2;
2309 minPMEnodes=nnodes/2;
2313 /* PME : PP load is in the range 0..1, let's test around the guess */
2314 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2315 minPMEnodes = floor(0.7*guessPMEnodes);
2316 maxPMEnodes = ceil(1.6*guessPMEnodes);
2317 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2322 /* Determine the npme range based on user input */
2323 maxPMEnodes = floor(maxPMEfraction*nnodes);
2324 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2325 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2326 if (maxPMEnodes != minPMEnodes)
2328 fprintf(stdout, "- %d ", maxPMEnodes);
2330 fprintf(stdout, "PME-only nodes.\n Note that the automatic number of PME-only nodes and no separate PME nodes are always tested.\n");
2339 /* Get the commands we need to set up the runs from environment variables */
2340 get_program_paths(bThreads, &cmd_mpirun, cmd_np, &cmd_mdrun, repeats);
2342 /* Print some header info to file */
2344 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2346 fprintf(fp, "%s for Gromacs %s\n", ShortProgram(),GromacsVersion());
2349 fprintf(fp, "Number of nodes : %d\n", nnodes);
2350 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2351 if ( strcmp(procstring[0], "none") != 0)
2353 fprintf(fp, "Passing # of nodes via : %s\n", procstring[0]);
2357 fprintf(fp, "Not setting number of nodes in system call\n");
2362 fprintf(fp, "Number of threads : %d\n", nnodes);
2365 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2366 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2367 fprintf(fp, "Benchmark steps : ");
2368 fprintf(fp, gmx_large_int_pfmt, bench_nsteps);
2370 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2373 fprintf(fp, "Checkpoint time step : ");
2374 fprintf(fp, gmx_large_int_pfmt, cpt_steps);
2377 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2379 if (new_sim_nsteps >= 0)
2382 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so",NFILE,fnm));
2383 fprintf(stderr, gmx_large_int_pfmt, new_sim_nsteps+cpt_steps);
2384 fprintf(stderr, " steps.\n");
2385 fprintf(fp, "Simulation steps : ");
2386 fprintf(fp, gmx_large_int_pfmt, new_sim_nsteps);
2390 fprintf(fp, "Repeats for each test : %d\n", repeats);
2392 if (npme_fixed >= -1)
2394 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2397 fprintf(fp, "Input file : %s\n", opt2fn("-s",NFILE,fnm));
2398 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2400 /* Allocate memory for the inputinfo struct: */
2402 info->nr_inputfiles = ntprs;
2403 for (i=0; i<ntprs; i++)
2405 snew(info->rcoulomb , ntprs);
2406 snew(info->rvdw , ntprs);
2407 snew(info->rlist , ntprs);
2408 snew(info->rlistlong, ntprs);
2409 snew(info->nkx , ntprs);
2410 snew(info->nky , ntprs);
2411 snew(info->nkz , ntprs);
2412 snew(info->fsx , ntprs);
2413 snew(info->fsy , ntprs);
2414 snew(info->fsz , ntprs);
2416 /* Make alternative tpr files to test: */
2417 snew(tpr_names, ntprs);
2418 for (i=0; i<ntprs; i++)
2420 snew(tpr_names[i], STRLEN);
2423 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2424 * different grids could be found. */
2425 make_benchmark_tprs(opt2fn("-s",NFILE,fnm), tpr_names, bench_nsteps+presteps,
2426 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2428 /********************************************************************************/
2429 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2430 /********************************************************************************/
2431 snew(perfdata, ntprs);
2434 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2435 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2436 cmd_args_bench, fnm, NFILE, presteps, cpt_steps);
2438 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gettime()-seconds)/60.0);
2440 /* Analyse the results and give a suggestion for optimal settings: */
2441 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2442 repeats, info, &best_tpr, &best_npme);
2444 /* Take the best-performing tpr file and enlarge nsteps to original value */
2445 if ( bKeepTPR && !bOverwrite )
2447 simulation_tpr = opt2fn("-s",NFILE,fnm);
2451 simulation_tpr = opt2fn("-so",NFILE,fnm);
2452 modify_PMEsettings(bOverwrite? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2453 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2456 /* Let's get rid of the temporary benchmark input files */
2457 for (i=0; i<ntprs; i++)
2459 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2460 remove(tpr_names[i]);
2463 /* Now start the real simulation if the user requested it ... */
2464 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2465 cmd_args_launch, simulation_tpr, best_npme);
2469 /* ... or simply print the performance results to screen: */
2472 finalize(opt2fn("-p", NFILE, fnm));