scran_pca
Principal component analysis for single-cell data
Loading...
Searching...
No Matches
scran_pca::BlockedPcaOptions< EigenVector_ > Struct Template Reference

Options for blocked_pca(). More...

#include <blocked_pca.hpp>

Public Attributes

int number = 25
 
bool scale = false
 
bool transpose = true
 
scran_blocks::WeightPolicy block_weight_policy = scran_blocks::WeightPolicy::VARIABLE
 
scran_blocks::VariableWeightParameters variable_block_weight_parameters
 
bool components_from_residuals = true
 
bool realize_matrix = true
 
int num_threads = 1
 
irlba::Options< EigenVector_ > irlba_options
 

Detailed Description

template<typename EigenVector_ = Eigen::VectorXd>
struct scran_pca::BlockedPcaOptions< EigenVector_ >

Options for blocked_pca().

Template Parameters
EigenVector_A floating-point Eigen::Vector class.

Member Data Documentation

◆ block_weight_policy

template<typename EigenVector_ = Eigen::VectorXd>
scran_blocks::WeightPolicy scran_pca::BlockedPcaOptions< EigenVector_ >::block_weight_policy = scran_blocks::WeightPolicy::VARIABLE

Policy for weighting the contribution of blocks of different size.

The default of scran_blocks::WeightPolicy::VARIABLE is to define equal weights for blocks once they reach a certain size (see BlockedPcaOptions::variable_block_weight_parameters). For smaller blocks, the weight is linearly proportional to its size to avoid outsized contributions from very small blocks.

Other options include scran_blocks::WeightPolicy::EQUAL, where all blocks are equally weighted regardless of size; and scran_blocks::WeightPolicy::NONE, where the contribution of each block is proportional to its size.

◆ components_from_residuals

template<typename EigenVector_ = Eigen::VectorXd>
bool scran_pca::BlockedPcaOptions< EigenVector_ >::components_from_residuals = true

Whether to compute the principal components from the residuals. If false, only the rotation vector is computed from the residuals and the original expression values are projected onto the new axes. This avoids strong assumptions about the nature of the differences between blocks as discussed in blocked_pca().

◆ irlba_options

template<typename EigenVector_ = Eigen::VectorXd>
irlba::Options<EigenVector_> scran_pca::BlockedPcaOptions< EigenVector_ >::irlba_options

Further options to pass to irlba::compute().

◆ num_threads

template<typename EigenVector_ = Eigen::VectorXd>
int scran_pca::BlockedPcaOptions< EigenVector_ >::num_threads = 1

Number of threads to use. The parallelization scheme is determined by tatami::parallelize() and irlba::parallelize(). Note that the exact values returned by blocked_pca() will change slightly with different num_threads, due to (deterministic) differences in the order of floating-point summations.

◆ number

template<typename EigenVector_ = Eigen::VectorXd>
int scran_pca::BlockedPcaOptions< EigenVector_ >::number = 25

Number of the top principal components (PCs) to compute. Retaining more PCs will capture more biological signal at the cost of increasing noise and compute time. If this is greater than the maximum number of PCs (i.e., the smaller dimension of the input matrix), only the maximum number of PCs will be reported in the results.

◆ realize_matrix

template<typename EigenVector_ = Eigen::VectorXd>
bool scran_pca::BlockedPcaOptions< EigenVector_ >::realize_matrix = true

Whether to realize tatami::Matrix objects into an appropriate in-memory format before PCA. This is typically faster but increases memory usage.

◆ scale

template<typename EigenVector_ = Eigen::VectorXd>
bool scran_pca::BlockedPcaOptions< EigenVector_ >::scale = false

Should genes be scaled to unit variance? This ensures that each gene contributes equally to the PCA, favoring consistent variation across many genes rather than large variation in a few genes. In the presence of a blocking factor, each gene's variance is calculated as a weighted sum of the variances from each block. Genes with zero variance are ignored.

◆ transpose

template<typename EigenVector_ = Eigen::VectorXd>
bool scran_pca::BlockedPcaOptions< EigenVector_ >::transpose = true

Should the PC matrix be transposed on output? If true, the output matrix is column-major with cells in the columns, which is compatible with downstream libscran steps.

◆ variable_block_weight_parameters

template<typename EigenVector_ = Eigen::VectorXd>
scran_blocks::VariableWeightParameters scran_pca::BlockedPcaOptions< EigenVector_ >::variable_block_weight_parameters

Parameters for the variable block weights, including the threshold at which blocks are considered to be large enough to have equal weight. Only used when BlockedPcaOptions::block_weight_policy = scran_blocks::WeightPolicy::VARIABLE.


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