This function collapses isomiRs into different groups. It is a similar concept than how to work with gene isoforms. With this function, different changes can be put together into a single miRNA variant. For instance all sequences with variants at 3' end can be considered as different elements in the table or analysis having the following naming hsa-miR-124a-5p.iso.t3:AAA.

isoCounts(
  ids,
  ref = FALSE,
  iso5 = FALSE,
  iso3 = FALSE,
  add = FALSE,
  snv = FALSE,
  seed = FALSE,
  all = FALSE,
  minc = 1,
  mins = 1,
  merge_by = NULL
)

Arguments

ids

Object of class IsomirDataSeq.

ref

Differentiate reference miRNA from rest.

iso5

Differentiate trimming at 5 miRNA from rest.

iso3

Differentiate trimming at 3 miRNA from rest.

add

Differentiate additions miRNA from rest.

snv

Differentiate nt substitution miRNA from rest.

seed

Differentiate changes in 2-7 nts from rest.

all

Differentiate all isomiRs.

minc

Int minimum number of isomiR sequences to be included.

mins

Int minimum number of samples with number of sequences bigger than minc counts.

merge_by

Column in coldata to merge samples into a single column in counts. Useful to combine technical replicates.

Value

IsomirDataSeq object with new count table. The count matrix can be access with counts(ids).

Details

You can merge all isomiRs into miRNAs by calling the function only with the first parameter isoCounts(ids). You can get a table with isomiRs altogether and the reference miRBase sequences by calling the function with ref=TRUE. You can get a table with 5' trimming isomiRS, miRBase reference and the rest by calling with isoCounts(ids, ref=TRUE, iso5=TRUE). If you set up all parameters to TRUE, you will get a table for each different sequence mapping to a miRNA (i.e. all isomiRs).

Examples for the naming used for the isomiRs are at http://seqcluster.readthedocs.org/mirna_annotation.html#mirna-annotation.

Examples

data(mirData) ids <- isoCounts(mirData, ref=TRUE) head(counts(ids))
#> cc1 cc2 cc3 cc4 cc5 cc6 cc7 ct1 #> hsa-let-7a-2-3p 2 7 2 2 6 0 3 0 #> hsa-let-7a-2-3p;ref 3 6 2 7 5 6 4 0 #> hsa-let-7a-3p 653 593 543 530 594 335 574 219 #> hsa-let-7a-3p;ref 114 114 87 79 137 54 107 39 #> hsa-let-7a-5p 81905 73478 54347 80277 75900 45619 67502 65481 #> hsa-let-7a-5p;ref 235825 171354 149541 180654 168884 107430 153061 143030 #> ct2 ct3 ct4 ct5 ct6 ct7 #> hsa-let-7a-2-3p 0 2 6 6 0 0 #> hsa-let-7a-2-3p;ref 4 0 6 6 8 3 #> hsa-let-7a-3p 448 312 553 510 403 443 #> hsa-let-7a-3p;ref 48 31 82 112 49 76 #> hsa-let-7a-5p 58497 40218 52069 56641 71781 59651 #> hsa-let-7a-5p;ref 163569 114028 123454 133092 158909 140272
# taking into account isomiRs and reference sequence. ids <- isoCounts(mirData, ref=TRUE, minc=10, mins=6) head(counts(ids))
#> cc1 cc2 cc3 cc4 cc5 cc6 cc7 ct1 #> hsa-let-7a-3p 653 593 543 530 594 335 574 219 #> hsa-let-7a-3p;ref 114 114 87 79 137 54 107 39 #> hsa-let-7a-5p 81905 73478 54347 80277 75900 45619 67502 65481 #> hsa-let-7a-5p;ref 235825 171354 149541 180654 168884 107430 153061 143030 #> hsa-let-7b-3p 1033 1729 668 1351 1741 644 1326 335 #> hsa-let-7b-3p;ref 4 220 11 34 143 53 173 3 #> ct2 ct3 ct4 ct5 ct6 ct7 #> hsa-let-7a-3p 448 312 553 510 403 443 #> hsa-let-7a-3p;ref 48 31 82 112 49 76 #> hsa-let-7a-5p 58497 40218 52069 56641 71781 59651 #> hsa-let-7a-5p;ref 163569 114028 123454 133092 158909 140272 #> hsa-let-7b-3p 808 483 1058 1051 1253 890 #> hsa-let-7b-3p;ref 123 11 91 143 178 98