Add fatal errors for VV and twin-range MTS
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 13 May 2014 15:52:50 +0000 (17:52 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 22 May 2014 07:13:46 +0000 (09:13 +0200)
Michael never implemented the multiple-time stepping with the VV
integrator family and constraints (see code that calls
combine_forces() from update_coords() in src/mdlib/update.c). Probably
that means the multiple-time-step regime was not tested with VV
either. Strictly speaking, these new fatal errors have scope that is
wider than is clearly warranted, but it is not clear the
no-constraints VV path was only ever as bad as the broken leap-frog
path is (see #1400).

I suspect VV+constraints will work with the incoming fixes for
leap-frog, but until someone wants to use it (and why would they?),
then I'm not going to test that it works as well as it does with
leap-frog.

Change-Id: Ib61d0fb7661bca2101c04423a6af1744420c06ab

src/kernel/md.c
src/kernel/readir.c
src/mdlib/md_support.c

index 2f35c6d0b854dd49361458988bfac96cc9a4068e..7091591483648a50ec4bc5882dc721f4c424f520 100644 (file)
@@ -1217,6 +1217,10 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
              * step to combine the long-range forces on these steps.
              * For nstcalclr=1 this is not done, since the forces would have been added
              * directly to the short-range forces already.
+             *
+             * TODO Remove various aspects of VV+twin-range in master
+             * branch, because VV integrators did not ever support
+             * twin-range multiple time stepping with constraints.
              */
             bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
 
index 1474cb4aa4eb2a833fe4d04da5d1e04acd6a69d4..77ebe78672abdc998ec2ed8fff61574dd085036e 100644 (file)
@@ -1240,6 +1240,12 @@ void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
         }
     }
 
+    if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
+    {
+        sprintf(warn_buf, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
+        warning_error(wi, warn_buf);
+    }
+
     /* IMPLICIT SOLVENT */
     if (ir->coulombtype == eelGB_NOTUSED)
     {
index b1b73bab5992810ef85239c30e1a10975fd4f23e..ba800e56b6d24135bd43668edf8adaa253ce5951 100644 (file)
@@ -50,6 +50,7 @@
 #include "nrnb.h"
 #include "md_logging.h"
 #include "md_support.h"
+#include "names.h"
 
 /* Is the signal in one simulation independent of other simulations? */
 gmx_bool gs_simlocal[eglsNR] = { TRUE, FALSE, FALSE, TRUE };
@@ -741,6 +742,11 @@ void check_ir_old_tpx_versions(t_commrec *cr, FILE *fplog,
                             "nstdhdl", &ir->fepvals->nstdhdl);
         }
     }
+
+    if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
+    {
+        gmx_fatal(FARGS, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
+    }
 }
 
 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,