Converted runner.c to C++
[alexxy/gromacs.git] / src / programs / mdrun / runner.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <algorithm>
43
44 #include <assert.h>
45 #include <signal.h>
46 #include <stdlib.h>
47 #include <string.h>
48 #ifdef HAVE_UNISTD_H
49 #include <unistd.h>
50 #endif
51
52 #include "typedefs.h"
53 #include "copyrite.h"
54 #include "force.h"
55 #include "mdrun.h"
56 #include "md_logging.h"
57 #include "md_support.h"
58 #include "network.h"
59 #include "names.h"
60 #include "disre.h"
61 #include "orires.h"
62 #include "pme.h"
63 #include "mdatoms.h"
64 #include "repl_ex.h"
65 #include "deform.h"
66 #include "qmmm.h"
67 #include "domdec.h"
68 #include "coulomb.h"
69 #include "constr.h"
70 #include "mvdata.h"
71 #include "checkpoint.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "sighandler.h"
74 #include "txtdump.h"
75 #include "gmx_detect_hardware.h"
76 #include "gmx_omp_nthreads.h"
77 #include "gromacs/gmxpreprocess/calc_verletbuf.h"
78 #include "membed.h"
79 #include "gmx_thread_affinity.h"
80 #include "inputrec.h"
81 #include "main.h"
82
83 #include "gromacs/essentialdynamics/edsam.h"
84 #include "gromacs/fileio/tpxio.h"
85 #include "gromacs/math/vec.h"
86 #include "gromacs/mdlib/nbnxn_search.h"
87 #include "gromacs/mdlib/nbnxn_consts.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/pulling/pull_rotation.h"
91 #include "gromacs/swap/swapcoords.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/gmxmpi.h"
95 #include "gromacs/utility/smalloc.h"
96
97 #ifdef GMX_FAHCORE
98 #include "corewrap.h"
99 #endif
100
101 #include "gpu_utils.h"
102 #include "nbnxn_cuda_data_mgmt.h"
103
104 typedef struct {
105     gmx_integrator_t *func;
106 } gmx_intp_t;
107
108 /* The array should match the eI array in include/types/enums.h */
109 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}};
110
111 gmx_int64_t         deform_init_init_step_tpx;
112 matrix              deform_init_box_tpx;
113 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
114
115
116 #ifdef GMX_THREAD_MPI
117 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
118  * the number of threads will get lowered.
119  */
120 #define MIN_ATOMS_PER_MPI_THREAD    90
121 #define MIN_ATOMS_PER_GPU           900
122
123 struct mdrunner_arglist
124 {
125     gmx_hw_opt_t    hw_opt;
126     FILE           *fplog;
127     t_commrec      *cr;
128     int             nfile;
129     const t_filenm *fnm;
130     output_env_t    oenv;
131     gmx_bool        bVerbose;
132     gmx_bool        bCompact;
133     int             nstglobalcomm;
134     ivec            ddxyz;
135     int             dd_node_order;
136     real            rdd;
137     real            rconstr;
138     const char     *dddlb_opt;
139     real            dlb_scale;
140     const char     *ddcsx;
141     const char     *ddcsy;
142     const char     *ddcsz;
143     const char     *nbpu_opt;
144     int             nstlist_cmdline;
145     gmx_int64_t     nsteps_cmdline;
146     int             nstepout;
147     int             resetstep;
148     int             nmultisim;
149     int             repl_ex_nst;
150     int             repl_ex_nex;
151     int             repl_ex_seed;
152     real            pforce;
153     real            cpt_period;
154     real            max_hours;
155     const char     *deviceOptions;
156     int             imdport;
157     unsigned long   Flags;
158 };
159
160
161 /* The function used for spawning threads. Extracts the mdrunner()
162    arguments from its one argument and calls mdrunner(), after making
163    a commrec. */
164 static void mdrunner_start_fn(void *arg)
165 {
166     struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
167     struct mdrunner_arglist  mc  = *mda; /* copy the arg list to make sure
168                                             that it's thread-local. This doesn't
169                                             copy pointed-to items, of course,
170                                             but those are all const. */
171     t_commrec *cr;                       /* we need a local version of this */
172     FILE      *fplog = NULL;
173     t_filenm  *fnm;
174
175     fnm = dup_tfn(mc.nfile, mc.fnm);
176
177     cr = reinitialize_commrec_for_this_thread(mc.cr);
178
179     if (MASTER(cr))
180     {
181         fplog = mc.fplog;
182     }
183
184     mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
185              mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
186              mc.ddxyz, mc.dd_node_order, mc.rdd,
187              mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
188              mc.ddcsx, mc.ddcsy, mc.ddcsz,
189              mc.nbpu_opt, mc.nstlist_cmdline,
190              mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
191              mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
192              mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.imdport, mc.Flags);
193 }
194
195 /* called by mdrunner() to start a specific number of threads (including
196    the main thread) for thread-parallel runs. This in turn calls mdrunner()
197    for each thread.
198    All options besides nthreads are the same as for mdrunner(). */
199 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
200                                          FILE *fplog, t_commrec *cr, int nfile,
201                                          const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
202                                          gmx_bool bCompact, int nstglobalcomm,
203                                          ivec ddxyz, int dd_node_order, real rdd, real rconstr,
204                                          const char *dddlb_opt, real dlb_scale,
205                                          const char *ddcsx, const char *ddcsy, const char *ddcsz,
206                                          const char *nbpu_opt, int nstlist_cmdline,
207                                          gmx_int64_t nsteps_cmdline,
208                                          int nstepout, int resetstep,
209                                          int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
210                                          real pforce, real cpt_period, real max_hours,
211                                          const char *deviceOptions, unsigned long Flags)
212 {
213     int                      ret;
214     struct mdrunner_arglist *mda;
215     t_commrec               *crn; /* the new commrec */
216     t_filenm                *fnmn;
217
218     /* first check whether we even need to start tMPI */
219     if (hw_opt->nthreads_tmpi < 2)
220     {
221         return cr;
222     }
223
224     /* a few small, one-time, almost unavoidable memory leaks: */
225     snew(mda, 1);
226     fnmn = dup_tfn(nfile, fnm);
227
228     /* fill the data structure to pass as void pointer to thread start fn */
229     /* hw_opt contains pointers, which should all be NULL at this stage */
230     mda->hw_opt          = *hw_opt;
231     mda->fplog           = fplog;
232     mda->cr              = cr;
233     mda->nfile           = nfile;
234     mda->fnm             = fnmn;
235     mda->oenv            = oenv;
236     mda->bVerbose        = bVerbose;
237     mda->bCompact        = bCompact;
238     mda->nstglobalcomm   = nstglobalcomm;
239     mda->ddxyz[XX]       = ddxyz[XX];
240     mda->ddxyz[YY]       = ddxyz[YY];
241     mda->ddxyz[ZZ]       = ddxyz[ZZ];
242     mda->dd_node_order   = dd_node_order;
243     mda->rdd             = rdd;
244     mda->rconstr         = rconstr;
245     mda->dddlb_opt       = dddlb_opt;
246     mda->dlb_scale       = dlb_scale;
247     mda->ddcsx           = ddcsx;
248     mda->ddcsy           = ddcsy;
249     mda->ddcsz           = ddcsz;
250     mda->nbpu_opt        = nbpu_opt;
251     mda->nstlist_cmdline = nstlist_cmdline;
252     mda->nsteps_cmdline  = nsteps_cmdline;
253     mda->nstepout        = nstepout;
254     mda->resetstep       = resetstep;
255     mda->nmultisim       = nmultisim;
256     mda->repl_ex_nst     = repl_ex_nst;
257     mda->repl_ex_nex     = repl_ex_nex;
258     mda->repl_ex_seed    = repl_ex_seed;
259     mda->pforce          = pforce;
260     mda->cpt_period      = cpt_period;
261     mda->max_hours       = max_hours;
262     mda->deviceOptions   = deviceOptions;
263     mda->Flags           = Flags;
264
265     /* now spawn new threads that start mdrunner_start_fn(), while
266        the main thread returns, we set thread affinity later */
267     ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
268                        mdrunner_start_fn, (void*)(mda) );
269     if (ret != TMPI_SUCCESS)
270     {
271         return NULL;
272     }
273
274     crn = reinitialize_commrec_for_this_thread(cr);
275     return crn;
276 }
277
278
279 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
280                                         const gmx_hw_opt_t  *hw_opt,
281                                         int                  nthreads_tot,
282                                         int                  ngpu)
283 {
284     int nthreads_tmpi;
285
286     /* There are no separate PME nodes here, as we ensured in
287      * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
288      * and a conditional ensures we would not have ended up here.
289      * Note that separate PME nodes might be switched on later.
290      */
291     if (ngpu > 0)
292     {
293         nthreads_tmpi = ngpu;
294         if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
295         {
296             nthreads_tmpi = nthreads_tot;
297         }
298     }
299     else if (hw_opt->nthreads_omp > 0)
300     {
301         /* Here we could oversubscribe, when we do, we issue a warning later */
302         nthreads_tmpi = std::max(1, nthreads_tot/hw_opt->nthreads_omp);
303     }
304     else
305     {
306         /* TODO choose nthreads_omp based on hardware topology
307            when we have a hardware topology detection library */
308         /* In general, when running up to 4 threads, OpenMP should be faster.
309          * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
310          * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
311          * even on two CPUs it's usually faster (but with many OpenMP threads
312          * it could be faster not to use HT, currently we always use HT).
313          * On Nehalem/Westmere we want to avoid running 16 threads over
314          * two CPUs with HT, so we need a limit<16; thus we use 12.
315          * A reasonable limit for Intel Sandy and Ivy bridge,
316          * not knowing the topology, is 16 threads.
317          */
318         const int nthreads_omp_always_faster             =  4;
319         const int nthreads_omp_always_faster_Nehalem     = 12;
320         const int nthreads_omp_always_faster_SandyBridge = 16;
321         gmx_bool  bIntel_Family6;
322
323         bIntel_Family6 =
324             (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
325              gmx_cpuid_family(hwinfo->cpuid_info) == 6);
326
327         if (nthreads_tot <= nthreads_omp_always_faster ||
328             (bIntel_Family6 &&
329              ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
330               (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
331         {
332             /* Use pure OpenMP parallelization */
333             nthreads_tmpi = 1;
334         }
335         else
336         {
337             /* Don't use OpenMP parallelization */
338             nthreads_tmpi = nthreads_tot;
339         }
340     }
341
342     return nthreads_tmpi;
343 }
344
345
346 /* Get the number of threads to use for thread-MPI based on how many
347  * were requested, which algorithms we're using,
348  * and how many particles there are.
349  * At the point we have already called check_and_update_hw_opt.
350  * Thus all options should be internally consistent and consistent
351  * with the hardware, except that ntmpi could be larger than #GPU.
352  */
353 static int get_nthreads_mpi(const gmx_hw_info_t *hwinfo,
354                             gmx_hw_opt_t *hw_opt,
355                             t_inputrec *inputrec, gmx_mtop_t *mtop,
356                             const t_commrec *cr,
357                             FILE *fplog)
358 {
359     int      nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
360     int      min_atoms_per_mpi_thread;
361     gmx_bool bCanUseGPU;
362
363     if (hw_opt->nthreads_tmpi > 0)
364     {
365         /* Trivial, return right away */
366         return hw_opt->nthreads_tmpi;
367     }
368
369     nthreads_hw = hwinfo->nthreads_hw_avail;
370
371     /* How many total (#tMPI*#OpenMP) threads can we start? */
372     if (hw_opt->nthreads_tot > 0)
373     {
374         nthreads_tot_max = hw_opt->nthreads_tot;
375     }
376     else
377     {
378         nthreads_tot_max = nthreads_hw;
379     }
380
381     bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET &&
382                   hwinfo->gpu_info.ncuda_dev_compatible > 0);
383     if (bCanUseGPU)
384     {
385         ngpu = hwinfo->gpu_info.ncuda_dev_compatible;
386     }
387     else
388     {
389         ngpu = 0;
390     }
391
392     if (inputrec->cutoff_scheme == ecutsGROUP)
393     {
394         /* We checked this before, but it doesn't hurt to do it once more */
395         assert(hw_opt->nthreads_omp == 1);
396     }
397
398     nthreads_tmpi =
399         get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
400
401     if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
402     {
403         /* Dims/steps are divided over the nodes iso splitting the atoms */
404         min_atoms_per_mpi_thread = 0;
405     }
406     else
407     {
408         if (bCanUseGPU)
409         {
410             min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
411         }
412         else
413         {
414             min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
415         }
416     }
417
418     /* Check if an algorithm does not support parallel simulation.  */
419     if (nthreads_tmpi != 1 &&
420         ( inputrec->eI == eiLBFGS ||
421           inputrec->coulombtype == eelEWALD ) )
422     {
423         nthreads_tmpi = 1;
424
425         md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
426         if (hw_opt->nthreads_tmpi > nthreads_tmpi)
427         {
428             gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
429         }
430     }
431     else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
432     {
433         /* the thread number was chosen automatically, but there are too many
434            threads (too few atoms per thread) */
435         nthreads_new = std::max(1, mtop->natoms/min_atoms_per_mpi_thread);
436
437         /* Avoid partial use of Hyper-Threading */
438         if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
439             nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
440         {
441             nthreads_new = nthreads_hw/2;
442         }
443
444         /* Avoid large prime numbers in the thread count */
445         if (nthreads_new >= 6)
446         {
447             /* Use only 6,8,10 with additional factors of 2 */
448             int fac;
449
450             fac = 2;
451             while (3*fac*2 <= nthreads_new)
452             {
453                 fac *= 2;
454             }
455
456             nthreads_new = (nthreads_new/fac)*fac;
457         }
458         else
459         {
460             /* Avoid 5 */
461             if (nthreads_new == 5)
462             {
463                 nthreads_new = 4;
464             }
465         }
466
467         nthreads_tmpi = nthreads_new;
468
469         fprintf(stderr, "\n");
470         fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
471         fprintf(stderr, "      only starting %d thread-MPI threads.\n", nthreads_tmpi);
472         fprintf(stderr, "      You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
473     }
474
475     return nthreads_tmpi;
476 }
477 #endif /* GMX_THREAD_MPI */
478
479
480 /* We determine the extra cost of the non-bonded kernels compared to
481  * a reference nstlist value of 10 (which is the default in grompp).
482  */
483 static const int    nbnxnReferenceNstlist = 10;
484 /* The values to try when switching  */
485 const int           nstlist_try[] = { 20, 25, 40 };
486 #define NNSTL  sizeof(nstlist_try)/sizeof(nstlist_try[0])
487 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
488  * but never more than listfac_max.
489  * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
490  * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
491  * Note that both CPU and GPU factors are conservative. Performance should
492  * not go down due to this tuning, except with a relatively slow GPU.
493  * On the other hand, at medium/high parallelization or with fast GPUs
494  * nstlist will not be increased enough to reach optimal performance.
495  */
496 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
497 static const float  nbnxn_cpu_listfac_ok    = 1.05;
498 static const float  nbnxn_cpu_listfac_max   = 1.09;
499 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
500 static const float  nbnxn_gpu_listfac_ok    = 1.20;
501 static const float  nbnxn_gpu_listfac_max   = 1.30;
502
503 /* Try to increase nstlist when using the Verlet cut-off scheme */
504 static void increase_nstlist(FILE *fp, t_commrec *cr,
505                              t_inputrec *ir, int nstlist_cmdline,
506                              const gmx_mtop_t *mtop, matrix box,
507                              gmx_bool bGPU)
508 {
509     float                  listfac_ok, listfac_max;
510     int                    nstlist_orig, nstlist_prev;
511     verletbuf_list_setup_t ls;
512     real                   rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
513     real                   rlist_new, rlist_prev;
514     size_t                 nstlist_ind = 0;
515     t_state                state_tmp;
516     gmx_bool               bBox, bDD, bCont;
517     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";
518     const char            *nve_err  = "Can not increase nstlist because an NVE ensemble is used";
519     const char            *vbd_err  = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
520     const char            *box_err  = "Can not increase nstlist because the box is too small";
521     const char            *dd_err   = "Can not increase nstlist because of domain decomposition limitations";
522     char                   buf[STRLEN];
523     const float            oneThird = 1.0f / 3.0f;
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 the reference value for nstlist (10).
610      */
611     nstlist_prev = ir->nstlist;
612     ir->nstlist  = nbnxnReferenceNstlist;
613     calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
614                             &rlistWithReferenceNstlist);
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  = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
621     rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - 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     const 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*std::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 != NULL && 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     t_inputrec               *inputrec;
1098     t_state                  *state = NULL;
1099     matrix                    box;
1100     gmx_ddbox_t               ddbox = {0};
1101     int                       npme_major, npme_minor;
1102     t_nrnb                   *nrnb;
1103     gmx_mtop_t               *mtop          = NULL;
1104     t_mdatoms                *mdatoms       = NULL;
1105     t_forcerec               *fr            = NULL;
1106     t_fcdata                 *fcd           = NULL;
1107     real                      ewaldcoeff_q  = 0;
1108     real                      ewaldcoeff_lj = 0;
1109     gmx_pme_t                *pmedata       = NULL;
1110     gmx_vsite_t              *vsite         = NULL;
1111     gmx_constr_t              constr;
1112     int                       nChargePerturbed = -1, nTypePerturbed = 0, status;
1113     gmx_wallcycle_t           wcycle;
1114     gmx_bool                  bReadEkin;
1115     gmx_walltime_accounting_t walltime_accounting = NULL;
1116     int                       rc;
1117     gmx_int64_t               reset_counters;
1118     gmx_edsam_t               ed           = NULL;
1119     int                       nthreads_pme = 1;
1120     int                       nthreads_pp  = 1;
1121     gmx_membed_t              membed       = NULL;
1122     gmx_hw_info_t            *hwinfo       = NULL;
1123     /* The master rank decides early on bUseGPU and broadcasts this later */
1124     gmx_bool                  bUseGPU      = FALSE;
1125
1126     /* CAUTION: threads may be started later on in this function, so
1127        cr doesn't reflect the final parallel state right now */
1128     snew(inputrec, 1);
1129     snew(mtop, 1);
1130
1131     if (Flags & MD_APPENDFILES)
1132     {
1133         fplog = NULL;
1134     }
1135
1136     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1137     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1138
1139     /* Detect hardware, gather information. This is an operation that is
1140      * global for this process (MPI rank). */
1141     hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
1142
1143
1144     snew(state, 1);
1145     if (SIMMASTER(cr))
1146     {
1147         /* Read (nearly) all data required for the simulation */
1148         read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
1149
1150         if (inputrec->cutoff_scheme != ecutsVERLET &&
1151             ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1152         {
1153             convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
1154         }
1155
1156         if (inputrec->cutoff_scheme == ecutsVERLET)
1157         {
1158             /* Here the master rank decides if all ranks will use GPUs */
1159             bUseGPU = (hwinfo->gpu_info.ncuda_dev_compatible > 0 ||
1160                        getenv("GMX_EMULATE_GPU") != NULL);
1161
1162             /* TODO add GPU kernels for this and replace this check by:
1163              * (bUseGPU && (ir->vdwtype == evdwPME &&
1164              *               ir->ljpme_combination_rule == eljpmeLB))
1165              * update the message text and the content of nbnxn_acceleration_supported.
1166              */
1167             if (bUseGPU &&
1168                 !nbnxn_acceleration_supported(fplog, cr, inputrec, bUseGPU))
1169             {
1170                 /* Fallback message printed by nbnxn_acceleration_supported */
1171                 if (bForceUseGPU)
1172                 {
1173                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
1174                 }
1175                 bUseGPU = FALSE;
1176             }
1177
1178             prepare_verlet_scheme(fplog, cr,
1179                                   inputrec, nstlist_cmdline, mtop, state->box,
1180                                   bUseGPU);
1181         }
1182         else
1183         {
1184             if (nstlist_cmdline > 0)
1185             {
1186                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
1187             }
1188
1189             if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
1190             {
1191                 md_print_warn(cr, fplog,
1192                               "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1193                               "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1194                               "      (for quick performance testing you can use the -testverlet option)\n");
1195             }
1196
1197             if (bForceUseGPU)
1198             {
1199                 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1200             }
1201
1202 #ifdef GMX_TARGET_BGQ
1203             md_print_warn(cr, fplog,
1204                           "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
1205                           "      BlueGene/Q. You will observe better performance from using the\n"
1206                           "      Verlet cut-off scheme.\n");
1207 #endif
1208         }
1209     }
1210
1211     /* Check and update the hardware options for internal consistency */
1212     check_and_update_hw_opt_1(hw_opt, SIMMASTER(cr));
1213
1214     /* Early check for externally set process affinity. */
1215     gmx_check_thread_affinity_set(fplog, cr,
1216                                   hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1217     if (SIMMASTER(cr))
1218     {
1219
1220 #ifdef GMX_THREAD_MPI
1221         if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1222         {
1223             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1224         }
1225 #endif
1226
1227         if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1228             cr->npmenodes <= 0)
1229         {
1230             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");
1231         }
1232     }
1233
1234 #ifdef GMX_THREAD_MPI
1235     if (SIMMASTER(cr))
1236     {
1237         /* Since the master knows the cut-off scheme, update hw_opt for this.
1238          * This is done later for normal MPI and also once more with tMPI
1239          * for all tMPI ranks.
1240          */
1241         check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1242
1243         /* NOW the threads will be started: */
1244         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1245                                                  hw_opt,
1246                                                  inputrec, mtop,
1247                                                  cr, fplog);
1248         if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1249         {
1250             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1251         }
1252
1253         if (hw_opt->nthreads_tmpi > 1)
1254         {
1255             t_commrec *cr_old       = cr;
1256             /* now start the threads. */
1257             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1258                                         oenv, bVerbose, bCompact, nstglobalcomm,
1259                                         ddxyz, dd_node_order, rdd, rconstr,
1260                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1261                                         nbpu_opt, nstlist_cmdline,
1262                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
1263                                         repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1264                                         cpt_period, max_hours, deviceOptions,
1265                                         Flags);
1266             /* the main thread continues here with a new cr. We don't deallocate
1267                the old cr because other threads may still be reading it. */
1268             if (cr == NULL)
1269             {
1270                 gmx_comm("Failed to spawn threads");
1271             }
1272         }
1273     }
1274 #endif
1275     /* END OF CAUTION: cr is now reliable */
1276
1277     /* g_membed initialisation *
1278      * Because we change the mtop, init_membed is called before the init_parallel *
1279      * (in case we ever want to make it run in parallel) */
1280     if (opt2bSet("-membed", nfile, fnm))
1281     {
1282         if (MASTER(cr))
1283         {
1284             fprintf(stderr, "Initializing membed");
1285         }
1286         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1287     }
1288
1289     if (PAR(cr))
1290     {
1291         /* now broadcast everything to the non-master nodes/threads: */
1292         init_parallel(cr, inputrec, mtop);
1293     }
1294     if (fplog != NULL)
1295     {
1296         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1297     }
1298
1299     /* now make sure the state is initialized and propagated */
1300     set_state_entries(state, inputrec);
1301
1302     /* A parallel command line option consistency check that we can
1303        only do after any threads have started. */
1304     if (!PAR(cr) &&
1305         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1306     {
1307         gmx_fatal(FARGS,
1308                   "The -dd or -npme option request a parallel simulation, "
1309 #ifndef GMX_MPI
1310                   "but %s was compiled without threads or MPI enabled"
1311 #else
1312 #ifdef GMX_THREAD_MPI
1313                   "but the number of threads (option -nt) is 1"
1314 #else
1315                   "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1316 #endif
1317 #endif
1318                   , ShortProgram()
1319                   );
1320     }
1321
1322     if ((Flags & MD_RERUN) &&
1323         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1324     {
1325         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1326     }
1327
1328     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
1329     {
1330         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
1331     }
1332
1333     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1334     {
1335         if (cr->npmenodes > 0)
1336         {
1337             gmx_fatal_collective(FARGS, cr, NULL,
1338                                  "PME nodes are requested, but the system does not use PME electrostatics or LJ-PME");
1339         }
1340
1341         cr->npmenodes = 0;
1342     }
1343
1344 #ifdef GMX_FAHCORE
1345     if (MASTER(cr))
1346     {
1347         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1348     }
1349 #endif
1350
1351     /* NMR restraints must be initialized before load_checkpoint,
1352      * since with time averaging the history is added to t_state.
1353      * For proper consistency check we therefore need to extend
1354      * t_state here.
1355      * So the PME-only nodes (if present) will also initialize
1356      * the distance restraints.
1357      */
1358     snew(fcd, 1);
1359
1360     /* This needs to be called before read_checkpoint to extend the state */
1361     init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
1362
1363     init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
1364                 state);
1365
1366     if (DEFORM(*inputrec))
1367     {
1368         /* Store the deform reference box before reading the checkpoint */
1369         if (SIMMASTER(cr))
1370         {
1371             copy_mat(state->box, box);
1372         }
1373         if (PAR(cr))
1374         {
1375             gmx_bcast(sizeof(box), box, cr);
1376         }
1377         /* Because we do not have the update struct available yet
1378          * in which the reference values should be stored,
1379          * we store them temporarily in static variables.
1380          * This should be thread safe, since they are only written once
1381          * and with identical values.
1382          */
1383         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1384         deform_init_init_step_tpx = inputrec->init_step;
1385         copy_mat(box, deform_init_box_tpx);
1386         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1387     }
1388
1389     if (opt2bSet("-cpi", nfile, fnm))
1390     {
1391         /* Check if checkpoint file exists before doing continuation.
1392          * This way we can use identical input options for the first and subsequent runs...
1393          */
1394         if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1395         {
1396             load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1397                             cr, ddxyz,
1398                             inputrec, state, &bReadEkin,
1399                             (Flags & MD_APPENDFILES),
1400                             (Flags & MD_APPENDFILESSET));
1401
1402             if (bReadEkin)
1403             {
1404                 Flags |= MD_READ_EKIN;
1405             }
1406         }
1407     }
1408
1409     if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1410 #ifdef GMX_THREAD_MPI
1411         /* With thread MPI only the master node/thread exists in mdrun.c,
1412          * therefore non-master nodes need to open the "seppot" log file here.
1413          */
1414         || (!MASTER(cr) && (Flags & MD_SEPPOT))
1415 #endif
1416         )
1417     {
1418         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1419                      Flags, &fplog);
1420     }
1421
1422     /* override nsteps with value from cmdline */
1423     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1424
1425     if (SIMMASTER(cr))
1426     {
1427         copy_mat(state->box, box);
1428     }
1429
1430     if (PAR(cr))
1431     {
1432         gmx_bcast(sizeof(box), box, cr);
1433     }
1434
1435     /* Essential dynamics */
1436     if (opt2bSet("-ei", nfile, fnm))
1437     {
1438         /* Open input and output files, allocate space for ED data structure */
1439         ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1440     }
1441
1442     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1443                      inputrec->eI == eiNM))
1444     {
1445         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1446                                            dddlb_opt, dlb_scale,
1447                                            ddcsx, ddcsy, ddcsz,
1448                                            mtop, inputrec,
1449                                            box, state->x,
1450                                            &ddbox, &npme_major, &npme_minor);
1451
1452         make_dd_communicators(fplog, cr, dd_node_order);
1453
1454         /* Set overallocation to avoid frequent reallocation of arrays */
1455         set_over_alloc_dd(TRUE);
1456     }
1457     else
1458     {
1459         /* PME, if used, is done on all nodes with 1D decomposition */
1460         cr->npmenodes = 0;
1461         cr->duty      = (DUTY_PP | DUTY_PME);
1462         npme_major    = 1;
1463         npme_minor    = 1;
1464
1465         if (inputrec->ePBC == epbcSCREW)
1466         {
1467             gmx_fatal(FARGS,
1468                       "pbc=%s is only implemented with domain decomposition",
1469                       epbc_names[inputrec->ePBC]);
1470         }
1471     }
1472
1473     if (PAR(cr))
1474     {
1475         /* After possible communicator splitting in make_dd_communicators.
1476          * we can set up the intra/inter node communication.
1477          */
1478         gmx_setup_nodecomm(fplog, cr);
1479     }
1480
1481     /* Initialize per-physical-node MPI process/thread ID and counters. */
1482     gmx_init_intranode_counters(cr);
1483
1484 #ifdef GMX_MPI
1485     md_print_info(cr, fplog, "Using %d MPI %s\n",
1486                   cr->nnodes,
1487 #ifdef GMX_THREAD_MPI
1488                   cr->nnodes == 1 ? "thread" : "threads"
1489 #else
1490                   cr->nnodes == 1 ? "process" : "processes"
1491 #endif
1492                   );
1493     fflush(stderr);
1494 #endif
1495
1496     /* Check and update hw_opt for the cut-off scheme */
1497     check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1498
1499     gmx_omp_nthreads_init(fplog, cr,
1500                           hwinfo->nthreads_hw_avail,
1501                           hw_opt->nthreads_omp,
1502                           hw_opt->nthreads_omp_pme,
1503                           (cr->duty & DUTY_PP) == 0,
1504                           inputrec->cutoff_scheme == ecutsVERLET);
1505
1506     if (PAR(cr))
1507     {
1508         /* The master rank decided on the use of GPUs,
1509          * broadcast this information to all ranks.
1510          */
1511         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
1512     }
1513
1514     if (bUseGPU)
1515     {
1516         if (cr->npmenodes == -1)
1517         {
1518             /* Don't automatically use PME-only nodes with GPUs */
1519             cr->npmenodes = 0;
1520         }
1521
1522         /* Select GPU id's to use */
1523         gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1524                            &hw_opt->gpu_opt);
1525     }
1526     else
1527     {
1528         /* Ignore (potentially) manually selected GPUs */
1529         hw_opt->gpu_opt.ncuda_dev_use = 0;
1530     }
1531
1532     /* check consistency across ranks of things like SIMD
1533      * support and number of GPUs selected */
1534     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1535
1536     if (DOMAINDECOMP(cr))
1537     {
1538         /* When we share GPUs over ranks, we need to know this for the DLB */
1539         dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1540     }
1541
1542     /* getting number of PP/PME threads
1543        PME: env variable should be read only on one node to make sure it is
1544        identical everywhere;
1545      */
1546     /* TODO nthreads_pp is only used for pinning threads.
1547      * This is a temporary solution until we have a hw topology library.
1548      */
1549     nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
1550     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1551
1552     wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1553
1554     if (PAR(cr))
1555     {
1556         /* Master synchronizes its value of reset_counters with all nodes
1557          * including PME only nodes */
1558         reset_counters = wcycle_get_reset_counters(wcycle);
1559         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1560         wcycle_set_reset_counters(wcycle, reset_counters);
1561     }
1562
1563     snew(nrnb, 1);
1564     if (cr->duty & DUTY_PP)
1565     {
1566         bcast_state(cr, state);
1567
1568         /* Initiate forcerecord */
1569         fr          = mk_forcerec();
1570         fr->hwinfo  = hwinfo;
1571         fr->gpu_opt = &hw_opt->gpu_opt;
1572         init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1573                       opt2fn("-table", nfile, fnm),
1574                       opt2fn("-tabletf", nfile, fnm),
1575                       opt2fn("-tablep", nfile, fnm),
1576                       opt2fn("-tableb", nfile, fnm),
1577                       nbpu_opt,
1578                       FALSE,
1579                       pforce);
1580
1581         /* version for PCA_NOT_READ_NODE (see md.c) */
1582         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1583            "nofile","nofile","nofile","nofile",FALSE,pforce);
1584          */
1585         fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1586
1587         /* Initialize QM-MM */
1588         if (fr->bQMMM)
1589         {
1590             init_QMMMrec(cr, mtop, inputrec, fr);
1591         }
1592
1593         /* Initialize the mdatoms structure.
1594          * mdatoms is not filled with atom data,
1595          * as this can not be done now with domain decomposition.
1596          */
1597         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1598
1599         /* Initialize the virtual site communication */
1600         vsite = init_vsite(mtop, cr, FALSE);
1601
1602         calc_shifts(box, fr->shift_vec);
1603
1604         /* With periodic molecules the charge groups should be whole at start up
1605          * and the virtual sites should not be far from their proper positions.
1606          */
1607         if (!inputrec->bContinuation && MASTER(cr) &&
1608             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1609         {
1610             /* Make molecules whole at start of run */
1611             if (fr->ePBC != epbcNONE)
1612             {
1613                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1614             }
1615             if (vsite)
1616             {
1617                 /* Correct initial vsite positions are required
1618                  * for the initial distribution in the domain decomposition
1619                  * and for the initial shell prediction.
1620                  */
1621                 construct_vsites_mtop(vsite, mtop, state->x);
1622             }
1623         }
1624
1625         if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1626         {
1627             ewaldcoeff_q  = fr->ewaldcoeff_q;
1628             ewaldcoeff_lj = fr->ewaldcoeff_lj;
1629             pmedata       = &fr->pmedata;
1630         }
1631         else
1632         {
1633             pmedata = NULL;
1634         }
1635     }
1636     else
1637     {
1638         /* This is a PME only node */
1639
1640         /* We don't need the state */
1641         done_state(state);
1642
1643         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1644         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1645         snew(pmedata, 1);
1646     }
1647
1648     if (hw_opt->thread_affinity != threadaffOFF)
1649     {
1650         /* Before setting affinity, check whether the affinity has changed
1651          * - which indicates that probably the OpenMP library has changed it
1652          * since we first checked).
1653          */
1654         gmx_check_thread_affinity_set(fplog, cr,
1655                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1656
1657         /* Set the CPU affinity */
1658         gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1659     }
1660
1661     /* Initiate PME if necessary,
1662      * either on all nodes or on dedicated PME nodes only. */
1663     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1664     {
1665         if (mdatoms)
1666         {
1667             nChargePerturbed = mdatoms->nChargePerturbed;
1668             if (EVDW_PME(inputrec->vdwtype))
1669             {
1670                 nTypePerturbed   = mdatoms->nTypePerturbed;
1671             }
1672         }
1673         if (cr->npmenodes > 0)
1674         {
1675             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1676             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1677             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1678         }
1679
1680         if (cr->duty & DUTY_PME)
1681         {
1682             status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1683                                   mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1684                                   (Flags & MD_REPRODUCIBLE), nthreads_pme);
1685             if (status != 0)
1686             {
1687                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1688             }
1689         }
1690     }
1691
1692
1693     if (integrator[inputrec->eI].func == do_md)
1694     {
1695         /* Turn on signal handling on all nodes */
1696         /*
1697          * (A user signal from the PME nodes (if any)
1698          * is communicated to the PP nodes.
1699          */
1700         signal_handler_install();
1701     }
1702
1703     if (cr->duty & DUTY_PP)
1704     {
1705         /* Assumes uniform use of the number of OpenMP threads */
1706         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1707
1708         if (inputrec->ePull != epullNO)
1709         {
1710             /* Initialize pull code */
1711             init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1712                       EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1713         }
1714
1715         if (inputrec->bRot)
1716         {
1717             /* Initialize enforced rotation code */
1718             init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1719                      bVerbose, Flags);
1720         }
1721
1722         if (inputrec->eSwapCoords != eswapNO)
1723         {
1724             /* Initialize ion swapping code */
1725             init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1726                             mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1727         }
1728
1729         constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1730
1731         if (DOMAINDECOMP(cr))
1732         {
1733             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1734             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1735                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1736
1737             set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1738
1739             setup_dd_grid(fplog, cr->dd);
1740         }
1741
1742         /* Now do whatever the user wants us to do (how flexible...) */
1743         integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1744                                       oenv, bVerbose, bCompact,
1745                                       nstglobalcomm,
1746                                       vsite, constr,
1747                                       nstepout, inputrec, mtop,
1748                                       fcd, state,
1749                                       mdatoms, nrnb, wcycle, ed, fr,
1750                                       repl_ex_nst, repl_ex_nex, repl_ex_seed,
1751                                       membed,
1752                                       cpt_period, max_hours,
1753                                       deviceOptions,
1754                                       imdport,
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         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1772         /* do PME only */
1773         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1774         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1775     }
1776
1777     wallcycle_stop(wcycle, ewcRUN);
1778
1779     /* Finish up, write some stuff
1780      * if rerunMD, don't write last frame again
1781      */
1782     finish_run(fplog, cr,
1783                inputrec, nrnb, wcycle, walltime_accounting,
1784                fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1785                nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1786                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1787
1788
1789     /* Free GPU memory and context */
1790     free_gpu_resources(fr, cr);
1791
1792     if (opt2bSet("-membed", nfile, fnm))
1793     {
1794         sfree(membed);
1795     }
1796
1797     gmx_hardware_info_free(hwinfo);
1798
1799     /* Does what it says */
1800     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1801     walltime_accounting_destroy(walltime_accounting);
1802
1803     /* Close logfile already here if we were appending to it */
1804     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1805     {
1806         gmx_log_close(fplog);
1807     }
1808
1809     rc = (int)gmx_get_stop_condition();
1810
1811 #ifdef GMX_THREAD_MPI
1812     /* we need to join all threads. The sub-threads join when they
1813        exit this function, but the master thread needs to be told to
1814        wait for that. */
1815     if (PAR(cr) && MASTER(cr))
1816     {
1817         tMPI_Finalize();
1818     }
1819 #endif
1820
1821     return rc;
1822 }