Add fatal error for Andersen+constraints+DD
authorMark Abraham <mark.j.abraham@gmail.com>
Mon, 24 Feb 2014 16:31:07 +0000 (17:31 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Sat, 1 Mar 2014 15:41:59 +0000 (16:41 +0100)
This combination produced a temperature that was 6.5 degrees higher
than the same with one domain, for a 1ns PME lysozyme in water
simulation.

Change-Id: I9f80276c47de955a5053bcabb6fe7c9bfdceaf0e

include/update.h
src/kernel/md.c
src/mdlib/update.c

index 93bbf0129edfbe6c2220e9fdb1bc75ffae8b1771..8fde948258b6b0c8aada62a1c3e41cabb19ee23f 100644 (file)
@@ -121,7 +121,7 @@ void update_coords(FILE             *fplog,
 /* Return TRUE if OK, FALSE in case of Shake Error */
 
 GMX_LIBMD_EXPORT
-extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_large_int_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr);
+extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_large_int_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr, gmx_bool bIsDomainDecomposition);
 
 GMX_LIBMD_EXPORT
 void update_constraints(FILE             *fplog,
index f802f55fa995e876b6787857ba34424e47162f09..2f35c6d0b854dd49361458988bfac96cc9a4068e 100644 (file)
@@ -1638,7 +1638,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
             if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
             {
                 gmx_bool bIfRandomize;
-                bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
+                bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr, DOMAINDECOMP(cr));
                 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
                 if (constr && bIfRandomize)
                 {
index bc18d8887f836bbe2375c3ebd913ab7dfef2a748..d1f851fbfcb9a1549c99ae4c85f9fc34a18b150e 100644 (file)
@@ -2083,11 +2083,17 @@ void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[]
             mv[XX], mv[YY], mv[ZZ]);
 }
 
-extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_large_int_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr)
+extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_large_int_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr, gmx_bool bIsDomainDecomposition)
 {
 
     int  i;
     real rate = (ir->delta_t)/ir->opts.tau_t[0];
+
+    if (ir->etc == etcANDERSEN && constr && bIsDomainDecomposition)
+    {
+        gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints and domain decomposition. Please consider the massive Andersen thermostat.");
+    }
+
     /* proceed with andersen if 1) it's fixed probability per
        particle andersen or 2) it's massive andersen and it's tau_t/dt */
     if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))