Fix constraints checks in grompp
authorMark Abraham <mark.j.abraham@gmail.com>
Sat, 29 Nov 2014 21:50:25 +0000 (22:50 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 15 Dec 2014 18:26:08 +0000 (19:26 +0100)
This adds F_CONSTRNC to the set of interactions that should trigger
checks on SHAKE and LINCS parameters, and permits the MTTK check to
test also for SETTLE for its warning. The latter will also be used to
prevent MTTK + any constraints in a subsequent patch.

Change-Id: I175a2f050623192ff35a4cf93605fcfdd909f41a

src/gromacs/gmxpreprocess/grompp.c
src/gromacs/gmxpreprocess/readir.c
src/gromacs/gmxpreprocess/readir.h

index 8a7c681315351031ca616b21b16f926420c07532..e449a07d02efc06d0eddb43598be31cc07cc79cc 100644 (file)
@@ -402,6 +402,8 @@ static void check_shells_inputrec(gmx_mtop_t *mtop,
     }
 }
 
+/* TODO Decide whether this function can be consolidated with
+ * gmx_mtop_ftype_count */
 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
 {
     int nint, mb;
@@ -624,11 +626,21 @@ new_status(const char *topfile, const char *topppfile, const char *confin,
                     nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
             warning(wi, buf);
         }
+
+        /* Do more checks, mostly related to constraints */
         if (bVerbose)
         {
             fprintf(stderr, "double-checking input for internal consistency...\n");
         }
-        double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
+        {
+            int bHasNormalConstraints = 0 < (nint_ftype(sys, molinfo, F_CONSTR) +
+                                             nint_ftype(sys, molinfo, F_CONSTRNC));
+            int bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, molinfo, F_SETTLE);
+            double_check(ir, state->box,
+                         bHasNormalConstraints,
+                         bHasAnyConstraints,
+                         wi);
+        }
     }
 
     if (bGenVel)
index 37819fb4b19641c6419e03c8a997ff2e52d18f9c..391707a6473d955592bd646393bbfe25f539335e 100644 (file)
@@ -4281,7 +4281,10 @@ void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
     check_disre(sys);
 }
 
-void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
+void double_check(t_inputrec *ir, matrix box,
+                  gmx_bool bHasNormalConstraints,
+                  gmx_bool bHasAnyConstraints,
+                  warninp_t wi)
 {
     real        min_size;
     gmx_bool    bTWIN;
@@ -4294,7 +4297,7 @@ void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
         warning_error(wi, ptr);
     }
 
-    if (bConstr && ir->eConstrAlg == econtSHAKE)
+    if (bHasNormalConstraints && ir->eConstrAlg == econtSHAKE)
     {
         if (ir->shake_tol <= 0.0)
         {
@@ -4317,7 +4320,7 @@ void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
         }
     }
 
-    if ( (ir->eConstrAlg == econtLINCS) && bConstr)
+    if ( (ir->eConstrAlg == econtLINCS) && bHasNormalConstraints)
     {
         /* If we have Lincs constraints: */
         if (ir->eI == eiMD && ir->etc == etcNO &&
@@ -4338,7 +4341,7 @@ void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
         }
     }
 
-    if (bConstr && ir->epc == epcMTTK)
+    if (bHasAnyConstraints && ir->epc == epcMTTK)
     {
         warning_note(wi, "MTTK with constraints is deprecated, and will be removed in GROMACS 5.1");
     }
index e20453e20042a3cfbd8f2aafb223bfbea18dce0c..dcfd3fc2fcfd1132769eeaec72aba3d2607ea535 100644 (file)
@@ -90,7 +90,9 @@ void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
 int search_string(const char *s, int ng, char *gn[]);
 /* Returns the index of string s in the index groups */
 
-void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr,
+void double_check(t_inputrec *ir, matrix box,
+                  gmx_bool bHasNormalConstraints,
+                  gmx_bool bHasAnyConstraints,
                   warninp_t wi);
 /* Do more checks */