Title: | Open-source toolkit to analyse data from xCELLigence System (RTCA) |
---|---|
Description: | Import, analyze and visualize data from Roche(R) xCELLigence RTCA systems. The package imports real-time cell electrical impedance data into R. As an alternative to commercial software shipped along the system, the Bioconductor package RTCA provides several unique transformation (normalization) strategies and various visualization tools. |
Authors: | Jitao David Zhang |
Maintainer: | Jitao David Zhang <[email protected]> |
License: | LGPL-3 |
Version: | 1.59.0 |
Built: | 2024-10-31 04:28:06 UTC |
Source: | https://github.com/bioc/RTCA |
Functions to manipulate indices or names of microtitre plates
alphaNames(row = 8, column = 12, order=c("column","row")) repairAlphaName(x) alphaNames2Pos(x) rowcol2pos(row = 1, column=1, plateFormat=c("96","384"))
alphaNames(row = 8, column = 12, order=c("column","row")) repairAlphaName(x) alphaNames2Pos(x) rowcol2pos(row = 1, column=1, plateFormat=c("96","384"))
row |
integer, row index, 1,...,8 for 96-well plates |
column |
integer, column index, 1,...,12 for 96-well plates |
x |
character, Well alpha name, in the form of [A-Z][0-9][0-9], like 'A01' |
order |
character, should the alpha names returned in a row-first or column-first order? |
plateFormat |
integer, the microtitre format, either 96 or 384 |
alphaNames
returns so-called alpha well names in the form of
[A-H][0-9][0-9] (i.e., A01, C03, D11, H12) for microtitre plates. The
order of returned alphaNames is controlled by the option order
,
which can be set either as col
or row
repairAlphaName
attempts to fix incomplete alpha well names. Now it
is mainly used to fix well names missing the leading 0 of numeric
index, like A1.
alphaName2Pos
returns the row and column number of the given
alpha well name, in the form of two-column data frame with row
and col as colnames.
rowcol2pos
returns the row-wise position index of given row and
column index.
See details
Jitao David Zhang [email protected]
wells <- alphaNames() repairAlphaName("A1") alphaNames2Pos(c("A01","B02","C03","H12")) rowcol2pos(3,1)
wells <- alphaNames() repairAlphaName("A1") alphaNames2Pos(c("A01","B02","C03","H12")) rowcol2pos(3,1)
Combine a list of RTCA objects
combineRTCA(list)
combineRTCA(list)
list |
A list of |
The current implementation requires all the objects have exactly the same time-points recorded (or at least of same length).
The combined RTCA
object has an obligatory column in the
phenoData
‘Plate’ (upper-case!), which matches the names of the
RTCA
list. When the list
has no names, the
‘Plate’ field is filled with integer index starting from 1.
A new RTCA
object
Special attention should be given to the cases where the list
parameter partially has names. In this case all items without name
will be assigned to a ‘Plate’ field of empty string
(“”). Therefore it is advised either to assign names to all
items of the list, or leave them all off.
Jitao David Zhang [email protected]
## An artificial example require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xSub1 <- x[,1:3] xSub2 <- x[,4:ncol(x)] xComb <- combineRTCA(list(sub1=xSub1, sub2=xSub2)) identical(exprs(x), exprs(xComb)) pData(xComb)$Plate ## in case of nameless list pData(combineRTCA(list(xSub1, xSub2)))$Plate ## partial names pData(combineRTCA(list(a=xSub1, xSub2)))$Plate
## An artificial example require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xSub1 <- x[,1:3] xSub2 <- x[,4:ncol(x)] xComb <- combineRTCA(list(sub1=xSub1, sub2=xSub2)) identical(exprs(x), exprs(xComb)) pData(xComb)$Plate ## in case of nameless list pData(combineRTCA(list(xSub1, xSub2)))$Plate ## partial names pData(combineRTCA(list(a=xSub1, xSub2)))$Plate
A convenience function to plot sample wells with control wells on an
E-plate in RTCA system. To use the function the phenoData field
of the RTCA
object must contain a field named “GeneSymbol”.
controlView(rtca, genesymbol = c("Allstar", "COPB2", "GFP", "mock", "PLK1", "WEE1"), cols, ylim, smooth = FALSE, group = TRUE, ylab = "Normalized cell index", xlab = "Time interval (hour)", drawsd = TRUE, normline = TRUE, ncol = 1, legendpos = "topleft", pData.column="GeneSymbol",...)
controlView(rtca, genesymbol = c("Allstar", "COPB2", "GFP", "mock", "PLK1", "WEE1"), cols, ylim, smooth = FALSE, group = TRUE, ylab = "Normalized cell index", xlab = "Time interval (hour)", drawsd = TRUE, normline = TRUE, ncol = 1, legendpos = "topleft", pData.column="GeneSymbol",...)
rtca |
An object of |
genesymbol |
character, gene symbols to be plotted. |
cols |
character, colors used by the provided gene symbols |
ylim |
y-axis lim |
smooth |
logical, whether the |
group |
logical. If ‘group’ is set to |
ylab |
y axis label |
xlab |
x axis label |
drawsd |
logical, should the error bar be drawn to represent standard deviation? |
normline |
logical, should the base-time indicated by a line? See
|
ncol |
integer, legend column number |
legendpos |
character, legend position |
pData.column |
The column which the |
... |
other parameters passed to the |
The function is often called to draw sample and control in one plot.
NULL
, the function is called for its side effect
Jitao David Zhang [email protected]
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") pfile <- system.file("extdata/testOutputPhenoData.csv", package="RTCA") pData <- read.csv(pfile, sep="\t", row.names="Well") metaData <- data.frame(labelDescription=c( "Rack number", "siRNA catalogue number", "siRNA gene symbol", "siRNA EntrezGene ID", "siRNA targeting accession" )) phData <- new("AnnotatedDataFrame", data=pData, varMetadata=metaData) x <- parseRTCA(ofile, phenoData=phData) controlView(x, genesymbol=c("mock","COPB2","PLK1"),ylim=c(0,2))
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") pfile <- system.file("extdata/testOutputPhenoData.csv", package="RTCA") pData <- read.csv(pfile, sep="\t", row.names="Well") metaData <- data.frame(labelDescription=c( "Rack number", "siRNA catalogue number", "siRNA gene symbol", "siRNA EntrezGene ID", "siRNA targeting accession" )) phData <- new("AnnotatedDataFrame", data=pData, varMetadata=metaData) x <- parseRTCA(ofile, phenoData=phData) controlView(x, genesymbol=c("mock","COPB2","PLK1"),ylim=c(0,2))
Derivative transform of RTCA object, returning the change rate of cell impedance
derivativeTransform(object)
derivativeTransform(object)
object |
An object of |
The first derivative of the cell impedance curve measured by RTCA. The derivative of the last time point is estimated by that of the next to last point.
An RTCA
object populated with derivative values
Jitao David Zhang [email protected]
smoothTransform
and interpolationTransform
for smoothing and interpolating the RTCA
data. rgrTransform
calculates relative growth rate,
which calls derivativeTransform
.
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xDeriv <- derivativeTransform(x)
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xDeriv <- derivativeTransform(x)
The functions implement easy interface to certain tasks of factor. See datails for explaination
factor2numeric(x) relevels(x, refs)
factor2numeric(x) relevels(x, refs)
x |
A vector of factor |
refs |
A vector of character, reference vector to give the orderof levels |
relevels
re-arrange the order of levels by the given
character refs
. Alternatively user could use
factor(...,levels=refs)
to achieve a similar effect, however
the relevels
enables also partial list. The missing
levels in refs
will be ordered to the last.
factor2numeric
converts factor of numerics into their
numeric form.
A vector of factor
Jitao David Zhang [email protected]
## factor2numeric numFac <- factor(c(3.5, 2.5, 2.5,3.5, 1)) numFac levels(numFac) factor2numeric(numFac) class(factor2numeric(numFac)) ## relevels relevels(numFac, c("3.5", "1", "2.5")) relevels(numFac, c("3.5", "2.5"))
## factor2numeric numFac <- factor(c(3.5, 2.5, 2.5,3.5, 1)) numFac levels(numFac) factor2numeric(numFac) class(factor2numeric(numFac)) ## relevels relevels(numFac, c("3.5", "1", "2.5")) relevels(numFac, c("3.5", "2.5"))
Interpolate RTCA data
interpolationTransform(object, interval=0.01, method=c("linear","constant","fmm","periodic","natural", "monoH.FC"))
interpolationTransform(object, interval=0.01, method=c("linear","constant","fmm","periodic","natural", "monoH.FC"))
object |
An |
... |
other parameters, |
interval |
numeric, the interval between interpolated points, set to 0.01 by default |
method |
character, specifying the method for interpolation,
“linear” by default (for linear interpolation). Allowed options are: “linear” and
“constant” for |
Since most RTCA experiements record the experiments in the irregular
time-series, sometimes however it is desired to have regular
intervals. interpolationTransform
interpolate between data
points to estimate results of regular intervals.
Two classes of interpolations are supported by now: linear (using
approx
) and cubic spline
(spline
) interpolation. By default linear
interpolation is used.
An interpolated object of RTCA
.
Jitao David Zhang [email protected]
rgrTransform
stands for relative growth rate
transformation, ratioTransform
for ratio normalization
adopted by Roche commercial software. smoothTransform
to
smooth the RTCA readout.
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xInter <- interpolationTransform(x)
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xInter <- interpolationTransform(x)
Get index for the nearest time point to the given one. Called internally in many time-point related functions.
nearestTimeIndex(rtca, time)
nearestTimeIndex(rtca, time)
rtca |
An object of |
time |
numeric, a time point |
The function finds the time point with minimum absolute difference to the given time and returns its index.
An integer, the index of the nearest time point
Jitao David Zhang [email protected]
timepoints
to return all time points of an
RTCA
object.
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) x xIndex <- nearestTimeIndex(x, 25) timepoints(x)[xIndex]
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) x xIndex <- nearestTimeIndex(x, 25) timepoints(x)[xIndex]
The function parses RTCA output file into RTCA
object
parseRTCA(file, dec = ".", phenoData, maskWell, ...)
parseRTCA(file, dec = ".", phenoData, maskWell, ...)
file |
character, name of the RTCA output file |
dec |
decimal sign of the file |
phenoData |
phenoData |
maskWell |
character, either names or regular expression pattern(s) for well(s) to mask |
... |
other parameters passed to |
A csv-like format file can be exported from the RTCA device, which can
be fed into this function to set up an instance of
RTCA
object.
In the /extdata/ directory of the package, such a file is provided as an example. The first line contains the experiment ID, which is followed by a matrix of recorded data in the tabular form. The first and second column records the time-interval in the unit of hour and hour-minute-second format respectively. The rest columns then record the read-out (‘Cell-Index’, or ‘CI’) of the device, with each well a role.
phenoData
allows user to annotate the wells.Its usage mimicks
the ExpressionSet
object in the Biobase
package.
maskWell
allows to mask wells in case, for example, they are known to be
contaminated. The values can be either a vector of well names, or a
regular expression pattern for wells to be masked. To learn
regular expression patterns see grep
.
An object of RTCA-class
Jitao David Zhang [email protected]
http://www.roche-applied-science.com/proddata/gpip/3_8_9_1_1_1.html
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") pfile <- system.file("extdata/testOutputPhenoData.csv", package="RTCA") pData <- read.csv(pfile, sep="\t", row.names="Well") metaData <- data.frame(labelDescription=c( "Rack number", "siRNA catalogue number", "siRNA gene symbol", "siRNA EntrezGene ID", "siRNA targeting accession" )) phData <- new("AnnotatedDataFrame", data=pData, varMetadata=metaData) x <- parseRTCA(ofile, phenoData=phData) print(x) ## mask wells, e.g. due to unusual values x.skip <- parseRTCA(ofile, phenoData=phData, maskWell=c("D09")) x.skip.multiWells <- parseRTCA(ofile, phenoData=phData, maskWell=c("A01", "B01", "C02")) ## skip the last row x.skip.pattern <- parseRTCA(ofile, phenoData=phData, maskWell=c("H[0-9]{2}")) ## check the number of masked wells noMasked <- function(x) sum(apply(x, 2, function(x) all(is.na(x)))) noMasked(exprs(x)) noMasked(exprs(x.skip)) noMasked(exprs(x.skip.multiWells)) noMasked(exprs(x.skip.pattern))
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") pfile <- system.file("extdata/testOutputPhenoData.csv", package="RTCA") pData <- read.csv(pfile, sep="\t", row.names="Well") metaData <- data.frame(labelDescription=c( "Rack number", "siRNA catalogue number", "siRNA gene symbol", "siRNA EntrezGene ID", "siRNA targeting accession" )) phData <- new("AnnotatedDataFrame", data=pData, varMetadata=metaData) x <- parseRTCA(ofile, phenoData=phData) print(x) ## mask wells, e.g. due to unusual values x.skip <- parseRTCA(ofile, phenoData=phData, maskWell=c("D09")) x.skip.multiWells <- parseRTCA(ofile, phenoData=phData, maskWell=c("A01", "B01", "C02")) ## skip the last row x.skip.pattern <- parseRTCA(ofile, phenoData=phData, maskWell=c("H[0-9]{2}")) ## check the number of masked wells noMasked <- function(x) sum(apply(x, 2, function(x) all(is.na(x)))) noMasked(exprs(x)) noMasked(exprs(x.skip)) noMasked(exprs(x.skip.multiWells)) noMasked(exprs(x.skip.pattern))
Plots a E-plate in RTCA assays in one plot to convey an overview of the plate
plateView(rtca, ylim, titles,...)
plateView(rtca, ylim, titles,...)
rtca |
An object of |
ylim |
ylab lim |
titles |
Titles of sub-figures representing each well. If
missing, the function seeks whether a Well column is
available in the pData of the RTCA object, and if so, its value will
be used. If not, the sample names (by |
... |
Other parameters passed to the |
For now the function only supports the visualization of a 96-well E-plate.
The plate view plot draws lines indicating cell index (or its
transformations) in a birdview. When ...
are not specified,
default color, line style and width are used. col
,lty
and lwd
can be a vector, and if needed they will be expanded
to have the same length as wells.
NULL
, the function is called for the side effect
Jitao David Zhang [email protected]
RTCA
for data structure, plot
for
the basic plot function.
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") rtca <- parseRTCA(ofile) ## Not run automatically, because of 'margin too large' ## plateView(rtca) ## plateView(rtca, lty=2) ## plateView(rtca, col=rep(1:8, each=12)) rtca.skip <- parseRTCA(ofile, maskWell="H[0-9]{2}") ## plateView(rtca.skip)
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") rtca <- parseRTCA(ofile) ## Not run automatically, because of 'margin too large' ## plateView(rtca) ## plateView(rtca, lty=2) ## plateView(rtca, col=rep(1:8, each=12)) rtca.skip <- parseRTCA(ofile, maskWell="H[0-9]{2}") ## plateView(rtca.skip)
Plot the mean and deviation of rows/columns of a RTCA E-plate, to provide hints of potential row/column effect of the plate
plotGridEffect(rtca, mode = c("column", "row"), xlab = "time point", ylab = "readout", legend = TRUE, col, ...)
plotGridEffect(rtca, mode = c("column", "row"), xlab = "time point", ylab = "readout", legend = TRUE, col, ...)
rtca |
An object of |
mode |
character, either “column” or “row”, to choose which effect to depict |
xlab |
x-axis label |
ylab |
y-axis label |
legend |
logical, whether the legend should be added |
col |
Color of the curves |
... |
Further parameters passed to |
The error bars depicts the standard deviations
NULL
, the funciton is called for its side effect
Jitao David Zhang
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) plotGridEffect(x)
require(RTCA) ofile <- system.file("extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) plotGridEffect(x)
Performs ratio transformation (normalisation) of RTCA data, as recommended by the producer Roche.
ratioTransform(object, time)
ratioTransform(object, time)
object |
An object of |
time |
numeric, the time point used to normalize the whole series of data |
The xCelligence software provided by Roche performs ratio transform implicitly by dividing the time-series impedance measurement by the value of a selected time point (so-called 'base-time'), for instance 5 hours after compound transfection, in each cell. The aim of this transformation was to scale (normalize) the data of different wells, since the normalized values of all wells are uniformly 1 at the base-time.
However, this method is vulnerable to arbitrary selection of the time point chosen to normalize. It may be helpful to try several base-time values before comparing normalized results.
See derivativeTransform
and rgrTransform
for other normalization (scaling) possibilities.
An object of RTCA
, populated with normalized
value. The normalized values of all wells are uniformly 1 at the base-time.
Jitao David Zhang [email protected]
smoothTransform
and interpolationTransform
for smoothing and interpolating the RTCA
data. rgrTransform
calculates relative growth rate,
derivativeTransform
calculates derivatitve. The later two
methods are not sensative to the selection of base-time point.
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xNorm <- ratioTransform(x, 35)
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xNorm <- ratioTransform(x, 35)
Transform RTCA data into relative growth rate
rgrTransform(object, smooth)
rgrTransform(object, smooth)
object |
An object of |
smooth |
logical, should the object be smooth transformed after
the |
TODO: relative growth rate
An object of RTCA
populated with relative growth
rate instead of input data
Jitao David Zhang <[email protected]>
TODO: reference
derivativeTransform
for first derivative. ratioTransform
for ratio normalization
adopted by Roche commercial software. smoothTransform
and interpolationTransform
for other transformation possibilities.
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xRgr <- rgrTransform(x)
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xRgr <- rgrTransform(x)
RTCA object
Objects can be created by calls of the form new("RTCA", assayData,
phenoData, featureData, experimentData, annotation, exprs,
...)
. However, it is more common to be constructed by
parseRTCA
function by reading in RTCA output data directly.
expID
:Object of class "character"
, experiment ID
timeline
:Object of class "RTCAtimeline"
,
recording action track along the time line
assayData
:Object of class "AssayData"
, assay
data inherited from ExpressionSet-class
phenoData
:Object of class
"AnnotatedDataFrame"
, pheno data of the assay, annotating
the wells
featureData
:Object of class
"AnnotatedDataFrame"
, feature data of the assay, preserved
for time-line recording by the package
experimentData
:Object of class "MIAME"
, idle
annotation
:Object of class "character"
, idle
.__classVersion__
:Object of class "Versions"
,idle
Class ExpressionSet-class
, directly.
Class eSet-class
, by class "ExpressionSet", distance 2.
Class VersionedBiobase-class
, by class "ExpressionSet", distance 3.
Class Versioned-class
, by class "ExpressionSet", distance 4.
signature(object = "RTCA", time = "numeric",
action = "character")
: add action at the specified time, passed
to the RTCAtimeline
slot
signature(object = "RTCA", time =
"numeric")
: get action at the specified time, passed
to the RTCAtimeline
slot
signature(x = "RTCA")
: plot RTCA
signature(object = "RTCA", time = "numeric")
:
remove action at the specified time, passed
to the RTCAtimeline
slot
signature(object = "RTCA")
: print method
codesignature(object = "RTCA"): get Experiment ID
codesignature(object = "RTCA", value = "ANY"): set Experiment ID
signature(x = "RTCA")
: deprecated
signature(object = "RTCA")
: get the RTCAtimeline
slot
signature(object = "RTCA")
: assign the
RTCAtimeline
slot
signature(object = "RTCA")
: get the
recording time points in a vector
signature(object = "RTCA")
: assign the
recording time points
signature(object = "RTCA", time =
"numeric", action = "character")
: update the action at the specified time, passed
to the RTCAtimeline
slot
signature(x = "RTCA", y)
: plot the RTCA running
plot with matplot
. y
is interpretated as the
indices of the columns to be plotted, and will be expanded to all
the columns in case it is missing.
Jitao David Zhang [email protected]
https://www.roche-applied-science.com/sis/xcelligence/index.jsp?id=xcect_000000 introduces xCelligence system.
http://www.roche-applied-science.com/proddata/gpip/3_8_9_1_1_1.html for brief introduction into RTCA
new("RTCA", expID="testExp01")
new("RTCA", expID="testExp01")
Time line of actions performed by the xCelligence device, supporting CRUD manipulations (create, read, update and delete).
Objects can be created by calls of the form
new("RTCAtimeline")
. However, it is more common to be called
implicitly by creating an instance of RTCA
object.
actionTrack
:Object of class "data.frame"
,
records action track in the form of two-column
data.frame
. The two columns must have the names
‘time’ and ‘action’.
timeUnit
:Object of class "character"
,
recording the unit of time points stored in the actionTrack
slot.
startTime
:Object of class "POSIXct"
, the
absolute time when the measurement started (at the time point
‘0’)
signature(object = "RTCAtimeline", time =
"numeric", action = "character")
: add action at the specified time
signature(object = "RTCAtimeline")
: get
the action track in the form of data.frame
signature(object = "RTCAtimeline", value
= "data.frame")
: assign the action track
signature(object = "RTCAtimeline", time = "numeric")
: get
action at the specified time
signature(object = "RTCAtimeline")
: order
the action track by the time
signature(object = "RTCAtimeline")
: undo all
editing of the object and reset it to the initial state
signature(object = "RTCAtimeline", time = "numeric")
: remove
the action at the specified time
signature(object = "RTCAtimeline")
: return the
time unit used by the actiont track
signature(object = "RTCAtimeline", value = "character")
: assign
the time unit used by the actiont track
signature(object = "RTCAtimeline")
: return
the starting POSIXct time of the experiment
signature(object = "RTCAtimeline", value =
"character")
: assign
the starting POSIXct time of the experiment
Jitao David Zhang [email protected]
http://www.xcelligence.roche.com/ introduces xCelligence system.
http://www.roche-applied-science.com/proddata/gpip/3_8_9_1_1_1.html for brief introduction into RTCA
tl <- new("RTCAtimeline") show(tl)
tl <- new("RTCAtimeline") show(tl)
Subset (slice) RTCA object with starting- and ending-time
sliceRTCA(x, start, end)
sliceRTCA(x, start, end)
x |
An object of |
start |
numeric, start time |
end |
numeric, end time |
In case the exact starting- or ending-time is not matched, the nearest time point will be used to subset.
An object of RTCA
Jitao David Zhang [email protected]
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) subx <- sliceRTCA(x, 20, 50)
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) subx <- sliceRTCA(x, 20, 50)
Smoothing the RTCA cell impedance measurement
smoothTransform(object, ...)
smoothTransform(object, ...)
object |
An object of |
... |
Parameters passed to |
smoothTransform
smooths the RTCA cell impedance measurement by
calling the function smooth.spline
. This feature
can be useful for visualiation purposes and in conjuction with other transformations.
An RTCA
object populated with smoothed values
ratioTransform
performs ratio transformation recommended
by the machine provider. interpolationTransform
for interpolating the RTCA
data. derivativeTransform
returns cell impedance change
rates and rgrTransform
calculates relative growth rate.
Jitao David Zhang [email protected]
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xSmooth <- smoothTransform(x)
require(RTCA) ofile <- system.file("/extdata/testOutput.csv", package="RTCA") x <- parseRTCA(ofile) xSmooth <- smoothTransform(x)
Import output files from Spectramax spectrophotometer (plate reader) into the list format compatible with the cellHTS2 package.
spectramaxImport(file, encoding="latin1")
spectramaxImport(file, encoding="latin1")
file |
A Spectramax file |
encoding |
File character encoding, by default “latin1” |
The function imports output files from Spectramax plate reader, with which single-channel cell-based assays could be performed. Such assay includes WST-1 viability assay, which can be used to validate RTCA assay results.
A list of two items: one data frame (no name) and one character vector (txt). The data frame contains following columns:
well |
Well indices ([A-Z][0-9][0-9] format) on the microtitre plate |
val |
Value of each well |
The character vector txt contains a copy of the file contents.
Jitao David Zhang <[email protected]>
cellHTS2 package documentation.
wstFiles <- dir(system.file("extdata", package="RTCA"), pattern="^WST.*csv$", full.names=TRUE) spectramaxImport(wstFiles[1]) ## NOT RUN ## spectramaxImport also supports multiple files, in which case the ## result is a list of individual lists spectramaxImport(wstFiles) ## END NOT RUN
wstFiles <- dir(system.file("extdata", package="RTCA"), pattern="^WST.*csv$", full.names=TRUE) spectramaxImport(wstFiles[1]) ## NOT RUN ## spectramaxImport also supports multiple files, in which case the ## result is a list of individual lists spectramaxImport(wstFiles) ## END NOT RUN