Skip to content

Commit

Permalink
Added wrap transformations (#48)
Browse files Browse the repository at this point in the history
* Added coordinate wrapping with tests.
  • Loading branch information
Estefania Barreto-Ojeda authored Aug 2, 2021
1 parent 8df8b78 commit 8a52295
Show file tree
Hide file tree
Showing 5 changed files with 316 additions and 16 deletions.
5 changes: 3 additions & 2 deletions membrane_curvature/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -11,8 +11,9 @@
#

# Add imports here
from membrane_curvature.surface import normalized_grid, derive_surface, get_z_surface
from membrane_curvature.curvature import mean_curvature, gaussian_curvature
from .surface import normalized_grid, derive_surface, get_z_surface
from .curvature import mean_curvature, gaussian_curvature
from .base import MembraneCurvature

# Handle versioneer
from ._version import get_versions
Expand Down
42 changes: 36 additions & 6 deletions membrane_curvature/base.py
Original file line number Diff line number Diff line change
Expand Up @@ -35,8 +35,8 @@ class MembraneCurvature(AnalysisBase):
select : str or iterable of str, optional.
The selection string of an atom selection to use as a
reference to derive a surface.
pbc : bool, optional
Apply periodic boundary conditions.
wrap : bool, optional
Apply coordinate wrapping to pack atoms into the primary unit cell.
n_x_bins : int, optional, default: '100'
Number of bins in grid in the x dimension.
n_y_bins : int, optional, default: '100'
Expand Down Expand Up @@ -68,8 +68,22 @@ class MembraneCurvature(AnalysisBase):
Each array has shape (`n_x_bins`, `n_y_bins`)
See also
--------
`MDAnalysis.transformations.wrap
<https://docs.mdanalysis.org/1.0.0/documentation_pages/transformations/wrap.html>`_
Notes
-----
Use `wrap=True` to translates the atoms of your `mda.Universe` back
in the unit cell. Use `wrap=False` for processed trajectories where
rotational/translational fit is performed.
For more details on when to use `wrap=True`, check the `Usage
<https://membrane-curvature.readthedocs.io/en/latest/source/pages/Usage.html>`_
page.
The derived surface and calculated curvatures are available in the
:attr:`results` attributes.
Expand All @@ -95,20 +109,22 @@ class MembraneCurvature(AnalysisBase):
mean_curvature = mc.results.average_mean_curvature
gaussian_curvature = mc.results.average_gaussian_curvature
The respective 2D curvature plots can be obtained using the `matplotlib` package for
data visualization via `imshow`. We recommend using the `gaussian` interpolation.
The respective 2D curvature plots can be obtained using the `matplotlib`
package for data visualization via `imshow`. We recommend using the
`gaussian` interpolation.
"""

def __init__(self, universe, select='all',
n_x_bins=100, n_y_bins=100,
x_range=None,
y_range=None,
pbc=True, **kwargs):
wrap=True, **kwargs):

super().__init__(universe.universe.trajectory, **kwargs)
self.ag = universe.select_atoms(select)
self.pbc = pbc
self.wrap = wrap
self.n_x_bins = n_x_bins
self.n_y_bins = n_y_bins
self.x_range = x_range if x_range else (0, universe.dimensions[0])
Expand All @@ -128,6 +144,17 @@ def __init__(self, universe, select='all',
warnings.warn(msg)
logger.warn(msg)

# Apply wrapping coordinates
if not self.wrap:
# Warning
msg = (" `wrap == False` may result in inaccurate calculation "
"of membrane curvature. Surfaces will be derived from "
"a reduced number of atoms. \n "
" Ignore this warning if your trajectory has "
" rotational/translational fit rotations! ")
warnings.warn(msg)
logger.warn(msg)

def _prepare(self):
# Initialize empty np.array with results
self.results.z_surface = np.full((self.n_frames,
Expand All @@ -141,6 +168,9 @@ def _prepare(self):
self.n_y_bins), np.nan)

def _single_frame(self):
# Apply wrapping coordinates
if self.wrap:
self.ag.wrap()
# Populate a slice with np.arrays of surface, mean, and gaussian per frame
self.results.z_surface[self._frame_index] = get_z_surface(self.ag.positions,
n_x_bins=self.n_x_bins,
Expand Down
2 changes: 1 addition & 1 deletion membrane_curvature/data/test_po4_small.gro
Original file line number Diff line number Diff line change
Expand Up @@ -9,4 +9,4 @@ Test file 9 lipids in grid
2708POPC PO4 7 0.000 2.000 12.000
2713POPC PO4 8 1.000 2.000 12.000
2999POPC PO4 9 2.000 2.000 12.000
3.00000 3.00000 3.00000
3.00000 3.00000 30.00000
9 changes: 8 additions & 1 deletion membrane_curvature/surface.py
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,7 @@
"""

import numpy as np
import warnings


def derive_surface(atoms, n_cells_x, n_cells_y, max_width_x, max_width_y):
Expand Down Expand Up @@ -87,11 +88,17 @@ def get_z_surface(coordinates, n_x_bins=10, n_y_bins=10, x_range=(0, 100), y_ran
for l, m, z in zip(cell_x_floor, cell_y_floor, z_coords):

try:
if l < 0 or m < 0:
msg = ("Atom outside grid boundaries. Skipping atom.")
warnings.warn(msg)
continue

grid_z_coordinates[l, m] += z
grid_norm_unit[l, m] += 1

except IndexError:
print("Atom outside grid boundaries. Skipping atom.")
msg = ("Atom outside grid boundaries. Skipping atom.")
warnings.warn(msg)

z_surface = normalized_grid(grid_z_coordinates, grid_norm_unit)

Expand Down
Loading

0 comments on commit 8a52295

Please sign in to comment.