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