Evolutionary Analysis Using Likelihood#

Specifying substitution models#

The available pre-defined substitution models#

In cases where code takes a substitution model as an argument, you can use the value under “Abbreviation” as a string.

from cogent3 import available_models

available_models()
Specify a model using 'Abbreviation' (case sensitive).
Model TypeAbbreviationDescription
nucleotideBHBarry and Hartigan Discrete Time substitution model Barry and Hartigan 1987. Biometrics 43: 261–76.
nucleotideDTDiscrete Time substitution model (non-stationary, non-reversible). motif_length=2 makes this a dinucleotide model, motif_length=3 a trinucleotide model.
nucleotideGNGeneral Markov Nucleotide (non-stationary, non-reversible). Kaehler, Yap, Zhang, Huttley, 2015, Sys Biol 64 (2): 281–93
nucleotidessGNstrand-symmetric general Markov nucleotide (non-stationary, non-reversible). Kaehler, 2017, Journal of Theoretical Biology 420: 144–51
nucleotideK80Kimura 1980
nucleotideJC69Jukes and Cantor's 1969 model
nucleotideGTRGeneral Time Reversible nucleotide substitution model.
nucleotideTN93Tamura and Nei 1993 model
nucleotideHKY85Hasegawa, Kishino and Yano 1985 model
nucleotideF81Felsenstein's 1981 model
codonCNFGTRConditional nucleotide frequency codon substitution model, GTR variant (with params analagous to the nucleotide GTR model). Yap, Lindsay, Easteal and Huttley, 2010, Mol Biol Evol 27: 726-734
codonCNFHKYConditional nucleotide frequency codon substitution model, HKY variant (with kappa, the ratio of transitions to transversions) Yap, Lindsay, Easteal and Huttley, 2010, Mol Biol Evol 27: 726-734
codonMG94HKYMuse and Gaut 1994 codon substitution model, HKY variant (with kappa, the ratio of transitions to transversions) Muse and Gaut, 1994, Mol Biol Evol, 11, 715-24
codonMG94GTRMuse and Gaut 1994 codon substitution model, GTR variant (with params analagous to the nucleotide GTR model) Muse and Gaut, 1994, Mol Biol Evol, 11, 715-24
codonGY94Goldman and Yang 1994 codon substitution model. N Goldman and Z Yang, 1994, Mol Biol Evol, 11(5):725-36.
codonY98Yang's 1998 substitution model, a derivative of the GY94. Z Yang, 1998, Mol Biol Evol, 15(5):568-73
codonH04GHuttley 2004 CpG substitution model. Includes a term for substitutions to or from CpG's. GA Huttley, 2004, Mol Biol Evol, 21(9):1760-8
codonH04GKHuttley 2004 CpG substitution model. Includes a term for transition substitutions to or from CpG's. GA Huttley, 2004, Mol Biol Evol, 21(9):1760-8
codonH04GGKHuttley 2004 CpG substitution model. Includes a general term for substitutions to or from CpG's and an adjustment for CpG transitions. GA Huttley, 2004, Mol Biol Evol, 21(9):1760-8
codonGNCGeneral Nucleotide Codon, a non-reversible codon model. Kaehler, Yap, Huttley, 2017, Gen Biol Evol 9(1): 134–49
proteinDSO78Dayhoff et al 1978 empirical protein model Dayhoff, MO, Schwartz RM, and Orcutt, BC. 1978 A model of evolutionary change in proteins. Pp. 345-352. Atlas of protein sequence and structure, Vol 5, Suppl. 3. National Biomedical Research Foundation, Washington D. C Matrix imported from PAML dayhoff.dat file
proteinJTT92Jones, Taylor and Thornton 1992 empirical protein model Jones DT, Taylor WR, Thornton JM. The rapid generation of mutation data matrices from protein sequences. Comput Appl Biosci. 1992 Jun;8(3):275-82. Matrix imported from PAML jones.dat file
proteinAH96Adachi and Hasegawa 1996 empirical model for mitochondrial proteins. Adachi J, Hasegawa M. Model of amino acid substitution in proteins encoded by mitochondrial DNA. J Mol Evol. 1996 Apr;42(4):459-68. Matrix imported from PAML mtREV24.dat file
proteinAH96_mtmammalsAdachi and Hasegawa 1996 empirical model for mammalian mitochondrial proteins. Adachi J, Hasegawa M. Model of amino acid substitution in proteins encoded by mitochondrial DNA. J Mol Evol. 1996 Apr;42(4):459-68. Matrix imported from PAML mtmam.dat file
proteinWG01Whelan and Goldman 2001 empirical model for globular proteins. Whelan S, Goldman N. A general empirical model of protein evolution derived from multiple protein families using a maximum-likelihood approach. Mol Biol Evol. 2001 May;18(5):691-9. Matrix imported from PAML wag.dat file

25 rows x 3 columns

Getting a substitution model with get_model()#

from cogent3.evolve.models import get_model

hky = get_model("HKY85")
print(hky)
TimeReversibleNucleotide(name='HKY85'; params=['kappa']; num_motifs=4; motifs=['T', 'C', 'A', 'G']))

Rate heterogeneity models#

We illustrate this for the gamma distributed case using examples of the canned models displayed above. Creating rate heterogeneity variants of the canned models can be done by using optional arguments that get passed to the substitution model class.

For nucleotide#

We specify a general time reversible nucleotide model with gamma distributed rate heterogeneity.

from cogent3.evolve.models import get_model

sub_mod = get_model("GTR", with_rate=True, distribution="gamma")
print(sub_mod)
TimeReversibleNucleotide(name='GTR'; params=['A/C', 'A/G', 'A/T', 'C/G', 'C/T']; num_motifs=4; motifs=['T', 'C', 'A', 'G']))

For codon#

We specify a conditional nucleotide frequency codon model with nucleotide general time reversible parameters and a parameter for the ratio of nonsynonymous to synonymous substitutions (omega) with gamma distributed rate heterogeneity.

from cogent3.evolve.models import get_model

sub_mod = get_model("CNFGTR", with_rate=True, distribution="gamma")
print(sub_mod)
TimeReversibleCodon(name='CNFGTR'; params=['A/C', 'A/G', 'A/T', 'C/G', 'C/T', 'omega']; num_motifs=61; motifs=['TTT', 'TTC', 'TTA', 'TTG', 'TCT', 'TCC', 'TCA', 'TCG', 'TAT', 'TAC', 'TGT', 'TGC', 'TGG', 'CTT', 'CTC', 'CTA', 'CTG', 'CCT', 'CCC', 'CCA', 'CCG', 'CAT', 'CAC', 'CAA', 'CAG', 'CGT', 'CGC', 'CGA', 'CGG', 'ATT', 'ATC', 'ATA', 'ATG', 'ACT', 'ACC', 'ACA', 'ACG', 'AAT', 'AAC', 'AAA', 'AAG', 'AGT', 'AGC', 'AGA', 'AGG', 'GTT', 'GTC', 'GTA', 'GTG', 'GCT', 'GCC', 'GCA', 'GCG', 'GAT', 'GAC', 'GAA', 'GAG', 'GGT', 'GGC', 'GGA', 'GGG']))

For protein#

We specify a Jones, Taylor and Thornton 1992 empirical protein substitution model with gamma distributed rate heterogeneity.

from cogent3.evolve.models import get_model

sub_mod = get_model("JTT92", with_rate=True, distribution="gamma")
print(sub_mod)
Empirical(name='JTT92'; num_motifs=20; motifs=['A', 'C', 'D', 'E', 'F', 'G', 'H', 'I', 'K', 'L', 'M', 'N', 'P', 'Q', 'R', 'S', 'T', 'V', 'W', 'Y']))

Making a likelihood function#

You start by specifying a substitution model and use that to construct a likelihood function for a specific tree.

from cogent3 import make_tree
from cogent3.evolve.models import get_model

sub_mod = get_model("F81")
tree = make_tree("(a,b,(c,d))")
lf = sub_mod.make_likelihood_function(tree)

Providing an alignment to a likelihood function#

You need to load an alignment and then provide it a likelihood function. I construct very simple trees and alignments for this example.

from cogent3 import make_aligned_seqs, make_tree
from cogent3.evolve.models import get_model

sub_mod = get_model("F81")
tree = make_tree("(a,b,(c,d))")
lf = sub_mod.make_likelihood_function(tree)
aln = make_aligned_seqs(
    [("a", "ACGT"), ("b", "AC-T"), ("c", "ACGT"), ("d", "AC-T")]
)
lf.set_alignment(aln)

Scoping parameters on trees – time heterogeneous models#

For many evolutionary analyses, it’s desirable to allow different branches on a tree to have different values of a parameter. We show this for a simple codon model case here where we want the great apes (the clade that includes human and orangutan) to have a different value of the ratio of nonsynonymous to synonymous substitutions. This parameter is identified in the precanned CNFGTR model as omega.

from cogent3 import load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
print(tree.ascii_art())
          /-Galago
         |
-root----|--HowlerMon
         |
         |          /-Rhesus
          \edge.3--|
                   |          /-Orangutan
                    \edge.2--|
                             |          /-Gorilla
                              \edge.1--|
                                       |          /-Human
                                        \edge.0--|
                                                  \-Chimpanzee
sm = get_model("CNFGTR")
lf = sm.make_likelihood_function(tree, digits=2)
lf.set_param_rule(
    "omega",
    tip_names=["Human", "Orangutan"],
    outgroup_name="Galago",
    clade=True,
    init=0.5,
)

We’ve set an initial value for this clade so that the edges affected by this rule are evident below.

lf

CNFGTR

number of free parameters = 78

Global params
A/CA/GA/TC/GC/T
1.001.001.001.001.00
Edge params
edgeparentlengthomega
Galagoroot1.001.00
HowlerMonroot1.001.00
Rhesusedge.31.001.00
Orangutanedge.21.000.50
Gorillaedge.11.000.50
Humanedge.01.000.50
Chimpanzeeedge.01.000.50
edge.0edge.11.000.50
edge.1edge.21.000.50
edge.2edge.31.001.00
edge.3root1.001.00
Motif params
AAAAACAAGAATACAACCACGACTAGAAGCAGGAGTATA
0.020.020.020.020.020.020.020.020.020.020.020.020.02
continuation
ATCATGATTCAACACCAGCATCCACCCCCGCCTCGACGC
0.020.020.020.020.020.020.020.020.020.020.020.020.02
continuation
CGGCGTCTACTCCTGCTTGAAGACGAGGATGCAGCCGCG
0.020.020.020.020.020.020.020.020.020.020.020.020.02
continuation
GCTGGAGGCGGGGGTGTAGTCGTGGTTTACTATTCATCC
0.020.020.020.020.020.020.020.020.020.020.020.020.02
continuation
TCGTCTTGCTGGTGTTTATTCTTGTTT
0.020.020.020.020.020.020.020.020.02

A more extensive description of capabilities is in Allowing substitution model parameters to differ between branches.

Specifying a parameter as constant#

This means the parameter will not be modified during likelihood maximisation. We show this here by making the omega parameter constant at the value 1 – essentially the condition of selective neutrality.

from cogent3 import load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
sm = get_model("CNFGTR")
lf = sm.make_likelihood_function(tree, digits=2)
lf.set_param_rule("omega", is_constant=True)

Providing a starting value for a parameter#

This can be useful to improve performance, the closer you are to the maximum likelihood estimator the quicker optimisation will be.

from cogent3 import load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
sm = get_model("CNFGTR")
lf = sm.make_likelihood_function(tree, digits=2)
lf.set_param_rule("omega", init=0.1)

Setting parameter bounds for optimisation#

This can be useful for stopping optimisers from getting stuck in a bad part of parameter space. The following is for omega in a codon model. I’m also providing an initial guess for the parameter (init=0.1) as well as a lower bound. An initial guess that is close to the maximum likelihood estimate will speed up optimisation.

from cogent3 import load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
sm = get_model("CNFGTR")
lf = sm.make_likelihood_function(tree, digits=2)
lf.set_param_rule("omega", init=0.1, lower=1e-9, upper=20.0)

Setting an upper bound for branch length#

If the branch length estimates seem too large, setting just an upper bound can be sensible. This will apply to all edges on the tree.

from cogent3 import load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
sm = get_model("F81")
lf = sm.make_likelihood_function(tree)
lf.set_param_rule("length", upper=1.0)

Note

If, after optimising, the branch lengths equal to the upper value you set then the function has not been fully maximised and you should consider adjusting the boundary again.

Specifying rate heterogeneity functions#

We extend the simple gamma distributed rate heterogeneity case for nucleotides from above to construction of the actual likelihood function. We do this for 4 bins and constraint the bin probabilities to be equal.

from cogent3 import load_tree
from cogent3.evolve.models import get_model

sm = get_model("GTR", with_rate=True, distribution="gamma")
tree = load_tree("data/primate_brca1.tree")
lf = sm.make_likelihood_function(tree, bins=4, digits=2)
lf.set_param_rule("bprobs", is_constant=True)

For more detailed discussion of defining and using these models see Analysis of rate heterogeneity.

Specifying Phylo-HMMs#

from cogent3 import load_tree
from cogent3.evolve.models import get_model

sm = get_model("GTR", with_rate=True, distribution="gamma")
tree = load_tree("data/primate_brca1.tree")
lf = sm.make_likelihood_function(tree, bins=4, sites_independent=False, digits=2)
lf.set_param_rule("bprobs", is_constant=True)

For more detailed discussion of defining and using these models see Evaluate process heterogeneity using a Hidden Markov Model.

Fitting likelihood functions - Choice of optimisers#

There are 2 types of optimiser: simulated annealing, a global optimiser; and Powell, a local optimiser. The simulated annealing method is slow compared to Powell and in general Powell is an adequate choice. I setup a simple nucleotide model to illustrate these.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
aln = load_aligned_seqs("data/primate_brca1.fasta")
sm = get_model("F81")
lf = sm.make_likelihood_function(tree, digits=3, space=2)
lf.set_alignment(aln)
lf.optimise(show_progress=False)

The default is to use Powell. For Powell, it’s recommended to set the max_restarts argument since this provides a mechanism for Powell to attempt restarting the optimisation from a slightly different spot which can help in overcoming local maxima.

lf.optimise(local=True, max_restarts=5, show_progress=False)

We might want to do crude simulated annealing following by more rigorous Powell. To do this we first need to use the global optimiser, setting local=False setting a large value for global_tolerance.

lf.optimise(local=False, global_tolerance=1.0, show_progress=False)

Followed by a standard call to optimise().

lf.optimise(show_progress=False, max_restarts=5, tolerance=1e-8)

How to check your optimisation was successful#

There is no guarantee that an optimised function has achieved a global maximum. We can, however, be sure that a maximum was achieved by validating that the optimiser stopped because the specified tolerance condition was met, rather than exceeding the maximum number of evaluations. The latter number is set to ensure optimisation doesn’t proceed endlessly. If the optimiser exited because this limit was exceeded you can be sure that the function has not been successfully optimised.

We can monitor this situation using the limit_action argument to optimise. Providing the value raise causes an exception to be raised if this condition occurs, as shown below. Providing warn (default) instead will cause a warning message to be printed to screen but execution will continue. The value ignore hides any such message.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
aln = load_aligned_seqs("data/primate_brca1.fasta")
sm = get_model("F81")
lf = sm.make_likelihood_function(tree, digits=3, space=2)
lf.set_alignment(aln)
try:
    lf.optimise(
        show_progress=False,
        limit_action="raise",
        max_evaluations=10,
        return_calculator=True,
    )
except ArithmeticError as err:
    print(err)
FORCED EXIT from optimiser after 10 evaluations

Note

We recommend using limit_action='raise' and catching the ArithmeticError error explicitly (as demonstrated above). You really shouldn’t be using results from such an optimisation run.

Overview of the fitted likelihood function#

In Jupyter, the likelihood function object presents a representation of the main object features.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

sm = get_model("GTR")
tree = load_tree("data/primate_brca1.tree")
lf = sm.make_likelihood_function(tree)
aln = load_aligned_seqs("data/primate_brca1.fasta")
lf.set_alignment(aln)
lf.optimise(local=True, show_progress=False)
lf

GTR

log-likelihood = -6992.7690

number of free parameters = 16

Global params
A/CA/GA/TC/GC/T
1.23165.25340.95852.31595.9700
Edge params
edgeparentlength
Galagoroot0.1731
HowlerMonroot0.0449
Rhesusedge.30.0216
Orangutanedge.20.0077
Gorillaedge.10.0025
Humanedge.00.0061
Chimpanzeeedge.00.0028
edge.0edge.10.0000
edge.1edge.20.0034
edge.2edge.30.0120
edge.3root0.0076
Motif params
ACGT
0.37570.17420.20950.2406

Log likelihood and number of free parameters#

Reusing the optimised lf object from above, we can get the log-likelihood and the number of free parameters.

lnL = lf.lnL
lnL
-6992.768994254335
nfp = lf.nfp
nfp
16

Warning

The number of free parameters (nfp) refers only to the number of parameters that were modifiable by the optimiser. Typically, the degrees-of-freedom of a likelihood ratio test statistic is computed as the difference in nfp between models. This will not be correct for models in which a boundary conditions exist (rate heterogeneity models where a parameter value boundary is set between bins).

Aikake Information Criterion#

Reusing the optimised lf object from above.

lf.get_aic()
14017.53798850867

We can also get the second-order AIC.

lf.get_aic(second_order=True)
14017.732482609492

Bayesian Information Criterion#

Reusing the optimised lf object from above.

lf.get_bic()
np.float64(14112.61578431146)

Getting maximum likelihood estimates#

Reusing the optimised lf object from above.

One at a time#

We get the statistics out individually. We get the length for the Human edge and the exchangeability parameter A/G.

a_g = lf.get_param_value("A/G")
a_g
np.float64(5.253426081232493)
human = lf.get_param_value("length", "Human")
human
np.float64(0.006060124932976551)
Just the motif probabilities#
mprobs = lf.get_motif_probs()
mprobs
TCAG
0.24060.17420.37570.2095
As tables#
tables = lf.get_statistics(with_motif_probs=True, with_titles=True)
tables[0]  # just displaying the first
global params
A/CA/GA/TC/GC/T
1.23165.25340.95852.31595.9700

1 rows x 5 columns

Testing Hypotheses - Using Likelihood Ratio Tests#

We test the molecular clock hypothesis for human and chimpanzee lineages. The null has these two branches constrained to be equal.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
aln = load_aligned_seqs("data/primate_brca1.fasta")
sm = get_model("F81")
lf = sm.make_likelihood_function(tree, digits=3, space=2)
lf.set_alignment(aln)
lf.set_param_rule(
    "length",
    tip_names=["Human", "Chimpanzee"],
    outgroup_name="Galago",
    clade=True,
    is_independent=False,
)
lf.set_name("Null Hypothesis")
lf.optimise(local=True, show_progress=False)
null_lnL = lf.lnL
null_nfp = lf.nfp
lf

Null Hypothesis

log-likelihood = -7177.4403

number of free parameters = 10

Edge params
edgeparentlength
Galagoroot0.167
HowlerMonroot0.044
Rhesusedge.30.021
Orangutanedge.20.008
Gorillaedge.10.002
Humanedge.00.004
Chimpanzeeedge.00.004
edge.0edge.10.000
edge.1edge.20.003
edge.2edge.30.012
edge.3root0.009
Motif params
ACGT
0.3760.1740.2090.241

The alternate allows the human and chimpanzee branches to differ by just setting all lengths to be independent.

lf.set_param_rule("length", is_independent=True)
lf.set_name("Alt Hypothesis")
lf.optimise(local=True, show_progress=False)
alt_lnL = lf.lnL
alt_nfp = lf.nfp
lf

Alt Hypothesis

log-likelihood = -7175.7756

number of free parameters = 11

Edge params
edgeparentlength
Galagoroot0.167
HowlerMonroot0.044
Rhesusedge.30.021
Orangutanedge.20.008
Gorillaedge.10.002
Humanedge.00.006
Chimpanzeeedge.00.003
edge.0edge.10.000
edge.1edge.20.003
edge.2edge.30.012
edge.3root0.009
Motif params
ACGT
0.3760.1740.2090.241

We import the function for computing the probability of a chi-square test statistic, compute the likelihood ratio test statistic, degrees of freedom and the corresponding probability.

from scipy.stats.distributions import chi2

LR = 2 * (alt_lnL - null_lnL)  # the likelihood ratio statistic
df = alt_nfp - null_nfp  # the test degrees of freedom
p = chi2.sf(LR, df)
print(f"LR={LR:.4f} ; df={df}; p={df:.4f}")
LR=3.3294 ; df=1; p=1.0000

Testing Hypotheses - By parametric bootstrapping#

If we can’t rely on the asymptotic behaviour of the LRT, e.g. due to small alignment length, we can use a parametric bootstrap. Convenience functions for that are described in more detail here Performing a parametric bootstrap.

In general, however, this capability derives from the ability of any defined evolve likelihood function to simulate an alignment. This property is provided as simulate_alignment method on likelihood function objects.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
aln = load_aligned_seqs("data/primate_brca1.fasta")

sm = get_model("F81")
lf = sm.make_likelihood_function(tree, digits=3, space=2)
lf.set_alignment(aln)
lf.set_param_rule(
    "length",
    tip_names=["Human", "Chimpanzee"],
    outgroup_name="Galago",
    clade=True,
    is_independent=False,
)
lf.set_name("Null Hypothesis")
lf.optimise(local=True, show_progress=False)
sim_aln = lf.simulate_alignment()
sim_aln[:60]
0
ChimpanzeeAATAGAACACTAGTATTGCAAACGAACCTACTTAAATTGAGTTAAATGGAACATTAAAAA
Galago....A.......A....A.....A............GA.......T........A.....
HowlerMon............A.......................G........T..............
Rhesus......C.....A.......................G........T..............
Orangutan............................................................
Gorilla............................................................
Human............................................................

7 x 60 dna alignment

Determining confidence intervals on MLEs#

The profile method is used to calculate a confidence interval for a named parameter. We show it here for a global substitution model exchangeability parameter (kappa, the ratio of transition to transversion rates) and for an edge specific parameter (just the human branch length).

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
aln = load_aligned_seqs("data/primate_brca1.fasta")
sm = get_model("HKY85")
lf = sm.make_likelihood_function(tree)
lf.set_alignment(aln)
lf.optimise(local=True, show_progress=False)
kappa_lo, kappa_mle, kappa_hi = lf.get_param_interval("kappa")
print(f"lo={kappa_lo:.2f} ; mle={kappa_mle:.2f} ; hi={kappa_hi:.2f}")
human_lo, human_mle, human_hi = lf.get_param_interval("length", "Human")
print(f"lo={human_lo:.2f} ; mle={human_mle:.2f} ; hi={human_hi:.2f}")
lo=3.78 ; mle=4.44 ; hi=5.22
lo=0.00 ; mle=0.01 ; hi=0.01

Saving results#

The best approach is to use the json string from the to_json() method. The saved data can be later reloaded using cogent3.util.deserialise.deserialise_object(). The json data contains the alignment, tree topology, substitution model, parameter values, etc..

To illustrate this, I create a very simple likelihood function. The json variable below is just a string that can be saved to disk.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

aln = make_aligned_seqs(data=dict(a="ACGG", b="ATAG", c="ATGG"))
tree = make_tree(tip_names=aln.names)
sm = get_model("F81")
lf = sm.make_likelihood_function(tree)
lf.set_alignment(aln)
json = lf.to_json()
json[:60]  # just truncating the displayed string
'{"model": {"alphabet": {"motifset": ["T", "C", "A", "G"], "g'

We deserialise the object from the string.

from cogent3.util.deserialise import deserialise_object

newlf = deserialise_object(json)
newlf

F81

log-likelihood = -14.2727

number of free parameters = 3

Edge params
edgeparentlength
aroot1.0000
broot1.0000
croot1.0000
Motif params
ACGT
0.33330.08330.41670.1667

Reconstructing ancestral sequences#

We first fit a likelihood function.

from cogent3 import load_aligned_seqs, load_tree
from cogent3.evolve.models import get_model

tree = load_tree("data/primate_brca1.tree")
aln = load_aligned_seqs("data/primate_brca1.fasta")
sm = get_model("F81")
lf = sm.make_likelihood_function(tree, digits=3, space=2)
lf.set_alignment(aln)
lf.optimise(show_progress=False)

We then get the most likely ancestral sequences.

ancestors = lf.likely_ancestral_seqs()
ancestors[:60]
0
edge.0TGTGGCACAAATACTCATGCCAGCTCATTACAGCATGAGAACAGTTTATTACTCACTAAA
edge.1............................................................
edge.2............................................................
edge.3...............................................G............
root...............................................G............

5 x 60 dna alignment

Or we can get the posterior probabilities (returned as a DictArray) of sequence states at each node.

ancestral_probs = lf.reconstruct_ancestral_seqs()
ancestral_probs["root"][:5]
TCAG
00.18160.00000.00000.0000
10.00000.00000.00000.1561
20.18160.00000.00000.0000
30.00000.00000.00000.1561
40.00000.00000.00000.1561

There’s nothing that improves performance quite like being close to the maximum likelihood values. So using the set_param_rule method to provide good starting values can be very useful. As this can be difficult to do one easy way is to build simpler models that are nested within the one you’re interested in. Fitting those models and then relaxing constraints until you’re at the parameterisation of interest can markedly improve optimisation speed.