This function adds region information as an IRanges
object, START
and END
information, to a
DNAStringSet
or an AAStringSet
and puts them into the
metadata
information.
This information can be used to restrict the distance calculation to
specific regions of the DNAStringSet
or the AAStringSet
.
addregion2string(seq, region = NULL, append = TRUE)
DNAStringSet
or AAStringSet
[mandatory]
IRanges
object [mandatory]
indicate if region should be appended or overwritten [default: TRUE]
An object of class DNAStringSet
or AAStringSet
## load example sequence data
data(iupac, package="MSA2dist")
iupac.aa <- iupac |> cds2aa(shorten = TRUE)
## create region
region1 <- IRanges::IRanges(start=c(1,41), end=c(20,50))
## add region
iupac.aa <- iupac.aa |> addregion2string(region=region1)
#(iupac.aa |> slot("metadata"))$region
iupac.aa |> region()
#> IRanges object with 2 ranges and 0 metadata columns:
#> start end width
#> <integer> <integer> <integer>
#> [1] 1 20 20
#> [2] 41 50 10
## append region
region2 <- IRanges::IRanges(start=c(21), end=c(30))
iupac.aa <- iupac.aa |> addregion2string(region=region2)
#(iupac.aa |> slot("metadata"))$region
iupac.aa |> region()
#> IRanges object with 2 ranges and 0 metadata columns:
#> start end width
#> <integer> <integer> <integer>
#> [1] 1 30 30
#> [2] 41 50 10
## overwrite region
iupac.aa <- iupac.aa |> addregion2string(region=region2, append=FALSE)
#(iupac.aa |> slot("metadata"))$region
iupac.aa |> region()
#> IRanges object with 1 range and 0 metadata columns:
#> start end width
#> <integer> <integer> <integer>
#> [1] 21 30 10
## reduce by region
#iupac.aa.region <- iupac.aa |> string2region(region=
# (iupac.aa |> slot("metadata"))$region)
iupac.aa.region <- iupac.aa |> string2region(region=
region(iupac.aa))
#iupac.aa.region |> slot("metadata")
iupac.aa.region |> region()
#> IRanges object with 1 range and 0 metadata columns:
#> start end width
#> <integer> <integer> <integer>
#> [1] 21 30 10