scran
C++ library for basic single-cell RNA-seq analyses
Loading...
Searching...
No Matches
Namespaces | Classes | Enumerations | Functions
scran Namespace Reference

Functions for single-cell RNA-seq analyses. More...

Namespaces

namespace  differential_analysis
 Utilities for differential analysis.
 
namespace  igraph
 Utilities for working with igraph in libscran.
 
namespace  quick_grouped_size_factors
 Quickly compute grouped size factors.
 

Classes

class  AggregateAcrossCells
 Aggregate expression values across cells. More...
 
class  BuildSnnGraph
 Build a shared nearest-neighbor graph with cells as nodes. More...
 
class  CenterSizeFactors
 Center size factors prior to scaling normalization. More...
 
class  ChooseHvgs
 Choose highly variable genes for downstream analyses. More...
 
struct  ChooseOutlierFilters
 Define outlier filters from the median and MAD. More...
 
class  ChoosePseudoCount
 Choose a pseudo-count for log-transformation. More...
 
class  ClusterSnnGraphLeiden
 Leiden clustering on a shared nearest-neighbor graph. More...
 
class  ClusterSnnGraphMultiLevel
 Multi-level clustering on a shared nearest-neighbor graph. More...
 
class  ClusterSnnGraphWalktrap
 Walktrap clustering on a shared nearest-neighbor graph. More...
 
class  ComputeMedianMad
 Compute the median and MAD from an array of values. More...
 
class  DownsampleByNeighbors
 Downsample a dataset based on its neighbors. More...
 
class  FilterCells
 Filter out low-quality cells. More...
 
class  FitVarianceTrend
 Fit a mean-variance trend to log-count data. More...
 
class  GroupedSizeFactors
 Compute grouped size factors to handle composition bias. More...
 
class  HypergeometricTail
 Compute hypergeometric tail probabilities. More...
 
class  LogNormCounts
 Compute log-normalized expression values. More...
 
class  MedianSizeFactors
 Compute median-based size factors to handle composition bias. More...
 
class  ModelGeneVariances
 Compute and model the per-gene variances in log-expression data. More...
 
class  MultiBatchPca
 Compute PCA after adjusting for differences between batch sizes. More...
 
class  PairwiseEffects
 Compute pairwise effect size between groups of cells. More...
 
class  PerCellAdtQcMetrics
 Compute per-cell quality control metrics from an ADT count matrix. More...
 
class  PerCellCrisprQcMetrics
 Compute per-cell quality control metrics from a CRISPR guide count matrix. More...
 
class  PerCellQcMetrics
 Compute a variety of per-cell quality control metrics from a count matrix. More...
 
class  PerCellRnaQcMetrics
 Compute typical per-cell quality control metrics from an RNA count matrix. More...
 
class  ResidualPca
 Compute PCA after regressing out an uninteresting factor. More...
 
class  SanitizeSizeFactors
 Sanitize invalid size factors. More...
 
class  ScaleByNeighbors
 Scale multi-modal embeddings to adjust for differences in variance. More...
 
class  ScoreFeatureSet
 Compute per-cell scores for a given feature set. More...
 
class  ScoreMarkers
 Score each gene as a candidate marker for each group of cells. More...
 
class  SimplePca
 Perform a simple PCA on a gene-cell matrix. More...
 
struct  SizeFactorValidity
 Validity of size factors. More...
 
class  SuggestAdtQcFilters
 Create filters to identify low-quality cells from ADT-derived QC metrics. More...
 
class  SuggestCrisprQcFilters
 Create filters to identify low-quality cells from CRISPR-derived QC metrics. More...
 
class  SuggestRnaQcFilters
 Create filters to identify low-quality cells from RNA-derived QC metrics. More...
 
class  SummarizeEffects
 Summarize pairwise effects into summary statistics per group. More...
 
struct  VariableBlockWeightParameters
 Parameters for variable_block_weight(). More...
 

Enumerations

enum class  WeightPolicy : char { NONE , VARIABLE , EQUAL }
 

Functions

template<typename T >
SizeFactorValidity validate_size_factors (size_t n, const T *size_factors)
 
template<bool check_nan_ = true, typename Stat_ , typename Output_ >
void average_vectors (size_t n, std::vector< Stat_ * > in, Output_ *out)
 
template<bool check_nan_ = true, typename Output_ = double, typename Stat_ >
std::vector< Output_ > average_vectors (size_t n, std::vector< Stat_ * > in)
 
template<bool check_nan_ = true, typename Stat_ , typename Weight_ , typename Output_ >
void average_vectors_weighted (size_t n, std::vector< Stat_ * > in, const Weight_ *w, Output_ *out)
 
template<bool check_nan_ = true, typename Output_ = double, typename Stat_ , typename Weight_ >
std::vector< Output_ > average_vectors_weighted (size_t n, std::vector< Stat_ * > in, const Weight_ *w)
 
template<typename Id_ >
size_t count_ids (size_t length, const Id_ *ids)
 
template<typename Output_ = int, bool allow_zeros_ = false, typename Id_ >
std::vector< Output_ > tabulate_ids (size_t length, const Id_ *ids, bool allow_zeros=false)
 
double variable_block_weight (double s, const VariableBlockWeightParameters &params)
 
template<typename Size_ >
std::vector< double > compute_block_weights (const std::vector< Size_ > &sizes, WeightPolicy policy, const VariableBlockWeightParameters &param)
 
template<bool retain, class Vector , typename Subset >
Vector subset_vector (const Vector &vec, const Subset *sub)
 
template<typename T >
std::vector< T * > vector_to_pointers (std::vector< std::vector< T > > &input)
 
template<typename T >
std::vector< const T * > vector_to_pointers (const std::vector< std::vector< T > > &input)
 
template<typename T >
std::vector< std::vector< T * > > vector_to_pointers (std::vector< std::vector< std::vector< T > > > &input)
 

Detailed Description

Functions for single-cell RNA-seq analyses.

Enumeration Type Documentation

◆ WeightPolicy

enum class scran::WeightPolicy : char
strong

Policy to use for weighting blocks based on their size, i.e., the number of cells in each block. This controls the calculation of weighted averages across blocks.

  • NONE: no weighting is performed. Larger blocks will contribute more to the weighted average.
  • EQUAL: each block receives equal weight, regardless of its size. Equivalent to averaging across blocks without weights.
  • VARIABLE: each batch is weighted using the logic in variable_block_weight(). This penalizes small blocks with unreliable statistics while equally weighting all large blocks.

Function Documentation

◆ validate_size_factors()

template<typename T >
SizeFactorValidity scran::validate_size_factors ( size_t  n,
const T *  size_factors 
)

Check whether there are any invalid size factors. Size factors are only technically valid if they are finite and positive.

Template Parameters
TFloating-point type for the size factors.
Parameters
nNumber of size factors.
[in]size_factorsPointer to an array of size factors of length n.
Returns
Validation results, indicating whether any zero or non-finite size factors exist.

◆ average_vectors() [1/2]

template<bool check_nan_ = true, typename Stat_ , typename Output_ >
void scran::average_vectors ( size_t  n,
std::vector< Stat_ * >  in,
Output_ *  out 
)

Average parallel elements across multiple arrays.

Template Parameters
check_nan_Whether to check for NaNs. If true, NaNs are ignored in the average calculations for each element, at the cost of some efficiency.
Stat_Type of the input statistic, typically floating point.
Output_Floating-point output type.
Parameters
nLength of each array.
[in]inVector of pointers to input arrays of length n.
[out]outPointer to an output array of length n. On completion, out is filled with the average of all arrays in in. Specifically, each element of out is set to the average of the corresponding elements across all in arrays.

◆ average_vectors() [2/2]

template<bool check_nan_ = true, typename Output_ = double, typename Stat_ >
std::vector< Output_ > scran::average_vectors ( size_t  n,
std::vector< Stat_ * >  in 
)

Average parallel elements across multiple arrays.

Template Parameters
check_nan_Whether to check for NaNs, see average_vectors().
OutputFloating-point output type.
StatType of the input statistic, typically floating point.
Parameters
nLength of each array.
[in]inVector of pointers to input arrays of length n.
Returns
A vector of length n is returned, containing the average of all arrays in in.

◆ average_vectors_weighted() [1/2]

template<bool check_nan_ = true, typename Stat_ , typename Weight_ , typename Output_ >
void scran::average_vectors_weighted ( size_t  n,
std::vector< Stat_ * >  in,
const Weight_ *  w,
Output_ *  out 
)

Compute a weighted average of parallel elements across multiple arrays.

Template Parameters
check_nan_Whether to check for NaNs. If true, NaNs are ignored in the average calculations for each element, at the cost of some efficiency.
Stat_Type of the input statistic, typically floating point.
Weight_Type of the weight, typically floating point.
Output_Floating-point output type.
Parameters
nLength of each array.
[in]inVector of pointers to input arrays of length n.
[in]wPointer to an array of length equal to in.size(), containing the weight to use for each input array. Weights should be non-negative and finite.
[out]outPointer to an output array of length n. On output, out is filled with the weighted average of all arrays in in. Specifically, each element of out is set to the weighted average of the corresponding elements across all in arrays.

◆ average_vectors_weighted() [2/2]

template<bool check_nan_ = true, typename Output_ = double, typename Stat_ , typename Weight_ >
std::vector< Output_ > scran::average_vectors_weighted ( size_t  n,
std::vector< Stat_ * >  in,
const Weight_ *  w 
)

Compute a weighted average of parallel elements across multiple arrays.

Template Parameters
check_nan_Whether to check for NaNs, see average_vectors_weighted() for details.
Output_Floating-point output type.
Weight_Type of the weight, typically floating point.
Stat_Type of the input statistic, typically floating point.
Parameters
nLength of each array.
[in]inVector of pointers to input arrays of the same length.
[in]wPointer to an array of length equal to in.size(), containing the weight to use for each input array. Weights should be non-negative and finite.
Returns
A vector is returned containing with the average of all arrays in in.

◆ count_ids()

template<typename Id_ >
size_t scran::count_ids ( size_t  length,
const Id_ *  ids 
)

Count the number of unique 0-based IDs, e.g., for block or group assignments. All IDs are assumed to be integers in [0, x) where x is the return value of this function.

Template Parameters
Id_Integer type for the IDs.
Parameters
lengthLength of the array in ids.
[in]idsPointer to an array containing 0-based IDs of some kind.
Returns
Number of IDs, or 0 if length = 0.

◆ tabulate_ids()

template<typename Output_ = int, bool allow_zeros_ = false, typename Id_ >
std::vector< Output_ > scran::tabulate_ids ( size_t  length,
const Id_ *  ids,
bool  allow_zeros = false 
)

Count the frequency of 0-based IDs, e.g., for block or group assignments. All IDs are assumed to be integers in [0, x) where x is the number of unique IDs.

Template Parameters
Output_Numeric type for the output frequencies.
Id_Integer type for the IDs.
Parameters
lengthLength of the array in ids.
[in]idsPointer to an array containing 0-based IDs of some kind.
allow_zerosWhether to throw an error if frequencies of zero are detected.
Returns
Vector of length equal to the number of IDs, containing the frequency of each ID. An error is raised if an ID has zero frequency and allow_zeros = false.

◆ variable_block_weight()

double scran::variable_block_weight ( double  s,
const VariableBlockWeightParameters params 
)
inline

Weight each block of cells for use in computing a weighted average across blocks. The weight for each block is calcualted from the size of that block.

  • If the block is empty smaller than some lower bound, it has zero weight.
  • If the block is greater than some upper bound, it has weight of 1.
  • Otherwise, the block has weight proportional to its size, increasing linearly from 0 to 1 between the two bounds.

Blocks that are "large enough" are considered to be equally trustworthy and receive the same weight, ensuring that each block contributes equally to the weighted average. By comparison, very small blocks receive lower weight as their statistics are generally less stable. If both cap and block_size are zero, the weight is also set to zero.

Parameters
sSize of the block, in terms of the number of cells in that block.
paramsParameters for the weight calculation, consisting of the lower and upper bounds.
Returns
Weight of the block, to use for computing a weighted average across blocks.

◆ compute_block_weights()

template<typename Size_ >
std::vector< double > scran::compute_block_weights ( const std::vector< Size_ > &  sizes,
WeightPolicy  policy,
const VariableBlockWeightParameters param 
)

Compute block weights for multiple blocks based on their size and the weighting policy. For variable weights, this function will call variable_block_weight() for each block.

Template Parameters
Size_Numeric type for the block size.
Parameters
sizesVector of block sizes.
policyPolicy for weighting blocks of different sizes.
paramParameters for the variable block weights.
Returns
Vector of block weights.

◆ subset_vector()

template<bool retain, class Vector , typename Subset >
Vector scran::subset_vector ( const Vector &  vec,
const Subset *  sub 
)

Subset a vector to retain/discard elements.

Template Parameters
retainShould the non-zero elements in sub be retained (true) or discarded (false)?
VectorA vector class that has a size() method, a [] operator, and a constructor that initializes the vector to a specified length.
SubsetInteger/boolean type for the subsetting vector.
Parameters
vecA vector of arbitrary (copy-able) elements.
subPointer to an array of integer/boolean elements of length equal to vec. Non-zero values indicate that an element should be retained (if retain=true) or discarded (otherwise).
Returns
A vector containing the desired subset of elements in vec.

◆ vector_to_pointers() [1/3]

template<typename T >
std::vector< T * > scran::vector_to_pointers ( std::vector< std::vector< T > > &  input)

Extract a vector of pointers from a vector of vectors. This is a convenient utility as many scran functions accept the former but can return the latter.

Template Parameters
TType of data.
Parameters
inputVector of vector of values.
Returns
A vector of pointers to each inner vector in input.

◆ vector_to_pointers() [2/3]

template<typename T >
std::vector< const T * > scran::vector_to_pointers ( const std::vector< std::vector< T > > &  input)

Extract a vector of const pointers from a vector of vectors.

Template Parameters
TType of data.
Parameters
inputVector of vector of values.
Returns
A vector of pointers to each inner vector in input.

◆ vector_to_pointers() [3/3]

template<typename T >
std::vector< std::vector< T * > > scran::vector_to_pointers ( std::vector< std::vector< std::vector< T > > > &  input)

Extract a vector of vector of pointers from a vector of vectors of vectors..

Template Parameters
TType of data.
Parameters
inputVector of vector of vector of values.
Returns
A vector of vector of pointers to each inner vector in input.