DESeqDataSetFromSlidingWindows: create DESeq data object

View source: R/preprocess.R

DESeqDataSetFromSlidingWindowsR Documentation

create DESeq data object

Description

create DESeq data object from sliding window counts, phenotype data and annotation data

Usage

DESeqDataSetFromSlidingWindows(
  countData,
  colData,
  annotObj,
  design,
  tidy = FALSE,
  ignoreRank = FALSE,
  start0based = TRUE
)

Arguments

countData

data.frame or matrix, sliding window count data

colData

DataFrame or data.frame, phenotype data, see DESeqDataSet

annotObj

data.frame or character, can either be a data.frame or a file name, see details

design

formula or matrix, design of the experiment, see DESeqDataSet

tidy

logical, If TRUE, first column is of countData is treated as rownames (defalt: FALSE), see DESeqDataSet

ignoreRank

logical, ignore rank, see DESeqDataSet

start0based

logical, TRUE (default) or FALSE. If TRUE, then the start positions in annotObj is considered to be 0-based

Details

If annotObj is a file name, the input file MUST be <TAB> separated, and supports reading in .gz files.
If annotObj is a data.frame, colnames(annotObj) MUST not be empty.
This function checks for the following columns after reading in the file or on data.frame:

  • chromosome: chromosome name

  • unique_id: unique id of the window, rownames(object) must match this column

  • begin: window start co-ordinate, see parameter start0based

  • end: window end co-ordinate

  • strand: strand

  • gene_id: gene id

  • gene_name: gene name

  • gene_type: gene type annotation

  • gene_region: gene region

  • Nr_of_region: number of the current region

  • Total_nr_of_region: total number of regions

  • window_number: window number

This function creates a DESeqDataSet using supplied countData, phenotype data and annotation data. The chromosomal locations and annotations of the sliding windows (parsed from annotObj) can be accessed from the returned object using: rowRanges(object)

Value

DESeq object

Examples


data("SLBP_K562_w50s20")
slbpDat <- counts(SLBP_K562_w50s20)
phenoDat <- DataFrame(conditions=as.factor(c(rep('IP',2),'SMI')),
row.names = colnames(slbpDat))
phenoDat$conditions <- relevel(phenoDat$conditions,ref='SMI')
annotDat <- as.data.frame(rowRanges(SLBP_K562_w50s20))
# by default chromsome column is 'seqnames'
# and begin co-ordinate column is 'start'
# rename these columns
colnames(annotDat)[1:2] <- c('chromosome','begin')
slbpDds <- DESeqDataSetFromSlidingWindows(countData = slbpDat,
colData = phenoDat,annotObj = annotDat,design=~conditions)


EMBL-Hentze-group/DEWSeq documentation built on Oct. 17, 2023, 10:41 p.m.