Title: | Analysis of patient-derived xenograft (PDX) data |
---|---|
Description: | The Xeva package provides efficient and powerful functions for patient-drived xenograft (PDX) based pharmacogenomic data analysis. This package contains a set of functions to perform analysis of patient-derived xenograft data. This package was developed by the BHKLab, for further information please see our documentation. |
Authors: | Arvind Mer [aut], Benjamin Haibe-Kains [aut, cre] |
Maintainer: | Benjamin Haibe-Kains <[email protected]> |
License: | GPL-3 |
Version: | 1.23.0 |
Built: | 2024-10-31 06:37:04 UTC |
Source: | https://github.com/bioc/Xeva |
area between curves Computes the area between two time-volume curves.
ABC( contr.time = NULL, contr.volume = NULL, treat.time = NULL, treat.volume = NULL )
ABC( contr.time = NULL, contr.volume = NULL, treat.time = NULL, treat.volume = NULL )
contr.time |
Time vector for control. |
contr.volume |
Volume vector for control. |
treat.time |
Time vector for treatment. |
treat.volume |
Volume vector for treatment. |
Returns batch response object.
contr.time <- treat.time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) contr.volume<- contr.time * tan(60*pi/180) treat.volume<- treat.time * tan(15*pi/180) abc <- ABC(contr.time, contr.volume, treat.time, treat.volume) par(pty="s") xylimit <- range(c(contr.time, contr.volume, treat.time, treat.volume)) plot(contr.time, contr.volume, type = "b", xlim = xylimit, ylim = xylimit) lines(treat.time, treat.volume, type = "b") polygon(c(treat.time, rev(treat.time)), c(contr.volume, rev(treat.volume)), col = "#fa9fb5", border = NA)
contr.time <- treat.time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) contr.volume<- contr.time * tan(60*pi/180) treat.volume<- treat.time * tan(15*pi/180) abc <- ABC(contr.time, contr.volume, treat.time, treat.volume) par(pty="s") xylimit <- range(c(contr.time, contr.volume, treat.time, treat.volume)) plot(contr.time, contr.volume, type = "b", xlim = xylimit, ylim = xylimit) lines(treat.time, treat.volume, type = "b") polygon(c(treat.time, rev(treat.time)), c(contr.volume, rev(treat.volume)), col = "#fa9fb5", border = NA)
Add a new experimental design in the expDesign
slot.
addExperimentalDesign( object, treatment = NULL, control = NULL, batch.id = NULL, replace = FALSE ) ## S4 method for signature 'XevaSet' addExperimentalDesign( object, treatment = NULL, control = NULL, batch.id = NULL, replace = FALSE )
addExperimentalDesign( object, treatment = NULL, control = NULL, batch.id = NULL, replace = FALSE ) ## S4 method for signature 'XevaSet' addExperimentalDesign( object, treatment = NULL, control = NULL, batch.id = NULL, replace = FALSE )
object |
The |
treatment |
The |
control |
The |
batch.id |
The |
replace |
If |
Returns Xeva
dataset with new experimental design added.
data(brca) brca <- addExperimentalDesign(object=brca, treatment=c("X.6047.LL71"), control=c("X.6047.uned"), batch.id="new.batch", replace=FALSE)
data(brca) brca <- addExperimentalDesign(object=brca, treatment=c("X.6047.LL71"), control=c("X.6047.uned"), batch.id="new.batch", replace=FALSE)
compute angle Computes the angle between two time-volume curves.
angle( contr.time = NULL, contr.volume = NULL, treat.time = NULL, treat.volume = NULL, degree = TRUE )
angle( contr.time = NULL, contr.volume = NULL, treat.time = NULL, treat.volume = NULL, degree = TRUE )
contr.time |
Time vector for control. |
contr.volume |
Volume vector for control. |
treat.time |
Time vector for treatment. |
treat.volume |
Volume vector for treatment. |
degree |
Default |
Returns batch response object.
contr.time <- treat.time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) contr.volume<- contr.time * tan(60*pi/180) treat.volume<- treat.time * tan(15*pi/180) ang <- angle(contr.time, contr.volume, treat.time, treat.volume) print(ang) par(pty="s") xylimit <- range(c(contr.time, contr.volume, treat.time, treat.volume)) plot(contr.time, contr.volume, type = "b", xlim = xylimit, ylim = xylimit) lines(treat.time, treat.volume, type = "b") abline(lm(contr.volume~contr.time)) abline(lm(treat.volume~treat.time))
contr.time <- treat.time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) contr.volume<- contr.time * tan(60*pi/180) treat.volume<- treat.time * tan(15*pi/180) ang <- angle(contr.time, contr.volume, treat.time, treat.volume) print(ang) par(pty="s") xylimit <- range(c(contr.time, contr.volume, treat.time, treat.volume)) plot(contr.time, contr.volume, type = "b", xlim = xylimit, ylim = xylimit) lines(treat.time, treat.volume, type = "b") abline(lm(contr.volume~contr.time)) abline(lm(treat.volume~treat.time))
AUC
Returns area under the curvearea under the curve
AUC
Returns area under the curve
AUC(time, volume)
AUC(time, volume)
time |
A |
volume |
First |
Returns angle
and slope
object.
time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) volume1<- time * tan(30*pi/180) volume2<- time * tan(45*pi/180) auc1 <- AUC(time, volume1) auc2 <- AUC(time, volume2) par(pty="s") xylimit <- range(c(time, volume1, volume2)) plot(time, volume1, type = "b", xlim = xylimit, ylim = xylimit) lines(time, volume2, type = "b") abline(lm(volume1~time)) abline(lm(volume2~time))
time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) volume1<- time * tan(30*pi/180) volume2<- time * tan(45*pi/180) auc1 <- AUC(time, volume1) auc2 <- AUC(time, volume2) par(pty="s") xylimit <- range(c(time, volume1, volume2)) plot(time, volume1, type = "b", xlim = xylimit, ylim = xylimit) lines(time, volume2, type = "b") abline(lm(volume1~time)) abline(lm(volume2~time))
Get batch information from a Xeva dataset.
batchInfo( object, batch = NULL, model.id = NULL, model.id.type = c("any", "control", "treatment") ) ## S4 method for signature 'XevaSet' batchInfo( object, batch = NULL, model.id = NULL, model.id.type = c("any", "control", "treatment") )
batchInfo( object, batch = NULL, model.id = NULL, model.id.type = c("any", "control", "treatment") ) ## S4 method for signature 'XevaSet' batchInfo( object, batch = NULL, model.id = NULL, model.id.type = c("any", "control", "treatment") )
object |
The Xeva object from which batch information is obtained. |
batch |
Name of the batch. Default |
model.id |
Model ID for which need to be searched in the batches. Default |
model.id.type |
Type of the model ID in a batch. See the Details section below. |
By default this function will return the names of all the batches present in the
dataset. If a batch specified, it will return the experiment design (control
and treatment model IDs) of that particular batch. If model.id
is specified,
it will return the names of all the batches where this particuler model.id
is present.
If both batch
and model.id
are specified, batch
will take precedent.
For model.id.type
, the default value 'any'
will return all the batch IDs
where the given model ID is present in any arm (ie. control or treatment) of the
batch. It can also be restricted to look only for treatment (or control) arm by
specifying the type.
A Vector
with batch names.
data(brca) ##to get all the batch names batch.name <- batchInfo(brca) ##to get a specific batch batch.design <- batchInfo(brca, batch=batch.name[1]) ##to get all the batches where a model.id is present batchInfo(brca, model.id="X.6047.uned")
data(brca) ##to get all the batch names batch.name <- batchInfo(brca) ##to get a specific batch batch.design <- batchInfo(brca, batch=batch.name[1]) ##to get all the batches where a model.id is present batchInfo(brca, model.id="X.6047.uned")
A Xeva object containing only breast cancer PDXs from the PDXE dataset For details about PDX-MI, see: Gao et al. High-throughput screening using patient-derived tumor xenografts to predict clinical trial drug response. Nature medicine, 21(11):1318, 2015.
data(brca)
data(brca)
An object of class XevaSet
of length 1.
https://www.nature.com/articles/nm.3954?draft=journal
A constructor to create XevaSet. Only objects returned by this constructor are expected to work with the XevaSet methods.
createXevaSet( name, model = data.frame(), drug = data.frame(), experiment = data.frame(), expDesign = list(), modelSensitivity = data.frame(), batchSensitivity = data.frame(), molecularProfiles = list(), modToBiobaseMap = data.frame() )
createXevaSet( name, model = data.frame(), drug = data.frame(), experiment = data.frame(), expDesign = list(), modelSensitivity = data.frame(), batchSensitivity = data.frame(), molecularProfiles = list(), modToBiobaseMap = data.frame() )
name |
A |
model |
A |
drug |
A |
experiment |
A |
expDesign |
A list containing name of the batch, control and treatment model.id |
modelSensitivity |
A |
batchSensitivity |
A |
molecularProfiles |
A |
modToBiobaseMap |
A |
This function creates a XevaSet object. It takes different model infromation and genomic data as input. For detailed discription of all varaibles please see Xeva vignette section "Creating new Xeva object"
Returns Xeva object
## read raw data files containg PDX experiment information and genomic data model = read.csv(system.file("extdata", "model.csv", package = "Xeva")) drug = read.csv(system.file("extdata", "drug.csv", package = "Xeva")) experiment= read.csv(system.file("extdata", "experiments.csv", package = "Xeva")) expDesign=readRDS(system.file("extdata", "batch_list.rds", package = "Xeva")) RNASeq=readRDS(system.file("extdata", "rnaseq.rds", package = "Xeva")) modToBiobaseMap=read.csv(system.file("extdata", "modelToExpressionMap.csv", package = "Xeva")) ## create Xeva object xeva.set = createXevaSet(name="example xevaSet", model=model, drug=drug, experiment=experiment, expDesign=expDesign, molecularProfiles=list(RNASeq = RNASeq), modToBiobaseMap = modToBiobaseMap) print(xeva.set)
## read raw data files containg PDX experiment information and genomic data model = read.csv(system.file("extdata", "model.csv", package = "Xeva")) drug = read.csv(system.file("extdata", "drug.csv", package = "Xeva")) experiment= read.csv(system.file("extdata", "experiments.csv", package = "Xeva")) expDesign=readRDS(system.file("extdata", "batch_list.rds", package = "Xeva")) RNASeq=readRDS(system.file("extdata", "rnaseq.rds", package = "Xeva")) modToBiobaseMap=read.csv(system.file("extdata", "modelToExpressionMap.csv", package = "Xeva")) ## create Xeva object xeva.set = createXevaSet(name="example xevaSet", model=model, drug=drug, experiment=experiment, expDesign=expDesign, molecularProfiles=list(RNASeq = RNASeq), modToBiobaseMap = modToBiobaseMap) print(xeva.set)
plot data for dose in model.id
dosePlot( object, model.id, max.time = NULL, treatment.only = FALSE, vol.normal = FALSE, concurrent.time = FALSE, point.shape = 21, point.size = 3, line.size = 4, point.color = "#878787", line.color = "#bababa", fill.col = c("#f5f5f5", "#E55100"), modify.x.axis = FALSE )
dosePlot( object, model.id, max.time = NULL, treatment.only = FALSE, vol.normal = FALSE, concurrent.time = FALSE, point.shape = 21, point.size = 3, line.size = 4, point.color = "#878787", line.color = "#bababa", fill.col = c("#f5f5f5", "#E55100"), modify.x.axis = FALSE )
object |
Xeva object. |
model.id |
one or multiple model.id |
max.time |
Maximum time point of the plot. Default |
treatment.only |
Default |
vol.normal |
Default |
concurrent.time |
Default |
point.shape |
shape of the point |
point.size |
size of the point |
line.size |
size of the line |
point.color |
color for point |
line.color |
color for line |
fill.col |
a vector with color to fill |
modify.x.axis |
Default |
A ggplot2 plot
data(brca) dosePlot(brca, model.id=c("X.6047.LJ16","X.6047.LJ16.trab"), fill.col=c("#f5f5f5", "#993404"))
data(brca) dosePlot(brca, model.id=c("X.6047.LJ16","X.6047.LJ16.trab"), fill.col=c("#f5f5f5", "#993404"))
This function allows you to see the available XevaSet object and download them for use with this package. The XevaSet have been extensively curated and organised within a XevaSet class, enabling use with all the analysis tools provided in Xeva.
downloadXevaSet( name = NULL, saveDir = file.path(".", "XevaSet"), XevaSetFileName = NULL, verbose = TRUE )
downloadXevaSet( name = NULL, saveDir = file.path(".", "XevaSet"), XevaSetFileName = NULL, verbose = TRUE )
name |
Character string, the name of the XevaSet to download. |
saveDir |
|
XevaSetFileName |
|
verbose |
|
A data.frame if name is NULL, showing all the available XevaSet objects. If name is specified, it will download the dataset from our server
downloadXevaSet() ##to download a dataset #library(Xeva) #PDXE_BRCA = downloadXevaSet(name="PDXE_BRCA", saveDir="XevaSet")
downloadXevaSet() ##to download a dataset #library(Xeva) #PDXE_BRCA = downloadXevaSet(name="PDXE_BRCA", saveDir="XevaSet")
Get drug information Get the drug information slot from a XevaSet object.
drugInform(object) ## S4 method for signature 'XevaSet' drugInform(object)
drugInform(object) ## S4 method for signature 'XevaSet' drugInform(object)
object |
The |
A data.frame
with the drug annotations.
data(brca) head(drugInform(brca))
data(brca) head(drugInform(brca))
Given a Xeva object and drug name, this function will return sensitivity values for all the genes/features.
drugSensitivitySig( object, drug, mDataType = NULL, molData = NULL, features = NULL, model.ids = NULL, model2bidMap = NULL, sensitivity.measure = "slope", fit = c("lm", "CI", "pearson", "spearman", NA), standardize = c("SD", "rescale", "none"), nthread = 1, tissue = NULL, verbose = TRUE )
drugSensitivitySig( object, drug, mDataType = NULL, molData = NULL, features = NULL, model.ids = NULL, model2bidMap = NULL, sensitivity.measure = "slope", fit = c("lm", "CI", "pearson", "spearman", NA), standardize = c("SD", "rescale", "none"), nthread = 1, tissue = NULL, verbose = TRUE )
object |
The |
drug |
Name of the drug. |
mDataType |
Molecular data type. |
molData |
External data matrix. Rows as features and columns as samples. |
features |
Set which molecular data features to use. Default |
model.ids |
Set which |
model2bidMap |
A |
sensitivity.measure |
Name of the sensitivity measure. |
fit |
Association method to use, can be 'lm', 'CI', 'pearson' or 'spearman'. If 'NA' only the data will be return. Default |
standardize |
Default |
nthread |
number of threads |
tissue |
tissue type. Default |
verbose |
Default |
Method to compute association can be specified by fit
. It can be one of the:
"lm" for linear models
"CI" for concordance index
"pearson" for Pearson correlation
"spearman" for Spearman correlation
If fit is set to NA, processed data (an ExpressionSet) will be returned.
A matrix of values can be directly passed to molData.
In case where a model.id
maps to multiple biobase.id
s, the first biobase.id
in the data.frame
will be used.
A data.frame
with features and values.
data(brca) senSig <- drugSensitivitySig(object=brca, drug="tamoxifen", mDataType="RNASeq", features=c(1,2,3,4,5), sensitivity.measure="slope", fit = "lm") ## example to compute the Pearson correlation between gene expression and PDX response senSig <- drugSensitivitySig(object=brca, drug="tamoxifen", mDataType="RNASeq", features=c(1,2,3,4,5), sensitivity.measure="slope", fit = "pearson")
data(brca) senSig <- drugSensitivitySig(object=brca, drug="tamoxifen", mDataType="RNASeq", features=c(1,2,3,4,5), sensitivity.measure="slope", fit = "lm") ## example to compute the Pearson correlation between gene expression and PDX response senSig <- drugSensitivitySig(object=brca, drug="tamoxifen", mDataType="RNASeq", features=c(1,2,3,4,5), sensitivity.measure="slope", fit = "pearson")
For a given model.id
, getExperiment
will
getExperiment( object, model.id = NULL, batch = NULL, patient.id = NULL, drug = NULL, control.name = NULL, treatment.only = FALSE, max.time = NULL, vol.normal = FALSE, log.volume = FALSE, return.list = FALSE, impute.value = FALSE, concurrent.time = FALSE ) ## S4 method for signature 'XevaSet' getExperiment( object, model.id = NULL, batch = NULL, patient.id = NULL, drug = NULL, control.name = NULL, treatment.only = FALSE, max.time = NULL, vol.normal = FALSE, log.volume = FALSE, return.list = FALSE, impute.value = FALSE, concurrent.time = FALSE )
getExperiment( object, model.id = NULL, batch = NULL, patient.id = NULL, drug = NULL, control.name = NULL, treatment.only = FALSE, max.time = NULL, vol.normal = FALSE, log.volume = FALSE, return.list = FALSE, impute.value = FALSE, concurrent.time = FALSE ) ## S4 method for signature 'XevaSet' getExperiment( object, model.id = NULL, batch = NULL, patient.id = NULL, drug = NULL, control.name = NULL, treatment.only = FALSE, max.time = NULL, vol.normal = FALSE, log.volume = FALSE, return.list = FALSE, impute.value = FALSE, concurrent.time = FALSE )
object |
The |
model.id |
The |
batch |
Batch name from the |
patient.id |
Patient id from the |
drug |
Name of the drug. |
control.name |
Name of drug used as control. Default |
treatment.only |
Default |
max.time |
Maximum time for data. |
vol.normal |
If TRUE it will normalize the volume. Default |
log.volume |
If TRUE log of the volume will be used. Default |
return.list |
Default |
impute.value |
Default |
concurrent.time |
Default |
a data.fram
will all the the values stored in experiment slot
data(brca) getExperiment(brca, model.id="X.6047.uned", treatment.only=TRUE) getExperiment(brca, model.id=c("X.6047.uned", "X.6047.pael"), treatment.only=TRUE) getExperiment(brca, batch="X-6047.paclitaxel", treatment.only=TRUE) ed <- list(batch.name="myBatch", treatment=c("X.6047.LJ16","X.6047.LJ16.trab"), control=c("X.6047.uned")) getExperiment(brca, batch=ed)
data(brca) getExperiment(brca, model.id="X.6047.uned", treatment.only=TRUE) getExperiment(brca, model.id=c("X.6047.uned", "X.6047.pael"), treatment.only=TRUE) getExperiment(brca, batch="X-6047.paclitaxel", treatment.only=TRUE) ed <- list(batch.name="myBatch", treatment=c("X.6047.LJ16","X.6047.LJ16.trab"), control=c("X.6047.uned")) getExperiment(brca, batch=ed)
This function serves to get molecular profiles from a XevaSet
object.
getMolecularProfiles(object, data.type)
getMolecularProfiles(object, data.type)
object |
The |
data.type |
|
An ExpressionSet
where sample names are the biobase.id
of the model.
data(brca) brca.RNA <- getMolecularProfiles(brca, data.type="RNASeq")
data(brca) brca.RNA <- getMolecularProfiles(brca, data.type="RNASeq")
Comput the linear mixed model (lmm) statistics for a PDX batch
lmm(data)
lmm(data)
data |
a data.frame containg a batch data |
The input data.frame (data) must contain these columns: model.id, volume, time, exp.type
Returns a fit object
data(repdx) data <- getExperiment(repdx, batch = "P1")$model lmm(data)
data(repdx) data <- getExperiment(repdx, batch = "P1")$model lmm(data)
modelInfo Generic Generic for modelInfo method
modelInfo(object, mDataType = NULL) ## S4 method for signature 'XevaSet' modelInfo(object, mDataType = NULL)
modelInfo(object, mDataType = NULL) ## S4 method for signature 'XevaSet' modelInfo(object, mDataType = NULL)
object |
Xeva object |
mDataType |
Molecular data type. |
A data.frame
with the model annotations.
data(brca) mid <- modelInfo(brca) head(mid)
data(brca) mid <- modelInfo(brca) head(mid)
mRECIST
Returns the mRECIST for given volume response.
mRECIST(time, volume, min.time = 10, return.detail = FALSE)
mRECIST(time, volume, min.time = 10, return.detail = FALSE)
time |
Value of best response. |
volume |
Value of best average response. |
min.time |
Minimum time after which tumor volume will be considered. |
return.detail |
Default |
Returns the mRECIST.
time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) volume<- c(250.8, 320.4, 402.3, 382.6, 384, 445.9, 460.2, 546.8, 554.3, 617.9) mRECIST(time, volume, min.time=10, return.detail=FALSE)
time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) volume<- c(250.8, 320.4, 402.3, 382.6, 384, 445.9, 460.2, 546.8, 554.3, 617.9) mRECIST(time, volume, min.time=10, return.detail=FALSE)
A dataset containing PDX models minimal information (PDX-MI) standard and corresponding Xeva variable.
data(PDXMI)
data(PDXMI)
An object of class data.frame
with 45 rows and 4 columns.
For details about PDX-MI, see:
Meehan, Terrence F., et al. "PDX-MI: minimal information for patient-derived tumor xenograft models." Cancer research 77.21 (2017): e62-e66.
http://cancerres.aacrjournals.org/lookup/doi/10.1158/0008-5472.CAN-17-0582
plotmRECIST
plots the mRECIST matrix obtained from summarizeResponse
.
plotmRECIST( mat, control.name = NA, control.col = "#238b45", drug.col = "black", colPalette = NULL, name = "Drug & Models", sort = TRUE, row_fontsize = 12, col_fontsize = 12, draw_plot = TRUE )
plotmRECIST( mat, control.name = NA, control.col = "#238b45", drug.col = "black", colPalette = NULL, name = "Drug & Models", sort = TRUE, row_fontsize = 12, col_fontsize = 12, draw_plot = TRUE )
mat |
The mRECIST matrix where rows are drugs and columns are patients. |
control.name |
Name of the control. |
control.col |
Color of the control. |
drug.col |
Color of the drug names. |
colPalette |
Color palette for mRECIST values. |
name |
Title of the plot. |
sort |
If matrix should be sorted before plotting. |
row_fontsize |
Size of the row name font. |
col_fontsize |
Size of the column name font. |
draw_plot |
Default |
mRECIST plot.
data(brca) brca.mr <- summarizeResponse(brca, response.measure = "mRECIST", group.by="patient.id") plotmRECIST(as.matrix(brca.mr), control.name = "untreated")
data(brca) brca.mr <- summarizeResponse(brca, response.measure = "mRECIST", group.by="patient.id") plotmRECIST(as.matrix(brca.mr), control.name = "untreated")
Plot data for a batch.id, experiment design or model.id
plotPDX( object, batch = NULL, patient.id = NULL, drug = NULL, model.id = NULL, model.color = NULL, control.name = NULL, max.time = NULL, treatment.only = FALSE, vol.normal = FALSE, impute.value = TRUE, concurrent.time = FALSE, control.col = "#e41a1c", treatment.col = "#377eb8", title = "", xlab = "Time", ylab = "Volume", log.y = FALSE, SE.plot = c("all", "none", "errorbar", "ribbon"), aspect.ratio = c(1, NULL), minor.line.size = 0.5, major.line.size = 0.7 ) plotBatch( object, batch = NULL, patient.id = NULL, drug = NULL, control.name = NULL, max.time = NULL, treatment.only = FALSE, vol.normal = FALSE, impute.value = TRUE, concurrent.time = FALSE, control.col = "#6baed6", treatment.col = "#fc8d59", title = "", xlab = "Time", ylab = "Volume", log.y = FALSE, SE.plot = c("all", "none", "errorbar", "ribbon"), aspect.ratio = c(1, NULL), minor.line.size = 0.5, major.line.size = 0.7 )
plotPDX( object, batch = NULL, patient.id = NULL, drug = NULL, model.id = NULL, model.color = NULL, control.name = NULL, max.time = NULL, treatment.only = FALSE, vol.normal = FALSE, impute.value = TRUE, concurrent.time = FALSE, control.col = "#e41a1c", treatment.col = "#377eb8", title = "", xlab = "Time", ylab = "Volume", log.y = FALSE, SE.plot = c("all", "none", "errorbar", "ribbon"), aspect.ratio = c(1, NULL), minor.line.size = 0.5, major.line.size = 0.7 ) plotBatch( object, batch = NULL, patient.id = NULL, drug = NULL, control.name = NULL, max.time = NULL, treatment.only = FALSE, vol.normal = FALSE, impute.value = TRUE, concurrent.time = FALSE, control.col = "#6baed6", treatment.col = "#fc8d59", title = "", xlab = "Time", ylab = "Volume", log.y = FALSE, SE.plot = c("all", "none", "errorbar", "ribbon"), aspect.ratio = c(1, NULL), minor.line.size = 0.5, major.line.size = 0.7 )
object |
Xeva object. |
batch |
Batch name or experiment design list. |
patient.id |
Patient id from the |
drug |
Name of the drug. Default |
model.id |
One or multiple model.id. Default |
model.color |
Color for |
control.name |
Name of the control sample. |
max.time |
Maximum time point of the plot. Default |
treatment.only |
Default |
vol.normal |
Default |
impute.value |
Default |
concurrent.time |
Default |
control.col |
Color for control plots. |
treatment.col |
Color for treatment plots. |
title |
Title of the plot. |
xlab |
Title of the x-axis. |
ylab |
Title of the y-axis. |
log.y |
Default |
SE.plot |
Plot type. Default |
aspect.ratio |
Default |
minor.line.size |
Line size for minor lines. Default |
major.line.size |
Line size for major lines. Default |
A ggplot2 plot with control and treatment batch data.
data(brca) plotPDX(brca, model.id=c("X.6047.LJ16","X.6047.LJ16.trab")) plotPDX(brca, batch="X-1004.BGJ398", vol.normal=TRUE) expDesign <- list(batch.name="myBatch", treatment=c("X.6047.LJ16","X.6047.LJ16.trab"), control=c("X.6047.uned")) plotBatch(brca, batch=expDesign, vol.normal=TRUE) plotBatch(brca, batch=expDesign, vol.normal=FALSE, SE.plot = "errorbar")
data(brca) plotPDX(brca, model.id=c("X.6047.LJ16","X.6047.LJ16.trab")) plotPDX(brca, batch="X-1004.BGJ398", vol.normal=TRUE) expDesign <- list(batch.name="myBatch", treatment=c("X.6047.LJ16","X.6047.LJ16.trab"), control=c("X.6047.uned")) plotBatch(brca, batch=expDesign, vol.normal=TRUE) plotBatch(brca, batch=expDesign, vol.normal=FALSE, SE.plot = "errorbar")
Print the batch response
## S3 method for class 'batchResponse' print(x, ...)
## S3 method for class 'batchResponse' print(x, ...)
x |
batchResponse object |
... |
Other arguments |
prints the batchResponse
Print the model response
## S3 method for class 'modelResponse' print(x, ...)
## S3 method for class 'modelResponse' print(x, ...)
x |
modelResponse object |
... |
Other arguments |
prints the modelResponse
Print the pdx batch
## S3 method for class 'pdxBatch' print(x, ...)
## S3 method for class 'pdxBatch' print(x, ...)
x |
pdxBatch object |
... |
Other arguments |
prints pdxBatch
A Xeva object containing anonymous PDX data with replicates. Each batch has 5 replicates.
data(repdx)
data(repdx)
An object of class XevaSet
of length 1.
response
Computes the drug response of an individual PDX model or batch.
response( object, model.id = NULL, batch = NULL, res.measure = c("mRECIST", "slope", "AUC", "angle", "abc", "TGI", "lmm"), treatment.only = FALSE, max.time = NULL, impute.value = TRUE, min.time = 10, concurrent.time = TRUE, vol.normal = FALSE, log.volume = FALSE, verbose = TRUE )
response( object, model.id = NULL, batch = NULL, res.measure = c("mRECIST", "slope", "AUC", "angle", "abc", "TGI", "lmm"), treatment.only = FALSE, max.time = NULL, impute.value = TRUE, min.time = 10, concurrent.time = TRUE, vol.normal = FALSE, log.volume = FALSE, verbose = TRUE )
object |
Xeva object. |
model.id |
|
batch |
|
res.measure |
Drug response measure. See |
treatment.only |
Default |
max.time |
Maximum time for data. |
impute.value |
Default |
min.time |
Default 10 days. Used for mRECIST computation. |
concurrent.time |
Default |
vol.normal |
If TRUE it will normalize the volume. Default |
log.volume |
If TRUE log of the volume will be used for response calculation. Default |
verbose |
Default |
At present the following response measures are implemented
mRECIST Computes mRECIST for individual PDX models
slope Computes slope of the fitted individual PDX curves
AUC Computes area under a PDX curve for individual PDX models
angle Computes angle between treatment and control PDX curves
abc Computes area between the treatment and control PDX curves
TGI Computes tumor growth inhibition using treatment and control PDX curves
lmm Computes linear mixed model (lmm) statistics for a PDX batch
Returns model or batch drug response object.
data(brca) response(brca, model.id="X.1004.BG98", res.measure="mRECIST") response(brca, batch="X-6047.paclitaxel", res.measure="angle") ed <- list(batch.name="myBatch", treatment=c("X.6047.LJ16","X.6047.LJ16.trab"), control=c("X.6047.uned")) response(brca, batch=ed, res.measure="angle")
data(brca) response(brca, model.id="X.1004.BG98", res.measure="mRECIST") response(brca, batch="X-6047.paclitaxel", res.measure="angle") ed <- list(batch.name="myBatch", treatment=c("X.6047.LJ16","X.6047.LJ16.trab"), control=c("X.6047.uned")) response(brca, batch=ed, res.measure="angle")
To select model IDs based on drug name and/or tissue type.
selectModelIds(object, drug = NULL, drug.match.exact = TRUE, tissue = NULL) ## S4 method for signature 'XevaSet' selectModelIds(object, drug = NULL, drug.match.exact = TRUE, tissue = NULL)
selectModelIds(object, drug = NULL, drug.match.exact = TRUE, tissue = NULL) ## S4 method for signature 'XevaSet' selectModelIds(object, drug = NULL, drug.match.exact = TRUE, tissue = NULL)
object |
The |
drug |
Name of the |
drug.match.exact |
Default |
tissue |
Tumor type. Default |
A vector
with the matched model.id
s.
data(brca) df = selectModelIds(brca, drug="trastuzumab", drug.match.exact=TRUE, tissue="BRCA") head(df) df2 = selectModelIds(brca, drug="trastuzumab", drug.match.exact=FALSE) head(df2)
data(brca) df = selectModelIds(brca, drug="trastuzumab", drug.match.exact=TRUE, tissue="BRCA") head(df) df2 = selectModelIds(brca, drug="trastuzumab", drug.match.exact=FALSE) head(df2)
Given a Xeva object, it will return a data.frame
detailing sensitivity information.
sensitivity(object, type = c("model", "batch"), sensitivity.measure = NULL)
sensitivity(object, type = c("model", "batch"), sensitivity.measure = NULL)
object |
The |
type |
Sensitivity type (either model or batch). |
sensitivity.measure |
Name of the |
A data.frame
with model or batch ID and sensitivity values.
data(brca) head(sensitivity(brca, type="batch")) head(sensitivity(brca, type="model"))
data(brca) head(sensitivity(brca, type="batch")) head(sensitivity(brca, type="model"))
setResponse
sets response of all PDXs in an Xeva object.
setResponse( object, res.measure = c("mRECIST", "slope", "AUC", "angle", "abc", "TGI", "lmm"), min.time = 10, treatment.only = FALSE, max.time = NULL, vol.normal = FALSE, impute.value = TRUE, concurrent.time = TRUE, log.volume = FALSE, verbose = TRUE )
setResponse( object, res.measure = c("mRECIST", "slope", "AUC", "angle", "abc", "TGI", "lmm"), min.time = 10, treatment.only = FALSE, max.time = NULL, vol.normal = FALSE, impute.value = TRUE, concurrent.time = TRUE, log.volume = FALSE, verbose = TRUE )
object |
Xeva object. |
res.measure |
Response measure, multiple measures are allowed. See |
min.time |
Minimum number of days for mRECIST computation. Default 10 days. |
treatment.only |
Default |
max.time |
Maximum number of days to consider for analysis. Data byond this will be discarded. Default |
vol.normal |
If TRUE it will will normalize the volume. Default |
impute.value |
Default |
concurrent.time |
Default |
log.volume |
If TRUE log of the volume will be used for response calculation. Default |
verbose |
Default |
At present fellowing response measure are implemented
mRECIST Computes mRECIST for indivial PDX model
slope Computes slope of the fitted indivial PDX curve
AUC Computes area under a PDX curve for indivial PDX model
angle Computes angle between treatment and control PDX curves
abc Computes area between the treatment and control PDX curves
TGI Computes tumor growth inhibition using treatment and control PDX curves
lmm Computes linear mixed model (lmm) statistics for a PDX batch
Returns updated Xeva object.
data(brca) brca <- setResponse(brca, res.measure = c("mRECIST"), verbose=FALSE)
data(brca) brca <- setResponse(brca, res.measure = c("mRECIST"), verbose=FALSE)
slope
returns the slope for given time and volume data.
slope(time, volume, degree = TRUE)
slope(time, volume, degree = TRUE)
time |
A |
volume |
A |
degree |
Default |
Returns the slope and a fit
object.
time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) volume<- c(250.8, 320.4, 402.3, 382.6, 384, 445.9, 460.2, 546.8, 554.3, 617.9) sl <- slope(time, volume) par(pty="s") xylimit <- range(c(time, volume)) plot(time, volume, type = "b", xlim = xylimit, ylim = xylimit) abline(lm(volume~time))
time <- c(0, 3, 7, 11, 18, 22, 26, 30, 32, 35) volume<- c(250.8, 320.4, 402.3, 382.6, 384, 445.9, 460.2, 546.8, 554.3, 617.9) sl <- slope(time, volume) par(pty="s") xylimit <- range(c(time, volume)) plot(time, volume, type = "b", xlim = xylimit, ylim = xylimit) abline(lm(volume~time))
Subset Xeva object.
subsetXeva(object, ids, id.name, keep.batch = TRUE)
subsetXeva(object, ids, id.name, keep.batch = TRUE)
object |
The |
ids |
IDs to be selected for. |
id.name |
Names of the IDs. |
keep.batch |
Default |
New Xeva object.
data(brca) print(brca) df <- subsetXeva(brca, ids = c("X-1004", "X-1008", "X-1286"), id.name="patient.id", keep.batch=TRUE) print(df)
data(brca) print(brca) df <- subsetXeva(brca, ids = c("X-1004", "X-1008", "X-1286"), id.name="patient.id", keep.batch=TRUE) print(df)
This function serves to get molecular profiles from a XevaSet
object.
summarizeMolecularProfiles( object, drug, mDataType, tissue = NULL, sensitivity.measure = NULL, unique.model = TRUE, batch = NULL )
summarizeMolecularProfiles( object, drug, mDataType, tissue = NULL, sensitivity.measure = NULL, unique.model = TRUE, batch = NULL )
object |
The |
drug |
Name of the drug. |
mDataType |
|
tissue |
Default |
sensitivity.measure |
Default |
unique.model |
Default |
batch |
Name of the batch. Default |
If a sequencing sample belongs to multiple models, summarizeMolecularProfiles
will create a separate column for each model.
All models without molecular data will be removed from the output ExpressionSet
.
An ExpressionSet
where sample names are model.id
and sensitivity measures will be presented in pData
.
data(brca) pacRNA <- summarizeMolecularProfiles(brca, drug="paclitaxel", mDataType="RNASeq", tissue= "BRCA", sensitivity.measure="mRECIST") print(pacRNA)
data(brca) pacRNA <- summarizeMolecularProfiles(brca, drug="paclitaxel", mDataType="RNASeq", tissue= "BRCA", sensitivity.measure="mRECIST") print(pacRNA)
This function summarizes the drug response information of PDXs.
summarizeResponse( object, response.measure = "mRECIST", model.id = NULL, batch.id = NULL, group.by = "patient.id", summary.stat = c(";", "mean", "median"), tissue = NULL )
summarizeResponse( object, response.measure = "mRECIST", model.id = NULL, batch.id = NULL, group.by = "patient.id", summary.stat = c(";", "mean", "median"), tissue = NULL )
object |
The |
response.measure |
|
model.id |
The |
batch.id |
A |
group.by |
Default |
summary.stat |
Dictates which summary method to use if multiple IDs are found. |
tissue |
Name of the tissue. Default |
There can be two types of drug response measure.
Per model response: One response value for each Model, eg. mRECIST_recomputed
for each model.
Per batch response: One response value for each Batch, eg. angle
between treatment and control groups.
For the per model response
output, columns will be model.id
(or group.by
).
For the per batch response
output, the group.by
value can be "batch.name"
.
A matrix
with rows as drug names, column as group.by
. Each cell contains response.measure
for the pair.
data(brca) brca.mR <- summarizeResponse(brca, response.measure = "mRECIST", group.by="patient.id")
data(brca) brca.mR <- summarizeResponse(brca, response.measure = "mRECIST", group.by="patient.id")
tumor growth inhibition (TGI) Computes the tumor growth inhibition (TGI) between two time-volume curves
TGI(contr.volume, treat.volume)
TGI(contr.volume, treat.volume)
contr.volume |
Volume vector for control |
treat.volume |
Volume vector for treatment |
Returns batch response object
contr.volume <- c(1.35, 6.57, 13.94, 20.39, 32.2, 39.26, 46.9, 53.91) treat.volume <- c(0.4, 1.26, 2.59, 3.62, 5.77, 6.67, 7.47, 8.98, 9.29, 9.44) TGI(contr.volume, treat.volume)
contr.volume <- c(1.35, 6.57, 13.94, 20.39, 32.2, 39.26, 46.9, 53.91) treat.volume <- c(0.4, 1.26, 2.59, 3.62, 5.77, 6.67, 7.47, 8.98, 9.29, 9.44) TGI(contr.volume, treat.volume)
waterfall plot Creates waterfall plot for a given drug.
waterfall( object, res.measure, drug = NULL, group.by = NULL, summary.stat = c(";", "mean", "median"), tissue = NULL, model.id = NULL, model.type = NULL, type.color = "#cc4c02", legend.name = NULL, yname = NULL, title = NULL, sort = TRUE )
waterfall( object, res.measure, drug = NULL, group.by = NULL, summary.stat = c(";", "mean", "median"), tissue = NULL, model.id = NULL, model.type = NULL, type.color = "#cc4c02", legend.name = NULL, yname = NULL, title = NULL, sort = TRUE )
object |
The |
res.measure |
PDX model drug response measure |
drug |
Name of the drug |
group.by |
Group drug response data |
summary.stat |
How to summarize multiple values |
tissue |
Tissue type |
model.id |
Indicates which |
model.type |
Type of model, such as mutated or wild type |
type.color |
A list with colors used for each type in the legend |
legend.name |
Name of the legend |
yname |
Name for the y-axis |
title |
Title of the plot |
sort |
Default |
waterfall plot in ggplot2
data(brca) waterfall(brca, drug="binimetinib", res.measure="best.avg.response_published") ## example with model.type where we color the models by TP53 mutation type mut <- summarizeMolecularProfiles(brca,drug = "binimetinib", mDataType="mutation") model.type <- Biobase::exprs(mut)["TP53", ] waterfall(brca, drug="binimetinib", res.measure="best.avg.response_published", tissue="BRCA", model.id=names(model.type), model.type= model.type)
data(brca) waterfall(brca, drug="binimetinib", res.measure="best.avg.response_published") ## example with model.type where we color the models by TP53 mutation type mut <- summarizeMolecularProfiles(brca,drug = "binimetinib", mDataType="mutation") model.type <- Biobase::exprs(mut)["TP53", ] waterfall(brca, drug="binimetinib", res.measure="best.avg.response_published", tissue="BRCA", model.id=names(model.type), model.type= model.type)