View source: R/Script_DROPLET_01_CREATE_MARVEL_OBJECT.R
| CreateMarvelObject.10x | R Documentation | 
Creates an S3 object named Marvel for downstream analysis, specifically for droplet-based RNA-sequencing data.
CreateMarvelObject.10x( gene.norm.matrix = NULL, gene.norm.pheno = NULL, gene.norm.feature = NULL, gene.count.matrix = NULL, gene.count.pheno = NULL, gene.count.feature = NULL, sj.count.matrix = NULL, sj.count.pheno = NULL, sj.count.feature = NULL, pca = NULL, gtf = NULL )
| gene.norm.matrix | Sparse matrix. UMI-collapsed, normalised, non-log2-transformed gene expression matrix. | 
| gene.norm.pheno | Data frame. Sample metadata for annotating  | 
| gene.norm.feature | Data frame. Gene metadata for annotating  | 
| gene.count.matrix | Sparse matrix. UMI-collapsed, non-normalised (raw counts), non-log2-transformed gene expression matrix. | 
| gene.count.pheno | Data frame. Sample metadata for annotating  | 
| gene.count.feature | Data frame. Gene metadata for annotating  | 
| sj.count.matrix | Sparse matrix. UMI-collapsed, non-normalised (raw counts), non-log2-transformed splice junction expression matrix. | 
| sj.count.pheno | Data frame. Sample metadata for annotating  | 
| sj.count.feature | Data frame. Splice junction metadata for annotating  | 
| pca | Data frame. Coordinates of PCA/tSNE/UMAP. | 
| gtf | Data frame. GTF used in cellranger. Will be used for annotating splice junctions downstream. | 
An object of class S3.
# Retrieve, observe format of pre-saved input files
marvel.demo.10x.raw <- readRDS(system.file("extdata/data",
                               "marvel.demo.10x.raw.rds",
                               package="MARVEL")
                               )
# Gene expression (Normalised)
    # Matrix
    df.gene.norm <- marvel.demo.10x.raw$gene.norm.matrix
    df.gene.norm[1:5, 1:5]
    # phenoData
    df.gene.norm.pheno <- marvel.demo.10x.raw$sample.metadata
    head(df.gene.norm.pheno)
    # featureData
    df.gene.norm.feature <- data.frame("gene_short_name"=rownames(df.gene.norm),
                                       stringsAsFactors=FALSE
                                       )
    head(df.gene.norm.feature)
# Gene expression (Counts)
    # Matrix
    df.gene.count <- marvel.demo.10x.raw$gene.count.matrix
    df.gene.count[1:5, 1:5]
    # phenoData
    df.gene.count.pheno <- data.frame("cell.id"=colnames(df.gene.count),
                                       stringsAsFactors=FALSE
                                       )
    head(df.gene.count.pheno)
    # featureData
    df.gene.count.feature <- data.frame("gene_short_name"=rownames(df.gene.count),
                                       stringsAsFactors=FALSE
                                       )
    head(df.gene.count.feature)
# SJ (Counts)
    # Matrix
    df.sj.count <- marvel.demo.10x.raw$sj.count.matrix
    df.sj.count[1:5, 1:5]
    # phenoData
    df.sj.count.pheno <- data.frame("cell.id"=colnames(df.sj.count),
                                     stringsAsFactors=FALSE
                                     )
    head(df.sj.count.pheno)
    # featureData
    df.sj.count.feature <- data.frame("coord.intron"=rownames(df.sj.count),
                                       stringsAsFactors=FALSE
                                       )
    head(df.sj.count.feature)
# tSNE coordinates
df.coord <- marvel.demo.10x.raw$pca
head(df.coord)
# GTF
gtf <- marvel.demo.10x.raw$gtf
head(gtf)
# Create MARVEL object
marvel.demo.10x <- CreateMarvelObject.10x(gene.norm.matrix=df.gene.norm,
                     gene.norm.pheno=df.gene.norm.pheno,
                     gene.norm.feature=df.gene.norm.feature,
                     gene.count.matrix=df.gene.count,
                     gene.count.pheno=df.gene.count.pheno,
                     gene.count.feature=df.gene.count.feature,
                     sj.count.matrix=df.sj.count,
                     sj.count.pheno=df.sj.count.pheno,
                     sj.count.feature=df.sj.count.feature,
                     pca=df.coord,
                     gtf=gtf
                     )
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.