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