Home | Trees | Indices | Help |
|
---|
|
object --+ | ProfileHMM
Describes a protein profile Hidden Markov Model. Optional parameters:
Instance Methods | |||
|
|||
|
|||
|
|||
Chain |
|
||
|
|||
list |
|
||
float |
|
||
ProfileHMMSegment |
|
||
|
|||
Structure |
|
||
|
|||
|
|||
Inherited from |
Static Methods | |||
|
Properties | |
A3MAlignment |
alignment Source multiple alignment |
list of HMMLayer |
all_layers A list of layers including start and start_insertion |
AbstractSequence |
consensus Consensus sequence |
SecondaryStructure |
dssp DSSP (calculated) secondary structure |
str |
dssp_solvent Solvent accessibility values |
effective_matches Number of effective matches (NEFF) |
|
bool | emission_pseudocounts |
State |
end Final state (at the end layer) |
EVDParameters |
evd Extreme-value distribution parameters (EVD) |
str |
family Alternative entry ID (FAM) |
bool |
has_structure True if this profile contains structural data |
str |
id Profile entry ID (FILE) |
HMMLayersCollection |
layers List-like access to the HMM's layers |
ProfileLength |
length Profile length |
float |
logbase Base of the logarithm used for score scaling |
str |
name Profile name (NAME) |
bool | pseudocounts |
SecondaryStructure |
psipred PSIPRED (predicted) secondary structure |
collection of Residue |
residues List of representative residues, attached to each layer |
float |
scale Score scaling factor |
ScoreUnits member |
score_units Current score units |
State |
start Start state (at the start layer) |
State |
start_insertion Insertion state at the start layer |
bool | transition_pseudocounts |
str |
version Format version number (HHsearch) |
Inherited from |
Method Details |
x.__init__(...) initializes x; see help(type(x)) for signature
|
Extract the structural information from the HMM. |
Convert emission and transition scores to the specified units.
|
De-serialize an HMM from a file.
|
Extract the emission scores of all match states in the profile. The metric of the emission scores returned depends on the current hmm.score_units setting - you may need to call hmm.convert_scores() to adjust the hmm to your particular needs.
|
Compute the Log-sum-of-odds score between the emission tables of self and other (Soeding 2004). If no observable Match state is found at a given layer, the Insertion state is used instead.
Note: This is not a full implementation of the formula since only emission vectors are involved in the computation and any transition probabilities are ignored. |
Extract a sub-segment of the profile.
|
Serialize this HMM to a file.
|
Extract the structural information from the HMM. |
Dump the profile in HHM format.
|
Property Details |
alignmentSource multiple alignment
|
all_layersA list of layers including start and start_insertion
|
consensusConsensus sequence
|
dsspDSSP (calculated) secondary structure
|
dssp_solventSolvent accessibility values
|
effective_matchesNumber of effective matches (NEFF)
|
emission_pseudocounts
|
endFinal state (at the end layer)
|
evdExtreme-value distribution parameters (EVD)
|
familyAlternative entry ID (FAM)
|
has_structureTrue if this profile contains structural data
|
idProfile entry ID (FILE)
|
layersList-like access to the HMM's layers
|
lengthProfile length
|
logbaseBase of the logarithm used for score scaling
|
nameProfile name (NAME)
|
pseudocounts
|
psipredPSIPRED (predicted) secondary structure
|
residuesList of representative residues, attached to each layer
|
scaleScore scaling factor
|
score_unitsCurrent score units
|
startStart state (at the start layer)
|
start_insertionInsertion state at the start layer
|
transition_pseudocounts
|
versionFormat version number (HHsearch)
|
Home | Trees | Indices | Help |
|
---|
Generated by Epydoc 3.0.1 on Tue Jul 4 20:19:04 2017 | http://epydoc.sourceforge.net |