{ "cells": [ { "cell_type": "markdown", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "# Evolutionary Analysis Using Likelihood" ] }, { "cell_type": "markdown", "metadata": { "toc-hr-collapsed": false }, "source": [ "## Specifying substitution models" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The available pre-defined substitution models\n", "\n", "In cases where code takes a substitution model as an argument, you can use the value under \"Abbreviation\" as a string." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Specify a model using 'Abbreviation' (case sensitive).
Model TypeAbbreviationDescription
nucleotideJC69Jukes and Cantor's 1969 model
nucleotideK80Kimura 1980
nucleotideF81Felsenstein's 1981 model
nucleotideHKY85Hasegawa, Kishino and Yanamo 1985 model
nucleotideTN93Tamura and Nei 1993 model
nucleotideGTRGeneral Time Reversible nucleotide substitution model.
nucleotidessGNstrand-symmetric general Markov nucleotide (non-stationary, non-reversible). Kaehler, 2017, Journal of Theoretical Biology 420: 144–51
nucleotideGNGeneral Markov Nucleotide (non-stationary, non-reversible). Kaehler, Yap, Zhang, Huttley, 2015, Sys Biol 64 (2): 281–93
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.
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
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
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
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
\n", "

\n", "25 rows x 3 columns

" ], "text/plain": [ "Specify a model using 'Abbreviation' (case sensitive).\n", "================================================================================================================================================================================================================================================================================================================================================\n", "Model Type Abbreviation Description\n", "------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "nucleotide JC69 Jukes and Cantor's 1969 model\n", "nucleotide K80 Kimura 1980\n", "nucleotide F81 Felsenstein's 1981 model\n", "nucleotide HKY85 Hasegawa, Kishino and Yanamo 1985 model\n", "nucleotide TN93 Tamura and Nei 1993 model\n", "nucleotide GTR General Time Reversible nucleotide substitution model.\n", "nucleotide ssGN strand-symmetric general Markov nucleotide (non-stationary, non-reversible). Kaehler, 2017, Journal of Theoretical Biology 420: 144–51\n", "nucleotide GN General Markov Nucleotide (non-stationary, non-reversible). Kaehler, Yap, Zhang, Huttley, 2015, Sys Biol 64 (2): 281–93\n", "nucleotide BH Barry and Hartigan Discrete Time substitution model Barry and Hartigan 1987. Biometrics 43: 261–76.\n", "nucleotide DT Discrete Time substitution model (non-stationary, non-reversible). motif_length=2 makes this a dinucleotide model, motif_length=3 a trinucleotide model.\n", " codon CNFGTR Conditional 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\n", " codon CNFHKY Conditional 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\n", " codon MG94HKY Muse 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\n", " codon MG94GTR Muse 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\n", " codon GY94 Goldman and Yang 1994 codon substitution model. N Goldman and Z Yang, 1994, Mol Biol Evol, 11(5):725-36.\n", " codon Y98 Yang's 1998 substitution model, a derivative of the GY94. Z Yang, 1998, Mol Biol Evol, 15(5):568-73\n", " codon H04G Huttley 2004 CpG substitution model. Includes a term for substitutions to or from CpG's. GA Huttley, 2004, Mol Biol Evol, 21(9):1760-8\n", " codon H04GK Huttley 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\n", " codon H04GGK Huttley 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\n", " codon GNC General Nucleotide Codon, a non-reversible codon model. Kaehler, Yap, Huttley, 2017, Gen Biol Evol 9(1): 134–49\n", " protein DSO78 Dayhoff 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\n", " protein AH96 Adachi 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\n", " protein AH96_mtmammals Adachi 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\n", " protein JTT92 Jones, 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\n", " protein WG01 Whelan 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\n", "------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "\n", "25 rows x 3 columns" ] }, "execution_count": 1, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cogent3 import available_models\n", "\n", "available_models()" ] }, { "cell_type": "markdown", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "### Getting a substitution model with ``get_model()``" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "TimeReversibleNucleotide ( name = 'HKY85'; type = 'None'; params = ['kappa']; number of motifs = 4; motifs = ['T', 'C', 'A', 'G'])\n", "\n" ] } ], "source": [ "from cogent3.evolve.models import get_model\n", "\n", "hky = get_model(\"HKY85\")\n", "print(hky)" ] }, { "cell_type": "markdown", "metadata": { "toc-hr-collapsed": true }, "source": [ "### Rate heterogeneity models\n", "\n", "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." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### For nucleotide\n", "\n", "We specify a general time reversible nucleotide model with gamma distributed rate heterogeneity." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "TimeReversibleNucleotide ( name = 'GTR'; type = 'None'; params = ['A/C', 'A/G', 'A/T', 'C/G', 'C/T']; number of motifs = 4; motifs = ['T', 'C', 'A', 'G'])\n", "\n" ] } ], "source": [ "from cogent3.evolve.models import get_model\n", "\n", "sub_mod = get_model(\"GTR\", with_rate=True, distribution='gamma')\n", "\n", "print(sub_mod)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### For codon\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "TimeReversibleCodon ( name = 'CNFGTR'; type = 'None'; params = ['A/C', 'A/G', 'A/T', 'C/G', 'C/T', 'omega']; number of 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'])\n", "\n" ] } ], "source": [ "from cogent3.evolve.models import get_model\n", "\n", "sub_mod = get_model(\"CNFGTR\", with_rate=True, distribution='gamma')\n", "\n", "print(sub_mod)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### For protein\n", "\n", "We specify a Jones, Taylor and Thornton 1992 empirical protein substitution model with gamma distributed rate heterogeneity." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "Empirical ( name = 'JTT92'; type = 'None'; number of motifs = 20; motifs = ['A', 'C', 'D', 'E', 'F', 'G', 'H', 'I', 'K', 'L', 'M', 'N', 'P', 'Q', 'R', 'S', 'T', 'V', 'W', 'Y'])\n", "\n" ] } ], "source": [ "from cogent3.evolve.models import get_model\n", "\n", "sub_mod = get_model(\"JTT92\", with_rate=True, distribution='gamma')\n", "\n", "print(sub_mod)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Making a likelihood function\n", "\n", "You start by specifying a substitution model and use that to construct a likelihood function for a specific tree." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "from cogent3 import make_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "sub_mod = get_model(\"F81\")\n", "tree = make_tree('(a,b,(c,d))')\n", "\n", "lf = sub_mod.make_likelihood_function(tree)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Providing an alignment to a likelihood function\n", "\n", "You need to load an alignment and then provide it a likelihood function. I construct very simple trees and alignments for this example." ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "from cogent3 import make_tree, make_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "sub_mod = get_model(\"F81\")\n", "tree = make_tree('(a,b,(c,d))')\n", "lf = sub_mod.make_likelihood_function(tree)\n", "aln = make_aligned_seqs([('a', 'ACGT'), ('b', 'AC-T'), ('c', 'ACGT'),\n", " ('d', 'AC-T')])\n", "lf.set_alignment(aln)" ] }, { "cell_type": "markdown", "metadata": { "toc-hr-collapsed": true }, "source": [ "### Scoping parameters on trees -- time heterogeneous models\n", "\n", "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``." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ " /-Galago\n", " |\n", "-root----|--HowlerMon\n", " |\n", " | /-Rhesus\n", " \\edge.3--|\n", " | /-Orangutan\n", " \\edge.2--|\n", " | /-Gorilla\n", " \\edge.1--|\n", " | /-Human\n", " \\edge.0--|\n", " \\-Chimpanzee\n" ] } ], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "\n", "print(tree.ascii_art())" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [], "source": [ "sm = get_model(\"CNFGTR\")\n", "lf = sm.make_likelihood_function(tree, digits=2)\n", "lf.set_param_rule('omega', tip_names=['Human', 'Orangutan'], outgroup_name='Galago', clade=True, init=0.5)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We've set an *initial* value for this clade so that the edges affected by this rule are evident below." ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

Likelihood function statistics

\n", "\n", "

number of free parameters = 78

\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Global params
A/CA/GA/TC/GC/T
1.001.001.001.001.00
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
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
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Motif params
AAAAACAAGAATACAACCACGACTAGAAGCAGGAGTATA
0.020.020.020.020.020.020.020.020.020.020.020.020.02
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
ATCATGATTCAACACCAGCATCCACCCCCGCCTCGACGC
0.020.020.020.020.020.020.020.020.020.020.020.020.02
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
CGGCGTCTACTCCTGCTTGAAGACGAGGATGCAGCCGCG
0.020.020.020.020.020.020.020.020.020.020.020.020.02
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
GCTGGAGGCGGGGGTGTAGTCGTGGTTTACTATTCATCC
0.020.020.020.020.020.020.020.020.020.020.020.020.02
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
TCGTCTTGCTGGTGTTTATTCTTGTTT
0.020.020.020.020.020.020.020.020.02
\n" ], "text/plain": [ "Likelihood function statistics\n", "number of free parameters = 78\n", "====================================\n", " A/C A/G A/T C/G C/T\n", "------------------------------------\n", "1.00 1.00 1.00 1.00 1.00\n", "------------------------------------\n", "=======================================\n", " edge parent length omega\n", "---------------------------------------\n", " Galago root 1.00 1.00\n", " HowlerMon root 1.00 1.00\n", " Rhesus edge.3 1.00 1.00\n", " Orangutan edge.2 1.00 0.50\n", " Gorilla edge.1 1.00 0.50\n", " Human edge.0 1.00 0.50\n", "Chimpanzee edge.0 1.00 0.50\n", " edge.0 edge.1 1.00 0.50\n", " edge.1 edge.2 1.00 0.50\n", " edge.2 edge.3 1.00 1.00\n", " edge.3 root 1.00 1.00\n", "---------------------------------------\n", "============================================================================\n", " AAA AAC AAG AAT ACA ACC ACG ACT AGA AGC\n", "----------------------------------------------------------------------------\n", "0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02\n", "----------------------------------------------------------------------------\n", "\n", "continued: \n", "============================================================================\n", " AGG AGT ATA ATC ATG ATT CAA CAC CAG CAT\n", "----------------------------------------------------------------------------\n", "0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02\n", "----------------------------------------------------------------------------\n", "\n", "continued: \n", "============================================================================\n", " CCA CCC CCG CCT CGA CGC CGG CGT CTA CTC\n", "----------------------------------------------------------------------------\n", "0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02\n", "----------------------------------------------------------------------------\n", "\n", "continued: \n", "============================================================================\n", " CTG CTT GAA GAC GAG GAT GCA GCC GCG GCT\n", "----------------------------------------------------------------------------\n", "0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02\n", "----------------------------------------------------------------------------\n", "\n", "continued: \n", "============================================================================\n", " GGA GGC GGG GGT GTA GTC GTG GTT TAC TAT\n", "----------------------------------------------------------------------------\n", "0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02\n", "----------------------------------------------------------------------------\n", "\n", "continued: \n", "============================================================================\n", " TCA TCC TCG TCT TGC TGG TGT TTA TTC TTG\n", "----------------------------------------------------------------------------\n", "0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02 0.02\n", "----------------------------------------------------------------------------\n", "\n", "continued: \n", "====\n", " TTT\n", "----\n", "0.02\n", "----" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "lf" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "A more extensive description of capabilities is in :ref:`scope-params-on-trees`." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Specifying a parameter as constant\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "sm = get_model(\"CNFGTR\")\n", "lf = sm.make_likelihood_function(tree, digits=2)\n", "lf.set_param_rule('omega', is_constant=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Providing a starting value for a parameter\n", "\n", "This can be useful to improve performance, the closer you are to the maximum likelihood estimator the quicker optimisation will be." ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "sm = get_model(\"CNFGTR\")\n", "lf = sm.make_likelihood_function(tree, digits=2)\n", "lf.set_param_rule('omega', init=0.1)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Setting parameter bounds for optimisation\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "sm = get_model(\"CNFGTR\")\n", "lf = sm.make_likelihood_function(tree, digits=2)\n", "lf.set_param_rule('omega', init=0.1, lower=1e-9, upper=20.0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Setting an upper bound for branch length\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree)\n", "lf.set_param_rule('length', upper=1.0)" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ ".. 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." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Specifying rate heterogeneity functions\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "sm = get_model(\"GTR\", with_rate=True, distribution='gamma')\n", "tree = load_tree('../data/primate_brca1.tree')\n", "lf = sm.make_likelihood_function(tree, bins=4, digits=2)\n", "lf.set_param_rule('bprobs', is_constant=True)" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "For more detailed discussion of defining and using these models see :ref:`rate-heterogeneity`." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Specifying Phylo-HMMs" ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree\n", "from cogent3.evolve.models import get_model\n", "\n", "sm = get_model(\"GTR\", with_rate=True, distribution='gamma')\n", "tree = load_tree('../data/primate_brca1.tree')\n", "lf = sm.make_likelihood_function(tree, bins=4, sites_independent=False,\n", " digits=2)\n", "lf.set_param_rule('bprobs', is_constant=True)" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "For more detailed discussion of defining and using these models see :ref:`rate-heterogeneity-hmm`." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Fitting likelihood functions - Choice of optimisers\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree, digits=3, space=2)\n", "lf.set_alignment(aln)\n", "lf.optimise(show_progress=False)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [], "source": [ "lf.optimise(local=True, max_restarts=5, show_progress=False)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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``." ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [], "source": [ "lf.optimise(local=False, global_tolerance=1.0, show_progress=False)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Followed by a standard call to `optimise()`." ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [], "source": [ "lf.optimise(show_progress=False, max_restarts=5, tolerance=1e-8)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### How to check your optimisation was successful\n", "\n", "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.\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 21, "metadata": { "raises-exception": "false" }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "FORCED EXIT from optimiser after 10 evaluations\n" ] } ], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree, digits=3, space=2)\n", "lf.set_alignment(aln)\n", "try:\n", " lf.optimise(show_progress=False, limit_action='raise',\n", " max_evaluations=10, return_calculator=True)\n", "except ArithmeticError as err:\n", " print(err)" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ ".. 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." ] }, { "cell_type": "markdown", "metadata": { "toc-hr-collapsed": true }, "source": [ "### Overview of the fitted likelihood function\n", "\n", "In Jupyter, the likelihood function object presents a representation of the main object features." ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

Likelihood function statistics

\n", "

log-likelihood = -6992.7690

\n", "

number of free parameters = 16

\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Global params
A/CA/GA/TC/GC/T
1.23165.25340.95852.31595.9700
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
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
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Motif params
ACGT
0.37570.17420.20950.2406
\n" ], "text/plain": [ "Likelihood function statistics\n", "log-likelihood = -6992.7690\n", "number of free parameters = 16\n", "==============================================\n", " A/C A/G A/T C/G C/T\n", "----------------------------------------------\n", "1.2316 5.2534 0.9585 2.3159 5.9700\n", "----------------------------------------------\n", "==============================\n", " edge parent length\n", "------------------------------\n", " Galago root 0.1731\n", " HowlerMon root 0.0449\n", " Rhesus edge.3 0.0216\n", " Orangutan edge.2 0.0077\n", " Gorilla edge.1 0.0025\n", " Human edge.0 0.0061\n", "Chimpanzee edge.0 0.0028\n", " edge.0 edge.1 0.0000\n", " edge.1 edge.2 0.0034\n", " edge.2 edge.3 0.0120\n", " edge.3 root 0.0076\n", "------------------------------\n", "====================================\n", " A C G T\n", "------------------------------------\n", "0.3757 0.1742 0.2095 0.2406\n", "------------------------------------" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "sm = get_model(\"GTR\")\n", "tree = load_tree('../data/primate_brca1.tree')\n", "lf = sm.make_likelihood_function(tree)\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "lf.set_alignment(aln)\n", "lf.optimise(local=True, show_progress=False)\n", "lf" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Log likelihood and number of free parameters" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Reusing the optimised `lf` object from above, we can get the log-likelihood and the number of free parameters." ] }, { "cell_type": "code", "execution_count": 23, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "-6992.768994254359" ] }, "execution_count": 23, "metadata": {}, "output_type": "execute_result" } ], "source": [ "lnL = lf.lnL\n", "lnL" ] }, { "cell_type": "code", "execution_count": 24, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "16" ] }, "execution_count": 24, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nfp = lf.nfp\n", "nfp" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ ".. 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)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Aikake Information Criterion\n", "\n", "Reusing the optimised `lf` object from above." ] }, { "cell_type": "code", "execution_count": 25, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "14017.537988508719" ] }, "execution_count": 25, "metadata": {}, "output_type": "execute_result" } ], "source": [ "lf.get_aic()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can also get the second-order AIC." ] }, { "cell_type": "code", "execution_count": 26, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "14017.732482609541" ] }, "execution_count": 26, "metadata": {}, "output_type": "execute_result" } ], "source": [ "lf.get_aic(second_order=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Bayesian Information Criterion\n", "\n", "Reusing the optimised `lf` object from above." ] }, { "cell_type": "code", "execution_count": 27, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "14112.615784311509" ] }, "execution_count": 27, "metadata": {}, "output_type": "execute_result" } ], "source": [ "lf.get_bic()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Getting maximum likelihood estimates\n", "\n", "Reusing the optimised `lf` object from above.\n", "\n", "##### One at a time\n", "\n", "We get the statistics out individually. We get the ``length`` for the Human edge and the exchangeability parameter ``A/G``." ] }, { "cell_type": "code", "execution_count": 28, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "5.253427222418745" ] }, "execution_count": 28, "metadata": {}, "output_type": "execute_result" } ], "source": [ "a_g = lf.get_param_value('A/G')\n", "a_g" ] }, { "cell_type": "code", "execution_count": 29, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "0.006060124812429323" ] }, "execution_count": 29, "metadata": {}, "output_type": "execute_result" } ], "source": [ "human = lf.get_param_value('length', 'Human')\n", "human" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "##### Just the motif probabilities" ] }, { "cell_type": "code", "execution_count": 30, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
TCAG
0.2410.1740.3760.209
\n" ], "text/plain": [ "====================================\n", " T C A G\n", "------------------------------------\n", "0.2406 0.1742 0.3757 0.2095\n", "------------------------------------" ] }, "execution_count": 30, "metadata": {}, "output_type": "execute_result" } ], "source": [ "mprobs = lf.get_motif_probs()\n", "mprobs" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "##### As tables" ] }, { "cell_type": "code", "execution_count": 31, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
global params
A/CA/GA/TC/GC/T
1.23165.25340.95852.31595.9700
\n", "

\n", "1 rows x 5 columns

" ], "text/plain": [ "global params\n", "==============================================\n", " A/C A/G A/T C/G C/T\n", "----------------------------------------------\n", "1.2316 5.2534 0.9585 2.3159 5.9700\n", "----------------------------------------------\n", "\n", "1 rows x 5 columns" ] }, "execution_count": 31, "metadata": {}, "output_type": "execute_result" } ], "source": [ "tables = lf.get_statistics(with_motif_probs=True, with_titles=True)\n", "tables[0] # just displaying the first" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Testing Hypotheses - Using Likelihood Ratio Tests\n", "\n", "We test the molecular clock hypothesis for human and chimpanzee lineages. The null has these two branches constrained to be equal." ] }, { "cell_type": "code", "execution_count": 32, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

Null Hypothesis

\n", "

log-likelihood = -7177.4403

\n", "

number of free parameters = 10

\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
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
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Motif params
ACGT
0.3760.1740.2090.241
\n" ], "text/plain": [ "Null Hypothesis\n", "log-likelihood = -7177.4403\n", "number of free parameters = 10\n", "==========================\n", " edge parent length\n", "--------------------------\n", " Galago root 0.167\n", " HowlerMon root 0.044\n", " Rhesus edge.3 0.021\n", " Orangutan edge.2 0.008\n", " Gorilla edge.1 0.002\n", " Human edge.0 0.004\n", "Chimpanzee edge.0 0.004\n", " edge.0 edge.1 0.000\n", " edge.1 edge.2 0.003\n", " edge.2 edge.3 0.012\n", " edge.3 root 0.009\n", "--------------------------\n", "==========================\n", " A C G T\n", "--------------------------\n", "0.376 0.174 0.209 0.241\n", "--------------------------" ] }, "execution_count": 32, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree, digits=3, space=2)\n", "lf.set_alignment(aln)\n", "lf.set_param_rule('length', tip_names=['Human', 'Chimpanzee'],\n", " outgroup_name='Galago', clade=True, is_independent=False)\n", "lf.set_name('Null Hypothesis')\n", "lf.optimise(local=True, show_progress=False)\n", "null_lnL = lf.lnL\n", "null_nfp = lf.nfp\n", "lf" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The alternate allows the human and chimpanzee branches to differ by just setting all lengths to be independent." ] }, { "cell_type": "code", "execution_count": 33, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

Alt Hypothesis

\n", "

log-likelihood = -7175.7756

\n", "

number of free parameters = 11

\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
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
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Motif params
ACGT
0.3760.1740.2090.241
\n" ], "text/plain": [ "Alt Hypothesis\n", "log-likelihood = -7175.7756\n", "number of free parameters = 11\n", "==========================\n", " edge parent length\n", "--------------------------\n", " Galago root 0.167\n", " HowlerMon root 0.044\n", " Rhesus edge.3 0.021\n", " Orangutan edge.2 0.008\n", " Gorilla edge.1 0.002\n", " Human edge.0 0.006\n", "Chimpanzee edge.0 0.003\n", " edge.0 edge.1 0.000\n", " edge.1 edge.2 0.003\n", " edge.2 edge.3 0.012\n", " edge.3 root 0.009\n", "--------------------------\n", "==========================\n", " A C G T\n", "--------------------------\n", "0.376 0.174 0.209 0.241\n", "--------------------------" ] }, "execution_count": 33, "metadata": {}, "output_type": "execute_result" } ], "source": [ "lf.set_param_rule('length', is_independent=True)\n", "lf.set_name('Alt Hypothesis')\n", "lf.optimise(local=True, show_progress=False)\n", "alt_lnL = lf.lnL\n", "alt_nfp = lf.nfp\n", "lf" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 34, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "LR=3.3294 ; df=1; p=1.0000\n" ] } ], "source": [ "from cogent3.maths.stats import chisqprob\n", "\n", "LR = 2 * (alt_lnL - null_lnL) # the likelihood ratio statistic\n", "df = (alt_nfp - null_nfp) # the test degrees of freedom\n", "p = chisqprob(LR, df)\n", "print(f'LR={LR:.4f} ; df={df}; p={df:.4f}')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Testing Hypotheses - By parametric bootstrapping" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "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 :ref:`parametric-bootstrap`.\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": 35, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
0
RhesusACATGGCAGCGTCGACAAGGATTCATTATAGTTCTTAACTTAGAGACATGAACAGAGATG
Chimpanzee...................................A........................
Galago..CA....A........T.............G..C.T..C.....GAC....G.......
Gorilla............................................................
HowlerMon.........................................................T..
Human.....................................................T......
Orangutan............................................................
\n", "

7 x 60 dna alignment

\n", "" ], "text/plain": [ "7 x 60 dna alignment: Chimpanzee[ACATGGCAGCG...], Galago[ACCAGGCAACG...], Gorilla[ACATGGCAGCG...], ..." ] }, "execution_count": 35, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree, digits=3, space=2)\n", "lf.set_alignment(aln)\n", "lf.set_param_rule('length', tip_names=['Human', 'Chimpanzee'],\n", " outgroup_name='Galago', clade=True, is_independent=False)\n", "lf.set_name('Null Hypothesis')\n", "lf.optimise(local=True, show_progress=False)\n", "sim_aln = lf.simulate_alignment()\n", "sim_aln[:60]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Determining confidence intervals on MLEs\n", "\n", "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)." ] }, { "cell_type": "code", "execution_count": 36, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "lo=3.78 ; mle=4.44 ; hi=5.22\n", "lo=0.00 ; mle=0.01 ; hi=0.01\n" ] } ], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "sm = get_model(\"HKY85\")\n", "lf = sm.make_likelihood_function(tree)\n", "lf.set_alignment(aln)\n", "lf.optimise(local=True, show_progress=False)\n", "kappa_lo, kappa_mle, kappa_hi = lf.get_param_interval('kappa')\n", "print(f\"lo={kappa_lo:.2f} ; mle={kappa_mle:.2f} ; hi={kappa_hi:.2f}\")\n", "human_lo, human_mle, human_hi = lf.get_param_interval('length', 'Human')\n", "print(f\"lo={human_lo:.2f} ; mle={human_mle:.2f} ; hi={human_hi:.2f}\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Saving results\n", "\n", "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..\n", "\n", "To illustrate this, I create a very simple likelihood function. The ``json`` variable below is just a string that can be saved to disk." ] }, { "cell_type": "code", "execution_count": 37, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "'{\"model\": {\"alphabet\": {\"motifset\": [\"T\", \"C\", \"A\", \"G\"], \"g'" ] }, "execution_count": 37, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "aln = make_aligned_seqs(data=dict(a=\"ACGG\", b=\"ATAG\", c=\"ATGG\"))\n", "tree = make_tree(tip_names=aln.names)\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree)\n", "lf.set_alignment(aln)\n", "json = lf.to_json()\n", "json[:60] # just truncating the displayed string" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We deserialise the object from the string." ] }, { "cell_type": "code", "execution_count": 38, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

unnamed

\n", "

log-likelihood = -14.2727

\n", "

number of free parameters = 3

\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Edge params
edgeparentlength
aroot1.0000
broot1.0000
croot1.0000
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
Motif params
ACGT
0.33330.08330.41670.1667
\n" ], "text/plain": [ "unnamed\n", "log-likelihood = -14.2727\n", "number of free parameters = 3\n", "========================\n", "edge parent length\n", "------------------------\n", " a root 1.0000\n", " b root 1.0000\n", " c root 1.0000\n", "------------------------\n", "====================================\n", " A C G T\n", "------------------------------------\n", "0.3333 0.0833 0.4167 0.1667\n", "------------------------------------" ] }, "execution_count": 38, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cogent3.util.deserialise import deserialise_object\n", "\n", "newlf = deserialise_object(json)\n", "newlf" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Reconstructing ancestral sequences\n", "\n", "We first fit a likelihood function." ] }, { "cell_type": "code", "execution_count": 39, "metadata": {}, "outputs": [], "source": [ "from cogent3 import load_tree, load_aligned_seqs\n", "from cogent3.evolve.models import get_model\n", "\n", "tree = load_tree('../data/primate_brca1.tree')\n", "aln = load_aligned_seqs('../data/primate_brca1.fasta')\n", "sm = get_model(\"F81\")\n", "lf = sm.make_likelihood_function(tree, digits=3, space=2)\n", "lf.set_alignment(aln)\n", "lf.optimise(show_progress=False)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We then get the most likely ancestral sequences." ] }, { "cell_type": "code", "execution_count": 40, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
0
rootTGTGGCACAAATACTCATGCCAGCTCATTACAGCATGAGAACAGTTTGTTACTCACTAAA
edge.0...............................................A............
edge.1...............................................A............
edge.2...............................................A............
edge.3............................................................
\n", "

5 x 60 dna alignment

\n", "" ], "text/plain": [ "5 x 60 dna alignment: edge.0[TGTGGCACAAA...], edge.1[TGTGGCACAAA...], edge.2[TGTGGCACAAA...], ..." ] }, "execution_count": 40, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ancestors = lf.likely_ancestral_seqs()\n", "ancestors[:60]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Or we can get the posterior probabilities (returned as a ``DictArray``) of sequence states at each node." ] }, { "cell_type": "code", "execution_count": 41, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
TCAG
00.1820.0000.0000.000
10.0000.0000.0000.156
20.1820.0000.0000.000
30.0000.0000.0000.156
40.0000.0000.0000.156
\n" ], "text/plain": [ "=========================================\n", " T C A G\n", "-----------------------------------------\n", "0 0.1816 0.0000 0.0000 0.0000\n", "1 0.0000 0.0000 0.0000 0.1561\n", "2 0.1816 0.0000 0.0000 0.0000\n", "3 0.0000 0.0000 0.0000 0.1561\n", "4 0.0000 0.0000 0.0000 0.1561\n", "-----------------------------------------" ] }, "execution_count": 41, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ancestral_probs = lf.reconstruct_ancestral_seqs()\n", "ancestral_probs[\"root\"][:5]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Tips for improved performance -- sequentially build the fitting\n", "\n", "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." ] } ], "metadata": { "celltoolbar": "Raw Cell Format", "kernelspec": { "display_name": "Python [conda env:c3dev] *", "language": "python", "name": "conda-env-c3dev-py" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.8.1" }, "pycharm": { "stem_cell": { "cell_type": "raw", "metadata": { "collapsed": false }, "source": [] } }, "widgets": { "application/vnd.jupyter.widget-state+json": { "state": {}, "version_major": 2, "version_minor": 0 } } }, "nbformat": 4, "nbformat_minor": 4 }