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