Bio.PDB.SASA module

Calculation of solvent accessible surface areas for Bio.PDB entities.

Uses the “rolling ball” algorithm developed by Shrake & Rupley algorithm, which uses a sphere (of equal radius to a solvent molecule) to probe the surface of the molecule.

Reference:

Shrake, A; Rupley, JA. (1973). J Mol Biol “Environment and exposure to solvent of protein atoms. Lysozyme and insulin”.

class Bio.PDB.SASA.ShrakeRupley(probe_radius=1.40, n_points=100, radii_dict=None)

Bases: object

Calculates SASAs using the Shrake-Rupley algorithm.

__init__(probe_radius=1.40, n_points=100, radii_dict=None)

Initialize the class.

Parameters:
  • probe_radius (float) – radius of the probe in A. Default is 1.40, roughly the radius of a water molecule.

  • n_points (int) – resolution of the surface of each atom. Default is 100. A higher number of points results in more precise measurements, but slows down the calculation.

  • radii_dict (dict) – user-provided dictionary of atomic radii to use in the calculation. Values will replace/complement those in the default ATOMIC_RADII dictionary.

>>> sr = ShrakeRupley()
>>> sr = ShrakeRupley(n_points=960)
>>> sr = ShrakeRupley(radii_dict={"O": 3.1415})
compute(entity, level='A')

Calculate surface accessibility surface area for an entity.

The resulting atomic surface accessibility values are attached to the .sasa attribute of each entity (or atom), depending on the level. For example, if level=”R”, all residues will have a .sasa attribute. Atoms will always be assigned a .sasa attribute with their individual values.

Parameters:
  • entity (Bio.PDB.Entity) – input entity.

  • level – the level at which ASA values are assigned, which can be one of “A” (Atom), “R” (Residue), “C” (Chain), “M” (Model), or “S” (Structure). The ASA value of an entity is the sum of all ASA values of its children. Defaults to “A”.

>>> from Bio.PDB import PDBParser
>>> from Bio.PDB.SASA import ShrakeRupley
>>> p = PDBParser(QUIET=1)
>>> # This assumes you have a local copy of 1LCD.pdb in a directory called "PDB"
>>> struct = p.get_structure("1LCD", "PDB/1LCD.pdb")
>>> sr = ShrakeRupley()
>>> sr.compute(struct, level="S")
>>> print(round(struct.sasa, 2))
7053.43
>>> print(round(struct[0]["A"][11]["OE1"].sasa, 2))
9.64