Title: | Hierarchical Epitope pROtein biNding |
---|---|
Description: | HERON is a software package for analyzing peptide binding array data. In addition to identifying significant binding probes, HERON also provides functions for finding epitopes (string of consecutive peptides within a protein). HERON also calculates significance on the probe, epitope, and protein level by employing meta p-value methods. HERON is designed for obtaining calls on the sample level and calculates fractions of hits for different conditions. |
Authors: | Sean McIlwain [aut, cre] , Irene Ong [aut] |
Maintainer: | Sean McIlwain <[email protected]> |
License: | GPL (>= 3) |
Version: | 1.5.0 |
Built: | 2024-10-30 07:27:18 UTC |
Source: | https://github.com/bioc/HERON |
Add Sequence Annotations for Epitopes
addSequenceAnnotations(eds)
addSequenceAnnotations(eds)
eds |
HERONEpitopeDataSet with probe_meta in metadata() |
HERONEpitopeDataSet with the rowData() set with sequence annotations
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_pr_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_pr_res) segments_res <- findEpitopeSegments(calls_res, "unique") epval_res <- calcEpitopePValues(calls_res, segments_res) epval_res <- addSequenceAnnotations(epval_res)
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_pr_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_pr_res) segments_res <- findEpitopeSegments(calls_res, "unique") epval_res <- calcEpitopePValues(calls_res, segments_res) epval_res <- addSequenceAnnotations(epval_res)
Calculate p-values using the "exprs" assay
calcCombPValues( obj, colData_in = NULL, d_sd_shift = NA, d_abs_shift = NA, d_paired = FALSE, g_sd_shift = 0, use = "tz", p_adjust_method = "BH" )
calcCombPValues( obj, colData_in = NULL, d_sd_shift = NA, d_abs_shift = NA, d_paired = FALSE, g_sd_shift = 0, use = "tz", p_adjust_method = "BH" )
obj |
HERONSequenceDataSet or HERONProbeDataSet |
colData_in |
optional column DataFrame (default: NULL => colData(obj))) |
d_sd_shift |
standard deviation shift for differential test |
d_abs_shift |
absolute shift for differential test |
d_paired |
run paired analysis |
g_sd_shift |
standard deviation shift for global test |
use |
use global-test ("z"), differential-test using t.test ("t"), differential-test using wilcox ("w"), or both global and differential ("tz") |
p_adjust_method |
method for adjusting p-values |
HERONSequenceDataSet/HERONProbeDataSet with the pvalue assay added
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan)
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan)
Calculate epitope-level p-values
calcEpitopePValues( probe_pds, epitope_ids, metap_method = "wmax1", p_adjust_method = "BH" )
calcEpitopePValues( probe_pds, epitope_ids, metap_method = "wmax1", p_adjust_method = "BH" )
probe_pds |
HERONProbeDataSet with the "pvalue" assay |
epitope_ids |
vector of epitope ids |
metap_method |
meta p-value method to use (see below) |
p_adjust_method |
what p.adjust method to use. |
The meta p-value methods supported by calcEpitopePValues
are:
min_bonf*,
min*,
max*,
fischer/sumlog,
hmp/harmonicmeanp,
wilkinsons_min1/tippets,
wilkinsons_min2/wmin2,
wilkinsons_min3,
wilkinsons_min4,
wilkinsons_min5,
wilkinsons_max1/wmax1,
wilkinsons_max2/wmax2,
and cct.
When choosing a p-value method, keep in mind that the epitope p-value should be one that requires most of the probe p-values to be small (e.g. *wmax1*) Other p-value methods such as the*cct* and the *hmp* have been shown to be more accurate with p-value that have dependencies.
HERONEpitopeDataSet with "pvalue" and "padj" assays
[stats::p.adjust()] for p_adjust_parameter.
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_pr_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_pr_res) segments_res <- findEpitopeSegments(calls_res, "unique") epval_res <- calcEpitopePValues(calls_res, segments_res)
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_pr_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_pr_res) segments_res <- findEpitopeSegments(calls_res, "unique") epval_res <- calcEpitopePValues(calls_res, segments_res)
Calculate Probe p-values using a differential paired t-test
calcProbePValuesTPaired( probe_mat, colData_in, sd_shift = NA, abs_shift = NA, debug = FALSE )
calcProbePValuesTPaired( probe_mat, colData_in, sd_shift = NA, abs_shift = NA, debug = FALSE )
probe_mat |
numeric matrix or data.frame of values |
colData_in |
design data.frame |
sd_shift |
standard deviation shift to use when calculating p-values. Either sd_shift or abs_shift should be set |
abs_shift |
absolute shift to use when calculating p-values. |
debug |
print debugging information |
matrix of p-values on the post columns defined in the colData matrix. Attributes of the matrix are:
pars - data.frame parameters used in the paired t-test for each row (e.g. df, sd)
mapping - data.frame of mapping used for pre-post column calculation diff_mat - data.frame containing the post-pre differences for each sample (column) and probe (row)
data(heffron2021_wuhan) colData_wu <- colData(heffron2021_wuhan) pre_idx = which(colData_wu$visit == "pre") ## Make some samples paired colData_post = colData_wu[colData_wu$visit == "post",] new_ids = rownames(colData_post)[seq_len(5)] colData_wu$ptid[pre_idx[seq_len(5)]] = new_ids exprs <- assay(heffron2021_wuhan, "exprs") pval_res <- calcProbePValuesTPaired(exprs, colData_wu)
data(heffron2021_wuhan) colData_wu <- colData(heffron2021_wuhan) pre_idx = which(colData_wu$visit == "pre") ## Make some samples paired colData_post = colData_wu[colData_wu$visit == "post",] new_ids = rownames(colData_post)[seq_len(5)] colData_wu$ptid[pre_idx[seq_len(5)]] = new_ids exprs <- assay(heffron2021_wuhan, "exprs") pval_res <- calcProbePValuesTPaired(exprs, colData_wu)
Calculate Probe p-values using a differential unpaired t-test
calcProbePValuesTUnpaired(probe_mat, colData_in, sd_shift = NA, abs_shift = NA)
calcProbePValuesTUnpaired(probe_mat, colData_in, sd_shift = NA, abs_shift = NA)
probe_mat |
numeric matrix or data.frame of values |
colData_in |
design data.frame |
sd_shift |
standard deviation shift to use when calculating p-values Either sd_shift or abs_shift should be set |
abs_shift |
absolute shift to use when calculating p-values |
matrix of p-values on the post columns defined in the colData matrix
data(heffron2021_wuhan) colData_wu <- colData(heffron2021_wuhan) pval_res <- calcProbePValuesTUnpaired(assay(heffron2021_wuhan), colData_wu)
data(heffron2021_wuhan) colData_wu <- colData(heffron2021_wuhan) pval_res <- calcProbePValuesTUnpaired(assay(heffron2021_wuhan), colData_wu)
Calculate Probe p-values using a two-sample wilcoxon test
calcProbePValuesWUnpaired(probe_mat, colData_in, exact = NULL, abs_shift = 0)
calcProbePValuesWUnpaired(probe_mat, colData_in, exact = NULL, abs_shift = 0)
probe_mat |
numeric matrix or data.frame of values |
colData_in |
design data.frame |
exact |
a logical indicating whether an exact p-value should be computed (see wilcox.test for details) |
abs_shift |
absolute shift to use when calculating p-values |
matrix of p-values on the post columns defined in the colData matrix
data(heffron2021_wuhan) colData_wu <- colData(heffron2021_wuhan) pval_res <- calcProbePValuesWUnpaired(assay(heffron2021_wuhan), colData_wu)
data(heffron2021_wuhan) colData_wu <- colData(heffron2021_wuhan) pval_res <- calcProbePValuesWUnpaired(assay(heffron2021_wuhan), colData_wu)
Calculate protein-level p-values
calcProteinPValues(epitope_ds, metap_method = "wmin1", p_adjust_method = "BH")
calcProteinPValues(epitope_ds, metap_method = "wmin1", p_adjust_method = "BH")
epitope_ds |
HERONEpitopeDataSet with the "pvalue" assay |
metap_method |
meta p-value method to use |
p_adjust_method |
p.adjust method to use |
see calcEpitopePValues for a list of meta p-value methods supported by HERON. the protein should be one that requires at least one of the epitope p-values to be small (e.g. wmax1).
HERONProteinDataSet with the "pvalue" and "padj" assays
[stats::p.adjust()] for p_adjust_parameter.
[calcEpitopePValues()] for meta p-value methods
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_pr_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_pr_res) segments_res <- findEpitopeSegments(calls_res, "unique") epval_res <- calcEpitopePValues(calls_res, segments_res) ppval_res <- calcProteinPValues(epval_res)
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_pr_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_pr_res) segments_res <- findEpitopeSegments(calls_res, "unique") epval_res <- calcEpitopePValues(calls_res, segments_res) ppval_res <- calcProteinPValues(epval_res)
Concatenate sequences together based upon their start positions. Assumes the probe sequences have an overlap.
catSequences(positions, sequences)
catSequences(positions, sequences)
positions |
start positions of probes in protein |
sequences |
probe sequences of probes |
concatenated sequence (character)
positions <- c(1,2) sequences <- c("MSGSASFEGGVFSPYL", "SGSASFEGGVFSPYLT") catSequences(positions, sequences)
positions <- c(1,2) sequences <- c("MSGSASFEGGVFSPYL", "SGSASFEGGVFSPYLT") catSequences(positions, sequences)
Convert HERONSequenceDataSet to HERONProbeDataSet
convertSequenceDSToProbeDS(seq_ds, probe_meta)
convertSequenceDSToProbeDS(seq_ds, probe_meta)
seq_ds |
a HERONSequenceDataSet object |
probe_meta |
optional data.frame with the PROBE_SEQUENCE, PROBE_ID columns the probe meta data frame can be provided within the metadata()$probe_meta or as a argument to the function. The argument supersedes the metadata list. |
HERONProbeDataSet
data(heffron2021_wuhan) probe_ds <- convertSequenceDSToProbeDS(heffron2021_wuhan) probe_meta <- metadata(heffron2021_wuhan)$probe_meta probe_ds <- convertSequenceDSToProbeDS(heffron2021_wuhan, probe_meta)
data(heffron2021_wuhan) probe_ds <- convertSequenceDSToProbeDS(heffron2021_wuhan) probe_meta <- metadata(heffron2021_wuhan)$probe_meta probe_ds <- convertSequenceDSToProbeDS(heffron2021_wuhan, probe_meta)
This function will find blocks of consecutive probes within the passed probe parameter
findBlocksProbeT( probes, protein_tiling, proteins = getProteinLabel(probes), starts = getProteinStart(probes) )
findBlocksProbeT( probes, protein_tiling, proteins = getProteinLabel(probes), starts = getProteinStart(probes) )
probes |
vector of probe identifiers of the format c(Prot1;1, ... Prot1;10) |
protein_tiling |
tiling of the associated proteins |
proteins |
associated proteins to probes (cache speed up) |
starts |
associated starts from probes (cache speed up) |
data.frame with the Protein, Start, Stop, and Number.Of.Probes columns
findBlocksProbeT(c("A;1", "A;2", "A;3", "B;2", "B;3", "C;10", "A;5", "A;6"))
findBlocksProbeT(c("A;1", "A;2", "A;3", "B;2", "B;3", "C;10", "A;5", "A;6"))
Find consecutive probes
findBlocksT(prot_df, protein_tiling)
findBlocksT(prot_df, protein_tiling)
prot_df |
data.frame with the Protein and Starting position of the probe |
protein_tiling |
tiling for information for each protein |
data.frame with the Protein, Start, Stop, and Number.Of.Probes columns
probes = c("A;1","A;2","A;3", "A;5","A;6", "A;8") prot_df = data.frame( Protein = getProteinLabel(probes), Pos = getProteinStart(probes) ) findBlocksT(prot_df)
probes = c("A;1","A;2","A;3", "A;5","A;6", "A;8") prot_df = data.frame( Protein = getProteinLabel(probes), Pos = getProteinStart(probes) ) findBlocksT(prot_df)
Find Epitopes from probe stats and calls.
findEpitopeSegments( PDS_obj, segment_method = "unique", segment_score_type = "binary", segment_dist_method = "hamming", segment_cutoff = "silhouette" )
findEpitopeSegments( PDS_obj, segment_method = "unique", segment_score_type = "binary", segment_dist_method = "hamming", segment_cutoff = "silhouette" )
PDS_obj |
HERONProbeDataSet with pvalues and calls in the assay |
segment_method |
which epitope finding method to use (binary or zscore, applies for hclust or skater) |
segment_score_type |
which type of scoring to use for probes |
segment_dist_method |
what kind of distance score method to use |
segment_cutoff |
for clustering methods, what cutoff to use (either numeric value or 'silhouette') |
a vector of epitope identifiers or segments found
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) segments_res <- findEpitopeSegments(pr_calls_res)
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) segments_res <- findEpitopeSegments(pr_calls_res)
Create EpitopeID from protein, first and last probes
getEpitopeID(protein, start, stop)
getEpitopeID(protein, start, stop)
protein |
vector of proteins |
start |
vector of first probe protein start positions |
stop |
vector of last probe protein start positions |
vector of epitope ids
getEpitopeID("A", 1, 2)
getEpitopeID("A", 1, 2)
Get probe ids from a vector of epitope ids
getEpitopeIDsToProbeIDs(epitope_ids, tiling = 1)
getEpitopeIDsToProbeIDs(epitope_ids, tiling = 1)
epitope_ids |
vector of epitope identifiers |
tiling |
tling of probes across proteins |
data.frame of epitope_to_probe mappings
getEpitopeIDsToProbeIDs(c("A_1_5","C_8_12"))
getEpitopeIDsToProbeIDs(c("A_1_5","C_8_12"))
Get the vector of probes from an epitope id
getEpitopeProbeIDs(epitope_id, tiling = 1)
getEpitopeProbeIDs(epitope_id, tiling = 1)
epitope_id |
EpitopeID to obtain probes from |
tiling |
Tiling of the probes across the protein (default 1) |
vector of probe_ids that are contained within the epitope
getEpitopeProbeIDs("A_1_5")
getEpitopeProbeIDs("A_1_5")
Format of EpitopeID is A_B_C, where A is the protein label B is the protein start position of the first probe in the epitope and C is the protein start position of the last probe in the epitope.
getEpitopeProtein(epitope_ids)
getEpitopeProtein(epitope_ids)
epitope_ids |
vector of epitope identifier character strings |
vector of protein labels
getEpitopeProtein("Prot1_1_5")
getEpitopeProtein("Prot1_1_5")
Obtain first probe's protein start position from Epitope ID
getEpitopeStart(epitope_ids)
getEpitopeStart(epitope_ids)
epitope_ids |
vector of epitope ids |
vector of integers indicating first probe start positions in the epitope(s)
getEpitopeStart("Prot1_1_5")
getEpitopeStart("Prot1_1_5")
Obtain last probe's protein start position from EpitopeID
getEpitopeStop(epitope_ids)
getEpitopeStop(epitope_ids)
epitope_ids |
vector of epitope ids |
vector of integers indicating the last probe protein start position
getEpitopeStop("Prot1_1_5")
getEpitopeStop("Prot1_1_5")
Calculates the number of samples (K), the frequency of samples (F), and the percentage of samples (P) called. If the colData DataFrame contains a condition column with at least two conditions, then a K, F, and P is calculated for each condition and the results are reported as separate columns.
getKofN(obj)
getKofN(obj)
obj |
HERON Dataset with a "calls" assay |
DataFrame with K (#calls), F (fraction calls), P (
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) getKofN(pr_calls_res)
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) getKofN(pr_calls_res)
Get Protein Label from Probe
getProteinLabel(probes)
getProteinLabel(probes)
probes |
vector of probes (i.e. c("A;1", "A;2")) |
vector of strings indicating the protein associated with the respective probes
getProteinLabel("A;1") getProteinLabel("B;2") getProteinLabel(c("A;1","B;2"))
getProteinLabel("A;1") getProteinLabel("B;2") getProteinLabel(c("A;1","B;2"))
Get the amino-acid starting position of the probe within the protein.
getProteinStart(probes)
getProteinStart(probes)
probes |
vector of probes (i.e. c("A;1", "A;2")) |
starting locations of the probes with their associated proteins
getProteinStart("A;1") getProteinStart("B;2") getProteinStart(c("A;1","B;2"))
getProteinStart("A;1") getProteinStart("B;2") getProteinStart(c("A;1","B;2"))
Given a set of probes, estimate the tiling of the probes across the protein. Usually, you will want to calculate this on all the probes available in the dataset.
getProteinTiling(probes, return.vector = TRUE)
getProteinTiling(probes, return.vector = TRUE)
probes |
vector of probes (i.e. A;1, A;2) |
return.vector |
Return result as vector or return as data.frame |
For each protein, the estimating tiling (spacing) of the probes across the amino acid sequence.
getProteinTiling(c("A;1","A;2","A;3", "B;2","B;3", "C;1", "C;3"))
getProteinTiling(c("A;1","A;2","A;3", "B;2","B;3", "C;1", "C;3"))
A subset of data from the paper https://www.ncbi.nlm.nih.gov/pmc/articles/PMC8245122/ publication.
data(heffron2021_wuhan)
data(heffron2021_wuhan)
## 'heffron2021_wuhan' A HERONSequenceDataSet with and "exprs" assay DataFrame with 1945 rows and 60 columns. Each column is a pre-processed binding signal from a serum sample peptide array set for the SARS-CoV-2. The matrix is a subset of the full matrix and contains sequences from the membrane, envelope, surface (spike), and nucleocapsid proteins.
The metadata()$probe_meta is a data frame with 1945 rows and 6 columns. The columns are POSITION - starting position of probe within protein, PROBE_SEQUENCE - amino acid sequence of probe, SEQ_ID - protein identifier SEQ_NAME - name of protein, PROBE_ID - combination of protein identifier and starting position, e.g. prot1;5.
The colData() is a DataFrame with 60 rows and 2 columns. The columns are SampleName - name of the sample, visit - either pre or post, ptid - subject id, and condition - all COVID
HERONSequenceDataSet
<https://github.com/Ong-Research/UW_Adult_Covid-19>
HERONEpitopeDataSet
is a subclass of SummarizedExperiment
used to hold assay information on the epitope-level
HERONEpitopeDataSet(pvalue, ...)
HERONEpitopeDataSet(pvalue, ...)
pvalue |
calculate epitope p-value matrix |
... |
arguments provided to |
HERONEpitopeDataSet object
pval <- matrix(runif(100),ncol=4) HERONEpitopeDataSet(pvalue = pval)
pval <- matrix(runif(100),ncol=4) HERONEpitopeDataSet(pvalue = pval)
HERONProbeDataSet
is a subclass of RangedSummarizedExperiment
used to hold assay information on the probe level
HERONProbeDataSet(...)
HERONProbeDataSet(...)
... |
arguments provided to |
HERONProbeDataSet object
pds <- HERONProbeDataSet()
pds <- HERONProbeDataSet()
HERONProteinDataSet
is a subclass of SummarizedExperiment
used to hold assay information on the protein-level
HERONProteinDataSet(pvalue, ...)
HERONProteinDataSet(pvalue, ...)
pvalue |
calculated protein p-value matrix |
... |
arguments provided to |
HERONProteinDataSet object
pval <- matrix(runif(100), ncol=4) HERONProteinDataSet(pvalue = pval)
pval <- matrix(runif(100), ncol=4) HERONProteinDataSet(pvalue = pval)
HERONSequenceDataSet
is a subclass of SummarizedExperiment
,
used to store the expression values, intermediate calculations, and
results of a differential binding code on the seeuqnce-level.
HERONSequenceDataSet(exprs, ...)
HERONSequenceDataSet(exprs, ...)
exprs |
binding values with rows as sequences and columns as samples |
... |
arguments provided to metadata can contain a probe DataFrame, that maps sequences (column PROBE_SEQUENCE) to probe identifiers ( column PROBE_ID) |
HERONSequenceDataSet object
exprs <- matrix(seq_len(100),ncol=4) colnames(exprs) <- c("C1", "C2", "C3", "C4") sds <- HERONSequenceDataSet(exprs = exprs)
exprs <- matrix(seq_len(100),ncol=4) colnames(exprs) <- c("C1", "C2", "C3", "C4") sds <- HERONSequenceDataSet(exprs = exprs)
log2 transform the "exprs" assay
log2Transform(se)
log2Transform(se)
se |
SummarizedExperiment with "exprs" assay |
SummarizedExperiment with "exprs" assay log2 transformed
data(heffron2021_wuhan) assay(heffron2021_wuhan, "exprs") <- 2^assay(heffron2021_wuhan, "exprs") res <- log2Transform(heffron2021_wuhan)
data(heffron2021_wuhan) assay(heffron2021_wuhan, "exprs") <- 2^assay(heffron2021_wuhan, "exprs") res <- log2Transform(heffron2021_wuhan)
Make Epitope Calls
makeEpitopeCalls(epi_ds, padj_cutoff = 0.05, one_hit_filter = TRUE)
makeEpitopeCalls(epi_ds, padj_cutoff = 0.05, one_hit_filter = TRUE)
epi_ds |
HERONEpitopeDataSet with pvalue assay |
padj_cutoff |
p-value cutoff to use |
one_hit_filter |
filter one hit epitopes? |
HERONEpitopeDataSet with calls assay added
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) epi_segments_uniq_res <- findEpitopeSegments( PDS_obj = pr_calls_res, segment_method = "unique" ) epi_padj_uniq <- calcEpitopePValues( probe_pds = pr_calls_res, epitope_ids = epi_segments_uniq_res, metap_method = "wilkinsons_max1" ) makeEpitopeCalls(epi_padj_uniq)
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) epi_segments_uniq_res <- findEpitopeSegments( PDS_obj = pr_calls_res, segment_method = "unique" ) epi_padj_uniq <- calcEpitopePValues( probe_pds = pr_calls_res, epitope_ids = epi_segments_uniq_res, metap_method = "wilkinsons_max1" ) makeEpitopeCalls(epi_padj_uniq)
makeProbeCalls
returns call information on a HERONProbeDataSet
using the "padj" assay
makeProbeCalls(pds, padj_cutoff = 0.05, one_hit_filter = TRUE)
makeProbeCalls(pds, padj_cutoff = 0.05, one_hit_filter = TRUE)
pds |
HERONProbeDataSet with the "padj" assay |
padj_cutoff |
cutoff to use |
one_hit_filter |
filter out one-hit probes? |
HERONProbeDataSet with the "calls" assay added
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_probe_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_probe_res)
data(heffron2021_wuhan) pval_seq_res <- calcCombPValues(heffron2021_wuhan) pval_probe_res <- convertSequenceDSToProbeDS(pval_seq_res) calls_res <- makeProbeCalls(pval_probe_res)
Make Protein-level Calls
makeProteinCalls(prot_ds, padj_cutoff = 0.05, one_hit_filter = FALSE)
makeProteinCalls(prot_ds, padj_cutoff = 0.05, one_hit_filter = FALSE)
prot_ds |
HERONProteinDataSet with the "padj" assay |
padj_cutoff |
cutoff to use |
one_hit_filter |
use the one-hit filter? |
HERONProteinDataSet with the "calls" assay added
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) epi_segments_uniq_res <- findEpitopeSegments( PDS_obj = pr_calls_res, segment_method = "unique" ) epi_padj_uniq <- calcEpitopePValues( probe_pds = pr_calls_res, epitope_ids = epi_segments_uniq_res, metap_method = "wilkinsons_max1" ) prot_padj_uniq <- calcProteinPValues( epitope_ds = epi_padj_uniq, metap_method = "tippetts" ) prot_calls <- makeProteinCalls(prot_padj_uniq)
data(heffron2021_wuhan) seq_pval_res <- calcCombPValues(heffron2021_wuhan) pr_pval_res <- convertSequenceDSToProbeDS(seq_pval_res) pr_calls_res <- makeProbeCalls(pr_pval_res) epi_segments_uniq_res <- findEpitopeSegments( PDS_obj = pr_calls_res, segment_method = "unique" ) epi_padj_uniq <- calcEpitopePValues( probe_pds = pr_calls_res, epitope_ids = epi_segments_uniq_res, metap_method = "wilkinsons_max1" ) prot_padj_uniq <- calcProteinPValues( epitope_ds = epi_padj_uniq, metap_method = "tippetts" ) prot_calls <- makeProteinCalls(prot_padj_uniq)
Cap vector at minimum/maximum values
min_max(val, min.value, max.value)
min_max(val, min.value, max.value)
val |
vector of values to cap |
min.value |
minimum value |
max.value |
maximum value |
vector of capped values
min_max(10, 1, 5)
min_max(10, 1, 5)
Find One-hit epitopes
oneHitEpitopes(sample_epitopes)
oneHitEpitopes(sample_epitopes)
sample_epitopes |
logical epitope matrix from makeCalls |
vector of one-hit, one-probe epitopes
hit_mat = data.frame( row.names = c("A_1_1","A_2_2","A_3_3","A_4_4"), sample1 = c(TRUE, FALSE, FALSE, TRUE), sample2 = c(TRUE, TRUE, FALSE, FALSE), sample3 = c(TRUE, TRUE, FALSE, FALSE) ) oneHitEpitopes(hit_mat)
hit_mat = data.frame( row.names = c("A_1_1","A_2_2","A_3_3","A_4_4"), sample1 = c(TRUE, FALSE, FALSE, TRUE), sample2 = c(TRUE, TRUE, FALSE, FALSE), sample3 = c(TRUE, TRUE, FALSE, FALSE) ) oneHitEpitopes(hit_mat)
Find one hit probes
oneHitProbes(sample_probes)
oneHitProbes(sample_probes)
sample_probes |
logical probe matrix from makeCalls |
vector of probes that are one-hits
hit_mat <- data.frame( row.names = c("A;1","A;2","A;3","A;4"), sample1 = c(TRUE, FALSE, FALSE, TRUE), sample2 = c(TRUE, TRUE, FALSE, FALSE), sample3 = c(TRUE, TRUE, FALSE, FALSE) ) oneHitProbes(hit_mat)
hit_mat <- data.frame( row.names = c("A;1","A;2","A;3","A;4"), sample1 = c(TRUE, FALSE, FALSE, TRUE), sample2 = c(TRUE, TRUE, FALSE, FALSE), sample3 = c(TRUE, TRUE, FALSE, FALSE) ) oneHitProbes(hit_mat)
Indicate which epitopes are just one probe.
oneProbeEpitopes(epitope_ids)
oneProbeEpitopes(epitope_ids)
epitope_ids |
vector of epitope ids |
vector of logical indicating epitopes that are one probe
oneProbeEpitopes(c("A_1_1", "B_1_1","C_1_2"))
oneProbeEpitopes(c("A_1_1", "B_1_1","C_1_2"))
Find probe hits with a consecutive probe or another sample
probeHitSupported(hit_mat)
probeHitSupported(hit_mat)
hit_mat |
matrix of logical values that indicate a hit with a TRUE value |
matrix of logical values indicate that the TRUE hit is supported by a consecutive probe hit in the sample sample or the within another sample
Convert p-value matrix to a z-score matrix
pvalue_to_zscore(mat.in, one.sided = TRUE, log.p = FALSE, inf.zscore = 16)
pvalue_to_zscore(mat.in, one.sided = TRUE, log.p = FALSE, inf.zscore = 16)
mat.in |
matrix of p-values |
one.sided |
p-values one-sided |
log.p |
are p-values log transformed? |
inf.zscore |
infinite z-scores are capped to this value |
matrix of z-scores
mat <- matrix(runif(100), nrow=10) rownames(mat) <- paste0("A;",seq_len(nrow(mat))) pvalue_to_zscore(mat)
mat <- matrix(runif(100), nrow=10) rownames(mat) <- paste0("A;",seq_len(nrow(mat))) pvalue_to_zscore(mat)
Normalize the exprs assay using quantile normalization
quantileNormalize(se)
quantileNormalize(se)
se |
SummarizedExperiment with exprs assay |
SummarizedExperiment with exprs assay normalized
data(heffron2021_wuhan) seq_ds_qn <- quantileNormalize(heffron2021_wuhan)
data(heffron2021_wuhan) seq_ds_qn <- quantileNormalize(heffron2021_wuhan)
Smooth probes across protein tiling
smoothProbeDS(probe_ds, w = 2, eps = 1e-06)
smoothProbeDS(probe_ds, w = 2, eps = 1e-06)
probe_ds |
HERONProbeDataSet to smooth |
w |
smoothing width, probes +/- w/2 before and after are used |
eps |
error tolerance |
HERONProbeDataSet with smoothed data in exprs object
data(heffron2021_wuhan) probe_ds <- convertSequenceDSToProbeDS(heffron2021_wuhan) smoothed_ds <- smoothProbeDS(probe_ds)
data(heffron2021_wuhan) probe_ds <- convertSequenceDSToProbeDS(heffron2021_wuhan) smoothed_ds <- smoothProbeDS(probe_ds)