Fixed non-working flat-bottom posres code.
authorJochen Hub <jhub@gwdg.de>
Tue, 10 Dec 2013 09:17:07 +0000 (10:17 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 9 Jan 2014 19:58:31 +0000 (20:58 +0100)
Change-Id: I905ba495981370bd9346909a707a61e10a852907

src/gromacs/gmxlib/bondfree.c
src/gromacs/legacyheaders/types/idef.h
src/gromacs/mdlib/sim_util.c

index 058398eafb3bcb8128a2a71f9d2ab2979bb63950..40989ed89c011d0e3c3939114c2da732f1bfbf65 100644 (file)
@@ -3796,7 +3796,7 @@ static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
 {
     return
         (interaction_function[ftype].flags & IF_BOND) &&
-        !(ftype == F_CONNBONDS || ftype == F_POSRES) &&
+        !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
         (ftype < F_GB12 || ftype > F_GB14);
 }
 
index 0883e1dcc5aba20164d26dfa77fc578b7b15d49e..3cf75bac3ac32aa412a77c8b7ea092ec22c289b2 100644 (file)
@@ -148,7 +148,7 @@ enum {
     F_NRE               /* This number is for the total number of energies     */
 };
 
-#define IS_RESTRAINT_TYPE(ifunc) (((ifunc == F_POSRES) || (ifunc == F_DISRES) || (ifunc == F_RESTRBONDS) || (ifunc == F_DISRESVIOL) || (ifunc == F_ORIRES) || (ifunc == F_ORIRESDEV) || (ifunc == F_ANGRES) || (ifunc == F_ANGRESZ) || (ifunc == F_DIHRES)))
+#define IS_RESTRAINT_TYPE(ifunc) (((ifunc == F_POSRES) || (ifunc == F_FBPOSRES) || (ifunc == F_DISRES) || (ifunc == F_RESTRBONDS) || (ifunc == F_DISRESVIOL) || (ifunc == F_ORIRES) || (ifunc == F_ORIRESDEV) || (ifunc == F_ANGRES) || (ifunc == F_ANGRESZ) || (ifunc == F_DIHRES)))
 
 /* A macro for checking if ftype is an explicit pair-listed LJ or COULOMB
  * interaction type:
index 86243251c7afe67a1e44df91ea2ac598e22a535f..b76cfcc948ca5e23e6c1dd2bf6578904bfe06cca 100644 (file)
@@ -351,6 +351,27 @@ static void posres_wrapper(FILE *fplog,
     }
 }
 
+static void fbposres_wrapper(t_inputrec *ir,
+                             t_nrnb *nrnb,
+                             gmx_localtop_t *top,
+                             matrix box, rvec x[],
+                             gmx_enerdata_t *enerd,
+                             t_forcerec *fr)
+{
+    t_pbc pbc;
+    real  v;
+
+    /* Flat-bottomed position restraints always require full pbc */
+    set_pbc(&pbc, ir->ePBC, box);
+    v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
+                 top->idef.iparams_fbposres,
+                 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
+                 ir->ePBC == epbcNONE ? NULL : &pbc,
+                 fr->rc_scaling, fr->ePBC, fr->posres_com);
+    enerd->term[F_FBPOSRES] += v;
+    inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
+}
+
 static void pull_potential_wrapper(FILE *fplog,
                                    gmx_bool bSepDVDL,
                                    t_commrec *cr,
@@ -1158,6 +1179,11 @@ void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
                        enerd, lambda, fr);
     }
 
+    if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
+    {
+        fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
+    }
+
     /* Compute the bonded and non-bonded energies and optionally forces */
     do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
                       cr, nrnb, wcycle, mdatoms,
@@ -1733,18 +1759,7 @@ void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
 
     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
     {
-        /* Flat-bottomed position restraints always require full pbc */
-        if (!(bStateChanged && bDoAdressWF))
-        {
-            set_pbc(&pbc, inputrec->ePBC, box);
-        }
-        v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
-                     top->idef.iparams_fbposres,
-                     (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
-                     inputrec->ePBC == epbcNONE ? NULL : &pbc,
-                     fr->rc_scaling, fr->ePBC, fr->posres_com);
-        enerd->term[F_FBPOSRES] += v;
-        inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
+        fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
     }
 
     /* Compute the bonded and non-bonded energies and optionally forces */