| Title: | Parameter Optimization for Flow Cytometry Data Transformation |
|---|---|
| Description: | Profile maximum likelihood estimation of parameters for flow cytometry data transformations. |
| Authors: | Greg Finak <[email protected]>, Juan Manuel-Perez <[email protected]>, Raphael Gottardo <[email protected]> |
| Maintainer: | Greg Finak <[email protected]> |
| License: | Artistic-2.0 |
| Version: | 1.65.0 |
| Built: | 2026-05-30 08:16:07 UTC |
| Source: | https://github.com/bioc/flowTrans |
Maximum likelihood estimation of parameters for common flow cytometry data transformations.
| Package: | flowTrans |
| Type: | Package |
| Version: | 0.6.0 |
| Date: | 2010-03-25 |
| License: | Artistic-2.0 |
| LazyLoad: | yes |
| biocViews: | Bioinformatics, FlowCytometry |
| Depends: | flowCore, methods |
| Imports: | flowCore, methods, stats, flowViz, flowClust |
| Collate: | logicleTransformWrapper.R mclMultivArcSinh.R mclMultivBiexp.R mclMultivLinLog.R mclMultivBoxCox.R arcsinhTransformWrapper.R biexponentialTransformWrapper.R boxcoxTransform.R flowTrans.R linLog.R |
| Packaged: | 2010-03-25 16:49:03 UTC; finak |
| Built: | R 2.10.0; ; 2009-11-30 16:49:05 UTC; unix |
Index:
flowTrans Optimizing transformations for flow cytometry
data
Greg Finak <[email protected]>, Juan Manuel-Perez <[email protected]>, Raphael Gottardo <[email protected]>
Maintainer: Greg Finak <[email protected]>
Finak G, Perez J M, Weng A, Gottardo R. Optimizing Transformations for Flow Cytometry. (Submitted)
#Load some data data(GvHD) #transform the first sample, forward and side scatter. result<-flowTrans(GvHD[[1]],"mclMultivArcSinh",colnames(GvHD[[1]])[1:2],n2f=FALSE,parameters.only=FALSE); plot(result$result); summary(result);#Load some data data(GvHD) #transform the first sample, forward and side scatter. result<-flowTrans(GvHD[[1]],"mclMultivArcSinh",colnames(GvHD[[1]])[1:2],n2f=FALSE,parameters.only=FALSE); plot(result$result); summary(result);
Extracts the transformation parameters from a flowTransResult object
extractParams(x, dims = NULL)extractParams(x, dims = NULL)
x |
An object of type |
dims |
A character vector specifying the dimensions for which to extract transformation parameters. |
A list of length length(dims). Each element contains a vector of parameters for transforming the particular dimension.The names of the list elements correspond to the names of the dimensions.
flowTrans estimates common transformation parameters.
Greg Finak <[email protected]>, Raphael Gottardo <[email protected]>
Finak G, Perez JM, Weng A, Gottardo R. Optimizing Data Transformation for Flow Cytometry. (Submitted)
summary,
flowTrans,
flowTransResult
m <- t(matrix(rnorm(10000),2)) colnames(m) <- c("A","B") m <- flowFrame(m) res <- flowTrans(m,"mclMultivBoxCox", c("A","B"), n2f=FALSE, parameters.only=FALSE) summary(res) extractParams(res)m <- t(matrix(rnorm(10000),2)) colnames(m) <- c("A","B") m <- flowFrame(m) res <- flowTrans(m,"mclMultivBoxCox", c("A","B"), n2f=FALSE, parameters.only=FALSE) summary(res) extractParams(res)
This is the principal function in the package. The function takes data,
a transformation function name, dimension names, and optional
preselection parameter as input, and runs the maximum likelihood
optimization on the data, performs the transformation, and returns the
transformed data together with a list of transformation paramters for
each dimension. The optimization is set to maximize the likelihood of
the parameters given the data, such that the transformed data is
normally distributed. An optional parameter n2f=[T,F] allows an
automated preselection of an approximately bivariate normal population
via the norm2Filter function. The transformation parameters are
then optimized for this preselected region, and finally applied globally
to the entire range of data. The optional argument
parameters.only=[T,F] specifies whether to return the parameters
or the transformed data.
flowTrans(dat, fun, dims, n2f,parameters.only)flowTrans(dat, fun, dims, n2f,parameters.only)
dat |
The data to be transformed, should be a |
fun |
A character string naming the transformation function to be
applied. Can be one of: |
dims |
A character vector identifying the dimensions to be transformed. |
n2f |
An optional |
parameters.only |
A logical specifying whether to return only the parameters, and
not the transformed data.
|
The transformation functions are multivariate, common parameter
transformations. The implementation utilizes a look up table to call
optimization routines for the different transformations as well as
optional flowCore and flowClust transformation
implementations to optimize parameters and transform data from a
single interface. Parameters are optimized such that the transformed
data is multivariate-normal.
if parameters.only = FALSE then:
flowTransResult |
A |
If parameters.only=TRUE the returned value will be a vector of
common transformation parameters:
if mclMultivArcSinh is called, the returned vector will
contain a, b, and c=0;
if mclMultivBiexp is called, the returned vector will
contain a, b, c, d, w;
if mclMultivBoxCox is called, the returned vector will
contain theta;
if mclMultivLinLog is called, the returned vector will
contain theta.
Greg Finak <[email protected]>, Raphael Gottardo <[email protected]>
Finak G, Perez J M, Weng A, Gottardo R. Optimizing Transformations for Flow Cytometry. (Submitted)
data(GvHD) result <- flowTrans(GvHD[[1]], "mclMultivArcSinh", colnames(GvHD[[1]])[1:2], n2f=FALSE, parameters.only=FALSE) par(mfrow=c(1,2)) contour(GvHD[[1]][,1:2], main="Untransformed FSC vs SSC") contour(result$result[,1:2], main="Transformed FSC vs SSC") summary(result)data(GvHD) result <- flowTrans(GvHD[[1]], "mclMultivArcSinh", colnames(GvHD[[1]])[1:2], n2f=FALSE, parameters.only=FALSE) par(mfrow=c(1,2)) contour(GvHD[[1]][,1:2], main="Untransformed FSC vs SSC") contour(result$result[,1:2], main="Transformed FSC vs SSC") summary(result)
Transforms a flowFrame and optimizes the parameters for the selected transformation.
Transforms the flowFrame using the function fun over the dimensions dims, with optional bivaraite normal prefiltering. Parameters are optimized to make the transformed data multivariate normal depending on the selected transformation
Class that contains the results of a call to flowTrans on a flowFrame
Objects can be created by calls of the form new("flowTransResult", ...).
The object contains three slots: the transformed flowFrame, a transform object from which the parameters can be extracted, and a dims slot containing a character vector of the dimensions that have been transformed.
result:Object of class "flowFrame". The transformed data.
trans:Object of class "transform". The transform applied to the flowFrame. Parameters can be extracted via summary(trans).
dims:Object of class "character". The names of the dimensions that have been transformed.
signature(object = "flowTransResult"): Summarize the transformed data.
Greg Finak <[email protected]>, Raphael Gottardo <[email protected]>.
Finak G, Perez JM, Weng A, Gottardo R. Optimizing Data Transformation for Flow Cytometry.
m <- t(matrix(rnorm(10000),2)) colnames(m) <- c("A","B") m <- flowFrame(m) res <- flowTrans(m,"mclMultivBoxCox", c("A","B"), n2f=FALSE, parameters.only=FALSE) summary(res) extractParams(res)m <- t(matrix(rnorm(10000),2)) colnames(m) <- c("A","B") m <- flowFrame(m) res <- flowTrans(m,"mclMultivBoxCox", c("A","B"), n2f=FALSE, parameters.only=FALSE) summary(res) extractParams(res)
Print a summary of the contents of a flowTransResult object.
Prints a summary of the flowTransResult object
Print a summary of the contents of a flowTransResult object, including the type of transformation, the transformed data frame, and the transformed dimensions and transformation parameters.
summary.flowTransResult(object, ...)summary.flowTransResult(object, ...)
object |
An object of |
... |
Additional parameters passed through. |
No return value
Greg Finak <[email protected]>, Raphael Gottardo <[email protected]>
Finak G, Perez JM, Weng A, Gottardo R. Optimizing Transformations for Flow Cytometry. (Submitted)
extractParams, flowTransResult
m <- t(matrix(rnorm(10000),2)) colnames(m) <- c("A","B") m <- flowFrame(m) res <- flowTrans(m, "mclMultivBoxCox", c("A","B"), n2f=FALSE, parameters.only=FALSE) summary(res)m <- t(matrix(rnorm(10000),2)) colnames(m) <- c("A","B") m <- flowFrame(m) res <- flowTrans(m, "mclMultivBoxCox", c("A","B"), n2f=FALSE, parameters.only=FALSE) summary(res)