docking_py package

Subpackages

Submodules

docking_py.cli module

Console script for docking_py.

docking_py.cli.main()

Console script for docking_py.

docking_py.docking module

Include the Docking class

class docking_py.docking.Docking(name, lig_pdb=None, rec_pdb=None, lig_pdbqt=None, rec_pdbqt=None, log_level=20)

Bases: object

Docking encapsulation class.

This class can be used to launch vina, smina, qvina and qvinaw.

Parameters:
  • name (str) – generic name of the system
  • lig_pdb (str, optional) – path of the ligand coordinate file (.pdb)
  • rec_pdb (str, optional) – path of the receptor coordinate file (.pdb)
  • lig_pdbqt (str, optional) – path of the ligand coordinate file (.pdbqt)
  • rec_pdbqt (str, optional) – path of the receptor coordinate file (.pdbqt)
  • dock_pdb (str, optional) – path of the docking ligand coordinate file (.pdb)
  • dock_log (str, optional) – path of the docking log file (.log)
align_receptor(ref_pdb, chain_ref=['A'], chain_rec=['A'])

Align self.rec_pdb to ref_pdb.

Example:
>>> pdb_manip.show_log()
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> dock_4yob = Docking(name='4yob')
>>> dock_4yob.extract_receptor(os.path.join(TEST_PATH, '4yob.pdb'),        TEST_OUT, {'res_name': pdb_manip.PROTEIN_RES})        #doctest: +ELLIPSIS
Succeed to read file ...4yob.pdb ,  916 atoms found
Succeed to save file ...4yob_rec.pdb
>>> dock_4yob.align_receptor(os.path.join(TEST_PATH, '1hsg.pdb'))
Succeed to read file .../4yob_rec.pdb ,  760 atoms found
Succeed to read file .../1hsg.pdb ,  1686 atoms found
PQITLWKRPIVTIKIGGQLKEALLNTGADDTVFEEVNLPGRWKPKLIGGIGGFVKVRQYDQVPIEICGHKVIGTVLVGPT
******|**|**************|*******|**||********|*******|*******| *******|*********
PQITLWQRPLVTIKIGGQLKEALLDTGADDTVLEEMSLPGRWKPKMIGGIGGFIKVRQYDQILIEICGHKAIGTVLVGPT
<BLANKLINE>
PTNVIGRNLMTQIGCTLNF
*|*|*****|*********
PVNIIGRNLLTQIGCTLNF
<BLANKLINE>
Succeed to save file ...4yob_rec.pdb
>>> coor_holo = pdb_manip.Coor(os.path.join(TEST_PATH, '1hsg.pdb'))        #doctest: +ELLIPSIS
Succeed to read file ...1hsg.pdb ,  1686 atoms found
>>> coor_rec = pdb_manip.Coor(dock_4yob.rec_pdb) #doctest: +ELLIPSIS
Succeed to read file ...4yob_rec.pdb ,  760 atoms found
>>> rmsd = coor_rec.compute_rmsd_to(coor_holo,        selec_dict={'name': ['CA'], 'chain':['A']})
>>> print('RMSD after alignement is {:.2f} Å'.format(rmsd))
RMSD after alignement is 1.50 Å
compute_dock_rmsd(ref_lig_pdb, selec_dict={})

Compute RMSD from docking pdb to ref_lig_pdb. By default use all atoms for RMSD calculation. To use only Calpha atoms define selec_dict={'name':['CA']}.

Parameters:
  • ref_lig_pdb (str) – PDB reference file
  • selec_dict (dict, optional, default={}) – Selection for RMSD calculation
Returns:

RMSD list

Return type:

list

display()

Display defined attribute of the Docking object.

dock_log
dock_pdb
dock_xml
extract_affinity()

Extract affinity from the docking .log file.

Returns:Affinity and RMSD informations as a dictionnary
Return type:dict
extract_autodock_pdb_affinity(out_pdb, reorder=True)

Extract pdb models from the the autodock log files.

CPU version

extract_autodock_pdb_affinity2(out_pdb, reorder=True)

Extract pdb models from the the autodock log files. Makes use of the xml generated by the gpu version.

GPU version

extract_lig_rec_pdb(coor_in, folder_out, rec_select_dict, lig_select_dict)
  • Extract receptor and ligand coordinates from a coor file
  • remove alternative location
  • Keep only amino acid residues
  • Save both coordinates and add it in the object
Parameters:
  • pdb_id (str) – PDB ID
  • rec_chain (list of str) – Chain(s) of the receptor
  • lig_chain (list of str) – Chain(s) of the ligand

Object field(s) changed:

  • self.rec_pdb
  • self.lig_pdb
Example:
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> dock_1hsg = Docking(name='1hsg')
>>> dock_1hsg.extract_lig_rec_pdb(os.path.join(TEST_PATH, '1hsg.pdb'),        TEST_OUT, {'res_name': pdb_manip.PROTEIN_RES}, {'res_name': 'MK1'})        #doctest: +ELLIPSIS
Succeed to read file ...1hsg.pdb ,  1686 atoms found
Succeed to save file ...1hsg_rec.pdb
Succeed to save file ...1hsg_input_lig.pdb
>>> coor_lig = pdb_manip.Coor(dock_1hsg.lig_pdb) #doctest: +ELLIPSIS
Succeed to read file ...1hsg_input_lig.pdb ,  45 atoms found
>>> coor_rec = pdb_manip.Coor(dock_1hsg.rec_pdb) #doctest: +ELLIPSIS
Succeed to read file ...1hsg_rec.pdb ,  1514 atoms found
extract_ligand(coor_in, folder_out, lig_select_dict)
  • Extract ligand coordinates
  • remove alternative location
  • Save coordinates and add it in the object

Object field(s) changed:

  • self.lig_pdb
Example:
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> dock_1hsg = Docking(name='1hsg')
>>> dock_1hsg.extract_ligand(os.path.join(TEST_PATH, '1hsg.pdb'),        TEST_OUT, {'res_name': 'MK1'})        #doctest: +ELLIPSIS
Succeed to read file ...1hsg.pdb ,  1686 atoms found
Succeed to save file ...1hsg_lig.pdb
>>> coor_lig = pdb_manip.Coor(dock_1hsg.lig_pdb) #doctest: +ELLIPSIS
Succeed to read file ...1hsg_lig.pdb ,  45 atoms found
extract_receptor(coor_in, folder_out, rec_select_dict)
  • Extract receptor coordinates
  • remove alternative location
  • Keep only amino acid residues
  • align structure on ref
  • Save coordinates and add it in the object
Parameters:
  • pdb_id (str) – PDB ID
  • ref_pdb (str) – Reference coordinates file
  • rec_chain (list of str) – Chain(s) of the receptor
  • rec_chain – Chain(s) of the reference file

Object field(s) changed:

  • self.rec_pdb
Example:
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> dock_1hsg = Docking(name='1hsg')
>>> dock_1hsg.extract_receptor(os.path.join(TEST_PATH, '1hsg.pdb'),        TEST_OUT, {'res_name': pdb_manip.PROTEIN_RES})        #doctest: +ELLIPSIS
Succeed to read file ...1hsg.pdb ,  1686 atoms found
Succeed to save file ...1hsg_rec.pdb
>>> coor_rec = pdb_manip.Coor(dock_1hsg.rec_pdb) #doctest: +ELLIPSIS
Succeed to read file ...1hsg_rec.pdb ,  1514 atoms found
get_gridfld()

Get gridfld from the .gpf file.

gpf
lig_pdb
lig_pdbqt
log_to_pdb(out_pdb)

Read autodock log, extract ligand model coordinates and extract affinities.

Parameters:out_pdb (str) – output pdb file

Warning

Difference between GPU and CPU version of autodock logs. Torsional Free Energy is not computed with GPU version.

prepare_grid(out_folder, gpf_out_prefix=None, spacing=0.375, grid_npts=None, center=None, check_file_out=True)

Grid preparation

Launch the prepare_gpf4.py command from MGLToolsPackage. And autogrid4.

prepare_ligand(lig_pdbqt=None, rigid=False, center=False, random_rot=False, check_file_out=True)

Ligand preparation to pdbqt format using the prepare_ligand4.py command. Can center the ligand, could be usefull with autodock (issues when x,y,z > 100 Å).

Parameters:
  • lig_pdbqt (str, optional, default=None) – output name
  • rigid (bool, optional, default=False) – Flag to define if ligand is rigid
  • center (bool, optional, default=False) – Flag to define if ligand have to centered
  • check_file_out (bool, optional, default=True) – flag to check or not if file has already been created. If the file is present then the command break.

Object requirement(s):

  • self.lig_pdb

Object field(s) changed:

  • self.lig_pdbqt
Example:
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> coor_1hsg = pdb_manip.Coor(os.path.join(TEST_PATH, '1hsg.pdb'))        #doctest: +ELLIPSIS
Succeed to read file ...tests/input/1hsg.pdb ,  1686 atoms found
>>> lig_coor = coor_1hsg.select_part_dict(            selec_dict={'res_name': 'MK1'})
>>> lig_atom_num = lig_coor.num
>>> print('Ligand has {} atoms'.format(lig_atom_num))
Ligand has 45 atoms
>>> out_lig = os.path.join(TEST_OUT,'lig.pdb')
>>> lig_coor.write_pdb(out_lig) #doctest: +ELLIPSIS
Succeed to save file .../lig.pdb
>>> test_dock = Docking('test', lig_pdb=out_lig)
>>> test_dock.prepare_ligand() #doctest: +ELLIPSIS
python2... .../prepare_ligand4.py -l lig.pdb -B none -A hydrogens -o lig.pdbqt
>>> coor_lig = pdb_manip.Coor(test_dock.lig_pdbqt)        #doctest: +ELLIPSIS
File name doesn't finish with .pdb read it as .pdb anyway
Succeed to read file .../lig.pdbqt ,  50 atoms found
>>> test_dock.display() #doctest: +ELLIPSIS
name         : test
lig_pdb      : .../lig.pdb
lig_pdbqt    : .../lig.pdbqt
ref_lig_pdb  : .../lig.pdb
prepare_receptor(rec_pdbqt=None, check_file_out=True)

Receptor preparation to pdbqt format using the prepare_receptor4.py command.

Parameters:
  • rec_pdbqt (str, optional, default=None) – output name
  • check_file_out (bool, optional, default=True) – flag to check or not if file has already been created. If the file is present then the command break.

Object requirement(s):

  • self.rec_pdb

Object field(s) changed:

  • self.rec_pdbqt
Example:
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> coor_1hsg = pdb_manip.Coor(os.path.join(TEST_PATH, '1hsg.pdb'))        #doctest: +ELLIPSIS
Succeed to read file .../1hsg.pdb ,  1686 atoms found
>>> # Keep only amino acid
>>> rec_coor = coor_1hsg.select_part_dict(selec_dict={'res_name': pdb_manip.PROTEIN_RES})
>>> out_rec = os.path.join(TEST_OUT,'rec.pdb')
>>> rec_coor.write_pdb(out_rec) #doctest: +ELLIPSIS
Succeed to save file .../rec.pdb
>>> rec_atom_num = rec_coor.num
>>> print('Receptor has {} atoms'.format(rec_atom_num))
Receptor has 1514 atoms
>>> test_dock = Docking('test', rec_pdb=out_rec)
>>> test_dock.prepare_receptor() #doctest: +ELLIPSIS
python2... .../prepare_receptor4.py -r .../rec.pdb -A checkhydrogens -o .../rec.pdbqt
>>> coor_rec = pdb_manip.Coor(test_dock.rec_pdbqt)        #doctest: +ELLIPSIS
File name doesn't finish with .pdb read it as .pdb anyway
Succeed to read file .../rec.pdbqt ,  1844 atoms found
>>> test_dock.display() #doctest: +ELLIPSIS
name         : test
rec_pdb      : .../rec.pdb
rec_pdbqt    : .../rec.pdbqt
random_rot_ligand()
  • Do a random rotation on ligand

Object field(s) changed:

  • self.lig_pdb
Example:
>>> TEST_OUT = str(getfixture('tmpdir'))
>>> dock_1hsg = Docking(name='1hsg')
>>> dock_1hsg.extract_ligand(os.path.join(TEST_PATH, '1hsg.pdb'),        TEST_OUT, {'res_name': 'MK1'})        #doctest: +ELLIPSIS
Succeed to read file ...1hsg.pdb ,  1686 atoms found
Succeed to save file ...1hsg_lig.pdb
>>> coor_lig = pdb_manip.Coor(dock_1hsg.lig_pdb) #doctest: +ELLIPSIS
Succeed to read file ...1hsg_lig.pdb ,  45 atoms found
>>> com_before = coor_lig.center_of_mass()
>>> dock_1hsg.random_rot_ligand()
Succeed to read file ...1hsg_lig.pdb ,  45 atoms found
Succeed to save file ...1hsg_lig.pdb
>>> coor_lig = pdb_manip.Coor(dock_1hsg.lig_pdb) #doctest: +ELLIPSIS
Succeed to read file ...1hsg_lig.pdb ,  45 atoms found
>>> com_after = coor_lig.center_of_mass()
>>> print('Same center of mass after rotation :{}'.format(com_before==com_after))
Same center of mass after rotation :[False False False]
..warning:
The function overwrite lig_pdb coordinates.
rec_com()

Get center of mass of the receptor pdb file.

rec_grid(buffer_space=30, spacing=1.0)

Compute grid from the receptor pdb file.

rec_pdb
rec_pdbqt
ref_lig_pdb
run_autodock(out_folder, dock_out_prefix=None, dock_log=None, dock_pdb=None, nrun=10, check_file_out=True, GPU=True)

Run autodock with cpu or gpu if available

run_autodock_cpu(out_folder, dock_out_prefix=None, dock_log=None, dock_pdb=None, dock_xml=None, dpf_out=None, nrun=10, param_list=[], check_file_out=True)
  1. Launch the prepare_dpf4.py command from MGLToolsPackage.
  2. Launch autodock4

This requires a gpf and associated map files, and ligand pdbqt It creates pose pdb + xml + dlg + dpf and smina like _log.txt files

run_autodock_docking(out_pdb, log=None, prepare_grid=True, num_modes=100, center=None, spacing=0.375, grid_size=None, grid_max_points=None, check_file_out=True, GPU=True)

Run docking using autodock.

Parameters:
  • out_pdb (str) – PDB output name
  • log (str, optional, default=None) – Log ouput name
  • prepare_grid (bool, optional, default=True) – perform grid setup
  • num_modes (int, optional, default=100) – maximum number of binding modes to generate
  • center (list, optional, default=None) – coordinate of the center (x, y, z, Angstroms)
  • grid_size (list, optional, default=None) – size in the docking box (x, y, z, Angstroms)
  • grid_max_points (int, optional, default=None) – max number of grid points per dimension (256 for GPU)
  • check_file_out (bool, optional, default=True) – flag to check or not if file has already been created. If the file is present then the command break.

Object requirement(s):

  • self.lig_pdbqt
  • self.rec_pdbqt

Object field(s) changed:

  • self.dock_pdb
  • self.dock_log
Example:
run_autodock_gpu(out_folder, dock_out_prefix=None, dock_log=None, dock_pdb=None, dock_xml=None, nrun=10, check_file_out=True)

Autodock GPU arguments:

mandatory: -ffile ./input/1stp/derived/1stp_protein.maps.fld -lfile ./input/1stp/derived/1stp_ligand.pdbqt

opyional: -nrun # LGA runs 1 -nev # Score evaluations (max.) per LGA run 2500000 -ngen # Generations (max.) per LGA run 27000 -lsmet Local-search method sw (Solis-Wets) -lsit # Local-search iterations (max.) 300 -psize Population size 150 -mrat Mutation rate 2 (%) -crat Crossover rate 80 (%) -lsrat Local-search rate 6 (%) -trat Tournament (selection) rate 60 (%) -resnam Name for docking output log “docking” -hsym Handle symmetry in RMSD calc. 1

This requires a gpf and associated map files, and ligand pdbqt It creates pose pdb + xml + dlg and smina like _log.txt files

Warning

Difference between GPU and CPU version of autodock logs. Torsional Free Energy is not computed with GPU version.

run_docking(out_pdb, log=None, dock_bin='vina', num_modes=100, energy_range=10, exhaustiveness=16, cpu=None, seed=None, autobox=False, center=None, grid_size=None, min_rmsd_filter=None, scoring=None, check_file_out=True)

Run docking using vina, qvina, qvinaw or smina.

Parameters:
  • out_pdb (str) – PDB output name
  • log (str, optional, default=None) – Log ouput name
  • dock_bin (str, optional, default='vina') – Docking software name (‘vina’, ‘qvina’, ‘qvinaw’, ‘smina’)
  • num_modes (int, optional, default=100) – maximum number of binding modes to generate
  • energy_range (int, optional, default=10) – maximum energy difference between the best binding mode and the worst one displayed (kcal/mol)
  • exhaustiveness (int, optional, default=16) – exhaustiveness of the global search (roughly proportional to time): 1+
  • cpu (int, optional, default=None) – the number of CPUs to use (the default is to try to detect the number of CPUs or, failing that, use 1)
  • seed (int, optional, default=None) – explicit random seed
  • autobox (bool, optional, default=False) – Flag to use ligand to define the docking box
  • center (list, optional, default=None) – coordinate of the center (x, y, z, Angstroms)
  • grid_size (list, optional, default=None) – size in the docking box (x, y, z, Angstroms)
  • check_file_out (bool, optional, default=True) – flag to check or not if file has already been created. If the file is present then the command break.

Object requirement(s):

  • self.lig_pdbqt
  • self.rec_pdbqt

Object field(s) changed:

  • self.dock_pdb
  • self.dock_log
Example:
view_dock(ref_pdb=None)

Return a nglview object to view the object coordinates in a jupyter notebook with the module nglview.

MDAnalysis module is required.

write_out_affinities(fn, affinities)

Save affinities in a file

docking_py.docking.set_log_level(level=20)

setup log verbose level

docking_py.docking.show_log()

To use only with Doctest !!! Redirect logger output to sys.stdout

Module contents

Top-level package for Docking Python.