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