Fix membed with partial revert of 29943f
authorErik Lindahl <erik@kth.se>
Sat, 9 Jul 2016 21:32:58 +0000 (23:32 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 11 Jul 2016 14:50:12 +0000 (16:50 +0200)
The membrane embedding algorithm must be initialized before
we call init_forcerec(), so it cannot trivially be moved into
do_md(). This has to be cleaned up anyway for release-2017
since we will remove the group scheme be then, but for now
this fix will allow us have the method working in release-2016.

Fixes #1998.

Change-Id: I8e16a67b8df16999793649d79eebca051b8b474f

src/gromacs/mdlib/integrator.h
src/gromacs/mdlib/minimize.cpp
src/gromacs/mdlib/tpi.cpp
src/programs/mdrun/md.cpp
src/programs/mdrun/runner.cpp

index d19a861b3a8656d2c611d957253096bbc128d449..aea0ec37ab5f7d9ea6f7fc0c970547fd6a113aaa 100644 (file)
@@ -87,6 +87,7 @@ namespace gmx
  * \param[in] repl_ex_nst         How often we do replica exchange (in steps)
  * \param[in] repl_ex_nex         How many replicas we have
  * \param[in] repl_ex_seed        The seed for Monte Carlo swaps
+ * \param[in] membed              Membrane embedding data structure
  * \param[in] cpt_period          How often to checkpoint the simulation
  * \param[in] max_hours           Maximume length of the simulation (wall time)
  * \param[in] imdport             Interactive MD port (socket)
@@ -107,6 +108,7 @@ typedef double integrator_t (FILE *fplog, t_commrec *cr,
                              gmx_edsam_t ed,
                              t_forcerec *fr,
                              int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
+                             gmx_membed_t gmx_unused * membed,
                              real cpt_period, real max_hours,
                              int imdport,
                              unsigned long Flags,
index 55200ce0ab8986d9817bad13d36f2965e2521dbf..3890e24e3734a167e3e6c9efab91a204945dcbc1 100644 (file)
@@ -1015,6 +1015,7 @@ namespace gmx
                            gmx_edsam_t ed,
                            t_forcerec *fr,
                            int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
+                           gmx_membed_t gmx_unused *membed,
                            real cpt_period, real max_hours,
                            int imdport,
                            unsigned long Flags,
@@ -1034,6 +1035,7 @@ double do_cg(FILE *fplog, t_commrec *cr,
              gmx_edsam_t gmx_unused ed,
              t_forcerec *fr,
              int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
+             gmx_membed_t gmx_unused *membed,
              real gmx_unused cpt_period, real gmx_unused max_hours,
              int imdport,
              unsigned long gmx_unused Flags,
@@ -1680,6 +1682,7 @@ double do_lbfgs(FILE *fplog, t_commrec *cr,
                 gmx_edsam_t gmx_unused ed,
                 t_forcerec *fr,
                 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
+                gmx_membed_t gmx_unused *membed,
                 real gmx_unused cpt_period, real gmx_unused max_hours,
                 int imdport,
                 unsigned long gmx_unused Flags,
@@ -2507,6 +2510,7 @@ double do_steep(FILE *fplog, t_commrec *cr,
                 gmx_edsam_t gmx_unused  ed,
                 t_forcerec *fr,
                 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
+                gmx_membed_t gmx_unused *membed,
                 real gmx_unused cpt_period, real gmx_unused max_hours,
                 int imdport,
                 unsigned long gmx_unused Flags,
@@ -2774,6 +2778,7 @@ double do_nm(FILE *fplog, t_commrec *cr,
              gmx_edsam_t  gmx_unused ed,
              t_forcerec *fr,
              int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
+             gmx_membed_t gmx_unused *membed,
              real gmx_unused cpt_period, real gmx_unused max_hours,
              int imdport,
              unsigned long gmx_unused Flags,
index 23bda0206af1351cbf6ac3946ea80c32e7e61cb4..aa1210a59fdc9641dc89ea5424e16311e0162ebb 100644 (file)
@@ -159,6 +159,7 @@ double do_tpi(FILE *fplog, t_commrec *cr,
               gmx_edsam_t gmx_unused ed,
               t_forcerec *fr,
               int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
+              gmx_membed_t gmx_unused *membed,
               real gmx_unused cpt_period, real gmx_unused max_hours,
               int gmx_unused imdport,
               unsigned long gmx_unused Flags,
index 1ac10d0e721cce2e28af672cf28a62805955a567..2a49aba4a37b3d20df5bb7b49f472fbf1f96000c 100644 (file)
@@ -219,6 +219,7 @@ double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
                   t_nrnb *nrnb, gmx_wallcycle_t wcycle,
                   gmx_edsam_t ed, t_forcerec *fr,
                   int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
+                  gmx_membed_t *membed,
                   real cpt_period, real max_hours,
                   int imdport,
                   unsigned long Flags,
@@ -283,7 +284,6 @@ double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
 
     /* Interactive MD */
     gmx_bool          bIMDstep = FALSE;
-    gmx_membed_t     *membed   = NULL;
 
 #ifdef GMX_FAHCORE
     /* Temporary addition for FAHCORE checkpointing */
@@ -340,18 +340,6 @@ double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
     }
     groups = &top_global->groups;
 
-    if (opt2bSet("-membed", nfile, fnm))
-    {
-        if (MASTER(cr))
-        {
-            fprintf(stderr, "Initializing membed");
-        }
-        /* Note that membed cannot work in parallel because mtop is
-         * changed here. Fix this if we ever want to make it run with
-         * multiple ranks. */
-        membed = init_membed(fplog, nfile, fnm, top_global, ir, state_global, cr, &cpt_period);
-    }
-
     if (ir->eSwapCoords != eswapNO)
     {
         /* Initialize ion swapping code */
@@ -1849,11 +1837,6 @@ double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
         finish_swapcoords(ir->swap);
     }
 
-    if (membed != nullptr)
-    {
-        free_membed(membed);
-    }
-
     /* IMD cleanup, if bIMD is TRUE. */
     IMD_finalize(ir->bIMD, ir->imd);
 
index ddc5895c956b11c1abdb083777dc6fce8ed9fa30..6fb9d5f5707f52af4c545f5a144880af27d4053b 100644 (file)
 
 #include "deform.h"
 #include "md.h"
+#include "membed.h"
 #include "repl_ex.h"
 #include "resource-division.h"
 
@@ -709,6 +710,7 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
     gmx_int64_t               reset_counters;
     gmx_edsam_t               ed           = NULL;
     int                       nthreads_pme = 1;
+    gmx_membed_t *            membed       = NULL;
     gmx_hw_info_t            *hwinfo       = NULL;
     /* The master rank decides early on bUseGPU and broadcasts this later */
     gmx_bool                  bUseGPU            = FALSE;
@@ -1142,6 +1144,19 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
         wcycle_set_reset_counters(wcycle, reset_counters);
     }
 
+    // Membrane embedding must be initialized before we call init_forcerec()
+    if (opt2bSet("-membed", nfile, fnm))
+    {
+        if (MASTER(cr))
+        {
+            fprintf(stderr, "Initializing membed");
+        }
+        /* Note that membed cannot work in parallel because mtop is
+         * changed here. Fix this if we ever want to make it run with
+         * multiple ranks. */
+        membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
+    }
+
     snew(nrnb, 1);
     if (cr->duty & DUTY_PP)
     {
@@ -1317,6 +1332,7 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
                                      fcd, state,
                                      mdatoms, nrnb, wcycle, ed, fr,
                                      repl_ex_nst, repl_ex_nex, repl_ex_seed,
+                                     membed,
                                      cpt_period, max_hours,
                                      imdport,
                                      Flags,
@@ -1355,6 +1371,11 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
     /* Free GPU memory and context */
     free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
 
+    if (membed != nullptr)
+    {
+        free_membed(membed);
+    }
+
     gmx_hardware_info_free(hwinfo);
 
     /* Does what it says */