Coverage for biobb_mem/mdanalysis_biobb/mda_hole.py: 63%
91 statements
« prev ^ index » next coverage.py v7.10.6, created at 2025-09-08 09:07 +0000
« prev ^ index » next coverage.py v7.10.6, created at 2025-09-08 09:07 +0000
1#!/usr/bin/env python3
3"""Module containing the MDAnalysis HOLE class and the command line interface."""
4import re
5import os
6import numpy as np
7import pandas as pd
8from biobb_common.generic.biobb_object import BiobbObject
9from biobb_common.tools.file_utils import launchlogger
10import MDAnalysis as mda
11from mdahole2.analysis import HoleAnalysis
14class MDAHole(BiobbObject):
15 """
16 | biobb_mem MDAHole
17 | Wrapper of the MDAnalysis HOLE module for analyzing ion channel pores or transporter pathways.
18 | MDAnalysis HOLE provides an interface to the HOLE suite of tools to analyze pore dimensions and properties along a channel or transporter pathway. The parameter names and defaults follow the `MDAnalysis HOLE <https://www.mdanalysis.org/mdahole2/api.html>`_ implementation.
20 Args:
21 input_top_path (str): Path to the input structure or topology file. File type: input. `Sample file <https://github.com/bioexcel/biobb_mem/raw/main/biobb_mem/test/data/A01JD/A01JD.pdb>`_. Accepted formats: crd (edam:3878), gro (edam:2033), mdcrd (edam:3878), mol2 (edam:3816), pdb (edam:1476), pdbqt (edam:1476), prmtop (edam:3881), psf (edam:3882), top (edam:3881), tpr (edam:2333), xml (edam:2332), xyz (edam:3887).
22 input_traj_path (str): Path to the input trajectory to be processed. File type: input. `Sample file <https://github.com/bioexcel/biobb_mem/raw/main/biobb_mem/test/data/A01JD/A01JD.xtc>`_. Accepted formats: arc (edam:2333), crd (edam:3878), dcd (edam:3878), ent (edam:1476), gro (edam:2033), inpcrd (edam:3878), mdcrd (edam:3878), mol2 (edam:3816), nc (edam:3650), pdb (edam:1476), pdbqt (edam:1476), restrt (edam:3886), tng (edam:3876), trr (edam:3910), xtc (edam:3875), xyz (edam:3887).
23 output_hole_path (str): Path to the output HOLE analysis results. File type: output. `Sample file <https://github.com/bioexcel/biobb_mem/raw/main/biobb_mem/test/reference/mdanalysis_biobb/hole.vmd>`_. Accepted formats: vmd (edam:format_2330).
24 output_csv_path (str): Path to the output CSV file containing the radius and coordinates of the pore. File type: output. `Sample file <https://github.com/bioexcel/biobb_mem/raw/main/biobb_mem/test/reference/mdanalysis_biobb/hole_profile.csv>`_. Accepted formats: csv (edam:format_3752).
25 properties (dic - Python dictionary object containing the tool parameters, not input/output files):
26 * **start** (*int*) - (None) Starting frame for slicing.
27 * **stop** (*int*) - (None) Ending frame for slicing.
28 * **steps** (*int*) - (None) Step for slicing.
29 * **executable** (*str*) - ("hole") Path to the HOLE executable.
30 * **select** (*str*) - ("protein") The selection string to create an atom selection that the HOLE analysis is applied to.
31 * **cpoint** (*list*) - (None) Coordinates of a point inside the pore (Å). If None, tries to guess based on the geometry.
32 * **cvect** (*list*) - (None) Search direction vector. If None, tries to guess based on the geometry.
33 * **sample** (*float*) - (0.2) Distance of sample points in Å. This value determines how many points in the pore profile are calculated.
34 * **end_radius** (*float*) - (22) Radius in Å, which is considered to be the end of the pore.
35 * **dot_density** (*int*) - (15) [5~35] Density of facets for generating a 3D pore representation.
36 * **remove_tmp** (*bool*) - (True) [WF property] Remove temporal files.
37 * **restart** (*bool*) - (False) [WF property] Do not execute if output files exist.
38 * **sandbox_path** (*str*) - ("./") [WF property] Parent path to the sandbox directory.
40 Examples:
41 This is a use example of how to use the building block from Python::
43 from biobb_mem.mdanalysis_biobb.mda_hole import mda_hole
44 prop = {
45 'select': 'protein',
46 'executable': 'hole'
47 }
48 mda_hole(input_top_path='/path/to/myTopology.pdb',
49 input_traj_path='/path/to/myTrajectory.xtc',
50 output_hole_path='/path/to/hole_analysis.csv',
51 output_hole_path='/path/to/hole_profile.csv',
52 properties=prop)
54 Info:
55 * wrapped_software:
56 * name: MDAnalysis
57 * version: 2.7.0
58 * license: GNU
59 * ontology:
60 * name: EDAM
61 * schema: http://edamontology.org/EDAM.owl
62 """
64 def __init__(self, input_top_path, input_traj_path, output_hole_path, output_csv_path=None,
65 properties=None, **kwargs) -> None:
66 properties = properties or {}
68 # Call parent class constructor
69 super().__init__(properties)
70 self.locals_var_dict = locals().copy()
72 # Input/Output files
73 self.io_dict = {
74 "in": {"input_top_path": input_top_path, "input_traj_path": input_traj_path},
75 "out": {"output_hole_path": output_hole_path, "output_csv_path": output_csv_path}
76 }
78 # Properties specific for MDAHole
79 self.start = properties.get('start', None)
80 self.stop = properties.get('stop', None)
81 self.steps = properties.get('steps', None)
82 self.executable = properties.get('executable', 'hole')
83 self.select = properties.get('select', 'protein')
84 self.cpoint = properties.get('cpoint', None)
85 self.cvect = properties.get('cvect', None)
86 self.sample = properties.get('sample', 0.2)
87 self.end_radius = properties.get('end_radius', 22)
88 self.dot_density = properties.get('dot_density', 15)
89 self.properties = properties
91 # Check the properties
92 self.check_properties(properties)
93 self.check_arguments()
95 @launchlogger
96 def launch(self) -> int:
97 """Execute the :class:`MDAHole <mdanalysis_biobb.mda_hole.MDAHole>` object."""
99 # Setup Biobb
100 if self.check_restart():
101 return 0
102 self.stage_files()
104 # Load the universe
105 u = mda.Universe(self.stage_io_dict["in"]["input_top_path"],
106 self.stage_io_dict["in"]["input_traj_path"])
107 # save current directory and move to temporary
108 cwd = os.getcwd()
109 os.chdir(self.stage_io_dict.get("unique_dir"))
110 # Create HoleAnalysis object
111 hole = HoleAnalysis(
112 universe=u,
113 select=self.select,
114 cpoint=self.cpoint,
115 cvect=self.cvect,
116 sample=self.sample,
117 executable=self.executable,
118 end_radius=self.end_radius
119 )
120 # Run the analysis with step parameter
121 hole.run(
122 start=self.start,
123 stop=self.stop,
124 step=self.steps
125 )
126 # Save the results to a CSV file
127 all_frames = []
128 for frame in hole.results.profiles.keys():
129 rxn_coord = hole.results.profiles[frame].rxn_coord
130 radius = hole.results.profiles[frame].radius
131 df_frame = pd.DataFrame({'Frame': frame, 'Pore Coordinate': rxn_coord, 'Radius': radius})
132 all_frames.append(df_frame)
133 # Concatenate all frames into a single DataFrame
134 df_all_frames = pd.concat(all_frames, ignore_index=True)
135 df_all_frames.to_csv(self.stage_io_dict["out"]["output_csv_path"], index=False)
137 hole.create_vmd_surface(
138 self.stage_io_dict["out"]["output_hole_path"],
139 dot_density=self.dot_density
140 )
141 hole.delete_temporary_files()
142 # move back to original directory
143 os.chdir(cwd)
144 # Copy files to host
145 self.copy_to_host()
146 # remove temporary folder(s)
147 self.remove_tmp_files()
149 self.check_arguments(output_files_created=True, raise_exception=False)
151 return self.return_code
154def mda_hole(input_top_path: str,
155 input_traj_path: str,
156 output_hole_path: str,
157 output_csv_path: str,
158 properties: dict = None,
159 **kwargs) -> int:
160 """Execute the :class:`MDAHole <mdanalysis_biobb.mda_hole.MDAHole>` class and
161 execute the :meth:`launch() <mdanalysis_biobb.mda_hole.MDAHole.launch>` method."""
162 return MDAHole(**dict(locals())).launch()
165mda_hole.__doc__ = MDAHole.__doc__
166main = MDAHole.get_main(mda_hole, "Analyze ion channel pores or transporter pathways.")
169def display_hole(input_top_path: str, input_traj_path: str,
170 output_hole_path: str = 'hole.vmd',
171 frame: int = 0, opacity: float = 0.9):
172 """
173 Visualize a channel using NGLView from a VMD file.
175 Args:
176 input_top_path (str): Path to the input topology file.
177 output_hole_path (str, optional): Path to the VMD file containing the channel data. Default is 'hole.vmd'.
178 frame (int, optional): Frame index to visualize. Default is 0.
179 opacity (float, optional): Opacity of the visualization. Default is 0.9.
180 Returns:
181 nglview.NGLWidget: NGLView widget for visualizing the channel.
182 """
184 try:
185 import nglview as nv
186 except ImportError:
187 raise ImportError('Please install the nglview package to visualize the channel.')
189 # Read the VMD file and parse triangles
190 with open(output_hole_path, 'r') as f:
191 lines = f.readlines()
193 # Find lines with triangle coordinates
194 trinorms = []
195 for i, line in enumerate(lines):
196 if i > 3 and 'set triangle' in line:
197 vmd_set = re.sub(r'set triangles\(\d+\)', '', line) # Remove set triangles(i)
198 vmd_set = re.sub(r'\{(\s*-?\d[^\s]*)(\s*-?\d[^\s]*)(\s*-?\d[^}]*)\}', r'[\1,\2,\3]', vmd_set) # Convert { x y z } to [x,y,z]
199 vmd_set = vmd_set.replace('{', '[').replace('}', ']') # Convert { to [ and } to ]
200 vmd_set = re.sub(r'\]\s*\[', '], [', vmd_set) # Add commas between brackets
201 vmd_set = eval(vmd_set.strip()) # Evaluate string as list
202 # different hole colors
203 trinorms.append(vmd_set)
204 # Create a list of positions, colors, and normals
205 colors = np.array([[1, 0, 0], # red
206 [0, 1, 0], # green
207 [0, 0, 1]]) # blue
208 poss, cols, nors = [], [], []
209 for i, color in enumerate(colors):
210 if len(trinorms[frame][i]) > 0:
211 col_dat = np.array(trinorms[frame][i])
212 poss.append(col_dat[:, :3, :].flatten()) # 3 first elements are positions
213 cols.append((np.zeros(col_dat.shape[0]*18).reshape(-1, 3) + color).flatten()) # 3 colors for each vertex
214 nors.append(col_dat[:, 3:, :].flatten()) # 3 last elements are normals
215 poss = np.concatenate(poss)
216 cols = np.concatenate(cols)
217 nors = np.concatenate(nors)
218 # Create NGLView widget
219 u = mda.Universe(input_top_path, input_traj_path)
220 view = nv.show_mdanalysis(u)
221 # view.clear_representations()
222 mesh = ('mesh', poss, cols)
223 view._add_shape([mesh], name='my_shape')
224 view.update_representation(component=1, repr_index=0, opacity=opacity)
225 return view
228if __name__ == '__main__':
229 main()