9947077818c4fe4f292ead1eb9180198d8f9bb10
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_tune_pme.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "config.h"
38
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
42 #include <ctime>
43
44 #include <algorithm>
45 #include <string>
46
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
50
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/smalloc.h"
75 #include "gromacs/utility/stringutil.h"
76
77 /* Enum for situations that can occur during log file parsing, the
78  * corresponding string entries can be found in do_the_tests() in
79  * const char* ParseLog[] */
80 /* TODO clean up CamelCasing of these enum names */
81 enum {
82     eParselogOK,
83     eParselogNotFound,
84     eParselogNoPerfData,
85     eParselogTerm,
86     eParselogResetProblem,
87     eParselogNoDDGrid,
88     eParselogTPXVersion,
89     eParselogNotParallel,
90     eParselogLargePrimeFactor,
91     eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs,
92     eParselogGpuProblem,
93     eParselogFatal,
94     eParselogNr
95 };
96
97
98 typedef struct
99 {
100     int     nPMEnodes;    /* number of PME-only nodes used in this test */
101     int     nx, ny, nz;   /* DD grid */
102     int     guessPME;     /* if nPMEnodes == -1, this is the guessed number of PME nodes */
103     double *Gcycles;      /* This can contain more than one value if doing multiple tests */
104     double  Gcycles_Av;
105     float  *ns_per_day;
106     float   ns_per_day_Av;
107     float  *PME_f_load;     /* PME mesh/force load average*/
108     float   PME_f_load_Av;  /* Average average ;) ... */
109     char   *mdrun_cmd_line; /* Mdrun command line used for this test */
110 } t_perf;
111
112
113 typedef struct
114 {
115     int             nr_inputfiles;   /* The number of tpr and mdp input files */
116     gmx_int64_t     orig_sim_steps;  /* Number of steps to be done in the real simulation */
117     gmx_int64_t     orig_init_step;  /* Init step for the real simulation */
118     real           *rcoulomb;        /* The coulomb radii [0...nr_inputfiles] */
119     real           *rvdw;            /* The vdW radii */
120     real           *rlist;           /* Neighbourlist cutoff radius */
121     int            *nkx, *nky, *nkz;
122     real           *fsx, *fsy, *fsz; /* Fourierspacing in x,y,z dimension */
123 } t_inputinfo;
124
125
126 static void sep_line(FILE *fp)
127 {
128     fprintf(fp, "\n------------------------------------------------------------\n");
129 }
130
131
132 /* Wrapper for system calls */
133 static int gmx_system_call(char *command)
134 {
135     return ( system(command) );
136 }
137
138
139 /* Check if string starts with substring */
140 static gmx_bool str_starts(const char *string, const char *substring)
141 {
142     return ( std::strncmp(string, substring, std::strlen(substring)) == 0);
143 }
144
145
146 static void cleandata(t_perf *perfdata, int test_nr)
147 {
148     perfdata->Gcycles[test_nr]    = 0.0;
149     perfdata->ns_per_day[test_nr] = 0.0;
150     perfdata->PME_f_load[test_nr] = 0.0;
151 }
152
153
154 static void remove_if_exists(const char *fn)
155 {
156     if (gmx_fexist(fn))
157     {
158         fprintf(stdout, "Deleting %s\n", fn);
159         remove(fn);
160     }
161 }
162
163
164 static void finalize(const char *fn_out)
165 {
166     char  buf[STRLEN];
167     FILE *fp;
168
169
170     fp = fopen(fn_out, "r");
171     fprintf(stdout, "\n\n");
172
173     while (fgets(buf, STRLEN-1, fp) != nullptr)
174     {
175         fprintf(stdout, "%s", buf);
176     }
177     fclose(fp);
178     fprintf(stdout, "\n\n");
179 }
180
181
182 enum {
183     eFoundNothing, eFoundDDStr, eFoundAccountingStr, eFoundCycleStr
184 };
185
186 static int parse_logfile(const char *logfile, const char *errfile,
187                          t_perf *perfdata, int test_nr, int presteps, gmx_int64_t cpt_steps,
188                          int nnodes)
189 {
190     FILE           *fp;
191     char            line[STRLEN], dumstring[STRLEN], dumstring2[STRLEN];
192     const char      matchstrdd[]  = "Domain decomposition grid";
193     const char      matchstrcr[]  = "resetting all time and cycle counters";
194     const char      matchstrbal[] = "Average PME mesh/force load:";
195     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";
196     const char      errSIG[]      = "signal, stopping at the next";
197     int             iFound;
198     float           dum1, dum2, dum3, dum4;
199     int             ndum;
200     int             npme;
201     gmx_int64_t     resetsteps     = -1;
202     gmx_bool        bFoundResetStr = FALSE;
203     gmx_bool        bResetChecked  = FALSE;
204
205
206     if (!gmx_fexist(logfile))
207     {
208         fprintf(stderr, "WARNING: Could not find logfile %s.\n", logfile);
209         cleandata(perfdata, test_nr);
210         return eParselogNotFound;
211     }
212
213     fp = fopen(logfile, "r");
214     perfdata->PME_f_load[test_nr] = -1.0;
215     perfdata->guessPME            = -1;
216
217     iFound = eFoundNothing;
218     if (1 == nnodes)
219     {
220         iFound = eFoundDDStr; /* Skip some case statements */
221     }
222
223     while (fgets(line, STRLEN, fp) != nullptr)
224     {
225         /* Remove leading spaces */
226         ltrim(line);
227
228         /* Check for TERM and INT signals from user: */
229         if (std::strstr(line, errSIG) != nullptr)
230         {
231             fclose(fp);
232             cleandata(perfdata, test_nr);
233             return eParselogTerm;
234         }
235
236         /* Check whether cycle resetting  worked */
237         if (presteps > 0 && !bFoundResetStr)
238         {
239             if (std::strstr(line, matchstrcr) != nullptr)
240             {
241                 sprintf(dumstring, "step %s", "%" GMX_SCNd64);
242                 sscanf(line, dumstring, &resetsteps);
243                 bFoundResetStr = TRUE;
244                 if (resetsteps == presteps+cpt_steps)
245                 {
246                     bResetChecked = TRUE;
247                 }
248                 else
249                 {
250                     sprintf(dumstring, "%" GMX_PRId64, resetsteps);
251                     sprintf(dumstring2, "%" GMX_PRId64, presteps+cpt_steps);
252                     fprintf(stderr, "WARNING: Time step counters were reset at step %s,\n"
253                             "         though they were supposed to be reset at step %s!\n",
254                             dumstring, dumstring2);
255                 }
256             }
257         }
258
259         /* Look for strings that appear in a certain order in the log file: */
260         switch (iFound)
261         {
262             case eFoundNothing:
263                 /* Look for domain decomp grid and separate PME nodes: */
264                 if (str_starts(line, matchstrdd))
265                 {
266                     sscanf(line, "Domain decomposition grid %d x %d x %d, separate PME ranks %d",
267                            &(perfdata->nx), &(perfdata->ny), &(perfdata->nz), &npme);
268                     if (perfdata->nPMEnodes == -1)
269                     {
270                         perfdata->guessPME = npme;
271                     }
272                     else if (perfdata->nPMEnodes != npme)
273                     {
274                         gmx_fatal(FARGS, "PME ranks from command line and output file are not identical");
275                     }
276                     iFound = eFoundDDStr;
277                 }
278                 /* Catch a few errors that might have occurred: */
279                 else if (str_starts(line, "There is no domain decomposition for"))
280                 {
281                     fclose(fp);
282                     return eParselogNoDDGrid;
283                 }
284                 else if (str_starts(line, "The number of ranks you selected"))
285                 {
286                     fclose(fp);
287                     return eParselogLargePrimeFactor;
288                 }
289                 else if (str_starts(line, "reading tpx file"))
290                 {
291                     fclose(fp);
292                     return eParselogTPXVersion;
293                 }
294                 else if (str_starts(line, "The -dd or -npme option request a parallel simulation"))
295                 {
296                     fclose(fp);
297                     return eParselogNotParallel;
298                 }
299                 break;
300             case eFoundDDStr:
301                 /* Even after the "Domain decomposition grid" string was found,
302                  * it could be that mdrun had to quit due to some error. */
303                 if (str_starts(line, "Incorrect launch configuration: mismatching number of"))
304                 {
305                     fclose(fp);
306                     return eParselogMismatchOfNumberOfPPRanksAndAvailableGPUs;
307                 }
308                 else if (str_starts(line, "Some of the requested GPUs do not exist"))
309                 {
310                     fclose(fp);
311                     return eParselogGpuProblem;
312                 }
313                 /* Look for PME mesh/force balance (not necessarily present, though) */
314                 else if (str_starts(line, matchstrbal))
315                 {
316                     sscanf(&line[std::strlen(matchstrbal)], "%f", &(perfdata->PME_f_load[test_nr]));
317                 }
318                 /* Look for matchstring */
319                 else if (str_starts(line, matchstring))
320                 {
321                     iFound = eFoundAccountingStr;
322                 }
323                 break;
324             case eFoundAccountingStr:
325                 /* Already found matchstring - look for cycle data */
326                 if (str_starts(line, "Total  "))
327                 {
328                     sscanf(line, "Total %*f %lf", &(perfdata->Gcycles[test_nr]));
329                     iFound = eFoundCycleStr;
330                 }
331                 break;
332             case eFoundCycleStr:
333                 /* Already found cycle data - look for remaining performance info and return */
334                 if (str_starts(line, "Performance:"))
335                 {
336                     ndum = sscanf(line, "%s %f %f %f %f", dumstring, &dum1, &dum2, &dum3, &dum4);
337                     /* (ns/day) is the second last entry, depending on whether GMX_DETAILED_PERF_STATS was set in print_perf(), nrnb.c */
338                     perfdata->ns_per_day[test_nr] = (ndum == 5) ? dum3 : dum1;
339                     fclose(fp);
340                     if (bResetChecked || presteps == 0)
341                     {
342                         return eParselogOK;
343                     }
344                     else
345                     {
346                         return eParselogResetProblem;
347                     }
348                 }
349                 break;
350         }
351     } /* while */
352
353     /* Close the log file */
354     fclose(fp);
355
356     /* Check why there is no performance data in the log file.
357      * Did a fatal errors occur? */
358     if (gmx_fexist(errfile))
359     {
360         fp = fopen(errfile, "r");
361         while (fgets(line, STRLEN, fp) != nullptr)
362         {
363             if (str_starts(line, "Fatal error:") )
364             {
365                 if (fgets(line, STRLEN, fp) != nullptr)
366                 {
367                     fprintf(stderr, "\nWARNING: An error occurred during this benchmark:\n"
368                             "%s\n", line);
369                 }
370                 fclose(fp);
371                 cleandata(perfdata, test_nr);
372                 return eParselogFatal;
373             }
374         }
375         fclose(fp);
376     }
377     else
378     {
379         fprintf(stderr, "WARNING: Could not find stderr file %s.\n", errfile);
380     }
381
382     /* Giving up ... we could not find out why there is no performance data in
383      * the log file. */
384     fprintf(stdout, "No performance data in log file.\n");
385     cleandata(perfdata, test_nr);
386
387     return eParselogNoPerfData;
388 }
389
390
391 static gmx_bool analyze_data(
392         FILE         *fp,
393         const char   *fn,
394         t_perf      **perfdata,
395         int           nnodes,
396         int           ntprs,
397         int           ntests,
398         int           nrepeats,
399         t_inputinfo  *info,
400         int          *index_tpr,    /* OUT: Nr of mdp file with best settings */
401         int          *npme_optimal) /* OUT: Optimal number of PME nodes */
402 {
403     int      i, j, k;
404     int      line  = 0, line_win = -1;
405     int      k_win = -1, i_win = -1, winPME;
406     double   s     = 0.0; /* standard deviation */
407     t_perf  *pd;
408     char     strbuf[STRLEN];
409     char     str_PME_f_load[13];
410     gmx_bool bCanUseOrigTPR;
411     gmx_bool bRefinedCoul, bRefinedVdW, bRefinedGrid;
412
413
414     if (nrepeats > 1)
415     {
416         sep_line(fp);
417         fprintf(fp, "Summary of successful runs:\n");
418         fprintf(fp, "Line tpr PME ranks  Gcycles Av.     Std.dev.       ns/day        PME/f");
419         if (nnodes > 1)
420         {
421             fprintf(fp, "    DD grid");
422         }
423         fprintf(fp, "\n");
424     }
425
426
427     for (k = 0; k < ntprs; k++)
428     {
429         for (i = 0; i < ntests; i++)
430         {
431             /* Select the right dataset: */
432             pd = &(perfdata[k][i]);
433
434             pd->Gcycles_Av    = 0.0;
435             pd->PME_f_load_Av = 0.0;
436             pd->ns_per_day_Av = 0.0;
437
438             if (pd->nPMEnodes == -1)
439             {
440                 sprintf(strbuf, "(%3d)", pd->guessPME);
441             }
442             else
443             {
444                 sprintf(strbuf, "     ");
445             }
446
447             /* Get the average run time of a setting */
448             for (j = 0; j < nrepeats; j++)
449             {
450                 pd->Gcycles_Av    += pd->Gcycles[j];
451                 pd->PME_f_load_Av += pd->PME_f_load[j];
452             }
453             pd->Gcycles_Av    /= nrepeats;
454             pd->PME_f_load_Av /= nrepeats;
455
456             for (j = 0; j < nrepeats; j++)
457             {
458                 if (pd->ns_per_day[j] > 0.0)
459                 {
460                     pd->ns_per_day_Av += pd->ns_per_day[j];
461                 }
462                 else
463                 {
464                     /* Somehow the performance number was not aquired for this run,
465                      * therefor set the average to some negative value: */
466                     pd->ns_per_day_Av = -1.0f*nrepeats;
467                     break;
468                 }
469             }
470             pd->ns_per_day_Av /= nrepeats;
471
472             /* Nicer output: */
473             if (pd->PME_f_load_Av > 0.0)
474             {
475                 sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load_Av);
476             }
477             else
478             {
479                 sprintf(str_PME_f_load, "%s", "         -  ");
480             }
481
482
483             /* We assume we had a successful run if both averages are positive */
484             if (pd->Gcycles_Av > 0.0 && pd->ns_per_day_Av > 0.0)
485             {
486                 /* Output statistics if repeats were done */
487                 if (nrepeats > 1)
488                 {
489                     /* Calculate the standard deviation */
490                     s = 0.0;
491                     for (j = 0; j < nrepeats; j++)
492                     {
493                         s += gmx::square( pd->Gcycles[j] - pd->Gcycles_Av );
494                     }
495                     s /= (nrepeats - 1);
496                     s  = std::sqrt(s);
497
498                     fprintf(fp, "%4d %3d %4d%s %12.3f %12.3f %12.3f %s",
499                             line, k, pd->nPMEnodes, strbuf, pd->Gcycles_Av, s,
500                             pd->ns_per_day_Av, str_PME_f_load);
501                     if (nnodes > 1)
502                     {
503                         fprintf(fp, "  %3d %3d %3d", pd->nx, pd->ny, pd->nz);
504                     }
505                     fprintf(fp, "\n");
506                 }
507                 /* Store the index of the best run found so far in 'winner': */
508                 if ( (k_win == -1) || (pd->Gcycles_Av < perfdata[k_win][i_win].Gcycles_Av) )
509                 {
510                     k_win    = k;
511                     i_win    = i;
512                     line_win = line;
513                 }
514                 line++;
515             }
516         }
517     }
518
519     if (k_win == -1)
520     {
521         gmx_fatal(FARGS, "None of the runs was successful! Check %s for problems.", fn);
522     }
523
524     sep_line(fp);
525
526     winPME = perfdata[k_win][i_win].nPMEnodes;
527
528     if (1 == ntests)
529     {
530         /* We stuck to a fixed number of PME-only nodes */
531         sprintf(strbuf, "settings No. %d", k_win);
532     }
533     else
534     {
535         /* We have optimized the number of PME-only nodes */
536         if (winPME == -1)
537         {
538             sprintf(strbuf, "%s", "the automatic number of PME ranks");
539         }
540         else
541         {
542             sprintf(strbuf, "%d PME ranks", winPME);
543         }
544     }
545     fprintf(fp, "Best performance was achieved with %s", strbuf);
546     if ((nrepeats > 1) && (ntests > 1))
547     {
548         fprintf(fp, " (see line %d)", line_win);
549     }
550     fprintf(fp, "\n");
551
552     /* Only mention settings if they were modified: */
553     bRefinedCoul = !gmx_within_tol(info->rcoulomb[k_win], info->rcoulomb[0], GMX_REAL_EPS);
554     bRefinedVdW  = !gmx_within_tol(info->rvdw[k_win], info->rvdw[0], GMX_REAL_EPS);
555     bRefinedGrid = !(info->nkx[k_win] == info->nkx[0] &&
556                      info->nky[k_win] == info->nky[0] &&
557                      info->nkz[k_win] == info->nkz[0]);
558
559     if (bRefinedCoul || bRefinedVdW || bRefinedGrid)
560     {
561         fprintf(fp, "Optimized PME settings:\n");
562         bCanUseOrigTPR = FALSE;
563     }
564     else
565     {
566         bCanUseOrigTPR = TRUE;
567     }
568
569     if (bRefinedCoul)
570     {
571         fprintf(fp, "   New Coulomb radius: %f nm (was %f nm)\n", info->rcoulomb[k_win], info->rcoulomb[0]);
572     }
573
574     if (bRefinedVdW)
575     {
576         fprintf(fp, "   New Van der Waals radius: %f nm (was %f nm)\n", info->rvdw[k_win], info->rvdw[0]);
577     }
578
579     if (bRefinedGrid)
580     {
581         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],
582                 info->nkx[0], info->nky[0], info->nkz[0]);
583     }
584
585     if (bCanUseOrigTPR && ntprs > 1)
586     {
587         fprintf(fp, "and original PME settings.\n");
588     }
589
590     fflush(fp);
591
592     /* Return the index of the mdp file that showed the highest performance
593      * and the optimal number of PME nodes */
594     *index_tpr    = k_win;
595     *npme_optimal = winPME;
596
597     return bCanUseOrigTPR;
598 }
599
600
601 /* Get the commands we need to set up the runs from environment variables */
602 static void get_program_paths(gmx_bool bThreads, char *cmd_mpirun[], char *cmd_mdrun[])
603 {
604     char      *cp;
605     const char def_mpirun[]   = "mpirun";
606
607     const char empty_mpirun[] = "";
608
609     /* Get the commands we need to set up the runs from environment variables */
610     if (!bThreads)
611     {
612         if ( (cp = getenv("MPIRUN")) != nullptr)
613         {
614             *cmd_mpirun = gmx_strdup(cp);
615         }
616         else
617         {
618             *cmd_mpirun = gmx_strdup(def_mpirun);
619         }
620     }
621     else
622     {
623         *cmd_mpirun = gmx_strdup(empty_mpirun);
624     }
625
626     if (*cmd_mdrun == nullptr)
627     {
628         /* The use of MDRUN is deprecated, but made available in 5.1
629            for backward compatibility. It may be removed in a future
630            version. */
631         if ( (cp = getenv("MDRUN" )) != nullptr)
632         {
633             *cmd_mdrun = gmx_strdup(cp);
634         }
635         else
636         {
637             gmx_fatal(FARGS, "The way to call mdrun must be set in the -mdrun command-line flag.");
638         }
639     }
640 }
641
642 /* Check that the commands will run mdrun (perhaps via mpirun) by
643  * running a very quick test simulation. Requires MPI environment or
644  * GPU support to be available if applicable. */
645 /* TODO implement feature to parse the log file to get the list of
646    compatible GPUs from mdrun, if the user of gmx tune-pme has not
647    given one. */
648 static void check_mdrun_works(gmx_bool    bThreads,
649                               const char *cmd_mpirun,
650                               const char *cmd_np,
651                               const char *cmd_mdrun,
652                               gmx_bool    bNeedGpuSupport)
653 {
654     char      *command = nullptr;
655     char      *cp;
656     char       line[STRLEN];
657     FILE      *fp;
658     const char filename[]     = "benchtest.log";
659
660     /* This string should always be identical to the one in copyrite.c,
661      * gmx_print_version_info() in the GMX_MPI section */
662     const char match_mpi[]     = "MPI library:        MPI";
663     const char match_mdrun[]   = "Executable: ";
664     const char match_nogpu[]   = "GPU support:        disabled";
665     gmx_bool   bMdrun          = FALSE;
666     gmx_bool   bMPI            = FALSE;
667     gmx_bool   bHaveGpuSupport = TRUE;
668
669     /* Run a small test to see whether mpirun + mdrun work  */
670     fprintf(stdout, "Making sure that mdrun can be executed. ");
671     if (bThreads)
672     {
673         snew(command, std::strlen(cmd_mdrun) + std::strlen(cmd_np) + std::strlen(filename) + 50);
674         sprintf(command, "%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mdrun, cmd_np, filename);
675     }
676     else
677     {
678         snew(command, std::strlen(cmd_mpirun) + std::strlen(cmd_np) + std::strlen(cmd_mdrun) + std::strlen(filename) + 50);
679         sprintf(command, "%s%s%s -version -maxh 0.001 1> %s 2>&1", cmd_mpirun, cmd_np, cmd_mdrun, filename);
680     }
681     fprintf(stdout, "Trying '%s' ... ", command);
682     make_backup(filename);
683     gmx_system_call(command);
684
685     /* Check if we find the characteristic string in the output: */
686     if (!gmx_fexist(filename))
687     {
688         gmx_fatal(FARGS, "Output from test run could not be found.");
689     }
690
691     fp = fopen(filename, "r");
692     /* We need to scan the whole output file, since sometimes the queuing system
693      * also writes stuff to stdout/err */
694     while (!feof(fp) )
695     {
696         cp = fgets(line, STRLEN, fp);
697         if (cp != nullptr)
698         {
699             if (str_starts(line, match_mdrun) )
700             {
701                 bMdrun = TRUE;
702             }
703             if (str_starts(line, match_mpi) )
704             {
705                 bMPI = TRUE;
706             }
707             if (str_starts(line, match_nogpu) )
708             {
709                 bHaveGpuSupport = FALSE;
710             }
711         }
712     }
713     fclose(fp);
714
715     if (bThreads)
716     {
717         if (bMPI)
718         {
719             gmx_fatal(FARGS, "Need a threaded version of mdrun. This one\n"
720                       "(%s)\n"
721                       "seems to have been compiled with MPI instead.",
722                       cmd_mdrun);
723         }
724     }
725     else
726     {
727         if (bMdrun && !bMPI)
728         {
729             gmx_fatal(FARGS, "Need an MPI-enabled version of mdrun. This one\n"
730                       "(%s)\n"
731                       "seems to have been compiled without MPI support.",
732                       cmd_mdrun);
733         }
734     }
735
736     if (!bMdrun)
737     {
738         gmx_fatal(FARGS, "Cannot execute mdrun. Please check %s for problems!",
739                   filename);
740     }
741
742     if (bNeedGpuSupport && !bHaveGpuSupport)
743     {
744         gmx_fatal(FARGS, "The mdrun executable did not have the expected GPU support.");
745     }
746
747     fprintf(stdout, "passed.\n");
748
749     /* Clean up ... */
750     remove(filename);
751     sfree(command);
752 }
753
754 /* Handles the no-GPU case by emitting an empty string. */
755 static std::string make_gpu_id_command_line(const char *eligible_gpu_ids)
756 {
757     /* If the user has given no eligible GPU IDs, or we're trying the
758      * default behaviour, then there is nothing for tune_pme to give
759      * to mdrun -gpu_id */
760     if (eligible_gpu_ids != nullptr)
761     {
762         return gmx::formatString("-gpu_id %s", eligible_gpu_ids);
763     }
764
765
766     return std::string();
767 }
768
769 static void launch_simulation(
770         gmx_bool    bLaunch,          /* Should the simulation be launched? */
771         FILE       *fp,               /* General log file */
772         gmx_bool    bThreads,         /* whether to use threads */
773         char       *cmd_mpirun,       /* Command for mpirun */
774         char       *cmd_np,           /* Switch for -np or -ntmpi or empty */
775         char       *cmd_mdrun,        /* Command for mdrun */
776         char       *args_for_mdrun,   /* Arguments for mdrun */
777         const char *simulation_tpr,   /* This tpr will be simulated */
778         int         nPMEnodes,        /* Number of PME ranks to use */
779         const char *eligible_gpu_ids) /* Available GPU IDs for
780                                        * constructing mdrun command lines */
781 {
782     char  *command;
783
784
785     /* Make enough space for the system call command,
786      * (200 extra chars for -npme ... etc. options should suffice): */
787     snew(command, std::strlen(cmd_mpirun)+std::strlen(cmd_mdrun)+std::strlen(cmd_np)+std::strlen(args_for_mdrun)+std::strlen(simulation_tpr)+200);
788
789     auto cmd_gpu_ids = make_gpu_id_command_line(eligible_gpu_ids);
790
791     /* Note that the -passall options requires args_for_mdrun to be at the end
792      * of the command line string */
793     if (bThreads)
794     {
795         sprintf(command, "%s%s-npme %d -s %s %s %s",
796                 cmd_mdrun, cmd_np, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids.c_str());
797     }
798     else
799     {
800         sprintf(command, "%s%s%s -npme %d -s %s %s %s",
801                 cmd_mpirun, cmd_np, cmd_mdrun, nPMEnodes, simulation_tpr, args_for_mdrun, cmd_gpu_ids.c_str());
802     }
803
804     fprintf(fp, "%s this command line to launch the simulation:\n\n%s", bLaunch ? "Using" : "Please use", command);
805     sep_line(fp);
806     fflush(fp);
807
808     /* Now the real thing! */
809     if (bLaunch)
810     {
811         fprintf(stdout, "\nLaunching simulation with best parameters now.\nExecuting '%s'", command);
812         sep_line(stdout);
813         fflush(stdout);
814         gmx_system_call(command);
815     }
816 }
817
818
819 static void modify_PMEsettings(
820         gmx_int64_t     simsteps,    /* Set this value as number of time steps */
821         gmx_int64_t     init_step,   /* Set this value as init_step */
822         const char     *fn_best_tpr, /* tpr file with the best performance */
823         const char     *fn_sim_tpr)  /* name of tpr file to be launched */
824 {
825     t_state        state;
826     gmx_mtop_t     mtop;
827     char           buf[200];
828
829     t_inputrec     irInstance;
830     t_inputrec    *ir = &irInstance;
831     read_tpx_state(fn_best_tpr, ir, &state, &mtop);
832
833     /* Reset nsteps and init_step to the value of the input .tpr file */
834     ir->nsteps    = simsteps;
835     ir->init_step = init_step;
836
837     /* Write the tpr file which will be launched */
838     sprintf(buf, "Writing optimized simulation file %s with nsteps=%s.\n", fn_sim_tpr, "%" GMX_PRId64);
839     fprintf(stdout, buf, ir->nsteps);
840     fflush(stdout);
841     write_tpx_state(fn_sim_tpr, ir, &state, &mtop);
842 }
843
844 static gmx_bool can_scale_rvdw(int vdwtype)
845 {
846     return (evdwCUT == vdwtype ||
847             evdwPME == vdwtype);
848 }
849
850 #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
851
852 /* Make additional TPR files with more computational load for the
853  * direct space processors: */
854 static void make_benchmark_tprs(
855         const char     *fn_sim_tpr,      /* READ : User-provided tpr file                 */
856         char           *fn_bench_tprs[], /* WRITE: Names of benchmark tpr files           */
857         gmx_int64_t     benchsteps,      /* Number of time steps for benchmark runs       */
858         gmx_int64_t     statesteps,      /* Step counter in checkpoint file               */
859         real            rmin,            /* Minimal Coulomb radius                        */
860         real            rmax,            /* Maximal Coulomb radius                        */
861         real            bScaleRvdw,      /* Scale rvdw along with rcoulomb                */
862         int            *ntprs,           /* No. of TPRs to write, each with a different
863                                             rcoulomb and fourierspacing                   */
864         t_inputinfo    *info,            /* Contains information about mdp file options   */
865         FILE           *fp)              /* Write the output here                         */
866 {
867     int           i, j, d;
868     t_state       state;
869     gmx_mtop_t    mtop;
870     real          nlist_buffer;     /* Thickness of the buffer regions for PME-switch potentials */
871     char          buf[200];
872     rvec          box_size;
873     gmx_bool      bNote = FALSE;
874     real          add;              /* Add this to rcoul for the next test    */
875     real          fac = 1.0;        /* Scaling factor for Coulomb radius      */
876     real          fourierspacing;   /* Basic fourierspacing from tpr          */
877
878
879     sprintf(buf, "Making benchmark tpr file%s with %s time step%s",
880             *ntprs > 1 ? "s" : "", "%" GMX_PRId64, benchsteps > 1 ? "s" : "");
881     fprintf(stdout, buf, benchsteps);
882     if (statesteps > 0)
883     {
884         sprintf(buf, " (adding %s steps from checkpoint file)", "%" GMX_PRId64);
885         fprintf(stdout, buf, statesteps);
886         benchsteps += statesteps;
887     }
888     fprintf(stdout, ".\n");
889
890     t_inputrec  irInstance;
891     t_inputrec *ir = &irInstance;
892     read_tpx_state(fn_sim_tpr, ir, &state, &mtop);
893
894     /* Check if some kind of PME was chosen */
895     if (EEL_PME(ir->coulombtype) == FALSE)
896     {
897         gmx_fatal(FARGS, "Can only do optimizations for simulations with %s electrostatics.",
898                   EELTYPE(eelPME));
899     }
900
901     /* Check if rcoulomb == rlist, which is necessary for plain PME. */
902     if (  (ir->cutoff_scheme != ecutsVERLET) &&
903           (eelPME == ir->coulombtype) && !(ir->rcoulomb == ir->rlist))
904     {
905         gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to rlist (%f).",
906                   EELTYPE(eelPME), ir->rcoulomb, ir->rlist);
907     }
908     /* For other PME types, rcoulomb is allowed to be smaller than rlist */
909     else if (ir->rcoulomb > ir->rlist)
910     {
911         gmx_fatal(FARGS, "%s requires rcoulomb (%f) to be equal to or smaller than rlist (%f)",
912                   EELTYPE(ir->coulombtype), ir->rcoulomb, ir->rlist);
913     }
914
915     if (bScaleRvdw && ir->rvdw != ir->rcoulomb)
916     {
917         fprintf(stdout, "NOTE: input rvdw != rcoulomb, will not scale rvdw\n");
918         bScaleRvdw = FALSE;
919     }
920
921     /* Reduce the number of steps for the benchmarks */
922     info->orig_sim_steps = ir->nsteps;
923     ir->nsteps           = benchsteps;
924     /* We must not use init_step from the input tpr file for the benchmarks */
925     info->orig_init_step = ir->init_step;
926     ir->init_step        = 0;
927
928     /* For PME-switch potentials, keep the radial distance of the buffer region */
929     nlist_buffer   = ir->rlist - ir->rcoulomb;
930
931     /* Determine length of triclinic box vectors */
932     for (d = 0; d < DIM; d++)
933     {
934         box_size[d] = 0;
935         for (i = 0; i < DIM; i++)
936         {
937             box_size[d] += state.box[d][i]*state.box[d][i];
938         }
939         box_size[d] = std::sqrt(box_size[d]);
940     }
941
942     if (ir->fourier_spacing > 0)
943     {
944         info->fsx[0] = ir->fourier_spacing;
945         info->fsy[0] = ir->fourier_spacing;
946         info->fsz[0] = ir->fourier_spacing;
947     }
948     else
949     {
950         /* Reconstruct fourierspacing per dimension from the number of grid points and box size */
951         info->fsx[0] = box_size[XX]/ir->nkx;
952         info->fsy[0] = box_size[YY]/ir->nky;
953         info->fsz[0] = box_size[ZZ]/ir->nkz;
954     }
955
956     /* If no value for the fourierspacing was provided on the command line, we
957      * use the reconstruction from the tpr file */
958     if (ir->fourier_spacing > 0)
959     {
960         /* Use the spacing from the tpr */
961         fourierspacing = ir->fourier_spacing;
962     }
963     else
964     {
965         /* Use the maximum observed spacing */
966         fourierspacing = std::max(std::max(info->fsx[0], info->fsy[0]), info->fsz[0]);
967     }
968
969     fprintf(stdout, "Calculating PME grid points on the basis of a fourierspacing of %f nm\n", fourierspacing);
970
971     /* For performance comparisons the number of particles is useful to have */
972     fprintf(fp, "   Number of particles  : %d\n", mtop.natoms);
973
974     /* Print information about settings of which some are potentially modified: */
975     fprintf(fp, "   Coulomb type         : %s\n", EELTYPE(ir->coulombtype));
976     fprintf(fp, "   Grid spacing x y z   : %f %f %f\n",
977             box_size[XX]/ir->nkx, box_size[YY]/ir->nky, box_size[ZZ]/ir->nkz);
978     fprintf(fp, "   Van der Waals type   : %s\n", EVDWTYPE(ir->vdwtype));
979     if (ir_vdw_switched(ir))
980     {
981         fprintf(fp, "   rvdw_switch          : %f nm\n", ir->rvdw_switch);
982     }
983     if (EPME_SWITCHED(ir->coulombtype))
984     {
985         fprintf(fp, "   rlist                : %f nm\n", ir->rlist);
986     }
987
988     /* Print a descriptive line about the tpr settings tested */
989     fprintf(fp, "\nWill try these real/reciprocal workload settings:\n");
990     fprintf(fp, " No.   scaling  rcoulomb");
991     fprintf(fp, "  nkx  nky  nkz");
992     fprintf(fp, "   spacing");
993     if (can_scale_rvdw(ir->vdwtype))
994     {
995         fprintf(fp, "      rvdw");
996     }
997     if (EPME_SWITCHED(ir->coulombtype))
998     {
999         fprintf(fp, "     rlist");
1000     }
1001     fprintf(fp, "  tpr file\n");
1002
1003     /* Loop to create the requested number of tpr input files */
1004     for (j = 0; j < *ntprs; j++)
1005     {
1006         /* The first .tpr is the provided one, just need to modify nsteps,
1007          * so skip the following block */
1008         if (j != 0)
1009         {
1010             /* Determine which Coulomb radii rc to use in the benchmarks */
1011             add = (rmax-rmin)/(*ntprs-1);
1012             if (gmx_within_tol(rmin, info->rcoulomb[0], GMX_REAL_EPS))
1013             {
1014                 ir->rcoulomb = rmin + j*add;
1015             }
1016             else if (gmx_within_tol(rmax, info->rcoulomb[0], GMX_REAL_EPS))
1017             {
1018                 ir->rcoulomb = rmin + (j-1)*add;
1019             }
1020             else
1021             {
1022                 /* rmin != rcoul != rmax, ergo test between rmin and rmax */
1023                 add          = (rmax-rmin)/(*ntprs-2);
1024                 ir->rcoulomb = rmin + (j-1)*add;
1025             }
1026
1027             /* Determine the scaling factor fac */
1028             fac = ir->rcoulomb/info->rcoulomb[0];
1029
1030             /* Scale the Fourier grid spacing */
1031             ir->nkx = ir->nky = ir->nkz = 0;
1032             calcFftGrid(nullptr, state.box, fourierspacing*fac, minimalPmeGridSize(ir->pme_order),
1033                         &ir->nkx, &ir->nky, &ir->nkz);
1034
1035             /* Adjust other radii since various conditions need to be fulfilled */
1036             if (eelPME == ir->coulombtype)
1037             {
1038                 /* plain PME, rcoulomb must be equal to rlist TODO only in the group scheme? */
1039                 ir->rlist = ir->rcoulomb;
1040             }
1041             else
1042             {
1043                 /* rlist must be >= rcoulomb, we keep the size of the buffer region */
1044                 ir->rlist = ir->rcoulomb + nlist_buffer;
1045             }
1046
1047             if (bScaleRvdw && can_scale_rvdw(ir->vdwtype))
1048             {
1049                 if (ecutsVERLET == ir->cutoff_scheme ||
1050                     evdwPME == ir->vdwtype)
1051                 {
1052                     /* With either the Verlet cutoff-scheme or LJ-PME,
1053                        the van der Waals radius must always equal the
1054                        Coulomb radius */
1055                     ir->rvdw = ir->rcoulomb;
1056                 }
1057                 else
1058                 {
1059                     /* For vdw cutoff, rvdw >= rlist */
1060                     ir->rvdw = std::max(info->rvdw[0], ir->rlist);
1061                 }
1062             }
1063         } /* end of "if (j != 0)" */
1064
1065         /* for j==0: Save the original settings
1066          * for j >0: Save modified radii and Fourier grids */
1067         info->rcoulomb[j]  = ir->rcoulomb;
1068         info->rvdw[j]      = ir->rvdw;
1069         info->nkx[j]       = ir->nkx;
1070         info->nky[j]       = ir->nky;
1071         info->nkz[j]       = ir->nkz;
1072         info->rlist[j]     = ir->rlist;
1073         info->fsx[j]       = fac*fourierspacing;
1074         info->fsy[j]       = fac*fourierspacing;
1075         info->fsz[j]       = fac*fourierspacing;
1076
1077         /* Write the benchmark tpr file */
1078         std::strncpy(fn_bench_tprs[j], fn_sim_tpr, std::strlen(fn_sim_tpr)-std::strlen(".tpr"));
1079         sprintf(buf, "_bench%.2d.tpr", j);
1080         std::strcat(fn_bench_tprs[j], buf);
1081         fprintf(stdout, "Writing benchmark tpr %s with nsteps=", fn_bench_tprs[j]);
1082         fprintf(stdout, "%" GMX_PRId64, ir->nsteps);
1083         if (j > 0)
1084         {
1085             fprintf(stdout, ", scaling factor %f\n", fac);
1086         }
1087         else
1088         {
1089             fprintf(stdout, ", unmodified settings\n");
1090         }
1091
1092         write_tpx_state(fn_bench_tprs[j], ir, &state, &mtop);
1093
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))
1099         {
1100             fprintf(fp, "%10f", ir->rvdw);
1101         }
1102         if (EPME_SWITCHED(ir->coulombtype))
1103         {
1104             fprintf(fp, "%10f", ir->rlist);
1105         }
1106         fprintf(fp, "  %-14s\n", fn_bench_tprs[j]);
1107
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) )
1111         {
1112             bNote = TRUE;
1113         }
1114     }
1115     if (bNote)
1116     {
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");
1120     }
1121     fflush(stdout);
1122     fflush(fp);
1123 }
1124
1125
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)
1130 {
1131     char        numstring[STRLEN];
1132     char        newfilename[STRLEN];
1133     const char *fn = nullptr;
1134     int         i;
1135     const char *opt;
1136
1137
1138     fprintf(stdout, "Cleaning up, deleting benchmark temp files ...\n");
1139
1140     for (i = 0; i < nfile; i++)
1141     {
1142         opt = (char *)fnm[i].opt;
1143         if (std::strcmp(opt, "-p") == 0)
1144         {
1145             /* do nothing; keep this file */
1146             ;
1147         }
1148         else if (std::strcmp(opt, "-bg") == 0)
1149         {
1150             /* Give the log file a nice name so one can later see which parameters were used */
1151             numstring[0] = '\0';
1152             if (nr > 0)
1153             {
1154                 sprintf(numstring, "_%d", nr);
1155             }
1156             sprintf(newfilename, "%s_no%d_np%d_npme%d%s", opt2fn("-bg", nfile, fnm), k, nnodes, nPMEnodes, numstring);
1157             if (gmx_fexist(opt2fn("-bg", nfile, fnm)))
1158             {
1159                 fprintf(stdout, "renaming log file to %s\n", newfilename);
1160                 make_backup(newfilename);
1161                 rename(opt2fn("-bg", nfile, fnm), newfilename);
1162             }
1163         }
1164         else if (std::strcmp(opt, "-err") == 0)
1165         {
1166             /* This file contains the output of stderr. We want to keep it in
1167              * cases where there have been problems. */
1168             fn           = opt2fn(opt, nfile, fnm);
1169             numstring[0] = '\0';
1170             if (nr > 0)
1171             {
1172                 sprintf(numstring, "_%d", nr);
1173             }
1174             sprintf(newfilename, "%s_no%d_np%d_npme%d%s", fn, k, nnodes, nPMEnodes, numstring);
1175             if (gmx_fexist(fn))
1176             {
1177                 if (bKeepStderr)
1178                 {
1179                     fprintf(stdout, "Saving stderr output in %s\n", newfilename);
1180                     make_backup(newfilename);
1181                     rename(fn, newfilename);
1182                 }
1183                 else
1184                 {
1185                     fprintf(stdout, "Deleting %s\n", fn);
1186                     remove(fn);
1187                 }
1188             }
1189         }
1190         /* Delete the files which are created for each benchmark run: (options -b*) */
1191         else if ( (0 == std::strncmp(opt, "-b", 2)) && (opt2bSet(opt, nfile, fnm) || !is_optional(&fnm[i])) )
1192         {
1193             remove_if_exists(opt2fn(opt, nfile, fnm));
1194         }
1195     }
1196 }
1197
1198
1199 enum {
1200     eNpmeAuto, eNpmeAll, eNpmeReduced, eNpmeSubset, eNpmeNr
1201 };
1202
1203 /* Create a list of numbers of PME nodes to test */
1204 static void make_npme_list(
1205         const char *npmevalues_opt, /* Make a complete list with all
1206                                      * possibilities or a short list that keeps only
1207                                      * reasonable numbers of PME nodes                  */
1208         int        *nentries,       /* Number of entries we put in the nPMEnodes list   */
1209         int        *nPMEnodes[],    /* Each entry contains the value for -npme          */
1210         int         nnodes,         /* Total number of nodes to do the tests on         */
1211         int         minPMEnodes,    /* Minimum number of PME nodes                      */
1212         int         maxPMEnodes)    /* Maximum number of PME nodes                      */
1213 {
1214     int i, npme, npp;
1215     int min_factor = 1;   /* We request that npp and npme have this minimal
1216                            * largest common factor (depends on npp)           */
1217     int nlistmax;         /* Max. list size                                   */
1218     int nlist;            /* Actual number of entries in list                 */
1219     int eNPME = 0;
1220
1221
1222     /* Do we need to check all possible values for -npme or is a reduced list enough? */
1223     if (!std::strcmp(npmevalues_opt, "all") )
1224     {
1225         eNPME = eNpmeAll;
1226     }
1227     else if (!std::strcmp(npmevalues_opt, "subset") )
1228     {
1229         eNPME = eNpmeSubset;
1230     }
1231     else /* "auto" or "range" */
1232     {
1233         if (nnodes <= 64)
1234         {
1235             eNPME = eNpmeAll;
1236         }
1237         else if (nnodes < 128)
1238         {
1239             eNPME = eNpmeReduced;
1240         }
1241         else
1242         {
1243             eNPME = eNpmeSubset;
1244         }
1245     }
1246
1247     /* Calculate how many entries we could possibly have (in case of -npme all) */
1248     if (nnodes > 2)
1249     {
1250         nlistmax = maxPMEnodes - minPMEnodes + 3;
1251         if (0 == minPMEnodes)
1252         {
1253             nlistmax--;
1254         }
1255     }
1256     else
1257     {
1258         nlistmax = 1;
1259     }
1260
1261     /* Now make the actual list which is at most of size nlist */
1262     snew(*nPMEnodes, nlistmax);
1263     nlist = 0; /* start counting again, now the real entries in the list */
1264     for (i = 0; i < nlistmax - 2; i++)
1265     {
1266         npme = maxPMEnodes - i;
1267         npp  = nnodes-npme;
1268         switch (eNPME)
1269         {
1270             case eNpmeAll:
1271                 min_factor = 1;
1272                 break;
1273             case eNpmeReduced:
1274                 min_factor = 2;
1275                 break;
1276             case eNpmeSubset:
1277                 /* For 2d PME we want a common largest factor of at least the cube
1278                  * root of the number of PP nodes */
1279                 min_factor = static_cast<int>(std::cbrt(npp));
1280                 break;
1281             default:
1282                 gmx_fatal(FARGS, "Unknown option for eNPME in make_npme_list");
1283                 break;
1284         }
1285         if (gmx_greatest_common_divisor(npp, npme) >= min_factor)
1286         {
1287             (*nPMEnodes)[nlist] = npme;
1288             nlist++;
1289         }
1290     }
1291     /* We always test 0 PME nodes and the automatic number */
1292     *nentries             = nlist + 2;
1293     (*nPMEnodes)[nlist  ] =  0;
1294     (*nPMEnodes)[nlist+1] = -1;
1295
1296     fprintf(stderr, "Will try the following %d different values for -npme:\n", *nentries);
1297     for (i = 0; i < *nentries-1; i++)
1298     {
1299         fprintf(stderr, "%d, ", (*nPMEnodes)[i]);
1300     }
1301     fprintf(stderr, "and %d (auto).\n", (*nPMEnodes)[*nentries-1]);
1302 }
1303
1304
1305 /* Allocate memory to store the performance data */
1306 static void init_perfdata(t_perf *perfdata[], int ntprs, int datasets, int repeats)
1307 {
1308     int i, j, k;
1309
1310
1311     for (k = 0; k < ntprs; k++)
1312     {
1313         snew(perfdata[k], datasets);
1314         for (i = 0; i < datasets; i++)
1315         {
1316             for (j = 0; j < repeats; j++)
1317             {
1318                 snew(perfdata[k][i].Gcycles, repeats);
1319                 snew(perfdata[k][i].ns_per_day, repeats);
1320                 snew(perfdata[k][i].PME_f_load, repeats);
1321             }
1322         }
1323     }
1324 }
1325
1326
1327 /* Check for errors on mdrun -h */
1328 static void make_sure_it_runs(char *mdrun_cmd_line, int length, FILE *fp,
1329                               const t_filenm *fnm, int nfile)
1330 {
1331     char       *command, *msg;
1332     int         ret;
1333
1334     snew(command, length +  15);
1335     snew(msg, length + 500);
1336
1337     fprintf(stdout, "Making sure the benchmarks can be executed by running just 1 step...\n");
1338     sprintf(command, "%s -nsteps 1 -quiet", mdrun_cmd_line);
1339     fprintf(stdout, "Executing '%s' ...\n", command);
1340     ret = gmx_system_call(command);
1341
1342     if (0 != ret)
1343     {
1344         /* To prevent confusion, do not again issue a gmx_fatal here since we already
1345          * get the error message from mdrun itself */
1346         sprintf(msg,
1347                 "Cannot run the first benchmark simulation! Please check the error message of\n"
1348                 "mdrun for the source of the problem. Did you provide a command line\n"
1349                 "argument that neither gmx tune_pme nor mdrun understands? If you're\n"
1350                 "sure your command line should work, you can bypass this check with \n"
1351                 "gmx tune_pme -nocheck. The failing command was:\n"
1352                 "\n%s\n\n", command);
1353
1354         fprintf(stderr, "%s", msg);
1355         sep_line(fp);
1356         fprintf(fp, "%s", msg);
1357
1358         exit(ret);
1359     }
1360     fprintf(stdout, "Benchmarks can be executed!\n");
1361
1362     /* Clean up the benchmark output files we just created */
1363     fprintf(stdout, "Cleaning up ...\n");
1364     remove_if_exists(opt2fn("-bc", nfile, fnm));
1365     remove_if_exists(opt2fn("-be", nfile, fnm));
1366     remove_if_exists(opt2fn("-bcpo", nfile, fnm));
1367     remove_if_exists(opt2fn("-bg", nfile, fnm));
1368     remove_if_exists(opt2fn("-bo", nfile, fnm));
1369     remove_if_exists(opt2fn("-bx", nfile, fnm));
1370
1371     sfree(command);
1372     sfree(msg    );
1373 }
1374
1375 static void do_the_tests(
1376         FILE           *fp,               /* General tune_pme output file           */
1377         char          **tpr_names,        /* Filenames of the input files to test   */
1378         int             maxPMEnodes,      /* Max fraction of nodes to use for PME   */
1379         int             minPMEnodes,      /* Min fraction of nodes to use for PME   */
1380         int             npme_fixed,       /* If >= -1, test fixed number of PME
1381                                            * nodes only                             */
1382         const char     *npmevalues_opt,   /* Which -npme values should be tested    */
1383         t_perf        **perfdata,         /* Here the performace data is stored     */
1384         int            *pmeentries,       /* Entries in the nPMEnodes list          */
1385         int             repeats,          /* Repeat each test this often            */
1386         int             nnodes,           /* Total number of nodes = nPP + nPME     */
1387         int             nr_tprs,          /* Total number of tpr files to test      */
1388         gmx_bool        bThreads,         /* Threads or MPI?                        */
1389         char           *cmd_mpirun,       /* mpirun command string                  */
1390         char           *cmd_np,           /* "-np", "-n", whatever mpirun needs     */
1391         char           *cmd_mdrun,        /* mdrun command string                   */
1392         char           *cmd_args_bench,   /* arguments for mdrun in a string        */
1393         const t_filenm *fnm,              /* List of filenames from command line    */
1394         int             nfile,            /* Number of files specified on the cmdl. */
1395         int             presteps,         /* DLB equilibration steps, is checked    */
1396         gmx_int64_t     cpt_steps,        /* Time step counter in the checkpoint    */
1397         gmx_bool        bCheck,           /* Check whether benchmark mdrun works    */
1398         const char     *eligible_gpu_ids) /* GPU IDs for
1399                                            * constructing mdrun command lines */
1400 {
1401     int      i, nr, k, ret, count = 0, totaltests;
1402     int     *nPMEnodes = nullptr;
1403     t_perf  *pd        = nullptr;
1404     int      cmdline_length;
1405     char    *command, *cmd_stub;
1406     char     buf[STRLEN];
1407     gmx_bool bResetProblem = FALSE;
1408     gmx_bool bFirst        = TRUE;
1409
1410     /* This string array corresponds to the eParselog enum type at the start
1411      * of this file */
1412     const char* ParseLog[] = {
1413         "OK.",
1414         "Logfile not found!",
1415         "No timings, logfile truncated?",
1416         "Run was terminated.",
1417         "Counters were not reset properly.",
1418         "No DD grid found for these settings.",
1419         "TPX version conflict!",
1420         "mdrun was not started in parallel!",
1421         "Number of PP ranks has a prime factor that is too large.",
1422         "The number of PP ranks did not suit the number of GPUs.",
1423         "Some GPUs were not detected or are incompatible.",
1424         "An error occurred."
1425     };
1426     char        str_PME_f_load[13];
1427
1428
1429     /* Allocate space for the mdrun command line. 100 extra characters should
1430        be more than enough for the -npme etcetera arguments */
1431     cmdline_length =  std::strlen(cmd_mpirun)
1432         + std::strlen(cmd_np)
1433         + std::strlen(cmd_mdrun)
1434         + std::strlen(cmd_args_bench)
1435         + std::strlen(tpr_names[0]) + 100;
1436     snew(command, cmdline_length);
1437     snew(cmd_stub, cmdline_length);
1438
1439     /* Construct the part of the command line that stays the same for all tests: */
1440     if (bThreads)
1441     {
1442         sprintf(cmd_stub, "%s%s", cmd_mdrun, cmd_np);
1443     }
1444     else
1445     {
1446         sprintf(cmd_stub, "%s%s%s ", cmd_mpirun, cmd_np, cmd_mdrun);
1447     }
1448
1449     /* Create a list of numbers of PME nodes to test */
1450     if (npme_fixed < -1)
1451     {
1452         make_npme_list(npmevalues_opt, pmeentries, &nPMEnodes,
1453                        nnodes, minPMEnodes, maxPMEnodes);
1454     }
1455     else
1456     {
1457         *pmeentries  = 1;
1458         snew(nPMEnodes, 1);
1459         nPMEnodes[0] = npme_fixed;
1460         fprintf(stderr, "Will use a fixed number of %d PME-only ranks.\n", nPMEnodes[0]);
1461     }
1462
1463     if (0 == repeats)
1464     {
1465         fprintf(fp, "\nNo benchmarks done since number of repeats (-r) is 0.\n");
1466         gmx_ffclose(fp);
1467         finalize(opt2fn("-p", nfile, fnm));
1468         exit(0);
1469     }
1470
1471     /* Allocate one dataset for each tpr input file: */
1472     init_perfdata(perfdata, nr_tprs, *pmeentries, repeats);
1473
1474     /*****************************************/
1475     /* Main loop over all tpr files to test: */
1476     /*****************************************/
1477     totaltests = nr_tprs*(*pmeentries)*repeats;
1478     for (k = 0; k < nr_tprs; k++)
1479     {
1480         fprintf(fp, "\nIndividual timings for input file %d (%s):\n", k, tpr_names[k]);
1481         fprintf(fp, "PME ranks      Gcycles       ns/day        PME/f    Remark\n");
1482         /* Loop over various numbers of PME nodes: */
1483         for (i = 0; i < *pmeentries; i++)
1484         {
1485             pd = &perfdata[k][i];
1486
1487             auto cmd_gpu_ids = make_gpu_id_command_line(eligible_gpu_ids);
1488
1489             /* Loop over the repeats for each scenario: */
1490             for (nr = 0; nr < repeats; nr++)
1491             {
1492                 pd->nPMEnodes = nPMEnodes[i];
1493
1494                 /* Add -npme and -s to the command line and save it. Note that
1495                  * the -passall (if set) options requires cmd_args_bench to be
1496                  * at the end of the command line string */
1497                 snew(pd->mdrun_cmd_line, cmdline_length);
1498                 sprintf(pd->mdrun_cmd_line, "%s-npme %d -s %s %s %s",
1499                         cmd_stub, pd->nPMEnodes, tpr_names[k], cmd_args_bench, cmd_gpu_ids.c_str());
1500
1501                 /* To prevent that all benchmarks fail due to a show-stopper argument
1502                  * on the mdrun command line, we make a quick check first.
1503                  * This check can be turned off in cases where the automatically chosen
1504                  * number of PME-only ranks leads to a number of PP ranks for which no
1505                  * decomposition can be found (e.g. for large prime numbers) */
1506                 if (bFirst && bCheck)
1507                 {
1508                     /* TODO This check is really for a functional
1509                      * .tpr, and if we need it, it should take place
1510                      * for every .tpr, and the logic for it should be
1511                      * immediately inside the loop over k, not in
1512                      * this inner loop. */
1513                     char *temporary_cmd_line;
1514
1515                     snew(temporary_cmd_line, cmdline_length);
1516                     /* TODO -npme 0 is more likely to succeed at low
1517                        parallelism than the default of -npme -1, but
1518                        is more likely to fail at the scaling limit
1519                        when the PP domains may be too small. "mpirun
1520                        -np 1 mdrun" is probably a reasonable thing to
1521                        do for this check, but it'll be easier to
1522                        implement that after some refactoring of how
1523                        the number of MPI ranks is managed. */
1524                     sprintf(temporary_cmd_line, "%s-npme 0 -nb cpu -s %s %s",
1525                             cmd_stub, tpr_names[k], cmd_args_bench);
1526                     make_sure_it_runs(temporary_cmd_line, cmdline_length, fp, fnm, nfile);
1527                 }
1528                 bFirst = FALSE;
1529
1530                 /* Do a benchmark simulation: */
1531                 if (repeats > 1)
1532                 {
1533                     sprintf(buf, ", pass %d/%d", nr+1, repeats);
1534                 }
1535                 else
1536                 {
1537                     buf[0] = '\0';
1538                 }
1539                 fprintf(stdout, "\n=== Progress %2.0f%%, tpr %d/%d, run %d/%d%s:\n",
1540                         (100.0*count)/totaltests,
1541                         k+1, nr_tprs, i+1, *pmeentries, buf);
1542                 make_backup(opt2fn("-err", nfile, fnm));
1543                 sprintf(command, "%s 1> /dev/null 2>%s", pd->mdrun_cmd_line, opt2fn("-err", nfile, fnm));
1544                 fprintf(stdout, "%s\n", pd->mdrun_cmd_line);
1545                 gmx_system_call(command);
1546
1547                 /* Collect the performance data from the log file; also check stderr
1548                  * for fatal errors */
1549                 ret = parse_logfile(opt2fn("-bg", nfile, fnm), opt2fn("-err", nfile, fnm),
1550                                     pd, nr, presteps, cpt_steps, nnodes);
1551                 if ((presteps > 0) && (ret == eParselogResetProblem))
1552                 {
1553                     bResetProblem = TRUE;
1554                 }
1555
1556                 if (-1 == pd->nPMEnodes)
1557                 {
1558                     sprintf(buf, "(%3d)", pd->guessPME);
1559                 }
1560                 else
1561                 {
1562                     sprintf(buf, "     ");
1563                 }
1564
1565                 /* Nicer output */
1566                 if (pd->PME_f_load[nr] > 0.0)
1567                 {
1568                     sprintf(str_PME_f_load, "%12.3f", pd->PME_f_load[nr]);
1569                 }
1570                 else
1571                 {
1572                     sprintf(str_PME_f_load, "%s", "         -  ");
1573                 }
1574
1575                 /* Write the data we got to disk */
1576                 fprintf(fp, "%4d%s %12.3f %12.3f %s    %s", pd->nPMEnodes,
1577                         buf, pd->Gcycles[nr], pd->ns_per_day[nr], str_PME_f_load, ParseLog[ret]);
1578                 if (!(ret == eParselogOK || ret == eParselogNoDDGrid || ret == eParselogNotFound) )
1579                 {
1580                     fprintf(fp, " Check %s file for problems.", ret == eParselogFatal ? "err" : "log");
1581                 }
1582                 fprintf(fp, "\n");
1583                 fflush(fp);
1584                 count++;
1585
1586                 /* Do some cleaning up and delete the files we do not need any more */
1587                 cleanup(fnm, nfile, k, nnodes, pd->nPMEnodes, nr, ret == eParselogFatal);
1588
1589                 /* If the first run with this number of processors already failed, do not try again: */
1590                 if (pd->Gcycles[0] <= 0.0 && repeats > 1)
1591                 {
1592                     fprintf(stdout, "Skipping remaining passes of unsuccessful setting, see log file for details.\n");
1593                     count += repeats-(nr+1);
1594                     break;
1595                 }
1596             } /* end of repeats loop */
1597         }     /* end of -npme loop */
1598     }         /* end of tpr file loop */
1599
1600     if (bResetProblem)
1601     {
1602         sep_line(fp);
1603         fprintf(fp, "WARNING: The cycle and time step counters could not be reset properly. ");
1604         sep_line(fp);
1605     }
1606     sfree(command);
1607     sfree(cmd_stub);
1608 }
1609
1610
1611 static void check_input(
1612         int             nnodes,
1613         int             repeats,
1614         int            *ntprs,
1615         real           *rmin,
1616         real            rcoulomb,
1617         real           *rmax,
1618         real            maxPMEfraction,
1619         real            minPMEfraction,
1620         int             npme_fixed,
1621         gmx_int64_t     bench_nsteps,
1622         const t_filenm *fnm,
1623         int             nfile,
1624         int             sim_part,
1625         int             presteps,
1626         int             npargs,
1627         t_pargs        *pa)
1628 {
1629     int old;
1630
1631
1632     /* Make sure the input file exists */
1633     if (!gmx_fexist(opt2fn("-s", nfile, fnm)))
1634     {
1635         gmx_fatal(FARGS, "File %s not found.", opt2fn("-s", nfile, fnm));
1636     }
1637
1638     /* Make sure that the checkpoint file is not overwritten during benchmarking */
1639     if ( (0 == std::strcmp(opt2fn("-cpi", nfile, fnm), opt2fn("-bcpo", nfile, fnm)) ) && (sim_part > 1) )
1640     {
1641         gmx_fatal(FARGS, "Checkpoint input (-cpi) and benchmark checkpoint output (-bcpo) files must not be identical.\n"
1642                   "The checkpoint input file must not be overwritten during the benchmarks.\n");
1643     }
1644
1645     /* Make sure that repeats is >= 0 (if == 0, only write tpr files) */
1646     if (repeats < 0)
1647     {
1648         gmx_fatal(FARGS, "Number of repeats < 0!");
1649     }
1650
1651     /* Check number of nodes */
1652     if (nnodes < 1)
1653     {
1654         gmx_fatal(FARGS, "Number of ranks/threads must be a positive integer.");
1655     }
1656
1657     /* Automatically choose -ntpr if not set */
1658     if (*ntprs < 1)
1659     {
1660         if (nnodes < 16)
1661         {
1662             *ntprs = 1;
1663         }
1664         else
1665         {
1666             *ntprs = 3;
1667             /* Set a reasonable scaling factor for rcoulomb */
1668             if (*rmax <= 0)
1669             {
1670                 *rmax = rcoulomb * 1.2;
1671             }
1672         }
1673         fprintf(stderr, "Will test %d tpr file%s.\n", *ntprs, *ntprs == 1 ? "" : "s");
1674     }
1675     else
1676     {
1677         if (1 == *ntprs)
1678         {
1679             fprintf(stderr, "Note: Choose ntpr>1 to shift PME load between real and reciprocal space.\n");
1680         }
1681     }
1682
1683     /* Make shure that rmin <= rcoulomb <= rmax */
1684     if (*rmin <= 0)
1685     {
1686         *rmin = rcoulomb;
1687     }
1688     if (*rmax <= 0)
1689     {
1690         *rmax = rcoulomb;
1691     }
1692     if (!(*rmin <= *rmax) )
1693     {
1694         gmx_fatal(FARGS, "Please choose the Coulomb radii such that rmin <= rmax.\n"
1695                   "rmin = %g, rmax = %g, actual rcoul from .tpr file = %g\n", *rmin, *rmax, rcoulomb);
1696     }
1697     /* Add test scenarios if rmin or rmax were set */
1698     if (*ntprs <= 2)
1699     {
1700         if (!gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1701         {
1702             (*ntprs)++;
1703             fprintf(stderr, "NOTE: Setting -rmin to %g changed -ntpr to %d\n",
1704                     *rmin, *ntprs);
1705         }
1706         if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && (*ntprs == 1) )
1707         {
1708             (*ntprs)++;
1709             fprintf(stderr, "NOTE: Setting -rmax to %g changed -ntpr to %d\n",
1710                     *rmax, *ntprs);
1711         }
1712     }
1713     old = *ntprs;
1714     /* If one of rmin, rmax is set, we need 2 tpr files at minimum */
1715     if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) || !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1716     {
1717         *ntprs = std::max(*ntprs, 2);
1718     }
1719
1720     /* If both rmin, rmax are set, we need 3 tpr files at minimum */
1721     if (!gmx_within_tol(*rmax, rcoulomb, GMX_REAL_EPS) && !gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) )
1722     {
1723         *ntprs = std::max(*ntprs, 3);
1724     }
1725
1726     if (old != *ntprs)
1727     {
1728         fprintf(stderr, "NOTE: Your rmin, rmax setting changed -ntpr to %d\n", *ntprs);
1729     }
1730
1731     if (*ntprs > 1)
1732     {
1733         if (gmx_within_tol(*rmin, rcoulomb, GMX_REAL_EPS) && gmx_within_tol(rcoulomb, *rmax, GMX_REAL_EPS)) /* We have just a single rc */
1734         {
1735             fprintf(stderr, "WARNING: Resetting -ntpr to 1 since no Coulomb radius scaling is requested.\n"
1736                     "Please set rmin < rmax to test Coulomb radii in the [rmin, rmax] interval\n"
1737                     "with correspondingly adjusted PME grid settings\n");
1738             *ntprs = 1;
1739         }
1740     }
1741
1742     /* Check whether max and min fraction are within required values */
1743     if (maxPMEfraction > 0.5 || maxPMEfraction < 0)
1744     {
1745         gmx_fatal(FARGS, "-max must be between 0 and 0.5");
1746     }
1747     if (minPMEfraction > 0.5 || minPMEfraction < 0)
1748     {
1749         gmx_fatal(FARGS, "-min must be between 0 and 0.5");
1750     }
1751     if (maxPMEfraction < minPMEfraction)
1752     {
1753         gmx_fatal(FARGS, "-max must be larger or equal to -min");
1754     }
1755
1756     /* Check whether the number of steps - if it was set - has a reasonable value */
1757     if (bench_nsteps < 0)
1758     {
1759         gmx_fatal(FARGS, "Number of steps must be positive.");
1760     }
1761
1762     if (bench_nsteps > 10000 || bench_nsteps < 100)
1763     {
1764         fprintf(stderr, "WARNING: steps=");
1765         fprintf(stderr, "%" GMX_PRId64, bench_nsteps);
1766         fprintf(stderr, ". Are you sure you want to perform so %s steps for each benchmark?\n", (bench_nsteps < 100) ? "few" : "many");
1767     }
1768
1769     if (presteps < 0)
1770     {
1771         gmx_fatal(FARGS, "Cannot have a negative number of presteps.\n");
1772     }
1773
1774     /* Check for rcoulomb scaling if more than one .tpr file is tested */
1775     if (*ntprs > 1)
1776     {
1777         if (*rmin/rcoulomb < 0.75 || *rmax/rcoulomb > 1.25)
1778         {
1779             fprintf(stderr, "WARNING: Applying extreme scaling factor. I hope you know what you are doing.\n");
1780         }
1781     }
1782
1783     /* If a fixed number of PME nodes is set we do rcoulomb and PME gird tuning
1784      * only. We need to check whether the requested number of PME-only nodes
1785      * makes sense. */
1786     if (npme_fixed > -1)
1787     {
1788         /* No more than 50% of all nodes can be assigned as PME-only nodes. */
1789         if (2*npme_fixed > nnodes)
1790         {
1791             gmx_fatal(FARGS, "Cannot have more than %d PME-only ranks for a total of %d ranks (you chose %d).\n",
1792                       nnodes/2, nnodes, npme_fixed);
1793         }
1794         if ((npme_fixed > 0) && (5*npme_fixed < nnodes))
1795         {
1796             fprintf(stderr, "WARNING: Only %g percent of the ranks are assigned as PME-only ranks.\n",
1797                     (100.0*npme_fixed)/nnodes);
1798         }
1799         if (opt2parg_bSet("-min", npargs, pa) || opt2parg_bSet("-max", npargs, pa))
1800         {
1801             fprintf(stderr, "NOTE: The -min, -max, and -npme options have no effect when a\n"
1802                     "      fixed number of PME-only ranks is requested with -fix.\n");
1803         }
1804     }
1805 }
1806
1807
1808 /* Returns TRUE when "opt" is needed at launch time */
1809 static gmx_bool is_launch_file(char *opt, gmx_bool bSet)
1810 {
1811     if (0 == std::strncmp(opt, "-swap", 5))
1812     {
1813         return bSet;
1814     }
1815
1816     /* Apart from the input .tpr and the output log files we need all options that
1817      * were set on the command line and that do not start with -b */
1818     if    (0 == std::strncmp(opt, "-b", 2) || 0 == std::strncmp(opt, "-s", 2)
1819            || 0 == std::strncmp(opt, "-err", 4) || 0 == std::strncmp(opt, "-p", 2) )
1820     {
1821         return FALSE;
1822     }
1823
1824     return bSet;
1825 }
1826
1827
1828 /* Returns TRUE when "opt" defines a file which is needed for the benchmarks runs */
1829 static gmx_bool is_bench_file(char *opt, gmx_bool bSet, gmx_bool bOptional, gmx_bool bIsOutput)
1830 {
1831     /* Apart from the input .tpr, all files starting with "-b" are for
1832      * _b_enchmark files exclusively */
1833     if (0 == std::strncmp(opt, "-s", 2))
1834     {
1835         return FALSE;
1836     }
1837
1838     if (0 == std::strncmp(opt, "-b", 2) || 0 == std::strncmp(opt, "-s", 2))
1839     {
1840         if (!bOptional || bSet)
1841         {
1842             return TRUE;
1843         }
1844         else
1845         {
1846             return FALSE;
1847         }
1848     }
1849     else
1850     {
1851         if (bIsOutput)
1852         {
1853             return FALSE;
1854         }
1855         else
1856         {
1857             if (bSet) /* These are additional input files like -cpi -ei */
1858             {
1859                 return TRUE;
1860             }
1861             else
1862             {
1863                 return FALSE;
1864             }
1865         }
1866     }
1867 }
1868
1869
1870 /* Adds 'buf' to 'str' */
1871 static void add_to_string(char **str, const char *buf)
1872 {
1873     int len;
1874
1875
1876     len = std::strlen(*str) + std::strlen(buf) + 1;
1877     srenew(*str, len);
1878     std::strcat(*str, buf);
1879 }
1880
1881
1882 /* Create the command line for the benchmark as well as for the real run */
1883 static void create_command_line_snippets(
1884         gmx_bool  bAppendFiles,
1885         gmx_bool  bKeepAndNumCPT,
1886         gmx_bool  bResetHWay,
1887         int       presteps,
1888         int       nfile,
1889         t_filenm  fnm[],
1890         char     *cmd_args_bench[],  /* command line arguments for benchmark runs */
1891         char     *cmd_args_launch[], /* command line arguments for simulation run */
1892         char      extra_args[],      /* Add this to the end of the command line */
1893         char     *deffnm)            /* Default file names, or NULL if not set */
1894 {
1895     int         i;
1896     char       *opt;
1897     const char *name;
1898     char        strbuf[STRLEN];
1899
1900
1901     /* strlen needs at least '\0' as a string: */
1902     snew(*cmd_args_bench, 1);
1903     snew(*cmd_args_launch, 1);
1904     *cmd_args_launch[0] = '\0';
1905     *cmd_args_bench[0]  = '\0';
1906
1907
1908     /*******************************************/
1909     /* 1. Process other command line arguments */
1910     /*******************************************/
1911     if (presteps > 0)
1912     {
1913         /* Add equilibration steps to benchmark options */
1914         sprintf(strbuf, "-resetstep %d ", presteps);
1915         add_to_string(cmd_args_bench, strbuf);
1916     }
1917     /* These switches take effect only at launch time */
1918     if (deffnm)
1919     {
1920         sprintf(strbuf, "-deffnm %s ", deffnm);
1921         add_to_string(cmd_args_launch, strbuf);
1922     }
1923     if (FALSE == bAppendFiles)
1924     {
1925         add_to_string(cmd_args_launch, "-noappend ");
1926     }
1927     if (bKeepAndNumCPT)
1928     {
1929         add_to_string(cmd_args_launch, "-cpnum ");
1930     }
1931     if (bResetHWay)
1932     {
1933         add_to_string(cmd_args_launch, "-resethway ");
1934     }
1935
1936     /********************/
1937     /* 2. Process files */
1938     /********************/
1939     for (i = 0; i < nfile; i++)
1940     {
1941         opt  = (char *)fnm[i].opt;
1942         name = opt2fn(opt, nfile, fnm);
1943
1944         /* Strbuf contains the options, now let's sort out where we need that */
1945         sprintf(strbuf, "%s %s ", opt, name);
1946
1947         if (is_bench_file(opt, opt2bSet(opt, nfile, fnm), is_optional(&fnm[i]), is_output(&fnm[i])) )
1948         {
1949             /* All options starting with -b* need the 'b' removed,
1950              * therefore overwrite strbuf */
1951             if (0 == std::strncmp(opt, "-b", 2))
1952             {
1953                 sprintf(strbuf, "-%s %s ", &opt[2], name);
1954             }
1955
1956             add_to_string(cmd_args_bench, strbuf);
1957         }
1958
1959         if (is_launch_file(opt, opt2bSet(opt, nfile, fnm)) )
1960         {
1961             add_to_string(cmd_args_launch, strbuf);
1962         }
1963     }
1964
1965     add_to_string(cmd_args_bench, extra_args);
1966     add_to_string(cmd_args_launch, extra_args);
1967 }
1968
1969
1970 /* Set option opt */
1971 static void setopt(const char *opt, int nfile, t_filenm fnm[])
1972 {
1973     int i;
1974
1975     for (i = 0; (i < nfile); i++)
1976     {
1977         if (std::strcmp(opt, fnm[i].opt) == 0)
1978         {
1979             fnm[i].flag |= ffSET;
1980         }
1981     }
1982 }
1983
1984
1985 /* This routine inspects the tpr file and ...
1986  * 1. checks for output files that get triggered by a tpr option. These output
1987  *    files are marked as 'set' to allow for a proper cleanup after each
1988  *    tuning run.
1989  * 2. returns the PME:PP load ratio
1990  * 3. returns rcoulomb from the tpr */
1991 static float inspect_tpr(int nfile, t_filenm fnm[], real *rcoulomb)
1992 {
1993     gmx_bool     bTpi;      /* Is test particle insertion requested?          */
1994     gmx_bool     bFree;     /* Is a free energy simulation requested?         */
1995     gmx_bool     bNM;       /* Is a normal mode analysis requested?           */
1996     gmx_bool     bSwap;     /* Is water/ion position swapping requested?      */
1997     t_state      state;
1998     gmx_mtop_t   mtop;
1999
2000
2001     /* Check tpr file for options that trigger extra output files */
2002     t_inputrec  irInstance;
2003     t_inputrec *ir = &irInstance;
2004     read_tpx_state(opt2fn("-s", nfile, fnm), ir, &state, &mtop);
2005     bFree = (efepNO  != ir->efep );
2006     bNM   = (eiNM    == ir->eI   );
2007     bSwap = (eswapNO != ir->eSwapCoords);
2008     bTpi  = EI_TPI(ir->eI);
2009
2010     /* Set these output files on the tuning command-line */
2011     if (ir->bPull)
2012     {
2013         setopt("-pf", nfile, fnm);
2014         setopt("-px", nfile, fnm);
2015     }
2016     if (bFree)
2017     {
2018         setopt("-dhdl", nfile, fnm);
2019     }
2020     if (bTpi)
2021     {
2022         setopt("-tpi", nfile, fnm);
2023         setopt("-tpid", nfile, fnm);
2024     }
2025     if (bNM)
2026     {
2027         setopt("-mtx", nfile, fnm);
2028     }
2029     if (bSwap)
2030     {
2031         setopt("-swap", nfile, fnm);
2032     }
2033
2034     *rcoulomb = ir->rcoulomb;
2035
2036     /* Return the estimate for the number of PME nodes */
2037     float npme = pme_load_estimate(&mtop, ir, state.box);
2038     return npme;
2039 }
2040
2041
2042 static void couple_files_options(int nfile, t_filenm fnm[])
2043 {
2044     int      i;
2045     gmx_bool bSet, bBench;
2046     char    *opt;
2047     char     buf[20];
2048
2049
2050     for (i = 0; i < nfile; i++)
2051     {
2052         opt    = (char *)fnm[i].opt;
2053         bSet   = ((fnm[i].flag & ffSET) != 0);
2054         bBench = (0 == std::strncmp(opt, "-b", 2));
2055
2056         /* Check optional files */
2057         /* If e.g. -eo is set, then -beo also needs to be set */
2058         if (is_optional(&fnm[i]) && bSet && !bBench)
2059         {
2060             sprintf(buf, "-b%s", &opt[1]);
2061             setopt(buf, nfile, fnm);
2062         }
2063         /* If -beo is set, then -eo also needs to be! */
2064         if (is_optional(&fnm[i]) && bSet && bBench)
2065         {
2066             sprintf(buf, "-%s", &opt[2]);
2067             setopt(buf, nfile, fnm);
2068         }
2069     }
2070 }
2071
2072
2073 #define BENCHSTEPS (1000)
2074
2075 int gmx_tune_pme(int argc, char *argv[])
2076 {
2077     const char     *desc[] = {
2078         "For a given number [TT]-np[tt] or [TT]-ntmpi[tt] of ranks, [THISMODULE] systematically",
2079         "times [gmx-mdrun] with various numbers of PME-only ranks and determines",
2080         "which setting is fastest. It will also test whether performance can",
2081         "be enhanced by shifting load from the reciprocal to the real space",
2082         "part of the Ewald sum. ",
2083         "Simply pass your [REF].tpr[ref] file to [THISMODULE] together with other options",
2084         "for [gmx-mdrun] as needed.[PAR]",
2085         "[THISMODULE] needs to call [gmx-mdrun] and so requires that you",
2086         "specify how to call mdrun with the argument to the [TT]-mdrun[tt]",
2087         "parameter. Depending how you have built GROMACS, values such as",
2088         "'gmx mdrun', 'gmx_d mdrun', or 'mdrun_mpi' might be needed.[PAR]",
2089         "The program that runs MPI programs can be set in the environment variable",
2090         "MPIRUN (defaults to 'mpirun'). Note that for certain MPI frameworks,",
2091         "you need to provide a machine- or hostfile. This can also be passed",
2092         "via the MPIRUN variable, e.g.[PAR]",
2093         "[TT]export MPIRUN=\"/usr/local/mpirun -machinefile hosts\"[tt]",
2094         "Note that in such cases it is normally necessary to compile",
2095         "and/or run [THISMODULE] without MPI support, so that it can call",
2096         "the MPIRUN program.[PAR]",
2097         "Before doing the actual benchmark runs, [THISMODULE] will do a quick",
2098         "check whether [gmx-mdrun] works as expected with the provided parallel settings",
2099         "if the [TT]-check[tt] option is activated (the default).",
2100         "Please call [THISMODULE] with the normal options you would pass to",
2101         "[gmx-mdrun] and add [TT]-np[tt] for the number of ranks to perform the",
2102         "tests on, or [TT]-ntmpi[tt] for the number of threads. You can also add [TT]-r[tt]",
2103         "to repeat each test several times to get better statistics. [PAR]",
2104         "[THISMODULE] can test various real space / reciprocal space workloads",
2105         "for you. With [TT]-ntpr[tt] you control how many extra [REF].tpr[ref] files will be",
2106         "written with enlarged cutoffs and smaller Fourier grids respectively.",
2107         "Typically, the first test (number 0) will be with the settings from the input",
2108         "[REF].tpr[ref] file; the last test (number [TT]ntpr[tt]) will have the Coulomb cutoff",
2109         "specified by [TT]-rmax[tt] with a somewhat smaller PME grid at the same time. ",
2110         "In this last test, the Fourier spacing is multiplied with [TT]rmax[tt]/rcoulomb. ",
2111         "The remaining [REF].tpr[ref] files will have equally-spaced Coulomb radii (and Fourier ",
2112         "spacings) between these extremes. [BB]Note[bb] that you can set [TT]-ntpr[tt] to 1",
2113         "if you just seek the optimal number of PME-only ranks; in that case",
2114         "your input [REF].tpr[ref] file will remain unchanged.[PAR]",
2115         "For the benchmark runs, the default of 1000 time steps should suffice for most",
2116         "MD systems. The dynamic load balancing needs about 100 time steps",
2117         "to adapt to local load imbalances, therefore the time step counters",
2118         "are by default reset after 100 steps. For large systems (>1M atoms), as well as ",
2119         "for a higher accuracy of the measurements, you should set [TT]-resetstep[tt] to a higher value.",
2120         "From the 'DD' load imbalance entries in the md.log output file you",
2121         "can tell after how many steps the load is sufficiently balanced. Example call:[PAR]",
2122         "[TT]gmx tune_pme -np 64 -s protein.tpr -launch[tt][PAR]",
2123         "After calling [gmx-mdrun] several times, detailed performance information",
2124         "is available in the output file [TT]perf.out[tt].",
2125         "[BB]Note[bb] that during the benchmarks, a couple of temporary files are written",
2126         "(options [TT]-b*[tt]), these will be automatically deleted after each test.[PAR]",
2127         "If you want the simulation to be started automatically with the",
2128         "optimized parameters, use the command line option [TT]-launch[tt].[PAR]",
2129         "Basic support for GPU-enabled [TT]mdrun[tt] exists. Give a string containing the IDs",
2130         "of the GPUs that you wish to use in the optimization in the [TT]-gpu_id[tt]",
2131         "command-line argument. This works exactly like [TT]mdrun -gpu_id[tt], does not imply a mapping,",
2132         "and merely declares the eligible set of GPU devices. [TT]gmx-tune_pme[tt] will construct calls to",
2133         "mdrun that use this set appropriately. [TT]gmx-tune_pme[tt] does not support",
2134         "[TT]-gputasks[tt].[PAR]",
2135     };
2136
2137     int             nnodes         = 1;
2138     int             repeats        = 2;
2139     int             pmeentries     = 0; /* How many values for -npme do we actually test for each tpr file */
2140     real            maxPMEfraction = 0.50;
2141     real            minPMEfraction = 0.25;
2142     int             maxPMEnodes, minPMEnodes;
2143     float           guessPMEratio;                    /* guessed PME:PP ratio based on the tpr file */
2144     float           guessPMEnodes;
2145     int             npme_fixed     = -2;              /* If >= -1, use only this number
2146                                                        * of PME-only nodes                */
2147     int             ntprs          = 0;
2148     real            rmin           = 0.0, rmax = 0.0; /* min and max value for rcoulomb if scaling is requested */
2149     real            rcoulomb       = -1.0;            /* Coulomb radius as set in .tpr file */
2150     gmx_bool        bScaleRvdw     = TRUE;
2151     gmx_int64_t     bench_nsteps   = BENCHSTEPS;
2152     gmx_int64_t     new_sim_nsteps = -1;   /* -1 indicates: not set by the user */
2153     gmx_int64_t     cpt_steps      = 0;    /* Step counter in .cpt input file   */
2154     int             presteps       = 1500; /* Do a full cycle reset after presteps steps */
2155     gmx_bool        bOverwrite     = FALSE, bKeepTPR;
2156     gmx_bool        bLaunch        = FALSE;
2157     char           *ExtraArgs      = nullptr;
2158     char          **tpr_names      = nullptr;
2159     const char     *simulation_tpr = nullptr;
2160     char           *deffnm         = nullptr;
2161     int             best_npme, best_tpr;
2162     int             sim_part = 1; /* For benchmarks with checkpoint files */
2163     char            bbuf[STRLEN];
2164
2165
2166     /* Default program names if nothing else is found */
2167     char               *cmd_mpirun = nullptr, *cmd_mdrun = nullptr;
2168     char               *cmd_args_bench, *cmd_args_launch;
2169     char               *cmd_np           = nullptr;
2170
2171     /* IDs of GPUs that are eligible for computation */
2172     char               *eligible_gpu_ids = nullptr;
2173
2174     t_perf            **perfdata = nullptr;
2175     t_inputinfo        *info;
2176     int                 i;
2177     FILE               *fp;
2178
2179     /* Print out how long the tuning took */
2180     double          seconds;
2181
2182     static t_filenm fnm[] = {
2183         /* tune_pme */
2184         { efOUT, "-p",      "perf",     ffWRITE },
2185         { efLOG, "-err",    "bencherr", ffWRITE },
2186         { efTPR, "-so",     "tuned",    ffWRITE },
2187         /* mdrun: */
2188         { efTPR, nullptr,      nullptr,       ffREAD },
2189         { efTRN, "-o",      nullptr,       ffWRITE },
2190         { efCOMPRESSED, "-x", nullptr,     ffOPTWR },
2191         { efCPT, "-cpi",    nullptr,       ffOPTRD },
2192         { efCPT, "-cpo",    nullptr,       ffOPTWR },
2193         { efSTO, "-c",      "confout",  ffWRITE },
2194         { efEDR, "-e",      "ener",     ffWRITE },
2195         { efLOG, "-g",      "md",       ffWRITE },
2196         { efXVG, "-dhdl",   "dhdl",     ffOPTWR },
2197         { efXVG, "-field",  "field",    ffOPTWR },
2198         { efXVG, "-table",  "table",    ffOPTRD },
2199         { efXVG, "-tablep", "tablep",   ffOPTRD },
2200         { efXVG, "-tableb", "table",    ffOPTRD },
2201         { efTRX, "-rerun",  "rerun",    ffOPTRD },
2202         { efXVG, "-tpi",    "tpi",      ffOPTWR },
2203         { efXVG, "-tpid",   "tpidist",  ffOPTWR },
2204         { efEDI, "-ei",     "sam",      ffOPTRD },
2205         { efXVG, "-eo",     "edsam",    ffOPTWR },
2206         { efXVG, "-devout", "deviatie", ffOPTWR },
2207         { efXVG, "-runav",  "runaver",  ffOPTWR },
2208         { efXVG, "-px",     "pullx",    ffOPTWR },
2209         { efXVG, "-pf",     "pullf",    ffOPTWR },
2210         { efXVG, "-ro",     "rotation", ffOPTWR },
2211         { efLOG, "-ra",     "rotangles", ffOPTWR },
2212         { efLOG, "-rs",     "rotslabs", ffOPTWR },
2213         { efLOG, "-rt",     "rottorque", ffOPTWR },
2214         { efMTX, "-mtx",    "nm",       ffOPTWR },
2215         { efXVG, "-swap",   "swapions", ffOPTWR },
2216         /* Output files that are deleted after each benchmark run */
2217         { efTRN, "-bo",     "bench",    ffWRITE },
2218         { efXTC, "-bx",     "bench",    ffWRITE },
2219         { efCPT, "-bcpo",   "bench",    ffWRITE },
2220         { efSTO, "-bc",     "bench",    ffWRITE },
2221         { efEDR, "-be",     "bench",    ffWRITE },
2222         { efLOG, "-bg",     "bench",    ffWRITE },
2223         { efXVG, "-beo",    "benchedo", ffOPTWR },
2224         { efXVG, "-bdhdl",  "benchdhdl", ffOPTWR },
2225         { efXVG, "-bfield", "benchfld", ffOPTWR },
2226         { efXVG, "-btpi",   "benchtpi", ffOPTWR },
2227         { efXVG, "-btpid",  "benchtpid", ffOPTWR },
2228         { efXVG, "-bdevout", "benchdev", ffOPTWR },
2229         { efXVG, "-brunav", "benchrnav", ffOPTWR },
2230         { efXVG, "-bpx",    "benchpx",  ffOPTWR },
2231         { efXVG, "-bpf",    "benchpf",  ffOPTWR },
2232         { efXVG, "-bro",    "benchrot", ffOPTWR },
2233         { efLOG, "-bra",    "benchrota", ffOPTWR },
2234         { efLOG, "-brs",    "benchrots", ffOPTWR },
2235         { efLOG, "-brt",    "benchrott", ffOPTWR },
2236         { efMTX, "-bmtx",   "benchn",   ffOPTWR },
2237         { efNDX, "-bdn",    "bench",    ffOPTWR },
2238         { efXVG, "-bswap",  "benchswp", ffOPTWR }
2239     };
2240
2241     gmx_bool        bThreads     = FALSE;
2242
2243     int             nthreads = 1;
2244
2245     const char     *procstring[] =
2246     { nullptr, "np", "n", "none", nullptr };
2247     const char     *npmevalues_opt[] =
2248     { nullptr, "auto", "all", "subset", nullptr };
2249
2250     gmx_bool          bAppendFiles          = TRUE;
2251     gmx_bool          bKeepAndNumCPT        = FALSE;
2252     gmx_bool          bResetCountersHalfWay = FALSE;
2253     gmx_bool          bBenchmark            = TRUE;
2254     gmx_bool          bCheck                = TRUE;
2255
2256     gmx_output_env_t *oenv = nullptr;
2257
2258     t_pargs           pa[] = {
2259         /***********************/
2260         /* tune_pme options: */
2261         /***********************/
2262         { "-mdrun",    FALSE, etSTR,  {&cmd_mdrun},
2263           "Command line to run a simulation, e.g. 'gmx mdrun' or 'mdrun_mpi'" },
2264         { "-np",       FALSE, etINT,  {&nnodes},
2265           "Number of ranks to run the tests on (must be > 2 for separate PME ranks)" },
2266         { "-npstring", FALSE, etENUM, {procstring},
2267           "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)" },
2268         { "-ntmpi",    FALSE, etINT,  {&nthreads},
2269           "Number of MPI-threads to run the tests on (turns MPI & mpirun off)"},
2270         { "-r",        FALSE, etINT,  {&repeats},
2271           "Repeat each test this often" },
2272         { "-max",      FALSE, etREAL, {&maxPMEfraction},
2273           "Max fraction of PME ranks to test with" },
2274         { "-min",      FALSE, etREAL, {&minPMEfraction},
2275           "Min fraction of PME ranks to test with" },
2276         { "-npme",     FALSE, etENUM, {npmevalues_opt},
2277           "Within -min and -max, benchmark all possible values for [TT]-npme[tt], or just a reasonable subset. "
2278           "Auto neglects -min and -max and chooses reasonable values around a guess for npme derived from the .tpr"},
2279         { "-fix",      FALSE, etINT,  {&npme_fixed},
2280           "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."},
2281         { "-rmax",     FALSE, etREAL, {&rmax},
2282           "If >0, maximal rcoulomb for -ntpr>1 (rcoulomb upscaling results in fourier grid downscaling)" },
2283         { "-rmin",     FALSE, etREAL, {&rmin},
2284           "If >0, minimal rcoulomb for -ntpr>1" },
2285         { "-scalevdw",  FALSE, etBOOL, {&bScaleRvdw},
2286           "Scale rvdw along with rcoulomb"},
2287         { "-ntpr",     FALSE, etINT,  {&ntprs},
2288           "Number of [REF].tpr[ref] files to benchmark. Create this many files with different rcoulomb scaling factors depending on -rmin and -rmax. "
2289           "If < 1, automatically choose the number of [REF].tpr[ref] files to test" },
2290         { "-steps",    FALSE, etINT64, {&bench_nsteps},
2291           "Take timings for this many steps in the benchmark runs" },
2292         { "-resetstep", FALSE, etINT,  {&presteps},
2293           "Let dlb equilibrate this many steps before timings are taken (reset cycle counters after this many steps)" },
2294         { "-nsteps",   FALSE, etINT64, {&new_sim_nsteps},
2295           "If non-negative, perform this many steps in the real run (overwrites nsteps from [REF].tpr[ref], add [REF].cpt[ref] steps)" },
2296         { "-launch",   FALSE, etBOOL, {&bLaunch},
2297           "Launch the real simulation after optimization" },
2298         { "-bench",    FALSE, etBOOL, {&bBenchmark},
2299           "Run the benchmarks or just create the input [REF].tpr[ref] files?" },
2300         { "-check",    FALSE, etBOOL, {&bCheck},
2301           "Before the benchmark runs, check whether mdrun works in parallel" },
2302         { "-gpu_id",   FALSE, etSTR,  {&eligible_gpu_ids},
2303           "List of unique GPU device IDs that are eligible for use" },
2304         /******************/
2305         /* mdrun options: */
2306         /******************/
2307         /* We let tune_pme parse and understand these options, because we need to
2308          * prevent that they appear on the mdrun command line for the benchmarks */
2309         { "-append",   FALSE, etBOOL, {&bAppendFiles},
2310           "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names (for launch only)" },
2311         { "-cpnum",    FALSE, etBOOL, {&bKeepAndNumCPT},
2312           "Keep and number checkpoint files (launch only)" },
2313         { "-deffnm",   FALSE, etSTR,  {&deffnm},
2314           "Set the default filenames (launch only)" },
2315         { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
2316           "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt] (launch only)" }
2317     };
2318
2319 #define NFILE asize(fnm)
2320
2321     seconds = gmx_gettime();
2322
2323     if (!parse_common_args(&argc, argv, PCA_NOEXIT_ON_ARGS,
2324                            NFILE, fnm, asize(pa), pa, asize(desc), desc,
2325                            0, nullptr, &oenv))
2326     {
2327         return 0;
2328     }
2329
2330     // procstring[0] is used inside two different conditionals further down
2331     GMX_RELEASE_ASSERT(procstring[0] != nullptr, "Options inconsistency; procstring[0] is NULL");
2332
2333     /* Store the remaining unparsed command line entries in a string which
2334      * is then attached to the mdrun command line */
2335     snew(ExtraArgs, 1);
2336     ExtraArgs[0] = '\0';
2337     for (i = 1; i < argc; i++) /* argc will now be 1 if everything was understood */
2338     {
2339         add_to_string(&ExtraArgs, argv[i]);
2340         add_to_string(&ExtraArgs, " ");
2341     }
2342
2343     if (opt2parg_bSet("-ntmpi", asize(pa), pa))
2344     {
2345         bThreads = TRUE;
2346         if (opt2parg_bSet("-npstring", asize(pa), pa))
2347         {
2348             fprintf(stderr, "WARNING: -npstring has no effect when using threads.\n");
2349         }
2350
2351         if (nnodes > 1)
2352         {
2353             gmx_fatal(FARGS, "Can't run multi-threaded MPI simulation yet!");
2354         }
2355         /* and now we just set this; a bit of an ugly hack*/
2356         nnodes = nthreads;
2357     }
2358     /* Check for PME:PP ratio and whether tpr triggers additional output files */
2359     guessPMEratio = inspect_tpr(NFILE, fnm, &rcoulomb);
2360
2361     /* Automatically set -beo options if -eo is set etc. */
2362     couple_files_options(NFILE, fnm);
2363
2364     /* Construct the command line arguments for benchmark runs
2365      * as well as for the simulation run */
2366     if (bThreads)
2367     {
2368         sprintf(bbuf, " -ntmpi %d ", nthreads);
2369     }
2370     else
2371     {
2372         /* This string will be used for MPI runs and will appear after the
2373          * mpirun command. */
2374         if (std::strcmp(procstring[0], "none") != 0)
2375         {
2376             sprintf(bbuf, " -%s %d ", procstring[0], nnodes);
2377         }
2378         else
2379         {
2380             sprintf(bbuf, " ");
2381         }
2382     }
2383
2384     cmd_np = bbuf;
2385
2386     create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
2387                                  NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs, deffnm);
2388
2389     /* Prepare to use checkpoint file if requested */
2390     sim_part = 1;
2391     if (opt2bSet("-cpi", NFILE, fnm))
2392     {
2393         const char *filename = opt2fn("-cpi", NFILE, fnm);
2394         int         cpt_sim_part;
2395         read_checkpoint_part_and_step(filename,
2396                                       &cpt_sim_part, &cpt_steps);
2397         if (cpt_sim_part == 0)
2398         {
2399             gmx_fatal(FARGS, "Checkpoint file %s could not be read!", filename);
2400         }
2401         /* tune_pme will run the next part of the simulation */
2402         sim_part = cpt_sim_part + 1;
2403     }
2404
2405     /* Open performance output file and write header info */
2406     fp = gmx_ffopen(opt2fn("-p", NFILE, fnm), "w");
2407
2408     /* Make a quick consistency check of command line parameters */
2409     check_input(nnodes, repeats, &ntprs, &rmin, rcoulomb, &rmax,
2410                 maxPMEfraction, minPMEfraction, npme_fixed,
2411                 bench_nsteps, fnm, NFILE, sim_part, presteps,
2412                 asize(pa), pa);
2413
2414     /* Determine the maximum and minimum number of PME nodes to test,
2415      * the actual list of settings is build in do_the_tests(). */
2416     if ((nnodes > 2) && (npme_fixed < -1))
2417     {
2418         if (0 == std::strcmp(npmevalues_opt[0], "auto"))
2419         {
2420             /* Determine the npme range automatically based on the PME:PP load guess */
2421             if (guessPMEratio > 1.0)
2422             {
2423                 /* More PME than PP work, probably we do not need separate PME nodes at all! */
2424                 maxPMEnodes = nnodes/2;
2425                 minPMEnodes = nnodes/2;
2426             }
2427             else
2428             {
2429                 /* PME : PP load is in the range 0..1, let's test around the guess */
2430                 guessPMEnodes = static_cast<int>(nnodes/(1.0 + 1.0/guessPMEratio));
2431                 minPMEnodes   = static_cast<int>(std::floor(0.7*guessPMEnodes));
2432                 maxPMEnodes   = static_cast<int>(std::ceil(1.6*guessPMEnodes));
2433                 maxPMEnodes   = std::min(maxPMEnodes, nnodes/2);
2434             }
2435         }
2436         else
2437         {
2438             /* Determine the npme range based on user input */
2439             maxPMEnodes = static_cast<int>(std::floor(maxPMEfraction*nnodes));
2440             minPMEnodes = std::max(static_cast<int>(std::floor(minPMEfraction*nnodes)), 0);
2441             fprintf(stdout, "Will try runs with %d ", minPMEnodes);
2442             if (maxPMEnodes != minPMEnodes)
2443             {
2444                 fprintf(stdout, "- %d ", maxPMEnodes);
2445             }
2446             fprintf(stdout, "PME-only ranks.\n  Note that the automatic number of PME-only ranks and no separate PME ranks are always tested.\n");
2447         }
2448     }
2449     else
2450     {
2451         maxPMEnodes = 0;
2452         minPMEnodes = 0;
2453     }
2454
2455     /* Get the commands we need to set up the runs from environment variables */
2456     get_program_paths(bThreads, &cmd_mpirun, &cmd_mdrun);
2457     if (bBenchmark && repeats > 0)
2458     {
2459         check_mdrun_works(bThreads, cmd_mpirun, cmd_np, cmd_mdrun, nullptr != eligible_gpu_ids);
2460     }
2461
2462     /* Print some header info to file */
2463     sep_line(fp);
2464     fprintf(fp, "\n      P E R F O R M A N C E   R E S U L T S\n");
2465     sep_line(fp);
2466     fprintf(fp, "%s for GROMACS %s\n", output_env_get_program_display_name(oenv),
2467             gmx_version());
2468     if (!bThreads)
2469     {
2470         fprintf(fp, "Number of ranks         : %d\n", nnodes);
2471         fprintf(fp, "The mpirun command is   : %s\n", cmd_mpirun);
2472         if (std::strcmp(procstring[0], "none") != 0)
2473         {
2474             fprintf(fp, "Passing # of ranks via  : -%s\n", procstring[0]);
2475         }
2476         else
2477         {
2478             fprintf(fp, "Not setting number of ranks in system call\n");
2479         }
2480     }
2481     else
2482     {
2483         fprintf(fp, "Number of threads       : %d\n", nnodes);
2484     }
2485
2486     fprintf(fp, "The mdrun  command is   : %s\n", cmd_mdrun);
2487     fprintf(fp, "mdrun args benchmarks   : %s\n", cmd_args_bench);
2488     fprintf(fp, "Benchmark steps         : ");
2489     fprintf(fp, "%" GMX_PRId64, bench_nsteps);
2490     fprintf(fp, "\n");
2491     fprintf(fp, "dlb equilibration steps : %d\n", presteps);
2492     if (sim_part > 1)
2493     {
2494         fprintf(fp, "Checkpoint time step    : ");
2495         fprintf(fp, "%" GMX_PRId64, cpt_steps);
2496         fprintf(fp, "\n");
2497     }
2498     fprintf(fp, "mdrun args at launchtime: %s\n", cmd_args_launch);
2499
2500     if (new_sim_nsteps >= 0)
2501     {
2502         bOverwrite = TRUE;
2503         fprintf(stderr, "Note: Simulation input file %s will have ", opt2fn("-so", NFILE, fnm));
2504         fprintf(stderr, "%" GMX_PRId64, new_sim_nsteps+cpt_steps);
2505         fprintf(stderr, " steps.\n");
2506         fprintf(fp, "Simulation steps        : ");
2507         fprintf(fp, "%" GMX_PRId64, new_sim_nsteps);
2508         fprintf(fp, "\n");
2509     }
2510     if (repeats > 1)
2511     {
2512         fprintf(fp, "Repeats for each test   : %d\n", repeats);
2513     }
2514
2515     if (npme_fixed >= -1)
2516     {
2517         fprintf(fp, "Fixing -npme at         : %d\n", npme_fixed);
2518     }
2519
2520     fprintf(fp, "Input file              : %s\n", opt2fn("-s", NFILE, fnm));
2521     fprintf(fp, "   PME/PP load estimate : %g\n", guessPMEratio);
2522
2523     /* Allocate memory for the inputinfo struct: */
2524     snew(info, 1);
2525     info->nr_inputfiles = ntprs;
2526     for (i = 0; i < ntprs; i++)
2527     {
2528         snew(info->rcoulomb, ntprs);
2529         snew(info->rvdw, ntprs);
2530         snew(info->rlist, ntprs);
2531         snew(info->nkx, ntprs);
2532         snew(info->nky, ntprs);
2533         snew(info->nkz, ntprs);
2534         snew(info->fsx, ntprs);
2535         snew(info->fsy, ntprs);
2536         snew(info->fsz, ntprs);
2537     }
2538     /* Make alternative tpr files to test: */
2539     snew(tpr_names, ntprs);
2540     for (i = 0; i < ntprs; i++)
2541     {
2542         snew(tpr_names[i], STRLEN);
2543     }
2544
2545     /* It can be that ntprs is reduced by make_benchmark_tprs if not enough
2546      * different grids could be found. */
2547     make_benchmark_tprs(opt2fn("-s", NFILE, fnm), tpr_names, bench_nsteps+presteps,
2548                         cpt_steps, rmin, rmax, bScaleRvdw, &ntprs, info, fp);
2549
2550     /********************************************************************************/
2551     /* Main loop over all scenarios we need to test: tpr files, PME nodes, repeats  */
2552     /********************************************************************************/
2553     snew(perfdata, ntprs);
2554     if (bBenchmark)
2555     {
2556         GMX_RELEASE_ASSERT(npmevalues_opt[0] != nullptr, "Options inconsistency; npmevalues_opt[0] is NULL");
2557         do_the_tests(fp, tpr_names, maxPMEnodes, minPMEnodes, npme_fixed, npmevalues_opt[0], perfdata, &pmeentries,
2558                      repeats, nnodes, ntprs, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2559                      cmd_args_bench, fnm, NFILE, presteps, cpt_steps, bCheck, eligible_gpu_ids);
2560
2561         fprintf(fp, "\nTuning took%8.1f minutes.\n", (gmx_gettime()-seconds)/60.0);
2562
2563         /* Analyse the results and give a suggestion for optimal settings: */
2564         bKeepTPR = analyze_data(fp, opt2fn("-p", NFILE, fnm), perfdata, nnodes, ntprs, pmeentries,
2565                                 repeats, info, &best_tpr, &best_npme);
2566
2567         /* Take the best-performing tpr file and enlarge nsteps to original value */
2568         if (bKeepTPR && !bOverwrite)
2569         {
2570             simulation_tpr = opt2fn("-s", NFILE, fnm);
2571         }
2572         else
2573         {
2574             simulation_tpr = opt2fn("-so", NFILE, fnm);
2575             modify_PMEsettings(bOverwrite ? (new_sim_nsteps+cpt_steps) : info->orig_sim_steps,
2576                                info->orig_init_step, tpr_names[best_tpr], simulation_tpr);
2577         }
2578
2579         /* Let's get rid of the temporary benchmark input files */
2580         for (i = 0; i < ntprs; i++)
2581         {
2582             fprintf(stdout, "Deleting temporary benchmark input file %s\n", tpr_names[i]);
2583             remove(tpr_names[i]);
2584         }
2585
2586         /* Now start the real simulation if the user requested it ... */
2587         launch_simulation(bLaunch, fp, bThreads, cmd_mpirun, cmd_np, cmd_mdrun,
2588                           cmd_args_launch, simulation_tpr, best_npme, eligible_gpu_ids);
2589     }
2590     gmx_ffclose(fp);
2591
2592     /* ... or simply print the performance results to screen: */
2593     if (!bLaunch)
2594     {
2595         finalize(opt2fn("-p", NFILE, fnm));
2596     }
2597
2598     return 0;
2599 }