dnastring2dist: dnastring2dist

Description Usage Arguments Value Author(s) See Also Examples

View source: R/dnastring2dist.R

Description

This function calculates pairwise distances for all combinations of a DNAStringSet.

Usage

1
2
3
4
5
6
7
8
9
dnastring2dist(
  dna,
  model = "IUPAC",
  threads = 1,
  score = NULL,
  mask = NULL,
  region = NULL,
  ...
)

Arguments

dna

DNAStringSet [mandatory]

model

specify model either "IUPAC" or any model from ape::dist.dna [default: IUPAC]

threads

number of parallel threads [default: 1]

score

score matrix use score matrix to calculate distances [default: NULL]

mask

IRanges object indicating masked sites [default: NULL]

region

IRanges object indicating region to use for dist calculation. Default is null, meaning all sites are used [default: NULL]

...

other ape::dist.dna parameters (see dist.dna)

Value

A data.frame of pairwise distance values distSTRING and sites used sitesUsed

Author(s)

Kristian K Ullrich

See Also

dist.dna

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
## load example sequence data
data("hiv", package="distSTRING")
#dnastring2dist(hiv, model="IUPAC")
hiv |> dnastring2dist(model="IUPAC")
#dnastring2dist(hiv, model="K80")
hiv |> dnastring2dist(model="K80")
data("woodmouse", package="ape")
#dnastring2dist(dnabin2dnastring(woodmouse), score=iupacMatrix())
woodmouse |> dnabin2dnastring() |> dnastring2dist()
#dnastring2dist(hiv, model = "IUPAC", threads = 2)
hiv |> dnastring2dist(model = "IUPAC", threads = 2)
## create mask
mask1 <- IRanges::IRanges(start=c(1,61,121), end=c(30,90,150))
## use mask
hiv |> dnastring2dist(model="IUPAC", mask=mask1)
## use region
region1 <- IRanges::IRanges(start=c(1,139), end=c(75,225))
hiv |> dnastring2dist(model="IUPAC", region=region1)
## use mask and region
hiv |> dnastring2dist(model="IUPAC", mask=mask1, region=region1)

kullrich/distSTRING documentation built on Dec. 21, 2021, 8:42 a.m.