Coverage for biobb_gromacs/gromacs/mdrun.py: 57%
125 statements
« prev ^ index » next coverage.py v7.9.1, created at 2025-06-25 09:23 +0000
« prev ^ index » next coverage.py v7.9.1, created at 2025-06-25 09:23 +0000
1#!/usr/bin/env python3
3"""Module containing the MDrun class and the command line interface."""
4import argparse
5from typing import Optional
6from biobb_common.generic.biobb_object import BiobbObject
7from biobb_common.configuration import settings
8from biobb_common.tools import file_utils as fu
9from biobb_common.tools.file_utils import launchlogger
10from biobb_gromacs.gromacs.common import get_gromacs_version
13class Mdrun(BiobbObject):
14 """
15 | biobb_gromacs Mdrun
16 | Wrapper of the `GROMACS mdrun <http://manual.gromacs.org/current/onlinehelp/gmx-mdrun.html>`_ module.
17 | 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.
19 Args:
20 input_tpr_path (str): Path to the portable binary run input file TPR. File type: input. `Sample file <https://github.com/bioexcel/biobb_gromacs/raw/master/biobb_gromacs/test/data/gromacs/mdrun.tpr>`_. Accepted formats: tpr (edam:format_2333).
21 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).
22 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).
23 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_mdrun.log>`_. Accepted formats: log (edam:format_2330).
24 output_trr_path (str) (Optional): 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 input_cpt_path (str) (Optional): Path to the input GROMACS checkpoint file CPT. File type: input. Accepted formats: cpt (edam:format_2333).
26 output_xtc_path (str) (Optional): Path to the GROMACS compressed trajectory file XTC. File type: output. Accepted formats: xtc (edam:format_3875).
27 output_cpt_path (str) (Optional): Path to the output GROMACS checkpoint file CPT. File type: output. Accepted formats: cpt (edam:format_2333).
28 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).
29 properties (dict - Python dictionary object containing the tool parameters, not input/output files):
30 * **mpi_bin** (*str*) - (None) Path to the MPI runner. Usually "mpirun" or "srun".
31 * **mpi_np** (*int*) - (0) [0~1000|1] Number of MPI processes. Usually an integer bigger than 1.
32 * **mpi_flags** (*str*) - (None) Path to the MPI hostlist file.
33 * **checkpoint_time** (*int*) - (15) [0~1000|1] Checkpoint writing interval in minutes. Only enabled if an output_cpt_path is provided.
34 * **num_threads** (*int*) - (0) [0~1000|1] Let GROMACS guess. The number of threads that are going to be used.
35 * **num_threads_mpi** (*int*) - (0) [0~1000|1] Let GROMACS guess. The number of GROMACS MPI threads that are going to be used.
36 * **num_threads_omp** (*int*) - (0) [0~1000|1] Let GROMACS guess. The number of GROMACS OPENMP threads that are going to be used.
37 * **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.
38 * **use_gpu** (*bool*) - (False) Use settings appropriate for GPU. Adds: -nb gpu -pme gpu
39 * **gpu_id** (*str*) - (None) list of unique GPU device IDs available to use.
40 * **gpu_tasks** (*str*) - (None) list of GPU device IDs, mapping each PP task on each node to a device.
41 * **gmx_lib** (*str*) - (None) Path set GROMACS GMXLIB environment variable.
42 * **binary_path** (*str*) - ("gmx") Path to the GROMACS executable binary.
43 * **remove_tmp** (*bool*) - (True) [WF property] Remove temporal files.
44 * **restart** (*bool*) - (False) [WF property] Do not execute if output files exist.
45 * **sandbox_path** (*str*) - ("./") [WF property] Parent path to the sandbox directory.
46 * **container_path** (*str*) - (None) Path to the binary executable of your container.
47 * **container_image** (*str*) - (None) Container Image identifier.
48 * **container_volume_path** (*str*) - ("/data") Path to an internal directory in the container.
49 * **container_working_dir** (*str*) - (None) Path to the internal CWD in the container.
50 * **container_user_id** (*str*) - (None) User number id to be mapped inside the container.
51 * **container_shell_path** (*str*) - ("/bin/bash") Path to the binary executable of the container shell.
53 Examples:
54 This is a use example of how to use the building block from Python::
56 from biobb_gromacs.gromacs.mdrun import mdrun
57 prop = { 'num_threads': 0,
58 'binary_path': 'gmx' }
59 mdrun(input_tpr_path='/path/to/myPortableBinaryRunInputFile.tpr',
60 output_trr_path='/path/to/newTrajectory.trr',
61 output_gro_path='/path/to/newStructure.gro',
62 output_edr_path='/path/to/newEnergy.edr',
63 output_log_path='/path/to/newSimulationLog.log',
64 properties=prop)
66 Info:
67 * wrapped_software:
68 * name: GROMACS Mdrun
69 * version: 2024.5
70 * license: LGPL 2.1
71 * multinode: mpi
72 * ontology:
73 * name: EDAM
74 * schema: http://edamontology.org/EDAM.owl
75 """
77 def __init__(self, input_tpr_path: str, output_gro_path: str, output_edr_path: str,
78 output_log_path: str, output_trr_path: Optional[str] = None, input_cpt_path: Optional[str] = None,
79 output_xtc_path: Optional[str] = None, output_cpt_path: Optional[str] = None,
80 output_dhdl_path: Optional[str] = None, properties: Optional[dict] = None, **kwargs) -> None:
81 properties = properties or {}
83 # Call parent class constructor
84 super().__init__(properties)
85 self.locals_var_dict = locals().copy()
87 # Input/Output files
88 self.io_dict = {
89 "in": {"input_tpr_path": input_tpr_path, "input_cpt_path": input_cpt_path},
90 "out": {"output_trr_path": output_trr_path, "output_gro_path": output_gro_path,
91 "output_edr_path": output_edr_path, "output_log_path": output_log_path,
92 "output_xtc_path": output_xtc_path, "output_cpt_path": output_cpt_path,
93 "output_dhdl_path": output_dhdl_path}
94 }
96 # Properties specific for BB
97 # general mpi properties
98 self.mpi_bin = properties.get('mpi_bin')
99 self.mpi_np = properties.get('mpi_np')
100 self.mpi_flags = properties.get('mpi_flags')
101 # gromacs cpu mpi/openmp properties
102 self.num_threads = str(properties.get('num_threads', ''))
103 self.num_threads_mpi = str(properties.get('num_threads_mpi', ''))
104 self.num_threads_omp = str(properties.get('num_threads_omp', ''))
105 self.num_threads_omp_pme = str(properties.get('num_threads_omp_pme', ''))
106 # gromacs gpus
107 self.use_gpu = properties.get('use_gpu', False) # Adds: -nb gpu -pme gpu
108 self.gpu_id = str(properties.get('gpu_id', ''))
109 self.gpu_tasks = str(properties.get('gpu_tasks', ''))
110 # gromacs
111 self.checkpoint_time = properties.get('checkpoint_time')
113 # Properties common in all GROMACS BB
114 self.gmx_lib = properties.get('gmx_lib', None)
115 self.binary_path: str = properties.get('binary_path', 'gmx')
116 self.gmx_nobackup = properties.get('gmx_nobackup', True)
117 self.gmx_nocopyright = properties.get('gmx_nocopyright', True)
118 if self.gmx_nobackup:
119 self.binary_path += ' -nobackup'
120 if self.gmx_nocopyright:
121 self.binary_path += ' -nocopyright'
122 if (not self.mpi_bin) and (not self.container_path):
123 self.gmx_version = get_gromacs_version(self.binary_path)
125 # Check the properties
126 self.check_properties(properties)
127 self.check_arguments()
129 @launchlogger
130 def launch(self) -> int:
131 """Execute the :class:`Mdrun <gromacs.mdrun.Mdrun>` object."""
133 # Setup Biobb
134 if self.check_restart():
135 return 0
137 # Optional output files (if not added mrun will create them using a generic name)
138 if not self.stage_io_dict["out"].get("output_trr_path"):
139 self.stage_io_dict["out"]["output_trr_path"] = fu.create_name(prefix=self.prefix, step=self.step, name='trajectory.trr')
140 self.tmp_files.append(self.stage_io_dict["out"]["output_trr_path"])
142 self.stage_files()
144 self.cmd = [self.binary_path, 'mdrun',
145 '-o', self.stage_io_dict["out"]["output_trr_path"],
146 '-s', self.stage_io_dict["in"]["input_tpr_path"],
147 '-c', self.stage_io_dict["out"]["output_gro_path"],
148 '-e', self.stage_io_dict["out"]["output_edr_path"],
149 '-g', self.stage_io_dict["out"]["output_log_path"]]
151 if self.stage_io_dict["in"].get("input_cpt_path"):
152 self.cmd.append('-cpi')
153 self.cmd.append(self.stage_io_dict["in"]["input_cpt_path"])
154 if self.stage_io_dict["out"].get("output_xtc_path"):
155 self.cmd.append('-x')
156 self.cmd.append(self.stage_io_dict["out"]["output_xtc_path"])
157 else:
158 self.tmp_files.append('traj_comp.xtc')
159 if self.stage_io_dict["out"].get("output_cpt_path"):
160 self.cmd.append('-cpo')
161 self.cmd.append(self.stage_io_dict["out"]["output_cpt_path"])
162 if self.checkpoint_time:
163 self.cmd.append('-cpt')
164 self.cmd.append(str(self.checkpoint_time))
165 if self.stage_io_dict["out"].get("output_dhdl_path"):
166 self.cmd.append('-dhdl')
167 self.cmd.append(self.stage_io_dict["out"]["output_dhdl_path"])
169 # general mpi properties
170 if self.mpi_bin:
171 mpi_cmd = [self.mpi_bin]
172 if self.mpi_np:
173 mpi_cmd.append('-n')
174 mpi_cmd.append(str(self.mpi_np))
175 if self.mpi_flags:
176 mpi_cmd.extend(self.mpi_flags)
177 self.cmd = mpi_cmd + self.cmd
179 # gromacs cpu mpi/openmp properties
180 if self.num_threads:
181 fu.log(f'User added number of gmx threads: {self.num_threads}', self.out_log)
182 self.cmd.append('-nt')
183 self.cmd.append(self.num_threads)
184 if self.num_threads_mpi:
185 fu.log(f'User added number of gmx mpi threads: {self.num_threads_mpi}', self.out_log)
186 self.cmd.append('-ntmpi')
187 self.cmd.append(self.num_threads_mpi)
188 if self.num_threads_omp:
189 fu.log(f'User added number of gmx omp threads: {self.num_threads_omp}', self.out_log)
190 self.cmd.append('-ntomp')
191 self.cmd.append(self.num_threads_omp)
192 if self.num_threads_omp_pme:
193 fu.log(f'User added number of gmx omp_pme threads: {self.num_threads_omp_pme}', self.out_log)
194 self.cmd.append('-ntomp_pme')
195 self.cmd.append(self.num_threads_omp_pme)
196 # GMX gpu properties
197 if self.use_gpu:
198 fu.log('Adding GPU specific settings adds: -nb gpu -pme gpu', self.out_log)
199 self.cmd += ["-nb", "gpu", "-pme", "gpu"]
200 if self.gpu_id:
201 fu.log(f'list of unique GPU device IDs available to use: {self.gpu_id}', self.out_log)
202 self.cmd.append('-gpu_id')
203 self.cmd.append(self.gpu_id)
204 if self.gpu_tasks:
205 fu.log(f'list of GPU device IDs, mapping each PP task on each node to a device: {self.gpu_tasks}', self.out_log)
206 self.cmd.append('-gputasks')
207 self.cmd.append(self.gpu_tasks)
209 if self.gmx_lib:
210 self.env_vars_dict['GMXLIB'] = self.gmx_lib
212 # Run Biobb block
213 self.run_biobb()
215 # Copy files to host
216 self.copy_to_host()
218 # Remove temporal files
219 # self.tmp_files.append(self.stage_io_dict.get("unique_dir", ""))
220 self.remove_tmp_files()
222 self.check_arguments(output_files_created=True, raise_exception=False)
223 return self.return_code
226def mdrun(input_tpr_path: str, output_gro_path: str, output_edr_path: str,
227 output_log_path: str, output_trr_path: Optional[str] = None, input_cpt_path: Optional[str] = None,
228 output_xtc_path: Optional[str] = None, output_cpt_path: Optional[str] = None,
229 output_dhdl_path: Optional[str] = None, properties: Optional[dict] = None, **kwargs) -> int:
230 """Create :class:`Mdrun <gromacs.mdrun.Mdrun>` class and
231 execute the :meth:`launch() <gromacs.mdrun.Mdrun.launch>` method."""
233 return Mdrun(input_tpr_path=input_tpr_path, output_trr_path=output_trr_path,
234 output_gro_path=output_gro_path, output_edr_path=output_edr_path,
235 output_log_path=output_log_path, input_cpt_path=input_cpt_path,
236 output_xtc_path=output_xtc_path, output_cpt_path=output_cpt_path,
237 output_dhdl_path=output_dhdl_path, properties=properties,
238 **kwargs).launch()
241mdrun.__doc__ = Mdrun.__doc__
244def main():
245 """Command line execution of this building block. Please check the command line documentation."""
246 parser = argparse.ArgumentParser(description="Wrapper for the GROMACS mdrun module.",
247 formatter_class=lambda prog: argparse.RawTextHelpFormatter(prog, width=99999))
248 parser.add_argument('-c', '--config', required=False, help="This file can be a YAML file, JSON file or JSON string")
250 # Specific args of each building block
251 required_args = parser.add_argument_group('required arguments')
252 required_args.add_argument('--input_tpr_path', required=True)
253 required_args.add_argument('--output_gro_path', required=True)
254 required_args.add_argument('--output_edr_path', required=True)
255 required_args.add_argument('--output_log_path', required=True)
256 parser.add_argument('--output_trr_path', required=False)
257 parser.add_argument('--input_cpt_path', required=False)
258 parser.add_argument('--output_xtc_path', required=False)
259 parser.add_argument('--output_cpt_path', required=False)
260 parser.add_argument('--output_dhdl_path', required=False)
262 args = parser.parse_args()
263 config = args.config if args.config else None
264 properties = settings.ConfReader(config=config).get_prop_dic()
266 # Specific call of each building block
267 mdrun(input_tpr_path=args.input_tpr_path, output_trr_path=args.output_trr_path,
268 output_gro_path=args.output_gro_path, output_edr_path=args.output_edr_path,
269 output_log_path=args.output_log_path, input_cpt_path=args.input_cpt_path,
270 output_xtc_path=args.output_xtc_path, output_cpt_path=args.output_cpt_path,
271 output_dhdl_path=args.output_dhdl_path, properties=properties)
274if __name__ == '__main__':
275 main()