Skip to main content

Calculate solvent accessible surface areas of proteins

Project description

.. currentmodule:: freesasa

Introduction
============

This package provides Python bindings for the `FreeSASA C Library
<https://github.com/mittinatten/freesasa>`_.

It can be installed using

.. code::

pip install freesasa

Binaries are available for Python 2.7, 3.4, 3.5 and 3.6 for Mac OS X
and Windows, in addition to the source distribution.


Basic calculations
------------------

Using defaults everywhere a simple calculation can be carried out as
follows (assuming the file ``1ubq.pdb`` is available)

.. code:: python

import freesasa

structure = freesasa.Structure("1ubq.pdb")
result = freesasa.calc(structure)
area_classes = freesasa.classifyResults(result, structure)

print "Total : %.2f A2" % result.totalArea()
for key in area_classes:
print key, ": %.2f A2" % area_classes[key]

Which would give the following output

.. code::

Total : 4804.06 A2
Polar : 2504.22 A2
Apolar : 2299.84 A2

The following does a high precision L&R calculation

.. code:: python

result = freesasa.calc(structure,
freesasa.Parameters({'algorithm' : freesasa.LeeRichards,
'n-slices' : 100}))

Using the results from a calculation we can also integrate SASA over a selection of
atoms, using a subset of the Pymol `selection syntax`_:

.. _selection syntax: http://freesasa.github.io/doxygen/Selection.html

.. code:: python

selections = freesasa.selectArea(('alanine, resn ala', 'r1_10, resi 1-10'),
structure, result)
for key in selections:
print key, ": %.2f A2" % selections[key]

which gives the output

.. code::

alanine : 120.08 A2
r1_10 : 634.31 A2

Customizing atom classification
-------------------------------

This uses the NACCESS parameters (the file ``naccess.config`` is
available in the ``share/`` directory of the repository).

.. code:: python

classifier = freesasa.Classifier("naccess.config")
structure = freesasa.Structure("1ubq.pdb", classifier)
result = freesasa.calc(structure)
area_classes = freesasa.classifyResults(result, structure, classifier)

Classification can be customized also by extending the :py:class:`.Classifier`
interface. The code below is an illustration of a classifier that
classes nitrogens separately, and assigns radii based on element only
(and crudely).

.. code:: python

import freesasa
import re

class DerivedClassifier(Classifier):
def classify(self, residueName, atomName):
if re.match('\s*N', atomName):
return 'Nitrogen'
return 'Not-nitrogen'

def radius(self, residueName, atomName):
if re.match('\s*N',atomName): # Nitrogen
return 1.6
if re.match('\s*C',atomName): # Carbon
return 1.7
if re.match('\s*O',atomName): # Oxygen
return 1.4
if re.match('\s*S',atomName): # Sulfur
return 1.8
return 0; # everything else (Hydrogen, etc)

classifier = DerivedClassifier()

# use the DerivedClassifier to calculate atom radii
structure = freesasa.Structure("1ubq.pdb", classifier)
result = freesasa.calc(structure)

# use the DerivedClassifier to classify atoms
area_classes = freesasa.classifyResults(result,structure,classifier)

Of course, this example is somewhat contrived, if we only want the
integrated area of nitrogen atoms, the simpler choice would be

.. code:: python

selection = freesasa.selectArea('nitrogen, symbol n', structure, result)


However, extending :py:class:`.Classifier`, as illustrated above, allows
classification to arbitrary complexity and also lets us redefine the
radii used in the calculation.

Bio.PDB
-------

FreeSASA can also calculate the SASA of a ``Bio.PDB`` structure from BioPython

.. code:: python

from Bio.PDB import PDBParser
parser = PDBParser()
structure = parser.get_structure("Ubiquitin", "1ubq.pdb")
result, sasa_classes = freesasa.calcBioPDB(structure)

If one needs more control over the analysis the structure can be
converted to a :py:class:`.Structure` using :py:func:`.structureFromBioPDB()`
and the calculation can be performed the normal way using this
structure.

Project details


Download files

Download the file for your platform. If you're not sure which to choose, learn more about installing packages.

Source Distribution

freesasa-2.0.3.tar.gz (177.8 kB view hashes)

Uploaded Source

Built Distributions

freesasa-2.0.3-py3.6-macosx-10.12-x86_64.egg (136.9 kB view hashes)

Uploaded Source

freesasa-2.0.3-py3.6-linux-x86_64.egg (487.9 kB view hashes)

Uploaded Source

freesasa-2.0.3-py3.5-macosx-10.12-x86_64.egg (138.4 kB view hashes)

Uploaded Source

freesasa-2.0.3-py3.5-linux-x86_64.egg (469.6 kB view hashes)

Uploaded Source

freesasa-2.0.3-py3.4-macosx-10.12-x86_64.egg (137.7 kB view hashes)

Uploaded Source

freesasa-2.0.3-py3.4-linux-x86_64.egg (474.6 kB view hashes)

Uploaded Source

freesasa-2.0.3-py2.7-macosx-10.12-x86_64.egg (137.2 kB view hashes)

Uploaded Source

freesasa-2.0.3-py2.7-linux-x86_64.egg (447.8 kB view hashes)

Uploaded Source

freesasa-2.0.3-cp36-cp36m-win_amd64.whl (240.9 kB view hashes)

Uploaded CPython 3.6m Windows x86-64

freesasa-2.0.3-cp36-cp36m-win32.whl (224.4 kB view hashes)

Uploaded CPython 3.6m Windows x86

freesasa-2.0.3-cp36-cp36m-macosx_10_12_x86_64.whl (263.9 kB view hashes)

Uploaded CPython 3.6m macOS 10.12+ x86-64

freesasa-2.0.3-cp35-cp35m-win_amd64.whl (239.1 kB view hashes)

Uploaded CPython 3.5m Windows x86-64

freesasa-2.0.3-cp35-cp35m-win32.whl (222.9 kB view hashes)

Uploaded CPython 3.5m Windows x86

freesasa-2.0.3-cp35-cp35m-macosx_10_12_x86_64.whl (265.5 kB view hashes)

Uploaded CPython 3.5m macOS 10.12+ x86-64

freesasa-2.0.3-cp34-cp34m-win_amd64.whl (240.7 kB view hashes)

Uploaded CPython 3.4m Windows x86-64

freesasa-2.0.3-cp34-cp34m-win32.whl (225.7 kB view hashes)

Uploaded CPython 3.4m Windows x86

freesasa-2.0.3-cp34-cp34m-macosx_10_12_x86_64.whl (264.7 kB view hashes)

Uploaded CPython 3.4m macOS 10.12+ x86-64

freesasa-2.0.3-cp27-cp27m-win_amd64.whl (241.2 kB view hashes)

Uploaded CPython 2.7m Windows x86-64

freesasa-2.0.3-cp27-cp27m-win32.whl (223.7 kB view hashes)

Uploaded CPython 2.7m Windows x86

freesasa-2.0.3-cp27-cp27m-macosx_10_12_x86_64.whl (264.3 kB view hashes)

Uploaded CPython 2.7m macOS 10.12+ x86-64

Supported by

AWS AWS Cloud computing and Security Sponsor Datadog Datadog Monitoring Fastly Fastly CDN Google Google Download Analytics Microsoft Microsoft PSF Sponsor Pingdom Pingdom Monitoring Sentry Sentry Error logging StatusPage StatusPage Status page