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