scran
C++ library for basic single-cell RNA-seq analyses
Loading...
Searching...
No Matches
Classes | Public Member Functions | List of all members
scran::ModelGeneVariances Class Reference

Compute and model the per-gene variances in log-expression data. More...

#include <ModelGeneVariances.hpp>

Classes

struct  BlockResults
 Results of variance modelling with blocks. More...
 
struct  Defaults
 Default parameters for variance modelling. More...
 
struct  Results
 Results of variance modelling without blocks. More...
 

Public Member Functions

ModelGeneVariancesset_span (double s=FitVarianceTrend::Defaults::span)
 
ModelGeneVariancesset_minimum_mean (double m=FitVarianceTrend::Defaults::minimum_mean)
 
ModelGeneVariancesset_use_fixed_width (bool u=FitVarianceTrend::Defaults::use_fixed_width)
 
ModelGeneVariancesset_fixed_width (double f=FitVarianceTrend::Defaults::fixed_width)
 
ModelGeneVariancesset_minimum_window_count (int c=FitVarianceTrend::Defaults::minimum_window_count)
 
ModelGeneVariancesset_block_weight_policy (WeightPolicy w=Defaults::block_weight_policy)
 
ModelGeneVariancesset_variable_block_weight_parameters (VariableBlockWeightParameters v=Defaults::variable_block_weight_parameters)
 
ModelGeneVariancesset_compute_average (bool a=Defaults::compute_average)
 
ModelGeneVariancesset_num_threads (int n=Defaults::num_threads)
 
template<typename Value_ , typename Index_ , typename Stat_ >
void run (const tatami::Matrix< Value_, Index_ > *mat, Stat_ *means, Stat_ *variances, Stat_ *fitted, Stat_ *residuals) const
 
template<typename Value_ , typename Index_ , typename Block_ , typename Stat_ >
void run_blocked (const tatami::Matrix< Value_, Index_ > *mat, const Block_ *block, std::vector< Stat_ * > means, std::vector< Stat_ * > variances, std::vector< Stat_ * > fitted, std::vector< Stat_ * > residuals, Stat_ *ave_means, Stat_ *ave_variances, Stat_ *ave_fitted, Stat_ *ave_residuals) const
 
template<typename Value_ , typename Index_ , typename Block_ , typename Stat_ >
void run_blocked (const tatami::Matrix< Value_, Index_ > *mat, const Block_ *block, std::vector< Stat_ * > means, std::vector< Stat_ * > variances, std::vector< Stat_ * > fitted, std::vector< Stat_ * > residuals) const
 
template<typename Value_ , typename Index_ >
Results run (const tatami::Matrix< Value_, Index_ > *mat) const
 
template<typename Value_ , typename Index_ , typename Block_ >
BlockResults run_blocked (const tatami::Matrix< Value_, Index_ > *mat, const Block_ *block) const
 

Detailed Description

Compute and model the per-gene variances in log-expression data.

This scans through a log-transformed normalized expression matrix (e.g., from LogNormCounts) and computes per-feature means and variances. It then fits a trend to the variances with respect to the means using FitVarianceTrend. We assume that most genes at any given abundance are not highly variable, such that the fitted value of the trend is interpreted as the "uninteresting" variance - this is mostly attributed to technical variation like sequencing noise, but can also represent constitutive biological noise like transcriptional bursting. Under this assumption, the residual can be treated as a quantification of biologically interesting variation, and can be used to identify relevant features for downstream analyses.

Member Function Documentation

◆ set_span()

ModelGeneVariances & scran::ModelGeneVariances::set_span ( double  s = FitVarianceTrend::Defaults::span)
inline
Parameters
sSee FitVarianceTrend::set_span() for more details.
Returns
A reference to this ModelGeneVariances object.

◆ set_minimum_mean()

ModelGeneVariances & scran::ModelGeneVariances::set_minimum_mean ( double  m = FitVarianceTrend::Defaults::minimum_mean)
inline
Parameters
mSee FitVarianceTrend::set_minimum_mean() for more details.
Returns
A reference to this ModelGeneVariances object.

◆ set_use_fixed_width()

ModelGeneVariances & scran::ModelGeneVariances::set_use_fixed_width ( bool  u = FitVarianceTrend::Defaults::use_fixed_width)
inline
Parameters
uSee FitVarianceTrend::set_use_fixed_width() for more details.
Returns
A reference to this ModelGeneVariances object.

◆ set_fixed_width()

ModelGeneVariances & scran::ModelGeneVariances::set_fixed_width ( double  f = FitVarianceTrend::Defaults::fixed_width)
inline
Parameters
fSee FitVarianceTrend::set_fixed_width() for more details.
Returns
A reference to this ModelGeneVariances object.

◆ set_minimum_window_count()

ModelGeneVariances & scran::ModelGeneVariances::set_minimum_window_count ( int  c = FitVarianceTrend::Defaults::minimum_window_count)
inline
Parameters
cSee FitVarianceTrend::set_minimum_window_count() for more details.
Returns
A reference to this ModelGeneVariances object.

◆ set_block_weight_policy()

ModelGeneVariances & scran::ModelGeneVariances::set_block_weight_policy ( WeightPolicy  w = Defaults::block_weight_policy)
inline
Parameters
wWeighting policy to use for averaging statistics across blocks.
Returns
A reference to this ModelGeneVariances instance.

◆ set_variable_block_weight_parameters()

ModelGeneVariances & scran::ModelGeneVariances::set_variable_block_weight_parameters ( VariableBlockWeightParameters  v = Defaults::variable_block_weight_parameters)
inline
Parameters
vParameters for the variable block weights, see variable_block_weight() for more details. Only used when the block weight policy is set to WeightPolicy::VARIABLE.
Returns
A reference to this ModelGeneVariances instance.

◆ set_compute_average()

ModelGeneVariances & scran::ModelGeneVariances::set_compute_average ( bool  a = Defaults::compute_average)
inline
Parameters
aWhether to compute the average of each statistic across blocks. Note that this only affects the run_blocked() method that returns a BlockResults object.
Returns
A reference to this ModelGeneVariances object.

◆ set_num_threads()

ModelGeneVariances & scran::ModelGeneVariances::set_num_threads ( int  n = Defaults::num_threads)
inline
Parameters
nNumber of threads to use.
Returns
A reference to this ModelGeneVariances object.

◆ run() [1/2]

template<typename Value_ , typename Index_ , typename Stat_ >
void scran::ModelGeneVariances::run ( const tatami::Matrix< Value_, Index_ > *  mat,
Stat_ *  means,
Stat_ *  variances,
Stat_ *  fitted,
Stat_ *  residuals 
) const
inline

Compute and model the per-feature variances from a log-expression matrix. This returns the mean and variance for each feature, as well as the fitted value and residuals from the mean-variance trend fitted across features.

Template Parameters
Value_Data type of the matrix.
Index_Integer type for the row/column indices.
Stat_Floating-point type for the output statistics.
Parameters
matPointer to a feature-by-cells tatami matrix containing log-expression values.
[out]meansPointer to an output array of length equal to the number of rows in mat, used to store the mean of each feature.
[out]variancesPointer to an output array of length equal to the number of rows in mat, used to store the variance of each feature.
[out]fittedPointer to an output array of length equal to the number of rows in mat, used to store the fitted value of the trend.
[out]residualsPointer to an output array of length equal to the number of rows in mat, used to store the residual from the trend for each feature.

◆ run_blocked() [1/3]

template<typename Value_ , typename Index_ , typename Block_ , typename Stat_ >
void scran::ModelGeneVariances::run_blocked ( const tatami::Matrix< Value_, Index_ > *  mat,
const Block_ *  block,
std::vector< Stat_ * >  means,
std::vector< Stat_ * >  variances,
std::vector< Stat_ * >  fitted,
std::vector< Stat_ * >  residuals,
Stat_ *  ave_means,
Stat_ *  ave_variances,
Stat_ *  ave_fitted,
Stat_ *  ave_residuals 
) const
inline

Compute and model the per-feature variances from a log-expression matrix with blocking. The mean and variance of each gene is computed separately for all cells in each block, and a separate trend is fitted to each block to obtain residuals. This ensures that sample and batch effects do not confound the variance estimates.

We also compute the average of each statistic across blocks, using the weighting strategy described in weight_block(). The average residual is particularly useful for feature selection with ChooseHVGs.

Template Parameters
Value_Data type of the matrix.
Index_Integer type for the row/column indices.
Block_Integer type to hold the block IDs.
Stat_Floating-point type for the output statistics.
Parameters
matPointer to a feature-by-cells tatami matrix containing log-expression values.
[in]blockPointer to an array of length equal to the number of cells, containing a 0-based block ID for each cell - see tabulate_ids() for more details. This can also be a nullptr, in which case all cells are assumed to belong to the same block.
[out]meansVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the mean of each feature in the corresponding block of cells.
[out]variancesVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the variance of each feature in the corresponding block of cells.
[out]fittedVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the fitted value of the trend for each feature in the corresponding block of cells.
[out]residualsVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the residual from the trend for each feature in the corresponding block of cells.
[out]ave_meansPointer to an array of length equal to the number of rows in mat, storing the average mean across blocks for each gene. If nullptr, the average calculation is skipped.
[out]ave_variancesPointer to an array of length equal to the number of rows in mat, storing the average variance across blocks for each gene. If nullptr, the average calculation is skipped.
[out]ave_fittedPointer to an array of length equal to the number of rows in mat, storing the average fitted value across blocks for each gene. If nullptr, the average calculation is skipped.
[out]ave_residualsPointer to an array of length equal to the number of rows in mat, storing the average residual across blocks for each gene. If nullptr, the average calculation is skipped.

◆ run_blocked() [2/3]

template<typename Value_ , typename Index_ , typename Block_ , typename Stat_ >
void scran::ModelGeneVariances::run_blocked ( const tatami::Matrix< Value_, Index_ > *  mat,
const Block_ *  block,
std::vector< Stat_ * >  means,
std::vector< Stat_ * >  variances,
std::vector< Stat_ * >  fitted,
std::vector< Stat_ * >  residuals 
) const
inline

Compute and model the per-feature variances from a log-expression matrix with blocking. This overload omits the calculation of the averaged statistics across blocks.

Template Parameters
Value_Data type of the matrix.
Index_Integer type for the row/column indices.
Block_Integer type to hold the block IDs.
Stat_Floating-point type for the output statistics.
Parameters
matPointer to a feature-by-cells tatami matrix containing log-expression values.
[in]blockPointer to an array of length equal to the number of cells, containing block identifiers - see tabulate_ids() for more details. This can also be a nullptr, in which case all cells are assumed to belong to the same block.
[out]meansVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the mean of each feature in the corresponding block of cells.
[out]variancesVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the variance of each feature in the corresponding block of cells.
[out]fittedVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the fitted value of the trend for each feature in the corresponding block of cells.
[out]residualsVector of length equal to the number of blocks, containing pointers to output arrays of length equal to the number of rows in mat. Each vector stores the residual from the trend for each feature in the corresponding block of cells.

◆ run() [2/2]

template<typename Value_ , typename Index_ >
Results scran::ModelGeneVariances::run ( const tatami::Matrix< Value_, Index_ > *  mat) const
inline

Compute and model the per-feature variances from a log-expression matrix.

Template Parameters
Value_Data type of the matrix.
Index_Integer type for the row/column indices.
Parameters
matPointer to a feature-by-cells tatami matrix containing log-expression values.
Returns
A Results object containing the results of the variance modelling.

◆ run_blocked() [3/3]

template<typename Value_ , typename Index_ , typename Block_ >
BlockResults scran::ModelGeneVariances::run_blocked ( const tatami::Matrix< Value_, Index_ > *  mat,
const Block_ *  block 
) const
inline

Compute and model the per-feature variances from a log-expression matrix with blocking, see run_blocked() for details.

Template Parameters
Value_Data type of the matrix.
Index_Integer type for the row/column indices.
Block_Integer type, to hold the block IDs.
Parameters
matPointer to a feature-by-cells tatami matrix containing log-expression values.
[in]blockPointer to an array of block identifiers, see run_blocked() for details.
Returns
A BlockResults object containing the results of the variance modelling in each block. An average for each statistic is also computed by default unless set_compute_average() is set to false.

The documentation for this class was generated from the following file: