| writeCls | R Documentation |
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.
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
)
eset |
An object of the |
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
|
write.add.fData.file |
Logical, whether additional featureData should
be written into a file named |
write.add.pData.file |
Logical, whether additional phenoData should be
written into a file named |
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).
Functions are used for their side effects.
writeCls(): writeCls
Jitao David Zhang <jitao_david.zhang@roche.com>
http://www.broadinstitute.org/gsea/doc/GSEAUserGuideTEXT.htm
See dfFactor for possible values of the
sample.group.col option.
See readGctCls for importing functions.
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="")))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.