Merge release-5-0 into release-5-1
[alexxy/gromacs.git] / src / programs / mdrun / runner.cpp
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,2015, 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
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <assert.h>
43 #include <signal.h>
44 #include <stdlib.h>
45 #include <string.h>
46
47 #include <algorithm>
48
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/essentialdynamics/edsam.h"
51 #include "gromacs/ewald/pme.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
54 #include "gromacs/legacyheaders/checkpoint.h"
55 #include "gromacs/legacyheaders/constr.h"
56 #include "gromacs/legacyheaders/copyrite.h"
57 #include "gromacs/legacyheaders/disre.h"
58 #include "gromacs/legacyheaders/force.h"
59 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
60 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
61 #include "gromacs/legacyheaders/gmx_thread_affinity.h"
62 #include "gromacs/legacyheaders/inputrec.h"
63 #include "gromacs/legacyheaders/main.h"
64 #include "gromacs/legacyheaders/md_logging.h"
65 #include "gromacs/legacyheaders/md_support.h"
66 #include "gromacs/legacyheaders/mdatoms.h"
67 #include "gromacs/legacyheaders/mdrun.h"
68 #include "gromacs/legacyheaders/names.h"
69 #include "gromacs/legacyheaders/network.h"
70 #include "gromacs/legacyheaders/oenv.h"
71 #include "gromacs/legacyheaders/orires.h"
72 #include "gromacs/legacyheaders/qmmm.h"
73 #include "gromacs/legacyheaders/sighandler.h"
74 #include "gromacs/legacyheaders/txtdump.h"
75 #include "gromacs/legacyheaders/typedefs.h"
76 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
77 #include "gromacs/math/vec.h"
78 #include "gromacs/mdlib/calc_verletbuf.h"
79 #include "gromacs/mdlib/nbnxn_search.h"
80 #include "gromacs/pbcutil/pbc.h"
81 #include "gromacs/pulling/pull.h"
82 #include "gromacs/pulling/pull_rotation.h"
83 #include "gromacs/swap/swapcoords.h"
84 #include "gromacs/timing/wallcycle.h"
85 #include "gromacs/topology/mtop_util.h"
86 #include "gromacs/utility/gmxassert.h"
87 #include "gromacs/utility/gmxmpi.h"
88 #include "gromacs/utility/smalloc.h"
89
90 #include "deform.h"
91 #include "membed.h"
92 #include "repl_ex.h"
93 #include "resource-division.h"
94
95 #ifdef GMX_FAHCORE
96 #include "corewrap.h"
97 #endif
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 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
113  * the number of threads will get lowered.
114  */
115 #define MIN_ATOMS_PER_MPI_THREAD    90
116 #define MIN_ATOMS_PER_GPU           900
117
118 struct mdrunner_arglist
119 {
120     gmx_hw_opt_t    hw_opt;
121     FILE           *fplog;
122     t_commrec      *cr;
123     int             nfile;
124     const t_filenm *fnm;
125     output_env_t    oenv;
126     gmx_bool        bVerbose;
127     gmx_bool        bCompact;
128     int             nstglobalcomm;
129     ivec            ddxyz;
130     int             dd_node_order;
131     real            rdd;
132     real            rconstr;
133     const char     *dddlb_opt;
134     real            dlb_scale;
135     const char     *ddcsx;
136     const char     *ddcsy;
137     const char     *ddcsz;
138     const char     *nbpu_opt;
139     int             nstlist_cmdline;
140     gmx_int64_t     nsteps_cmdline;
141     int             nstepout;
142     int             resetstep;
143     int             nmultisim;
144     int             repl_ex_nst;
145     int             repl_ex_nex;
146     int             repl_ex_seed;
147     real            pforce;
148     real            cpt_period;
149     real            max_hours;
150     int             imdport;
151     unsigned long   Flags;
152 };
153
154
155 /* The function used for spawning threads. Extracts the mdrunner()
156    arguments from its one argument and calls mdrunner(), after making
157    a commrec. */
158 static void mdrunner_start_fn(void *arg)
159 {
160     struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
161     struct mdrunner_arglist  mc  = *mda; /* copy the arg list to make sure
162                                             that it's thread-local. This doesn't
163                                             copy pointed-to items, of course,
164                                             but those are all const. */
165     t_commrec *cr;                       /* we need a local version of this */
166     FILE      *fplog = NULL;
167     t_filenm  *fnm;
168
169     fnm = dup_tfn(mc.nfile, mc.fnm);
170
171     cr = reinitialize_commrec_for_this_thread(mc.cr);
172
173     if (MASTER(cr))
174     {
175         fplog = mc.fplog;
176     }
177
178     mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
179              mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
180              mc.ddxyz, mc.dd_node_order, mc.rdd,
181              mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
182              mc.ddcsx, mc.ddcsy, mc.ddcsz,
183              mc.nbpu_opt, mc.nstlist_cmdline,
184              mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
185              mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
186              mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
187 }
188
189 /* called by mdrunner() to start a specific number of threads (including
190    the main thread) for thread-parallel runs. This in turn calls mdrunner()
191    for each thread.
192    All options besides nthreads are the same as for mdrunner(). */
193 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
194                                          FILE *fplog, t_commrec *cr, int nfile,
195                                          const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
196                                          gmx_bool bCompact, int nstglobalcomm,
197                                          ivec ddxyz, int dd_node_order, real rdd, real rconstr,
198                                          const char *dddlb_opt, real dlb_scale,
199                                          const char *ddcsx, const char *ddcsy, const char *ddcsz,
200                                          const char *nbpu_opt, int nstlist_cmdline,
201                                          gmx_int64_t nsteps_cmdline,
202                                          int nstepout, int resetstep,
203                                          int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
204                                          real pforce, real cpt_period, real max_hours,
205                                          unsigned long Flags)
206 {
207     int                      ret;
208     struct mdrunner_arglist *mda;
209     t_commrec               *crn; /* the new commrec */
210     t_filenm                *fnmn;
211
212     /* first check whether we even need to start tMPI */
213     if (hw_opt->nthreads_tmpi < 2)
214     {
215         return cr;
216     }
217
218     /* a few small, one-time, almost unavoidable memory leaks: */
219     snew(mda, 1);
220     fnmn = dup_tfn(nfile, fnm);
221
222     /* fill the data structure to pass as void pointer to thread start fn */
223     /* hw_opt contains pointers, which should all be NULL at this stage */
224     mda->hw_opt          = *hw_opt;
225     mda->fplog           = fplog;
226     mda->cr              = cr;
227     mda->nfile           = nfile;
228     mda->fnm             = fnmn;
229     mda->oenv            = oenv;
230     mda->bVerbose        = bVerbose;
231     mda->bCompact        = bCompact;
232     mda->nstglobalcomm   = nstglobalcomm;
233     mda->ddxyz[XX]       = ddxyz[XX];
234     mda->ddxyz[YY]       = ddxyz[YY];
235     mda->ddxyz[ZZ]       = ddxyz[ZZ];
236     mda->dd_node_order   = dd_node_order;
237     mda->rdd             = rdd;
238     mda->rconstr         = rconstr;
239     mda->dddlb_opt       = dddlb_opt;
240     mda->dlb_scale       = dlb_scale;
241     mda->ddcsx           = ddcsx;
242     mda->ddcsy           = ddcsy;
243     mda->ddcsz           = ddcsz;
244     mda->nbpu_opt        = nbpu_opt;
245     mda->nstlist_cmdline = nstlist_cmdline;
246     mda->nsteps_cmdline  = nsteps_cmdline;
247     mda->nstepout        = nstepout;
248     mda->resetstep       = resetstep;
249     mda->nmultisim       = nmultisim;
250     mda->repl_ex_nst     = repl_ex_nst;
251     mda->repl_ex_nex     = repl_ex_nex;
252     mda->repl_ex_seed    = repl_ex_seed;
253     mda->pforce          = pforce;
254     mda->cpt_period      = cpt_period;
255     mda->max_hours       = max_hours;
256     mda->Flags           = Flags;
257
258     /* now spawn new threads that start mdrunner_start_fn(), while
259        the main thread returns, we set thread affinity later */
260     ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
261                        mdrunner_start_fn, (void*)(mda) );
262     if (ret != TMPI_SUCCESS)
263     {
264         return NULL;
265     }
266
267     crn = reinitialize_commrec_for_this_thread(cr);
268     return crn;
269 }
270
271 #endif /* GMX_THREAD_MPI */
272
273
274 /* We determine the extra cost of the non-bonded kernels compared to
275  * a reference nstlist value of 10 (which is the default in grompp).
276  */
277 static const int    nbnxnReferenceNstlist = 10;
278 /* The values to try when switching  */
279 const int           nstlist_try[] = { 20, 25, 40 };
280 #define NNSTL  sizeof(nstlist_try)/sizeof(nstlist_try[0])
281 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
282  * but never more than listfac_max.
283  * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
284  * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
285  * Note that both CPU and GPU factors are conservative. Performance should
286  * not go down due to this tuning, except with a relatively slow GPU.
287  * On the other hand, at medium/high parallelization or with fast GPUs
288  * nstlist will not be increased enough to reach optimal performance.
289  */
290 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
291 static const float  nbnxn_cpu_listfac_ok    = 1.05;
292 static const float  nbnxn_cpu_listfac_max   = 1.09;
293 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
294 static const float  nbnxn_gpu_listfac_ok    = 1.20;
295 static const float  nbnxn_gpu_listfac_max   = 1.30;
296
297 /* Try to increase nstlist when using the Verlet cut-off scheme */
298 static void increase_nstlist(FILE *fp, t_commrec *cr,
299                              t_inputrec *ir, int nstlist_cmdline,
300                              const gmx_mtop_t *mtop, matrix box,
301                              gmx_bool bGPU)
302 {
303     float                  listfac_ok, listfac_max;
304     int                    nstlist_orig, nstlist_prev;
305     verletbuf_list_setup_t ls;
306     real                   rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
307     real                   rlist_new, rlist_prev;
308     size_t                 nstlist_ind = 0;
309     t_state                state_tmp;
310     gmx_bool               bBox, bDD, bCont;
311     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";
312     const char            *nve_err  = "Can not increase nstlist because an NVE ensemble is used";
313     const char            *vbd_err  = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
314     const char            *box_err  = "Can not increase nstlist because the box is too small";
315     const char            *dd_err   = "Can not increase nstlist because of domain decomposition limitations";
316     char                   buf[STRLEN];
317     const float            oneThird = 1.0f / 3.0f;
318
319     if (nstlist_cmdline <= 0)
320     {
321         if (ir->nstlist == 1)
322         {
323             /* The user probably set nstlist=1 for a reason,
324              * don't mess with the settings.
325              */
326             return;
327         }
328
329         if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
330         {
331             fprintf(fp, nstl_gpu, ir->nstlist);
332         }
333         nstlist_ind = 0;
334         while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
335         {
336             nstlist_ind++;
337         }
338         if (nstlist_ind == NNSTL)
339         {
340             /* There are no larger nstlist value to try */
341             return;
342         }
343     }
344
345     if (EI_MD(ir->eI) && ir->etc == etcNO)
346     {
347         if (MASTER(cr))
348         {
349             fprintf(stderr, "%s\n", nve_err);
350         }
351         if (fp != NULL)
352         {
353             fprintf(fp, "%s\n", nve_err);
354         }
355
356         return;
357     }
358
359     if (ir->verletbuf_tol == 0 && bGPU)
360     {
361         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");
362     }
363
364     if (ir->verletbuf_tol < 0)
365     {
366         if (MASTER(cr))
367         {
368             fprintf(stderr, "%s\n", vbd_err);
369         }
370         if (fp != NULL)
371         {
372             fprintf(fp, "%s\n", vbd_err);
373         }
374
375         return;
376     }
377
378     if (bGPU)
379     {
380         listfac_ok  = nbnxn_gpu_listfac_ok;
381         listfac_max = nbnxn_gpu_listfac_max;
382     }
383     else
384     {
385         listfac_ok  = nbnxn_cpu_listfac_ok;
386         listfac_max = nbnxn_cpu_listfac_max;
387     }
388
389     nstlist_orig = ir->nstlist;
390     if (nstlist_cmdline > 0)
391     {
392         if (fp)
393         {
394             sprintf(buf, "Getting nstlist=%d from command line option",
395                     nstlist_cmdline);
396         }
397         ir->nstlist = nstlist_cmdline;
398     }
399
400     verletbuf_get_list_setup(TRUE, bGPU, &ls);
401
402     /* Allow rlist to make the list a given factor larger than the list
403      * would be with the reference value for nstlist (10).
404      */
405     nstlist_prev = ir->nstlist;
406     ir->nstlist  = nbnxnReferenceNstlist;
407     calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
408                             &rlistWithReferenceNstlist);
409     ir->nstlist  = nstlist_prev;
410
411     /* Determine the pair list size increase due to zero interactions */
412     rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
413                                               mtop->natoms/det(box));
414     rlist_ok  = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
415     rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
416     if (debug)
417     {
418         fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
419                 rlist_inc, rlist_ok, rlist_max);
420     }
421
422     nstlist_prev = nstlist_orig;
423     rlist_prev   = ir->rlist;
424     do
425     {
426         if (nstlist_cmdline <= 0)
427         {
428             ir->nstlist = nstlist_try[nstlist_ind];
429         }
430
431         /* Set the pair-list buffer size in ir */
432         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
433
434         /* Does rlist fit in the box? */
435         bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
436         bDD  = TRUE;
437         if (bBox && DOMAINDECOMP(cr))
438         {
439             /* Check if rlist fits in the domain decomposition */
440             if (inputrec2nboundeddim(ir) < DIM)
441             {
442                 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
443             }
444             copy_mat(box, state_tmp.box);
445             bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
446         }
447
448         if (debug)
449         {
450             fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
451                     ir->nstlist, rlist_new, bBox, bDD);
452         }
453
454         bCont = FALSE;
455
456         if (nstlist_cmdline <= 0)
457         {
458             if (bBox && bDD && rlist_new <= rlist_max)
459             {
460                 /* Increase nstlist */
461                 nstlist_prev = ir->nstlist;
462                 rlist_prev   = rlist_new;
463                 bCont        = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
464             }
465             else
466             {
467                 /* Stick with the previous nstlist */
468                 ir->nstlist = nstlist_prev;
469                 rlist_new   = rlist_prev;
470                 bBox        = TRUE;
471                 bDD         = TRUE;
472             }
473         }
474
475         nstlist_ind++;
476     }
477     while (bCont);
478
479     if (!bBox || !bDD)
480     {
481         gmx_warning(!bBox ? box_err : dd_err);
482         if (fp != NULL)
483         {
484             fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
485         }
486         ir->nstlist = nstlist_orig;
487     }
488     else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
489     {
490         sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
491                 nstlist_orig, ir->nstlist,
492                 ir->rlist, rlist_new);
493         if (MASTER(cr))
494         {
495             fprintf(stderr, "%s\n\n", buf);
496         }
497         if (fp != NULL)
498         {
499             fprintf(fp, "%s\n\n", buf);
500         }
501         ir->rlist     = rlist_new;
502         ir->rlistlong = rlist_new;
503     }
504 }
505
506 static void prepare_verlet_scheme(FILE                           *fplog,
507                                   t_commrec                      *cr,
508                                   t_inputrec                     *ir,
509                                   int                             nstlist_cmdline,
510                                   const gmx_mtop_t               *mtop,
511                                   matrix                          box,
512                                   gmx_bool                        bUseGPU)
513 {
514     /* For NVE simulations, we will retain the initial list buffer */
515     if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
516     {
517         /* Update the Verlet buffer size for the current run setup */
518         verletbuf_list_setup_t ls;
519         real                   rlist_new;
520
521         /* Here we assume SIMD-enabled kernels are being used. But as currently
522          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
523          * and 4x2 gives a larger buffer than 4x4, this is ok.
524          */
525         verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
526
527         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
528
529         if (rlist_new != ir->rlist)
530         {
531             if (fplog != NULL)
532             {
533                 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
534                         ir->rlist, rlist_new,
535                         ls.cluster_size_i, ls.cluster_size_j);
536             }
537             ir->rlist     = rlist_new;
538             ir->rlistlong = rlist_new;
539         }
540     }
541
542     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
543     {
544         gmx_fatal(FARGS, "Can not set nstlist without %s",
545                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
546     }
547
548     if (EI_DYNAMICS(ir->eI))
549     {
550         /* Set or try nstlist values */
551         increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
552     }
553 }
554
555 /* Override the value in inputrec with value passed on the command line (if any) */
556 static void override_nsteps_cmdline(FILE            *fplog,
557                                     gmx_int64_t      nsteps_cmdline,
558                                     t_inputrec      *ir,
559                                     const t_commrec *cr)
560 {
561     assert(ir);
562     assert(cr);
563
564     /* override with anything else than the default -2 */
565     if (nsteps_cmdline > -2)
566     {
567         char sbuf_steps[STEPSTRSIZE];
568         char sbuf_msg[STRLEN];
569
570         ir->nsteps = nsteps_cmdline;
571         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
572         {
573             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
574                     gmx_step_str(nsteps_cmdline, sbuf_steps),
575                     fabs(nsteps_cmdline*ir->delta_t));
576         }
577         else
578         {
579             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
580                     gmx_step_str(nsteps_cmdline, sbuf_steps));
581         }
582
583         md_print_warn(cr, fplog, "%s\n", sbuf_msg);
584     }
585     else if (nsteps_cmdline < -2)
586     {
587         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
588                   nsteps_cmdline);
589     }
590     /* Do nothing if nsteps_cmdline == -2 */
591 }
592
593 int mdrunner(gmx_hw_opt_t *hw_opt,
594              FILE *fplog, t_commrec *cr, int nfile,
595              const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
596              gmx_bool bCompact, int nstglobalcomm,
597              ivec ddxyz, int dd_node_order, real rdd, real rconstr,
598              const char *dddlb_opt, real dlb_scale,
599              const char *ddcsx, const char *ddcsy, const char *ddcsz,
600              const char *nbpu_opt, int nstlist_cmdline,
601              gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
602              int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
603              int repl_ex_seed, real pforce, real cpt_period, real max_hours,
604              int imdport, unsigned long Flags)
605 {
606     gmx_bool                  bForceUseGPU, bTryUseGPU, bRerunMD;
607     t_inputrec               *inputrec;
608     t_state                  *state = NULL;
609     matrix                    box;
610     gmx_ddbox_t               ddbox = {0};
611     int                       npme_major, npme_minor;
612     t_nrnb                   *nrnb;
613     gmx_mtop_t               *mtop          = NULL;
614     t_mdatoms                *mdatoms       = NULL;
615     t_forcerec               *fr            = NULL;
616     t_fcdata                 *fcd           = NULL;
617     real                      ewaldcoeff_q  = 0;
618     real                      ewaldcoeff_lj = 0;
619     struct gmx_pme_t        **pmedata       = NULL;
620     gmx_vsite_t              *vsite         = NULL;
621     gmx_constr_t              constr;
622     int                       nChargePerturbed = -1, nTypePerturbed = 0, status;
623     gmx_wallcycle_t           wcycle;
624     gmx_bool                  bReadEkin;
625     gmx_walltime_accounting_t walltime_accounting = NULL;
626     int                       rc;
627     gmx_int64_t               reset_counters;
628     gmx_edsam_t               ed           = NULL;
629     int                       nthreads_pme = 1;
630     int                       nthreads_pp  = 1;
631     gmx_membed_t              membed       = NULL;
632     gmx_hw_info_t            *hwinfo       = NULL;
633     /* The master rank decides early on bUseGPU and broadcasts this later */
634     gmx_bool                  bUseGPU      = FALSE;
635
636     /* CAUTION: threads may be started later on in this function, so
637        cr doesn't reflect the final parallel state right now */
638     snew(inputrec, 1);
639     snew(mtop, 1);
640
641     if (Flags & MD_APPENDFILES)
642     {
643         fplog = NULL;
644     }
645
646     bRerunMD     = (Flags & MD_RERUN);
647     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
648     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
649
650     /* Detect hardware, gather information. This is an operation that is
651      * global for this process (MPI rank). */
652     hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
653
654     gmx_print_detected_hardware(fplog, cr, hwinfo);
655
656     if (fplog != NULL)
657     {
658         /* Print references after all software/hardware printing */
659         please_cite(fplog, "Abraham2015");
660         please_cite(fplog, "Pall2015");
661         please_cite(fplog, "Pronk2013");
662         please_cite(fplog, "Hess2008b");
663         please_cite(fplog, "Spoel2005a");
664         please_cite(fplog, "Lindahl2001a");
665         please_cite(fplog, "Berendsen95a");
666     }
667
668     snew(state, 1);
669     if (SIMMASTER(cr))
670     {
671         /* Read (nearly) all data required for the simulation */
672         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, NULL, mtop);
673
674         if (inputrec->cutoff_scheme == ecutsVERLET)
675         {
676             /* Here the master rank decides if all ranks will use GPUs */
677             bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
678                        getenv("GMX_EMULATE_GPU") != NULL);
679
680             /* TODO add GPU kernels for this and replace this check by:
681              * (bUseGPU && (ir->vdwtype == evdwPME &&
682              *               ir->ljpme_combination_rule == eljpmeLB))
683              * update the message text and the content of nbnxn_acceleration_supported.
684              */
685             if (bUseGPU &&
686                 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
687             {
688                 /* Fallback message printed by nbnxn_acceleration_supported */
689                 if (bForceUseGPU)
690                 {
691                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
692                 }
693                 bUseGPU = FALSE;
694             }
695
696             prepare_verlet_scheme(fplog, cr,
697                                   inputrec, nstlist_cmdline, mtop, state->box,
698                                   bUseGPU);
699         }
700         else
701         {
702             if (nstlist_cmdline > 0)
703             {
704                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
705             }
706
707             if (hwinfo->gpu_info.n_dev_compatible > 0)
708             {
709                 md_print_warn(cr, fplog,
710                               "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
711                               "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
712             }
713
714             if (bForceUseGPU)
715             {
716                 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
717             }
718
719 #ifdef GMX_TARGET_BGQ
720             md_print_warn(cr, fplog,
721                           "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
722                           "      BlueGene/Q. You will observe better performance from using the\n"
723                           "      Verlet cut-off scheme.\n");
724 #endif
725         }
726
727         if (inputrec->eI == eiSD2)
728         {
729             md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
730                           "it is slower than integrator %s and is slightly less accurate\n"
731                           "with constraints. Use the %s integrator.",
732                           ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
733         }
734     }
735
736     /* Check and update the hardware options for internal consistency */
737     check_and_update_hw_opt_1(hw_opt, cr);
738
739     /* Early check for externally set process affinity. */
740     gmx_check_thread_affinity_set(fplog, cr,
741                                   hw_opt, hwinfo->nthreads_hw_avail, FALSE);
742
743 #ifdef GMX_THREAD_MPI
744     if (SIMMASTER(cr))
745     {
746         if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
747         {
748             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
749         }
750
751         /* Since the master knows the cut-off scheme, update hw_opt for this.
752          * This is done later for normal MPI and also once more with tMPI
753          * for all tMPI ranks.
754          */
755         check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
756
757         /* NOW the threads will be started: */
758         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
759                                                  hw_opt,
760                                                  inputrec, mtop,
761                                                  cr, fplog, bUseGPU);
762
763         if (hw_opt->nthreads_tmpi > 1)
764         {
765             t_commrec *cr_old       = cr;
766             /* now start the threads. */
767             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
768                                         oenv, bVerbose, bCompact, nstglobalcomm,
769                                         ddxyz, dd_node_order, rdd, rconstr,
770                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
771                                         nbpu_opt, nstlist_cmdline,
772                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
773                                         repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
774                                         cpt_period, max_hours,
775                                         Flags);
776             /* the main thread continues here with a new cr. We don't deallocate
777                the old cr because other threads may still be reading it. */
778             if (cr == NULL)
779             {
780                 gmx_comm("Failed to spawn threads");
781             }
782         }
783     }
784 #endif
785     /* END OF CAUTION: cr is now reliable */
786
787     /* g_membed initialisation *
788      * Because we change the mtop, init_membed is called before the init_parallel *
789      * (in case we ever want to make it run in parallel) */
790     if (opt2bSet("-membed", nfile, fnm))
791     {
792         if (MASTER(cr))
793         {
794             fprintf(stderr, "Initializing membed");
795         }
796         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
797     }
798
799     if (PAR(cr))
800     {
801         /* now broadcast everything to the non-master nodes/threads: */
802         init_parallel(cr, inputrec, mtop);
803
804         /* The master rank decided on the use of GPUs,
805          * broadcast this information to all ranks.
806          */
807         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
808     }
809
810     if (fplog != NULL)
811     {
812         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
813         fprintf(fplog, "\n");
814     }
815
816     /* now make sure the state is initialized and propagated */
817     set_state_entries(state, inputrec);
818
819     /* A parallel command line option consistency check that we can
820        only do after any threads have started. */
821     if (!PAR(cr) &&
822         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
823     {
824         gmx_fatal(FARGS,
825                   "The -dd or -npme option request a parallel simulation, "
826 #ifndef GMX_MPI
827                   "but %s was compiled without threads or MPI enabled"
828 #else
829 #ifdef GMX_THREAD_MPI
830                   "but the number of threads (option -nt) is 1"
831 #else
832                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
833 #endif
834 #endif
835                   , output_env_get_program_display_name(oenv)
836                   );
837     }
838
839     if (bRerunMD &&
840         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
841     {
842         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
843     }
844
845     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
846     {
847         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
848     }
849
850     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
851     {
852         if (cr->npmenodes > 0)
853         {
854             gmx_fatal_collective(FARGS, cr, NULL,
855                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
856         }
857
858         cr->npmenodes = 0;
859     }
860
861     if (bUseGPU && cr->npmenodes < 0)
862     {
863         /* With GPUs we don't automatically use PME-only ranks. PME ranks can
864          * improve performance with many threads per GPU, since our OpenMP
865          * scaling is bad, but it's difficult to automate the setup.
866          */
867         cr->npmenodes = 0;
868     }
869
870 #ifdef GMX_FAHCORE
871     if (MASTER(cr))
872     {
873         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
874     }
875 #endif
876
877     /* NMR restraints must be initialized before load_checkpoint,
878      * since with time averaging the history is added to t_state.
879      * For proper consistency check we therefore need to extend
880      * t_state here.
881      * So the PME-only nodes (if present) will also initialize
882      * the distance restraints.
883      */
884     snew(fcd, 1);
885
886     /* This needs to be called before read_checkpoint to extend the state */
887     init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
888
889     init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
890                 state);
891
892     if (DEFORM(*inputrec))
893     {
894         /* Store the deform reference box before reading the checkpoint */
895         if (SIMMASTER(cr))
896         {
897             copy_mat(state->box, box);
898         }
899         if (PAR(cr))
900         {
901             gmx_bcast(sizeof(box), box, cr);
902         }
903         /* Because we do not have the update struct available yet
904          * in which the reference values should be stored,
905          * we store them temporarily in static variables.
906          * This should be thread safe, since they are only written once
907          * and with identical values.
908          */
909         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
910         deform_init_init_step_tpx = inputrec->init_step;
911         copy_mat(box, deform_init_box_tpx);
912         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
913     }
914
915     if (opt2bSet("-cpi", nfile, fnm))
916     {
917         /* Check if checkpoint file exists before doing continuation.
918          * This way we can use identical input options for the first and subsequent runs...
919          */
920         if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
921         {
922             load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
923                             cr, ddxyz,
924                             inputrec, state, &bReadEkin,
925                             (Flags & MD_APPENDFILES),
926                             (Flags & MD_APPENDFILESSET));
927
928             if (bReadEkin)
929             {
930                 Flags |= MD_READ_EKIN;
931             }
932         }
933     }
934
935     if (MASTER(cr) && (Flags & MD_APPENDFILES))
936     {
937         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
938                      Flags, &fplog);
939     }
940
941     /* override nsteps with value from cmdline */
942     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
943
944     if (SIMMASTER(cr))
945     {
946         copy_mat(state->box, box);
947     }
948
949     if (PAR(cr))
950     {
951         gmx_bcast(sizeof(box), box, cr);
952     }
953
954     /* Essential dynamics */
955     if (opt2bSet("-ei", nfile, fnm))
956     {
957         /* Open input and output files, allocate space for ED data structure */
958         ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
959     }
960
961     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
962                      inputrec->eI == eiNM))
963     {
964         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
965                                            dddlb_opt, dlb_scale,
966                                            ddcsx, ddcsy, ddcsz,
967                                            mtop, inputrec,
968                                            box, state->x,
969                                            &ddbox, &npme_major, &npme_minor);
970
971         make_dd_communicators(fplog, cr, dd_node_order);
972
973         /* Set overallocation to avoid frequent reallocation of arrays */
974         set_over_alloc_dd(TRUE);
975     }
976     else
977     {
978         /* PME, if used, is done on all nodes with 1D decomposition */
979         cr->npmenodes = 0;
980         cr->duty      = (DUTY_PP | DUTY_PME);
981         npme_major    = 1;
982         npme_minor    = 1;
983
984         if (inputrec->ePBC == epbcSCREW)
985         {
986             gmx_fatal(FARGS,
987                       "pbc=%s is only implemented with domain decomposition",
988                       epbc_names[inputrec->ePBC]);
989         }
990     }
991
992     if (PAR(cr))
993     {
994         /* After possible communicator splitting in make_dd_communicators.
995          * we can set up the intra/inter node communication.
996          */
997         gmx_setup_nodecomm(fplog, cr);
998     }
999
1000     /* Initialize per-physical-node MPI process/thread ID and counters. */
1001     gmx_init_intranode_counters(cr);
1002 #ifdef GMX_MPI
1003     if (MULTISIM(cr))
1004     {
1005         md_print_info(cr, fplog,
1006                       "This is simulation %d out of %d running as a composite GROMACS\n"
1007                       "multi-simulation job. Setup for this simulation:\n\n",
1008                       cr->ms->sim, cr->ms->nsim);
1009     }
1010     md_print_info(cr, fplog, "Using %d MPI %s\n",
1011                   cr->nnodes,
1012 #ifdef GMX_THREAD_MPI
1013                   cr->nnodes == 1 ? "thread" : "threads"
1014 #else
1015                   cr->nnodes == 1 ? "process" : "processes"
1016 #endif
1017                   );
1018     fflush(stderr);
1019 #endif
1020
1021     /* Check and update hw_opt for the cut-off scheme */
1022     check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1023
1024     /* Check and update hw_opt for the number of MPI ranks */
1025     check_and_update_hw_opt_3(hw_opt);
1026
1027     gmx_omp_nthreads_init(fplog, cr,
1028                           hwinfo->nthreads_hw_avail,
1029                           hw_opt->nthreads_omp,
1030                           hw_opt->nthreads_omp_pme,
1031                           (cr->duty & DUTY_PP) == 0,
1032                           inputrec->cutoff_scheme == ecutsVERLET);
1033
1034 #ifndef NDEBUG
1035     if (integrator[inputrec->eI].func != do_tpi &&
1036         inputrec->cutoff_scheme == ecutsVERLET)
1037     {
1038         gmx_feenableexcept();
1039     }
1040 #endif
1041
1042     if (bUseGPU)
1043     {
1044         /* Select GPU id's to use */
1045         gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1046                            &hw_opt->gpu_opt);
1047     }
1048     else
1049     {
1050         /* Ignore (potentially) manually selected GPUs */
1051         hw_opt->gpu_opt.n_dev_use = 0;
1052     }
1053
1054     /* check consistency across ranks of things like SIMD
1055      * support and number of GPUs selected */
1056     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1057
1058     /* Now that we know the setup is consistent, check for efficiency */
1059     check_resource_division_efficiency(hwinfo, hw_opt, Flags & MD_NTOMPSET,
1060                                        cr, fplog);
1061
1062     if (DOMAINDECOMP(cr))
1063     {
1064         /* When we share GPUs over ranks, we need to know this for the DLB */
1065         dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1066     }
1067
1068     /* getting number of PP/PME threads
1069        PME: env variable should be read only on one node to make sure it is
1070        identical everywhere;
1071      */
1072     /* TODO nthreads_pp is only used for pinning threads.
1073      * This is a temporary solution until we have a hw topology library.
1074      */
1075     nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
1076     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1077
1078     wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1079
1080     if (PAR(cr))
1081     {
1082         /* Master synchronizes its value of reset_counters with all nodes
1083          * including PME only nodes */
1084         reset_counters = wcycle_get_reset_counters(wcycle);
1085         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1086         wcycle_set_reset_counters(wcycle, reset_counters);
1087     }
1088
1089     snew(nrnb, 1);
1090     if (cr->duty & DUTY_PP)
1091     {
1092         bcast_state(cr, state);
1093
1094         /* Initiate forcerecord */
1095         fr          = mk_forcerec();
1096         fr->hwinfo  = hwinfo;
1097         fr->gpu_opt = &hw_opt->gpu_opt;
1098         init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1099                       opt2fn("-table", nfile, fnm),
1100                       opt2fn("-tabletf", nfile, fnm),
1101                       opt2fn("-tablep", nfile, fnm),
1102                       opt2fn("-tableb", nfile, fnm),
1103                       nbpu_opt,
1104                       FALSE,
1105                       pforce);
1106
1107         /* version for PCA_NOT_READ_NODE (see md.c) */
1108         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1109            "nofile","nofile","nofile","nofile",FALSE,pforce);
1110          */
1111
1112         /* Initialize QM-MM */
1113         if (fr->bQMMM)
1114         {
1115             init_QMMMrec(cr, mtop, inputrec, fr);
1116         }
1117
1118         /* Initialize the mdatoms structure.
1119          * mdatoms is not filled with atom data,
1120          * as this can not be done now with domain decomposition.
1121          */
1122         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1123
1124         /* Initialize the virtual site communication */
1125         vsite = init_vsite(mtop, cr, FALSE);
1126
1127         calc_shifts(box, fr->shift_vec);
1128
1129         /* With periodic molecules the charge groups should be whole at start up
1130          * and the virtual sites should not be far from their proper positions.
1131          */
1132         if (!inputrec->bContinuation && MASTER(cr) &&
1133             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1134         {
1135             /* Make molecules whole at start of run */
1136             if (fr->ePBC != epbcNONE)
1137             {
1138                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1139             }
1140             if (vsite)
1141             {
1142                 /* Correct initial vsite positions are required
1143                  * for the initial distribution in the domain decomposition
1144                  * and for the initial shell prediction.
1145                  */
1146                 construct_vsites_mtop(vsite, mtop, state->x);
1147             }
1148         }
1149
1150         if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1151         {
1152             ewaldcoeff_q  = fr->ewaldcoeff_q;
1153             ewaldcoeff_lj = fr->ewaldcoeff_lj;
1154             pmedata       = &fr->pmedata;
1155         }
1156         else
1157         {
1158             pmedata = NULL;
1159         }
1160     }
1161     else
1162     {
1163         /* This is a PME only node */
1164
1165         /* We don't need the state */
1166         done_state(state);
1167
1168         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1169         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1170         snew(pmedata, 1);
1171     }
1172
1173     if (hw_opt->thread_affinity != threadaffOFF)
1174     {
1175         /* Before setting affinity, check whether the affinity has changed
1176          * - which indicates that probably the OpenMP library has changed it
1177          * since we first checked).
1178          */
1179         gmx_check_thread_affinity_set(fplog, cr,
1180                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1181
1182         /* Set the CPU affinity */
1183         gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1184     }
1185
1186     /* Initiate PME if necessary,
1187      * either on all nodes or on dedicated PME nodes only. */
1188     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1189     {
1190         if (mdatoms)
1191         {
1192             nChargePerturbed = mdatoms->nChargePerturbed;
1193             if (EVDW_PME(inputrec->vdwtype))
1194             {
1195                 nTypePerturbed   = mdatoms->nTypePerturbed;
1196             }
1197         }
1198         if (cr->npmenodes > 0)
1199         {
1200             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1201             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1202             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1203         }
1204
1205         if (cr->duty & DUTY_PME)
1206         {
1207             status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1208                                   mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1209                                   (Flags & MD_REPRODUCIBLE), nthreads_pme);
1210             if (status != 0)
1211             {
1212                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1213             }
1214         }
1215     }
1216
1217
1218     if (integrator[inputrec->eI].func == do_md)
1219     {
1220         /* Turn on signal handling on all nodes */
1221         /*
1222          * (A user signal from the PME nodes (if any)
1223          * is communicated to the PP nodes.
1224          */
1225         signal_handler_install();
1226     }
1227
1228     if (cr->duty & DUTY_PP)
1229     {
1230         /* Assumes uniform use of the number of OpenMP threads */
1231         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1232
1233         if (inputrec->bPull)
1234         {
1235             /* Initialize pull code */
1236             inputrec->pull_work =
1237                 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1238                           mtop, cr, oenv, inputrec->fepvals->init_lambda,
1239                           EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1240         }
1241
1242         if (inputrec->bRot)
1243         {
1244             /* Initialize enforced rotation code */
1245             init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1246                      bVerbose, Flags);
1247         }
1248
1249         if (inputrec->eSwapCoords != eswapNO)
1250         {
1251             /* Initialize ion swapping code */
1252             init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1253                             mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1254         }
1255
1256         constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1257
1258         if (DOMAINDECOMP(cr))
1259         {
1260             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1261             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1262                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1263
1264             set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1265
1266             setup_dd_grid(fplog, cr->dd);
1267         }
1268
1269         /* Now do whatever the user wants us to do (how flexible...) */
1270         integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1271                                       oenv, bVerbose, bCompact,
1272                                       nstglobalcomm,
1273                                       vsite, constr,
1274                                       nstepout, inputrec, mtop,
1275                                       fcd, state,
1276                                       mdatoms, nrnb, wcycle, ed, fr,
1277                                       repl_ex_nst, repl_ex_nex, repl_ex_seed,
1278                                       membed,
1279                                       cpt_period, max_hours,
1280                                       imdport,
1281                                       Flags,
1282                                       walltime_accounting);
1283
1284         if (inputrec->bPull)
1285         {
1286             finish_pull(inputrec->pull_work);
1287         }
1288
1289         if (inputrec->bRot)
1290         {
1291             finish_rot(inputrec->rot);
1292         }
1293
1294     }
1295     else
1296     {
1297         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1298         /* do PME only */
1299         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1300         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1301     }
1302
1303     wallcycle_stop(wcycle, ewcRUN);
1304
1305     /* Finish up, write some stuff
1306      * if rerunMD, don't write last frame again
1307      */
1308     finish_run(fplog, cr,
1309                inputrec, nrnb, wcycle, walltime_accounting,
1310                fr ? fr->nbv : NULL,
1311                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1312
1313
1314     /* Free GPU memory and context */
1315     free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1316
1317     if (opt2bSet("-membed", nfile, fnm))
1318     {
1319         sfree(membed);
1320     }
1321
1322     gmx_hardware_info_free(hwinfo);
1323
1324     /* Does what it says */
1325     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1326     walltime_accounting_destroy(walltime_accounting);
1327
1328     /* Close logfile already here if we were appending to it */
1329     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1330     {
1331         gmx_log_close(fplog);
1332     }
1333
1334     rc = (int)gmx_get_stop_condition();
1335
1336     done_ed(&ed);
1337
1338 #ifdef GMX_THREAD_MPI
1339     /* we need to join all threads. The sub-threads join when they
1340        exit this function, but the master thread needs to be told to
1341        wait for that. */
1342     if (PAR(cr) && MASTER(cr))
1343     {
1344         tMPI_Finalize();
1345     }
1346 #endif
1347
1348     return rc;
1349 }