From: Berk Hess Date: Wed, 30 Jul 2014 12:59:17 +0000 (+0200) Subject: Added checks for TPI + Verlet scheme X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=ebfd9140ea5bfd364ba7f586d2316de59df21456;p=alexxy%2Fgromacs.git Added checks for TPI + Verlet scheme Change-Id: I040f4cc9dbb8ff25ed97d18cca0b4a2c1e7014b0 --- diff --git a/src/gromacs/gmxpreprocess/readir.c b/src/gromacs/gmxpreprocess/readir.c index 9e2dbf0311..9d28fe0252 100644 --- a/src/gromacs/gmxpreprocess/readir.c +++ b/src/gromacs/gmxpreprocess/readir.c @@ -575,6 +575,8 @@ void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts, CHECK(ir->nstlist <= 0); sprintf(err_buf, "TPI does not work with full electrostatics other than PME"); CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype)); + sprintf(err_buf, "TPI does not work (yet) with the Verlet cut-off scheme"); + CHECK(ir->cutoff_scheme == ecutsVERLET); } /* SHAKE / LINCS */ diff --git a/src/gromacs/mdlib/tpi.c b/src/gromacs/mdlib/tpi.c index 72a6901f9e..8eaa88b740 100644 --- a/src/gromacs/mdlib/tpi.c +++ b/src/gromacs/mdlib/tpi.c @@ -167,6 +167,11 @@ double do_tpi(FILE *fplog, t_commrec *cr, real bU_bin_limit = 50; real bU_logV_bin_limit = bU_bin_limit + 10; + if (inputrec->cutoff_scheme == ecutsVERLET) + { + gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme"); + } + nnodes = cr->nnodes; top = gmx_mtop_generate_local_top(top_global, inputrec);