18 from libc.stdio cimport FILE, fopen, fclose
19 from libc.stdlib cimport free, realloc, malloc
20 from libc.string cimport memcpy
21 from cfreesasa cimport *
24 ShrakeRupley =
'ShrakeRupley'
27 LeeRichards =
'LeeRichards'
36 silent = FREESASA_V_SILENT
39 nowarnings = FREESASA_V_NOWARNINGS
42 normal = FREESASA_V_NORMAL
45 debug = FREESASA_V_DEBUG
49 'algorithm' : LeeRichards,
50 'probe-radius' : freesasa_default_parameters.probe_radius,
51 'n-points' : freesasa_default_parameters.shrake_rupley_n_points,
52 'n-slices' : freesasa_default_parameters.lee_richards_n_slices,
53 'n-threads' : freesasa_default_parameters.n_threads
61 cdef freesasa_parameters _c_param
70 self.
_c_param = freesasa_default_parameters
72 if 'algorithm' in param: self.
setAlgorithm(param[
'algorithm'])
73 if 'probe-radius' in param: self.
setProbeRadius(param[
'probe-radius'])
74 if 'n-points' in param: self.
setNPoints(param[
'n-points'])
75 if 'n-slices' in param: self.
setNSlices(param[
'n-slices'])
76 if 'n-threads' in param: self.
setNThreads(param[
'n-threads'])
79 if not key
in defaultParameters:
80 unknownKeys.append(key)
81 if len(unknownKeys) > 0:
82 raise AssertionError(
'Key(s): ',unknownKeys,
', unknown');
89 if alg == ShrakeRupley:
90 self._c_param.alg = FREESASA_SHRAKE_RUPLEY
91 elif alg == LeeRichards:
92 self._c_param.alg = FREESASA_LEE_RICHARDS
94 raise AssertionError(
"Algorithm '%s' is unknown" % alg)
100 if self._c_param.alg == FREESASA_SHRAKE_RUPLEY:
102 if self._c_param.alg == FREESASA_LEE_RICHARDS:
104 raise Exception(
"No algorithm specified, shouldn't be possible")
111 self._c_param.probe_radius = r
116 return self._c_param.probe_radius
123 self._c_param.shrake_rupley_n_points = n
128 return self._c_param.shrake_rupley_n_points
135 self._c_param.lee_richards_n_slices = n
140 return self._c_param.lee_richards_n_slices
148 self._c_param.n_threads = n
154 return self._c_param.n_threads
157 def _get_address(self, size_t ptr2ptr):
158 cdef freesasa_parameters **p = <freesasa_parameters**> ptr2ptr
167 cdef freesasa_result* _c_result
182 return self._c_result.n_atoms
192 return self._c_result.total
201 assert(i < self._c_result.n_atoms)
202 return self._c_result.sasa[i]
204 def _get_address(self, size_t ptr2ptr):
205 cdef freesasa_result **p = <freesasa_result**> ptr2ptr
223 cdef freesasa_classifier* _c_classifier
238 if fileName
is not None:
239 config = fopen(fileName,
'r')
241 raise IOError(
"File '%s' could not be opened." % fileName)
245 raise Exception(
"Error parsing configuration in '%s'." % fileName)
263 classIndex = self._c_classifier.sasa_class(residueName,atomName,self.
_c_classifier)
264 if classIndex
is FREESASA_WARN:
266 return self._c_classifier.class2str(classIndex,self.
_c_classifier)
278 return self._c_classifier.radius(residueName,atomName,self.
_c_classifier)
280 def _get_address(self, size_t ptr2ptr):
281 cdef freesasa_classifier **p = <freesasa_classifier**> ptr2ptr
294 cdef freesasa_structure* _c_structure
298 defaultOptions = {
'hetatm' :
False,
300 'join-models' :
False,
301 'skip-unknown' :
False,
302 'halt-at-unknown' :
False
321 def __init__(self,fileName=None,classifier=None,
322 options = defaultOptions):
330 input = fopen(fileName,
'r')
332 raise IOError(
"File '%s' could not be opened." % fileName)
333 structure_options = Structure._get_structure_options(options)
338 raise Exception(
"Error reading '%s'." % fileName)
343 if (classifier
is not None):
368 def addAtom(self, atomName, residueName, residueNumber, chainLabel, x, y, z):
369 if (type(residueNumber)
is str):
370 resnum = residueNumber
371 elif (type(residueNumber)
is int):
372 resnum =
"%d" % residueNumber
374 raise Exception(
"Residue-number invalid, must be either string or number")
375 cdef const char *label = chainLabel
377 residueName, resnum, label[0],
379 assert(ret != FREESASA_FAIL)
402 assert len(radiusArray) == n
403 cdef double *r = <double *>malloc(sizeof(double)*n)
404 assert(r
is not NULL)
406 r[i] = radiusArray[i]
423 assert(i >= 0
and i < self.
nAtoms())
426 assert(r
is not NULL)
434 assert(i >= 0
and i < self.
nAtoms())
443 assert(i >= 0
and i < self.
nAtoms())
452 assert(i >= 0
and i < self.
nAtoms())
461 assert(i >= 0
and i < self.
nAtoms())
473 assert(i >= 0
and i < self.
nAtoms())
476 return [coord[3*i], coord[3*i+1], coord[3*i+2]]
479 def _get_structure_options(param):
483 knownOptions = {
'hetatm',
'hydrogen',
'join-models',
'separate-models',
484 'separate-chains',
'skip-unknown',
'halt-at-unknown'}
487 if not key
in knownOptions:
488 unknownOptions.append(key)
489 if len(unknownOptions) > 0:
490 raise AssertionError(
"Option(s): ",unknownOptions,
" unknown.");
493 if 'hetatm' in param
and param[
'hetatm']:
494 options |= FREESASA_INCLUDE_HETATM
495 if 'hydrogen' in param
and param[
'hydrogen']:
496 options |= FREESASA_INCLUDE_HYDROGEN
497 if 'join-models' in param
and param[
'join-models']:
498 options |= FREESASA_JOIN_MODELS
499 if 'separate-models' in param
and param[
'separate-models']:
500 options |= FREESASA_SEPARATE_MODELS
501 if 'separate-chains' in param
and param[
'separate-chains']:
502 options |= FREESASA_SEPARATE_CHAINS
503 if 'skip-unknown' in param
and param[
'skip-unknown']:
504 options |= FREESASA_SKIP_UNKNOWN
505 if 'halt-at-unknown' in param
and param[
'halt-at-unknown']:
506 options |= FREESASA_HALT_AT_UNKNOWN
509 def _get_address(self, size_t ptr2ptr):
510 cdef freesasa_structure **p = <freesasa_structure**> ptr2ptr
513 def _set_address(self, size_t ptr2ptr):
514 cdef freesasa_structure **p = <freesasa_structure**> ptr2ptr
524 defaultStructureArrayOptions = {
527 'separate-chains' :
True,
528 'separate-models' :
False
550 options = defaultStructureArrayOptions,
552 assert fileName
is not None
554 assert((
'separate-chains' in options
and options[
'separate-chains']
is True)
555 or (
'separate-models' in options
and options[
'separate-models']
is True))
556 structure_options = Structure._get_structure_options(options)
558 input = fopen(fileName,
'r')
560 raise IOError(
"File '%s' could not be opened." % fileName)
565 raise Exception(
"Problems reading structures in '%s'." % fileName)
569 structures[-1]._set_address(<size_t> &sArray[i])
570 if classifier
is not None:
571 structures[-1].setRadiiWithClassifier(classifier)
581 def calc(structure,parameters=None):
582 cdef const freesasa_parameters *p = NULL
583 cdef const freesasa_structure *s = NULL
584 if parameters
is not None: parameters._get_address(<size_t>&p)
585 structure._get_address(<size_t>&s)
588 if result._c_result
is NULL:
589 raise Exception(
"Error calculating SASA.")
600 if classifier
is None:
603 for i
in range(0,structure.nAtoms()):
604 name = classifier.classify(structure.residueName(i),structure.atomName(i))
607 ret[name] += result.atomArea(i)
621 cdef freesasa_structure *s
622 cdef freesasa_result *r
624 cdef char *name = <char*>malloc(FREESASA_MAX_SELECTION_NAME+1);
625 structure._get_address(<size_t> &s)
626 result._get_address(<size_t> &r)
630 if ret == FREESASA_FAIL:
631 raise Exception(
"Error parsing '%s'" % cmd)
640 assert(verbosity
in [silent, nowarnings, normal])
664 if (classifier
is None):
666 optbitfield = Structure._get_structure_options(options)
668 atoms = bioPDBStructure.get_atoms()
672 hetflag, resseq, icode = r.get_id()
674 if (hetflag
is not ' ' and not (optbitfield & FREESASA_INCLUDE_HETATM)):
680 if (classifier.classify(r.get_resname(), a.get_fullname())
is 'Unknown'):
681 if (optbitfield & FREESASA_SKIP_UNKNOWN):
683 if (optbitfield & FREESASA_HALT_AT_UNKNOWN):
684 raise Exception(
"Halting at unknown atom")
686 structure.addAtom(a.get_fullname(), r.get_resname(), resseq, c.get_id(),
689 structure.setRadiiWithClassifier(classifier)
711 def calcBioPDB(bioPDBStructure, parameters = Parameters(),
712 classifier =
None, options = Structure.defaultOptions):
714 result =
calc(structure, parameters)
716 return result, sasa_classes
def structureFromBioPDB
Create a freesasa structure from a Bio.PDB structure.
def setProbeRadius(self, r)
Set probe radius.
const double * freesasa_structure_coord_array(const freesasa_structure *structure)
Get array of coordinates.
Assigns class and radius to atom by residue and atom name.
Stores parameter values to be used by calculation.
def calc
Calculate SASA of Structure.
def probeRadius(self)
Get probe radius.
def getVerbosity()
Get global verbosity.
def __dealloc__(self)
The destructor.
def selectArea(commands, structure, result)
Sum SASA result over a selection of atoms.
def atomArea(self, i)
SASA for a given atom.
def __init__
Initializes Parameters object.
def nSlices(self)
Get the number of slices per atom in Lee & Richards algorithm.
def structureArray
Create array of structures from PDB file.
def radius(self, residueName, atomName)
Radius of atom.
const char * freesasa_structure_atom_res_number(const freesasa_structure *structure, int i)
Get residue number.
Represents a protein structure, including its atomic radii.
def addAtom(self, atomName, residueName, residueNumber, chainLabel, x, y, z)
Add atom to structure.
char freesasa_structure_atom_chain(const freesasa_structure *structure, int i)
Get chain label.
freesasa_verbosity freesasa_get_verbosity(void)
Get the current verbosity level.
freesasa_structure * freesasa_structure_new(void)
Allocate empty structure.
freesasa_result * freesasa_calc_structure(const freesasa_structure *structure, const freesasa_parameters *parameters)
Calculates SASA based on a given structure.
int freesasa_structure_n(const freesasa_structure *structure)
Get number of atoms.
def setRadiiWithClassifier(self, classifier)
Assign radii to atoms in structure using a classifier.
def radius(self, i)
Radius of atom.
def nAtoms(self)
Number of atoms.
def __cinit__(self)
The constructor.
def calcBioPDB
Calc SASA from Bio.PDB structure.
def atomName(self, i)
Get atom name.
def setRadii(self, radiusArray)
Set atomic radii from an array.
def nThreads(self)
Get the number of threads to use in calculations.
const double * freesasa_structure_radius(const freesasa_structure *structure)
Returns a pointer to an array of the radii of each atom.
def classifyResults
Break SASA result down into classes.
def coord(self, i)
Get coordinates of given atom.
int freesasa_structure_add_atom(freesasa_structure *structure, const char *atom_name, const char *residue_name, const char *residue_number, char chain_label, double x, double y, double z)
Add individual atom to structure using default behavior.
def nAtoms(self)
Number of atoms in the results.
def residueNumber(self, i)
Get residue number for given atom.
def algorithm(self)
Get algorithm.
void freesasa_classifier_free(freesasa_classifier *classifier)
Frees a classifier object.
def setVerbosity(verbosity)
Set global verbosity.
int freesasa_select_area(const char *command, char *name, double *area, const freesasa_structure *structure, const freesasa_result *result)
Get area of a selection.
int freesasa_set_verbosity(freesasa_verbosity v)
Set the global verbosity level.
def __cinit__
Constructor.
freesasa_classifier * freesasa_classifier_from_file(FILE *file)
Generate a classifier from a config-file.
def __dealloc__(self)
The destructor.
def nPoints(self)
Get number of test points in Shrake & Rupley algorithm.
def setNPoints(self, n)
Set number of test points in Shrake & Rupley algorithm.
const char * freesasa_structure_atom_res_name(const freesasa_structure *structure, int i)
Get residue name.
def setAlgorithm(self, alg)
Set algorithm.
Stores results from SASA calculation.
def residueName(self, i)
Get residue name of given atom.
void freesasa_structure_set_radius(freesasa_structure *structure, const double *radii)
Override the radii set when the structure was initialized.
def setNThreads(self, n)
Set the number of threads to use in calculations.
freesasa_structure ** freesasa_structure_array(FILE *pdb, int *n, const freesasa_classifier *classifier, int options)
Init array of structures from PDB.
def classify(self, residueName, atomName)
Class of atom.
def totalArea(self)
Total SASA.
def chainLabel(self, i)
Get chain label for given atom.
freesasa_structure * freesasa_structure_from_pdb(FILE *pdb, const freesasa_classifier *classifier, int options)
Init structure with coordinates from pdb-file.
def setNSlices(self, n)
Set the number of slices per atom in Lee & Richards algorithm.
const char * freesasa_structure_atom_name(const freesasa_structure *structure, int i)
Get atom name.
void freesasa_structure_free(freesasa_structure *structure)
Free structure.
def __dealloc__(self)
The destructor.
void freesasa_result_free(freesasa_result *result)
Frees a freesasa_result object.