2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018,2019, 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.
47 #ifdef HAVE_SYS_TIME_H
51 #include "gromacs/commandline/pargs.h"
52 #include "gromacs/ewald/pme.h"
53 #include "gromacs/fft/calcgrid.h"
54 #include "gromacs/fileio/checkpoint.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/gmxana/gmx_ana.h"
57 #include "gromacs/math/utilities.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/mdlib/perf_est.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/state.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/taskassignment/usergpuids.h"
66 #include "gromacs/timing/walltime_accounting.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/utility/arraysize.h"
69 #include "gromacs/utility/baseversion.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/futil.h"
73 #include "gromacs/utility/gmxassert.h"
74 #include "gromacs/utility/path.h"
75 #include "gromacs/utility/smalloc.h"
76 #include "gromacs/utility/stringutil.h"
78 /* Enum for situations that can occur during log file parsing, the
79 * corresponding string entries can be found in do_the_tests() in
80 * const char* ParseLog[] */
81 /* TODO clean up CamelCasing of these enum names */
87 eParselogResetProblem
,
91 eParselogLargePrimeFactor
,
92 eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs
,
101 int nPMEnodes
; /* number of PME-only nodes used in this test */
102 int nx
, ny
, nz
; /* DD grid */
103 int guessPME
; /* if nPMEnodes == -1, this is the guessed number of PME nodes */
104 double *Gcycles
; /* This can contain more than one value if doing multiple tests */
108 float *PME_f_load
; /* PME mesh/force load average*/
109 float PME_f_load_Av
; /* Average average ;) ... */
110 char *mdrun_cmd_line
; /* Mdrun command line used for this test */
116 int nr_inputfiles
; /* The number of tpr and mdp input files */
117 int64_t orig_sim_steps
; /* Number of steps to be done in the real simulation */
118 int64_t orig_init_step
; /* Init step for the real simulation */
119 real
*rcoulomb
; /* The coulomb radii [0...nr_inputfiles] */
120 real
*rvdw
; /* The vdW radii */
121 real
*rlist
; /* Neighbourlist cutoff radius */
122 int *nkx
, *nky
, *nkz
;
123 real
*fsx
, *fsy
, *fsz
; /* Fourierspacing in x,y,z dimension */
127 static void sep_line(FILE *fp
)
129 fprintf(fp
, "\n------------------------------------------------------------\n");
133 /* Wrapper for system calls */
134 static int gmx_system_call(char *command
)
136 return ( system(command
) );
140 /* Check if string starts with substring */
141 static gmx_bool
str_starts(const char *string
, const char *substring
)
143 return ( std::strncmp(string
, substring
, std::strlen(substring
)) == 0);
147 static void cleandata(t_perf
*perfdata
, int test_nr
)
149 perfdata
->Gcycles
[test_nr
] = 0.0;
150 perfdata
->ns_per_day
[test_nr
] = 0.0;
151 perfdata
->PME_f_load
[test_nr
] = 0.0;
155 static void remove_if_exists(const char *fn
)
159 fprintf(stdout
, "Deleting %s\n", fn
);
165 static void finalize(const char *fn_out
)
171 fp
= fopen(fn_out
, "r");
172 fprintf(stdout
, "\n\n");
174 while (fgets(buf
, STRLEN
-1, fp
) != nullptr)
176 fprintf(stdout
, "%s", buf
);
179 fprintf(stdout
, "\n\n");
184 eFoundNothing
, eFoundDDStr
, eFoundAccountingStr
, eFoundCycleStr
187 static int parse_logfile(const char *logfile
, const char *errfile
,
188 t_perf
*perfdata
, int test_nr
, int presteps
, int64_t cpt_steps
,
192 char line
[STRLEN
], dumstring
[STRLEN
], dumstring2
[STRLEN
];
193 const char matchstrdd
[] = "Domain decomposition grid";
194 const char matchstrcr
[] = "resetting all time and cycle counters";
195 const char matchstrbal
[] = "Average PME mesh/force load:";
196 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";
197 const char errSIG
[] = "signal, stopping at the next";
199 float dum1
, dum2
, dum3
, dum4
;
202 int64_t resetsteps
= -1;
203 gmx_bool bFoundResetStr
= FALSE
;
204 gmx_bool bResetChecked
= FALSE
;
207 if (!gmx_fexist(logfile
))
209 fprintf(stderr
, "WARNING: Could not find logfile %s.\n", logfile
);
210 cleandata(perfdata
, test_nr
);
211 return eParselogNotFound
;
214 fp
= fopen(logfile
, "r");
215 perfdata
->PME_f_load
[test_nr
] = -1.0;
216 perfdata
->guessPME
= -1;
218 iFound
= eFoundNothing
;
221 iFound
= eFoundDDStr
; /* Skip some case statements */
224 while (fgets(line
, STRLEN
, fp
) != nullptr)
226 /* Remove leading spaces */
229 /* Check for TERM and INT signals from user: */
230 if (std::strstr(line
, errSIG
) != nullptr)
233 cleandata(perfdata
, test_nr
);
234 return eParselogTerm
;
237 /* Check whether cycle resetting worked */
238 if (presteps
> 0 && !bFoundResetStr
)
240 if (std::strstr(line
, matchstrcr
) != nullptr)
242 sprintf(dumstring
, "step %s", "%" SCNd64
);
243 sscanf(line
, dumstring
, &resetsteps
);
244 bFoundResetStr
= TRUE
;
245 if (resetsteps
== presteps
+cpt_steps
)
247 bResetChecked
= TRUE
;
251 sprintf(dumstring
, "%" PRId64
, resetsteps
);
252 sprintf(dumstring2
, "%" PRId64
, presteps
+cpt_steps
);
253 fprintf(stderr
, "WARNING: Time step counters were reset at step %s,\n"
254 " though they were supposed to be reset at step %s!\n",
255 dumstring
, dumstring2
);
260 /* Look for strings that appear in a certain order in the log file: */
264 /* Look for domain decomp grid and separate PME nodes: */
265 if (str_starts(line
, matchstrdd
))
267 sscanf(line
, "Domain decomposition grid %d x %d x %d, separate PME ranks %d",
268 &(perfdata
->nx
), &(perfdata
->ny
), &(perfdata
->nz
), &npme
);
269 if (perfdata
->nPMEnodes
== -1)
271 perfdata
->guessPME
= npme
;
273 else if (perfdata
->nPMEnodes
!= npme
)
275 gmx_fatal(FARGS
, "PME ranks from command line and output file are not identical");
277 iFound
= eFoundDDStr
;
279 /* Catch a few errors that might have occurred: */
280 else if (str_starts(line
, "There is no domain decomposition for"))
283 return eParselogNoDDGrid
;
285 else if (str_starts(line
, "The number of ranks you selected"))
288 return eParselogLargePrimeFactor
;
290 else if (str_starts(line
, "reading tpx file"))
293 return eParselogTPXVersion
;
295 else if (str_starts(line
, "The -dd or -npme option request a parallel simulation"))
298 return eParselogNotParallel
;
302 /* Even after the "Domain decomposition grid" string was found,
303 * it could be that mdrun had to quit due to some error. */
304 if (str_starts(line
, "Incorrect launch configuration: mismatching number of"))
307 return eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs
;
309 else if (str_starts(line
, "Some of the requested GPUs do not exist"))
312 return eParselogGpuProblem
;
314 /* Look for PME mesh/force balance (not necessarily present, though) */
315 else if (str_starts(line
, matchstrbal
))
317 sscanf(&line
[std::strlen(matchstrbal
)], "%f", &(perfdata
->PME_f_load
[test_nr
]));
319 /* Look for matchstring */
320 else if (str_starts(line
, matchstring
))
322 iFound
= eFoundAccountingStr
;
325 case eFoundAccountingStr
:
326 /* Already found matchstring - look for cycle data */
327 if (str_starts(line
, "Total "))
329 sscanf(line
, "Total %*f %lf", &(perfdata
->Gcycles
[test_nr
]));
330 iFound
= eFoundCycleStr
;
334 /* Already found cycle data - look for remaining performance info and return */
335 if (str_starts(line
, "Performance:"))
337 ndum
= sscanf(line
, "%s %f %f %f %f", dumstring
, &dum1
, &dum2
, &dum3
, &dum4
);
338 /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
339 perfdata
->ns_per_day
[test_nr
] = (ndum
== 5) ? dum3
: dum1
;
341 if (bResetChecked
|| presteps
== 0)
347 return eParselogResetProblem
;
354 /* Close the log file */
357 /* Check why there is no performance data in the log file.
358 * Did a fatal errors occur? */
359 if (gmx_fexist(errfile
))
361 fp
= fopen(errfile
, "r");
362 while (fgets(line
, STRLEN
, fp
) != nullptr)
364 if (str_starts(line
, "Fatal error:") )
366 if (fgets(line
, STRLEN
, fp
) != nullptr)
368 fprintf(stderr
, "\nWARNING: An error occurred during this benchmark:\n"
372 cleandata(perfdata
, test_nr
);
373 return eParselogFatal
;
380 fprintf(stderr
, "WARNING: Could not find stderr file %s.\n", errfile
);
383 /* Giving up ... we could not find out why there is no performance data in
385 fprintf(stdout
, "No performance data in log file.\n");
386 cleandata(perfdata
, test_nr
);
388 return eParselogNoPerfData
;
392 static gmx_bool
analyze_data(
401 int *index_tpr
, /* OUT: Nr of mdp file with best settings */
402 int *npme_optimal
) /* OUT: Optimal number of PME nodes */
405 int line
= 0, line_win
= -1;
406 int k_win
= -1, i_win
= -1, winPME
;
407 double s
= 0.0; /* standard deviation */
410 char str_PME_f_load
[13];
411 gmx_bool bCanUseOrigTPR
;
412 gmx_bool bRefinedCoul
, bRefinedVdW
, bRefinedGrid
;
418 fprintf(fp
, "Summary of successful runs:\n");
419 fprintf(fp
, "Line tpr PME ranks Gcycles Av. Std.dev. ns/day PME/f");
422 fprintf(fp
, " DD grid");
428 for (k
= 0; k
< ntprs
; k
++)
430 for (i
= 0; i
< ntests
; i
++)
432 /* Select the right dataset: */
433 pd
= &(perfdata
[k
][i
]);
435 pd
->Gcycles_Av
= 0.0;
436 pd
->PME_f_load_Av
= 0.0;
437 pd
->ns_per_day_Av
= 0.0;
439 if (pd
->nPMEnodes
== -1)
441 sprintf(strbuf
, "(%3d)", pd
->guessPME
);
445 sprintf(strbuf
, " ");
448 /* Get the average run time of a setting */
449 for (j
= 0; j
< nrepeats
; j
++)
451 pd
->Gcycles_Av
+= pd
->Gcycles
[j
];
452 pd
->PME_f_load_Av
+= pd
->PME_f_load
[j
];
454 pd
->Gcycles_Av
/= nrepeats
;
455 pd
->PME_f_load_Av
/= nrepeats
;
457 for (j
= 0; j
< nrepeats
; j
++)
459 if (pd
->ns_per_day
[j
] > 0.0)
461 pd
->ns_per_day_Av
+= pd
->ns_per_day
[j
];
465 /* Somehow the performance number was not aquired for this run,
466 * therefor set the average to some negative value: */
467 pd
->ns_per_day_Av
= -1.0f
*nrepeats
;
471 pd
->ns_per_day_Av
/= nrepeats
;
474 if (pd
->PME_f_load_Av
> 0.0)
476 sprintf(str_PME_f_load
, "%12.3f", pd
->PME_f_load_Av
);
480 sprintf(str_PME_f_load
, "%s", " - ");
484 /* We assume we had a successful run if both averages are positive */
485 if (pd
->Gcycles_Av
> 0.0 && pd
->ns_per_day_Av
> 0.0)
487 /* Output statistics if repeats were done */
490 /* Calculate the standard deviation */
492 for (j
= 0; j
< nrepeats
; j
++)
494 s
+= gmx::square( pd
->Gcycles
[j
] - pd
->Gcycles_Av
);
499 fprintf(fp
, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
500 line
, k
, pd
->nPMEnodes
, strbuf
, pd
->Gcycles_Av
, s
,
501 pd
->ns_per_day_Av
, str_PME_f_load
);
504 fprintf(fp
, " %3d %3d %3d", pd
->nx
, pd
->ny
, pd
->nz
);
508 /* Store the index of the best run found so far in 'winner': */
509 if ( (k_win
== -1) || (pd
->Gcycles_Av
< perfdata
[k_win
][i_win
].Gcycles_Av
) )
522 gmx_fatal(FARGS
, "None of the runs was successful! Check %s for problems.", fn
);
527 winPME
= perfdata
[k_win
][i_win
].nPMEnodes
;
531 /* We stuck to a fixed number of PME-only nodes */
532 sprintf(strbuf
, "settings No. %d", k_win
);
536 /* We have optimized the number of PME-only nodes */
539 sprintf(strbuf
, "%s", "the automatic number of PME ranks");
543 sprintf(strbuf
, "%d PME ranks", winPME
);
546 fprintf(fp
, "Best performance was achieved with %s", strbuf
);
547 if ((nrepeats
> 1) && (ntests
> 1))
549 fprintf(fp
, " (see line %d)", line_win
);
553 /* Only mention settings if they were modified: */
554 bRefinedCoul
= !gmx_within_tol(info
->rcoulomb
[k_win
], info
->rcoulomb
[0], GMX_REAL_EPS
);
555 bRefinedVdW
= !gmx_within_tol(info
->rvdw
[k_win
], info
->rvdw
[0], GMX_REAL_EPS
);
556 bRefinedGrid
= !(info
->nkx
[k_win
] == info
->nkx
[0] &&
557 info
->nky
[k_win
] == info
->nky
[0] &&
558 info
->nkz
[k_win
] == info
->nkz
[0]);
560 if (bRefinedCoul
|| bRefinedVdW
|| bRefinedGrid
)
562 fprintf(fp
, "Optimized PME settings:\n");
563 bCanUseOrigTPR
= FALSE
;
567 bCanUseOrigTPR
= TRUE
;
572 fprintf(fp
, " New Coulomb radius: %f nm (was %f nm)\n", info
->rcoulomb
[k_win
], info
->rcoulomb
[0]);
577 fprintf(fp
, " New Van der Waals radius: %f nm (was %f nm)\n", info
->rvdw
[k_win
], info
->rvdw
[0]);
582 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
],
583 info
->nkx
[0], info
->nky
[0], info
->nkz
[0]);
586 if (bCanUseOrigTPR
&& ntprs
> 1)
588 fprintf(fp
, "and original PME settings.\n");
593 /* Return the index of the mdp file that showed the highest performance
594 * and the optimal number of PME nodes */
596 *npme_optimal
= winPME
;
598 return bCanUseOrigTPR
;
602 /* Get the commands we need to set up the runs from environment variables */
603 static void get_program_paths(gmx_bool bThreads
, char *cmd_mpirun
[], char *cmd_mdrun
[])
606 const char def_mpirun
[] = "mpirun";
608 const char empty_mpirun
[] = "";
610 /* Get the commands we need to set up the runs from environment variables */
613 if ( (cp
= getenv("MPIRUN")) != nullptr)
615 *cmd_mpirun
= gmx_strdup(cp
);
619 *cmd_mpirun
= gmx_strdup(def_mpirun
);
624 *cmd_mpirun
= gmx_strdup(empty_mpirun
);
627 if (*cmd_mdrun
== nullptr)
629 /* The use of MDRUN is deprecated, but made available in 5.1
630 for backward compatibility. It may be removed in a future
632 if ( (cp
= getenv("MDRUN" )) != nullptr)
634 *cmd_mdrun
= gmx_strdup(cp
);
638 gmx_fatal(FARGS
, "The way to call mdrun must be set in the -mdrun command-line flag.");
643 /* Check that the commands will run mdrun (perhaps via mpirun) by
644 * running a very quick test simulation. Requires MPI environment or
645 * GPU support to be available if applicable. */
646 /* TODO implement feature to parse the log file to get the list of
647 compatible GPUs from mdrun, if the user of gmx tune-pme has not
649 static void check_mdrun_works(gmx_bool bThreads
,
650 const char *cmd_mpirun
,
652 const char *cmd_mdrun
,
653 gmx_bool bNeedGpuSupport
)
655 char *command
= nullptr;
659 const char filename
[] = "benchtest.log";
661 /* This string should always be identical to the one in copyrite.c,
662 * gmx_print_version_info() in the GMX_MPI section */
663 const char match_mpi
[] = "MPI library: MPI";
664 const char match_mdrun
[] = "Executable: ";
665 const char match_nogpu
[] = "GPU support: disabled";
666 gmx_bool bMdrun
= FALSE
;
667 gmx_bool bMPI
= FALSE
;
668 gmx_bool bHaveGpuSupport
= TRUE
;
670 /* Run a small test to see whether mpirun + mdrun work */
671 fprintf(stdout
, "Making sure that mdrun can be executed. ");
674 snew(command
, std::strlen(cmd_mdrun
) + std::strlen(cmd_np
) + std::strlen(filename
) + 50);
675 sprintf(command
, "%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mdrun
, cmd_np
, filename
);
679 snew(command
, std::strlen(cmd_mpirun
) + std::strlen(cmd_np
) + std::strlen(cmd_mdrun
) + std::strlen(filename
) + 50);
680 sprintf(command
, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun
, cmd_np
, cmd_mdrun
, filename
);
682 fprintf(stdout
, "Trying '%s' ... ", command
);
683 make_backup(filename
);
684 gmx_system_call(command
);
686 /* Check if we find the characteristic string in the output: */
687 if (!gmx_fexist(filename
))
689 gmx_fatal(FARGS
, "Output from test run could not be found.");
692 fp
= fopen(filename
, "r");
693 /* We need to scan the whole output file, since sometimes the queuing system
694 * also writes stuff to stdout/err */
697 cp
= fgets(line
, STRLEN
, fp
);
700 if (str_starts(line
, match_mdrun
) )
704 if (str_starts(line
, match_mpi
) )
708 if (str_starts(line
, match_nogpu
) )
710 bHaveGpuSupport
= FALSE
;
720 gmx_fatal(FARGS
, "Need a threaded version of mdrun. This one\n"
722 "seems to have been compiled with MPI instead.",
730 gmx_fatal(FARGS
, "Need an MPI-enabled version of mdrun. This one\n"
732 "seems to have been compiled without MPI support.",
739 gmx_fatal(FARGS
, "Cannot execute mdrun. Please check %s for problems!",
743 if (bNeedGpuSupport
&& !bHaveGpuSupport
)
745 gmx_fatal(FARGS
, "The mdrun executable did not have the expected GPU support.");
748 fprintf(stdout
, "passed.\n");
755 /* Handles the no-GPU case by emitting an empty string. */
756 static std::string
make_gpu_id_command_line(const char *eligible_gpu_ids
)
758 /* If the user has given no eligible GPU IDs, or we're trying the
759 * default behaviour, then there is nothing for tune_pme to give
760 * to mdrun -gpu_id */
761 if (eligible_gpu_ids
!= nullptr)
763 return gmx::formatString("-gpu_id %s", eligible_gpu_ids
);
767 return std::string();
770 static void launch_simulation(
771 gmx_bool bLaunch
, /* Should the simulation be launched? */
772 FILE *fp
, /* General log file */
773 gmx_bool bThreads
, /* whether to use threads */
774 char *cmd_mpirun
, /* Command for mpirun */
775 char *cmd_np
, /* Switch for -np or -ntmpi or empty */
776 char *cmd_mdrun
, /* Command for mdrun */
777 char *args_for_mdrun
, /* Arguments for mdrun */
778 const char *simulation_tpr
, /* This tpr will be simulated */
779 int nPMEnodes
, /* Number of PME ranks to use */
780 const char *eligible_gpu_ids
) /* Available GPU IDs for
781 * constructing mdrun command lines */
786 /* Make enough space for the system call command,
787 * (200 extra chars for -npme ... etc. options should suffice): */
788 snew(command
, std::strlen(cmd_mpirun
)+std::strlen(cmd_mdrun
)+std::strlen(cmd_np
)+std::strlen(args_for_mdrun
)+std::strlen(simulation_tpr
)+200);
790 auto cmd_gpu_ids
= make_gpu_id_command_line(eligible_gpu_ids
);
792 /* Note that the -passall options requires args_for_mdrun to be at the end
793 * of the command line string */
796 sprintf(command
, "%s%s-npme %d -s %s %s %s",
797 cmd_mdrun
, cmd_np
, nPMEnodes
, simulation_tpr
, args_for_mdrun
, cmd_gpu_ids
.c_str());
801 sprintf(command
, "%s%s%s -npme %d -s %s %s %s",
802 cmd_mpirun
, cmd_np
, cmd_mdrun
, nPMEnodes
, simulation_tpr
, args_for_mdrun
, cmd_gpu_ids
.c_str());
805 fprintf(fp
, "%s this command line to launch the simulation:\n\n%s", bLaunch
? "Using" : "Please use", command
);
809 /* Now the real thing! */
812 fprintf(stdout
, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command
);
815 gmx_system_call(command
);
820 static void modify_PMEsettings(
821 int64_t simsteps
, /* Set this value as number of time steps */
822 int64_t init_step
, /* Set this value as init_step */
823 const char *fn_best_tpr
, /* tpr file with the best performance */
824 const char *fn_sim_tpr
) /* name of tpr file to be launched */
830 t_inputrec irInstance
;
831 t_inputrec
*ir
= &irInstance
;
832 read_tpx_state(fn_best_tpr
, ir
, &state
, &mtop
);
834 /* Reset nsteps and init_step to the value of the input .tpr file */
835 ir
->nsteps
= simsteps
;
836 ir
->init_step
= init_step
;
838 /* Write the tpr file which will be launched */
839 sprintf(buf
, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr
, "%" PRId64
);
840 fprintf(stdout
, buf
, ir
->nsteps
);
842 write_tpx_state(fn_sim_tpr
, ir
, &state
, &mtop
);
845 static gmx_bool
can_scale_rvdw(int vdwtype
)
847 return (evdwCUT
== vdwtype
||
851 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
853 /* Make additional TPR files with more computational load for the
854 * direct space processors: */
855 static void make_benchmark_tprs(
856 const char *fn_sim_tpr
, /* READ : User-provided tpr file */
857 char *fn_bench_tprs
[], /* WRITE: Names of benchmark tpr files */
858 int64_t benchsteps
, /* Number of time steps for benchmark runs */
859 int64_t statesteps
, /* Step counter in checkpoint file */
860 real rmin
, /* Minimal Coulomb radius */
861 real rmax
, /* Maximal Coulomb radius */
862 bool bScaleRvdw
, /* Scale rvdw along with rcoulomb */
863 const int *ntprs
, /* No. of TPRs to write, each with a different
864 rcoulomb and fourierspacing */
865 t_inputinfo
*info
, /* Contains information about mdp file options */
866 FILE *fp
) /* Write the output here */
871 real nlist_buffer
; /* Thickness of the buffer regions for PME-switch potentials */
874 gmx_bool bNote
= FALSE
;
875 real add
; /* Add this to rcoul for the next test */
876 real fac
= 1.0; /* Scaling factor for Coulomb radius */
877 real fourierspacing
; /* Basic fourierspacing from tpr */
880 sprintf(buf
, "Making benchmark tpr file%s with %s time step%s",
881 *ntprs
> 1 ? "s" : "", "%" PRId64
, benchsteps
> 1 ? "s" : "");
882 fprintf(stdout
, buf
, benchsteps
);
885 sprintf(buf
, " (adding %s steps from checkpoint file)", "%" PRId64
);
886 fprintf(stdout
, buf
, statesteps
);
887 benchsteps
+= statesteps
;
889 fprintf(stdout
, ".\n");
891 t_inputrec irInstance
;
892 t_inputrec
*ir
= &irInstance
;
893 read_tpx_state(fn_sim_tpr
, ir
, &state
, &mtop
);
895 /* Check if some kind of PME was chosen */
896 if (EEL_PME(ir
->coulombtype
) == FALSE
)
898 gmx_fatal(FARGS
, "Can only do optimizations for simulations with %s electrostatics.",
902 /* Check if rcoulomb == rlist, which is necessary for plain PME. */
903 if ( (ir
->cutoff_scheme
!= ecutsVERLET
) &&
904 (eelPME
== ir
->coulombtype
) && !(ir
->rcoulomb
== ir
->rlist
))
906 gmx_fatal(FARGS
, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
907 EELTYPE(eelPME
), ir
->rcoulomb
, ir
->rlist
);
909 /* For other PME types, rcoulomb is allowed to be smaller than rlist */
910 else if (ir
->rcoulomb
> ir
->rlist
)
912 gmx_fatal(FARGS
, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
913 EELTYPE(ir
->coulombtype
), ir
->rcoulomb
, ir
->rlist
);
916 if (bScaleRvdw
&& ir
->rvdw
!= ir
->rcoulomb
)
918 fprintf(stdout
, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
922 /* Reduce the number of steps for the benchmarks */
923 info
->orig_sim_steps
= ir
->nsteps
;
924 ir
->nsteps
= benchsteps
;
925 /* We must not use init_step from the input tpr file for the benchmarks */
926 info
->orig_init_step
= ir
->init_step
;
929 /* For PME-switch potentials, keep the radial distance of the buffer region */
930 nlist_buffer
= ir
->rlist
- ir
->rcoulomb
;
932 /* Determine length of triclinic box vectors */
933 for (d
= 0; d
< DIM
; d
++)
936 for (i
= 0; i
< DIM
; i
++)
938 box_size
[d
] += state
.box
[d
][i
]*state
.box
[d
][i
];
940 box_size
[d
] = std::sqrt(box_size
[d
]);
943 if (ir
->fourier_spacing
> 0)
945 info
->fsx
[0] = ir
->fourier_spacing
;
946 info
->fsy
[0] = ir
->fourier_spacing
;
947 info
->fsz
[0] = ir
->fourier_spacing
;
951 /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
952 info
->fsx
[0] = box_size
[XX
]/ir
->nkx
;
953 info
->fsy
[0] = box_size
[YY
]/ir
->nky
;
954 info
->fsz
[0] = box_size
[ZZ
]/ir
->nkz
;
957 /* If no value for the fourierspacing was provided on the command line, we
958 * use the reconstruction from the tpr file */
959 if (ir
->fourier_spacing
> 0)
961 /* Use the spacing from the tpr */
962 fourierspacing
= ir
->fourier_spacing
;
966 /* Use the maximum observed spacing */
967 fourierspacing
= std::max(std::max(info
->fsx
[0], info
->fsy
[0]), info
->fsz
[0]);
970 fprintf(stdout
, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing
);
972 /* For performance comparisons the number of particles is useful to have */
973 fprintf(fp
, " Number of particles : %d\n", mtop
.natoms
);
975 /* Print information about settings of which some are potentially modified: */
976 fprintf(fp
, " Coulomb type : %s\n", EELTYPE(ir
->coulombtype
));
977 fprintf(fp
, " Grid spacing x y z : %f %f %f\n",
978 box_size
[XX
]/ir
->nkx
, box_size
[YY
]/ir
->nky
, box_size
[ZZ
]/ir
->nkz
);
979 fprintf(fp
, " Van der Waals type : %s\n", EVDWTYPE(ir
->vdwtype
));
980 if (ir_vdw_switched(ir
))
982 fprintf(fp
, " rvdw_switch : %f nm\n", ir
->rvdw_switch
);
984 if (EPME_SWITCHED(ir
->coulombtype
))
986 fprintf(fp
, " rlist : %f nm\n", ir
->rlist
);
989 /* Print a descriptive line about the tpr settings tested */
990 fprintf(fp
, "\nWill try these real/reciprocal workload settings:\n");
991 fprintf(fp
, " No. scaling rcoulomb");
992 fprintf(fp
, " nkx nky nkz");
993 fprintf(fp
, " spacing");
994 if (can_scale_rvdw(ir
->vdwtype
))
996 fprintf(fp
, " rvdw");
998 if (EPME_SWITCHED(ir
->coulombtype
))
1000 fprintf(fp
, " rlist");
1002 fprintf(fp
, " tpr file\n");
1004 /* Loop to create the requested number of tpr input files */
1005 for (j
= 0; j
< *ntprs
; j
++)
1007 /* The first .tpr is the provided one, just need to modify nsteps,
1008 * so skip the following block */
1011 /* Determine which Coulomb radii rc to use in the benchmarks */
1012 add
= (rmax
-rmin
)/(*ntprs
-1);
1013 if (gmx_within_tol(rmin
, info
->rcoulomb
[0], GMX_REAL_EPS
))
1015 ir
->rcoulomb
= rmin
+ j
*add
;
1017 else if (gmx_within_tol(rmax
, info
->rcoulomb
[0], GMX_REAL_EPS
))
1019 ir
->rcoulomb
= rmin
+ (j
-1)*add
;
1023 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
1024 add
= (rmax
-rmin
)/(*ntprs
-2);
1025 ir
->rcoulomb
= rmin
+ (j
-1)*add
;
1028 /* Determine the scaling factor fac */
1029 fac
= ir
->rcoulomb
/info
->rcoulomb
[0];
1031 /* Scale the Fourier grid spacing */
1032 ir
->nkx
= ir
->nky
= ir
->nkz
= 0;
1033 calcFftGrid(nullptr, state
.box
, fourierspacing
*fac
, minimalPmeGridSize(ir
->pme_order
),
1034 &ir
->nkx
, &ir
->nky
, &ir
->nkz
);
1036 /* Adjust other radii since various conditions need to be fulfilled */
1037 if (eelPME
== ir
->coulombtype
)
1039 /* plain PME, rcoulomb must be equal to rlist TODO only in the group scheme? */
1040 ir
->rlist
= ir
->rcoulomb
;
1044 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1045 ir
->rlist
= ir
->rcoulomb
+ nlist_buffer
;
1048 if (bScaleRvdw
&& can_scale_rvdw(ir
->vdwtype
))
1050 if (ecutsVERLET
== ir
->cutoff_scheme
||
1051 evdwPME
== ir
->vdwtype
)
1053 /* With either the Verlet cutoff-scheme or LJ-PME,
1054 the van der Waals radius must always equal the
1056 ir
->rvdw
= ir
->rcoulomb
;
1060 /* For vdw cutoff, rvdw >= rlist */
1061 ir
->rvdw
= std::max(info
->rvdw
[0], ir
->rlist
);
1064 } /* end of "if (j != 0)" */
1066 /* for j==0: Save the original settings
1067 * for j >0: Save modified radii and Fourier grids */
1068 info
->rcoulomb
[j
] = ir
->rcoulomb
;
1069 info
->rvdw
[j
] = ir
->rvdw
;
1070 info
->nkx
[j
] = ir
->nkx
;
1071 info
->nky
[j
] = ir
->nky
;
1072 info
->nkz
[j
] = ir
->nkz
;
1073 info
->rlist
[j
] = ir
->rlist
;
1074 info
->fsx
[j
] = fac
*fourierspacing
;
1075 info
->fsy
[j
] = fac
*fourierspacing
;
1076 info
->fsz
[j
] = fac
*fourierspacing
;
1078 /* Write the benchmark tpr file */
1079 fn_bench_tprs
[j
] = gmx_strdup(gmx::Path::concatenateBeforeExtension(fn_sim_tpr
, gmx::formatString("_bench%.2d", j
)).c_str());
1081 fprintf(stdout
, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs
[j
]);
1082 fprintf(stdout
, "%" PRId64
, ir
->nsteps
);
1085 fprintf(stdout
, ", scaling factor %f\n", fac
);
1089 fprintf(stdout
, ", unmodified settings\n");
1092 write_tpx_state(fn_bench_tprs
[j
], ir
, &state
, &mtop
);
1094 /* Write information about modified tpr settings to log file */
1095 fprintf(fp
, "%4d%10f%10f", j
, fac
, ir
->rcoulomb
);
1096 fprintf(fp
, "%5d%5d%5d", ir
->nkx
, ir
->nky
, ir
->nkz
);
1097 fprintf(fp
, " %9f ", info
->fsx
[j
]);
1098 if (can_scale_rvdw(ir
->vdwtype
))
1100 fprintf(fp
, "%10f", ir
->rvdw
);
1102 if (EPME_SWITCHED(ir
->coulombtype
))
1104 fprintf(fp
, "%10f", ir
->rlist
);
1106 fprintf(fp
, " %-14s\n", fn_bench_tprs
[j
]);
1108 /* Make it clear to the user that some additional settings were modified */
1109 if (!gmx_within_tol(ir
->rvdw
, info
->rvdw
[0], GMX_REAL_EPS
)
1110 || !gmx_within_tol(ir
->rlist
, info
->rlist
[0], GMX_REAL_EPS
) )
1117 fprintf(fp
, "\nNote that in addition to the Coulomb radius and the Fourier grid\n"
1118 "other input settings were also changed (see table above).\n"
1119 "Please check if the modified settings are appropriate.\n");
1126 /* Rename the files we want to keep to some meaningful filename and
1127 * delete the rest */
1128 static void cleanup(const t_filenm
*fnm
, int nfile
, int k
, int nnodes
,
1129 int nPMEnodes
, int nr
, gmx_bool bKeepStderr
)
1131 char numstring
[STRLEN
];
1132 const char *fn
= nullptr;
1137 fprintf(stdout
, "Cleaning up, deleting benchmark temp files ...\n");
1139 for (i
= 0; i
< nfile
; i
++)
1141 opt
= const_cast<char *>(fnm
[i
].opt
);
1142 if (std::strcmp(opt
, "-p") == 0)
1144 /* do nothing; keep this file */
1147 else if (std::strcmp(opt
, "-bg") == 0)
1149 /* Give the log file a nice name so one can later see which parameters were used */
1150 numstring
[0] = '\0';
1153 sprintf(numstring
, "_%d", nr
);
1155 std::string newfilename
= gmx::formatString("%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile
, fnm
), k
, nnodes
, nPMEnodes
, numstring
);
1156 if (gmx_fexist(opt2fn("-bg", nfile
, fnm
)))
1158 fprintf(stdout
, "renaming log file to %s\n", newfilename
.c_str());
1159 make_backup(newfilename
);
1160 rename(opt2fn("-bg", nfile
, fnm
), newfilename
.c_str());
1163 else if (std::strcmp(opt
, "-err") == 0)
1165 /* This file contains the output of stderr. We want to keep it in
1166 * cases where there have been problems. */
1167 fn
= opt2fn(opt
, nfile
, fnm
);
1168 numstring
[0] = '\0';
1171 sprintf(numstring
, "_%d", nr
);
1173 std::string newfilename
= gmx::formatString("%s_no%d_np%d_npme%d%s", fn
, k
, nnodes
, nPMEnodes
, numstring
);
1178 fprintf(stdout
, "Saving stderr output in %s\n", newfilename
.c_str());
1179 make_backup(newfilename
);
1180 rename(fn
, newfilename
.c_str());
1184 fprintf(stdout
, "Deleting %s\n", fn
);
1189 /* Delete the files which are created for each benchmark run: (options -b*) */
1190 else if ( (0 == std::strncmp(opt
, "-b", 2)) && (opt2bSet(opt
, nfile
, fnm
) || !is_optional(&fnm
[i
])) )
1192 remove_if_exists(opt2fn(opt
, nfile
, fnm
));
1199 eNpmeAuto
, eNpmeAll
, eNpmeReduced
, eNpmeSubset
, eNpmeNr
1202 /* Create a list of numbers of PME nodes to test */
1203 static void make_npme_list(
1204 const char *npmevalues_opt
, /* Make a complete list with all
1205 * possibilities or a short list that keeps only
1206 * reasonable numbers of PME nodes */
1207 int *nentries
, /* Number of entries we put in the nPMEnodes list */
1208 int *nPMEnodes
[], /* Each entry contains the value for -npme */
1209 int nnodes
, /* Total number of nodes to do the tests on */
1210 int minPMEnodes
, /* Minimum number of PME nodes */
1211 int maxPMEnodes
) /* Maximum number of PME nodes */
1214 int min_factor
= 1; /* We request that npp and npme have this minimal
1215 * largest common factor (depends on npp) */
1216 int nlistmax
; /* Max. list size */
1217 int nlist
; /* Actual number of entries in list */
1221 /* Do we need to check all possible values for -npme or is a reduced list enough? */
1222 if (!std::strcmp(npmevalues_opt
, "all") )
1226 else if (!std::strcmp(npmevalues_opt
, "subset") )
1228 eNPME
= eNpmeSubset
;
1230 else /* "auto" or "range" */
1236 else if (nnodes
< 128)
1238 eNPME
= eNpmeReduced
;
1242 eNPME
= eNpmeSubset
;
1246 /* Calculate how many entries we could possibly have (in case of -npme all) */
1249 nlistmax
= maxPMEnodes
- minPMEnodes
+ 3;
1250 if (0 == minPMEnodes
)
1260 /* Now make the actual list which is at most of size nlist */
1261 snew(*nPMEnodes
, nlistmax
);
1262 nlist
= 0; /* start counting again, now the real entries in the list */
1263 for (i
= 0; i
< nlistmax
- 2; i
++)
1265 npme
= maxPMEnodes
- i
;
1276 /* For 2d PME we want a common largest factor of at least the cube
1277 * root of the number of PP nodes */
1278 min_factor
= static_cast<int>(std::cbrt(npp
));
1281 gmx_fatal(FARGS
, "Unknown option for eNPME in make_npme_list");
1283 if (gmx_greatest_common_divisor(npp
, npme
) >= min_factor
)
1285 (*nPMEnodes
)[nlist
] = npme
;
1289 /* We always test 0 PME nodes and the automatic number */
1290 *nentries
= nlist
+ 2;
1291 (*nPMEnodes
)[nlist
] = 0;
1292 (*nPMEnodes
)[nlist
+1] = -1;
1294 fprintf(stderr
, "Will try the following %d different values for -npme:\n", *nentries
);
1295 for (i
= 0; i
< *nentries
-1; i
++)
1297 fprintf(stderr
, "%d, ", (*nPMEnodes
)[i
]);
1299 fprintf(stderr
, "and %d (auto).\n", (*nPMEnodes
)[*nentries
-1]);
1303 /* Allocate memory to store the performance data */
1304 static void init_perfdata(t_perf
*perfdata
[], int ntprs
, int datasets
, int repeats
)
1309 for (k
= 0; k
< ntprs
; k
++)
1311 snew(perfdata
[k
], datasets
);
1312 for (i
= 0; i
< datasets
; i
++)
1314 for (j
= 0; j
< repeats
; j
++)
1316 snew(perfdata
[k
][i
].Gcycles
, repeats
);
1317 snew(perfdata
[k
][i
].ns_per_day
, repeats
);
1318 snew(perfdata
[k
][i
].PME_f_load
, repeats
);
1325 /* Check for errors on mdrun -h */
1326 static void make_sure_it_runs(char *mdrun_cmd_line
, int length
, FILE *fp
,
1327 const t_filenm
*fnm
, int nfile
)
1329 char *command
, *msg
;
1332 snew(command
, length
+ 15);
1333 snew(msg
, length
+ 500);
1335 fprintf(stdout
, "Making sure the benchmarks can be executed by running just 1 step...\n");
1336 sprintf(command
, "%s -nsteps 1 -quiet", mdrun_cmd_line
);
1337 fprintf(stdout
, "Executing '%s' ...\n", command
);
1338 ret
= gmx_system_call(command
);
1342 /* To prevent confusion, do not again issue a gmx_fatal here since we already
1343 * get the error message from mdrun itself */
1345 "Cannot run the first benchmark simulation! Please check the error message of\n"
1346 "mdrun for the source of the problem. Did you provide a command line\n"
1347 "argument that neither gmx tune_pme nor mdrun understands? If you're\n"
1348 "sure your command line should work, you can bypass this check with \n"
1349 "gmx tune_pme -nocheck. The failing command was:\n"
1350 "\n%s\n\n", command
);
1352 fprintf(stderr
, "%s", msg
);
1354 fprintf(fp
, "%s", msg
);
1358 fprintf(stdout
, "Benchmarks can be executed!\n");
1360 /* Clean up the benchmark output files we just created */
1361 fprintf(stdout
, "Cleaning up ...\n");
1362 remove_if_exists(opt2fn("-bc", nfile
, fnm
));
1363 remove_if_exists(opt2fn("-be", nfile
, fnm
));
1364 remove_if_exists(opt2fn("-bcpo", nfile
, fnm
));
1365 remove_if_exists(opt2fn("-bg", nfile
, fnm
));
1366 remove_if_exists(opt2fn("-bo", nfile
, fnm
));
1367 remove_if_exists(opt2fn("-bx", nfile
, fnm
));
1373 static void do_the_tests(
1374 FILE *fp
, /* General tune_pme output file */
1375 char **tpr_names
, /* Filenames of the input files to test */
1376 int maxPMEnodes
, /* Max fraction of nodes to use for PME */
1377 int minPMEnodes
, /* Min fraction of nodes to use for PME */
1378 int npme_fixed
, /* If >= -1, test fixed number of PME
1380 const char *npmevalues_opt
, /* Which -npme values should be tested */
1381 t_perf
**perfdata
, /* Here the performace data is stored */
1382 int *pmeentries
, /* Entries in the nPMEnodes list */
1383 int repeats
, /* Repeat each test this often */
1384 int nnodes
, /* Total number of nodes = nPP + nPME */
1385 int nr_tprs
, /* Total number of tpr files to test */
1386 gmx_bool bThreads
, /* Threads or MPI? */
1387 char *cmd_mpirun
, /* mpirun command string */
1388 char *cmd_np
, /* "-np", "-n", whatever mpirun needs */
1389 char *cmd_mdrun
, /* mdrun command string */
1390 char *cmd_args_bench
, /* arguments for mdrun in a string */
1391 const t_filenm
*fnm
, /* List of filenames from command line */
1392 int nfile
, /* Number of files specified on the cmdl. */
1393 int presteps
, /* DLB equilibration steps, is checked */
1394 int64_t cpt_steps
, /* Time step counter in the checkpoint */
1395 gmx_bool bCheck
, /* Check whether benchmark mdrun works */
1396 const char *eligible_gpu_ids
) /* GPU IDs for
1397 * constructing mdrun command lines */
1399 int i
, nr
, k
, ret
, count
= 0, totaltests
;
1400 int *nPMEnodes
= nullptr;
1401 t_perf
*pd
= nullptr;
1403 char *command
, *cmd_stub
;
1405 gmx_bool bResetProblem
= FALSE
;
1406 gmx_bool bFirst
= TRUE
;
1408 /* This string array corresponds to the eParselog enum type at the start
1410 const char* ParseLog
[] = {
1412 "Logfile not found!",
1413 "No timings, logfile truncated?",
1414 "Run was terminated.",
1415 "Counters were not reset properly.",
1416 "No DD grid found for these settings.",
1417 "TPX version conflict!",
1418 "mdrun was not started in parallel!",
1419 "Number of PP ranks has a prime factor that is too large.",
1420 "The number of PP ranks did not suit the number of GPUs.",
1421 "Some GPUs were not detected or are incompatible.",
1422 "An error occurred."
1424 char str_PME_f_load
[13];
1427 /* Allocate space for the mdrun command line. 100 extra characters should
1428 be more than enough for the -npme etcetera arguments */
1429 cmdline_length
= std::strlen(cmd_mpirun
)
1430 + std::strlen(cmd_np
)
1431 + std::strlen(cmd_mdrun
)
1432 + std::strlen(cmd_args_bench
)
1433 + std::strlen(tpr_names
[0]) + 100;
1434 snew(command
, cmdline_length
);
1435 snew(cmd_stub
, cmdline_length
);
1437 /* Construct the part of the command line that stays the same for all tests: */
1440 sprintf(cmd_stub
, "%s%s", cmd_mdrun
, cmd_np
);
1444 sprintf(cmd_stub
, "%s%s%s ", cmd_mpirun
, cmd_np
, cmd_mdrun
);
1447 /* Create a list of numbers of PME nodes to test */
1448 if (npme_fixed
< -1)
1450 make_npme_list(npmevalues_opt
, pmeentries
, &nPMEnodes
,
1451 nnodes
, minPMEnodes
, maxPMEnodes
);
1457 nPMEnodes
[0] = npme_fixed
;
1458 fprintf(stderr
, "Will use a fixed number of %d PME-only ranks.\n", nPMEnodes
[0]);
1463 fprintf(fp
, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1465 finalize(opt2fn("-p", nfile
, fnm
));
1469 /* Allocate one dataset for each tpr input file: */
1470 init_perfdata(perfdata
, nr_tprs
, *pmeentries
, repeats
);
1472 /*****************************************/
1473 /* Main loop over all tpr files to test: */
1474 /*****************************************/
1475 totaltests
= nr_tprs
*(*pmeentries
)*repeats
;
1476 for (k
= 0; k
< nr_tprs
; k
++)
1478 fprintf(fp
, "\nIndividual timings for input file %d (%s):\n", k
, tpr_names
[k
]);
1479 fprintf(fp
, "PME ranks Gcycles ns/day PME/f Remark\n");
1480 /* Loop over various numbers of PME nodes: */
1481 for (i
= 0; i
< *pmeentries
; i
++)
1483 pd
= &perfdata
[k
][i
];
1485 auto cmd_gpu_ids
= make_gpu_id_command_line(eligible_gpu_ids
);
1487 /* Loop over the repeats for each scenario: */
1488 for (nr
= 0; nr
< repeats
; nr
++)
1490 pd
->nPMEnodes
= nPMEnodes
[i
];
1492 /* Add -npme and -s to the command line and save it. Note that
1493 * the -passall (if set) options requires cmd_args_bench to be
1494 * at the end of the command line string */
1495 snew(pd
->mdrun_cmd_line
, cmdline_length
);
1496 sprintf(pd
->mdrun_cmd_line
, "%s-npme %d -s %s %s %s",
1497 cmd_stub
, pd
->nPMEnodes
, tpr_names
[k
], cmd_args_bench
, cmd_gpu_ids
.c_str());
1499 /* To prevent that all benchmarks fail due to a show-stopper argument
1500 * on the mdrun command line, we make a quick check first.
1501 * This check can be turned off in cases where the automatically chosen
1502 * number of PME-only ranks leads to a number of PP ranks for which no
1503 * decomposition can be found (e.g. for large prime numbers) */
1504 if (bFirst
&& bCheck
)
1506 /* TODO This check is really for a functional
1507 * .tpr, and if we need it, it should take place
1508 * for every .tpr, and the logic for it should be
1509 * immediately inside the loop over k, not in
1510 * this inner loop. */
1511 char *temporary_cmd_line
;
1513 snew(temporary_cmd_line
, cmdline_length
);
1514 /* TODO -npme 0 is more likely to succeed at low
1515 parallelism than the default of -npme -1, but
1516 is more likely to fail at the scaling limit
1517 when the PP domains may be too small. "mpirun
1518 -np 1 mdrun" is probably a reasonable thing to
1519 do for this check, but it'll be easier to
1520 implement that after some refactoring of how
1521 the number of MPI ranks is managed. */
1522 sprintf(temporary_cmd_line
, "%s-npme 0 -nb cpu -s %s %s",
1523 cmd_stub
, tpr_names
[k
], cmd_args_bench
);
1524 make_sure_it_runs(temporary_cmd_line
, cmdline_length
, fp
, fnm
, nfile
);
1528 /* Do a benchmark simulation: */
1531 sprintf(buf
, ", pass %d/%d", nr
+1, repeats
);
1537 fprintf(stdout
, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1538 (100.0*count
)/totaltests
,
1539 k
+1, nr_tprs
, i
+1, *pmeentries
, buf
);
1540 make_backup(opt2fn("-err", nfile
, fnm
));
1541 sprintf(command
, "%s 1> /dev/null 2>%s", pd
->mdrun_cmd_line
, opt2fn("-err", nfile
, fnm
));
1542 fprintf(stdout
, "%s\n", pd
->mdrun_cmd_line
);
1543 gmx_system_call(command
);
1545 /* Collect the performance data from the log file; also check stderr
1546 * for fatal errors */
1547 ret
= parse_logfile(opt2fn("-bg", nfile
, fnm
), opt2fn("-err", nfile
, fnm
),
1548 pd
, nr
, presteps
, cpt_steps
, nnodes
);
1549 if ((presteps
> 0) && (ret
== eParselogResetProblem
))
1551 bResetProblem
= TRUE
;
1554 if (-1 == pd
->nPMEnodes
)
1556 sprintf(buf
, "(%3d)", pd
->guessPME
);
1564 if (pd
->PME_f_load
[nr
] > 0.0)
1566 sprintf(str_PME_f_load
, "%12.3f", pd
->PME_f_load
[nr
]);
1570 sprintf(str_PME_f_load
, "%s", " - ");
1573 /* Write the data we got to disk */
1574 fprintf(fp
, "%4d%s %12.3f %12.3f %s %s", pd
->nPMEnodes
,
1575 buf
, pd
->Gcycles
[nr
], pd
->ns_per_day
[nr
], str_PME_f_load
, ParseLog
[ret
]);
1576 if (!(ret
== eParselogOK
|| ret
== eParselogNoDDGrid
|| ret
== eParselogNotFound
) )
1578 fprintf(fp
, " Check %s file for problems.", ret
== eParselogFatal
? "err" : "log");
1584 /* Do some cleaning up and delete the files we do not need any more */
1585 cleanup(fnm
, nfile
, k
, nnodes
, pd
->nPMEnodes
, nr
, ret
== eParselogFatal
);
1587 /* If the first run with this number of processors already failed, do not try again: */
1588 if (pd
->Gcycles
[0] <= 0.0 && repeats
> 1)
1590 fprintf(stdout
, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1591 count
+= repeats
-(nr
+1);
1594 } /* end of repeats loop */
1595 } /* end of -npme loop */
1596 } /* end of tpr file loop */
1601 fprintf(fp
, "WARNING: The cycle and time step counters could not be reset properly. ");
1609 static void check_input(
1616 real maxPMEfraction
,
1617 real minPMEfraction
,
1619 int64_t bench_nsteps
,
1620 const t_filenm
*fnm
,
1630 /* Make sure the input file exists */
1631 if (!gmx_fexist(opt2fn("-s", nfile
, fnm
)))
1633 gmx_fatal(FARGS
, "File %s not found.", opt2fn("-s", nfile
, fnm
));
1636 /* Make sure that the checkpoint file is not overwritten during benchmarking */
1637 if ( (0 == std::strcmp(opt2fn("-cpi", nfile
, fnm
), opt2fn("-bcpo", nfile
, fnm
)) ) && (sim_part
> 1) )
1639 gmx_fatal(FARGS
, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1640 "The checkpoint input file must not be overwritten during the benchmarks.\n");
1643 /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1646 gmx_fatal(FARGS
, "Number of repeats < 0!");
1649 /* Check number of nodes */
1652 gmx_fatal(FARGS
, "Number of ranks/threads must be a positive integer.");
1655 /* Automatically choose -ntpr if not set */
1665 /* Set a reasonable scaling factor for rcoulomb */
1668 *rmax
= rcoulomb
* 1.2;
1671 fprintf(stderr
, "Will test %d tpr file%s.\n", *ntprs
, *ntprs
== 1 ? "" : "s");
1677 fprintf(stderr
, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1681 /* Make shure that rmin <= rcoulomb <= rmax */
1690 if (!(*rmin
<= *rmax
) )
1692 gmx_fatal(FARGS
, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1693 "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin
, *rmax
, rcoulomb
);
1695 /* Add test scenarios if rmin or rmax were set */
1698 if (!gmx_within_tol(*rmin
, rcoulomb
, GMX_REAL_EPS
) && (*ntprs
== 1) )
1701 fprintf(stderr
, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1704 if (!gmx_within_tol(*rmax
, rcoulomb
, GMX_REAL_EPS
) && (*ntprs
== 1) )
1707 fprintf(stderr
, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1712 /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1713 if (!gmx_within_tol(*rmax
, rcoulomb
, GMX_REAL_EPS
) || !gmx_within_tol(*rmin
, rcoulomb
, GMX_REAL_EPS
) )
1715 *ntprs
= std::max(*ntprs
, 2);
1718 /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1719 if (!gmx_within_tol(*rmax
, rcoulomb
, GMX_REAL_EPS
) && !gmx_within_tol(*rmin
, rcoulomb
, GMX_REAL_EPS
) )
1721 *ntprs
= std::max(*ntprs
, 3);
1726 fprintf(stderr
, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs
);
1731 if (gmx_within_tol(*rmin
, rcoulomb
, GMX_REAL_EPS
) && gmx_within_tol(rcoulomb
, *rmax
, GMX_REAL_EPS
)) /* We have just a single rc */
1733 fprintf(stderr
, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1734 "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1735 "with correspondingly adjusted PME grid settings\n");
1740 /* Check whether max and min fraction are within required values */
1741 if (maxPMEfraction
> 0.5 || maxPMEfraction
< 0)
1743 gmx_fatal(FARGS
, "-max must be between 0 and 0.5");
1745 if (minPMEfraction
> 0.5 || minPMEfraction
< 0)
1747 gmx_fatal(FARGS
, "-min must be between 0 and 0.5");
1749 if (maxPMEfraction
< minPMEfraction
)
1751 gmx_fatal(FARGS
, "-max must be larger or equal to -min");
1754 /* Check whether the number of steps - if it was set - has a reasonable value */
1755 if (bench_nsteps
< 0)
1757 gmx_fatal(FARGS
, "Number of steps must be positive.");
1760 if (bench_nsteps
> 10000 || bench_nsteps
< 100)
1762 fprintf(stderr
, "WARNING: steps=");
1763 fprintf(stderr
, "%" PRId64
, bench_nsteps
);
1764 fprintf(stderr
, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps
< 100) ? "few" : "many");
1769 gmx_fatal(FARGS
, "Cannot have a negative number of presteps.\n");
1772 /* Check for rcoulomb scaling if more than one .tpr file is tested */
1775 if (*rmin
/rcoulomb
< 0.75 || *rmax
/rcoulomb
> 1.25)
1777 fprintf(stderr
, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1781 /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1782 * only. We need to check whether the requested number of PME-only nodes
1784 if (npme_fixed
> -1)
1786 /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1787 if (2*npme_fixed
> nnodes
)
1789 gmx_fatal(FARGS
, "Cannot have more than %d PME-only ranks for a total of %d ranks (you chose %d).\n",
1790 nnodes
/2, nnodes
, npme_fixed
);
1792 if ((npme_fixed
> 0) && (5*npme_fixed
< nnodes
))
1794 fprintf(stderr
, "WARNING: Only %g percent of the ranks are assigned as PME-only ranks.\n",
1795 (100.0*npme_fixed
)/nnodes
);
1797 if (opt2parg_bSet("-min", npargs
, pa
) || opt2parg_bSet("-max", npargs
, pa
))
1799 fprintf(stderr
, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1800 " fixed number of PME-only ranks is requested with -fix.\n");
1806 /* Returns TRUE when "opt" is needed at launch time */
1807 static gmx_bool
is_launch_file(char *opt
, gmx_bool bSet
)
1809 if (0 == std::strncmp(opt
, "-swap", 5))
1814 /* Apart from the input .tpr and the output log files we need all options that
1815 * were set on the command line and that do not start with -b */
1816 if (0 == std::strncmp(opt
, "-b", 2) || 0 == std::strncmp(opt
, "-s", 2)
1817 || 0 == std::strncmp(opt
, "-err", 4) || 0 == std::strncmp(opt
, "-p", 2) )
1826 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1827 static gmx_bool
is_bench_file(char *opt
, gmx_bool bSet
, gmx_bool bOptional
, gmx_bool bIsOutput
)
1829 /* Apart from the input .tpr, all files starting with "-b" are for
1830 * _b_enchmark files exclusively */
1831 if (0 == std::strncmp(opt
, "-s", 2))
1836 if (0 == std::strncmp(opt
, "-b", 2) || 0 == std::strncmp(opt
, "-s", 2))
1838 return !bOptional
|| bSet
;
1848 return bSet
; /* These are additional input files like -cpi -ei */
1854 /* Adds 'buf' to 'str' */
1855 static void add_to_string(char **str
, const char *buf
)
1860 len
= std::strlen(*str
) + std::strlen(buf
) + 1;
1862 std::strcat(*str
, buf
);
1866 /* Create the command line for the benchmark as well as for the real run */
1867 static void create_command_line_snippets(
1868 gmx_bool bAppendFiles
,
1869 gmx_bool bKeepAndNumCPT
,
1870 gmx_bool bResetHWay
,
1874 char *cmd_args_bench
[], /* command line arguments for benchmark runs */
1875 char *cmd_args_launch
[], /* command line arguments for simulation run */
1876 char extra_args
[], /* Add this to the end of the command line */
1877 char *deffnm
) /* Default file names, or NULL if not set */
1882 char strbuf
[STRLEN
];
1885 /* strlen needs at least '\0' as a string: */
1886 snew(*cmd_args_bench
, 1);
1887 snew(*cmd_args_launch
, 1);
1888 *cmd_args_launch
[0] = '\0';
1889 *cmd_args_bench
[0] = '\0';
1892 /*******************************************/
1893 /* 1. Process other command line arguments */
1894 /*******************************************/
1897 /* Add equilibration steps to benchmark options */
1898 sprintf(strbuf
, "-resetstep %d ", presteps
);
1899 add_to_string(cmd_args_bench
, strbuf
);
1901 /* These switches take effect only at launch time */
1904 sprintf(strbuf
, "-deffnm %s ", deffnm
);
1905 add_to_string(cmd_args_launch
, strbuf
);
1909 add_to_string(cmd_args_launch
, "-noappend ");
1913 add_to_string(cmd_args_launch
, "-cpnum ");
1917 add_to_string(cmd_args_launch
, "-resethway ");
1920 /********************/
1921 /* 2. Process files */
1922 /********************/
1923 for (i
= 0; i
< nfile
; i
++)
1925 opt
= const_cast<char *>(fnm
[i
].opt
);
1926 name
= opt2fn(opt
, nfile
, fnm
);
1928 /* Strbuf contains the options, now let's sort out where we need that */
1929 sprintf(strbuf
, "%s %s ", opt
, name
);
1931 if (is_bench_file(opt
, opt2bSet(opt
, nfile
, fnm
), is_optional(&fnm
[i
]), is_output(&fnm
[i
])) )
1933 /* All options starting with -b* need the 'b' removed,
1934 * therefore overwrite strbuf */
1935 if (0 == std::strncmp(opt
, "-b", 2))
1937 sprintf(strbuf
, "-%s %s ", &opt
[2], name
);
1940 add_to_string(cmd_args_bench
, strbuf
);
1943 if (is_launch_file(opt
, opt2bSet(opt
, nfile
, fnm
)) )
1945 add_to_string(cmd_args_launch
, strbuf
);
1949 add_to_string(cmd_args_bench
, extra_args
);
1950 add_to_string(cmd_args_launch
, extra_args
);
1954 /* Set option opt */
1955 static void setopt(const char *opt
, int nfile
, t_filenm fnm
[])
1959 for (i
= 0; (i
< nfile
); i
++)
1961 if (std::strcmp(opt
, fnm
[i
].opt
) == 0)
1963 fnm
[i
].flag
|= ffSET
;
1969 /* This routine inspects the tpr file and ...
1970 * 1. checks for output files that get triggered by a tpr option. These output
1971 * files are marked as 'set' to allow for a proper cleanup after each
1973 * 2. returns the PME:PP load ratio
1974 * 3. returns rcoulomb from the tpr */
1975 static float inspect_tpr(int nfile
, t_filenm fnm
[], real
*rcoulomb
)
1977 gmx_bool bTpi
; /* Is test particle insertion requested? */
1978 gmx_bool bFree
; /* Is a free energy simulation requested? */
1979 gmx_bool bNM
; /* Is a normal mode analysis requested? */
1980 gmx_bool bSwap
; /* Is water/ion position swapping requested? */
1985 /* Check tpr file for options that trigger extra output files */
1986 t_inputrec irInstance
;
1987 t_inputrec
*ir
= &irInstance
;
1988 read_tpx_state(opt2fn("-s", nfile
, fnm
), ir
, &state
, &mtop
);
1989 bFree
= (efepNO
!= ir
->efep
);
1990 bNM
= (eiNM
== ir
->eI
);
1991 bSwap
= (eswapNO
!= ir
->eSwapCoords
);
1992 bTpi
= EI_TPI(ir
->eI
);
1994 /* Set these output files on the tuning command-line */
1997 setopt("-pf", nfile
, fnm
);
1998 setopt("-px", nfile
, fnm
);
2002 setopt("-dhdl", nfile
, fnm
);
2006 setopt("-tpi", nfile
, fnm
);
2007 setopt("-tpid", nfile
, fnm
);
2011 setopt("-mtx", nfile
, fnm
);
2015 setopt("-swap", nfile
, fnm
);
2018 *rcoulomb
= ir
->rcoulomb
;
2020 /* Return the estimate for the number of PME nodes */
2021 float npme
= pme_load_estimate(&mtop
, ir
, state
.box
);
2026 static void couple_files_options(int nfile
, t_filenm fnm
[])
2029 gmx_bool bSet
, bBench
;
2034 for (i
= 0; i
< nfile
; i
++)
2036 opt
= const_cast<char *>(fnm
[i
].opt
);
2037 bSet
= ((fnm
[i
].flag
& ffSET
) != 0);
2038 bBench
= (0 == std::strncmp(opt
, "-b", 2));
2040 /* Check optional files */
2041 /* If e.g. -eo is set, then -beo also needs to be set */
2042 if (is_optional(&fnm
[i
]) && bSet
&& !bBench
)
2044 sprintf(buf
, "-b%s", &opt
[1]);
2045 setopt(buf
, nfile
, fnm
);
2047 /* If -beo is set, then -eo also needs to be! */
2048 if (is_optional(&fnm
[i
]) && bSet
&& bBench
)
2050 sprintf(buf
, "-%s", &opt
[2]);
2051 setopt(buf
, nfile
, fnm
);
2057 #define BENCHSTEPS (1000)
2059 int gmx_tune_pme(int argc
, char *argv
[])
2061 const char *desc
[] = {
2062 "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of ranks, [THISMODULE] systematically",
2063 "times [gmx-mdrun] with various numbers of PME-only ranks and determines",
2064 "which setting is fastest. It will also test whether performance can",
2065 "be enhanced by shifting load from the reciprocal to the real space",
2066 "part of the Ewald sum. ",
2067 "Simply pass your [REF].tpr[ref] file to [THISMODULE] together with other options",
2068 "for [gmx-mdrun] as needed.[PAR]",
2069 "[THISMODULE] needs to call [gmx-mdrun] and so requires that you",
2070 "specify how to call mdrun with the argument to the [TT]-mdrun[tt]",
2071 "parameter. Depending how you have built GROMACS, values such as",
2072 "'gmx mdrun', 'gmx_d mdrun', or 'mdrun_mpi' might be needed.[PAR]",
2073 "The program that runs MPI programs can be set in the environment variable",
2074 "MPIRUN (defaults to 'mpirun'). Note that for certain MPI frameworks,",
2075 "you need to provide a machine- or hostfile. This can also be passed",
2076 "via the MPIRUN variable, e.g.[PAR]",
2077 "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt]",
2078 "Note that in such cases it is normally necessary to compile",
2079 "and/or run [THISMODULE] without MPI support, so that it can call",
2080 "the MPIRUN program.[PAR]",
2081 "Before doing the actual benchmark runs, [THISMODULE] will do a quick",
2082 "check whether [gmx-mdrun] works as expected with the provided parallel settings",
2083 "if the [TT]-check[tt] option is activated (the default).",
2084 "Please call [THISMODULE] with the normal options you would pass to",
2085 "[gmx-mdrun] and add [TT]-np[tt] for the number of ranks to perform the",
2086 "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2087 "to repeat each test several times to get better statistics. [PAR]",
2088 "[THISMODULE] can test various real space / reciprocal space workloads",
2089 "for you. With [TT]-ntpr[tt] you control how many extra [REF].tpr[ref] files will be",
2090 "written with enlarged cutoffs and smaller Fourier grids respectively.",
2091 "Typically, the first test (number 0) will be with the settings from the input",
2092 "[REF].tpr[ref] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2093 "specified by [TT]-rmax[tt] with a somewhat smaller PME grid at the same time. ",
2094 "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2095 "The remaining [REF].tpr[ref] files will have equally-spaced Coulomb radii (and Fourier ",
2096 "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2097 "if you just seek the optimal number of PME-only ranks; in that case",
2098 "your input [REF].tpr[ref] file will remain unchanged.[PAR]",
2099 "For the benchmark runs, the default of 1000 time steps should suffice for most",
2100 "MD systems. The dynamic load balancing needs about 100 time steps",
2101 "to adapt to local load imbalances, therefore the time step counters",
2102 "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2103 "for a higher accuracy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2104 "From the 'DD' load imbalance entries in the md.log output file you",
2105 "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]",
2106 "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2107 "After calling [gmx-mdrun] several times, detailed performance information",
2108 "is available in the output file [TT]perf.out[tt].",
2109 "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2110 "(options [TT]-b*[tt]), these will be automatically deleted after each test.[PAR]",
2111 "If you want the simulation to be started automatically with the",
2112 "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2113 "Basic support for GPU-enabled [TT]mdrun[tt] exists. Give a string containing the IDs",
2114 "of the GPUs that you wish to use in the optimization in the [TT]-gpu_id[tt]",
2115 "command-line argument. This works exactly like [TT]mdrun -gpu_id[tt], does not imply a mapping,",
2116 "and merely declares the eligible set of GPU devices. [TT]gmx-tune_pme[tt] will construct calls to",
2117 "mdrun that use this set appropriately. [TT]gmx-tune_pme[tt] does not support",
2118 "[TT]-gputasks[tt].[PAR]",
2123 int pmeentries
= 0; /* How many values for -npme do we actually test for each tpr file */
2124 real maxPMEfraction
= 0.50;
2125 real minPMEfraction
= 0.25;
2126 int maxPMEnodes
, minPMEnodes
;
2127 float guessPMEratio
; /* guessed PME:PP ratio based on the tpr file */
2128 float guessPMEnodes
;
2129 int npme_fixed
= -2; /* If >= -1, use only this number
2130 * of PME-only nodes */
2132 real rmin
= 0.0, rmax
= 0.0; /* min and max value for rcoulomb if scaling is requested */
2133 real rcoulomb
= -1.0; /* Coulomb radius as set in .tpr file */
2134 gmx_bool bScaleRvdw
= TRUE
;
2135 int64_t bench_nsteps
= BENCHSTEPS
;
2136 int64_t new_sim_nsteps
= -1; /* -1 indicates: not set by the user */
2137 int64_t cpt_steps
= 0; /* Step counter in .cpt input file */
2138 int presteps
= 1500; /* Do a full cycle reset after presteps steps */
2139 gmx_bool bOverwrite
= FALSE
, bKeepTPR
;
2140 gmx_bool bLaunch
= FALSE
;
2141 char *ExtraArgs
= nullptr;
2142 char **tpr_names
= nullptr;
2143 const char *simulation_tpr
= nullptr;
2144 char *deffnm
= nullptr;
2145 int best_npme
, best_tpr
;
2146 int sim_part
= 1; /* For benchmarks with checkpoint files */
2150 /* Default program names if nothing else is found */
2151 char *cmd_mpirun
= nullptr, *cmd_mdrun
= nullptr;
2152 char *cmd_args_bench
, *cmd_args_launch
;
2153 char *cmd_np
= nullptr;
2155 /* IDs of GPUs that are eligible for computation */
2156 char *eligible_gpu_ids
= nullptr;
2158 t_perf
**perfdata
= nullptr;
2163 /* Print out how long the tuning took */
2166 static t_filenm fnm
[] = {
2168 { efOUT
, "-p", "perf", ffWRITE
},
2169 { efLOG
, "-err", "bencherr", ffWRITE
},
2170 { efTPR
, "-so", "tuned", ffWRITE
},
2172 { efTPR
, "-s", nullptr, ffREAD
},
2173 { efTRN
, "-o", nullptr, ffWRITE
},
2174 { efCOMPRESSED
, "-x", nullptr, ffOPTWR
},
2175 { efCPT
, "-cpi", nullptr, ffOPTRD
},
2176 { efCPT
, "-cpo", nullptr, ffOPTWR
},
2177 { efSTO
, "-c", "confout", ffWRITE
},
2178 { efEDR
, "-e", "ener", ffWRITE
},
2179 { efLOG
, "-g", "md", ffWRITE
},
2180 { efXVG
, "-dhdl", "dhdl", ffOPTWR
},
2181 { efXVG
, "-field", "field", ffOPTWR
},
2182 { efXVG
, "-table", "table", ffOPTRD
},
2183 { efXVG
, "-tablep", "tablep", ffOPTRD
},
2184 { efXVG
, "-tableb", "table", ffOPTRD
},
2185 { efTRX
, "-rerun", "rerun", ffOPTRD
},
2186 { efXVG
, "-tpi", "tpi", ffOPTWR
},
2187 { efXVG
, "-tpid", "tpidist", ffOPTWR
},
2188 { efEDI
, "-ei", "sam", ffOPTRD
},
2189 { efXVG
, "-eo", "edsam", ffOPTWR
},
2190 { efXVG
, "-devout", "deviatie", ffOPTWR
},
2191 { efXVG
, "-runav", "runaver", ffOPTWR
},
2192 { efXVG
, "-px", "pullx", ffOPTWR
},
2193 { efXVG
, "-pf", "pullf", ffOPTWR
},
2194 { efXVG
, "-ro", "rotation", ffOPTWR
},
2195 { efLOG
, "-ra", "rotangles", ffOPTWR
},
2196 { efLOG
, "-rs", "rotslabs", ffOPTWR
},
2197 { efLOG
, "-rt", "rottorque", ffOPTWR
},
2198 { efMTX
, "-mtx", "nm", ffOPTWR
},
2199 { efXVG
, "-swap", "swapions", ffOPTWR
},
2200 /* Output files that are deleted after each benchmark run */
2201 { efTRN
, "-bo", "bench", ffWRITE
},
2202 { efXTC
, "-bx", "bench", ffWRITE
},
2203 { efCPT
, "-bcpo", "bench", ffWRITE
},
2204 { efSTO
, "-bc", "bench", ffWRITE
},
2205 { efEDR
, "-be", "bench", ffWRITE
},
2206 { efLOG
, "-bg", "bench", ffWRITE
},
2207 { efXVG
, "-beo", "benchedo", ffOPTWR
},
2208 { efXVG
, "-bdhdl", "benchdhdl", ffOPTWR
},
2209 { efXVG
, "-bfield", "benchfld", ffOPTWR
},
2210 { efXVG
, "-btpi", "benchtpi", ffOPTWR
},
2211 { efXVG
, "-btpid", "benchtpid", ffOPTWR
},
2212 { efXVG
, "-bdevout", "benchdev", ffOPTWR
},
2213 { efXVG
, "-brunav", "benchrnav", ffOPTWR
},
2214 { efXVG
, "-bpx", "benchpx", ffOPTWR
},
2215 { efXVG
, "-bpf", "benchpf", ffOPTWR
},
2216 { efXVG
, "-bro", "benchrot", ffOPTWR
},
2217 { efLOG
, "-bra", "benchrota", ffOPTWR
},
2218 { efLOG
, "-brs", "benchrots", ffOPTWR
},
2219 { efLOG
, "-brt", "benchrott", ffOPTWR
},
2220 { efMTX
, "-bmtx", "benchn", ffOPTWR
},
2221 { efNDX
, "-bdn", "bench", ffOPTWR
},
2222 { efXVG
, "-bswap", "benchswp", ffOPTWR
}
2225 gmx_bool bThreads
= FALSE
;
2229 const char *procstring
[] =
2230 { nullptr, "np", "n", "none", nullptr };
2231 const char *npmevalues_opt
[] =
2232 { nullptr, "auto", "all", "subset", nullptr };
2234 gmx_bool bAppendFiles
= TRUE
;
2235 gmx_bool bKeepAndNumCPT
= FALSE
;
2236 gmx_bool bResetCountersHalfWay
= FALSE
;
2237 gmx_bool bBenchmark
= TRUE
;
2238 gmx_bool bCheck
= TRUE
;
2240 gmx_output_env_t
*oenv
= nullptr;
2243 /***********************/
2244 /* tune_pme options: */
2245 /***********************/
2246 { "-mdrun", FALSE
, etSTR
, {&cmd_mdrun
},
2247 "Command line to run a simulation, e.g. 'gmx mdrun' or 'mdrun_mpi'" },
2248 { "-np", FALSE
, etINT
, {&nnodes
},
2249 "Number of ranks to run the tests on (must be > 2 for separate PME ranks)" },
2250 { "-npstring", FALSE
, etENUM
, {procstring
},
2251 "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)" },
2252 { "-ntmpi", FALSE
, etINT
, {&nthreads
},
2253 "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2254 { "-r", FALSE
, etINT
, {&repeats
},
2255 "Repeat each test this often" },
2256 { "-max", FALSE
, etREAL
, {&maxPMEfraction
},
2257 "Max fraction of PME ranks to test with" },
2258 { "-min", FALSE
, etREAL
, {&minPMEfraction
},
2259 "Min fraction of PME ranks to test with" },
2260 { "-npme", FALSE
, etENUM
, {npmevalues_opt
},
2261 "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2262 "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2263 { "-fix", FALSE
, etINT
, {&npme_fixed
},
2264 "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."},
2265 { "-rmax", FALSE
, etREAL
, {&rmax
},
2266 "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2267 { "-rmin", FALSE
, etREAL
, {&rmin
},
2268 "If >0, minimal rcoulomb for -ntpr>1" },
2269 { "-scalevdw", FALSE
, etBOOL
, {&bScaleRvdw
},
2270 "Scale rvdw along with rcoulomb"},
2271 { "-ntpr", FALSE
, etINT
, {&ntprs
},
2272 "Number of [REF].tpr[ref] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2273 "If < 1, automatically choose the number of [REF].tpr[ref] files to test" },
2274 { "-steps", FALSE
, etINT64
, {&bench_nsteps
},
2275 "Take timings for this many steps in the benchmark runs" },
2276 { "-resetstep", FALSE
, etINT
, {&presteps
},
2277 "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2278 { "-nsteps", FALSE
, etINT64
, {&new_sim_nsteps
},
2279 "If non-negative, perform this many steps in the real run (overwrites nsteps from [REF].tpr[ref], add [REF].cpt[ref] steps)" },
2280 { "-launch", FALSE
, etBOOL
, {&bLaunch
},
2281 "Launch the real simulation after optimization" },
2282 { "-bench", FALSE
, etBOOL
, {&bBenchmark
},
2283 "Run the benchmarks or just create the input [REF].tpr[ref] files?" },
2284 { "-check", FALSE
, etBOOL
, {&bCheck
},
2285 "Before the benchmark runs, check whether mdrun works in parallel" },
2286 { "-gpu_id", FALSE
, etSTR
, {&eligible_gpu_ids
},
2287 "List of unique GPU device IDs that are eligible for use" },
2288 /******************/
2289 /* mdrun options: */
2290 /******************/
2291 /* We let tune_pme parse and understand these options, because we need to
2292 * prevent that they appear on the mdrun command line for the benchmarks */
2293 { "-append", FALSE
, etBOOL
, {&bAppendFiles
},
2294 "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2295 { "-cpnum", FALSE
, etBOOL
, {&bKeepAndNumCPT
},
2296 "Keep and number checkpoint files (launch only)" },
2297 { "-deffnm", FALSE
, etSTR
, {&deffnm
},
2298 "Set the default filenames (launch only)" },
2299 { "-resethway", FALSE
, etBOOL
, {&bResetCountersHalfWay
},
2300 "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2303 #define NFILE asize(fnm)
2305 seconds
= gmx_gettime();
2307 if (!parse_common_args(&argc
, argv
, PCA_NOEXIT_ON_ARGS
,
2308 NFILE
, fnm
, asize(pa
), pa
, asize(desc
), desc
,
2314 // procstring[0]Â is used inside two different conditionals further down
2315 GMX_RELEASE_ASSERT(procstring
[0] != nullptr, "Options inconsistency; procstring[0]Â is NULL");
2317 /* Store the remaining unparsed command line entries in a string which
2318 * is then attached to the mdrun command line */
2320 ExtraArgs
[0] = '\0';
2321 for (i
= 1; i
< argc
; i
++) /* argc will now be 1 if everything was understood */
2323 add_to_string(&ExtraArgs
, argv
[i
]);
2324 add_to_string(&ExtraArgs
, " ");
2327 if (opt2parg_bSet("-ntmpi", asize(pa
), pa
))
2330 if (opt2parg_bSet("-npstring", asize(pa
), pa
))
2332 fprintf(stderr
, "WARNING: -npstring has no effect when using threads.\n");
2337 gmx_fatal(FARGS
, "Can't run multi-threaded MPI simulation yet!");
2339 /* and now we just set this; a bit of an ugly hack*/
2342 /* Check for PME:PP ratio and whether tpr triggers additional output files */
2343 guessPMEratio
= inspect_tpr(NFILE
, fnm
, &rcoulomb
);
2345 /* Automatically set -beo options if -eo is set etc. */
2346 couple_files_options(NFILE
, fnm
);
2348 /* Construct the command line arguments for benchmark runs
2349 * as well as for the simulation run */
2352 sprintf(bbuf
, " -ntmpi %d ", nthreads
);
2356 /* This string will be used for MPI runs and will appear after the
2357 * mpirun command. */
2358 if (std::strcmp(procstring
[0], "none") != 0)
2360 sprintf(bbuf
, " -%s %d ", procstring
[0], nnodes
);
2370 create_command_line_snippets(bAppendFiles
, bKeepAndNumCPT
, bResetCountersHalfWay
, presteps
,
2371 NFILE
, fnm
, &cmd_args_bench
, &cmd_args_launch
, ExtraArgs
, deffnm
);
2373 /* Prepare to use checkpoint file if requested */
2375 if (opt2bSet("-cpi", NFILE
, fnm
))
2377 const char *filename
= opt2fn("-cpi", NFILE
, fnm
);
2379 read_checkpoint_part_and_step(filename
,
2380 &cpt_sim_part
, &cpt_steps
);
2381 if (cpt_sim_part
== 0)
2383 gmx_fatal(FARGS
, "Checkpoint file %s could not be read!", filename
);
2385 /* tune_pme will run the next part of the simulation */
2386 sim_part
= cpt_sim_part
+ 1;
2389 /* Open performance output file and write header info */
2390 fp
= gmx_ffopen(opt2fn("-p", NFILE
, fnm
), "w");
2392 /* Make a quick consistency check of command line parameters */
2393 check_input(nnodes
, repeats
, &ntprs
, &rmin
, rcoulomb
, &rmax
,
2394 maxPMEfraction
, minPMEfraction
, npme_fixed
,
2395 bench_nsteps
, fnm
, NFILE
, sim_part
, presteps
,
2398 /* Determine the maximum and minimum number of PME nodes to test,
2399 * the actual list of settings is build in do_the_tests(). */
2400 if ((nnodes
> 2) && (npme_fixed
< -1))
2402 if (0 == std::strcmp(npmevalues_opt
[0], "auto"))
2404 /* Determine the npme range automatically based on the PME:PP load guess */
2405 if (guessPMEratio
> 1.0)
2407 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2408 maxPMEnodes
= nnodes
/2;
2409 minPMEnodes
= nnodes
/2;
2413 /* PME : PP load is in the range 0..1, let's test around the guess */
2414 guessPMEnodes
= static_cast<int>(nnodes
/(1.0 + 1.0/guessPMEratio
));
2415 minPMEnodes
= static_cast<int>(std::floor(0.7*guessPMEnodes
));
2416 maxPMEnodes
= static_cast<int>(std::ceil(1.6*guessPMEnodes
));
2417 maxPMEnodes
= std::min(maxPMEnodes
, nnodes
/2);
2422 /* Determine the npme range based on user input */
2423 maxPMEnodes
= static_cast<int>(std::floor(maxPMEfraction
*nnodes
));
2424 minPMEnodes
= std::max(static_cast<int>(std::floor(minPMEfraction
*nnodes
)), 0);
2425 fprintf(stdout
, "Will try runs with %d ", minPMEnodes
);
2426 if (maxPMEnodes
!= minPMEnodes
)
2428 fprintf(stdout
, "- %d ", maxPMEnodes
);
2430 fprintf(stdout
, "PME-only ranks.\n Note that the automatic number of PME-only ranks and no separate PME ranks are always tested.\n");
2439 /* Get the commands we need to set up the runs from environment variables */
2440 get_program_paths(bThreads
, &cmd_mpirun
, &cmd_mdrun
);
2441 if (bBenchmark
&& repeats
> 0)
2443 check_mdrun_works(bThreads
, cmd_mpirun
, cmd_np
, cmd_mdrun
, nullptr != eligible_gpu_ids
);
2446 /* Print some header info to file */
2448 fprintf(fp
, "\n P E R F O R M A N C E R E S U L T S\n");
2450 fprintf(fp
, "%s for GROMACS %s\n", output_env_get_program_display_name(oenv
),
2454 fprintf(fp
, "Number of ranks : %d\n", nnodes
);
2455 fprintf(fp
, "The mpirun command is : %s\n", cmd_mpirun
);
2456 if (std::strcmp(procstring
[0], "none") != 0)
2458 fprintf(fp
, "Passing # of ranks via : -%s\n", procstring
[0]);
2462 fprintf(fp
, "Not setting number of ranks in system call\n");
2467 fprintf(fp
, "Number of threads : %d\n", nnodes
);
2470 fprintf(fp
, "The mdrun command is : %s\n", cmd_mdrun
);
2471 fprintf(fp
, "mdrun args benchmarks : %s\n", cmd_args_bench
);
2472 fprintf(fp
, "Benchmark steps : ");
2473 fprintf(fp
, "%" PRId64
, bench_nsteps
);
2475 fprintf(fp
, "dlb equilibration steps : %d\n", presteps
);
2478 fprintf(fp
, "Checkpoint time step : ");
2479 fprintf(fp
, "%" PRId64
, cpt_steps
);
2482 fprintf(fp
, "mdrun args at launchtime: %s\n", cmd_args_launch
);
2484 if (new_sim_nsteps
>= 0)
2487 fprintf(stderr
, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE
, fnm
));
2488 fprintf(stderr
, "%" PRId64
, new_sim_nsteps
+cpt_steps
);
2489 fprintf(stderr
, " steps.\n");
2490 fprintf(fp
, "Simulation steps : ");
2491 fprintf(fp
, "%" PRId64
, new_sim_nsteps
);
2496 fprintf(fp
, "Repeats for each test : %d\n", repeats
);
2499 if (npme_fixed
>= -1)
2501 fprintf(fp
, "Fixing -npme at : %d\n", npme_fixed
);
2504 fprintf(fp
, "Input file : %s\n", opt2fn("-s", NFILE
, fnm
));
2505 fprintf(fp
, " PME/PP load estimate : %g\n", guessPMEratio
);
2507 /* Allocate memory for the inputinfo struct: */
2509 info
->nr_inputfiles
= ntprs
;
2510 for (i
= 0; i
< ntprs
; i
++)
2512 snew(info
->rcoulomb
, ntprs
);
2513 snew(info
->rvdw
, ntprs
);
2514 snew(info
->rlist
, ntprs
);
2515 snew(info
->nkx
, ntprs
);
2516 snew(info
->nky
, ntprs
);
2517 snew(info
->nkz
, ntprs
);
2518 snew(info
->fsx
, ntprs
);
2519 snew(info
->fsy
, ntprs
);
2520 snew(info
->fsz
, ntprs
);
2522 /* Make alternative tpr files to test: */
2523 snew(tpr_names
, ntprs
);
2524 for (i
= 0; i
< ntprs
; i
++)
2526 snew(tpr_names
[i
], STRLEN
);
2529 /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2530 * different grids could be found. */
2531 make_benchmark_tprs(opt2fn("-s", NFILE
, fnm
), tpr_names
, bench_nsteps
+presteps
,
2532 cpt_steps
, rmin
, rmax
, bScaleRvdw
, &ntprs
, info
, fp
);
2534 /********************************************************************************/
2535 /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats */
2536 /********************************************************************************/
2537 snew(perfdata
, ntprs
);
2540 GMX_RELEASE_ASSERT(npmevalues_opt
[0] != nullptr, "Options inconsistency; npmevalues_opt[0] is NULL");
2541 do_the_tests(fp
, tpr_names
, maxPMEnodes
, minPMEnodes
, npme_fixed
, npmevalues_opt
[0], perfdata
, &pmeentries
,
2542 repeats
, nnodes
, ntprs
, bThreads
, cmd_mpirun
, cmd_np
, cmd_mdrun
,
2543 cmd_args_bench
, fnm
, NFILE
, presteps
, cpt_steps
, bCheck
, eligible_gpu_ids
);
2545 fprintf(fp
, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds
)/60.0);
2547 /* Analyse the results and give a suggestion for optimal settings: */
2548 bKeepTPR
= analyze_data(fp
, opt2fn("-p", NFILE
, fnm
), perfdata
, nnodes
, ntprs
, pmeentries
,
2549 repeats
, info
, &best_tpr
, &best_npme
);
2551 /* Take the best-performing tpr file and enlarge nsteps to original value */
2552 if (bKeepTPR
&& !bOverwrite
)
2554 simulation_tpr
= opt2fn("-s", NFILE
, fnm
);
2558 simulation_tpr
= opt2fn("-so", NFILE
, fnm
);
2559 modify_PMEsettings(bOverwrite
? (new_sim_nsteps
+cpt_steps
) : info
->orig_sim_steps
,
2560 info
->orig_init_step
, tpr_names
[best_tpr
], simulation_tpr
);
2563 /* Let's get rid of the temporary benchmark input files */
2564 for (i
= 0; i
< ntprs
; i
++)
2566 fprintf(stdout
, "Deleting temporary benchmark input file %s\n", tpr_names
[i
]);
2567 remove(tpr_names
[i
]);
2570 /* Now start the real simulation if the user requested it ... */
2571 launch_simulation(bLaunch
, fp
, bThreads
, cmd_mpirun
, cmd_np
, cmd_mdrun
,
2572 cmd_args_launch
, simulation_tpr
, best_npme
, eligible_gpu_ids
);
2576 /* ... or simply print the performance results to screen: */
2579 finalize(opt2fn("-p", NFILE
, fnm
));