Title: | Compositional omics model based visual integration |
---|---|
Description: | This explorative ordination method combines quasi-likelihood estimation, compositional regression models and latent variable models for integrative visualization of several omics datasets. Both unconstrained and constrained integration are available. The results are shown as interpretable, compositional multiplots. |
Authors: | Stijn Hawinkel [cre, aut] |
Maintainer: | Stijn Hawinkel <[email protected]> |
License: | GPL-2 |
Version: | 1.19.0 |
Built: | 2024-10-30 05:33:54 UTC |
Source: | https://github.com/bioc/combi |
Add a link on a compositional plot
addLink( DIplot, links, Views, samples, variable = NULL, Dims = c(1, 2), addLabel = FALSE, labPos = NULL, projColour = "grey", latentSize = 0.25 )
addLink( DIplot, links, Views, samples, variable = NULL, Dims = c(1, 2), addLabel = FALSE, labPos = NULL, projColour = "grey", latentSize = 0.25 )
DIplot |
A list with ggplot object where the links are to be added, and data frames with coordinates (obtained by setting plot(..., returnCoords = TRUE)) |
links |
A matrix containing either feature names (two column matrix) or approximate coordinates (four column matrix) |
Views |
Indices or names of the views for which the links should be added |
samples |
Sample names or approximate sample coordinates |
variable |
Name of variable in environmental gradient for which link should be plotted |
Dims |
vector of length 2 referring to the model dimensions |
addLabel |
A boolean, should arrow with label be plotted? |
labPos |
The position of the label, as a numeric vector of length 2 |
projColour |
The colour of the projection, as character string |
latentSize |
Size of the line from the origin to the latent variable dot |
A ggplot object with the links added
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) ## End(Not run) load(system.file("extdata", "zhangFits.RData", package = "combi")) Plot = plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", returnCoords = TRUE) addLink(Plot, links = cbind("OTU0565b3","OTUa14fb5"), Views = 1, samples = c(1,1))
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) ## End(Not run) load(system.file("extdata", "zhangFits.RData", package = "combi")) Plot = plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", returnCoords = TRUE) addLink(Plot, links = cbind("OTU0565b3","OTUa14fb5"), Views = 1, samples = c(1,1))
Array multiplication
arrayMult(centralMat, outerMat, ncols = ncol(outerMat))
arrayMult(centralMat, outerMat, ncols = ncol(outerMat))
centralMat |
an nxp matrix |
outerMat |
an nxd matrix |
ncols |
an integer, the number of columns of outerMat |
an nxpxd matrix, the stacked matrices of centralMat multiplied to every column of outermat
A function to build a centering matrix based on a dataframe
buildCentMat(object)
buildCentMat(object)
object |
an modelDI object or dataframe |
a centering matrix consisting of ones and zeroes, or a list with components
centMat |
a centering matrix consisting of ones and zeroes |
datFrame |
The dataframe with factors with one level removed |
Build the composition matrix for a certain dimension m dimensions
buildCompMat( colMat, paramEsts, latentVar, m, norm = TRUE, id = seq_len(m), subtractMax = TRUE )
buildCompMat( colMat, paramEsts, latentVar, m, norm = TRUE, id = seq_len(m), subtractMax = TRUE )
colMat |
The nxp independence model composition matrix |
paramEsts |
The matrix of feature parameter estimates |
latentVar |
The matrix of latent variables |
m |
the required dimension |
norm |
a boolean, should the composition matrix be normalized? |
id |
The vector of dimensions to consider |
subtractMax |
A boolean, should the maximum be substracted from every composition prior to exponentiation? Recommended for numerical stability |
A matrix with compositions in the rows
Build confounder design matrices with and without intercepts
buildConfMat(confounders)
buildConfMat(confounders)
confounders |
A dataframe of confounding variables #' For the preliminary trimming, we do not include an intercept, but we do include all the levels of the factors using contrasts=FALSE: we want to do the trimming in every subgroup, so no hidden reference levels For the filtering we just use a model with an intercept and treatment coding, here the interest is only in adjusting the offset |
a list with components
confModelMatTrim |
A confounder matrix without intercept, with all levels of factors present. This will be used to trim out taxa that have zero abundances in any subgroup defined by confounders |
confModelMat |
A confounder matrix with intercept, and with reference levels for factors absent. This will be used to fit the model to modify the independence model, and may include continuous variables |
A function to build the covariate matrix of the constraints
buildCovMat(datFrame)
buildCovMat(datFrame)
datFrame |
the dataframe with which the covariate matrix is to be built In this case we will 1) Include dummy's for every level of the categorical variable, and force them to sum to zero. This is needed for plotting and required for reference level indepenent normalization. 2) Exclude an intercept. The density function f() will provide this already. |
a list with components
covModelMat |
The model matrix |
datFrame |
The dataframe used to construct the model matrix |
Prepare an empty Jacobian matrix, with useful entries prefilled. In case of distribution "gaussian", it returns the lhs matrix of the linear system for finding the feature paramters
buildEmptyJac( n, m, lower, distribution = "quasi", normal = FALSE, nLambda1s = 1, centMat = NULL, weights = 1 )
buildEmptyJac( n, m, lower, distribution = "quasi", normal = FALSE, nLambda1s = 1, centMat = NULL, weights = 1 )
n |
the number of parameters |
m |
the dimension |
lower |
the current parameter estimates |
distribution |
A character string, the distributional assumption for the data |
normal |
a boolean, are normalization restrictions in place? |
nLambda1s |
The number of centering restrictions |
centMat |
The centering matrix |
weights |
Vector of feature weights |
an empty jacobian matrix, or the lhs of the system of estimating equations
Build an offset matrix from an marginal model object
buildMarginalOffset(indepModel, invLink)
buildMarginalOffset(indepModel, invLink)
indepModel |
The fitted marginal model, a list |
invLink |
The inverse link function |
an offset matrix of the size of the data
A function to build the mu matrix
buildMu(offSet, latentVar, paramEsts, distribution, paramMatrix = FALSE)
buildMu(offSet, latentVar, paramEsts, distribution, paramMatrix = FALSE)
offSet |
the offset matrix |
latentVar , paramEsts , distribution
|
Latent variables, parameter estimates and distribution type |
paramMatrix |
A boolean, are feature parameters provided as matrix |
The mean matrix
Build the marginal mu matrix
buildMuMargins(x, otherMargin, col)
buildMuMargins(x, otherMargin, col)
x |
The marginal parameters begin estimated |
otherMargin |
The parameters of the other margin |
col |
A logical, are the column parameters being estimated? |
a matrix of means
Build a marginal offset matrix given a model
buildOffsetModel(modelObj, View, distributions, compositional)
buildOffsetModel(modelObj, View, distributions, compositional)
modelObj |
a modelDI object |
View |
The view for which to build the offset |
distributions , compositional
|
belong to the view |
The offset matrix
Check for alias structures in a dataframe, and throw an error when one is found
checkAlias(datFrame, covariatesNames)
checkAlias(datFrame, covariatesNames)
datFrame |
the data frame to be checked for alias structure |
covariatesNames |
The names of the variables to be considered |
Throws an error when an alias structure is detected, returns invisible otherwise
Quickly check if the mean variance trend provides a good fit
checkMeanVarTrend(data, meanVarFit = "spline", returnTrend = FALSE, ...)
checkMeanVarTrend(data, meanVarFit = "spline", returnTrend = FALSE, ...)
data |
Data in any acceptable format (see details ?combi) |
meanVarFit |
The type of mean variance fit, either "cubic" or "spline" |
returnTrend |
A boolean, should the estimated trend be returned (TRUE) or only plotted (FALSE)? |
... |
passed on to the estMeanVarTrend() function |
A plot object
data(Zhang) par(mfrow = c(1,2)) lapply(list("microbiome" = zhangMicrobio, "metabolome" = zhangMetabo), checkMeanVarTrend) par(mfrow = c(1,1))
data(Zhang) par(mfrow = c(1,2)) lapply(list("microbiome" = zhangMicrobio, "metabolome" = zhangMetabo), checkMeanVarTrend) par(mfrow = c(1,1))
Check for monotonicity in compositional datasets fro given dimensions
checkMonotonicity(modelObj, Dim)
checkMonotonicity(modelObj, Dim)
modelObj |
The combi fit |
Dim |
The dimensions considered |
A boolean matrix indicating monotonicity for every feature
Perform model-based data integration
combi( data, M = 2L, covariates = NULL, distributions, compositional, maxIt = 300L, tol = 0.001, verbose = FALSE, prevCutOff = 0.95, minFraction = 0.1, logTransformGaussian = TRUE, confounders = NULL, compositionalConf = rep(FALSE, length(data)), nleq.control = list(maxit = 1000L, cndtol = 1e-16), record = TRUE, weights = NULL, fTol = 1e-05, meanVarFit = "spline", maxFeats = 2000, dispFreq = 10L, allowMissingness = FALSE, biasReduction = TRUE, maxItFeat = 20L, initPower = 1 )
combi( data, M = 2L, covariates = NULL, distributions, compositional, maxIt = 300L, tol = 0.001, verbose = FALSE, prevCutOff = 0.95, minFraction = 0.1, logTransformGaussian = TRUE, confounders = NULL, compositionalConf = rep(FALSE, length(data)), nleq.control = list(maxit = 1000L, cndtol = 1e-16), record = TRUE, weights = NULL, fTol = 1e-05, meanVarFit = "spline", maxFeats = 2000, dispFreq = 10L, allowMissingness = FALSE, biasReduction = TRUE, maxItFeat = 20L, initPower = 1 )
data |
A list of data objects with the same number of samples. See details. |
M |
the required dimension of the fit, a non-negative integer |
covariates |
a dataframe of n samples with sample-specific variables. |
distributions |
a character vector describing which distributional assumption should be used. See details. |
compositional |
A logical vector with the same length as "data", indicating if the datasets should be treated as compositional |
maxIt |
an integer, the maximum number of iterations |
tol |
A small scalar, the convergence tolerance |
verbose |
Logical. Should verbose output be printed to the console? |
prevCutOff |
a scalar, the prevalance cutoff for the trimming. |
minFraction |
a scalar, each taxon's total abundance should equal at least the number of samples n times minFraction, otherwise it is trimmed. |
logTransformGaussian |
A boolean, should the gaussian data be logtransformed, i.e. are they log-normal? |
confounders |
A dataframe or a list of dataframes with the same length as data. In the former case the same dataframe is used for conditioning, In the latter case each view has its own conditioning variables (or NULL). |
compositionalConf |
A logical vector with the same length as "data", indicating if the datasets should be treated as compositional when correcting for confounders. Numerical problems may occur when set to TRUE |
nleq.control |
A list of arguments to the nleqslv function |
record |
A boolean, should intermediate estimates be stored? Can be useful to check convergence |
weights |
A character string, either 'marginal' or 'uniform', indicating rrhow the feature parameters should be weighted in the normalization |
fTol |
The tolerance for solving the estimating equations |
meanVarFit |
The type of mean variance fit, see details |
maxFeats |
The maximal number of features for a Newton-Raphson procedure to be feasible |
dispFreq |
An integer, the period after which the variances should be reestimated |
allowMissingness |
A boolean, should NA values be allowed? |
biasReduction |
A boolean, should bias reduction be applied to allow for confounder correction in groups with all zeroes? Not guaranteed to work |
maxItFeat |
Integers, the maximum allowed number of iterations in the estimation of the feature parameters |
initPower |
The power to be applied to the residual matrix used to calculate the starting value. Must be positive; can be tweaked in case of numerical problems (i.e. infinite values returned by nleqslv) |
Data can be provided as raw matrices with features in the columns, or as phyloseq, SummarizedExperiment or ExpressionSet objects. Estimation of independence model and view wise parameters can be parametrized. See ?BiocParallel::bplapply and ?BiocParallel::register. meanVarFit = "spline" yields a cubic spline fit for the abundance-variance trend, "cubic" gives a third degree polynomial. Both converge to the diagonal line with slope 1 for small means. Distribution can be either "quasi" for quasi likelihood or "gaussian" for Gaussian data
An object of the "combi" class, containing all information on the data integration and fitting procedure
data(Zhang) #The method works on several datasets at once, and simply is not very fast. #Hence the "Not run" statement ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) ## End(Not run)
data(Zhang) #The method works on several datasets at once, and simply is not very fast. #Hence the "Not run" statement ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) ## End(Not run)
Plot the convegrence of the different parameter estimates in a line plot
convPlot( model, latent = is.null(View), nVars = Inf, Dim = 1L, View = NULL, size = 0.125 )
convPlot( model, latent = is.null(View), nVars = Inf, Dim = 1L, View = NULL, size = 0.125 )
model |
A fitted modelDI object |
latent |
A boolean, should latent variable trajectory be plotted |
nVars |
An integer, the number of variables to plot. By default all are plotted |
Dim |
An integer, the dimension to be plotted |
View |
An integer or character string, indicating the view to be plotted (if latent = FALSE) |
size |
The line size (see ?geom_path) |
A ggplot object containing the convergence plot
## Not run: data(Zhang) #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) ## End(Not run) load(system.file("extdata", "zhangFits.RData", package = "combi")) convPlot(microMetaboInt) convPlot(microMetaboInt, Dim = 2) convPlot(microMetaboInt, View = "microbiome")
## Not run: data(Zhang) #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) ## End(Not run) load(system.file("extdata", "zhangFits.RData", package = "combi")) convPlot(microMetaboInt) convPlot(microMetaboInt, Dim = 2) convPlot(microMetaboInt, View = "microbiome")
The score function to estimate the latent variables
deriv2LagrangianFeatures( x, data, distribution, offSet, latentVars, numVar, paramEstsLower, mm, Jac, meanVarTrend, weights, compositional, indepModel, ... )
deriv2LagrangianFeatures( x, data, distribution, offSet, latentVars, numVar, paramEstsLower, mm, Jac, meanVarTrend, weights, compositional, indepModel, ... )
x |
parameter estimates |
data |
A list of data matrices |
distribution , compositional , meanVarTrend , offSet , numVar
|
Characteristics of the view |
latentVars |
A vector of latent variables |
paramEstsLower |
lower dimension estimates |
mm |
the current dimension |
Jac |
a prefab jacobian |
weights |
The normalization weights |
indepModel |
the independence model |
... |
Additional arguments passed on to the score and jacobian functions |
A vector of length n, the evaluation of the score functions of the latent variables
The jacobian function to estimate the latent variables
deriv2LagrangianLatentVars( x, data, distributions, offsets, paramEsts, paramMats, numVars, latentVarsLower, n, m, Jac, numSets, meanVarTrends, links, varPosts, indepModels, compositional, ... )
deriv2LagrangianLatentVars( x, data, distributions, offsets, paramEsts, paramMats, numVars, latentVarsLower, n, m, Jac, numSets, meanVarTrends, links, varPosts, indepModels, compositional, ... )
x |
The current estimates of the latent variables |
distributions , links , compositional , data , meanVarTrends , offsets , numVars , numSets , paramMats , paramEsts , varPosts , indepModels
|
Characteristics of the views |
latentVarsLower |
The parameter estimates of the lower dimensions |
n , m
|
integers, number of samples and dimensions |
Jac |
an empty jacobian matrix |
... |
arguments to the jacobian function, currently ignored |
A vector of length n, the evaluation of the score functions of the latent variables
The score function to estimate the latent variables
deriv2LagrangianLatentVarsConstr( x, data, distributions, offsets, paramEsts, paramMats, numVars, latentVarsLower, nn, m, Jac, numSets, meanVarTrends, links, numCov, covMat, nLambda1s, varPosts, compositional, indepModels, ... )
deriv2LagrangianLatentVarsConstr( x, data, distributions, offsets, paramEsts, paramMats, numVars, latentVarsLower, nn, m, Jac, numSets, meanVarTrends, links, numCov, covMat, nLambda1s, varPosts, compositional, indepModels, ... )
x |
The current estimates of the latent variables |
distributions , data , links , compositional , meanVarTrends , offsets , numVars , paramMats , paramEsts
|
Characteristics of the view |
latentVarsLower |
The parameter estimates of the lower dimensions |
nn |
number of samples |
m , numSets , varPosts , indepModels
|
other arguments |
Jac |
an empty jacobian matrix |
numCov |
The number of covariates |
covMat |
the covariates matrix |
nLambda1s |
The number of centering restrictions |
... |
arguments to the jacobian function, currently ignored |
A vector of length nn, the evaluation of the score functions of the latent variables
The score function to estimate the feature parameters
derivLagrangianFeatures( x, data, distribution, offSet, latentVars, numVar, paramEstsLower, mm, indepModel, meanVarTrend, weights, compositional, ... )
derivLagrangianFeatures( x, data, distribution, offSet, latentVars, numVar, paramEstsLower, mm, indepModel, meanVarTrend, weights, compositional, ... )
x |
current parameter estimates |
data |
A list of data matrices |
distribution , compositional , meanVarTrend , offSet , numVar , indepModel , paramEstsLower
|
Characteristics of the view |
latentVars |
A vector of latent variables |
mm |
the current dimension |
weights |
The normalization weights |
... |
arguments to the jacobian function, currently ignored |
A vector with the evaluation of the score functions of the feature parameters
The score function to estimate the latent variables
derivLagrangianLatentVars( x, data, distributions, offsets, paramEsts, paramMats, numVars, n, m, numSets, meanVarTrends, links, varPosts, latentVarsLower, compositional, indepModels, ... )
derivLagrangianLatentVars( x, data, distributions, offsets, paramEsts, paramMats, numVars, n, m, numSets, meanVarTrends, links, varPosts, latentVarsLower, compositional, indepModels, ... )
x |
The current estimates of the latent variables |
n |
The number of samples |
m |
The dimensions |
numSets |
The number of views |
latentVarsLower |
The parameter estimates of the lower dimensions |
compositional , links , indepModels , meanVarTrends , numVars , distributions , data , offsets , varPosts , paramMats , paramEsts
|
Lists of inforamtion on all the views |
... |
arguments to the jacobian function, currently ignored |
A vector of length n, the evaluation of the score functions of the latent variables
The score function to estimate the latent variables
derivLagrangianLatentVarsConstr( x, data, distributions, offsets, paramEsts, numVars, latentVarsLower, n, m, numSets, meanVarTrends, links, covMat, numCov, centMat, nLambda1s, varPosts, compositional, indepModels, paramMats, ... )
derivLagrangianLatentVarsConstr( x, data, distributions, offsets, paramEsts, numVars, latentVarsLower, n, m, numSets, meanVarTrends, links, covMat, numCov, centMat, nLambda1s, varPosts, compositional, indepModels, paramMats, ... )
x |
The current estimates of the latent variables |
latentVarsLower |
The parameter estimates of the lower dimensions |
n |
The number of samples |
m |
The dimensions |
numSets |
The number of views |
covMat |
The covariance matrix |
numCov |
The number of covariates |
centMat |
A centering matrix |
nLambda1s |
The number of dummy variables |
compositional , links , indepModels , meanVarTrends , numVars , distributions , data , offsets , varPosts , paramMats , paramEsts
|
Lists of information on all the views |
... |
arguments to the jacobian function, currently ignored |
A vector of length n, the evaluation of the score functions of the latent variables
Estimate the feature parameters
estFeatureParameters( paramEsts, lambdasParams, seqSets, data, distributions, offsets, nCores, m, JacFeatures, meanVarTrends, latentVars, numVars, control, weights, compositional, indepModels, fTol, allowMissingness, maxItFeat, ... )
estFeatureParameters( paramEsts, lambdasParams, seqSets, data, distributions, offsets, nCores, m, JacFeatures, meanVarTrends, latentVars, numVars, control, weights, compositional, indepModels, fTol, allowMissingness, maxItFeat, ... )
paramEsts |
Current list of parameter estimates for the different views |
lambdasParams |
The lagrange multipliers |
seqSets |
A vector with view indices |
data |
A list of data matrices |
distributions |
A character vector describing the distributions |
offsets |
A list of offset matrices |
nCores |
The number of cores to use in multithreading |
m |
The dimension |
JacFeatures |
An empty Jacobian matrix |
meanVarTrends |
The mean-variance trends of the different views |
latentVars |
A vector of latent variables |
numVars |
The number of variables |
control |
A list of control arguments for the nleqslv function |
weights |
The normalization weights |
compositional |
A list of booleans indicating compositionality |
indepModels |
A list of independence model |
fTol |
A convergence tolerance |
allowMissingness |
A boolean indicating whether missing values are allowed |
maxItFeat |
An integer, the maximum number of iterations |
... |
Additional arguments passed on to the score and jacobian functions |
A vector with estimates of the feature parameters
Estimate the independence model belonging to one view
estIndepModel( data, distribution, compositional, maxIt, tol, link, invLink, meanVarFit, newtonRaphson, dispFreq, ... )
estIndepModel( data, distribution, compositional, maxIt, tol, link, invLink, meanVarFit, newtonRaphson, dispFreq, ... )
data |
a list of data matrices with the same number of samples n in the rows. Also phyloseq objects are acceptable |
distribution |
a character string describing which distributional assumption should be used. |
compositional |
A logical indicating if the dataset should be treated as compositional |
maxIt |
an integer, the maximum number of iterations |
tol |
A small scalar, the convergence tolerance |
link , invLink
|
link and inverse link function |
meanVarFit |
mean variance model |
newtonRaphson |
a boolean, should newton-raphson be used |
dispFreq |
An integer, frequency of dispersion estimation |
... |
passed on to the estOff() function |
A list with elements
rowOff |
The row offsets |
colOff |
The column offsets |
converged |
A logical flag, indicating whether the fit converged |
iter |
An integer, the number of iterations |
Estimate the latent variables
estLatentVars(latentVars, lambdasLatent, constrained, fTol, ...)
estLatentVars(latentVars, lambdasLatent, constrained, fTol, ...)
latentVars |
A vector of latent variables |
lambdasLatent |
A vector of Lagrange multipliers |
constrained |
A boolean, is the ordination constrained? |
fTol |
The convergence tolerance |
... |
additional arguments passed on to score and jacobian functions |
A vector of length n, the estimates of the latent variables
Estimate a column-wise mean-variance trend
estMeanVarTrend( data, meanMat, baseAbundances, libSizes, plot = FALSE, meanVarFit, degree = 2L, constraint = "none", ... )
estMeanVarTrend( data, meanMat, baseAbundances, libSizes, plot = FALSE, meanVarFit, degree = 2L, constraint = "none", ... )
data |
the data matrix with n rows |
meanMat |
the estimated mean matrix |
baseAbundances |
The baseline abundances |
libSizes |
Library sizes |
plot |
A boolean, should the trend be plotted? |
meanVarFit |
A character string describing the type of trend to be fitted: either "spline" or "cubic" |
degree |
The degree of the spline |
constraint |
Constraint to the spline |
... |
additional arguments passed on to the plot() function |
A list with components
meanVarTrend |
An smoothed trend function, that can map a mean on a variance |
meanVarTrendDeriv |
A derivative function of this |
Estimate the row/column parameters of the independence model
estOff( data, distribution, rowOff, colOff, meanVarTrend, col, newtonRaphson, libSizes, ... )
estOff( data, distribution, rowOff, colOff, meanVarTrend, col, newtonRaphson, libSizes, ... )
data |
a list of data matrices with the same number of samples n in the rows. Also phyloseq objects are acceptable |
distribution |
a character string describing which distributional assumption should be used. |
rowOff , colOff
|
current row and column offset estimates |
meanVarTrend |
The estimated mean-variance trend |
col |
A logical, should column offsets be estimated |
newtonRaphson |
A boolean, should Newton-Raphson be used to solve the estimating equations |
libSizes |
The library sizes, used to evaluate the mean-variance trend |
... |
passed onto nleqslv |
The estimated marginal parameters
Extract coordinates from fitted object
extractCoords(modelObj, Dim)
extractCoords(modelObj, Dim)
modelObj |
The fitted model |
Dim |
the required dimensions |
A list with components (matrices with two columns)
latentData |
The latent variables |
featureData |
The feature parameters |
varData |
The variables |
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) ## End(Not run) #Load the fits load(system.file("extdata", "zhangFits.RData", package = "combi")) extractCoords(microMetaboInt, Dim = c(1,2))
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) ## End(Not run) #Load the fits load(system.file("extdata", "zhangFits.RData", package = "combi")) extractCoords(microMetaboInt, Dim = c(1,2))
Helper function to extract data matrix from phyloseq, expressionset objects etc. Also filers out all zero rows
extractData(data, logTransformGaussian = TRUE)
extractData(data, logTransformGaussian = TRUE)
data |
The list of data objects, either matrix, phyloseq or ExpressionSet objects |
logTransformGaussian |
A boolean, should array data be logtransformed |
the raw data matrices, samples in the rows
data(Zhang) matrixList = extractData(list("microbiome" = zhangMicrobio, "metabolome" = zhangMetabo))
data(Zhang) matrixList = extractData(list("microbiome" = zhangMicrobio, "metabolome" = zhangMetabo))
A function to extract a data matrix from a number of objects
extractMat(Y, ...) ## S4 method for signature 'ExpressionSet' extractMat(Y, logTransformGaussian, ...) ## S4 method for signature 'SummarizedExperiment' extractMat(Y, ...) ## S4 method for signature 'matrix' extractMat(Y, ...)
extractMat(Y, ...) ## S4 method for signature 'ExpressionSet' extractMat(Y, logTransformGaussian, ...) ## S4 method for signature 'SummarizedExperiment' extractMat(Y, ...) ## S4 method for signature 'matrix' extractMat(Y, ...)
Y |
a phyloseq or eSet object, or another object, or a raw data matrix |
... |
additional arguments for the extractor function |
logTransformGaussian |
A boolean, should array data be logtransformed |
A data matrix with samples in the rows and features in the columns
Filter out the effect of known confounders
filterConfounders( confMat, data, distribution, link, invLink, compositional, control, meanVarTrend, offSet, numVar, marginModel, biasReduction, allowMissingness )
filterConfounders( confMat, data, distribution, link, invLink, compositional, control, meanVarTrend, offSet, numVar, marginModel, biasReduction, allowMissingness )
confMat |
A confounder design matrix |
data |
data matrix |
distribution , link , invLink , compositional , meanVarTrend , offSet , numVar , marginModel
|
Characteristics of the view |
control |
A list of control elements to the nleqslv function |
biasReduction |
A boolean, should bias reduction be applied |
allowMissingness |
A boolean, are missing values allowed? |
Parameter estimates accounting for the effects of the confounders
Extract the influence on the estimation of the latent variable
getInflLatentVar(score, InvJac, i)
getInflLatentVar(score, InvJac, i)
score |
The score matrix |
InvJac |
The inverse Jacobian |
i |
the sample index |
The influence of all observations on the i-th latent variable
Gram schimdt orhtogonalize a with respect to b, and normalize
gramSchmidtOrth(a, b, weights = 1, norm = TRUE)
gramSchmidtOrth(a, b, weights = 1, norm = TRUE)
a |
the vector to be orthogonalized |
b |
the vector to be orthogonalized to |
weights |
weights vector |
norm |
a boolean, should the result be normalized? |
The orthogonalized vector
Functions to indent the plot to include the entire labels
indentPlot(plt, xInd = 0, yInd = 0)
indentPlot(plt, xInd = 0, yInd = 0)
plt |
a ggplot object |
xInd |
a scalar or a vector of length 2, specifying the indentation left and right of the plot to allow for the labels to be printed entirely |
yInd |
a a scalar or a vector of length 2, specifying the indentation top and bottom of the plot to allow for the labels to be printed entirely |
a ggplot object, squared
A ggplot line plot showing the influences
inflPlot( modelObj, plotType = ifelse(length(modelObj$data) <= 2, "pointplot", "boxplot"), pointFun = "sum", lineSize = 0.07, Dim = 1, samples = seq_len(nrow(if (is.null(modelObj$covariates)) modelObj$latentVars else modelObj$alphas)), ... )
inflPlot( modelObj, plotType = ifelse(length(modelObj$data) <= 2, "pointplot", "boxplot"), pointFun = "sum", lineSize = 0.07, Dim = 1, samples = seq_len(nrow(if (is.null(modelObj$covariates)) modelObj$latentVars else modelObj$alphas)), ... )
modelObj |
The fitted data integration |
plotType |
The type of plot requested, see details |
pointFun |
The function to calculate the summary measure to be plotted |
lineSize |
The line size |
Dim |
The dimension required |
samples |
index vector of which samples to be plotted |
... |
additional arguments passed on to the influence() function |
The options for plotType are: "pointPlot": Dot plot of total influence per view and sample, "boxplot": plot boxplot of influence of all observations per view and sample, "boxplotSingle": boxplot of log absolute total influence per view, "lineplot": line plot of total influence per view and sample. In the pointplot, dots crosses represent parameter estimates
A ggplot object
data(Zhang) #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) load(system.file("extdata", "zhangFits.RData", package = "combi")) inflPlot(microMetaboInt) #Constrained inflPlot(microMetaboIntConstr)
data(Zhang) #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) load(system.file("extdata", "zhangFits.RData", package = "combi")) inflPlot(microMetaboInt) #Constrained inflPlot(microMetaboIntConstr)
Evaluate the influence function
## S3 method for class 'combi' influence(modelObj, samples = is.null(View), Dim = 1, View = NULL)
## S3 method for class 'combi' influence(modelObj, samples = is.null(View), Dim = 1, View = NULL)
modelObj |
The model object |
samples |
A boolean, should we look at sample variables? Throws an error otherwise |
Dim , View
|
Integers, the dimension and views required |
Especially the influence of the different views on the latent variable or gradient estimation may be of interest. The influence values are not all calculated. Rather, the score values and inverse jacobian are returned so they can easily be calculated
A list with components
score |
The evaluation of the score function |
InvJac |
The inverted jacobian matrix |
Jacobian when estimating confounder variables
jacConfounders( confMat, data, distribution, x, meanVarTrend, offSet, CompMat, libSizes, allowMissingness )
jacConfounders( confMat, data, distribution, x, meanVarTrend, offSet, CompMat, libSizes, allowMissingness )
data , confMat , meanVarTrend
|
Characteristics of the views |
distribution , offSet
|
distribution and offset of the view |
x |
the parameter estimates |
libSizes , CompMat
|
Library sizes and relative abunance |
allowMissingness |
a boolean, should missing values be allowed |
the jacobian matrix
Jacobian for conditioning under compositionality
jacConfoundersComp( x, confMat, data, meanVarTrend, marginModel, allowMissingness, biasReduction, subtractMax = TRUE )
jacConfoundersComp( x, confMat, data, meanVarTrend, marginModel, allowMissingness, biasReduction, subtractMax = TRUE )
x |
the parameter estimates |
confMat , data , meanVarTrend
|
arguments belonging to views |
marginModel , biasReduction , subtractMax
|
The marginal mode, and booleans indicating bias reduction and maximum subtraction |
allowMissingness |
a boolean, should missing values be allowed |
the jacobian matrix
Evaluate the jacobian for estimating the feature parameters for one view
jacFeatures( latentVars, data, distribution, paramEsts, meanVarTrend, offSet, compositional, indepModel, m, paramEstsLower, allowMissingness, ... )
jacFeatures( latentVars, data, distribution, paramEsts, meanVarTrend, offSet, compositional, indepModel, m, paramEstsLower, allowMissingness, ... )
latentVars |
A vector of latent variables |
data |
A list of data matrices |
distribution , compositional , meanVarTrend , offSet , paramEsts , paramEstsLower , indepModel
|
Characteristics of each view |
m |
dimension |
allowMissingness |
a boolean, are missing values allowed? |
... |
Additional arguments passed on to the score and jacobian functions |
The jacobian matrix
Evaluate the jacobian for estimating the latent variable for one view
jacLatentVars( latentVar, data, distribution, paramEsts, paramMats, offSet, meanVarTrend, n, varPosts, mm, indepModel, latentVarsLower, compositional, allowMissingness, ... )
jacLatentVars( latentVar, data, distribution, paramEsts, paramMats, offSet, meanVarTrend, n, varPosts, mm, indepModel, latentVarsLower, compositional, allowMissingness, ... )
latentVar |
the latent variable estimates |
distribution , data , varPosts , compositional , meanVarTrend , offSet , paramEsts , paramMats , indepModel
|
Characteristics of each view |
n |
the number of samples |
mm |
the current dimension |
latentVarsLower |
the lower dimensional latent variables |
allowMissingness |
a boolean, should missing values be allowed |
... |
additional arguments passed on to score and jacobian functions |
The diagonal of the jacobian matrix
Evaluate the jacobian for estimating the latent variable for one view for constrained ordination
jacLatentVarsConstr( latentVar, data, distribution, paramEsts, offSet, meanVarTrend, numCov, covMat, varPosts, compositional, mm, indepModel, latentVarsLower, ... )
jacLatentVarsConstr( latentVar, data, distribution, paramEsts, offSet, meanVarTrend, numCov, covMat, varPosts, compositional, mm, indepModel, latentVarsLower, ... )
latentVar |
current latent variable estimates |
distribution , compositional , meanVarTrend , offSet , paramEsts , indepModel , varPosts , data
|
Characteristics of each view |
numCov |
the number of covariates |
covMat |
the covariates matrix |
mm |
the dimension |
latentVarsLower |
latent variable estimates of lower dimensions |
... |
additional arguments passed on to score and jacobian functions |
The jacobian matrix
Make multiplots of the data integration object
## S3 method for class 'combi' plot( x, ..., Dim = c(1, 2), samDf = NULL, samShape = NULL, samCol = NULL, featurePlot = "threshold", featNum = 15L, samColValues = NULL, manExpFactorTaxa = 0.975, featSize = switch(featurePlot, threshold = 2.5, points = samSize * 0.7, density = 0.35), crossSize = 4, manExpFactorVar = 0.975, varNum = nrow(x$alphas), varSize = 2.5, samSize = 1.75, featCols = c("darkblue", "darkgreen", "grey10", "turquoise4", "blue", "green", "grey", "cornflowerblue", "lightgreen", "grey75"), strokeSize = 0.05, warnMonotonicity = FALSE, returnCoords = FALSE, squarePlot = TRUE, featAlpha = 0.5, featShape = 8, xInd = 0, yInd = 0, checkOverlap = FALSE, shapeValues = (21:(21 + length(unique(samDf[[samShape]])))) )
## S3 method for class 'combi' plot( x, ..., Dim = c(1, 2), samDf = NULL, samShape = NULL, samCol = NULL, featurePlot = "threshold", featNum = 15L, samColValues = NULL, manExpFactorTaxa = 0.975, featSize = switch(featurePlot, threshold = 2.5, points = samSize * 0.7, density = 0.35), crossSize = 4, manExpFactorVar = 0.975, varNum = nrow(x$alphas), varSize = 2.5, samSize = 1.75, featCols = c("darkblue", "darkgreen", "grey10", "turquoise4", "blue", "green", "grey", "cornflowerblue", "lightgreen", "grey75"), strokeSize = 0.05, warnMonotonicity = FALSE, returnCoords = FALSE, squarePlot = TRUE, featAlpha = 0.5, featShape = 8, xInd = 0, yInd = 0, checkOverlap = FALSE, shapeValues = (21:(21 + length(unique(samDf[[samShape]])))) )
x |
the fit |
... |
additional arguments, currently ignored |
Dim |
the dimensions to be plotted |
samDf |
a dataframe of sample variables |
samShape |
A variable name from samDf used to shape the samples |
samCol |
A variable name from samDf used to colour the samples |
featurePlot |
A character string, either "threshold", "points" or "density". See details |
featNum , varNum
|
The number of features and variables to plot |
samColValues |
Colours for the samples |
manExpFactorTaxa , manExpFactorVar
|
Expansion factors for taxa and variables, normally calculated natively |
featSize , crossSize , varSize , samSize , strokeSize
|
Size parameters for the features (text, dots or density contour lines), central cross, variable labels, sample dots, sample strokes and feature contour lines |
featCols |
Colours for the features |
warnMonotonicity |
A boolean, should a warning be thrown when the feature proportions of compositional views do not all vary monotonically with all latent variables? |
returnCoords |
A boolean, should coordinates be returned, e.g. for use in thrird party software |
squarePlot |
A boolean, should the axes be square? Strongly recommended |
featAlpha |
Controls the transparacny of the features |
featShape |
Shape of feature dots when featurePlot = "points" |
xInd , yInd
|
x and y indentations |
checkOverlap |
A boolean, should overlapping labels be omitted? |
shapeValues |
the shapes, as numeric values |
It is usually impossible to plot all features with their labels. Therefore, he default option of the 'featurePlot' parameter is "threshold", whereby only the 'featNum" features furthest away from the origin are shown. Alternatively, the "points" or "density" options are available to plot all features as a point or density cloud, but without labels.
A ggplot object containing the plot
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) ## End(Not run) #Load the fits load(system.file("extdata", "zhangFits.RData", package = "combi")) plot(microMetaboInt) plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX") #Plot all features as points or density plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", featurePlot = "points") plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", featurePlot = "density") #Constrained plot(microMetaboIntConstr, samDf = zhangMetavars, samCol = "ABX")
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) ## End(Not run) #Load the fits load(system.file("extdata", "zhangFits.RData", package = "combi")) plot(microMetaboInt) plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX") #Plot all features as points or density plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", featurePlot = "points") plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", featurePlot = "density") #Constrained plot(microMetaboIntConstr, samDf = zhangMetavars, samCol = "ABX")
Horner's method to evaluate a polynomial, copied from the polynom package. the most efficient way
polyHorner(coefs, x)
polyHorner(coefs, x)
coefs |
the polynomial coefficients |
x |
the input values for the polynomial function |
the evaluated polynomial
A custom spline prediction function, extending linearly with a slope such that prediction never drops below first bisectant
predictSpline( fit, newdata, linX, coefsQuad, deriv = 0L, meanVarFit, minFit, new.knots, degree )
predictSpline( fit, newdata, linX, coefsQuad, deriv = 0L, meanVarFit, minFit, new.knots, degree )
fit |
The existing spline fit |
newdata |
points in which the spline needs to be evaluated |
linX |
The x at which the fit becomes linear and intersects the diagonal line |
coefsQuad |
parameters of a quadratic fit |
deriv |
An integer. Which derivative is required? |
meanVarFit |
A character string, indicating which type of mean variance fit is being used |
minFit |
The lower bound of the cubic fit |
new.knots |
The knots at which the spline is to be evaluated |
degree |
The degree of the polynomial fit |
The evaluation of the spline, i.e. the predicted variance
prepare the jacobian matrix
prepareJacMat(mu, data, meanVarTrend, CompMat, libSizes)
prepareJacMat(mu, data, meanVarTrend, CompMat, libSizes)
mu |
the mean matrix |
data |
the count matrix |
meanVarTrend |
The mean variance trend |
CompMat |
The compisition matrix |
libSizes |
The library sizes |
the matrix which can be summed over
prepare the jacobian for the latent variabels compostional
prepareJacMatComp(mu, paramEsts, CompMat0, meanVarTrend, data, libSizes)
prepareJacMatComp(mu, paramEsts, CompMat0, meanVarTrend, data, libSizes)
mu |
the mean matrix |
paramEsts |
Current parameter estimates |
CompMat0 |
The compisition matrix |
meanVarTrend |
The mean variance trend |
data |
the count matrix |
libSizes |
The library sizesv |
The empty jacobian matrix with entries maximally filled out
Prepare a helper matrix for score function evaluation under quasi-likelihood
prepareScoreMat(data, mu, meanVarTrend, CompMat, libSizes)
prepareScoreMat(data, mu, meanVarTrend, CompMat, libSizes)
data |
the count matrix |
mu |
the mean matrix |
meanVarTrend |
The mean variance trend |
CompMat |
The compisition matrix |
libSizes |
The library sizes |
The helper matrix
Print an overview of a fitted combi x
## S3 method for class 'combi' print(x, ...)
## S3 method for class 'combi' print(x, ...)
x |
a fitted combi x |
... |
Further arguments, currently ignored |
An overview of the number of dimensions, views and parameters, type of ordination and importance parameters
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) ## End(Not run) #Load the fits load(system.file("extdata", "zhangFits.RData", package = "combi")) print(microMetaboInt) print(microMetaboIntConstr) #Or simply microMetaboInt
data(Zhang) ## Not run: #Unconstrained microMetaboInt = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, verbose = TRUE) #Constrained microMetaboIntConstr = combi( list("microbiome" = zhangMicrobio, "metabolomics" = zhangMetabo), distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE, covariates = zhangMetavars, verbose = TRUE) ## End(Not run) #Load the fits load(system.file("extdata", "zhangFits.RData", package = "combi")) print(microMetaboInt) print(microMetaboIntConstr) #Or simply microMetaboInt
The jacobian for column offset estimation
quasiJacIndep(x, data, otherMargin, meanVarTrend, col, libSizes, ...)
quasiJacIndep(x, data, otherMargin, meanVarTrend, col, libSizes, ...)
x |
the initial guess for the current margin |
data |
the data matrix |
otherMargin |
The other margin |
meanVarTrend |
the function describing the mean-variance trend |
col |
A logical, is the column being estimated? |
libSizes |
The library sizes |
... |
passed on to prepareJacMat |
the jacobian matrix
Quasi score equations for column offset parameters of sequence count data
quasiScoreIndep(x, data, otherMargin, meanVarTrend, col, libSizes, ...)
quasiScoreIndep(x, data, otherMargin, meanVarTrend, col, libSizes, ...)
x |
the initial guess for the current margin |
data |
the data matrix |
otherMargin |
The other margin |
meanVarTrend |
the function describing the mean-variance trend |
col |
A logical, is the column being estimated? |
libSizes |
The library sizes |
... |
passed on to prepareJacMat |
the evaluated estimating equation
A function to efficiently row multiply a matrix and a vector
rowMultiply(matrix, vector)
rowMultiply(matrix, vector)
matrix |
a numeric matrix of dimension a-by-b |
vector |
a numeric vector of length b t(t(matrix)*vector) but then faster |
Memory intensive but that does not matter with given matrix sizes
a matrix, row multplied by the vector
A helper function to rescale coordinates
scaleCoords(featCoords, latentData, manExpFactorTaxa, featNum = NULL)
scaleCoords(featCoords, latentData, manExpFactorTaxa, featNum = NULL)
featCoords |
the feature coordinates to be rescaled |
latentData |
latent variables |
manExpFactorTaxa |
an expansion factor |
featNum |
the number of features to retain |
The rescaled feature coordinates
Score functions for confounder variables
scoreConfounders( x, data, distribution, offSet, confMat, meanVarTrend, allowMissingness, libSizes, CompMat )
scoreConfounders( x, data, distribution, offSet, confMat, meanVarTrend, allowMissingness, libSizes, CompMat )
x |
the parameter estimates |
data , distribution , offSet , confMat , meanVarTrend
|
Characteristics of the views |
allowMissingness |
a boolean, should missing values be allowed |
libSizes , CompMat
|
Library sizes and relative abunance |
The evaluation of the estimating equations
Score equations for conditioning under compositionality
scoreConfoundersComp( x, confMat, data, meanVarTrend, marginModel, biasReduction, allowMissingness, subtractMax = TRUE )
scoreConfoundersComp( x, confMat, data, meanVarTrend, marginModel, biasReduction, allowMissingness, subtractMax = TRUE )
x |
Confounder parameter estimates |
confMat |
confounder matrix |
data |
data |
meanVarTrend |
mean variance trend |
marginModel |
marginal models |
biasReduction |
A boolean, should a bias reduced estimation be applied? |
allowMissingness |
A boolean, are missing values allowed |
subtractMax |
A boolean, shuold the maximum be subtracted before softmax transformation? Recommended for numerical stability |
The evaluation of the estimating equations
Evaluate the score functions for the estimation of the feature parameters for a single dataset
scoreFeatureParams( x, data, distribution, offSet, latentVar, meanVarTrend, mm, indepModel, compositional, paramEstsLower, allowMissingness, ... )
scoreFeatureParams( x, data, distribution, offSet, latentVar, meanVarTrend, mm, indepModel, compositional, paramEstsLower, allowMissingness, ... )
x |
the parameter estimates |
data , distribution , offSet , meanVarTrend , indepModel , compositional , paramEstsLower
|
Characteristics of the views |
latentVar |
the latent variables |
mm |
the dimension |
allowMissingness |
a boolean, should missing values be allowed |
... |
Additional arguments passed on to the score and jacobian functions |
A vector with evaluated score function
Evaluate the score functions for the estimation of the latent variables for a single dataset
scoreLatentVars( data, distribution, paramEsts, paramMats, offSet, latentVar, meanVarTrend, constrained = FALSE, covMat = NULL, varPosts, compositional, indepModel, mm, latentVarsLower, allowMissingness, ... )
scoreLatentVars( data, distribution, paramEsts, paramMats, offSet, latentVar, meanVarTrend, constrained = FALSE, covMat = NULL, varPosts, compositional, indepModel, mm, latentVarsLower, allowMissingness, ... )
data , distribution , offSet , meanVarTrend , indepModel , varPosts , paramEsts , paramMats , compositional
|
Characteristics of the views |
latentVar |
the latent variable estimates |
constrained |
a boolean, is this a constrained analysis |
covMat |
a matrix of constraining covariates |
mm |
the current dimension |
latentVarsLower |
the lower dimensional latent variables |
allowMissingness |
a boolean, should missing values be allowed |
... |
additional arguments passed on to score and jacobian functions |
A vector of length n, with evaluated score function
A small auxiliary function for the indices of the lagrange multipliers
seqM(y, normal = TRUE, nLambda1s = 1)
seqM(y, normal = TRUE, nLambda1s = 1)
y |
an integer, the current dimension |
normal |
a logical, is there a normalization restriction? |
nLambda1s |
the number of centering restrictions |
a vector containing the ranks of the current lagrangian multipliers
Trim based on confounders to avoid taxa with only zero counts
trimOnConfounders(confounders, data, prevCutOff, minFraction, n)
trimOnConfounders(confounders, data, prevCutOff, minFraction, n)
confounders |
a nxt confounder matrix |
data |
the data matrix |
prevCutOff |
a scalar between 0 and 1, the prevalence cut off |
minFraction |
a scalar between 0 and 1, each taxon's total abundance should equal at least the number of samples n times minFraction, otherwise it is trimmed |
n |
the number of samples Should be called prior to fitting the independence model |
A trimmed data matrix nxp'
Metabolome of mice that underwent Pulsed Antibiotic Treatment (PAT) and controls
data(Zhang)
data(Zhang)
SummarizedExperiment with metabolome data
The metabolome data as a SummarizedExperiment object
Baseline covariates of PAT mice and healthy controls
data(Zhang)
data(Zhang)
A dataframe with baseline sample variables
The metadata on the mice
Microbiome of mice that underwent Pulsed Antibiotic Treatment (PAT) and controls. The data were extracted from the source https://www.ibdmdb.org/, and then only the samples matching between microbiome and metabolome were retained.
data(Zhang)
data(Zhang)
A phyloseq object containing microbiome data
The microbiome dataset pruned for matches with the metabolome object