Contents

1 Introduction

SynExtend is a package of tools for working with objects of class Synteny built from the package DECIPHER’s FindSynteny() function.

Synteny maps provide a powerful tool for quantifying and visualizing where pairs of genomes share order. Typically these maps are built from predictions of orthologous pairs, where groups of pairs that provide contiguous and sequential blocks in their respective genomes are deemed a ‘syntenic block’. That designation of synteny can then used to further interrogate the predicted orthologs themselves, or query topics like genomic rearrangements or ancestor genome reconstruction.

FindSynteny takes a different approach, finding exactly matched shared k-mers and determining where shared k-mers, or blocks of proximate shared k-mers are significant. Combining the information generated by FindSynteny with locations of genomic features allows us to simply mark where features are linked by syntenic k-mers. These linked features represent potential orthologous pairs, and can be easily evaluated on the basis of the k-mers that they share, or alignment.

2 Package Structure

Currently SynExtend contains one set of functions for performig orthology predictions, as well as a rearrangement estimation function that is currently under construction.

2.1 Installation

  1. Install the latest version of R using CRAN.
  2. Install SynExtend in R by running the following commands:
if (!requireNamespace("BiocManager",
                      quietly = TRUE)) {
  install.packages("BiocManager")
}
BiocManager::install("SynExtend")

3 Usage

Using the FindSynteny function in DECIPHER builds an object of class Synteny. In this tutorial, a prebuilt DECIPHER database is used. For database construction see ?Seqs2DB in DECIPHER. This example starts with a database containing three endosymbiont genomes that were chosen to keep package data to a minimum.

library(SynExtend)
## Loading required package: DECIPHER
## Loading required package: Biostrings
## Loading required package: BiocGenerics
## Loading required package: generics
## 
## Attaching package: 'generics'
## The following objects are masked from 'package:base':
## 
##     as.difftime, as.factor, as.ordered, intersect, is.element, setdiff,
##     setequal, union
## 
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:stats':
## 
##     IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
## 
##     Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
##     as.data.frame, basename, cbind, colnames, dirname, do.call,
##     duplicated, eval, evalq, get, grep, grepl, is.unsorted, lapply,
##     mapply, match, mget, order, paste, pmax, pmax.int, pmin, pmin.int,
##     rank, rbind, rownames, sapply, saveRDS, table, tapply, unique,
##     unsplit, which.max, which.min
## Loading required package: S4Vectors
## Loading required package: stats4
## 
## Attaching package: 'S4Vectors'
## The following object is masked from 'package:utils':
## 
##     findMatches
## The following objects are masked from 'package:base':
## 
##     I, expand.grid, unname
## Loading required package: IRanges
## Loading required package: XVector
## Loading required package: Seqinfo
## 
## Attaching package: 'Biostrings'
## The following object is masked from 'package:base':
## 
##     strsplit
## 
## Attaching package: 'SynExtend'
## The following object is masked from 'package:stats':
## 
##     dendrapply
DBPATH <- system.file("extdata",
                      "Endosymbionts_v05a.sqlite",
                      package = "SynExtend")

Syn <- FindSynteny(dbFile = DBPATH)
## ================================================================================
## 
## Time difference of 25.67 secs

Synteny maps represent where genomes share order. Simply printing a synteny object to the console displays a gross level view of the data inside. Objects of class Synteny can also be plotted to provide clear visual representations of the data inside. The genomes used in this example are distantly related and fairly dissimilar.

Syn
##           1         2         3         4
## 1     1 seq  53% hits 8.9% hits  94% hits
## 2   1 block     1 seq 7.1% hits  95% hits
## 3 96 blocks 79 blocks     1 seq 9.3% hits
## 4  2 blocks  2 blocks 99 blocks     1 seq
pairs(Syn)

Data present inside objects of class Synteny can also be accessed relatively easily. The object itself is functionally a matrix of lists, with data describing exactly matched k-mers present in the upper triangle, and data describing blocks of chained k-mers in the lower triangle. For more information see ?FindSynteny in the package DECIPHER.

print(head(Syn[[1, 2]]))
##      index1 index2 strand width start1 start2 frame1 frame2
## [1,]      1      1      1   271 372826 791111      0      0
## [2,]      1      1      1   363 373097 790839      0      0
## [3,]      1      1      1  1283 373461 790476      3      3
## [4,]      1      1      1   197 374745 789193      0      0
## [5,]      1      1      1     4 374943 788996      3      1
## [6,]      1      1      1    43 374947 788992      1      2
print(head(Syn[[2, 1]]))
##      index1 index2 strand  score start1 start2   end1   end2 first_hit last_hit
## [1,]      1      1      1 475595 372826 373207 790687 791111         1      327

The above printed objects show the data for the comparison between the first and second genome in our database.

To take advantage of these synteny maps, we can then overlay the gene calls for each genome present on top of our map.

Next, GFF annotations for the associated genomes are parsed to provide gene calls in a use-able format. GFFs are not the only possible source of appropriate gene calls, but they are the source that was used during package construction and testing. Parsed GFFs can be constructed with gffToDataFrame, for full functionality, or GFFs can be imported via rtracklater::import() for limited functionality. GeneCalls for both the PairSummaries and NucleotideOverlap functions must be named list, and those names must match dimnames(Syn)[[1]].

# generating genecalls with local data:
GC <- gffToDataFrame(GFF = system.file("extdata",
                                       "GCF_023585825.1_ASM2358582v1_genomic.gff.gz",
                                       package = "SynExtend"),
                     Verbose = TRUE)
## ================================================================================
## Time difference of 7.460854 secs
# in an effort to be space conscious, not all original gffs are kept within this package
GeneCalls <- get(data("Endosymbionts_GeneCalls", package = "SynExtend"))

SynExtend’s gffToDataFrame function will directly import gff files into a usable format, and includes other extracted information.

print(head(GeneCalls[[1]]))
## DataFrame with 6 rows and 11 columns
##       Index    Strand     Start      Stop        Type                 ID
##   <integer> <integer> <integer> <integer> <character>        <character>
## 1         1         0       170       493        gene gene-M9394_RS00005
## 2         1         0       574      1935        gene gene-M9394_RS00010
## 3         1         0      1995      2522        gene gene-M9394_RS00015
## 4         1         0      2581      3321        gene gene-M9394_RS00020
## 5         1         1      4383      5396        gene gene-M9394_RS00025
## 6         1         0      5682      6377        gene gene-M9394_RS00030
##           Range                Product    Coding Translation_Table
##   <IRangesList>            <character> <logical>       <character>
## 1       170-493 cell envelope integr..      TRUE                11
## 2      574-1935 Tol-Pal system beta ..      TRUE                11
## 3     1995-2522 peptidoglycan-associ..      TRUE                11
## 4     2581-3321 tol-pal system prote..      TRUE                11
## 5     4383-5396 6-phosphogluconolact..      TRUE                11
## 6     5682-6377 2%2C3-diphosphoglyce..      TRUE                11
##          Contig
##     <character>
## 1 NZ_CP097751.1
## 2 NZ_CP097751.1
## 3 NZ_CP097751.1
## 4 NZ_CP097751.1
## 5 NZ_CP097751.1
## 6 NZ_CP097751.1

Raw GFF imports are also acceptable, but prevent alignments in amino acid space with PairSummaries().

X01 <- rtracklayer::import(system.file("extdata",
                                       "GCF_023585825.1_ASM2358582v1_genomic.gff.gz",
                                       package = "SynExtend"))
class(X01)
print(X01)

SynExtend’s primary functions provide a way to identify where pairs of genes are explicitly linked by syntenic hits, and then summarize those links. The first step is just identifying those links.

Links <- NucleotideOverlap(SyntenyObject = Syn,
                           GeneCalls = GeneCalls,
                           Verbose = TRUE)
## 
## Reconciling genecalls.
## ================================================================================
## Finding connected features.
## ================================================================================
## Time difference of 0.479465 secs

The Links object generated by NucleotideOverlap is a raw representation of positions on the synteny map where shared k-mers link genes between paired genomes. As such, it is analagous in shape to objects of class Synteny. This raw object is unlikely to be useful to most users, but has been left exposed to ensure that this data remains accessible should a user desire to have access to it.

class(Links)
## [1] "LinkedPairs"
print(Links)
##            1          2         3         4
## 1   1 Contig  378 Pairs 111 Pairs 711 Pairs
## 2  415 Kmers   1 Contig  81 Pairs 711 Pairs
## 3  268 Kmers  233 Kmers  1 Contig 106 Pairs
## 4 1728 Kmers 2000 Kmers 292 Kmers  1 Contig

This raw data can be processed to provide a straightforward summary of predicted pairs.

LinkedPairs1 <- PairSummaries(SyntenyLinks = Links,
                              DBPATH = DBPATH,
                              PIDs = FALSE,
                              Verbose = TRUE)
## 
## Preparing overhead data.
## Overhead complete.
## Collecting pairs.
## ================================================================================
## Time difference of 2.965107 secs

The object LinkedPairs1 is a data.frame where each row is populated by information about a predicted orthologous pair. By default PairSummaries uses a simple model to determine whether the k-mers that link a pair of genes are likely to provide an erroneous link. When set to Model = "Global", is is simply a prediction of whether the involved nucleotides are likely to describe a pair of genomic features whose alignment would result in a PID that falls within a random distribution. This model is effective if somewhat permissive, but is significantly faster than performing many pairwise alignments.

print(head(LinkedPairs1))
##        p1      p2 ExactMatch Consensus TotalKmers MaxKmer p1FeatureLength
## 1 1_1_312 2_1_661        261 0.9945852          1     261             849
## 2 1_1_313 2_1_660        927 1.0000000          1     927             927
## 3 1_1_314 2_1_659        588 1.0000000          1     588             588
## 4 1_1_315 2_1_658       1656 1.0000000          2    1414            1656
## 5 1_1_316 2_1_657        885 1.0000000          1     885             885
## 6 1_1_317 2_1_656         74 1.0000000          1      74              74
##   p2FeatureLength Adjacent     TetDist PIDType PredictedPID
## 1             880        1 0.011881038      AA    0.9226646
## 2             927        2 0.005927733      AA    0.9998543
## 3             588        2 0.004834918      AA    0.9988983
## 4            1656        2 0.002963690      AA    0.9999549
## 5             885        2 0.004535147      AA    0.9998307
## 6              74        2 0.000000000      NT    0.8953289

PairSummaries includes arguments that allow for aligning all pairs that are predicted, via PIDs = TRUE, while IgnoreDefaultStringSet = FALSE indicates that alignments should be performed in nucleotide or amino acid space as is appropriate for the linked sequences. Setting IgnoreDefaultStringSet = TRUE will force all alignments into nucleotide space.

As of SynExtend v 1.3.13, the functions ExtractBy and DisjointSet have been added to provide users with direct tools to work with PairSummaries objects.

SingleLinkageClusters <- DisjointSet(Pairs = LinkedPairs1,
                                     Verbose = TRUE)
## 
## Assigning initial root:
## ================================================================================
## 
## Time difference of 0.01638627 secs
## 
## Assigning final root:
## ================================================================================
## 
## Time difference of 0.01547384 secs
## 
## Assigning single linkage clusters.
## Assignments complete.
## 
## Time difference of 0.03682733 secs
# extract the first 10 clusters
Sets <- ExtractBy(x = LinkedPairs1,
                  y = DBPATH,
                  z = SingleLinkageClusters[1:10],
                  Verbose = TRUE)
## 
## Extracting Sequences:
## ================================================================================
## 
## Arranging Sequences:
## ================================================================================
## 
## Time difference of 0.3287015 secs
head(Sets)
## [[1]]
## DNAStringSet object of length 3:
##     width seq                                               names               
## [1]   849 ATGCCTAAATATATACAAAAAAA...TGAAACTAATTATACAAAATTAA 1_1_312
## [2]   880 ATGCCTAAATATATACAAAAAAA...TGAAACTAATTATACAAAATTAA 2_1_661
## [3]   849 ATGCCTAAATATATACAAAAAAA...TGAAACTAATTATACAAAATTAA 4_1_656
## 
## [[2]]
## DNAStringSet object of length 4:
##     width seq                                               names               
## [1]   927 ATGAATCGTGCATATAATTTATT...TTGAAGCTATTTTAAAAGAGTAA 1_1_313
## [2]   927 ATGAATCGTGCATATAATTTATT...TTGAAGCTATTTTAAAAGAGTAA 2_1_660
## [3]   927 ATGAATCGTGCATATAATTTATT...TTGAAGCTATTTTAAAAGAGTAA 4_1_657
## [4]  1503 ATGTTTAATTGGTTTAACCGTAA...TAGAGGCTTTACTTAATGATTGA 3_1_288
## 
## [[3]]
## DNAStringSet object of length 3:
##     width seq                                               names               
## [1]   588 ATGAAAAATAAATTACATTTAAG...ATCTCTTTTTAAGAAGATCATAG 1_1_314
## [2]   588 ATGAAAAATAAATTACATTTAAG...ATCTCTTTTTAAGAAGATCATAG 2_1_659
## [3]   588 GTGAAAAATAAATTACATTTAGG...ATCTTTTTTTAAGAAGATCATAG 4_1_658
## 
## [[4]]
## DNAStringSet object of length 3:
##     width seq                                               names               
## [1]  1656 ATGAAAAATATAAATCCTAGTCG...GTTATAAATCTTGGTGTGATTGA 1_1_315
## [2]  1656 ATGAAAAATATAAATCCTAGTCG...GTTATAAATCTTGGTGTGATTGA 2_1_658
## [3]  1656 ATGAAAAATATAAATCCTAATCG...GTTACAAATCCTGGTGTGGTTAA 4_1_659
## 
## [[5]]
## DNAStringSet object of length 3:
##     width seq                                               names               
## [1]   885 ATGCCAGTTCGGATATTAAATAA...GGCTGAATCATTATGTACAATAA 1_1_316
## [2]   885 ATGCCAGTTCGGATATTAAATAA...GGCTGAATCATTATGTACAATAA 2_1_657
## [3]   915 ATGCCAGTTCAGGTATTAAACAG...TAAAACATTTAAGTTCTTTCTGA 4_1_660
## 
## [[6]]
## DNAStringSet object of length 3:
##     width seq                                               names               
## [1]    74 AGGCTTGTAGCTCAGGCGGTTAG...TGGTTCAATTCCACTCAGGCCTA 1_1_317
## [2]    74 AGGCTTGTAGCTCAGGCGGTTAG...TGGTTCAATTCCACTCAGGCCTA 2_1_656
## [3]    74 AGGCTTGTAGCTCAGGCGGTTAG...TGGTTCAATTCCACTCAGGCCTA 4_1_661

Session Info:

sessionInfo()
## R version 4.6.0 RC (2026-04-17 r89917)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
## 
## 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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] SynExtend_1.24.0    DECIPHER_3.8.0      Biostrings_2.80.0  
##  [4] Seqinfo_1.2.0       XVector_0.52.0      IRanges_2.46.0     
##  [7] S4Vectors_0.50.0    BiocGenerics_0.58.0 generics_0.1.4     
## [10] BiocStyle_2.40.0   
## 
## loaded via a namespace (and not attached):
##  [1] bit_4.6.0           jsonlite_2.0.0      compiler_4.6.0     
##  [4] BiocManager_1.30.27 crayon_1.5.3        Rcpp_1.1.1-1.1     
##  [7] tinytex_0.59        blob_1.3.0          magick_2.9.1       
## [10] parallel_4.6.0      jquerylib_0.1.4     yaml_2.3.12        
## [13] fastmap_1.2.0       R6_2.6.1            knitr_1.51         
## [16] bookdown_0.46       DBI_1.3.0           bslib_0.10.0       
## [19] rlang_1.2.0         cachem_1.1.0        xfun_0.57          
## [22] sass_0.4.10         bit64_4.8.0         otel_0.2.0         
## [25] RSQLite_2.4.6       memoise_2.0.1       cli_3.6.6          
## [28] magrittr_2.0.5      digest_0.6.39       lifecycle_1.0.5    
## [31] vctrs_0.7.3         evaluate_1.0.5      rmarkdown_2.31     
## [34] tools_4.6.0         pkgconfig_2.0.3     htmltools_0.5.9