RNAlib-2.4.14
Global Partition Function and Equilibrium Probabilities

Variations of the global partition function algorithm. More...

Detailed Description

Variations of the global partition function algorithm.

We provide implementations of the global partition function algorithm for

+ Collaboration diagram for Global Partition Function and Equilibrium Probabilities:

Modules

 Computing Partition Functions of a Distance Based Partitioning
 Compute the partition function and stochastically sample secondary structures for a partitioning of the secondary structure space according to the base pair distance to two fixed reference structures.
 
 Deprecated Interface for Global Partition Function Computation
 

Files

file  part_func.h
 Partition function implementations.
 

Data Structures

struct  vrna_dimer_pf_s
 Data structure returned by vrna_pf_dimer() More...
 

Functions

void vrna_pf_dimer_probs (double FAB, double FA, double FB, vrna_ep_t *prAB, const vrna_ep_t *prA, const vrna_ep_t *prB, int Alength, const vrna_exp_param_t *exp_params)
 Compute Boltzmann probabilities of dimerization without homodimers. More...
 
double vrna_pr_structure (vrna_fold_compound_t *fc, const char *structure)
 Compute the equilibrium probability of a particular secondary structure. More...
 
vrna_ep_tvrna_plist_from_probs (vrna_fold_compound_t *vc, double cut_off)
 Create a vrna_ep_t from base pair probability matrix. More...
 

Base pair related probability computations

double vrna_mean_bp_distance_pr (int length, FLT_OR_DBL *pr)
 Get the mean base pair distance in the thermodynamic ensemble from a probability matrix. More...
 
double vrna_mean_bp_distance (vrna_fold_compound_t *vc)
 Get the mean base pair distance in the thermodynamic ensemble. More...
 
double vrna_ensemble_defect (vrna_fold_compound_t *fc, const char *structure)
 Compute the Ensemble Defect for a given target structure. More...
 
vrna_ep_tvrna_stack_prob (vrna_fold_compound_t *vc, double cutoff)
 Compute stacking probabilities. More...
 
double * vrna_positional_entropy (vrna_fold_compound_t *fc)
 Compute a vector of positional entropies. More...
 

Basic global partition function interface

float vrna_pf (vrna_fold_compound_t *vc, char *structure)
 Compute the partition function $Q$ for a given RNA sequence, or sequence alignment. More...
 
vrna_dimer_pf_t vrna_pf_dimer (vrna_fold_compound_t *vc, char *structure)
 Calculate partition function and base pair probabilities of nucleic acid/nucleic acid dimers. More...
 

Simplified global partition function computation using sequence(s) or multiple sequence alignment(s)

float vrna_pf_fold (const char *sequence, char *structure, vrna_ep_t **pl)
 Compute Partition function $Q$ (and base pair probabilities) for an RNA sequence using a comparative method. More...
 
float vrna_pf_circfold (const char *sequence, char *structure, vrna_ep_t **pl)
 Compute Partition function $Q$ (and base pair probabilities) for a circular RNA sequences using a comparative method. More...
 
float vrna_pf_alifold (const char **sequences, char *structure, vrna_ep_t **pl)
 Compute Partition function $Q$ (and base pair probabilities) for an RNA sequence alignment using a comparative method. More...
 
float vrna_pf_circalifold (const char **sequences, char *structure, vrna_ep_t **pl)
 Compute Partition function $Q$ (and base pair probabilities) for an alignment of circular RNA sequences using a comparative method. More...
 
vrna_dimer_pf_t vrna_pf_co_fold (const char *seq, char *structure, vrna_ep_t **pl)
 Calculate partition function and base pair probabilities of nucleic acid/nucleic acid dimers. More...
 

Data Structure Documentation

struct vrna_dimer_pf_s

Data structure returned by vrna_pf_dimer()

Data Fields

double F0AB
 Null model without DuplexInit.
 
double FAB
 all states with DuplexInit correction
 
double FcAB
 true hybrid states only
 
double FA
 monomer A
 
double FB
 monomer B
 

Function Documentation

double vrna_mean_bp_distance_pr ( int  length,
FLT_OR_DBL pr 
)

#include <ViennaRNA/equilibrium_probs.h>

Get the mean base pair distance in the thermodynamic ensemble from a probability matrix.

\[ <d> = \sum_{a,b} p_a p_b d(S_a,S_b) \]

this can be computed from the pair probs $ p_{ij} $ as

\[ <d> = \sum_{ij} p_{ij}(1-p_{ij}) \]

Parameters
lengthThe length of the sequence
prThe matrix containing the base pair probabilities
Returns
The mean pair distance of the structure ensemble
double vrna_mean_bp_distance ( vrna_fold_compound_t vc)

#include <ViennaRNA/equilibrium_probs.h>

Get the mean base pair distance in the thermodynamic ensemble.

\[ <d> = \sum_{a,b} p_a p_b d(S_a,S_b) \]

this can be computed from the pair probs $p_{ij}$ as

\[ <d> = \sum_{ij} p_{ij}(1-p_{ij}) \]

Parameters
vcThe fold compound data structure
Returns
The mean pair distance of the structure ensemble
SWIG Wrapper Notes:
This function is attached as method mean_bp_distance() to objects of type fold_compound
double vrna_ensemble_defect ( vrna_fold_compound_t fc,
const char *  structure 
)

#include <ViennaRNA/equilibrium_probs.h>

Compute the Ensemble Defect for a given target structure.

Given a target structure $s$, compute the average dissimilarity of a randomly drawn structure from the ensemble, i.e.:

\[ ED(s) = 1 - \frac{1}{n} \sum_{ij, (i,j) \in s} p_{ij} - \frac{1}{n} \sum_{i}(1 - s_i)q_i \]

with sequence length $n$, the probability $p_{ij}$ of a base pair $(i,j)$, the probability $q_i = 1 - \sum_j p_{ij}$ of nucleotide $i$ being unpaired, and the indicator variable $s_i = 1$ if $\exists (i,j) \in s$, and $s_i = 0$ otherwise.

Precondition
The vrna_fold_compound_t input parameter fc must contain a valid base pair probability matrix. This means that partition function and base pair probabilities must have been computed using fc before execution of this function!
See also
vrna_pf(), vrna_pairing_probs()
Parameters
fcA fold_compound with pre-computed base pair probabilities
structureA target structure in dot-bracket notation
Returns
The ensemble defect with respect to the target structure, or -1. upon failure, e.g. pre-conditions are not met
SWIG Wrapper Notes:
This function is attached as method ensemble_defect() to objects of type fold_compound
vrna_ep_t* vrna_stack_prob ( vrna_fold_compound_t vc,
double  cutoff 
)

#include <ViennaRNA/equilibrium_probs.h>

Compute stacking probabilities.

For each possible base pair $(i,j)$, compute the probability of a stack $(i,j)$, $(i+1, j-1)$.

Parameters
vcThe fold compound data structure with precomputed base pair probabilities
cutoffA cutoff value that limits the output to stacks with $ p > \textrm{cutoff} $.
Returns
A list of stacks with enclosing base pair $(i,j)$ and probabiltiy $ p $
void vrna_pf_dimer_probs ( double  FAB,
double  FA,
double  FB,
vrna_ep_t prAB,
const vrna_ep_t prA,
const vrna_ep_t prB,
int  Alength,
const vrna_exp_param_t exp_params 
)

#include <ViennaRNA/equilibrium_probs.h>

Compute Boltzmann probabilities of dimerization without homodimers.

Given the pair probabilities and free energies (in the null model) for a dimer AB and the two constituent monomers A and B, compute the conditional pair probabilities given that a dimer AB actually forms. Null model pair probabilities are given as a list as produced by vrna_plist_from_probs(), the dimer probabilities 'prAB' are modified in place.

Parameters
FABfree energy of dimer AB
FAfree energy of monomer A
FBfree energy of monomer B
prABpair probabilities for dimer
prApair probabilities monomer
prBpair probabilities monomer
AlengthLength of molecule A
exp_paramsThe precomputed Boltzmann factors
double vrna_pr_structure ( vrna_fold_compound_t fc,
const char *  structure 
)

#include <ViennaRNA/equilibrium_probs.h>

Compute the equilibrium probability of a particular secondary structure.

The probability $p(s)$ of a particular secondary structure $s$ can be computed as

\[ p(s) = \frac{exp(-\beta E(s)}{Z} \]

from the structures free energy $E(s)$ and the partition function

\[ Z = \sum_s exp(-\beta E(s)),\quad\mathrm{with}\quad\beta = \frac{1}{RT} \]

where $R$ is the gas constant and $T$ the thermodynamic temperature.

Precondition
The fold compound fc must have went through a call to vrna_pf() to fill the dynamic programming matrices with the corresponding partition function.
Parameters
fcThe fold compound data structure with precomputed partition function
structureThe secondary structure to compute the probability for in dot-bracket notation
Returns
The probability of the input structure (range $[0:1]$)
SWIG Wrapper Notes:
This function is attached as method pr_structure() to objects of type fold_compound
float vrna_pf ( vrna_fold_compound_t vc,
char *  structure 
)

#include <ViennaRNA/part_func.h>

Compute the partition function $Q$ for a given RNA sequence, or sequence alignment.

If structure is not a NULL pointer on input, it contains on return a string consisting of the letters " . , | { } ( ) " denoting bases that are essentially unpaired, weakly paired, strongly paired without preference, weakly upstream (downstream) paired, or strongly up- (down-)stream paired bases, respectively. If the model's compute_bpp is set to 0 base pairing probabilities will not be computed (saving CPU time), otherwise after calculations took place pr will contain the probability that bases i and j pair.

Note
This function is polymorphic. It accepts vrna_fold_compound_t of type VRNA_FC_TYPE_SINGLE, and VRNA_FC_TYPE_COMPARATIVE.
This function may return INF / 100. in case of contradicting constraints or numerical over-/underflow. In the latter case, a corresponding warning will be issued to stdout.
See also
vrna_fold_compound_t, vrna_fold_compound(), vrna_pf_fold(), vrna_pf_circfold(), vrna_fold_compound_comparative(), vrna_pf_alifold(), vrna_pf_circalifold(), vrna_db_from_probs(), vrna_exp_params(), vrna_aln_pinfo()
Parameters
[in,out]vcThe fold compound data structure
[in,out]structureA pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
Returns
The ensemble free energy $G = -RT \cdot \log(Q) $ in kcal/mol
SWIG Wrapper Notes:
This function is attached as method pf() to objects of type fold_compound
vrna_dimer_pf_t vrna_pf_dimer ( vrna_fold_compound_t vc,
char *  structure 
)

#include <ViennaRNA/part_func.h>

Calculate partition function and base pair probabilities of nucleic acid/nucleic acid dimers.

This is the cofold partition function folding.

Note
This function may return INF / 100. for the FA, FB, FAB, F0AB members of the output data structure in case of contradicting constraints or numerical over-/underflow. In the latter case, a corresponding warning will be issued to stdout.
See also
vrna_fold_compound() for how to retrieve the necessary data structure
Parameters
vcthe fold compound data structure
structureWill hold the structure or constraints
Returns
vrna_dimer_pf_t structure containing a set of energies needed for concentration computations.
SWIG Wrapper Notes:
This function is attached as method pf_dimer() to objects of type fold_compound
float vrna_pf_fold ( const char *  sequence,
char *  structure,
vrna_ep_t **  pl 
)

#include <ViennaRNA/part_func.h>

Compute Partition function $Q$ (and base pair probabilities) for an RNA sequence using a comparative method.

This simplified interface to vrna_pf() computes the partition function and, if required, base pair probabilities for an RNA sequence using default options. Memory required for dynamic programming (DP) matrices will be allocated and free'd on-the-fly. Hence, after return of this function, the recursively filled matrices are not available any more for any post-processing.

Note
In case you want to use the filled DP matrices for any subsequent post-processing step, or you require other conditions than specified by the default model details, use vrna_pf(), and the data structure vrna_fold_compound_t instead.
See also
vrna_pf_circfold(), vrna_pf(), vrna_fold_compound(), vrna_fold_compound_t
Parameters
sequenceRNA sequence
structureA pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
plA pointer to a list of vrna_ep_t to store pairing probabilities (Maybe NULL)
Returns
The ensemble free energy $G = -RT \cdot \log(Q) $ in kcal/mol
float vrna_pf_circfold ( const char *  sequence,
char *  structure,
vrna_ep_t **  pl 
)

#include <ViennaRNA/part_func.h>

Compute Partition function $Q$ (and base pair probabilities) for a circular RNA sequences using a comparative method.

This simplified interface to vrna_pf() computes the partition function and, if required, base pair probabilities for a circular RNA sequence using default options. Memory required for dynamic programming (DP) matrices will be allocated and free'd on-the-fly. Hence, after return of this function, the recursively filled matrices are not available any more for any post-processing.

Note
In case you want to use the filled DP matrices for any subsequent post-processing step, or you require other conditions than specified by the default model details, use vrna_pf(), and the data structure vrna_fold_compound_t instead.

Folding of circular RNA sequences is handled as a post-processing step of the forward recursions. See [10] for further details.

See also
vrna_pf_fold(), vrna_pf(), vrna_fold_compound(), vrna_fold_compound_t
Parameters
sequenceA circular RNA sequence
structureA pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
plA pointer to a list of vrna_ep_t to store pairing probabilities (Maybe NULL)
Returns
The ensemble free energy $G = -RT \cdot \log(Q) $ in kcal/mol
float vrna_pf_alifold ( const char **  sequences,
char *  structure,
vrna_ep_t **  pl 
)

#include <ViennaRNA/part_func.h>

Compute Partition function $Q$ (and base pair probabilities) for an RNA sequence alignment using a comparative method.

This simplified interface to vrna_pf() computes the partition function and, if required, base pair probabilities for an RNA sequence alignment using default options. Memory required for dynamic programming (DP) matrices will be allocated and free'd on-the-fly. Hence, after return of this function, the recursively filled matrices are not available any more for any post-processing.

Note
In case you want to use the filled DP matrices for any subsequent post-processing step, or you require other conditions than specified by the default model details, use vrna_pf(), and the data structure vrna_fold_compound_t instead.
See also
vrna_pf_circalifold(), vrna_pf(), vrna_fold_compound_comparative(), vrna_fold_compound_t
Parameters
sequencesRNA sequence alignment
structureA pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
plA pointer to a list of vrna_ep_t to store pairing probabilities (Maybe NULL)
Returns
The ensemble free energy $G = -RT \cdot \log(Q) $ in kcal/mol
float vrna_pf_circalifold ( const char **  sequences,
char *  structure,
vrna_ep_t **  pl 
)

#include <ViennaRNA/part_func.h>

Compute Partition function $Q$ (and base pair probabilities) for an alignment of circular RNA sequences using a comparative method.

This simplified interface to vrna_pf() computes the partition function and, if required, base pair probabilities for an RNA sequence alignment using default options. Memory required for dynamic programming (DP) matrices will be allocated and free'd on-the-fly. Hence, after return of this function, the recursively filled matrices are not available any more for any post-processing.

Note
In case you want to use the filled DP matrices for any subsequent post-processing step, or you require other conditions than specified by the default model details, use vrna_pf(), and the data structure vrna_fold_compound_t instead.

Folding of circular RNA sequences is handled as a post-processing step of the forward recursions. See [10] for further details.

See also
vrna_pf_alifold(), vrna_pf(), vrna_fold_compound_comparative(), vrna_fold_compound_t
Parameters
sequencesSequence alignment of circular RNAs
structureA pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
plA pointer to a list of vrna_ep_t to store pairing probabilities (Maybe NULL)
Returns
The ensemble free energy $G = -RT \cdot \log(Q) $ in kcal/mol
vrna_ep_t* vrna_plist_from_probs ( vrna_fold_compound_t vc,
double  cut_off 
)

#include <ViennaRNA/utils/structures.h>

Create a vrna_ep_t from base pair probability matrix.

The probability matrix provided via the vrna_fold_compound_t is parsed and all pair probabilities above the given threshold are used to create an entry in the plist

The end of the plist is marked by sequence positions i as well as j equal to 0. This condition should be used to stop looping over its entries

Parameters
[in]vcThe fold compound
[in]cut_offThe cutoff value
Returns
A pointer to the plist that is to be created
double * vrna_positional_entropy ( vrna_fold_compound_t fc)

#include <ViennaRNA/equilibrium_probs.h>

Compute a vector of positional entropies.

This function computes the positional entropies from base pair probabilities as

\[ S(i) = - \sum_i p_{ij} \log(p_{ij}) - q_i \log(q_i) \]

with unpaired probabilities $ q_i = 1 - \sum_j p_{ij} $.

Low entropy regions have little structural flexibility and the reliability of the predicted structure is high. High entropy implies many structural alternatives. While these alternatives may be functionally important, they make structure prediction more difficult and thus less reliable.

Precondition
This function requires pre-computed base pair probabilities! Thus, vrna_pf() must be called beforehand.
Parameters
fcA fold_compound with pre-computed base pair probabilities
Returns
A 1-based vector of positional entropies $ S(i) $. (position 0 contains the sequence length)
SWIG Wrapper Notes:
This function is attached as method positional_entropy() to objects of type fold_compound
vrna_dimer_pf_t vrna_pf_co_fold ( const char *  seq,
char *  structure,
vrna_ep_t **  pl 
)

#include <ViennaRNA/part_func.h>

Calculate partition function and base pair probabilities of nucleic acid/nucleic acid dimers.

This simplified interface to vrna_pf_dimer() computes the partition function and, if required, base pair probabilities for an RNA-RNA interaction using default options. Memory required for dynamic programming (DP) matrices will be allocated and free'd on-the-fly. Hence, after return of this function, the recursively filled matrices are not available any more for any post-processing.

Note
In case you want to use the filled DP matrices for any subsequent post-processing step, or you require other conditions than specified by the default model details, use vrna_pf_dimer(), and the data structure vrna_fold_compound_t instead.
See also
vrna_pf_dimer()
Parameters
seqTwo concatenated RNA sequences with a delimiting '&' in between
structureA pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
plA pointer to a list of vrna_ep_t to store pairing probabilities (Maybe NULL)
Returns
vrna_dimer_pf_t structure containing a set of energies needed for concentration computations.