View source: R/addregion2string.R
addregion2string | R Documentation |
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)
seq |
|
region |
|
append |
indicate if region should be appended or overwritten [default: TRUE] |
An object of class DNAStringSet
or AAStringSet
Kristian K Ullrich
addmask2string
,
addpop2string
,
addpos2string
## 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()
## 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()
## overwrite region
iupac.aa <- iupac.aa |> addregion2string(region=region2, append=FALSE)
#(iupac.aa |> slot("metadata"))$region
iupac.aa |> region()
## 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()
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.