Contents

library(comapr)
library(GenomicRanges)
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> 
#> 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, append,
#>     as.data.frame, basename, cbind, colnames, dirname, do.call,
#>     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
#>     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
#>     pmin.int, rank, rbind, rownames, sapply, setdiff, sort, table,
#>     tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: S4Vectors
#> 
#> Attaching package: 'S4Vectors'
#> The following objects are masked from 'package:base':
#> 
#>     I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
library(BiocParallel)

1 Installation

Install via BiocManager as follow:

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

BiocManager::install("comapr")

2 Introduction

comapr lets you interrogate the genotyping results for a sequence of markers across the chromosome and detect meiotic crossovers from genotype shifts. In this document, we demonstrate how genetic distances are calculated from genotyping results for a group of samples using functions available from comapr.

3 The demo data set of genotyping results

The package includes a small data object that contains 70 markers and their genotyping results for 22 samples. The 22 samples are the BC1F1 progenies generated by

Therefore, the genotype shifts detected in BC1F1 samples represent crossovers that have happened during meiosis for the F1 parents.

BiocParallel::register(SerialParam())
BiocParallel::bpparam()
#> class: SerialParam
#>   bpisup: FALSE; bpnworkers: 1; bptasks: 0; bpjobname: BPJOB
#>   bplog: FALSE; bpthreshold: INFO; bpstopOnError: TRUE
#>   bpRNGseed: ; bptimeout: NA; bpprogressbar: FALSE
#>   bpexportglobals: FALSE; bpexportvariables: FALSE; bpforceGC: FALSE; bpfallback: FALSE
#>   bplogdir: NA
#>   bpresultdir: NA
data(snp_geno_gr)
data(parents_geno)

3.1 Format the genotype

In order to detect crossovers from the 22 samples’ genotype results (BC1F1 samples), we need to format the result into Homo_ref, Homo_alt and Het encodings. This can be simply done through comparing with the parents’ genotypes.

Here, the genotype noted as “Fail” will be converted to NA. Please also note that we supplied “Homo_ref” as one of the fail options because “Homo_ref” is not in the possible genotypes of markers in BC1F1 samples.

corrected_geno <- correctGT(gt_matrix = mcols(snp_geno_gr),
                            ref = parents_geno$ref,
                            alt = parents_geno$alt,
                            fail = "Fail",
                            wrong_label = "Homo_ref")

mcols(snp_geno_gr) <- corrected_geno

The corrected_geno matrix

head(mcols(snp_geno_gr)[,1:5])
#> DataFrame with 6 rows and 5 columns
#>           X92         X93         X94         X95         X96
#>   <character> <character> <character> <character> <character>
#> 1    Homo_alt    Homo_alt    Homo_alt    Homo_alt    Homo_alt
#> 2    Homo_alt    Homo_alt    Homo_alt    Homo_alt    Homo_alt
#> 3    Homo_alt    Homo_alt    Homo_alt    Homo_alt    Homo_alt
#> 4    Homo_alt    Homo_alt    Homo_alt    Homo_alt    Homo_alt
#> 5          NA    Homo_alt    Homo_alt    Homo_alt    Homo_alt
#> 6         Het    Homo_alt    Homo_alt    Homo_alt    Homo_alt

Note that there are missing values in this resulting matrix that can be resulted from:

  • No genotype data from the genotyping result
  • The genotype is wrong and it is converted to NA

3.2 Find outlier samples/markers

In this step, we try to identify markers that have NA genotype across many samples or samples that have a lot markers failed for removal. We use the countGT function for find bad markers/samples.

genotype_counts <- countGT(mcols(snp_geno_gr))

genotype_counts$plot

The number of markes and samples are saved in a list returned by countGT

genotype_counts$n_markers
#>  [1] 30 33 33 33 33 33 31 31 32 33 33 33 32 33 33 33 32 32 33 32 33 33
genotype_counts$n_samples
#>  [1] 22 22 22 22 18 22 22 22 22 22 22 21 22 22 22 22 22 22 22 22 22 22 22 22 22
#> [26] 22 22 22 16 22 22 21 22

We now filter out markers/samples by using function filterGT. min_markers specifies at least how many markers a sample needs to be kept. Likewise for min_samples.

A printed message contains information about how much markers or samples have been filtered.

corrected_geno <- filterGT(snp_geno_gr,
                           min_markers = 30,
                           min_samples = 2)
#> filter out 0 marker(s)
#> filter out 0 sample(s)

3.3 Find sample duplicates

Sample duplicates are identified by finding samples that share exactly same genotypes across all available markers. findDupSamples can be applied and a threshold value is provided and used as a cut-off on the percentage of same genotype markers the duplicated samples should share.

dups <- findDupSamples(mcols(corrected_geno),
                       threshold = 0.99)
dups
#>      X98  
#> [1,] "X98"
#> [2,] "X99"

Now we remove the duplicated samples.

mcols(corrected_geno) <- mcols(corrected_geno)[,colnames(mcols(corrected_geno))!="X98"]
#corrected_geno

3.4 Count crossovers

Crossovers are detected and counted through examining the patterns of genotypes along the chromosome. When there is a shift from one genotype block to another, a crossover is observed. This is done through calling countCOs function which returns a GRange object with crossover counts for the list of marker intervals.

The crossover count values in the columns can be non-integer when one observed crossover can not be determined to be completely distributed to the marker interval in the corresponding row. The observed crossover is then distributed to the adjacent intervals proportionally to their interval base pair sizes.

marker_gr_cos <- countCOs(corrected_geno)

marker_gr_cos[1:5,1:5]
#> GRanges object with 5 ranges and 5 metadata columns:
#>       seqnames            ranges strand |       X92       X93       X94
#>          <Rle>         <IRanges>  <Rle> | <numeric> <numeric> <numeric>
#>   [1]        1  6655965-21638463      * | 0.9358742         0         0
#>   [2]        1 21638465-22665059      * | 0.0641258         0         0
#>   [3]        1 22665061-34590735      * | 0.0000000         0         0
#>   [4]        1 35033642-38996025      * | 0.0000000         0         0
#>   [5]        2   4248665-5348752      * | 0.0000000         0         0
#>             X95       X96
#>       <numeric> <numeric>
#>   [1]         0         0
#>   [2]         0         0
#>   [3]         0         0
#>   [4]         0         0
#>   [5]         0         0
#>   -------
#>   seqinfo: 4 sequences from an unspecified genome; no seqlengths

3.5 Calculate genetic distance

The genetic distances of marker intervals are calcuated based on the crossover rates via applying mapping function, Kosambi or Haldane by calling the calGeneticDist function. The returned genetic distances are in unit of centiMorgan.

dist_gr <- calGeneticDist(marker_gr_cos,
                          mapping_fun = "k")

dist_gr[1:5,]
#> GRanges object with 5 ranges and 1 metadata column:
#>       seqnames            ranges strand | kosambi_cm
#>          <Rle>         <IRanges>  <Rle> |  <numeric>
#>   [1]        1  6655965-21638463      * |   14.36279
#>   [2]        1 21638465-22665059      * |    5.08472
#>   [3]        1 22665061-34590735      * |    9.64156
#>   [4]        1 35033642-38996025      * |    9.64156
#>   [5]        2   4248665-5348752      * |    4.77638
#>   -------
#>   seqinfo: 4 sequences from an unspecified genome; no seqlengths

3.6 Calculate genetic distance for equally binned intervals

Alternatively, instead of returning the genetic distances in supplied marker intervals, we can specify a bin_size which tells calGeneticDist to return calculated genetic distances for equally sized chromosome bins.

dist_bin_gr <- calGeneticDist(marker_gr_cos,bin_size = 1e6)

dist_bin_gr[1:5,]
#> GRanges object with 5 ranges and 1 metadata column:
#>       seqnames          ranges strand | kosambi_cm
#>          <Rle>       <IRanges>  <Rle> |  <numeric>
#>   [1]        1        1-998753      * |          0
#>   [2]        1  998754-1997506      * |          0
#>   [3]        1 1997507-2996258      * |          0
#>   [4]        1 2996259-3995011      * |          0
#>   [5]        1 3995012-4993763      * |          0
#>   -------
#>   seqinfo: 4 sequences from mm10 genome

3.7 Total genetic distances

With genetic distances calculated, we can do a sum of all genetic distances across all marker intervals. We can see that we got the same total genetic distances for marker based intervals and the equally binned intervals.

sum(dist_bin_gr$kosambi_cm)
#> [1] 168.0926
sum(dist_gr$kosambi_cm)
#> [1] 168.0926

3.8 Plotting genetic distance for each bin

comapr also includes functions for visulising genetic distances of marker intervals or binned intervals.

plotGeneticDist(dist_bin_gr,chr = "1")

3.9 Plot cumulative distances

We can also plot the cumulative genetic distances of certain chromosomes

plotGeneticDist(dist_bin_gr,chr=c("1"),cumulative = TRUE)

Multiple chromosomes

plotGeneticDist(dist_bin_gr,cumulative = TRUE)

3.10 Whole Genome Genetic distances Plot

comapr implements a whole genome plot function too, that takes all chromosomes available in the result and plot a cumulative genetic distances by cumulatively summing all intervals across all chromosomes.

plotWholeGenome(dist_bin_gr)

4 Sessioninfo

sessionInfo()
#> R version 4.2.0 RC (2022-04-19 r82224)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.4 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.15-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] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#> [1] BiocParallel_1.30.0  GenomicRanges_1.48.0 GenomeInfoDb_1.32.0 
#> [4] IRanges_2.30.0       S4Vectors_0.34.0     BiocGenerics_0.42.0 
#> [7] comapr_1.0.0         BiocStyle_2.24.0    
#> 
#> loaded via a namespace (and not attached):
#>   [1] colorspace_2.0-3            rjson_0.2.21               
#>   [3] ellipsis_0.3.2              circlize_0.4.14            
#>   [5] biovizBase_1.44.0           htmlTable_2.4.0            
#>   [7] XVector_0.36.0              GlobalOptions_0.1.2        
#>   [9] base64enc_0.1-3             dichromat_2.0-0            
#>  [11] rstudioapi_0.13             farver_2.1.0               
#>  [13] bit64_4.0.5                 AnnotationDbi_1.58.0       
#>  [15] fansi_1.0.3                 xml2_1.3.3                 
#>  [17] codetools_0.2-18            splines_4.2.0              
#>  [19] cachem_1.0.6                knitr_1.38                 
#>  [21] Formula_1.2-4               jsonlite_1.8.0             
#>  [23] Rsamtools_2.12.0            cluster_2.1.3              
#>  [25] dbplyr_2.1.1                png_0.1-7                  
#>  [27] BiocManager_1.30.17         compiler_4.2.0             
#>  [29] httr_1.4.2                  backports_1.4.1            
#>  [31] assertthat_0.2.1            Matrix_1.4-1               
#>  [33] fastmap_1.1.0               lazyeval_0.2.2             
#>  [35] cli_3.3.0                   htmltools_0.5.2            
#>  [37] prettyunits_1.1.1           tools_4.2.0                
#>  [39] gtable_0.3.0                glue_1.6.2                 
#>  [41] GenomeInfoDbData_1.2.8      reshape2_1.4.4             
#>  [43] dplyr_1.0.8                 rappdirs_0.3.3             
#>  [45] Rcpp_1.0.8.3                Biobase_2.56.0             
#>  [47] jquerylib_0.1.4             vctrs_0.4.1                
#>  [49] Biostrings_2.64.0           rtracklayer_1.56.0         
#>  [51] iterators_1.0.14            xfun_0.30                  
#>  [53] stringr_1.4.0               lifecycle_1.0.1            
#>  [55] ensembldb_2.20.0            restfulr_0.0.13            
#>  [57] XML_3.99-0.9                zlibbioc_1.42.0            
#>  [59] scales_1.2.0                BSgenome_1.64.0            
#>  [61] VariantAnnotation_1.42.0    ProtGenerics_1.28.0        
#>  [63] hms_1.1.1                   MatrixGenerics_1.8.0       
#>  [65] parallel_4.2.0              SummarizedExperiment_1.26.0
#>  [67] AnnotationFilter_1.20.0     RColorBrewer_1.1-3         
#>  [69] yaml_2.3.5                  curl_4.3.2                 
#>  [71] memoise_2.0.1               gridExtra_2.3              
#>  [73] ggplot2_3.3.5               sass_0.4.1                 
#>  [75] biomaRt_2.52.0              rpart_4.1.16               
#>  [77] latticeExtra_0.6-29         stringi_1.7.6              
#>  [79] RSQLite_2.2.12              highr_0.9                  
#>  [81] BiocIO_1.6.0                foreach_1.5.2              
#>  [83] checkmate_2.1.0             GenomicFeatures_1.48.0     
#>  [85] filelock_1.0.2              shape_1.4.6                
#>  [87] rlang_1.0.2                 pkgconfig_2.0.3            
#>  [89] matrixStats_0.62.0          bitops_1.0-7               
#>  [91] evaluate_0.15               lattice_0.20-45            
#>  [93] purrr_0.3.4                 labeling_0.4.2             
#>  [95] GenomicAlignments_1.32.0    htmlwidgets_1.5.4          
#>  [97] bit_4.0.4                   tidyselect_1.1.2           
#>  [99] plyr_1.8.7                  magrittr_2.0.3             
#> [101] bookdown_0.26               R6_2.5.1                   
#> [103] magick_2.7.3                generics_0.1.2             
#> [105] Hmisc_4.7-0                 DelayedArray_0.22.0        
#> [107] DBI_1.1.2                   withr_2.5.0                
#> [109] pillar_1.7.0                foreign_0.8-82             
#> [111] survival_3.3-1              KEGGREST_1.36.0            
#> [113] RCurl_1.98-1.6              nnet_7.3-17                
#> [115] tibble_3.1.6                crayon_1.5.1               
#> [117] utf8_1.2.2                  plotly_4.10.0              
#> [119] BiocFileCache_2.4.0         rmarkdown_2.14             
#> [121] jpeg_0.1-9                  progress_1.2.2             
#> [123] grid_4.2.0                  data.table_1.14.2          
#> [125] blob_1.2.3                  digest_0.6.29              
#> [127] tidyr_1.2.0                 munsell_0.5.0              
#> [129] viridisLite_0.4.0           Gviz_1.40.0                
#> [131] bslib_0.3.1