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