Title: | Classes and Functions to Serve as the Basis for Other 'Gx' Packages |
---|---|
Description: | A collection of functions and classes which serve as the foundation for our lab's suite of R packages, such as 'PharmacoGx' and 'RadioGx'. This package was created to abstract shared functionality from other lab package releases to increase ease of maintainability and reduce code repetition in current and future 'Gx' suite programs. Major features include a 'CoreSet' class, from which 'RadioSet' and 'PharmacoSet' are derived, along with get and set methods for each respective slot. Additional functions related to fitting and plotting dose response curves, quantifying statistical correlation and calculating area under the curve (AUC) or survival fraction (SF) are included. For more details please see the included documentation, as well as: Smirnov, P., Safikhani, Z., El-Hachem, N., Wang, D., She, A., Olsen, C., Freeman, M., Selby, H., Gendoo, D., Grossman, P., Beck, A., Aerts, H., Lupien, M., Goldenberg, A. (2015) <doi:10.1093/bioinformatics/btv723>. Manem, V., Labie, M., Smirnov, P., Kofia, V., Freeman, M., Koritzinksy, M., Abazeed, M., Haibe-Kains, B., Bratman, S. (2018) <doi:10.1101/449793>. |
Authors: | Jermiah Joseph [aut], Petr Smirnov [aut], Ian Smith [aut], Christopher Eeles [aut], Feifei Li [aut], Benjamin Haibe-Kains [aut, cre] |
Maintainer: | Benjamin Haibe-Kains <[email protected]> |
License: | GPL (>= 3) |
Version: | 2.11.0 |
Built: | 2024-10-30 05:21:22 UTC |
Source: | https://github.com/bioc/CoreGx |
This is used to pass through unevaluated R expressions into subset and
[
, where they will be evaluated in the correct context.
.(...)
.(...)
... |
|
call
An R call object containing the quoted expression.
.(sample_line1 == 'A2058')
.(sample_line1 == 'A2058')
Convert a LongTable assay into a BumpyMatrix object
.assayToBumpyMatrix(LT, assay, rows, cols, sparse = TRUE)
.assayToBumpyMatrix(LT, assay, rows, cols, sparse = TRUE)
LT |
|
assay |
|
rows |
|
cols |
|
sparse |
|
BumpyMatrix
containing the data from assay
.
stats::optim
L-BFGS-B with fall-back grid/pattern search
if convergence is not achieved.Function to fit curve via stats::optim
.fitCurve2( par, x, y, fn, loss, lower = -Inf, upper = Inf, precision = 1e-04, density = c(2, 10, 5), step = 0.5/density, ..., loss_args = list(), span = 1, optim_only = FALSE, control = list(factr = 1e-08, ndeps = rep(1e-04, times = length(par)), trace = 0) )
.fitCurve2( par, x, y, fn, loss, lower = -Inf, upper = Inf, precision = 1e-04, density = c(2, 10, 5), step = 0.5/density, ..., loss_args = list(), span = 1, optim_only = FALSE, control = list(factr = 1e-08, ndeps = rep(1e-04, times = length(par)), trace = 0) )
par |
|
x |
|
y |
|
fn |
|
loss |
|
lower |
|
upper |
|
precision |
|
density |
|
step |
initial step size for pattern search. |
... |
|
loss_args |
|
span |
|
optim_only |
|
control |
|
TODO
numeric
Vector of optimal parameters for fn
fit against y
on the values of x
.
## Not run: # Four parameter hill curve equation hillEqn <- function(x, Emin, Emax, EC50, lambda) { (Emin + Emax * (x / EC50)^lambda) / (1 + (x / EC50)^lambda) } # Make some dummy data doses <- rev(1000 / (2^(1:20))) lambda <- 1 Emin <- 1 Emax <- 0.1 EC50 <- median(doses) response <- hillEqn(doses, Emin=Emin, lambda=lambda, Emax=Emax, EC50=EC50) nresponse <- response + rnorm(length(response), sd=sd(response)*0.1) # add noise # 3-parameter optimization 3par <- .fitCurve2( par=c(Emax, EC50, lambda), x=doses, y=nresponse, fn=hillEqn, Emin=Emin, # set this as constant in the function being optimized (via ...) loss=.normal_loss, loss_args=list(trunc=FALSE, n=1, scale=0.07), upper=c(1, max(doses), 6), lower=c(0, min(doses), 0) ) # 2-parameter optimization 2par <- .fitCurve2( par=c(Emax, EC50), x=doses, y=nresponse, fn=hillEqn, Emin=Emin, # set this as constant in the function being optimized (via ...) lambda=1, loss=.normal_loss, loss_args=list(trunc=FALSE, n=1, scale=0.07), upper=c(1, max(doses)), lower=c(0, min(doses)) ) ## End(Not run)
## Not run: # Four parameter hill curve equation hillEqn <- function(x, Emin, Emax, EC50, lambda) { (Emin + Emax * (x / EC50)^lambda) / (1 + (x / EC50)^lambda) } # Make some dummy data doses <- rev(1000 / (2^(1:20))) lambda <- 1 Emin <- 1 Emax <- 0.1 EC50 <- median(doses) response <- hillEqn(doses, Emin=Emin, lambda=lambda, Emax=Emax, EC50=EC50) nresponse <- response + rnorm(length(response), sd=sd(response)*0.1) # add noise # 3-parameter optimization 3par <- .fitCurve2( par=c(Emax, EC50, lambda), x=doses, y=nresponse, fn=hillEqn, Emin=Emin, # set this as constant in the function being optimized (via ...) loss=.normal_loss, loss_args=list(trunc=FALSE, n=1, scale=0.07), upper=c(1, max(doses), 6), lower=c(0, min(doses), 0) ) # 2-parameter optimization 2par <- .fitCurve2( par=c(Emax, EC50), x=doses, y=nresponse, fn=hillEqn, Emin=Emin, # set this as constant in the function being optimized (via ...) lambda=1, loss=.normal_loss, loss_args=list(trunc=FALSE, n=1, scale=0.07), upper=c(1, max(doses)), lower=c(0, min(doses)) ) ## End(Not run)
Convert LongTable to gDR Style SummarizedExperiment
.longTableToSummarizedExperiment(LT, assay_names)
.longTableToSummarizedExperiment(LT, assay_names)
LT |
|
assay_names |
|
SummarizedExperiment
object with all assay from LT
as
BumpyMatrix
es.
Single bracket subsetting for a LongTable object. See subset for more details.
## S4 method for signature 'LongTable,ANY,ANY,ANY' x[i, j, assays = assayNames(x), ..., drop = FALSE]
## S4 method for signature 'LongTable,ANY,ANY,ANY' x[i, j, assays = assayNames(x), ..., drop = FALSE]
x |
|
i |
|
j |
|
assays |
|
... |
Included to ensure drop can only be set by name. |
drop |
|
This function is endomorphic, it always returns a LongTable object.
A LongTable
containing only the data specified in the function
parameters.
# Character merckLongTable['ABT-888', 'CAOV3'] # Numeric merckLongTable[1, c(1, 2)] # Logical merckLongTable[, colData(merckLongTable)$sampleid == 'A2058'] # Call merckLongTable[ .(drug1id == 'Dasatinib' & drug2id != '5-FU'), .(sampleid == 'A2058'), ]
# Character merckLongTable['ABT-888', 'CAOV3'] # Numeric merckLongTable[1, c(1, 2)] # Logical merckLongTable[, colData(merckLongTable)$sampleid == 'A2058'] # Call merckLongTable[ .(drug1id == 'Dasatinib' & drug2id != '5-FU'), .(sampleid == 'A2058'), ]
[[<-
Method for LongTable ClassJust a wrapper around assay<- for convenience. See
?'assay<-,LongTable,character-method'.
## S4 replacement method for signature 'LongTable,ANY,ANY' x[[i]] <- value
## S4 replacement method for signature 'LongTable,ANY,ANY' x[[i]] <- value
x |
A |
i |
The name of the assay to update, must be in |
value |
A |
A LongTable
object with the assay i
updated using value
.
merckLongTable[['sensitivity']] <- merckLongTable[['sensitivity']]
merckLongTable[['sensitivity']] <- merckLongTable[['sensitivity']]
Select an assay from a LongTable object
## S4 method for signature 'LongTable' x$name
## S4 method for signature 'LongTable' x$name
x |
A |
name |
|
data.frame
The assay object.
merckLongTable$sensitivity
merckLongTable$sensitivity
Update an assay from a LongTable object
## S4 replacement method for signature 'LongTable' x$name <- value
## S4 replacement method for signature 'LongTable' x$name <- value
x |
A |
name |
|
value |
A |
Updates the assay name
in x
with value
, returning an invisible
NULL.
merckLongTable$sensitivity <- merckLongTable$sensitivity
merckLongTable$sensitivity <- merckLongTable$sensitivity
data.table
object.Compute a group-by operation over a data.table
in a functional, pipe
compatible format.
## S4 method for signature 'data.table' aggregate( x, by, ..., subset = TRUE, nthread = 1, progress = TRUE, BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
## S4 method for signature 'data.table' aggregate( x, by, ..., subset = TRUE, nthread = 1, progress = TRUE, BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
x |
|
by |
|
... |
|
subset |
|
nthread |
|
progress |
|
BPPARAM |
|
enlist |
|
moreArgs |
|
This S4 method override the default aggregate
method for a data.frame
and as such you need to call aggregate.data.frame
directly to get the
original S3 method for a data.table
.
Arguments in ...
are substituted and wrapped in a list, which is passed
through to the j argument of [.data.table
internally. The function currently
tries to build informative column names for unnamed arguments in ...
by
appending the name of each function call with the name of its first argument,
which is assumed to be the column name being aggregated over. If an argument
to ...
is named, that will be the column name of its value in the resulting
data.table
.
The primary use case for enlist=FALSE
is to allow computation of dependent
aggregations, where the output from a previous aggregation is required in a
subsequent one. For this case, wrap your call in {
and assign intermediate
results to variables, returning the final results as a list where each list
item will become a column in the final table with the corresponding name.
Name inference is disabled for this case, since it is assumed you will name
the returned list items appropriately.
A major advantage over multiple calls to aggregate
is that
the overhead of parallelization is paid only once even for complex multi-step
computations like fitting a model, capturing its paramters, and making
predictions using it. It also allows capturing arbitrarily complex calls
which can be recomputed later using the
update,TreatmentResponseExperiment-method
A potential disadvantage is increased RAM usage per
thread due to storing intermediate values in variables, as well as any
memory allocation overhead associate therewith.
data.table
of aggregated results with an aggregations
attribute
capturing metadata about the last aggregation performed on the table.
LongTable
or inheriting classCompute a group-by operation over a LongTable
object or it's inhering
classes.
## S4 method for signature 'LongTable' aggregate( x, assay, by, ..., subset = TRUE, nthread = 1, progress = TRUE, BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
## S4 method for signature 'LongTable' aggregate( x, assay, by, ..., subset = TRUE, nthread = 1, progress = TRUE, BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
x |
|
assay |
|
by |
|
... |
|
subset |
|
nthread |
|
progress |
|
BPPARAM |
|
enlist |
|
moreArgs |
|
Arguments in ...
are substituted and wrapped in a list, which is passed
through to the j argument of [.data.table
internally. The function currently
tries to build informative column names for unnamed arguments in ...
by
appending the name of each function call with the name of its first argument,
which is assumed to be the column name being aggregated over. If an argument
to ...
is named, that will be the column name of its value in the resulting
data.table
.
The primary use case for enlist=FALSE
is to allow computation of dependent
aggregations, where the output from a previous aggregation is required in a
subsequent one. For this case, wrap your call in {
and assign intermediate
results to variables, returning the final results as a list where each list
item will become a column in the final table with the corresponding name.
Name inference is disabled for this case, since it is assumed you will name
the returned list items appropriately.
A major advantage over multiple calls to aggregate
is that
the overhead of parallelization is paid only once even for complex multi-step
computations like fitting a model, capturing its paramters, and making
predictions using it. It also allows capturing arbitrarily complex calls
which can be recomputed later using the
update,TreatmentResponseExperiment-method
A potential disadvantage is increased RAM usage per
thread due to storing intermediate values in variables, as well as any
memory allocation overhead associate therewith.
data.table
of aggregation results.
data.table::[.data.table
, BiocParallel::bplapply
Functional API for data.table aggregation which allows capture of associated aggregate calls so they can be recomputed later.
aggregate2( x, by, ..., nthread = 1, progress = interactive(), BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
aggregate2( x, by, ..., nthread = 1, progress = interactive(), BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
x |
|
by |
|
... |
|
nthread |
|
progress |
|
BPPARAM |
|
enlist |
|
moreArgs |
|
Arguments in ...
are substituted and wrapped in a list, which is passed
through to the j argument of [.data.table
internally. The function currently
tries to build informative column names for unnamed arguments in ...
by
appending the name of each function call with the name of its first argument,
which is assumed to be the column name being aggregated over. If an argument
to ...
is named, that will be the column name of its value in the resulting
data.table
.
The primary use case for enlist=FALSE
is to allow computation of dependent
aggregations, where the output from a previous aggregation is required in a
subsequent one. For this case, wrap your call in {
and assign intermediate
results to variables, returning the final results as a list where each list
item will become a column in the final table with the corresponding name.
Name inference is disabled for this case, since it is assumed you will name
the returned list items appropriately.
A major advantage over multiple calls to aggregate
is that
the overhead of parallelization is paid only once even for complex multi-step
computations like fitting a model, capturing its paramters, and making
predictions using it. It also allows capturing arbitrarily complex calls
which can be recomputed later using the
update,TreatmentResponseExperiment-method
A potential disadvantage is increased RAM usage per
thread due to storing intermediate values in variables, as well as any
memory allocation overhead associate therewith.
data.table
of aggregation results.
data.table::[.data.table
, BiocParallel::bplapply
This function calculates an Adaptive Matthews Correlation Coefficient (AMCC) for two vectors of values of the same length. It assumes the entries in the two vectors are paired. The Adaptive Matthews Correlation Coefficient for two vectors of values is defined as the Maximum Matthews Coefficient over all possible binary splits of the ranks of the two vectors. In this way, it calculates the best possible agreement of a binary classifier on the two vectors of data. If the AMCC is low, then it is impossible to find any binary classification of the two vectors with a high degree of concordance.
amcc(x, y, step.prct = 0, min.cat = 3, nperm = 1000, nthread = 1, ...)
amcc(x, y, step.prct = 0, min.cat = 3, nperm = 1000, nthread = 1, ...)
x , y
|
Two paired vectors of values. Could be replicates of observations for the same experiments for example. |
step.prct |
Instead of testing all possible splits of the data, it is possible to test steps of a percentage size of the total number of ranks in x/y. If this variable is 0, function defaults to testing all possible splits. |
min.cat |
The minimum number of members per category. Classifications with less members fitting into both categories will not be considered. |
nperm |
The number of perumatation to use for estimating significance. If 0, then no p-value is calculated. |
nthread |
Number of threads to parallize over. Both the AMCC calculation and the permutation testing is done in parallel. |
... |
Additional arguments |
Returns a list with two elements. $amcc contains the highest 'mcc' value over all the splits, the p value, as well as the rank at which the split was done.
x <- c(1,2,3,4,5,6,7) y <- c(1,3,5,4,2,7,6) amcc(x,y, min.cat=2)
x <- c(1,2,3,4,5,6,7) y <- c(1,3,5,4,2,7,6) amcc(x,y, min.cat=2)
LongTable
to a TreatmentResponseExperiment
Coerce a LongTable into a data.table
.
Currently only supports coercing to data.table or data.frame
Coerce a data.table with the proper configuration attributes back to a LongTable
from |
A |
The data in object
, as the child-class
TreatmentResponseExperiment
.
A data.table
with the data from a LongTable.
data.table
containing the data from the LongTable, with the
‘longTableDataMapper’ attribute containg the metadata needed to reverse
the coercing operation.
LongTable
object configured with the longTableDataMapper
data.table
with long format of data in from
data.frame
with long format of data in from
.
SummarizedExperiment
with each assay as a BumpyMatrix
A TREDataMapper
object.
data(clevelandSmall_cSet) TRE <- as(treatmentResponse(clevelandSmall_cSet), "TreatmentResponseExperiment") TRE as(merckLongTable, 'data.table') dataTable <- as(merckLongTable, 'data.table') print(attr(dataTable, 'longTableDataMapper')) # Method doesn't work without this as(dataTable, 'LongTable') SE <- molecularProfilesSlot(clevelandSmall_cSet)[[1]] as(SE, 'data.table') SE <- molecularProfilesSlot(clevelandSmall_cSet)[[1]] as(SE, 'data.frame')
data(clevelandSmall_cSet) TRE <- as(treatmentResponse(clevelandSmall_cSet), "TreatmentResponseExperiment") TRE as(merckLongTable, 'data.table') dataTable <- as(merckLongTable, 'data.table') print(attr(dataTable, 'longTableDataMapper')) # Method doesn't work without this as(dataTable, 'LongTable') SE <- molecularProfilesSlot(clevelandSmall_cSet)[[1]] as(SE, 'data.table') SE <- molecularProfilesSlot(clevelandSmall_cSet)[[1]] as(SE, 'data.frame')
Coerce a data.table with the proper configuration attributes back to a LongTable
as.long.table(x)
as.long.table(x)
x |
A |
LongTable
object configured with the longTableDataMapper
dataTable <- as(merckLongTable, 'data.table') print(attr(dataTable, 'longTableDataMapper')) # Method doesn't work without this as.long.table(dataTable)
dataTable <- as(merckLongTable, 'data.table') print(attr(dataTable, 'longTableDataMapper')) # Method doesn't work without this as.long.table(dataTable)
Generic to access the assay columns of a rectangular object.
assayCols(object, ...)
assayCols(object, ...)
object |
|
... |
Allow new arguments to this generic. |
Depends on the implemented method.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Retrieve and assayIndex
assayIndex(x, ...)
assayIndex(x, ...)
x |
An |
... |
|
An object representing the "assayIndex" of an S4
object.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Retrieve a set of assayKeys
assayKeys(x, ...)
assayKeys(x, ...)
x |
An |
... |
|
An object representing the "assayKeys" of an S4
object.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Prevents modification of objects labelled with the "immutable" S3-class by intercepting assignment during S3-method dispatch and returning an error.
\method{subset}{immutable}(object, ...) <- value ## S3 replacement method for class 'immutable' object[...] <- value ## S3 replacement method for class 'immutable' object[[...]] <- value ## S3 replacement method for class 'immutable' object$... <- value ## S3 replacement method for class 'immutable' names(x) <- value ## S3 replacement method for class 'immutable' dimnames(x) <- value \method{colnames}{immutable}(x) <- value \method{rownames}{immutable}(x) <- value
\method{subset}{immutable}(object, ...) <- value ## S3 replacement method for class 'immutable' object[...] <- value ## S3 replacement method for class 'immutable' object[[...]] <- value ## S3 replacement method for class 'immutable' object$... <- value ## S3 replacement method for class 'immutable' names(x) <- value ## S3 replacement method for class 'immutable' dimnames(x) <- value \method{colnames}{immutable}(x) <- value \method{rownames}{immutable}(x) <- value
object , x
|
An R object inherting from the "immutable" S3-class. |
... |
Catch subset arguments for various dimensions. |
value |
Not used. |
None, throws an error.
immutable_df <- immutable(data.frame(a=1:5, b=letters[1:5])) # return immutable data.frame immutable_df[1:4, ] # return immutable vector immutable_df$a
immutable_df <- immutable(data.frame(a=1:5, b=letters[1:5])) # return immutable data.frame immutable_df[1:4, ] # return immutable vector immutable_df$a
S4
object.Build an assay table with an S4
object.
buildComboProfiles(object, ...)
buildComboProfiles(object, ...)
object |
|
... |
Allow new arguments to be defined for this generic. |
data.table
.
"This is a generic method!"
"This is a generic method!"
Build an assay table with selected assay profiles for drug combinations
## S4 method for signature 'LongTable' buildComboProfiles(object, profiles)
## S4 method for signature 'LongTable' buildComboProfiles(object, profiles)
object |
|
profiles |
|
A data.table
containing fields
treatment1id
, treatment1dose
, treatment2id
, treatment2dose
, sampleid
,
which are used as keys to keep track of profiles,
along with columns of selected profiles from their assays.
Each *_1
is the monothearpy profile of treatment 1 in the combination,
and the same rule applies to treatment 2.
## Not run: combo_profile_1 <- buildComboProfiles(tre, c("auc", "SCORE")) combo_profile_2 <- buildComboProfiles(tre, c("HS", "EC50", "E_inf", "ZIP")) ## End(Not run)
## Not run: combo_profile_1 <- buildComboProfiles(tre, c("auc", "SCORE")) combo_profile_2 <- buildComboProfiles(tre, c("HS", "EC50", "E_inf", "ZIP")) ## End(Not run)
Build a LongTable object
buildLongTable(from, ...)
buildLongTable(from, ...)
from |
What to build the LongTable from? |
... |
|
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
LongTable Create a LongTable object from a single .csv file
## S4 method for signature 'character' buildLongTable(from, rowDataCols, colDataCols, assayCols)
## S4 method for signature 'character' buildLongTable(from, rowDataCols, colDataCols, assayCols)
from |
|
rowDataCols |
|
colDataCols |
|
assayCols |
|
A LongTable
object containing one or more assays, indexed by
rowID and colID.
Create a LongTable object from a single data.table or data.frame object.
## S4 method for signature 'data.frame' buildLongTable(from, rowDataCols, colDataCols, assayCols)
## S4 method for signature 'data.frame' buildLongTable(from, rowDataCols, colDataCols, assayCols)
from |
|
rowDataCols |
|
colDataCols |
|
assayCols |
|
A LongTable
object containing one or more assays, indexed by
rowID and colID.
Create a LongTable object from a list containing file paths, data.frames and data.tables.
## S4 method for signature 'list' buildLongTable(from, rowDataCols, colDataCols, assayCols)
## S4 method for signature 'list' buildLongTable(from, rowDataCols, colDataCols, assayCols)
from |
|
rowDataCols |
|
colDataCols |
|
assayCols |
|
A LongTable
object constructed with the data in from
.
## Not run: assayList <- assays(merckLongTable, withDimnames=TRUE) rowDataCols <- list(rowIDs(merckLongTable), rowMeta(merckLongTable)) colDataCols <- list(colIDs(merckLongTable), colMeta(merckLongTable)) assayCols <- assayCols(merckLongTable) longTable <- buildLongTable(from=assayList, rowDataCols, colDataCols, assayCols) ## End(Not run)
## Not run: assayList <- assays(merckLongTable, withDimnames=TRUE) rowDataCols <- list(rowIDs(merckLongTable), rowMeta(merckLongTable)) colDataCols <- list(colIDs(merckLongTable), colMeta(merckLongTable)) assayCols <- assayCols(merckLongTable) longTable <- buildLongTable(from=assayList, rowDataCols, colDataCols, assayCols) ## End(Not run)
Ensures that c
and append
to an "immutable" class object return an
immutable class object.
## S3 method for class 'immutable' c(x, ...)
## S3 method for class 'immutable' c(x, ...)
x |
An R object inheriting from the "immutable" S3-clas |
... |
Objects to concatenate to |
x with one or more values appended to it.
cardinality
relationships between a group
of columns (your identifiers) and all other columns.Search a data.frame for 1:cardinality
relationships between a group
of columns (your identifiers) and all other columns.
checkColumnCardinality(df, group, cardinality = 1, ...)
checkColumnCardinality(df, group, cardinality = 1, ...)
df |
A |
group |
A |
cardinality |
The cardinality of to search for (i.e., 1: |
... |
Fall through arguments to data.table:: |
A character
vector with the names of the columns with
cardinality of 1:cardinality
with the columns listed in group
.
df <- rawdata(exampleDataMapper) checkColumnCardinality(df, group='treatmentid')
df <- rawdata(exampleDataMapper) checkColumnCardinality(df, group='treatmentid')
This function checks the structure of a PharamcoSet, ensuring that the correct annotations are in place and all the required slots are filled so that matching of samples and drugs can be properly done across different types of data and with other studies.
checkCsetStructure(object, plotDist = FALSE, result.dir = tempdir())
checkCsetStructure(object, plotDist = FALSE, result.dir = tempdir())
object |
A |
plotDist |
Should the function also plot the distribution of molecular data? |
result.dir |
The path to the directory for saving the plots as a string.
Defaults to this R sessions |
Prints out messages whenever describing the errors found in the structure of the cSet object passed in.
checkCsetStructure(clevelandSmall_cSet)
checkCsetStructure(clevelandSmall_cSet)
This dataset is just a dummy object derived from the Cleveland_mut RadioSet in the RadioGx R package. It's contents should not be interpreted and it is only present to test the functions in this package and provide examples
data(clevelandSmall_cSet)
data(clevelandSmall_cSet)
CoreSet object
Lamb et al. The Connectivity Map: using gene-expression signatures to connect small molecules, genes, and disease. Science, 2006.
colData
out of the rawdata
slot using
the assigned colDataMap
metadata.Convenience method to subset the colData
out of the rawdata
slot using
the assigned colDataMap
metadata.
## S4 method for signature 'LongTableDataMapper' colData(x, key = TRUE)
## S4 method for signature 'LongTableDataMapper' colData(x, key = TRUE)
x |
|
key |
|
data.table
The colData
as specified in the colDataMap
slot.
colData
out of the rawdata
slot using
the assigned colDataMap
metadata.Convenience method to subset the colData
out of the rawdata
slot using
the assigned colDataMap
metadata.
## S4 method for signature 'TREDataMapper' colData(x, key = TRUE)
## S4 method for signature 'TREDataMapper' colData(x, key = TRUE)
x |
|
key |
|
data.table
The colData
as specified in the colDataMap
slot.
Generic to access the row identifiers for an object.
colIDs(object, ...)
colIDs(object, ...)
object |
|
... |
ALlow new arguments to this generic |
Depends on the implemented method.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Useful for converting a regular function into a function amenable to
optimization via stats::optim
, which requires all free parameters be
passed as a single vector par
.
collect_fn_params(fn)
collect_fn_params(fn)
fn |
|
Takes a function of the form f(x, ...), where ... is any number of additional
function parameters (bot not literal ...
!) and parses it to a function of
the form f(par, x) where par
is a vector of values for ... in
the same order as the arguments appear in fn
.
function
A new non-primitive function where the first argument is
par
, which takes a vector of parameters being optimized, and the
second argument is the old first argument to fn
(usually x
since this
is the independent variable to optimize the function over).
Generic to access the column identifiers for a rectangular object.
colMeta(object, ...)
colMeta(object, ...)
object |
|
... |
ALlow new arguments to this generic |
Depends on impemented method.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
A function for finding the connectivity between two signatures, using either the GSEA method based on the KS statistic, or the gwc method based on a weighted spearman statistic. The GSEA analysis is implemented in the piano package.
connectivityScore( x, y, method = c("fgsea", "gwc"), nperm = 10000, nthread = 1, gwc.method = c("spearman", "pearson"), ... )
connectivityScore( x, y, method = c("fgsea", "gwc"), nperm = 10000, nthread = 1, gwc.method = c("spearman", "pearson"), ... )
x |
A |
y |
A |
method |
|
nperm |
|
nthread |
|
gwc.method |
|
... |
Additional arguments passed down to gsea and gwc functions |
numeric
a numeric vector with the score and the p-value associated
with it
F. Pozzi, T. Di Matteo, T. Aste, 'Exponential smoothing weighted correlations', The European Physical Journal B, Vol. 85, No 6, 2012. DOI: 10.1140/epjb/e2012-20697-x
Varemo, L., Nielsen, J. and Nookaew, I. (2013) Enriching the gene set analysis of genome-wide data by incorporating directionality of gene expression and combining statistical hypotheses and methods. Nucleic Acids Research. 41 (8), 4378-4391. doi: 10.1093/nar/gkt111
xValue <- c(1,5,23,4,8,9,2,19,11,12,13) xSig <- c(0.01, 0.001, .97, 0.01,0.01,0.28,0.7,0.01,0.01,0.01,0.01) yValue <- c(1,5,10,4,8,19,22,19,11,12,13) ySig <- c(0.01, 0.001, .97,0.01, 0.01,0.78,0.9,0.01,0.01,0.01,0.01) xx <- cbind(xValue, xSig) yy <- cbind(yValue, ySig) rownames(xx) <- rownames(yy) <- c('1','2','3','4','5','6','7','8','9','10','11') data.cor <- connectivityScore(xx,yy,method='gwc', gwc.method='spearman', nperm=300)
xValue <- c(1,5,23,4,8,9,2,19,11,12,13) xSig <- c(0.01, 0.001, .97, 0.01,0.01,0.28,0.7,0.01,0.01,0.01,0.01) yValue <- c(1,5,10,4,8,19,22,19,11,12,13) ySig <- c(0.01, 0.001, .97,0.01, 0.01,0.78,0.9,0.01,0.01,0.01,0.01) xx <- cbind(xValue, xSig) yy <- cbind(yValue, ySig) rownames(xx) <- rownames(yy) <- c('1','2','3','4','5','6','7','8','9','10','11') data.cor <- connectivityScore(xx,yy,method='gwc', gwc.method='spearman', nperm=300)
deprecated
or defunct
methods in the CoreGx
R package.List of deprecated
or defunct
methods in the CoreGx
R package.
CoreSet
: The CoreSet
constructor is being updated to have a new API. This
API is currently available via the CoreSet2
constructor. In Bioconductor
3.16, the old constructor will be renamed CoreSet2
and the new constructor
will be renamed CoreSet
.
buildLongTable
: This function no longer works as building a LongTable
or
TreatmentResponseExperiment
now uses a DataMapper
and the metaConstruct
method. See vignette("LongTable")
for a detailed description of how
to create a LongTable object.
A constructor that simplifies the process of creating CoreSets, as well as creates empty objects for data not provided to the constructor. Only objects returned by this constructor are expected to work with the CoreSet methods.
CoreSet( name, molecularProfiles = list(), sample = data.frame(), sensitivityInfo = data.frame(), sensitivityRaw = array(dim = c(0, 0, 0)), sensitivityProfiles = matrix(), sensitivityN = matrix(nrow = 0, ncol = 0), perturbationN = array(NA, dim = c(0, 0, 0)), curationSample = data.frame(), curationTissue = data.frame(), curationTreatment = data.frame(), treatment = data.frame(), datasetType = c("sensitivity", "perturbation", "both"), verify = TRUE, ... )
CoreSet( name, molecularProfiles = list(), sample = data.frame(), sensitivityInfo = data.frame(), sensitivityRaw = array(dim = c(0, 0, 0)), sensitivityProfiles = matrix(), sensitivityN = matrix(nrow = 0, ncol = 0), perturbationN = array(NA, dim = c(0, 0, 0)), curationSample = data.frame(), curationTissue = data.frame(), curationTreatment = data.frame(), treatment = data.frame(), datasetType = c("sensitivity", "perturbation", "both"), verify = TRUE, ... )
name |
A |
molecularProfiles |
A |
sample |
A |
sensitivityInfo |
A |
sensitivityRaw |
A 3 Dimensional |
sensitivityProfiles |
|
sensitivityN , perturbationN
|
A |
curationSample , curationTissue , curationTreatment
|
A |
treatment |
A |
datasetType |
A |
verify |
|
... |
Catch and parse any renamed constructor arguments. |
Parameters to this function have been renamed!
cell is now sample
drug is now treatment
An object of class CoreSet
data(clevelandSmall_cSet) clevelandSmall_cSet
data(clevelandSmall_cSet) clevelandSmall_cSet
CoreSet
Documentation for the various setters and getters which allow manipulation
of data in the slots of a CoreSet
object.
## S4 method for signature 'CoreSet' annotation(object) ## S4 replacement method for signature 'CoreSet,list' annotation(object) <- value ## S4 method for signature 'CoreSet' dateCreated(object) ## S4 replacement method for signature 'CoreSet,character' dateCreated(object) <- value ## S4 method for signature 'CoreSet' name(object) ## S4 replacement method for signature 'CoreSet' name(object) <- value ## S4 method for signature 'CoreSet' sampleInfo(object) ## S4 replacement method for signature 'CoreSet,data.frame' sampleInfo(object) <- value ## S4 method for signature 'CoreSet' sampleNames(object) ## S4 replacement method for signature 'CoreSet,character' sampleNames(object) <- value ## S4 method for signature 'CoreSet' treatmentInfo(object) ## S4 replacement method for signature 'CoreSet,data.frame' treatmentInfo(object) <- value ## S4 method for signature 'CoreSet' treatmentNames(object) ## S4 replacement method for signature 'CoreSet,character' treatmentNames(object) <- value ## S4 method for signature 'CoreSet' curation(object) ## S4 replacement method for signature 'CoreSet,list' curation(object) <- value ## S4 method for signature 'CoreSet' datasetType(object) ## S4 replacement method for signature 'CoreSet,character' datasetType(object) <- value ## S4 method for signature 'CoreSet' molecularProfiles(object, mDataType, assay) ## S4 replacement method for signature 'CoreSet,character,character,matrix' molecularProfiles(object, mDataType, assay) <- value ## S4 replacement method for signature 'CoreSet,character,missing,matrix' molecularProfiles(object, mDataType, assay) <- value ## S4 replacement method for signature 'CoreSet,missing,missing,list_OR_MAE' molecularProfiles(object, mDataType, assay) <- value ## S4 method for signature 'CoreSet' featureInfo(object, mDataType) ## S4 replacement method for signature 'CoreSet,character,data.frame' featureInfo(object, mDataType) <- value ## S4 method for signature 'CoreSet,character' phenoInfo(object, mDataType) ## S4 replacement method for signature 'CoreSet,character,data.frame' phenoInfo(object, mDataType) <- value ## S4 method for signature 'CoreSet,character' fNames(object, mDataType) ## S4 replacement method for signature 'CoreSet,character,character' fNames(object, mDataType) <- value ## S4 method for signature 'CoreSet' mDataNames(object) ## S4 replacement method for signature 'CoreSet' mDataNames(object) <- value ## S4 method for signature 'CoreSet' molecularProfilesSlot(object) ## S4 replacement method for signature 'CoreSet,list_OR_MAE' molecularProfilesSlot(object) <- value ## S4 method for signature 'CoreSet' sensitivityInfo(object, dimension, ...) ## S4 replacement method for signature 'CoreSet,data.frame' sensitivityInfo(object, dimension, ...) <- value ## S4 method for signature 'CoreSet' sensitivityMeasures(object) ## S4 replacement method for signature 'CoreSet,character' sensitivityMeasures(object) <- value ## S4 method for signature 'CoreSet' sensitivityProfiles(object) ## S4 replacement method for signature 'CoreSet,data.frame' sensitivityProfiles(object) <- value ## S4 method for signature 'CoreSet' sensitivityRaw(object) ## S4 replacement method for signature 'CoreSet,array' sensitivityRaw(object) <- value ## S4 method for signature 'CoreSet' treatmentResponse(object) ## S4 replacement method for signature 'CoreSet,list_OR_LongTable' treatmentResponse(object) <- value ## S4 method for signature 'CoreSet' sensNumber(object) ## S4 replacement method for signature 'CoreSet,matrix' sensNumber(object) <- value ## S4 method for signature 'CoreSet' pertNumber(object) ## S4 replacement method for signature 'CoreSet,array' pertNumber(object) <- value
## S4 method for signature 'CoreSet' annotation(object) ## S4 replacement method for signature 'CoreSet,list' annotation(object) <- value ## S4 method for signature 'CoreSet' dateCreated(object) ## S4 replacement method for signature 'CoreSet,character' dateCreated(object) <- value ## S4 method for signature 'CoreSet' name(object) ## S4 replacement method for signature 'CoreSet' name(object) <- value ## S4 method for signature 'CoreSet' sampleInfo(object) ## S4 replacement method for signature 'CoreSet,data.frame' sampleInfo(object) <- value ## S4 method for signature 'CoreSet' sampleNames(object) ## S4 replacement method for signature 'CoreSet,character' sampleNames(object) <- value ## S4 method for signature 'CoreSet' treatmentInfo(object) ## S4 replacement method for signature 'CoreSet,data.frame' treatmentInfo(object) <- value ## S4 method for signature 'CoreSet' treatmentNames(object) ## S4 replacement method for signature 'CoreSet,character' treatmentNames(object) <- value ## S4 method for signature 'CoreSet' curation(object) ## S4 replacement method for signature 'CoreSet,list' curation(object) <- value ## S4 method for signature 'CoreSet' datasetType(object) ## S4 replacement method for signature 'CoreSet,character' datasetType(object) <- value ## S4 method for signature 'CoreSet' molecularProfiles(object, mDataType, assay) ## S4 replacement method for signature 'CoreSet,character,character,matrix' molecularProfiles(object, mDataType, assay) <- value ## S4 replacement method for signature 'CoreSet,character,missing,matrix' molecularProfiles(object, mDataType, assay) <- value ## S4 replacement method for signature 'CoreSet,missing,missing,list_OR_MAE' molecularProfiles(object, mDataType, assay) <- value ## S4 method for signature 'CoreSet' featureInfo(object, mDataType) ## S4 replacement method for signature 'CoreSet,character,data.frame' featureInfo(object, mDataType) <- value ## S4 method for signature 'CoreSet,character' phenoInfo(object, mDataType) ## S4 replacement method for signature 'CoreSet,character,data.frame' phenoInfo(object, mDataType) <- value ## S4 method for signature 'CoreSet,character' fNames(object, mDataType) ## S4 replacement method for signature 'CoreSet,character,character' fNames(object, mDataType) <- value ## S4 method for signature 'CoreSet' mDataNames(object) ## S4 replacement method for signature 'CoreSet' mDataNames(object) <- value ## S4 method for signature 'CoreSet' molecularProfilesSlot(object) ## S4 replacement method for signature 'CoreSet,list_OR_MAE' molecularProfilesSlot(object) <- value ## S4 method for signature 'CoreSet' sensitivityInfo(object, dimension, ...) ## S4 replacement method for signature 'CoreSet,data.frame' sensitivityInfo(object, dimension, ...) <- value ## S4 method for signature 'CoreSet' sensitivityMeasures(object) ## S4 replacement method for signature 'CoreSet,character' sensitivityMeasures(object) <- value ## S4 method for signature 'CoreSet' sensitivityProfiles(object) ## S4 replacement method for signature 'CoreSet,data.frame' sensitivityProfiles(object) <- value ## S4 method for signature 'CoreSet' sensitivityRaw(object) ## S4 replacement method for signature 'CoreSet,array' sensitivityRaw(object) <- value ## S4 method for signature 'CoreSet' treatmentResponse(object) ## S4 replacement method for signature 'CoreSet,list_OR_LongTable' treatmentResponse(object) <- value ## S4 method for signature 'CoreSet' sensNumber(object) ## S4 replacement method for signature 'CoreSet,matrix' sensNumber(object) <- value ## S4 method for signature 'CoreSet' pertNumber(object) ## S4 replacement method for signature 'CoreSet,array' pertNumber(object) <- value
object |
A |
value |
See details. |
mDataType |
|
assay |
|
dimension |
See details. |
... |
See details. |
annotation: A list
of CoreSet annotations with items: 'name',
the name of the object; 'dateCreated', date the object was created; 'sessionInfo',
the sessionInfo()
when the object was created; 'call', the R constructor call;
and 'version', the object version.
annotation<-: Setter method for the annotation slot. Arguments:
value: a list
of annotations to update the CoreSet with.
dateCreated: character(1)
The date the CoreSet
object was
created, as returned by the date()
function.
dateCreated<-: Update the 'dateCreated' item in the annotation
slot of
a CoreSet
object. Arguments:
value: A character(1)
vector, as returned by the date()
function.
name: character(1)
The name of the CoreSet
, retreived from
the @annotation
slot.
name<-: Update the @annotation$name
value in a CoreSet
object.
value: character(1)
The name of the CoreSet
object.
cellInfo: data.frame
Metadata for all sample in a CoreSet
object.
sampleInfo<-: assign updated sample annotations to the CoreSet
object.
Arguments:
value: a data.frame
object.
sampleNames: character
Retrieve the rownames of the data.frame
in
the sample
slot from a CoreSet object.
sampleNames<-: assign new rownames to the sampleInfo data.frame
for
a CoreSet object.
Arguments:
value: character
vector of rownames for the sampleInfo(object)
data.frame
.
treatmentInfo: data.frame
Metadata for all treatments in a CoreSet
object. Arguments:
object: CoreSet
An object to retrieve treatment metadata from.
treatmentInfo<-: CoreSet
object with updated treatment metadata.
object. Arguments:
object: CoreSet
An object to set treatment metadata for.
value: data.frame
A new table of treatment metadata for object
.
treatmentNames: character
Names for all treatments in a CoreSet
object. Arguments:
object: CoreSet
An object to retrieve treatment names from.
treatmentNames<-: CoreSet
Object with updates treatment names.
object. Arguments:
object: CoreSet
An object to set treatment names from.
value: character
A character vector of updated treatment names.
curation: A list
of curated mappings between identifiers in the
CoreSet object and the original data publication. Contains two data.frame
s, 'sample' with sample ids and
'tissue' with tissue ids.
curation<-: Update the curation
slot of a CoreSet object. Arugments:
value: A list
of data.frame
s, one for each type of curated
identifier. For a CoreSet
object the slot should contain tissue and
sample id data.frame
s.
datasetType: character(1)
The type treatment response in the
sensitivity
slot. Valid values are 'sensitivity', 'perturbation' or 'both'.
datasetType<-: Update the datasetType slot of a CoreSet object. Arguments:
value: A character(1)
vector with one of 'sensitivity', 'perturbation'
or 'both'
molecularProfiles: matrix()
Retrieve an assay in a
SummarizedExperiment
from the molecularProfiles
slot of a CoreSet
object with the specified mDataType
. Valid mDataType
arguments can be
found with mDataNames(object)
. Exclude mDataType
and assay
to
access the entire slot. Arguments:
assay: Optional character(1)
vector specifying an assay in the
SummarizedExperiment
of the molecularProfiles
slot of the
CoreSet
object for the specified mDataType
. If excluded,
defaults to modifying the first assay in the SummarizedExperiment
for
the given mDataType
.
molecularProfiles<-: Update an assay in a SummarizedExperiment
from
the molecularProfiles
slot of a CoreSet object with the specified
mDataType
. Valid mDataType
arguments can be found with
mDataNames(object)
. Omit mDataType
and assay
to update the slot.
assay: Optional character(1)
vector specifying an assay in the
SummarizedExperiment
of the molecularProfiles
slot of the
CoreSet
object for the specified mDataType
. If excluded,
defaults to modifying the first assay in the SummarizedExperiment
for
the given mDataType
.
value: A matrix
of values to assign to the assay
slot of the
SummarizedExperiment
for the selected mDataType
. The rownames and
column names must match the associated SummarizedExperiment
.
featureInfo: Retrieve a DataFrame
of feature metadata for the specified
mDataType
from the molecularProfiles
slot of a CoreSet
object. More
specifically, retrieve the @rowData
slot from the SummarizedExperiment
from the @molecularProfiles
of a CoreSet
object with the name
mDataType
.
featureInfo<-: Update the featureInfo(object, mDataType)
DataFrame
with new feature metadata. Arguments:
value: A data.frame
or DataFrame
with updated feature metadata for
the specified molecular profile in the molecularProfiles
slot of a
CoreSet
object.
phenoInfo: Return the @colData
slot from the SummarizedExperiment
of
mDataType
, containing sample-level metadata, from a CoreSet
object.
phenoInfo<-: Update the @colData
slot of the SummarizedExperiment
of mDataType
in the @molecularProfiles
slot of a CoreSet
object.
This updates the sample-level metadata in-place.
value: A data.frame
or DataFrame
object where rows are samples
and columns are sample metadata.
fNames: character()
The features names from the rowData
slot of a
SummarizedExperiment
of mDataType
within a CoreSet
object.
fNames: Updates the rownames of the feature metadata (i.e., rowData
)
for a SummarizedExperiment
of mDataType
within a CoreSet
object.
value: character()
A character vector of new features names for the
rowData
of the SummarizedExperiment
of mDataType
in the
@molecularProfiles
slot of a CoreSet
object. Must be the same
length as nrow(featureInfo(object, mDataType))
,
the number of rows in the feature metadata.
mDataNames: character
Retrieve the names of the molecular data types
available in the molecularProfiles
slot of a CoreSet
object. These
are the options which can be used in the mDataType
parameter of various
molecularProfiles
slot accessors methods.
mDataNames: Update the molecular data type names of the
molecularProfiles
slot of a CoreSet object. Arguments:
value: character
vector of molecular datatype names, with length
equal to length(molecularProfilesSlot(object))
.
molecularProfilesSlot: Return the contents of the @molecularProfiles
slot of a CoreSet
object. This will either be a list
or
MultiAssayExperiment
of SummarizedExperiment
s.
molecularProfilesSlot<-: Update the contents of the @molecularProfiles
slot of a CoreSet
object. Arguemnts:
value: A list
or MultiAssayExperiment
of SummarizedExperiment
s. The
list
and assays
should be named for the molecular datatype in each
SummarizedExperiment
.
dimension
: Optional character(1)
One of 'treatment', 'sample' or
'assay' to retrieve rowData
, colData
or the 'assay_metadata' assay from
the CoreSet
@sensitvity
LongTable
object, respectively. Ignored with
warning if @treatmentResponse
is not a LongTable
object.
...
: Additional arguments to the rowData
or colData
.
LongTable
methods. Only used if the sensitivity slot contains a
LongTable
object instead of a list
and the dimension
argument is
specified.
sensitivityInfo: DataFrame
or data.frame
of sensitivity treatment combo
by sample metadata for the CoreSet
object. When the dimension
parameter is used, it allows retrieval of the dimension specific metadata
from the LongTable
object in @treatmentResponse
of a CoreSet object.
sensitivityInfo<-: Update the @treatmentResponse
slot metadata for a
CoreSet
object. When used without the dimension
argument is behaves
similar to the old CoreSet implementation, where the @treatmentResponse
slot
contained a list with a $info
data.frame
item. When the dimension
arugment is used, more complicated assignments can occur where 'sample'
modifies the @sensitvity
LongTable
colData, 'treatment' the rowData and
'assay' the 'assay_metadata' assay.
Arguments:
value: A data.frame
of treatment response experiment metadata,
documenting experiment level metadata (mapping to treatments and samples). If
the @treatmentResponse
slot doesn't contain a LongTable
and dimension
is
not specified, you can only modify existing columns as returned by
sensitivityInfo(object)
.
sensitivityMeaures: Get the 'sensitivityMeasures' available in a CoreSet
object. Each measure reprents some summary of sample sensitivity to a given
treatment, such as ic50, ec50, AUC, AAC, etc. The results are returned as a
character
vector with all available metrics for the PSet object.
sensitivityMeaures: Update the sensitivity meaure in a CoreSet
object. Thesee values are the column names of the 'profiles' assay and
represent various compued sensitviity metrics such as ic50, ec50, AUC, AAC,
etc.
value: A character
vector of new sensitivity measure names, the
then length of the character vector must matcht he number of columns of the
'profiles' assay, excluding metadata and key columns.
sensitivityProfiles: Return the sensitivity profile summaries from the sensitivity slot. This data.frame cotanins vaarious sensitivity summary metrics, such as ic50, amax, EC50, aac, HS, etc as columns, with rows as treatment by sample experiments.
sensitivityProfiles<-: Update the sensitivity profile summaries the
sensitivity slot. Arguments:
-value: A data.frame
the the same number of rows as as returned by
sensitivityProfiles(object)
, but potentially modified columns, such as the
computation of additional summary metrics.
sensitivityRaw: Access the raw sensitiity measurents for a CoreSet
object. A 3D array
where rows are experiment_ids, columns are doses
and the third dimension is metric, either 'Dose' for the doses used or
'Viability' for the sample viability at that dose.
sensitvityRaw<-: Update the raw dose and viability data in a CoreSet
.
value: A 3D array
object where rows are experiment_ids, columns are
replicates and pages are c('Dose', 'Viability'), with the corresponding
dose or viability measurement for that experiment_id and replicate.
sensNumber: Return a count of viability observations in a CoreSet
object for each treatment-combo by sample combination.
sensNumber<-: Update the 'n' item, which holds a matrix with a count
of treatment by sample-line experiment counts, in the list
in @treatmentResponse
slot of a CoreSet
object. Will error when @sensitviity
contains
a LongTable
object, since the counts are computed on the fly. Arguments:
value: A matrix
where rows are samples and columns are treatments, with a
count of the number of experiments for each combination as the values.
pertNumber: array
Summary of available perturbation experiments
from in a CoreSet
object. Returns a 3D array
with the number of
perturbation experiments per treatment and sample, and data type.
pertNumber<-: Update the @perturbation$n
value in a CoreSet
object,
which stores a summary of the available perturbation experiments. Arguments:
value: A new 3D array
with the number of perturbation experiments per
treatment and sample, and data type
Accessors: See details.
Setters: An updated CoreSet
object, returned invisibly.
data(clevelandSmall_cSet) ## @annotation annotation(clevelandSmall_cSet) annotation(clevelandSmall_cSet) <- annotation(clevelandSmall_cSet) dateCreated(clevelandSmall_cSet) ## dateCreated dateCreated(clevelandSmall_cSet) <- date() name(clevelandSmall_cSet) name(clevelandSmall_cSet) <- 'new_name' sampleInfo(clevelandSmall_cSet) <- sampleInfo(clevelandSmall_cSet) sampleNames(clevelandSmall_cSet) sampleNames(clevelandSmall_cSet) <- sampleNames(clevelandSmall_cSet) treatmentInfo(clevelandSmall_cSet) treatmentInfo(clevelandSmall_cSet) <- treatmentInfo(clevelandSmall_cSet) treatmentNames(clevelandSmall_cSet) treatmentNames(clevelandSmall_cSet) <- treatmentNames(clevelandSmall_cSet) ## curation curation(clevelandSmall_cSet) curation(clevelandSmall_cSet) <- curation(clevelandSmall_cSet) datasetType(clevelandSmall_cSet) datasetType(clevelandSmall_cSet) <- 'both' # No assay specified molecularProfiles(clevelandSmall_cSet, 'rna') <- molecularProfiles(clevelandSmall_cSet, 'rna') # Specific assay molecularProfiles(clevelandSmall_cSet, 'rna', 'exprs') <- molecularProfiles(clevelandSmall_cSet, 'rna', 'exprs') # Replace the whole slot molecularProfiles(clevelandSmall_cSet) <- molecularProfiles(clevelandSmall_cSet) featureInfo(clevelandSmall_cSet, 'rna') featureInfo(clevelandSmall_cSet, 'rna') <- featureInfo(clevelandSmall_cSet, 'rna') phenoInfo(clevelandSmall_cSet, 'rna') phenoInfo(clevelandSmall_cSet, 'rna') <- phenoInfo(clevelandSmall_cSet, 'rna') fNames(clevelandSmall_cSet, 'rna') fNames(clevelandSmall_cSet, 'rna') <- fNames(clevelandSmall_cSet, 'rna') mDataNames(clevelandSmall_cSet) mDataNames(clevelandSmall_cSet) <- mDataNames(clevelandSmall_cSet) molecularProfilesSlot(clevelandSmall_cSet) molecularProfilesSlot(clevelandSmall_cSet) <- molecularProfilesSlot(clevelandSmall_cSet) sensitivityInfo(clevelandSmall_cSet) sensitivityInfo(clevelandSmall_cSet) <- sensitivityInfo(clevelandSmall_cSet) sensitivityMeasures(clevelandSmall_cSet) <- sensitivityMeasures(clevelandSmall_cSet) sensitivityMeasures(clevelandSmall_cSet) <- sensitivityMeasures(clevelandSmall_cSet) sensitivityProfiles(clevelandSmall_cSet) sensitivityProfiles(clevelandSmall_cSet) <- sensitivityProfiles(clevelandSmall_cSet) head(sensitivityRaw(clevelandSmall_cSet)) sensitivityRaw(clevelandSmall_cSet) <- sensitivityRaw(clevelandSmall_cSet) treatmentResponse(clevelandSmall_cSet) treatmentResponse(clevelandSmall_cSet) <- treatmentResponse(clevelandSmall_cSet) sensNumber(clevelandSmall_cSet) sensNumber(clevelandSmall_cSet) <- sensNumber(clevelandSmall_cSet) pertNumber(clevelandSmall_cSet) pertNumber(clevelandSmall_cSet) <- pertNumber(clevelandSmall_cSet)
data(clevelandSmall_cSet) ## @annotation annotation(clevelandSmall_cSet) annotation(clevelandSmall_cSet) <- annotation(clevelandSmall_cSet) dateCreated(clevelandSmall_cSet) ## dateCreated dateCreated(clevelandSmall_cSet) <- date() name(clevelandSmall_cSet) name(clevelandSmall_cSet) <- 'new_name' sampleInfo(clevelandSmall_cSet) <- sampleInfo(clevelandSmall_cSet) sampleNames(clevelandSmall_cSet) sampleNames(clevelandSmall_cSet) <- sampleNames(clevelandSmall_cSet) treatmentInfo(clevelandSmall_cSet) treatmentInfo(clevelandSmall_cSet) <- treatmentInfo(clevelandSmall_cSet) treatmentNames(clevelandSmall_cSet) treatmentNames(clevelandSmall_cSet) <- treatmentNames(clevelandSmall_cSet) ## curation curation(clevelandSmall_cSet) curation(clevelandSmall_cSet) <- curation(clevelandSmall_cSet) datasetType(clevelandSmall_cSet) datasetType(clevelandSmall_cSet) <- 'both' # No assay specified molecularProfiles(clevelandSmall_cSet, 'rna') <- molecularProfiles(clevelandSmall_cSet, 'rna') # Specific assay molecularProfiles(clevelandSmall_cSet, 'rna', 'exprs') <- molecularProfiles(clevelandSmall_cSet, 'rna', 'exprs') # Replace the whole slot molecularProfiles(clevelandSmall_cSet) <- molecularProfiles(clevelandSmall_cSet) featureInfo(clevelandSmall_cSet, 'rna') featureInfo(clevelandSmall_cSet, 'rna') <- featureInfo(clevelandSmall_cSet, 'rna') phenoInfo(clevelandSmall_cSet, 'rna') phenoInfo(clevelandSmall_cSet, 'rna') <- phenoInfo(clevelandSmall_cSet, 'rna') fNames(clevelandSmall_cSet, 'rna') fNames(clevelandSmall_cSet, 'rna') <- fNames(clevelandSmall_cSet, 'rna') mDataNames(clevelandSmall_cSet) mDataNames(clevelandSmall_cSet) <- mDataNames(clevelandSmall_cSet) molecularProfilesSlot(clevelandSmall_cSet) molecularProfilesSlot(clevelandSmall_cSet) <- molecularProfilesSlot(clevelandSmall_cSet) sensitivityInfo(clevelandSmall_cSet) sensitivityInfo(clevelandSmall_cSet) <- sensitivityInfo(clevelandSmall_cSet) sensitivityMeasures(clevelandSmall_cSet) <- sensitivityMeasures(clevelandSmall_cSet) sensitivityMeasures(clevelandSmall_cSet) <- sensitivityMeasures(clevelandSmall_cSet) sensitivityProfiles(clevelandSmall_cSet) sensitivityProfiles(clevelandSmall_cSet) <- sensitivityProfiles(clevelandSmall_cSet) head(sensitivityRaw(clevelandSmall_cSet)) sensitivityRaw(clevelandSmall_cSet) <- sensitivityRaw(clevelandSmall_cSet) treatmentResponse(clevelandSmall_cSet) treatmentResponse(clevelandSmall_cSet) <- treatmentResponse(clevelandSmall_cSet) sensNumber(clevelandSmall_cSet) sensNumber(clevelandSmall_cSet) <- sensNumber(clevelandSmall_cSet) pertNumber(clevelandSmall_cSet) pertNumber(clevelandSmall_cSet) <- pertNumber(clevelandSmall_cSet)
CoreSet - A generic data container for molecular profiles and treatment response data
The CoreSet (cSet) class was developed as a superclass for pSets in the PharmacoGx and RadioGx packages to contain the data generated in screens of cancer sample lines for their genetic profile and sensitivities to therapy (Pharmacological or Radiation). This class is meant to be a superclass which is contained within the PharmacoSet (pSet) and RadioSet (rSet) objects exported by PharmacoGx and RadioGx. The format of the data is similar for both pSets and rSets, allowing much of the code to be abstracted into the CoreSet super-class. However, the models involved with quantifying sampleular response to Pharmacological and Radiation therapy are widely different, and extension of the cSet class allows the packages to apply the correct model for the given data.
annotation
See Slots section.
molecularProfiles
See Slots section.
sample
See Slots section.
treatment
See Slots section.
treatmentResponse
See Slots section.
perturbation
See Slots section.
curation
See Slots section.
datasetType
See Slots section.
annotation: A list
of annotation data about the CoreSet
,
including the $name
and the session information for how the object
was created, detailing the exact versions of R and all the packages used.
molecularProfiles: A list
or MultiAssayExperiment
containing
CoreSet
object.
sample: A data.frame
containg the annotations for all the samples
profiled in the data set, across all molecular data types and
treatment response experiments.
treatment: A data.frame
containing the annotations for all treatments
in the dataset, including the mandatory 'treatmentid' column to uniquely
identify each treatment.
treatmentResponse: A list
or LongTable
containing all the data for the
treatment response experiment, including $info
, a data.frame
containing the experimental info, $raw
a 3D array
containing raw data,
$profiles
, a data.frame
containing sensitivity profiles
statistics, and $n
, a data.frame
detailing the number of
experiments for each sample-drug/radiationInfo pair
perturbation: list
containing $n
, a data.frame
summarizing the available perturbation data. This slot is currently
being deprecated.
curation: list
containing mappings for treatment
,
sample
and tissue
names used in the data set to universal
identifiers used between different CoreSet
objects.
datasetType: character
string of 'sensitivity',
'perturbation', or both detailing what type of data can be found in the
CoreSet
, for proper processing of the data
CoreSet
object.Documentation for utility methods for a CoreSet
object, such as
set operations like subset and intersect. See @details for information
on different types of methods and their implementations.
## S4 method for signature 'CoreSet' subsetBySample(x, samples) ## S4 method for signature 'CoreSet' subsetByTreatment(x, treatments) ## S4 method for signature 'CoreSet' subsetByFeature(x, features, mDataTypes)
## S4 method for signature 'CoreSet' subsetBySample(x, samples) ## S4 method for signature 'CoreSet' subsetByTreatment(x, treatments) ## S4 method for signature 'CoreSet' subsetByFeature(x, features, mDataTypes)
x |
A |
samples |
|
treatments |
|
features |
|
mDataTypes |
|
subsetBySample: Subset a CoreSet
object by sample identifier.
value: a CoreSet
object containing only samples
.
subsetByTreatment: Subset a CoreSet
object by treatment identifier.
value: a CoreSet
object containing only treatments
.
subsetByFeature: Subset a CoreSet
object by molecular feature
identifier.
value: a CoreSet
object containing only features
.
See details.
data(clevelandSmall_cSet) ## subset methods ### subsetBySample samples <- sampleInfo(clevelandSmall_cSet)$sampleid[seq_len(10)] clevelandSmall_cSet_sub <- subsetBySample(clevelandSmall_cSet, samples) ## subset methods ### subsetByTreatment #treatments <- treatmentInfo(clevelandSmall_cSet)$treatmentid[seq_len(10)] #clevelandSmall_cSet_sub <- subsetByTreatment(clevelandSmall_cSet, treatments) ## subset methods ### subsetByFeature features <- fNames(clevelandSmall_cSet, 'rna')[seq_len(5)] clevelandSmall_cSet_sub <- subsetByFeature(clevelandSmall_cSet, features, 'rna')
data(clevelandSmall_cSet) ## subset methods ### subsetBySample samples <- sampleInfo(clevelandSmall_cSet)$sampleid[seq_len(10)] clevelandSmall_cSet_sub <- subsetBySample(clevelandSmall_cSet, samples) ## subset methods ### subsetByTreatment #treatments <- treatmentInfo(clevelandSmall_cSet)$treatmentid[seq_len(10)] #clevelandSmall_cSet_sub <- subsetByTreatment(clevelandSmall_cSet, treatments) ## subset methods ### subsetByFeature features <- fNames(clevelandSmall_cSet, 'rna')[seq_len(5)] clevelandSmall_cSet_sub <- subsetByFeature(clevelandSmall_cSet, features, 'rna')
New implementation of the CoreSet constructor to support MAE and TRE. This
constructor will be swapped with the original CoreSet
constructor as
part of an overhaul of the CoreSet class structure.
CoreSet2( name = "emptySet", treatment = data.frame(), sample = data.frame(), molecularProfiles = MultiAssayExperiment(), treatmentResponse = LongTable(), datasetType = "sensitivity", perturbation = list(n = array(dim = 3), info = "No perturbation data!"), curation = list(sample = data.frame(), treatment = data.frame()) )
CoreSet2( name = "emptySet", treatment = data.frame(), sample = data.frame(), molecularProfiles = MultiAssayExperiment(), treatmentResponse = LongTable(), datasetType = "sensitivity", perturbation = list(n = array(dim = 3), info = "No perturbation data!"), curation = list(sample = data.frame(), treatment = data.frame()) )
name |
A |
treatment |
A |
sample |
A |
molecularProfiles |
A |
treatmentResponse |
A |
datasetType |
A deprecated slot in a |
perturbation |
A deprecated slot in a |
curation |
A |
A CoreSet
object storing standardized and curated treatment
response and multiomic profile data associated with a given publication.
data(clevelandSmall_cSet) clevelandSmall_cSet
data(clevelandSmall_cSet) clevelandSmall_cSet
Computes the cosine similarity and significance using permutation test. This
function uses random numbers, to ensure reproducibility please call
set.seed()
before running the function.
cosinePerm( x, y, nperm = 1000, alternative = c("two.sided", "less", "greater"), include.perm = FALSE, nthread = 1, ... )
cosinePerm( x, y, nperm = 1000, alternative = c("two.sided", "less", "greater"), include.perm = FALSE, nthread = 1, ... )
x |
|
y |
|
nperm |
|
alternative |
|
include.perm |
|
nthread |
|
... |
A |
A list
estimate of the cosine similarity, p-value and
estimates after random permutations (null distribution) in include.perm is
set to 'TRUE'
x <- factor(c(1,2,1,2,1)) y <- factor(c(2,2,1,1,1)) cosinePerm(x, y)
x <- factor(c(1,2,1,2,1)) y <- factor(c(2,2,1,1,1)) cosinePerm(x, y)
DataMapper
object.Documentation for the various setters and getters which allow manipulation
of data in the slots of a DataMapper
object.
## S4 method for signature 'DataMapper' rawdata(object) ## S4 replacement method for signature 'DataMapper,ANY' rawdata(object) <- value
## S4 method for signature 'DataMapper' rawdata(object) ## S4 replacement method for signature 'DataMapper,ANY' rawdata(object) <- value
object |
A |
value |
A |
rawdata: Get the raw data slot from a DataMapper
object. Returns
a list-like containing one or more raw data inputs to the
DataMapper
object.
rawdata: Set the raw data slot from a DataMapper
object.
value: The list
-like object to set for the rawdata slot. Note: this
currently only supports data.frame
or data.table
objects.
Accessors: See details
Setters: An update DataMapper
object, returned invisibly.
Other DataMapper-accessors:
LongTableDataMapper-accessors
,
TREDataMapper-accessors
This object will be used as a way to abstract away data preprocessing.
rawdata: A list-like object containing one or more pieces of raw data
that will be processed and mapped to the slots of an S4
object.
metadata: A List
of object level metadata.
Drop parameters from a function and replace them with constants inside the function body.
drop_fn_params(fn, args)
drop_fn_params(fn, args)
fn |
|
args |
|
function
A new non-primitize function with the parameters named in
args
deleted and their values fixed with the values from args
in the
function body.
Perform aggregation over an S4 object, but return an object of the same class.
endoaggregate(x, ...)
endoaggregate(x, ...)
x |
An |
... |
|
An object with the same class as x
.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
LongTable
or
inheriting classCompute a group-by operation over a LongTable
object or its inhering
classes.
## S4 method for signature 'LongTable' endoaggregate( x, ..., assay, target = assay, by, subset = TRUE, nthread = 1, progress = TRUE, BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
## S4 method for signature 'LongTable' endoaggregate( x, ..., assay, target = assay, by, subset = TRUE, nthread = 1, progress = TRUE, BPPARAM = NULL, enlist = TRUE, moreArgs = list() )
x |
|
... |
|
assay |
|
target |
|
by |
|
subset |
|
nthread |
|
progress |
|
BPPARAM |
|
enlist |
|
moreArgs |
|
Arguments in ...
are substituted and wrapped in a list, which is passed
through to the j argument of [.data.table
internally. The function currently
tries to build informative column names for unnamed arguments in ...
by
appending the name of each function call with the name of its first argument,
which is assumed to be the column name being aggregated over. If an argument
to ...
is named, that will be the column name of its value in the resulting
data.table
.
The primary use case for enlist=FALSE
is to allow computation of dependent
aggregations, where the output from a previous aggregation is required in a
subsequent one. For this case, wrap your call in {
and assign intermediate
results to variables, returning the final results as a list where each list
item will become a column in the final table with the corresponding name.
Name inference is disabled for this case, since it is assumed you will name
the returned list items appropriately.
A major advantage over multiple calls to aggregate
is that
the overhead of parallelization is paid only once even for complex multi-step
computations like fitting a model, capturing its paramters, and making
predictions using it. It also allows capturing arbitrarily complex calls
which can be recomputed later using the
update,TreatmentResponseExperiment-method
A potential disadvantage is increased RAM usage per
thread due to storing intermediate values in variables, as well as any
memory allocation overhead associate therewith.
Object with the same class as x
, with the aggregation results
assigned to target
, using strategy
if target
is an existing assay in
x
.
data.table::[.data.table
, BiocParallel::bplapply
A dummy LongTableDataMapper object to be used in package examples.
data(exampleDataMapper)
data(exampleDataMapper)
LongTableDataMapper object
Internal slot for storing metadata relevant to the internal operation of an S4 object.
getIntern(object, x, ...)
getIntern(object, x, ...)
object |
|
x |
|
... |
Allow new parmeters to be defined for this generic. |
Warning: This method is intended for developer use and can be ignored by users.
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Generic for Guessing the Mapping Between Some Raw Data and an S4 Object
guessMapping(object, ...)
guessMapping(object, ...)
object |
An |
... |
Allow new arguments to be defined for this generic. |
A list
with mapping guesses as items.
"Generics shouldn't need examples!"
"Generics shouldn't need examples!"
Checks for columns which are uniquely identified by a group of identifiers.
This should be used to help identify the columns required to uniquely
identify the rows, columns, assays and metadata of a DataMapper
class
object.
## S4 method for signature 'LongTableDataMapper' guessMapping(object, groups, subset, data = FALSE)
## S4 method for signature 'LongTableDataMapper' guessMapping(object, groups, subset, data = FALSE)
object |
A |
groups |
A |
subset |
A |
data |
A |
Any unmapped columns will be added to the end of the returned list
in an
item called unmapped.
The function automatically guesses metadata by checking if any columns have only a single value. This is returned as an additional item in the list.
A list
, where each item is named for the associated groups
item
the guess is for. The character vector in each item are columns which are
uniquely identified by the identifiers from that group.
guessMapping(exampleDataMapper, groups=list(rows='treatmentid', cols='sampleid'), subset=FALSE)
guessMapping(exampleDataMapper, groups=list(rows='treatmentid', cols='sampleid'), subset=FALSE)
Calculate the gwc score between two vectors, using either a weighted spearman or pearson correlation
gwc( x1, p1, x2, p2, method.cor = c("pearson", "spearman"), nperm = 10000, truncate.p = 1e-16, ... )
gwc( x1, p1, x2, p2, method.cor = c("pearson", "spearman"), nperm = 10000, truncate.p = 1e-16, ... )
x1 |
|
p1 |
|
x2 |
|
p2 |
|
method.cor |
|
nperm |
|
truncate.p |
|
... |
Other passed down to internal functions |
numeric
a vector of two values, the correlation and associated p-value.
data(clevelandSmall_cSet) x <- molecularProfiles(clevelandSmall_cSet,'rna')[,1] y <- molecularProfiles(clevelandSmall_cSet,'rna')[,2] x_p <- rep(0.05, times=length(x)) y_p <- rep(0.05, times=length(y)) names(x_p) <- names(x) names(y_p) <- names(y) gwc(x,x_p,y,y_p, nperm=100)
data(clevelandSmall_cSet) x <- molecularProfiles(clevelandSmall_cSet,'rna')[,1] y <- molecularProfiles(clevelandSmall_cSet,'rna')[,2] x_p <- rep(0.05, times=length(x)) y_p <- rep(0.05, times=length(y)) names(x_p) <- names(x) names(y_p) <- names(y) gwc(x,x_p,y,y_p, nperm=100)
Generic to access the unique id columns in an S4 object used to
idCols(object, ...)
idCols(object, ...)
object |
An |
... |
Allow new arguments to this generic. |
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
This method should allow any S3 object in R to become immutable by
intercepting [<-
, [[<-
, $<-
and c
during S3-method dispatch and
returning an error.
Reverse with call to the mutable
function.
immutable(object) is.immutable(object) ## S3 method for class 'immutable' print(x, ...) show.immutable(x)
immutable(object) is.immutable(object) ## S3 method for class 'immutable' print(x, ...) show.immutable(x)
object , x
|
Any R object which uses S3 method dispatch |
... |
Fallthrough arguments to |
The motivation for this class was to create pseudo-private slots in an R S4 object by preventing mutation of those slots outside of the accessors written for the class. It should behave as expected for R object which operate with 'copy-on-modify' semantics, including most base R functions and S3 objects.
An environment was not suitable for this case due to the 'copy-by-reference' semantics, such that normal R assignment, which users assume makes a copy of the object, actually references the same environment in both the original and copy of the object.
WARNING: This implementation is unable to intercept modifications to a
data.table
via the set*
group of methods. This is because these methods
are not S3 generics and therefore no mechanism exists for hooking into them
to extend their functionality. In general, this helper class will only work
for objects with an S3 interface.
The object
with "immutable" prepended to its class attribute.
logical(1)
Does the object inherit from the "immutable" S3-class?
None, invisible(NULL)
assignment-immutable
, setOps-immutable
immutable_list <- immutable(as.list(1:5)) class(immutable_list) # errors during assignment operations tryCatch({ immutable_list$new <- 1 }, error=print) immutable_list <- immutable(as.list(1:5)) is.immutable(immutable_list)
immutable_list <- immutable(as.list(1:5)) class(immutable_list) # errors during assignment operations tryCatch({ immutable_list$new <- 1 }, error=print) immutable_list <- immutable(as.list(1:5)) is.immutable(immutable_list)
stats::optim
.Functions compatible with optim
have the parameter named par
as their
first formal argument where each value is a respective free parameter to
be optimized.
is_optim_compatible(fn)
is_optim_compatible(fn)
fn |
|
logical(1)
TRUE
if the first value of formalArg(fn)
is "par",
otherwise FALSE
.
Get the types of all items in a list
is.items(list, ..., FUN = is)
is.items(list, ..., FUN = is)
list |
A |
... |
|
FUN |
|
logical
A vector indicating if the list item is the specified
type.
list <- list(c(1,2,3), c('a','b','c')) is.items(list, 'character')
list <- list(c(1,2,3), c('a','b','c')) is.items(list, 'character')
MultiAssayExperiment
lapply
lapply method for MultiAssayExperiment
## S4 method for signature 'MultiAssayExperiment' lapply(X, FUN, ...)
## S4 method for signature 'MultiAssayExperiment' lapply(X, FUN, ...)
X |
A |
FUN |
A function to be applied to each |
... |
Fall through parameters to |
A MultiAssayExperiment
object, modifed such that
experiments(X) <- endoapply(experiments(X), FUN, ...)
.s
A class union to allow multiple types in a CoreSet slot
LongTable constructor method
LongTable( rowData, rowIDs, colData, colIDs, assays, assayIDs, metadata = list(), keep.rownames = FALSE )
LongTable( rowData, rowIDs, colData, colIDs, assays, assayIDs, metadata = list(), keep.rownames = FALSE )
rowData |
|
rowIDs |
|
colData |
|
colIDs |
|
assays |
|
assayIDs |
|
metadata |
|
keep.rownames |
|
A LongTable
object containing the data for a treatment response
experiment and configured according to the rowIDs and colIDs arguments.
"See vignette('The LongTable Class', package='CoreGx')"
"See vignette('The LongTable Class', package='CoreGx')"
LongTable
Documentation for the various setters and getters which allow manipulation
of data in the slots of a LongTable
object.
Accessors: See details.
Setters: An updated LongTable
object, returned invisibly.
data(merckLongTable)
data(merckLongTable)
LongTableDataMapper
class, which maps from one or
more raw experimental data files to the slots of a LongTable
object.Constructor for the LongTableDataMapper
class, which maps from one or
more raw experimental data files to the slots of a LongTable
object.
LongTableDataMapper( rawdata = data.frame(), rowDataMap = list(character(), character()), colDataMap = list(character(), character()), assayMap = list(list(character(), character())), metadataMap = list(character()) )
LongTableDataMapper( rawdata = data.frame(), rowDataMap = list(character(), character()), colDataMap = list(character(), character()), assayMap = list(list(character(), character())), metadataMap = list(character()) )
rawdata |
A |
rowDataMap |
A list-like object containing two |
colDataMap |
A list-like object containing two |
assayMap |
A list-like where each item is a |
metadataMap |
A list-like where each item is a |
The guessMapping
method can be used to test hypotheses about the
cardinality of one or more sets of identifier columns. This is helpful
to determine the id columns for rowDataMap
and colDataMap
, as well
as identify columns mapping to assays
or metadata
.
To attach metadata not associated with rawdata
, please use the metadata
assignment method on your LongTableDataMapper
. This metadata will be
merged with any metadata from metadataMap
and added to the LongTable
which this object ultimately constructs.
A LongTable
object, with columns mapped to it's slots according
to the various maps in the LongTableDataMapper
object.
data(exampleDataMapper) exampleDataMapper
data(exampleDataMapper) exampleDataMapper
LongTableDataMapper
object.Documentation for the various setters and getters which allow manipulation
of data in the slots of a LongTableDataMapper
object.
## S4 replacement method for signature 'LongTableDataMapper,list' rawdata(object) <- value ## S4 method for signature 'LongTableDataMapper' rowDataMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' rowDataMap(object) <- value ## S4 method for signature 'LongTableDataMapper' colDataMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' colDataMap(object) <- value ## S4 method for signature 'LongTableDataMapper' assayMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' assayMap(object) <- value ## S4 method for signature 'LongTableDataMapper' metadataMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' metadataMap(object) <- value
## S4 replacement method for signature 'LongTableDataMapper,list' rawdata(object) <- value ## S4 method for signature 'LongTableDataMapper' rowDataMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' rowDataMap(object) <- value ## S4 method for signature 'LongTableDataMapper' colDataMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' colDataMap(object) <- value ## S4 method for signature 'LongTableDataMapper' assayMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' assayMap(object) <- value ## S4 method for signature 'LongTableDataMapper' metadataMap(object) ## S4 replacement method for signature 'LongTableDataMapper,list_OR_List' metadataMap(object) <- value
object |
A |
value |
See details. |
rawdata: Get the raw data slot from a LongTableDataMapper
object. Returns
a list-like containing one or more raw data inputs to the
LongTableDataMapper
object.
rawdata: Set the raw data slot from a LongTableDataMapper
object.
value: The list
-like object to set for the rawdata slot. Note: this
currently only supports data.frame
or data.table
objects.
rowDataMap: list
of two character
vectors, the first are the
columns required to uniquely identify each row of a LongTableDataMapper
and the
second any additional row-level metadata. If the character vectors
have names, the resulting columns are automatically renamed to the
item name of the specified column.
rowDataMap<-: Update the @rowDataMap
slot of a LongTableDataMapper
object,
returning an invisible NULL. Arguments:
value: A list
or List
where the first item is the names of the
identifier columns – columns needed to uniquely identify each row in
rowData – and the second item is the metadata associated with those
the identifier columns, but not required to uniquely identify rows in
the object rowData.
colDataMap: list
of two character
vectors, the first are the
columns required to uniquely identify each row of a LongTableDataMapper
and the
second any additional col-level metadata. If the character vectors
have names, the resulting columns are automatically renamed to the
item name of the specified column.
colDataMap<-: Update the @colDataMap
slot of a LongTableDataMapper
object,
returning an invisible NULL. Arguments:
value: A list
or List
where the first item is the names of the
identifier columns – columns needed to uniquely identify each row in
colData – and the second item is the metadata associated with those
the identifier columns, but not required to uniquely identify rows in
the object rowData.
assayMap: A list
of character vectors. The name of each list item
will be the assay in a LongTableDataMapper
object that the columns in the
character
vector will be assigned to. Column renaming occurs automatically
when the character vectors have names (from the value to the name).
assayMap<-: Updates the @assayMap
slot of a LongTableDataMapper
object,
returning an invisible NULL. Arguments:
value: A list
of character vectors, where the name of each list
item is the name of an assay and the values of each character vector
specify the columns mapping to the assay in the S4
object the
LongTableDataMapper
constructs.
metadataMap: A list
of character
vectors. Each item is an element
of the constructed objects @metadata
slot.
metadataMap<-: Updates LongTableDataMapper
object in-place, then returns an
invisible(NULL)
. Arguments:
value: A list
of character
vectors. The name of each list item
is the name of the item in the @metadata
slot of the LongTableDataMapper
object
created when metaConstruct
is called on the DataMapper
, and a
character vector specifies the columns of @rawdata
to assign to each item.
Accessors: See details
Setters: An update LongTableDataMapper
object, returned invisibly.
Other DataMapper-accessors:
DataMapper-accessors
,
TREDataMapper-accessors
rowDataMap(exampleDataMapper) rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) assayMap(exampleDataMapper) <- list(sensitivity=c(viability1='viability')) metadataMap(exampleDataMapper) metadataMap(exampleDataMapper) <- list(object_metadata=c('metadata'))
rowDataMap(exampleDataMapper) rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) assayMap(exampleDataMapper) <- list(sensitivity=c(viability1='viability')) metadataMap(exampleDataMapper) metadataMap(exampleDataMapper) <- list(object_metadata=c('metadata'))
LongTable
ObjectA Class for Mapping Between Raw Data and an LongTable
Object
## S4 method for signature 'LongTableDataMapper' show(object)
## S4 method for signature 'LongTableDataMapper' show(object)
object |
A |
invisible
Prints to console.
show(LongTableDataMapper)
: Show method for LongTableDataMapper.
Determines how the object is displayed in the console.
rawdata
See Slots section.
rowDataMap
See Slots section.
colDataMap
See Slots section.
assayMap
See Slots section.
metadataMap
See Slots section.
rowDataMap: A list-like object containing two character
vectors.
The first is column names in rawdata
needed to uniquely identify each row,
the second is additional columns which map to rows, but are not required to
uniquely identify them. Rows should be drugs.
colDataMap: A list-like object containing two character
vectors.
The first is column names in rawdata
needed to uniquely identify each
column, the second is additional columns which map to rows, but are not
required to uniquely identify them. Columns should be samples.
assayMap A list-like where each item is a list
with two elements
specifying an assay, the first being the identifier columns in rawdata
needed to uniquely identify each row an assay, and the second a list of
rawdata
columns to be mapped to that assay. The names of assayMap
will be the names of the assays in the LongTable
that is created when
calling metaConstruct
on this DataMapper
object.
metadataMap: A list-like where each item is a character
vector of
rawdata
column names to assign to the @metadata
of the LongTable
,
where the name of that assay is the name of the list item. If names are
omitted, assays will be numbered by their index in the list.
rawdata: A list-like object containing one or more pieces of raw data
that will be processed and mapped to the slots of an S4
object.
metadata: A List
of object level metadata.
show(exampleDataMapper)
show(exampleDataMapper)
stats::optim
.Takes a non-primitive R function and refactors it to be compatible with optimization via stats::optim
.
make_optim_function(fn, ...)
make_optim_function(fn, ...)
fn |
|
... |
Arguments to |
drop_fn_params
, collect_fn_params
The function computes a Matthews correlation coefficient for two factors provided to the function. It assumes each factor is a factor of class labels, and the enteries are paired in order of the vectors.
mcc( x, y, nperm = 1000, nthread = 1, alternative = c("two.sided", "less", "greater"), ... )
mcc( x, y, nperm = 1000, nthread = 1, alternative = c("two.sided", "less", "greater"), ... )
x , y
|
|
nperm |
|
nthread |
|
alternative |
indicates the alternative hypothesis and must be one of ‘"two.sided"’, ‘"greater"’ or ‘"less"’. You can specify just the initial letter. ‘"greater"’ corresponds to positive association, ‘"less"’ to negative association. |
... |
|
Please note: we recommend you call set.seed() before using this function to ensure the reproducibility of your results. Write down the seed number or save it in a script if you intend to use the results in a publication.
A list with the MCC as the $estimate, and p value as $p.value
x <- factor(c(1,2,1,2,3,1)) y <- factor(c(2,1,1,1,2,2)) mcc(x,y)
x <- factor(c(1,2,1,2,3,1)) y <- factor(c(2,1,1,1,2,2)) mcc(x,y)
This is a LongTable object created from some drug combination data provided to our lab by Merck.
data(merckLongTable)
data(merckLongTable)
LongTable object
TODO:: Include a reference
S4
object.Merge assays with an S4
object.
mergeAssays(object, ...)
mergeAssays(object, ...)
object |
|
... |
Allow new arguments to be defined for this generic. |
A modified version of object
.
"This is a generic method!"
"This is a generic method!"
LongTable
or inheriting classEndomorphically merge assays within a LongTable
or inheriting class
## S4 method for signature 'LongTable' mergeAssays(object, x, y, target = x, ..., metadata = FALSE)
## S4 method for signature 'LongTable' mergeAssays(object, x, y, target = x, ..., metadata = FALSE)
object |
A |
x |
|
y |
|
target |
|
... |
Fallthrough arguments to merge.data.table to specify the join type. Use this to specify which columns to merge on. If excluded, defaults to by=assayKeys(objecty, y). |
metadata |
|
A copy of object
with assays x
and y
merged and assigned to
target
.
Christopher Eeles
S4
container object.This method is intended to abstract away complex constructor arguments
and data preprocessing steps needed to transform raw data, such as that
produced in a treatment-response or next-gen sequencing experiment, and
automate building of the appropriate S4
container object. This is
is intended to allow mapping between different experimental designs,
in the form of an S4
configuration object, and various S4 class
containers in the Bioconductor community and beyond.
metaConstruct(mapper, ...) ## S4 method for signature 'LongTableDataMapper' metaConstruct(mapper) ## S4 method for signature 'TREDataMapper' metaConstruct(mapper)
metaConstruct(mapper, ...) ## S4 method for signature 'LongTableDataMapper' metaConstruct(mapper) ## S4 method for signature 'TREDataMapper' metaConstruct(mapper)
mapper |
An |
... |
Allow new arguments to be defined for this generic. |
An S4
object for which the class corresponds to the type of
the build configuration object passed to this method.
A LongTable
object, as specified in the mapper.
A TreatmentResponseExperiment
object, as specified in the mapper.
data(exampleDataMapper) rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) <- list(sensitivity=list(c("treatmentid", "sampleid"), c('viability'))) metadataMap(exampleDataMapper) <- list(experiment_metadata=c('metadata')) longTable <- metaConstruct(exampleDataMapper) longTable data(exampleDataMapper) exampleDataMapper <- as(exampleDataMapper, "TREDataMapper") rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) <- list(sensitivity=list(c("treatmentid", "sampleid"), c('viability'))) metadataMap(exampleDataMapper) <- list(experiment_metadata=c('metadata')) tre <- metaConstruct(exampleDataMapper) tre
data(exampleDataMapper) rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) <- list(sensitivity=list(c("treatmentid", "sampleid"), c('viability'))) metadataMap(exampleDataMapper) <- list(experiment_metadata=c('metadata')) longTable <- metaConstruct(exampleDataMapper) longTable data(exampleDataMapper) exampleDataMapper <- as(exampleDataMapper, "TREDataMapper") rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) <- list(sensitivity=list(c("treatmentid", "sampleid"), c('viability'))) metadataMap(exampleDataMapper) <- list(experiment_metadata=c('metadata')) tre <- metaConstruct(exampleDataMapper) tre
LongTable
objectGetter method for the metadata slot of a LongTable
object
## S4 method for signature 'LongTable' metadata(x)
## S4 method for signature 'LongTable' metadata(x)
x |
The |
list
The contents of the metadata
slot of the LongTable
object.
LongTable
objectSetter method for the metadata slot of a LongTable
object
## S4 replacement method for signature 'LongTable' metadata(x) <- value
## S4 replacement method for signature 'LongTable' metadata(x) <- value
x |
|
value |
|
LongTable
A copy of the LongTable
object with the value
in
the metadata slot.
Remove the "immutable" S3-class from an R object, allowing it to be modified normally again.
mutable(object)
mutable(object)
object |
An R object inheriting from the "immutable" class. |
The object
with the "immutable" class stripped from it.
immut_list <- immutable(list()) mutable(immut_list)
immut_list <- immutable(list()) mutable(immut_list)
This is a TreatmentResponseExperiment
object containing a subset of
NCI-ALMANAC drug combination screening data,
with 2347 unique treatment combinations on 10 cancer cell lines selected.
data(nci_TRE_small)
data(nci_TRE_small)
TreatmentResponseExperiment object
Susan L. Holbeck, Richard Camalier, James A. Crowell, Jeevan Prasaad Govindharajulu, Melinda Hollingshead, Lawrence W. Anderson, Eric Polley, Larry Rubinstein, Apurva Srivastava, Deborah Wilsker, Jerry M. Collins, James H. Doroshow; The National Cancer Institute ALMANAC: A Comprehensive Screening Resource for the Detection of Anticancer Drug Pairs with Enhanced Therapeutic Activity. Cancer Res 1 July 2017; 77 (13): 3564–3576. https://doi.org/10.1158/0008-5472.CAN-17-0489
A helper method to find the best multithreading configuration for your computer
optimizeCoreGx(sample_data, set = FALSE, report = !set)
optimizeCoreGx(sample_data, set = FALSE, report = !set)
sample_data |
|
set |
|
report |
|
If set=TRUE
, modifies data.table
threads via setDTthreads()
, otherwise
displays a message indicating the optimal number of threads.
If report=TRUE
, also returns a data.frame
of the benchmark results.
data(merckLongTable) optimizeCoreGx(merckLongTable)
data(merckLongTable) optimizeCoreGx(merckLongTable)
This method allows integer indexes used to maintain referential integrity internal to an S4 object to be reset. This is useful particularly after subsetting an object, as certain indexes may no longer be present in the object data. Reindexing removes gaps integer indexes and ensures that the smallest contiguous integer values are used in an objects indexes.
reindex(object, ...)
reindex(object, ...)
object |
|
... |
|
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
After subsetting a LongTable, it is possible that values of rowKey or colKey could no longer be present in the object. As a result there the indexes will no longer be contiguous integers. This method will calcualte a new set of rowKey and colKey values such that integer indexes are the smallest set of contiguous integers possible for the data.
## S4 method for signature 'LongTable' reindex(object)
## S4 method for signature 'LongTable' reindex(object)
object |
The |
A copy of the LongTable
with all keys as the smallest set of
contiguous integers possible given the current data.
rowData
out of the rawdata
slot using
the assigned rowDataMap
metadata.Convenience method to subset the rowData
out of the rawdata
slot using
the assigned rowDataMap
metadata.
## S4 method for signature 'LongTableDataMapper' rowData(x, key = TRUE)
## S4 method for signature 'LongTableDataMapper' rowData(x, key = TRUE)
x |
|
key |
|
data.table
The rowData
as specified in the rowDataMap
slot.
rowData
out of the rawdata
slot using
the assigned rowDataMap
metadata.Convenience method to subset the rowData
out of the rawdata
slot using
the assigned rowDataMap
metadata.
## S4 method for signature 'TREDataMapper' rowData(x, key = TRUE)
## S4 method for signature 'TREDataMapper' rowData(x, key = TRUE)
x |
|
key |
|
data.table
The rowData
as specified in the rowDataMap
slot.
Generic to access the row identifiers from
rowIDs(object, ...)
rowIDs(object, ...)
object |
|
... |
Allow new arguments to this generic. |
Depends on the implemented method.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Generic to access the row identifiers from
rowMeta(object, ...)
rowMeta(object, ...)
object |
|
... |
Allow new arguments to this generic. |
Depends on the implemented method.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Generic function to get the annotations for a treatment response experiment from an S4 class
sensitivityInfo(object, ...)
sensitivityInfo(object, ...)
object |
An |
... |
Allow new arguments to be defined for this generic. |
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Generic function to get the annotations for a treatment response experiment from an S4 class.
sensitivityInfo(object, ...) <- value
sensitivityInfo(object, ...) <- value
object |
An |
... |
Allow new arguments to be defined for this generic. |
value |
The new treatment response experiment annotations. |
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Get the names of the sensitivity summary metrics available in an S4 object.
sensitivityMeasures(object, ...)
sensitivityMeasures(object, ...)
object |
An |
... |
Fallthrough arguements for defining new methods |
Depends on the implemented method
sensitivityMeasures(clevelandSmall_cSet)
sensitivityMeasures(clevelandSmall_cSet)
Set the names of the sensitivity summary metrics available in an S4 object.
sensitivityMeasures(object, ...) <- value
sensitivityMeasures(object, ...) <- value
object |
An |
... |
Allow new methods to be defined for this generic. |
value |
A set of names for sensitivity measures to use to update the object with. |
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
A generic for sensitivityProfiles getter method
sensitivityProfiles(object, ...)
sensitivityProfiles(object, ...)
object |
The |
... |
|
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
A generic for the sensitivityProfiles replacement method
sensitivityProfiles(object, ...) <- value
sensitivityProfiles(object, ...) <- value
object |
An |
... |
Fallthrough arguments for defining new methods |
value |
An object with the new sensitivity profiles. If a matrix object is passed in, converted to data.frame before assignment |
Updated CoreSet
Generic function to get the raw data array for a treatment response experiment from an S4 class.
sensitivityRaw(object, ...)
sensitivityRaw(object, ...)
object |
An |
... |
|
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Generic function to set the raw data array for a treatment response experiment in an S4 class.
sensitivityRaw(object, ...) <- value
sensitivityRaw(object, ...) <- value
object |
An |
... |
|
value |
An object containing dose and viability metrics to update the object with. |
Depends on the implemented method
Convert the sensitivity slot in an object inheriting from a CoreSet from a list to a LongTable.
sensitivitySlotToLongTable(object, ...)
sensitivitySlotToLongTable(object, ...)
object |
|
... |
Allow new arguments to be defined on this generic. |
A LongTable
object containing the data in the sensitivity slot.
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Subset an immutable object, returning another immutable object.
subset.immutable(x, ...) ## S3 method for class 'immutable' x[...] ## S3 method for class 'immutable' x[[...]] ## S3 method for class 'immutable' x$...
subset.immutable(x, ...) ## S3 method for class 'immutable' x[...] ## S3 method for class 'immutable' x[[...]] ## S3 method for class 'immutable' x$...
x |
An R object inheriting from the "immutable" S3-class. |
... |
Catch any additional parameters. Lets objects with arbitrary dimensions be made immutable. |
An immutable subset of x
.
immut_mat <- immutable(matrix(1:100, 10, 10)) immut_mat[1:5, 1:5]
immut_mat <- immutable(matrix(1:100, 10, 10)) immut_mat[1:5, 1:5]
Show a CoreSet
## S4 method for signature 'CoreSet' show(object)
## S4 method for signature 'CoreSet' show(object)
object |
|
Prints the CoreSet object to the output stream, and returns invisible NULL.
show(clevelandSmall_cSet)
show(clevelandSmall_cSet)
Show method for the LongTable class
## S4 method for signature 'LongTable' show(object)
## S4 method for signature 'LongTable' show(object)
object |
A |
invisible
Prints to console.
show(merckLongTable)
show(merckLongTable)
Signature
class object, as returned by
drugSensitivitysig
or radSensitivtySig
functions available in
PharmacoGx
and RadioGx
, respectively.Get the annotations for a Signature
class object, as returned by
drugSensitivitysig
or radSensitivtySig
functions available in
PharmacoGx
and RadioGx
, respectively.
showSigAnnot(object, ...)
showSigAnnot(object, ...)
object |
A |
... |
Allow definition of new arguments to this generic |
NULL Prints the signature annotations to console
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Allows use of the colData and rowData data.table
objects to query based on
rowID and colID, which is then used to subset all assay data.table
s stored
in the assays
slot.
This function is endomorphic, it always returns a LongTable object.
## S4 method for signature 'LongTable' subset(x, i, j, assays = assayNames(x), reindex = TRUE)
## S4 method for signature 'LongTable' subset(x, i, j, assays = assayNames(x), reindex = TRUE)
x |
|
i |
|
j |
|
assays |
|
reindex |
|
LongTable
A new LongTable
object subset based on the specified
parameters.
# Character subset(merckLongTable, 'ABT-888', 'CAOV3') # Numeric subset(merckLongTable, 1, c(1, 2)) # Logical subset(merckLongTable, , colData(merckLongTable)$sampleid == 'A2058') # Call subset(merckLongTable, drug1id == 'Dasatinib' & drug2id != '5-FU', sampleid == 'A2058')
# Character subset(merckLongTable, 'ABT-888', 'CAOV3') # Numeric subset(merckLongTable, 1, c(1, 2)) # Logical subset(merckLongTable, , colData(merckLongTable)$sampleid == 'A2058') # Call subset(merckLongTable, drug1id == 'Dasatinib' & drug2id != '5-FU', sampleid == 'A2058')
Summarize molecular profile data such that there is a single entry for each sample line/treatment combination
summarizeMolecularProfiles(object, ...)
summarizeMolecularProfiles(object, ...)
object |
An |
... |
Allow definition of new arguments to this generic |
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Summarize across replicates for a sensitivity dose-response experiment
summarizeSensitivityProfiles(object, ...)
summarizeSensitivityProfiles(object, ...)
object |
An |
... |
Allow definition of new arguments to this generic |
Depends on the implemented method
print("Generics shouldn't need examples?")
print("Generics shouldn't need examples?")
Builds a TreatmentResponseExperiment
object from rectangular
objects. The rowData
argument should contain row level metadata, while
the colData
argument should contain column level metadata, for the
experimental assays
in the assays
list. The rowIDs
and colIDs
lists are used to configure
the internal keys mapping rows or columns to rows in the assays. Each list
should contain at minimum one character vector, specifying which columns
in rowData
or colData
are required to uniquely identify each row. An
optional second character vector can be included, specifying any metadata
columns for either dimension. These should contain information about each
row but NOT be required to uniquely identify a row in the colData
or
rowData
objects. Additional metadata can be attached to a
TreatmentResponseExperiment
by passing a list to the metadata argument.
TreatmentResponseExperiment( rowData, rowIDs, colData, colIDs, assays, assayIDs, metadata = list(), keep.rownames = FALSE )
TreatmentResponseExperiment( rowData, rowIDs, colData, colIDs, assays, assayIDs, metadata = list(), keep.rownames = FALSE )
rowData |
|
rowIDs |
|
colData |
|
colIDs |
|
assays |
A |
assayIDs |
|
metadata |
A |
keep.rownames |
|
For now this class is simply a wrapper around a LongTable
class. In the
future we plan to refactor CoreGx such that the LongTable
class is in a
separate pacakge. We can then specialize the implementation of
TreatmentResponseExperiment
to better capture the biomedical nature of
this object.
A TreatmentResponseExperiment
object containing the data for a
treatment response experiment configured according to the rowIDs and
colIDs arguments.
TREDataMapper
class, which maps from one or
more raw experimental data files to the slots of a LongTable
object.Constructor for the TREDataMapper
class, which maps from one or
more raw experimental data files to the slots of a LongTable
object.
TREDataMapper( rawdata = data.frame(), rowDataMap = list(character(), character()), colDataMap = list(character(), character()), assayMap = list(list(character(), character())), metadataMap = list(character()) )
TREDataMapper( rawdata = data.frame(), rowDataMap = list(character(), character()), colDataMap = list(character(), character()), assayMap = list(list(character(), character())), metadataMap = list(character()) )
rawdata |
A |
rowDataMap |
A list-like object containing two |
colDataMap |
A list-like object containing two |
assayMap |
A list-like where each item is a |
metadataMap |
A list-like where each item is a |
The guessMapping
method can be used to test hypotheses about the
cardinality of one or more sets of identifier columns. This is helpful
to determine the id columns for rowDataMap
and colDataMap
, as well
as identify columns mapping to assays
or metadata
.
To attach metadata not associated with rawdata
, please use the metadata
assignment method on your TREDataMapper
. This metadata will be
merge with any metadata from metadataMap
and added to the LongTable
which this object ultimately constructs.
A TREDataMapper
object, with columns mapped to it's slots according
to the various maps in the LongTableDataMapper
object.
TREDataMapper
object.Documentation for the various setters and getters which allow manipulation
of data in the slots of a TREDataMapper
object.
## S4 replacement method for signature 'TREDataMapper,list' rawdata(object) <- value ## S4 method for signature 'TREDataMapper' rowDataMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' rowDataMap(object) <- value ## S4 method for signature 'TREDataMapper' colDataMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' colDataMap(object) <- value ## S4 method for signature 'TREDataMapper' assayMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' assayMap(object) <- value ## S4 method for signature 'TREDataMapper' metadataMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' metadataMap(object) <- value
## S4 replacement method for signature 'TREDataMapper,list' rawdata(object) <- value ## S4 method for signature 'TREDataMapper' rowDataMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' rowDataMap(object) <- value ## S4 method for signature 'TREDataMapper' colDataMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' colDataMap(object) <- value ## S4 method for signature 'TREDataMapper' assayMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' assayMap(object) <- value ## S4 method for signature 'TREDataMapper' metadataMap(object) ## S4 replacement method for signature 'TREDataMapper,list_OR_List' metadataMap(object) <- value
object |
A |
value |
See details. |
rawdata: Get the raw data slot from a TREDataMapper
object. Returns
a list-like containing one or more raw data inputs to the
TREDataMapper
object.
rawdata: Set the raw data slot from a TREDataMapper
object.
value: The list
-like object to set for the rawdata slot. Note: this
currently only supports data.frame
or data.table
objects.
rowDataMap: list
of two character
vectors, the first are the
columns required to uniquely identify each row of a TREDataMapper
and the
second any additional row-level metadata. If the character vectors
have names, the resulting columns are automatically renamed to the
item name of the specified column.
rowDataMap<-: Update the @rowDataMap
slot of a TREDataMapper
object,
returning an invisible NULL. Arguments:
value: A list
or List
where the first item is the names of the
identifier columns – columns needed to uniquely identify each row in
rowData – and the second item is the metadata associated with those
the identifier columns, but not required to uniquely identify rows in
the object rowData.
colDataMap: list
of two character
vectors, the first are the
columns required to uniquely identify each row of a TREDataMapper
and the
second any additional col-level metadata. If the character vectors
have names, the resulting columns are automatically renamed to the
item name of the specified column.
colDataMap<-: Update the @colDataMap
slot of a TREDataMapper
object,
returning an invisible NULL. Arguments:
value: A list
or List
where the first item is the names of the
identifier columns – columns needed to uniquely identify each row in
colData – and the second item is the metadata associated with those
the identifier columns, but not required to uniquely identify rows in
the object rowData.
assayMap: A list
of character vectors. The name of each list item
will be the assay in a LongTableDataMapper
object that the columns in the
character
vector will be assigned to. Column renaming occurs automatically
when the character vectors have names (from the value to the name).
assayMap<-: Updates the @assayMap
slot of a TREDataMapper
object,
returning an invisible NULL. Arguments:
value: A list
of character vectors, where the name of each list
item is the name of an assay and the values of each character vector
specify the columns mapping to the assay in the S4
object the
TREDataMapper
constructs.
metadataMap: A list
of character
vectors. Each item is an element
of the constructed objects @metadata
slot.
metadataMap<-: Updates TREDataMapper
object in-place, then returns an
invisible(NULL)
. Arguments:
value: A list
of character
vectors. The name of each list item
is the name of the item in the @metadata
slot of the TREDataMapper
object
created when metaConstruct
is called on the DataMapper
, and a
character vector specifies the columns of @rawdata
to assign to each item.
Accessors: See details
Setters: An update TREDataMapper
object, returned invisibly.
Other DataMapper-accessors:
DataMapper-accessors
,
LongTableDataMapper-accessors
rowDataMap(exampleDataMapper) rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) assayMap(exampleDataMapper) <- list(sensitivity=c(viability1='viability')) metadataMap(exampleDataMapper) metadataMap(exampleDataMapper) <- list(object_metadata=c('metadata'))
rowDataMap(exampleDataMapper) rowDataMap(exampleDataMapper) <- list(c('treatmentid'), c()) colDataMap(exampleDataMapper) colDataMap(exampleDataMapper) <- list(c('sampleid'), c()) assayMap(exampleDataMapper) assayMap(exampleDataMapper) <- list(sensitivity=c(viability1='viability')) metadataMap(exampleDataMapper) metadataMap(exampleDataMapper) <- list(object_metadata=c('metadata'))
TreatmentResponseExperiment
ObjectA Class for Mapping Between Raw Data and an TreatmentResponseExperiment
Object
rawdata
See Slots section.
rowDataMap
See Slots section.
colDataMap
See Slots section.
assayMap
See Slots section.
metadataMap
See Slots section.
rowDataMap: A list-like object containing two character
vectors.
The first is column names in rawdata
needed to uniquely identify each row,
the second is additional columns which map to rows, but are not required to
uniquely identify them. Rows should be drugs.
colDataMap: A list-like object containing two character
vectors.
The first is column names in rawdata
needed to uniquely identify each
column, the second is additional columns which map to rows, but are not
required to uniquely identify them. Columns should be samples.
assayMap A list-like where each item is a list
with two elements
specifying an assay, the first being the identifier columns in rawdata
needed to uniquely identify each row an assay, and the second a list of
rawdata
columns to be mapped to that assay. The names of assayMap
will be the names of the assays in the LongTable
that is created when
calling metaConstruct
on this DataMapper
object.
metadataMap: A list-like where each item is a character
vector of
rawdata
column names to assign to the @metadata
of the LongTable
,
where the name of that assay is the name of the list item. If names are
omitted, assays will be numbered by their index in the list.
rawdata: A list-like object containing one or more pieces of raw data
that will be processed and mapped to the slots of an S4
object.
metadata: A List
of object level metadata.
CoreSet
class after changes in it struture or APIUpdate the CoreSet
class after changes in it struture or API
## S4 method for signature 'CoreSet' updateObject(object, verify = FALSE, verbose = FALSE)
## S4 method for signature 'CoreSet' updateObject(object, verify = FALSE, verbose = FALSE)
object |
A |
verify |
A |
verbose |
TRUE or FALSE, indicating whether information about the update should be reported |
CoreSet
with update class structure.
LongTable
class after changes in it struture or APIUpdate the LongTable
class after changes in it struture or API
## S4 method for signature 'LongTable' updateObject(object, verify = FALSE, verbose = FALSE)
## S4 method for signature 'LongTable' updateObject(object, verify = FALSE, verbose = FALSE)
object |
A |
verify |
A |
verbose |
TRUE or FALSE, indicating whether information about the update should be reported |
LongTable
with update class structure.