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