Contents

1 Introduction

This R package provides methods for genetic finemapping in inbred mice by taking advantage of their very high homozygosity rate (>95%).

Method prio allows to select strain combinations which best refine a specified genetic region. E.g. if a crossing experiment with two inbred mouse strains ‘strain1’ and ‘strain2’ resulted in a QTL, the outputted strain combinations can be used to refine the respective region in further crossing experiments and to select candidate genes.

2 Installation

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

3 Loading package

library(MouseFM)

4 Example function calls

Available mouse strains

avail_strains()
#>              id        strain
#> 1  129P2_OlaHsd  129P2/OlaHsd
#> 2   129S1_SvImJ   129S1/SvImJ
#> 3  129S5SvEvBrd 129S5/SvEvBrd
#> 4           A_J           A/J
#> 5         AKR_J         AKR/J
#> 6       BALB_cJ       BALB/cJ
#> 7          BTBR          BTBR
#> 8       BUB_BnJ       BUB/BnJ
#> 9       C3H_HeH       C3H/HeH
#> 10      C3H_HeJ       C3H/HeJ
#> 11    C57BL_10J     C57BL/10J
#> 12     C57BL_6J      C57BL/6J
#> 13    C57BL_6NJ     C57BL/6NJ
#> 14    C57BR_cdJ     C57BR/cdJ
#> 15       C57L_J        C57L/J
#> 16        C58_J         C58/J
#> 17     CAST_EiJ      CAST/EiJ
#> 18        CBA_J         CBA/J
#> 19       DBA_1J        DBA/1J
#> 20       DBA_2J        DBA/2J
#> 21       FVB_NJ        FVB/NJ
#> 22        I_LnJ         I/LnJ
#> 23       KK_HiJ        KK/HiJ
#> 24    LEWES_EiJ     LEWES/EiJ
#> 25         LP_J          LP/J
#> 26     MOLF_EiJ      MOLF/EiJ
#> 27   NOD_ShiLtJ    NOD/ShiLtJ
#> 28     NZB_B1NJ      NZB/B1NJ
#> 29    NZO_HlLtJ     NZO/HlLtJ
#> 30     NZW_LacJ      NZW/LacJ
#> 31      PWK_PhJ       PWK/PhJ
#> 32         RF_J          RF/J
#> 33      SEA_GnJ       SEA/GnJ
#> 34    SPRET_EiJ     SPRET/EiJ
#> 35        ST_bJ         ST/bJ
#> 36      WSB_EiJ       WSB/EiJ
#> 37  ZALENDE_EiJ   ZALENDE/EiJ

Prioritize additional mouse strains for a given region which was identified in a crossing experiment with strain1 C57BL_6J and strain2 AKR_J.

df = prio("chr1", start=5000000, end=6000000, strain1="C57BL_6J", strain2="AKR_J")
#> Query chr1:5,000,000-6,000,000
#> Calculate reduction factors...
#> Set size 1: 35 combinations
#> Set size 1: continue with 20 of 35 strains
#> Set size 2: 190 combinations
#> Set size 3: 1,140 combinations

View meta information

comment(df)
#> NULL

Extract the combinations with the best refinement

get_top(df$reduction, n_top=3)
#>    strain1 strain2              combination      mean       min       max n
#> 8 C57BL_6J   AKR_J C3H_HeH,DBA_1J,SPRET_EiJ 0.8068057 0.7467057 0.9926794 3
#> 7 C57BL_6J   AKR_J C3H_HeH,DBA_2J,SPRET_EiJ 0.8068057 0.7467057 0.9926794 3
#> 6 C57BL_6J   AKR_J C3H_HeJ,DBA_1J,SPRET_EiJ 0.8068057 0.7467057 0.9926794 3

Create plots

plots = vis_reduction_factors(df$genotypes, df$reduction, 2)
plots[[1]]

plots[[2]]

5 Output of Session Info

The output of sessionInfo() on the system on which this document was compiled:

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] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] MouseFM_1.6.0    BiocStyle_2.24.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] Biobase_2.56.0         httr_1.4.2             tidyr_1.2.0           
#>  [4] sass_0.4.1             bit64_4.0.5            jsonlite_1.8.0        
#>  [7] gtools_3.9.2           bslib_0.3.1            assertthat_0.2.1      
#> [10] highr_0.9              BiocManager_1.30.17    stats4_4.2.0          
#> [13] BiocFileCache_2.4.0    blob_1.2.3             GenomeInfoDbData_1.2.8
#> [16] yaml_2.3.5             progress_1.2.2         pillar_1.7.0          
#> [19] RSQLite_2.2.12         rlist_0.4.6.2          glue_1.6.2            
#> [22] digest_0.6.29          GenomicRanges_1.48.0   XVector_0.36.0        
#> [25] colorspace_2.0-3       plyr_1.8.7             htmltools_0.5.2       
#> [28] XML_3.99-0.9           pkgconfig_2.0.3        biomaRt_2.52.0        
#> [31] magick_2.7.3           bookdown_0.26          zlibbioc_1.42.0       
#> [34] purrr_0.3.4            scales_1.2.0           tibble_3.1.6          
#> [37] KEGGREST_1.36.0        farver_2.1.0           generics_0.1.2        
#> [40] IRanges_2.30.0         ggplot2_3.3.5          ellipsis_0.3.2        
#> [43] cachem_1.0.6           BiocGenerics_0.42.0    cli_3.3.0             
#> [46] magrittr_2.0.3         crayon_1.5.1           memoise_2.0.1         
#> [49] evaluate_0.15          fansi_1.0.3            xml2_1.3.3            
#> [52] tools_4.2.0            data.table_1.14.2      prettyunits_1.1.1     
#> [55] hms_1.1.1              lifecycle_1.0.1        stringr_1.4.0         
#> [58] S4Vectors_0.34.0       munsell_0.5.0          AnnotationDbi_1.58.0  
#> [61] Biostrings_2.64.0      compiler_4.2.0         jquerylib_0.1.4       
#> [64] GenomeInfoDb_1.32.0    rlang_1.0.2            grid_4.2.0            
#> [67] RCurl_1.98-1.6         rappdirs_0.3.3         bitops_1.0-7          
#> [70] rmarkdown_2.14         gtable_0.3.0           DBI_1.1.2             
#> [73] curl_4.3.2             reshape2_1.4.4         R6_2.5.1              
#> [76] knitr_1.38             dplyr_1.0.8            fastmap_1.1.0         
#> [79] bit_4.0.4              utf8_1.2.2             filelock_1.0.2        
#> [82] stringi_1.7.6          Rcpp_1.0.8.3           vctrs_0.4.1           
#> [85] png_0.1-7              dbplyr_2.1.1           tidyselect_1.1.2      
#> [88] xfun_0.30