Spelling fixes
[gromacs.git] / src / gromacs / gmxana / gmx_tune_pme.c
blobcfb6927ce93ec241f96ace209ade515492cdfab5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 #include "gmxpre.h"
37 #include "config.h"
39 #include <stdlib.h>
40 #include <time.h>
42 #ifdef HAVE_SYS_TIME_H
43 #include <sys/time.h>
44 #endif
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/tpxio.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/legacyheaders/calcgrid.h"
50 #include "gromacs/legacyheaders/checkpoint.h"
51 #include "gromacs/legacyheaders/inputrec.h"
52 #include "gromacs/legacyheaders/macros.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/perf_est.h"
55 #include "gromacs/legacyheaders/readinp.h"
56 #include "gromacs/legacyheaders/typedefs.h"
57 #include "gromacs/legacyheaders/types/commrec.h"
58 #include "gromacs/math/utilities.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/timing/walltime_accounting.h"
61 #include "gromacs/utility/baseversion.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/smalloc.h"
66 /* Enum for situations that can occur during log file parsing, the
67 * corresponding string entries can be found in do_the_tests() in
68 * const char* ParseLog[] */
69 /* TODO clean up CamelCasing of these enum names */
70 enum {
71 eParselogOK,
72 eParselogNotFound,
73 eParselogNoPerfData,
74 eParselogTerm,
75 eParselogResetProblem,
76 eParselogNoDDGrid,
77 eParselogTPXVersion,
78 eParselogNotParallel,
79 eParselogLargePrimeFactor,
80 eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs,
81 eParselogGpuProblem,
82 eParselogFatal,
83 eParselogNr
87 typedef struct
89 int nPMEnodes; /* number of PME-only nodes used in this test */
90 int nx, ny, nz; /* DD grid */
91 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
92 double *Gcycles; /* This can contain more than one value if doing multiple tests */
93 double Gcycles_Av;
94 float *ns_per_day;
95 float ns_per_day_Av;
96 float *PME_f_load; /* PME mesh/force load average*/
97 float PME_f_load_Av; /* Average average ;) ... */
98 char *mdrun_cmd_line; /* Mdrun command line used for this test */
99 } t_perf;
102 typedef struct
104 int nr_inputfiles; /* The number of tpr and mdp input files */
105 gmx_int64_t orig_sim_steps; /* Number of steps to be done in the real simulation */
106 gmx_int64_t orig_init_step; /* Init step for the real simulation */
107 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
108 real *rvdw; /* The vdW radii */
109 real *rlist; /* Neighbourlist cutoff radius */
110 real *rlistlong;
111 int *nkx, *nky, *nkz;
112 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
113 } t_inputinfo;
116 static void sep_line(FILE *fp)
118 fprintf(fp, "\n------------------------------------------------------------\n");
122 /* Wrapper for system calls */
123 static int gmx_system_call(char *command)
125 return ( system(command) );
129 /* Check if string starts with substring */
130 static gmx_bool str_starts(const char *string, const char *substring)
132 return ( strncmp(string, substring, strlen(substring)) == 0);
136 static void cleandata(t_perf *perfdata, int test_nr)
138 perfdata->Gcycles[test_nr] = 0.0;
139 perfdata->ns_per_day[test_nr] = 0.0;
140 perfdata->PME_f_load[test_nr] = 0.0;
142 return;
146 static void remove_if_exists(const char *fn)
148 if (gmx_fexist(fn))
150 fprintf(stdout, "Deleting %s\n", fn);
151 remove(fn);
156 static void finalize(const char *fn_out)
158 char buf[STRLEN];
159 FILE *fp;
162 fp = fopen(fn_out, "r");
163 fprintf(stdout, "\n\n");
165 while (fgets(buf, STRLEN-1, fp) != NULL)
167 fprintf(stdout, "%s", buf);
169 fclose(fp);
170 fprintf(stdout, "\n\n");
174 enum {
175 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_int64_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 float dum1, dum2, dum3, dum4;
191 int ndum;
192 int npme;
193 gmx_int64_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)
212 iFound = eFoundDDStr; /* Skip some case statements */
215 while (fgets(line, STRLEN, fp) != NULL)
217 /* Remove leading spaces */
218 ltrim(line);
220 /* Check for TERM and INT signals from user: */
221 if (strstr(line, errSIG) != NULL)
223 fclose(fp);
224 cleandata(perfdata, test_nr);
225 return eParselogTerm;
228 /* Check whether cycle resetting worked */
229 if (presteps > 0 && !bFoundResetStr)
231 if (strstr(line, matchstrcr) != NULL)
233 sprintf(dumstring, "step %s", "%"GMX_SCNd64);
234 sscanf(line, dumstring, &resetsteps);
235 bFoundResetStr = TRUE;
236 if (resetsteps == presteps+cpt_steps)
238 bResetChecked = TRUE;
240 else
242 sprintf(dumstring, "%"GMX_PRId64, resetsteps);
243 sprintf(dumstring2, "%"GMX_PRId64, presteps+cpt_steps);
244 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
245 " though they were supposed to be reset at step %s!\n",
246 dumstring, dumstring2);
251 /* Look for strings that appear in a certain order in the log file: */
252 switch (iFound)
254 case eFoundNothing:
255 /* Look for domain decomp grid and separate PME nodes: */
256 if (str_starts(line, matchstrdd))
258 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME ranks %d",
259 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
260 if (perfdata->nPMEnodes == -1)
262 perfdata->guessPME = npme;
264 else if (perfdata->nPMEnodes != npme)
266 gmx_fatal(FARGS, "PME ranks from command line and output file are not identical");
268 iFound = eFoundDDStr;
270 /* Catch a few errors that might have occurred: */
271 else if (str_starts(line, "There is no domain decomposition for"))
273 fclose(fp);
274 return eParselogNoDDGrid;
276 else if (str_starts(line, "The number of ranks you selected"))
278 fclose(fp);
279 return eParselogLargePrimeFactor;
281 else if (str_starts(line, "reading tpx file"))
283 fclose(fp);
284 return eParselogTPXVersion;
286 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
288 fclose(fp);
289 return eParselogNotParallel;
291 break;
292 case eFoundDDStr:
293 /* Even after the "Domain decomposition grid" string was found,
294 * it could be that mdrun had to quit due to some error. */
295 if (str_starts(line, "Incorrect launch configuration: mismatching number of"))
297 fclose(fp);
298 return eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs;
300 else if (str_starts(line, "Some of the requested GPUs do not exist"))
302 fclose(fp);
303 return eParselogGpuProblem;
305 /* Look for PME mesh/force balance (not necessarily present, though) */
306 else if (str_starts(line, matchstrbal))
308 sscanf(&line[strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
310 /* Look for matchstring */
311 else if (str_starts(line, matchstring))
313 iFound = eFoundAccountingStr;
315 break;
316 case eFoundAccountingStr:
317 /* Already found matchstring - look for cycle data */
318 if (str_starts(line, "Total "))
320 sscanf(line, "Total %*f %lf", &(perfdata->Gcycles[test_nr]));
321 iFound = eFoundCycleStr;
323 break;
324 case eFoundCycleStr:
325 /* Already found cycle data - look for remaining performance info and return */
326 if (str_starts(line, "Performance:"))
328 ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
329 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
330 perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
331 fclose(fp);
332 if (bResetChecked || presteps == 0)
334 return eParselogOK;
336 else
338 return eParselogResetProblem;
341 break;
343 } /* while */
345 /* Close the log file */
346 fclose(fp);
348 /* Check why there is no performance data in the log file.
349 * Did a fatal errors occur? */
350 if (gmx_fexist(errfile))
352 fp = fopen(errfile, "r");
353 while (fgets(line, STRLEN, fp) != NULL)
355 if (str_starts(line, "Fatal error:") )
357 if (fgets(line, STRLEN, fp) != NULL)
359 fprintf(stderr, "\nWARNING: An error occurred during this benchmark:\n"
360 "%s\n", line);
362 fclose(fp);
363 cleandata(perfdata, test_nr);
364 return eParselogFatal;
367 fclose(fp);
369 else
371 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
374 /* Giving up ... we could not find out why there is no performance data in
375 * the log file. */
376 fprintf(stdout, "No performance data in log file.\n");
377 cleandata(perfdata, test_nr);
379 return eParselogNoPerfData;
383 static gmx_bool analyze_data(
384 FILE *fp,
385 const char *fn,
386 t_perf **perfdata,
387 int nnodes,
388 int ntprs,
389 int ntests,
390 int nrepeats,
391 t_inputinfo *info,
392 int *index_tpr, /* OUT: Nr of mdp file with best settings */
393 int *npme_optimal) /* OUT: Optimal number of PME nodes */
395 int i, j, k;
396 int line = 0, line_win = -1;
397 int k_win = -1, i_win = -1, winPME;
398 double s = 0.0; /* standard deviation */
399 t_perf *pd;
400 char strbuf[STRLEN];
401 char str_PME_f_load[13];
402 gmx_bool bCanUseOrigTPR;
403 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
406 if (nrepeats > 1)
408 sep_line(fp);
409 fprintf(fp, "Summary of successful runs:\n");
410 fprintf(fp, "Line tpr PME ranks Gcycles Av. Std.dev. ns/day PME/f");
411 if (nnodes > 1)
413 fprintf(fp, " DD grid");
415 fprintf(fp, "\n");
419 for (k = 0; k < ntprs; k++)
421 for (i = 0; i < ntests; i++)
423 /* Select the right dataset: */
424 pd = &(perfdata[k][i]);
426 pd->Gcycles_Av = 0.0;
427 pd->PME_f_load_Av = 0.0;
428 pd->ns_per_day_Av = 0.0;
430 if (pd->nPMEnodes == -1)
432 sprintf(strbuf, "(%3d)", pd->guessPME);
434 else
436 sprintf(strbuf, " ");
439 /* Get the average run time of a setting */
440 for (j = 0; j < nrepeats; j++)
442 pd->Gcycles_Av += pd->Gcycles[j];
443 pd->PME_f_load_Av += pd->PME_f_load[j];
445 pd->Gcycles_Av /= nrepeats;
446 pd->PME_f_load_Av /= nrepeats;
448 for (j = 0; j < nrepeats; j++)
450 if (pd->ns_per_day[j] > 0.0)
452 pd->ns_per_day_Av += pd->ns_per_day[j];
454 else
456 /* Somehow the performance number was not aquired for this run,
457 * therefor set the average to some negative value: */
458 pd->ns_per_day_Av = -1.0f*nrepeats;
459 break;
462 pd->ns_per_day_Av /= nrepeats;
464 /* Nicer output: */
465 if (pd->PME_f_load_Av > 0.0)
467 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
469 else
471 sprintf(str_PME_f_load, "%s", " - ");
475 /* We assume we had a successful run if both averages are positive */
476 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
478 /* Output statistics if repeats were done */
479 if (nrepeats > 1)
481 /* Calculate the standard deviation */
482 s = 0.0;
483 for (j = 0; j < nrepeats; j++)
485 s += pow( pd->Gcycles[j] - pd->Gcycles_Av, 2 );
487 s /= (nrepeats - 1);
488 s = sqrt(s);
490 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
491 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
492 pd->ns_per_day_Av, str_PME_f_load);
493 if (nnodes > 1)
495 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
497 fprintf(fp, "\n");
499 /* Store the index of the best run found so far in 'winner': */
500 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
502 k_win = k;
503 i_win = i;
504 line_win = line;
506 line++;
511 if (k_win == -1)
513 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
516 sep_line(fp);
518 winPME = perfdata[k_win][i_win].nPMEnodes;
520 if (1 == ntests)
522 /* We stuck to a fixed number of PME-only nodes */
523 sprintf(strbuf, "settings No. %d", k_win);
525 else
527 /* We have optimized the number of PME-only nodes */
528 if (winPME == -1)
530 sprintf(strbuf, "%s", "the automatic number of PME ranks");
532 else
534 sprintf(strbuf, "%d PME ranks", winPME);
537 fprintf(fp, "Best performance was achieved with %s", strbuf);
538 if ((nrepeats > 1) && (ntests > 1))
540 fprintf(fp, " (see line %d)", line_win);
542 fprintf(fp, "\n");
544 /* Only mention settings if they were modified: */
545 bRefinedCoul = !gmx_within_tol(info->rcoulomb[k_win], info->rcoulomb[0], GMX_REAL_EPS);
546 bRefinedVdW = !gmx_within_tol(info->rvdw[k_win], info->rvdw[0], GMX_REAL_EPS);
547 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
548 info->nky[k_win] == info->nky[0] &&
549 info->nkz[k_win] == info->nkz[0]);
551 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
553 fprintf(fp, "Optimized PME settings:\n");
554 bCanUseOrigTPR = FALSE;
556 else
558 bCanUseOrigTPR = TRUE;
561 if (bRefinedCoul)
563 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
566 if (bRefinedVdW)
568 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
571 if (bRefinedGrid)
573 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],
574 info->nkx[0], info->nky[0], info->nkz[0]);
577 if (bCanUseOrigTPR && ntprs > 1)
579 fprintf(fp, "and original PME settings.\n");
582 fflush(fp);
584 /* Return the index of the mdp file that showed the highest performance
585 * and the optimal number of PME nodes */
586 *index_tpr = k_win;
587 *npme_optimal = winPME;
589 return bCanUseOrigTPR;
593 /* Get the commands we need to set up the runs from environment variables */
594 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
596 char *cp;
597 FILE *fp;
598 const char def_mpirun[] = "mpirun";
600 const char empty_mpirun[] = "";
602 /* Get the commands we need to set up the runs from environment variables */
603 if (!bThreads)
605 if ( (cp = getenv("MPIRUN")) != NULL)
607 *cmd_mpirun = gmx_strdup(cp);
609 else
611 *cmd_mpirun = gmx_strdup(def_mpirun);
614 else
616 *cmd_mpirun = gmx_strdup(empty_mpirun);
619 if (*cmd_mdrun == NULL)
621 /* The use of MDRUN is deprecated, but made available in 5.1
622 for backward compatibility. It may be removed in a future
623 version. */
624 if ( (cp = getenv("MDRUN" )) != NULL)
626 *cmd_mdrun = gmx_strdup(cp);
628 else
630 gmx_fatal(FARGS, "The way to call mdrun must be set in the -mdrun command-line flag.");
635 /* Check that the commands will run mdrun (perhaps via mpirun) by
636 * running a very quick test simulation. Requires MPI environment or
637 * GPU support to be available if applicable. */
638 /* TODO implement feature to parse the log file to get the list of
639 compatible GPUs from mdrun, if the user of gmx tune-pme has not
640 given one. */
641 static void check_mdrun_works(gmx_bool bThreads,
642 const char *cmd_mpirun,
643 const char *cmd_np,
644 const char *cmd_mdrun,
645 gmx_bool bNeedGpuSupport)
647 char *command = NULL;
648 char *cp;
649 char line[STRLEN];
650 FILE *fp;
651 const char filename[] = "benchtest.log";
653 /* This string should always be identical to the one in copyrite.c,
654 * gmx_print_version_info() in the defined(GMX_MPI) section */
655 const char match_mpi[] = "MPI library: MPI";
656 const char match_mdrun[] = "Executable: ";
657 const char match_gpu[] = "GPU support: enabled";
658 gmx_bool bMdrun = FALSE;
659 gmx_bool bMPI = FALSE;
660 gmx_bool bHaveGpuSupport = FALSE;
662 /* Run a small test to see whether mpirun + mdrun work */
663 fprintf(stdout, "Making sure that mdrun can be executed. ");
664 if (bThreads)
666 snew(command, strlen(cmd_mdrun) + strlen(cmd_np) + strlen(filename) + 50);
667 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
669 else
671 snew(command, strlen(cmd_mpirun) + strlen(cmd_np) + strlen(cmd_mdrun) + strlen(filename) + 50);
672 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
674 fprintf(stdout, "Trying '%s' ... ", command);
675 make_backup(filename);
676 gmx_system_call(command);
678 /* Check if we find the characteristic string in the output: */
679 if (!gmx_fexist(filename))
681 gmx_fatal(FARGS, "Output from test run could not be found.");
684 fp = fopen(filename, "r");
685 /* We need to scan the whole output file, since sometimes the queuing system
686 * also writes stuff to stdout/err */
687 while (!feof(fp) )
689 cp = fgets(line, STRLEN, fp);
690 if (cp != NULL)
692 if (str_starts(line, match_mdrun) )
694 bMdrun = TRUE;
696 if (str_starts(line, match_mpi) )
698 bMPI = TRUE;
700 if (str_starts(line, match_gpu) )
702 bHaveGpuSupport = TRUE;
706 fclose(fp);
708 if (bThreads)
710 if (bMPI)
712 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
713 "(%s)\n"
714 "seems to have been compiled with MPI instead.",
715 cmd_mdrun);
718 else
720 if (bMdrun && !bMPI)
722 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
723 "(%s)\n"
724 "seems to have been compiled without MPI support.",
725 cmd_mdrun);
729 if (!bMdrun)
731 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
732 filename);
735 if (bNeedGpuSupport && !bHaveGpuSupport)
737 gmx_fatal(FARGS, "The mdrun executable did not have the expected GPU support.");
740 fprintf(stdout, "passed.\n");
742 /* Clean up ... */
743 remove(filename);
744 sfree(command);
747 /*! \brief Helper struct so we can parse the string with eligible GPU
748 IDs outside do_the_tests. */
749 typedef struct eligible_gpu_ids
751 int n; /**< Length of ids */
752 int *ids; /**< Array of length n. NULL if no GPUs in use */
753 } t_eligible_gpu_ids;
755 /* Handles the no-GPU case by emitting an empty string. */
756 static char *make_gpu_id_command_line(int numRanks, int numPmeRanks, const t_eligible_gpu_ids *gpu_ids)
758 char *command_line, *flag = "-gpu_id ", *ptr;
759 int flag_length;
761 /* Reserve enough room for the option name, enough single-digit
762 GPU ids (since that is currently all that is possible to use
763 with mdrun), and a terminating NULL. */
764 flag_length = strlen(flag);
765 snew(command_line, flag_length + numRanks + 1);
766 ptr = command_line;
768 /* If the user has given no eligible GPU IDs, or we're trying the
769 * default behaviour, then there is nothing for g_tune_pme to give
770 * to mdrun -gpu_id */
771 if (gpu_ids->n > 0 && numPmeRanks > -1)
773 int numPpRanks, max_num_ranks_for_each_GPU;
774 int gpu_id, rank;
776 /* Write the option flag */
777 strcpy(ptr, flag);
778 ptr += flag_length;
780 numPpRanks = numRanks - numPmeRanks;
781 max_num_ranks_for_each_GPU = numPpRanks / gpu_ids->n;
782 if (max_num_ranks_for_each_GPU * gpu_ids->n != numPpRanks)
784 /* Some GPUs will receive more work than others, which
785 * we choose to be those with the lowest indices */
786 max_num_ranks_for_each_GPU++;
789 /* Loop over all eligible GPU ids */
790 for (gpu_id = 0, rank = 0; gpu_id < gpu_ids->n; gpu_id++)
792 int rank_for_this_GPU;
793 /* Loop over all PP ranks for GPU with ID gpu_id, building the
794 assignment string. */
795 for (rank_for_this_GPU = 0;
796 rank_for_this_GPU < max_num_ranks_for_each_GPU && rank < numPpRanks;
797 rank++, rank_for_this_GPU++)
799 *ptr = '0' + gpu_ids->ids[gpu_id];
800 ptr++;
804 *ptr = '\0';
806 return command_line;
809 static void launch_simulation(
810 gmx_bool bLaunch, /* Should the simulation be launched? */
811 FILE *fp, /* General log file */
812 gmx_bool bThreads, /* whether to use threads */
813 char *cmd_mpirun, /* Command for mpirun */
814 char *cmd_np, /* Switch for -np or -ntmpi or empty */
815 char *cmd_mdrun, /* Command for mdrun */
816 char *args_for_mdrun, /* Arguments for mdrun */
817 const char *simulation_tpr, /* This tpr will be simulated */
818 int nnodes, /* Number of ranks to use */
819 int nPMEnodes, /* Number of PME ranks to use */
820 const t_eligible_gpu_ids *gpu_ids) /* Struct containing GPU IDs for
821 * constructing mdrun command lines */
823 char *command, *cmd_gpu_ids;
826 /* Make enough space for the system call command,
827 * (200 extra chars for -npme ... etc. options should suffice): */
828 snew(command, strlen(cmd_mpirun)+strlen(cmd_mdrun)+strlen(cmd_np)+strlen(args_for_mdrun)+strlen(simulation_tpr)+200);
830 cmd_gpu_ids = make_gpu_id_command_line(nnodes, nPMEnodes, gpu_ids);
832 /* Note that the -passall options requires args_for_mdrun to be at the end
833 * of the command line string */
834 if (bThreads)
836 sprintf(command, "%s%s-npme %d -s %s %s %s",
837 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids);
839 else
841 sprintf(command, "%s%s%s -npme %d -s %s %s %s",
842 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids);
845 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
846 sep_line(fp);
847 fflush(fp);
849 /* Now the real thing! */
850 if (bLaunch)
852 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
853 sep_line(stdout);
854 fflush(stdout);
855 gmx_system_call(command);
860 static void modify_PMEsettings(
861 gmx_int64_t simsteps, /* Set this value as number of time steps */
862 gmx_int64_t init_step, /* Set this value as init_step */
863 const char *fn_best_tpr, /* tpr file with the best performance */
864 const char *fn_sim_tpr) /* name of tpr file to be launched */
866 t_inputrec *ir;
867 t_state state;
868 gmx_mtop_t mtop;
869 char buf[200];
871 snew(ir, 1);
872 read_tpx_state(fn_best_tpr, ir, &state, NULL, &mtop);
874 /* Reset nsteps and init_step to the value of the input .tpr file */
875 ir->nsteps = simsteps;
876 ir->init_step = init_step;
878 /* Write the tpr file which will be launched */
879 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%"GMX_PRId64);
880 fprintf(stdout, buf, ir->nsteps);
881 fflush(stdout);
882 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
884 sfree(ir);
887 static gmx_bool can_scale_rvdw(int vdwtype)
889 return (evdwCUT == vdwtype ||
890 evdwPME == vdwtype);
893 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
895 /* Make additional TPR files with more computational load for the
896 * direct space processors: */
897 static void make_benchmark_tprs(
898 const char *fn_sim_tpr, /* READ : User-provided tpr file */
899 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
900 gmx_int64_t benchsteps, /* Number of time steps for benchmark runs */
901 gmx_int64_t statesteps, /* Step counter in checkpoint file */
902 real rmin, /* Minimal Coulomb radius */
903 real rmax, /* Maximal Coulomb radius */
904 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
905 int *ntprs, /* No. of TPRs to write, each with a different
906 rcoulomb and fourierspacing */
907 t_inputinfo *info, /* Contains information about mdp file options */
908 FILE *fp) /* Write the output here */
910 int i, j, d;
911 t_inputrec *ir;
912 t_state state;
913 gmx_mtop_t mtop;
914 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
915 char buf[200];
916 rvec box_size;
917 gmx_bool bNote = FALSE;
918 real add; /* Add this to rcoul for the next test */
919 real fac = 1.0; /* Scaling factor for Coulomb radius */
920 real fourierspacing; /* Basic fourierspacing from tpr */
923 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
924 *ntprs > 1 ? "s" : "", "%"GMX_PRId64, benchsteps > 1 ? "s" : "");
925 fprintf(stdout, buf, benchsteps);
926 if (statesteps > 0)
928 sprintf(buf, " (adding %s steps from checkpoint file)", "%"GMX_PRId64);
929 fprintf(stdout, buf, statesteps);
930 benchsteps += statesteps;
932 fprintf(stdout, ".\n");
935 snew(ir, 1);
936 read_tpx_state(fn_sim_tpr, ir, &state, NULL, &mtop);
938 /* Check if some kind of PME was chosen */
939 if (EEL_PME(ir->coulombtype) == FALSE)
941 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
942 EELTYPE(eelPME));
945 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
946 if ( (ir->cutoff_scheme != ecutsVERLET) &&
947 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
949 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
950 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
952 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
953 else if (ir->rcoulomb > ir->rlist)
955 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
956 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
959 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
961 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
962 bScaleRvdw = FALSE;
965 /* Reduce the number of steps for the benchmarks */
966 info->orig_sim_steps = ir->nsteps;
967 ir->nsteps = benchsteps;
968 /* We must not use init_step from the input tpr file for the benchmarks */
969 info->orig_init_step = ir->init_step;
970 ir->init_step = 0;
972 /* For PME-switch potentials, keep the radial distance of the buffer region */
973 nlist_buffer = ir->rlist - ir->rcoulomb;
975 /* Determine length of triclinic box vectors */
976 for (d = 0; d < DIM; d++)
978 box_size[d] = 0;
979 for (i = 0; i < DIM; i++)
981 box_size[d] += state.box[d][i]*state.box[d][i];
983 box_size[d] = sqrt(box_size[d]);
986 if (ir->fourier_spacing > 0)
988 info->fsx[0] = ir->fourier_spacing;
989 info->fsy[0] = ir->fourier_spacing;
990 info->fsz[0] = ir->fourier_spacing;
992 else
994 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
995 info->fsx[0] = box_size[XX]/ir->nkx;
996 info->fsy[0] = box_size[YY]/ir->nky;
997 info->fsz[0] = box_size[ZZ]/ir->nkz;
1000 /* If no value for the fourierspacing was provided on the command line, we
1001 * use the reconstruction from the tpr file */
1002 if (ir->fourier_spacing > 0)
1004 /* Use the spacing from the tpr */
1005 fourierspacing = ir->fourier_spacing;
1007 else
1009 /* Use the maximum observed spacing */
1010 fourierspacing = max(max(info->fsx[0], info->fsy[0]), info->fsz[0]);
1013 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
1015 /* For performance comparisons the number of particles is useful to have */
1016 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
1018 /* Print information about settings of which some are potentially modified: */
1019 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
1020 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
1021 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
1022 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
1023 if (ir_vdw_switched(ir))
1025 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
1027 if (EPME_SWITCHED(ir->coulombtype))
1029 fprintf(fp, " rlist : %f nm\n", ir->rlist);
1031 if (ir->rlistlong != max_cutoff(ir->rvdw, ir->rcoulomb))
1033 fprintf(fp, " rlistlong : %f nm\n", ir->rlistlong);
1036 /* Print a descriptive line about the tpr settings tested */
1037 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
1038 fprintf(fp, " No. scaling rcoulomb");
1039 fprintf(fp, " nkx nky nkz");
1040 fprintf(fp, " spacing");
1041 if (can_scale_rvdw(ir->vdwtype))
1043 fprintf(fp, " rvdw");
1045 if (EPME_SWITCHED(ir->coulombtype))
1047 fprintf(fp, " rlist");
1049 if (ir->rlistlong != max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb)) )
1051 fprintf(fp, " rlistlong");
1053 fprintf(fp, " tpr file\n");
1055 /* Loop to create the requested number of tpr input files */
1056 for (j = 0; j < *ntprs; j++)
1058 /* The first .tpr is the provided one, just need to modify nsteps,
1059 * so skip the following block */
1060 if (j != 0)
1062 /* Determine which Coulomb radii rc to use in the benchmarks */
1063 add = (rmax-rmin)/(*ntprs-1);
1064 if (gmx_within_tol(rmin, info->rcoulomb[0], GMX_REAL_EPS))
1066 ir->rcoulomb = rmin + j*add;
1068 else if (gmx_within_tol(rmax, info->rcoulomb[0], GMX_REAL_EPS))
1070 ir->rcoulomb = rmin + (j-1)*add;
1072 else
1074 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
1075 add = (rmax-rmin)/(*ntprs-2);
1076 ir->rcoulomb = rmin + (j-1)*add;
1079 /* Determine the scaling factor fac */
1080 fac = ir->rcoulomb/info->rcoulomb[0];
1082 /* Scale the Fourier grid spacing */
1083 ir->nkx = ir->nky = ir->nkz = 0;
1084 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
1086 /* Adjust other radii since various conditions need to be fulfilled */
1087 if (eelPME == ir->coulombtype)
1089 /* plain PME, rcoulomb must be equal to rlist */
1090 ir->rlist = ir->rcoulomb;
1092 else
1094 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1095 ir->rlist = ir->rcoulomb + nlist_buffer;
1098 if (bScaleRvdw && can_scale_rvdw(ir->vdwtype))
1100 if (ecutsVERLET == ir->cutoff_scheme ||
1101 evdwPME == ir->vdwtype)
1103 /* With either the Verlet cutoff-scheme or LJ-PME,
1104 the van der Waals radius must always equal the
1105 Coulomb radius */
1106 ir->rvdw = ir->rcoulomb;
1108 else
1110 /* For vdw cutoff, rvdw >= rlist */
1111 ir->rvdw = max(info->rvdw[0], ir->rlist);
1115 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1117 } /* end of "if (j != 0)" */
1119 /* for j==0: Save the original settings
1120 * for j >0: Save modified radii and Fourier grids */
1121 info->rcoulomb[j] = ir->rcoulomb;
1122 info->rvdw[j] = ir->rvdw;
1123 info->nkx[j] = ir->nkx;
1124 info->nky[j] = ir->nky;
1125 info->nkz[j] = ir->nkz;
1126 info->rlist[j] = ir->rlist;
1127 info->rlistlong[j] = ir->rlistlong;
1128 info->fsx[j] = fac*fourierspacing;
1129 info->fsy[j] = fac*fourierspacing;
1130 info->fsz[j] = fac*fourierspacing;
1132 /* Write the benchmark tpr file */
1133 strncpy(fn_bench_tprs[j], fn_sim_tpr, strlen(fn_sim_tpr)-strlen(".tpr"));
1134 sprintf(buf, "_bench%.2d.tpr", j);
1135 strcat(fn_bench_tprs[j], buf);
1136 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1137 fprintf(stdout, "%"GMX_PRId64, ir->nsteps);
1138 if (j > 0)
1140 fprintf(stdout, ", scaling factor %f\n", fac);
1142 else
1144 fprintf(stdout, ", unmodified settings\n");
1147 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1149 /* Write information about modified tpr settings to log file */
1150 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1151 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1152 fprintf(fp, " %9f ", info->fsx[j]);
1153 if (can_scale_rvdw(ir->vdwtype))
1155 fprintf(fp, "%10f", ir->rvdw);
1157 if (EPME_SWITCHED(ir->coulombtype))
1159 fprintf(fp, "%10f", ir->rlist);
1161 if (info->rlistlong[0] != max_cutoff(info->rlist[0], max_cutoff(info->rvdw[0], info->rcoulomb[0])) )
1163 fprintf(fp, "%10f", ir->rlistlong);
1165 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1167 /* Make it clear to the user that some additional settings were modified */
1168 if (!gmx_within_tol(ir->rvdw, info->rvdw[0], GMX_REAL_EPS)
1169 || !gmx_within_tol(ir->rlistlong, info->rlistlong[0], GMX_REAL_EPS) )
1171 bNote = TRUE;
1174 if (bNote)
1176 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1177 "other input settings were also changed (see table above).\n"
1178 "Please check if the modified settings are appropriate.\n");
1180 fflush(stdout);
1181 fflush(fp);
1182 sfree(ir);
1186 /* Rename the files we want to keep to some meaningful filename and
1187 * delete the rest */
1188 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1189 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1191 char numstring[STRLEN];
1192 char newfilename[STRLEN];
1193 const char *fn = NULL;
1194 int i;
1195 const char *opt;
1198 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1200 for (i = 0; i < nfile; i++)
1202 opt = (char *)fnm[i].opt;
1203 if (strcmp(opt, "-p") == 0)
1205 /* do nothing; keep this file */
1208 else if (strcmp(opt, "-bg") == 0)
1210 /* Give the log file a nice name so one can later see which parameters were used */
1211 numstring[0] = '\0';
1212 if (nr > 0)
1214 sprintf(numstring, "_%d", nr);
1216 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1217 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1219 fprintf(stdout, "renaming log file to %s\n", newfilename);
1220 make_backup(newfilename);
1221 rename(opt2fn("-bg", nfile, fnm), newfilename);
1224 else if (strcmp(opt, "-err") == 0)
1226 /* This file contains the output of stderr. We want to keep it in
1227 * cases where there have been problems. */
1228 fn = opt2fn(opt, nfile, fnm);
1229 numstring[0] = '\0';
1230 if (nr > 0)
1232 sprintf(numstring, "_%d", nr);
1234 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1235 if (gmx_fexist(fn))
1237 if (bKeepStderr)
1239 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1240 make_backup(newfilename);
1241 rename(fn, newfilename);
1243 else
1245 fprintf(stdout, "Deleting %s\n", fn);
1246 remove(fn);
1250 /* Delete the files which are created for each benchmark run: (options -b*) */
1251 else if ( (0 == strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1253 remove_if_exists(opt2fn(opt, nfile, fnm));
1259 enum {
1260 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1263 /* Create a list of numbers of PME nodes to test */
1264 static void make_npme_list(
1265 const char *npmevalues_opt, /* Make a complete list with all
1266 * possibilities or a short list that keeps only
1267 * reasonable numbers of PME nodes */
1268 int *nentries, /* Number of entries we put in the nPMEnodes list */
1269 int *nPMEnodes[], /* Each entry contains the value for -npme */
1270 int nnodes, /* Total number of nodes to do the tests on */
1271 int minPMEnodes, /* Minimum number of PME nodes */
1272 int maxPMEnodes) /* Maximum number of PME nodes */
1274 int i, npme, npp;
1275 int min_factor = 1; /* We request that npp and npme have this minimal
1276 * largest common factor (depends on npp) */
1277 int nlistmax; /* Max. list size */
1278 int nlist; /* Actual number of entries in list */
1279 int eNPME = 0;
1282 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1283 if (0 == strcmp(npmevalues_opt, "all") )
1285 eNPME = eNpmeAll;
1287 else if (0 == strcmp(npmevalues_opt, "subset") )
1289 eNPME = eNpmeSubset;
1291 else /* "auto" or "range" */
1293 if (nnodes <= 64)
1295 eNPME = eNpmeAll;
1297 else if (nnodes < 128)
1299 eNPME = eNpmeReduced;
1301 else
1303 eNPME = eNpmeSubset;
1307 /* Calculate how many entries we could possibly have (in case of -npme all) */
1308 if (nnodes > 2)
1310 nlistmax = maxPMEnodes - minPMEnodes + 3;
1311 if (0 == minPMEnodes)
1313 nlistmax--;
1316 else
1318 nlistmax = 1;
1321 /* Now make the actual list which is at most of size nlist */
1322 snew(*nPMEnodes, nlistmax);
1323 nlist = 0; /* start counting again, now the real entries in the list */
1324 for (i = 0; i < nlistmax - 2; i++)
1326 npme = maxPMEnodes - i;
1327 npp = nnodes-npme;
1328 switch (eNPME)
1330 case eNpmeAll:
1331 min_factor = 1;
1332 break;
1333 case eNpmeReduced:
1334 min_factor = 2;
1335 break;
1336 case eNpmeSubset:
1337 /* For 2d PME we want a common largest factor of at least the cube
1338 * root of the number of PP nodes */
1339 min_factor = (int) pow(npp, 1.0/3.0);
1340 break;
1341 default:
1342 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1343 break;
1345 if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1347 (*nPMEnodes)[nlist] = npme;
1348 nlist++;
1351 /* We always test 0 PME nodes and the automatic number */
1352 *nentries = nlist + 2;
1353 (*nPMEnodes)[nlist ] = 0;
1354 (*nPMEnodes)[nlist+1] = -1;
1356 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1357 for (i = 0; i < *nentries-1; i++)
1359 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1361 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1365 /* Allocate memory to store the performance data */
1366 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1368 int i, j, k;
1371 for (k = 0; k < ntprs; k++)
1373 snew(perfdata[k], datasets);
1374 for (i = 0; i < datasets; i++)
1376 for (j = 0; j < repeats; j++)
1378 snew(perfdata[k][i].Gcycles, repeats);
1379 snew(perfdata[k][i].ns_per_day, repeats);
1380 snew(perfdata[k][i].PME_f_load, repeats);
1387 /* Check for errors on mdrun -h */
1388 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp,
1389 const t_filenm *fnm, int nfile)
1391 const char *fn = NULL;
1392 char *command, *msg;
1393 int ret;
1396 snew(command, length + 15);
1397 snew(msg, length + 500);
1399 fprintf(stdout, "Making sure the benchmarks can be executed by running just 1 step...\n");
1400 sprintf(command, "%s -nsteps 1 -quiet", mdrun_cmd_line);
1401 fprintf(stdout, "Executing '%s' ...\n", command);
1402 ret = gmx_system_call(command);
1404 if (0 != ret)
1406 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1407 * get the error message from mdrun itself */
1408 sprintf(msg,
1409 "Cannot run the first benchmark simulation! Please check the error message of\n"
1410 "mdrun for the source of the problem. Did you provide a command line\n"
1411 "argument that neither gmx tune_pme nor mdrun understands? If you're\n"
1412 "sure your command line should work, you can bypass this check with \n"
1413 "gmx tune_pme -nocheck. The failing command was:\n"
1414 "\n%s\n\n", command);
1416 fprintf(stderr, "%s", msg);
1417 sep_line(fp);
1418 fprintf(fp, "%s", msg);
1420 exit(ret);
1422 fprintf(stdout, "Benchmarks can be executed!\n");
1424 /* Clean up the benchmark output files we just created */
1425 fprintf(stdout, "Cleaning up ...\n");
1426 remove_if_exists(opt2fn("-bc", nfile, fnm));
1427 remove_if_exists(opt2fn("-be", nfile, fnm));
1428 remove_if_exists(opt2fn("-bcpo", nfile, fnm));
1429 remove_if_exists(opt2fn("-bg", nfile, fnm));
1430 remove_if_exists(opt2fn("-bo", nfile, fnm));
1431 remove_if_exists(opt2fn("-bx", nfile, fnm));
1433 sfree(command);
1434 sfree(msg );
1437 static void do_the_tests(
1438 FILE *fp, /* General g_tune_pme output file */
1439 char **tpr_names, /* Filenames of the input files to test */
1440 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1441 int minPMEnodes, /* Min fraction of nodes to use for PME */
1442 int npme_fixed, /* If >= -1, test fixed number of PME
1443 * nodes only */
1444 const char *npmevalues_opt, /* Which -npme values should be tested */
1445 t_perf **perfdata, /* Here the performace data is stored */
1446 int *pmeentries, /* Entries in the nPMEnodes list */
1447 int repeats, /* Repeat each test this often */
1448 int nnodes, /* Total number of nodes = nPP + nPME */
1449 int nr_tprs, /* Total number of tpr files to test */
1450 gmx_bool bThreads, /* Threads or MPI? */
1451 char *cmd_mpirun, /* mpirun command string */
1452 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1453 char *cmd_mdrun, /* mdrun command string */
1454 char *cmd_args_bench, /* arguments for mdrun in a string */
1455 const t_filenm *fnm, /* List of filenames from command line */
1456 int nfile, /* Number of files specified on the cmdl. */
1457 int presteps, /* DLB equilibration steps, is checked */
1458 gmx_int64_t cpt_steps, /* Time step counter in the checkpoint */
1459 gmx_bool bCheck, /* Check whether benchmark mdrun works */
1460 const t_eligible_gpu_ids *gpu_ids) /* Struct containing GPU IDs for
1461 * constructing mdrun command lines */
1463 int i, nr, k, ret, count = 0, totaltests;
1464 int *nPMEnodes = NULL;
1465 t_perf *pd = NULL;
1466 int cmdline_length;
1467 char *command, *cmd_stub;
1468 char buf[STRLEN];
1469 gmx_bool bResetProblem = FALSE;
1470 gmx_bool bFirst = TRUE;
1471 gmx_bool bUsingGpus = 0 < gpu_ids->n;
1473 /* This string array corresponds to the eParselog enum type at the start
1474 * of this file */
1475 const char* ParseLog[] = {
1476 "OK.",
1477 "Logfile not found!",
1478 "No timings, logfile truncated?",
1479 "Run was terminated.",
1480 "Counters were not reset properly.",
1481 "No DD grid found for these settings.",
1482 "TPX version conflict!",
1483 "mdrun was not started in parallel!",
1484 "Number of PP ranks has a prime factor that is too large.",
1485 "The number of PP ranks did not suit the number of GPUs.",
1486 "Some GPUs were not detected or are incompatible.",
1487 "An error occurred."
1489 char str_PME_f_load[13];
1492 /* Allocate space for the mdrun command line. 100 extra characters should
1493 be more than enough for the -npme etcetera arguments */
1494 cmdline_length = strlen(cmd_mpirun)
1495 + strlen(cmd_np)
1496 + strlen(cmd_mdrun)
1497 + strlen(cmd_args_bench)
1498 + strlen(tpr_names[0]) + 100;
1499 snew(command, cmdline_length);
1500 snew(cmd_stub, cmdline_length);
1502 /* Construct the part of the command line that stays the same for all tests: */
1503 if (bThreads)
1505 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1507 else
1509 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1512 /* Create a list of numbers of PME nodes to test */
1513 if (npme_fixed < -1)
1515 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1516 nnodes, minPMEnodes, maxPMEnodes);
1518 else
1520 *pmeentries = 1;
1521 snew(nPMEnodes, 1);
1522 nPMEnodes[0] = npme_fixed;
1523 fprintf(stderr, "Will use a fixed number of %d PME-only ranks.\n", nPMEnodes[0]);
1526 if (0 == repeats)
1528 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1529 gmx_ffclose(fp);
1530 finalize(opt2fn("-p", nfile, fnm));
1531 exit(0);
1534 /* Allocate one dataset for each tpr input file: */
1535 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1537 /*****************************************/
1538 /* Main loop over all tpr files to test: */
1539 /*****************************************/
1540 totaltests = nr_tprs*(*pmeentries)*repeats;
1541 for (k = 0; k < nr_tprs; k++)
1543 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1544 fprintf(fp, "PME ranks Gcycles ns/day PME/f Remark\n");
1545 /* Loop over various numbers of PME nodes: */
1546 for (i = 0; i < *pmeentries; i++)
1548 char *cmd_gpu_ids = NULL;
1550 pd = &perfdata[k][i];
1552 cmd_gpu_ids = make_gpu_id_command_line(nnodes, nPMEnodes[i], gpu_ids);
1554 /* Loop over the repeats for each scenario: */
1555 for (nr = 0; nr < repeats; nr++)
1557 pd->nPMEnodes = nPMEnodes[i];
1559 /* Add -npme and -s to the command line and save it. Note that
1560 * the -passall (if set) options requires cmd_args_bench to be
1561 * at the end of the command line string */
1562 snew(pd->mdrun_cmd_line, cmdline_length);
1563 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s %s",
1564 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench, cmd_gpu_ids);
1566 /* To prevent that all benchmarks fail due to a show-stopper argument
1567 * on the mdrun command line, we make a quick check first.
1568 * This check can be turned off in cases where the automatically chosen
1569 * number of PME-only ranks leads to a number of PP ranks for which no
1570 * decomposition can be found (e.g. for large prime numbers) */
1571 if (bFirst && bCheck)
1573 /* TODO This check is really for a functional
1574 * .tpr, and if we need it, it should take place
1575 * for every .tpr, and the logic for it should be
1576 * immediately inside the loop over k, not in
1577 * this inner loop. */
1578 char *temporary_cmd_line;
1580 snew(temporary_cmd_line, cmdline_length);
1581 /* TODO -npme 0 is more likely to succeed at low
1582 parallelism than the default of -npme -1, but
1583 is more likely to fail at the scaling limit
1584 when the PP domains may be too small. "mpirun
1585 -np 1 mdrun" is probably a reasonable thing to
1586 do for this check, but it'll be easier to
1587 implement that after some refactoring of how
1588 the number of MPI ranks is managed. */
1589 sprintf(temporary_cmd_line, "%s-npme 0 -nb cpu -s %s %s",
1590 cmd_stub, tpr_names[k], cmd_args_bench);
1591 make_sure_it_runs(temporary_cmd_line, cmdline_length, fp, fnm, nfile);
1593 bFirst = FALSE;
1595 /* Do a benchmark simulation: */
1596 if (repeats > 1)
1598 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1600 else
1602 buf[0] = '\0';
1604 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1605 (100.0*count)/totaltests,
1606 k+1, nr_tprs, i+1, *pmeentries, buf);
1607 make_backup(opt2fn("-err", nfile, fnm));
1608 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1609 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1610 gmx_system_call(command);
1612 /* Collect the performance data from the log file; also check stderr
1613 * for fatal errors */
1614 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1615 pd, nr, presteps, cpt_steps, nnodes);
1616 if ((presteps > 0) && (ret == eParselogResetProblem))
1618 bResetProblem = TRUE;
1621 if (-1 == pd->nPMEnodes)
1623 sprintf(buf, "(%3d)", pd->guessPME);
1625 else
1627 sprintf(buf, " ");
1630 /* Nicer output */
1631 if (pd->PME_f_load[nr] > 0.0)
1633 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1635 else
1637 sprintf(str_PME_f_load, "%s", " - ");
1640 /* Write the data we got to disk */
1641 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1642 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1643 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1645 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1647 fprintf(fp, "\n");
1648 fflush(fp);
1649 count++;
1651 /* Do some cleaning up and delete the files we do not need any more */
1652 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1654 /* If the first run with this number of processors already failed, do not try again: */
1655 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1657 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1658 count += repeats-(nr+1);
1659 break;
1661 } /* end of repeats loop */
1662 sfree(cmd_gpu_ids);
1663 } /* end of -npme loop */
1664 } /* end of tpr file loop */
1666 if (bResetProblem)
1668 sep_line(fp);
1669 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1670 sep_line(fp);
1672 sfree(command);
1673 sfree(cmd_stub);
1677 static void check_input(
1678 int nnodes,
1679 int repeats,
1680 int *ntprs,
1681 real *rmin,
1682 real rcoulomb,
1683 real *rmax,
1684 real maxPMEfraction,
1685 real minPMEfraction,
1686 int npme_fixed,
1687 gmx_int64_t bench_nsteps,
1688 const t_filenm *fnm,
1689 int nfile,
1690 int sim_part,
1691 int presteps,
1692 int npargs,
1693 t_pargs *pa)
1695 int old;
1698 /* Make sure the input file exists */
1699 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1701 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1704 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1705 if ( (0 == strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1707 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1708 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1711 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1712 if (repeats < 0)
1714 gmx_fatal(FARGS, "Number of repeats < 0!");
1717 /* Check number of nodes */
1718 if (nnodes < 1)
1720 gmx_fatal(FARGS, "Number of ranks/threads must be a positive integer.");
1723 /* Automatically choose -ntpr if not set */
1724 if (*ntprs < 1)
1726 if (nnodes < 16)
1728 *ntprs = 1;
1730 else
1732 *ntprs = 3;
1733 /* Set a reasonable scaling factor for rcoulomb */
1734 if (*rmax <= 0)
1736 *rmax = rcoulomb * 1.2;
1739 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1741 else
1743 if (1 == *ntprs)
1745 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1749 /* Make shure that rmin <= rcoulomb <= rmax */
1750 if (*rmin <= 0)
1752 *rmin = rcoulomb;
1754 if (*rmax <= 0)
1756 *rmax = rcoulomb;
1758 if (!(*rmin <= *rmax) )
1760 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1761 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1763 /* Add test scenarios if rmin or rmax were set */
1764 if (*ntprs <= 2)
1766 if (!gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1768 (*ntprs)++;
1769 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1770 *rmin, *ntprs);
1772 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1774 (*ntprs)++;
1775 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1776 *rmax, *ntprs);
1779 old = *ntprs;
1780 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1781 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) || !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1783 *ntprs = max(*ntprs, 2);
1786 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1787 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1789 *ntprs = max(*ntprs, 3);
1792 if (old != *ntprs)
1794 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1797 if (*ntprs > 1)
1799 if (gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && gmx_within_tol(rcoulomb, *rmax, GMX_REAL_EPS)) /* We have just a single rc */
1801 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1802 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1803 "with correspondingly adjusted PME grid settings\n");
1804 *ntprs = 1;
1808 /* Check whether max and min fraction are within required values */
1809 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1811 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1813 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1815 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1817 if (maxPMEfraction < minPMEfraction)
1819 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1822 /* Check whether the number of steps - if it was set - has a reasonable value */
1823 if (bench_nsteps < 0)
1825 gmx_fatal(FARGS, "Number of steps must be positive.");
1828 if (bench_nsteps > 10000 || bench_nsteps < 100)
1830 fprintf(stderr, "WARNING: steps=");
1831 fprintf(stderr, "%"GMX_PRId64, bench_nsteps);
1832 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1835 if (presteps < 0)
1837 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1840 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1841 if (*ntprs > 1)
1843 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1845 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1849 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1850 * only. We need to check whether the requested number of PME-only nodes
1851 * makes sense. */
1852 if (npme_fixed > -1)
1854 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1855 if (2*npme_fixed > nnodes)
1857 gmx_fatal(FARGS, "Cannot have more than %d PME-only ranks for a total of %d ranks (you chose %d).\n",
1858 nnodes/2, nnodes, npme_fixed);
1860 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1862 fprintf(stderr, "WARNING: Only %g percent of the ranks are assigned as PME-only ranks.\n",
1863 100.0*((real)npme_fixed / (real)nnodes));
1865 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1867 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1868 " fixed number of PME-only ranks is requested with -fix.\n");
1874 /* Returns TRUE when "opt" is needed at launch time */
1875 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1877 if (0 == strncmp(opt, "-swap", 5))
1879 return bSet;
1882 /* Apart from the input .tpr and the output log files we need all options that
1883 * were set on the command line and that do not start with -b */
1884 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2)
1885 || 0 == strncmp(opt, "-err", 4) || 0 == strncmp(opt, "-p", 2) )
1887 return FALSE;
1890 return bSet;
1894 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1895 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1897 /* Apart from the input .tpr, all files starting with "-b" are for
1898 * _b_enchmark files exclusively */
1899 if (0 == strncmp(opt, "-s", 2))
1901 return FALSE;
1904 if (0 == strncmp(opt, "-b", 2) || 0 == strncmp(opt, "-s", 2))
1906 if (!bOptional || bSet)
1908 return TRUE;
1910 else
1912 return FALSE;
1915 else
1917 if (bIsOutput)
1919 return FALSE;
1921 else
1923 if (bSet) /* These are additional input files like -cpi -ei */
1925 return TRUE;
1927 else
1929 return FALSE;
1936 /* Adds 'buf' to 'str' */
1937 static void add_to_string(char **str, char *buf)
1939 int len;
1942 len = strlen(*str) + strlen(buf) + 1;
1943 srenew(*str, len);
1944 strcat(*str, buf);
1948 /* Create the command line for the benchmark as well as for the real run */
1949 static void create_command_line_snippets(
1950 gmx_bool bAppendFiles,
1951 gmx_bool bKeepAndNumCPT,
1952 gmx_bool bResetHWay,
1953 int presteps,
1954 int nfile,
1955 t_filenm fnm[],
1956 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1957 char *cmd_args_launch[], /* command line arguments for simulation run */
1958 char extra_args[], /* Add this to the end of the command line */
1959 char *deffnm) /* Default file names, or NULL if not set */
1961 int i;
1962 char *opt;
1963 const char *name;
1964 char strbuf[STRLEN];
1967 /* strlen needs at least '\0' as a string: */
1968 snew(*cmd_args_bench, 1);
1969 snew(*cmd_args_launch, 1);
1970 *cmd_args_launch[0] = '\0';
1971 *cmd_args_bench[0] = '\0';
1974 /*******************************************/
1975 /* 1. Process other command line arguments */
1976 /*******************************************/
1977 if (presteps > 0)
1979 /* Add equilibration steps to benchmark options */
1980 sprintf(strbuf, "-resetstep %d ", presteps);
1981 add_to_string(cmd_args_bench, strbuf);
1983 /* These switches take effect only at launch time */
1984 if (deffnm)
1986 sprintf(strbuf, "-deffnm %s ", deffnm);
1987 add_to_string(cmd_args_launch, strbuf);
1989 if (FALSE == bAppendFiles)
1991 add_to_string(cmd_args_launch, "-noappend ");
1993 if (bKeepAndNumCPT)
1995 add_to_string(cmd_args_launch, "-cpnum ");
1997 if (bResetHWay)
1999 add_to_string(cmd_args_launch, "-resethway ");
2002 /********************/
2003 /* 2. Process files */
2004 /********************/
2005 for (i = 0; i < nfile; i++)
2007 opt = (char *)fnm[i].opt;
2008 name = opt2fn(opt, nfile, fnm);
2010 /* Strbuf contains the options, now let's sort out where we need that */
2011 sprintf(strbuf, "%s %s ", opt, name);
2013 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
2015 /* All options starting with -b* need the 'b' removed,
2016 * therefore overwrite strbuf */
2017 if (0 == strncmp(opt, "-b", 2))
2019 sprintf(strbuf, "-%s %s ", &opt[2], name);
2022 add_to_string(cmd_args_bench, strbuf);
2025 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
2027 add_to_string(cmd_args_launch, strbuf);
2031 add_to_string(cmd_args_bench, extra_args);
2032 add_to_string(cmd_args_launch, extra_args);
2036 /* Set option opt */
2037 static void setopt(const char *opt, int nfile, t_filenm fnm[])
2039 int i;
2041 for (i = 0; (i < nfile); i++)
2043 if (strcmp(opt, fnm[i].opt) == 0)
2045 fnm[i].flag |= ffSET;
2051 /* This routine inspects the tpr file and ...
2052 * 1. checks for output files that get triggered by a tpr option. These output
2053 * files are marked as 'set' to allow for a proper cleanup after each
2054 * tuning run.
2055 * 2. returns the PME:PP load ratio
2056 * 3. returns rcoulomb from the tpr */
2057 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
2059 gmx_bool bTpi; /* Is test particle insertion requested? */
2060 gmx_bool bFree; /* Is a free energy simulation requested? */
2061 gmx_bool bNM; /* Is a normal mode analysis requested? */
2062 gmx_bool bSwap; /* Is water/ion position swapping requested? */
2063 t_inputrec ir;
2064 t_state state;
2065 gmx_mtop_t mtop;
2068 /* Check tpr file for options that trigger extra output files */
2069 read_tpx_state(opt2fn("-s", nfile, fnm), &ir, &state, NULL, &mtop);
2070 bFree = (efepNO != ir.efep );
2071 bNM = (eiNM == ir.eI );
2072 bSwap = (eswapNO != ir.eSwapCoords);
2073 bTpi = EI_TPI(ir.eI);
2075 /* Set these output files on the tuning command-line */
2076 if (ir.bPull)
2078 setopt("-pf", nfile, fnm);
2079 setopt("-px", nfile, fnm);
2081 if (bFree)
2083 setopt("-dhdl", nfile, fnm);
2085 if (bTpi)
2087 setopt("-tpi", nfile, fnm);
2088 setopt("-tpid", nfile, fnm);
2090 if (bNM)
2092 setopt("-mtx", nfile, fnm);
2094 if (bSwap)
2096 setopt("-swap", nfile, fnm);
2099 *rcoulomb = ir.rcoulomb;
2101 /* Return the estimate for the number of PME nodes */
2102 return pme_load_estimate(&mtop, &ir, state.box);
2106 static void couple_files_options(int nfile, t_filenm fnm[])
2108 int i;
2109 gmx_bool bSet, bBench;
2110 char *opt;
2111 char buf[20];
2114 for (i = 0; i < nfile; i++)
2116 opt = (char *)fnm[i].opt;
2117 bSet = ((fnm[i].flag & ffSET) != 0);
2118 bBench = (0 == strncmp(opt, "-b", 2));
2120 /* Check optional files */
2121 /* If e.g. -eo is set, then -beo also needs to be set */
2122 if (is_optional(&fnm[i]) && bSet && !bBench)
2124 sprintf(buf, "-b%s", &opt[1]);
2125 setopt(buf, nfile, fnm);
2127 /* If -beo is set, then -eo also needs to be! */
2128 if (is_optional(&fnm[i]) && bSet && bBench)
2130 sprintf(buf, "-%s", &opt[2]);
2131 setopt(buf, nfile, fnm);
2137 #define BENCHSTEPS (1000)
2139 int gmx_tune_pme(int argc, char *argv[])
2141 const char *desc[] = {
2142 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of ranks, [THISMODULE] systematically",
2143 "times [gmx-mdrun] with various numbers of PME-only ranks and determines",
2144 "which setting is fastest. It will also test whether performance can",
2145 "be enhanced by shifting load from the reciprocal to the real space",
2146 "part of the Ewald sum. ",
2147 "Simply pass your [REF].tpr[ref] file to [THISMODULE] together with other options",
2148 "for [gmx-mdrun] as needed.[PAR]",
2149 "[THISMODULE] needs to call [gmx-mdrun] and so requires that you",
2150 "specify how to call mdrun with the argument to the [TT]-mdrun[tt]",
2151 "parameter. Depending how you have built GROMACS, values such as",
2152 "'gmx mdrun', 'gmx_d mdrun', or 'mdrun_mpi' might be needed.[PAR]",
2153 "The program that runs MPI programs can be set in the environment variable",
2154 "MPIRUN (defaults to 'mpirun'). Note that for certain MPI frameworks,",
2155 "you need to provide a machine- or hostfile. This can also be passed",
2156 "via the MPIRUN variable, e.g.[PAR]",
2157 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt]",
2158 "Note that in such cases it is normally necessary to compile",
2159 "and/or run [THISMODULE] without MPI support, so that it can call",
2160 "the MPIRUN program.[PAR]",
2161 "Before doing the actual benchmark runs, [THISMODULE] will do a quick",
2162 "check whether [gmx-mdrun] works as expected with the provided parallel settings",
2163 "if the [TT]-check[tt] option is activated (the default).",
2164 "Please call [THISMODULE] with the normal options you would pass to",
2165 "[gmx-mdrun] and add [TT]-np[tt] for the number of ranks to perform the",
2166 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2167 "to repeat each test several times to get better statistics. [PAR]",
2168 "[THISMODULE] can test various real space / reciprocal space workloads",
2169 "for you. With [TT]-ntpr[tt] you control how many extra [REF].tpr[ref] files will be",
2170 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2171 "Typically, the first test (number 0) will be with the settings from the input",
2172 "[REF].tpr[ref] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2173 "specified by [TT]-rmax[tt] with a somewhat smaller PME grid at the same time. ",
2174 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2175 "The remaining [REF].tpr[ref] files will have equally-spaced Coulomb radii (and Fourier "
2176 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2177 "if you just seek the optimal number of PME-only ranks; in that case",
2178 "your input [REF].tpr[ref] file will remain unchanged.[PAR]",
2179 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2180 "MD systems. The dynamic load balancing needs about 100 time steps",
2181 "to adapt to local load imbalances, therefore the time step counters",
2182 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2183 "for a higher accuracy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2184 "From the 'DD' load imbalance entries in the md.log output file you",
2185 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2186 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2187 "After calling [gmx-mdrun] several times, detailed performance information",
2188 "is available in the output file [TT]perf.out[tt].",
2189 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2190 "(options [TT]-b*[tt]), these will be automatically deleted after each test.[PAR]",
2191 "If you want the simulation to be started automatically with the",
2192 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2193 "Basic support for GPU-enabled [TT]mdrun[tt] exists. Give a string containing the IDs",
2194 "of the GPUs that you wish to use in the optimization in the [TT]-gpu_id[tt]",
2195 "command-line argument. Unlike [TT]mdrun -gpu_id[tt], this does not imply a mapping",
2196 "but merely the eligible set. [TT]g_tune_pme[tt] will construct calls to",
2197 "mdrun that use this set appropriately, assuming that PP ranks with low indices",
2198 "should map to GPUs with low indices, and increasing both monotonically",
2199 "over the respective sets.[PAR]",
2202 int nnodes = 1;
2203 int repeats = 2;
2204 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2205 real maxPMEfraction = 0.50;
2206 real minPMEfraction = 0.25;
2207 int maxPMEnodes, minPMEnodes;
2208 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2209 float guessPMEnodes;
2210 int npme_fixed = -2; /* If >= -1, use only this number
2211 * of PME-only nodes */
2212 int ntprs = 0;
2213 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2214 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2215 gmx_bool bScaleRvdw = TRUE;
2216 gmx_int64_t bench_nsteps = BENCHSTEPS;
2217 gmx_int64_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2218 gmx_int64_t cpt_steps = 0; /* Step counter in .cpt input file */
2219 int presteps = 100; /* Do a full cycle reset after presteps steps */
2220 gmx_bool bOverwrite = FALSE, bKeepTPR;
2221 gmx_bool bLaunch = FALSE;
2222 char *ExtraArgs = NULL;
2223 char **tpr_names = NULL;
2224 const char *simulation_tpr = NULL;
2225 char *deffnm = NULL;
2226 int best_npme, best_tpr;
2227 int sim_part = 1; /* For benchmarks with checkpoint files */
2228 char bbuf[STRLEN];
2231 /* Default program names if nothing else is found */
2232 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2233 char *cmd_args_bench, *cmd_args_launch;
2234 char *cmd_np = NULL;
2236 /* IDs of GPUs that are eligible for computation */
2237 char *eligible_gpu_ids = NULL;
2238 t_eligible_gpu_ids *gpu_ids = NULL;
2240 t_perf **perfdata = NULL;
2241 t_inputinfo *info;
2242 int i;
2243 FILE *fp;
2244 t_commrec *cr;
2246 /* Print out how long the tuning took */
2247 double seconds;
2249 static t_filenm fnm[] = {
2250 /* g_tune_pme */
2251 { efOUT, "-p", "perf", ffWRITE },
2252 { efLOG, "-err", "bencherr", ffWRITE },
2253 { efTPR, "-so", "tuned", ffWRITE },
2254 /* mdrun: */
2255 { efTPR, NULL, NULL, ffREAD },
2256 { efTRN, "-o", NULL, ffWRITE },
2257 { efCOMPRESSED, "-x", NULL, ffOPTWR },
2258 { efCPT, "-cpi", NULL, ffOPTRD },
2259 { efCPT, "-cpo", NULL, ffOPTWR },
2260 { efSTO, "-c", "confout", ffWRITE },
2261 { efEDR, "-e", "ener", ffWRITE },
2262 { efLOG, "-g", "md", ffWRITE },
2263 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2264 { efXVG, "-field", "field", ffOPTWR },
2265 { efXVG, "-table", "table", ffOPTRD },
2266 { efXVG, "-tabletf", "tabletf", ffOPTRD },
2267 { efXVG, "-tablep", "tablep", ffOPTRD },
2268 { efXVG, "-tableb", "table", ffOPTRD },
2269 { efTRX, "-rerun", "rerun", ffOPTRD },
2270 { efXVG, "-tpi", "tpi", ffOPTWR },
2271 { efXVG, "-tpid", "tpidist", ffOPTWR },
2272 { efEDI, "-ei", "sam", ffOPTRD },
2273 { efXVG, "-eo", "edsam", ffOPTWR },
2274 { efXVG, "-devout", "deviatie", ffOPTWR },
2275 { efXVG, "-runav", "runaver", ffOPTWR },
2276 { efXVG, "-px", "pullx", ffOPTWR },
2277 { efXVG, "-pf", "pullf", ffOPTWR },
2278 { efXVG, "-ro", "rotation", ffOPTWR },
2279 { efLOG, "-ra", "rotangles", ffOPTWR },
2280 { efLOG, "-rs", "rotslabs", ffOPTWR },
2281 { efLOG, "-rt", "rottorque", ffOPTWR },
2282 { efMTX, "-mtx", "nm", ffOPTWR },
2283 { efNDX, "-dn", "dipole", ffOPTWR },
2284 { efXVG, "-swap", "swapions", ffOPTWR },
2285 /* Output files that are deleted after each benchmark run */
2286 { efTRN, "-bo", "bench", ffWRITE },
2287 { efXTC, "-bx", "bench", ffWRITE },
2288 { efCPT, "-bcpo", "bench", ffWRITE },
2289 { efSTO, "-bc", "bench", ffWRITE },
2290 { efEDR, "-be", "bench", ffWRITE },
2291 { efLOG, "-bg", "bench", ffWRITE },
2292 { efXVG, "-beo", "benchedo", ffOPTWR },
2293 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2294 { efXVG, "-bfield", "benchfld", ffOPTWR },
2295 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2296 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2297 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2298 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2299 { efXVG, "-bpx", "benchpx", ffOPTWR },
2300 { efXVG, "-bpf", "benchpf", ffOPTWR },
2301 { efXVG, "-bro", "benchrot", ffOPTWR },
2302 { efLOG, "-bra", "benchrota", ffOPTWR },
2303 { efLOG, "-brs", "benchrots", ffOPTWR },
2304 { efLOG, "-brt", "benchrott", ffOPTWR },
2305 { efMTX, "-bmtx", "benchn", ffOPTWR },
2306 { efNDX, "-bdn", "bench", ffOPTWR },
2307 { efXVG, "-bswap", "benchswp", ffOPTWR }
2310 gmx_bool bThreads = FALSE;
2312 int nthreads = 1;
2314 const char *procstring[] =
2315 { NULL, "np", "n", "none", NULL };
2316 const char *npmevalues_opt[] =
2317 { NULL, "auto", "all", "subset", NULL };
2319 gmx_bool bAppendFiles = TRUE;
2320 gmx_bool bKeepAndNumCPT = FALSE;
2321 gmx_bool bResetCountersHalfWay = FALSE;
2322 gmx_bool bBenchmark = TRUE;
2323 gmx_bool bCheck = TRUE;
2325 output_env_t oenv = NULL;
2327 t_pargs pa[] = {
2328 /***********************/
2329 /* g_tune_pme options: */
2330 /***********************/
2331 { "-mdrun", FALSE, etSTR, {&cmd_mdrun},
2332 "Command line to run a simulation, e.g. 'gmx mdrun' or 'mdrun_mpi'" },
2333 { "-np", FALSE, etINT, {&nnodes},
2334 "Number of ranks to run the tests on (must be > 2 for separate PME ranks)" },
2335 { "-npstring", FALSE, etENUM, {procstring},
2336 "Name of the [TT]$MPIRUN[tt] option that specifies the number of ranks to use ('np', or 'n'; use 'none' if there is no such option)" },
2337 { "-ntmpi", FALSE, etINT, {&nthreads},
2338 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2339 { "-r", FALSE, etINT, {&repeats},
2340 "Repeat each test this often" },
2341 { "-max", FALSE, etREAL, {&maxPMEfraction},
2342 "Max fraction of PME ranks to test with" },
2343 { "-min", FALSE, etREAL, {&minPMEfraction},
2344 "Min fraction of PME ranks to test with" },
2345 { "-npme", FALSE, etENUM, {npmevalues_opt},
2346 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2347 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2348 { "-fix", FALSE, etINT, {&npme_fixed},
2349 "If >= -1, do not vary the number of PME-only ranks, instead use this fixed value and only vary rcoulomb and the PME grid spacing."},
2350 { "-rmax", FALSE, etREAL, {&rmax},
2351 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2352 { "-rmin", FALSE, etREAL, {&rmin},
2353 "If >0, minimal rcoulomb for -ntpr>1" },
2354 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2355 "Scale rvdw along with rcoulomb"},
2356 { "-ntpr", FALSE, etINT, {&ntprs},
2357 "Number of [REF].tpr[ref] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2358 "If < 1, automatically choose the number of [REF].tpr[ref] files to test" },
2359 { "-steps", FALSE, etINT64, {&bench_nsteps},
2360 "Take timings for this many steps in the benchmark runs" },
2361 { "-resetstep", FALSE, etINT, {&presteps},
2362 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2363 { "-nsteps", FALSE, etINT64, {&new_sim_nsteps},
2364 "If non-negative, perform this many steps in the real run (overwrites nsteps from [REF].tpr[ref], add [REF].cpt[ref] steps)" },
2365 { "-launch", FALSE, etBOOL, {&bLaunch},
2366 "Launch the real simulation after optimization" },
2367 { "-bench", FALSE, etBOOL, {&bBenchmark},
2368 "Run the benchmarks or just create the input [REF].tpr[ref] files?" },
2369 { "-check", FALSE, etBOOL, {&bCheck},
2370 "Before the benchmark runs, check whether mdrun works in parallel" },
2371 { "-gpu_id", FALSE, etSTR, {&eligible_gpu_ids},
2372 "List of GPU device id-s that are eligible for use (unlike mdrun, does not imply any mapping)" },
2373 /******************/
2374 /* mdrun options: */
2375 /******************/
2376 /* We let g_tune_pme parse and understand these options, because we need to
2377 * prevent that they appear on the mdrun command line for the benchmarks */
2378 { "-append", FALSE, etBOOL, {&bAppendFiles},
2379 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2380 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2381 "Keep and number checkpoint files (launch only)" },
2382 { "-deffnm", FALSE, etSTR, {&deffnm},
2383 "Set the default filenames (launch only)" },
2384 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2385 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2388 #define NFILE asize(fnm)
2390 seconds = gmx_gettime();
2392 if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2393 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2394 0, NULL, &oenv))
2396 return 0;
2399 /* Store the remaining unparsed command line entries in a string which
2400 * is then attached to the mdrun command line */
2401 snew(ExtraArgs, 1);
2402 ExtraArgs[0] = '\0';
2403 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2405 add_to_string(&ExtraArgs, argv[i]);
2406 add_to_string(&ExtraArgs, " ");
2409 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2411 bThreads = TRUE;
2412 if (opt2parg_bSet("-npstring", asize(pa), pa))
2414 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2417 if (nnodes > 1)
2419 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2421 /* and now we just set this; a bit of an ugly hack*/
2422 nnodes = nthreads;
2424 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2425 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2427 /* Automatically set -beo options if -eo is set etc. */
2428 couple_files_options(NFILE, fnm);
2430 /* Construct the command line arguments for benchmark runs
2431 * as well as for the simulation run */
2432 if (bThreads)
2434 sprintf(bbuf, " -ntmpi %d ", nthreads);
2436 else
2438 /* This string will be used for MPI runs and will appear after the
2439 * mpirun command. */
2440 if (strcmp(procstring[0], "none") != 0)
2442 sprintf(bbuf, " -%s %d ", procstring[0], nnodes);
2444 else
2446 sprintf(bbuf, " ");
2450 cmd_np = bbuf;
2452 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2453 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs, deffnm);
2455 /* Prepare to use checkpoint file if requested */
2456 sim_part = 1;
2457 if (opt2bSet("-cpi", NFILE, fnm))
2459 const char *filename = opt2fn("-cpi", NFILE, fnm);
2460 int cpt_sim_part;
2461 read_checkpoint_part_and_step(filename,
2462 &cpt_sim_part, &cpt_steps);
2463 if (cpt_sim_part == 0)
2465 gmx_fatal(FARGS, "Checkpoint file %s could not be read!", filename);
2467 /* tune_pme will run the next part of the simulation */
2468 sim_part = cpt_sim_part + 1;
2471 /* Open performance output file and write header info */
2472 fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2474 /* Make a quick consistency check of command line parameters */
2475 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2476 maxPMEfraction, minPMEfraction, npme_fixed,
2477 bench_nsteps, fnm, NFILE, sim_part, presteps,
2478 asize(pa), pa);
2479 /* Check any GPU IDs passed make sense, and fill the data structure for them */
2480 snew(gpu_ids, 1);
2481 parse_digits_from_plain_string(eligible_gpu_ids, &gpu_ids->n, &gpu_ids->ids);
2483 /* Determine the maximum and minimum number of PME nodes to test,
2484 * the actual list of settings is build in do_the_tests(). */
2485 if ((nnodes > 2) && (npme_fixed < -1))
2487 if (0 == strcmp(npmevalues_opt[0], "auto"))
2489 /* Determine the npme range automatically based on the PME:PP load guess */
2490 if (guessPMEratio > 1.0)
2492 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2493 maxPMEnodes = nnodes/2;
2494 minPMEnodes = nnodes/2;
2496 else
2498 /* PME : PP load is in the range 0..1, let's test around the guess */
2499 guessPMEnodes = nnodes/(1.0 + 1.0/guessPMEratio);
2500 minPMEnodes = floor(0.7*guessPMEnodes);
2501 maxPMEnodes = ceil(1.6*guessPMEnodes);
2502 maxPMEnodes = min(maxPMEnodes, nnodes/2);
2505 else
2507 /* Determine the npme range based on user input */
2508 maxPMEnodes = floor(maxPMEfraction*nnodes);
2509 minPMEnodes = max(floor(minPMEfraction*nnodes), 0);
2510 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2511 if (maxPMEnodes != minPMEnodes)
2513 fprintf(stdout, "- %d ", maxPMEnodes);
2515 fprintf(stdout, "PME-only ranks.\n Note that the automatic number of PME-only ranks and no separate PME ranks are always tested.\n");
2518 else
2520 maxPMEnodes = 0;
2521 minPMEnodes = 0;
2524 /* Get the commands we need to set up the runs from environment variables */
2525 get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2526 if (bBenchmark && repeats > 0)
2528 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun, NULL != eligible_gpu_ids);
2531 /* Print some header info to file */
2532 sep_line(fp);
2533 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2534 sep_line(fp);
2535 fprintf(fp, "%s for GROMACS %s\n", output_env_get_program_display_name(oenv),
2536 gmx_version());
2537 if (!bThreads)
2539 fprintf(fp, "Number of ranks : %d\n", nnodes);
2540 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2541 if (strcmp(procstring[0], "none") != 0)
2543 fprintf(fp, "Passing # of ranks via : -%s\n", procstring[0]);
2545 else
2547 fprintf(fp, "Not setting number of ranks in system call\n");
2550 else
2552 fprintf(fp, "Number of threads : %d\n", nnodes);
2555 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2556 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2557 fprintf(fp, "Benchmark steps : ");
2558 fprintf(fp, "%"GMX_PRId64, bench_nsteps);
2559 fprintf(fp, "\n");
2560 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2561 if (sim_part > 1)
2563 fprintf(fp, "Checkpoint time step : ");
2564 fprintf(fp, "%"GMX_PRId64, cpt_steps);
2565 fprintf(fp, "\n");
2567 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2569 if (new_sim_nsteps >= 0)
2571 bOverwrite = TRUE;
2572 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2573 fprintf(stderr, "%"GMX_PRId64, new_sim_nsteps+cpt_steps);
2574 fprintf(stderr, " steps.\n");
2575 fprintf(fp, "Simulation steps : ");
2576 fprintf(fp, "%"GMX_PRId64, new_sim_nsteps);
2577 fprintf(fp, "\n");
2579 if (repeats > 1)
2581 fprintf(fp, "Repeats for each test : %d\n", repeats);
2584 if (npme_fixed >= -1)
2586 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2589 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2590 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2592 /* Allocate memory for the inputinfo struct: */
2593 snew(info, 1);
2594 info->nr_inputfiles = ntprs;
2595 for (i = 0; i < ntprs; i++)
2597 snew(info->rcoulomb, ntprs);
2598 snew(info->rvdw, ntprs);
2599 snew(info->rlist, ntprs);
2600 snew(info->rlistlong, ntprs);
2601 snew(info->nkx, ntprs);
2602 snew(info->nky, ntprs);
2603 snew(info->nkz, ntprs);
2604 snew(info->fsx, ntprs);
2605 snew(info->fsy, ntprs);
2606 snew(info->fsz, ntprs);
2608 /* Make alternative tpr files to test: */
2609 snew(tpr_names, ntprs);
2610 for (i = 0; i < ntprs; i++)
2612 snew(tpr_names[i], STRLEN);
2615 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2616 * different grids could be found. */
2617 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2618 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2620 /********************************************************************************/
2621 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2622 /********************************************************************************/
2623 snew(perfdata, ntprs);
2624 if (bBenchmark)
2626 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2627 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2628 cmd_args_bench, fnm, NFILE, presteps, cpt_steps, bCheck, gpu_ids);
2630 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2632 /* Analyse the results and give a suggestion for optimal settings: */
2633 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2634 repeats, info, &best_tpr, &best_npme);
2636 /* Take the best-performing tpr file and enlarge nsteps to original value */
2637 if (bKeepTPR && !bOverwrite)
2639 simulation_tpr = opt2fn("-s", NFILE, fnm);
2641 else
2643 simulation_tpr = opt2fn("-so", NFILE, fnm);
2644 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2645 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2648 /* Let's get rid of the temporary benchmark input files */
2649 for (i = 0; i < ntprs; i++)
2651 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2652 remove(tpr_names[i]);
2655 /* Now start the real simulation if the user requested it ... */
2656 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2657 cmd_args_launch, simulation_tpr, nnodes, best_npme, gpu_ids);
2659 gmx_ffclose(fp);
2661 /* ... or simply print the performance results to screen: */
2662 if (!bLaunch)
2664 finalize(opt2fn("-p", NFILE, fnm));
2667 return 0;