gromacs package

Submodules

gromacs.gmx_cluster module

Module containing the GMX Cluster class and the command line interface.

class gromacs.gmx_cluster.GMXCluster(input_structure_path, input_traj_path, output_pdb_path, input_index_path=None, output_cluster_log_path=None, output_rmsd_cluster_xpm_path=None, output_rmsd_dist_xvg_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXCluster
Wrapper of the GROMACS cluster module for clustering structures from a given GROMACS compatible trajectory.
GROMACS cluster can cluster structures using several different methods. Distances between structures can be determined from a trajectory. RMS deviation after fitting or RMS deviation of atom-pair distances can be used to define the distance between structures.
Parameters:
  • input_structure_path (str) – Path to the input structure file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_traj_path (str) –

    Path to the GROMACS trajectory file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_pdb_path (str) –

    Path to the output cluster file. File type: output. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • output_cluster_log_path (str) (Optional) –

    Path to the output log file. File type: output. Sample file. Accepted formats: log (edam:format_2330).

  • output_rmsd_cluster_xpm_path (str) (Optional) –

    Path to the output X PixMap compatible matrix file. File type: output. Sample file. Accepted formats: xpm (edam:format_3599).

  • output_rmsd_dist_xvg_path (str) (Optional) –

    Path to xvgr/xmgr file. File type: output. Sample file. Accepted formats: xvg (edam:format_2330).

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

    • fit_selection (str) - (“System”) Group where the fitting will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups).

    • output_selection (str) - (“System”) Group that is going to be written in the output trajectory. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups).

    • dista (bool) - (False) Use RMSD of distances instead of RMS deviation.

    • nofit (bool) - (False) Do not use least squares fitting before RMSD calculation.

    • method (str) - (“linkage”) Method for cluster determination. Values: linkage (Add a structure to a cluster when its distance to any element of the cluster is less than cutoff), jarvis-patrick (Add a structure to a cluster when this structure and a structure in the cluster have each other as neighbors and they have a least P neighbors in common), monte-carlo (Reorder the RMSD matrix using Monte Carlo such that the order of the frames is using the smallest possible increments), diagonalization (Diagonalize the RMSD matrix), gromos (Count number of neighbors using cut-off and take structure with largest number of neighbors with all its neighbors as cluster and eliminate it from the pool of clusters).

    • cutoff (float) - (0.1) [0~10|0.1] RMSD cut-off (nm) for two structures to be neighbor.

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_cluster import gmx_cluster
prop = {
    'fit_selection': 'System',
    'output_selection': 'System',
    'method': 'linkage'
}
gmx_cluster(input_structure_path='/path/to/myStructure.tpr',
            input_traj_path='/path/to/myTrajectory.trr',
            input_index_path='/path/to/myIndex.ndx',
            output_pdb_path='/path/to/newStructure.pdb',
            properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXCluster gromacs.gmx_cluster.GMXCluster object.

gromacs.gmx_cluster.gmx_cluster(input_structure_path: str, input_traj_path: str, output_pdb_path: str, input_index_path: str | None = None, output_cluster_log_path: str | None = None, output_rmsd_cluster_xpm_path: str | None = None, output_rmsd_dist_xvg_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXCluster class and execute the launch() method.

gromacs.gmx_cluster.main()[source]

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

gromacs.gmx_rms module

Module containing the GMX Rms class and the command line interface.

class gromacs.gmx_rms.GMXRms(input_structure_path, input_traj_path, output_xvg_path, input_index_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXRms
Wrapper of the GROMACS rms module for performing a Root Mean Square deviation (RMSd) analysis from a given GROMACS compatible trajectory.
GROMACS rms compares two structures by computing the root mean square deviation (RMSD), the size-independent rho similarity parameter (rho) or the scaled rho (rhosc).
Parameters:
  • input_structure_path (str) –

    Path to the input structure file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_traj_path (str) –

    Path to the GROMACS trajectory file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_xvg_path (str) –

    Path to the XVG output file. File type: output. Sample file. Accepted formats: xvg (edam:format_2030).

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

    • xvg (str) - (“none”) XVG plot formatting. Values: xmgrace, xmgr, none.

    • selection (str) - (“System”) Group where the rms will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_rms import gmx_rms
prop = {
    'xvg': 'xmgr',
    'selection': 'Water_and_ions'
}
gmx_rms(input_structure_path='/path/to/myStructure.tpr',
        input_traj_path='/path/to/myTrajectory.trr',
        output_xvg_path='/path/to/newXVG.xvg',
        input_index_path='/path/to/myIndex.ndx',
        properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXRms gromacs.gmx_rms.GMXRms object.

gromacs.gmx_rms.gmx_rms(input_structure_path: str, input_traj_path: str, output_xvg_path: str, input_index_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXRms class and execute the launch() method.

gromacs.gmx_rms.main()[source]

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

gromacs.gmx_energy module

Module containing the GMX Energy class and the command line interface.

class gromacs.gmx_energy.GMXEnergy(input_energy_path, output_xvg_path, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXEnergy
Wrapper of the GROMACS energy module for extracting energy components from a given GROMACS energy file.
GROMACS energy extracts energy components from an energy file. The user is prompted to interactively select the desired energy terms.
Parameters:
  • input_energy_path (str) –

    Path to the input EDR file. File type: input. Sample file. Accepted formats: edr (edam:format_2330).

  • output_xvg_path (str) –

    Path to the XVG output file. File type: output. Sample file. Accepted formats: xvg (edam:format_2030).

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

    • xvg (str) - (“none”) XVG plot formatting. Values: xmgrace, xmgr, none.

    • terms (list) - ([“Potential”]) Energy terms. Values: Angle, Proper-Dih., Improper-Dih., LJ-14, Coulomb-14, LJ-(SR), Coulomb-(SR), Coul.-recip., Position-Rest., Potential, Kinetic-En., Total-Energy, Temperature, Pressure, Constr.-rmsd, Box-X, Box-Y, Box-Z, Volume, Density, pV, Enthalpy, Vir-XX, Vir-XY, Vir-XZ, Vir-YX, Vir-YY, Vir-YZ, Vir-ZX, Vir-ZY, Vir-ZZ, Pres-XX, Pres-XY, Pres-XZ, Pres-YX, Pres-YY, Pres-YZ, Pres-ZX, Pres-ZY, Pres-ZZ, #Surf*SurfTen, Box-Vel-XX, Box-Vel-YY, Box-Vel-ZZ, Mu-X, Mu-Y, Mu-Z, T-Protein, T-non-Protein, Lamb-Protein, Lamb-non-Protein.

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_energy import gmx_energy
prop = {
    'xvg': 'xmgr',
    'terms': ['Potential', 'Pressure']
}
gmx_energy(input_energy_path='/path/to/myEnergyFile.edr',
            output_xvg_path='/path/to/newXVG.xvg',
            properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

create_instructions_file()[source]

Creates an input file using the properties file settings

launch() int[source]

Execute the GMXEnergy gromacs.gmx_energy.GMXEnergy object.

gromacs.gmx_energy.gmx_energy(input_energy_path: str, output_xvg_path: str, properties: dict | None = None, **kwargs) int[source]

Execute the GMXEnergy class and execute the launch() method.

gromacs.gmx_energy.main()[source]

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

gromacs.gmx_rgyr module

Module containing the GMX Rgyr class and the command line interface.

class gromacs.gmx_rgyr.GMXRgyr(input_structure_path, input_traj_path, output_xvg_path, input_index_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXRgyr
Wrapper of the GROMACS gyrate module for computing the radius of gyration (Rgyr) of a molecule about the x-, y- and z-axes, as a function of time, from a given GROMACS compatible trajectory.
GROMACS gyrate computes the radius of gyration of a molecule and the radii of gyration about the x-, y- and z-axes, as a function of time. The atoms are explicitly mass weighted.
Parameters:
  • input_structure_path (str) –

    Path to the input structure file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_traj_path (str) –

    Path to the GROMACS trajectory file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_xvg_path (str) –

    Path to the XVG output file. File type: output. Sample file. Accepted formats: xvg (edam:format_2030).

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

    • xvg (str) - (“none”) XVG plot formatting. Values: xmgrace, xmgr, none.

    • selection (str) - (“System”) Group where the rgyr will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_rgyr import gmx_rgyr
prop = {
    'xvg': 'xmgr',
    'selection': 'Water_and_ions'
}
gmx_rgyr(input_structure_path='/path/to/myStructure.tpr',
        input_traj_path='/path/to/myTrajectory.trr',
        output_xvg_path='/path/to/newXVG.xvg',
        input_index_path='/path/to/myIndex.ndx',
        properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXRgyr gromacs.gmx_rgyr.GMXRgyr object.

gromacs.gmx_rgyr.gmx_rgyr(input_structure_path: str, input_traj_path: str, output_xvg_path: str, input_index_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXRgyr class and execute the launch() method.

gromacs.gmx_rgyr.main()[source]

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

gromacs.gmx_image module

Module containing the GMX TrjConvStr class and the command line interface.

class gromacs.gmx_image.GMXImage(input_traj_path, input_top_path, output_traj_path, input_index_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXImage
Wrapper of the GROMACS trjconv module for correcting periodicity (image) from a given GROMACS compatible trajectory file.
GROMACS trjconv module can convert trajectory files in many ways. See the GROMACS trjconv official documentation for further information.
Parameters:
  • input_traj_path (str) –

    Path to the GROMACS trajectory file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_top_path (str) –

    Path to the GROMACS input topology file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_traj_path (str) –

    Path to the output file. File type: output. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

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

    • fit_selection (str) - (“System”) Group where the fitting will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • center_selection (str) - (“System”) Group where the trjconv will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • cluster_selection (str) - (“System”) Group assigned to be the cluster, onto which all atoms are wrapped around the box, such that they are closest to the center of mass of the cluster, which is iteratively updated. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • output_selection (str) - (“System”) Group that is going to be written in the output trajectory. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • pbc (str) - (“mol”) PBC treatment (see help text for full description). Values: none (No PBC treatment), mol (Puts the center of mass of molecules in the box), res (Puts the center of mass of residues in the box), atom (Puts all the atoms in the box), nojump (Checks if atoms jump across the box and then puts them back), cluster (Clusters all the atoms in the selected index such that they are all closest to the center of mass of the cluster which is iteratively updated), whole (Only makes broken molecules whole).

    • center (bool) - (True) Center atoms in box.

    • ur (str) - (“compact”) Unit-cell representation. Values: rect (It’s the ordinary brick shape), tric (It’s the triclinic unit cell), compact (Puts all atoms at the closest distance from the center of the box).

    • fit (str) - (“none”) Fit molecule to ref structure in the structure file. Values: none, rot+trans, rotxy+transxy, translation, transxy, progressive.

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_image import gmx_image
prop = {
    'fit_selection': 'System',
    'center_selection': 'Water_and_ions',
    'output_selection': 'System',
    'pbc': 'mol'
}
gmx_image(input_traj_path='/path/to/myTrajectory.trr',
            input_top_path='/path/to/myTopology.tpr',
            output_traj_path='/path/to/newTrajectory.xtc',
            input_index_path='/path/to/myIndex.ndx',
            properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXImage gromacs.gmx_image.GMXImage object.

gromacs.gmx_image.gmx_image(input_traj_path: str, input_top_path: str, output_traj_path: str, input_index_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXImage class and execute the launch() method.

gromacs.gmx_image.main()[source]

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

gromacs.gmx_trjconv_str module

Module containing the GMX TrjConvStr class and the command line interface.

class gromacs.gmx_trjconv_str.GMXTrjConvStr(input_structure_path, input_top_path, output_str_path, input_index_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXTrjConvStr
Wrapper of the GROMACS trjconv module for converting between GROMACS compatible structure file formats and/or extracting a selection of atoms.
GROMACS trjconv module can convert trajectory files in many ways. See the GROMACS trjconv official documentation for further information.
Parameters:
  • input_structure_path (str) –

    Path to the input structure file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_top_path (str) –

    Path to the GROMACS input topology file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_str_path (str) –

    Path to the output file. File type: output. Sample file. Accepted formats: pdb (edam:format_1476), xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), tng (edam:format_3876).

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

    • selection (str) - (“System”) Group where the trjconv will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • pbc (str) - (“mol”) PBC treatment (see help text for full description). Values: none (No PBC treatment), mol (Puts the center of mass of molecules in the box), res (Puts the center of mass of residues in the box), atom (Puts all the atoms in the box), nojump (Checks if atoms jump across the box and then puts them back), cluster (Clusters all the atoms in the selected index such that they are all closest to the center of mass of the cluster which is iteratively updated), whole (Only makes broken molecules whole).

    • center (bool) - (True) Center atoms in box.

    • ur (str) - (“compact”) Unit-cell representation. Values: rect (It’s the ordinary brick shape), tric (It’s the triclinic unit cell), compact (Puts all atoms at the closest distance from the center of the box).

    • fit (str) - (“none”) Fit molecule to ref structure in the structure file. Values: none, rot+trans, rotxy+transxy, translation, transxy, progressive.

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_trjconv_str import gmx_trjconv_str
prop = {
    'selection': 'System'
}
gmx_trjconv_str(input_structure_path='/path/to/myStructure.trr',
                input_top_path='/path/to/myTopology.tpr',
                output_str_path='/path/to/newStructure.pdb',
                input_index_path='/path/to/myIndex.ndx',
                properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXTrjConvStr gromacs.gmx_trjconv_str.GMXTrjConvStr object.

gromacs.gmx_trjconv_str.gmx_trjconv_str(input_structure_path: str, input_top_path: str, output_str_path: str, input_index_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXTrjConvStr class and execute the launch() method.

gromacs.gmx_trjconv_str.main()[source]

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

gromacs.gmx_trjconv_str_ens module

Module containing the GMX TrjConvStr class and the command line interface.

class gromacs.gmx_trjconv_str_ens.GMXTrjConvStrEns(input_traj_path, input_top_path, output_str_ens_path, input_index_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXTrjConvStrEns
Wrapper of the GROMACS trjconv module for extracting an ensemble of frames containing a selection of atoms from GROMACS compatible trajectory files.
GROMACS trjconv module can convert trajectory files in many ways. See the GROMACS trjconv official documentation for further information.
Parameters:
  • input_traj_path (str) –

    Path to the GROMACS trajectory file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_top_path (str) –

    Path to the GROMACS input topology file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_str_ens_path (str) –

    Path to the output file. File type: output. Sample file. Accepted formats: zip (edam:format_3987).

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

    • selection (str) - (“System”) Group where the trjconv will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • skip (int) - (1) [0~10000|1] Only write every nr-th frame.

    • start (int) - (0) [0~10000|1] Time of first frame to read from trajectory (default unit ps).

    • end (int) - (0) [0~10000|1] Time of last frame to read from trajectory (default unit ps).

    • dt (int) - (0) [0~10000|1] Only write frame when t MOD dt = first time (ps).

    • output_name (str) - (“output”) File name for ensemble of output files.

    • output_type (str) - (“pdb”) File type for ensemble of output files. Values: gro (Contains a molecular structure in Gromos87 format), g96 (Can be a GROMOS-96 initial/final configuration file or a coordinate trajectory file or a combination of both), pdb (Molecular structure files in the protein databank file format).

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_trjconv_str_ens import gmx_trjconv_str_ens
prop = {
    'selection': 'System',
    'start': 0,
    'end': 10,
    'dt': 1
}
gmx_trjconv_str_ens(input_traj_path='/path/to/myStructure.trr',
                    input_top_path='/path/to/myTopology.tpr',
                    output_str_ens_path='/path/to/newStructureEnsemble.zip',
                    input_index_path='/path/to/myIndex.ndx',
                    properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXTrjConvStrEns gromacs.gmx_trjconv_str_ens.GMXTrjConvStrEns object.

gromacs.gmx_trjconv_str_ens.gmx_trjconv_str_ens(input_traj_path: str, input_top_path: str, output_str_ens_path: str, input_index_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXTrjConvStrEns class and execute the launch() method.

gromacs.gmx_trjconv_str_ens.main()[source]

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

gromacs.gmx_trjconv_trj module

Module containing the GMX TrjConvStr class and the command line interface.

class gromacs.gmx_trjconv_trj.GMXTrjConvTrj(input_traj_path, output_traj_path, input_index_path=None, input_top_path=None, properties=None, **kwargs)[source]

Bases: BiobbObject

biobb_analysis GMXTrjConvTrj
Wrapper of the GROMACS trjconv module for converting between GROMACS compatible trajectory file formats and/or extracts a selection of atoms.
GROMACS trjconv module can convert trajectory files in many ways. See the GROMACS trjconv official documentation for further information.
Parameters:
  • input_traj_path (str) –

    Path to the GROMACS trajectory file. File type: input. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

  • input_top_path (str) (Optional) –

    Path to the GROMACS input topology file. File type: input. Sample file. Accepted formats: tpr (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), brk (edam:format_2033), ent (edam:format_1476).

  • input_index_path (str) (Optional) –

    Path to the GROMACS index file. File type: input. Sample file. Accepted formats: ndx (edam:format_2033).

  • output_traj_path (str) –

    Path to the output file. File type: output. Sample file. Accepted formats: xtc (edam:format_3875), trr (edam:format_3910), cpt (edam:format_2333), gro (edam:format_2033), g96 (edam:format_2033), pdb (edam:format_1476), tng (edam:format_3876).

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

    • selection (str) - (“System”) Group where the trjconv will be performed. If input_index_path provided, check the file for the accepted values. Values: System (all atoms in the system), Protein (all protein atoms), Protein-H (protein atoms excluding hydrogens), C-alpha (C-alpha atoms), Backbone (protein backbone atoms: N; C-alpha and C), MainChain (protein main chain atoms: N; C-alpha; C and O; including oxygens in C-terminus), MainChain+Cb (protein main chain atoms including C-beta), MainChain+H (protein main chain atoms including backbone amide hydrogens and hydrogens on the N-terminus), SideChain (protein side chain atoms: that is all atoms except N; C-alpha; C; O; backbone amide hydrogens and oxygens in C-terminus and hydrogens on the N-terminus), SideChain-H (protein side chain atoms excluding all hydrogens), Prot-Masses (protein atoms excluding dummy masses), non-Protein (all non-protein atoms), Water (water molecules), SOL (water molecules), non-Water (anything not covered by the Water group), Ion (any name matching an Ion entry in residuetypes.dat), NA (all NA atoms), CL (all CL atoms), Water_and_ions (combination of the Water and Ions groups), DNA (all DNA atoms), RNA (all RNA atoms), Protein_DNA (all Protein-DNA complex atoms), Protein_RNA (all Protein-RNA complex atoms), Protein_DNA_RNA (all Protein-DNA-RNA complex atoms), DNA_RNA (all DNA-RNA complex atoms).

    • start (int) - (0) [0~10000|1] Time of first frame to read from trajectory (default unit ps).

    • end (int) - (0) [0~10000|1] Time of last frame to read from trajectory (default unit ps).

    • dt (int) - (0) [0~10000|1] Only write frame when t MOD dt = first time (ps).

    • binary_path (str) - (“gmx”) Path to the GROMACS executable binary.

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

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

    • container_path (str) - (None) Container path definition.

    • container_image (str) - (‘gromacs/gromacs:2022.2’) Container image definition.

    • container_volume_path (str) - (‘/tmp’) Container volume path definition.

    • container_working_dir (str) - (None) Container working directory definition.

    • container_user_id (str) - (None) Container user_id definition.

    • container_shell_path (str) - (‘/bin/bash’) Path to default shell inside the container.

Examples

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

from biobb_analysis.gromacs.gmx_trjconv_trj import gmx_trjconv_trj
prop = {
    'selection': 'System',
    'start': 0,
    'end': 0
}
gmx_trjconv_trj(input_traj_path='/path/to/myStructure.trr',
                output_traj_path='/path/to/newTrajectory.xtc',
                input_index_path='/path/to/myIndex.ndx',
                properties=prop)
Info:
check_data_params(out_log, err_log)[source]

Checks all the input/output paths and parameters

launch() int[source]

Execute the GMXTrjConvTrj gromacs.gmx_trjconv_trj.GMXTrjConvTrj object.

gromacs.gmx_trjconv_trj.gmx_trjconv_trj(input_traj_path: str, output_traj_path: str, input_index_path: str | None = None, input_top_path: str | None = None, properties: dict | None = None, **kwargs) int[source]

Execute the GMXTrjConvTrj class and execute the launch() method.

gromacs.gmx_trjconv_trj.main()[source]

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