Implemented LJ-PME nbnxn kernels
[alexxy/gromacs.git] / src / programs / mdrun / runner.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40 #include <signal.h>
41 #include <stdlib.h>
42 #ifdef HAVE_UNISTD_H
43 #include <unistd.h>
44 #endif
45 #include <string.h>
46 #include <assert.h>
47
48 #include "typedefs.h"
49 #include "smalloc.h"
50 #include "sysstuff.h"
51 #include "copyrite.h"
52 #include "force.h"
53 #include "mdrun.h"
54 #include "md_logging.h"
55 #include "md_support.h"
56 #include "network.h"
57 #include "names.h"
58 #include "disre.h"
59 #include "orires.h"
60 #include "pme.h"
61 #include "mdatoms.h"
62 #include "repl_ex.h"
63 #include "qmmm.h"
64 #include "domdec.h"
65 #include "coulomb.h"
66 #include "constr.h"
67 #include "mvdata.h"
68 #include "checkpoint.h"
69 #include "mtop_util.h"
70 #include "sighandler.h"
71 #include "txtdump.h"
72 #include "gmx_detect_hardware.h"
73 #include "gmx_omp_nthreads.h"
74 #include "calc_verletbuf.h"
75 #include "gmx_fatal_collective.h"
76 #include "membed.h"
77 #include "macros.h"
78 #include "gmx_thread_affinity.h"
79 #include "inputrec.h"
80
81 #include "gromacs/fileio/tpxio.h"
82 #include "gromacs/mdlib/nbnxn_search.h"
83 #include "gromacs/mdlib/nbnxn_consts.h"
84 #include "gromacs/timing/wallcycle.h"
85 #include "gromacs/utility/gmxmpi.h"
86 #include "gromacs/utility/gmxomp.h"
87 #include "gromacs/swap/swapcoords.h"
88 #include "gromacs/essentialdynamics/edsam.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/pulling/pull_rotation.h"
91
92 #ifdef GMX_FAHCORE
93 #include "corewrap.h"
94 #endif
95
96 #include "gpu_utils.h"
97 #include "nbnxn_cuda_data_mgmt.h"
98
99 typedef struct {
100     gmx_integrator_t *func;
101 } gmx_intp_t;
102
103 /* The array should match the eI array in include/types/enums.h */
104 const gmx_intp_t    integrator[eiNR] = { {do_md}, {do_steep}, {do_cg}, {do_md}, {do_md}, {do_nm}, {do_lbfgs}, {do_tpi}, {do_tpi}, {do_md}, {do_md}, {do_md}};
105
106 gmx_int64_t         deform_init_init_step_tpx;
107 matrix              deform_init_box_tpx;
108 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
109
110
111 #ifdef GMX_THREAD_MPI
112 struct mdrunner_arglist
113 {
114     gmx_hw_opt_t    hw_opt;
115     FILE           *fplog;
116     t_commrec      *cr;
117     int             nfile;
118     const t_filenm *fnm;
119     output_env_t    oenv;
120     gmx_bool        bVerbose;
121     gmx_bool        bCompact;
122     int             nstglobalcomm;
123     ivec            ddxyz;
124     int             dd_node_order;
125     real            rdd;
126     real            rconstr;
127     const char     *dddlb_opt;
128     real            dlb_scale;
129     const char     *ddcsx;
130     const char     *ddcsy;
131     const char     *ddcsz;
132     const char     *nbpu_opt;
133     int             nstlist_cmdline;
134     gmx_int64_t     nsteps_cmdline;
135     int             nstepout;
136     int             resetstep;
137     int             nmultisim;
138     int             repl_ex_nst;
139     int             repl_ex_nex;
140     int             repl_ex_seed;
141     real            pforce;
142     real            cpt_period;
143     real            max_hours;
144     const char     *deviceOptions;
145     unsigned long   Flags;
146     int             ret; /* return value */
147 };
148
149
150 /* The function used for spawning threads. Extracts the mdrunner()
151    arguments from its one argument and calls mdrunner(), after making
152    a commrec. */
153 static void mdrunner_start_fn(void *arg)
154 {
155     struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
156     struct mdrunner_arglist  mc  = *mda; /* copy the arg list to make sure
157                                             that it's thread-local. This doesn't
158                                             copy pointed-to items, of course,
159                                             but those are all const. */
160     t_commrec *cr;                       /* we need a local version of this */
161     FILE      *fplog = NULL;
162     t_filenm  *fnm;
163
164     fnm = dup_tfn(mc.nfile, mc.fnm);
165
166     cr = reinitialize_commrec_for_this_thread(mc.cr);
167
168     if (MASTER(cr))
169     {
170         fplog = mc.fplog;
171     }
172
173     mda->ret = mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
174                         mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
175                         mc.ddxyz, mc.dd_node_order, mc.rdd,
176                         mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
177                         mc.ddcsx, mc.ddcsy, mc.ddcsz,
178                         mc.nbpu_opt, mc.nstlist_cmdline,
179                         mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
180                         mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
181                         mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
182 }
183
184 /* called by mdrunner() to start a specific number of threads (including
185    the main thread) for thread-parallel runs. This in turn calls mdrunner()
186    for each thread.
187    All options besides nthreads are the same as for mdrunner(). */
188 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
189                                          FILE *fplog, t_commrec *cr, int nfile,
190                                          const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
191                                          gmx_bool bCompact, int nstglobalcomm,
192                                          ivec ddxyz, int dd_node_order, real rdd, real rconstr,
193                                          const char *dddlb_opt, real dlb_scale,
194                                          const char *ddcsx, const char *ddcsy, const char *ddcsz,
195                                          const char *nbpu_opt, int nstlist_cmdline,
196                                          gmx_int64_t nsteps_cmdline,
197                                          int nstepout, int resetstep,
198                                          int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
199                                          real pforce, real cpt_period, real max_hours,
200                                          const char *deviceOptions, unsigned long Flags)
201 {
202     int                      ret;
203     struct mdrunner_arglist *mda;
204     t_commrec               *crn; /* the new commrec */
205     t_filenm                *fnmn;
206
207     /* first check whether we even need to start tMPI */
208     if (hw_opt->nthreads_tmpi < 2)
209     {
210         return cr;
211     }
212
213     /* a few small, one-time, almost unavoidable memory leaks: */
214     snew(mda, 1);
215     fnmn = dup_tfn(nfile, fnm);
216
217     /* fill the data structure to pass as void pointer to thread start fn */
218     /* hw_opt contains pointers, which should all be NULL at this stage */
219     mda->hw_opt          = *hw_opt;
220     mda->fplog           = fplog;
221     mda->cr              = cr;
222     mda->nfile           = nfile;
223     mda->fnm             = fnmn;
224     mda->oenv            = oenv;
225     mda->bVerbose        = bVerbose;
226     mda->bCompact        = bCompact;
227     mda->nstglobalcomm   = nstglobalcomm;
228     mda->ddxyz[XX]       = ddxyz[XX];
229     mda->ddxyz[YY]       = ddxyz[YY];
230     mda->ddxyz[ZZ]       = ddxyz[ZZ];
231     mda->dd_node_order   = dd_node_order;
232     mda->rdd             = rdd;
233     mda->rconstr         = rconstr;
234     mda->dddlb_opt       = dddlb_opt;
235     mda->dlb_scale       = dlb_scale;
236     mda->ddcsx           = ddcsx;
237     mda->ddcsy           = ddcsy;
238     mda->ddcsz           = ddcsz;
239     mda->nbpu_opt        = nbpu_opt;
240     mda->nstlist_cmdline = nstlist_cmdline;
241     mda->nsteps_cmdline  = nsteps_cmdline;
242     mda->nstepout        = nstepout;
243     mda->resetstep       = resetstep;
244     mda->nmultisim       = nmultisim;
245     mda->repl_ex_nst     = repl_ex_nst;
246     mda->repl_ex_nex     = repl_ex_nex;
247     mda->repl_ex_seed    = repl_ex_seed;
248     mda->pforce          = pforce;
249     mda->cpt_period      = cpt_period;
250     mda->max_hours       = max_hours;
251     mda->deviceOptions   = deviceOptions;
252     mda->Flags           = Flags;
253
254     /* now spawn new threads that start mdrunner_start_fn(), while
255        the main thread returns, we set thread affinity later */
256     ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
257                        mdrunner_start_fn, (void*)(mda) );
258     if (ret != TMPI_SUCCESS)
259     {
260         return NULL;
261     }
262
263     crn = reinitialize_commrec_for_this_thread(cr);
264     return crn;
265 }
266
267
268 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
269                                         const gmx_hw_opt_t  *hw_opt,
270                                         int                  nthreads_tot,
271                                         int                  ngpu)
272 {
273     int nthreads_tmpi;
274
275     /* There are no separate PME nodes here, as we ensured in
276      * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
277      * and a conditional ensures we would not have ended up here.
278      * Note that separate PME nodes might be switched on later.
279      */
280     if (ngpu > 0)
281     {
282         nthreads_tmpi = ngpu;
283         if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
284         {
285             nthreads_tmpi = nthreads_tot;
286         }
287     }
288     else if (hw_opt->nthreads_omp > 0)
289     {
290         /* Here we could oversubscribe, when we do, we issue a warning later */
291         nthreads_tmpi = max(1, nthreads_tot/hw_opt->nthreads_omp);
292     }
293     else
294     {
295         /* TODO choose nthreads_omp based on hardware topology
296            when we have a hardware topology detection library */
297         /* In general, when running up to 4 threads, OpenMP should be faster.
298          * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
299          * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
300          * even on two CPUs it's usually faster (but with many OpenMP threads
301          * it could be faster not to use HT, currently we always use HT).
302          * On Nehalem/Westmere we want to avoid running 16 threads over
303          * two CPUs with HT, so we need a limit<16; thus we use 12.
304          * A reasonable limit for Intel Sandy and Ivy bridge,
305          * not knowing the topology, is 16 threads.
306          */
307         const int nthreads_omp_always_faster             =  4;
308         const int nthreads_omp_always_faster_Nehalem     = 12;
309         const int nthreads_omp_always_faster_SandyBridge = 16;
310         const int first_model_Nehalem                    = 0x1A;
311         const int first_model_SandyBridge                = 0x2A;
312         gmx_bool  bIntel_Family6;
313
314         bIntel_Family6 =
315             (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
316              gmx_cpuid_family(hwinfo->cpuid_info) == 6);
317
318         if (nthreads_tot <= nthreads_omp_always_faster ||
319             (bIntel_Family6 &&
320              ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
321               (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
322         {
323             /* Use pure OpenMP parallelization */
324             nthreads_tmpi = 1;
325         }
326         else
327         {
328             /* Don't use OpenMP parallelization */
329             nthreads_tmpi = nthreads_tot;
330         }
331     }
332
333     return nthreads_tmpi;
334 }
335
336
337 /* Get the number of threads to use for thread-MPI based on how many
338  * were requested, which algorithms we're using,
339  * and how many particles there are.
340  * At the point we have already called check_and_update_hw_opt.
341  * Thus all options should be internally consistent and consistent
342  * with the hardware, except that ntmpi could be larger than #GPU.
343  */
344 static int get_nthreads_mpi(const gmx_hw_info_t *hwinfo,
345                             gmx_hw_opt_t *hw_opt,
346                             t_inputrec *inputrec, gmx_mtop_t *mtop,
347                             const t_commrec *cr,
348                             FILE *fplog)
349 {
350     int      nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
351     int      min_atoms_per_mpi_thread;
352     char    *env;
353     char     sbuf[STRLEN];
354     gmx_bool bCanUseGPU;
355
356     if (hw_opt->nthreads_tmpi > 0)
357     {
358         /* Trivial, return right away */
359         return hw_opt->nthreads_tmpi;
360     }
361
362     nthreads_hw = hwinfo->nthreads_hw_avail;
363
364     /* How many total (#tMPI*#OpenMP) threads can we start? */
365     if (hw_opt->nthreads_tot > 0)
366     {
367         nthreads_tot_max = hw_opt->nthreads_tot;
368     }
369     else
370     {
371         nthreads_tot_max = nthreads_hw;
372     }
373
374     bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET &&
375                   hwinfo->gpu_info.ncuda_dev_compatible > 0);
376     if (bCanUseGPU)
377     {
378         ngpu = hwinfo->gpu_info.ncuda_dev_compatible;
379     }
380     else
381     {
382         ngpu = 0;
383     }
384
385     if (inputrec->cutoff_scheme == ecutsGROUP)
386     {
387         /* We checked this before, but it doesn't hurt to do it once more */
388         assert(hw_opt->nthreads_omp == 1);
389     }
390
391     nthreads_tmpi =
392         get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
393
394     if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
395     {
396         /* Dims/steps are divided over the nodes iso splitting the atoms */
397         min_atoms_per_mpi_thread = 0;
398     }
399     else
400     {
401         if (bCanUseGPU)
402         {
403             min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
404         }
405         else
406         {
407             min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
408         }
409     }
410
411     /* Check if an algorithm does not support parallel simulation.  */
412     if (nthreads_tmpi != 1 &&
413         ( inputrec->eI == eiLBFGS ||
414           inputrec->coulombtype == eelEWALD ) )
415     {
416         nthreads_tmpi = 1;
417
418         md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
419         if (hw_opt->nthreads_tmpi > nthreads_tmpi)
420         {
421             gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
422         }
423     }
424     else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
425     {
426         /* the thread number was chosen automatically, but there are too many
427            threads (too few atoms per thread) */
428         nthreads_new = max(1, mtop->natoms/min_atoms_per_mpi_thread);
429
430         /* Avoid partial use of Hyper-Threading */
431         if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
432             nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
433         {
434             nthreads_new = nthreads_hw/2;
435         }
436
437         /* Avoid large prime numbers in the thread count */
438         if (nthreads_new >= 6)
439         {
440             /* Use only 6,8,10 with additional factors of 2 */
441             int fac;
442
443             fac = 2;
444             while (3*fac*2 <= nthreads_new)
445             {
446                 fac *= 2;
447             }
448
449             nthreads_new = (nthreads_new/fac)*fac;
450         }
451         else
452         {
453             /* Avoid 5 */
454             if (nthreads_new == 5)
455             {
456                 nthreads_new = 4;
457             }
458         }
459
460         nthreads_tmpi = nthreads_new;
461
462         fprintf(stderr, "\n");
463         fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
464         fprintf(stderr, "      only starting %d thread-MPI threads.\n", nthreads_tmpi);
465         fprintf(stderr, "      You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
466     }
467
468     return nthreads_tmpi;
469 }
470 #endif /* GMX_THREAD_MPI */
471
472
473 /* We determine the extra cost of the non-bonded kernels compared to
474  * a reference nstlist value of 10 (which is the default in grompp).
475  */
476 static const int    nbnxn_reference_nstlist = 10;
477 /* The values to try when switching  */
478 const int           nstlist_try[] = { 20, 25, 40 };
479 #define NNSTL  sizeof(nstlist_try)/sizeof(nstlist_try[0])
480 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
481  * but never more than listfac_max.
482  * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
483  * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
484  * Note that both CPU and GPU factors are conservative. Performance should
485  * not go down due to this tuning, except with a relatively slow GPU.
486  * On the other hand, at medium/high parallelization or with fast GPUs
487  * nstlist will not be increased enough to reach optimal performance.
488  */
489 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
490 static const float  nbnxn_cpu_listfac_ok    = 1.05;
491 static const float  nbnxn_cpu_listfac_max   = 1.09;
492 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
493 static const float  nbnxn_gpu_listfac_ok    = 1.20;
494 static const float  nbnxn_gpu_listfac_max   = 1.30;
495
496 /* Try to increase nstlist when using the Verlet cut-off scheme */
497 static void increase_nstlist(FILE *fp, t_commrec *cr,
498                              t_inputrec *ir, int nstlist_cmdline,
499                              const gmx_mtop_t *mtop, matrix box,
500                              gmx_bool bGPU)
501 {
502     float                  listfac_ok, listfac_max;
503     int                    nstlist_orig, nstlist_prev;
504     verletbuf_list_setup_t ls;
505     real                   rlist_nstlist10, rlist_inc, rlist_ok, rlist_max;
506     real                   rlist_new, rlist_prev;
507     int                    nstlist_ind = 0;
508     t_state                state_tmp;
509     gmx_bool               bBox, bDD, bCont;
510     const char            *nstl_gpu = "\nFor optimal performance with a GPU nstlist (now %d) should be larger.\nThe optimum depends on your CPU and GPU resources.\nYou might want to try several nstlist values.\n";
511     const char            *nve_err  = "Can not increase nstlist because an NVE ensemble is used";
512     const char            *vbd_err  = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
513     const char            *box_err  = "Can not increase nstlist because the box is too small";
514     const char            *dd_err   = "Can not increase nstlist because of domain decomposition limitations";
515     char                   buf[STRLEN];
516
517     if (nstlist_cmdline <= 0)
518     {
519         if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
520         {
521             fprintf(fp, nstl_gpu, ir->nstlist);
522         }
523         nstlist_ind = 0;
524         while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
525         {
526             nstlist_ind++;
527         }
528         if (nstlist_ind == NNSTL)
529         {
530             /* There are no larger nstlist value to try */
531             return;
532         }
533     }
534
535     if (EI_MD(ir->eI) && ir->etc == etcNO)
536     {
537         if (MASTER(cr))
538         {
539             fprintf(stderr, "%s\n", nve_err);
540         }
541         if (fp != NULL)
542         {
543             fprintf(fp, "%s\n", nve_err);
544         }
545
546         return;
547     }
548
549     if (ir->verletbuf_tol == 0 && bGPU)
550     {
551         gmx_fatal(FARGS, "You are using an old tpr file with a GPU, please generate a new tpr file with an up to date version of grompp");
552     }
553
554     if (ir->verletbuf_tol < 0)
555     {
556         if (MASTER(cr))
557         {
558             fprintf(stderr, "%s\n", vbd_err);
559         }
560         if (fp != NULL)
561         {
562             fprintf(fp, "%s\n", vbd_err);
563         }
564
565         return;
566     }
567
568     if (bGPU)
569     {
570         listfac_ok  = nbnxn_gpu_listfac_ok;
571         listfac_max = nbnxn_gpu_listfac_max;
572     }
573     else
574     {
575         listfac_ok  = nbnxn_cpu_listfac_ok;
576         listfac_max = nbnxn_cpu_listfac_max;
577     }
578
579     nstlist_orig = ir->nstlist;
580     if (nstlist_cmdline > 0)
581     {
582         if (fp)
583         {
584             sprintf(buf, "Getting nstlist=%d from command line option",
585                     nstlist_cmdline);
586         }
587         ir->nstlist = nstlist_cmdline;
588     }
589
590     verletbuf_get_list_setup(bGPU, &ls);
591
592     /* Allow rlist to make the list a given factor larger than the list
593      * would be with nstlist=10.
594      */
595     nstlist_prev = ir->nstlist;
596     ir->nstlist  = 10;
597     calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
598                             &rlist_nstlist10);
599     ir->nstlist  = nstlist_prev;
600
601     /* Determine the pair list size increase due to zero interactions */
602     rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
603                                               mtop->natoms/det(box));
604     rlist_ok  = (rlist_nstlist10 + rlist_inc)*pow(listfac_ok, 1.0/3.0) - rlist_inc;
605     rlist_max = (rlist_nstlist10 + rlist_inc)*pow(listfac_max, 1.0/3.0) - rlist_inc;
606     if (debug)
607     {
608         fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
609                 rlist_inc, rlist_ok, rlist_max);
610     }
611
612     nstlist_prev = nstlist_orig;
613     rlist_prev   = ir->rlist;
614     do
615     {
616         if (nstlist_cmdline <= 0)
617         {
618             ir->nstlist = nstlist_try[nstlist_ind];
619         }
620
621         /* Set the pair-list buffer size in ir */
622         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
623
624         /* Does rlist fit in the box? */
625         bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
626         bDD  = TRUE;
627         if (bBox && DOMAINDECOMP(cr))
628         {
629             /* Check if rlist fits in the domain decomposition */
630             if (inputrec2nboundeddim(ir) < DIM)
631             {
632                 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
633             }
634             copy_mat(box, state_tmp.box);
635             bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
636         }
637
638         if (debug)
639         {
640             fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
641                     ir->nstlist, rlist_new, bBox, bDD);
642         }
643
644         bCont = FALSE;
645
646         if (nstlist_cmdline <= 0)
647         {
648             if (bBox && bDD && rlist_new <= rlist_max)
649             {
650                 /* Increase nstlist */
651                 nstlist_prev = ir->nstlist;
652                 rlist_prev   = rlist_new;
653                 bCont        = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
654             }
655             else
656             {
657                 /* Stick with the previous nstlist */
658                 ir->nstlist = nstlist_prev;
659                 rlist_new   = rlist_prev;
660                 bBox        = TRUE;
661                 bDD         = TRUE;
662             }
663         }
664
665         nstlist_ind++;
666     }
667     while (bCont);
668
669     if (!bBox || !bDD)
670     {
671         gmx_warning(!bBox ? box_err : dd_err);
672         if (fp != NULL)
673         {
674             fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
675         }
676         ir->nstlist = nstlist_orig;
677     }
678     else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
679     {
680         sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
681                 nstlist_orig, ir->nstlist,
682                 ir->rlist, rlist_new);
683         if (MASTER(cr))
684         {
685             fprintf(stderr, "%s\n\n", buf);
686         }
687         if (fp != NULL)
688         {
689             fprintf(fp, "%s\n\n", buf);
690         }
691         ir->rlist     = rlist_new;
692         ir->rlistlong = rlist_new;
693     }
694 }
695
696 static void prepare_verlet_scheme(FILE                           *fplog,
697                                   t_commrec                      *cr,
698                                   t_inputrec                     *ir,
699                                   int                             nstlist_cmdline,
700                                   const gmx_mtop_t               *mtop,
701                                   matrix                          box,
702                                   gmx_bool                        bUseGPU)
703 {
704     /* For NVE simulations, we will retain the initial list buffer */
705     if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
706     {
707         /* Update the Verlet buffer size for the current run setup */
708         verletbuf_list_setup_t ls;
709         real                   rlist_new;
710
711         /* Here we assume SIMD-enabled kernels are being used. But as currently
712          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
713          * and 4x2 gives a larger buffer than 4x4, this is ok.
714          */
715         verletbuf_get_list_setup(bUseGPU, &ls);
716
717         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
718
719         if (rlist_new != ir->rlist)
720         {
721             if (fplog != NULL)
722             {
723                 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
724                         ir->rlist, rlist_new,
725                         ls.cluster_size_i, ls.cluster_size_j);
726             }
727             ir->rlist     = rlist_new;
728             ir->rlistlong = rlist_new;
729         }
730     }
731
732     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
733     {
734         gmx_fatal(FARGS, "Can not set nstlist without %s",
735                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
736     }
737
738     if (EI_DYNAMICS(ir->eI))
739     {
740         /* Set or try nstlist values */
741         increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
742     }
743 }
744
745 static void convert_to_verlet_scheme(FILE *fplog,
746                                      t_inputrec *ir,
747                                      gmx_mtop_t *mtop, real box_vol)
748 {
749     char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
750
751     md_print_warn(NULL, fplog, "%s\n", conv_mesg);
752
753     ir->cutoff_scheme = ecutsVERLET;
754     ir->verletbuf_tol = 0.005;
755
756     if (ir->rcoulomb != ir->rvdw)
757     {
758         gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
759     }
760
761     if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
762     {
763         gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
764     }
765     else if (ir_vdw_switched(ir) || ir_coulomb_switched(ir))
766     {
767         if (ir_vdw_switched(ir) && ir->vdw_modifier == eintmodNONE)
768         {
769             ir->vdwtype = evdwCUT;
770
771             switch (ir->vdwtype)
772             {
773                 case evdwSHIFT:  ir->vdw_modifier = eintmodFORCESWITCH; break;
774                 case evdwSWITCH: ir->vdw_modifier = eintmodPOTSWITCH; break;
775                 default: gmx_fatal(FARGS, "The Verlet scheme does not support Van der Waals interactions of type '%s'", evdw_names[ir->vdwtype]);
776             }
777         }
778         if (ir_coulomb_switched(ir) && ir->coulomb_modifier == eintmodNONE)
779         {
780             if (EEL_FULL(ir->coulombtype))
781             {
782                 /* With full electrostatic only PME can be switched */
783                 ir->coulombtype      = eelPME;
784                 ir->coulomb_modifier = eintmodPOTSHIFT;
785             }
786             else
787             {
788                 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
789                 ir->coulombtype      = eelRF;
790                 ir->epsilon_rf       = 0.0;
791                 ir->coulomb_modifier = eintmodPOTSHIFT;
792             }
793         }
794
795         /* We set the pair energy error tolerance to a small number.
796          * Note that this is only for testing. For production the user
797          * should think about this and set the mdp options.
798          */
799         ir->verletbuf_tol = 1e-4;
800     }
801
802     if (inputrec2nboundeddim(ir) != 3)
803     {
804         gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
805     }
806
807     if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
808     {
809         gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
810     }
811
812     if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
813     {
814         verletbuf_list_setup_t ls;
815
816         verletbuf_get_list_setup(FALSE, &ls);
817         calc_verlet_buffer_size(mtop, box_vol, ir, -1, &ls, NULL, &ir->rlist);
818     }
819     else
820     {
821         real rlist_fac;
822
823         if (EI_MD(ir->eI))
824         {
825             rlist_fac       = 1 + verlet_buffer_ratio_NVE_T0;
826         }
827         else
828         {
829             rlist_fac       = 1 + verlet_buffer_ratio_nodynamics;
830         }
831         ir->verletbuf_tol   = -1;
832         ir->rlist           = rlist_fac*max(ir->rvdw, ir->rcoulomb);
833     }
834
835     gmx_mtop_remove_chargegroups(mtop);
836 }
837
838 static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
839 {
840     fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
841             hw_opt->nthreads_tot,
842             hw_opt->nthreads_tmpi,
843             hw_opt->nthreads_omp,
844             hw_opt->nthreads_omp_pme,
845             hw_opt->gpu_opt.gpu_id != NULL ? hw_opt->gpu_opt.gpu_id : "");
846 }
847
848 /* Checks we can do when we don't (yet) know the cut-off scheme */
849 static void check_and_update_hw_opt_1(gmx_hw_opt_t *hw_opt,
850                                       gmx_bool      bIsSimMaster)
851 {
852     gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
853
854 #ifndef GMX_THREAD_MPI
855     if (hw_opt->nthreads_tot > 0)
856     {
857         gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
858     }
859     if (hw_opt->nthreads_tmpi > 0)
860     {
861         gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
862     }
863 #endif
864
865 #ifndef GMX_OPENMP
866     if (hw_opt->nthreads_omp > 1)
867     {
868         gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but Gromacs was compiled without OpenMP support");
869     }
870     hw_opt->nthreads_omp = 1;
871 #endif
872
873     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
874     {
875         /* We have the same number of OpenMP threads for PP and PME processes,
876          * thus we can perform several consistency checks.
877          */
878         if (hw_opt->nthreads_tmpi > 0 &&
879             hw_opt->nthreads_omp > 0 &&
880             hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
881         {
882             gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
883                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
884         }
885
886         if (hw_opt->nthreads_tmpi > 0 &&
887             hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
888         {
889             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
890                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
891         }
892
893         if (hw_opt->nthreads_omp > 0 &&
894             hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
895         {
896             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
897                       hw_opt->nthreads_tot, hw_opt->nthreads_omp);
898         }
899
900         if (hw_opt->nthreads_tmpi > 0 &&
901             hw_opt->nthreads_omp <= 0)
902         {
903             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
904         }
905     }
906
907 #ifndef GMX_OPENMP
908     if (hw_opt->nthreads_omp > 1)
909     {
910         gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
911     }
912 #endif
913
914     if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
915     {
916         gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
917     }
918
919     if (hw_opt->nthreads_tot == 1)
920     {
921         hw_opt->nthreads_tmpi = 1;
922
923         if (hw_opt->nthreads_omp > 1)
924         {
925             gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
926                       hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
927         }
928         hw_opt->nthreads_omp = 1;
929     }
930
931     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
932     {
933         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
934     }
935
936     /* Parse GPU IDs, if provided.
937      * We check consistency with the tMPI thread count later.
938      */
939     gmx_parse_gpu_ids(&hw_opt->gpu_opt);
940
941 #ifdef GMX_THREAD_MPI
942     if (hw_opt->gpu_opt.ncuda_dev_use > 0 && hw_opt->nthreads_tmpi == 0)
943     {
944         /* Set the number of MPI threads equal to the number of GPUs */
945         hw_opt->nthreads_tmpi = hw_opt->gpu_opt.ncuda_dev_use;
946
947         if (hw_opt->nthreads_tot > 0 &&
948             hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
949         {
950             /* We have more GPUs than total threads requested.
951              * We choose to (later) generate a mismatch error,
952              * instead of launching more threads than requested.
953              */
954             hw_opt->nthreads_tmpi = hw_opt->nthreads_tot;
955         }
956     }
957 #endif
958
959     if (debug)
960     {
961         print_hw_opt(debug, hw_opt);
962     }
963 }
964
965 /* Checks we can do when we know the cut-off scheme */
966 static void check_and_update_hw_opt_2(gmx_hw_opt_t *hw_opt,
967                                       int           cutoff_scheme)
968 {
969     if (cutoff_scheme == ecutsGROUP)
970     {
971         /* We only have OpenMP support for PME only nodes */
972         if (hw_opt->nthreads_omp > 1)
973         {
974             gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
975                       ecutscheme_names[cutoff_scheme],
976                       ecutscheme_names[ecutsVERLET]);
977         }
978         hw_opt->nthreads_omp = 1;
979     }
980
981     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
982     {
983         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
984     }
985
986     if (debug)
987     {
988         print_hw_opt(debug, hw_opt);
989     }
990 }
991
992
993 /* Override the value in inputrec with value passed on the command line (if any) */
994 static void override_nsteps_cmdline(FILE            *fplog,
995                                     gmx_int64_t      nsteps_cmdline,
996                                     t_inputrec      *ir,
997                                     const t_commrec *cr)
998 {
999     char sbuf[STEPSTRSIZE];
1000
1001     assert(ir);
1002     assert(cr);
1003
1004     /* override with anything else than the default -2 */
1005     if (nsteps_cmdline > -2)
1006     {
1007         char stmp[STRLEN];
1008
1009         ir->nsteps = nsteps_cmdline;
1010         if (EI_DYNAMICS(ir->eI))
1011         {
1012             sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
1013                     gmx_step_str(nsteps_cmdline, sbuf),
1014                     nsteps_cmdline*ir->delta_t);
1015         }
1016         else
1017         {
1018             sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
1019                     gmx_step_str(nsteps_cmdline, sbuf));
1020         }
1021
1022         md_print_warn(cr, fplog, "%s\n", stmp);
1023     }
1024 }
1025
1026 /* Frees GPU memory and destroys the CUDA context.
1027  *
1028  * Note that this function needs to be called even if GPUs are not used
1029  * in this run because the PME ranks have no knowledge of whether GPUs
1030  * are used or not, but all ranks need to enter the barrier below.
1031  */
1032 static void free_gpu_resources(const t_forcerec *fr,
1033                                const t_commrec  *cr)
1034 {
1035     gmx_bool bIsPPrankUsingGPU;
1036     char     gpu_err_str[STRLEN];
1037
1038     bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU;
1039
1040     if (bIsPPrankUsingGPU)
1041     {
1042         /* free nbnxn data in GPU memory */
1043         nbnxn_cuda_free(fr->nbv->cu_nbv);
1044
1045         /* With tMPI we need to wait for all ranks to finish deallocation before
1046          * destroying the context in free_gpu() as some ranks may be sharing
1047          * GPU and context.
1048          * Note: as only PP ranks need to free GPU resources, so it is safe to
1049          * not call the barrier on PME ranks.
1050          */
1051 #ifdef GMX_THREAD_MPI
1052         if (PAR(cr))
1053         {
1054             gmx_barrier(cr);
1055         }
1056 #endif  /* GMX_THREAD_MPI */
1057
1058         /* uninitialize GPU (by destroying the context) */
1059         if (!free_gpu(gpu_err_str))
1060         {
1061             gmx_warning("On node %d failed to free GPU #%d: %s",
1062                         cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1063         }
1064     }
1065 }
1066
1067 int mdrunner(gmx_hw_opt_t *hw_opt,
1068              FILE *fplog, t_commrec *cr, int nfile,
1069              const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
1070              gmx_bool bCompact, int nstglobalcomm,
1071              ivec ddxyz, int dd_node_order, real rdd, real rconstr,
1072              const char *dddlb_opt, real dlb_scale,
1073              const char *ddcsx, const char *ddcsy, const char *ddcsz,
1074              const char *nbpu_opt, int nstlist_cmdline,
1075              gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
1076              int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
1077              int repl_ex_seed, real pforce, real cpt_period, real max_hours,
1078              const char *deviceOptions, unsigned long Flags)
1079 {
1080     gmx_bool                  bForceUseGPU, bTryUseGPU;
1081     double                    nodetime = 0, realtime;
1082     t_inputrec               *inputrec;
1083     t_state                  *state = NULL;
1084     matrix                    box;
1085     gmx_ddbox_t               ddbox = {0};
1086     int                       npme_major, npme_minor;
1087     real                      tmpr1, tmpr2;
1088     t_nrnb                   *nrnb;
1089     gmx_mtop_t               *mtop          = NULL;
1090     t_mdatoms                *mdatoms       = NULL;
1091     t_forcerec               *fr            = NULL;
1092     t_fcdata                 *fcd           = NULL;
1093     real                      ewaldcoeff_q  = 0;
1094     real                      ewaldcoeff_lj = 0;
1095     gmx_pme_t                *pmedata       = NULL;
1096     gmx_vsite_t              *vsite         = NULL;
1097     gmx_constr_t              constr;
1098     int                       i, m, nChargePerturbed = -1, nTypePerturbed = 0, status, nalloc;
1099     char                     *gro;
1100     gmx_wallcycle_t           wcycle;
1101     gmx_bool                  bReadRNG, bReadEkin;
1102     int                       list;
1103     gmx_walltime_accounting_t walltime_accounting = NULL;
1104     int                       rc;
1105     gmx_int64_t               reset_counters;
1106     gmx_edsam_t               ed           = NULL;
1107     t_commrec                *cr_old       = cr;
1108     int                       nthreads_pme = 1;
1109     int                       nthreads_pp  = 1;
1110     gmx_membed_t              membed       = NULL;
1111     gmx_hw_info_t            *hwinfo       = NULL;
1112     /* The master rank decides early on bUseGPU and broadcasts this later */
1113     gmx_bool                  bUseGPU      = FALSE;
1114
1115     /* CAUTION: threads may be started later on in this function, so
1116        cr doesn't reflect the final parallel state right now */
1117     snew(inputrec, 1);
1118     snew(mtop, 1);
1119
1120     if (Flags & MD_APPENDFILES)
1121     {
1122         fplog = NULL;
1123     }
1124
1125     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1126     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1127
1128     /* Detect hardware, gather information. This is an operation that is
1129      * global for this process (MPI rank). */
1130     hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
1131
1132
1133     snew(state, 1);
1134     if (SIMMASTER(cr))
1135     {
1136         /* Read (nearly) all data required for the simulation */
1137         read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
1138
1139         if (inputrec->cutoff_scheme != ecutsVERLET &&
1140             ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1141         {
1142             convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
1143         }
1144
1145         if (inputrec->cutoff_scheme == ecutsVERLET)
1146         {
1147             /* Here the master rank decides if all ranks will use GPUs */
1148             bUseGPU = (hwinfo->gpu_info.ncuda_dev_compatible > 0 ||
1149                        getenv("GMX_EMULATE_GPU") != NULL);
1150
1151             /* TODO add GPU kernels for this and replace this check by:
1152              * (bUseGPU && (ir->vdwtype == evdwPME &&
1153              *               ir->ljpme_combination_rule == eljpmeLB))
1154              * update the message text and the content of nbnxn_acceleration_supported.
1155              */
1156             if (bUseGPU &&
1157                 !nbnxn_acceleration_supported(fplog, cr, inputrec, bUseGPU))
1158             {
1159                 /* Fallback message printed by nbnxn_acceleration_supported */
1160                 if (bForceUseGPU)
1161                 {
1162                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
1163                 }
1164                 bUseGPU = FALSE;
1165             }
1166
1167             prepare_verlet_scheme(fplog, cr,
1168                                   inputrec, nstlist_cmdline, mtop, state->box,
1169                                   bUseGPU);
1170         }
1171         else
1172         {
1173             if (nstlist_cmdline > 0)
1174             {
1175                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
1176             }
1177
1178             if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
1179             {
1180                 md_print_warn(cr, fplog,
1181                               "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1182                               "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1183                               "      (for quick performance testing you can use the -testverlet option)\n");
1184             }
1185
1186             if (bForceUseGPU)
1187             {
1188                 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1189             }
1190
1191 #ifdef GMX_TARGET_BGQ
1192             md_print_warn(cr, fplog,
1193                           "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
1194                           "      BlueGene/Q. You will observe better performance from using the\n"
1195                           "      Verlet cut-off scheme.\n");
1196 #endif
1197         }
1198     }
1199
1200     /* Check for externally set OpenMP affinity and turn off internal
1201      * pinning if any is found. We need to do this check early to tell
1202      * thread-MPI whether it should do pinning when spawning threads.
1203      * TODO: the above no longer holds, we should move these checks down
1204      */
1205     gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1206
1207     /* Check and update the hardware options for internal consistency */
1208     check_and_update_hw_opt_1(hw_opt, SIMMASTER(cr));
1209
1210     if (SIMMASTER(cr))
1211     {
1212 #ifdef GMX_THREAD_MPI
1213         /* Early check for externally set process affinity.
1214          * With thread-MPI this is needed as pinning might get turned off,
1215          * which needs to be known before starting thread-MPI.
1216          * With thread-MPI hw_opt is processed here on the master rank
1217          * and passed to the other ranks later, so we only do this on master.
1218          */
1219         gmx_check_thread_affinity_set(fplog,
1220                                       NULL,
1221                                       hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1222 #endif
1223
1224 #ifdef GMX_THREAD_MPI
1225         if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1226         {
1227             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1228         }
1229 #endif
1230
1231         if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1232             cr->npmenodes <= 0)
1233         {
1234             gmx_fatal(FARGS, "You need to explicitly specify the number of PME nodes (-npme) when using different number of OpenMP threads for PP and PME nodes");
1235         }
1236     }
1237
1238 #ifdef GMX_THREAD_MPI
1239     if (SIMMASTER(cr))
1240     {
1241         /* Since the master knows the cut-off scheme, update hw_opt for this.
1242          * This is done later for normal MPI and also once more with tMPI
1243          * for all tMPI ranks.
1244          */
1245         check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1246
1247         /* NOW the threads will be started: */
1248         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1249                                                  hw_opt,
1250                                                  inputrec, mtop,
1251                                                  cr, fplog);
1252         if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1253         {
1254             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1255         }
1256
1257         if (hw_opt->nthreads_tmpi > 1)
1258         {
1259             /* now start the threads. */
1260             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1261                                         oenv, bVerbose, bCompact, nstglobalcomm,
1262                                         ddxyz, dd_node_order, rdd, rconstr,
1263                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1264                                         nbpu_opt, nstlist_cmdline,
1265                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
1266                                         repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1267                                         cpt_period, max_hours, deviceOptions,
1268                                         Flags);
1269             /* the main thread continues here with a new cr. We don't deallocate
1270                the old cr because other threads may still be reading it. */
1271             if (cr == NULL)
1272             {
1273                 gmx_comm("Failed to spawn threads");
1274             }
1275         }
1276     }
1277 #endif
1278     /* END OF CAUTION: cr is now reliable */
1279
1280     /* g_membed initialisation *
1281      * Because we change the mtop, init_membed is called before the init_parallel *
1282      * (in case we ever want to make it run in parallel) */
1283     if (opt2bSet("-membed", nfile, fnm))
1284     {
1285         if (MASTER(cr))
1286         {
1287             fprintf(stderr, "Initializing membed");
1288         }
1289         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1290     }
1291
1292     if (PAR(cr))
1293     {
1294         /* now broadcast everything to the non-master nodes/threads: */
1295         init_parallel(cr, inputrec, mtop);
1296     }
1297     if (fplog != NULL)
1298     {
1299         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1300     }
1301
1302     /* now make sure the state is initialized and propagated */
1303     set_state_entries(state, inputrec, cr->nnodes);
1304
1305     /* A parallel command line option consistency check that we can
1306        only do after any threads have started. */
1307     if (!PAR(cr) &&
1308         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1309     {
1310         gmx_fatal(FARGS,
1311                   "The -dd or -npme option request a parallel simulation, "
1312 #ifndef GMX_MPI
1313                   "but %s was compiled without threads or MPI enabled"
1314 #else
1315 #ifdef GMX_THREAD_MPI
1316                   "but the number of threads (option -nt) is 1"
1317 #else
1318                   "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1319 #endif
1320 #endif
1321                   , ShortProgram()
1322                   );
1323     }
1324
1325     if ((Flags & MD_RERUN) &&
1326         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1327     {
1328         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1329     }
1330
1331     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
1332     {
1333         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
1334     }
1335
1336     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1337     {
1338         if (cr->npmenodes > 0)
1339         {
1340             gmx_fatal_collective(FARGS, cr, NULL,
1341                                  "PME nodes are requested, but the system does not use PME electrostatics or LJ-PME");
1342         }
1343
1344         cr->npmenodes = 0;
1345     }
1346
1347 #ifdef GMX_FAHCORE
1348     if (MASTER(cr))
1349     {
1350         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1351     }
1352 #endif
1353
1354     /* NMR restraints must be initialized before load_checkpoint,
1355      * since with time averaging the history is added to t_state.
1356      * For proper consistency check we therefore need to extend
1357      * t_state here.
1358      * So the PME-only nodes (if present) will also initialize
1359      * the distance restraints.
1360      */
1361     snew(fcd, 1);
1362
1363     /* This needs to be called before read_checkpoint to extend the state */
1364     init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
1365
1366     init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
1367                 state);
1368
1369     if (DEFORM(*inputrec))
1370     {
1371         /* Store the deform reference box before reading the checkpoint */
1372         if (SIMMASTER(cr))
1373         {
1374             copy_mat(state->box, box);
1375         }
1376         if (PAR(cr))
1377         {
1378             gmx_bcast(sizeof(box), box, cr);
1379         }
1380         /* Because we do not have the update struct available yet
1381          * in which the reference values should be stored,
1382          * we store them temporarily in static variables.
1383          * This should be thread safe, since they are only written once
1384          * and with identical values.
1385          */
1386         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1387         deform_init_init_step_tpx = inputrec->init_step;
1388         copy_mat(box, deform_init_box_tpx);
1389         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1390     }
1391
1392     if (opt2bSet("-cpi", nfile, fnm))
1393     {
1394         /* Check if checkpoint file exists before doing continuation.
1395          * This way we can use identical input options for the first and subsequent runs...
1396          */
1397         if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1398         {
1399             load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1400                             cr, ddxyz,
1401                             inputrec, state, &bReadRNG, &bReadEkin,
1402                             (Flags & MD_APPENDFILES),
1403                             (Flags & MD_APPENDFILESSET));
1404
1405             if (bReadRNG)
1406             {
1407                 Flags |= MD_READ_RNG;
1408             }
1409             if (bReadEkin)
1410             {
1411                 Flags |= MD_READ_EKIN;
1412             }
1413         }
1414     }
1415
1416     if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1417 #ifdef GMX_THREAD_MPI
1418         /* With thread MPI only the master node/thread exists in mdrun.c,
1419          * therefore non-master nodes need to open the "seppot" log file here.
1420          */
1421         || (!MASTER(cr) && (Flags & MD_SEPPOT))
1422 #endif
1423         )
1424     {
1425         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1426                      Flags, &fplog);
1427     }
1428
1429     /* override nsteps with value from cmdline */
1430     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1431
1432     if (SIMMASTER(cr))
1433     {
1434         copy_mat(state->box, box);
1435     }
1436
1437     if (PAR(cr))
1438     {
1439         gmx_bcast(sizeof(box), box, cr);
1440     }
1441
1442     /* Essential dynamics */
1443     if (opt2bSet("-ei", nfile, fnm))
1444     {
1445         /* Open input and output files, allocate space for ED data structure */
1446         ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1447     }
1448
1449     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1450                      inputrec->eI == eiNM))
1451     {
1452         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1453                                            dddlb_opt, dlb_scale,
1454                                            ddcsx, ddcsy, ddcsz,
1455                                            mtop, inputrec,
1456                                            box, state->x,
1457                                            &ddbox, &npme_major, &npme_minor);
1458
1459         make_dd_communicators(fplog, cr, dd_node_order);
1460
1461         /* Set overallocation to avoid frequent reallocation of arrays */
1462         set_over_alloc_dd(TRUE);
1463     }
1464     else
1465     {
1466         /* PME, if used, is done on all nodes with 1D decomposition */
1467         cr->npmenodes = 0;
1468         cr->duty      = (DUTY_PP | DUTY_PME);
1469         npme_major    = 1;
1470         npme_minor    = 1;
1471
1472         if (inputrec->ePBC == epbcSCREW)
1473         {
1474             gmx_fatal(FARGS,
1475                       "pbc=%s is only implemented with domain decomposition",
1476                       epbc_names[inputrec->ePBC]);
1477         }
1478     }
1479
1480     if (PAR(cr))
1481     {
1482         /* After possible communicator splitting in make_dd_communicators.
1483          * we can set up the intra/inter node communication.
1484          */
1485         gmx_setup_nodecomm(fplog, cr);
1486     }
1487
1488     /* Initialize per-physical-node MPI process/thread ID and counters. */
1489     gmx_init_intranode_counters(cr);
1490
1491 #ifdef GMX_MPI
1492     md_print_info(cr, fplog, "Using %d MPI %s\n",
1493                   cr->nnodes,
1494 #ifdef GMX_THREAD_MPI
1495                   cr->nnodes == 1 ? "thread" : "threads"
1496 #else
1497                   cr->nnodes == 1 ? "process" : "processes"
1498 #endif
1499                   );
1500     fflush(stderr);
1501 #endif
1502
1503     /* Check and update hw_opt for the cut-off scheme */
1504     check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1505
1506     gmx_omp_nthreads_init(fplog, cr,
1507                           hwinfo->nthreads_hw_avail,
1508                           hw_opt->nthreads_omp,
1509                           hw_opt->nthreads_omp_pme,
1510                           (cr->duty & DUTY_PP) == 0,
1511                           inputrec->cutoff_scheme == ecutsVERLET);
1512
1513     if (PAR(cr))
1514     {
1515         /* The master rank decided on the use of GPUs,
1516          * broadcast this information to all ranks.
1517          */
1518         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
1519     }
1520
1521     if (bUseGPU)
1522     {
1523         if (cr->npmenodes == -1)
1524         {
1525             /* Don't automatically use PME-only nodes with GPUs */
1526             cr->npmenodes = 0;
1527         }
1528
1529         /* Select GPU id's to use */
1530         gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1531                            &hw_opt->gpu_opt);
1532     }
1533     else
1534     {
1535         /* Ignore (potentially) manually selected GPUs */
1536         hw_opt->gpu_opt.ncuda_dev_use = 0;
1537     }
1538
1539     /* check consistency across ranks of things like SIMD
1540      * support and number of GPUs selected */
1541     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1542
1543     if (DOMAINDECOMP(cr))
1544     {
1545         /* When we share GPUs over ranks, we need to know this for the DLB */
1546         dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1547     }
1548
1549     /* getting number of PP/PME threads
1550        PME: env variable should be read only on one node to make sure it is
1551        identical everywhere;
1552      */
1553     /* TODO nthreads_pp is only used for pinning threads.
1554      * This is a temporary solution until we have a hw topology library.
1555      */
1556     nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
1557     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1558
1559     wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1560
1561     if (PAR(cr))
1562     {
1563         /* Master synchronizes its value of reset_counters with all nodes
1564          * including PME only nodes */
1565         reset_counters = wcycle_get_reset_counters(wcycle);
1566         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1567         wcycle_set_reset_counters(wcycle, reset_counters);
1568     }
1569
1570     snew(nrnb, 1);
1571     if (cr->duty & DUTY_PP)
1572     {
1573         bcast_state(cr, state);
1574
1575         /* Initiate forcerecord */
1576         fr          = mk_forcerec();
1577         fr->hwinfo  = hwinfo;
1578         fr->gpu_opt = &hw_opt->gpu_opt;
1579         init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1580                       opt2fn("-table", nfile, fnm),
1581                       opt2fn("-tabletf", nfile, fnm),
1582                       opt2fn("-tablep", nfile, fnm),
1583                       opt2fn("-tableb", nfile, fnm),
1584                       nbpu_opt,
1585                       FALSE,
1586                       pforce);
1587
1588         /* version for PCA_NOT_READ_NODE (see md.c) */
1589         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1590            "nofile","nofile","nofile","nofile",FALSE,pforce);
1591          */
1592         fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1593
1594         /* Initialize QM-MM */
1595         if (fr->bQMMM)
1596         {
1597             init_QMMMrec(cr, mtop, inputrec, fr);
1598         }
1599
1600         /* Initialize the mdatoms structure.
1601          * mdatoms is not filled with atom data,
1602          * as this can not be done now with domain decomposition.
1603          */
1604         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1605
1606         if (mdatoms->nPerturbed > 0 && inputrec->cutoff_scheme == ecutsVERLET)
1607         {
1608             gmx_fatal(FARGS, "The Verlet cut-off scheme does not (yet) support free-energy calculations with perturbed atoms, only perturbed interactions. This will be implemented soon. Use the group scheme for now.");
1609         }
1610
1611         /* Initialize the virtual site communication */
1612         vsite = init_vsite(mtop, cr, FALSE);
1613
1614         calc_shifts(box, fr->shift_vec);
1615
1616         /* With periodic molecules the charge groups should be whole at start up
1617          * and the virtual sites should not be far from their proper positions.
1618          */
1619         if (!inputrec->bContinuation && MASTER(cr) &&
1620             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1621         {
1622             /* Make molecules whole at start of run */
1623             if (fr->ePBC != epbcNONE)
1624             {
1625                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1626             }
1627             if (vsite)
1628             {
1629                 /* Correct initial vsite positions are required
1630                  * for the initial distribution in the domain decomposition
1631                  * and for the initial shell prediction.
1632                  */
1633                 construct_vsites_mtop(vsite, mtop, state->x);
1634             }
1635         }
1636
1637         if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1638         {
1639             ewaldcoeff_q  = fr->ewaldcoeff_q;
1640             ewaldcoeff_lj = fr->ewaldcoeff_lj;
1641             pmedata       = &fr->pmedata;
1642         }
1643         else
1644         {
1645             pmedata = NULL;
1646         }
1647     }
1648     else
1649     {
1650         /* This is a PME only node */
1651
1652         /* We don't need the state */
1653         done_state(state);
1654
1655         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1656         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1657         snew(pmedata, 1);
1658     }
1659
1660     if (hw_opt->thread_affinity != threadaffOFF)
1661     {
1662         /* Before setting affinity, check whether the affinity has changed
1663          * - which indicates that probably the OpenMP library has changed it
1664          * since we first checked).
1665          */
1666         gmx_check_thread_affinity_set(fplog, cr,
1667                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1668
1669         /* Set the CPU affinity */
1670         gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1671     }
1672
1673     /* Initiate PME if necessary,
1674      * either on all nodes or on dedicated PME nodes only. */
1675     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1676     {
1677         if (mdatoms)
1678         {
1679             nChargePerturbed = mdatoms->nChargePerturbed;
1680             if (EVDW_PME(inputrec->vdwtype))
1681             {
1682                 nTypePerturbed   = mdatoms->nTypePerturbed;
1683             }
1684         }
1685         if (cr->npmenodes > 0)
1686         {
1687             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1688             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1689             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1690         }
1691
1692         if (cr->duty & DUTY_PME)
1693         {
1694             status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1695                                   mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1696                                   (Flags & MD_REPRODUCIBLE), nthreads_pme);
1697             if (status != 0)
1698             {
1699                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1700             }
1701         }
1702     }
1703
1704
1705     if (integrator[inputrec->eI].func == do_md)
1706     {
1707         /* Turn on signal handling on all nodes */
1708         /*
1709          * (A user signal from the PME nodes (if any)
1710          * is communicated to the PP nodes.
1711          */
1712         signal_handler_install();
1713     }
1714
1715     if (cr->duty & DUTY_PP)
1716     {
1717         /* Assumes uniform use of the number of OpenMP threads */
1718         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1719
1720         if (inputrec->ePull != epullNO)
1721         {
1722             /* Initialize pull code */
1723             init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1724                       EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1725         }
1726
1727         if (inputrec->bRot)
1728         {
1729             /* Initialize enforced rotation code */
1730             init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1731                      bVerbose, Flags);
1732         }
1733
1734         if (inputrec->eSwapCoords != eswapNO)
1735         {
1736             /* Initialize ion swapping code */
1737             init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1738                             mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1739         }
1740
1741         constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1742
1743         if (DOMAINDECOMP(cr))
1744         {
1745             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1746                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1747
1748             set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1749
1750             setup_dd_grid(fplog, cr->dd);
1751         }
1752
1753         /* Now do whatever the user wants us to do (how flexible...) */
1754         integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1755                                       oenv, bVerbose, bCompact,
1756                                       nstglobalcomm,
1757                                       vsite, constr,
1758                                       nstepout, inputrec, mtop,
1759                                       fcd, state,
1760                                       mdatoms, nrnb, wcycle, ed, fr,
1761                                       repl_ex_nst, repl_ex_nex, repl_ex_seed,
1762                                       membed,
1763                                       cpt_period, max_hours,
1764                                       deviceOptions,
1765                                       Flags,
1766                                       walltime_accounting);
1767
1768         if (inputrec->ePull != epullNO)
1769         {
1770             finish_pull(inputrec->pull);
1771         }
1772
1773         if (inputrec->bRot)
1774         {
1775             finish_rot(inputrec->rot);
1776         }
1777
1778     }
1779     else
1780     {
1781         /* do PME only */
1782         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1783         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1784     }
1785
1786     wallcycle_stop(wcycle, ewcRUN);
1787
1788     /* Finish up, write some stuff
1789      * if rerunMD, don't write last frame again
1790      */
1791     finish_run(fplog, cr,
1792                inputrec, nrnb, wcycle, walltime_accounting,
1793                fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1794                nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1795                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1796
1797
1798     /* Free GPU memory and context */
1799     free_gpu_resources(fr, cr);
1800
1801     if (opt2bSet("-membed", nfile, fnm))
1802     {
1803         sfree(membed);
1804     }
1805
1806     gmx_hardware_info_free(hwinfo);
1807
1808     /* Does what it says */
1809     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", walltime_accounting);
1810     walltime_accounting_destroy(walltime_accounting);
1811
1812     /* Close logfile already here if we were appending to it */
1813     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1814     {
1815         gmx_log_close(fplog);
1816     }
1817
1818     rc = (int)gmx_get_stop_condition();
1819
1820 #ifdef GMX_THREAD_MPI
1821     /* we need to join all threads. The sub-threads join when they
1822        exit this function, but the master thread needs to be told to
1823        wait for that. */
1824     if (PAR(cr) && MASTER(cr))
1825     {
1826         tMPI_Finalize();
1827     }
1828 #endif
1829
1830     return rc;
1831 }