Package 'Rtpca'

Title: Thermal proximity co-aggregation with R
Description: R package for performing thermal proximity co-aggregation analysis with thermal proteome profiling datasets to analyse protein complex assembly and (differential) protein-protein interactions across conditions.
Authors: Nils Kurzawa [aut, cre], André Mateus [aut], Mikhail M. Savitski [aut]
Maintainer: Nils Kurzawa <[email protected]>
License: GPL-3
Version: 1.15.0
Built: 2024-07-24 05:30:50 UTC
Source: https://github.com/bioc/Rtpca

Help Index


Extract CommonFeatures

Description

Extract CommonFeatures

Usage

## S4 method for signature 'tpcaResult'
CommonFeatures(object)

Arguments

object

and object of class tpcaResult

Value

a vector of common features across replicates

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)

ppi_anno <- tibble(
    x = "2",
    y = "3",
    combined_score = 700,
    pair = "2:3")

tpcaObj <- runTPCA(
    objList = mat_list,
    complexAnno = NULL,
    ppiAnno = ppi_anno
)

CommonFeatures(tpcaObj)

Extract ComplexAnnotation

Description

Extract ComplexAnnotation

Usage

## S4 method for signature 'tpcaResult'
ComplexAnnotation(object)

Arguments

object

and object of class tpcaResult

Value

a data frame containing the complex annotation

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
ComplexAnnotation(tpcaObj)

Extract ComplexBackgroundDistributionList

Description

Extract ComplexBackgroundDistributionList

Usage

## S4 method for signature 'tpcaResult'
ComplexBackgroundDistributionList(object)

Arguments

object

and object of class tpcaResult

Value

a list of data frames containing distances of random complexes with different number of members

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
ComplexBackgroundDistributionList(tpcaObj)

Extract ComplexRocTable

Description

Extract ComplexRocTable

Usage

## S4 method for signature 'tpcaResult'
ComplexRocTable(object)

Arguments

object

and object of class tpcaResult

Value

a data frame containing a complex analysis roc table

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
ComplexRocTable(tpcaObj)

Extract ContrastCondName

Description

Extract ContrastCondName

Usage

## S4 method for signature 'tpcaResult'
ContrastCondName(object)

Arguments

object

and object of class tpcaResult

Value

a character string describing the contrast condition

Examples

tpcaObj <- new("tpcaResult")
ContrastCondName(tpcaObj)

Extract ContrastDistMat

Description

Extract ContrastDistMat

Usage

## S4 method for signature 'tpcaResult'
ContrastDistMat(object)

Arguments

object

an object of class tpcaResult

Value

a matrix containing the constrast distance matrix of all pariwise protein-protein melting curve distances computed from a TPP experiment

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
ContrastDistMat(tpcaObj)

Extract ContrastList

Description

Extract ContrastList

Usage

## S4 method for signature 'tpcaResult'
ContrastList(object)

Arguments

object

an object of class tpcaResult

Value

an object list containing TPP data

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
ContrastList(tpcaObj)

Create distance matrix of all vs all protein melting profiles

Description

Create distance matrix of all vs all protein melting profiles

Usage

createDistMat(
  objList,
  rownameCol = NULL,
  summaryMethodStr = "median",
  distMethodStr = "euclidean"
)

Arguments

objList

list of objects suitable for the analysis, currently allowed classes of objects are: matrices, data.frames, tibbles and ExpressionSets

rownameCol

in case the input objects are tibbles this parameter takes in the name (character) of the column specifying protein names or ids

summaryMethodStr

character string indicating a method to use to summarize measurements across replicates, default is "median", other options are c("mean", "rbind")

distMethodStr

method to use within dist function, default is 'euclidean'

Value

a distance matrix of all pairwise protein melting profiles

Examples

library(Biobase)

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

colnames(m1) <- paste0("X", 1:4)
colnames(m2) <- paste0("X", 1:4)
colnames(m3) <- paste0("X", 1:4)

mat_list <- list(
    m1, m2, m3
)

createDistMat(mat_list)

expr1 <- ExpressionSet(m1)
expr2 <- ExpressionSet(m2)
expr3 <- ExpressionSet(m3)

exprSet_list <- list(
    expr1, expr2, expr3
)

createDistMat(exprSet_list)

Extract CtrlCondName

Description

Extract CtrlCondName

Usage

## S4 method for signature 'tpcaResult'
CtrlCondName(object)

Arguments

object

and object of class tpcaResult

Value

a character string describing the control condition

Examples

tpcaObj <- new("tpcaResult")
CtrlCondName(tpcaObj)

Extract diffTpcaResultTable

Description

Extract diffTpcaResultTable

Usage

## S4 method for signature 'tpcaResult'
diffTpcaResultTable(object)

Arguments

object

an object of class tpcaResult

Value

a data frame containing the results from a diffTpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
diffTpcaResultTable(tpcaObj)

Extract DistMat

Description

Extract DistMat

Usage

## S4 method for signature 'tpcaResult'
DistMat(object)

Arguments

object

an object of class tpcaResult

Value

a matrix containing the distance matrix of all pariwise protein-protein melting curve distances computed from a TPP experiment

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
DistMat(tpcaObj)

Extract DistMethod

Description

Extract DistMethod

Usage

## S4 method for signature 'tpcaResult'
DistMethod(object)

Arguments

object

and object of class tpcaResult

Value

a character string of the dist method

Examples

tpcaObj <- new("tpcaResult")
DistMethod(tpcaObj)

Extract ObjList

Description

Extract ObjList

Usage

## S4 method for signature 'tpcaResult'
ObjList(object)

Arguments

object

an object of class tpcaResult

Value

an object list containing TPP data

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
ObjList(tpcaObj)

Data frame of eukaryotic protein-protein interactions inferred from annotated protein complexes by Ori et al. and StringDB interations with a combined score of at least 900

Description

data frame assigning proteins to (in)directly interacting proteins within protein complexes

Usage

data("ori_et_al_complex_ppis")

Format

data frame with columns complex_name, x, y, pair (unique pair id)

References

Ori et al. (2016), Genome Biology, 17, 47; Jensen et al. (2009), Nucleic Acids Research, 37, D412–D416

Examples

data("ori_et_al_complex_ppis")

Data frame of annotated protein complexes by Ori et al.

Description

data frame assigning proteins to annotated protein complexes

Usage

data("ori_et_al_complexes_df")

Format

data frame with columns ensembl_id, protein and id (complex identifier)

References

Ori et al. (2016), Genome Biology, 17, 47

Examples

data("ori_et_al_complexes_df")

Plot Complex ROC curve

Description

Plots a ROC curve representing how well a given TPP dataset recovers annotated proteins complexes. The ROC curve is generated based on the supplied protein complex annotation specificity is assessed by comparing the given complex annotation to random permutations of that table, i.e. proteins randomly grouped together.

Usage

plotComplexRoc(tpcaObj, computeAUC = FALSE)

Arguments

tpcaObj

tpcaResult object

computeAUC

logical parameter indicating whether area under the ROC should be computed and indicated in the lower right corner of the plot

Value

ggplot object of a receiver operating curve (ROC)

Examples

rocTab = data.frame(
    TPR = c(0, 0.1, 0.2, 0.4, 0.5, 0.7, 0.9, 1),
    FPR = c(0, 0.05, 0.1, 0.2, 0.5, 0.7, 0.9, 1)
)

tpcaTest <- new(
    "tpcaResult",
    ComplexRocTable = rocTab)

plotComplexRoc(tpcaTest)

Plot differential TPCA analysis results

Description

Plot differential TPCA analysis results

Usage

plotDiffTpcaVolcano(
  tpcaObj,
  alpha = 0.1,
  setXLim = FALSE,
  xlimit = c(-0.75, 0.75)
)

Arguments

tpcaObj

a tpcaObj after having performed a differential analysis, see runDiffTPCA

alpha

significance level / FDR at which null hypothesis should be rejected

setXLim

logical determining whether x-axis limits should be set according to xlimit

xlimit

numeric vector with two elements determing the x-axis limits, only is implemented if setXLim is set to TRUE

Value

ggplot displaying a volcano plot

Examples

library(dplyr)
library(Biobase)

m1 <- matrix(1:28, ncol = 4)
m2 <- matrix(2:25, ncol = 4)
m3 <- matrix(c(2:10, 1:19), ncol = 4)

rownames(m1) <- 1:7
rownames(m2) <- 3:8
rownames(m3) <- 2:8

mat_list <- list(
    m1, m2, m3
)

c1 <- matrix(29:2, ncol = 4)
c2 <- matrix(26:3, ncol = 4)
c3 <- matrix(c(11:3, 20:2), ncol = 4)

rownames(c1) <- 1:7
rownames(c2) <- 3:8
rownames(c3) <- 2:8

contrast_list <- list(
    c1, c2, c3
)

ppi_anno <- tibble(
    x = c("3", "3"),
    y = c("5", "7"),
    pair = c("3:5", "3:7"))

ref_df <- tibble(
    pair = c("3:5", "3:7"),
    valueC2 = c(4, 8)
)

diff_tpca <- runDiffTPCA(
  mat_list, contrast_list, ppiAnno = ppi_anno)

plotDiffTpcaVolcano(diff_tpca)

Plot thermal profile of protein pairs

Description

Plot thermal profile of protein pairs

Usage

plotPPiProfiles(tpcaObj, pair, splinesDf = 4)

Arguments

tpcaObj

a tpcaObj after having performed a differential analysis, see runDiffTPCA

pair

character vector of one or more protein names

splinesDf

numeric, degree of freedom of the spline fit to the melting curves

Value

ggplot displaying the thermal profile of a protein pair across conditions

Examples

library(Biobase)

set.seed(12)
m1 <- matrix(rnorm(50), ncol = 10)
m2 <- matrix(rnorm(50), ncol = 10)

rownames(m1) <- letters[1:5]
rownames(m2) <- letters[1:5]

colnames(m1) <- paste("fc", 1:10, sep = "_")
colnames(m2) <- paste("fc", 1:10, sep = "_")

pheno <- data.frame(
    temperature = seq(37, 67, length.out = 10))
rownames(pheno) <- paste("fc", 1:10, sep = "_")

eset1 <- ExpressionSet(
    assayData = m1,
    phenoData = AnnotatedDataFrame(pheno)
)

eset2 <- ExpressionSet(
    assayData = m2,
    phenoData = AnnotatedDataFrame(pheno)
)

tpcaObj <- new("tpcaResult",
               ObjList = list(eset1),

ContrastList = list(eset2),
                CtrlCondName = "control",
                ContrastCondName = "treatment")

plotPPiProfiles(tpcaObj, pair = c("b", "d"))

Plot PPI ROC curve

Description

Plot PPI ROC curve

Usage

plotPPiRoc(tpcaObj, computeAUC = FALSE)

Arguments

tpcaObj

tpcaResult object

computeAUC

logical parameter indicating whether area under the ROC should be computed and indicated in the lower right corner of the plot

Value

ggplot object of a receiver operating curve (ROC)

Examples

rocTab = data.frame(
    TPR = c(0, 0.1, 0.2, 0.4, 0.5, 0.7, 0.9, 1),
    FPR = c(0, 0.05, 0.1, 0.2, 0.5, 0.7, 0.9, 1)
)

tpcaTest <- new(
    "tpcaResult",
    PPiRocTable = rocTab)

plotPPiRoc(tpcaTest)

Plot TPCA analysis results

Description

Plot TPCA analysis results

Usage

plotTpcaVolcano(tpcaObj, alpha = 0.1)

Arguments

tpcaObj

a tpcaObj after having performed a differential analysis, see runDiffTPCA

alpha

significance level / FDR at which null hypothesis should be rejected

Value

ggplot displaying a volcano plot

Examples

library(dplyr)
library(Biobase)

m1 <- matrix(1:28, ncol = 4)
m2 <- matrix(2:25, ncol = 4)
m3 <- matrix(c(2:10, 1:19), ncol = 4)

rownames(m1) <- 1:7
rownames(m2) <- 3:8
rownames(m3) <- 2:8

mat_list <- list(
    m1, m2, m3
)

complex_anno <- tibble(
    protein = c("3", "4", "5", 
       "4", "5", "6", "7"),
    id = c(rep("1", 3), rep("2", 4)),
    count = c(rep(3, 3), rep(4, 4)))

tpca_result <- runTPCA(
  mat_list, complexAnno = complex_anno)

plotTpcaVolcano(tpca_result)

Extract PPiAnnotation

Description

Extract PPiAnnotation

Usage

## S4 method for signature 'tpcaResult'
PPiAnnotation(object)

Arguments

object

and object of class tpcaResult

Value

a data frame containing the results from a tpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
PPiAnnotation(tpcaObj)

Extract PPiRocTable

Description

Extract PPiRocTable

Usage

## S4 method for signature 'tpcaResult'
PPiRocTable(object)

Arguments

object

an object of class tpcaResult

Value

a data frame containing the results from a tpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
PPiRocTable(tpcaObj)

Extract PPiRocTableAnno

Description

Extract PPiRocTableAnno

Usage

## S4 method for signature 'tpcaResult'
PPiRocTableAnno(object)

Arguments

object

an object of class tpcaResult

Value

a data frame containing annotation information for PPiRocTable

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
PPiRocTableAnno(tpcaObj)

Run differential TPCA analysis

Description

Run differential TPCA analysis

Usage

runDiffTPCA(
  objList,
  contrastList,
  ctrlCondName = "control",
  contrastCondName = "treatment",
  ppiAnno = NULL,
  complexAnno = NULL,
  rownameCol = NULL,
  summaryMethodStr = "median",
  distMethodStr = "euclidean",
  n = 10000,
  p_adj_method = "BH"
)

Arguments

objList

input list of objects, e.g. ExpressionSets retrieved after TPP data import or matrices or data frames

contrastList

input list of objects for comparison at e.g. different treatment condtion, same file formats work as for objList

ctrlCondName

character string indicating the name of the control condition, default is "control"

contrastCondName

character string indicating the name of the contrast condition, default is "treatment"

ppiAnno

data frame annotation known protein-protein interactions (PPI) to test

complexAnno

data frame annotating known protein complexes of interest to test

rownameCol

in case the input objects are tibbles this parameter takes in the name (character) of the column specifying protein names or ids

summaryMethodStr

character string indicating a method to use to summarize measurements across replicates, default is "median", other options are c("mean", "rbind")

distMethodStr

method to use within dist function, default is 'euclidean'

n

number of random protein pair draws to obtain empirical p-value, default is 10000

p_adj_method

method to be used for multiple testing adjusment, default is "BH"

Value

an object of class tpcaResult with the following slots: 1) ObjList: containing the supplied list of objects

Examples

library(dplyr)
library(Biobase)

m1 <- matrix(1:28, ncol = 4)
m2 <- matrix(2:25, ncol = 4)
m3 <- matrix(c(2:10, 1:19), ncol = 4)

rownames(m1) <- 1:7
rownames(m2) <- 3:8
rownames(m3) <- 2:8

mat_list <- list(
    m1, m2, m3
)

c1 <- matrix(29:2, ncol = 4)
c2 <- matrix(26:3, ncol = 4)
c3 <- matrix(c(11:3, 20:2), ncol = 4)

rownames(c1) <- 1:7
rownames(c2) <- 3:8
rownames(c3) <- 2:8

contrast_list <- list(
    c1, c2, c3
)

ppi_anno <- tibble(
    x = c("3", "3"),
    y = c("5", "7"),
    pair = c("3:5", "3:7"))

ref_df <- tibble(
    pair = c("3:5", "3:7"),
    valueC2 = c(4, 8)
)

diff_tpca <- Rtpca:::runDiffTPCA(
  mat_list, contrast_list, ppiAnno = ppi_anno)

Run the TPCA analysis

Description

Run the TPCA analysis

Usage

runTPCA(
  objList,
  complexAnno = NULL,
  ppiAnno = NULL,
  rownameCol = NULL,
  summaryMethodStr = "median",
  distMethodStr = "euclidean",
  doRocAnalysis = TRUE,
  minCount = 3,
  nSamp = 10000,
  p_adj_method = "BH"
)

Arguments

objList

inout list of objects, e.g. ExpressionSets retrieved after TPP data import or matrices or data frames

complexAnno

data frame annotating known protein complexes of interest to test

ppiAnno

data frame annotation known protein-protein interactions (PPI) to test

rownameCol

in case the input objects are tibbles this parameter takes in the name (character) of the column specifying protein names or ids

summaryMethodStr

character string indicating a method to use to summarize measurements across replicates, default is "median", other options are c("mean", "rbind")

distMethodStr

method to use within dist function, default is 'euclidean'

doRocAnalysis

logical indicating whether a ROC analysis should be performed which can be used to assess the predictive power of the dataset for protein-protein interactions / protein complexes based on distanc between melting curves of protein interactions partners

minCount

integer indicating how many subunits of a complex should be qunatified to inlucde it into the analysis, default is 3

nSamp

integer indicating the number of random samples which should be performed to estimate empirical null distributions, default is 10000

p_adj_method

character string indicating a valid method to be used for multiple testing adjusment, default is "BH" which makes p.adjust use benjamini-hochberg, for additional options check ?p.adjust

Value

an object of class tpcaResult with the following slots: 1) ObjList: containing the supplied list of objects

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

colnames(m1) <- paste0("X", 1:4)
colnames(m2) <- paste0("X", 1:4)
colnames(m3) <- paste0("X", 1:4)

mat_list <- list(
    m1, m2, m3
)

ppi_anno <- tibble(
    x = "2",
    y = "3",
    combined_score = 700,
    pair = "2:3")

runTPCA(
    objList = mat_list,
    complexAnno = NULL,
    ppiAnno = ppi_anno
)

Set CommonFeatures

Description

Set CommonFeatures

Usage

## S4 method for signature 'tpcaResult'
SetCommonFeatures(object, commonFeatures)

Arguments

object

and object of class tpcaResult

commonFeatures

a vector of characters indicating the common features across replicates

Value

a vector of common features across replicates

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)

tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetCommonFeatures(tpcaObj, c("2", "3"))

Set ComplexAnnotation

Description

Set ComplexAnnotation

Usage

## S4 method for signature 'tpcaResult'
SetComplexAnnotation(object, df)

Arguments

object

an object of class tpcaResult

df

data frame containing complex annotation

Value

an object of class tpcaResult

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetComplexAnnotation(tpcaObj, data.frame(id = "complex1"))

Set ComplexBackgroundDistributionList

Description

Set ComplexBackgroundDistributionList

Usage

## S4 method for signature 'tpcaResult'
SetComplexBackgroundDistributionList(object, lt)

Arguments

object

an object of class tpcaResult

lt

a list of data frames containing distances of random complexes with different number of members

Value

an object of class tpcaResult

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetComplexBackgroundDistributionList(tpcaObj, 
  list('3' = data.frame(pair = "A:B")))

Set ComplexRocTable

Description

Set ComplexRocTable

Usage

## S4 method for signature 'tpcaResult'
SetComplexRocTable(object, df)

Arguments

object

and object of class tpcaResult

df

data.frame containg ComplexRocTable to set

Value

a data frame containing the results from a tpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetComplexRocTable(tpcaObj, data.frame(FPR = 1, TPR = 1))

Set ContrastCondName

Description

Set ContrastCondName

Usage

## S4 method for signature 'tpcaResult'
SetContrastCondName(object, name)

Arguments

object

an object of class tpcaResult

name

a character string describing the contrast condition

Value

an object of class tpcaResult

Examples

tpcaObj <- new("tpcaResult")
SetContrastCondName(tpcaObj, "DMSO")

Set ContrastDistMat

Description

Set ContrastDistMat

Usage

## S4 method for signature 'tpcaResult'
SetContrastDistMat(object, mat)

Arguments

object

and object of class tpcaResult

mat

matrix containg contrast distance matrix to set

Value

a data frame containing the results from a tpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetContrastDistMat(tpcaObj, matrix(c(0, 1, 0, 1), ncol = 2))

Set CtrlCondName

Description

Set CtrlCondName

Usage

## S4 method for signature 'tpcaResult'
SetCtrlCondName(object, name)

Arguments

object

an object of class tpcaResult

name

a character string describing the control condition

Value

an object of class tpcaResult

Examples

tpcaObj <- new("tpcaResult")
SetCtrlCondName(tpcaObj, "DMSO")

Set diffTpcaResultTable

Description

Set diffTpcaResultTable

Usage

## S4 method for signature 'tpcaResult'
SetDiffTpcaResultTable(object, df)

Arguments

object

an object of class tpcaResult

df

a data frame containing the results from a differential tpca analysis

Value

an object of class tpcaResult

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetDiffTpcaResultTable(tpcaObj, data.frame(pair = "A:B"))

Set DistMat

Description

Set DistMat

Usage

## S4 method for signature 'tpcaResult'
SetDistMat(object, mat)

Arguments

object

an object of class tpcaResult

mat

matrix containg distance matrix to set

Value

a data frame containing the results from a tpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetDistMat(tpcaObj, matrix(c(0, 1, 0, 1), ncol = 2))

Set distMethod

Description

Set distMethod

Usage

## S4 method for signature 'tpcaResult'
SetDistMethod(object, method)

Arguments

object

an object of class tpcaResult

method

character string of dist method

Value

an object of class tpcaResult

Examples

tpcaObj <- new("tpcaResult")
SetDistMethod(tpcaObj, "euclidean")

Set PPiRocTable

Description

Set PPiRocTable

Usage

## S4 method for signature 'tpcaResult'
SetPPiRocTable(object, df)

Arguments

object

an object of class tpcaResult

df

data.frame containg PPiRocTable to set

Value

an object of class tpcaResult

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetPPiRocTable(tpcaObj, data.frame(FPR = 1, TPR = 1))

Set PPiRocTableAnno

Description

Set PPiRocTableAnno

Usage

## S4 method for signature 'tpcaResult'
SetPPiRocTableAnno(object, df)

Arguments

object

an object of class tpcaResult

df

data.frame containg PPiRocTable annotation to set

Value

an object of class tpcaResult

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetPPiRocTableAnno(tpcaObj, data.frame(pair = "A:B"))

Set summaryMethod

Description

Set summaryMethod

Usage

## S4 method for signature 'tpcaResult'
SetSummaryMethod(object, method)

Arguments

object

an object of class tpcaResult

method

character string of summary method

Value

an object of class tpcaResult

Examples

tpcaObj <- new("tpcaResult")
SetSummaryMethod(tpcaObj, "median")

Set TpcaResultTable

Description

Set TpcaResultTable

Usage

## S4 method for signature 'tpcaResult'
SetTpcaResultTable(object, df)

Arguments

object

an object of class tpcaResult

df

a data frame containing the results from a tpca analysis

Value

an object of class tpcaResult

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
SetTpcaResultTable(tpcaObj, data.frame(pair = "A:B"))

Data frame of annotated human protein-protein interactions retrieved from stringDB with a combined interaction score equal or higher than 700

Description

data frame assigning proteins to interacting proteins

Usage

data("string_ppi_df")

Format

data frame with columns x, y (gene symbol of interactors), combined_score, pair (unique pair id)

References

Jensen et al. (2009), Nucleic Acids Research, 37, D412–D416

Examples

data("string_ppi_df")

Extract SummaryMethod

Description

Extract SummaryMethod

Usage

## S4 method for signature 'tpcaResult'
SummaryMethod(object)

Arguments

object

and object of class tpcaResult

Value

a character string of the summary method

Examples

tpcaObj <- new("tpcaResult")
SummaryMethod(tpcaObj)

S4 TPCA Result Class

Description

S4 TPCA Result Class

Value

an object of class tpcaResult with the following slots: 1) ObjList: containing the supplied list of objects (e.g. a list of Expression Sets summarizing a TPP experiment) 2) ContrastList: containing the supplied list of constrast objects (if supplied for performance of a differential Rtpca analysis) 3) CtrlCondName: character string indicating the control condition, e.g. "control" 4) ContrastCondName: character string indicating the contrast condition, e.g. "drug treatment" 5) DistMat: a matrix containing all pairwise protein-protein distances obtained from comparing their melting curves in the control condition 6) ContrastDistMat: a matrix containing all pairwise protein-protein distances obtained from comparing their melting curves in the contrast condition 7) CommonFeatures: a vector containing the features (proteins) found in common between control and contrast condition 8) ComplexAnnotation: a data frame supplied by the user annotating protein to protein complexes 9) ComplexBackgroundDistributionList: a list of distances drawn for random groups of proteins with different number of members 10) PPiAnnotation: a data frame supplied by the user annotating protein-protein interactions 11) PPiRocTable: data frame containing false positive rate and true positive rate based on ranking the TPCA analysis results by euclidean distance of melting curves of protein pairs, annotated PPIs are considered true positives 12) PPiRocTableAnno: annotation to PPiRocTable 13) ComplexRocTable: data frame containing false positive rate and true positive rate based on ranking the TPCA analysis results by euclidean distance of melting curves of proteins within annotated complexes, annotated complexes are considered true positives, proteins in randomly permuted complex annotations are considered false positives 14) summaryMethod: character string of summarization method used to summarize data across replicates 15) distMethod: character string of distance method used to compare melting curves of proteins 16) tpcaResultTable: data frame containing the results from a tpca analysis 17) diffTpcaResultTable: data frame containing the results from a differential tpca analysis

Slots

ObjList

list.

ContrastList

list.

CtrlCondName

character.

ContrastCondName

character.

DistMat

matrix.

ContrastDistMat

matrix

CommonFeatures

vector.

ComplexAnnotation

data.frame.

ComplexBackgroundDistributionList

list.

PPiAnnotation

data.frame.

PPiRocTable

data.frame

PPiRocTableAnno

data.frame

ComplexRocTable

data.frame

summaryMethod

character.

distMethod

character.

tpcaResultTable

data.frame.

diffTpcaResultTable

data.frame.

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)

Extract tpcaResultTable

Description

Extract tpcaResultTable

Usage

## S4 method for signature 'tpcaResult'
tpcaResultTable(object)

Arguments

object

an object of class tpcaResult

Value

a data frame containing the results from a tpca analysis

Examples

m1 <- matrix(1:12, ncol = 4)
m2 <- matrix(2:13, ncol = 4)
m3 <- matrix(c(2:10, 1:7), ncol = 4)

rownames(m1) <- 1:3
rownames(m2) <- 2:4
rownames(m3) <- 2:5

mat_list <- list(
    m1, m2, m3
)
tpcaObj <- new("tpcaResult", ObjList = mat_list)
tpcaResultTable(tpcaObj)