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