Add membrane embedding docs
authorErik Lindahl <erik@kth.se>
Sat, 13 Jun 2015 08:22:28 +0000 (10:22 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 16 Jun 2015 20:37:14 +0000 (22:37 +0200)
These were killed by mistake in 86601a237df when the
separate g_membed program was merged into mdrun.

Fixes #1724.

Change-Id: If97dbf13c269e0711cda9798fac03225920c0d6f

docs/user-guide/mdrun-features.rst
src/programs/mdrun/mdrun.cpp
src/programs/mdrun/membed.c

index 5a2e5bf38e82d0e04820f14a913a4cc49fd04e82..e8724511c97483cb75600d8c7174accef5af4373 100644 (file)
@@ -131,3 +131,104 @@ and executes 100 steps. `mdrun -maxh 2.5` will terminate the
 simulation shortly before 2.5 hours elapse, which can be useful when
 running under cluster queues (as long as the queueing system does not
 ever suspend the simulation).
+
+Running a membrane protein embedding simulation
+-----------------------------------------------
+
+This is a module to help embed a membrane protein into an equilibrated
+lipid bilayer at a position and orientation specified by the user. The
+main advantage is that it is possible to use very complex lipid bilayers
+with a number of different components that have been relaxed for a
+long time in a previous simulation. In theory that could be accomplished
+with a procedure similar to genbox, but since lipids are much larger
+than water molecules it will lead to a large vacuum layer between the
+protein and membrane if we remove all molecules where any atom is
+overlapping. Instead, this module works by first artifically shrinking
+the protein in the xy-plane, then it removes lipids that overlap with
+a much smaller core, after which we gradually push the protein atoms
+back to their initial positions, while using normal dynamics for the
+rest of the system so lipids adapt to the protein.
+
+To use membrane embedding, start by building a lipid bilayer that is
+just-so-slightly larger in the xy-plane than what you expect to need
+in the end, and make sure you have enough water outside the membrane
+to accomodate globular domains. Place the protein in the same coordinate
+file (and topology) as your lipid bilayer, and make sure it is in the
+orientation and position you want right in the middle of the bilayer.
+
+The first settings have to be entered in the mdp file that controls
+your simulation. You need an energy group corresponding to your
+protein, this group should be frozen (all dimensions), and we should
+exclude all interactions inside the protein to avoid problems when it
+is distorted. For instance:
+
+::
+
+    integrator     = md
+    energygrps     = Protein
+    freezegrps     = Protein
+    freezedim      = Y Y Y
+    energygrp_excl = Protein Protein
+
+You will also need a number of settings for the actual membrane
+embedding process. These are entered as similar name and value pairs,
+but in the separate text data file ``embed.dat`` that you provide as
+the argument to the ``-membed`` option (we refer to the below
+when explaining the process). The embedding works in for stages:
+
+1. The protein is resized around its center of mass by a factor
+   ``xy`` in the xy-plane (the bilayer plane), and a factor ``z``
+   along the z-axis (normal to the bilayer). If the height of the
+   protein is the same or smaller than the thickness of the
+   membrane, a z-fraction larger than 1.0 can prevent the protein
+   from being enveloped by the lipids.
+
+2. All lipid and solvent molecules overlapping with the resized
+   protein are removed. All interactions inside the protein are
+   turned off to prevent numerical issues for small values of the
+   scaling fractions.
+
+3. A single md step is performed, where atoms in the rest of the
+   system are moved.
+
+4. The resize factors are adjusted by the small amounts
+   (1-xy)/nxy and (1-z)/nz, where ``nxy`` and ``nz`` are the
+   number of iterations to use.  The resize factor for the xy-plane
+   is adjusted first. The resize factor for the z-direction is not
+   changed until the xy factor is 1.0 (after ``nxy`` iterations).
+
+5. Steps 3 and 4 are repeated until the protein has again reached
+   its original size, i.e. after nxy+nz iterations. After the
+   embedding you might still want to perform a short relaxation.
+
+Parameters that can be specified in ``embed.dat``, with default
+values that will be used if the setting is omitted:
+
+- ``xyinit`` (0.5) Resize factor for the protein in the xy
+  dimension before starting embedding.
+
+- ``xyend`` (1.0) Final resize factor in the xy dimension.
+
+- ``zinit`` (1.0) Resize factor for the protein in the z
+   dimension before starting embedding.
+
+- ``zend`` (1.0) Final resize faction in the z dimension.
+
+- ``nxy`` (1000) Number of iteration for the xy dimension.
+
+- ``nz`` (0) Number of iterations for the z dimension.
+
+- ``rad`` (0.22) Probe radius to check for overlap between
+  the group to embed and the membrane.
+
+- ``pieces`` (1) Perform piecewise resize. Select parts of the group
+  to insert and resize these with respect to their own geometrical center.
+
+- ``asymmetry`` (no) Allow asymmetric insertion, i.e. the number of
+  lipids removed from the upper and lower leaflet will not be checked.
+
+- ``ndiff`` (0) Number of lipids that will additionally be removed
+  from the lower (negative number) or upper (positive number)
+  membrane leaflet.
+
+- ``maxwarn`` (0) Largest number of membed warnings allowed.
index f924821d93c1c7cc20bb9ccacce4ec9c8641de03..b5fd731267dbafca004db2b71abaf79d924e48e6 100644 (file)
@@ -146,9 +146,11 @@ int gmx_mdrun(int argc, char *argv[])
         "investigation are: polarizability.",
         "[PAR]",
         "The option [TT]-membed[tt] does what used to be g_membed, i.e. embed",
-        "a protein into a membrane. The data file should contain the options",
-        "that where passed to g_membed before. The [TT]-mn[tt] and [TT]-mp[tt]",
-        "both apply to this as well.",
+        "a protein into a membrane. This module requires a number of settings",
+        "that are provided in a data file that is the argument of this option.",
+        "For more details in membrane embedding, see the documentation in the",
+        "user guide. The options [TT]-mn[tt] and [TT]-mp[tt] are used to provide",
+        "the index and topology files used for the embedding.",
         "[PAR]",
         "The option [TT]-pforce[tt] is useful when you suspect a simulation",
         "crashes due to too large forces. With this option coordinates and",
index fae5d4ddc3b9f919dc8df0b17c9821d91d387766..64c3680ffa6e37d4b11f220218f2d623401b56fd 100644 (file)
@@ -1142,14 +1142,15 @@ gmx_membed_t init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop_
         {
             warn++;
             fprintf(stderr, "\nWarning %d;\nThe number of steps used to grow the z-coordinate of %s (%d)"
-                    " is probably too small.\nIncrease -nz or maxwarn.\n\n", warn, ins, it_z);
+                    " is probably too small.\nIncrease -nz or the maxwarn setting in the membed input file.\n\n", warn, ins, it_z);
         }
 
         if (it_xy+it_z > inputrec->nsteps)
         {
             warn++;
             fprintf(stderr, "\nWarning %d:\nThe number of growth steps (-nxy + -nz) is larger than the "
-                    "number of steps in the tpr.\n\n", warn);
+                    "number of steps in the tpr.\n"
+                    "(increase maxwarn in the membed input file to override)\n\n", warn);
         }
 
         fr_id = -1;
@@ -1225,13 +1226,13 @@ gmx_membed_t init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop_
             warn++;
             fprintf(stderr, "\nWarning %d:\nThe xy-area is very small compared to the area of the protein.\n"
                     "This might cause pressure problems during the growth phase. Just try with\n"
-                    "current setup (-maxwarn + 1), but if pressure problems occur, lower the\n"
-                    "compressibility in the mdp-file or use no pressure coupling at all.\n\n", warn);
+                    "current setup and increase 'maxwarn' in your membed settings file, but lower the\n"
+                    "compressibility in the mdp-file or disable pressure coupling if problems occur.\n\n", warn);
         }
 
         if (warn > maxwarn)
         {
-            gmx_fatal(FARGS, "Too many warnings.\n");
+            gmx_fatal(FARGS, "Too many warnings (override by setting maxwarn in the membed input file)\n");
         }
 
         printf("The estimated area of the protein in the membrane is %.3f nm^2\n", prot_area);
@@ -1294,7 +1295,7 @@ gmx_membed_t init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop_
             warn++;
             fprintf(stderr, "\nWarning %d:\nTrying to remove a larger lipid area than the estimated "
                     "protein area\nTry making the -xyinit resize factor smaller or increase "
-                    "maxwarn.\n\n", warn);
+                    "maxwarn in the membed input file.\n\n", warn);
         }
 
         /*remove all lipids and waters overlapping and update all important structures*/
@@ -1310,8 +1311,8 @@ gmx_membed_t init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop_
 
         if (warn > maxwarn)
         {
-            gmx_fatal(FARGS, "Too many warnings.\nIf you are sure these warnings are harmless, "
-                      "you can increase -maxwarn");
+            gmx_fatal(FARGS, "Too many warnings.\nIf you are sure these warnings are harmless,\n"
+                      "you can increase the maxwarn setting in the membed input file.");
         }
 
         if (ftp2bSet(efTOP, nfile, fnm))