FreeSASA  1.1
Open source SASA calculations
View on GitHub
Public Member Functions | Static Public Attributes
freesasa.Structure Class Reference

Represents a protein structure, including its atomic radii. More...

Public Member Functions

def __init__
 Constructor. More...
 
def addAtom (self, atomName, residueName, residueNumber, chainLabel, x, y, z)
 Add atom to structure. More...
 
def setRadiiWithClassifier (self, classifier)
 Assign radii to atoms in structure using a classifier. More...
 
def setRadii (self, radiusArray)
 Set atomic radii from an array. More...
 
def nAtoms (self)
 Number of atoms. More...
 
def radius (self, i)
 Radius of atom. More...
 
def atomName (self, i)
 Get atom name. More...
 
def residueName (self, i)
 Get residue name of given atom. More...
 
def residueNumber (self, i)
 Get residue number for given atom. More...
 
def chainLabel (self, i)
 Get chain label for given atom. More...
 
def coord (self, i)
 Get coordinates of given atom. More...
 
def __dealloc__ (self)
 The destructor.
 

Static Public Attributes

dictionary defaultOptions
 By default ignore HETATM, Hydrogens, only use first model. More...
 

Detailed Description

Represents a protein structure, including its atomic radii.

Initialized from PDB-file. Calculates atomic radii using default classifier, or custom one provided as argument to initalizer

Wraps the C struct freesasa_structure.

Since it is intended to be a static structure the word 'get' is omitted in the getter-functions.

Definition at line 293 of file freesasa.pyx.

Constructor & Destructor Documentation

def freesasa.Structure.__init__ (   self,
  fileName = None,
  classifier = None,
  options = defaultOptions 
)

Constructor.

If PDB file is provided, the structure will be constructed based on the file. If not, this simply initializes an empty structure and the other arguments are ignored. In this case atoms will have to be added manually using addAtom().

Parameters
fileName(str) PDB file (if None empty structure generated).
classifierAn optional Classifier to calculate atomic radii, uses default if none provided
optionsspecify which atoms and models to include
Exceptions
IOErrorProblem opening/reading file.
ExceptionProblem parsing PDB file or calculating atomic radii.
ExceptionIf option 'halt-at-unknown' selected and unknown atom encountered.

Definition at line 322 of file freesasa.pyx.

Member Function Documentation

def freesasa.Structure.addAtom (   self,
  atomName,
  residueName,
  residueNumber,
  chainLabel,
  x,
  y,
  z 
)

Add atom to structure.

This function is meant to be used if the structure was not initialized from a PDB. Default radii will be assigned to each atom. This can be overriden by calling setRadiiWithClassifier() afterwards.

There are no restraints on string lengths for the arguments, but the atom won't be added if the default classifier doesn't recognize the atom and also cannot deduce its element from the atom name.

Parameters
atomName(str) atom name (e.g. "CA")
residueName(str) residue name (e.g. "ALA")
residueNumber(str or int) residue number (e.g. '12') or integer. Some PDBs have residue-numbers that aren't regular numbers. Therefore treated as a string primarily.
chainLabel(str) 1-character string with chain label (e.g. 'A')
x,y,z(float) coordinates
Exceptions
ExceptionResidue-number invalid

Definition at line 368 of file freesasa.pyx.

def freesasa.Structure.setRadiiWithClassifier (   self,
  classifier 
)

Assign radii to atoms in structure using a classifier.

Parameters
classifierA classifier to use to calculate radii
Exceptions
AssertionErrorif structure not properly initialized

Definition at line 385 of file freesasa.pyx.

def freesasa.Structure.setRadii (   self,
  radiusArray 
)

Set atomic radii from an array.

Parameters
radiusArrayArray of atomic radii in Ångström, should have nAtoms() elements.
Exceptions
AssertionErrorif radiusArray has wrong dimension, structure not properly initialized, or if the array contains negative radii (not properly classified?)

Definition at line 399 of file freesasa.pyx.

def freesasa.Structure.nAtoms (   self)

Number of atoms.

Returns
Number of atoms
Exceptions
AssertionErrorif not properly initialized

Definition at line 414 of file freesasa.pyx.

def freesasa.Structure.radius (   self,
  i 
)

Radius of atom.

Parameters
i(int) Index of atom.
Returns
Radius in Å.
Exceptions
AssertionErrorif index out of bounds, object not properly initalized.

Definition at line 422 of file freesasa.pyx.

def freesasa.Structure.atomName (   self,
  i 
)

Get atom name.

Parameters
i(int) Atom index.
Returns
Atom name as 4-character string.
Exceptions
AssertionErrorif index out of range or Structure not properly initialized.

Definition at line 433 of file freesasa.pyx.

def freesasa.Structure.residueName (   self,
  i 
)

Get residue name of given atom.

Parameters
i(int) Atom index.
Returns
Residue name as 3-character string.
Exceptions
AssertionErrorif index out of range or Structure not properly initialized

Definition at line 442 of file freesasa.pyx.

def freesasa.Structure.residueNumber (   self,
  i 
)

Get residue number for given atom.

Parameters
i(int) Atom index.
Returns
Residue number as 4-character string
Exceptions
AssertionErrorif index out of range or Structure not properly initialized

Definition at line 451 of file freesasa.pyx.

def freesasa.Structure.chainLabel (   self,
  i 
)

Get chain label for given atom.

Parameters
i(int) Atom index.
Returns
Chain label as 1-character string
Exceptions
AssertionErrorif index out of range or Structure not properly initialized

Definition at line 460 of file freesasa.pyx.

def freesasa.Structure.coord (   self,
  i 
)

Get coordinates of given atom.

Parameters
i(int) Atom index.
Returns
array of x, y, and z coordinates
Exceptions
AssertionErrorif index out of range or Structure not properly initialized

Definition at line 472 of file freesasa.pyx.

Field Documentation

dictionary freesasa.Structure.defaultOptions
static
Initial value:
1 = { 'hetatm' : False,
2  'hydrogen' : False,
3  'join-models' : False,
4  'skip-unknown' : False,
5  'halt-at-unknown' : False
6  }

By default ignore HETATM, Hydrogens, only use first model.

For unknown atoms try to guess the radius, if this fails, assign radius 0 (to allow changing the radius later).

Definition at line 298 of file freesasa.pyx.


The documentation for this class was generated from the following file: