Coverage for biobb_mem/mdanalysis_biobb/mda_hole.py: 57%

104 statements  

« prev     ^ index     » next       coverage.py v7.6.11, created at 2025-02-10 11:25 +0000

1#!/usr/bin/env python3 

2 

3"""Module containing the MDAnalysis HOLE class and the command line interface.""" 

4import re 

5import os 

6import argparse 

7import numpy as np 

8import pandas as pd 

9from biobb_common.generic.biobb_object import BiobbObject 

10from biobb_common.configuration import settings 

11from biobb_common.tools.file_utils import launchlogger 

12import MDAnalysis as mda 

13from mdahole2.analysis import HoleAnalysis 

14 

15 

16class MDAHole(BiobbObject): 

17 """ 

18 | biobb_mem MDAHole 

19 | Wrapper of the MDAnalysis HOLE module for analyzing ion channel pores or transporter pathways. 

20 | 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. 

21 

22 Args: 

23 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). 

24 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). 

25 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). 

26 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. 

27 properties (dic - Python dictionary object containing the tool parameters, not input/output files): 

28 * **start** (*int*) - (None) Starting frame for slicing. 

29 * **stop** (*int*) - (None) Ending frame for slicing. 

30 * **steps** (*int*) - (None) Step for slicing. 

31 * **executable** (*str*) - ("hole") Path to the HOLE executable. 

32 * **select** (*str*) - ("protein") The selection string to create an atom selection that the HOLE analysis is applied to. 

33 * **cpoint** (*list*) - (None) Coordinates of a point inside the pore (Å). If None, tries to guess based on the geometry. 

34 * **cvect** (*list*) - (None) Search direction vector. If None, tries to guess based on the geometry. 

35 * **sample** (*float*) - (0.2) Distance of sample points in Å. This value determines how many points in the pore profile are calculated. 

36 * **end_radius** (*float*) - (22) Radius in Å, which is considered to be the end of the pore. 

37 * **dot_density** (*int*) - (15) [5~35] Density of facets for generating a 3D pore representation. 

38 * **remove_tmp** (*bool*) - (True) [WF property] Remove temporal files. 

39 * **restart** (*bool*) - (False) [WF property] Do not execute if output files exist. 

40 * **sandbox_path** (*str*) - ("./") [WF property] Parent path to the sandbox directory. 

41 

42 Examples: 

43 This is a use example of how to use the building block from Python:: 

44 

45 from biobb_mem.mdanalysis_biobb.mda_hole import mda_hole 

46 prop = { 

47 'select': 'protein', 

48 'executable': 'hole' 

49 } 

50 mda_hole(input_top_path='/path/to/myTopology.pdb', 

51 input_traj_path='/path/to/myTrajectory.xtc', 

52 output_hole_path='/path/to/hole_analysis.csv', 

53 output_hole_path='/path/to/hole_profile.csv', 

54 properties=prop) 

55 

56 Info: 

57 * wrapped_software: 

58 * name: MDAnalysis 

59 * version: 2.7.0 

60 * license: GNU 

61 * ontology: 

62 * name: EDAM 

63 * schema: http://edamontology.org/EDAM.owl 

64 """ 

65 

66 def __init__(self, input_top_path, input_traj_path, output_hole_path, output_csv_path=None, 

67 properties=None, **kwargs) -> None: 

68 properties = properties or {} 

69 

70 # Call parent class constructor 

71 super().__init__(properties) 

72 self.locals_var_dict = locals().copy() 

73 

74 # Input/Output files 

75 self.io_dict = { 

76 "in": {"input_top_path": input_top_path, "input_traj_path": input_traj_path}, 

77 "out": {"output_hole_path": output_hole_path, "output_csv_path": output_csv_path} 

78 } 

79 

80 # Properties specific for MDAHole 

81 self.start = properties.get('start', None) 

82 self.stop = properties.get('stop', None) 

83 self.steps = properties.get('steps', None) 

84 self.executable = properties.get('executable', 'hole') 

85 self.select = properties.get('select', 'protein') 

86 self.cpoint = properties.get('cpoint', None) 

87 self.cvect = properties.get('cvect', None) 

88 self.sample = properties.get('sample', 0.2) 

89 self.end_radius = properties.get('end_radius', 22) 

90 self.dot_density = properties.get('dot_density', 15) 

91 self.properties = properties 

92 

93 # Check the properties 

94 self.check_properties(properties) 

95 self.check_arguments() 

96 

97 @launchlogger 

98 def launch(self) -> int: 

99 """Execute the :class:`MDAHole <mdanalysis_biobb.mda_hole.MDAHole>` class.""" 

100 

101 # Setup Biobb 

102 if self.check_restart(): 

103 return 0 

104 self.stage_files() 

105 

106 # Load the universe 

107 u = mda.Universe(self.stage_io_dict["in"]["input_top_path"], 

108 self.stage_io_dict["in"]["input_traj_path"]) 

109 # save current directory and move to temporary 

110 cwd = os.getcwd() 

111 os.chdir(self.stage_io_dict.get("unique_dir")) 

112 # Create HoleAnalysis object 

113 hole = HoleAnalysis( 

114 universe=u, 

115 select=self.select, 

116 cpoint=self.cpoint, 

117 cvect=self.cvect, 

118 sample=self.sample, 

119 executable=self.executable, 

120 end_radius=self.end_radius 

121 ) 

122 # Run the analysis with step parameter 

123 hole.run( 

124 start=self.start, 

125 stop=self.stop, 

126 step=self.steps 

127 ) 

128 # Save the results to a CSV file 

129 all_frames = [] 

130 for frame in hole.results.profiles.keys(): 

131 rxn_coord = hole.results.profiles[frame].rxn_coord 

132 radius = hole.results.profiles[frame].radius 

133 df_frame = pd.DataFrame({'Frame': frame, 'Pore Coordinate': rxn_coord, 'Radius': radius}) 

134 all_frames.append(df_frame) 

135 # Concatenate all frames into a single DataFrame 

136 df_all_frames = pd.concat(all_frames, ignore_index=True) 

137 df_all_frames.to_csv(self.stage_io_dict["out"]["output_csv_path"], index=False) 

138 

139 hole.create_vmd_surface( 

140 self.stage_io_dict["out"]["output_hole_path"], 

141 dot_density=self.dot_density 

142 ) 

143 hole.delete_temporary_files() 

144 # move back to original directory 

145 os.chdir(cwd) 

146 # Copy files to host 

147 self.copy_to_host() 

148 # remove temporary folder(s) 

149 self.tmp_files.extend([ 

150 self.stage_io_dict.get("unique_dir") 

151 ]) 

152 self.remove_tmp_files() 

153 

154 self.check_arguments(output_files_created=True, raise_exception=False) 

155 

156 return self.return_code 

157 

158 

159def mda_hole(input_top_path: str, input_traj_path: str, output_hole_path: str, output_csv_path: str, properties: dict = None, **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 

163 return MDAHole(input_top_path=input_top_path, 

164 input_traj_path=input_traj_path, 

165 output_hole_path=output_hole_path, 

166 output_csv_path=output_csv_path, 

167 properties=properties, **kwargs).launch() 

168 

169 

170def main(): 

171 """Command line execution of this building block. Please check the command line documentation.""" 

172 parser = argparse.ArgumentParser(description="Analyze ion channel pores or transporter pathways.", formatter_class=lambda prog: argparse.RawTextHelpFormatter(prog, width=99999)) 

173 parser.add_argument('--config', required=False, help='Configuration file') 

174 

175 # Specific args of each building block 

176 required_args = parser.add_argument_group('required arguments') 

177 required_args.add_argument('--input_top_path', required=True, help='Path to the input structure or topology file. Accepted formats: crd, gro, mdcrd, mol2, pdb, pdbqt, prmtop, psf, top, tpr, xml, xyz.') 

178 required_args.add_argument('--input_traj_path', required=True, help='Path to the input trajectory to be processed. Accepted formats: arc, crd, dcd, ent, gro, inpcrd, mdcrd, mol2, nc, pdb, pdbqt, restrt, tng, trr, xtc, xyz.') 

179 required_args.add_argument('--output_hole_path', required=True, help='Path to the output HOLE analysis results. Accepted formats: vmd.') 

180 required_args.add_argument('--output_csv_path', required=True, help='Path to the output CSV file containing the radius and coordinates of the pore. Accepted formats: csv.') 

181 

182 args = parser.parse_args() 

183 args.config = args.config or "{}" 

184 properties = settings.ConfReader(config=args.config).get_prop_dic() 

185 

186 # Specific call of each building block 

187 mda_hole(input_top_path=args.input_top_path, 

188 input_traj_path=args.input_traj_path, 

189 output_hole_path=args.output_hole_path, 

190 output_csv_path=args.output_csv_path, 

191 properties=properties) 

192 

193 

194def display_hole(input_top_path: str, input_traj_path: str, 

195 output_hole_path: str = 'hole.vmd', 

196 frame: int = 0, opacity: float = 0.9): 

197 """ 

198 Visualize a channel using NGLView from a VMD file. 

199 

200 Args: 

201 input_top_path (str): Path to the input topology file. 

202 output_hole_path (str, optional): Path to the VMD file containing the channel data. Default is 'hole.vmd'. 

203 frame (int, optional): Frame index to visualize. Default is 0. 

204 opacity (float, optional): Opacity of the visualization. Default is 0.9. 

205 Returns: 

206 nglview.NGLWidget: NGLView widget for visualizing the channel. 

207 """ 

208 

209 try: 

210 import nglview as nv 

211 except ImportError: 

212 raise ImportError('Please install the nglview package to visualize the channel.') 

213 

214 # Read the VMD file and parse triangles 

215 with open(output_hole_path, 'r') as f: 

216 lines = f.readlines() 

217 

218 # Find lines with triangle coordinates 

219 trinorms = [] 

220 for i, line in enumerate(lines): 

221 if i > 3 and 'set triangle' in line: 

222 vmd_set = re.sub(r'set triangles\(\d+\)', '', line) # Remove set triangles(i) 

223 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] 

224 vmd_set = vmd_set.replace('{', '[').replace('}', ']') # Convert { to [ and } to ] 

225 vmd_set = re.sub(r'\]\s*\[', '], [', vmd_set) # Add commas between brackets 

226 vmd_set = eval(vmd_set.strip()) # Evaluate string as list 

227 # different hole colors 

228 trinorms.append(vmd_set) 

229 # Create a list of positions, colors, and normals 

230 colors = np.array([[1, 0, 0], # red 

231 [0, 1, 0], # green 

232 [0, 0, 1]]) # blue 

233 poss, cols, nors = [], [], [] 

234 for i, color in enumerate(colors): 

235 if len(trinorms[frame][i]) > 0: 

236 col_dat = np.array(trinorms[frame][i]) 

237 poss.append(col_dat[:, :3, :].flatten()) # 3 first elements are positions 

238 cols.append((np.zeros(col_dat.shape[0]*18).reshape(-1, 3) + color).flatten()) # 3 colors for each vertex 

239 nors.append(col_dat[:, 3:, :].flatten()) # 3 last elements are normals 

240 poss = np.concatenate(poss) 

241 cols = np.concatenate(cols) 

242 nors = np.concatenate(nors) 

243 # Create NGLView widget 

244 u = mda.Universe(input_top_path, input_traj_path) 

245 view = nv.show_mdanalysis(u) 

246 # view.clear_representations() 

247 mesh = ('mesh', poss, cols) 

248 view._add_shape([mesh], name='my_shape') 

249 view.update_representation(component=1, repr_index=0, opacity=opacity) 

250 return view 

251 

252 

253if __name__ == '__main__': 

254 main()