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