fitPoisthNorm | R Documentation |
Poisson threshold model based normalization-log2 transformation for single slide or for multiple slides
fitPoisthNorm(object, ...) ## S4 method for signature 'NanoStringGeoMxSet' fitPoisthNorm( object, split = FALSE, ROIs_high = NULL, features_high = NULL, features_all = NULL, sizefact_start = NULL, sizefact_BG = NULL, threshold_mean = NULL, preci2 = 10000, iterations = 2, prior_type = c("contrast", "equal"), sizefactrec = TRUE, size_scale = c("sum", "first"), sizescalebythreshold = FALSE, covrob = FALSE, preci1con = 1/25, cutoff = 15, confac = 1, calhes = FALSE ) ## S4 method for signature 'matrix' fitPoisthNorm( object, probenum = rep(1, NROW(object)), features_high, features_all, sizefact_start, sizefact_BG, threshold_mean, preci2 = 10000, iterations = 2, prior_type = c("contrast", "equal"), sizefactrec = TRUE, size_scale = c("sum", "first"), sizescalebythreshold = FALSE, covrob = FALSE, preci1con = 1/25, cutoff = 15, confac = 1, calhes = FALSE )
object |
count matrix with features in rows and samples in columns |
... |
additional argument list that might be used |
split |
indicator variable on whether it is for multiple slides (Yes, TRUE; No, FALSE) |
ROIs_high |
ROIs with high expressions defined based on featfact and featfact |
features_high |
subset of features which are well above the background |
features_all |
full feature vector to apply the normalization on |
sizefact_start |
initial value for size factors |
sizefact_BG |
size factor for background |
threshold_mean |
average threshold level |
preci2 |
precision for threshold, default=10000 |
iterations |
iteration number, default=2, the first iteration using the features_high to construct the prior for parameters then refit the model on all features. precision matrix for threshold: preci2 |
prior_type |
prior type for preci1, "equal" or "contrast", default="contrast" |
sizefactrec |
XXXX, default = TRUE |
size_scale |
method to scale the sizefact, sum(sizefact)=1 when size_scale="sum", sizefact[1]=1 when size_scale="first" |
sizescalebythreshold |
XXXX, default = FALSE |
covrob |
whether to use robust covariance in calculating the prior precision matrix 1, default = FALSE |
preci1con |
The user input constant term in specifying precision matrix 1, default=1/25 |
cutoff |
term in calculating precision matrix 1, default=15 |
confac |
The user input factor for contrast in precision matrix 1, default=1 |
calhes |
The user input whether to calculate hessian: calhes, default=FALSE |
probenum |
a vector of numbers of probes in each gene |
if split is FALSE, a valid GeoMx S4 object including the following items
para0_norm, matrix of estimated parameters for iter=1, features in columns and parameters(log2 expression, threshold) in rows, in featureData.
para_norm, matrix of estimated parameters for iter=2, features in columns and parameters(log2 expression, threshold) in rows, in featureData.
normmat0, matrix of log2 expression for iter=1, features in columns and log2 expression in rows, in assay slot.
normmat, matrix of log2 expression for iter=2, features in columns and log2 expression in rows, in assay lot.
sizefact_norm, estimated sizefact, in phenoData.
sizefact0_norm, estimated sizefact in iter=1, in phenoData.
preci1, precision matrix 1, in experimentData.
conv0, vector of convergence for iter=1, 0 converged, 1 not converged, in featureData
conv, vector of convergence for iter=2, 0 converged, 1 not converged, in featureData
features_high, same as the input features_high, in featureData
features_all, same as the input features_all, in featureData
if split is TRUE, a valid GeoMx S4 object with the following items appended.
threshold0, matrix of estimated threshold for iter=1, features in columns and threshold for different slides in rows, in featureData.
threshold, matrix of estimated threshold for iter=2, features in columns and threshold for different slides in rows, in featureData.
normmat0_sp, matrix of log2 expression for iter=1, features in columns and log2 expression in rows, in assay slot.
normmat_sp, matrix of log2 expression for iter=2, features in columns and log2 expression in rows, in assay slot.
sizefact_norm_sp, estimated sizefact, in phenoData
sizefact0_norm_sp, estimated sizefact in iter=1, in phenoData
preci1, precision matrix 1, in experimentData
conv0_sp_XX, vector of convergence for each unique slide value for iter=1, 0 converged, 1 not converged, in featureData for each unique slide.
conv_sp_XX, vector of convergence for each unique slide value for iter=2, 0 converged, 1 not converged, in featureData for each unique slide.
features_high_sp, same as the input features_high, in featureData.
features_all_sp, same as the input features_all, in featureData.
a list of following items
para0, matrix of estimated parameters for iter=1, features in columns and parameters(log2 expression, threshold) in rows.
para, matrix of estimated parameters for iter=2, features in columns and parameters(log2 expression, threshold) in rows.
normmat0, matrix of log2 expression for iter=1, features in columns and log2 expression in rows.
normmat, matrix of log2 expression for iter=2, features in columns and log2 expression in rows.
sizefact, estimated sizefact
sizefact0, estimated sizefact in iter=1
preci1, precision matrix 1
Im0, Information matrix of parameters in iter=1
Im, Information matrix of parameters in iter=2
conv0, vector of convergence for iter=1, 0 converged, 1 not converged
conv, vector of convergence for iter=2, 0 converged, 1 not converged
features_high, same as the input features_high
features_all, same as the input features_all
library(Biobase) library(dplyr) data(demoData) demoData <- fitPoisBG(demoData, size_scale = "sum") demoData <- aggreprobe(demoData, use = "cor") demoData <- BGScoreTest(demoData) thmean <- 1 * mean(fData(demoData)$featfact, na.rm = TRUE) demo_pos <- demoData[which(!fData(demoData)$CodeClass == "Negative"), ] demo_neg <- demoData[which(fData(demoData)$CodeClass == "Negative"), ] sc1_scores <- fData(demo_pos)[, "scores"] names(sc1_scores) <- fData(demo_pos)[, "TargetName"] features_high <- ((sc1_scores > quantile(sc1_scores, probs = 0.4)) & (sc1_scores < quantile(sc1_scores, probs = 0.95))) |> which() |> names() set.seed(123) features_high <- sample(features_high, 100) demoData <- fitNBth(demoData, features_high = features_high, sizefact_BG = demo_neg$sizefact, threshold_start = thmean, iterations = 5, start_para = c(200, 1), lower_sizefact = 0, lower_threshold = 100, tol = 1e-8) ROIs_high <- sampleNames(demoData)[which((quantile(fData(demoData)[["para"]][, 1], probs = 0.90, na.rm = TRUE) - notes(demoData)[["threshold"]]) * demoData$sizefact_fitNBth > 2)] features_all <- rownames(demo_pos) thmean <- mean(fData(demo_neg)[["featfact"]]) demoData <- fitPoisthNorm( object = demoData, split = FALSE, ROIs_high = ROIs_high, features_high = features_high, features_all = features_all, sizefact_start = demoData[, ROIs_high][["sizefact_fitNBth"]], sizefact_BG = demoData[, ROIs_high][["sizefact"]], threshold_mean = thmean, preci2 = 10000, prior_type = "contrast", covrob = FALSE, preci1con = 1 / 25 )
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.