Remove continuation from convert_tpr
[gromacs.git] / src / gromacs / gmxana / gmx_tune_pme.cpp
blob93df67f2503d56645b9b639e869466a2205f8d86
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 <cmath>
40 #include <cstdlib>
41 #include <cstring>
42 #include <ctime>
44 #include <algorithm>
46 #ifdef HAVE_SYS_TIME_H
47 #include <sys/time.h>
48 #endif
50 #include "gromacs/commandline/pargs.h"
51 #include "gromacs/fft/calcgrid.h"
52 #include "gromacs/fileio/checkpoint.h"
53 #include "gromacs/fileio/readinp.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/gmxana/gmx_ana.h"
56 #include "gromacs/math/utilities.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/perf_est.h"
59 #include "gromacs/mdrunutility/mdmodules.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/timing/walltime_accounting.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/utility/arraysize.h"
67 #include "gromacs/utility/baseversion.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/smalloc.h"
74 /* Enum for situations that can occur during log file parsing, the
75 * corresponding string entries can be found in do_the_tests() in
76 * const char* ParseLog[] */
77 /* TODO clean up CamelCasing of these enum names */
78 enum {
79 eParselogOK,
80 eParselogNotFound,
81 eParselogNoPerfData,
82 eParselogTerm,
83 eParselogResetProblem,
84 eParselogNoDDGrid,
85 eParselogTPXVersion,
86 eParselogNotParallel,
87 eParselogLargePrimeFactor,
88 eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs,
89 eParselogGpuProblem,
90 eParselogFatal,
91 eParselogNr
95 typedef struct
97 int nPMEnodes; /* number of PME-only nodes used in this test */
98 int nx, ny, nz; /* DD grid */
99 int guessPME; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
100 double *Gcycles; /* This can contain more than one value if doing multiple tests */
101 double Gcycles_Av;
102 float *ns_per_day;
103 float ns_per_day_Av;
104 float *PME_f_load; /* PME mesh/force load average*/
105 float PME_f_load_Av; /* Average average ;) ... */
106 char *mdrun_cmd_line; /* Mdrun command line used for this test */
107 } t_perf;
110 typedef struct
112 int nr_inputfiles; /* The number of tpr and mdp input files */
113 gmx_int64_t orig_sim_steps; /* Number of steps to be done in the real simulation */
114 gmx_int64_t orig_init_step; /* Init step for the real simulation */
115 real *rcoulomb; /* The coulomb radii [0...nr_inputfiles] */
116 real *rvdw; /* The vdW radii */
117 real *rlist; /* Neighbourlist cutoff radius */
118 int *nkx, *nky, *nkz;
119 real *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
120 } t_inputinfo;
123 static void sep_line(FILE *fp)
125 fprintf(fp, "\n------------------------------------------------------------\n");
129 /* Wrapper for system calls */
130 static int gmx_system_call(char *command)
132 return ( system(command) );
136 /* Check if string starts with substring */
137 static gmx_bool str_starts(const char *string, const char *substring)
139 return ( std::strncmp(string, substring, std::strlen(substring)) == 0);
143 static void cleandata(t_perf *perfdata, int test_nr)
145 perfdata->Gcycles[test_nr] = 0.0;
146 perfdata->ns_per_day[test_nr] = 0.0;
147 perfdata->PME_f_load[test_nr] = 0.0;
149 return;
153 static void remove_if_exists(const char *fn)
155 if (gmx_fexist(fn))
157 fprintf(stdout, "Deleting %s\n", fn);
158 remove(fn);
163 static void finalize(const char *fn_out)
165 char buf[STRLEN];
166 FILE *fp;
169 fp = fopen(fn_out, "r");
170 fprintf(stdout, "\n\n");
172 while (fgets(buf, STRLEN-1, fp) != NULL)
174 fprintf(stdout, "%s", buf);
176 fclose(fp);
177 fprintf(stdout, "\n\n");
181 enum {
182 eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr
185 static int parse_logfile(const char *logfile, const char *errfile,
186 t_perf *perfdata, int test_nr, int presteps, gmx_int64_t cpt_steps,
187 int nnodes)
189 FILE *fp;
190 char line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
191 const char matchstrdd[] = "Domain decomposition grid";
192 const char matchstrcr[] = "resetting all time and cycle counters";
193 const char matchstrbal[] = "Average PME mesh/force load:";
194 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";
195 const char errSIG[] = "signal, stopping at the next";
196 int iFound;
197 float dum1, dum2, dum3, dum4;
198 int ndum;
199 int npme;
200 gmx_int64_t resetsteps = -1;
201 gmx_bool bFoundResetStr = FALSE;
202 gmx_bool bResetChecked = FALSE;
205 if (!gmx_fexist(logfile))
207 fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
208 cleandata(perfdata, test_nr);
209 return eParselogNotFound;
212 fp = fopen(logfile, "r");
213 perfdata->PME_f_load[test_nr] = -1.0;
214 perfdata->guessPME = -1;
216 iFound = eFoundNothing;
217 if (1 == nnodes)
219 iFound = eFoundDDStr; /* Skip some case statements */
222 while (fgets(line, STRLEN, fp) != NULL)
224 /* Remove leading spaces */
225 ltrim(line);
227 /* Check for TERM and INT signals from user: */
228 if (std::strstr(line, errSIG) != NULL)
230 fclose(fp);
231 cleandata(perfdata, test_nr);
232 return eParselogTerm;
235 /* Check whether cycle resetting worked */
236 if (presteps > 0 && !bFoundResetStr)
238 if (std::strstr(line, matchstrcr) != NULL)
240 sprintf(dumstring, "step %s", "%" GMX_SCNd64);
241 sscanf(line, dumstring, &resetsteps);
242 bFoundResetStr = TRUE;
243 if (resetsteps == presteps+cpt_steps)
245 bResetChecked = TRUE;
247 else
249 sprintf(dumstring, "%" GMX_PRId64, resetsteps);
250 sprintf(dumstring2, "%" GMX_PRId64, presteps+cpt_steps);
251 fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
252 " though they were supposed to be reset at step %s!\n",
253 dumstring, dumstring2);
258 /* Look for strings that appear in a certain order in the log file: */
259 switch (iFound)
261 case eFoundNothing:
262 /* Look for domain decomp grid and separate PME nodes: */
263 if (str_starts(line, matchstrdd))
265 sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME ranks %d",
266 &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
267 if (perfdata->nPMEnodes == -1)
269 perfdata->guessPME = npme;
271 else if (perfdata->nPMEnodes != npme)
273 gmx_fatal(FARGS, "PME ranks from command line and output file are not identical");
275 iFound = eFoundDDStr;
277 /* Catch a few errors that might have occurred: */
278 else if (str_starts(line, "There is no domain decomposition for"))
280 fclose(fp);
281 return eParselogNoDDGrid;
283 else if (str_starts(line, "The number of ranks you selected"))
285 fclose(fp);
286 return eParselogLargePrimeFactor;
288 else if (str_starts(line, "reading tpx file"))
290 fclose(fp);
291 return eParselogTPXVersion;
293 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
295 fclose(fp);
296 return eParselogNotParallel;
298 break;
299 case eFoundDDStr:
300 /* Even after the "Domain decomposition grid" string was found,
301 * it could be that mdrun had to quit due to some error. */
302 if (str_starts(line, "Incorrect launch configuration: mismatching number of"))
304 fclose(fp);
305 return eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs;
307 else if (str_starts(line, "Some of the requested GPUs do not exist"))
309 fclose(fp);
310 return eParselogGpuProblem;
312 /* Look for PME mesh/force balance (not necessarily present, though) */
313 else if (str_starts(line, matchstrbal))
315 sscanf(&line[std::strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
317 /* Look for matchstring */
318 else if (str_starts(line, matchstring))
320 iFound = eFoundAccountingStr;
322 break;
323 case eFoundAccountingStr:
324 /* Already found matchstring - look for cycle data */
325 if (str_starts(line, "Total "))
327 sscanf(line, "Total %*f %lf", &(perfdata->Gcycles[test_nr]));
328 iFound = eFoundCycleStr;
330 break;
331 case eFoundCycleStr:
332 /* Already found cycle data - look for remaining performance info and return */
333 if (str_starts(line, "Performance:"))
335 ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
336 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
337 perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
338 fclose(fp);
339 if (bResetChecked || presteps == 0)
341 return eParselogOK;
343 else
345 return eParselogResetProblem;
348 break;
350 } /* while */
352 /* Close the log file */
353 fclose(fp);
355 /* Check why there is no performance data in the log file.
356 * Did a fatal errors occur? */
357 if (gmx_fexist(errfile))
359 fp = fopen(errfile, "r");
360 while (fgets(line, STRLEN, fp) != NULL)
362 if (str_starts(line, "Fatal error:") )
364 if (fgets(line, STRLEN, fp) != NULL)
366 fprintf(stderr, "\nWARNING: An error occurred during this benchmark:\n"
367 "%s\n", line);
369 fclose(fp);
370 cleandata(perfdata, test_nr);
371 return eParselogFatal;
374 fclose(fp);
376 else
378 fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
381 /* Giving up ... we could not find out why there is no performance data in
382 * the log file. */
383 fprintf(stdout, "No performance data in log file.\n");
384 cleandata(perfdata, test_nr);
386 return eParselogNoPerfData;
390 static gmx_bool analyze_data(
391 FILE *fp,
392 const char *fn,
393 t_perf **perfdata,
394 int nnodes,
395 int ntprs,
396 int ntests,
397 int nrepeats,
398 t_inputinfo *info,
399 int *index_tpr, /* OUT: Nr of mdp file with best settings */
400 int *npme_optimal) /* OUT: Optimal number of PME nodes */
402 int i, j, k;
403 int line = 0, line_win = -1;
404 int k_win = -1, i_win = -1, winPME;
405 double s = 0.0; /* standard deviation */
406 t_perf *pd;
407 char strbuf[STRLEN];
408 char str_PME_f_load[13];
409 gmx_bool bCanUseOrigTPR;
410 gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
413 if (nrepeats > 1)
415 sep_line(fp);
416 fprintf(fp, "Summary of successful runs:\n");
417 fprintf(fp, "Line tpr PME ranks Gcycles Av. Std.dev. ns/day PME/f");
418 if (nnodes > 1)
420 fprintf(fp, " DD grid");
422 fprintf(fp, "\n");
426 for (k = 0; k < ntprs; k++)
428 for (i = 0; i < ntests; i++)
430 /* Select the right dataset: */
431 pd = &(perfdata[k][i]);
433 pd->Gcycles_Av = 0.0;
434 pd->PME_f_load_Av = 0.0;
435 pd->ns_per_day_Av = 0.0;
437 if (pd->nPMEnodes == -1)
439 sprintf(strbuf, "(%3d)", pd->guessPME);
441 else
443 sprintf(strbuf, " ");
446 /* Get the average run time of a setting */
447 for (j = 0; j < nrepeats; j++)
449 pd->Gcycles_Av += pd->Gcycles[j];
450 pd->PME_f_load_Av += pd->PME_f_load[j];
452 pd->Gcycles_Av /= nrepeats;
453 pd->PME_f_load_Av /= nrepeats;
455 for (j = 0; j < nrepeats; j++)
457 if (pd->ns_per_day[j] > 0.0)
459 pd->ns_per_day_Av += pd->ns_per_day[j];
461 else
463 /* Somehow the performance number was not aquired for this run,
464 * therefor set the average to some negative value: */
465 pd->ns_per_day_Av = -1.0f*nrepeats;
466 break;
469 pd->ns_per_day_Av /= nrepeats;
471 /* Nicer output: */
472 if (pd->PME_f_load_Av > 0.0)
474 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
476 else
478 sprintf(str_PME_f_load, "%s", " - ");
482 /* We assume we had a successful run if both averages are positive */
483 if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
485 /* Output statistics if repeats were done */
486 if (nrepeats > 1)
488 /* Calculate the standard deviation */
489 s = 0.0;
490 for (j = 0; j < nrepeats; j++)
492 s += gmx::square( pd->Gcycles[j] - pd->Gcycles_Av );
494 s /= (nrepeats - 1);
495 s = std::sqrt(s);
497 fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
498 line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
499 pd->ns_per_day_Av, str_PME_f_load);
500 if (nnodes > 1)
502 fprintf(fp, " %3d %3d %3d", pd->nx, pd->ny, pd->nz);
504 fprintf(fp, "\n");
506 /* Store the index of the best run found so far in 'winner': */
507 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
509 k_win = k;
510 i_win = i;
511 line_win = line;
513 line++;
518 if (k_win == -1)
520 gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
523 sep_line(fp);
525 winPME = perfdata[k_win][i_win].nPMEnodes;
527 if (1 == ntests)
529 /* We stuck to a fixed number of PME-only nodes */
530 sprintf(strbuf, "settings No. %d", k_win);
532 else
534 /* We have optimized the number of PME-only nodes */
535 if (winPME == -1)
537 sprintf(strbuf, "%s", "the automatic number of PME ranks");
539 else
541 sprintf(strbuf, "%d PME ranks", winPME);
544 fprintf(fp, "Best performance was achieved with %s", strbuf);
545 if ((nrepeats > 1) && (ntests > 1))
547 fprintf(fp, " (see line %d)", line_win);
549 fprintf(fp, "\n");
551 /* Only mention settings if they were modified: */
552 bRefinedCoul = !gmx_within_tol(info->rcoulomb[k_win], info->rcoulomb[0], GMX_REAL_EPS);
553 bRefinedVdW = !gmx_within_tol(info->rvdw[k_win], info->rvdw[0], GMX_REAL_EPS);
554 bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
555 info->nky[k_win] == info->nky[0] &&
556 info->nkz[k_win] == info->nkz[0]);
558 if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
560 fprintf(fp, "Optimized PME settings:\n");
561 bCanUseOrigTPR = FALSE;
563 else
565 bCanUseOrigTPR = TRUE;
568 if (bRefinedCoul)
570 fprintf(fp, " New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
573 if (bRefinedVdW)
575 fprintf(fp, " New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
578 if (bRefinedGrid)
580 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],
581 info->nkx[0], info->nky[0], info->nkz[0]);
584 if (bCanUseOrigTPR && ntprs > 1)
586 fprintf(fp, "and original PME settings.\n");
589 fflush(fp);
591 /* Return the index of the mdp file that showed the highest performance
592 * and the optimal number of PME nodes */
593 *index_tpr = k_win;
594 *npme_optimal = winPME;
596 return bCanUseOrigTPR;
600 /* Get the commands we need to set up the runs from environment variables */
601 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
603 char *cp;
604 const char def_mpirun[] = "mpirun";
606 const char empty_mpirun[] = "";
608 /* Get the commands we need to set up the runs from environment variables */
609 if (!bThreads)
611 if ( (cp = getenv("MPIRUN")) != NULL)
613 *cmd_mpirun = gmx_strdup(cp);
615 else
617 *cmd_mpirun = gmx_strdup(def_mpirun);
620 else
622 *cmd_mpirun = gmx_strdup(empty_mpirun);
625 if (*cmd_mdrun == NULL)
627 /* The use of MDRUN is deprecated, but made available in 5.1
628 for backward compatibility. It may be removed in a future
629 version. */
630 if ( (cp = getenv("MDRUN" )) != NULL)
632 *cmd_mdrun = gmx_strdup(cp);
634 else
636 gmx_fatal(FARGS, "The way to call mdrun must be set in the -mdrun command-line flag.");
641 /* Check that the commands will run mdrun (perhaps via mpirun) by
642 * running a very quick test simulation. Requires MPI environment or
643 * GPU support to be available if applicable. */
644 /* TODO implement feature to parse the log file to get the list of
645 compatible GPUs from mdrun, if the user of gmx tune-pme has not
646 given one. */
647 static void check_mdrun_works(gmx_bool bThreads,
648 const char *cmd_mpirun,
649 const char *cmd_np,
650 const char *cmd_mdrun,
651 gmx_bool bNeedGpuSupport)
653 char *command = NULL;
654 char *cp;
655 char line[STRLEN];
656 FILE *fp;
657 const char filename[] = "benchtest.log";
659 /* This string should always be identical to the one in copyrite.c,
660 * gmx_print_version_info() in the GMX_MPI section */
661 const char match_mpi[] = "MPI library: MPI";
662 const char match_mdrun[] = "Executable: ";
663 const char match_gpu[] = "GPU support: enabled";
664 gmx_bool bMdrun = FALSE;
665 gmx_bool bMPI = FALSE;
666 gmx_bool bHaveGpuSupport = FALSE;
668 /* Run a small test to see whether mpirun + mdrun work */
669 fprintf(stdout, "Making sure that mdrun can be executed. ");
670 if (bThreads)
672 snew(command, std::strlen(cmd_mdrun) + std::strlen(cmd_np) + std::strlen(filename) + 50);
673 sprintf(command, "%s%s-version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
675 else
677 snew(command, std::strlen(cmd_mpirun) + std::strlen(cmd_np) + std::strlen(cmd_mdrun) + std::strlen(filename) + 50);
678 sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
680 fprintf(stdout, "Trying '%s' ... ", command);
681 make_backup(filename);
682 gmx_system_call(command);
684 /* Check if we find the characteristic string in the output: */
685 if (!gmx_fexist(filename))
687 gmx_fatal(FARGS, "Output from test run could not be found.");
690 fp = fopen(filename, "r");
691 /* We need to scan the whole output file, since sometimes the queuing system
692 * also writes stuff to stdout/err */
693 while (!feof(fp) )
695 cp = fgets(line, STRLEN, fp);
696 if (cp != NULL)
698 if (str_starts(line, match_mdrun) )
700 bMdrun = TRUE;
702 if (str_starts(line, match_mpi) )
704 bMPI = TRUE;
706 if (str_starts(line, match_gpu) )
708 bHaveGpuSupport = TRUE;
712 fclose(fp);
714 if (bThreads)
716 if (bMPI)
718 gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
719 "(%s)\n"
720 "seems to have been compiled with MPI instead.",
721 cmd_mdrun);
724 else
726 if (bMdrun && !bMPI)
728 gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
729 "(%s)\n"
730 "seems to have been compiled without MPI support.",
731 cmd_mdrun);
735 if (!bMdrun)
737 gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
738 filename);
741 if (bNeedGpuSupport && !bHaveGpuSupport)
743 gmx_fatal(FARGS, "The mdrun executable did not have the expected GPU support.");
746 fprintf(stdout, "passed.\n");
748 /* Clean up ... */
749 remove(filename);
750 sfree(command);
753 /*! \brief Helper struct so we can parse the string with eligible GPU
754 IDs outside do_the_tests. */
755 typedef struct eligible_gpu_ids
757 int n; /**< Length of ids */
758 int *ids; /**< Array of length n. NULL if no GPUs in use */
759 } t_eligible_gpu_ids;
761 /* Handles the no-GPU case by emitting an empty string. */
762 static char *make_gpu_id_command_line(int numRanks, int numPmeRanks, const t_eligible_gpu_ids *gpu_ids)
764 char *command_line, *ptr;
765 const char *flag = "-gpu_id ";
766 int flag_length;
768 /* Reserve enough room for the option name, enough single-digit
769 GPU ids (since that is currently all that is possible to use
770 with mdrun), and a terminating NULL. */
771 flag_length = std::strlen(flag);
772 snew(command_line, flag_length + numRanks + 1);
773 ptr = command_line;
775 /* If the user has given no eligible GPU IDs, or we're trying the
776 * default behaviour, then there is nothing for g_tune_pme to give
777 * to mdrun -gpu_id */
778 if (gpu_ids->n > 0 && numPmeRanks > -1)
780 int numPpRanks, max_num_ranks_for_each_GPU;
781 int gpu_id, rank;
783 /* Write the option flag */
784 std::strcpy(ptr, flag);
785 ptr += flag_length;
787 numPpRanks = numRanks - numPmeRanks;
788 max_num_ranks_for_each_GPU = numPpRanks / gpu_ids->n;
789 if (max_num_ranks_for_each_GPU * gpu_ids->n != numPpRanks)
791 /* Some GPUs will receive more work than others, which
792 * we choose to be those with the lowest indices */
793 max_num_ranks_for_each_GPU++;
796 /* Loop over all eligible GPU ids */
797 for (gpu_id = 0, rank = 0; gpu_id < gpu_ids->n; gpu_id++)
799 int rank_for_this_GPU;
800 /* Loop over all PP ranks for GPU with ID gpu_id, building the
801 assignment string. */
802 for (rank_for_this_GPU = 0;
803 rank_for_this_GPU < max_num_ranks_for_each_GPU && rank < numPpRanks;
804 rank++, rank_for_this_GPU++)
806 *ptr = '0' + gpu_ids->ids[gpu_id];
807 ptr++;
811 *ptr = '\0';
813 return command_line;
816 static void launch_simulation(
817 gmx_bool bLaunch, /* Should the simulation be launched? */
818 FILE *fp, /* General log file */
819 gmx_bool bThreads, /* whether to use threads */
820 char *cmd_mpirun, /* Command for mpirun */
821 char *cmd_np, /* Switch for -np or -ntmpi or empty */
822 char *cmd_mdrun, /* Command for mdrun */
823 char *args_for_mdrun, /* Arguments for mdrun */
824 const char *simulation_tpr, /* This tpr will be simulated */
825 int nnodes, /* Number of ranks to use */
826 int nPMEnodes, /* Number of PME ranks to use */
827 const t_eligible_gpu_ids *gpu_ids) /* Struct containing GPU IDs for
828 * constructing mdrun command lines */
830 char *command, *cmd_gpu_ids;
833 /* Make enough space for the system call command,
834 * (200 extra chars for -npme ... etc. options should suffice): */
835 snew(command, std::strlen(cmd_mpirun)+std::strlen(cmd_mdrun)+std::strlen(cmd_np)+std::strlen(args_for_mdrun)+std::strlen(simulation_tpr)+200);
837 cmd_gpu_ids = make_gpu_id_command_line(nnodes, nPMEnodes, gpu_ids);
839 /* Note that the -passall options requires args_for_mdrun to be at the end
840 * of the command line string */
841 if (bThreads)
843 sprintf(command, "%s%s-npme %d -s %s %s %s",
844 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids);
846 else
848 sprintf(command, "%s%s%s -npme %d -s %s %s %s",
849 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids);
852 fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
853 sep_line(fp);
854 fflush(fp);
856 /* Now the real thing! */
857 if (bLaunch)
859 fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
860 sep_line(stdout);
861 fflush(stdout);
862 gmx_system_call(command);
867 static void modify_PMEsettings(
868 gmx_int64_t simsteps, /* Set this value as number of time steps */
869 gmx_int64_t init_step, /* Set this value as init_step */
870 const char *fn_best_tpr, /* tpr file with the best performance */
871 const char *fn_sim_tpr) /* name of tpr file to be launched */
873 t_inputrec *ir;
874 t_state state;
875 gmx_mtop_t mtop;
876 char buf[200];
878 gmx::MDModules mdModules;
879 ir = mdModules.inputrec();
880 read_tpx_state(fn_best_tpr, ir, &state, &mtop);
882 /* Reset nsteps and init_step to the value of the input .tpr file */
883 ir->nsteps = simsteps;
884 ir->init_step = init_step;
886 /* Write the tpr file which will be launched */
887 sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%" GMX_PRId64);
888 fprintf(stdout, buf, ir->nsteps);
889 fflush(stdout);
890 write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
893 static gmx_bool can_scale_rvdw(int vdwtype)
895 return (evdwCUT == vdwtype ||
896 evdwPME == vdwtype);
899 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
901 /* Make additional TPR files with more computational load for the
902 * direct space processors: */
903 static void make_benchmark_tprs(
904 const char *fn_sim_tpr, /* READ : User-provided tpr file */
905 char *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files */
906 gmx_int64_t benchsteps, /* Number of time steps for benchmark runs */
907 gmx_int64_t statesteps, /* Step counter in checkpoint file */
908 real rmin, /* Minimal Coulomb radius */
909 real rmax, /* Maximal Coulomb radius */
910 real bScaleRvdw, /* Scale rvdw along with rcoulomb */
911 int *ntprs, /* No. of TPRs to write, each with a different
912 rcoulomb and fourierspacing */
913 t_inputinfo *info, /* Contains information about mdp file options */
914 FILE *fp) /* Write the output here */
916 int i, j, d;
917 t_inputrec *ir;
918 t_state state;
919 gmx_mtop_t mtop;
920 real nlist_buffer; /* Thickness of the buffer regions for PME-switch potentials */
921 char buf[200];
922 rvec box_size;
923 gmx_bool bNote = FALSE;
924 real add; /* Add this to rcoul for the next test */
925 real fac = 1.0; /* Scaling factor for Coulomb radius */
926 real fourierspacing; /* Basic fourierspacing from tpr */
929 sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
930 *ntprs > 1 ? "s" : "", "%" GMX_PRId64, benchsteps > 1 ? "s" : "");
931 fprintf(stdout, buf, benchsteps);
932 if (statesteps > 0)
934 sprintf(buf, " (adding %s steps from checkpoint file)", "%" GMX_PRId64);
935 fprintf(stdout, buf, statesteps);
936 benchsteps += statesteps;
938 fprintf(stdout, ".\n");
940 gmx::MDModules mdModules;
941 ir = mdModules.inputrec();
942 read_tpx_state(fn_sim_tpr, ir, &state, &mtop);
944 /* Check if some kind of PME was chosen */
945 if (EEL_PME(ir->coulombtype) == FALSE)
947 gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
948 EELTYPE(eelPME));
951 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
952 if ( (ir->cutoff_scheme != ecutsVERLET) &&
953 (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
955 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
956 EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
958 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
959 else if (ir->rcoulomb > ir->rlist)
961 gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
962 EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
965 if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
967 fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
968 bScaleRvdw = FALSE;
971 /* Reduce the number of steps for the benchmarks */
972 info->orig_sim_steps = ir->nsteps;
973 ir->nsteps = benchsteps;
974 /* We must not use init_step from the input tpr file for the benchmarks */
975 info->orig_init_step = ir->init_step;
976 ir->init_step = 0;
978 /* For PME-switch potentials, keep the radial distance of the buffer region */
979 nlist_buffer = ir->rlist - ir->rcoulomb;
981 /* Determine length of triclinic box vectors */
982 for (d = 0; d < DIM; d++)
984 box_size[d] = 0;
985 for (i = 0; i < DIM; i++)
987 box_size[d] += state.box[d][i]*state.box[d][i];
989 box_size[d] = std::sqrt(box_size[d]);
992 if (ir->fourier_spacing > 0)
994 info->fsx[0] = ir->fourier_spacing;
995 info->fsy[0] = ir->fourier_spacing;
996 info->fsz[0] = ir->fourier_spacing;
998 else
1000 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
1001 info->fsx[0] = box_size[XX]/ir->nkx;
1002 info->fsy[0] = box_size[YY]/ir->nky;
1003 info->fsz[0] = box_size[ZZ]/ir->nkz;
1006 /* If no value for the fourierspacing was provided on the command line, we
1007 * use the reconstruction from the tpr file */
1008 if (ir->fourier_spacing > 0)
1010 /* Use the spacing from the tpr */
1011 fourierspacing = ir->fourier_spacing;
1013 else
1015 /* Use the maximum observed spacing */
1016 fourierspacing = std::max(std::max(info->fsx[0], info->fsy[0]), info->fsz[0]);
1019 fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
1021 /* For performance comparisons the number of particles is useful to have */
1022 fprintf(fp, " Number of particles : %d\n", mtop.natoms);
1024 /* Print information about settings of which some are potentially modified: */
1025 fprintf(fp, " Coulomb type : %s\n", EELTYPE(ir->coulombtype));
1026 fprintf(fp, " Grid spacing x y z : %f %f %f\n",
1027 box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
1028 fprintf(fp, " Van der Waals type : %s\n", EVDWTYPE(ir->vdwtype));
1029 if (ir_vdw_switched(ir))
1031 fprintf(fp, " rvdw_switch : %f nm\n", ir->rvdw_switch);
1033 if (EPME_SWITCHED(ir->coulombtype))
1035 fprintf(fp, " rlist : %f nm\n", ir->rlist);
1038 /* Print a descriptive line about the tpr settings tested */
1039 fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
1040 fprintf(fp, " No. scaling rcoulomb");
1041 fprintf(fp, " nkx nky nkz");
1042 fprintf(fp, " spacing");
1043 if (can_scale_rvdw(ir->vdwtype))
1045 fprintf(fp, " rvdw");
1047 if (EPME_SWITCHED(ir->coulombtype))
1049 fprintf(fp, " rlist");
1051 fprintf(fp, " tpr file\n");
1053 /* Loop to create the requested number of tpr input files */
1054 for (j = 0; j < *ntprs; j++)
1056 /* The first .tpr is the provided one, just need to modify nsteps,
1057 * so skip the following block */
1058 if (j != 0)
1060 /* Determine which Coulomb radii rc to use in the benchmarks */
1061 add = (rmax-rmin)/(*ntprs-1);
1062 if (gmx_within_tol(rmin, info->rcoulomb[0], GMX_REAL_EPS))
1064 ir->rcoulomb = rmin + j*add;
1066 else if (gmx_within_tol(rmax, info->rcoulomb[0], GMX_REAL_EPS))
1068 ir->rcoulomb = rmin + (j-1)*add;
1070 else
1072 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
1073 add = (rmax-rmin)/(*ntprs-2);
1074 ir->rcoulomb = rmin + (j-1)*add;
1077 /* Determine the scaling factor fac */
1078 fac = ir->rcoulomb/info->rcoulomb[0];
1080 /* Scale the Fourier grid spacing */
1081 ir->nkx = ir->nky = ir->nkz = 0;
1082 calc_grid(NULL, state.box, fourierspacing*fac, &ir->nkx, &ir->nky, &ir->nkz);
1084 /* Adjust other radii since various conditions need to be fulfilled */
1085 if (eelPME == ir->coulombtype)
1087 /* plain PME, rcoulomb must be equal to rlist TODO only in the group scheme? */
1088 ir->rlist = ir->rcoulomb;
1090 else
1092 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1093 ir->rlist = ir->rcoulomb + nlist_buffer;
1096 if (bScaleRvdw && can_scale_rvdw(ir->vdwtype))
1098 if (ecutsVERLET == ir->cutoff_scheme ||
1099 evdwPME == ir->vdwtype)
1101 /* With either the Verlet cutoff-scheme or LJ-PME,
1102 the van der Waals radius must always equal the
1103 Coulomb radius */
1104 ir->rvdw = ir->rcoulomb;
1106 else
1108 /* For vdw cutoff, rvdw >= rlist */
1109 ir->rvdw = std::max(info->rvdw[0], ir->rlist);
1112 } /* end of "if (j != 0)" */
1114 /* for j==0: Save the original settings
1115 * for j >0: Save modified radii and Fourier grids */
1116 info->rcoulomb[j] = ir->rcoulomb;
1117 info->rvdw[j] = ir->rvdw;
1118 info->nkx[j] = ir->nkx;
1119 info->nky[j] = ir->nky;
1120 info->nkz[j] = ir->nkz;
1121 info->rlist[j] = ir->rlist;
1122 info->fsx[j] = fac*fourierspacing;
1123 info->fsy[j] = fac*fourierspacing;
1124 info->fsz[j] = fac*fourierspacing;
1126 /* Write the benchmark tpr file */
1127 std::strncpy(fn_bench_tprs[j], fn_sim_tpr, std::strlen(fn_sim_tpr)-std::strlen(".tpr"));
1128 sprintf(buf, "_bench%.2d.tpr", j);
1129 std::strcat(fn_bench_tprs[j], buf);
1130 fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1131 fprintf(stdout, "%" GMX_PRId64, ir->nsteps);
1132 if (j > 0)
1134 fprintf(stdout, ", scaling factor %f\n", fac);
1136 else
1138 fprintf(stdout, ", unmodified settings\n");
1141 write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1143 /* Write information about modified tpr settings to log file */
1144 fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb);
1145 fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz);
1146 fprintf(fp, " %9f ", info->fsx[j]);
1147 if (can_scale_rvdw(ir->vdwtype))
1149 fprintf(fp, "%10f", ir->rvdw);
1151 if (EPME_SWITCHED(ir->coulombtype))
1153 fprintf(fp, "%10f", ir->rlist);
1155 fprintf(fp, " %-14s\n", fn_bench_tprs[j]);
1157 /* Make it clear to the user that some additional settings were modified */
1158 if (!gmx_within_tol(ir->rvdw, info->rvdw[0], GMX_REAL_EPS)
1159 || !gmx_within_tol(ir->rlist, info->rlist[0], GMX_REAL_EPS) )
1161 bNote = TRUE;
1164 if (bNote)
1166 fprintf(fp, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1167 "other input settings were also changed (see table above).\n"
1168 "Please check if the modified settings are appropriate.\n");
1170 fflush(stdout);
1171 fflush(fp);
1175 /* Rename the files we want to keep to some meaningful filename and
1176 * delete the rest */
1177 static void cleanup(const t_filenm *fnm, int nfile, int k, int nnodes,
1178 int nPMEnodes, int nr, gmx_bool bKeepStderr)
1180 char numstring[STRLEN];
1181 char newfilename[STRLEN];
1182 const char *fn = NULL;
1183 int i;
1184 const char *opt;
1187 fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1189 for (i = 0; i < nfile; i++)
1191 opt = (char *)fnm[i].opt;
1192 if (std::strcmp(opt, "-p") == 0)
1194 /* do nothing; keep this file */
1197 else if (std::strcmp(opt, "-bg") == 0)
1199 /* Give the log file a nice name so one can later see which parameters were used */
1200 numstring[0] = '\0';
1201 if (nr > 0)
1203 sprintf(numstring, "_%d", nr);
1205 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1206 if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1208 fprintf(stdout, "renaming log file to %s\n", newfilename);
1209 make_backup(newfilename);
1210 rename(opt2fn("-bg", nfile, fnm), newfilename);
1213 else if (std::strcmp(opt, "-err") == 0)
1215 /* This file contains the output of stderr. We want to keep it in
1216 * cases where there have been problems. */
1217 fn = opt2fn(opt, nfile, fnm);
1218 numstring[0] = '\0';
1219 if (nr > 0)
1221 sprintf(numstring, "_%d", nr);
1223 sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1224 if (gmx_fexist(fn))
1226 if (bKeepStderr)
1228 fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1229 make_backup(newfilename);
1230 rename(fn, newfilename);
1232 else
1234 fprintf(stdout, "Deleting %s\n", fn);
1235 remove(fn);
1239 /* Delete the files which are created for each benchmark run: (options -b*) */
1240 else if ( (0 == std::strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1242 remove_if_exists(opt2fn(opt, nfile, fnm));
1248 enum {
1249 eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1252 /* Create a list of numbers of PME nodes to test */
1253 static void make_npme_list(
1254 const char *npmevalues_opt, /* Make a complete list with all
1255 * possibilities or a short list that keeps only
1256 * reasonable numbers of PME nodes */
1257 int *nentries, /* Number of entries we put in the nPMEnodes list */
1258 int *nPMEnodes[], /* Each entry contains the value for -npme */
1259 int nnodes, /* Total number of nodes to do the tests on */
1260 int minPMEnodes, /* Minimum number of PME nodes */
1261 int maxPMEnodes) /* Maximum number of PME nodes */
1263 int i, npme, npp;
1264 int min_factor = 1; /* We request that npp and npme have this minimal
1265 * largest common factor (depends on npp) */
1266 int nlistmax; /* Max. list size */
1267 int nlist; /* Actual number of entries in list */
1268 int eNPME = 0;
1271 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1272 if (!std::strcmp(npmevalues_opt, "all") )
1274 eNPME = eNpmeAll;
1276 else if (!std::strcmp(npmevalues_opt, "subset") )
1278 eNPME = eNpmeSubset;
1280 else /* "auto" or "range" */
1282 if (nnodes <= 64)
1284 eNPME = eNpmeAll;
1286 else if (nnodes < 128)
1288 eNPME = eNpmeReduced;
1290 else
1292 eNPME = eNpmeSubset;
1296 /* Calculate how many entries we could possibly have (in case of -npme all) */
1297 if (nnodes > 2)
1299 nlistmax = maxPMEnodes - minPMEnodes + 3;
1300 if (0 == minPMEnodes)
1302 nlistmax--;
1305 else
1307 nlistmax = 1;
1310 /* Now make the actual list which is at most of size nlist */
1311 snew(*nPMEnodes, nlistmax);
1312 nlist = 0; /* start counting again, now the real entries in the list */
1313 for (i = 0; i < nlistmax - 2; i++)
1315 npme = maxPMEnodes - i;
1316 npp = nnodes-npme;
1317 switch (eNPME)
1319 case eNpmeAll:
1320 min_factor = 1;
1321 break;
1322 case eNpmeReduced:
1323 min_factor = 2;
1324 break;
1325 case eNpmeSubset:
1326 /* For 2d PME we want a common largest factor of at least the cube
1327 * root of the number of PP nodes */
1328 min_factor = static_cast<int>(std::cbrt(npp));
1329 break;
1330 default:
1331 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1332 break;
1334 if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1336 (*nPMEnodes)[nlist] = npme;
1337 nlist++;
1340 /* We always test 0 PME nodes and the automatic number */
1341 *nentries = nlist + 2;
1342 (*nPMEnodes)[nlist ] = 0;
1343 (*nPMEnodes)[nlist+1] = -1;
1345 fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1346 for (i = 0; i < *nentries-1; i++)
1348 fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1350 fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1354 /* Allocate memory to store the performance data */
1355 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1357 int i, j, k;
1360 for (k = 0; k < ntprs; k++)
1362 snew(perfdata[k], datasets);
1363 for (i = 0; i < datasets; i++)
1365 for (j = 0; j < repeats; j++)
1367 snew(perfdata[k][i].Gcycles, repeats);
1368 snew(perfdata[k][i].ns_per_day, repeats);
1369 snew(perfdata[k][i].PME_f_load, repeats);
1376 /* Check for errors on mdrun -h */
1377 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp,
1378 const t_filenm *fnm, int nfile)
1380 char *command, *msg;
1381 int ret;
1383 snew(command, length + 15);
1384 snew(msg, length + 500);
1386 fprintf(stdout, "Making sure the benchmarks can be executed by running just 1 step...\n");
1387 sprintf(command, "%s -nsteps 1 -quiet", mdrun_cmd_line);
1388 fprintf(stdout, "Executing '%s' ...\n", command);
1389 ret = gmx_system_call(command);
1391 if (0 != ret)
1393 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1394 * get the error message from mdrun itself */
1395 sprintf(msg,
1396 "Cannot run the first benchmark simulation! Please check the error message of\n"
1397 "mdrun for the source of the problem. Did you provide a command line\n"
1398 "argument that neither gmx tune_pme nor mdrun understands? If you're\n"
1399 "sure your command line should work, you can bypass this check with \n"
1400 "gmx tune_pme -nocheck. The failing command was:\n"
1401 "\n%s\n\n", command);
1403 fprintf(stderr, "%s", msg);
1404 sep_line(fp);
1405 fprintf(fp, "%s", msg);
1407 exit(ret);
1409 fprintf(stdout, "Benchmarks can be executed!\n");
1411 /* Clean up the benchmark output files we just created */
1412 fprintf(stdout, "Cleaning up ...\n");
1413 remove_if_exists(opt2fn("-bc", nfile, fnm));
1414 remove_if_exists(opt2fn("-be", nfile, fnm));
1415 remove_if_exists(opt2fn("-bcpo", nfile, fnm));
1416 remove_if_exists(opt2fn("-bg", nfile, fnm));
1417 remove_if_exists(opt2fn("-bo", nfile, fnm));
1418 remove_if_exists(opt2fn("-bx", nfile, fnm));
1420 sfree(command);
1421 sfree(msg );
1424 static void do_the_tests(
1425 FILE *fp, /* General g_tune_pme output file */
1426 char **tpr_names, /* Filenames of the input files to test */
1427 int maxPMEnodes, /* Max fraction of nodes to use for PME */
1428 int minPMEnodes, /* Min fraction of nodes to use for PME */
1429 int npme_fixed, /* If >= -1, test fixed number of PME
1430 * nodes only */
1431 const char *npmevalues_opt, /* Which -npme values should be tested */
1432 t_perf **perfdata, /* Here the performace data is stored */
1433 int *pmeentries, /* Entries in the nPMEnodes list */
1434 int repeats, /* Repeat each test this often */
1435 int nnodes, /* Total number of nodes = nPP + nPME */
1436 int nr_tprs, /* Total number of tpr files to test */
1437 gmx_bool bThreads, /* Threads or MPI? */
1438 char *cmd_mpirun, /* mpirun command string */
1439 char *cmd_np, /* "-np", "-n", whatever mpirun needs */
1440 char *cmd_mdrun, /* mdrun command string */
1441 char *cmd_args_bench, /* arguments for mdrun in a string */
1442 const t_filenm *fnm, /* List of filenames from command line */
1443 int nfile, /* Number of files specified on the cmdl. */
1444 int presteps, /* DLB equilibration steps, is checked */
1445 gmx_int64_t cpt_steps, /* Time step counter in the checkpoint */
1446 gmx_bool bCheck, /* Check whether benchmark mdrun works */
1447 const t_eligible_gpu_ids *gpu_ids) /* Struct containing GPU IDs for
1448 * constructing mdrun command lines */
1450 int i, nr, k, ret, count = 0, totaltests;
1451 int *nPMEnodes = NULL;
1452 t_perf *pd = NULL;
1453 int cmdline_length;
1454 char *command, *cmd_stub;
1455 char buf[STRLEN];
1456 gmx_bool bResetProblem = FALSE;
1457 gmx_bool bFirst = TRUE;
1459 /* This string array corresponds to the eParselog enum type at the start
1460 * of this file */
1461 const char* ParseLog[] = {
1462 "OK.",
1463 "Logfile not found!",
1464 "No timings, logfile truncated?",
1465 "Run was terminated.",
1466 "Counters were not reset properly.",
1467 "No DD grid found for these settings.",
1468 "TPX version conflict!",
1469 "mdrun was not started in parallel!",
1470 "Number of PP ranks has a prime factor that is too large.",
1471 "The number of PP ranks did not suit the number of GPUs.",
1472 "Some GPUs were not detected or are incompatible.",
1473 "An error occurred."
1475 char str_PME_f_load[13];
1478 /* Allocate space for the mdrun command line. 100 extra characters should
1479 be more than enough for the -npme etcetera arguments */
1480 cmdline_length = std::strlen(cmd_mpirun)
1481 + std::strlen(cmd_np)
1482 + std::strlen(cmd_mdrun)
1483 + std::strlen(cmd_args_bench)
1484 + std::strlen(tpr_names[0]) + 100;
1485 snew(command, cmdline_length);
1486 snew(cmd_stub, cmdline_length);
1488 /* Construct the part of the command line that stays the same for all tests: */
1489 if (bThreads)
1491 sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1493 else
1495 sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1498 /* Create a list of numbers of PME nodes to test */
1499 if (npme_fixed < -1)
1501 make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1502 nnodes, minPMEnodes, maxPMEnodes);
1504 else
1506 *pmeentries = 1;
1507 snew(nPMEnodes, 1);
1508 nPMEnodes[0] = npme_fixed;
1509 fprintf(stderr, "Will use a fixed number of %d PME-only ranks.\n", nPMEnodes[0]);
1512 if (0 == repeats)
1514 fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1515 gmx_ffclose(fp);
1516 finalize(opt2fn("-p", nfile, fnm));
1517 exit(0);
1520 /* Allocate one dataset for each tpr input file: */
1521 init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1523 /*****************************************/
1524 /* Main loop over all tpr files to test: */
1525 /*****************************************/
1526 totaltests = nr_tprs*(*pmeentries)*repeats;
1527 for (k = 0; k < nr_tprs; k++)
1529 fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1530 fprintf(fp, "PME ranks Gcycles ns/day PME/f Remark\n");
1531 /* Loop over various numbers of PME nodes: */
1532 for (i = 0; i < *pmeentries; i++)
1534 char *cmd_gpu_ids = NULL;
1536 pd = &perfdata[k][i];
1538 cmd_gpu_ids = make_gpu_id_command_line(nnodes, nPMEnodes[i], gpu_ids);
1540 /* Loop over the repeats for each scenario: */
1541 for (nr = 0; nr < repeats; nr++)
1543 pd->nPMEnodes = nPMEnodes[i];
1545 /* Add -npme and -s to the command line and save it. Note that
1546 * the -passall (if set) options requires cmd_args_bench to be
1547 * at the end of the command line string */
1548 snew(pd->mdrun_cmd_line, cmdline_length);
1549 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s %s",
1550 cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench, cmd_gpu_ids);
1552 /* To prevent that all benchmarks fail due to a show-stopper argument
1553 * on the mdrun command line, we make a quick check first.
1554 * This check can be turned off in cases where the automatically chosen
1555 * number of PME-only ranks leads to a number of PP ranks for which no
1556 * decomposition can be found (e.g. for large prime numbers) */
1557 if (bFirst && bCheck)
1559 /* TODO This check is really for a functional
1560 * .tpr, and if we need it, it should take place
1561 * for every .tpr, and the logic for it should be
1562 * immediately inside the loop over k, not in
1563 * this inner loop. */
1564 char *temporary_cmd_line;
1566 snew(temporary_cmd_line, cmdline_length);
1567 /* TODO -npme 0 is more likely to succeed at low
1568 parallelism than the default of -npme -1, but
1569 is more likely to fail at the scaling limit
1570 when the PP domains may be too small. "mpirun
1571 -np 1 mdrun" is probably a reasonable thing to
1572 do for this check, but it'll be easier to
1573 implement that after some refactoring of how
1574 the number of MPI ranks is managed. */
1575 sprintf(temporary_cmd_line, "%s-npme 0 -nb cpu -s %s %s",
1576 cmd_stub, tpr_names[k], cmd_args_bench);
1577 make_sure_it_runs(temporary_cmd_line, cmdline_length, fp, fnm, nfile);
1579 bFirst = FALSE;
1581 /* Do a benchmark simulation: */
1582 if (repeats > 1)
1584 sprintf(buf, ", pass %d/%d", nr+1, repeats);
1586 else
1588 buf[0] = '\0';
1590 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1591 (100.0*count)/totaltests,
1592 k+1, nr_tprs, i+1, *pmeentries, buf);
1593 make_backup(opt2fn("-err", nfile, fnm));
1594 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1595 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1596 gmx_system_call(command);
1598 /* Collect the performance data from the log file; also check stderr
1599 * for fatal errors */
1600 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1601 pd, nr, presteps, cpt_steps, nnodes);
1602 if ((presteps > 0) && (ret == eParselogResetProblem))
1604 bResetProblem = TRUE;
1607 if (-1 == pd->nPMEnodes)
1609 sprintf(buf, "(%3d)", pd->guessPME);
1611 else
1613 sprintf(buf, " ");
1616 /* Nicer output */
1617 if (pd->PME_f_load[nr] > 0.0)
1619 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1621 else
1623 sprintf(str_PME_f_load, "%s", " - ");
1626 /* Write the data we got to disk */
1627 fprintf(fp, "%4d%s %12.3f %12.3f %s %s", pd->nPMEnodes,
1628 buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1629 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1631 fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1633 fprintf(fp, "\n");
1634 fflush(fp);
1635 count++;
1637 /* Do some cleaning up and delete the files we do not need any more */
1638 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1640 /* If the first run with this number of processors already failed, do not try again: */
1641 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1643 fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1644 count += repeats-(nr+1);
1645 break;
1647 } /* end of repeats loop */
1648 sfree(cmd_gpu_ids);
1649 } /* end of -npme loop */
1650 } /* end of tpr file loop */
1652 if (bResetProblem)
1654 sep_line(fp);
1655 fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1656 sep_line(fp);
1658 sfree(command);
1659 sfree(cmd_stub);
1663 static void check_input(
1664 int nnodes,
1665 int repeats,
1666 int *ntprs,
1667 real *rmin,
1668 real rcoulomb,
1669 real *rmax,
1670 real maxPMEfraction,
1671 real minPMEfraction,
1672 int npme_fixed,
1673 gmx_int64_t bench_nsteps,
1674 const t_filenm *fnm,
1675 int nfile,
1676 int sim_part,
1677 int presteps,
1678 int npargs,
1679 t_pargs *pa)
1681 int old;
1684 /* Make sure the input file exists */
1685 if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1687 gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1690 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1691 if ( (0 == std::strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1693 gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1694 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1697 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1698 if (repeats < 0)
1700 gmx_fatal(FARGS, "Number of repeats < 0!");
1703 /* Check number of nodes */
1704 if (nnodes < 1)
1706 gmx_fatal(FARGS, "Number of ranks/threads must be a positive integer.");
1709 /* Automatically choose -ntpr if not set */
1710 if (*ntprs < 1)
1712 if (nnodes < 16)
1714 *ntprs = 1;
1716 else
1718 *ntprs = 3;
1719 /* Set a reasonable scaling factor for rcoulomb */
1720 if (*rmax <= 0)
1722 *rmax = rcoulomb * 1.2;
1725 fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1727 else
1729 if (1 == *ntprs)
1731 fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1735 /* Make shure that rmin <= rcoulomb <= rmax */
1736 if (*rmin <= 0)
1738 *rmin = rcoulomb;
1740 if (*rmax <= 0)
1742 *rmax = rcoulomb;
1744 if (!(*rmin <= *rmax) )
1746 gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1747 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1749 /* Add test scenarios if rmin or rmax were set */
1750 if (*ntprs <= 2)
1752 if (!gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1754 (*ntprs)++;
1755 fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1756 *rmin, *ntprs);
1758 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1760 (*ntprs)++;
1761 fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1762 *rmax, *ntprs);
1765 old = *ntprs;
1766 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1767 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) || !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1769 *ntprs = std::max(*ntprs, 2);
1772 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1773 if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1775 *ntprs = std::max(*ntprs, 3);
1778 if (old != *ntprs)
1780 fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1783 if (*ntprs > 1)
1785 if (gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && gmx_within_tol(rcoulomb, *rmax, GMX_REAL_EPS)) /* We have just a single rc */
1787 fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1788 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1789 "with correspondingly adjusted PME grid settings\n");
1790 *ntprs = 1;
1794 /* Check whether max and min fraction are within required values */
1795 if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1797 gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1799 if (minPMEfraction > 0.5 || minPMEfraction < 0)
1801 gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1803 if (maxPMEfraction < minPMEfraction)
1805 gmx_fatal(FARGS, "-max must be larger or equal to -min");
1808 /* Check whether the number of steps - if it was set - has a reasonable value */
1809 if (bench_nsteps < 0)
1811 gmx_fatal(FARGS, "Number of steps must be positive.");
1814 if (bench_nsteps > 10000 || bench_nsteps < 100)
1816 fprintf(stderr, "WARNING: steps=");
1817 fprintf(stderr, "%" GMX_PRId64, bench_nsteps);
1818 fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1821 if (presteps < 0)
1823 gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1826 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1827 if (*ntprs > 1)
1829 if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1831 fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1835 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1836 * only. We need to check whether the requested number of PME-only nodes
1837 * makes sense. */
1838 if (npme_fixed > -1)
1840 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1841 if (2*npme_fixed > nnodes)
1843 gmx_fatal(FARGS, "Cannot have more than %d PME-only ranks for a total of %d ranks (you chose %d).\n",
1844 nnodes/2, nnodes, npme_fixed);
1846 if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1848 fprintf(stderr, "WARNING: Only %g percent of the ranks are assigned as PME-only ranks.\n",
1849 (100.0*npme_fixed)/nnodes);
1851 if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1853 fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1854 " fixed number of PME-only ranks is requested with -fix.\n");
1860 /* Returns TRUE when "opt" is needed at launch time */
1861 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1863 if (0 == std::strncmp(opt, "-swap", 5))
1865 return bSet;
1868 /* Apart from the input .tpr and the output log files we need all options that
1869 * were set on the command line and that do not start with -b */
1870 if (0 == std::strncmp(opt, "-b", 2) || 0 == std::strncmp(opt, "-s", 2)
1871 || 0 == std::strncmp(opt, "-err", 4) || 0 == std::strncmp(opt, "-p", 2) )
1873 return FALSE;
1876 return bSet;
1880 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1881 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1883 /* Apart from the input .tpr, all files starting with "-b" are for
1884 * _b_enchmark files exclusively */
1885 if (0 == std::strncmp(opt, "-s", 2))
1887 return FALSE;
1890 if (0 == std::strncmp(opt, "-b", 2) || 0 == std::strncmp(opt, "-s", 2))
1892 if (!bOptional || bSet)
1894 return TRUE;
1896 else
1898 return FALSE;
1901 else
1903 if (bIsOutput)
1905 return FALSE;
1907 else
1909 if (bSet) /* These are additional input files like -cpi -ei */
1911 return TRUE;
1913 else
1915 return FALSE;
1922 /* Adds 'buf' to 'str' */
1923 static void add_to_string(char **str, const char *buf)
1925 int len;
1928 len = std::strlen(*str) + std::strlen(buf) + 1;
1929 srenew(*str, len);
1930 std::strcat(*str, buf);
1934 /* Create the command line for the benchmark as well as for the real run */
1935 static void create_command_line_snippets(
1936 gmx_bool bAppendFiles,
1937 gmx_bool bKeepAndNumCPT,
1938 gmx_bool bResetHWay,
1939 int presteps,
1940 int nfile,
1941 t_filenm fnm[],
1942 char *cmd_args_bench[], /* command line arguments for benchmark runs */
1943 char *cmd_args_launch[], /* command line arguments for simulation run */
1944 char extra_args[], /* Add this to the end of the command line */
1945 char *deffnm) /* Default file names, or NULL if not set */
1947 int i;
1948 char *opt;
1949 const char *name;
1950 char strbuf[STRLEN];
1953 /* strlen needs at least '\0' as a string: */
1954 snew(*cmd_args_bench, 1);
1955 snew(*cmd_args_launch, 1);
1956 *cmd_args_launch[0] = '\0';
1957 *cmd_args_bench[0] = '\0';
1960 /*******************************************/
1961 /* 1. Process other command line arguments */
1962 /*******************************************/
1963 if (presteps > 0)
1965 /* Add equilibration steps to benchmark options */
1966 sprintf(strbuf, "-resetstep %d ", presteps);
1967 add_to_string(cmd_args_bench, strbuf);
1969 /* These switches take effect only at launch time */
1970 if (deffnm)
1972 sprintf(strbuf, "-deffnm %s ", deffnm);
1973 add_to_string(cmd_args_launch, strbuf);
1975 if (FALSE == bAppendFiles)
1977 add_to_string(cmd_args_launch, "-noappend ");
1979 if (bKeepAndNumCPT)
1981 add_to_string(cmd_args_launch, "-cpnum ");
1983 if (bResetHWay)
1985 add_to_string(cmd_args_launch, "-resethway ");
1988 /********************/
1989 /* 2. Process files */
1990 /********************/
1991 for (i = 0; i < nfile; i++)
1993 opt = (char *)fnm[i].opt;
1994 name = opt2fn(opt, nfile, fnm);
1996 /* Strbuf contains the options, now let's sort out where we need that */
1997 sprintf(strbuf, "%s %s ", opt, name);
1999 if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
2001 /* All options starting with -b* need the 'b' removed,
2002 * therefore overwrite strbuf */
2003 if (0 == std::strncmp(opt, "-b", 2))
2005 sprintf(strbuf, "-%s %s ", &opt[2], name);
2008 add_to_string(cmd_args_bench, strbuf);
2011 if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
2013 add_to_string(cmd_args_launch, strbuf);
2017 add_to_string(cmd_args_bench, extra_args);
2018 add_to_string(cmd_args_launch, extra_args);
2022 /* Set option opt */
2023 static void setopt(const char *opt, int nfile, t_filenm fnm[])
2025 int i;
2027 for (i = 0; (i < nfile); i++)
2029 if (std::strcmp(opt, fnm[i].opt) == 0)
2031 fnm[i].flag |= ffSET;
2037 /* This routine inspects the tpr file and ...
2038 * 1. checks for output files that get triggered by a tpr option. These output
2039 * files are marked as 'set' to allow for a proper cleanup after each
2040 * tuning run.
2041 * 2. returns the PME:PP load ratio
2042 * 3. returns rcoulomb from the tpr */
2043 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
2045 gmx_bool bTpi; /* Is test particle insertion requested? */
2046 gmx_bool bFree; /* Is a free energy simulation requested? */
2047 gmx_bool bNM; /* Is a normal mode analysis requested? */
2048 gmx_bool bSwap; /* Is water/ion position swapping requested? */
2049 t_inputrec *ir;
2050 t_state state;
2051 gmx_mtop_t mtop;
2054 /* Check tpr file for options that trigger extra output files */
2055 gmx::MDModules mdModules;
2056 ir = mdModules.inputrec();
2057 read_tpx_state(opt2fn("-s", nfile, fnm), ir, &state, &mtop);
2058 bFree = (efepNO != ir->efep );
2059 bNM = (eiNM == ir->eI );
2060 bSwap = (eswapNO != ir->eSwapCoords);
2061 bTpi = EI_TPI(ir->eI);
2063 /* Set these output files on the tuning command-line */
2064 if (ir->bPull)
2066 setopt("-pf", nfile, fnm);
2067 setopt("-px", nfile, fnm);
2069 if (bFree)
2071 setopt("-dhdl", nfile, fnm);
2073 if (bTpi)
2075 setopt("-tpi", nfile, fnm);
2076 setopt("-tpid", nfile, fnm);
2078 if (bNM)
2080 setopt("-mtx", nfile, fnm);
2082 if (bSwap)
2084 setopt("-swap", nfile, fnm);
2087 *rcoulomb = ir->rcoulomb;
2089 /* Return the estimate for the number of PME nodes */
2090 float npme = pme_load_estimate(&mtop, ir, state.box);
2091 return npme;
2095 static void couple_files_options(int nfile, t_filenm fnm[])
2097 int i;
2098 gmx_bool bSet, bBench;
2099 char *opt;
2100 char buf[20];
2103 for (i = 0; i < nfile; i++)
2105 opt = (char *)fnm[i].opt;
2106 bSet = ((fnm[i].flag & ffSET) != 0);
2107 bBench = (0 == std::strncmp(opt, "-b", 2));
2109 /* Check optional files */
2110 /* If e.g. -eo is set, then -beo also needs to be set */
2111 if (is_optional(&fnm[i]) && bSet && !bBench)
2113 sprintf(buf, "-b%s", &opt[1]);
2114 setopt(buf, nfile, fnm);
2116 /* If -beo is set, then -eo also needs to be! */
2117 if (is_optional(&fnm[i]) && bSet && bBench)
2119 sprintf(buf, "-%s", &opt[2]);
2120 setopt(buf, nfile, fnm);
2126 #define BENCHSTEPS (1000)
2128 int gmx_tune_pme(int argc, char *argv[])
2130 const char *desc[] = {
2131 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of ranks, [THISMODULE] systematically",
2132 "times [gmx-mdrun] with various numbers of PME-only ranks and determines",
2133 "which setting is fastest. It will also test whether performance can",
2134 "be enhanced by shifting load from the reciprocal to the real space",
2135 "part of the Ewald sum. ",
2136 "Simply pass your [REF].tpr[ref] file to [THISMODULE] together with other options",
2137 "for [gmx-mdrun] as needed.[PAR]",
2138 "[THISMODULE] needs to call [gmx-mdrun] and so requires that you",
2139 "specify how to call mdrun with the argument to the [TT]-mdrun[tt]",
2140 "parameter. Depending how you have built GROMACS, values such as",
2141 "'gmx mdrun', 'gmx_d mdrun', or 'mdrun_mpi' might be needed.[PAR]",
2142 "The program that runs MPI programs can be set in the environment variable",
2143 "MPIRUN (defaults to 'mpirun'). Note that for certain MPI frameworks,",
2144 "you need to provide a machine- or hostfile. This can also be passed",
2145 "via the MPIRUN variable, e.g.[PAR]",
2146 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt]",
2147 "Note that in such cases it is normally necessary to compile",
2148 "and/or run [THISMODULE] without MPI support, so that it can call",
2149 "the MPIRUN program.[PAR]",
2150 "Before doing the actual benchmark runs, [THISMODULE] will do a quick",
2151 "check whether [gmx-mdrun] works as expected with the provided parallel settings",
2152 "if the [TT]-check[tt] option is activated (the default).",
2153 "Please call [THISMODULE] with the normal options you would pass to",
2154 "[gmx-mdrun] and add [TT]-np[tt] for the number of ranks to perform the",
2155 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2156 "to repeat each test several times to get better statistics. [PAR]",
2157 "[THISMODULE] can test various real space / reciprocal space workloads",
2158 "for you. With [TT]-ntpr[tt] you control how many extra [REF].tpr[ref] files will be",
2159 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2160 "Typically, the first test (number 0) will be with the settings from the input",
2161 "[REF].tpr[ref] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2162 "specified by [TT]-rmax[tt] with a somewhat smaller PME grid at the same time. ",
2163 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2164 "The remaining [REF].tpr[ref] files will have equally-spaced Coulomb radii (and Fourier "
2165 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2166 "if you just seek the optimal number of PME-only ranks; in that case",
2167 "your input [REF].tpr[ref] file will remain unchanged.[PAR]",
2168 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2169 "MD systems. The dynamic load balancing needs about 100 time steps",
2170 "to adapt to local load imbalances, therefore the time step counters",
2171 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2172 "for a higher accuracy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2173 "From the 'DD' load imbalance entries in the md.log output file you",
2174 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]"
2175 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2176 "After calling [gmx-mdrun] several times, detailed performance information",
2177 "is available in the output file [TT]perf.out[tt].",
2178 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2179 "(options [TT]-b*[tt]), these will be automatically deleted after each test.[PAR]",
2180 "If you want the simulation to be started automatically with the",
2181 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2182 "Basic support for GPU-enabled [TT]mdrun[tt] exists. Give a string containing the IDs",
2183 "of the GPUs that you wish to use in the optimization in the [TT]-gpu_id[tt]",
2184 "command-line argument. Unlike [TT]mdrun -gpu_id[tt], this does not imply a mapping",
2185 "but merely the eligible set. [TT]g_tune_pme[tt] will construct calls to",
2186 "mdrun that use this set appropriately, assuming that PP ranks with low indices",
2187 "should map to GPUs with low indices, and increasing both monotonically",
2188 "over the respective sets.[PAR]",
2191 int nnodes = 1;
2192 int repeats = 2;
2193 int pmeentries = 0; /* How many values for -npme do we actually test for each tpr file */
2194 real maxPMEfraction = 0.50;
2195 real minPMEfraction = 0.25;
2196 int maxPMEnodes, minPMEnodes;
2197 float guessPMEratio; /* guessed PME:PP ratio based on the tpr file */
2198 float guessPMEnodes;
2199 int npme_fixed = -2; /* If >= -1, use only this number
2200 * of PME-only nodes */
2201 int ntprs = 0;
2202 real rmin = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2203 real rcoulomb = -1.0; /* Coulomb radius as set in .tpr file */
2204 gmx_bool bScaleRvdw = TRUE;
2205 gmx_int64_t bench_nsteps = BENCHSTEPS;
2206 gmx_int64_t new_sim_nsteps = -1; /* -1 indicates: not set by the user */
2207 gmx_int64_t cpt_steps = 0; /* Step counter in .cpt input file */
2208 int presteps = 100; /* Do a full cycle reset after presteps steps */
2209 gmx_bool bOverwrite = FALSE, bKeepTPR;
2210 gmx_bool bLaunch = FALSE;
2211 char *ExtraArgs = NULL;
2212 char **tpr_names = NULL;
2213 const char *simulation_tpr = NULL;
2214 char *deffnm = NULL;
2215 int best_npme, best_tpr;
2216 int sim_part = 1; /* For benchmarks with checkpoint files */
2217 char bbuf[STRLEN];
2220 /* Default program names if nothing else is found */
2221 char *cmd_mpirun = NULL, *cmd_mdrun = NULL;
2222 char *cmd_args_bench, *cmd_args_launch;
2223 char *cmd_np = NULL;
2225 /* IDs of GPUs that are eligible for computation */
2226 char *eligible_gpu_ids = NULL;
2227 t_eligible_gpu_ids *gpu_ids = NULL;
2229 t_perf **perfdata = NULL;
2230 t_inputinfo *info;
2231 int i;
2232 FILE *fp;
2234 /* Print out how long the tuning took */
2235 double seconds;
2237 static t_filenm fnm[] = {
2238 /* g_tune_pme */
2239 { efOUT, "-p", "perf", ffWRITE },
2240 { efLOG, "-err", "bencherr", ffWRITE },
2241 { efTPR, "-so", "tuned", ffWRITE },
2242 /* mdrun: */
2243 { efTPR, NULL, NULL, ffREAD },
2244 { efTRN, "-o", NULL, ffWRITE },
2245 { efCOMPRESSED, "-x", NULL, ffOPTWR },
2246 { efCPT, "-cpi", NULL, ffOPTRD },
2247 { efCPT, "-cpo", NULL, ffOPTWR },
2248 { efSTO, "-c", "confout", ffWRITE },
2249 { efEDR, "-e", "ener", ffWRITE },
2250 { efLOG, "-g", "md", ffWRITE },
2251 { efXVG, "-dhdl", "dhdl", ffOPTWR },
2252 { efXVG, "-field", "field", ffOPTWR },
2253 { efXVG, "-table", "table", ffOPTRD },
2254 { efXVG, "-tablep", "tablep", ffOPTRD },
2255 { efXVG, "-tableb", "table", ffOPTRD },
2256 { efTRX, "-rerun", "rerun", ffOPTRD },
2257 { efXVG, "-tpi", "tpi", ffOPTWR },
2258 { efXVG, "-tpid", "tpidist", ffOPTWR },
2259 { efEDI, "-ei", "sam", ffOPTRD },
2260 { efXVG, "-eo", "edsam", ffOPTWR },
2261 { efXVG, "-devout", "deviatie", ffOPTWR },
2262 { efXVG, "-runav", "runaver", ffOPTWR },
2263 { efXVG, "-px", "pullx", ffOPTWR },
2264 { efXVG, "-pf", "pullf", ffOPTWR },
2265 { efXVG, "-ro", "rotation", ffOPTWR },
2266 { efLOG, "-ra", "rotangles", ffOPTWR },
2267 { efLOG, "-rs", "rotslabs", ffOPTWR },
2268 { efLOG, "-rt", "rottorque", ffOPTWR },
2269 { efMTX, "-mtx", "nm", ffOPTWR },
2270 { efXVG, "-swap", "swapions", ffOPTWR },
2271 /* Output files that are deleted after each benchmark run */
2272 { efTRN, "-bo", "bench", ffWRITE },
2273 { efXTC, "-bx", "bench", ffWRITE },
2274 { efCPT, "-bcpo", "bench", ffWRITE },
2275 { efSTO, "-bc", "bench", ffWRITE },
2276 { efEDR, "-be", "bench", ffWRITE },
2277 { efLOG, "-bg", "bench", ffWRITE },
2278 { efXVG, "-beo", "benchedo", ffOPTWR },
2279 { efXVG, "-bdhdl", "benchdhdl", ffOPTWR },
2280 { efXVG, "-bfield", "benchfld", ffOPTWR },
2281 { efXVG, "-btpi", "benchtpi", ffOPTWR },
2282 { efXVG, "-btpid", "benchtpid", ffOPTWR },
2283 { efXVG, "-bdevout", "benchdev", ffOPTWR },
2284 { efXVG, "-brunav", "benchrnav", ffOPTWR },
2285 { efXVG, "-bpx", "benchpx", ffOPTWR },
2286 { efXVG, "-bpf", "benchpf", ffOPTWR },
2287 { efXVG, "-bro", "benchrot", ffOPTWR },
2288 { efLOG, "-bra", "benchrota", ffOPTWR },
2289 { efLOG, "-brs", "benchrots", ffOPTWR },
2290 { efLOG, "-brt", "benchrott", ffOPTWR },
2291 { efMTX, "-bmtx", "benchn", ffOPTWR },
2292 { efNDX, "-bdn", "bench", ffOPTWR },
2293 { efXVG, "-bswap", "benchswp", ffOPTWR }
2296 gmx_bool bThreads = FALSE;
2298 int nthreads = 1;
2300 const char *procstring[] =
2301 { NULL, "np", "n", "none", NULL };
2302 const char *npmevalues_opt[] =
2303 { NULL, "auto", "all", "subset", NULL };
2305 gmx_bool bAppendFiles = TRUE;
2306 gmx_bool bKeepAndNumCPT = FALSE;
2307 gmx_bool bResetCountersHalfWay = FALSE;
2308 gmx_bool bBenchmark = TRUE;
2309 gmx_bool bCheck = TRUE;
2311 gmx_output_env_t *oenv = NULL;
2313 t_pargs pa[] = {
2314 /***********************/
2315 /* g_tune_pme options: */
2316 /***********************/
2317 { "-mdrun", FALSE, etSTR, {&cmd_mdrun},
2318 "Command line to run a simulation, e.g. 'gmx mdrun' or 'mdrun_mpi'" },
2319 { "-np", FALSE, etINT, {&nnodes},
2320 "Number of ranks to run the tests on (must be > 2 for separate PME ranks)" },
2321 { "-npstring", FALSE, etENUM, {procstring},
2322 "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)" },
2323 { "-ntmpi", FALSE, etINT, {&nthreads},
2324 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2325 { "-r", FALSE, etINT, {&repeats},
2326 "Repeat each test this often" },
2327 { "-max", FALSE, etREAL, {&maxPMEfraction},
2328 "Max fraction of PME ranks to test with" },
2329 { "-min", FALSE, etREAL, {&minPMEfraction},
2330 "Min fraction of PME ranks to test with" },
2331 { "-npme", FALSE, etENUM, {npmevalues_opt},
2332 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2333 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2334 { "-fix", FALSE, etINT, {&npme_fixed},
2335 "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."},
2336 { "-rmax", FALSE, etREAL, {&rmax},
2337 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2338 { "-rmin", FALSE, etREAL, {&rmin},
2339 "If >0, minimal rcoulomb for -ntpr>1" },
2340 { "-scalevdw", FALSE, etBOOL, {&bScaleRvdw},
2341 "Scale rvdw along with rcoulomb"},
2342 { "-ntpr", FALSE, etINT, {&ntprs},
2343 "Number of [REF].tpr[ref] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2344 "If < 1, automatically choose the number of [REF].tpr[ref] files to test" },
2345 { "-steps", FALSE, etINT64, {&bench_nsteps},
2346 "Take timings for this many steps in the benchmark runs" },
2347 { "-resetstep", FALSE, etINT, {&presteps},
2348 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2349 { "-nsteps", FALSE, etINT64, {&new_sim_nsteps},
2350 "If non-negative, perform this many steps in the real run (overwrites nsteps from [REF].tpr[ref], add [REF].cpt[ref] steps)" },
2351 { "-launch", FALSE, etBOOL, {&bLaunch},
2352 "Launch the real simulation after optimization" },
2353 { "-bench", FALSE, etBOOL, {&bBenchmark},
2354 "Run the benchmarks or just create the input [REF].tpr[ref] files?" },
2355 { "-check", FALSE, etBOOL, {&bCheck},
2356 "Before the benchmark runs, check whether mdrun works in parallel" },
2357 { "-gpu_id", FALSE, etSTR, {&eligible_gpu_ids},
2358 "List of GPU device id-s that are eligible for use (unlike mdrun, does not imply any mapping)" },
2359 /******************/
2360 /* mdrun options: */
2361 /******************/
2362 /* We let g_tune_pme parse and understand these options, because we need to
2363 * prevent that they appear on the mdrun command line for the benchmarks */
2364 { "-append", FALSE, etBOOL, {&bAppendFiles},
2365 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2366 { "-cpnum", FALSE, etBOOL, {&bKeepAndNumCPT},
2367 "Keep and number checkpoint files (launch only)" },
2368 { "-deffnm", FALSE, etSTR, {&deffnm},
2369 "Set the default filenames (launch only)" },
2370 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2371 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2374 #define NFILE asize(fnm)
2376 seconds = gmx_gettime();
2378 if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2379 NFILE, fnm, asize(pa), pa, asize(desc), desc,
2380 0, NULL, &oenv))
2382 return 0;
2385 // procstring[0] is used inside two different conditionals further down
2386 GMX_RELEASE_ASSERT(procstring[0] != NULL, "Options inconsistency; procstring[0] is NULL");
2388 /* Store the remaining unparsed command line entries in a string which
2389 * is then attached to the mdrun command line */
2390 snew(ExtraArgs, 1);
2391 ExtraArgs[0] = '\0';
2392 for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2394 add_to_string(&ExtraArgs, argv[i]);
2395 add_to_string(&ExtraArgs, " ");
2398 if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2400 bThreads = TRUE;
2401 if (opt2parg_bSet("-npstring", asize(pa), pa))
2403 fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2406 if (nnodes > 1)
2408 gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2410 /* and now we just set this; a bit of an ugly hack*/
2411 nnodes = nthreads;
2413 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2414 guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2416 /* Automatically set -beo options if -eo is set etc. */
2417 couple_files_options(NFILE, fnm);
2419 /* Construct the command line arguments for benchmark runs
2420 * as well as for the simulation run */
2421 if (bThreads)
2423 sprintf(bbuf, " -ntmpi %d ", nthreads);
2425 else
2427 /* This string will be used for MPI runs and will appear after the
2428 * mpirun command. */
2429 if (std::strcmp(procstring[0], "none") != 0)
2431 sprintf(bbuf, " -%s %d ", procstring[0], nnodes);
2433 else
2435 sprintf(bbuf, " ");
2439 cmd_np = bbuf;
2441 create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2442 NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs, deffnm);
2444 /* Prepare to use checkpoint file if requested */
2445 sim_part = 1;
2446 if (opt2bSet("-cpi", NFILE, fnm))
2448 const char *filename = opt2fn("-cpi", NFILE, fnm);
2449 int cpt_sim_part;
2450 read_checkpoint_part_and_step(filename,
2451 &cpt_sim_part, &cpt_steps);
2452 if (cpt_sim_part == 0)
2454 gmx_fatal(FARGS, "Checkpoint file %s could not be read!", filename);
2456 /* tune_pme will run the next part of the simulation */
2457 sim_part = cpt_sim_part + 1;
2460 /* Open performance output file and write header info */
2461 fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2463 /* Make a quick consistency check of command line parameters */
2464 check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2465 maxPMEfraction, minPMEfraction, npme_fixed,
2466 bench_nsteps, fnm, NFILE, sim_part, presteps,
2467 asize(pa), pa);
2468 /* Check any GPU IDs passed make sense, and fill the data structure for them */
2469 snew(gpu_ids, 1);
2470 parse_digits_from_string(eligible_gpu_ids, &gpu_ids->n, &gpu_ids->ids);
2472 /* Determine the maximum and minimum number of PME nodes to test,
2473 * the actual list of settings is build in do_the_tests(). */
2474 if ((nnodes > 2) && (npme_fixed < -1))
2476 if (0 == std::strcmp(npmevalues_opt[0], "auto"))
2478 /* Determine the npme range automatically based on the PME:PP load guess */
2479 if (guessPMEratio > 1.0)
2481 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2482 maxPMEnodes = nnodes/2;
2483 minPMEnodes = nnodes/2;
2485 else
2487 /* PME : PP load is in the range 0..1, let's test around the guess */
2488 guessPMEnodes = static_cast<int>(nnodes/(1.0 + 1.0/guessPMEratio));
2489 minPMEnodes = static_cast<int>(std::floor(0.7*guessPMEnodes));
2490 maxPMEnodes = static_cast<int>(std::ceil(1.6*guessPMEnodes));
2491 maxPMEnodes = std::min(maxPMEnodes, nnodes/2);
2494 else
2496 /* Determine the npme range based on user input */
2497 maxPMEnodes = static_cast<int>(std::floor(maxPMEfraction*nnodes));
2498 minPMEnodes = std::max(static_cast<int>(std::floor(minPMEfraction*nnodes)), 0);
2499 fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2500 if (maxPMEnodes != minPMEnodes)
2502 fprintf(stdout, "- %d ", maxPMEnodes);
2504 fprintf(stdout, "PME-only ranks.\n Note that the automatic number of PME-only ranks and no separate PME ranks are always tested.\n");
2507 else
2509 maxPMEnodes = 0;
2510 minPMEnodes = 0;
2513 /* Get the commands we need to set up the runs from environment variables */
2514 get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2515 if (bBenchmark && repeats > 0)
2517 check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun, NULL != eligible_gpu_ids);
2520 /* Print some header info to file */
2521 sep_line(fp);
2522 fprintf(fp, "\n P E R F O R M A N C E R E S U L T S\n");
2523 sep_line(fp);
2524 fprintf(fp, "%s for GROMACS %s\n", output_env_get_program_display_name(oenv),
2525 gmx_version());
2526 if (!bThreads)
2528 fprintf(fp, "Number of ranks : %d\n", nnodes);
2529 fprintf(fp, "The mpirun command is : %s\n", cmd_mpirun);
2530 if (std::strcmp(procstring[0], "none") != 0)
2532 fprintf(fp, "Passing # of ranks via : -%s\n", procstring[0]);
2534 else
2536 fprintf(fp, "Not setting number of ranks in system call\n");
2539 else
2541 fprintf(fp, "Number of threads : %d\n", nnodes);
2544 fprintf(fp, "The mdrun command is : %s\n", cmd_mdrun);
2545 fprintf(fp, "mdrun args benchmarks : %s\n", cmd_args_bench);
2546 fprintf(fp, "Benchmark steps : ");
2547 fprintf(fp, "%" GMX_PRId64, bench_nsteps);
2548 fprintf(fp, "\n");
2549 fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2550 if (sim_part > 1)
2552 fprintf(fp, "Checkpoint time step : ");
2553 fprintf(fp, "%" GMX_PRId64, cpt_steps);
2554 fprintf(fp, "\n");
2556 fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2558 if (new_sim_nsteps >= 0)
2560 bOverwrite = TRUE;
2561 fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2562 fprintf(stderr, "%" GMX_PRId64, new_sim_nsteps+cpt_steps);
2563 fprintf(stderr, " steps.\n");
2564 fprintf(fp, "Simulation steps : ");
2565 fprintf(fp, "%" GMX_PRId64, new_sim_nsteps);
2566 fprintf(fp, "\n");
2568 if (repeats > 1)
2570 fprintf(fp, "Repeats for each test : %d\n", repeats);
2573 if (npme_fixed >= -1)
2575 fprintf(fp, "Fixing -npme at : %d\n", npme_fixed);
2578 fprintf(fp, "Input file : %s\n", opt2fn("-s", NFILE, fnm));
2579 fprintf(fp, " PME/PP load estimate : %g\n", guessPMEratio);
2581 /* Allocate memory for the inputinfo struct: */
2582 snew(info, 1);
2583 info->nr_inputfiles = ntprs;
2584 for (i = 0; i < ntprs; i++)
2586 snew(info->rcoulomb, ntprs);
2587 snew(info->rvdw, ntprs);
2588 snew(info->rlist, ntprs);
2589 snew(info->nkx, ntprs);
2590 snew(info->nky, ntprs);
2591 snew(info->nkz, ntprs);
2592 snew(info->fsx, ntprs);
2593 snew(info->fsy, ntprs);
2594 snew(info->fsz, ntprs);
2596 /* Make alternative tpr files to test: */
2597 snew(tpr_names, ntprs);
2598 for (i = 0; i < ntprs; i++)
2600 snew(tpr_names[i], STRLEN);
2603 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2604 * different grids could be found. */
2605 make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2606 cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2608 /********************************************************************************/
2609 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2610 /********************************************************************************/
2611 snew(perfdata, ntprs);
2612 if (bBenchmark)
2614 GMX_RELEASE_ASSERT(npmevalues_opt[0] != NULL, "Options inconsistency; npmevalues_opt[0] is NULL");
2615 do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2616 repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2617 cmd_args_bench, fnm, NFILE, presteps, cpt_steps, bCheck, gpu_ids);
2619 fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2621 /* Analyse the results and give a suggestion for optimal settings: */
2622 bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2623 repeats, info, &best_tpr, &best_npme);
2625 /* Take the best-performing tpr file and enlarge nsteps to original value */
2626 if (bKeepTPR && !bOverwrite)
2628 simulation_tpr = opt2fn("-s", NFILE, fnm);
2630 else
2632 simulation_tpr = opt2fn("-so", NFILE, fnm);
2633 modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2634 info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2637 /* Let's get rid of the temporary benchmark input files */
2638 for (i = 0; i < ntprs; i++)
2640 fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2641 remove(tpr_names[i]);
2644 /* Now start the real simulation if the user requested it ... */
2645 launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2646 cmd_args_launch, simulation_tpr, nnodes, best_npme, gpu_ids);
2648 gmx_ffclose(fp);
2650 /* ... or simply print the performance results to screen: */
2651 if (!bLaunch)
2653 finalize(opt2fn("-p", NFILE, fnm));
2656 return 0;