Global Partition Function and Equilibrium Probabilities
Variations of the global partition function algorithm.
We provide implementations of the global partition function algorithm for
Single sequences,
Multiple sequence alignments (MSA), and
RNA-RNA hybrids
Basic global partition function interface
-
FLT_OR_DBL vrna_pf(vrna_fold_compound_t *fc, 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.
- SWIG Wrapper Notes:
This function is attached as method
pf()
to objects of typefold_compound
. See, e.g.RNA.fold_compound.pf()
in the Python API.
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()
Note
This function is polymorphic. It accepts vrna_fold_compound_t of type VRNA_FC_TYPE_SINGLE, and VRNA_FC_TYPE_COMPARATIVE. Also, 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
.- Parameters
fc – [inout] The fold compound data structure
structure – [inout] A 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
-
vrna_dimer_pf_t vrna_pf_dimer(vrna_fold_compound_t *fc, 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.
- SWIG Wrapper Notes:
This function is attached as method
pf_dimer()
to objects of typefold_compound
. See, e.g.RNA.fold_compound.pf_dimer()
in the Python API.
See also
vrna_fold_compound() for how to retrieve the necessary data structure
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 tostdout
.- Parameters
fc – the fold compound data structure
structure – Will hold the structure or constraints
- Returns
vrna_dimer_pf_t structure containing a set of energies needed for concentration computations.
-
FLT_OR_DBL *vrna_pf_substrands(vrna_fold_compound_t *fc, size_t complex_size)
- #include <ViennaRNA/part_func.h>
-
FLT_OR_DBL vrna_pf_add(FLT_OR_DBL dG1, FLT_OR_DBL dG2, double kT)
- #include <ViennaRNA/part_func.h>
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)
- #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.
- Parameters
sequence – RNA sequence
structure – A pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
pl – A 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.
Folding of circular RNA sequences is handled as a post-processing step of the forward recursions. See Hofacker and Stadler [2006] for further details.
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.
- Parameters
sequence – A circular RNA sequence
structure – A pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
pl – A 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.
- Parameters
sequences – RNA sequence alignment
structure – A pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
pl – A 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.
Folding of circular RNA sequences is handled as a post-processing step of the forward recursions. See Hofacker and Stadler [2006] for further details.
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.
- Parameters
sequences – Sequence alignment of circular RNAs
structure – A pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
pl – A 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_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.
See also
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.
- Parameters
seq – Two concatenated RNA sequences with a delimiting ‘&’ in between
structure – A pointer to the character array where position-wise pairing propensity will be stored. (Maybe NULL)
pl – A 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.
Functions
-
vrna_ep_t *vrna_plist_from_probs(vrna_fold_compound_t *fc, 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
fc – [in] The fold compound
cut_off – [in] The cutoff value
- Returns
A pointer to the plist that is to be created
-
struct vrna_dimer_pf_s
- #include <ViennaRNA/part_func.h>
Data structure returned by vrna_pf_dimer()
-
struct vrna_multimer_pf_s
-
FLT_OR_DBL vrna_pf(vrna_fold_compound_t *fc, char *structure)