writeGctCls: Export ExpressionSet as Gct/Cls files

View source: R/io_gct_cls.R

writeClsR Documentation

Export ExpressionSet as Gct/Cls files

Description

Gct/Cls file formats are required by the Gene Set Enrichment Analysis (GSEA) tool. Functions writeGct and writeCls exports file of two formats respectively, and writeGctCls calls the two function internally to write two files.

Usage

writeCls(eset, file = stdout(), sample.group.col = "group")

writeGctCls(
  eset,
  file.base,
  feat.name,
  feat.desc,
  sample.group.col,
  write.add.fData.file = TRUE,
  write.add.pData.file = TRUE
)

Arguments

eset

An object of the eSet class, for example an ExpressionSet object

file

Name of the Gct/Cls file. If left missing, the file is printed on the standard output.

sample.group.col

Integer, character or a factor vector of the same length as the sample number, indicating classes (groups) of samples. See details.

file.base

For writeGctCls, the base name of the two files: the suffix (.gct and .cls) will be appended

feat.name

Integer or character, indicating which column of the featureData should be used as feature descriptions. If the value is missing, the Description column of the Gct file will be left blank

feat.desc

Integer or character, indicating which column of the featureData should be used as feature names; if missing, results of the featureNames function will be used as identifiers in the Gct file. See details.

write.add.fData.file

Logical, whether additional featureData should be written into a file named ${file.base}.add.fData.txt

write.add.pData.file

Logical, whether additional phenoData should be written into a file named ${file.base}.add.pData.txt

Details

The feat.name option specifies what identifiers should be used for features (probesets). When the value is missing, featureNames is called to provide feature identifiers.

In contrast, the sample.group.col cannot be missing: since cls files encode groups (classes) of samples, and if sample.group.col was missing, it is usually impossible to get class information from sampleNames.

Internally writeCls calls dfFactor function to determine factor of samples. Therefore sample.group.col is to a certain degree generic: it can be a character string or integer index of the pData(eset) data matrix, or a factor vector of the same length as ncol(eset).

Value

Functions are used for their side effects.

Functions

  • writeCls(): writeCls

Author(s)

Jitao David Zhang <jitao_david.zhang@roche.com>

References

http://www.broadinstitute.org/gsea/doc/GSEAUserGuideTEXT.htm

See Also

See dfFactor for possible values of the sample.group.col option.

See readGctCls for importing functions.

Examples


data(sample.ExpressionSet, package="Biobase")
writeGct(sample.ExpressionSet[1:5, 1:4], file=stdout())
writeCls(sample.ExpressionSet, file=stdout(), sample.group.col="type")

tmpfile <- tempfile()
writeGctCls(sample.ExpressionSet, file.base=tmpfile, sample.group.col="type")
readLines(paste(tmpfile, ".cls",sep=""))
unlink(c(paste(tmpfile, ".cls", sep=""), paste(tmpfile, ".gct", sep="")))


bedapub/ribiosExpression documentation built on Sept. 2, 2023, 4:37 a.m.