pyrexMD.core

Hint

This module contains functions enabling interactive analyses. Its main parts are the iPlayer and iPlot classes, which allow the use of a trajectory viewer or a dynamic linking of the trajectory viewer and any 2D graph.

Example:

import MDAnalysis as mda
import pyrexMD.misc as misc
import pyrexMD.core as core
import pyrexMD.analysis.contacts as con

ref = mda.Universe("<pdb_file>")
mobile = mda.Universe("<tpr_file>", "<xtc_file>")

# show ref structure in trajectory viewer
tv = core.iPlayer(ref)
tv()

# check for formed bias contacts
bias = misc.read_file("path/to/bias/contacts", usecols=(0,1))
FRAMES, QBIAS, CM = con.get_Qbias(mobile, bias)

# interactive plot (ctrl-click into the plot to jump to frame)
ip = core.iPlot(u1, xdata=FRAMES, ydata=QBIAS, xlabel="frame", ylabel="Qbias")
ip()

Content:

pyrexMD.core.cp = [(0.12156862745098039, 0.4666666666666667, 0.7058823529411765), (1.0, 0.4980392156862745, 0.054901960784313725), (0.17254901960784313, 0.6274509803921569, 0.17254901960784313), (0.8392156862745098, 0.15294117647058825, 0.1568627450980392), (0.5803921568627451, 0.403921568627451, 0.7411764705882353), (0.5490196078431373, 0.33725490196078434, 0.29411764705882354), (0.8901960784313725, 0.4666666666666667, 0.7607843137254902), (0.4980392156862745, 0.4980392156862745, 0.4980392156862745), (0.7372549019607844, 0.7411764705882353, 0.13333333333333333), (0.09019607843137255, 0.7450980392156863, 0.8117647058823529)]

access seaborn colors via core.cp[0] - core.cp[9]

class pyrexMD.core.iColor(iPlayer_cls)[source]

Bases: object

get_color_scheme(color, universe=None)[source]

Returns color scheme.

Parameters
  • color (str) –

  • universe (universe, atomgrp) –

Returns

color_scheme (array)

HEX colors for each RES taken from iPlayer.universe (default) or universe.

Accepted colors (simple):
‘red’ ‘black’
‘orange’ ‘grey_80’
‘yellow’ ‘grey_60’
‘green’ ‘grey_40’
‘cyan’ ‘grey_20’
‘blue’ ‘white’
‘purple’
‘magenta’
Accepted colors (complex):
‘baw’: black and white
‘b/w’: black and white
‘rainbow’: red to magenta
set_color_scheme(color, ci=None)[source]

Applies color scheme. If ci=None, then target is determinded by <iPlayer_object>.widgets.Components.description.

Parameters
  • color (str) –

  • universe (universe, atomgrp) –

Returns

color_scheme (array)

HEX colors for each RES taken from iPlayer.universe (default) or universe.

Accepted colors (simple):
‘red’ ‘black’
‘orange’ ‘grey_80’
‘yellow’ ‘grey_60’
‘green’ ‘grey_40’
cyan’ ‘grey_20’
‘blue’ ‘white’
‘purple’
‘magenta’
Accepted colors (complex):
‘baw’: black and white
‘b/w’: black and white
‘rainbow’: red to magenta
class pyrexMD.core.iPlayer(universe=None)[source]

Bases: object

show_player(layout='default')[source]

Show trajectory viewer with GUI.

Parameters

layout (str) –

‘default’: default layout
’Marie’: special layout for special person

Alternative:

Execute show_player() method by calling the object.

Example

tv = core.iPlayer(<universe>) tv() # short version of tv.show_player()

sync_view()[source]

Alias for <iPlayer_object>.player.sync_view().

class pyrexMD.core.iPlot(universe=None, xdata=None, ydata=None, xlabel='X', ylabel='Y', title='', tu='ps', figsize=(8, 4.5), layout='default')[source]

Bases: pyrexMD.core.iPlayer

show_plot(xdata=None, ydata=None, xlabel='X', ylabel='Y', title='', tu='ps', figsize=(8, 4.5))[source]

Show trajectory viewer with GUI and interactive matplotlib plot. Interactive red bar can be moved by pressing any key.

Parameters
  • xdata (array) –

  • ydata (array) –

  • xlabel (str) –

  • ylabel (str) –

  • title (str) –

  • tu (str) –

    time unit of plot. Either ‘ps’, ‘ns’ or ‘frame’.
    If ‘ns’ is selected, time stamps of MDAnalysis universe are converted from ps to ns.
    Important to make the interactive red bar work properly.

  • figsize (tuple) –

Note

Alternatively execute show_plot() method by calling the object.

Example

ip = core.iPlot(<universe>)
ip() # short version of ip.show_plot()
class pyrexMD.core.iWidgets(iPlayer_cls)[source]

Bases: object

toggle_UI(a=None)[source]

Toggle UI with a=True/False. If a=None, then UI will switch to other state.

pyrexMD.gmx

Hint

This module contains modified GromacsWrapper functions for streamlining the interaction with GROMACS for system setups etc.

Example:

import pyrexMD.gmx as gmx
import pyrexMD.misc as misc

# create ref pdb:
pdb = "path/to/pdb"
ref = gmx.get_ref_structure(pdb, ff='amber99sb-ildn', water='tip3p', ignh=True)

# generate topology & box
gmx.pdb2gmx(f=ref, o="protein.gro", ff='amber99sb-ildn', water='tip3p', ignh=True)
gmx.editconf(f="protein.gro", o="box.gro", d=2.0, c=True, bt="cubic")

# generate solvent & ions
gmx.solvate(cp="box.gro", o="solvent.gro")
gmx.grompp(f="ions.mdp", o="ions.tpr",c="solvent.gro")
gmx.genion(s="ions.tpr", o="ions.gro", neutral=True, input="SOL")

# copy mdp files (min.mdp, nvt.mdp, npt.mdp, md.mdp) into working directory
misc.cp("path/to/mdp/files", ".")

# minimize
gmx.grompp(f="min.mdp", o="min.tpr", c="ions.gro")
gmx.mdrun(deffnm="min")

# NVT equilibration
gmx.grompp(f="nvt.mdp", o="nvt.tpr", c="min.gro", r="min.gro")
gmx.mdrun(deffnm="nvt")

# NPT equilibration
gmx.grompp(f="npt.mdp", o="npt.tpr", c="nvt.gro", r="nvt.gro", t="nvt.cpt")
gmx.mdrun(deffnm="npt")

# MD run
gmx.grompp(f="md.mdp", o="traj.tpr", c="npt.gro", t="npt.cpt")
gmx.mdrun(deffn="traj")

Content:

pyrexMD.gmx.clean_up(path='./', pattern='gmx_pattern', ignore=None, verbose=True)[source]

Clean up GROMACS backup files with the patterns #*# and .*offsets.npz

Note

pattern can be implemented into path variable directly.

Parameters
  • path (str) – directory path

  • pattern (None, str) –

    pattern
    None: check for files with path only
    str: check for files with joined path of path + pattern
    ”gmx_pattern”: remove files with the patterns #*# and .*offsets.npz

  • ignore (None, str, list of str, tuple of str) –

    one or multiple ignore pattern
    None: no ignore pattern
    str: single ignore pattern
    list of str, tuple of str: list/tuple with multiple ignore patterns

  • verbose (bool) – print message (‘removed file: … ‘)

pyrexMD.gmx.convert_TPR(s, o='default', odir='./', sel='protein', verbose=True, **kwargs)[source]

Modified function of gromacs.convert_tpr().

GROMACS info:

‘gmx convert-tpr’ can edit run input files in three ways: - modify tpr settings - create new tpr for a subset of original tpr - setting the charges of a specified group to zero(useful when doing free energy estimates using the LIE (Linear Interaction Energy) method.

Parameters
  • s (str) – tpr

  • o (str) –

    tpr
    ”default”: <s>.tpr -> <s>_<sel>.tpr

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • sel (str) –

    selection string(case independant)
    ”system” or “all”: all atoms
    ”protein”: protein atoms
    ”ca” or “calpha”: CA atoms
    ”bb” or “backbone”: backbone atoms

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments

cprint_color (str) –

Hint

Find more valid Keyword Args via

  • python -> gromacs.convert_tpr.help()

  • terminal -> gmx convert_tpr -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.convert_TPR2PDB(f, o='default', odir='./', verbose=True, **kwargs)
  • Modified function of gromacs.editconf()

  • Alias function of convert_TPR2PDB()

GROMACS info:

‘gmx editconf’ converts generic structure format to .pdb, .gro, or .g96.

Parameters
  • f (str) – tpr (gro g96 pdb brk ent esp)

  • o (str) –

    pdb (gro g96)
    ”default”: <f>.tpr -> <f>.pdb
    <f>.gro -> box.gro

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • box (int) –

    box vector lengths xyz (only 1 value for bt=”cubic”).
    requires center=False to work.

  • bt (str) – box type: cubic triclinic dodecahedron octahedron

  • c (bool) – center molecule in box

  • d (str) – distance between the solute and the box

  • cprint_color (str) –

  • log (bool) – save log file

  • log_overwrite (bool) –

    True: save log file as <logs/editconf.log>
    False: save log file as <logs/editconf_{i}.log> for i=1,…,999

Hint

Find more valid Keyword Args via

  • python -> gromacs.editconf.help()

  • terminal -> gmx editconf -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.create_complex(f=[], o='complex.gro', verbose=True)[source]

Create complex using multiple .gro files. First .gro file in <f> must have largest box.

Parameters
  • f (list) – list of gro files

  • o (str) – output structure file: gro

  • verbose (bool) –

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.create_positions_dat(box_vec=[25, 25, 25], nmol=range(100, 110, 10), verbose=True)[source]

Create folder with position.dat files. Each file contains <nmol> copies of the box center. Use in combination with “gmx insert-molecule module”

Parameters
  • box_vec (list) – box vectors

  • nmol (int, tuple, list, range) – number of molecules to add (replaces -nmol from insert-molecule cmd)

  • verbose (bool) –

Returns

file_dir (str)

realpath of folder containing positions.dat files

Example

>> gmx insert-molecules -f “box_25.gro” -ci “../model_ATP/atp_ini.pdb”
-o “complex.gro” -ip “positions_dat/positions_100.dat” -dr 5 5 5 -try 100
pyrexMD.gmx.editconf(f, o='default', odir='./', verbose=True, **kwargs)[source]
  • Modified function of gromacs.editconf()

  • Alias function of convert_TPR2PDB()

GROMACS info:

‘gmx editconf’ converts generic structure format to .pdb, .gro, or .g96.

Parameters
  • f (str) – tpr (gro g96 pdb brk ent esp)

  • o (str) –

    pdb (gro g96)
    ”default”: <f>.tpr -> <f>.pdb
    <f>.gro -> box.gro

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • box (int) –

    box vector lengths xyz (only 1 value for bt=”cubic”).
    requires center=False to work.

  • bt (str) – box type: cubic triclinic dodecahedron octahedron

  • c (bool) – center molecule in box

  • d (str) – distance between the solute and the box

  • cprint_color (str) –

  • log (bool) – save log file

  • log_overwrite (bool) –

    True: save log file as <logs/editconf.log>
    False: save log file as <logs/editconf_{i}.log> for i=1,…,999

Hint

Find more valid Keyword Args via

  • python -> gromacs.editconf.help()

  • terminal -> gmx editconf -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.extend_complex_topology(ligand_name, ligand_itp, ligand_prm, ligand_nmol, top='topol.top', top_out='topol_complex.top', verbose=True)[source]

Extend an existing .top file.

Parameters
  • ligand_name (str) – ligand name

  • ligand_itp (str) – ligand .itp file

  • ligand_prm (str) – ligand .prm file

  • ligand_nmol (int) – number of ligand molecules

  • top (str) – topology input file (will be extended)

  • top_out (str) – topology output file

  • verbose (bool) –

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.fix_TRAJ(tpr, xtc, o='default', odir='./', tu='ns', sel='protein', pbc='mol', center=True, verbose=True, **kwargs)[source]
Fix topology:
  • convert tpr to selection only

  • if o=”default”: save as <s>_<sel>.tpr

Fix trajectory:
  • selection only

  • tu ns

  • pbc mol

  • center

  • if o=”default”: save as <f>_<sel>.xtc

Parameters
  • tpr (str) – structure: tpr

  • xtc (str) – trajectory: xtc (trr cpt gro g96 pdb tng)

  • o (str, list) –

    filenames of new structure and new trajectory (both selection only)
    ”default” (str):
    <tpr>.tpr -> <xtc>_<sel>.tpr
    <xtc>.xtc -> <xtc>_<sel>.xtc
    [os, of] (list):
    os (str): filename of new structure (selection only)
    of (str): filename of new trajectory (selection only)

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • tu (str) – time unit: fs ps ns us ms s

  • sel (str) –

    selection string(case independant)
    ”system” or “all”: all atoms
    ”protein”: protein atoms
    ”ca” or “calpha”: CA atoms
    ”bb” or “backbone”: backbone atoms

  • pbc (str) –

    PBC treatment: none mol res atom nojump cluster whole
    ”mol”: puts the center of mass of molecules in the box. requires structure file s.
    ”nojump”: checks if atoms jump across the box and then puts them back (i.e. all molecules will remain whole)
    (see GROMACS help text for other descriptions)

  • center (bool) – center atoms in box

  • verbose (bool) – print/mute GROMACS messages

Hint

Find more valid Keyword Args via

  • python -> gromacs.trjconv.help()

  • terminal -> gmx trjconv -h

Returns

tpr_file (str)

realpath of new tpr file (selection only)

xtc_file (str)

realpath of new xtc file (selection only)

pyrexMD.gmx.genion(s, o, p='topol.top', input='SOL', pname='NA', nname='CL', conc=0.15, neutral=True, verbose=True, **kwargs)[source]

Modified fuction of gromacs.genion().

GROMACS info:

gmx genion randomly replaces solvent molecules with monoatomic ions. The group of solvent molecules should be continuous and all molecules should have the same number of atoms.

Parameters
  • s (str) – structure file: tpr

  • o (str) – output file: gro pdb (g96 brk ent esp)

  • p (str) – topology file: topol.top

  • input (str) –

    selection grp
    default: 13 (SOL)

  • pname (str) – positive ion name

  • nname (str) – negative ion name

  • conc (float) – add salt concentration (mol/liter) and rescale to box size

  • neutral (bool) – add enough ions to neutralize the system. These ions are added on top of those specified with -np/-nn or -conc

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • cprint_color (str) –

  • log (bool) – save log file

  • log_overwrite (bool) –

    True: save log file as <logs/genion.log>
    False: save log file as <logs/genion_{i}.log> for i=1,…,999

Hint

Find more valid Keyword Args via

  • python -> gromacs.genion.help()

  • terminal -> gmx genion -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.get_RMSD(ref, xtc, o='default', odir='./', tu='ns', sel=['bb', 'bb'], verbose=True, **kwargs)[source]

Modified function of gromacs.rms(). Calculate backbone RMSD.

Parameters
  • ref (str) – reference structure: pdb (tpr gro g96 brk ent)

  • xtc (str) – trajectory: xtc (trr cpt gro g96 pdb tng)

  • o (str) –

    xvg
    ”default”: rmsd.xvg

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • tu (str) – time unit: fs ps ns us ms s

  • sel (list) –

    list with selection strings (case independant)
    ”system” or “all”: all atoms
    ”protein”: protein atoms
    ”ca” or “calpha”: CA atoms
    ”bb” or “backbone”: backbone atoms

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments

cprint_color (str) –

Hint

Find valid Keyword Args via

  • python -> gromacs.rms.help()

  • terminal -> gmx rms -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.get_ref_structure(f, o='default', odir='./', ff='amber99sb-ildn', water='tip3p', ignh=True, verbose=True, **kwargs)[source]

Creates ref structure of input structure via pdb2gmx to fix atom count by applying force field.

Parameters
  • f (str) – input structure file: pdb gro tpr (g96 brk ent esp)

  • o (str) – output structure file: pdb gro trp (g96 brk ent esp)

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • ff (str) –

    force field (see <GROMACS_path>/top/<ff> for valid inputs)
    ”amber99sb-ildn” etc. for proteins
    ”amber14sb_OL15” etc. for RNAs

  • water (str) – water model

  • ignh (bool) – ignore hydrogen

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments

cprint_color (str) –

Hint

Find valid Keyword Args via

  • python -> gromacs.pdb2gmx.help()

  • terminal -> gmx pdb2gmx -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.get_template_1_EM(gmx_mpi=True)[source]
Parameters

gmx_mpi (bool) – use gmx_mpi instead of gmx

pyrexMD.gmx.get_template_2_NVT(gmx_mpi=True)[source]
Parameters

gmx_mpi (bool) – use gmx_mpi instead of gmx

pyrexMD.gmx.get_template_3_NPT(gmx_mpi=True)[source]
Parameters

gmx_mpi (bool) – use gmx_mpi instead of gmx

pyrexMD.gmx.get_template_4_MD(gmx_mpi=True)[source]
Parameters

gmx_mpi (bool) – use gmx_mpi instead of gmx

pyrexMD.gmx.grompp(f, o, c, p='topol.top', verbose=True, **kwargs)[source]

Modified function of gromacs.grompp().

Parameters
  • f (str) – input file: mdp

  • o (str) – output file: tpr

  • c (str) – structure file: gro tpr pdb (g96 brk ent esp)

  • p (str) – topology file: top

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • maxwarn (int) –

    Number of allowed warnings during input processing.
    Not for normal use and may generate unstable systems.

  • cprint_color (str) –

  • log (bool) – save log file

  • log_overwrite (bool) –

    True: save log file as <logs/grompp.log>
    False: save log file as <logs/grompp_{i}.log> for i=1,…,999

Hint

Find more valid Keyword Args via

  • python -> gromacs.grompp.help()

  • terminal -> gmx grompp -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.mdrun(verbose=True, **kwargs)[source]

Alias fuction of gromacs.mdrun().

Note

output can be printed only after mdrun has finished. To see realtime updates invoke the command using “!gmx_mpi mdrun <parameters>” instead

Parameters

verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • deffnm (str) – default filename

  • nsteps (int) – maximum number of steps (used instead of mdp setting)

Hint

Find more valid Keyword Args via

  • python -> gromacs.mdrun.help()

  • terminal -> gmx mdrun -h

pyrexMD.gmx.pdb2gmx(f, o='protein.gro', odir='./', ff='amber99sb-ildn', water='tip3p', ignh=True, verbose=True, **kwargs)[source]

Modified function of gromacs.pdb2gmx().

GROMACS info:

‘pdb2gmx’ reads a .pdb (or .gro) file, reads some database files, adds hydrogens to the molecules and generates coordinates in GROMACS (GROMOS), or optionally .pdb, format and a topology in GROMACS format. These files can subsequently be processed to generate a run input file.

Parameters
  • f (str) – input structure file: pdb gro tpr (g96 brk ent esp)

  • o (str) – output structure file: pdb gro (g96 brk ent esp)

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • ff (str) –

    force field (see <GROMACS_path>/top/<ff> for valid inputs)
    ”amber99sb-ildn” etc. for proteins
    ”amber14sb_OL15” etc. for RNAs

  • water (str) – water model

  • ignh (bool) – ignore hydrogen

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • p (str) – topology file: topol.top

  • i (str) – include file for topology: posre.itp

  • n (str) – index file: index.ndx

  • cprint_color (str) –

  • log (bool) – save log file

  • log_overwrite (bool) –

    True: save log file as <logs/pdb2gmx.log>
    False: save log file as <logs/pdb2gmx_{i}.log> for i=1,…,999

Hint

Find more valid Keyword Args via

  • python -> gromacs.pdb2gmx.help()

  • terminal -> gmx pdb2gmx -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.solvate(cp, cs='spc216.gro', o='solvent.gro', p='topol.top', verbose=True, **kwargs)[source]

Modified function of gromacs.solvate(). cp is usually “box.gro”

Parameters
  • cp (str) – structure file ~ solute: gro pdb tpr (g96 brk ent esp)

  • cs (str) –

    structure file ~ solvent: gro pdb tpr (g96 brk ent esp)
    Note: “spc216.gro” is used for all 3-point water models

  • o (str) –

    gro pdb (g96 brk ent esp)
    default case: save in same directory as cp
    special case: if o is relative/absolute path -> save there

  • p (str) – topology file: topol.top

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments
  • maxsol (int) –

    maximum number of solvent molecules to add if they fit in the box.
    0 (default): ignore this setting

  • cprint_color (str) –

  • log (bool) – save log file

  • log_overwrite (bool) –

    True: save log file as <logs/solvate.log>
    False: save log file as <logs/solvate_{i}.log> for i=1,…,999

Hint

Find more valid Keyword Args via

  • python -> gromacs.solvate.help()

  • terminal -> gmx solvate -h

Returns

ofile (str)

realpath of output file

pyrexMD.gmx.trjconv(s, f, o='default', odir='./', sel='protein', verbose=True, **kwargs)[source]

Modified function of gromacs.trjconv().

GROMACS info:

gmx trjconv can convert trajectory files in many ways

  • from one format to another

  • select a subset of atoms

  • center atoms in the box

  • fit atoms to reference structure

  • reduce the number of frames

  • etc.

Parameters
  • s (str) – structure: tpr (gro g96 pdb brk ent)

  • f (str) – trajectory: xtc (trr cpt gro g96 pdb tng)

  • o (str) – trajectory: | xtc (trr gro g96 pdb tng) | “default”: <f>.xtc -> <f>_<sel>.xtc

  • odir (str) –

    output directory
    special case: odir is ignored when o is relative/absolute path

  • sel (str) –

    selection string (case independant)
    ”system” or “all”: all atoms
    ”protein”: protein atoms
    ”ca” or “calpha”: CA atoms
    ”bb” or “backbone”: backbone atoms

  • verbose (bool) – print/mute GROMACS messages

Keyword Arguments

cprint_color (str) –

Hint

Find more valid Keyword Args via

  • python -> gromacs.trjconv.help()

  • terminal -> gmx trjconv -h

Returns

ofile (str)

realpath of output file

pyrexMD.rex

Hint

This module contains functions related to (contact-guided) Replica Exchange Molecular Dynamics, mainly for automating and speeding up the simulation setup.

Example:

import pyrexMD.misc as misc
import pyrexMD.rex as rex
import pyrexMD.topology as top

decoy_dir = "path/to/decoy/directory"

# create rex_i directories and assign decoys
rex.assign_best_decoys(decoy_dir)
rex_dirs = rex.get_REX_DIRS()

# check for consistent topology
rex.check_REX_PDBS(decoy_dir)

# copy mdp files (min.mdp, nvt.mdp, npt.mdp, rex.mdp) into working directory
misc.cp("path/to/mdp/files", ".")

# get parameters for fixed box size and solvent molecules
boxsize, maxsol = rex.WF_get_system_parameters(wdir="./rex_0_get_system_parameters/")

# create systems for each replica
rex.WF_REX_setup(rex_dirs=rex_dirs, boxsize=boxsize, maxsol=maxsol)

# minimize
rex.WF_REX_setup_energy_minimization(rex_dirs=rex_dirs, nsteps=10, verbose=False)

# add bias contacts (RES pairs defined in DCA_fin)
top.DCA_res2atom_mapping(ref_pdb=<ref_pdb>, DCA_fin=<file_path>, n_DCA=50, usecols=(0,1))
top.DCA_modify_topology(top_fin="topol.top", DCA_used_fin=<file_path> , k=10, save_as="topol_mod.top")

# prepare temperature distribution
rex.prep_REX_temps(T_0=280, n_REX=len(rex_dirs), k=0.006)

# create mdp and tpr files
rex.prep_REX_mdp(main_dir="./", n_REX=len(rex_dirs))
rex.prep_REX_tpr(main_dir="./", n_REX=len(rex_dirs))

# next: upload REX MD run files on HPC and execute production run

Content:

pyrexMD.rex.WF_REX_setup(rex_dirs, boxsize, maxsol, ff='amber99sb-ildn', water='tip3p', ignh=True, verbose=False, verbose_gmx=False)[source]

Workflow: REX setup (without energy minimization)

Parameters
  • rex_dirs (list) – list with rex_dirs (output of rex.get_REX_DIRS())

  • boxsize (float) – suggested boxsize parameter (output of gmx.WF_getParameter_boxsize())

  • maxsol (int) – suggested max solution parameter (output of gmx.WF_getParameter_max_solution())

  • ff (str) –

    force field (see <GROMACS_path>/top/<ff> for valid inputs)
    ”amber99sb-ildn” etc. for proteins
    ”amber14sb_OL15” etc. for RNAs

  • water (str) – water model

  • ignh (bool) – ignore hydrogen

  • verbose (bool) – show blue log messages (saved file as, saved log as)

  • verbose_gmx (bool) – show gmx module messages (requires verbose=True)

pyrexMD.rex.WF_REX_setup_energy_minimization(rex_dirs, nsteps=None, verbose=False)[source]

Workflow: REX energy minimization

Parameters
  • rex_dirs (list) – list with rex_dirs (output of rex.get_REX_DIRS())

  • nsteps (None,int) –

    maximum number of steps
    None: use .mdp option
    int: use instead of .mdp option

  • verbose (bool) – show/hide GROMACS output

pyrexMD.rex.WF_getParameter_boxsize(logfile='./logs/editconf.log', base=0.2, verbose=True)[source]

Read <logfile> and suggest a ‘fixed boxsize’ parameter for REX simulations.

Parameters
  • logfile (str) – path to <editconf.log> containing the line ‘system size: X Y Z’

  • base (float) – round up base for highest box dimension taken from <logfile>

  • verbose (bool) –

Returns

boxsize (float)

suggested boxsize parameter

pyrexMD.rex.WF_getParameter_maxsol(logfile='./logs/solvate.log', maxsol_reduce=50, verbose=True)[source]

Read <logfile> and suggest a ‘max solution’ parameter for REX simulations.

Parameters
  • logfile (str) – path to <editconf.log> containing the line ‘system size: X Y Z’

  • maxsol_reduce (int) –

    reduce max solution number taken from <logfile>
    -> guarantees fixed solution number for different start configurations

  • verbose (bool) –

Returns

maxsol (int)

suggested max solution parameter

pyrexMD.rex.WF_get_system_parameters(wdir, verbose=False, **kwargs)[source]

Get system parameters “boxsize” and “maxsol” for system in working directory wdir. Typically, wdir is supposed to be the folder “rex_0_get_system_parameters” which is created during the Workflow for setup of REX MD simulations during the execution of rex.assign_best_decoys().

Parameters
  • wdir (str) – working directory

  • verbose (bool) –

Keyword Arguments
  • d (int) – distance between model and system edges. Defaults to 2.

  • bt (str) – box type. Defaults to “cubic”.

Returns

boxsize (float)

suggested value for maximum boxsize

maxsol (int)

suggested value for maximum solvent molecules

pyrexMD.rex.apply_ff_best_decoys(best_decoys_dir, n_decoys=None, odir='./PDBID_best_decoys_ref', create_dir=True, verbose=False, **kwargs)[source]

Apply forcefield on best decoys and save as <filename>_ref.pdb.

Parameters
  • best_decoys_dir (str) –

    directory with

    • best decoys (output of cluster.rank_cluster_decoys() -> cluster.copy_cluster.decoys())

    • decoy_scores.log (output of cluster.log_cluster.decoys())

  • n_decoys (None, int) –

    None: use all decoys from best_decoys_dir
    int: use this number of decoys from best_decoys_dir

  • odir (str) –

    output directory
    Note: if “PDBID” in <odir> and no <pdbid> kwarg is passed -> find and replace “PDBID” in <odir> automatically based on filenames

  • create_dir (bool) –

  • verbose (bool) –

Keyword Arguments
  • pdbid (str) –

  • logfile (str) – logfile name (default: “decoy_scores.log”)

  • water (str) – water model (default: “tip3p”)

  • input (str) – forcefield number (default: “6”)

  • cprint_color (str) –

Returns

odir (str)

output directory with ref pdb files (forcefield is applied)

pyrexMD.rex.assign_best_decoys(best_decoys_dir, rex_dir='./', verbose=False, **kwargs)[source]

Assigns decoys based on ranking taken from <best_decoys_dir>/decoy_scores.log to each rex subdirectory rex_1, rex_2, … Also creates a ‘rex_0_get_system_parameters’ directory, which is a copy of rex_1.

Parameters
  • best_decoys_dir (str) –

    directory with

    • best decoys (output of cluster.rank_cluster_decoys() -> cluster.copy_cluster.decoys())

    • decoys_score.log (output of cluster.log_cluster.decoys())

  • rex_dir (str) – rex directory with folders rex_1, rex_2, …

  • verbose (bool) –

Keyword Arguments
  • cprint_color (str) –

  • realpath (bool) – prints realpaths

pyrexMD.rex.check2_REX_PDBS(REX_PDBS, ref_pdb=None, verbose=False)[source]

Used for debbuging if check_REX_PDBS() fails.

Tests if all REX PDBS have equal topologies (RES, ID, NAME arrays). Uses ref_pdb as template for tests if passed, otherwise uses first REX pdb as template.

Note

  • if target is known -> apply ff -> save as ref -> use as ref

  • if target is unknown -> use ref_pdb=None or use one of REX_PDBS -> apply ff -> use as ref

Parameters
  • REX_PDBS (list) – output of get_REX_PDBS()

  • ref_pdb (None, str) –

    reference pdb
    None: use ref_pdb as template for tests
    str: use first REX pdb as template for tests

  • sel (str) – selection string

  • ref_is_known (bool) –

    True: use ref_pdb as template for tests
    False: use first REX DPB as template for tests

  • verbose (bool) –

pyrexMD.rex.check_REX_PDBS(REX_PDBS, ref_pdb=None, sel='protein', verbose=True, **kwargs)[source]

Tests if all REX PDBS have equal topologies (RES, ID, NAME arrays). Uses ref_pdb as template for tests if passed, otherwise uses first REX pdb as template.

Note

  • if target is known -> apply ff -> save as ref -> use as ref

  • if target is unknown -> use ref_pdb=None or use one of REX_PDBS -> apply ff -> use as ref

  • use check2_REX_PDBS() for debugging if check_REX_PDBS() fails.

Parameters
  • REX_PDBS (list) – output of get_REX_PDBS()

  • ref_pdb (None, str) –

    reference pdb
    None: use ref_pdb as template for tests
    str: use first REX pdb as template for tests

  • sel (str) – selection string

  • ref_is_known (bool) –

    True: use ref_pdb as template for tests
    False: use first REX DPB as template for tests

  • verbose (bool) –

pyrexMD.rex.create_pull_group_ndx(ref_pdb, sel='name CA', save_as='pull_group.ndx')[source]

Create index file for pull group.

Parameters
  • ref_pdb (str) – PDB file (after applying force field)

  • sel (str) – selection string

  • save_as (str) –

pyrexMD.rex.create_special_group_ndx(ref_pdb, sel='name CA', save_as='special_group.ndx')[source]

Create index file for special group.

Parameters
  • ref_pdb (str) – PDB file (after applying force field)

  • sel (str) – selection string

  • save_as (str) –

pyrexMD.rex.get_REX_DIRS(main_dir='./', realpath=True)[source]

get list with REX DIRS.

Parameters
  • main_dir (str) – directory with folders rex_1, rex_2, etc.

  • realpath (bool) – return realpath

Returns

REX_DIRS (list)

list with REX directory paths

pyrexMD.rex.get_REX_PDBS(main_dir='./', realpath=True)[source]

get list with REX PDBS.

Parameters
  • main_dir (str) – directory with folders rex_1, rex_2, etc.

  • realpath (bool) – return realpath

Returns

REX_PDBS (list)

list with PDB paths within folders rex_1, rex_2, etc.

pyrexMD.rex.prep_REX_mdp(main_dir='./', template='rex.mdp', n_REX=None, verbose=True)[source]

Prepare REX mdp -> copy template and change tempereratures according to rex_temps.log

Parameters
  • main_dir (str) – main directory with rex_1, rex_2, etc.

  • template (str) – template file

  • n_REX (None/int) – number of replica

  • verbose (bool) –

pyrexMD.rex.prep_REX_temps(T_0=None, n_REX=None, k=None)[source]

Prepare REX temperatures.

Parameters
  • T_0 (None/float) – starting temperature (Kelvin)

  • n_REX (None/int) – number of replica

  • k (None/float) – temperature distrubtion scaling factor

pyrexMD.rex.prep_REX_tpr(main_dir='./', n_REX=None, verbose=False, **kwargs)[source]

Prepare REX tpr.

Parameters
  • main_dir (str) – main directory with rex_1, rex_2, etc.

  • n_REX (None/int) – number of replica

  • verbose (bool) –

Keyword Arguments
  • f (str) –

  • o (str) –

  • c (str) –

  • p (str) –

Note

Keyword Args are parameter of gmx.grompp(f,o,c,r,p)

pyrexMD.topology

Hint

This module contains functions for modifying universe topologies, e.g., align atoms/residues of two universes, get matching selection strings, include bias contacts.

pyrexMD.topology.DCA_modify_scoreFile(score_fin, shift_res, res_cols=(0, 1), score_col=2, **kwargs)[source]

Modify score file (MSA scores) by shifting residues.

Parameters
  • score_fin (str) – path to score file

  • shift_res (int) – shift residues by this value

  • outputFileName (str) –

    realpath to modified score file.
    if “PDBID_mod.score” is left as default: try to automatically detect
    pattern based on score_fin and add the “_mod” part into filename.

  • res_cols (tuple/list) – score_fin columns with residue numbers

  • score_col (tuple/list) – score_fin column with score/confidence

Keyword Arguments

save_as (str) –

“PDBID_mod.score”
if “PDBID_mod.score” (default): try to automatically detect pattern
based on score_fin and insert the “_mod” part into filename.

Returns

save_as (str)

realpath to modified score file

pyrexMD.topology.DCA_modify_topology(top_fin, DCA_used_fin, n_DCA=None, k=10, skiprows='auto', **kwargs)[source]

Modifies topology by using the contacts written in “DCA_used.txt” file. “DCA_used.txt” is supposed to have 4 columns: RESi, RESj, ATOMi, ATOMj.

Modify topology:

  • top_fin (topol.top file): use as template

  • DCA_used_fin (DCA_used.txt): use all contacts as restraints

  • modify bond section of new topology by adding contacts with force constant

Parameters
  • top_fin (str) – topology file (path)

  • DCA_used_fin (str) – DCA file (path)

  • n_DCA (None, int) –

    number of used DCA contacts
    None: search header of DCA_used_fin for entry with topn_DCA

  • k (int, float) – force coefficient of contact pairs

  • skiprows (int) –

    ignore header rows of DCA_used_fin
    -1 or “auto”: auto detect

Keyword Arguments
  • pdbid (str) – “auto” (default): detect PDBID based on ref_PDB path

  • save_as (str) –

    “PDBID_topol_mod.top” (default)
    detect and replace PDBID in save_as based on ref_PDB path
    if kwarg <pdbid> is “auto” (default).

pyrexMD.topology.DCA_res2atom_mapping(ref_pdb, DCA_fin, n_DCA, usecols, DCA_skiprows='auto', filter_DCA=True, save_log=True, **kwargs)[source]

Get contact mapping. Returns lists of matching RES pairs and ATOM pairs for contacts specified in DCA_fin file.

Note

  • maps only CA atoms if ref_pdb is protein.

  • maps only N1 or N3 atoms if ref_pdb is nucleic.

Parameters
  • ref_pdb (str) – reference PDB (path)

  • DCA_fin (str) – DCA file (path)

  • n_DCA (int) – number of DCA contacts which should be used

  • usecols (tuple, list) – columns containing the RES pairs in DCA_fin

  • DCA_skiprows (int) –

    ignore header rows of DCA_fin
    -1 or “auto”: auto detect

  • filter_DCA (bool) –

    True: ignore DCA pairs with abs(i-j) < 3
    False: use all DCA pairs w/o applying filter

  • save_log (bool) – create log file with contact mapping

Keyword Arguments
  • cprint_color (None, str) – colored print color

  • pdbid (str) – “auto” (default): detect PDBID based on ref_PDB path

  • default_dir (str) – “./”

  • save_as (str) –

    “PDBID_DCA_used.txt”
    detect and replace PDBID in kwarg <save_as> based on ref_PDB path
    if kwarg <pdbid> is “auto” (default).

Returns

RES_PAIR (list)

list with RES pair tuples (RESi, RESj)

ATOM_PAIR (list)

list with ATOM pair tuples (ATOMi, ATOMj)

pyrexMD.topology.align_resids(mobile, ref, norm=True, verbose=True, **kwargs)[source]

Align resids of mobile and ref by comparing universe.residues.resnames and shifting the resids of reference (usually smaller than mobile).

Parameters
  • mobile (universe) –

  • ref (universe) –

  • norm (bool) – apply topology.norm_resids()

  • verbose (bool) –

Keyword Arguments

cprint_color (None, str) – colored print color

pyrexMD.topology.check_matching_selection(mobile, ref, sel='protein', verbose=True)[source]

Gets matching selection strings for mobile and reference and then tests if matching selections have equal atoms and residues. Returns matching selection strings.

Parameters
  • mobile (universe) –

  • ref (universe) –

  • sel (str) – selection string

  • verbose (bool) –

    True: prints detailed report
    False: prints only selection strings and if atoms/residues match

pyrexMD.topology.check_residues(mobile, ref, sel1='all', sel2='all', verbose=True)[source]

Compares residues of mobile and reference universe using specific selections.

Parameters
  • mobile (universe) –

  • ref (universe) –

  • sel1 (str) – selection string of mobile

  • sel2 (str) – selection string of ref

  • verbose (bool) –

    True: prints detailed report
    False: prints only if residue ids and names match

pyrexMD.topology.dump_structure(u, frames, save_as, default_dir='./structures', sel='protein')[source]

Dump structures for a list of frames with the file name “save_as”. Automatically prepends frame number to the extension.

Parameters
  • u (universe, atomgrp) – universe containing the structure

  • frames (array, list) –

  • save_as (str) – save name or realpath to save file. frame number will be prepended to the extension

  • default_dir (str) –

  • sel (str) – selection string

Returns

dirpath (str)

realpath to directory with dumped structures

pyrexMD.topology.get_matching_selection(mobile, ref, sel='protein and name CA', norm=True, verbose=True)[source]

Get matching selection strings of mobile and reference after resids alignment.

Parameters
  • mobile (universe) –

  • ref (universe) –

  • sel (str) – selection string

  • norm (bool) – apply topology.norm_resids()

  • verbose (bool) –

Returns

sel1 (str)

matching selection string of mobile

sel2 (str)

matching selection string of ref

pyrexMD.topology.get_resids_shift(mobile, ref)[source]

Compares universe.residues.resnames between mobile and ref in order to get resids shift.

Parameters
  • mobile (universe) –

  • ref (universe) –

Returns

shift (int)

shift value of ref residues to match mobile residues

pyrexMD.topology.norm_and_align_universe(mobile, ref, verbose=True)[source]
  • Norm reference and mobile universe

  • Align mobile universe on reference universe (matching atom ids and res ids)

Parameters
  • mobile (universe, atomgrp) –

  • ref (universe, atomgrp) –

  • verbose (bool) –

pyrexMD.topology.norm_ids(u, info='', verbose=True)[source]

Modify existing MDA universe/atomgrp and normalize ids according to min(universe.atoms.ids) = 1

Parameters
  • u (universe, atomgrp) – structure

  • info (str) –

    additional info for print message
    ’reference’
    ’mobile’

  • verbose (bool) –

Example

>> ref = mda.Universe(<top>)
>> print(ref.atoms.ids)
[0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19]

>> norm_ids(ref, ‘reference’)
Norming atom ids of reference…

>> print(ref.atoms.ids)
[1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20]
pyrexMD.topology.norm_resids(u, info='', verbose=True)[source]

Modify existing MDA universe/atomgrp and normalize resids according to min(universe.residues.resids) = 1

Parameters
  • u (universe, atomgrp) – structure

  • info (str) –

    additional info for print message
    ’reference’
    ’topology’

Example

>> ref = mda.Universe(<top>)
>> print(ref.residues.resids)
[0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19]

>> norm_resids(u, ‘reference’)
Norming resids of reference…

>> print(ref.residues.resids)
[1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20]
pyrexMD.topology.norm_universe(u, *args, info='', verbose=True)[source]

Executes for each universe:

  • norm_ids(u, info)

  • norm_resids(u, info)

Parameters
  • u (universe, atomgrp, list, array) – universe or list/array with universes

  • info (str) –

    additional info for print message
    ’reference’
    ’mobile’

  • verbose (bool) –

pyrexMD.topology.parsePDB(fin, sel, norm=False, filter_rna=True, verbose=False)[source]

Reads PDB file and returns the columns for residue id, residue name, atom id and atom name as lists.

Parameters
  • fin (str) – PDB file

  • sel (str) – selection string

  • norm (bool) – apply norm_universe() before parsing.

  • filter_rna (bool) – filter selection only for N1 and N3 atoms based on nucleic residues

  • verbose (bool) –

Returns

RESID (list)

residue id column of PDB

RESNAME (list)

residue name column of PDB

ID (list)

atom id column of PDB

NAME (list)

atom name column of PDB

pyrexMD.topology.sel2selid(u, sel, norm=True, filter_rna=True)[source]

Converts selection string to selection id string.

Parameters
  • u (universe, str) – structure universe or path to structure pdb

  • sel (str) – selection string

  • norm (bool) – apply norm_universe()

  • filter_rna (bool) – filter selection only for N1 and N3 atoms based on nucleic residues

Returns

selid (str)

selection id string

Example

>> u = mda.Universe(“path/to/pdb”)
>> selid = sel2selid(u, sel=”protein and name CA”)
>> selid
‘index 4 18 37 58 77 94 118 137 159 171 178 193 199 210 221 228 260 274 288 294’
>> u.select_atoms(“protein and name CA”) == u.selectg_atoms(selid)
True
pyrexMD.topology.shift_resids(u, shift=None, verbose=True)[source]

Shift mda.universe.residues by <shift> value.

Parameters
  • u (universe) –

  • shift (None, int) – shift value

  • verbose (bool) –

pyrexMD.topology.true_select_atoms(u, sel='protein', ignh=True, norm=True)[source]

Create “true selection” (atomgroup copy in NEW universe) after applying

  • norm_ids()

  • norm_resids()

Note that atom indices are reassigned according to range(0, len(u.atoms))

Parameters
  • u (str, universe, atomgrp) – structure path (PDB id) or universe with structure

  • sel (str) – selection string

  • ignh (bool) – ignore hydrogen (mass < 1.2)

  • norm (bool) – apply topology.norm_universe()

Returns

a (atomgrp)

”true selection” ~ atomgroup copy in NEW universe

Note

mda.select_atoms() remembers information about the original universe configuration which is sometimes NOT WANTED.

Example

# Init universe
>> u = mda.Universe(<top>)

# select_atoms() behaviour
>> a = u.select_atoms(“protein and type C”)
>> a
<AtomGroup with 72 atoms>
>> a.residues.atoms
<AtomGroup with 964 atoms>

# true_select_atoms() behaviour
>> a = topology.true_select_atoms(u, sel=”protein and type C”)
>> a
<Universe with 72 atoms>
>> a.residue.atoms
<AtomGroup with 72 atoms>

Subpackages