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

Aggregate expression values across cells. More...

#include <AggregateAcrossCells.hpp>

Classes

struct  Combinations
 Unique combinations of factors. More...
 
struct  Defaults
 Default parameters for aggregation. More...
 
struct  Results
 Container for the aggregation results. More...
 

Public Member Functions

template<typename Data , typename Index , typename Factor , typename Sum , typename Detected >
void run (const tatami::Matrix< Data, Index > *input, const Factor *factor, std::vector< Sum * > sums, std::vector< Detected * > detected)
 
AggregateAcrossCellsset_compute_sums (bool c=Defaults::compute_sums)
 
AggregateAcrossCellsset_compute_detected (bool c=Defaults::compute_detected)
 
AggregateAcrossCellsset_num_threads (int n=Defaults::num_threads)
 
template<typename Sum = double, typename Detected = int, typename Data , typename Index , typename Factor >
Results< Sum, Detectedrun (const tatami::Matrix< Data, Index > *input, const Factor *factor)
 

Static Public Member Functions

template<typename Factor , typename Combined >
static Combinations< Factorcombine_factors (size_t n, std::vector< const Factor * > factors, Combined *combined)
 
template<typename Combined = int, typename Factor >
static std::pair< Combinations< Factor >, std::vector< Combined > > combine_factors (size_t n, std::vector< const Factor * > factors)
 

Detailed Description

Aggregate expression values across cells.

This class computes the sum of expression values for each grouping of cells, typically for the creation of pseudo-bulk expression profiles for cluster/sample combinations. Expression values are generally expected to be counts, though the same code can be trivially re-used to compute the average log-expression. We can also report the number of cells with detected (i.e., positive) expression values in each grouping.

Member Function Documentation

◆ combine_factors() [1/2]

static Combinations< Factor > scran::AggregateAcrossCells::combine_factors ( size_t  n,
std::vector< const Factor * >  factors,
Combined combined 
)
inlinestatic
Template Parameters
FactorFactor type.
CombinedType of the combined factor. May need to be different from Factor if the latter does not have enough unique levels.
Parameters
nNumber of observations (i.e., cells).
[in]factorsPointers to arrays of length n, each containing a different factor.
[out]combinedPointer to an array of length n, in which the combined factor is to be stored.
Returns
A Combinations object is returned containing the unique combinations of levels observed in factors. A combined factor is saved to combined; values are indices of combinations in the output Combinations object.

This function compresses multiple factors into a single combined factor, which can then be used as the factor in run(). In this manner, AggregateAcrossCells can easily handle aggregation across any number of factors.

◆ combine_factors() [2/2]

template<typename Combined = int, typename Factor >
static std::pair< Combinations< Factor >, std::vector< Combined > > scran::AggregateAcrossCells::combine_factors ( size_t  n,
std::vector< const Factor * >  factors 
)
inlinestatic
Template Parameters
FactorFactor type.
CombinedType of the combined factor. May need to be different from Factor if the latter does not have enough unique levels.
Parameters
nNumber of observations (i.e., cells).
[in]factorsPointers to arrays of length n, each containing a different factor.
Returns
A pair containing:
  • A Combinations object, containing the unique combinations of levels observed in factors.
  • A vector of length n containing the combined factor.

See the other combine_factors() method for details.

◆ run() [1/2]

void scran::AggregateAcrossCells::run ( const tatami::Matrix< Data, Index > *  input,
const Factor factor,
std::vector< Sum * >  sums,
std::vector< Detected * >  detected 
)
inline

Compute the per-gene sum and number of cells with detected expression for each level of a grouping factor.

Template Parameters
DataType of data in the input matrix, should be numeric.
IndexInteger type of index in the input matrix.
FactorInteger type of the factor.
SumType of the sum, usually the same as Data.
DetectedType for the number of detected cells, usually integer.
Parameters
inputThe input matrix where rows are features and columns are cells.
[in]factorPointer to an array of length equal to the number of columns of input, containing the factor level for each cell. All levels should be integers in $[0, N)$ where $N$ is the number of unique levels.
[out]sumsVector of length $N$ (see factor), containing pointers to arrays of length equal to the number of columns of input. These will be filled with the summed expression across all cells in the corresponding level for each gene. Alternatively, if the vector is of length 0, no sums will be computed.
[out]detectedVector of length $N$ (see factor), containing pointers to arrays of length equal to the number of columns of input. These will be filled with the number of cells with detected expression in the corresponding level for each gene. Alternatively, if the vector is of length 0, no numbers will be computed.

◆ set_compute_sums()

AggregateAcrossCells & scran::AggregateAcrossCells::set_compute_sums ( bool  c = Defaults::compute_sums)
inline
Parameters
cWhether to compute the sum within each factor level.
Returns
A reference to this AggregateAcrossCells object.

This function only affects run() when sums and detected are not supplied as inputs.

◆ set_compute_detected()

AggregateAcrossCells & scran::AggregateAcrossCells::set_compute_detected ( bool  c = Defaults::compute_detected)
inline
Parameters
cWhether to compute the number of detected cells within each factor level.
Returns
A reference to this AggregateAcrossCells object.

This function only affects run() when sums and detected are not supplied as inputs.

◆ set_num_threads()

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

◆ run() [2/2]

Results< Sum, Detected > scran::AggregateAcrossCells::run ( const tatami::Matrix< Data, Index > *  input,
const Factor factor 
)
inline
Template Parameters
SumType of the sum, should be numeric.
DetectedType for the number of detected cells, usually integer.
DataType of data in the input matrix, should be numeric.
IndexInteger type of index in the input matrix.
FactorInteger type of the factor.
Parameters
inputThe input matrix where rows are features and columns are cells.
[in]factorPointer to an array of length equal to the number of columns of input, containing the factor level for each cell. All levels should be integers in $[0, N)$ where $N$ is the number of unique levels.
Returns
A Results object is returned containing the summed expression and the number of detected cells within each factor level across all genes.

This function will respect any user-supplied setting of set_compute_sums() and set_compute_detected(). If either/both are false, the corresponding statistic(s) will not be computed and the corresponding vector in Results will be empty.


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