Coverage for biobb_gromacs/gromacs/grompp_mdrun.py: 73%

82 statements  

« prev     ^ index     » next       coverage.py v7.9.1, created at 2025-06-25 09:23 +0000

1#!/usr/bin/env python3 

2 

3"""Module containing the GromppMDrun class and the command line interface.""" 

4import argparse 

5from typing import Optional 

6from pathlib import Path 

7from biobb_common.generic.biobb_object import BiobbObject 

8from biobb_common.configuration import settings 

9from biobb_common.tools import file_utils as fu 

10from biobb_common.tools.file_utils import launchlogger 

11from biobb_gromacs.gromacs.grompp import grompp 

12from biobb_gromacs.gromacs.mdrun import mdrun 

13 

14 

15class GromppMdrun(BiobbObject): 

16 """ 

17 | biobb_gromacs GromppMdrun 

18 | Wrapper of the `GROMACS grompp <http://manual.gromacs.org/current/onlinehelp/gmx-grompp.html>`_ module and the `GROMACS mdrun <http://manual.gromacs.org/current/onlinehelp/gmx-mdrun.html>`_ module. 

19 | Grompp The GROMACS preprocessor module needs to be fed with the input system and the dynamics parameters to create a portable binary run input file TPR. The simulation parameters can be specified by two methods: 1.The predefined mdp settings defined at simulation_type property or 2.A custom mdp file defined at the input_mdp_path argument. These two methods are mutually exclusive. In both cases can be further modified by adding parameters to the mdp section in the yaml configuration file. The simulation parameter names and default values can be consulted in the `official MDP specification <http://manual.gromacs.org/current/user-guide/mdp-options.html>`_. MDRun is the main computational chemistry engine within GROMACS. It performs Molecular Dynamics simulations, but it can also perform Stochastic Dynamics, Energy Minimization, test particle insertion or (re)calculation of energies. 

20 

21 Args: 

22 input_gro_path (str): Path to the input GROMACS structure GRO file. File type: input. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/data/gromacs/grompp.gro>`_. Accepted formats: gro (edam:format_2033). 

23 input_top_zip_path (str): Path to the input GROMACS topology TOP and ITP files in zip format. File type: input. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/data/gromacs/grompp.zip>`_. Accepted formats: zip (edam:format_3987). 

24 output_trr_path (str): Path to the GROMACS uncompressed raw trajectory file TRR. File type: output. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/reference/gromacs/ref_mdrun.trr>`_. Accepted formats: trr (edam:format_3910). 

25 output_gro_path (str): Path to the output GROMACS structure GRO file. File type: output. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/reference/gromacs/ref_mdrun.gro>`_. Accepted formats: gro (edam:format_2033). 

26 output_edr_path (str): Path to the output GROMACS portable energy file EDR. File type: output. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/reference/gromacs/ref_mdrun.edr>`_. Accepted formats: edr (edam:format_2330). 

27 output_log_path (str): Path to the output GROMACS trajectory log file LOG. File type: output. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/reference/gromacs/ref_gmx_mdrun.log>`_. Accepted formats: log (edam:format_2330). 

28 input_cpt_path (str) (Optional): Path to the input GROMACS checkpoint file CPT. File type: input. Accepted formats: cpt (edam:format_2333). 

29 input_ndx_path (str) (Optional): Path to the input GROMACS index files NDX. File type: input. Accepted formats: ndx (edam:format_2033). 

30 input_mdp_path (str) (Optional): Path to the input GROMACS `MDP file <http://manual.gromacs.org/current/user-guide/mdp-options.html>`_. File type: input. Accepted formats: mdp (edam:format_2330). 

31 output_xtc_path (str) (Optional): Path to the GROMACS compressed trajectory file XTC. File type: output. Accepted formats: xtc (edam:format_3875). 

32 output_cpt_path (str) (Optional): Path to the output GROMACS checkpoint file CPT. File type: output. Accepted formats: cpt (edam:format_2333). 

33 output_dhdl_path (str) (Optional): Path to the output dhdl.xvg file only used when free energy calculation is turned on. File type: output. Accepted formats: xvg (edam:format_2033). 

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

35 * **mdp** (*dict*) - ({}) MDP options specification. 

36 * **simulation_type** (*str*) - ("minimization") Default options for the mdp file. Each creates a different mdp file. Values: `minimization <https://biobb-gromacs.readthedocs.io/en/latest/_static/mdp/minimization.mdp>`_ (Energy minimization using steepest descent algorithm is used), `nvt <https://biobb-gromacs.readthedocs.io/en/latest/_static/mdp/nvt.mdp>`_ (substance N Volume V and Temperature T are conserved), `npt <https://biobb-gromacs.readthedocs.io/en/latest/_static/mdp/npt.mdp>`_ (substance N pressure P and Temperature T are conserved), `free <https://biobb-gromacs.readthedocs.io/en/latest/_static/mdp/free.mdp>`_ (No design constraints applied; Free MD), `ions <https://biobb-gromacs.readthedocs.io/en/latest/_static/mdp/minimization.mdp>`_ (Synonym of minimization), index (Creates an empty mdp file). 

37 * **maxwarn** (*int*) - (10) [0~1000|1] Maximum number of allowed warnings. 

38 * **mpi_bin** (*str*) - (None) Path to the MPI runner. Usually "mpirun" or "srun". 

39 * **mpi_np** (*str*) - (None) Number of MPI processes. Usually an integer bigger than 1. 

40 * **mpi_hostlist** (*str*) - (None) Path to the MPI hostlist file. 

41 * **checkpoint_time** (*int*) - (15) [0~1000|1] Checkpoint writing interval in minutes. Only enabled if an output_cpt_path is provided. 

42 * **num_threads** (*int*) - (0) [0-1000|1] Let GROMACS guess. The number of threads that are going to be used. 

43 * **num_threads_mpi** (*int*) - (0) [0-1000|1] Let GROMACS guess. The number of GROMACS MPI threads that are going to be used. 

44 * **num_threads_omp** (*int*) - (0) [0-1000|1] Let GROMACS guess. The number of GROMACS OPENMP threads that are going to be used. 

45 * **num_threads_omp_pme** (*int*) - (0) [0-1000|1] Let GROMACS guess. The number of GROMACS OPENMP_PME threads that are going to be used. 

46 * **use_gpu** (*bool*) - (False) Use settings appropriate for GPU. Adds: -nb gpu -pme gpu 

47 * **gpu_id** (*str*) - (None) list of unique GPU device IDs available to use. 

48 * **gpu_tasks** (*str*) - (None) list of GPU device IDs, mapping each PP task on each node to a device. 

49 * **gmx_lib** (*str*) - (None) Path set GROMACS GMXLIB environment variable. 

50 * **binary_path** (*str*) - ("gmx") Path to the GROMACS executable binary. 

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

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

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

54 * **container_path** (*str*) - (None) Path to the binary executable of your container. 

55 * **container_image** (*str*) - ("gromacs/gromacs:latest") Container Image identifier. 

56 * **container_volume_path** (*str*) - ("/data") Path to an internal directory in the container. 

57 * **container_working_dir** (*str*) - (None) Path to the internal CWD in the container. 

58 * **container_user_id** (*str*) - (None) User number id to be mapped inside the container. 

59 * **container_shell_path** (*str*) - ("/bin/bash") Path to the binary executable of the container shell. 

60 

61 Examples: 

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

63 

64 from biobb_gromacs.gromacs.grompp_mdrun import grompp_mdrun 

65 prop = { 'num_threads': 0, 

66 'binary_path': 'gmx', 

67 'mdp': 

68 { 'simulation_type': 'minimization', 

69 'emtol':'500', 

70 'nsteps':'5000'} 

71 } 

72 grompp_mdrun(input_gro_path='/path/to/myStructure.gro', 

73 input_top_zip_path='/path/to/myTopology.zip', 

74 output_trr_path='/path/to/newTrajectory.trr', 

75 output_gro_path='/path/to/newStructure.gro', 

76 output_edr_path='/path/to/newEnergy.edr', 

77 output_log_path='/path/to/newSimulationLog.log', 

78 properties=prop) 

79 

80 Info: 

81 * wrapped_software: 

82 * name: GROMACS Grompp & MDRun 

83 * version: 2024.5 

84 * license: LGPL 2.1 

85 * ontology: 

86 * name: EDAM 

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

88 """ 

89 

90 def __init__(self, input_gro_path: str, input_top_zip_path: str, output_trr_path: str, 

91 output_gro_path: str, output_edr_path: str, output_log_path: str, 

92 input_cpt_path: Optional[str] = None, input_ndx_path: Optional[str] = None, input_mdp_path: Optional[str] = None, 

93 output_xtc_path: Optional[str] = None, output_cpt_path: Optional[str] = None, output_dhdl_path: Optional[str] = None, 

94 output_tpr_path: Optional[str] = None, properties: Optional[dict] = None, **kwargs) -> None: 

95 # Properties management 

96 properties = properties or {} 

97 

98 # Call parent class constructor 

99 super().__init__(properties) 

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

101 

102 grompp_properties_keys = ['mdp', 'maxwarn', 'simulation_type'] 

103 mdrun_properties_keys = ['mpi_bin', 'mpi_np', 'mpi_hostlist', 'checkpoint_time', 'num_threads', 'num_threads_mpi', 'num_threads_omp', 'num_threads_omp_pme', 'use_gpu', 'gpu_id', 'gpu_tasks', 'dev'] 

104 self.properties_grompp = {} 

105 self.properties_mdrun = {} 

106 if properties: 

107 self.global_log = properties.get('global_log', None) 

108 self.properties_grompp = properties.copy() 

109 for key in mdrun_properties_keys: 

110 self.properties_grompp.pop(key, None) 

111 self.properties_mdrun = properties.copy() 

112 for key in grompp_properties_keys: 

113 self.properties_mdrun.pop(key, None) 

114 

115 # Grompp arguments 

116 self.input_gro_path = input_gro_path 

117 self.input_top_zip_path = input_top_zip_path 

118 self.output_tpr_path = output_tpr_path 

119 if not self.output_tpr_path: 

120 self.output_tpr_path = str(Path(fu.create_unique_dir()).joinpath('internal.tpr')) 

121 self.input_cpt_path = input_cpt_path 

122 self.input_ndx_path = input_ndx_path 

123 self.input_mdp_path = input_mdp_path 

124 

125 # MDRun arguments 

126 self.input_tpr_path = self.output_tpr_path 

127 self.output_trr_path = output_trr_path 

128 self.output_gro_path = output_gro_path 

129 self.output_edr_path = output_edr_path 

130 self.output_log_path = output_log_path 

131 self.output_xtc_path = output_xtc_path 

132 self.output_cpt_path = output_cpt_path 

133 self.output_dhdl_path = output_dhdl_path 

134 

135 @launchlogger 

136 def launch(self) -> int: 

137 """Execute the :class:`GromppMdrun <gromacs.grompp_mdrun.GromppMdrun>` object.""" 

138 

139 fu.log('Calling Grompp class', self.out_log, self.global_log) 

140 grompp_return_code = grompp(input_gro_path=self.input_gro_path, input_top_zip_path=self.input_top_zip_path, 

141 output_tpr_path=self.output_tpr_path, input_cpt_path=self.input_cpt_path, # type: ignore 

142 input_ndx_path=self.input_ndx_path, input_mdp_path=self.input_mdp_path, 

143 properties=self.properties_grompp) 

144 fu.log(f'Grompp return code: {grompp_return_code}', self.out_log, self.global_log) 

145 

146 if not grompp_return_code: 

147 fu.log('Grompp return code is correct. Calling MDRun class', self.out_log, self.global_log) 

148 mdrun_return_code = mdrun(input_tpr_path=self.output_tpr_path, output_trr_path=self.output_trr_path, # type: ignore 

149 output_gro_path=self.output_gro_path, output_edr_path=self.output_edr_path, 

150 output_log_path=self.output_log_path, output_xtc_path=self.output_xtc_path, 

151 output_cpt_path=self.output_cpt_path, output_dhdl_path=self.output_dhdl_path, 

152 properties=self.properties_mdrun) 

153 fu.log(f'MDRun return code: {mdrun_return_code}', self.out_log, self.global_log) 

154 else: 

155 return 1 

156 

157 # Remove temporal files 

158 self.tmp_files.extend([ 

159 # self.stage_io_dict.get("unique_dir", ''), 

160 Path(str(self.output_tpr_path)).parent] 

161 ) 

162 self.remove_tmp_files() 

163 

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

165 return mdrun_return_code 

166 

167 

168def grompp_mdrun(input_gro_path: str, input_top_zip_path: str, output_trr_path: str, 

169 output_gro_path: str, output_edr_path: str, output_log_path: str, 

170 input_cpt_path: Optional[str] = None, input_ndx_path: Optional[str] = None, input_mdp_path: Optional[str] = None, 

171 output_xtc_path: Optional[str] = None, output_cpt_path: Optional[str] = None, output_dhdl_path: Optional[str] = None, 

172 output_tpr_path: Optional[str] = None, properties: Optional[dict] = None, **kwargs) -> int: 

173 return GromppMdrun(input_gro_path=input_gro_path, input_top_zip_path=input_top_zip_path, 

174 output_trr_path=output_trr_path, output_gro_path=output_gro_path, 

175 output_edr_path=output_edr_path, output_log_path=output_log_path, 

176 input_cpt_path=input_cpt_path, input_ndx_path=input_ndx_path, 

177 input_mdp_path=input_mdp_path, output_xtc_path=output_xtc_path, 

178 output_cpt_path=output_cpt_path, output_dhdl_path=output_dhdl_path, 

179 output_tpr_path=output_tpr_path, properties=properties, **kwargs).launch() 

180 

181 

182grompp_mdrun.__doc__ = GromppMdrun.__doc__ 

183 

184 

185def main(): 

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

187 parser = argparse.ArgumentParser(description="Wrapper for the GROMACS grompp_mdrun module.", 

188 formatter_class=lambda prog: argparse.RawTextHelpFormatter(prog, width=99999)) 

189 parser.add_argument('-c', '--config', required=False, help="This file can be a YAML file, JSON file or JSON string") 

190 

191 # Specific args of each building block 

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

193 required_args.add_argument('--input_gro_path', required=True) 

194 required_args.add_argument('--input_top_zip_path', required=True) 

195 required_args.add_argument('--output_trr_path', required=True) 

196 required_args.add_argument('--output_gro_path', required=True) 

197 required_args.add_argument('--output_edr_path', required=True) 

198 required_args.add_argument('--output_log_path', required=True) 

199 parser.add_argument('--input_cpt_path', required=False) 

200 parser.add_argument('--input_ndx_path', required=False) 

201 parser.add_argument('--input_mdp_path', required=False) 

202 parser.add_argument('--output_xtc_path', required=False) 

203 parser.add_argument('--output_cpt_path', required=False) 

204 parser.add_argument('--output_dhdl_path', required=False) 

205 parser.add_argument('--output_tpr_path', required=False) 

206 

207 args = parser.parse_args() 

208 config = args.config if args.config else None 

209 properties = settings.ConfReader(config=config).get_prop_dic() 

210 

211 # Specific call of each building block 

212 grompp_mdrun(input_gro_path=args.input_gro_path, input_top_zip_path=args.input_top_zip_path, 

213 output_trr_path=args.output_trr_path, output_gro_path=args.output_gro_path, 

214 output_edr_path=args.output_edr_path, output_log_path=args.output_log_path, 

215 input_cpt_path=args.input_cpt_path, input_ndx_path=args.input_ndx_path, 

216 input_mdp_path=args.input_mdp_path, output_xtc_path=args.output_xtc_path, 

217 output_cpt_path=args.output_cpt_path, output_dhdl_path=args.output_dhdl_path, 

218 output_tpr_path=args.output_tpr_path, properties=properties) 

219 

220 

221if __name__ == '__main__': 

222 main()