Last data update: 2014.03.03

R: Plot the features with the highest expression values
plotHighestExprsR Documentation

Plot the features with the highest expression values

Description

Plot the features with the highest expression values

Usage

plotHighestExprs(object, col_by_variable = "total_features", n = 50,
  drop_features = NULL, exprs_values = "counts")

Arguments

object

an SCESet object containing expression values and experimental information. Must have been appropriately prepared.

col_by_variable

variable name (must be a column name of pData(object)) to be used to assign colours to cell-level values.

n

numeric scalar giving the number of the most expressed features to show. Default value is 50.

drop_features

a character, logical or numeric vector indicating which features (e.g. genes, transcripts) to drop when producing the plot. For example, control genes might be dropped to focus attention on contribution from endogenous rather than synthetic genes.

exprs_values

which slot of the assayData in the object should be used to define expression? Valid options are "counts" (default), "tpm", "fpkm" and "exprs".

Details

Plot the percentage of counts accounted for by the top n most highly expressed features across the dataset.

Value

a ggplot plot object

Examples

data("sc_example_counts")
data("sc_example_cell_info")
pd <- new("AnnotatedDataFrame", data = sc_example_cell_info)
rownames(pd) <- pd$Cell
example_sceset <- newSCESet(countData = sc_example_counts, phenoData = pd)
example_sceset <- calculateQCMetrics(example_sceset, feature_controls = 1:500)
plotHighestExprs(example_sceset, col_by_variable="total_features")
plotHighestExprs(example_sceset, col_by_variable="Mutation_Status")

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(scater)
Loading required package: Biobase
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

Loading required package: ggplot2

Attaching package: 'scater'

The following object is masked from 'package:stats':

    filter

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/scater/plotHighestExprs.Rd_%03d_medium.png", width=480, height=480)
> ### Name: plotHighestExprs
> ### Title: Plot the features with the highest expression values
> ### Aliases: plotHighestExprs
> 
> ### ** Examples
> 
> data("sc_example_counts")
> data("sc_example_cell_info")
> pd <- new("AnnotatedDataFrame", data = sc_example_cell_info)
> rownames(pd) <- pd$Cell
> example_sceset <- newSCESet(countData = sc_example_counts, phenoData = pd)
> example_sceset <- calculateQCMetrics(example_sceset, feature_controls = 1:500)
> plotHighestExprs(example_sceset, col_by_variable="total_features")
> plotHighestExprs(example_sceset, col_by_variable="Mutation_Status")
> 
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>