Coverage for biobb_mem/lipyphilic_biobb/lpp_zpositions.py: 88%

56 statements  

« prev     ^ index     » next       coverage.py v7.10.7, created at 2025-10-03 16:03 +0000

1#!/usr/bin/env python3 

2 

3"""Module containing the Lipyphilic ZPositions class and the command line interface.""" 

4 

5from biobb_common.generic.biobb_object import BiobbObject 

6from biobb_common.tools.file_utils import launchlogger 

7import MDAnalysis as mda 

8from biobb_mem.lipyphilic_biobb.common import ignore_no_box 

9from lipyphilic.lib.z_positions import ZPositions 

10import pandas as pd 

11import numpy as np 

12 

13 

14class LPPZPositions(BiobbObject): 

15 """ 

16 | biobb_mem LPPZPositions 

17 | Wrapper of the LiPyphilic ZPositions module for calculating the z distance of lipids to the bilayer center. 

18 | LiPyphilic is a Python package for analyzing MD simulations of lipid bilayers. The parameter names and defaults are the same as the ones in the official `Lipyphilic documentation <https://lipyphilic.readthedocs.io/en/latest/reference/analysis/z_pos.html>`_. 

19 

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_positions_path (str): Path to the output z positions. File type: output. `Sample file <https://github.com/bioexcel/biobb_mem/raw/main/biobb_mem/test/reference/lipyphilic_biobb/zpositions.csv>`_. Accepted formats: csv (edam:format_3752). 

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

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

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

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

28 * **lipid_sel** (*str*) - ("all") Selection string for the lipids in a membrane. The selection should cover **all** residues in the membrane, including cholesterol. 

29 * **height_sel** (*str*) - ("all") Atom selection for the molecules for which the z position will be calculated. 

30 * **n_bins** (*int*) - (1) Number of bins in *x* and *y* to use to create a grid of membrane patches. Local membrane midpoints are computed for each patch, and lipids assigned a leaflet based on the distance to their local membrane midpoint. The default is `1`, which is equivalent to computing a single global midpoint. 

31 * **ignore_no_box** (*bool*) - (False) Ignore the absence of box information in the trajectory. If the trajectory does not contain box information, the box will be set to the minimum and maximum positions of the atoms in the trajectory. 

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

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

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

35 

36 Examples: 

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

38 

39 from biobb_mem.lipyphilic_biobb.lpp_zpositions import lpp_zpositions 

40 prop = { 

41 'lipid_sel': 'name GL1 GL2 ROH', 

42 } 

43 lpp_zpositions(input_top_path='/path/to/myTopology.tpr', 

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

45 output_positions_path='/path/to/zpositions.csv.csv', 

46 properties=prop) 

47 

48 Info: 

49 * wrapped_software: 

50 * name: LiPyphilic 

51 * version: 0.10.0 

52 * license: GPL-2.0 

53 * ontology: 

54 * name: EDAM 

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

56 

57 """ 

58 

59 def __init__(self, input_top_path, input_traj_path, output_positions_path, 

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

61 properties = properties or {} 

62 

63 # Call parent class constructor 

64 super().__init__(properties) 

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

66 

67 # Input/Output files 

68 self.io_dict = { 

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

70 "out": {"output_positions_path": output_positions_path} 

71 } 

72 

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

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

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

76 self.lipid_sel = properties.get('lipid_sel', 'all') 

77 self.height_sel = properties.get('height_sel', 'all') 

78 self.n_bins = properties.get('n_bins', 1) 

79 self.ignore_no_box = properties.get('ignore_no_box', False) 

80 # Properties specific for BB 

81 self.remove_tmp = properties.get('remove_tmp', True) 

82 self.restart = properties.get('restart', False) 

83 self.sandbox_path = properties.get('sandbox_path', "./") 

84 

85 # Check the properties 

86 self.check_properties(properties) 

87 self.check_arguments() 

88 

89 @launchlogger 

90 def launch(self) -> int: 

91 """Execute the :class:`LPPZPositions <lipyphilic_biobb.lpp_zpositions.LPPZPositions>` object.""" 

92 

93 # Setup Biobb 

94 if self.check_restart(): 

95 return 0 

96 self.stage_files() 

97 

98 # Load the trajectory 

99 u = mda.Universe(self.stage_io_dict["in"]["input_top_path"], self.stage_io_dict["in"]["input_traj_path"]) 

100 ignore_no_box(u, self.ignore_no_box, self.out_log, self.global_log) 

101 

102 # Create ZPositions object 

103 positions = ZPositions( 

104 universe=u, 

105 lipid_sel=self.lipid_sel, 

106 height_sel=self.height_sel, 

107 n_bins=self.n_bins 

108 ) 

109 # Run the analysis 

110 positions.run(start=self.start, stop=self.stop, step=self.steps) 

111 # Save the results 

112 frames = positions.z_positions.shape[1] 

113 resnames = np.repeat(positions._height_atoms.resnames, frames) 

114 resindices = np.tile(positions._height_atoms.resindices, frames) 

115 frame_numbers = np.repeat(np.arange(frames), positions._height_atoms.n_residues) 

116 

117 df = pd.DataFrame({ 

118 'resname': resnames, 

119 'resindex': resindices, 

120 'frame': frame_numbers, 

121 'zposition': positions.z_positions.T.flatten() 

122 }) 

123 

124 # Save the DataFrame to a CSV file 

125 df.to_csv(self.stage_io_dict["out"]["output_positions_path"], index=False) 

126 

127 # Copy files to host 

128 self.copy_to_host() 

129 self.remove_tmp_files() 

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

131 

132 return self.return_code 

133 

134 

135def lpp_zpositions(input_top_path: str, input_traj_path: str, output_positions_path: str = None, properties: dict = None, **kwargs) -> int: 

136 """Execute the :class:`LPPZPositions <lipyphilic_biobb.lpp_zpositions.LPPZPositions>` class and 

137 execute the :meth:`launch() <lipyphilic_biobb.lpp_zpositions.LPPZPositions.launch>` method.""" 

138 return LPPZPositions(**dict(locals())).launch() 

139 

140 

141lpp_zpositions.__doc__ = LPPZPositions.__doc__ 

142main = LPPZPositions.get_main(lpp_zpositions, "Calculate the z distance in of lipids to the bilayer center.") 

143 

144 

145def frame_df(output_positions_path): 

146 """ 

147 Processes a CSV file containing z-position data and calculates the mean positive, mean negative, 

148 thickness, and standard deviation of thickness for each frame. 

149 

150 Args: 

151 output_positions_path (str): Path to the CSV file containing z-position data. 

152 Returns: 

153 pandas.DataFrame: A DataFrame with the following columns: 

154 - mean_positive: Mean of positive z-positions for each frame. 

155 - mean_negative: Mean of negative z-positions for each frame. 

156 - thickness: Difference between mean_positive and mean_negative for each frame. 

157 - std_thickness: Standard deviation of the absolute z-positions for each frame. 

158 """ 

159 

160 df = pd.read_csv(output_positions_path) 

161 grouped = df.groupby('frame')['zposition'].agg( 

162 mean_positive=lambda x: x[x > 0].mean(), 

163 mean_negative=lambda x: x[x < 0].mean(), 

164 std_positive=lambda x: x[x > 0].std(), 

165 std_negative=lambda x: x[x < 0].std() 

166 ) 

167 grouped['thickness'] = grouped['mean_positive'] - grouped['mean_negative'] 

168 grouped['std_thickness'] = df.groupby('frame')['zposition'].apply(lambda x: x.abs().std()) 

169 return grouped 

170 

171 

172if __name__ == '__main__': 

173 main()