added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / src / tools / gmx_tune_pme.c
blob7f96afbf7aa52a6dd3f854e5366710849f418c93
1 /*
2 *
3 * This source code is part of
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
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
31 * And Hey:
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34 #ifdef HAVE_CONFIG_H
35 #include <config.h>
36 #endif
39 #include <time.h>
40 #ifdef HAVE_SYS_TIME_H
41 #include <sys/time.h>
42 #endif
46 #include "statutil.h"
47 #include "typedefs.h"
48 #include "smalloc.h"
49 #include "vec.h"
50 #include "copyrite.h"
51 #include "statutil.h"
52 #include "tpxio.h"
53 #include "string2.h"
54 #include "readinp.h"
55 #include "calcgrid.h"
56 #include "checkpoint.h"
57 #include "gmx_ana.h"
58 #include "names.h"
59 #include "perf_est.h"
63 /* Enum for situations that can occur during log file parsing, the
64 * corresponding string entries can be found in do_the_tests() in
65 * const char* ParseLog[] */
66 enum {
67 eParselogOK,
68 eParselogNotFound,
69 eParselogNoPerfData,
70 eParselogTerm,
71 eParselogResetProblem,
72 eParselogNoDDGrid,
73 eParselogTPXVersion,
74 eParselogNotParallel,
75 eParselogFatal,
76 eParselogNr
80 typedef struct
82 int nPMEnodes; /* number of PME-only nodes used in this test */
83 int nx, ny, nz; /* DD grid */
84 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
85 double *Gcycles; /* This can contain more than one value if doing multiple tests */
86 double Gcycles_Av;
87 float *ns_per_day;
88 float ns_per_day_Av;
89 float *PME_f_load; /* PME mesh/force load average*/
90 float PME_f_load_Av; /* Average average ;) ... */
91 char *mdrun_cmd_line; /* Mdrun command line used for this test */
92 } t_perf;
95 typedef struct
97 int nr_inputfiles; /* The number of tpr and mdp input files */
98 gmx_large_int_t orig_sim_steps; /* Number of steps to be done in the real simulation */
99 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
100 real *rvdw; /* The vdW radii */
101 real *rlist; /* Neighbourlist cutoff radius */
102 real *rlistlong;
103 int *nkx, *nky, *nkz;
104 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
105 } t_inputinfo;
108 static void sep_line(FILE *fp)
110 fprintf(fp, "\n------------------------------------------------------------\n");
114 /* Wrapper for system calls */
115 static int gmx_system_call(char *command)
117 #ifdef GMX_NO_SYSTEM
118 gmx_fatal(FARGS,"No calls to system(3) supported on this platform. Attempted to call:\n'%s'\n",command);
119 #else
120 return ( system(command) );
121 #endif
125 /* Check if string starts with substring */
126 static gmx_bool str_starts(const char *string, const char *substring)
128 return ( strncmp(string, substring, strlen(substring)) == 0);
132 static void cleandata(t_perf *perfdata, int test_nr)
134 perfdata->Gcycles[test_nr] = 0.0;
135 perfdata->ns_per_day[test_nr] = 0.0;
136 perfdata->PME_f_load[test_nr] = 0.0;
138 return;
142 static gmx_bool is_equal(real a, real b)
144 real diff, eps=1.0e-7;
147 diff = a - b;
149 if (diff < 0.0) diff = -diff;
151 if (diff < eps)
152 return TRUE;
153 else
154 return FALSE;
158 static void finalize(const char *fn_out)
160 char buf[STRLEN];
161 FILE *fp;
164 fp = fopen(fn_out,"r");
165 fprintf(stdout,"\n\n");
167 while( fgets(buf,STRLEN-1,fp) != NULL )
169 fprintf(stdout,"%s",buf);
171 fclose(fp);
172 fprintf(stdout,"\n\n");
176 enum {eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr};
178 static int parse_logfile(const char *logfile, const char *errfile,
179 t_perf *perfdata, int test_nr, int presteps, gmx_large_int_t cpt_steps,
180 int nnodes)
182 FILE *fp;
183 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
184 const char matchstrdd[]="Domain decomposition grid";
185 const char matchstrcr[]="resetting all time and cycle counters";
186 const char matchstrbal[]="Average PME mesh/force load:";
187 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";
188 const char errSIG[]="signal, stopping at the next";
189 int iFound;
190 int procs;
191 float dum1,dum2,dum3;
192 int npme;
193 gmx_large_int_t resetsteps=-1;
194 gmx_bool bFoundResetStr = FALSE;
195 gmx_bool bResetChecked = FALSE;
198 if (!gmx_fexist(logfile))
200 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
201 cleandata(perfdata, test_nr);
202 return eParselogNotFound;
205 fp = fopen(logfile, "r");
206 perfdata->PME_f_load[test_nr] = -1.0;
207 perfdata->guessPME = -1;
209 iFound = eFoundNothing;
210 if (1 == nnodes)
211 iFound = eFoundDDStr; /* Skip some case statements */
213 while (fgets(line, STRLEN, fp) != NULL)
215 /* Remove leading spaces */
216 ltrim(line);
218 /* Check for TERM and INT signals from user: */
219 if ( strstr(line, errSIG) != NULL )
221 fclose(fp);
222 cleandata(perfdata, test_nr);
223 return eParselogTerm;
226 /* Check whether cycle resetting worked */
227 if (presteps > 0 && !bFoundResetStr)
229 if (strstr(line, matchstrcr) != NULL)
231 sprintf(dumstring, "Step %s", gmx_large_int_pfmt);
232 sscanf(line, dumstring, &resetsteps);
233 bFoundResetStr = TRUE;
234 if (resetsteps == presteps+cpt_steps)
236 bResetChecked = TRUE;
238 else
240 sprintf(dumstring , gmx_large_int_pfmt, resetsteps);
241 sprintf(dumstring2, gmx_large_int_pfmt, presteps+cpt_steps);
242 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
243 " though they were supposed to be reset at step %s!\n",
244 dumstring, dumstring2);
249 /* Look for strings that appear in a certain order in the log file: */
250 switch(iFound)
252 case eFoundNothing:
253 /* Look for domain decomp grid and separate PME nodes: */
254 if (str_starts(line, matchstrdd))
256 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME nodes %d",
257 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
258 if (perfdata->nPMEnodes == -1)
259 perfdata->guessPME = npme;
260 else if (perfdata->nPMEnodes != npme)
261 gmx_fatal(FARGS, "PME nodes from command line and output file are not identical");
262 iFound = eFoundDDStr;
264 /* Catch a few errors that might have occured: */
265 else if (str_starts(line, "There is no domain decomposition for"))
267 fclose(fp);
268 return eParselogNoDDGrid;
270 else if (str_starts(line, "reading tpx file"))
272 fclose(fp);
273 return eParselogTPXVersion;
275 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
277 fclose(fp);
278 return eParselogNotParallel;
280 break;
281 case eFoundDDStr:
282 /* Look for PME mesh/force balance (not necessarily present, though) */
283 if (str_starts(line, matchstrbal))
284 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
285 /* Look for matchstring */
286 if (str_starts(line, matchstring))
287 iFound = eFoundAccountingStr;
288 break;
289 case eFoundAccountingStr:
290 /* Already found matchstring - look for cycle data */
291 if (str_starts(line, "Total "))
293 sscanf(line,"Total %d %lf",&procs,&(perfdata->Gcycles[test_nr]));
294 iFound = eFoundCycleStr;
296 break;
297 case eFoundCycleStr:
298 /* Already found cycle data - look for remaining performance info and return */
299 if (str_starts(line, "Performance:"))
301 sscanf(line,"%s %f %f %f %f", dumstring, &dum1, &dum2, &(perfdata->ns_per_day[test_nr]), &dum3);
302 fclose(fp);
303 if (bResetChecked || presteps == 0)
304 return eParselogOK;
305 else
306 return eParselogResetProblem;
308 break;
310 } /* while */
312 /* Close the log file */
313 fclose(fp);
315 /* Check why there is no performance data in the log file.
316 * Did a fatal errors occur? */
317 if (gmx_fexist(errfile))
319 fp = fopen(errfile, "r");
320 while (fgets(line, STRLEN, fp) != NULL)
322 if ( str_starts(line, "Fatal error:") )
324 if (fgets(line, STRLEN, fp) != NULL)
325 fprintf(stderr, "\nWARNING: An error occured during this benchmark:\n"
326 "%s\n", line);
327 fclose(fp);
328 cleandata(perfdata, test_nr);
329 return eParselogFatal;
332 fclose(fp);
334 else
336 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
339 /* Giving up ... we could not find out why there is no performance data in
340 * the log file. */
341 fprintf(stdout, "No performance data in log file.\n");
342 cleandata(perfdata, test_nr);
344 return eParselogNoPerfData;
348 static gmx_bool analyze_data(
349 FILE *fp,
350 const char *fn,
351 t_perf **perfdata,
352 int nnodes,
353 int ntprs,
354 int ntests,
355 int nrepeats,
356 t_inputinfo *info,
357 int *index_tpr, /* OUT: Nr of mdp file with best settings */
358 int *npme_optimal) /* OUT: Optimal number of PME nodes */
360 int i,j,k;
361 int line=0, line_win=-1;
362 int k_win=-1, i_win=-1, winPME;
363 double s=0.0; /* standard deviation */
364 t_perf *pd;
365 char strbuf[STRLEN];
366 char str_PME_f_load[13];
367 gmx_bool bCanUseOrigTPR;
368 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
371 if (nrepeats > 1)
373 sep_line(fp);
374 fprintf(fp, "Summary of successful runs:\n");
375 fprintf(fp, "Line tpr PME nodes Gcycles Av. Std.dev. ns/day PME/f");
376 if (nnodes > 1)
377 fprintf(fp, " DD grid");
378 fprintf(fp, "\n");
382 for (k=0; k<ntprs; k++)
384 for (i=0; i<ntests; i++)
386 /* Select the right dataset: */
387 pd = &(perfdata[k][i]);
389 pd->Gcycles_Av = 0.0;
390 pd->PME_f_load_Av = 0.0;
391 pd->ns_per_day_Av = 0.0;
393 if (pd->nPMEnodes == -1)
394 sprintf(strbuf, "(%3d)", pd->guessPME);
395 else
396 sprintf(strbuf, " ");
398 /* Get the average run time of a setting */
399 for (j=0; j<nrepeats; j++)
401 pd->Gcycles_Av += pd->Gcycles[j];
402 pd->PME_f_load_Av += pd->PME_f_load[j];
404 pd->Gcycles_Av /= nrepeats;
405 pd->PME_f_load_Av /= nrepeats;
407 for (j=0; j<nrepeats; j++)
409 if (pd->ns_per_day[j] > 0.0)
410 pd->ns_per_day_Av += pd->ns_per_day[j];
411 else
413 /* Somehow the performance number was not aquired for this run,
414 * therefor set the average to some negative value: */
415 pd->ns_per_day_Av = -1.0f*nrepeats;
416 break;
419 pd->ns_per_day_Av /= nrepeats;
421 /* Nicer output: */
422 if (pd->PME_f_load_Av > 0.0)
423 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
424 else
425 sprintf(str_PME_f_load, "%s", " - ");
428 /* We assume we had a successful run if both averages are positive */
429 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
431 /* Output statistics if repeats were done */
432 if (nrepeats > 1)
434 /* Calculate the standard deviation */
435 s = 0.0;
436 for (j=0; j<nrepeats; j++)
437 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
438 s /= (nrepeats - 1);
439 s = sqrt(s);
441 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
442 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
443 pd->ns_per_day_Av, str_PME_f_load);
444 if (nnodes > 1)
445 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
446 fprintf(fp, "\n");
448 /* Store the index of the best run found so far in 'winner': */
449 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
451 k_win = k;
452 i_win = i;
453 line_win = line;
455 line++;
460 if (k_win == -1)
461 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
463 sep_line(fp);
465 winPME = perfdata[k_win][i_win].nPMEnodes;
467 if (1 == ntests)
469 /* We stuck to a fixed number of PME-only nodes */
470 sprintf(strbuf, "settings No. %d", k_win);
472 else
474 /* We have optimized the number of PME-only nodes */
475 if (winPME == -1)
476 sprintf(strbuf, "%s", "the automatic number of PME nodes");
477 else
478 sprintf(strbuf, "%d PME nodes", winPME);
480 fprintf(fp, "Best performance was achieved with %s", strbuf);
481 if ((nrepeats > 1) && (ntests > 1))
482 fprintf(fp, " (see line %d)", line_win);
483 fprintf(fp, "\n");
485 /* Only mention settings if they were modified: */
486 bRefinedCoul = !is_equal(info->rcoulomb[k_win], info->rcoulomb[0]);
487 bRefinedVdW = !is_equal(info->rvdw[k_win] , info->rvdw[0] );
488 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
489 info->nky[k_win] == info->nky[0] &&
490 info->nkz[k_win] == info->nkz[0]);
492 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
494 fprintf(fp, "Optimized PME settings:\n");
495 bCanUseOrigTPR = FALSE;
497 else
499 bCanUseOrigTPR = TRUE;
502 if (bRefinedCoul)
503 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
505 if (bRefinedVdW)
506 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
508 if (bRefinedGrid)
509 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],
510 info->nkx[0], info->nky[0], info->nkz[0]);
512 if (bCanUseOrigTPR && ntprs > 1)
513 fprintf(fp, "and original PME settings.\n");
515 fflush(fp);
517 /* Return the index of the mdp file that showed the highest performance
518 * and the optimal number of PME nodes */
519 *index_tpr = k_win;
520 *npme_optimal = winPME;
522 return bCanUseOrigTPR;
526 /* Get the commands we need to set up the runs from environment variables */
527 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char cmd_np[],
528 char *cmd_mdrun[], int repeats)
530 char *command=NULL;
531 char *cp;
532 char *cp2;
533 char line[STRLEN];
534 FILE *fp;
535 const char def_mpirun[] = "mpirun";
536 const char def_mdrun[] = "mdrun";
537 const char filename[] = "benchtest.log";
538 const char match_mpi[] = "NNODES=";
539 const char match_mdrun[]= "Program: ";
540 const char empty_mpirun[] = "";
541 gmx_bool bMdrun = FALSE;
542 gmx_bool bMPI = FALSE;
545 /* Get the commands we need to set up the runs from environment variables */
546 if (!bThreads)
548 if ( (cp = getenv("MPIRUN")) != NULL)
549 *cmd_mpirun = strdup(cp);
550 else
551 *cmd_mpirun = strdup(def_mpirun);
553 else
555 *cmd_mpirun = strdup(empty_mpirun);
558 if ( (cp = getenv("MDRUN" )) != NULL )
559 *cmd_mdrun = strdup(cp);
560 else
561 *cmd_mdrun = strdup(def_mdrun);
564 /* If no simulations have to be performed, we are done here */
565 if (repeats <= 0)
566 return;
568 /* Run a small test to see whether mpirun + mdrun work */
569 fprintf(stdout, "Making sure that mdrun can be executed. ");
570 if (bThreads)
572 snew(command, strlen(*cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
573 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", *cmd_mdrun, cmd_np, filename);
575 else
577 snew(command, strlen(*cmd_mpirun) + strlen(cmd_np) + strlen(*cmd_mdrun) + strlen(filename) + 50);
578 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", *cmd_mpirun, cmd_np, *cmd_mdrun, filename);
580 fprintf(stdout, "Trying '%s' ... ", command);
581 make_backup(filename);
582 gmx_system_call(command);
584 /* Check if we find the characteristic string in the output: */
585 if (!gmx_fexist(filename))
586 gmx_fatal(FARGS, "Output from test run could not be found.");
588 fp = fopen(filename, "r");
589 /* We need to scan the whole output file, since sometimes the queuing system
590 * also writes stuff to stdout/err */
591 while ( !feof(fp) )
593 cp2=fgets(line, STRLEN, fp);
594 if (cp2!=NULL)
596 if ( str_starts(line, match_mdrun) )
597 bMdrun = TRUE;
598 if ( str_starts(line, match_mpi) )
599 bMPI = TRUE;
602 fclose(fp);
604 if (bThreads)
606 if (bMPI)
608 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
609 "(%s)\n"
610 "seems to have been compiled with MPI instead.",
611 *cmd_mdrun);
614 else
616 if (bMdrun && !bMPI)
618 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
619 "(%s)\n"
620 "seems to have been compiled without MPI support.",
621 *cmd_mdrun);
625 if (!bMdrun)
627 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
628 filename);
631 fprintf(stdout, "passed.\n");
633 /* Clean up ... */
634 remove(filename);
635 sfree(command);
639 static void launch_simulation(
640 gmx_bool bLaunch, /* Should the simulation be launched? */
641 FILE *fp, /* General log file */
642 gmx_bool bThreads, /* whether to use threads */
643 char *cmd_mpirun, /* Command for mpirun */
644 char *cmd_np, /* Switch for -np or -nt or empty */
645 char *cmd_mdrun, /* Command for mdrun */
646 char *args_for_mdrun, /* Arguments for mdrun */
647 const char *simulation_tpr, /* This tpr will be simulated */
648 int nnodes, /* Number of nodes to run on */
649 int nPMEnodes) /* Number of PME nodes to use */
651 char *command;
654 /* Make enough space for the system call command,
655 * (100 extra chars for -npme ... etc. options should suffice): */
656 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+100);
658 /* Note that the -passall options requires args_for_mdrun to be at the end
659 * of the command line string */
660 if (bThreads)
662 sprintf(command, "%s%s-npme %d -s %s %s",
663 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun);
665 else
667 sprintf(command, "%s%s%s -npme %d -s %s %s",
668 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun);
671 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch? "Using":"Please use", command);
672 sep_line(fp);
673 fflush(fp);
675 /* Now the real thing! */
676 if (bLaunch)
678 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
679 sep_line(stdout);
680 fflush(stdout);
681 gmx_system_call(command);
686 static void modify_PMEsettings(
687 gmx_large_int_t simsteps, /* Set this value as number of time steps */
688 const char *fn_best_tpr, /* tpr file with the best performance */
689 const char *fn_sim_tpr) /* name of tpr file to be launched */
691 t_inputrec *ir;
692 t_state state;
693 gmx_mtop_t mtop;
694 char buf[200];
696 snew(ir,1);
697 read_tpx_state(fn_best_tpr,ir,&state,NULL,&mtop);
699 /* Set nsteps to the right value */
700 ir->nsteps = simsteps;
702 /* Write the tpr file which will be launched */
703 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, gmx_large_int_pfmt);
704 fprintf(stdout,buf,ir->nsteps);
705 fflush(stdout);
706 write_tpx_state(fn_sim_tpr,ir,&state,&mtop);
708 sfree(ir);
712 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
714 /* Make additional TPR files with more computational load for the
715 * direct space processors: */
716 static void make_benchmark_tprs(
717 const char *fn_sim_tpr, /* READ : User-provided tpr file */
718 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
719 gmx_large_int_t benchsteps, /* Number of time steps for benchmark runs */
720 gmx_large_int_t statesteps, /* Step counter in checkpoint file */
721 real rmin, /* Minimal Coulomb radius */
722 real rmax, /* Maximal Coulomb radius */
723 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
724 int *ntprs, /* No. of TPRs to write, each with a different
725 rcoulomb and fourierspacing */
726 t_inputinfo *info, /* Contains information about mdp file options */
727 FILE *fp) /* Write the output here */
729 int i,j,d;
730 t_inputrec *ir;
731 t_state state;
732 gmx_mtop_t mtop;
733 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
734 char buf[200];
735 rvec box_size;
736 gmx_bool bNote = FALSE;
737 real add; /* Add this to rcoul for the next test */
738 real fac = 1.0; /* Scaling factor for Coulomb radius */
739 real fourierspacing; /* Basic fourierspacing from tpr */
742 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
743 *ntprs > 1? "s":"", gmx_large_int_pfmt, benchsteps>1?"s":"");
744 fprintf(stdout, buf, benchsteps);
745 if (statesteps > 0)
747 sprintf(buf, " (adding %s steps from checkpoint file)", gmx_large_int_pfmt);
748 fprintf(stdout, buf, statesteps);
749 benchsteps += statesteps;
751 fprintf(stdout, ".\n");
754 snew(ir,1);
755 read_tpx_state(fn_sim_tpr,ir,&state,NULL,&mtop);
757 /* Check if some kind of PME was chosen */
758 if (EEL_PME(ir->coulombtype) == FALSE)
759 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
760 EELTYPE(eelPME));
762 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
763 if ( (ir->cutoff_scheme != ecutsVERLET) &&
764 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
766 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
767 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
769 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
770 else if (ir->rcoulomb > ir->rlist)
772 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
773 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
776 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
778 fprintf(stdout,"NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
779 bScaleRvdw = FALSE;
782 /* Reduce the number of steps for the benchmarks */
783 info->orig_sim_steps = ir->nsteps;
784 ir->nsteps = benchsteps;
786 /* For PME-switch potentials, keep the radial distance of the buffer region */
787 nlist_buffer = ir->rlist - ir->rcoulomb;
789 /* Determine length of triclinic box vectors */
790 for(d=0; d<DIM; d++)
792 box_size[d] = 0;
793 for(i=0;i<DIM;i++)
794 box_size[d] += state.box[d][i]*state.box[d][i];
795 box_size[d] = sqrt(box_size[d]);
798 if (ir->fourier_spacing > 0)
800 info->fsx[0] = ir->fourier_spacing;
801 info->fsy[0] = ir->fourier_spacing;
802 info->fsz[0] = ir->fourier_spacing;
804 else
806 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
807 info->fsx[0] = box_size[XX]/ir->nkx;
808 info->fsy[0] = box_size[YY]/ir->nky;
809 info->fsz[0] = box_size[ZZ]/ir->nkz;
812 /* If no value for the fourierspacing was provided on the command line, we
813 * use the reconstruction from the tpr file */
814 if (ir->fourier_spacing > 0)
816 /* Use the spacing from the tpr */
817 fourierspacing = ir->fourier_spacing;
819 else
821 /* Use the maximum observed spacing */
822 fourierspacing = max(max(info->fsx[0],info->fsy[0]),info->fsz[0]);
825 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
827 /* For performance comparisons the number of particles is useful to have */
828 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
830 /* Print information about settings of which some are potentially modified: */
831 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
832 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
833 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
834 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
835 if (EVDW_SWITCHED(ir->vdwtype))
836 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
837 if (EPME_SWITCHED(ir->coulombtype))
838 fprintf(fp, " rlist : %f nm\n", ir->rlist);
839 if (ir->rlistlong != max_cutoff(ir->rvdw,ir->rcoulomb))
840 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
842 /* Print a descriptive line about the tpr settings tested */
843 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
844 fprintf(fp, " No. scaling rcoulomb");
845 fprintf(fp, " nkx nky nkz");
846 fprintf(fp, " spacing");
847 if (evdwCUT == ir->vdwtype)
848 fprintf(fp, " rvdw");
849 if (EPME_SWITCHED(ir->coulombtype))
850 fprintf(fp, " rlist");
851 if ( ir->rlistlong != max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb)) )
852 fprintf(fp, " rlistlong");
853 fprintf(fp, " tpr file\n");
855 /* Loop to create the requested number of tpr input files */
856 for (j = 0; j < *ntprs; j++)
858 /* The first .tpr is the provided one, just need to modify nsteps,
859 * so skip the following block */
860 if (j != 0)
862 /* Determine which Coulomb radii rc to use in the benchmarks */
863 add = (rmax-rmin)/(*ntprs-1);
864 if (is_equal(rmin,info->rcoulomb[0]))
866 ir->rcoulomb = rmin + j*add;
868 else if (is_equal(rmax,info->rcoulomb[0]))
870 ir->rcoulomb = rmin + (j-1)*add;
872 else
874 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
875 add = (rmax-rmin)/(*ntprs-2);
876 ir->rcoulomb = rmin + (j-1)*add;
879 /* Determine the scaling factor fac */
880 fac = ir->rcoulomb/info->rcoulomb[0];
882 /* Scale the Fourier grid spacing */
883 ir->nkx = ir->nky = ir->nkz = 0;
884 calc_grid(NULL,state.box,fourierspacing*fac,&ir->nkx,&ir->nky,&ir->nkz);
886 /* Adjust other radii since various conditions neet to be fulfilled */
887 if (eelPME == ir->coulombtype)
889 /* plain PME, rcoulomb must be equal to rlist */
890 ir->rlist = ir->rcoulomb;
892 else
894 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
895 ir->rlist = ir->rcoulomb + nlist_buffer;
898 if (bScaleRvdw && evdwCUT == ir->vdwtype)
900 /* For vdw cutoff, rvdw >= rlist */
901 ir->rvdw = max(info->rvdw[0], ir->rlist);
904 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
906 } /* end of "if (j != 0)" */
908 /* for j==0: Save the original settings
909 * for j >0: Save modified radii and Fourier grids */
910 info->rcoulomb[j] = ir->rcoulomb;
911 info->rvdw[j] = ir->rvdw;
912 info->nkx[j] = ir->nkx;
913 info->nky[j] = ir->nky;
914 info->nkz[j] = ir->nkz;
915 info->rlist[j] = ir->rlist;
916 info->rlistlong[j] = ir->rlistlong;
917 info->fsx[j] = fac*fourierspacing;
918 info->fsy[j] = fac*fourierspacing;
919 info->fsz[j] = fac*fourierspacing;
921 /* Write the benchmark tpr file */
922 strncpy(fn_bench_tprs[j],fn_sim_tpr,strlen(fn_sim_tpr)-strlen(".tpr"));
923 sprintf(buf, "_bench%.2d.tpr", j);
924 strcat(fn_bench_tprs[j], buf);
925 fprintf(stdout,"Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
926 fprintf(stdout, gmx_large_int_pfmt, ir->nsteps);
927 if (j > 0)
928 fprintf(stdout,", scaling factor %f\n", fac);
929 else
930 fprintf(stdout,", unmodified settings\n");
932 write_tpx_state(fn_bench_tprs[j],ir,&state,&mtop);
934 /* Write information about modified tpr settings to log file */
935 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
936 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
937 fprintf(fp, " %9f ", info->fsx[j]);
938 if (evdwCUT == ir->vdwtype)
939 fprintf(fp, "%10f", ir->rvdw);
940 if (EPME_SWITCHED(ir->coulombtype))
941 fprintf(fp, "%10f", ir->rlist);
942 if ( info->rlistlong[0] != max_cutoff(info->rlist[0],max_cutoff(info->rvdw[0],info->rcoulomb[0])) )
943 fprintf(fp, "%10f", ir->rlistlong);
944 fprintf(fp, " %-14s\n",fn_bench_tprs[j]);
946 /* Make it clear to the user that some additional settings were modified */
947 if ( !is_equal(ir->rvdw , info->rvdw[0])
948 || !is_equal(ir->rlistlong, info->rlistlong[0]) )
950 bNote = TRUE;
953 if (bNote)
954 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
955 "other input settings were also changed (see table above).\n"
956 "Please check if the modified settings are appropriate.\n");
957 fflush(stdout);
958 fflush(fp);
959 sfree(ir);
963 /* Rename the files we want to keep to some meaningful filename and
964 * delete the rest */
965 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
966 int nPMEnodes, int nr, gmx_bool bKeepStderr)
968 char numstring[STRLEN];
969 char newfilename[STRLEN];
970 const char *fn=NULL;
971 int i;
972 const char *opt;
975 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
977 for (i=0; i<nfile; i++)
979 opt = (char *)fnm[i].opt;
980 if ( strcmp(opt, "-p") == 0 )
982 /* do nothing; keep this file */
985 else if (strcmp(opt, "-bg") == 0)
987 /* Give the log file a nice name so one can later see which parameters were used */
988 numstring[0] = '\0';
989 if (nr > 0)
990 sprintf(numstring, "_%d", nr);
991 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg",nfile,fnm), k, nnodes, nPMEnodes, numstring);
992 if (gmx_fexist(opt2fn("-bg",nfile,fnm)))
994 fprintf(stdout, "renaming log file to %s\n", newfilename);
995 make_backup(newfilename);
996 rename(opt2fn("-bg",nfile,fnm), newfilename);
999 else if (strcmp(opt, "-err") == 0)
1001 /* This file contains the output of stderr. We want to keep it in
1002 * cases where there have been problems. */
1003 fn = opt2fn(opt, nfile, fnm);
1004 numstring[0] = '\0';
1005 if (nr > 0)
1006 sprintf(numstring, "_%d", nr);
1007 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1008 if (gmx_fexist(fn))
1010 if (bKeepStderr)
1012 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1013 make_backup(newfilename);
1014 rename(fn, newfilename);
1016 else
1018 fprintf(stdout, "Deleting %s\n", fn);
1019 remove(fn);
1023 /* Delete the files which are created for each benchmark run: (options -b*) */
1024 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt,nfile,fnm) || !is_optional(&fnm[i])) )
1026 fn = opt2fn(opt, nfile, fnm);
1027 if (gmx_fexist(fn))
1029 fprintf(stdout, "Deleting %s\n", fn);
1030 remove(fn);
1037 /* Returns the largest common factor of n1 and n2 */
1038 static int largest_common_factor(int n1, int n2)
1040 int factor, nmax;
1042 nmax = min(n1, n2);
1043 for (factor=nmax; factor > 0; factor--)
1045 if ( 0==(n1 % factor) && 0==(n2 % factor) )
1047 return(factor);
1050 return 0; /* one for the compiler */
1053 enum {eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr};
1055 /* Create a list of numbers of PME nodes to test */
1056 static void make_npme_list(
1057 const char *npmevalues_opt, /* Make a complete list with all
1058 * possibilities or a short list that keeps only
1059 * reasonable numbers of PME nodes */
1060 int *nentries, /* Number of entries we put in the nPMEnodes list */
1061 int *nPMEnodes[], /* Each entry contains the value for -npme */
1062 int nnodes, /* Total number of nodes to do the tests on */
1063 int minPMEnodes, /* Minimum number of PME nodes */
1064 int maxPMEnodes) /* Maximum number of PME nodes */
1066 int i,npme,npp;
1067 int min_factor=1; /* We request that npp and npme have this minimal
1068 * largest common factor (depends on npp) */
1069 int nlistmax; /* Max. list size */
1070 int nlist; /* Actual number of entries in list */
1071 int eNPME=0;
1074 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1075 if ( 0 == strcmp(npmevalues_opt, "all") )
1077 eNPME = eNpmeAll;
1079 else if ( 0 == strcmp(npmevalues_opt, "subset") )
1081 eNPME = eNpmeSubset;
1083 else /* "auto" or "range" */
1085 if (nnodes <= 64)
1086 eNPME = eNpmeAll;
1087 else if (nnodes < 128)
1088 eNPME = eNpmeReduced;
1089 else
1090 eNPME = eNpmeSubset;
1093 /* Calculate how many entries we could possibly have (in case of -npme all) */
1094 if (nnodes > 2)
1096 nlistmax = maxPMEnodes - minPMEnodes + 3;
1097 if (0 == minPMEnodes)
1098 nlistmax--;
1100 else
1101 nlistmax = 1;
1103 /* Now make the actual list which is at most of size nlist */
1104 snew(*nPMEnodes, nlistmax);
1105 nlist = 0; /* start counting again, now the real entries in the list */
1106 for (i = 0; i < nlistmax - 2; i++)
1108 npme = maxPMEnodes - i;
1109 npp = nnodes-npme;
1110 switch (eNPME)
1112 case eNpmeAll:
1113 min_factor = 1;
1114 break;
1115 case eNpmeReduced:
1116 min_factor = 2;
1117 break;
1118 case eNpmeSubset:
1119 /* For 2d PME we want a common largest factor of at least the cube
1120 * root of the number of PP nodes */
1121 min_factor = (int) pow(npp, 1.0/3.0);
1122 break;
1123 default:
1124 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1125 break;
1127 if (largest_common_factor(npp, npme) >= min_factor)
1129 (*nPMEnodes)[nlist] = npme;
1130 nlist++;
1133 /* We always test 0 PME nodes and the automatic number */
1134 *nentries = nlist + 2;
1135 (*nPMEnodes)[nlist ] = 0;
1136 (*nPMEnodes)[nlist+1] = -1;
1138 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1139 for (i=0; i<*nentries-1; i++)
1140 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1141 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1145 /* Allocate memory to store the performance data */
1146 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1148 int i, j, k;
1151 for (k=0; k<ntprs; k++)
1153 snew(perfdata[k], datasets);
1154 for (i=0; i<datasets; i++)
1156 for (j=0; j<repeats; j++)
1158 snew(perfdata[k][i].Gcycles , repeats);
1159 snew(perfdata[k][i].ns_per_day, repeats);
1160 snew(perfdata[k][i].PME_f_load, repeats);
1167 /* Check for errors on mdrun -h */
1168 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp)
1170 char *command, *msg;
1171 int ret;
1174 snew(command, length + 15);
1175 snew(msg , length + 500);
1177 fprintf(stdout, "Making shure the benchmarks can be executed ...\n");
1178 sprintf(command, "%s-h -quiet", mdrun_cmd_line);
1179 ret = gmx_system_call(command);
1181 if (0 != ret)
1183 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1184 * get the error message from mdrun itself */
1185 sprintf(msg, "Cannot run the benchmark simulations! Please check the error message of\n"
1186 "mdrun for the source of the problem. Did you provide a command line\n"
1187 "argument that neither g_tune_pme nor mdrun understands? Offending command:\n"
1188 "\n%s\n\n", command);
1190 fprintf(stderr, "%s", msg);
1191 sep_line(fp);
1192 fprintf(fp , "%s", msg);
1194 exit(ret);
1197 sfree(command);
1198 sfree(msg );
1202 static void do_the_tests(
1203 FILE *fp, /* General g_tune_pme output file */
1204 char **tpr_names, /* Filenames of the input files to test */
1205 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1206 int minPMEnodes, /* Min fraction of nodes to use for PME */
1207 int npme_fixed, /* If >= -1, test fixed number of PME
1208 * nodes only */
1209 const char *npmevalues_opt, /* Which -npme values should be tested */
1210 t_perf **perfdata, /* Here the performace data is stored */
1211 int *pmeentries, /* Entries in the nPMEnodes list */
1212 int repeats, /* Repeat each test this often */
1213 int nnodes, /* Total number of nodes = nPP + nPME */
1214 int nr_tprs, /* Total number of tpr files to test */
1215 gmx_bool bThreads, /* Threads or MPI? */
1216 char *cmd_mpirun, /* mpirun command string */
1217 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1218 char *cmd_mdrun, /* mdrun command string */
1219 char *cmd_args_bench, /* arguments for mdrun in a string */
1220 const t_filenm *fnm, /* List of filenames from command line */
1221 int nfile, /* Number of files specified on the cmdl. */
1222 int sim_part, /* For checkpointing */
1223 int presteps, /* DLB equilibration steps, is checked */
1224 gmx_large_int_t cpt_steps) /* Time step counter in the checkpoint */
1226 int i,nr,k,ret,count=0,totaltests;
1227 int *nPMEnodes=NULL;
1228 t_perf *pd=NULL;
1229 int cmdline_length;
1230 char *command, *cmd_stub;
1231 char buf[STRLEN];
1232 gmx_bool bResetProblem=FALSE;
1233 gmx_bool bFirst=TRUE;
1236 /* This string array corresponds to the eParselog enum type at the start
1237 * of this file */
1238 const char* ParseLog[] = {"OK.",
1239 "Logfile not found!",
1240 "No timings, logfile truncated?",
1241 "Run was terminated.",
1242 "Counters were not reset properly.",
1243 "No DD grid found for these settings.",
1244 "TPX version conflict!",
1245 "mdrun was not started in parallel!",
1246 "An error occured." };
1247 char str_PME_f_load[13];
1250 /* Allocate space for the mdrun command line. 100 extra characters should
1251 be more than enough for the -npme etcetera arguments */
1252 cmdline_length = strlen(cmd_mpirun)
1253 + strlen(cmd_np)
1254 + strlen(cmd_mdrun)
1255 + strlen(cmd_args_bench)
1256 + strlen(tpr_names[0]) + 100;
1257 snew(command , cmdline_length);
1258 snew(cmd_stub, cmdline_length);
1260 /* Construct the part of the command line that stays the same for all tests: */
1261 if (bThreads)
1263 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1265 else
1267 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1270 /* Create a list of numbers of PME nodes to test */
1271 if (npme_fixed < -1)
1273 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1274 nnodes, minPMEnodes, maxPMEnodes);
1276 else
1278 *pmeentries = 1;
1279 snew(nPMEnodes, 1);
1280 nPMEnodes[0] = npme_fixed;
1281 fprintf(stderr, "Will use a fixed number of %d PME-only nodes.\n", nPMEnodes[0]);
1284 if (0 == repeats)
1286 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1287 ffclose(fp);
1288 finalize(opt2fn("-p", nfile, fnm));
1289 exit(0);
1292 /* Allocate one dataset for each tpr input file: */
1293 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1295 /*****************************************/
1296 /* Main loop over all tpr files to test: */
1297 /*****************************************/
1298 totaltests = nr_tprs*(*pmeentries)*repeats;
1299 for (k=0; k<nr_tprs;k++)
1301 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1302 fprintf(fp, "PME nodes Gcycles ns/day PME/f Remark\n");
1303 /* Loop over various numbers of PME nodes: */
1304 for (i = 0; i < *pmeentries; i++)
1306 pd = &perfdata[k][i];
1308 /* Loop over the repeats for each scenario: */
1309 for (nr = 0; nr < repeats; nr++)
1311 pd->nPMEnodes = nPMEnodes[i];
1313 /* Add -npme and -s to the command line and save it. Note that
1314 * the -passall (if set) options requires cmd_args_bench to be
1315 * at the end of the command line string */
1316 snew(pd->mdrun_cmd_line, cmdline_length);
1317 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s",
1318 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench);
1320 /* To prevent that all benchmarks fail due to a show-stopper argument
1321 * on the mdrun command line, we make a quick check with mdrun -h first */
1322 if (bFirst)
1323 make_sure_it_runs(pd->mdrun_cmd_line, cmdline_length, fp);
1324 bFirst = FALSE;
1326 /* Do a benchmark simulation: */
1327 if (repeats > 1)
1328 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1329 else
1330 buf[0]='\0';
1331 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1332 (100.0*count)/totaltests,
1333 k+1, nr_tprs, i+1, *pmeentries, buf);
1334 make_backup(opt2fn("-err",nfile,fnm));
1335 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err",nfile,fnm));
1336 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1337 gmx_system_call(command);
1339 /* Collect the performance data from the log file; also check stderr
1340 * for fatal errors */
1341 ret = parse_logfile(opt2fn("-bg",nfile,fnm), opt2fn("-err",nfile,fnm),
1342 pd, nr, presteps, cpt_steps, nnodes);
1343 if ((presteps > 0) && (ret == eParselogResetProblem))
1344 bResetProblem = TRUE;
1346 if (-1 == pd->nPMEnodes)
1347 sprintf(buf, "(%3d)", pd->guessPME);
1348 else
1349 sprintf(buf, " ");
1351 /* Nicer output */
1352 if (pd->PME_f_load[nr] > 0.0)
1353 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1354 else
1355 sprintf(str_PME_f_load, "%s", " - ");
1357 /* Write the data we got to disk */
1358 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1359 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1360 if (! (ret==eParselogOK || ret==eParselogNoDDGrid || ret==eParselogNotFound) )
1361 fprintf(fp, " Check %s file for problems.", ret==eParselogFatal? "err":"log");
1362 fprintf(fp, "\n");
1363 fflush(fp);
1364 count++;
1366 /* Do some cleaning up and delete the files we do not need any more */
1367 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret==eParselogFatal);
1369 /* If the first run with this number of processors already failed, do not try again: */
1370 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1372 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1373 count += repeats-(nr+1);
1374 break;
1376 } /* end of repeats loop */
1377 } /* end of -npme loop */
1378 } /* end of tpr file loop */
1380 if (bResetProblem)
1382 sep_line(fp);
1383 fprintf(fp, "WARNING: The cycle and time step counters could not be reset\n"
1384 "properly. The reason could be that mpirun did not manage to\n"
1385 "export the environment variable GMX_RESET_COUNTER. You might\n"
1386 "have to give a special switch to mpirun for that.\n"
1387 "Alternatively, you can manually set GMX_RESET_COUNTER to the\n"
1388 "value normally provided by -presteps.");
1389 sep_line(fp);
1391 sfree(command);
1392 sfree(cmd_stub);
1396 static void check_input(
1397 int nnodes,
1398 int repeats,
1399 int *ntprs,
1400 real *rmin,
1401 real rcoulomb,
1402 real *rmax,
1403 real maxPMEfraction,
1404 real minPMEfraction,
1405 int npme_fixed,
1406 gmx_large_int_t bench_nsteps,
1407 const t_filenm *fnm,
1408 int nfile,
1409 int sim_part,
1410 int presteps,
1411 int npargs,
1412 t_pargs *pa)
1414 int old;
1417 /* Make sure the input file exists */
1418 if (!gmx_fexist(opt2fn("-s",nfile,fnm)))
1419 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s",nfile,fnm));
1421 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1422 if ( (0 == strcmp(opt2fn("-cpi",nfile,fnm), opt2fn("-bcpo",nfile,fnm)) ) && (sim_part > 1) )
1423 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1424 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1426 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1427 if (repeats < 0)
1428 gmx_fatal(FARGS, "Number of repeats < 0!");
1430 /* Check number of nodes */
1431 if (nnodes < 1)
1432 gmx_fatal(FARGS, "Number of nodes/threads must be a positive integer.");
1434 /* Automatically choose -ntpr if not set */
1435 if (*ntprs < 1)
1437 if (nnodes < 16)
1438 *ntprs = 1;
1439 else
1441 *ntprs = 3;
1442 /* Set a reasonable scaling factor for rcoulomb */
1443 if (*rmax <= 0)
1444 *rmax = rcoulomb * 1.2;
1446 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs==1?"":"s");
1448 else
1450 if (1 == *ntprs)
1451 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1454 /* Make shure that rmin <= rcoulomb <= rmax */
1455 if (*rmin <= 0) *rmin = rcoulomb;
1456 if (*rmax <= 0) *rmax = rcoulomb;
1457 if ( !(*rmin <= *rmax) )
1459 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1460 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1462 /* Add test scenarios if rmin or rmax were set */
1463 if (*ntprs <= 2)
1465 if ( !is_equal(*rmin, rcoulomb) && (*ntprs == 1) )
1467 (*ntprs)++;
1468 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1469 *rmin, *ntprs);
1471 if (!is_equal(*rmax, rcoulomb) && (*ntprs == 1) )
1473 (*ntprs)++;
1474 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1475 *rmax, *ntprs);
1478 old = *ntprs;
1479 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1480 if ( !is_equal(*rmax, rcoulomb) || !is_equal(*rmin, rcoulomb) )
1481 *ntprs = max(*ntprs, 2);
1483 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1484 if ( !is_equal(*rmax, rcoulomb) && !is_equal(*rmin, rcoulomb) )
1485 *ntprs = max(*ntprs, 3);
1487 if (old != *ntprs)
1489 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1492 if (*ntprs > 1)
1494 if (is_equal(*rmin,rcoulomb) && is_equal(rcoulomb,*rmax)) /* We have just a single rc */
1496 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1497 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1498 "with correspondingly adjusted PME grid settings\n");
1499 *ntprs = 1;
1503 /* Check whether max and min fraction are within required values */
1504 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1505 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1506 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1507 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1508 if (maxPMEfraction < minPMEfraction)
1509 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1511 /* Check whether the number of steps - if it was set - has a reasonable value */
1512 if (bench_nsteps < 0)
1513 gmx_fatal(FARGS, "Number of steps must be positive.");
1515 if (bench_nsteps > 10000 || bench_nsteps < 100)
1517 fprintf(stderr, "WARNING: steps=");
1518 fprintf(stderr, gmx_large_int_pfmt, bench_nsteps);
1519 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100)? "few" : "many");
1522 if (presteps < 0)
1524 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1527 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1528 if (*ntprs > 1)
1530 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1531 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1534 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1535 * only. We need to check whether the requested number of PME-only nodes
1536 * makes sense. */
1537 if (npme_fixed > -1)
1539 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1540 if (2*npme_fixed > nnodes)
1542 gmx_fatal(FARGS, "Cannot have more than %d PME-only nodes for a total of %d nodes (you chose %d).\n",
1543 nnodes/2, nnodes, npme_fixed);
1545 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1547 fprintf(stderr, "WARNING: Only %g percent of the nodes are assigned as PME-only nodes.\n",
1548 100.0*((real)npme_fixed / (real)nnodes));
1550 if (opt2parg_bSet("-min",npargs,pa) || opt2parg_bSet("-max",npargs,pa))
1552 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1553 " fixed number of PME-only nodes is requested with -fix.\n");
1559 /* Returns TRUE when "opt" is needed at launch time */
1560 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1562 /* Apart from the input .tpr we need all options that were set
1563 * on the command line and that do not start with -b */
1564 if (0 == strncmp(opt,"-b", 2) || 0 == strncmp(opt,"-s", 2))
1565 return FALSE;
1567 if (bSet)
1568 return TRUE;
1569 else
1570 return FALSE;
1574 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1575 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1577 /* Apart from the input .tpr, all files starting with "-b" are for
1578 * _b_enchmark files exclusively */
1579 if (0 == strncmp(opt,"-s", 2)) return FALSE;
1580 if (0 == strncmp(opt,"-b", 2) || 0 == strncmp(opt,"-s", 2))
1582 if (!bOptional || bSet)
1583 return TRUE;
1584 else
1585 return FALSE;
1587 else
1589 if (bIsOutput)
1590 return FALSE;
1591 else
1592 if (bSet) /* These are additional input files like -cpi -ei */
1593 return TRUE;
1594 else
1595 return FALSE;
1600 /* Adds 'buf' to 'str' */
1601 static void add_to_string(char **str, char *buf)
1603 int len;
1606 len = strlen(*str) + strlen(buf) + 1;
1607 srenew(*str, len);
1608 strcat(*str, buf);
1612 /* Create the command line for the benchmark as well as for the real run */
1613 static void create_command_line_snippets(
1614 gmx_bool bThreads,
1615 gmx_bool bAppendFiles,
1616 gmx_bool bKeepAndNumCPT,
1617 gmx_bool bResetHWay,
1618 int presteps,
1619 int nfile,
1620 t_filenm fnm[],
1621 int npargs,
1622 t_pargs *pa,
1623 const char *procstring, /* How to pass the number of processors to $MPIRUN */
1624 char *cmd_np[], /* Actual command line snippet, e.g. '-np <N>' */
1625 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1626 char *cmd_args_launch[], /* command line arguments for simulation run */
1627 char extra_args[]) /* Add this to the end of the command line */
1629 int i;
1630 char *opt;
1631 const char *name;
1632 char strbuf[STRLEN];
1635 /* strlen needs at least '\0' as a string: */
1636 snew(*cmd_args_bench ,1);
1637 snew(*cmd_args_launch,1);
1638 *cmd_args_launch[0]='\0';
1639 *cmd_args_bench[0] ='\0';
1642 /*******************************************/
1643 /* 1. Process other command line arguments */
1644 /*******************************************/
1645 if (presteps > 0)
1647 /* Add equilibration steps to benchmark options */
1648 sprintf(strbuf, "-resetstep %d ", presteps);
1649 add_to_string(cmd_args_bench, strbuf);
1651 /* These switches take effect only at launch time */
1652 if (FALSE == bAppendFiles)
1654 add_to_string(cmd_args_launch, "-noappend ");
1656 if (bKeepAndNumCPT)
1658 add_to_string(cmd_args_launch, "-cpnum ");
1660 if (bResetHWay)
1662 add_to_string(cmd_args_launch, "-resethway ");
1665 /********************/
1666 /* 2. Process files */
1667 /********************/
1668 for (i=0; i<nfile; i++)
1670 opt = (char *)fnm[i].opt;
1671 name = opt2fn(opt,nfile,fnm);
1673 /* Strbuf contains the options, now let's sort out where we need that */
1674 sprintf(strbuf, "%s %s ", opt, name);
1676 if ( is_bench_file(opt, opt2bSet(opt,nfile,fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1678 /* All options starting with -b* need the 'b' removed,
1679 * therefore overwrite strbuf */
1680 if (0 == strncmp(opt, "-b", 2))
1681 sprintf(strbuf, "-%s %s ", &opt[2], name);
1683 add_to_string(cmd_args_bench, strbuf);
1686 if ( is_launch_file(opt,opt2bSet(opt,nfile,fnm)) )
1687 add_to_string(cmd_args_launch, strbuf);
1690 add_to_string(cmd_args_bench , extra_args);
1691 add_to_string(cmd_args_launch, extra_args);
1695 /* Set option opt */
1696 static void setopt(const char *opt,int nfile,t_filenm fnm[])
1698 int i;
1700 for(i=0; (i<nfile); i++)
1701 if (strcmp(opt,fnm[i].opt)==0)
1702 fnm[i].flag |= ffSET;
1706 /* This routine inspects the tpr file and ...
1707 * 1. checks for output files that get triggered by a tpr option. These output
1708 * files are marked as 'set' to allow for a proper cleanup after each
1709 * tuning run.
1710 * 2. returns the PME:PP load ratio
1711 * 3. returns rcoulomb from the tpr */
1712 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1714 gmx_bool bPull; /* Is pulling requested in .tpr file? */
1715 gmx_bool bTpi; /* Is test particle insertion requested? */
1716 gmx_bool bFree; /* Is a free energy simulation requested? */
1717 gmx_bool bNM; /* Is a normal mode analysis requested? */
1718 t_inputrec ir;
1719 t_state state;
1720 gmx_mtop_t mtop;
1723 /* Check tpr file for options that trigger extra output files */
1724 read_tpx_state(opt2fn("-s",nfile,fnm),&ir,&state,NULL,&mtop);
1725 bPull = (epullNO != ir.ePull);
1726 bFree = (efepNO != ir.efep );
1727 bNM = (eiNM == ir.eI );
1728 bTpi = EI_TPI(ir.eI);
1730 /* Set these output files on the tuning command-line */
1731 if (bPull)
1733 setopt("-pf" , nfile, fnm);
1734 setopt("-px" , nfile, fnm);
1736 if (bFree)
1738 setopt("-dhdl", nfile, fnm);
1740 if (bTpi)
1742 setopt("-tpi" , nfile, fnm);
1743 setopt("-tpid", nfile, fnm);
1745 if (bNM)
1747 setopt("-mtx" , nfile, fnm);
1750 *rcoulomb = ir.rcoulomb;
1752 /* Return the estimate for the number of PME nodes */
1753 return pme_load_estimate(&mtop,&ir,state.box);
1757 static void couple_files_options(int nfile, t_filenm fnm[])
1759 int i;
1760 gmx_bool bSet,bBench;
1761 char *opt;
1762 char buf[20];
1765 for (i=0; i<nfile; i++)
1767 opt = (char *)fnm[i].opt;
1768 bSet = ((fnm[i].flag & ffSET) != 0);
1769 bBench = (0 == strncmp(opt,"-b", 2));
1771 /* Check optional files */
1772 /* If e.g. -eo is set, then -beo also needs to be set */
1773 if (is_optional(&fnm[i]) && bSet && !bBench)
1775 sprintf(buf, "-b%s", &opt[1]);
1776 setopt(buf,nfile,fnm);
1778 /* If -beo is set, then -eo also needs to be! */
1779 if (is_optional(&fnm[i]) && bSet && bBench)
1781 sprintf(buf, "-%s", &opt[2]);
1782 setopt(buf,nfile,fnm);
1788 static double gettime()
1790 #ifdef HAVE_GETTIMEOFDAY
1791 struct timeval t;
1792 double seconds;
1794 gettimeofday(&t,NULL);
1796 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
1798 return seconds;
1799 #else
1800 double seconds;
1802 seconds = time(NULL);
1804 return seconds;
1805 #endif
1809 #define BENCHSTEPS (1000)
1811 int gmx_tune_pme(int argc,char *argv[])
1813 const char *desc[] = {
1814 "For a given number [TT]-np[tt] or [TT]-nt[tt] of processors/threads, this program systematically",
1815 "times [TT]mdrun[tt] with various numbers of PME-only nodes and determines",
1816 "which setting is fastest. It will also test whether performance can",
1817 "be enhanced by shifting load from the reciprocal to the real space",
1818 "part of the Ewald sum. ",
1819 "Simply pass your [TT].tpr[tt] file to [TT]g_tune_pme[tt] together with other options",
1820 "for [TT]mdrun[tt] as needed.[PAR]",
1821 "Which executables are used can be set in the environment variables",
1822 "MPIRUN and MDRUN. If these are not present, 'mpirun' and 'mdrun'",
1823 "will be used as defaults. Note that for certain MPI frameworks you",
1824 "need to provide a machine- or hostfile. This can also be passed",
1825 "via the MPIRUN variable, e.g.[PAR]",
1826 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt][PAR]",
1827 "Please call [TT]g_tune_pme[tt] with the normal options you would pass to",
1828 "[TT]mdrun[tt] and add [TT]-np[tt] for the number of processors to perform the",
1829 "tests on, or [TT]-nt[tt] for the number of threads. You can also add [TT]-r[tt]",
1830 "to repeat each test several times to get better statistics. [PAR]",
1831 "[TT]g_tune_pme[tt] can test various real space / reciprocal space workloads",
1832 "for you. With [TT]-ntpr[tt] you control how many extra [TT].tpr[tt] files will be",
1833 "written with enlarged cutoffs and smaller Fourier grids respectively.",
1834 "Typically, the first test (number 0) will be with the settings from the input",
1835 "[TT].tpr[tt] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
1836 "specified by [TT]-rmax[tt] with a somwhat smaller PME grid at the same time. ",
1837 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
1838 "The remaining [TT].tpr[tt] files will have equally-spaced Coulomb radii (and Fourier "
1839 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
1840 "if you just seek the optimal number of PME-only nodes; in that case",
1841 "your input [TT].tpr[tt] file will remain unchanged.[PAR]",
1842 "For the benchmark runs, the default of 1000 time steps should suffice for most",
1843 "MD systems. The dynamic load balancing needs about 100 time steps",
1844 "to adapt to local load imbalances, therefore the time step counters",
1845 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
1846 "for a higher accuarcy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
1847 "From the 'DD' load imbalance entries in the md.log output file you",
1848 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
1849 "[TT]g_tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
1850 "After calling [TT]mdrun[tt] several times, detailed performance information",
1851 "is available in the output file [TT]perf.out.[tt] ",
1852 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
1853 "(options [TT]-b[tt]*), these will be automatically deleted after each test.[PAR]",
1854 "If you want the simulation to be started automatically with the",
1855 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
1858 int nnodes =1;
1859 int repeats=2;
1860 int pmeentries=0; /* How many values for -npme do we actually test for each tpr file */
1861 real maxPMEfraction=0.50;
1862 real minPMEfraction=0.25;
1863 int maxPMEnodes, minPMEnodes;
1864 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
1865 float guessPMEnodes;
1866 int npme_fixed=-2; /* If >= -1, use only this number
1867 * of PME-only nodes */
1868 int ntprs=0;
1869 real rmin=0.0,rmax=0.0; /* min and max value for rcoulomb if scaling is requested */
1870 real rcoulomb=-1.0; /* Coulomb radius as set in .tpr file */
1871 gmx_bool bScaleRvdw=TRUE;
1872 gmx_large_int_t bench_nsteps=BENCHSTEPS;
1873 gmx_large_int_t new_sim_nsteps=-1; /* -1 indicates: not set by the user */
1874 gmx_large_int_t cpt_steps=0; /* Step counter in .cpt input file */
1875 int presteps=100; /* Do a full cycle reset after presteps steps */
1876 gmx_bool bOverwrite=FALSE, bKeepTPR;
1877 gmx_bool bLaunch=FALSE;
1878 char *ExtraArgs=NULL;
1879 char **tpr_names=NULL;
1880 const char *simulation_tpr=NULL;
1881 int best_npme, best_tpr;
1882 int sim_part = 1; /* For benchmarks with checkpoint files */
1883 char bbuf[STRLEN];
1885 /* Default program names if nothing else is found */
1886 char *cmd_mpirun=NULL, *cmd_mdrun=NULL;
1887 char *cmd_args_bench, *cmd_args_launch;
1888 char *cmd_np=NULL;
1890 t_perf **perfdata=NULL;
1891 t_inputinfo *info;
1892 int i;
1893 FILE *fp;
1894 t_commrec *cr;
1896 /* Print out how long the tuning took */
1897 double seconds;
1899 static t_filenm fnm[] = {
1900 /* g_tune_pme */
1901 { efOUT, "-p", "perf", ffWRITE },
1902 { efLOG, "-err", "bencherr", ffWRITE },
1903 { efTPX, "-so", "tuned", ffWRITE },
1904 /* mdrun: */
1905 { efTPX, NULL, NULL, ffREAD },
1906 { efTRN, "-o", NULL, ffWRITE },
1907 { efXTC, "-x", NULL, ffOPTWR },
1908 { efCPT, "-cpi", NULL, ffOPTRD },
1909 { efCPT, "-cpo", NULL, ffOPTWR },
1910 { efSTO, "-c", "confout", ffWRITE },
1911 { efEDR, "-e", "ener", ffWRITE },
1912 { efLOG, "-g", "md", ffWRITE },
1913 { efXVG, "-dhdl", "dhdl", ffOPTWR },
1914 { efXVG, "-field", "field", ffOPTWR },
1915 { efXVG, "-table", "table", ffOPTRD },
1916 { efXVG, "-tabletf", "tabletf", ffOPTRD },
1917 { efXVG, "-tablep", "tablep", ffOPTRD },
1918 { efXVG, "-tableb", "table", ffOPTRD },
1919 { efTRX, "-rerun", "rerun", ffOPTRD },
1920 { efXVG, "-tpi", "tpi", ffOPTWR },
1921 { efXVG, "-tpid", "tpidist", ffOPTWR },
1922 { efEDI, "-ei", "sam", ffOPTRD },
1923 { efEDO, "-eo", "sam", ffOPTWR },
1924 { efGCT, "-j", "wham", ffOPTRD },
1925 { efGCT, "-jo", "bam", ffOPTWR },
1926 { efXVG, "-ffout", "gct", ffOPTWR },
1927 { efXVG, "-devout", "deviatie", ffOPTWR },
1928 { efXVG, "-runav", "runaver", ffOPTWR },
1929 { efXVG, "-px", "pullx", ffOPTWR },
1930 { efXVG, "-pf", "pullf", ffOPTWR },
1931 { efXVG, "-ro", "rotation", ffOPTWR },
1932 { efLOG, "-ra", "rotangles",ffOPTWR },
1933 { efLOG, "-rs", "rotslabs", ffOPTWR },
1934 { efLOG, "-rt", "rottorque",ffOPTWR },
1935 { efMTX, "-mtx", "nm", ffOPTWR },
1936 { efNDX, "-dn", "dipole", ffOPTWR },
1937 /* Output files that are deleted after each benchmark run */
1938 { efTRN, "-bo", "bench", ffWRITE },
1939 { efXTC, "-bx", "bench", ffWRITE },
1940 { efCPT, "-bcpo", "bench", ffWRITE },
1941 { efSTO, "-bc", "bench", ffWRITE },
1942 { efEDR, "-be", "bench", ffWRITE },
1943 { efLOG, "-bg", "bench", ffWRITE },
1944 { efEDO, "-beo", "bench", ffOPTWR },
1945 { efXVG, "-bdhdl", "benchdhdl",ffOPTWR },
1946 { efXVG, "-bfield", "benchfld" ,ffOPTWR },
1947 { efXVG, "-btpi", "benchtpi", ffOPTWR },
1948 { efXVG, "-btpid", "benchtpid",ffOPTWR },
1949 { efGCT, "-bjo", "bench", ffOPTWR },
1950 { efXVG, "-bffout", "benchgct", ffOPTWR },
1951 { efXVG, "-bdevout","benchdev", ffOPTWR },
1952 { efXVG, "-brunav", "benchrnav",ffOPTWR },
1953 { efXVG, "-bpx", "benchpx", ffOPTWR },
1954 { efXVG, "-bpf", "benchpf", ffOPTWR },
1955 { efXVG, "-bro", "benchrot", ffOPTWR },
1956 { efLOG, "-bra", "benchrota",ffOPTWR },
1957 { efLOG, "-brs", "benchrots",ffOPTWR },
1958 { efLOG, "-brt", "benchrott",ffOPTWR },
1959 { efMTX, "-bmtx", "benchn", ffOPTWR },
1960 { efNDX, "-bdn", "bench", ffOPTWR }
1963 gmx_bool bThreads = FALSE;
1965 int nthreads=1;
1967 const char *procstring[] =
1968 { NULL, "-np", "-n", "none", NULL };
1969 const char *npmevalues_opt[] =
1970 { NULL, "auto", "all", "subset", NULL };
1972 gmx_bool bAppendFiles=TRUE;
1973 gmx_bool bKeepAndNumCPT=FALSE;
1974 gmx_bool bResetCountersHalfWay=FALSE;
1975 gmx_bool bBenchmark=TRUE;
1977 output_env_t oenv=NULL;
1979 t_pargs pa[] = {
1980 /***********************/
1981 /* g_tune_pme options: */
1982 /***********************/
1983 { "-np", FALSE, etINT, {&nnodes},
1984 "Number of nodes to run the tests on (must be > 2 for separate PME nodes)" },
1985 { "-npstring", FALSE, etENUM, {procstring},
1986 "Specify the number of processors to [TT]$MPIRUN[tt] using this string"},
1987 { "-nt", FALSE, etINT, {&nthreads},
1988 "Number of threads to run the tests on (turns MPI & mpirun off)"},
1989 { "-r", FALSE, etINT, {&repeats},
1990 "Repeat each test this often" },
1991 { "-max", FALSE, etREAL, {&maxPMEfraction},
1992 "Max fraction of PME nodes to test with" },
1993 { "-min", FALSE, etREAL, {&minPMEfraction},
1994 "Min fraction of PME nodes to test with" },
1995 { "-npme", FALSE, etENUM, {npmevalues_opt},
1996 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
1997 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
1998 { "-fix", FALSE, etINT, {&npme_fixed},
1999 "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."},
2000 { "-rmax", FALSE, etREAL, {&rmax},
2001 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2002 { "-rmin", FALSE, etREAL, {&rmin},
2003 "If >0, minimal rcoulomb for -ntpr>1" },
2004 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2005 "Scale rvdw along with rcoulomb"},
2006 { "-ntpr", FALSE, etINT, {&ntprs},
2007 "Number of [TT].tpr[tt] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2008 "If < 1, automatically choose the number of [TT].tpr[tt] files to test" },
2009 { "-steps", FALSE, etGMX_LARGE_INT, {&bench_nsteps},
2010 "Take timings for this many steps in the benchmark runs" },
2011 { "-resetstep",FALSE, etINT, {&presteps},
2012 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2013 { "-simsteps", FALSE, etGMX_LARGE_INT, {&new_sim_nsteps},
2014 "If non-negative, perform this many steps in the real run (overwrites nsteps from [TT].tpr[tt], add [TT].cpt[tt] steps)" },
2015 { "-launch", FALSE, etBOOL, {&bLaunch},
2016 "Launch the real simulation after optimization" },
2017 { "-bench", FALSE, etBOOL, {&bBenchmark},
2018 "Run the benchmarks or just create the input [TT].tpr[tt] files?" },
2019 /******************/
2020 /* mdrun options: */
2021 /******************/
2022 /* We let g_tune_pme parse and understand these options, because we need to
2023 * prevent that they appear on the mdrun command line for the benchmarks */
2024 { "-append", FALSE, etBOOL, {&bAppendFiles},
2025 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2026 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2027 "Keep and number checkpoint files (launch only)" },
2028 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2029 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2033 #define NFILE asize(fnm)
2035 CopyRight(stderr,argv[0]);
2037 seconds = gettime();
2039 parse_common_args(&argc,argv,PCA_NOEXIT_ON_ARGS,
2040 NFILE,fnm,asize(pa),pa,asize(desc),desc,
2041 0,NULL,&oenv);
2043 /* Store the remaining unparsed command line entries in a string which
2044 * is then attached to the mdrun command line */
2045 snew(ExtraArgs, 1);
2046 ExtraArgs[0] = '\0';
2047 for (i=1; i<argc; i++) /* argc will now be 1 if everything was understood */
2049 add_to_string(&ExtraArgs, argv[i]);
2050 add_to_string(&ExtraArgs, " ");
2053 if (opt2parg_bSet("-nt",asize(pa),pa))
2055 bThreads=TRUE;
2056 if (opt2parg_bSet("-npstring",asize(pa),pa))
2057 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2059 if (nnodes > 1)
2060 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2061 /* and now we just set this; a bit of an ugly hack*/
2062 nnodes=nthreads;
2064 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2065 guessPMEratio = inspect_tpr(NFILE,fnm,&rcoulomb);
2067 /* Automatically set -beo options if -eo is set etc. */
2068 couple_files_options(NFILE,fnm);
2070 /* Construct the command line arguments for benchmark runs
2071 * as well as for the simulation run */
2072 if (bThreads)
2073 sprintf(bbuf," -nt %d ", nthreads);
2074 else
2075 sprintf(bbuf," -np %d ", nnodes);
2077 cmd_np = bbuf;
2079 create_command_line_snippets(bThreads,bAppendFiles,bKeepAndNumCPT,bResetCountersHalfWay,presteps,
2080 NFILE,fnm,asize(pa),pa,procstring[0],
2081 &cmd_np, &cmd_args_bench, &cmd_args_launch,
2082 ExtraArgs);
2084 /* Read in checkpoint file if requested */
2085 sim_part = 1;
2086 if(opt2bSet("-cpi",NFILE,fnm))
2088 snew(cr,1);
2089 cr->duty=DUTY_PP; /* makes the following routine happy */
2090 read_checkpoint_simulation_part(opt2fn("-cpi",NFILE,fnm),
2091 &sim_part,&cpt_steps,cr,
2092 FALSE,NFILE,fnm,NULL,NULL);
2093 sfree(cr);
2094 sim_part++;
2095 /* sim_part will now be 1 if no checkpoint file was found */
2096 if (sim_part<=1)
2097 gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi",NFILE,fnm));
2100 /* Open performance output file and write header info */
2101 fp = ffopen(opt2fn("-p",NFILE,fnm),"w");
2103 /* Make a quick consistency check of command line parameters */
2104 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2105 maxPMEfraction, minPMEfraction, npme_fixed,
2106 bench_nsteps, fnm, NFILE, sim_part, presteps,
2107 asize(pa),pa);
2109 /* Determine the maximum and minimum number of PME nodes to test,
2110 * the actual list of settings is build in do_the_tests(). */
2111 if ((nnodes > 2) && (npme_fixed < -1))
2113 if (0 == strcmp(npmevalues_opt[0], "auto"))
2115 /* Determine the npme range automatically based on the PME:PP load guess */
2116 if (guessPMEratio > 1.0)
2118 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2119 maxPMEnodes=nnodes/2;
2120 minPMEnodes=nnodes/2;
2122 else
2124 /* PME : PP load is in the range 0..1, let's test around the guess */
2125 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2126 minPMEnodes = floor(0.7*guessPMEnodes);
2127 maxPMEnodes = ceil(1.6*guessPMEnodes);
2128 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2131 else
2133 /* Determine the npme range based on user input */
2134 maxPMEnodes = floor(maxPMEfraction*nnodes);
2135 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2136 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2137 if (maxPMEnodes != minPMEnodes)
2138 fprintf(stdout, "- %d ", maxPMEnodes);
2139 fprintf(stdout, "PME-only nodes.\n Note that the automatic number of PME-only nodes and no separate PME nodes are always tested.\n");
2142 else
2144 maxPMEnodes = 0;
2145 minPMEnodes = 0;
2148 /* Get the commands we need to set up the runs from environment variables */
2149 get_program_paths(bThreads, &cmd_mpirun, cmd_np, &cmd_mdrun, repeats);
2151 /* Print some header info to file */
2152 sep_line(fp);
2153 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2154 sep_line(fp);
2155 fprintf(fp, "%s for Gromacs %s\n", ShortProgram(),GromacsVersion());
2156 if (!bThreads)
2158 fprintf(fp, "Number of nodes : %d\n", nnodes);
2159 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2160 if ( strcmp(procstring[0], "none") != 0)
2161 fprintf(fp, "Passing # of nodes via : %s\n", procstring[0]);
2162 else
2163 fprintf(fp, "Not setting number of nodes in system call\n");
2165 else
2166 fprintf(fp, "Number of threads : %d\n", nnodes);
2168 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2169 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2170 fprintf(fp, "Benchmark steps : ");
2171 fprintf(fp, gmx_large_int_pfmt, bench_nsteps);
2172 fprintf(fp, "\n");
2173 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2174 if (sim_part > 1)
2176 fprintf(fp, "Checkpoint time step : ");
2177 fprintf(fp, gmx_large_int_pfmt, cpt_steps);
2178 fprintf(fp, "\n");
2180 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2182 if (new_sim_nsteps >= 0)
2184 bOverwrite = TRUE;
2185 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so",NFILE,fnm));
2186 fprintf(stderr, gmx_large_int_pfmt, new_sim_nsteps+cpt_steps);
2187 fprintf(stderr, " steps.\n");
2188 fprintf(fp, "Simulation steps : ");
2189 fprintf(fp, gmx_large_int_pfmt, new_sim_nsteps);
2190 fprintf(fp, "\n");
2192 if (repeats > 1)
2193 fprintf(fp, "Repeats for each test : %d\n", repeats);
2195 if (npme_fixed >= -1)
2197 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2200 fprintf(fp, "Input file : %s\n", opt2fn("-s",NFILE,fnm));
2201 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2203 /* Allocate memory for the inputinfo struct: */
2204 snew(info, 1);
2205 info->nr_inputfiles = ntprs;
2206 for (i=0; i<ntprs; i++)
2208 snew(info->rcoulomb , ntprs);
2209 snew(info->rvdw , ntprs);
2210 snew(info->rlist , ntprs);
2211 snew(info->rlistlong, ntprs);
2212 snew(info->nkx , ntprs);
2213 snew(info->nky , ntprs);
2214 snew(info->nkz , ntprs);
2215 snew(info->fsx , ntprs);
2216 snew(info->fsy , ntprs);
2217 snew(info->fsz , ntprs);
2219 /* Make alternative tpr files to test: */
2220 snew(tpr_names, ntprs);
2221 for (i=0; i<ntprs; i++)
2222 snew(tpr_names[i], STRLEN);
2224 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2225 * different grids could be found. */
2226 make_benchmark_tprs(opt2fn("-s",NFILE,fnm), tpr_names, bench_nsteps+presteps,
2227 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2229 /********************************************************************************/
2230 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2231 /********************************************************************************/
2232 snew(perfdata, ntprs);
2233 if (bBenchmark)
2235 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2236 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2237 cmd_args_bench, fnm, NFILE, sim_part, presteps, cpt_steps);
2239 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gettime()-seconds)/60.0);
2241 /* Analyse the results and give a suggestion for optimal settings: */
2242 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2243 repeats, info, &best_tpr, &best_npme);
2245 /* Take the best-performing tpr file and enlarge nsteps to original value */
2246 if ( bKeepTPR && !bOverwrite )
2248 simulation_tpr = opt2fn("-s",NFILE,fnm);
2250 else
2252 simulation_tpr = opt2fn("-so",NFILE,fnm);
2253 modify_PMEsettings(bOverwrite? (new_sim_nsteps+cpt_steps) :
2254 info->orig_sim_steps, tpr_names[best_tpr],
2255 simulation_tpr);
2258 /* Let's get rid of the temporary benchmark input files */
2259 for (i=0; i<ntprs; i++)
2261 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2262 remove(tpr_names[i]);
2265 /* Now start the real simulation if the user requested it ... */
2266 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2267 cmd_args_launch, simulation_tpr, nnodes, best_npme);
2269 ffclose(fp);
2271 /* ... or simply print the performance results to screen: */
2272 if (!bLaunch)
2273 finalize(opt2fn("-p", NFILE, fnm));
2275 return 0;