The analysis of synteny (i.e., conserved gene content and order in a genomic segment across species) can help understand the trajectory of duplicated genes through evolution. In particular, synteny analyses are widely used to investigate the evolution of genes derived from whole-genome duplication (WGD) events, as they can reveal genomic rearrangements that happened following the duplication of all chromosomes. However, synteny analysis are typically performed in a pairwise manner, which can be difficult to explore and interpret when comparing several species. To understand global patterns of synteny, Zhao and Schranz (2017) proposed a network-based approach to analyze synteny. In synteny networks, genes in a given syntenic block are represented as nodes connected by an edge. Synteny networks have been used to explore, among others, global synteny patterns in mammalian and angiosperm genomes (Zhao and Schranz 2019), the evolution of MADS-box transcription factors (Zhao et al. 2017), and infer a microsynteny-based phylogeny for angiosperms (Zhao et al. 2021). syntenet is a package that can be used to infer synteny networks from protein sequences and perform downstream network analyses that include:
Network clustering;
Phylogenomic profiling, which consists in identifying which species contain which clusters. This analysis can reveal highly conserved synteny clusters and taxon-specific ones (e.g., family- and order-specific clusters);
Microsynteny-based phylogeny reconstruction with maximum likelihood, which can be achieved by inferring a phylogeny from a binary matrix of phylogenomic profiles with IQTREE.
if(!requireNamespace('BiocManager', quietly = TRUE))
install.packages('BiocManager')
BiocManager::install("syntenet")
# Load package after installation
library(syntenet)
For this vignette, we will use the proteomes and gene annotation of the algae species Ostreococcus lucimarinus and Ostreococcus sp RCC809, which were obtained from Pico-PLAZA 3.0 (Vandepoele et al. 2013).
# Protein sequences
data(proteomes)
head(proteomes)
#> $Olucimarinus
#> AAStringSet object of length 1901:
#> width seq names
#> [1] 911 MTTMADERASIARVSVVKYGAI...VQLYTYPGSTNDPNFLLKLA* OL01G00010
#> [2] 789 MGGRRCFCSRSSPVGVGAPAPA...PPQCGADIEAGSEPPPDKCG* OL01G00020
#> [3] 618 MTRAKDAIVVDDGNDDDDDDDD...RDASASLALALAFSSEESVV* OL01G00030
#> [4] 547 MPTKAQCWVVSYARVRDGASRS...TGSVSARASIFGEQASFRKA* OL01G00040
#> [5] 319 MFTASHTTSKVTLRARVATQPR...HNGMALWRETTPKDSLIPAL* OL01G00050
#> ... ... ...
#> [1897] 106 MAANDGETKLPEDGWIQPCFAC...RAIVDQVGGEHLKGSLMPIE* OL03G05910
#> [1898] 70 RMGIVKLATDGSVWVHSPIELD...QQWKDAYPGATLYACPGLKSK OL03G05920
#> [1899] 680 MDDAHDARWATTSARDGERARA...RSVGPSASDKILEALFPVAD* OL03G05930
#> [1900] 179 MRAVRERSKANLAARVKEEATR...ELERTRELFARARVRAYECI* OL03G05940
#> [1901] 83 MFVREARRAIPRFIKDPPQAFH...ESGDVRSVEGEVCGAVLVDE* OL03G05950
#>
#> $OspRCC809
#> AAStringSet object of length 1433:
#> width seq names
#> [1] 274 MASTTGSAARRVFVDVEKTVNG...DVLSLGQGSLSGESSSSDEE* ORCC809_01G00010
#> [2] 175 MDQMRAANAQRSYLLFFVLFFL...SRRLLGRLDSEHTDLHPSWR* ORCC809_01G00020
#> [3] 403 MTAPRVRASRRATATAAATVTA...LTERDLRYMEPKATIEEWMG* ORCC809_01G00030
#> [4] 217 MTIDADGDDTLAPHAPAHGEVS...LIRLRGVEKTPTVDPPPPPP* ORCC809_01G00040
#> [5] 1691 RIEADEKSLLVFGKESPVRTAC...VRMGNNVVTSRYASSESEEDV ORCC809_01G00050
#> ... ... ...
#> [1429] 428 MVDANATTQTFVLEAEQELRVE...DLPSNVLLVGNLKWLGEDGK* ORCC809_03G02980
#> [1430] 378 MSVPRTTLRRIPLGNARDVLVT...ETLKAIDAVHAQCRDPCIAT* ORCC809_03G02990
#> [1431] 1156 MRATSAPSIVSFVARVACLFVA...CAFGTSLASFVVERARRLEN* ORCC809_03G03000
#> [1432] 541 MAITVFLTDHGRRASALTFLVV...GFGVGAVKFMLAPEMVKSLA* ORCC809_03G03010
#> [1433] 289 MSLSSLRSFSRSISSAPGGRSC...EPEEPEPEEPEPEEPEEPEP* ORCC809_03G03020
# Annotation (ranges)
data(annotation)
head(annotation)
#> GRangesList object of length 2:
#> $Olucimarinus
#> GRanges object with 1903 ranges and 4 metadata columns:
#> seqnames ranges strand | type ID Name
#> <Rle> <IRanges> <Rle> | <factor> <character> <character>
#> [1] Chr_1 939-3671 - | gene OL01G00010 OL01G00010
#> [2] Chr_1 3907-6927 + | gene OL01G00020 OL01G00020
#> [3] Chr_1 7085-9160 + | gene OL01G00030 OL01G00030
#> [4] Chr_1 9830-11480 + | gene OL01G00040 OL01G00040
#> [5] Chr_1 11467-12599 - | gene OL01G00050 OL01G00050
#> ... ... ... ... . ... ... ...
#> [1899] Chr_3 977435-977752 - | gene OL03G05910 OL03G05910
#> [1900] Chr_3 978702-978911 - | gene OL03G05920 OL03G05920
#> [1901] Chr_3 979281-981320 - | gene OL03G05930 OL03G05930
#> [1902] Chr_3 981778-982314 + | gene OL03G05940 OL03G05940
#> [1903] Chr_3 982498-982746 + | gene OL03G05950 OL03G05950
#> gene_id
#> <character>
#> [1] OL01G00010
#> [2] OL01G00020
#> [3] OL01G00030
#> [4] OL01G00040
#> [5] OL01G00050
#> ... ...
#> [1899] OL03G05910
#> [1900] OL03G05920
#> [1901] OL03G05930
#> [1902] OL03G05940
#> [1903] OL03G05950
#> -------
#> seqinfo: 6 sequences from an unspecified genome; no seqlengths
#>
#> $OspRCC809
#> GRanges object with 1433 ranges and 4 metadata columns:
#> seqnames ranges strand | type ID
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> [1] chr_1 321-1142 - | gene ORCC809_01G00010
#> [2] chr_1 1463-2089 + | gene ORCC809_01G00020
#> [3] chr_1 2162-3370 - | gene ORCC809_01G00030
#> [4] chr_1 3774-4424 - | gene ORCC809_01G00040
#> [5] chr_1 4693-9924 - | gene ORCC809_01G00050
#> ... ... ... ... . ... ...
#> [1429] chr_3 504915-506198 - | gene ORCC809_03G02980
#> [1430] chr_3 506377-507510 + | gene ORCC809_03G02990
#> [1431] chr_3 507856-511323 + | gene ORCC809_03G03000
#> [1432] chr_3 511533-513155 - | gene ORCC809_03G03010
#> [1433] chr_3 513841-514707 + | gene ORCC809_03G03020
#> Name gene_id
#> <character> <character>
#> [1] ORCC809_01G00010 ORCC809_01G00010
#> [2] ORCC809_01G00020 ORCC809_01G00020
#> [3] ORCC809_01G00030 ORCC809_01G00030
#> [4] ORCC809_01G00040 ORCC809_01G00040
#> [5] ORCC809_01G00050 ORCC809_01G00050
#> ... ... ...
#> [1429] ORCC809_03G02980 ORCC809_03G02980
#> [1430] ORCC809_03G02990 ORCC809_03G02990
#> [1431] ORCC809_03G03000 ORCC809_03G03000
#> [1432] ORCC809_03G03010 ORCC809_03G03010
#> [1433] ORCC809_03G03020 ORCC809_03G03020
#> -------
#> seqinfo: 6 sequences from an unspecified genome; no seqlengths
To detect synteny and infer synteny networks, syntenet requires two objects as input:
AAStringSet
objects containing the translated sequences
of primary transcripts for each species.GRangesList
or CompressedGRangesList
object containing
the coordinates for the genes in seq.If you have whole-genome protein sequences in FASTA files, store all FASTA
files in the same directory and use the function fasta2AAStringSetlist()
to
read all FASTA files into a list of AAStringSet
objects.
Likewise, if you have gene annotation in GFF/GFF3/GTF files,
store all files in the same directory and use the function gff2GRangesList()
to read all GFF/GFF3/GTF files into a GRangesList object
.
For a demonstration, we will read example FASTA and GFF3 files stored in
subdirectories named sequences/ and annotation/, which are located
in the extdata/
directory of this package.
AAStringSet
objectsHere is how you can use fasta2AAStringSetlist()
to read FASTA files
in a directory as a list of AAStringSet
objects.
# Path to directory containing FASTA files
fasta_dir <- system.file("extdata", "sequences", package = "syntenet")
fasta_dir
#> [1] "/tmp/RtmpjZSSnQ/Rinst27476741df4097/syntenet/extdata/sequences"
dir(fasta_dir) # see the contents of the directory
#> [1] "Olucimarinus.fa.gz" "OspRCC809.fa.gz"
# Read all FASTA files in `fasta_dir`
aastringsetlist <- fasta2AAStringSetlist(fasta_dir)
aastringsetlist
#> $Olucimarinus
#> AAStringSet object of length 100:
#> width seq names
#> [1] 911 MTTMADERASIARVSVVKYGAI...DVQLYTYPGSTNDPNFLLKLA* OL01G00010
#> [2] 789 MGGRRCFCSRSSPVGVGAPAPA...FPPQCGADIEAGSEPPPDKCG* OL01G00020
#> [3] 618 MTRAKDAIVVDDGNDDDDDDDD...DRDASASLALALAFSSEESVV* OL01G00030
#> [4] 547 MPTKAQCWVVSYARVRDGASRS...VTGSVSARASIFGEQASFRKA* OL01G00040
#> [5] 319 MFTASHTTSKVTLRARVATQPR...LHNGMALWRETTPKDSLIPAL* OL01G00050
#> ... ... ...
#> [96] 476 MVPARNFLDGANAREVELDRVV...VMRKLREPDSVARLAGQTGVR* OL01G00960
#> [97] 771 MARHRGTRGGWNATTTEGGDGR...SIPDDGFDESSSVSASTIDGF* OL01G00970
#> [98] 494 MDSEFWGCVIPAGRAVRVEVAT...FIKSRKDLFTIDGAYVRLVKK* OL01G00980
#> [99] 264 VRAIVGATTRIQTRAPPRANHR...DWSFISDEFQDDASDSEVIDR* OL01G00990
#> [100] 565 MQLDAFRKATVKGVATRVGGAD...QLADLLRKNMGVPAKFIDAQN* OL01G01000
#>
#> $OspRCC809
#> AAStringSet object of length 100:
#> width seq names
#> [1] 274 MASTTGSAARRVFVDVEKTVNG...WDVLSLGQGSLSGESSSSDEE* ORCC809_01G00010
#> [2] 175 MDQMRAANAQRSYLLFFVLFFL...SSRRLLGRLDSEHTDLHPSWR* ORCC809_01G00020
#> [3] 403 MTAPRVRASRRATATAAATVTA...ALTERDLRYMEPKATIEEWMG* ORCC809_01G00030
#> [4] 217 MTIDADGDDTLAPHAPAHGEVS...SLIRLRGVEKTPTVDPPPPPP* ORCC809_01G00040
#> [5] 1691 RIEADEKSLLVFGKESPVRTAC...SVRMGNNVVTSRYASSESEEDV ORCC809_01G00050
#> ... ... ...
#> [96] 357 MSRGLADNWDDAEGYYCARIGE...TVNEALQHPFIVERIRTTAPN* ORCC809_01G00960
#> [97] 164 MAMDSFRSAPRSRRRVEATSRE...SKPVKPVREPVRMVEASTGAH* ORCC809_01G00970
#> [98] 85 MPEGTVFIGNIPYDATESSLTE...NLNAREYNGRQLRVDHAETMKG ORCC809_01G00980
#> [99] 229 MKGGGGASGAAASANGNGAVGG...PDQRAQVEYLRQLAAQQGMVR* ORCC809_01G00990
#> [100] 103 RKAGGERWEDSSLAEWPENDFR...EMAGKFIGNRPVKLRKSAWNER ORCC809_01G01000
And that’s it! Now you have a list of AAStringSet
objects.
GRangesList
objectHere is how you can use gff2GRangesList()
to read GFF/GFF3/GTF files
in a directory as a GRangesList
object.
# Path to directory containing FASTA files
gff_dir <- system.file("extdata", "annotation", package = "syntenet")
gff_dir
#> [1] "/tmp/RtmpjZSSnQ/Rinst27476741df4097/syntenet/extdata/annotation"
dir(gff_dir) # see the contents of the directory
#> [1] "Olucimarinus.gff3.gz" "OspRCC809.gff3.gz"
# Read all FASTA files in `fasta_dir`
grangeslist <- gff2GRangesList(gff_dir)
grangeslist
#> GRangesList object of length 2:
#> $Olucimarinus
#> GRanges object with 100 ranges and 7 metadata columns:
#> seqnames ranges strand | source type score
#> <Rle> <IRanges> <Rle> | <factor> <factor> <numeric>
#> [1] Chr_1 939-3671 - | rtracklayer gene NA
#> [2] Chr_1 3907-6927 + | rtracklayer gene NA
#> [3] Chr_1 7085-9160 + | rtracklayer gene NA
#> [4] Chr_1 9830-11480 + | rtracklayer gene NA
#> [5] Chr_1 11467-12599 - | rtracklayer gene NA
#> ... ... ... ... . ... ... ...
#> [96] Chr_1 170975-172402 + | rtracklayer gene NA
#> [97] Chr_1 172445-174757 - | rtracklayer gene NA
#> [98] Chr_1 175358-176839 + | rtracklayer gene NA
#> [99] Chr_1 176901-177692 - | rtracklayer gene NA
#> [100] Chr_1 177742-179436 - | rtracklayer gene NA
#> phase ID Name gene_id
#> <integer> <character> <character> <character>
#> [1] <NA> OL01G00010 OL01G00010 OL01G00010
#> [2] <NA> OL01G00020 OL01G00020 OL01G00020
#> [3] <NA> OL01G00030 OL01G00030 OL01G00030
#> [4] <NA> OL01G00040 OL01G00040 OL01G00040
#> [5] <NA> OL01G00050 OL01G00050 OL01G00050
#> ... ... ... ... ...
#> [96] <NA> OL01G00960 OL01G00960 OL01G00960
#> [97] <NA> OL01G00970 OL01G00970 OL01G00970
#> [98] <NA> OL01G00980 OL01G00980 OL01G00980
#> [99] <NA> OL01G00990 OL01G00990 OL01G00990
#> [100] <NA> OL01G01000 OL01G01000 OL01G01000
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
#>
#> $OspRCC809
#> GRanges object with 100 ranges and 7 metadata columns:
#> seqnames ranges strand | source type score
#> <Rle> <IRanges> <Rle> | <factor> <factor> <numeric>
#> [1] chr_1 321-1142 - | rtracklayer gene NA
#> [2] chr_1 1463-2089 + | rtracklayer gene NA
#> [3] chr_1 2162-3370 - | rtracklayer gene NA
#> [4] chr_1 3774-4424 - | rtracklayer gene NA
#> [5] chr_1 4693-9924 - | rtracklayer gene NA
#> ... ... ... ... . ... ... ...
#> [96] chr_1 165459-166529 - | rtracklayer gene NA
#> [97] chr_1 166654-167213 - | rtracklayer gene NA
#> [98] chr_1 167296-167550 + | rtracklayer gene NA
#> [99] chr_1 167542-168228 + | rtracklayer gene NA
#> [100] chr_1 168639-168947 - | rtracklayer gene NA
#> phase ID Name gene_id
#> <integer> <character> <character> <character>
#> [1] <NA> ORCC809_01G00010 ORCC809_01G00010 ORCC809_01G00010
#> [2] <NA> ORCC809_01G00020 ORCC809_01G00020 ORCC809_01G00020
#> [3] <NA> ORCC809_01G00030 ORCC809_01G00030 ORCC809_01G00030
#> [4] <NA> ORCC809_01G00040 ORCC809_01G00040 ORCC809_01G00040
#> [5] <NA> ORCC809_01G00050 ORCC809_01G00050 ORCC809_01G00050
#> ... ... ... ... ...
#> [96] <NA> ORCC809_01G00960 ORCC809_01G00960 ORCC809_01G00960
#> [97] <NA> ORCC809_01G00970 ORCC809_01G00970 ORCC809_01G00970
#> [98] <NA> ORCC809_01G00980 ORCC809_01G00980 ORCC809_01G00980
#> [99] <NA> ORCC809_01G00990 ORCC809_01G00990 ORCC809_01G00990
#> [100] <NA> ORCC809_01G01000 ORCC809_01G01000 ORCC809_01G01000
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
And now you have a GRangesList
object.
The first part of the pipeline consists in processing the data to make it
match a standard structure. However, before processing the data for synteny
detection, you must use the function check_input()
to check if your data can
enter the pipeline. This function checks the input data for 3
required conditions:
Names of seq list (i.e., names(seq)
) match
the names of annotation GRangesList
/CompressedGRangesList
(i.e., names(annotation)
)
For each species (list elements), the number of sequences in seq is not greater than the number of genes in annotation. This is a way to ensure users do not input the translated sequences for multiple isoforms of the same gene (generated by alternative splicing). Ideally, the number of sequences in seq should be equal to the number of genes in annotation, but this may not always stand true because of non-protein-coding genes.
For each species, sequence names (i.e., names(seq[[x]])
,
equivalent to FASTA headers) match gene names in annotation
.
Let’s check if the example data sets satisfy these 3 criteria:
check_input(proteomes, annotation)
#> [1] TRUE
As you can see, the data passed the checks. Now, let’s process them
with the function process_input()
.
This function processes the input sequences and annotation to:
Remove whitespace and anything after it in sequence names
(i.e., names(seq[[x]])
, which is equivalent to FASTA headers), if
there is any.
Remove period followed by number at the end of sequence names, which typically indicates different isoforms of the same gene (e.g., Arabidopsis thaliana’s transcript AT1G01010.1, which belongs to gene AT1G01010).
Add a unique species identifier to sequence names. The species identifier consists of the first 3-5 strings of the element name. For instance, if the first element of the seq list is named “Athaliana”, each sequence in it will have an identifier “Atha_” added to the beginning of each gene name (e.g., Atha_AT1G01010).
If sequences have an asterisk (*) representing stop codon, remove it.
Add a unique species identifier (same as above) to
gene and chromosome names of each element of the annotation
GRangesList
/CompressedGRangesList
.
Filter each element of the annotation
GRangesList
/CompressedGRangesList
to keep only seqnames,
ranges, and gene ID.
Let’s process our input data:
pdata <- process_input(proteomes, annotation)
# Looking at the processed data
pdata$seq
#> $Olucimarinus
#> AAStringSet object of length 1901:
#> width seq names
#> [1] 910 MTTMADERASIARVSVVKYGAI...DVQLYTYPGSTNDPNFLLKLA Olu_OL01G00010
#> [2] 788 MGGRRCFCSRSSPVGVGAPAPA...FPPQCGADIEAGSEPPPDKCG Olu_OL01G00020
#> [3] 617 MTRAKDAIVVDDGNDDDDDDDD...DRDASASLALALAFSSEESVV Olu_OL01G00030
#> [4] 546 MPTKAQCWVVSYARVRDGASRS...VTGSVSARASIFGEQASFRKA Olu_OL01G00040
#> [5] 318 MFTASHTTSKVTLRARVATQPR...LHNGMALWRETTPKDSLIPAL Olu_OL01G00050
#> ... ... ...
#> [1897] 105 MAANDGETKLPEDGWIQPCFAC...LRAIVDQVGGEHLKGSLMPIE Olu_OL03G05910
#> [1898] 69 RMGIVKLATDGSVWVHSPIELD...AQQWKDAYPGATLYACPGLKS Olu_OL03G05920
#> [1899] 679 MDDAHDARWATTSARDGERARA...ARSVGPSASDKILEALFPVAD Olu_OL03G05930
#> [1900] 178 MRAVRERSKANLAARVKEEATR...LELERTRELFARARVRAYECI Olu_OL03G05940
#> [1901] 82 MFVREARRAIPRFIKDPPQAFH...QESGDVRSVEGEVCGAVLVDE Olu_OL03G05950
#>
#> $OspRCC809
#> AAStringSet object of length 1433:
#> width seq names
#> [1] 273 MASTTGSAARRVFVDVEKTVNG...WDVLSLGQGSLSGESSSSDEE Osp_ORCC809_01G00010
#> [2] 174 MDQMRAANAQRSYLLFFVLFFL...SSRRLLGRLDSEHTDLHPSWR Osp_ORCC809_01G00020
#> [3] 402 MTAPRVRASRRATATAAATVTA...ALTERDLRYMEPKATIEEWMG Osp_ORCC809_01G00030
#> [4] 216 MTIDADGDDTLAPHAPAHGEVS...SLIRLRGVEKTPTVDPPPPPP Osp_ORCC809_01G00040
#> [5] 1690 RIEADEKSLLVFGKESPVRTAC...SVRMGNNVVTSRYASSESEED Osp_ORCC809_01G00050
#> ... ... ...
#> [1429] 427 MVDANATTQTFVLEAEQELRVE...GDLPSNVLLVGNLKWLGEDGK Osp_ORCC809_03G02980
#> [1430] 377 MSVPRTTLRRIPLGNARDVLVT...KETLKAIDAVHAQCRDPCIAT Osp_ORCC809_03G02990
#> [1431] 1155 MRATSAPSIVSFVARVACLFVA...ACAFGTSLASFVVERARRLEN Osp_ORCC809_03G03000
#> [1432] 540 MAITVFLTDHGRRASALTFLVV...PGFGVGAVKFMLAPEMVKSLA Osp_ORCC809_03G03010
#> [1433] 288 MSLSSLRSFSRSISSAPGGRSC...EEPEEPEPEEPEPEEPEEPEP Osp_ORCC809_03G03020
pdata$annotation
#> $Olucimarinus
#> GRanges object with 1903 ranges and 1 metadata column:
#> seqnames ranges strand | gene
#> <Rle> <IRanges> <Rle> | <character>
#> [1] Olu_Chr_1 939-3671 * | Olu_OL01G00010
#> [2] Olu_Chr_1 3907-6927 * | Olu_OL01G00020
#> [3] Olu_Chr_1 7085-9160 * | Olu_OL01G00030
#> [4] Olu_Chr_1 9830-11480 * | Olu_OL01G00040
#> [5] Olu_Chr_1 11467-12599 * | Olu_OL01G00050
#> ... ... ... ... . ...
#> [1899] Olu_Chr_3 977435-977752 * | Olu_OL03G05910
#> [1900] Olu_Chr_3 978702-978911 * | Olu_OL03G05920
#> [1901] Olu_Chr_3 979281-981320 * | Olu_OL03G05930
#> [1902] Olu_Chr_3 981778-982314 * | Olu_OL03G05940
#> [1903] Olu_Chr_3 982498-982746 * | Olu_OL03G05950
#> -------
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
#>
#> $OspRCC809
#> GRanges object with 1433 ranges and 1 metadata column:
#> seqnames ranges strand | gene
#> <Rle> <IRanges> <Rle> | <character>
#> [1] Osp_chr_1 321-1142 * | Osp_ORCC809_01G00010
#> [2] Osp_chr_1 1463-2089 * | Osp_ORCC809_01G00020
#> [3] Osp_chr_1 2162-3370 * | Osp_ORCC809_01G00030
#> [4] Osp_chr_1 3774-4424 * | Osp_ORCC809_01G00040
#> [5] Osp_chr_1 4693-9924 * | Osp_ORCC809_01G00050
#> ... ... ... ... . ...
#> [1429] Osp_chr_3 504915-506198 * | Osp_ORCC809_03G02980
#> [1430] Osp_chr_3 506377-507510 * | Osp_ORCC809_03G02990
#> [1431] Osp_chr_3 507856-511323 * | Osp_ORCC809_03G03000
#> [1432] Osp_chr_3 511533-513155 * | Osp_ORCC809_03G03010
#> [1433] Osp_chr_3 513841-514707 * | Osp_ORCC809_03G03020
#> -------
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
Now that we have our processed data, we can infer the synteny network. To detect synteny, we need the tabular output from BLASTp (Altschul et al. 1997) or similar programs. Here, we give you two options:
Running BLASTp, DIAMOND (Buchfink, Reuter, and Drost 2021) or similar programs
on the command line, saving the output in tabular format and reading
the output table as a data frame
(e.g., with the base R function read.table()
).
Using a helper function in syntenet
named run_diamond()
, which runs DIAMOND from the R session and
automatically parses its output to a data frame.
Here, we will use run_diamond()
to demonstrate how it works.
Needless to say, you need to have DIAMOND installed in your machine
and in your PATH to run this function. To check if you have DIAMOND installed,
use the function diamond_is_installed()
.
data(blast_list)
if(diamond_is_installed()) {
blast_list <- run_diamond(seq = pdata$seq)
}
The output of run_diamond()
is a list of data frames containing the tabular
output of pairwise BLAST-like searches. Let’s take a look at it.
# List names
names(blast_list)
#> [1] "Olucimarinus_Olucimarinus" "Olucimarinus_OspRCC809"
#> [3] "OspRCC809_Olucimarinus" "OspRCC809_OspRCC809"
# Inspect first data frame
head(blast_list$Olucimarinus_Olucimarinus)
#> query db perc_identity length mismatches gap_open qstart
#> 1 Olu_OL01G00010 Olu_OL01G00010 100 910 0 0 1
#> 3 Olu_OL01G00020 Olu_OL01G00020 100 788 0 0 1
#> 5 Olu_OL01G00030 Olu_OL01G00030 100 617 0 0 1
#> 6 Olu_OL01G00040 Olu_OL01G00040 100 546 0 0 1
#> 7 Olu_OL01G00050 Olu_OL01G00050 100 318 0 0 1
#> 8 Olu_OL01G00060 Olu_OL01G00060 100 361 0 0 1
#> qend tstart tend evalue bitscore
#> 1 910 1 910 0.00e+00 1623
#> 3 788 1 788 0.00e+00 1363
#> 5 617 1 617 0.00e+00 1046
#> 6 546 1 546 0.00e+00 1035
#> 7 318 1 318 2.98e-219 595
#> 8 361 1 361 1.72e-267 721
Now we can use this list of BLAST-like data frames to detect synteny. Here, we reimplemented the popular MCScanX algorithm (Wang et al. 2012), originally written in C++, using the Rcpp (Eddelbuettel and François 2011) framework for R and C++ integration. This means that syntenet comes with a native version of MCScanX, so you can run MCScanX in R without having to install it yourself. Amazing, huh?
To detect synteny and infer the synteny network, use the
function infer_syntenet()
. The output is a network represented as a
so-called edge list (i.e., a 2-column data frame with node 1 and node 2
in columns 1 and 2, respectively).
# Infer synteny network
net <- infer_syntenet(blast_list, pdata$annotation)
# Look at the output
head(net)
#> Anchor1 Anchor2
#> 1 Olu_OL01G00100 Osp_ORCC809_01G06480
#> 2 Olu_OL01G00130 Osp_ORCC809_01G06440
#> 3 Olu_OL01G00150 Osp_ORCC809_01G06420
#> 4 Olu_OL01G00160 Osp_ORCC809_01G06410
#> 5 Olu_OL01G00170 Osp_ORCC809_01G06400
#> 6 Olu_OL01G00180 Osp_ORCC809_01G06390
In a synteny network, each row of the edge list represents an anchor pair. In the edge list above, for example, the genes Olu_OL01G00100 and Osp_ORCC809_01G06480 are in the same syntenic block.
After inferring the synteny network, the first thing you would want to do is cluster your network and identify which phylogenetic groups are contained in each cluster. This is what we call phylogenomic profiling. This way, you can identify clade-specific clusters, and highly conserved clusters, for instance. Here, we will use an example network of BUSCO genes for 25 eudicot species, which was obtained from Zhao and Schranz (2019).
To obtain the phylogenomic profiles, you first need to cluster your network.
This can be done with cluster_network()
.1 Friendly tip: syntenet uses the Infomap algorithm
to cluster networks, which has been shown to have the best performance
(Zhao and Schranz 2019). This algorithm assigns each gene to a single cluster.
However, for some cases (e.g., detection of tandem arrays),
you may want to use an algorithm that allows community overlap
(i.e., a gene can be part of more than one cluster).
If this is your case, we recommend the clique percolation algorithm,
which is implemented in the R package
CliquePercolation (Lange 2021).
# Load example data and take a look at it
data(network)
head(network)
#> node1 node2
#> 1 cca_23646 Lang_109327134
#> 2 cca_23646 Lang_109328075
#> 3 cca_23646 Mnot_21394516
#> 4 cca_23646 Zjuj_107413994
#> 5 cca_23646 adu_Aradu.8SN53
#> 6 cca_23646 car_14082.1
# Cluster network
clusters <- cluster_network(network)
head(clusters)
#> Gene Cluster
#> 1 cca_23646 1
#> 2 cca_23668 2
#> 3 cca_32926 3
#> 4 cca_26186 4
#> 5 cca_24381 5
#> 6 cca_24396 6
Now that each gene has been assigned to a cluster, we can identify the phylogenomic profiles of each cluster. This function returns a list of 2 elements:
A matrix of phylogenomic profiles, with clusters in rows and species in columns.
An hclust
object containing the results of Ward’s clustering on a matrix
of Ruzicka distances obtained from the phylogenomic profiles matrix.
# Phylogenomic profiling
profiles <- phylogenomic_profile(clusters)
# Exploring the output
str(profiles)
#> List of 2
#> $ profile_matrix: int [1:2411, 1:25] 1 1 0 1 0 0 0 1 1 1 ...
#> ..- attr(*, "dimnames")=List of 2
#> .. ..$ : chr [1:2411] "335" "262" "1142" "1202" ...
#> .. ..$ : chr [1:25] "Lang" "Mnot" "Zjuj" "adu" ...
#> $ hclust :List of 7
#> ..$ merge : int [1:2410, 1:2] -18 -1125 -21 -1010 -36 -83 -1174 -115 -125 -127 ...
#> ..$ height : num [1:2410] 0 0 0 0 0 0 0 0 0 0 ...
#> ..$ order : int [1:2411] 335 262 1142 1202 1258 395 1393 663 599 595 ...
#> ..$ labels : chr [1:2411] "1" "2" "3" "4" ...
#> ..$ method : chr "ward.D"
#> ..$ call : language stats::hclust(d = dist_mat, method = "ward.D")
#> ..$ dist.method: chr "ruzicka"
#> ..- attr(*, "class")= chr "hclust"
head(profiles$profile_matrix)
#>
#> Lang Mnot Zjuj adu car cca fve gma hlu jcu lja lus mdo mes mtr pbr pmu
#> 335 1 1 0 1 0 1 0 1 0 0 0 0 1 0 1 0 1
#> 262 1 1 0 1 0 1 0 1 0 0 1 0 1 0 1 1 1
#> 1142 0 1 0 1 0 1 0 1 0 0 1 0 2 1 1 0 1
#> 1202 1 1 1 1 1 1 1 2 0 1 0 2 1 0 0 0 0
#> 1258 0 1 1 1 1 1 0 2 0 0 1 2 1 0 0 0 1
#> 395 0 1 1 1 0 1 1 1 0 0 0 1 1 1 1 1 1
#>
#> ppe ptr pvu rco roc tpr van vra
#> 335 1 1 0 0 1 1 0 0
#> 262 1 0 1 0 1 1 1 0
#> 1142 1 0 1 0 1 1 1 0
#> 1202 1 1 1 1 0 1 0 0
#> 1258 1 1 1 0 0 1 0 0
#> 395 1 1 0 0 1 1 0 0
As a plot is worth a thousand words (or numbers), you can visualize the
phylogenomic profiles as a heatmap with plot_profiles()
. Let’s also add some
species annotation.
# Create a vector of custom species order to plot
species_order <- c(
"vra", "van", "pvu", "gma", "cca", "tpr", "mtr", "adu", "lja",
"Lang", "car", "pmu", "ppe", "pbr", "mdo", "roc", "fve",
"Mnot", "Zjuj", "jcu", "mes", "rco", "lus", "ptr"
)
# Create a data frame of families for each species
species_annotation <- data.frame(
Species = species_order,
Family = c(rep("Fabaceae", 11), rep("Rosaceae", 6),
"Moraceae", "Ramnaceae", rep("Euphorbiaceae", 3),
"Linaceae", "Salicaceae")
)
head(species_annotation)
#> Species Family
#> 1 vra Fabaceae
#> 2 van Fabaceae
#> 3 pvu Fabaceae
#> 4 gma Fabaceae
#> 5 cca Fabaceae
#> 6 tpr Fabaceae
# Plot phylogenomic profiles
plot_profiles(
profiles, species_annotation,
cluster_species = species_order,
cluster_columns = TRUE
)
The heatmap is a nice way to observe patterns. For instance, you can see some Rosaceae-specific clusters, Fabaceae-specific clusters, and highly conserved ones as well.
If you want to explore in more details the group-specific clusters,
you can use the function find_GS_clusters()
. For that, you only need to
input the profiles matrix and a data frame of species annotation (i.e.,
species groups).
# Find group-specific clusters
gs_clusters <- find_GS_clusters(profiles$profile_matrix, species_annotation)
#> Could not find annotation for species:
#> hlu
head(gs_clusters)
#> Group Percentage Cluster
#> 2 Fabaceae 90.91 1156
#> 21 Fabaceae 81.82 1170
#> 5 Ramnaceae 100.00 1279
#> 22 Fabaceae 90.91 1305
#> 23 Fabaceae 81.82 1309
#> 24 Fabaceae 90.91 1310
# How many family-specific clusters are there?
nrow(gs_clusters)
#> [1] 394
As you can see, there are 394 family-specific clusters in the network. Let’s plot a heatmap of group-specific clusters only.
# Filter profiles matrix to only include group-specific clusters
pgs <- profiles
idx <- rownames(pgs$profile_matrix) %in% gs_clusters$Cluster
pgs$profile_matrix <- pgs$profile_matrix[idx, ]
# Plot heatmap
plot_profiles(
pgs, species_annotation,
cluster_species = species_order,
cluster_columns = TRUE
)
Pretty cool, huh? You can also visualize clusters as a network plot with
the function plot_network()
. For example, let’s visualize the
group-specific clusters.
# 1) Visualize a network of first 5 GS-clusters
id <- gs_clusters$Cluster[1:5]
plot_network(network, clusters, cluster_id = id)
# 2) Coloring nodes by family
genes <- unique(c(network$node1, network$node2))
gene_df <- data.frame(
Gene = genes,
Species = unlist(lapply(strsplit(genes, "_"), head, 1))
)
gene_df <- merge(gene_df, species_annotation)[, c("Gene", "Family")]
head(gene_df)
#> Gene Family
#> 1 Lang_109343937 Fabaceae
#> 2 Lang_109342839 Fabaceae
#> 3 Lang_109356231 Fabaceae
#> 4 Lang_109349826 Fabaceae
#> 5 Lang_109342812 Fabaceae
#> 6 Lang_109347788 Fabaceae
plot_network(network, clusters, cluster_id = id, color_by = gene_df)
# 3) Interactive visualization (zoom out and in to explore it)
plot_network(network, clusters, cluster_id = id,
interactive = TRUE, dim_interactive = c(500, 300))
Finally, you can use the information on presence/absence of clusters in each species to reconstruct a microsynteny-based phylogeny.
To do that, you first need to binarize the profiles matrix (0s and 1s
representing absence and presence, respectively) and transpose it. This can
be done with binarize_and_tranpose()
.
bt_mat <- binarize_and_transpose(profiles$profile_matrix)
# Looking at the first 5 rows and 5 columns of the matrix
bt_mat[1:5, 1:5]
#>
#> 335 262 1142 1202 1258
#> Lang 1 1 0 1 0
#> Mnot 1 1 1 1 1
#> Zjuj 0 0 0 1 1
#> adu 1 1 1 1 1
#> car 0 0 0 1 1
Now, you can use this transposed binary matrix as input to IQTREE (Minh et al. 2020) to infer a phylogeny. Here, once again, we give you 2 options:
Using profiles2phylip()
to write the transposed binary matrix to a
PHYLIP file and running IQTREE from the command line. As you are inferring
a phylogeny from the profiles matrix, you need to specify -st MORPH
to
indicate you have a binary alignment.
Using the helper function infer_microsynteny_phylogeny()
,
which allows you to run IQTREE from an R session. You need to have
IQTREE installed in your machine and in your PATH to run this function.
You can check if you have IQTREE installed with iqtree_is_installed()
.
For the sake of demonstration, we will infer a phylogeny with
infer_microsynteny_phylogeny()
using the profiles for BUSCO genes for
six legume species only. We will also remove non-variable sites (i.e.,
clusters that are present in all species or absent in all species).
However, we’re only using this filtered data set for speed issues.
For real-life applications, the best thing to do is to
include phylogenomic profiles for all genes (not only BUSCO genes).
# Leave only 6 legume species and P. mume as an outgroup for testing purposes
included <- c("gma", "pvu", "vra", "van", "cca", "pmu")
bt_mat <- bt_mat[rownames(bt_mat) %in% included, ]
# Remove non-variable sites
bt_mat <- bt_mat[, colSums(bt_mat) != length(included)]
if(iqtree_is_installed()) {
phylo <- infer_microsynteny_phylogeny(bt_mat, outgroup = "pmu",
threads = 1)
}
The output of infer_microsynteny_phylogeny()
is a character vector with paths
to the output files from IQTREE. Usually, you are interested in the file
ending in .treefile. This is the species tree in Newick format, and it can
be visualized with your favorite tree viewer. We strongly recommend using
the read.tree()
function from the Bioconductor package treeio (Wang et al. 2020) to
read the tree, and visualizing it with the ggtree
Bioc package (Yu et al. 2017).
For example, let’s visualize a microsynteny-based phylogeny for 123 angiosperm
species, whose phylogenomic profiles were obtained from Zhao et al. (2021).
data(angiosperm_phylogeny)
# Plotting the tree
library(ggtree)
ggtree(angiosperm_phylogeny) +
geom_tiplab(size = 3) +
xlim(0, 0.3)
This document was created under the following conditions:
sessionInfo()
#> R version 4.2.1 (2022-06-23)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.5 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_GB LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] ggtree_3.6.0 syntenet_1.0.0 BiocStyle_2.26.0
#>
#> loaded via a namespace (and not attached):
#> [1] nlme_3.1-160 bitops_1.0-7
#> [3] matrixStats_0.62.0 labdsv_2.0-1
#> [5] RColorBrewer_1.1-3 GenomeInfoDb_1.34.0
#> [7] tools_4.2.1 bslib_0.4.0
#> [9] utf8_1.2.2 R6_2.5.1
#> [11] lazyeval_0.2.2 DBI_1.1.3
#> [13] BiocGenerics_0.44.0 mgcv_1.8-41
#> [15] colorspace_2.0-3 tidyselect_1.2.0
#> [17] compiler_4.2.1 cli_3.4.1
#> [19] Biobase_2.58.0 intergraph_2.0-2
#> [21] network_1.18.0 DelayedArray_0.24.0
#> [23] labeling_0.4.2 rtracklayer_1.58.0
#> [25] bookdown_0.29 sass_0.4.2
#> [27] scales_1.2.1 yulab.utils_0.0.5
#> [29] stringr_1.4.1 digest_0.6.30
#> [31] Rsamtools_2.14.0 rmarkdown_2.17
#> [33] XVector_0.38.0 pkgconfig_2.0.3
#> [35] htmltools_0.5.3 MatrixGenerics_1.10.0
#> [37] highr_0.9 fastmap_1.1.0
#> [39] htmlwidgets_1.5.4 rlang_1.0.6
#> [41] gridGraphics_0.5-1 farver_2.1.1
#> [43] jquerylib_0.1.4 BiocIO_1.8.0
#> [45] generics_0.1.3 jsonlite_1.8.3
#> [47] statnet.common_4.7.0 BiocParallel_1.32.0
#> [49] dplyr_1.0.10 RCurl_1.98-1.9
#> [51] magrittr_2.0.3 ggplotify_0.1.0
#> [53] ggnetwork_0.5.10 GenomeInfoDbData_1.2.9
#> [55] patchwork_1.1.2 Matrix_1.5-1
#> [57] Rcpp_1.0.9 munsell_0.5.0
#> [59] S4Vectors_0.36.0 fansi_1.0.3
#> [61] ape_5.6-2 lifecycle_1.0.3
#> [63] stringi_1.7.8 yaml_2.3.6
#> [65] MASS_7.3-58.1 SummarizedExperiment_1.28.0
#> [67] zlibbioc_1.44.0 Rtsne_0.16
#> [69] grid_4.2.1 parallel_4.2.1
#> [71] crayon_1.5.2 lattice_0.20-45
#> [73] Biostrings_2.66.0 splines_4.2.1
#> [75] knitr_1.40 pillar_1.8.1
#> [77] igraph_1.3.5 GenomicRanges_1.50.0
#> [79] rjson_0.2.21 codetools_0.2-18
#> [81] stats4_4.2.1 XML_3.99-0.12
#> [83] glue_1.6.2 evaluate_0.17
#> [85] ggfun_0.0.7 BiocManager_1.30.19
#> [87] treeio_1.22.0 vctrs_0.5.0
#> [89] networkD3_0.4 purrr_0.3.5
#> [91] tidyr_1.2.1 gtable_0.3.1
#> [93] assertthat_0.2.1 cachem_1.0.6
#> [95] ggplot2_3.3.6 xfun_0.34
#> [97] restfulr_0.0.15 tidytree_0.4.1
#> [99] coda_0.19-4 tibble_3.1.8
#> [101] pheatmap_1.0.12 aplot_0.1.8
#> [103] GenomicAlignments_1.34.0 IRanges_2.32.0
#> [105] cluster_2.1.4
Altschul, Stephen F, Thomas L Madden, Alejandro A Schäffer, Jinghui Zhang, Zheng Zhang, Webb Miller, and David J Lipman. 1997. “Gapped Blast and Psi-Blast: A New Generation of Protein Database Search Programs.” Nucleic Acids Research 25 (17): 3389–3402.
Buchfink, Benjamin, Klaus Reuter, and Hajk-Georg Drost. 2021. “Sensitive Protein Alignments at Tree-of-Life Scale Using Diamond.” Nature Methods 18 (4): 366–68.
Eddelbuettel, Dirk, and Romain François. 2011. “Rcpp: Seamless R and C++ Integration.” Journal of Statistical Software 40: 1–18.
Lange, Jens. 2021. “CliquePercolation: An R Package for Conducting and Visualizing Results of the Clique Percolation Network Community Detection Algorithm.” Journal of Open Source Software 6 (62): 3210.
Minh, Bui Quang, Heiko A Schmidt, Olga Chernomor, Dominik Schrempf, Michael D Woodhams, Arndt Von Haeseler, and Robert Lanfear. 2020. “IQ-Tree 2: New Models and Efficient Methods for Phylogenetic Inference in the Genomic Era.” Molecular Biology and Evolution 37 (5): 1530–4.
Vandepoele, Klaas, Michiel Van Bel, Guilhem Richard, Sofie Van Landeghem, Bram Verhelst, Hervé Moreau, Yves Van de Peer, Nigel Grimsley, and Gwenael Piganeau. 2013. “Pico-Plaza, a Genome Database of Microbial Photosynthetic Eukaryotes.” Environmental Microbiology 15 (8): 2147–53.
Wang, Li-Gen, Tommy Tsan-Yuk Lam, Shuangbin Xu, Zehan Dai, Lang Zhou, Tingze Feng, Pingfan Guo, et al. 2020. “Treeio: An R Package for Phylogenetic Tree Input and Output with Richly Annotated and Associated Data.” Molecular Biology and Evolution 37 (2): 599–603.
Wang, Yupeng, Haibao Tang, Jeremy D DeBarry, Xu Tan, Jingping Li, Xiyin Wang, Tae-ho Lee, et al. 2012. “MCScanX: A Toolkit for Detection and Evolutionary Analysis of Gene Synteny and Collinearity.” Nucleic Acids Research 40 (7): e49–e49.
Yu, Guangchuang, David K Smith, Huachen Zhu, Yi Guan, and Tommy Tsan-Yuk Lam. 2017. “Ggtree: An R Package for Visualization and Annotation of Phylogenetic Trees with Their Covariates and Other Associated Data.” Methods in Ecology and Evolution 8 (1): 28–36.
Zhao, Tao, Rens Holmer, Suzanne de Bruijn, Gerco C Angenent, Harrold A van den Burg, and M Eric Schranz. 2017. “Phylogenomic Synteny Network Analysis of Mads-Box Transcription Factor Genes Reveals Lineage-Specific Transpositions, Ancient Tandem Duplications, and Deep Positional Conservation.” The Plant Cell 29 (6): 1278–92.
Zhao, Tao, and M Eric Schranz. 2017. “Network Approaches for Plant Phylogenomic Synteny Analysis.” Current Opinion in Plant Biology 36: 129–34.
———. 2019. “Network-Based Microsynteny Analysis Identifies Major Differences and Genomic Outliers in Mammalian and Angiosperm Genomes.” Proceedings of the National Academy of Sciences 116 (6): 2165–74.
Zhao, Tao, Arthur Zwaenepoel, Jia-Yu Xue, Shu-Min Kao, Zhen Li, M Eric Schranz, and Yves Van de Peer. 2021. “Whole-Genome Microsynteny-Based Phylogeny of Angiosperms.” Nature Communications 12 (1): 1–14.